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