MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
MemoryManager.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, lgp, i4b, i8b
4  use constantsmodule, only: dzero, done, &
5  dem3, dem6, dem9, dep3, dep6, dep9, &
11  tabright
12  use simvariablesmodule, only: errmsg
14  use memorytypemodule, only: memorytype
19  use tablemodule, only: tabletype, table_cr
21 
22  implicit none
23  private
24  public :: mem_allocate
25  public :: mem_checkin
26  public :: mem_reallocate
27  public :: mem_setptr
28  public :: mem_copyptr
29  public :: mem_reassignptr
30  public :: mem_deallocate
31  public :: mem_write_usage
32  public :: mem_da
33  public :: mem_set_print_option
34  public :: get_from_memorystore
35 
36  public :: get_mem_type
37  public :: get_mem_rank
38  public :: get_mem_elem_size
39  public :: get_mem_shape
40  public :: get_isize
41  public :: copy_dbl1d
42 
43  public :: memorystore
44  public :: mem_print_detailed
45 
47  type(tabletype), pointer :: memtab => null()
48  integer(I8B) :: nvalues_alogical = 0
49  integer(I8B) :: nvalues_astr = 0
50  integer(I8B) :: nvalues_aint = 0
51  integer(I8B) :: nvalues_adbl = 0
52  integer(I4B) :: iprmem = 0
53 
54  interface mem_allocate
55  module procedure &
57  allocate_str, &
59  allocate_int, &
63  allocate_dbl, &
68  end interface mem_allocate
69 
70  interface mem_checkin
71  module procedure &
72  checkin_int1d, &
73  checkin_int2d, &
74  checkin_dbl1d, &
75  checkin_dbl2d, &
77  end interface mem_checkin
78 
79  interface mem_reallocate
80  module procedure &
87  end interface mem_reallocate
88 
89  interface mem_setptr
90  module procedure &
92  setptr_int, &
93  setptr_int1d, &
94  setptr_int2d, &
95  setptr_int3d, &
96  setptr_dbl, &
97  setptr_dbl1d, &
98  setptr_dbl2d, &
99  setptr_dbl3d, &
100  setptr_str, &
101  setptr_str1d, &
103  end interface mem_setptr
104 
105  interface mem_copyptr
106  module procedure &
107  copyptr_int1d, &
108  copyptr_int2d, &
109  copyptr_dbl1d, &
111  end interface mem_copyptr
112 
113  interface mem_reassignptr
114  module procedure &
115  reassignptr_int, &
120  end interface mem_reassignptr
121 
122  interface mem_deallocate
123  module procedure &
125  deallocate_str, &
128  deallocate_int, &
132  deallocate_dbl, &
136  end interface mem_deallocate
137 
138 contains
139 
140  !> @ brief Get the variable memory type
141  !!
142  !! Returns any of 'LOGICAL', 'INTEGER', 'DOUBLE', 'STRING'.
143  !! returns 'UNKNOWN' when the variable is not found.
144  !<
145  subroutine get_mem_type(name, mem_path, var_type)
146  character(len=*), intent(in) :: name !< variable name
147  character(len=*), intent(in) :: mem_path !< path where the variable is stored
148  character(len=LENMEMTYPE), intent(out) :: var_type !< memory type
149  ! -- local
150  type(memorytype), pointer :: mt
151  logical(LGP) :: found
152  ! -- code
153  mt => null()
154  var_type = 'UNKNOWN'
155  call get_from_memorystore(name, mem_path, mt, found)
156  if (found) then
157  var_type = mt%memtype
158  end if
159  end subroutine get_mem_type
160 
161  !> @ brief Get the variable rank
162  !!
163  !! Returns rank = -1 when not found.
164  !<
165  subroutine get_mem_rank(name, mem_path, rank)
166  character(len=*), intent(in) :: name !< variable name
167  character(len=*), intent(in) :: mem_path !< mem_path
168  integer(I4B), intent(out) :: rank !< rank
169  ! -- local
170  type(memorytype), pointer :: mt => null()
171  logical(LGP) :: found
172  ! -- code
173  !
174  ! -- initialize rank to a value to communicate failure
175  rank = -1
176  !
177  ! -- get the entry from the memory manager
178  call get_from_memorystore(name, mem_path, mt, found)
179  !
180  ! -- set rank
181  if (found) then
182  if (associated(mt%logicalsclr)) rank = 0
183  if (associated(mt%intsclr)) rank = 0
184  if (associated(mt%dblsclr)) rank = 0
185  if (associated(mt%aint1d)) rank = 1
186  if (associated(mt%aint2d)) rank = 2
187  if (associated(mt%aint3d)) rank = 3
188  if (associated(mt%adbl1d)) rank = 1
189  if (associated(mt%adbl2d)) rank = 2
190  if (associated(mt%adbl3d)) rank = 3
191  if (associated(mt%strsclr)) rank = 0
192  if (associated(mt%astr1d)) rank = 1
193  if (associated(mt%acharstr1d)) rank = 1
194  end if
195  end subroutine get_mem_rank
196 
197  !> @ brief Get the memory size of a single element of the stored variable
198  !!
199  !! Memory size in bytes, returns size = -1 when not found. This is
200  !< also string length.
201  subroutine get_mem_elem_size(name, mem_path, size)
202  character(len=*), intent(in) :: name !< variable name
203  character(len=*), intent(in) :: mem_path !< path where the variable is stored
204  integer(I4B), intent(out) :: size !< size of the variable in bytes
205  ! -- local
206  type(memorytype), pointer :: mt => null()
207  logical(LGP) :: found
208  ! -- code
209  !
210  ! -- initialize size to a value to communicate failure
211  size = -1
212  !
213  ! -- get the entry from the memory manager
214  call get_from_memorystore(name, mem_path, mt, found)
215  !
216  ! -- set memory size
217  if (found) then
218  size = mt%element_size
219  end if
220  end subroutine get_mem_elem_size
221 
222  !> @ brief Get the variable memory shape
223  !!
224  !! Returns an integer array with the shape (Fortran ordering),
225  !! and set shape(1) = -1 when not found.
226  !<
227  subroutine get_mem_shape(name, mem_path, mem_shape)
228  character(len=*), intent(in) :: name !< variable name
229  character(len=*), intent(in) :: mem_path !< path where the variable is stored
230  integer(I4B), dimension(:), intent(out) :: mem_shape !< shape of the variable
231  ! -- local
232  type(memorytype), pointer :: mt => null()
233  logical(LGP) :: found
234  ! -- code
235  !
236  ! -- get the entry from the memory manager
237  call get_from_memorystore(name, mem_path, mt, found)
238  !
239  ! -- set shape
240  if (found) then
241  if (associated(mt%logicalsclr)) mem_shape = shape(mt%logicalsclr)
242  if (associated(mt%intsclr)) mem_shape = shape(mt%logicalsclr)
243  if (associated(mt%dblsclr)) mem_shape = shape(mt%dblsclr)
244  if (associated(mt%aint1d)) mem_shape = shape(mt%aint1d)
245  if (associated(mt%aint2d)) mem_shape = shape(mt%aint2d)
246  if (associated(mt%aint3d)) mem_shape = shape(mt%aint3d)
247  if (associated(mt%adbl1d)) mem_shape = shape(mt%adbl1d)
248  if (associated(mt%adbl2d)) mem_shape = shape(mt%adbl2d)
249  if (associated(mt%adbl3d)) mem_shape = shape(mt%adbl3d)
250  if (associated(mt%strsclr)) mem_shape = shape(mt%strsclr)
251  if (associated(mt%astr1d)) mem_shape = shape(mt%astr1d)
252  if (associated(mt%acharstr1d)) mem_shape = shape(mt%acharstr1d)
253  ! -- to communicate failure
254  else
255  mem_shape(1) = -1
256  end if
257  end subroutine get_mem_shape
258 
259  !> @ brief Get the number of elements for this variable
260  !!
261  !! Returns with isize = -1 when not found.
262  !! Return 1 for scalars.
263  !<
264  subroutine get_isize(name, mem_path, isize)
265  character(len=*), intent(in) :: name !< variable name
266  character(len=*), intent(in) :: mem_path !< path where the variable is stored
267  integer(I4B), intent(out) :: isize !< number of elements (flattened)
268  ! -- local
269  type(memorytype), pointer :: mt => null()
270  logical(LGP) :: found
271  logical(LGP) :: terminate
272  ! -- code
273  !
274  ! -- initialize isize to a value to communicate failure
275  isize = -1
276  !
277  ! -- don't exit program if variable not found
278  terminate = .false.
279  !
280  ! -- get the entry from the memory manager
281  call get_from_memorystore(name, mem_path, mt, found, terminate)
282  !
283  ! -- set isize
284  if (found) then
285  isize = mt%isize
286  end if
287  end subroutine get_isize
288 
289  !> @ brief Get a memory type entry from the memory list
290  !!
291  !! Default value for @par check is .true. which means that this
292  !! routine will kill the program when the memory entry cannot be found.
293  !<
294  subroutine get_from_memorystore(name, mem_path, mt, found, check)
295  character(len=*), intent(in) :: name !< variable name
296  character(len=*), intent(in) :: mem_path !< path where the variable is stored
297  type(memorytype), pointer, intent(inout) :: mt !< memory type entry
298  logical(LGP), intent(out) :: found !< set to .true. when found
299  logical(LGP), intent(in), optional :: check !< to suppress aborting the program when not found,
300  !! set check = .false.
301  ! -- local
302  logical(LGP) check_opt
303  ! -- code
304  mt => memorystore%get(name, mem_path)
305  found = associated(mt)
306 
307  check_opt = .true.
308  if (present(check)) then
309  check_opt = check
310  end if
311  if (check_opt) then
312  if (.not. found) then
313  errmsg = "Programming error in memory manager. Variable '"// &
314  trim(name)//"' in '"//trim(mem_path)//"' cannot be "// &
315  "assigned because it does not exist in memory manager."
316  call store_error(errmsg, terminate=.true.)
317  end if
318  end if
319  end subroutine get_from_memorystore
320 
321  !> @brief Issue allocation error message and stop program execution
322  !<
323  subroutine allocate_error(varname, mem_path, istat, isize)
324  character(len=*), intent(in) :: varname !< variable name
325  character(len=*), intent(in) :: mem_path !< path where the variable is stored
326  integer(I4B), intent(in) :: istat !< status code
327  integer(I4B), intent(in) :: isize !< size of allocation
328  ! -- local
329  character(len=20) :: csize
330  character(len=20) :: cstat
331  ! -- code
332  !
333  ! -- initialize character variables
334  write (csize, '(i0)') isize
335  write (cstat, '(i0)') istat
336  !
337  ! -- create error message
338  errmsg = "Error trying to allocate memory. Path '"//trim(mem_path)// &
339  "' variable name '"//trim(varname)//"' size '"//trim(csize)// &
340  "'. Error message is '"//trim(adjustl(errmsg))// &
341  "'. Status code is "//trim(cstat)//'.'
342  !
343  ! -- store error and stop program execution
344  call store_error(errmsg, terminate=.true.)
345  end subroutine allocate_error
346 
347  !> @brief Allocate a logical scalar
348  !<
349  subroutine allocate_logical(sclr, name, mem_path)
350  logical(LGP), pointer, intent(inout) :: sclr !< variable for allocation
351  character(len=*), intent(in) :: name !< variable name
352  character(len=*), intent(in) :: mem_path !< path where the variable is stored
353  ! -- local
354  integer(I4B) :: istat
355  type(memorytype), pointer :: mt
356  ! -- code
357  !
358  ! -- check variable name length
359  call mem_check_length(name, lenvarname, "variable")
360  !
361  ! -- allocate the logical scalar
362  allocate (sclr, stat=istat, errmsg=errmsg)
363  if (istat /= 0) then
364  call allocate_error(name, mem_path, istat, 1)
365  end if
366  !
367  ! -- update counter
369  !
370  ! -- allocate memory type
371  allocate (mt)
372  !
373  ! -- set memory type
374  mt%logicalsclr => sclr
375  mt%element_size = lgp
376  mt%isize = 1
377  mt%name = name
378  mt%path = mem_path
379  write (mt%memtype, "(a)") 'LOGICAL'
380  !
381  ! -- add memory type to the memory list
382  call memorystore%add(mt)
383  end subroutine allocate_logical
384 
385  !> @brief Allocate a character string
386  !<
387  subroutine allocate_str(sclr, ilen, name, mem_path)
388  integer(I4B), intent(in) :: ilen !< string length
389  character(len=ilen), pointer, intent(inout) :: sclr !< variable for allocation
390  character(len=*), intent(in) :: name !< variable name
391  character(len=*), intent(in) :: mem_path !< path where the variable is stored
392  ! -- local
393  integer(I4B) :: istat
394  type(memorytype), pointer :: mt
395  ! -- format
396  ! -- code
397  !
398  ! -- make sure ilen is greater than 0
399  if (ilen < 1) then
400  errmsg = 'Programming error in allocate_str. ILEN must be greater than 0.'
401  call store_error(errmsg, terminate=.true.)
402  end if
403  !
404  ! -- check variable name length
405  call mem_check_length(name, lenvarname, "variable")
406  !
407  ! -- allocate string
408  allocate (character(len=ilen) :: sclr, stat=istat, errmsg=errmsg)
409  if (istat /= 0) then
410  call allocate_error(name, mem_path, istat, 1)
411  end if
412  !
413  ! -- set sclr to a empty string
414  sclr = ' '
415  !
416  ! -- update counter
417  nvalues_astr = nvalues_astr + ilen
418  !
419  ! -- allocate memory type
420  allocate (mt)
421  !
422  ! -- set memory type
423  mt%strsclr => sclr
424  mt%element_size = ilen
425  mt%isize = 1
426  mt%name = name
427  mt%path = mem_path
428  write (mt%memtype, "(a,' LEN=',i0)") 'STRING', ilen
429  !
430  ! -- add defined length string to the memory manager list
431  call memorystore%add(mt)
432  end subroutine allocate_str
433 
434  !> @brief Allocate a 1-dimensional defined length string array
435  !<
436  subroutine allocate_str1d(astr1d, ilen, nrow, name, mem_path)
437  integer(I4B), intent(in) :: ilen !< string length
438  character(len=ilen), dimension(:), &
439  pointer, contiguous, intent(inout) :: astr1d !< variable for allocation
440  integer(I4B), intent(in) :: nrow !< number of strings in array
441  character(len=*), intent(in) :: name !< variable name
442  character(len=*), intent(in) :: mem_path !< path where the variable is stored
443  ! -- local variables
444  type(memorytype), pointer :: mt
445  character(len=ilen) :: string
446  integer(I4B) :: n
447  integer(I4B) :: istat
448  integer(I4B) :: isize
449  ! -- code
450  !
451  ! -- initialize string
452  string = ''
453  !
454  ! -- make sure ilen is greater than 0
455  if (ilen < 1) then
456  errmsg = 'Programming error in allocate_str1d. '// &
457  'ILEN must be greater than 0.'
458  call store_error(errmsg, terminate=.true.)
459  end if
460  !
461  ! -- check variable name length
462  call mem_check_length(name, lenvarname, "variable")
463  !
464  ! -- calculate isize
465  isize = nrow
466  !
467  ! -- allocate defined length string array
468  allocate (character(len=ilen) :: astr1d(nrow), stat=istat, errmsg=errmsg)
469  !
470  ! -- check for error condition
471  if (istat /= 0) then
472  call allocate_error(name, mem_path, istat, isize)
473  end if
474  !
475  ! -- fill deferred length string with empty string
476  do n = 1, nrow
477  astr1d(n) = string
478  end do
479  !
480  ! -- update counter
481  nvalues_astr = nvalues_astr + isize
482  !
483  ! -- allocate memory type
484  allocate (mt)
485  !
486  ! -- set memory type
487  ! this does not work with gfortran 11.3 and 12.1
488  ! so we have to disable the pointing to astr1d
489  ! mt%astr1d => astr1d
490  mt%element_size = ilen
491  mt%isize = isize
492  mt%name = name
493  mt%path = mem_path
494  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
495  !
496  ! -- add deferred length character array to the memory manager list
497  call memorystore%add(mt)
498  end subroutine allocate_str1d
499 
500  !> @brief Allocate a 1-dimensional array of deferred-length CharacterStringType
501  !<
502  subroutine allocate_charstr1d(acharstr1d, ilen, nrow, name, mem_path)
503  type(characterstringtype), dimension(:), &
504  pointer, contiguous, intent(inout) :: acharstr1d !< variable for allocation
505  integer(I4B), intent(in) :: ilen !< string length
506  integer(I4B), intent(in) :: nrow !< number of strings in array
507  character(len=*), intent(in) :: name !< variable name
508  character(len=*), intent(in) :: mem_path !< path where the variable is stored
509  ! -- local variables
510  character(len=ilen) :: string
511  type(memorytype), pointer :: mt
512  integer(I4B) :: n
513  integer(I4B) :: istat
514  integer(I4B) :: isize
515  ! -- code
516  !
517  ! -- initialize string
518  string = ''
519  !
520  ! -- check variable name length
521  call mem_check_length(name, lenvarname, "variable")
522  !
523  ! -- calculate isize
524  isize = nrow
525  !
526  ! -- allocate deferred length string array
527  allocate (acharstr1d(nrow), stat=istat, errmsg=errmsg)
528  !
529  ! -- check for error condition
530  if (istat /= 0) then
531  call allocate_error(name, mem_path, istat, isize)
532  end if
533  !
534  ! -- fill deferred length string with empty string
535  do n = 1, nrow
536  acharstr1d(n) = string
537  end do
538  !
539  ! -- update counter
540  nvalues_astr = nvalues_astr + isize
541  !
542  ! -- allocate memory type
543  allocate (mt)
544  !
545  ! -- set memory type
546  mt%acharstr1d => acharstr1d
547  mt%element_size = ilen
548  mt%isize = isize
549  mt%name = name
550  mt%path = mem_path
551  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
552  !
553  ! -- add deferred length character array to the memory manager list
554  call memorystore%add(mt)
555  end subroutine allocate_charstr1d
556 
557  !> @brief Allocate a integer scalar
558  !<
559  subroutine allocate_int(sclr, name, mem_path)
560  integer(I4B), pointer, intent(inout) :: sclr !< variable for allocation
561  character(len=*), intent(in) :: name !< variable name
562  character(len=*), intent(in) :: mem_path !< path where the variable is stored
563  ! -- local
564  type(memorytype), pointer :: mt
565  integer(I4B) :: istat
566  ! -- code
567  !
568  ! -- check variable name length
569  call mem_check_length(name, lenvarname, "variable")
570  !
571  ! -- allocate integer scalar
572  allocate (sclr, stat=istat, errmsg=errmsg)
573  if (istat /= 0) then
574  call allocate_error(name, mem_path, istat, 1)
575  end if
576  !
577  ! -- update counter
579  !
580  ! -- allocate memory type
581  allocate (mt)
582  !
583  ! -- set memory type
584  mt%intsclr => sclr
585  mt%element_size = i4b
586  mt%isize = 1
587  mt%name = name
588  mt%path = mem_path
589  write (mt%memtype, "(a)") 'INTEGER'
590  !
591  ! -- add memory type to the memory list
592  call memorystore%add(mt)
593  end subroutine allocate_int
594 
595  !> @brief Allocate a 1-dimensional integer array
596  !<
597  subroutine allocate_int1d(aint, nrow, name, mem_path)
598  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< variable for allocation
599  integer(I4B), intent(in) :: nrow !< integer array number of rows
600  character(len=*), intent(in) :: name !< variable name
601  character(len=*), intent(in) :: mem_path !< path where variable is stored
602  ! --local
603  type(memorytype), pointer :: mt
604  integer(I4B) :: istat
605  integer(I4B) :: isize
606  ! -- code
607  !
608  ! -- check variable name length
609  call mem_check_length(name, lenvarname, "variable")
610  !
611  ! -- set isize
612  isize = nrow
613  !
614  ! -- allocate integer array
615  allocate (aint(nrow), stat=istat, errmsg=errmsg)
616  if (istat /= 0) then
617  call allocate_error(name, mem_path, istat, isize)
618  end if
619  !
620  ! -- update counter
621  nvalues_aint = nvalues_aint + isize
622  !
623  ! -- allocate memory type
624  allocate (mt)
625  !
626  ! -- set memory type
627  mt%aint1d => aint
628  mt%element_size = i4b
629  mt%isize = isize
630  mt%name = name
631  mt%path = mem_path
632  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', isize
633  !
634  ! -- add memory type to the memory list
635  call memorystore%add(mt)
636  end subroutine allocate_int1d
637 
638  !> @brief Allocate a 2-dimensional integer array
639  !<
640  subroutine allocate_int2d(aint, ncol, nrow, name, mem_path)
641  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< variable for allocation
642  integer(I4B), intent(in) :: ncol !< number of columns
643  integer(I4B), intent(in) :: nrow !< number of rows
644  character(len=*), intent(in) :: name !< variable name
645  character(len=*), intent(in) :: mem_path !< path where variable is stored
646  ! -- local
647  type(memorytype), pointer :: mt
648  integer(I4B) :: istat
649  integer(I4B) :: isize
650  ! -- code
651  !
652  ! -- check the variable name length
653  call mem_check_length(name, lenvarname, "variable")
654  !
655  ! -- set isize
656  isize = ncol * nrow
657  !
658  ! -- allocate the integer array
659  allocate (aint(ncol, nrow), stat=istat, errmsg=errmsg)
660  if (istat /= 0) then
661  call allocate_error(name, mem_path, istat, isize)
662  end if
663  !
664  ! -- update the counter
665  nvalues_aint = nvalues_aint + isize
666  !
667  ! -- allocate memory type
668  allocate (mt)
669  !
670  ! -- set memory type
671  mt%aint2d => aint
672  mt%element_size = i4b
673  mt%isize = isize
674  mt%name = name
675  mt%path = mem_path
676  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
677  !
678  ! -- add memory type to the memory list
679  call memorystore%add(mt)
680  end subroutine allocate_int2d
681 
682  !> @brief Allocate a 3-dimensional integer array
683  !<
684  subroutine allocate_int3d(aint, ncol, nrow, nlay, name, mem_path)
685  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< variable for allocation
686  integer(I4B), intent(in) :: ncol !< number of columns
687  integer(I4B), intent(in) :: nrow !< number of rows
688  integer(I4B), intent(in) :: nlay !< number of layers
689  character(len=*), intent(in) :: name !< variable name
690  character(len=*), intent(in) :: mem_path !< path where variable is stored
691  ! -- local
692  type(memorytype), pointer :: mt
693  integer(I4B) :: istat
694  integer(I4B) :: isize
695  ! -- code
696  !
697  ! -- check variable name length
698  call mem_check_length(name, lenvarname, "variable")
699  !
700  ! -- set isize
701  isize = ncol * nrow * nlay
702  !
703  ! -- allocate integer array
704  allocate (aint(ncol, nrow, nlay), stat=istat, errmsg=errmsg)
705  if (istat /= 0) then
706  call allocate_error(name, mem_path, istat, isize)
707  end if
708  !
709  ! -- update counter
710  nvalues_aint = nvalues_aint + isize
711  !
712  ! -- allocate memory type
713  allocate (mt)
714  !
715  ! -- set memory type
716  mt%aint3d => aint
717  mt%element_size = i4b
718  mt%isize = isize
719  mt%name = name
720  mt%path = mem_path
721  write (mt%memtype, "(a,' (',i0,',',i0,',',i0,')')") 'INTEGER', ncol, &
722  nrow, nlay
723  !
724  ! -- add memory type to the memory list
725  call memorystore%add(mt)
726  end subroutine allocate_int3d
727 
728  !> @brief Allocate a real scalar
729  !<
730  subroutine allocate_dbl(sclr, name, mem_path)
731  real(DP), pointer, intent(inout) :: sclr !< variable for allocation
732  character(len=*), intent(in) :: name !< variable name
733  character(len=*), intent(in) :: mem_path !< path where variable is stored
734  ! -- local
735  type(memorytype), pointer :: mt
736  integer(I4B) :: istat
737  ! -- code
738  !
739  ! -- check variable name length
740  call mem_check_length(name, lenvarname, "variable")
741  !
742  ! -- allocate real scalar
743  allocate (sclr, stat=istat, errmsg=errmsg)
744  if (istat /= 0) then
745  call allocate_error(name, mem_path, istat, 1)
746  end if
747  !
748  ! -- update counter
750  !
751  ! -- allocate memory type
752  allocate (mt)
753  !
754  ! -- set memory type
755  mt%dblsclr => sclr
756  mt%element_size = dp
757  mt%isize = 1
758  mt%name = name
759  mt%path = mem_path
760  write (mt%memtype, "(a)") 'DOUBLE'
761  !
762  ! -- add memory type to the memory list
763  call memorystore%add(mt)
764  end subroutine allocate_dbl
765 
766  !> @brief Allocate a 1-dimensional real array
767  !<
768  subroutine allocate_dbl1d(adbl, nrow, name, mem_path)
769  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< variable for allocation
770  integer(I4B), intent(in) :: nrow !< number of rows
771  character(len=*), intent(in) :: name !< variable name
772  character(len=*), intent(in) :: mem_path !< path where variable is stored
773  ! -- local
774  type(memorytype), pointer :: mt
775  integer(I4B) :: istat
776  integer(I4B) :: isize
777  ! -- code
778  !
779  ! -- check the variable name length
780  call mem_check_length(name, lenvarname, "variable")
781  !
782  ! -- set isize
783  isize = nrow
784  !
785  ! -- allocate the real array
786  allocate (adbl(nrow), stat=istat, errmsg=errmsg)
787  if (istat /= 0) then
788  call allocate_error(name, mem_path, istat, isize)
789  end if
790  !
791  ! -- update counter
792  nvalues_adbl = nvalues_adbl + isize
793  !
794  ! -- allocate memory type
795  allocate (mt)
796  !
797  ! -- set memory type
798  mt%adbl1d => adbl
799  mt%element_size = dp
800  mt%isize = isize
801  mt%name = name
802  mt%path = mem_path
803  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', isize
804  !
805  ! -- add memory type to the memory list
806  call memorystore%add(mt)
807  end subroutine allocate_dbl1d
808 
809  !> @brief Allocate a 2-dimensional real array
810  !<
811  subroutine allocate_dbl2d(adbl, ncol, nrow, name, mem_path)
812  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< variable for allocation
813  integer(I4B), intent(in) :: ncol !< number of columns
814  integer(I4B), intent(in) :: nrow !< number of rows
815  character(len=*), intent(in) :: name !< variable name
816  character(len=*), intent(in) :: mem_path !< path where variable is stored
817  ! -- local
818  type(memorytype), pointer :: mt
819  integer(I4B) :: istat
820  integer(I4B) :: isize
821  ! -- code
822  !
823  ! -- check the variable name length
824  call mem_check_length(name, lenvarname, "variable")
825  !
826  ! -- set isize
827  isize = ncol * nrow
828  !
829  ! -- allocate the real array
830  allocate (adbl(ncol, nrow), stat=istat, errmsg=errmsg)
831  if (istat /= 0) then
832  call allocate_error(name, mem_path, istat, isize)
833  end if
834  !
835  ! -- update counter
836  nvalues_adbl = nvalues_adbl + isize
837  !
838  ! -- allocate memory type
839  allocate (mt)
840  !
841  ! -- set memory type
842  mt%adbl2d => adbl
843  mt%element_size = dp
844  mt%isize = isize
845  mt%name = name
846  mt%path = mem_path
847  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
848  !
849  ! -- add memory type to the memory list
850  call memorystore%add(mt)
851  end subroutine allocate_dbl2d
852 
853  !> @brief Allocate a 3-dimensional real array
854  !<
855  subroutine allocate_dbl3d(adbl, ncol, nrow, nlay, name, mem_path)
856  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< variable for allocation
857  integer(I4B), intent(in) :: ncol !< number of columns
858  integer(I4B), intent(in) :: nrow !< number of rows
859  integer(I4B), intent(in) :: nlay !< number of layers
860  character(len=*), intent(in) :: name !< variable name
861  character(len=*), intent(in) :: mem_path !< path where variable is stored
862  ! -- local
863  type(memorytype), pointer :: mt
864  integer(I4B) :: istat
865  integer(I4B) :: isize
866  ! -- code
867  !
868  ! -- check the variable name length
869  call mem_check_length(name, lenvarname, "variable")
870  !
871  ! -- set isize
872  isize = ncol * nrow * nlay
873  !
874  ! -- allocate the real array
875  allocate (adbl(ncol, nrow, nlay), stat=istat, errmsg=errmsg)
876  if (istat /= 0) then
877  call allocate_error(name, mem_path, istat, isize)
878  end if
879  !
880  ! -- update the counter
881  nvalues_adbl = nvalues_adbl + isize
882  !
883  ! -- allocate memory type
884  allocate (mt)
885  !
886  ! -- set memory type
887  mt%adbl3d => adbl
888  mt%element_size = dp
889  mt%isize = isize
890  mt%name = name
891  mt%path = mem_path
892  write (mt%memtype, "(a,' (',i0,',',i0,',',i0,')')") 'DOUBLE', ncol, &
893  nrow, nlay
894  !
895  ! -- add memory type to the memory list
896  call memorystore%add(mt)
897  end subroutine allocate_dbl3d
898 
899  !> @brief Check in an existing 1d integer array with a new address (name + path)
900  !<
901  subroutine checkin_int1d(aint, name, mem_path, name2, mem_path2)
902  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: aint !< the existing array
903  character(len=*), intent(in) :: name !< new variable name
904  character(len=*), intent(in) :: mem_path !< new path where variable is stored
905  character(len=*), intent(in) :: name2 !< existing variable name
906  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
907  ! --local
908  type(memorytype), pointer :: mt
909  integer(I4B) :: isize
910  ! -- code
911  !
912  ! -- check variable name length
913  call mem_check_length(name, lenvarname, "variable")
914  !
915  ! -- set isize
916  isize = size(aint)
917  !
918  ! -- allocate memory type
919  allocate (mt)
920  !
921  ! -- set memory type
922  mt%aint1d => aint
923  mt%element_size = i4b
924  mt%isize = isize
925  mt%name = name
926  mt%path = mem_path
927  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', isize
928  !
929  ! -- set master information
930  mt%master = .false.
931  mt%mastername = name2
932  mt%masterPath = mem_path2
933  !
934  ! -- add memory type to the memory list
935  call memorystore%add(mt)
936  end subroutine checkin_int1d
937 
938  !> @brief Check in an existing 2d integer array with a new address (name + path)
939  !<
940  subroutine checkin_int2d(aint2d, name, mem_path, name2, mem_path2)
941  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint2d !< the existing 2d array
942  character(len=*), intent(in) :: name !< new variable name
943  character(len=*), intent(in) :: mem_path !< new path where variable is stored
944  character(len=*), intent(in) :: name2 !< existing variable name
945  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
946  ! -- local
947  type(memorytype), pointer :: mt
948  integer(I4B) :: ncol, nrow, isize
949  ! -- code
950  !
951  ! -- check the variable name length
952  call mem_check_length(name, lenvarname, "variable")
953  !
954  ! -- set isize
955  ncol = size(aint2d, dim=1)
956  nrow = size(aint2d, dim=2)
957  isize = ncol * nrow
958  !
959  ! -- allocate memory type
960  allocate (mt)
961  !
962  ! -- set memory type
963  mt%aint2d => aint2d
964  mt%isize = isize
965  mt%name = name
966  mt%path = mem_path
967  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
968  !
969  ! -- set master information
970  mt%master = .false.
971  mt%mastername = name2
972  mt%masterPath = mem_path2
973  !
974  ! -- add memory type to the memory list
975  call memorystore%add(mt)
976  end subroutine checkin_int2d
977 
978  !> @brief Check in an existing 1d double precision array with a new address (name + path)
979  !<
980  subroutine checkin_dbl1d(adbl, name, mem_path, name2, mem_path2)
981  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< the existing array
982  character(len=*), intent(in) :: name !< new variable name
983  character(len=*), intent(in) :: mem_path !< new path where variable is stored
984  character(len=*), intent(in) :: name2 !< existing variable name
985  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
986  ! -- local
987  type(memorytype), pointer :: mt
988  integer(I4B) :: isize
989  ! -- code
990  !
991  ! -- check the variable name length
992  call mem_check_length(name, lenvarname, "variable")
993  !
994  ! -- set isize
995  isize = size(adbl)
996  !
997  ! -- allocate memory type
998  allocate (mt)
999  !
1000  ! -- set memory type
1001  mt%adbl1d => adbl
1002  mt%element_size = dp
1003  mt%isize = isize
1004  mt%name = name
1005  mt%path = mem_path
1006  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', isize
1007  !
1008  ! -- set master information
1009  mt%master = .false.
1010  mt%mastername = name2
1011  mt%masterPath = mem_path2
1012  !
1013  ! -- add memory type to the memory list
1014  call memorystore%add(mt)
1015  end subroutine checkin_dbl1d
1016 
1017  !> @brief Check in an existing 2d double precision array with a new address (name + path)
1018  !<
1019  subroutine checkin_dbl2d(adbl2d, name, mem_path, name2, mem_path2)
1020  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl2d !< the existing 2d array
1021  character(len=*), intent(in) :: name !< new variable name
1022  character(len=*), intent(in) :: mem_path !< new path where variable is stored
1023  character(len=*), intent(in) :: name2 !< existing variable name
1024  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
1025  ! -- local
1026  type(memorytype), pointer :: mt
1027  integer(I4B) :: ncol, nrow, isize
1028  ! -- code
1029  !
1030  ! -- check the variable name length
1031  call mem_check_length(name, lenvarname, "variable")
1032  !
1033  ! -- set isize
1034  ncol = size(adbl2d, dim=1)
1035  nrow = size(adbl2d, dim=2)
1036  isize = ncol * nrow
1037  !
1038  ! -- allocate memory type
1039  allocate (mt)
1040  !
1041  ! -- set memory type
1042  mt%adbl2d => adbl2d
1043  mt%isize = isize
1044  mt%name = name
1045  mt%path = mem_path
1046  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
1047  !
1048  ! -- set master information
1049  mt%master = .false.
1050  mt%mastername = name2
1051  mt%masterPath = mem_path2
1052  !
1053  ! -- add memory type to the memory list
1054  call memorystore%add(mt)
1055  end subroutine checkin_dbl2d
1056 
1057  !> @brief Check in an existing 1d CharacterStringType array with a new address (name + path)
1058  !<
1059  subroutine checkin_charstr1d(acharstr1d, ilen, name, mem_path, name2, mem_path2)
1060  type(characterstringtype), dimension(:), &
1061  pointer, contiguous, intent(inout) :: acharstr1d !< the existing array
1062  integer(I4B), intent(in) :: ilen
1063  character(len=*), intent(in) :: name !< new variable name
1064  character(len=*), intent(in) :: mem_path !< new path where variable is stored
1065  character(len=*), intent(in) :: name2 !< existing variable name
1066  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
1067  ! --local
1068  type(memorytype), pointer :: mt
1069  integer(I4B) :: isize
1070  ! -- code
1071  !
1072  ! -- check variable name length
1073  call mem_check_length(name, lenvarname, "variable")
1074  !
1075  ! -- set isize
1076  isize = size(acharstr1d)
1077  !
1078  ! -- allocate memory type
1079  allocate (mt)
1080  !
1081  ! -- set memory type
1082  mt%acharstr1d => acharstr1d
1083  mt%element_size = ilen
1084  mt%isize = isize
1085  mt%name = name
1086  mt%path = mem_path
1087  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, isize
1088  !
1089  ! -- set master information
1090  mt%master = .false.
1091  mt%mastername = name2
1092  mt%masterPath = mem_path2
1093  !
1094  ! -- add memory type to the memory list
1095  call memorystore%add(mt)
1096  end subroutine checkin_charstr1d
1097 
1098  !> @brief Reallocate a 1-dimensional defined length string array
1099  !<
1100  subroutine reallocate_str1d(astr, ilen, nrow, name, mem_path)
1101  integer(I4B), intent(in) :: ilen !< string length
1102  integer(I4B), intent(in) :: nrow !< number of rows
1103  character(len=ilen), dimension(:), pointer, contiguous, intent(inout) :: astr !< the reallocated string array
1104  character(len=*), intent(in) :: name !< variable name
1105  character(len=*), intent(in) :: mem_path !< path where variable is stored
1106  ! -- local
1107  type(memorytype), pointer :: mt
1108  logical(LGP) :: found
1109  character(len=ilen), dimension(:), allocatable :: astrtemp
1110  integer(I4B) :: istat
1111  integer(I4B) :: isize
1112  integer(I4B) :: isize_old
1113  integer(I4B) :: nrow_old
1114  integer(I4B) :: n
1115  !
1116  ! -- Find and assign mt
1117  call get_from_memorystore(name, mem_path, mt, found)
1118  !
1119  ! -- reallocate astr1d
1120  if (found) then
1121  isize_old = mt%isize
1122  if (isize_old > 0) then
1123  nrow_old = size(astr)
1124  else
1125  nrow_old = 0
1126  end if
1127  !
1128  ! -- calculate isize
1129  isize = nrow
1130  !
1131  ! -- allocate astrtemp
1132  allocate (astrtemp(nrow), stat=istat, errmsg=errmsg)
1133  if (istat /= 0) then
1134  call allocate_error(name, mem_path, istat, isize)
1135  end if
1136  !
1137  ! -- copy existing values
1138  do n = 1, nrow_old
1139  astrtemp(n) = astr(n)
1140  end do
1141  !
1142  ! -- fill new values with missing values
1143  do n = nrow_old + 1, nrow
1144  astrtemp(n) = ''
1145  end do
1146  !
1147  ! -- deallocate mt pointer, repoint, recalculate isize
1148  deallocate (astr)
1149  !
1150  ! -- allocate astr1d
1151  allocate (astr(nrow), stat=istat, errmsg=errmsg)
1152  if (istat /= 0) then
1153  call allocate_error(name, mem_path, istat, isize)
1154  end if
1155  !
1156  ! -- fill the reallocate character array
1157  do n = 1, nrow
1158  astr(n) = astrtemp(n)
1159  end do
1160  !
1161  ! -- deallocate temporary storage
1162  deallocate (astrtemp)
1163  !
1164  ! -- reset memory manager values
1165  mt%element_size = ilen
1166  mt%isize = isize
1167  mt%nrealloc = mt%nrealloc + 1
1168  mt%master = .true.
1169  nvalues_astr = nvalues_astr + isize - isize_old
1170  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
1171  else
1172  errmsg = "Programming error, variable '"//trim(name)//"' from '"// &
1173  trim(mem_path)//"' is not defined in the memory manager. Use "// &
1174  "mem_allocate instead."
1175  call store_error(errmsg, terminate=.true.)
1176  end if
1177  end subroutine reallocate_str1d
1178 
1179  !> @brief Reallocate a 1-dimensional deferred length string array
1180  !<
1181  subroutine reallocate_charstr1d(acharstr1d, ilen, nrow, name, mem_path)
1182  type(characterstringtype), dimension(:), pointer, contiguous, &
1183  intent(inout) :: acharstr1d !< the reallocated charstring array
1184  integer(I4B), intent(in) :: ilen !< string length
1185  integer(I4B), intent(in) :: nrow !< number of rows
1186  character(len=*), intent(in) :: name !< variable name
1187  character(len=*), intent(in) :: mem_path !< path where variable is stored
1188  ! -- local
1189  type(memorytype), pointer :: mt
1190  logical(LGP) :: found
1191  type(characterstringtype), dimension(:), allocatable :: astrtemp
1192  character(len=ilen) :: string
1193  integer(I4B) :: istat
1194  integer(I4B) :: isize
1195  integer(I4B) :: isize_old
1196  integer(I4B) :: nrow_old
1197  integer(I4B) :: n
1198  !
1199  ! -- Initialize string
1200  string = ''
1201  !
1202  ! -- Find and assign mt
1203  call get_from_memorystore(name, mem_path, mt, found)
1204  !
1205  ! -- reallocate astr1d
1206  if (found) then
1207  isize_old = mt%isize
1208  if (isize_old > 0) then
1209  nrow_old = size(acharstr1d)
1210  else
1211  nrow_old = 0
1212  end if
1213  !
1214  ! -- calculate isize
1215  isize = nrow
1216  !
1217  ! -- allocate astrtemp
1218  allocate (astrtemp(nrow), stat=istat, errmsg=errmsg)
1219  if (istat /= 0) then
1220  call allocate_error(name, mem_path, istat, isize)
1221  end if
1222  !
1223  ! -- copy existing values
1224  do n = 1, nrow_old
1225  astrtemp(n) = acharstr1d(n)
1226  call acharstr1d(n)%destroy()
1227  end do
1228  !
1229  ! -- fill new values with missing values
1230  do n = nrow_old + 1, nrow
1231  astrtemp(n) = string
1232  end do
1233  !
1234  ! -- deallocate mt pointer, repoint, recalculate isize
1235  deallocate (acharstr1d)
1236  !
1237  ! -- allocate astr1d
1238  allocate (acharstr1d(nrow), stat=istat, errmsg=errmsg)
1239  if (istat /= 0) then
1240  call allocate_error(name, mem_path, istat, isize)
1241  end if
1242  !
1243  ! -- fill the reallocated character array
1244  do n = 1, nrow
1245  acharstr1d(n) = astrtemp(n)
1246  call astrtemp(n)%destroy()
1247  end do
1248  !
1249  ! -- deallocate temporary storage
1250  deallocate (astrtemp)
1251  !
1252  ! -- reset memory manager values
1253  mt%acharstr1d => acharstr1d
1254  mt%element_size = ilen
1255  mt%isize = isize
1256  mt%nrealloc = mt%nrealloc + 1
1257  mt%master = .true.
1258  nvalues_astr = nvalues_astr + isize - isize_old
1259  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
1260  else
1261  errmsg = "Programming error, variable '"//trim(name)//"' from '"// &
1262  trim(mem_path)//"' is not defined in the memory manager. Use "// &
1263  "mem_allocate instead."
1264  call store_error(errmsg, terminate=.true.)
1265  end if
1266  end subroutine reallocate_charstr1d
1267 
1268  !> @brief Reallocate a 1-dimensional integer array
1269  !<
1270  subroutine reallocate_int1d(aint, nrow, name, mem_path)
1271  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< the reallocated integer array
1272  integer(I4B), intent(in) :: nrow !< number of rows
1273  character(len=*), intent(in) :: name !< variable name
1274  character(len=*), intent(in) :: mem_path !< path where variable is stored
1275  ! -- local
1276  type(memorytype), pointer :: mt
1277  logical(LGP) :: found
1278  integer(I4B) :: istat
1279  integer(I4B) :: isize
1280  integer(I4B) :: i
1281  integer(I4B) :: isizeold
1282  integer(I4B) :: ifill
1283  ! -- code
1284  !
1285  ! -- Find and assign mt
1286  call get_from_memorystore(name, mem_path, mt, found)
1287  !
1288  ! -- Allocate aint and then refill
1289  isize = nrow
1290  isizeold = size(mt%aint1d)
1291  ifill = min(isizeold, isize)
1292  allocate (aint(nrow), stat=istat, errmsg=errmsg)
1293  if (istat /= 0) then
1294  call allocate_error(name, mem_path, istat, isize)
1295  end if
1296  do i = 1, ifill
1297  aint(i) = mt%aint1d(i)
1298  end do
1299  !
1300  ! -- deallocate mt pointer, repoint, recalculate isize
1301  deallocate (mt%aint1d)
1302  mt%aint1d => aint
1303  mt%element_size = i4b
1304  mt%isize = isize
1305  mt%nrealloc = mt%nrealloc + 1
1306  mt%master = .true.
1307  nvalues_aint = nvalues_aint + isize - isizeold
1308  end subroutine reallocate_int1d
1309 
1310  !> @brief Reallocate a 2-dimensional integer array
1311  !<
1312  subroutine reallocate_int2d(aint, ncol, nrow, name, mem_path)
1313  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< the reallocated 2d integer array
1314  integer(I4B), intent(in) :: ncol !< number of columns
1315  integer(I4B), intent(in) :: nrow !< number of rows
1316  character(len=*), intent(in) :: name !< variable name
1317  character(len=*), intent(in) :: mem_path !< path where variable is stored
1318  ! -- local
1319  type(memorytype), pointer :: mt
1320  logical(LGP) :: found
1321  integer(I4B) :: istat
1322  integer(I4B), dimension(2) :: ishape
1323  integer(I4B) :: i
1324  integer(I4B) :: j
1325  integer(I4B) :: isize
1326  integer(I4B) :: isizeold
1327  ! -- code
1328  !
1329  ! -- Find and assign mt
1330  call get_from_memorystore(name, mem_path, mt, found)
1331  !
1332  ! -- Allocate aint and then refill
1333  ishape = shape(mt%aint2d)
1334  isize = nrow * ncol
1335  isizeold = ishape(1) * ishape(2)
1336  allocate (aint(ncol, nrow), stat=istat, errmsg=errmsg)
1337  if (istat /= 0) then
1338  call allocate_error(name, mem_path, istat, isize)
1339  end if
1340  do i = 1, ishape(2)
1341  do j = 1, ishape(1)
1342  aint(j, i) = mt%aint2d(j, i)
1343  end do
1344  end do
1345  !
1346  ! -- deallocate mt pointer, repoint, recalculate isize
1347  deallocate (mt%aint2d)
1348  mt%aint2d => aint
1349  mt%element_size = i4b
1350  mt%isize = isize
1351  mt%nrealloc = mt%nrealloc + 1
1352  mt%master = .true.
1353  nvalues_aint = nvalues_aint + isize - isizeold
1354  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
1355  end subroutine reallocate_int2d
1356 
1357  !> @brief Reallocate a 1-dimensional real array
1358  !<
1359  subroutine reallocate_dbl1d(adbl, nrow, name, mem_path)
1360  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< the reallocated 1d real array
1361  integer(I4B), intent(in) :: nrow !< number of rows
1362  character(len=*), intent(in) :: name !< variable name
1363  character(len=*), intent(in) :: mem_path !< path where variable is stored
1364  ! -- local
1365  type(memorytype), pointer :: mt
1366  integer(I4B) :: istat
1367  integer(I4B) :: isize
1368  integer(I4B) :: i
1369  integer(I4B) :: isizeold
1370  integer(I4B) :: ifill
1371  logical(LGP) :: found
1372  ! -- code
1373  !
1374  ! -- Find and assign mt
1375  call get_from_memorystore(name, mem_path, mt, found)
1376  !
1377  ! -- Allocate adbl and then refill
1378  isize = nrow
1379  isizeold = size(mt%adbl1d)
1380  ifill = min(isizeold, isize)
1381  allocate (adbl(nrow), stat=istat, errmsg=errmsg)
1382  if (istat /= 0) then
1383  call allocate_error(name, mem_path, istat, isize)
1384  end if
1385  do i = 1, ifill
1386  adbl(i) = mt%adbl1d(i)
1387  end do
1388  !
1389  ! -- deallocate mt pointer, repoint, recalculate isize
1390  deallocate (mt%adbl1d)
1391  mt%adbl1d => adbl
1392  mt%element_size = dp
1393  mt%isize = isize
1394  mt%nrealloc = mt%nrealloc + 1
1395  mt%master = .true.
1396  nvalues_adbl = nvalues_adbl + isize - isizeold
1397  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', isize
1398  end subroutine reallocate_dbl1d
1399 
1400  !> @brief Reallocate a 2-dimensional real array
1401  !<
1402  subroutine reallocate_dbl2d(adbl, ncol, nrow, name, mem_path)
1403  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< the reallocated 2d real array
1404  integer(I4B), intent(in) :: ncol !< number of columns
1405  integer(I4B), intent(in) :: nrow !< number of rows
1406  character(len=*), intent(in) :: name !< variable name
1407  character(len=*), intent(in) :: mem_path !< path where variable is stored
1408  ! -- local
1409  type(memorytype), pointer :: mt
1410  logical(LGP) :: found
1411  integer(I4B) :: istat
1412  integer(I4B), dimension(2) :: ishape
1413  integer(I4B) :: i
1414  integer(I4B) :: j
1415  integer(I4B) :: isize
1416  integer(I4B) :: isizeold
1417  ! -- code
1418  !
1419  ! -- Find and assign mt
1420  call get_from_memorystore(name, mem_path, mt, found)
1421  !
1422  ! -- Allocate adbl and then refill
1423  ishape = shape(mt%adbl2d)
1424  isize = nrow * ncol
1425  isizeold = ishape(1) * ishape(2)
1426  allocate (adbl(ncol, nrow), stat=istat, errmsg=errmsg)
1427  if (istat /= 0) then
1428  call allocate_error(name, mem_path, istat, isize)
1429  end if
1430  do i = 1, ishape(2)
1431  do j = 1, ishape(1)
1432  adbl(j, i) = mt%adbl2d(j, i)
1433  end do
1434  end do
1435  !
1436  ! -- deallocate mt pointer, repoint, recalculate isize
1437  deallocate (mt%adbl2d)
1438  mt%adbl2d => adbl
1439  mt%element_size = dp
1440  mt%isize = isize
1441  mt%nrealloc = mt%nrealloc + 1
1442  mt%master = .true.
1443  nvalues_adbl = nvalues_adbl + isize - isizeold
1444  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
1445  end subroutine reallocate_dbl2d
1446 
1447  !> @brief Set pointer to a logical scalar
1448  !<
1449  subroutine setptr_logical(sclr, name, mem_path)
1450  logical(LGP), pointer, intent(inout) :: sclr !< pointer to logical scalar
1451  character(len=*), intent(in) :: name !< variable name
1452  character(len=*), intent(in) :: mem_path !< path where variable is stored
1453  ! -- local
1454  type(memorytype), pointer :: mt
1455  logical(LGP) :: found
1456  ! -- code
1457  call get_from_memorystore(name, mem_path, mt, found)
1458  sclr => mt%logicalsclr
1459  end subroutine setptr_logical
1460 
1461  !> @brief Set pointer to integer scalar
1462  !<
1463  subroutine setptr_int(sclr, name, mem_path)
1464  integer(I4B), pointer, intent(inout) :: sclr !< pointer to integer scalar
1465  character(len=*), intent(in) :: name !< variable name
1466  character(len=*), intent(in) :: mem_path !< path where variable is stored
1467  ! -- local
1468  type(memorytype), pointer :: mt
1469  logical(LGP) :: found
1470  ! -- code
1471  call get_from_memorystore(name, mem_path, mt, found)
1472  sclr => mt%intsclr
1473  end subroutine setptr_int
1474 
1475  !> @brief Set pointer to 1d integer array
1476  !<
1477  subroutine setptr_int1d(aint, name, mem_path)
1478  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< pointer to 1d integer array
1479  character(len=*), intent(in) :: name !< variable name
1480  character(len=*), intent(in) :: mem_path !< path where variable is stored
1481  ! -- local
1482  type(memorytype), pointer :: mt
1483  logical(LGP) :: found
1484  ! -- code
1485  call get_from_memorystore(name, mem_path, mt, found)
1486  aint => mt%aint1d
1487  end subroutine setptr_int1d
1488 
1489  !> @brief Set pointer to 2d integer array
1490  !<
1491  subroutine setptr_int2d(aint, name, mem_path)
1492  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< pointer to 2d integer array
1493  character(len=*), intent(in) :: name !< variable name
1494  character(len=*), intent(in) :: mem_path !< path where variable is stored
1495  ! -- local
1496  type(memorytype), pointer :: mt
1497  logical(LGP) :: found
1498  ! -- code
1499  call get_from_memorystore(name, mem_path, mt, found)
1500  aint => mt%aint2d
1501  end subroutine setptr_int2d
1502 
1503  !> @brief Set pointer to 3d integer array
1504  !<
1505  subroutine setptr_int3d(aint, name, mem_path)
1506  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< pointer to 3d integer array
1507  character(len=*), intent(in) :: name !< variable name
1508  character(len=*), intent(in) :: mem_path !< path where variable is stored
1509  ! -- local
1510  type(memorytype), pointer :: mt
1511  logical(LGP) :: found
1512  ! -- code
1513  call get_from_memorystore(name, mem_path, mt, found)
1514  aint => mt%aint3d
1515  end subroutine setptr_int3d
1516 
1517  !> @brief Set pointer to a real scalar
1518  !<
1519  subroutine setptr_dbl(sclr, name, mem_path)
1520  real(DP), pointer, intent(inout) :: sclr !< pointer to a real scalar
1521  character(len=*), intent(in) :: name !< variable name
1522  character(len=*), intent(in) :: mem_path !< path where variable is stored
1523  ! -- local
1524  type(memorytype), pointer :: mt
1525  logical(LGP) :: found
1526  ! -- code
1527  call get_from_memorystore(name, mem_path, mt, found)
1528  sclr => mt%dblsclr
1529  end subroutine setptr_dbl
1530 
1531  !> @brief Set pointer to a 1d real array
1532  !<
1533  subroutine setptr_dbl1d(adbl, name, mem_path)
1534  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< pointer to 1d real array
1535  character(len=*), intent(in) :: name !< variable name
1536  character(len=*), intent(in) :: mem_path !< path where variable is stored
1537  ! -- local
1538  type(memorytype), pointer :: mt
1539  logical(LGP) :: found
1540  ! -- code
1541  call get_from_memorystore(name, mem_path, mt, found)
1542  adbl => mt%adbl1d
1543  end subroutine setptr_dbl1d
1544 
1545  !> @brief Set pointer to a 2d real array
1546  !<
1547  subroutine setptr_dbl2d(adbl, name, mem_path)
1548  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 2d real array
1549  character(len=*), intent(in) :: name !< variable name
1550  character(len=*), intent(in) :: mem_path !< path where variable is stored
1551  ! -- local
1552  type(memorytype), pointer :: mt
1553  logical(LGP) :: found
1554  ! -- code
1555  call get_from_memorystore(name, mem_path, mt, found)
1556  adbl => mt%adbl2d
1557  end subroutine setptr_dbl2d
1558 
1559  !> @brief Set pointer to a 3d real array
1560  !<
1561  subroutine setptr_dbl3d(adbl, name, mem_path)
1562  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 3d real array
1563  character(len=*), intent(in) :: name !< variable name
1564  character(len=*), intent(in) :: mem_path !< path where variable is stored
1565  ! -- local
1566  type(memorytype), pointer :: mt
1567  logical(LGP) :: found
1568  ! -- code
1569  call get_from_memorystore(name, mem_path, mt, found)
1570  adbl => mt%adbl3d
1571  end subroutine setptr_dbl3d
1572 
1573  !> @brief Set pointer to a string (scalar)
1574  !<
1575  subroutine setptr_str(asrt, name, mem_path)
1576  character(len=:), pointer :: asrt !< pointer to the character string
1577  character(len=*), intent(in) :: name !< variable name
1578  character(len=*), intent(in) :: mem_path !< path where variable is stored
1579  ! -- local
1580  type(memorytype), pointer :: mt
1581  logical(LGP) :: found
1582  ! -- code
1583  call get_from_memorystore(name, mem_path, mt, found)
1584  asrt => mt%strsclr
1585  end subroutine setptr_str
1586 
1587  !> @brief Set pointer to a fixed-length string array
1588  !<
1589  subroutine setptr_str1d(astr1d, name, mem_path)
1590  character(len=:), dimension(:), &
1591  pointer, contiguous, intent(inout) :: astr1d !< pointer to the string array
1592  character(len=*), intent(in) :: name !< variable name
1593  character(len=*), intent(in) :: mem_path !< path where variable is stored
1594  ! -- local
1595  type(memorytype), pointer :: mt
1596  logical(LGP) :: found
1597  ! -- code
1598  call get_from_memorystore(name, mem_path, mt, found)
1599  astr1d => mt%astr1d
1600  end subroutine setptr_str1d
1601 
1602  !> @brief Set pointer to an array of CharacterStringType
1603  !<
1604  subroutine setptr_charstr1d(acharstr1d, name, mem_path)
1605  type(characterstringtype), dimension(:), pointer, contiguous, &
1606  intent(inout) :: acharstr1d !< the reallocated charstring array
1607  character(len=*), intent(in) :: name !< variable name
1608  character(len=*), intent(in) :: mem_path !< path where variable is stored
1609  ! -- local
1610  type(memorytype), pointer :: mt
1611  logical(LGP) :: found
1612  ! -- code
1613  call get_from_memorystore(name, mem_path, mt, found)
1614  acharstr1d => mt%acharstr1d
1615  end subroutine setptr_charstr1d
1616 
1617  !> @brief Make a copy of a 1-dimensional integer array
1618  !<
1619  subroutine copyptr_int1d(aint, name, mem_path, mem_path_copy)
1620  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< returned copy of 1d integer array
1621  character(len=*), intent(in) :: name !< variable name
1622  character(len=*), intent(in) :: mem_path !< path where variable is stored
1623  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1624  !! if passed then the copy is added to the
1625  !! memory manager
1626  ! -- local
1627  type(memorytype), pointer :: mt
1628  logical(LGP) :: found
1629  integer(I4B) :: n
1630  ! -- code
1631  call get_from_memorystore(name, mem_path, mt, found)
1632  aint => null()
1633  ! -- check the copy into the memory manager
1634  if (present(mem_path_copy)) then
1635  call allocate_int1d(aint, size(mt%aint1d), mt%name, mem_path_copy)
1636  ! -- create a local copy
1637  else
1638  allocate (aint(size(mt%aint1d)))
1639  end if
1640  do n = 1, size(mt%aint1d)
1641  aint(n) = mt%aint1d(n)
1642  end do
1643  end subroutine copyptr_int1d
1644 
1645  !> @brief Make a copy of a 2-dimensional integer array
1646  !<
1647  subroutine copyptr_int2d(aint, name, mem_path, mem_path_copy)
1648  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< returned copy of 2d integer array
1649  character(len=*), intent(in) :: name !< variable name
1650  character(len=*), intent(in) :: mem_path !< path where variable is stored
1651  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1652  !! if passed then the copy is added to the
1653  !! memory manager
1654  ! -- local
1655  type(memorytype), pointer :: mt
1656  logical(LGP) :: found
1657  integer(I4B) :: i
1658  integer(I4B) :: j
1659  integer(I4B) :: ncol
1660  integer(I4B) :: nrow
1661  ! -- code
1662  call get_from_memorystore(name, mem_path, mt, found)
1663  aint => null()
1664  ncol = size(mt%aint2d, dim=1)
1665  nrow = size(mt%aint2d, dim=2)
1666  ! -- check the copy into the memory manager
1667  if (present(mem_path_copy)) then
1668  call allocate_int2d(aint, ncol, nrow, mt%name, mem_path_copy)
1669  ! -- create a local copy
1670  else
1671  allocate (aint(ncol, nrow))
1672  end if
1673  do i = 1, nrow
1674  do j = 1, ncol
1675  aint(j, i) = mt%aint2d(j, i)
1676  end do
1677  end do
1678  end subroutine copyptr_int2d
1679 
1680  !> @brief Make a copy of a 1-dimensional real array
1681  !<
1682  subroutine copyptr_dbl1d(adbl, name, mem_path, mem_path_copy)
1683  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< returned copy of 1d real array
1684  character(len=*), intent(in) :: name !< variable name
1685  character(len=*), intent(in) :: mem_path !< path where variable is stored
1686  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1687  !! if passed then the copy is added to the
1688  !! memory manager
1689  ! -- local
1690  type(memorytype), pointer :: mt
1691  logical(LGP) :: found
1692  integer(I4B) :: n
1693  ! -- code
1694  call get_from_memorystore(name, mem_path, mt, found)
1695  adbl => null()
1696  ! -- check the copy into the memory manager
1697  if (present(mem_path_copy)) then
1698  call allocate_dbl1d(adbl, size(mt%adbl1d), mt%name, mem_path_copy)
1699  ! -- create a local copy
1700  else
1701  allocate (adbl(size(mt%adbl1d)))
1702  end if
1703  do n = 1, size(mt%adbl1d)
1704  adbl(n) = mt%adbl1d(n)
1705  end do
1706  end subroutine copyptr_dbl1d
1707 
1708  !> @brief Make a copy of a 2-dimensional real array
1709  !<
1710  subroutine copyptr_dbl2d(adbl, name, mem_path, mem_path_copy)
1711  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< returned copy of 2d real array
1712  character(len=*), intent(in) :: name !< variable name
1713  character(len=*), intent(in) :: mem_path !< path where variable is stored
1714  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1715  !! if passed then the copy is added to the
1716  !! memory manager
1717  ! -- local
1718  type(memorytype), pointer :: mt
1719  logical(LGP) :: found
1720  integer(I4B) :: i
1721  integer(I4B) :: j
1722  integer(I4B) :: ncol
1723  integer(I4B) :: nrow
1724  ! -- code
1725  call get_from_memorystore(name, mem_path, mt, found)
1726  adbl => null()
1727  ncol = size(mt%adbl2d, dim=1)
1728  nrow = size(mt%adbl2d, dim=2)
1729  ! -- check the copy into the memory manager
1730  if (present(mem_path_copy)) then
1731  call allocate_dbl2d(adbl, ncol, nrow, mt%name, mem_path_copy)
1732  ! -- create a local copy
1733  else
1734  allocate (adbl(ncol, nrow))
1735  end if
1736  do i = 1, nrow
1737  do j = 1, ncol
1738  adbl(j, i) = mt%adbl2d(j, i)
1739  end do
1740  end do
1741  end subroutine copyptr_dbl2d
1742 
1743  !> @brief Copy values from a 1-dimensional real array in the memory
1744  !< manager to a passed 1-dimensional real array
1745  subroutine copy_dbl1d(adbl, name, mem_path)
1746  real(dp), dimension(:), intent(inout) :: adbl !< target array
1747  character(len=*), intent(in) :: name !< variable name
1748  character(len=*), intent(in) :: mem_path !< path where variable is stored
1749  ! -- local
1750  type(memorytype), pointer :: mt
1751  logical(LGP) :: found
1752  integer(I4B) :: n
1753  ! -- code
1754  call get_from_memorystore(name, mem_path, mt, found)
1755  do n = 1, size(mt%adbl1d)
1756  adbl(n) = mt%adbl1d(n)
1757  end do
1758  end subroutine copy_dbl1d
1759 
1760  !> @brief Set the pointer for an integer scalar to
1761  !< a target array already stored in the memory manager
1762  subroutine reassignptr_int(sclr, name, mem_path, name_target, mem_path_target)
1763  integer(I4B), pointer, intent(inout) :: sclr !< pointer to integer scalar
1764  character(len=*), intent(in) :: name !< variable name
1765  character(len=*), intent(in) :: mem_path !< path where variable is stored
1766  character(len=*), intent(in) :: name_target !< name of target variable
1767  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1768  ! -- local
1769  type(memorytype), pointer :: mt
1770  type(memorytype), pointer :: mt2
1771  logical(LGP) :: found
1772  ! -- code
1773  call get_from_memorystore(name, mem_path, mt, found)
1774  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1775  if (associated(sclr)) then
1777  deallocate (sclr)
1778  end if
1779  sclr => mt2%intsclr
1780  mt%intsclr => sclr
1781  mt%element_size = i4b
1782  mt%isize = 1
1783  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', mt%isize
1784  !
1785  ! -- set master information
1786  mt%master = .false.
1787  mt%mastername = name_target
1788  mt%masterPath = mem_path_target
1789  end subroutine reassignptr_int
1790 
1791  !> @brief Set the pointer for a 1-dimensional integer array to
1792  !< a target array already stored in the memory manager
1793  subroutine reassignptr_int1d(aint, name, mem_path, name_target, mem_path_target)
1794  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< pointer to 1d integer array
1795  character(len=*), intent(in) :: name !< variable name
1796  character(len=*), intent(in) :: mem_path !< path where variable is stored
1797  character(len=*), intent(in) :: name_target !< name of target variable
1798  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1799  ! -- local
1800  type(memorytype), pointer :: mt
1801  type(memorytype), pointer :: mt2
1802  logical(LGP) :: found
1803  ! -- code
1804  call get_from_memorystore(name, mem_path, mt, found)
1805  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1806  if (size(aint) > 0) then
1807  nvalues_aint = nvalues_aint - size(aint)
1808  deallocate (aint)
1809  end if
1810  aint => mt2%aint1d
1811  mt%aint1d => aint
1812  mt%element_size = i4b
1813  mt%isize = size(aint)
1814  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', mt%isize
1815  !
1816  ! -- set master information
1817  mt%master = .false.
1818  mt%mastername = name_target
1819  mt%masterPath = mem_path_target
1820  end subroutine reassignptr_int1d
1821 
1822  !> @brief Set the pointer for a 2-dimensional integer array to
1823  !< a target array already stored in the memory manager
1824  subroutine reassignptr_int2d(aint, name, mem_path, name_target, mem_path_target)
1825  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< pointer to 2d integer array
1826  character(len=*), intent(in) :: name !< variable name
1827  character(len=*), intent(in) :: mem_path !< path where variable is stored
1828  character(len=*), intent(in) :: name_target !< name of target variable
1829  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1830  ! -- local
1831  type(memorytype), pointer :: mt
1832  type(memorytype), pointer :: mt2
1833  logical(LGP) :: found
1834  integer(I4B) :: ncol
1835  integer(I4B) :: nrow
1836  ! -- code
1837  call get_from_memorystore(name, mem_path, mt, found)
1838  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1839  if (size(aint) > 0) then
1840  nvalues_aint = nvalues_aint - size(aint)
1841  deallocate (aint)
1842  end if
1843  aint => mt2%aint2d
1844  mt%aint2d => aint
1845  mt%element_size = i4b
1846  mt%isize = size(aint)
1847  ncol = size(aint, dim=1)
1848  nrow = size(aint, dim=2)
1849  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
1850  !
1851  ! -- set master information
1852  mt%master = .false.
1853  mt%mastername = name_target
1854  mt%masterPath = mem_path_target
1855  end subroutine reassignptr_int2d
1856 
1857  !> @brief Set the pointer for a 1-dimensional real array to
1858  !< a target array already stored in the memory manager
1859  subroutine reassignptr_dbl1d(adbl, name, mem_path, name_target, mem_path_target)
1860  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< pointer to 1d real array
1861  character(len=*), intent(in) :: name !< variable name
1862  character(len=*), intent(in) :: mem_path !< path where variable is stored
1863  character(len=*), intent(in) :: name_target !< name of target variable
1864  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1865  ! -- local
1866  type(memorytype), pointer :: mt
1867  type(memorytype), pointer :: mt2
1868  logical(LGP) :: found
1869  ! -- code
1870  call get_from_memorystore(name, mem_path, mt, found)
1871  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1872  if (size(adbl) > 0) then
1873  nvalues_adbl = nvalues_adbl - size(adbl)
1874  deallocate (adbl)
1875  end if
1876  adbl => mt2%adbl1d
1877  mt%adbl1d => adbl
1878  mt%element_size = dp
1879  mt%isize = size(adbl)
1880  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', mt%isize
1881  !
1882  ! -- set master information
1883  mt%master = .false.
1884  mt%mastername = name_target
1885  mt%masterPath = mem_path_target
1886  end subroutine reassignptr_dbl1d
1887 
1888  !> @brief Set the pointer for a 2-dimensional real array to
1889  !< a target array already stored in the memory manager
1890  subroutine reassignptr_dbl2d(adbl, name, mem_path, name_target, mem_path_target)
1891  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 2d real array
1892  character(len=*), intent(in) :: name !< variable name
1893  character(len=*), intent(in) :: mem_path !< path where variable is stored
1894  character(len=*), intent(in) :: name_target !< name of target variable
1895  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1896  ! -- local
1897  type(memorytype), pointer :: mt
1898  type(memorytype), pointer :: mt2
1899  logical(LGP) :: found
1900  integer(I4B) :: ncol
1901  integer(I4b) :: nrow
1902  ! -- code
1903  call get_from_memorystore(name, mem_path, mt, found)
1904  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1905  if (size(adbl) > 0) then
1906  nvalues_adbl = nvalues_adbl - size(adbl)
1907  deallocate (adbl)
1908  end if
1909  adbl => mt2%adbl2d
1910  mt%adbl2d => adbl
1911  mt%element_size = dp
1912  mt%isize = size(adbl)
1913  ncol = size(adbl, dim=1)
1914  nrow = size(adbl, dim=2)
1915  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
1916  !
1917  ! -- set master information
1918  mt%master = .false.
1919  mt%mastername = name_target
1920  mt%masterPath = mem_path_target
1921  end subroutine reassignptr_dbl2d
1922 
1923  !> @brief Deallocate a variable-length character string
1924  !<
1925  subroutine deallocate_str(sclr, name, mem_path)
1926  character(len=*), pointer, intent(inout) :: sclr !< pointer to string
1927  character(len=*), intent(in), optional :: name !< variable name
1928  character(len=*), intent(in), optional :: mem_path !< path where variable is stored
1929  ! -- local
1930  type(memorytype), pointer :: mt
1931  logical(LGP) :: found
1932  type(memorycontaineriteratortype), allocatable :: itr
1933  ! -- code
1934  found = .false.
1935  if (present(name) .and. present(mem_path)) then
1936  call get_from_memorystore(name, mem_path, mt, found)
1937  nullify (mt%strsclr)
1938  else
1939  itr = memorystore%iterator()
1940  do while (itr%has_next())
1941  call itr%next()
1942  mt => itr%value()
1943  if (associated(mt%strsclr, sclr)) then
1944  nullify (mt%strsclr)
1945  found = .true.
1946  exit
1947  end if
1948  end do
1949  end if
1950  if (.not. found) then
1951  call store_error('Programming error in deallocate_str.', terminate=.true.)
1952  else
1953  if (mt%master) then
1954  deallocate (sclr)
1955  else
1956  nullify (sclr)
1957  end if
1958  end if
1959  end subroutine deallocate_str
1960 
1961  !> @brief Deallocate an array of defined-length character strings
1962  !!
1963  !<
1964  subroutine deallocate_str1d(astr1d, name, mem_path)
1965  character(len=*), dimension(:), pointer, contiguous, intent(inout) :: astr1d !< array of strings
1966  character(len=*), optional, intent(in) :: name !< variable name
1967  character(len=*), optional, intent(in) :: mem_path !< path where variable is stored
1968  ! -- local
1969  type(memorytype), pointer :: mt
1970  logical(LGP) :: found
1971  type(memorycontaineriteratortype), allocatable :: itr
1972  ! -- code
1973  !
1974  ! -- process optional variables
1975  found = .false.
1976  if (present(name) .and. present(mem_path)) then
1977  call get_from_memorystore(name, mem_path, mt, found)
1978  nullify (mt%astr1d)
1979  else
1980  itr = memorystore%iterator()
1981  do while (itr%has_next())
1982  call itr%next()
1983  mt => itr%value()
1984  if (associated(mt%astr1d, astr1d)) then
1985  nullify (mt%astr1d)
1986  found = .true.
1987  exit
1988  end if
1989  end do
1990  end if
1991  if (.not. found .and. size(astr1d) > 0) then
1992  call store_error('programming error in deallocate_str1d', terminate=.true.)
1993  else
1994  if (mt%master) then
1995  deallocate (astr1d)
1996  else
1997  nullify (astr1d)
1998  end if
1999  end if
2000  end subroutine deallocate_str1d
2001 
2002  !> @brief Deallocate an array of deferred-length character strings
2003  !!
2004  !<
2005  subroutine deallocate_charstr1d(astr1d, name, mem_path)
2006  type(characterstringtype), dimension(:), pointer, contiguous, &
2007  intent(inout) :: astr1d !< array of strings
2008  character(len=*), optional, intent(in) :: name !< variable name
2009  character(len=*), optional, intent(in) :: mem_path !< path where variable is stored
2010  ! -- local
2011  type(memorytype), pointer :: mt
2012  logical(LGP) :: found
2013  type(memorycontaineriteratortype), allocatable :: itr
2014  integer(I4B) :: n
2015  ! -- code
2016  !
2017  ! -- process optional variables
2018  found = .false.
2019  if (present(name) .and. present(mem_path)) then
2020  call get_from_memorystore(name, mem_path, mt, found)
2021  nullify (mt%acharstr1d)
2022  else
2023  itr = memorystore%iterator()
2024  do while (itr%has_next())
2025  call itr%next()
2026  mt => itr%value()
2027  if (associated(mt%acharstr1d, astr1d)) then
2028  nullify (mt%acharstr1d)
2029  found = .true.
2030  exit
2031  end if
2032  end do
2033  end if
2034  if (.not. found .and. size(astr1d) > 0) then
2035  call store_error('programming error in deallocate_charstr1d', &
2036  terminate=.true.)
2037  else
2038  if (mt%master) then
2039  do n = 1, size(astr1d)
2040  call astr1d(n)%destroy()
2041  end do
2042  deallocate (astr1d)
2043  else
2044  nullify (astr1d)
2045  end if
2046  end if
2047  end subroutine deallocate_charstr1d
2048 
2049  !> @brief Deallocate a logical scalar
2050  !<
2051  subroutine deallocate_logical(sclr)
2052  logical(LGP), pointer, intent(inout) :: sclr !< logical scalar to deallocate
2053  ! -- local
2054  class(memorytype), pointer :: mt
2055  logical(LGP) :: found
2056  type(memorycontaineriteratortype), allocatable :: itr
2057  ! -- code
2058  found = .false.
2059  itr = memorystore%iterator()
2060  do while (itr%has_next())
2061  call itr%next()
2062  mt => itr%value()
2063  if (associated(mt%logicalsclr, sclr)) then
2064  nullify (mt%logicalsclr)
2065  found = .true.
2066  exit
2067  end if
2068  end do
2069  if (.not. found) then
2070  call store_error('programming error in deallocate_logical', &
2071  terminate=.true.)
2072  else
2073  if (mt%master) then
2074  deallocate (sclr)
2075  else
2076  nullify (sclr)
2077  end if
2078  end if
2079  end subroutine deallocate_logical
2080 
2081  !> @brief Deallocate a integer scalar
2082  !<
2083  subroutine deallocate_int(sclr)
2084  integer(I4B), pointer, intent(inout) :: sclr !< integer variable to deallocate
2085  ! -- local
2086  class(memorytype), pointer :: mt
2087  logical(LGP) :: found
2088  type(memorycontaineriteratortype), allocatable :: itr
2089  ! -- code
2090  found = .false.
2091  itr = memorystore%iterator()
2092  do while (itr%has_next())
2093  call itr%next()
2094  mt => itr%value()
2095  if (associated(mt%intsclr, sclr)) then
2096  nullify (mt%intsclr)
2097  found = .true.
2098  exit
2099  end if
2100  end do
2101  if (.not. found) then
2102  call store_error('Programming error in deallocate_int.', terminate=.true.)
2103  else
2104  if (mt%master) then
2105  deallocate (sclr)
2106  else
2107  nullify (sclr)
2108  end if
2109  end if
2110  end subroutine deallocate_int
2111 
2112  !> @brief Deallocate a real scalar
2113  !<
2114  subroutine deallocate_dbl(sclr)
2115  real(DP), pointer, intent(inout) :: sclr !< real variable to deallocate
2116  ! -- local
2117  class(memorytype), pointer :: mt
2118  logical(LGP) :: found
2119  type(memorycontaineriteratortype), allocatable :: itr
2120  ! -- code
2121  found = .false.
2122  itr = memorystore%iterator()
2123  do while (itr%has_next())
2124  call itr%next()
2125  mt => itr%value()
2126  if (associated(mt%dblsclr, sclr)) then
2127  nullify (mt%dblsclr)
2128  found = .true.
2129  exit
2130  end if
2131  end do
2132  if (.not. found) then
2133  call store_error('Programming error in deallocate_dbl.', terminate=.true.)
2134  else
2135  if (mt%master) then
2136  deallocate (sclr)
2137  else
2138  nullify (sclr)
2139  end if
2140  end if
2141  end subroutine deallocate_dbl
2142 
2143  !> @brief Deallocate a 1-dimensional integer array
2144  !<
2145  subroutine deallocate_int1d(aint, name, mem_path)
2146  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< 1d integer array to deallocate
2147  character(len=*), optional :: name !< variable name
2148  character(len=*), optional :: mem_path !< path where variable is stored
2149  ! -- local
2150  type(memorytype), pointer :: mt
2151  logical(LGP) :: found
2152  type(memorycontaineriteratortype), allocatable :: itr
2153  ! -- code
2154  !
2155  ! -- process optional variables
2156  found = .false.
2157  if (present(name) .and. present(mem_path)) then
2158  call get_from_memorystore(name, mem_path, mt, found)
2159  nullify (mt%aint1d)
2160  else
2161  itr = memorystore%iterator()
2162  do while (itr%has_next())
2163  call itr%next()
2164  mt => itr%value()
2165  if (associated(mt%aint1d, aint)) then
2166  nullify (mt%aint1d)
2167  found = .true.
2168  exit
2169  end if
2170  end do
2171  end if
2172  if (.not. found .and. size(aint) > 0) then
2173  call store_error('programming error in deallocate_int1d', terminate=.true.)
2174  else
2175  if (mt%master) then
2176  deallocate (aint)
2177  else
2178  nullify (aint)
2179  end if
2180  end if
2181  end subroutine deallocate_int1d
2182 
2183  !> @brief Deallocate a 2-dimensional integer array
2184  !<
2185  subroutine deallocate_int2d(aint, name, mem_path)
2186  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< 2d integer array to deallocate
2187  character(len=*), optional :: name !< variable name
2188  character(len=*), optional :: mem_path !< path where variable is stored
2189  ! -- local
2190  type(memorytype), pointer :: mt
2191  logical(LGP) :: found
2192  type(memorycontaineriteratortype), allocatable :: itr
2193  ! -- code
2194  !
2195  ! -- process optional variables
2196  found = .false.
2197  if (present(name) .and. present(mem_path)) then
2198  call get_from_memorystore(name, mem_path, mt, found)
2199  nullify (mt%aint2d)
2200  else
2201  itr = memorystore%iterator()
2202  do while (itr%has_next())
2203  call itr%next()
2204  mt => itr%value()
2205  if (associated(mt%aint2d, aint)) then
2206  nullify (mt%aint2d)
2207  found = .true.
2208  exit
2209  end if
2210  end do
2211  end if
2212  if (.not. found .and. size(aint) > 0) then
2213  call store_error('programming error in deallocate_int2d', terminate=.true.)
2214  else
2215  if (mt%master) then
2216  deallocate (aint)
2217  else
2218  nullify (aint)
2219  end if
2220  end if
2221  end subroutine deallocate_int2d
2222 
2223  !> @brief Deallocate a 3-dimensional integer array
2224  !<
2225  subroutine deallocate_int3d(aint, name, mem_path)
2226  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< 3d integer array to deallocate
2227  character(len=*), optional :: name !< variable name
2228  character(len=*), optional :: mem_path !< path where variable is stored
2229  ! -- local
2230  type(memorytype), pointer :: mt
2231  logical(LGP) :: found
2232  type(memorycontaineriteratortype), allocatable :: itr
2233  ! -- code
2234  !
2235  ! -- process optional variables
2236  found = .false.
2237  if (present(name) .and. present(mem_path)) then
2238  call get_from_memorystore(name, mem_path, mt, found)
2239  nullify (mt%aint3d)
2240  else
2241  itr = memorystore%iterator()
2242  do while (itr%has_next())
2243  call itr%next()
2244  mt => itr%value()
2245  if (associated(mt%aint3d, aint)) then
2246  nullify (mt%aint3d)
2247  found = .true.
2248  exit
2249  end if
2250  end do
2251  end if
2252  if (.not. found .and. size(aint) > 0) then
2253  call store_error('programming error in deallocate_int3d', terminate=.true.)
2254  else
2255  if (mt%master) then
2256  deallocate (aint)
2257  else
2258  nullify (aint)
2259  end if
2260  end if
2261  end subroutine deallocate_int3d
2262 
2263  !> @brief Deallocate a 1-dimensional real array
2264  !<
2265  subroutine deallocate_dbl1d(adbl, name, mem_path)
2266  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< 1d real array to deallocate
2267  character(len=*), optional :: name !< variable name
2268  character(len=*), optional :: mem_path !< path where variable is stored
2269  ! -- local
2270  type(memorytype), pointer :: mt
2271  logical(LGP) :: found
2272  type(memorycontaineriteratortype), allocatable :: itr
2273  ! -- code
2274  !
2275  ! -- process optional variables
2276  found = .false.
2277  if (present(name) .and. present(mem_path)) then
2278  call get_from_memorystore(name, mem_path, mt, found)
2279  nullify (mt%adbl1d)
2280  else
2281  itr = memorystore%iterator()
2282  do while (itr%has_next())
2283  call itr%next()
2284  mt => itr%value()
2285  if (associated(mt%adbl1d, adbl)) then
2286  nullify (mt%adbl1d)
2287  found = .true.
2288  exit
2289  end if
2290  end do
2291  end if
2292  if (.not. found .and. size(adbl) > 0) then
2293  call store_error('programming error in deallocate_dbl1d', terminate=.true.)
2294  else
2295  if (mt%master) then
2296  deallocate (adbl)
2297  else
2298  nullify (adbl)
2299  end if
2300  end if
2301  end subroutine deallocate_dbl1d
2302 
2303  !> @brief Deallocate a 2-dimensional real array
2304  !<
2305  subroutine deallocate_dbl2d(adbl, name, mem_path)
2306  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< 2d real array to deallocate
2307  character(len=*), optional :: name !< variable name
2308  character(len=*), optional :: mem_path !< path where variable is stored
2309  ! -- local
2310  type(memorytype), pointer :: mt
2311  logical(LGP) :: found
2312  type(memorycontaineriteratortype), allocatable :: itr
2313  ! -- code
2314  !
2315  ! -- process optional variables
2316  found = .false.
2317  if (present(name) .and. present(mem_path)) then
2318  call get_from_memorystore(name, mem_path, mt, found)
2319  nullify (mt%adbl2d)
2320  else
2321  itr = memorystore%iterator()
2322  do while (itr%has_next())
2323  call itr%next()
2324  mt => itr%value()
2325  if (associated(mt%adbl2d, adbl)) then
2326  nullify (mt%adbl2d)
2327  found = .true.
2328  exit
2329  end if
2330  end do
2331  end if
2332  if (.not. found .and. size(adbl) > 0) then
2333  call store_error('programming error in deallocate_dbl2d', terminate=.true.)
2334  else
2335  if (mt%master) then
2336  deallocate (adbl)
2337  else
2338  nullify (adbl)
2339  end if
2340  end if
2341  end subroutine deallocate_dbl2d
2342 
2343  !> @brief Deallocate a 3-dimensional real array
2344  !<
2345  subroutine deallocate_dbl3d(adbl, name, mem_path)
2346  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< 3d real array to deallocate
2347  character(len=*), optional :: name !< variable name
2348  character(len=*), optional :: mem_path !< path where variable is stored
2349  ! -- local
2350  type(memorytype), pointer :: mt
2351  logical(LGP) :: found
2352  type(memorycontaineriteratortype), allocatable :: itr
2353  ! -- code
2354  !
2355  ! -- process optional variables
2356  found = .false.
2357  if (present(name) .and. present(mem_path)) then
2358  call get_from_memorystore(name, mem_path, mt, found)
2359  nullify (mt%adbl3d)
2360  else
2361  itr = memorystore%iterator()
2362  do while (itr%has_next())
2363  call itr%next()
2364  mt => itr%value()
2365  if (associated(mt%adbl3d, adbl)) then
2366  nullify (mt%adbl3d)
2367  found = .true.
2368  exit
2369  end if
2370  end do
2371  end if
2372  if (.not. found .and. size(adbl) > 0) then
2373  call store_error('programming error in deallocate_dbl3d', terminate=.true.)
2374  else
2375  if (mt%master) then
2376  deallocate (adbl)
2377  else
2378  nullify (adbl)
2379  end if
2380  end if
2381  end subroutine deallocate_dbl3d
2382 
2383  !> @brief Set the memory print option
2384  !<
2385  subroutine mem_set_print_option(iout, keyword, error_msg)
2386  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2387  character(len=*), intent(in) :: keyword !< memory print option
2388  character(len=*), intent(inout) :: error_msg !< returned error message if keyword is not valid option
2389  ! -- local
2390  ! -- format
2391  ! -- code
2392  select case (keyword)
2393  case ('NONE')
2394  iprmem = 0
2395  write (iout, '(4x, a)') &
2396  'LIMITED MEMORY INFORMATION WILL BE WRITTEN.'
2397  case ('SUMMARY')
2398  iprmem = 1
2399  write (iout, '(4x, a)') &
2400  'A SUMMARY OF SIMULATION MEMORY INFORMATION WILL BE WRITTEN.'
2401  case ('ALL')
2402  iprmem = 2
2403  write (iout, '(4x, a)') &
2404  'ALL SIMULATION MEMORY INFORMATION WILL BE WRITTEN.'
2405  case default
2406  error_msg = "Unknown memory print option '"//trim(keyword)//"."
2407  end select
2408  end subroutine mem_set_print_option
2409 
2410  !> @brief Create a table if memory_print_option is 'SUMMARY'
2411  !<
2412  subroutine mem_summary_table(iout, nrows, cunits)
2413  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2414  integer(I4B), intent(in) :: nrows !< number of table rows
2415  character(len=*), intent(in) :: cunits !< memory units (bytes, kilobytes, megabytes, or gigabytes)
2416  ! -- local
2417  character(len=LINELENGTH) :: title
2418  character(len=LINELENGTH) :: text
2419  integer(I4B) :: nterms
2420  ! -- formats
2421  ! -- code
2422  nterms = 6
2423  !
2424  ! -- set up table title
2425  title = 'SUMMARY INFORMATION ON VARIABLES STORED IN THE MEMORY MANAGER, '// &
2426  'IN '//trim(cunits)
2427  !
2428  ! -- set up stage tableobj
2429  call table_cr(memtab, 'MEM SUM', title)
2430  call memtab%table_df(nrows, nterms, iout)
2431  !
2432  ! -- data type
2433  text = 'COMPONENT'
2434  call memtab%initialize_column(text, 20, alignment=tableft)
2435  !
2436  ! -- memory allocated for characters
2437  text = 'CHARACTER'
2438  call memtab%initialize_column(text, 15, alignment=tabcenter)
2439  !
2440  ! -- memory allocated for logical
2441  text = 'LOGICAL'
2442  call memtab%initialize_column(text, 15, alignment=tabcenter)
2443  !
2444  ! -- memory allocated for integers
2445  text = 'INTEGER'
2446  call memtab%initialize_column(text, 15, alignment=tabcenter)
2447  !
2448  ! -- memory allocated for reals
2449  text = 'REAL'
2450  call memtab%initialize_column(text, 15, alignment=tabcenter)
2451  !
2452  ! -- total memory allocated
2453  text = 'TOTAL'
2454  call memtab%initialize_column(text, 15, alignment=tabcenter)
2455  end subroutine mem_summary_table
2456 
2457  !> @brief Create a table if memory_print_option is 'ALL'
2458  !<
2459  subroutine mem_detailed_table(iout, nrows)
2460  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2461  integer(I4B), intent(in) :: nrows !< number of table rows
2462  ! -- local
2463  character(len=LINELENGTH) :: title
2464  character(len=LINELENGTH) :: text
2465  integer(I4B) :: nterms
2466  ! -- formats
2467  ! -- code
2468  nterms = 5
2469  !
2470  ! -- set up table title
2471  title = 'DETAILED INFORMATION ON VARIABLES STORED IN THE MEMORY MANAGER'
2472  !
2473  ! -- set up stage tableobj
2474  call table_cr(memtab, 'MEM DET', title)
2475  call memtab%table_df(nrows, nterms, iout)
2476  !
2477  ! -- origin
2478  text = 'ORIGIN'
2479  call memtab%initialize_column(text, lenmempath, alignment=tableft)
2480  !
2481  ! -- variable
2482  text = 'VARIABLE NAME'
2483  call memtab%initialize_column(text, lenvarname, alignment=tableft)
2484  !
2485  ! -- data type
2486  text = 'DATA TYPE'
2487  call memtab%initialize_column(text, 16, alignment=tableft)
2488  !
2489  ! -- size
2490  text = 'NUMBER OF ITEMS'
2491  call memtab%initialize_column(text, 20, alignment=tabright)
2492  !
2493  ! -- is it a pointer
2494  text = 'ASSOCIATED VARIABLE'
2495  call memtab%initialize_column(text, lenmemaddress, alignment=tableft)
2496  end subroutine mem_detailed_table
2497 
2498  !> @brief Write a row for the memory_print_option 'SUMMARY' table
2499  !<
2500  subroutine mem_summary_line(component, rchars, rlog, rint, rreal, bytes)
2501  character(len=*), intent(in) :: component !< character defining the program component (e.g. solution)
2502  real(DP), intent(in) :: rchars !< allocated size of characters (in common units)
2503  real(DP), intent(in) :: rlog !< allocated size of logical (in common units)
2504  real(DP), intent(in) :: rint !< allocated size of integer variables (in common units)
2505  real(DP), intent(in) :: rreal !< allocated size of real variables (in common units)
2506  real(DP), intent(in) :: bytes !< total allocated memory in memory manager (in common units)
2507  ! -- formats
2508  ! -- code
2509  !
2510  ! -- write line
2511  call memtab%add_term(component)
2512  call memtab%add_term(rchars)
2513  call memtab%add_term(rlog)
2514  call memtab%add_term(rint)
2515  call memtab%add_term(rreal)
2516  call memtab%add_term(bytes)
2517  end subroutine mem_summary_line
2518 
2519  !> @brief Determine appropriate memory unit and conversion factor
2520  !<
2521  subroutine mem_units(bytes, fact, cunits)
2522  ! -- dummy
2523  real(DP), intent(in) :: bytes !< total nr. of bytes
2524  real(DP), intent(inout) :: fact !< conversion factor
2525  character(len=*), intent(inout) :: cunits !< string with memory unit
2526  ! -- local
2527  ! -- formats
2528  ! -- code
2529  !
2530  ! -- initialize factor and unit string
2531  cunits = 'UNKNOWN'
2532  fact = done
2533  !
2534  ! -- set factor
2535  if (bytes < dep3) then
2536  fact = done
2537  cunits = 'BYTES'
2538  else if (bytes < dep6) then
2539  fact = dem3
2540  cunits = 'KILOBYTES'
2541  else if (bytes < dep9) then
2542  fact = dem6
2543  cunits = 'MEGABYTES'
2544  else
2545  fact = dem9
2546  cunits = 'GIGABYTES'
2547  end if
2548  end subroutine mem_units
2549 
2550  !> @brief Create and fill a table with the total allocated memory
2551  !< in the memory manager
2552  subroutine mem_summary_total(iout, bytes)
2553  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2554  real(DP), intent(in) :: bytes !< total number of bytes allocated in the memory manager
2555  ! -- local
2556  character(len=LINELENGTH) :: title
2557  character(len=LINELENGTH) :: text
2558  character(LEN=10) :: cunits
2559  integer(I4B) :: nterms
2560  integer(I4B) :: nrows
2561  real(DP) :: fact
2562  real(DP) :: smb
2563  ! -- formats
2564  ! -- code
2565  !
2566  ! -- calculate factor and memory units
2567  call mem_units(bytes, fact, cunits)
2568  !
2569  ! -- set table terms
2570  nterms = 2
2571  nrows = 6
2572  !
2573  ! -- set up table title
2574  title = 'MEMORY MANAGER TOTAL STORAGE BY DATA TYPE, IN '//trim(cunits)
2575  !
2576  ! -- set up stage tableobj
2577  call table_cr(memtab, 'MEM TOT', title)
2578  call memtab%table_df(nrows, nterms, iout)
2579  !
2580  ! -- data type
2581  text = 'DATA TYPE'
2582  call memtab%initialize_column(text, 15, alignment=tableft)
2583  !
2584  ! -- number of values
2585  text = 'ALLOCATED MEMORY'
2586  call memtab%initialize_column(text, 15, alignment=tabcenter)
2587  !
2588  ! -- write data
2589  !
2590  ! -- characters
2591  smb = real(nvalues_astr, dp) * fact
2592  call memtab%add_term('Character')
2593  call memtab%add_term(smb)
2594  !
2595  ! -- logicals
2596  smb = real(nvalues_alogical * lgp, dp) * fact
2597  call memtab%add_term('Logical')
2598  call memtab%add_term(smb)
2599  !
2600  ! -- integers
2601  smb = real(nvalues_aint * i4b, dp) * fact
2602  call memtab%add_term('Integer')
2603  call memtab%add_term(smb)
2604  !
2605  ! -- reals
2606  smb = real(nvalues_adbl * dp, dp) * fact
2607  call memtab%add_term('Real')
2608  call memtab%add_term(smb)
2609  !
2610  ! -- total memory usage
2611  call memtab%print_separator()
2612  smb = bytes * fact
2613  call memtab%add_term('Total')
2614  call memtab%add_term(smb)
2615  !
2616  ! -- Virtual memory
2617  smb = calc_virtual_mem() * fact
2618  call memtab%add_term('Virtual')
2619  call memtab%add_term(smb)
2620  !
2621  ! -- deallocate table
2622  call mem_cleanup_table()
2623  end subroutine mem_summary_total
2624 
2625  !> @brief Generic function to clean a memory manager table
2626  !<
2627  subroutine mem_cleanup_table()
2628  ! -- local
2629  ! -- formats
2630  ! -- code
2631  call memtab%table_da()
2632  deallocate (memtab)
2633  nullify (memtab)
2634  end subroutine mem_cleanup_table
2635 
2636  !> @brief Write memory manager memory usage based on the
2637  !! user-specified memory_print_option
2638  !!
2639  !! The total memory usage by data types (int, real, etc.)
2640  !! is written for every simulation.
2641  !<
2642  subroutine mem_write_usage(iout)
2643  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2644  ! -- local
2645  class(memorytype), pointer :: mt
2646  character(len=LENMEMADDRESS), allocatable, dimension(:) :: cunique
2647  ! character(len=LENMEMPATH) :: mem_path
2648  character(len=LENMEMPATH) :: context
2649  character(len=LENCOMPONENTNAME) :: component
2650  character(len=LENCOMPONENTNAME) :: subcomponent
2651  character(len=LENMEMADDRESS) :: context_component
2652  character(LEN=10) :: cunits
2653  type(memorycontaineriteratortype), allocatable :: itr
2654  integer(I4B) :: icomp
2655  integer(I4B) :: ilen
2656  integer(I8B) :: nchars
2657  integer(I8B) :: nlog
2658  integer(I8B) :: nint
2659  integer(I8B) :: nreal
2660  real(dp) :: simbytes
2661  real(dp) :: fact
2662  real(dp) :: rchars
2663  real(dp) :: rlog
2664  real(dp) :: rint
2665  real(dp) :: rreal
2666  real(dp) :: bytes
2667  ! -- formats
2668  ! -- code
2669  !
2670  ! -- Calculate simulation memory allocation
2671  simbytes = (nvalues_astr + &
2672  nvalues_alogical * lgp + &
2673  nvalues_aint * i4b + &
2674  nvalues_adbl * dp)
2675  simbytes = real(simbytes, dp)
2676  !
2677  ! -- calculate factor and memory units
2678  call mem_units(simbytes, fact, cunits)
2679  !
2680  ! -- Write summary table for simulation components
2681  if (iprmem == 1) then
2682  !
2683  ! -- Find unique names of simulation components
2684  call mem_unique_origins(cunique)
2685  call mem_summary_table(iout, size(cunique), cunits)
2686  do icomp = 1, size(cunique)
2687  nchars = 0
2688  nlog = 0
2689  nint = 0
2690  nreal = 0
2691  bytes = dzero
2692  ilen = len_trim(cunique(icomp))
2693  itr = memorystore%iterator()
2694  do while (itr%has_next())
2695  call itr%next()
2696  mt => itr%value()
2697  call split_mem_path(mt%path, component, subcomponent)
2698  context = get_mem_path_context(mt%path)
2699  context_component = trim(context)//component
2700  if (cunique(icomp) /= context_component(1:ilen)) cycle
2701  if (.not. mt%master) cycle
2702  if (mt%memtype(1:6) == 'STRING') then
2703  nchars = nchars + mt%isize * mt%element_size
2704  else if (mt%memtype(1:7) == 'LOGICAL') then
2705  nlog = nlog + mt%isize
2706  else if (mt%memtype(1:7) == 'INTEGER') then
2707  nint = nint + mt%isize
2708  else if (mt%memtype(1:6) == 'DOUBLE') then
2709  nreal = nreal + mt%isize
2710  end if
2711  end do
2712  !
2713  ! -- calculate size of each data type in bytes
2714  rchars = real(nchars, dp) * fact
2715  rlog = real(nlog * lgp, dp) * fact
2716  rint = real(nint * i4b, dp) * fact
2717  rreal = real(nreal * dp, dp) * fact
2718  !
2719  ! -- calculate total storage in bytes
2720  bytes = rchars + rlog + rint + rreal
2721  !
2722  ! -- write data
2723  call mem_summary_line(cunique(icomp), rchars, rlog, rint, rreal, bytes)
2724  end do
2725  call mem_cleanup_table()
2726  end if
2727  !
2728  ! -- Write table with all variables for iprmem == 2
2729  if (iprmem == 2) then
2730  call mem_print_detailed(iout)
2731  end if
2732  !
2733  ! -- Write total memory allocation
2734  call mem_summary_total(iout, simbytes)
2735  end subroutine mem_write_usage
2736 
2737  subroutine mem_print_detailed(iout)
2738  integer(I4B) :: iout
2739  ! local
2740  class(memorytype), pointer :: mt
2741  type(memorycontaineriteratortype), allocatable :: itr
2742 
2743  call mem_detailed_table(iout, memorystore%count())
2744  itr = memorystore%iterator()
2745  do while (itr%has_next())
2746  call itr%next()
2747  mt => itr%value()
2748  call mt%table_entry(memtab)
2749  end do
2750  call mem_cleanup_table()
2751 
2752  end subroutine mem_print_detailed
2753 
2754  !> @brief Sum up virtual memory, i.e. memory
2755  !< that is owned by other processes
2756  function calc_virtual_mem() result(vmem_size)
2757  real(dp) :: vmem_size
2758  ! local
2759  type(memorycontaineriteratortype), allocatable :: itr
2760  type(memorytype), pointer :: mt
2761 
2762  vmem_size = dzero
2763  itr = memorystore%iterator()
2764  do while (itr%has_next())
2765  call itr%next()
2766  mt => itr%value()
2767  if (index(mt%path, "__P") == 1) then
2768  vmem_size = mt%element_size * mt%isize + vmem_size
2769  end if
2770  end do
2771 
2772  end function calc_virtual_mem
2773 
2774  !> @brief Deallocate memory in the memory manager
2775  !<
2776  subroutine mem_da()
2777  ! -- modules
2778  use versionmodule, only: idevelopmode
2779  use inputoutputmodule, only: upcase
2780  ! -- local
2781  class(memorytype), pointer :: mt
2782  character(len=LINELENGTH) :: error_msg
2783  character(len=LENVARNAME) :: ucname
2784  type(memorycontaineriteratortype), allocatable :: itr
2785  ! -- code
2786  itr = memorystore%iterator()
2787  do while (itr%has_next())
2788  call itr%next()
2789  mt => itr%value()
2790  if (idevelopmode == 1) then
2791  !
2792  ! -- check if memory has been deallocated
2793  if (mt%mt_associated() .and. mt%element_size == -1) then
2794  error_msg = trim(adjustl(mt%path))//' '// &
2795  trim(adjustl(mt%name))//' has invalid element size'
2796  call store_error(trim(error_msg))
2797  end if
2798  !
2799  ! -- check if memory has been deallocated
2800  if (mt%mt_associated() .and. mt%isize > 0) then
2801  error_msg = trim(adjustl(mt%path))//' '// &
2802  trim(adjustl(mt%name))//' not deallocated'
2803  call store_error(trim(error_msg))
2804  end if
2805  !
2806  ! -- check case of varname
2807  ucname = mt%name
2808  call upcase(ucname)
2809  if (mt%name /= ucname) then
2810  error_msg = trim(adjustl(mt%path))//' '// &
2811  trim(adjustl(mt%name))//' not upper case'
2812  call store_error(trim(error_msg))
2813  end if
2814  end if
2815  !
2816  ! -- deallocate instance of memory type
2817  deallocate (mt)
2818  end do
2819  call memorystore%clear()
2820  if (count_errors() > 0) then
2821  call store_error('Could not clear memory list.', terminate=.true.)
2822  end if
2823  end subroutine mem_da
2824 
2825  !> @brief Create a array with unique first components from all memory paths.
2826  !! Only the first component of the memory path is evaluated.
2827  !<
2828  subroutine mem_unique_origins(cunique)
2829  ! -- modules
2831  ! -- dummy
2832  character(len=LENMEMADDRESS), allocatable, dimension(:), intent(inout) :: &
2833  cunique !< array with unique first components
2834  ! -- local
2835  class(memorytype), pointer :: mt
2836  character(len=LENMEMPATH) :: context
2837  character(len=LENCOMPONENTNAME) :: component
2838  character(len=LENCOMPONENTNAME) :: subcomponent
2839  character(len=LENMEMADDRESS) :: context_component
2840  type(memorycontaineriteratortype), allocatable :: itr
2841  integer(I4B) :: ipa
2842  ! -- code
2843  !
2844  ! -- initialize cunique
2845  allocate (cunique(0))
2846  !
2847  ! -- find unique origins
2848  itr = memorystore%iterator()
2849  do while (itr%has_next())
2850  call itr%next()
2851  mt => itr%value()
2852  call split_mem_path(mt%path, component, subcomponent)
2853  context = get_mem_path_context(mt%path)
2854  context_component = trim(context)//component
2855  ipa = ifind(cunique, context_component)
2856  if (ipa < 1) then
2857  call expandarray(cunique, 1)
2858  cunique(size(cunique)) = context_component
2859  end if
2860  end do
2861  end subroutine mem_unique_origins
2862 
2863 end module memorymanagermodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lencomponentname
maximum length of a component name
Definition: Constants.f90:18
@ tabcenter
centered table column
Definition: Constants.f90:172
@ tabright
right justified table column
Definition: Constants.f90:173
@ tableft
left justified table column
Definition: Constants.f90:171
integer(i4b), parameter lenmemaddress
maximum length of the full memory address, including variable name
Definition: Constants.f90:31
real(dp), parameter dep3
real constant 1000
Definition: Constants.f90:88
@ tabucstring
upper case string table data
Definition: Constants.f90:180
@ tabstring
string table data
Definition: Constants.f90:179
@ tabreal
real table data
Definition: Constants.f90:182
@ tabinteger
integer table data
Definition: Constants.f90:181
real(dp), parameter dep6
real constant 1000000
Definition: Constants.f90:89
integer(i4b), parameter lenmemseparator
maximum length of the memory path separator used, currently a '/'
Definition: Constants.f90:26
real(dp), parameter dep9
real constant 1e9
Definition: Constants.f90:90
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
integer(i4b), parameter, public lenmemtype
maximum length of a memory manager type
Definition: Constants.f90:62
real(dp), parameter dem3
real constant 1e-3
Definition: Constants.f90:106
real(dp), parameter dem6
real constant 1e-6
Definition: Constants.f90:109
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dem9
real constant 1e-9
Definition: Constants.f90:112
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine, public upcase(word)
Convert to upper case.
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function get_mem_path_context(mem_path)
Return the context from the memory path.
subroutine mem_check_length(name, max_length, description)
Generic routine to check the length of (parts of) the memory address.
subroutine strip_context_mem_path(mem_path, mem_path_no_context)
Remove the context from the memory path.
subroutine split_mem_path(mem_path, component, subcomponent)
Split the memory path into component(s)
subroutine reallocate_int2d(aint, ncol, nrow, name, mem_path)
Reallocate a 2-dimensional integer array.
subroutine allocate_logical(sclr, name, mem_path)
Allocate a logical scalar.
subroutine allocate_int2d(aint, ncol, nrow, name, mem_path)
Allocate a 2-dimensional integer array.
subroutine mem_summary_line(component, rchars, rlog, rint, rreal, bytes)
Write a row for the memory_print_option 'SUMMARY' table.
subroutine allocate_dbl3d(adbl, ncol, nrow, nlay, name, mem_path)
Allocate a 3-dimensional real array.
integer(i8b) nvalues_adbl
subroutine deallocate_str1d(astr1d, name, mem_path)
Deallocate an array of defined-length character strings.
subroutine deallocate_int(sclr)
Deallocate a integer scalar.
type(tabletype), pointer memtab
subroutine setptr_int3d(aint, name, mem_path)
Set pointer to 3d integer array.
subroutine mem_detailed_table(iout, nrows)
Create a table if memory_print_option is 'ALL'.
subroutine deallocate_charstr1d(astr1d, name, mem_path)
Deallocate an array of deferred-length character strings.
subroutine allocate_charstr1d(acharstr1d, ilen, nrow, name, mem_path)
Allocate a 1-dimensional array of deferred-length CharacterStringType.
subroutine reallocate_charstr1d(acharstr1d, ilen, nrow, name, mem_path)
Reallocate a 1-dimensional deferred length string array.
subroutine, public mem_write_usage(iout)
Write memory manager memory usage based on the user-specified memory_print_option.
subroutine reallocate_dbl2d(adbl, ncol, nrow, name, mem_path)
Reallocate a 2-dimensional real array.
integer(i8b) nvalues_aint
subroutine setptr_int(sclr, name, mem_path)
Set pointer to integer scalar.
subroutine checkin_charstr1d(acharstr1d, ilen, name, mem_path, name2, mem_path2)
Check in an existing 1d CharacterStringType array with a new address (name + path)
subroutine mem_cleanup_table()
Generic function to clean a memory manager table.
subroutine copyptr_int1d(aint, name, mem_path, mem_path_copy)
Make a copy of a 1-dimensional integer array.
subroutine deallocate_dbl1d(adbl, name, mem_path)
Deallocate a 1-dimensional real array.
subroutine allocate_str1d(astr1d, ilen, nrow, name, mem_path)
Allocate a 1-dimensional defined length string array.
subroutine allocate_str(sclr, ilen, name, mem_path)
Allocate a character string.
subroutine, public get_mem_type(name, mem_path, var_type)
@ brief Get the variable memory type
subroutine setptr_int2d(aint, name, mem_path)
Set pointer to 2d integer array.
subroutine, public mem_da()
Deallocate memory in the memory manager.
subroutine deallocate_int1d(aint, name, mem_path)
Deallocate a 1-dimensional integer array.
integer(i8b) nvalues_astr
subroutine deallocate_dbl(sclr)
Deallocate a real scalar.
subroutine copyptr_int2d(aint, name, mem_path, mem_path_copy)
Make a copy of a 2-dimensional integer array.
subroutine allocate_dbl(sclr, name, mem_path)
Allocate a real scalar.
real(dp) function calc_virtual_mem()
Sum up virtual memory, i.e. memory.
subroutine copyptr_dbl2d(adbl, name, mem_path, mem_path_copy)
Make a copy of a 2-dimensional real array.
subroutine reassignptr_int2d(aint, name, mem_path, name_target, mem_path_target)
Set the pointer for a 2-dimensional integer array to.
subroutine setptr_dbl(sclr, name, mem_path)
Set pointer to a real scalar.
subroutine deallocate_int3d(aint, name, mem_path)
Deallocate a 3-dimensional integer array.
subroutine deallocate_int2d(aint, name, mem_path)
Deallocate a 2-dimensional integer array.
subroutine reassignptr_dbl1d(adbl, name, mem_path, name_target, mem_path_target)
Set the pointer for a 1-dimensional real array to.
subroutine allocate_int3d(aint, ncol, nrow, nlay, name, mem_path)
Allocate a 3-dimensional integer array.
subroutine, public get_mem_shape(name, mem_path, mem_shape)
@ brief Get the variable memory shape
subroutine deallocate_dbl2d(adbl, name, mem_path)
Deallocate a 2-dimensional real array.
subroutine mem_units(bytes, fact, cunits)
Determine appropriate memory unit and conversion factor.
subroutine setptr_int1d(aint, name, mem_path)
Set pointer to 1d integer array.
subroutine setptr_str1d(astr1d, name, mem_path)
Set pointer to a fixed-length string array.
subroutine checkin_int2d(aint2d, name, mem_path, name2, mem_path2)
Check in an existing 2d integer array with a new address (name + path)
type(memorystoretype), public memorystore
subroutine reallocate_dbl1d(adbl, nrow, name, mem_path)
Reallocate a 1-dimensional real array.
subroutine setptr_dbl2d(adbl, name, mem_path)
Set pointer to a 2d real array.
subroutine deallocate_logical(sclr)
Deallocate a logical scalar.
subroutine checkin_dbl2d(adbl2d, name, mem_path, name2, mem_path2)
Check in an existing 2d double precision array with a new address (name + path)
subroutine deallocate_dbl3d(adbl, name, mem_path)
Deallocate a 3-dimensional real array.
subroutine reallocate_int1d(aint, nrow, name, mem_path)
Reallocate a 1-dimensional integer array.
subroutine, public get_from_memorystore(name, mem_path, mt, found, check)
@ brief Get a memory type entry from the memory list
subroutine, public mem_print_detailed(iout)
subroutine reallocate_str1d(astr, ilen, nrow, name, mem_path)
Reallocate a 1-dimensional defined length string array.
subroutine copyptr_dbl1d(adbl, name, mem_path, mem_path_copy)
Make a copy of a 1-dimensional real array.
subroutine reassignptr_dbl2d(adbl, name, mem_path, name_target, mem_path_target)
Set the pointer for a 2-dimensional real array to.
subroutine allocate_int1d(aint, nrow, name, mem_path)
Allocate a 1-dimensional integer array.
subroutine allocate_int(sclr, name, mem_path)
Allocate a integer scalar.
subroutine reassignptr_int1d(aint, name, mem_path, name_target, mem_path_target)
Set the pointer for a 1-dimensional integer array to.
subroutine mem_summary_total(iout, bytes)
Create and fill a table with the total allocated memory.
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
subroutine allocate_dbl2d(adbl, ncol, nrow, name, mem_path)
Allocate a 2-dimensional real array.
subroutine deallocate_str(sclr, name, mem_path)
Deallocate a variable-length character string.
subroutine, public get_mem_rank(name, mem_path, rank)
@ brief Get the variable rank
subroutine checkin_dbl1d(adbl, name, mem_path, name2, mem_path2)
Check in an existing 1d double precision array with a new address (name + path)
subroutine checkin_int1d(aint, name, mem_path, name2, mem_path2)
Check in an existing 1d integer array with a new address (name + path)
subroutine reassignptr_int(sclr, name, mem_path, name_target, mem_path_target)
Set the pointer for an integer scalar to.
subroutine, public copy_dbl1d(adbl, name, mem_path)
Copy values from a 1-dimensional real array in the memory.
subroutine setptr_dbl3d(adbl, name, mem_path)
Set pointer to a 3d real array.
subroutine, public get_mem_elem_size(name, mem_path, size)
@ brief Get the memory size of a single element of the stored variable
subroutine setptr_dbl1d(adbl, name, mem_path)
Set pointer to a 1d real array.
subroutine setptr_str(asrt, name, mem_path)
Set pointer to a string (scalar)
integer(i8b) nvalues_alogical
subroutine setptr_charstr1d(acharstr1d, name, mem_path)
Set pointer to an array of CharacterStringType.
subroutine mem_unique_origins(cunique)
Create a array with unique first components from all memory paths. Only the first component of the me...
subroutine mem_summary_table(iout, nrows, cunits)
Create a table if memory_print_option is 'SUMMARY'.
subroutine setptr_logical(sclr, name, mem_path)
Set pointer to a logical scalar.
subroutine, public mem_set_print_option(iout, keyword, error_msg)
Set the memory print option.
subroutine allocate_dbl1d(adbl, nrow, name, mem_path)
Allocate a 1-dimensional real array.
subroutine allocate_error(varname, mem_path, istat, isize)
Issue allocation error message and stop program execution.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
subroutine, public table_cr(this, name, title)
Definition: Table.f90:87
This module contains version information.
Definition: version.f90:7
integer(i4b), parameter idevelopmode
Definition: version.f90:19
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.