MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
mf6bmi.f90
Go to the documentation of this file.
1 !> @brief This module contains the MODFLOW 6 BMI
2 !!
3 !! This BMI interface matches the CSDMS standard, with a few modifications:
4 !!
5 !! - This interface will build into a shared library that can be called from other
6 !! executables and scripts, not necessarily written in Fortran. Therefore we have
7 !! omitted the type-boundness of the routines, since they cannot have the
8 !! bind(C,"...") attribute.
9 !! - MODFLOW has internal data arrays with rank > 1 that we would like to expose. An
10 !! example would be access to data in the BOUND array of GWF boundary packages (BndType).
11 !! The get_value_ptr calls below support this, returning a C-style pointer to the arrays
12 !! and methods have been added to query the variable's rank and shape.
13 !!
14 !! Note on style: BMI apparently uses underscores, we use underscores in some
15 !! places but camelcase in other. Since this is a dedicated BMI interface module,
16 !! we'll use underscores here as well.
17 !<
18 module mf6bmi
19  use mf6bmiutil
20  use mf6bmierror
21  use mf6coremodule
22  use tdismodule, only: kper, kstp
23  use iso_c_binding, only: c_int, c_char, c_double, c_null_char, c_loc, c_ptr, &
24  c_f_pointer
25  use kindmodule, only: dp, i4b, lgp
32  use memorytypemodule, only: memorytype
35  use inputoutputmodule, only: getunit
36  implicit none
37 
38  integer(c_int), bind(C, name="ISTDOUTTOFILE") :: istdout_to_file = 1 !< output control: =0 to screen, >0 to file
39  !DIR$ ATTRIBUTES DLLEXPORT :: istdout_to_file
40 
41 contains
42 
43  function bmi_get_component_name(name) result(bmi_status) &
44  bind(C, name="get_component_name")
45  !DIR$ ATTRIBUTES DLLEXPORT :: bmi_get_component_name
46  ! -- dummy variables
47  character(kind=c_char), intent(inout) :: name(bmi_lencomponentname)
48  integer(kind=c_int) :: bmi_status !< BMI status code
49  ! -- local variables
50 
51  name = string_to_char_array('MODFLOW 6', 9)
52  bmi_status = bmi_success
53 
54  end function bmi_get_component_name
55 
56  !> @brief Initialize the computational core
57  !!
58  !! It is required to have the MODFLOW 6 configuration file 'mfsim.nam'
59  !! available in the current working directory when calling this routine.
60  !!
61  !! NOTE: initialization should always be matched with a call to
62  !! bmi_finalize(). However, we currently do not support the reinitialization
63  !! of a model in the same memory space... You would have to create a new
64  !! process for that.
65  !<
66  function bmi_initialize() result(bmi_status) bind(C, name="initialize")
67  !DIR$ ATTRIBUTES DLLEXPORT :: bmi_initialize
68  ! -- dummy variables
69  integer(kind=c_int) :: bmi_status !< BMI status code
70  ! -- local variables
71 
72  if (istdout_to_file > 0) then
73  ! -- open stdout file mfsim.stdout
74  istdout = getunit()
75  !
76  ! -- set STDOUT to a physical file unit
77  open (unit=istdout, file=simstdout)
78  end if
79  !
80  ! -- initialize MODFLOW 6
81  call mf6initialize()
82 
83  bmi_status = bmi_success
84 
85  end function bmi_initialize
86 
87  !> @brief Perform a computational time step
88  !!
89  !! It will prepare the timestep, call the calculation routine
90  !! on all the solution groups in the simulation, and finalize
91  !! the timestep by printing out diagnostics and writing output.
92  !! It can be called in succession to perform multiple steps up
93  !! to the simulation's end time is reached.
94  !<
95  function bmi_update() result(bmi_status) bind(C, name="update")
96  !DIR$ ATTRIBUTES DLLEXPORT :: bmi_update
97  ! -- dummy variables
98  integer(kind=c_int) :: bmi_status !< BMI status code
99  ! -- local variables
100  logical :: hasconverged
101 
102  hasconverged = mf6update()
103 
104  bmi_status = bmi_success
105  end function bmi_update
106 
107  !> @brief Clean up the initialized simulation
108  !!
109  !! Performs teardown tasks for the initialized simulation, this
110  !! call should match the call to bmi_initialize()
111  !<
112  function bmi_finalize() result(bmi_status) bind(C, name="finalize")
113  !DIR$ ATTRIBUTES DLLEXPORT :: bmi_finalize
114  ! -- modules
115  use simvariablesmodule, only: iforcestop
116  ! -- dummy variables
117  integer(kind=c_int) :: bmi_status !< BMI status code
118 
119  ! we don't want a full stop() here, this disables it:
120  iforcestop = 0
121  call mf6finalize()
122 
123  bmi_status = bmi_success
124 
125  end function bmi_finalize
126 
127  !> @brief Get the start time of the simulation
128  !!
129  !! As MODFLOW currently does not have internal time, this will be
130  !! returning 0.0 for now. New version...
131  !<
132  function get_start_time(start_time) result(bmi_status) &
133  bind(C, name="get_start_time")
134  !DIR$ ATTRIBUTES DLLEXPORT :: get_start_time
135  ! -- dummy variables
136  real(kind=c_double), intent(out) :: start_time !< start time
137  integer(kind=c_int) :: bmi_status !< BMI status code
138 
139  start_time = 0.0_dp
140  bmi_status = bmi_success
141 
142  end function get_start_time
143 
144  !> @brief Get the end time of the simulation
145  !!
146  !! As MODFLOW does currently does not have internal time, this will be
147  !! equal to the total runtime.
148  !<
149  function get_end_time(end_time) result(bmi_status) bind(C, name="get_end_time")
150  !DIR$ ATTRIBUTES DLLEXPORT :: get_end_time
151  ! -- modules
152  use tdismodule, only: totalsimtime
153  ! -- dummy variables
154  real(kind=c_double), intent(out) :: end_time !< end time
155  integer(kind=c_int) :: bmi_status !< BMI status code
156 
157  end_time = totalsimtime
158  bmi_status = bmi_success
159 
160  end function get_end_time
161 
162  !> @brief Get the current time of the simulation
163  !!
164  !! As MODFLOW currently does not have internal time, this will be
165  !! equal to the time passed w.r.t. the start time of the simulation.
166  !<
167  function get_current_time(current_time) result(bmi_status) &
168  bind(C, name="get_current_time")
169  !DIR$ ATTRIBUTES DLLEXPORT :: get_current_time
170  ! -- modules
171  use tdismodule, only: totim
172  ! -- dummy variables
173  real(kind=c_double), intent(out) :: current_time !< current time
174  integer(kind=c_int) :: bmi_status !< BMI status code
175 
176  current_time = totim
177  bmi_status = bmi_success
178 
179  end function get_current_time
180 
181  !> @brief Get the time step for the simulation
182  !!
183  !! Note that the returned value may vary between and within stress periods,
184  !! depending on your time discretization settings in the TDIS package.
185  !<
186  function get_time_step(time_step) result(bmi_status) &
187  bind(C, name="get_time_step")
188  !DIR$ ATTRIBUTES DLLEXPORT :: get_time_step
189  ! -- modules
190  use tdismodule, only: delt
191  ! -- dummy variables
192  real(kind=c_double), intent(out) :: time_step !< current time step
193  integer(kind=c_int) :: bmi_status !< BMI status code
194 
195  time_step = delt
196  bmi_status = bmi_success
197 
198  end function get_time_step
199 
200  !> @brief Get the number of input variables in the simulation
201  !!
202  !! This concerns all variables stored in the memory manager
203  !<
204  function get_input_item_count(count) result(bmi_status) &
205  bind(C, name="get_input_item_count")
206  !DIR$ ATTRIBUTES DLLEXPORT :: get_input_item_count
207  ! -- dummy variables
208  integer(kind=c_int), intent(out) :: count !< the number of input variables
209  integer(kind=c_int) :: bmi_status !< BMI status code
210 
211  count = memorystore%count()
212 
213  bmi_status = bmi_success
214 
215  end function get_input_item_count
216 
217  !> @brief Get the number of output variables in the simulation
218  !!
219  !! This concerns all variables stored in the memory manager
220  !<
221  function get_output_item_count(count) result(bmi_status) &
222  bind(C, name="get_output_item_count")
223  !DIR$ ATTRIBUTES DLLEXPORT :: get_output_item_count
224  ! -- dummy variables
225  integer(kind=c_int), intent(out) :: count !< the number of output variables
226  integer(kind=c_int) :: bmi_status !< BMI status code
227 
228  count = memorystore%count()
229 
230  bmi_status = bmi_success
231 
232  end function get_output_item_count
233 
234  !> @brief Returns all input variables in the simulation
235  !!
236  !! This functions returns the full address for all variables in the
237  !! memory manager
238  !!
239  !! The array @p c_names should be pre-allocated of proper size:
240  !!
241  !! size = BMI_LENVARADDRESS * get_input_item_count()
242  !!
243  !! The strings will be written contiguously with stride equal to
244  !! BMI_LENVARADDRESS and nul-terminated where the trailing spaces start:
245  !!
246  !! c_names = 'variable_address_1\x00 ... variable_address_2\x00 ... ' etc.
247  !<
248  function get_input_var_names(c_names) result(bmi_status) &
249  bind(C, name="get_input_var_names")
250  !DIR$ ATTRIBUTES DLLEXPORT :: get_input_var_names
251  ! -- dummy variables
252  character(kind=c_char, len=1), intent(inout) :: c_names(*) !< array with memory paths for input variables
253  integer(kind=c_int) :: bmi_status !< BMI status code
254  ! -- local variables
255  integer(I4B) :: start, i
256  type(memorycontaineriteratortype), allocatable :: itr
257  type(memorytype), pointer :: mt => null()
258  character(len=LENMEMADDRESS) :: var_address
259 
260  start = 1
261  itr = memorystore%iterator()
262  do while (itr%has_next())
263  call itr%next()
264  mt => itr%value()
265  var_address = create_mem_address(mt%path, mt%name)
266  do i = 1, len(trim(var_address))
267  c_names(start + i - 1) = var_address(i:i)
268  end do
269  c_names(start + i) = c_null_char
270  start = start + bmi_lenvaraddress
271  end do
272 
273  bmi_status = bmi_success
274 
275  end function get_input_var_names
276 
277  !> @brief Returns all output variables in the simulation
278  !!
279  !! This function works analogously to get_input_var_names(),
280  !! and currently returns the same set of memory variables,
281  !! which is all of them!
282  !<
283  function get_output_var_names(c_names) result(bmi_status) &
284  bind(C, name="get_output_var_names")
285  !DIR$ ATTRIBUTES DLLEXPORT :: get_output_var_names
286  ! -- dummy variables
287  character(kind=c_char, len=1), intent(inout) :: c_names(*) !< array with memory paths for output variables
288  integer(kind=c_int) :: bmi_status !< BMI status code
289  ! -- local variables
290  integer(I4B) :: start, i
291  type(memorycontaineriteratortype), allocatable :: itr
292  type(memorytype), pointer :: mt => null()
293  character(len=LENMEMADDRESS) :: var_address
294 
295  start = 1
296  itr = memorystore%iterator()
297  do while (itr%has_next())
298  call itr%next()
299  mt => itr%value()
300  var_address = create_mem_address(mt%path, mt%name)
301  do i = 1, len(trim(var_address))
302  c_names(start + i - 1) = var_address(i:i)
303  end do
304  c_names(start + i) = c_null_char
305  start = start + bmi_lenvaraddress
306  end do
307 
308  bmi_status = bmi_success
309 
310  end function get_output_var_names
311 
312  !> @brief Get the size (in bytes) of a single element of a variable
313  !<
314  function get_var_itemsize(c_var_address, var_size) result(bmi_status) &
315  bind(C, name="get_var_itemsize")
316  !DIR$ ATTRIBUTES DLLEXPORT :: get_var_itemsize
317  ! -- dummy variables
318  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
319  integer(kind=c_int), intent(out) :: var_size !< size of the element in bytes
320  integer(kind=c_int) :: bmi_status !< BMI status code
321  ! -- local variables
322  character(len=LENMEMPATH) :: mem_path
323  character(len=LENVARNAME) :: var_name
324  logical(LGP) :: valid
325 
326  bmi_status = bmi_success
327 
328  call split_address(c_var_address, mem_path, var_name, valid)
329  if (.not. valid) then
330  bmi_status = bmi_failure
331  return
332  end if
333 
334  call get_mem_elem_size(var_name, mem_path, var_size)
335  if (var_size == -1) bmi_status = bmi_failure
336 
337  end function get_var_itemsize
338 
339  !> @brief Get size of the variable, in bytes
340  !<
341  function get_var_nbytes(c_var_address, var_nbytes) result(bmi_status) &
342  bind(C, name="get_var_nbytes")
343  !DIR$ ATTRIBUTES DLLEXPORT :: get_var_nbytes
344  ! -- dummy variables
345  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
346  integer(kind=c_int), intent(out) :: var_nbytes !< size in bytes
347  integer(kind=c_int) :: bmi_status !< BMI status code
348  ! -- local variables
349  integer(I4B) :: var_size, isize
350  character(len=LENMEMPATH) :: mem_path
351  character(len=LENVARNAME) :: var_name
352  logical(LGP) :: valid
353 
354  bmi_status = bmi_success
355 
356  call split_address(c_var_address, mem_path, var_name, valid)
357  if (.not. valid) then
358  bmi_status = bmi_failure
359  return
360  end if
361 
362  call get_mem_elem_size(var_name, mem_path, var_size)
363  if (var_size == -1) bmi_status = bmi_failure
364  call get_isize(var_name, mem_path, isize)
365  if (isize == -1) bmi_status = bmi_failure
366 
367  var_nbytes = var_size * isize
368 
369  end function get_var_nbytes
370 
371  !> @brief Copy the double precision values of a variable into the array
372  !!
373  !! The copied variable us located at @p c_var_address. The caller should
374  !! provide @p c_arr_ptr pointing to an array of the proper shape (the
375  !! BMI function get_var_shape() can be used to create it). Multi-dimensional
376  !! arrays are supported.
377  !<
378  function get_value_double(c_var_address, c_arr_ptr) result(bmi_status) &
379  bind(C, name="get_value_double")
380  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_double
381  ! -- modules
383  ! -- dummy variables
384  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
385  type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the double precision array
386  integer(kind=c_int) :: bmi_status !< BMI status code
387  ! -- local variables
388  character(len=LENMEMPATH) :: mem_path
389  character(len=LENVARNAME) :: var_name
390  logical(LGP) :: valid
391  integer(I4B) :: rank
392  real(dp), pointer :: src_ptr, tgt_ptr
393  real(dp), dimension(:), pointer, contiguous :: src1d_ptr, tgt1d_ptr
394  real(dp), dimension(:, :), pointer, contiguous :: src2d_ptr, tgt2d_ptr
395  real(dp), dimension(:, :, :), pointer, contiguous :: src3d_ptr, tgt3d_ptr
396  integer(I4B) :: i, j, k
397 
398  bmi_status = bmi_success
399 
400  call split_address(c_var_address, mem_path, var_name, valid)
401  if (.not. valid) then
402  bmi_status = bmi_failure
403  return
404  end if
405 
406  ! convert pointer and copy data from memory manager into
407  ! the passed array, using loops to avoid stack overflow
408  rank = -1
409  call get_mem_rank(var_name, mem_path, rank)
410 
411  if (rank == 0) then
412  call mem_setptr(src_ptr, var_name, mem_path)
413  call c_f_pointer(c_arr_ptr, tgt_ptr)
414  tgt_ptr = src_ptr
415  else if (rank == 1) then
416  call mem_setptr(src1d_ptr, var_name, mem_path)
417  call c_f_pointer(c_arr_ptr, tgt1d_ptr, shape(src1d_ptr))
418  do i = 1, size(tgt1d_ptr)
419  tgt1d_ptr(i) = src1d_ptr(i)
420  end do
421  else if (rank == 2) then
422  call mem_setptr(src2d_ptr, var_name, mem_path)
423  call c_f_pointer(c_arr_ptr, tgt2d_ptr, shape(src2d_ptr))
424  do j = 1, size(tgt2d_ptr, 2)
425  do i = 1, size(tgt2d_ptr, 1)
426  tgt2d_ptr(i, j) = src2d_ptr(i, j)
427  end do
428  end do
429  else if (rank == 3) then
430  call mem_setptr(src3d_ptr, var_name, mem_path)
431  call c_f_pointer(c_arr_ptr, tgt3d_ptr, shape(src3d_ptr))
432  do k = 1, size(tgt3d_ptr, 3)
433  do j = 1, size(tgt3d_ptr, 2)
434  do i = 1, size(tgt3d_ptr, 1)
435  tgt3d_ptr(i, j, k) = src3d_ptr(i, j, k)
436  end do
437  end do
438  end do
439  else
440  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
442  bmi_status = bmi_failure
443  return
444  end if
445 
446  end function get_value_double
447 
448  !> @brief Copy the integer values of a variable into the array
449  !!
450  !! The copied variable is located at @p c_var_address. The caller should
451  !! provide @p c_arr_ptr pointing to an array of the proper shape (the
452  !! BMI function get_var_shape() can be used to create it). Multi-dimensional
453  !! arrays are supported.
454  !<
455  function get_value_int(c_var_address, c_arr_ptr) result(bmi_status) &
456  bind(C, name="get_value_int")
457  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_int
458  ! -- modules
460  ! -- dummy variables
461  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
462  type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the integer array
463  integer(kind=c_int) :: bmi_status !< BMI status code
464  ! -- local variables
465  character(len=LENMEMPATH) :: mem_path
466  character(len=LENVARNAME) :: var_name
467  logical(LGP) :: valid
468  integer(I4B) :: rank
469  integer(I4B), pointer :: src_ptr, tgt_ptr
470  integer(I4B), dimension(:), pointer, contiguous :: src1d_ptr, tgt1d_ptr
471  integer(I4B), dimension(:, :), pointer, contiguous :: src2d_ptr, tgt2d_ptr
472  integer(I4B), dimension(:, :, :), pointer, contiguous :: src3d_ptr, tgt3d_ptr
473  integer(I4B) :: i, j, k
474 
475  bmi_status = bmi_success
476 
477  call split_address(c_var_address, mem_path, var_name, valid)
478  if (.not. valid) then
479  bmi_status = bmi_failure
480  return
481  end if
482 
483  ! convert pointer and copy data from memory manager into
484  ! the passed array, using loops to avoid stack overflow
485  rank = -1
486  call get_mem_rank(var_name, mem_path, rank)
487 
488  if (rank == 0) then
489  call mem_setptr(src_ptr, var_name, mem_path)
490  call c_f_pointer(c_arr_ptr, tgt_ptr)
491  tgt_ptr = src_ptr
492  else if (rank == 1) then
493  call mem_setptr(src1d_ptr, var_name, mem_path)
494  call c_f_pointer(c_arr_ptr, tgt1d_ptr, shape(src1d_ptr))
495  do i = 1, size(tgt1d_ptr)
496  tgt1d_ptr(i) = src1d_ptr(i)
497  end do
498  else if (rank == 2) then
499  call mem_setptr(src2d_ptr, var_name, mem_path)
500  call c_f_pointer(c_arr_ptr, tgt2d_ptr, shape(src2d_ptr))
501  do j = 1, size(tgt2d_ptr, 2)
502  do i = 1, size(tgt2d_ptr, 1)
503  tgt2d_ptr(i, j) = src2d_ptr(i, j)
504  end do
505  end do
506  else if (rank == 3) then
507  call mem_setptr(src3d_ptr, var_name, mem_path)
508  call c_f_pointer(c_arr_ptr, tgt3d_ptr, shape(src3d_ptr))
509  do k = 1, size(tgt3d_ptr, 3)
510  do j = 1, size(tgt3d_ptr, 2)
511  do i = 1, size(tgt3d_ptr, 1)
512  tgt3d_ptr(i, j, k) = src3d_ptr(i, j, k)
513  end do
514  end do
515  end do
516  else
517  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
519  bmi_status = bmi_failure
520  return
521  end if
522 
523  end function get_value_int
524 
525  !> @brief Copy the logical scalar value into the array
526  !!
527  !! The copied variable us located at @p c_var_address. The caller should
528  !! provide @p c_arr_ptr pointing to a scalar array with rank=0.
529  !<
530  function get_value_bool(c_var_address, c_arr_ptr) result(bmi_status) &
531  bind(C, name="get_value_bool")
532  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_bool
533  ! -- modules
535  ! -- dummy variables
536  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
537  type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the logical array
538  integer(kind=c_int) :: bmi_status !< BMI status code
539  ! -- local variables
540  character(len=LENMEMPATH) :: mem_path
541  character(len=LENVARNAME) :: var_name
542  logical(LGP) :: valid
543  integer(I4B) :: rank
544  logical(LGP), pointer :: src_ptr, tgt_ptr
545 
546  bmi_status = bmi_success
547 
548  call split_address(c_var_address, mem_path, var_name, valid)
549  if (.not. valid) then
550  bmi_status = bmi_failure
551  return
552  end if
553 
554  rank = -1
555  call get_mem_rank(var_name, mem_path, rank)
556 
557  ! convert pointer
558  if (rank == 0) then
559  call mem_setptr(src_ptr, var_name, mem_path)
560  call c_f_pointer(c_arr_ptr, tgt_ptr)
561  tgt_ptr = src_ptr
562  else
563  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
565  bmi_status = bmi_failure
566  return
567  end if
568 
569  end function get_value_bool
570 
571  !> @brief Copy the string(s) of a variable into the array
572  !!
573  !! The copied variable is located at @p c_var_address. The caller should
574  !! provide @p c_arr_ptr pointing to an array of the proper shape (the
575  !! BMI function get_var_shape() can be used to create it). For strings
576  !! currently scalars and 1d arrays (of CharacterStringType) are
577  !< supported
578  function get_value_string(c_var_address, c_arr_ptr) result(bmi_status) &
579  bind(C, name="get_value_string")
580  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_string
581  ! -- modules
582  ! -- dummy variables
583  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
584  type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the string array
585  integer(kind=c_int) :: bmi_status !< BMI status code
586  ! -- local variables
587  character(len=LENMEMPATH) :: mem_path
588  character(len=LENVARNAME) :: var_name
589  logical(LGP) :: valid
590  integer(I4B) :: rank
591  character(len=:), pointer :: srcstr
592  character(kind=c_char), pointer :: tgtstr(:)
593  type(characterstringtype), dimension(:), pointer, contiguous :: srccharstr1d
594  character(kind=c_char), pointer :: tgtstr1d(:, :)
595  character(:), allocatable :: tempstr
596  integer(I4B) :: i, ilen, isize
597 
598  bmi_status = bmi_success
599 
600  call split_address(c_var_address, mem_path, var_name, valid)
601  if (.not. valid) then
602  bmi_status = bmi_failure
603  return
604  end if
605 
606  ! single string, or array of strings (CharacterStringType)
607  rank = -1
608  call get_mem_rank(var_name, mem_path, rank)
609 
610  if (rank == 0) then
611  ! a string scalar
612  call mem_setptr(srcstr, var_name, mem_path)
613  call get_mem_elem_size(var_name, mem_path, ilen)
614  call c_f_pointer(c_arr_ptr, tgtstr, shape=[ilen + 1])
615 
616  tgtstr(1:len(srcstr) + 1) = string_to_char_array(srcstr, len(srcstr))
617 
618  else if (rank == 1) then
619  ! an array of strings
620  call mem_setptr(srccharstr1d, var_name, mem_path)
621  if (.not. associated(srccharstr1d)) then
622  write (bmi_last_error, fmt_general_err) 'string type not supported in API'
624  bmi_status = bmi_failure
625  return
626  end if
627 
628  ! create fortran pointer to C data array
629  call get_isize(var_name, mem_path, isize)
630  call get_mem_elem_size(var_name, mem_path, ilen)
631  call c_f_pointer(c_arr_ptr, tgtstr1d, shape=[ilen + 1, isize])
632 
633  ! allocate work array to handle CharacterStringType,
634  ! and copy the strings
635  allocate (character(ilen) :: tempstr)
636  do i = 1, isize
637  tempstr = srccharstr1d(i)
638  tgtstr1d(1:ilen + 1, i) = string_to_char_array(tempstr, ilen)
639  end do
640  deallocate (tempstr)
641  else
642  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
644  bmi_status = bmi_failure
645  return
646  end if
647 
648  end function get_value_string
649 
650  !> @brief Copy the value of a variable into the array
651  !!
652  !! The copied variable is located at @p c_var_address. The caller should
653  !! provide @p c_arr_ptr pointing to an array of the proper shape (the
654  !! BMI function get_var_shape() can be used to create it). Multi-dimensional
655  !! arrays are supported.
656  !<
657  function get_value(c_var_address, c_arr_ptr) result(bmi_status) &
658  bind(C, name="get_value")
659  !DIR$ ATTRIBUTES DLLEXPORT :: get_value
660  ! -- modules
661  use constantsmodule, only: lenmemtype
662  ! -- dummy variables
663  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
664  type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
665  integer(kind=c_int) :: bmi_status !< BMI status code
666  ! -- local variables
667  character(len=LENMEMPATH) :: mem_path
668  character(len=LENMEMTYPE) :: mem_type
669  character(len=LENVARNAME) :: var_name
670  logical(LGP) :: valid
671 
672  bmi_status = bmi_success
673 
674  call split_address(c_var_address, mem_path, var_name, valid)
675  if (.not. valid) then
676  bmi_status = bmi_failure
677  return
678  end if
679 
680  call get_mem_type(var_name, mem_path, mem_type)
681 
682  if (index(mem_type, "DOUBLE") /= 0) then
683  bmi_status = get_value_double(c_var_address, c_arr_ptr)
684  else if (index(mem_type, "INTEGER") /= 0) then
685  bmi_status = get_value_int(c_var_address, c_arr_ptr)
686  else if (index(mem_type, "LOGICAL") /= 0) then
687  bmi_status = get_value_bool(c_var_address, c_arr_ptr)
688  else if (index(mem_type, "STRING") /= 0) then
689  bmi_status = get_value_string(c_var_address, c_arr_ptr)
690  else
691  write (bmi_last_error, fmt_unsupported_type) trim(var_name)
693  bmi_status = bmi_failure
694  return
695  end if
696 
697  end function get_value
698 
699  !> @brief Get a pointer to an array
700  !!
701  !! The array is located at @p c_var_address. There is no copying of data involved.
702  !! Multi-dimensional arrays are supported and the get_var_rank() function
703  !! can be used to get the variable's dimensionality, and get_var_shape() for
704  !! its shape.
705  !<
706  function get_value_ptr(c_var_address, c_arr_ptr) result(bmi_status) &
707  bind(C, name="get_value_ptr")
708  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr
709  ! -- modules
710  use constantsmodule, only: lenmemtype
711  ! -- dummy variables
712  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
713  type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
714  integer(kind=c_int) :: bmi_status !< BMI status code
715  ! -- local variables
716  character(len=LENMEMPATH) :: mem_path
717  character(len=LENMEMTYPE) :: mem_type
718  character(len=LENVARNAME) :: var_name
719  logical(LGP) :: valid
720 
721  bmi_status = bmi_success
722 
723  call split_address(c_var_address, mem_path, var_name, valid)
724  if (.not. valid) then
725  bmi_status = bmi_failure
726  return
727  end if
728 
729  call get_mem_type(var_name, mem_path, mem_type)
730 
731  if (index(mem_type, "DOUBLE") /= 0) then
732  bmi_status = get_value_ptr_double(c_var_address, c_arr_ptr)
733  else if (index(mem_type, "INTEGER") /= 0) then
734  bmi_status = get_value_ptr_int(c_var_address, c_arr_ptr)
735  else if (index(mem_type, "LOGICAL") /= 0) then
736  bmi_status = get_value_ptr_bool(c_var_address, c_arr_ptr)
737  else
738  write (bmi_last_error, fmt_unsupported_type) trim(var_name)
740  bmi_status = bmi_failure
741  return
742  end if
743 
744  end function get_value_ptr
745 
746  !> @brief Get a pointer to the array of double precision numbers
747  !!
748  !! The array is located at @p c_var_address. There is no copying of data involved.
749  !! Multi-dimensional arrays are supported and the get_var_rank() function
750  !! can be used to get the variable's dimensionality, and get_var_shape() for
751  !! its shape.
752  !<
753  function get_value_ptr_double(c_var_address, c_arr_ptr) result(bmi_status) &
754  bind(C, name="get_value_ptr_double")
755  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr_double
756  ! -- dummy variables
757  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
758  type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
759  integer(kind=c_int) :: bmi_status !< BMI status code
760  ! -- local variables
761  character(len=LENMEMPATH) :: mem_path
762  character(len=LENVARNAME) :: var_name
763  logical(LGP) :: valid
764  real(dp), pointer :: scalar_ptr
765  real(dp), dimension(:), pointer, contiguous :: array_ptr
766  real(dp), dimension(:, :), pointer, contiguous :: array2d_ptr
767  real(dp), dimension(:, :, :), pointer, contiguous :: array3d_ptr
768  integer(I4B) :: rank
769 
770  bmi_status = bmi_success
771 
772  call split_address(c_var_address, mem_path, var_name, valid)
773  if (.not. valid) then
774  bmi_status = bmi_failure
775  return
776  end if
777 
778  rank = -1
779  call get_mem_rank(var_name, mem_path, rank)
780  if (rank == 0) then
781  call mem_setptr(scalar_ptr, var_name, mem_path)
782  c_arr_ptr = c_loc(scalar_ptr)
783  else if (rank == 1) then
784  call mem_setptr(array_ptr, var_name, mem_path)
785  c_arr_ptr = c_loc(array_ptr)
786  else if (rank == 2) then
787  call mem_setptr(array2d_ptr, var_name, mem_path)
788  c_arr_ptr = c_loc(array2d_ptr)
789  else if (rank == 3) then
790  call mem_setptr(array3d_ptr, var_name, mem_path)
791  c_arr_ptr = c_loc(array3d_ptr)
792  else
793  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
795  bmi_status = bmi_failure
796  return
797  end if
798 
799  end function get_value_ptr_double
800 
801  !> @brief Get a pointer to the array of integer numbers
802  !!
803  !! The array is located at @p c_var_address. There is no copying of data involved.
804  !! Multi-dimensional arrays are supported and the get_var_rank() function
805  !! can be used to get the variable's dimensionality.
806  !<
807  function get_value_ptr_int(c_var_address, c_arr_ptr) result(bmi_status) &
808  bind(C, name="get_value_ptr_int")
809  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr_int
810  ! -- dummy variables
811  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
812  type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
813  integer(kind=c_int) :: bmi_status !< BMI status code
814  ! -- local variables
815  character(len=LENMEMPATH) :: mem_path
816  character(len=LENVARNAME) :: var_name
817  logical(LGP) :: valid
818  integer(I4B) :: rank
819  integer(I4B), pointer :: scalar_ptr
820  integer(I4B), dimension(:), pointer, contiguous :: array_ptr
821  integer(I4B), dimension(:, :), pointer, contiguous :: array2d_ptr
822  integer(I4B), dimension(:, :, :), pointer, contiguous :: array3d_ptr
823 
824  bmi_status = bmi_success
825 
826  call split_address(c_var_address, mem_path, var_name, valid)
827  if (.not. valid) then
828  bmi_status = bmi_failure
829  return
830  end if
831 
832  rank = -1
833  call get_mem_rank(var_name, mem_path, rank)
834 
835  if (rank == 0) then
836  call mem_setptr(scalar_ptr, var_name, mem_path)
837  c_arr_ptr = c_loc(scalar_ptr)
838  else if (rank == 1) then
839  call mem_setptr(array_ptr, var_name, mem_path)
840  c_arr_ptr = c_loc(array_ptr)
841  else if (rank == 2) then
842  call mem_setptr(array2d_ptr, var_name, mem_path)
843  c_arr_ptr = c_loc(array2d_ptr)
844  else if (rank == 3) then
845  call mem_setptr(array3d_ptr, var_name, mem_path)
846  c_arr_ptr = c_loc(array3d_ptr)
847  else
848  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
850  bmi_status = bmi_failure
851  return
852  end if
853 
854  end function get_value_ptr_int
855 
856  !> @brief Get a pointer to the logical scalar value
857  !!
858  !! Only scalar values (with rank=0) are supported.
859  !<
860  function get_value_ptr_bool(c_var_address, c_arr_ptr) result(bmi_status) &
861  bind(C, name="get_value_ptr_bool")
862  !DIR$ ATTRIBUTES DLLEXPORT :: get_value_ptr_bool
863  ! -- dummy variables
864  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
865  type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
866  integer(kind=c_int) :: bmi_status !< BMI status code
867  ! -- local variables
868  character(len=LENMEMPATH) :: mem_path
869  character(len=LENVARNAME) :: var_name
870  logical(LGP) :: valid
871  logical(LGP), pointer :: scalar_ptr
872  integer(I4B) :: rank
873 
874  bmi_status = bmi_success
875 
876  call split_address(c_var_address, mem_path, var_name, valid)
877  if (.not. valid) then
878  bmi_status = bmi_failure
879  return
880  end if
881 
882  rank = -1
883  call get_mem_rank(var_name, mem_path, rank)
884  if (rank == 0) then
885  call mem_setptr(scalar_ptr, var_name, mem_path)
886  c_arr_ptr = c_loc(scalar_ptr)
887  else
888  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
890  bmi_status = bmi_failure
891  return
892  end if
893 
894  end function get_value_ptr_bool
895 
896  !> @brief Set new values for a given variable
897  !!
898  !! The array pointed to by @p c_arr_ptr can have rank equal to 0, 1, or 2
899  !! and should have a C-style layout, which is particularly important for
900  !! rank > 1.
901  !<
902  function set_value(c_var_address, c_arr_ptr) result(bmi_status) &
903  bind(C, name="set_value")
904  !DIR$ ATTRIBUTES DLLEXPORT :: set_value
905  ! -- modules
906  use constantsmodule, only: lenmemtype
907  ! -- dummy variables
908  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
909  type(c_ptr), intent(inout) :: c_arr_ptr !< pointer to the array
910  integer(kind=c_int) :: bmi_status !< BMI status code
911  ! -- local variables
912  character(len=LENMEMPATH) :: mem_path
913  character(len=LENMEMTYPE) :: mem_type
914  character(len=LENVARNAME) :: var_name
915  logical(LGP) :: valid
916 
917  bmi_status = bmi_success
918 
919  call split_address(c_var_address, mem_path, var_name, valid)
920  if (.not. valid) then
921  bmi_status = bmi_failure
922  return
923  end if
924 
925  call get_mem_type(var_name, mem_path, mem_type)
926 
927  if (index(mem_type, "DOUBLE") /= 0) then
928  bmi_status = set_value_double(c_var_address, c_arr_ptr)
929  else if (index(mem_type, "INTEGER") /= 0) then
930  bmi_status = set_value_int(c_var_address, c_arr_ptr)
931  else if (index(mem_type, "LOGICAL") /= 0) then
932  bmi_status = set_value_bool(c_var_address, c_arr_ptr)
933  else
934  write (bmi_last_error, fmt_unsupported_type) trim(var_name)
936  bmi_status = bmi_failure
937  return
938  end if
939 
940  end function set_value
941 
942  !> @brief Set new values for a variable of type double
943  !!
944  !! The array pointed to by @p c_arr_ptr can have rank equal to 0, 1, or 2
945  !! and should have a C-style layout, which is particularly important for
946  !! rank > 1.
947  !<
948  function set_value_double(c_var_address, c_arr_ptr) result(bmi_status) &
949  bind(C, name="set_value_double")
950  !DIR$ ATTRIBUTES DLLEXPORT :: set_value_double
951  ! -- modules
953  ! -- dummy variables
954  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
955  type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the double precision array
956  integer(kind=c_int) :: bmi_status !< BMI status code
957  ! -- local variables
958  character(len=LENMEMPATH) :: mem_path
959  character(len=LENVARNAME) :: var_name
960  logical(LGP) :: valid
961  integer(I4B) :: rank
962  real(dp), pointer :: src_ptr, tgt_ptr
963  real(dp), dimension(:), pointer, contiguous :: src1d_ptr, tgt1d_ptr
964  real(dp), dimension(:, :), pointer, contiguous :: src2d_ptr, tgt2d_ptr
965  integer(I4B) :: i, j
966  integer(I4B) :: status
967 
968  bmi_status = bmi_success
969 
970  call split_address(c_var_address, mem_path, var_name, valid)
971  if (.not. valid) then
972  bmi_status = bmi_failure
973  return
974  end if
975 
976  ! convert pointer and copy, using loops to avoid stack overflow
977  rank = -1
978  call get_mem_rank(var_name, mem_path, rank)
979 
980  if (rank == 0) then
981  call mem_setptr(tgt_ptr, var_name, mem_path)
982  call c_f_pointer(c_arr_ptr, src_ptr)
983  tgt_ptr = src_ptr
984  else if (rank == 1) then
985  call mem_setptr(tgt1d_ptr, var_name, mem_path)
986  call c_f_pointer(c_arr_ptr, src1d_ptr, shape(tgt1d_ptr))
987  do i = 1, size(tgt1d_ptr)
988  tgt1d_ptr(i) = src1d_ptr(i)
989  end do
990  else if (rank == 2) then
991  call mem_setptr(tgt2d_ptr, var_name, mem_path)
992  call c_f_pointer(c_arr_ptr, src2d_ptr, shape(tgt2d_ptr))
993  do j = 1, size(tgt2d_ptr, 2)
994  do i = 1, size(tgt2d_ptr, 1)
995  tgt2d_ptr(i, j) = src2d_ptr(i, j)
996  end do
997  end do
998  else
999  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
1001  bmi_status = bmi_failure
1002  return
1003  end if
1004 
1005  ! trigger event:
1006  call on_memory_set(var_name, mem_path, status)
1007  if (status /= 0) then
1008  ! something went terribly wrong here, aborting
1009  write (bmi_last_error, fmt_invalid_mem_access) trim(var_name)
1011  bmi_status = bmi_failure
1012  return
1013  end if
1014 
1015  end function set_value_double
1016 
1017  !> @brief Set new values for a variable of type integer
1018  !!
1019  !! The array pointed to by @p c_arr_ptr can have rank equal to 0, 1, or 2.
1020  !<
1021  function set_value_int(c_var_address, c_arr_ptr) result(bmi_status) &
1022  bind(C, name="set_value_int")
1023  !DIR$ ATTRIBUTES DLLEXPORT :: set_value_int
1024  ! -- modules
1026  ! -- dummy variables
1027  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
1028  type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the integer array
1029  integer(kind=c_int) :: bmi_status !< BMI status code
1030  ! -- local variables
1031  character(len=LENMEMPATH) :: mem_path
1032  character(len=LENVARNAME) :: var_name
1033  logical(LGP) :: valid
1034  integer(I4B) :: rank
1035  integer(I4B), pointer :: src_ptr, tgt_ptr
1036  integer(I4B), dimension(:), pointer, contiguous :: src1d_ptr, tgt1d_ptr
1037  integer(I4B), dimension(:, :), pointer, contiguous :: src2d_ptr, tgt2d_ptr
1038  integer(I4B) :: i, j
1039  integer(I4B) :: status
1040 
1041  bmi_status = bmi_success
1042 
1043  call split_address(c_var_address, mem_path, var_name, valid)
1044  if (.not. valid) then
1045  bmi_status = bmi_failure
1046  return
1047  end if
1048 
1049  ! convert pointer and copy, using loops to avoid stack overflow
1050  rank = -1
1051  call get_mem_rank(var_name, mem_path, rank)
1052 
1053  if (rank == 0) then
1054  call mem_setptr(tgt_ptr, var_name, mem_path)
1055  call c_f_pointer(c_arr_ptr, src_ptr)
1056  tgt_ptr = src_ptr
1057  else if (rank == 1) then
1058  call mem_setptr(tgt1d_ptr, var_name, mem_path)
1059  call c_f_pointer(c_arr_ptr, src1d_ptr, shape(tgt1d_ptr))
1060  do i = 1, size(tgt1d_ptr)
1061  tgt1d_ptr(i) = src1d_ptr(i)
1062  end do
1063  else if (rank == 2) then
1064  call mem_setptr(tgt2d_ptr, var_name, mem_path)
1065  call c_f_pointer(c_arr_ptr, src2d_ptr, shape(tgt2d_ptr))
1066  do j = 1, size(tgt2d_ptr, 2)
1067  do i = 1, size(tgt2d_ptr, 1)
1068  tgt2d_ptr(i, j) = src2d_ptr(i, j)
1069  end do
1070  end do
1071  else
1072  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
1074  bmi_status = bmi_failure
1075  return
1076  end if
1077 
1078  ! trigger event:
1079  call on_memory_set(var_name, mem_path, status)
1080  if (status /= 0) then
1081  ! something went terribly wrong here, aborting
1082  write (bmi_last_error, fmt_invalid_mem_access) trim(var_name)
1084  bmi_status = bmi_failure
1085  return
1086  end if
1087 
1088  end function set_value_int
1089 
1090  !> @brief Set new value for a logical scalar variable
1091  !!
1092  !! The array pointed to by @p c_arr_ptr must have a rank equal to 0.
1093  !<
1094  function set_value_bool(c_var_address, c_arr_ptr) result(bmi_status) &
1095  bind(C, name="set_value_bool")
1096  !DIR$ ATTRIBUTES DLLEXPORT :: set_value_bool
1097  ! -- modules
1099  ! -- dummy variables
1100  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
1101  type(c_ptr), intent(in) :: c_arr_ptr !< pointer to the logical array
1102  integer(kind=c_int) :: bmi_status !< BMI status code
1103  ! -- local variables
1104  character(len=LENMEMPATH) :: mem_path
1105  character(len=LENVARNAME) :: var_name
1106  logical(LGP) :: valid
1107  integer(I4B) :: rank
1108  logical(LGP), pointer :: src_ptr, tgt_ptr
1109  integer(I4B) :: status
1110 
1111  bmi_status = bmi_success
1112 
1113  call split_address(c_var_address, mem_path, var_name, valid)
1114  if (.not. valid) then
1115  bmi_status = bmi_failure
1116  return
1117  end if
1118 
1119  rank = -1
1120  call get_mem_rank(var_name, mem_path, rank)
1121 
1122  ! convert pointer
1123  if (rank == 0) then
1124  call mem_setptr(tgt_ptr, var_name, mem_path)
1125  call c_f_pointer(c_arr_ptr, src_ptr)
1126  tgt_ptr = src_ptr
1127  else
1128  write (bmi_last_error, fmt_unsupported_rank) trim(var_name)
1130  bmi_status = bmi_failure
1131  return
1132  end if
1133 
1134  ! trigger event:
1135  call on_memory_set(var_name, mem_path, status)
1136  if (status /= 0) then
1137  ! something went terribly wrong here, aborting
1138  write (bmi_last_error, fmt_invalid_mem_access) trim(var_name)
1140  bmi_status = bmi_failure
1141  return
1142  end if
1143 
1144  end function set_value_bool
1145 
1146  !> @brief Get the variable type as a string
1147  !!
1148  !! The type returned is that of a single element.
1149  !! When the variable cannot be found, the string
1150  !< 'UNKNOWN' is assigned.
1151  function get_var_type(c_var_address, c_var_type) result(bmi_status) &
1152  bind(C, name="get_var_type")
1153  !DIR$ ATTRIBUTES DLLEXPORT :: get_var_type
1154  ! -- modules
1155  use constantsmodule, only: lenmemtype
1156  ! -- dummy variables
1157  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
1158  character(kind=c_char), intent(out) :: c_var_type(bmi_lenvartype) !< variable type as a string
1159  integer(kind=c_int) :: bmi_status !< BMI status code
1160  ! -- local variables
1161  character(len=LENMEMPATH) :: mem_path
1162  character(len=LENVARNAME) :: var_name
1163  character(len=LENMEMTYPE) :: mem_type
1164  logical(LGP) :: valid
1165 
1166  bmi_status = bmi_success
1167 
1168  call split_address(c_var_address, mem_path, var_name, valid)
1169  if (.not. valid) then
1170  bmi_status = bmi_failure
1171  return
1172  end if
1173 
1174  call get_mem_type(var_name, mem_path, mem_type)
1175  c_var_type(1:len(trim(mem_type)) + 1) = &
1176  string_to_char_array(trim(mem_type), len(trim(mem_type)))
1177 
1178  if (mem_type == 'UNKNOWN') then
1179  write (bmi_last_error, fmt_general_err) 'unknown memory type'
1181  bmi_status = bmi_failure
1182  end if
1183 
1184  end function get_var_type
1185 
1186  !> @brief Get the variable rank (non-BMI)
1187  !!
1188  !! In order to support multi-dimensional arrays, this function gives
1189  !! access to the rank of the array.
1190  !<
1191  function get_var_rank(c_var_address, c_var_rank) result(bmi_status) &
1192  bind(C, name="get_var_rank")
1193  !DIR$ ATTRIBUTES DLLEXPORT :: get_var_rank
1194  ! -- dummy variables
1195  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
1196  integer(kind=c_int), intent(out) :: c_var_rank !< variable rank
1197  integer(kind=c_int) :: bmi_status !< BMI status code
1198  ! -- local variables
1199  character(len=LENMEMPATH) :: mem_path
1200  character(len=LENVARNAME) :: var_name
1201  logical(LGP) :: valid
1202 
1203  bmi_status = bmi_success
1204 
1205  call split_address(c_var_address, mem_path, var_name, valid)
1206  if (.not. valid) then
1207  bmi_status = bmi_failure
1208  return
1209  end if
1210 
1211  call get_mem_rank(var_name, mem_path, c_var_rank)
1212  if (c_var_rank == -1) then
1213  bmi_status = bmi_failure
1214  return
1215  end if
1216 
1217  end function get_var_rank
1218 
1219  !> @brief Get the shape of the array for the variable (non-BMI)
1220  !!
1221  !! The shape is an integer array with size equal to the rank of
1222  !! the variable (see get_var_rank()) and values that give the
1223  !! length of the array in each dimension. The target shape array
1224  !! @p c_var_shape should be allocated before calling this routine.
1225  !!
1226  !! Note that the returned shape representation will has been converted
1227  !! to C-style.
1228  !<
1229  function get_var_shape(c_var_address, c_var_shape) result(bmi_status) &
1230  bind(C, name="get_var_shape")
1231  !DIR$ ATTRIBUTES DLLEXPORT :: get_var_shape
1232  ! -- modules
1233  use constantsmodule, only: maxmemrank
1234  ! -- dummy variables
1235  character(kind=c_char), intent(in) :: c_var_address(*) !< memory address string of the variable
1236  integer(c_int), intent(inout) :: c_var_shape(*) !< 1D array with the variable's shape
1237  integer(kind=c_int) :: bmi_status !< BMI status code
1238  ! -- local variables
1239  integer(I4B), dimension(MAXMEMRANK) :: var_shape
1240  integer(I4B) :: var_rank
1241  character(len=LENMEMPATH) :: mem_path
1242  character(len=LENVARNAME) :: var_name
1243  logical(LGP) :: valid
1244 
1245  bmi_status = bmi_success
1246 
1247  call split_address(c_var_address, mem_path, var_name, valid)
1248  if (.not. valid) then
1249  bmi_status = bmi_failure
1250  return
1251  end if
1252 
1253  var_shape = 0
1254  var_rank = 0
1255  call get_mem_rank(var_name, mem_path, var_rank)
1256  call get_mem_shape(var_name, mem_path, var_shape)
1257  if (var_shape(1) == -1 .or. var_rank == -1) then
1258  bmi_status = bmi_failure
1259  return
1260  end if
1261 
1262  ! External calls to this BMI are assumed C style, so if the internal shape
1263  ! is (100,1) we get (100,1,undef) from the call get_mem_shape
1264  ! This we need to convert to C-style which should be (1,100).
1265  ! Hence, we reverse the array and drop undef:
1266  c_var_shape(1:var_rank) = var_shape(var_rank:1:-1)
1267 
1268  end function get_var_shape
1269 
1270 end module mf6bmi
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
integer(i4b), parameter, public maxmemrank
maximum memory manager length (up to 3-dimensional arrays)
Definition: Constants.f90:61
integer(i4b), parameter, public lenmemtype
maximum length of a memory manager type
Definition: Constants.f90:62
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
integer(i4b) function, public getunit()
Get a free unit number.
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmemaddress) function create_mem_address(mem_path, var_name)
returns the address string of the memory object
subroutine, public get_mem_type(name, mem_path, var_type)
@ brief Get the variable memory type
subroutine, public get_mem_shape(name, mem_path, mem_shape)
@ brief Get the variable memory shape
type(memorystoretype), public memorystore
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
subroutine, public get_mem_rank(name, mem_path, rank)
@ brief Get the variable rank
subroutine, public get_mem_elem_size(name, mem_path, size)
@ brief Get the memory size of a single element of the stored variable
subroutine, public on_memory_set(var_name, mem_path, status)
Triggers the calling of the side effect handler for this variable.
This module contains the MODFLOW 6 BMI.
Definition: mf6bmi.f90:18
integer(kind=c_int) function get_end_time(end_time)
Get the end time of the simulation.
Definition: mf6bmi.f90:150
integer(kind=c_int) function get_input_var_names(c_names)
Returns all input variables in the simulation.
Definition: mf6bmi.f90:250
integer(kind=c_int) function get_input_item_count(count)
Get the number of input variables in the simulation.
Definition: mf6bmi.f90:206
integer(kind=c_int) function get_value_ptr_double(c_var_address, c_arr_ptr)
Get a pointer to the array of double precision numbers.
Definition: mf6bmi.f90:755
integer(kind=c_int) function get_value_ptr_int(c_var_address, c_arr_ptr)
Get a pointer to the array of integer numbers.
Definition: mf6bmi.f90:809
integer(kind=c_int) function get_current_time(current_time)
Get the current time of the simulation.
Definition: mf6bmi.f90:169
integer(kind=c_int) function get_var_rank(c_var_address, c_var_rank)
Get the variable rank (non-BMI)
Definition: mf6bmi.f90:1193
integer(kind=c_int) function get_output_item_count(count)
Get the number of output variables in the simulation.
Definition: mf6bmi.f90:223
integer(kind=c_int) function get_value_ptr_bool(c_var_address, c_arr_ptr)
Get a pointer to the logical scalar value.
Definition: mf6bmi.f90:862
integer(kind=c_int) function set_value(c_var_address, c_arr_ptr)
Set new values for a given variable.
Definition: mf6bmi.f90:904
integer(kind=c_int) function get_time_step(time_step)
Get the time step for the simulation.
Definition: mf6bmi.f90:188
integer(kind=c_int) function get_var_nbytes(c_var_address, var_nbytes)
Get size of the variable, in bytes.
Definition: mf6bmi.f90:343
integer(kind=c_int) function set_value_double(c_var_address, c_arr_ptr)
Set new values for a variable of type double.
Definition: mf6bmi.f90:950
integer(kind=c_int) function get_start_time(start_time)
Get the start time of the simulation.
Definition: mf6bmi.f90:134
integer(kind=c_int) function set_value_int(c_var_address, c_arr_ptr)
Set new values for a variable of type integer.
Definition: mf6bmi.f90:1023
integer(kind=c_int) function get_var_itemsize(c_var_address, var_size)
Get the size (in bytes) of a single element of a variable.
Definition: mf6bmi.f90:316
integer(kind=c_int) function bmi_get_component_name(name)
Definition: mf6bmi.f90:45
integer(kind=c_int) function bmi_update()
Perform a computational time step.
Definition: mf6bmi.f90:96
integer(kind=c_int) function get_value_ptr(c_var_address, c_arr_ptr)
Get a pointer to an array.
Definition: mf6bmi.f90:708
integer(kind=c_int) function get_var_shape(c_var_address, c_var_shape)
Get the shape of the array for the variable (non-BMI)
Definition: mf6bmi.f90:1231
integer(kind=c_int) function bmi_initialize()
Initialize the computational core.
Definition: mf6bmi.f90:67
integer(kind=c_int) function bmi_finalize()
Clean up the initialized simulation.
Definition: mf6bmi.f90:113
integer(kind=c_int) function get_value_int(c_var_address, c_arr_ptr)
Copy the integer values of a variable into the array.
Definition: mf6bmi.f90:457
integer(kind=c_int) function get_output_var_names(c_names)
Returns all output variables in the simulation.
Definition: mf6bmi.f90:285
integer(kind=c_int) function get_var_type(c_var_address, c_var_type)
Get the variable type as a string.
Definition: mf6bmi.f90:1153
integer(kind=c_int) function get_value_double(c_var_address, c_arr_ptr)
Copy the double precision values of a variable into the array.
Definition: mf6bmi.f90:380
integer(kind=c_int) function set_value_bool(c_var_address, c_arr_ptr)
Set new value for a logical scalar variable.
Definition: mf6bmi.f90:1096
integer(kind=c_int) function get_value(c_var_address, c_arr_ptr)
Copy the value of a variable into the array.
Definition: mf6bmi.f90:659
integer(kind=c_int) function get_value_string(c_var_address, c_arr_ptr)
Copy the string(s) of a variable into the array.
Definition: mf6bmi.f90:580
integer(kind=c_int) function get_value_bool(c_var_address, c_arr_ptr)
Copy the logical scalar value into the array.
Definition: mf6bmi.f90:532
integer(c_int), bind(C, name="ISTDOUTTOFILE") istdout_to_file
output control: =0 to screen, >0 to file
Definition: mf6bmi.f90:38
Detailed error information for the BMI.
Definition: mf6bmiError.f90:6
character(len= *), parameter fmt_general_err
Definition: mf6bmiError.f90:22
integer, parameter bmi_failure
BMI status code for failure (taken from bmi.f90, CSDMS)
Definition: mf6bmiError.f90:12
character(len= *), parameter fmt_unsupported_type
Definition: mf6bmiError.f90:31
character(len=lenerrmessage) bmi_last_error
module variable containing the last error as a Fortran string
Definition: mf6bmiError.f90:20
subroutine report_bmi_error(err_msg)
Sets the last BMI error message and copies it to an exported C-string.
Definition: mf6bmiError.f90:47
character(len= *), parameter fmt_invalid_mem_access
Definition: mf6bmiError.f90:34
character(len= *), parameter fmt_unsupported_rank
Definition: mf6bmiError.f90:28
integer, parameter bmi_success
BMI status code for success (taken from bmi.f90, CSDMS)
Definition: mf6bmiError.f90:13
This module contains helper routines and parameters for the MODFLOW 6 BMI.
Definition: mf6bmiUtil.f90:4
integer(c_int), bind(C, name="BMI_LENVARTYPE") bmi_lenvartype
max. length for variable type C-strings
Definition: mf6bmiUtil.f90:26
integer(c_int), bind(C, name="BMI_LENCOMPONENTNAME") bmi_lencomponentname
component name length, i.e. 'MODFLOW 6'
Definition: mf6bmiUtil.f90:38
subroutine split_address(c_var_address, mem_path, var_name, success)
Split the variable address string.
Definition: mf6bmiUtil.f90:54
pure character(kind=c_char, len=1) function, dimension(length+1) string_to_char_array(string, length)
Convert Fortran string to C-style character string.
Definition: mf6bmiUtil.f90:149
integer(c_int), bind(C, name="BMI_LENVARADDRESS") bmi_lenvaraddress
max. length for the variable's address C-string
Definition: mf6bmiUtil.f90:34
Core MODFLOW 6 module.
Definition: mf6core.f90:8
logical(lgp) function mf6update()
Run a time step.
Definition: mf6core.f90:110
subroutine mf6initialize()
Initialize a simulation.
Definition: mf6core.f90:70
subroutine mf6finalize()
Finalize the simulation.
Definition: mf6core.f90:132
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) iforcestop
forced stop flag (1) forces a call to ustop(..) when the simulation has ended, (0) doesn't
character(len=linelength) simstdout
name of standard out file if screen output is piped to a file
integer(i4b) istdout
unit number for stdout
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
real(dp), pointer, public totalsimtime
time at end of simulation
Definition: tdis.f90:37
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
An iterator used to iterate through a MemoryContainer.