-
Notifications
You must be signed in to change notification settings - Fork 117
/
mf6bmi.f90
1262 lines (1107 loc) · 43.7 KB
/
mf6bmi.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!> @brief This module contains the MODFLOW 6 BMI
!!
!! This BMI interface matches the CSDMS standard, with a few modifications:
!!
!! - This interface will build into a shared library that can be called from other
!! executables and scripts, not necessarily written in Fortran. Therefore we have
!! omitted the type-boundness of the routines, since they cannot have the
!! bind(C,"...") attribute.
!! - MODFLOW has internal data arrays with rank > 1 that we would like to expose. An
!! example would be access to data in the BOUND array of GWF boundary packages (BndType).
!! The get_value_ptr calls below support this, returning a C-style pointer to the arrays
!! and methods have been added to query the variable's rank and shape.
!!
!! Note on style: BMI apparently uses underscores, we use underscores in some
!! places but camelcase in other. Since this is a dedicated BMI interface module,
!! we'll use underscores here as well.
!<
module mf6bmi
use mf6bmiUtil
use mf6bmiError
use Mf6CoreModule
use TdisModule, only: kper, kstp
use iso_c_binding, only: c_int, c_char, c_double, C_NULL_CHAR, c_loc, c_ptr, &
c_f_pointer
use KindModule, only: DP, I4B, LGP
use ConstantsModule, only: LENMEMPATH, LENVARNAME
use CharacterStringModule
use MemoryManagerModule, only: mem_setptr, get_mem_elem_size, get_isize, &
get_mem_rank, get_mem_shape, get_mem_type, &
memorylist, get_from_memorylist
use MemoryTypeModule, only: MemoryType
use MemoryHelperModule, only: create_mem_address
use SimVariablesModule, only: simstdout, istdout
use InputOutputModule, only: getunit
implicit none
integer(c_int), bind(C, name="ISTDOUTTOFILE") :: istdout_to_file = 1 !< output control: =0 to screen, >0 to file
!DIR$ ATTRIBUTES DLLEXPORT :: istdout_to_file
contains
function bmi_get_component_name(name) result(bmi_status) &
bind(C, name="get_component_name")
!DIR$ ATTRIBUTES DLLEXPORT :: bmi_get_component_name
! -- dummy variables
character(kind=c_char), intent(out) :: name(BMI_LENCOMPONENTNAME)
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
name = string_to_char_array('MODFLOW 6', 9)
bmi_status = BMI_SUCCESS
end function bmi_get_component_name
!> @brief Initialize the computational core
!!
!! It is required to have the MODFLOW 6 configuration file 'mfsim.nam'
!! available in the current working directory when calling this routine.
!!
!! NOTE: initialization should always be matched with a call to
!! bmi_finalize(). However, we currently do not support the reinitialization
!! of a model in the same memory space... You would have to create a new
!! process for that.
!<
function bmi_initialize() result(bmi_status) bind(C, name="initialize")
!DIR$ ATTRIBUTES DLLEXPORT :: bmi_initialize
! -- dummy variables
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
if (istdout_to_file > 0) then
! -- open stdout file mfsim.stdout
istdout = getunit()
!
! -- set STDOUT to a physical file unit
open (unit=istdout, file=simstdout)
end if
!
! -- initialize MODFLOW 6
call Mf6Initialize()
bmi_status = BMI_SUCCESS
end function bmi_initialize
!> @brief Perform a computational time step
!!
!! It will prepare the timestep, call the calculation routine
!! on all the solution groups in the simulation, and finalize
!! the timestep by printing out diagnostics and writing output.
!! It can be called in succession to perform multiple steps up
!! to the simulation's end time is reached.
!<
function bmi_update() result(bmi_status) bind(C, name="update")
!DIR$ ATTRIBUTES DLLEXPORT :: bmi_update
! -- dummy variables
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
logical :: hasConverged
hasConverged = Mf6Update()
bmi_status = BMI_SUCCESS
end function bmi_update
!> @brief Clean up the initialized simulation
!!
!! Performs teardown tasks for the initialized simulation, this
!! call should match the call to bmi_initialize()
!<
function bmi_finalize() result(bmi_status) bind(C, name="finalize")
!DIR$ ATTRIBUTES DLLEXPORT :: bmi_finalize
! -- modules
use SimVariablesModule, only: iforcestop
! -- dummy variables
integer(kind=c_int) :: bmi_status !< BMI status code
! we don't want a full stop() here, this disables it:
iforcestop = 0
call Mf6Finalize()
bmi_status = BMI_SUCCESS
end function bmi_finalize
!> @brief Get the start time of the simulation
!!
!! As MODFLOW currently does not have internal time, this will be
!! returning 0.0 for now. New version...
!<
function get_start_time(start_time) result(bmi_status) &
bind(C, name="get_start_time")
!DIR$ ATTRIBUTES DLLEXPORT :: get_start_time
! -- dummy variables
real(kind=c_double), intent(out) :: start_time !< start time
integer(kind=c_int) :: bmi_status !< BMI status code
start_time = 0.0_DP
bmi_status = BMI_SUCCESS
end function get_start_time
!> @brief Get the end time of the simulation
!!
!! As MODFLOW does currently does not have internal time, this will be
!! equal to the total runtime.
!<
function get_end_time(end_time) result(bmi_status) bind(C, name="get_end_time")
!DIR$ ATTRIBUTES DLLEXPORT :: get_end_time
! -- modules
use TdisModule, only: totalsimtime
! -- dummy variables
real(kind=c_double), intent(out) :: end_time !< end time
integer(kind=c_int) :: bmi_status !< BMI status code
end_time = totalsimtime
bmi_status = BMI_SUCCESS
end function get_end_time
!> @brief Get the current time of the simulation
!!
!! As MODFLOW currently does not have internal time, this will be
!! equal to the time passed w.r.t. the start time of the simulation.
!<
function get_current_time(current_time) result(bmi_status) &
bind(C, name="get_current_time")
!DIR$ ATTRIBUTES DLLEXPORT :: get_current_time
! -- modules
use TdisModule, only: totim
! -- dummy variables
real(kind=c_double), intent(out) :: current_time !< current time
integer(kind=c_int) :: bmi_status !< BMI status code
current_time = totim
bmi_status = BMI_SUCCESS
end function get_current_time
!> @brief Get the time step for the simulation
!!
!! Note that the returned value may vary between and within stress periods,
!! depending on your time discretization settings in the TDIS package.
!<
function get_time_step(time_step) result(bmi_status) &
bind(C, name="get_time_step")
!DIR$ ATTRIBUTES DLLEXPORT :: get_time_step
! -- modules
use TdisModule, only: delt
! -- dummy variables
real(kind=c_double), intent(out) :: time_step !< current time step
integer(kind=c_int) :: bmi_status !< BMI status code
time_step = delt
bmi_status = BMI_SUCCESS
end function get_time_step
!> @brief Get the number of input variables in the simulation
!!
!! This concerns all variables stored in the memory manager
!<
function get_input_item_count(count) result(bmi_status) &
bind(C, name="get_input_item_count")
!DIR$ ATTRIBUTES DLLEXPORT :: get_input_item_count
! -- dummy variables
integer(kind=c_int), intent(out) :: count !< the number of input variables
integer(kind=c_int) :: bmi_status !< BMI status code
count = memorylist%count()
bmi_status = BMI_SUCCESS
end function get_input_item_count
!> @brief Get the number of output variables in the simulation
!!
!! This concerns all variables stored in the memory manager
!<
function get_output_item_count(count) result(bmi_status) &
bind(C, name="get_output_item_count")
!DIR$ ATTRIBUTES DLLEXPORT :: get_output_item_count
! -- dummy variables
integer(kind=c_int), intent(out) :: count !< the number of output variables
integer(kind=c_int) :: bmi_status !< BMI status code
count = memorylist%count()
bmi_status = BMI_SUCCESS
end function get_output_item_count
!> @brief Returns all input variables in the simulation
!!
!! This functions returns the full address for all variables in the
!! memory manager
!!
!! The array @p c_names should be pre-allocated of proper size:
!!
!! size = BMI_LENVARADDRESS * get_input_item_count()
!!
!! The strings will be written contiguously with stride equal to
!! BMI_LENVARADDRESS and nul-terminated where the trailing spaces start:
!!
!! c_names = 'variable_address_1\x00 ... variable_address_2\x00 ... ' etc.
!<
function get_input_var_names(c_names) result(bmi_status) &
bind(C, name="get_input_var_names")
!DIR$ ATTRIBUTES DLLEXPORT :: get_input_var_names
! -- dummy variables
character(kind=c_char, len=1), intent(inout) :: c_names(*) !< array with memory paths for input variables
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
integer(I4B) :: imem, start, i
type(MemoryType), pointer :: mt => null()
character(len=LENMEMADDRESS) :: var_address
start = 1
do imem = 1, memorylist%count()
mt => memorylist%Get(imem)
var_address = create_mem_address(mt%path, mt%name)
do i = 1, len(trim(var_address))
c_names(start + i - 1) = var_address(i:i)
end do
c_names(start + i) = c_null_char
start = start + BMI_LENVARADDRESS
end do
bmi_status = BMI_SUCCESS
end function get_input_var_names
!> @brief Returns all output variables in the simulation
!!
!! This function works analogously to get_input_var_names(),
!! and currently returns the same set of memory variables,
!! which is all of them!
!<
function get_output_var_names(c_names) result(bmi_status) &
bind(C, name="get_output_var_names")
!DIR$ ATTRIBUTES DLLEXPORT :: get_output_var_names
! -- dummy variables
character(kind=c_char, len=1), intent(inout) :: c_names(*) !< array with memory paths for output variables
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
integer(I4B) :: imem, start, i
type(MemoryType), pointer :: mt => null()
character(len=LENMEMADDRESS) :: var_address
start = 1
do imem = 1, memorylist%count()
mt => memorylist%Get(imem)
var_address = create_mem_address(mt%path, mt%name)
do i = 1, len(trim(var_address))
c_names(start + i - 1) = var_address(i:i)
end do
c_names(start + i) = c_null_char
start = start + BMI_LENVARADDRESS
end do
bmi_status = BMI_SUCCESS
end function get_output_var_names
!> @brief Get the size (in bytes) of a single element of a variable
!<
function get_var_itemsize(c_var_address, var_size) result(bmi_status) &
bind(C, name="get_var_itemsize")
!DIR$ ATTRIBUTES DLLEXPORT :: get_var_itemsize
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
integer(kind=c_int), intent(out) :: var_size !< size of the element in bytes
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
call get_mem_elem_size(var_name, mem_path, var_size)
if (var_size == -1) bmi_status = BMI_FAILURE
end function get_var_itemsize
!> @brief Get size of the variable, in bytes
!<
function get_var_nbytes(c_var_address, var_nbytes) result(bmi_status) &
bind(C, name="get_var_nbytes")
!DIR$ ATTRIBUTES DLLEXPORT :: get_var_nbytes
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
integer(kind=c_int), intent(out) :: var_nbytes !< size in bytes
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
integer(I4B) :: var_size, isize
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
call get_mem_elem_size(var_name, mem_path, var_size)
if (var_size == -1) bmi_status = BMI_FAILURE
call get_isize(var_name, mem_path, isize)
if (isize == -1) bmi_status = BMI_FAILURE
var_nbytes = var_size * isize
end function get_var_nbytes
!> @brief Copy the double precision values of a variable into the array
!!
!! The copied variable us located at @p c_var_address. The caller should
!! provide @p c_arr_ptr pointing to an array of the proper shape (the
!! BMI function get_var_shape() can be used to create it). Multi-dimensional
!! arrays are supported.
!<
function get_value_double(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_double")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_double
! -- modules
use MemorySetHandlerModule, only: on_memory_set
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the double precision array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
integer(I4B) :: rank
real(DP), pointer :: src_ptr, tgt_ptr
real(DP), dimension(:), pointer, contiguous :: src1D_ptr, tgt1D_ptr
real(DP), dimension(:, :), pointer, contiguous :: src2D_ptr, tgt2D_ptr
real(DP), dimension(:, :, :), pointer, contiguous :: src3D_ptr, tgt3D_ptr
integer(I4B) :: i, j, k
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
! convert pointer and copy data from memory manager into
! the passed array, using loops to avoid stack overflow
rank = -1
call get_mem_rank(var_name, mem_path, rank)
if (rank == 0) then
call mem_setptr(src_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt_ptr)
tgt_ptr = src_ptr
else if (rank == 1) then
call mem_setptr(src1D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt1D_ptr, shape(src1D_ptr))
do i = 1, size(tgt1D_ptr)
tgt1D_ptr(i) = src1D_ptr(i)
end do
else if (rank == 2) then
call mem_setptr(src2D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt2D_ptr, shape(src2D_ptr))
do j = 1, size(tgt2D_ptr, 2)
do i = 1, size(tgt2D_ptr, 1)
tgt2D_ptr(i, j) = src2D_ptr(i, j)
end do
end do
else if (rank == 3) then
call mem_setptr(src3D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt3D_ptr, shape(src3D_ptr))
do k = 1, size(tgt3D_ptr, 3)
do j = 1, size(tgt3D_ptr, 2)
do i = 1, size(tgt3D_ptr, 1)
tgt3D_ptr(i, j, k) = src3D_ptr(i, j, k)
end do
end do
end do
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_double
!> @brief Copy the integer values of a variable into the array
!!
!! The copied variable is located at @p c_var_address. The caller should
!! provide @p c_arr_ptr pointing to an array of the proper shape (the
!! BMI function get_var_shape() can be used to create it). Multi-dimensional
!! arrays are supported.
!<
function get_value_int(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_int")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_int
! -- modules
use MemorySetHandlerModule, only: on_memory_set
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the integer array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
integer(I4B) :: rank
integer(I4B), pointer :: src_ptr, tgt_ptr
integer(I4B), dimension(:), pointer, contiguous :: src1D_ptr, tgt1D_ptr
integer(I4B), dimension(:, :), pointer, contiguous :: src2D_ptr, tgt2D_ptr
integer(I4B), dimension(:, :, :), pointer, contiguous :: src3D_ptr, tgt3D_ptr
integer(I4B) :: i, j, k
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
! convert pointer and copy data from memory manager into
! the passed array, using loops to avoid stack overflow
rank = -1
call get_mem_rank(var_name, mem_path, rank)
if (rank == 0) then
call mem_setptr(src_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt_ptr)
tgt_ptr = src_ptr
else if (rank == 1) then
call mem_setptr(src1D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt1D_ptr, shape(src1D_ptr))
do i = 1, size(tgt1D_ptr)
tgt1D_ptr(i) = src1D_ptr(i)
end do
else if (rank == 2) then
call mem_setptr(src2D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt2D_ptr, shape(src2D_ptr))
do j = 1, size(tgt2D_ptr, 2)
do i = 1, size(tgt2D_ptr, 1)
tgt2D_ptr(i, j) = src2D_ptr(i, j)
end do
end do
else if (rank == 3) then
call mem_setptr(src3D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt3D_ptr, shape(src3D_ptr))
do k = 1, size(tgt3D_ptr, 3)
do j = 1, size(tgt3D_ptr, 2)
do i = 1, size(tgt3D_ptr, 1)
tgt3D_ptr(i, j, k) = src3D_ptr(i, j, k)
end do
end do
end do
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_int
!> @brief Copy the logical scalar value into the array
!!
!! The copied variable us located at @p c_var_address. The caller should
!! provide @p c_arr_ptr pointing to a scalar array with rank=0.
!<
function get_value_bool(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_bool")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_bool
! -- modules
use MemorySetHandlerModule, only: on_memory_set
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the logical array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
integer(I4B) :: rank
logical(LGP), pointer :: src_ptr, tgt_ptr
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
rank = -1
call get_mem_rank(var_name, mem_path, rank)
! convert pointer
if (rank == 0) then
call mem_setptr(src_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, tgt_ptr)
tgt_ptr = src_ptr
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_bool
!> @brief Copy the string(s) of a variable into the array
!!
!! The copied variable is located at @p c_var_address. The caller should
!! provide @p c_arr_ptr pointing to an array of the proper shape (the
!! BMI function get_var_shape() can be used to create it). For strings
!! currently scalars and 1d arrays (of CharacterStringType) are
!< supported
function get_value_string(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_string")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_string
! -- modules
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the string array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
integer(I4B) :: rank
character(len=:), pointer :: srcstr
character(kind=c_char), pointer :: tgtstr(:)
type(CharacterStringType), dimension(:), pointer, contiguous :: srccharstr1d
character(kind=c_char), pointer :: tgtstr1d(:, :)
character(:), allocatable :: tempstr
integer(I4B) :: i, ilen, isize
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
! single string, or array of strings (CharacterStringType)
rank = -1
call get_mem_rank(var_name, mem_path, rank)
if (rank == 0) then
! a string scalar
call mem_setptr(srcstr, var_name, mem_path)
call get_mem_elem_size(var_name, mem_path, ilen)
call c_f_pointer(c_arr_ptr, tgtstr, shape=[ilen + 1])
tgtstr(1:len(srcstr) + 1) = string_to_char_array(srcstr, len(srcstr))
else if (rank == 1) then
! an array of strings
call mem_setptr(srccharstr1d, var_name, mem_path)
if (.not. associated(srccharstr1d)) then
write (bmi_last_error, fmt_general_err) 'string type not supported in API'
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
! create fortran pointer to C data array
call get_isize(var_name, mem_path, isize)
call get_mem_elem_size(var_name, mem_path, ilen)
call c_f_pointer(c_arr_ptr, tgtstr1d, shape=[ilen + 1, isize])
! allocate work array to handle CharacterStringType,
! and copy the strings
allocate (character(ilen) :: tempstr)
do i = 1, isize
tempstr = srccharstr1d(i)
tgtstr1d(1:ilen + 1, i) = string_to_char_array(tempstr, ilen)
end do
deallocate (tempstr)
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_string
!> @brief Copy the value of a variable into the array
!!
!! The copied variable is located at @p c_var_address. The caller should
!! provide @p c_arr_ptr pointing to an array of the proper shape (the
!! BMI function get_var_shape() can be used to create it). Multi-dimensional
!! arrays are supported.
!<
function get_value(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value
! -- modules
use ConstantsModule, only: LENMEMTYPE
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENMEMTYPE) :: mem_type
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
call get_mem_type(var_name, mem_path, mem_type)
if (index(mem_type, "DOUBLE") /= 0) then
bmi_status = get_value_double(c_var_address, c_arr_ptr)
else if (index(mem_type, "INTEGER") /= 0) then
bmi_status = get_value_int(c_var_address, c_arr_ptr)
else if (index(mem_type, "LOGICAL") /= 0) then
bmi_status = get_value_bool(c_var_address, c_arr_ptr)
else if (index(mem_type, "STRING") /= 0) then
bmi_status = get_value_string(c_var_address, c_arr_ptr)
else
write (bmi_last_error, fmt_unsupported_type) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value
!> @brief Get a pointer to an array
!!
!! The array is located at @p c_var_address. There is no copying of data involved.
!! Multi-dimensional arrays are supported and the get_var_rank() function
!! can be used to get the variable's dimensionality, and get_var_shape() for
!! its shape.
!<
function get_value_ptr(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_ptr")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr
! -- modules
use ConstantsModule, only: LENMEMTYPE
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENMEMTYPE) :: mem_type
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
call get_mem_type(var_name, mem_path, mem_type)
if (index(mem_type, "DOUBLE") /= 0) then
bmi_status = get_value_ptr_double(c_var_address, c_arr_ptr)
else if (index(mem_type, "INTEGER") /= 0) then
bmi_status = get_value_ptr_int(c_var_address, c_arr_ptr)
else if (index(mem_type, "LOGICAL") /= 0) then
bmi_status = get_value_ptr_bool(c_var_address, c_arr_ptr)
else
write (bmi_last_error, fmt_unsupported_type) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_ptr
!> @brief Get a pointer to the array of double precision numbers
!!
!! The array is located at @p c_var_address. There is no copying of data involved.
!! Multi-dimensional arrays are supported and the get_var_rank() function
!! can be used to get the variable's dimensionality, and get_var_shape() for
!! its shape.
!<
function get_value_ptr_double(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_ptr_double")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr_double
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
real(DP), pointer :: scalar_ptr
real(DP), dimension(:), pointer, contiguous :: array_ptr
real(DP), dimension(:, :), pointer, contiguous :: array2D_ptr
real(DP), dimension(:, :, :), pointer, contiguous :: array3D_ptr
integer(I4B) :: rank
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
rank = -1
call get_mem_rank(var_name, mem_path, rank)
if (rank == 0) then
call mem_setptr(scalar_ptr, var_name, mem_path)
c_arr_ptr = c_loc(scalar_ptr)
else if (rank == 1) then
call mem_setptr(array_ptr, var_name, mem_path)
c_arr_ptr = c_loc(array_ptr)
else if (rank == 2) then
call mem_setptr(array2D_ptr, var_name, mem_path)
c_arr_ptr = c_loc(array2D_ptr)
else if (rank == 3) then
call mem_setptr(array3D_ptr, var_name, mem_path)
c_arr_ptr = c_loc(array3D_ptr)
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_ptr_double
!> @brief Get a pointer to the array of integer numbers
!!
!! The array is located at @p c_var_address. There is no copying of data involved.
!! Multi-dimensional arrays are supported and the get_var_rank() function
!! can be used to get the variable's dimensionality.
!<
function get_value_ptr_int(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_ptr_int")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr_int
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
integer(I4B) :: rank
integer(I4B), pointer :: scalar_ptr
integer(I4B), dimension(:), pointer, contiguous :: array_ptr
integer(I4B), dimension(:, :), pointer, contiguous :: array2D_ptr
integer(I4B), dimension(:, :, :), pointer, contiguous :: array3D_ptr
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
rank = -1
call get_mem_rank(var_name, mem_path, rank)
if (rank == 0) then
call mem_setptr(scalar_ptr, var_name, mem_path)
c_arr_ptr = c_loc(scalar_ptr)
else if (rank == 1) then
call mem_setptr(array_ptr, var_name, mem_path)
c_arr_ptr = c_loc(array_ptr)
else if (rank == 2) then
call mem_setptr(array2D_ptr, var_name, mem_path)
c_arr_ptr = c_loc(array2D_ptr)
else if (rank == 3) then
call mem_setptr(array3D_ptr, var_name, mem_path)
c_arr_ptr = c_loc(array3D_ptr)
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_ptr_int
!> @brief Get a pointer to the logical scalar value
!!
!! Only scalar values (with rank=0) are supported.
!<
function get_value_ptr_bool(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="get_value_ptr_bool")
!DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr_bool
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
logical(LGP), pointer :: scalar_ptr
integer(I4B) :: rank
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
rank = -1
call get_mem_rank(var_name, mem_path, rank)
if (rank == 0) then
call mem_setptr(scalar_ptr, var_name, mem_path)
c_arr_ptr = c_loc(scalar_ptr)
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function get_value_ptr_bool
!> @brief Set new values for a given variable
!!
!! The array pointed to by @p c_arr_ptr can have rank equal to 0, 1, or 2
!! and should have a C-style layout, which is particularly important for
!! rank > 1.
!<
function set_value(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="set_value")
!DIR$ ATTRIBUTES DLLEXPORT :: set_value
! -- modules
use ConstantsModule, only: LENMEMTYPE
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENMEMTYPE) :: mem_type
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
call get_mem_type(var_name, mem_path, mem_type)
if (index(mem_type, "DOUBLE") /= 0) then
bmi_status = set_value_double(c_var_address, c_arr_ptr)
else if (index(mem_type, "INTEGER") /= 0) then
bmi_status = set_value_int(c_var_address, c_arr_ptr)
else if (index(mem_type, "LOGICAL") /= 0) then
bmi_status = set_value_bool(c_var_address, c_arr_ptr)
else
write (bmi_last_error, fmt_unsupported_type) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
end function set_value
!> @brief Set new values for a variable of type double
!!
!! The array pointed to by @p c_arr_ptr can have rank equal to 0, 1, or 2
!! and should have a C-style layout, which is particularly important for
!! rank > 1.
!<
function set_value_double(c_var_address, c_arr_ptr) result(bmi_status) &
bind(C, name="set_value_double")
!DIR$ ATTRIBUTES DLLEXPORT :: set_value_double
! -- modules
use MemorySetHandlerModule, only: on_memory_set
! -- dummy variables
character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the double precision array
integer(kind=c_int) :: bmi_status !< BMI status code
! -- local variables
character(len=LENMEMPATH) :: mem_path
character(len=LENVARNAME) :: var_name
logical(LGP) :: valid
integer(I4B) :: rank
real(DP), pointer :: src_ptr, tgt_ptr
real(DP), dimension(:), pointer, contiguous :: src1D_ptr, tgt1D_ptr
real(DP), dimension(:, :), pointer, contiguous :: src2D_ptr, tgt2D_ptr
integer(I4B) :: i, j
integer(I4B) :: status
bmi_status = BMI_SUCCESS
call split_address(c_var_address, mem_path, var_name, valid)
if (.not. valid) then
bmi_status = BMI_FAILURE
return
end if
! convert pointer and copy, using loops to avoid stack overflow
rank = -1
call get_mem_rank(var_name, mem_path, rank)
if (rank == 0) then
call mem_setptr(tgt_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, src_ptr)
tgt_ptr = src_ptr
else if (rank == 1) then
call mem_setptr(tgt1D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, src1D_ptr, shape(tgt1D_ptr))
do i = 1, size(tgt1D_ptr)
tgt1D_ptr(i) = src1D_ptr(i)
end do
else if (rank == 2) then
call mem_setptr(tgt2D_ptr, var_name, mem_path)
call c_f_pointer(c_arr_ptr, src2D_ptr, shape(tgt2D_ptr))
do j = 1, size(tgt2D_ptr, 2)
do i = 1, size(tgt2D_ptr, 1)
tgt2D_ptr(i, j) = src2D_ptr(i, j)
end do
end do
else
write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
call report_bmi_error(bmi_last_error)
bmi_status = BMI_FAILURE
return
end if
! trigger event:
call on_memory_set(var_name, mem_path, status)
if (status /= 0) then
! something went terribly wrong here, aborting