MODFLOW 6  version 6.6.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  end do
1227  !
1228  ! -- fill new values with missing values
1229  do n = nrow_old + 1, nrow
1230  astrtemp(n) = string
1231  end do
1232  !
1233  ! -- deallocate mt pointer, repoint, recalculate isize
1234  deallocate (acharstr1d)
1235  !
1236  ! -- allocate astr1d
1237  allocate (acharstr1d(nrow), stat=istat, errmsg=errmsg)
1238  if (istat /= 0) then
1239  call allocate_error(name, mem_path, istat, isize)
1240  end if
1241  !
1242  ! -- fill the reallocated character array
1243  do n = 1, nrow
1244  acharstr1d(n) = astrtemp(n)
1245  end do
1246  !
1247  ! -- deallocate temporary storage
1248  deallocate (astrtemp)
1249  !
1250  ! -- reset memory manager values
1251  mt%acharstr1d => acharstr1d
1252  mt%element_size = ilen
1253  mt%isize = isize
1254  mt%nrealloc = mt%nrealloc + 1
1255  mt%master = .true.
1256  nvalues_astr = nvalues_astr + isize - isize_old
1257  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
1258  else
1259  errmsg = "Programming error, variable '"//trim(name)//"' from '"// &
1260  trim(mem_path)//"' is not defined in the memory manager. Use "// &
1261  "mem_allocate instead."
1262  call store_error(errmsg, terminate=.true.)
1263  end if
1264  end subroutine reallocate_charstr1d
1265 
1266  !> @brief Reallocate a 1-dimensional integer array
1267  !<
1268  subroutine reallocate_int1d(aint, nrow, name, mem_path)
1269  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< the reallocated integer array
1270  integer(I4B), intent(in) :: nrow !< number of rows
1271  character(len=*), intent(in) :: name !< variable name
1272  character(len=*), intent(in) :: mem_path !< path where variable is stored
1273  ! -- local
1274  type(memorytype), pointer :: mt
1275  logical(LGP) :: found
1276  integer(I4B) :: istat
1277  integer(I4B) :: isize
1278  integer(I4B) :: i
1279  integer(I4B) :: isizeold
1280  integer(I4B) :: ifill
1281  ! -- code
1282  !
1283  ! -- Find and assign mt
1284  call get_from_memorystore(name, mem_path, mt, found)
1285  !
1286  ! -- Allocate aint and then refill
1287  isize = nrow
1288  isizeold = size(mt%aint1d)
1289  ifill = min(isizeold, isize)
1290  allocate (aint(nrow), stat=istat, errmsg=errmsg)
1291  if (istat /= 0) then
1292  call allocate_error(name, mem_path, istat, isize)
1293  end if
1294  do i = 1, ifill
1295  aint(i) = mt%aint1d(i)
1296  end do
1297  !
1298  ! -- deallocate mt pointer, repoint, recalculate isize
1299  deallocate (mt%aint1d)
1300  mt%aint1d => aint
1301  mt%element_size = i4b
1302  mt%isize = isize
1303  mt%nrealloc = mt%nrealloc + 1
1304  mt%master = .true.
1305  nvalues_aint = nvalues_aint + isize - isizeold
1306  end subroutine reallocate_int1d
1307 
1308  !> @brief Reallocate a 2-dimensional integer array
1309  !<
1310  subroutine reallocate_int2d(aint, ncol, nrow, name, mem_path)
1311  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< the reallocated 2d integer array
1312  integer(I4B), intent(in) :: ncol !< number of columns
1313  integer(I4B), intent(in) :: nrow !< number of rows
1314  character(len=*), intent(in) :: name !< variable name
1315  character(len=*), intent(in) :: mem_path !< path where variable is stored
1316  ! -- local
1317  type(memorytype), pointer :: mt
1318  logical(LGP) :: found
1319  integer(I4B) :: istat
1320  integer(I4B), dimension(2) :: ishape
1321  integer(I4B) :: i
1322  integer(I4B) :: j
1323  integer(I4B) :: isize
1324  integer(I4B) :: isizeold
1325  ! -- code
1326  !
1327  ! -- Find and assign mt
1328  call get_from_memorystore(name, mem_path, mt, found)
1329  !
1330  ! -- Allocate aint and then refill
1331  ishape = shape(mt%aint2d)
1332  isize = nrow * ncol
1333  isizeold = ishape(1) * ishape(2)
1334  allocate (aint(ncol, nrow), stat=istat, errmsg=errmsg)
1335  if (istat /= 0) then
1336  call allocate_error(name, mem_path, istat, isize)
1337  end if
1338  do i = 1, ishape(2)
1339  do j = 1, ishape(1)
1340  aint(j, i) = mt%aint2d(j, i)
1341  end do
1342  end do
1343  !
1344  ! -- deallocate mt pointer, repoint, recalculate isize
1345  deallocate (mt%aint2d)
1346  mt%aint2d => aint
1347  mt%element_size = i4b
1348  mt%isize = isize
1349  mt%nrealloc = mt%nrealloc + 1
1350  mt%master = .true.
1351  nvalues_aint = nvalues_aint + isize - isizeold
1352  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
1353  end subroutine reallocate_int2d
1354 
1355  !> @brief Reallocate a 1-dimensional real array
1356  !<
1357  subroutine reallocate_dbl1d(adbl, nrow, name, mem_path)
1358  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< the reallocated 1d real array
1359  integer(I4B), intent(in) :: nrow !< number of rows
1360  character(len=*), intent(in) :: name !< variable name
1361  character(len=*), intent(in) :: mem_path !< path where variable is stored
1362  ! -- local
1363  type(memorytype), pointer :: mt
1364  integer(I4B) :: istat
1365  integer(I4B) :: isize
1366  integer(I4B) :: i
1367  integer(I4B) :: isizeold
1368  integer(I4B) :: ifill
1369  logical(LGP) :: found
1370  ! -- code
1371  !
1372  ! -- Find and assign mt
1373  call get_from_memorystore(name, mem_path, mt, found)
1374  !
1375  ! -- Allocate adbl and then refill
1376  isize = nrow
1377  isizeold = size(mt%adbl1d)
1378  ifill = min(isizeold, isize)
1379  allocate (adbl(nrow), stat=istat, errmsg=errmsg)
1380  if (istat /= 0) then
1381  call allocate_error(name, mem_path, istat, isize)
1382  end if
1383  do i = 1, ifill
1384  adbl(i) = mt%adbl1d(i)
1385  end do
1386  !
1387  ! -- deallocate mt pointer, repoint, recalculate isize
1388  deallocate (mt%adbl1d)
1389  mt%adbl1d => adbl
1390  mt%element_size = dp
1391  mt%isize = isize
1392  mt%nrealloc = mt%nrealloc + 1
1393  mt%master = .true.
1394  nvalues_adbl = nvalues_adbl + isize - isizeold
1395  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', isize
1396  end subroutine reallocate_dbl1d
1397 
1398  !> @brief Reallocate a 2-dimensional real array
1399  !<
1400  subroutine reallocate_dbl2d(adbl, ncol, nrow, name, mem_path)
1401  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< the reallocated 2d real array
1402  integer(I4B), intent(in) :: ncol !< number of columns
1403  integer(I4B), intent(in) :: nrow !< number of rows
1404  character(len=*), intent(in) :: name !< variable name
1405  character(len=*), intent(in) :: mem_path !< path where variable is stored
1406  ! -- local
1407  type(memorytype), pointer :: mt
1408  logical(LGP) :: found
1409  integer(I4B) :: istat
1410  integer(I4B), dimension(2) :: ishape
1411  integer(I4B) :: i
1412  integer(I4B) :: j
1413  integer(I4B) :: isize
1414  integer(I4B) :: isizeold
1415  ! -- code
1416  !
1417  ! -- Find and assign mt
1418  call get_from_memorystore(name, mem_path, mt, found)
1419  !
1420  ! -- Allocate adbl and then refill
1421  ishape = shape(mt%adbl2d)
1422  isize = nrow * ncol
1423  isizeold = ishape(1) * ishape(2)
1424  allocate (adbl(ncol, nrow), stat=istat, errmsg=errmsg)
1425  if (istat /= 0) then
1426  call allocate_error(name, mem_path, istat, isize)
1427  end if
1428  do i = 1, ishape(2)
1429  do j = 1, ishape(1)
1430  adbl(j, i) = mt%adbl2d(j, i)
1431  end do
1432  end do
1433  !
1434  ! -- deallocate mt pointer, repoint, recalculate isize
1435  deallocate (mt%adbl2d)
1436  mt%adbl2d => adbl
1437  mt%element_size = dp
1438  mt%isize = isize
1439  mt%nrealloc = mt%nrealloc + 1
1440  mt%master = .true.
1441  nvalues_adbl = nvalues_adbl + isize - isizeold
1442  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
1443  end subroutine reallocate_dbl2d
1444 
1445  !> @brief Set pointer to a logical scalar
1446  !<
1447  subroutine setptr_logical(sclr, name, mem_path)
1448  logical(LGP), pointer, intent(inout) :: sclr !< pointer to logical scalar
1449  character(len=*), intent(in) :: name !< variable name
1450  character(len=*), intent(in) :: mem_path !< path where variable is stored
1451  ! -- local
1452  type(memorytype), pointer :: mt
1453  logical(LGP) :: found
1454  ! -- code
1455  call get_from_memorystore(name, mem_path, mt, found)
1456  sclr => mt%logicalsclr
1457  end subroutine setptr_logical
1458 
1459  !> @brief Set pointer to integer scalar
1460  !<
1461  subroutine setptr_int(sclr, name, mem_path)
1462  integer(I4B), pointer, intent(inout) :: sclr !< pointer to integer scalar
1463  character(len=*), intent(in) :: name !< variable name
1464  character(len=*), intent(in) :: mem_path !< path where variable is stored
1465  ! -- local
1466  type(memorytype), pointer :: mt
1467  logical(LGP) :: found
1468  ! -- code
1469  call get_from_memorystore(name, mem_path, mt, found)
1470  sclr => mt%intsclr
1471  end subroutine setptr_int
1472 
1473  !> @brief Set pointer to 1d integer array
1474  !<
1475  subroutine setptr_int1d(aint, name, mem_path)
1476  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< pointer to 1d integer array
1477  character(len=*), intent(in) :: name !< variable name
1478  character(len=*), intent(in) :: mem_path !< path where variable is stored
1479  ! -- local
1480  type(memorytype), pointer :: mt
1481  logical(LGP) :: found
1482  ! -- code
1483  call get_from_memorystore(name, mem_path, mt, found)
1484  aint => mt%aint1d
1485  end subroutine setptr_int1d
1486 
1487  !> @brief Set pointer to 2d integer array
1488  !<
1489  subroutine setptr_int2d(aint, name, mem_path)
1490  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< pointer to 2d integer array
1491  character(len=*), intent(in) :: name !< variable name
1492  character(len=*), intent(in) :: mem_path !< path where variable is stored
1493  ! -- local
1494  type(memorytype), pointer :: mt
1495  logical(LGP) :: found
1496  ! -- code
1497  call get_from_memorystore(name, mem_path, mt, found)
1498  aint => mt%aint2d
1499  end subroutine setptr_int2d
1500 
1501  !> @brief Set pointer to 3d integer array
1502  !<
1503  subroutine setptr_int3d(aint, name, mem_path)
1504  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< pointer to 3d integer array
1505  character(len=*), intent(in) :: name !< variable name
1506  character(len=*), intent(in) :: mem_path !< path where variable is stored
1507  ! -- local
1508  type(memorytype), pointer :: mt
1509  logical(LGP) :: found
1510  ! -- code
1511  call get_from_memorystore(name, mem_path, mt, found)
1512  aint => mt%aint3d
1513  end subroutine setptr_int3d
1514 
1515  !> @brief Set pointer to a real scalar
1516  !<
1517  subroutine setptr_dbl(sclr, name, mem_path)
1518  real(DP), pointer, intent(inout) :: sclr !< pointer to a real scalar
1519  character(len=*), intent(in) :: name !< variable name
1520  character(len=*), intent(in) :: mem_path !< path where variable is stored
1521  ! -- local
1522  type(memorytype), pointer :: mt
1523  logical(LGP) :: found
1524  ! -- code
1525  call get_from_memorystore(name, mem_path, mt, found)
1526  sclr => mt%dblsclr
1527  end subroutine setptr_dbl
1528 
1529  !> @brief Set pointer to a 1d real array
1530  !<
1531  subroutine setptr_dbl1d(adbl, name, mem_path)
1532  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< pointer to 1d real array
1533  character(len=*), intent(in) :: name !< variable name
1534  character(len=*), intent(in) :: mem_path !< path where variable is stored
1535  ! -- local
1536  type(memorytype), pointer :: mt
1537  logical(LGP) :: found
1538  ! -- code
1539  call get_from_memorystore(name, mem_path, mt, found)
1540  adbl => mt%adbl1d
1541  end subroutine setptr_dbl1d
1542 
1543  !> @brief Set pointer to a 2d real array
1544  !<
1545  subroutine setptr_dbl2d(adbl, name, mem_path)
1546  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 2d real array
1547  character(len=*), intent(in) :: name !< variable name
1548  character(len=*), intent(in) :: mem_path !< path where variable is stored
1549  ! -- local
1550  type(memorytype), pointer :: mt
1551  logical(LGP) :: found
1552  ! -- code
1553  call get_from_memorystore(name, mem_path, mt, found)
1554  adbl => mt%adbl2d
1555  end subroutine setptr_dbl2d
1556 
1557  !> @brief Set pointer to a 3d real array
1558  !<
1559  subroutine setptr_dbl3d(adbl, name, mem_path)
1560  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 3d real array
1561  character(len=*), intent(in) :: name !< variable name
1562  character(len=*), intent(in) :: mem_path !< path where variable is stored
1563  ! -- local
1564  type(memorytype), pointer :: mt
1565  logical(LGP) :: found
1566  ! -- code
1567  call get_from_memorystore(name, mem_path, mt, found)
1568  adbl => mt%adbl3d
1569  end subroutine setptr_dbl3d
1570 
1571  !> @brief Set pointer to a string (scalar)
1572  !<
1573  subroutine setptr_str(asrt, name, mem_path)
1574  character(len=:), pointer :: asrt !< pointer to the character string
1575  character(len=*), intent(in) :: name !< variable name
1576  character(len=*), intent(in) :: mem_path !< path where variable is stored
1577  ! -- local
1578  type(memorytype), pointer :: mt
1579  logical(LGP) :: found
1580  ! -- code
1581  call get_from_memorystore(name, mem_path, mt, found)
1582  asrt => mt%strsclr
1583  end subroutine setptr_str
1584 
1585  !> @brief Set pointer to a fixed-length string array
1586  !<
1587  subroutine setptr_str1d(astr1d, name, mem_path)
1588  character(len=:), dimension(:), &
1589  pointer, contiguous, intent(inout) :: astr1d !< pointer to the string array
1590  character(len=*), intent(in) :: name !< variable name
1591  character(len=*), intent(in) :: mem_path !< path where variable is stored
1592  ! -- local
1593  type(memorytype), pointer :: mt
1594  logical(LGP) :: found
1595  ! -- code
1596  call get_from_memorystore(name, mem_path, mt, found)
1597  astr1d => mt%astr1d
1598  end subroutine setptr_str1d
1599 
1600  !> @brief Set pointer to an array of CharacterStringType
1601  !<
1602  subroutine setptr_charstr1d(acharstr1d, name, mem_path)
1603  type(characterstringtype), dimension(:), pointer, contiguous, &
1604  intent(inout) :: acharstr1d !< the reallocated charstring array
1605  character(len=*), intent(in) :: name !< variable name
1606  character(len=*), intent(in) :: mem_path !< path where variable is stored
1607  ! -- local
1608  type(memorytype), pointer :: mt
1609  logical(LGP) :: found
1610  ! -- code
1611  call get_from_memorystore(name, mem_path, mt, found)
1612  acharstr1d => mt%acharstr1d
1613  end subroutine setptr_charstr1d
1614 
1615  !> @brief Make a copy of a 1-dimensional integer array
1616  !<
1617  subroutine copyptr_int1d(aint, name, mem_path, mem_path_copy)
1618  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< returned copy of 1d integer array
1619  character(len=*), intent(in) :: name !< variable name
1620  character(len=*), intent(in) :: mem_path !< path where variable is stored
1621  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1622  !! if passed then the copy is added to the
1623  !! memory manager
1624  ! -- local
1625  type(memorytype), pointer :: mt
1626  logical(LGP) :: found
1627  integer(I4B) :: n
1628  ! -- code
1629  call get_from_memorystore(name, mem_path, mt, found)
1630  aint => null()
1631  ! -- check the copy into the memory manager
1632  if (present(mem_path_copy)) then
1633  call allocate_int1d(aint, size(mt%aint1d), mt%name, mem_path_copy)
1634  ! -- create a local copy
1635  else
1636  allocate (aint(size(mt%aint1d)))
1637  end if
1638  do n = 1, size(mt%aint1d)
1639  aint(n) = mt%aint1d(n)
1640  end do
1641  end subroutine copyptr_int1d
1642 
1643  !> @brief Make a copy of a 2-dimensional integer array
1644  !<
1645  subroutine copyptr_int2d(aint, name, mem_path, mem_path_copy)
1646  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< returned copy of 2d integer array
1647  character(len=*), intent(in) :: name !< variable name
1648  character(len=*), intent(in) :: mem_path !< path where variable is stored
1649  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1650  !! if passed then the copy is added to the
1651  !! memory manager
1652  ! -- local
1653  type(memorytype), pointer :: mt
1654  logical(LGP) :: found
1655  integer(I4B) :: i
1656  integer(I4B) :: j
1657  integer(I4B) :: ncol
1658  integer(I4B) :: nrow
1659  ! -- code
1660  call get_from_memorystore(name, mem_path, mt, found)
1661  aint => null()
1662  ncol = size(mt%aint2d, dim=1)
1663  nrow = size(mt%aint2d, dim=2)
1664  ! -- check the copy into the memory manager
1665  if (present(mem_path_copy)) then
1666  call allocate_int2d(aint, ncol, nrow, mt%name, mem_path_copy)
1667  ! -- create a local copy
1668  else
1669  allocate (aint(ncol, nrow))
1670  end if
1671  do i = 1, nrow
1672  do j = 1, ncol
1673  aint(j, i) = mt%aint2d(j, i)
1674  end do
1675  end do
1676  end subroutine copyptr_int2d
1677 
1678  !> @brief Make a copy of a 1-dimensional real array
1679  !<
1680  subroutine copyptr_dbl1d(adbl, name, mem_path, mem_path_copy)
1681  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< returned copy of 1d real array
1682  character(len=*), intent(in) :: name !< variable name
1683  character(len=*), intent(in) :: mem_path !< path where variable is stored
1684  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1685  !! if passed then the copy is added to the
1686  !! memory manager
1687  ! -- local
1688  type(memorytype), pointer :: mt
1689  logical(LGP) :: found
1690  integer(I4B) :: n
1691  ! -- code
1692  call get_from_memorystore(name, mem_path, mt, found)
1693  adbl => null()
1694  ! -- check the copy into the memory manager
1695  if (present(mem_path_copy)) then
1696  call allocate_dbl1d(adbl, size(mt%adbl1d), mt%name, mem_path_copy)
1697  ! -- create a local copy
1698  else
1699  allocate (adbl(size(mt%adbl1d)))
1700  end if
1701  do n = 1, size(mt%adbl1d)
1702  adbl(n) = mt%adbl1d(n)
1703  end do
1704  end subroutine copyptr_dbl1d
1705 
1706  !> @brief Make a copy of a 2-dimensional real array
1707  !<
1708  subroutine copyptr_dbl2d(adbl, name, mem_path, mem_path_copy)
1709  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< returned copy of 2d real array
1710  character(len=*), intent(in) :: name !< variable name
1711  character(len=*), intent(in) :: mem_path !< path where variable is stored
1712  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1713  !! if passed then the copy is added to the
1714  !! memory manager
1715  ! -- local
1716  type(memorytype), pointer :: mt
1717  logical(LGP) :: found
1718  integer(I4B) :: i
1719  integer(I4B) :: j
1720  integer(I4B) :: ncol
1721  integer(I4B) :: nrow
1722  ! -- code
1723  call get_from_memorystore(name, mem_path, mt, found)
1724  adbl => null()
1725  ncol = size(mt%adbl2d, dim=1)
1726  nrow = size(mt%adbl2d, dim=2)
1727  ! -- check the copy into the memory manager
1728  if (present(mem_path_copy)) then
1729  call allocate_dbl2d(adbl, ncol, nrow, mt%name, mem_path_copy)
1730  ! -- create a local copy
1731  else
1732  allocate (adbl(ncol, nrow))
1733  end if
1734  do i = 1, nrow
1735  do j = 1, ncol
1736  adbl(j, i) = mt%adbl2d(j, i)
1737  end do
1738  end do
1739  end subroutine copyptr_dbl2d
1740 
1741  !> @brief Copy values from a 1-dimensional real array in the memory
1742  !< manager to a passed 1-dimensional real array
1743  subroutine copy_dbl1d(adbl, name, mem_path)
1744  real(dp), dimension(:), intent(inout) :: adbl !< target array
1745  character(len=*), intent(in) :: name !< variable name
1746  character(len=*), intent(in) :: mem_path !< path where variable is stored
1747  ! -- local
1748  type(memorytype), pointer :: mt
1749  logical(LGP) :: found
1750  integer(I4B) :: n
1751  ! -- code
1752  call get_from_memorystore(name, mem_path, mt, found)
1753  do n = 1, size(mt%adbl1d)
1754  adbl(n) = mt%adbl1d(n)
1755  end do
1756  end subroutine copy_dbl1d
1757 
1758  !> @brief Set the pointer for an integer scalar to
1759  !< a target array already stored in the memory manager
1760  subroutine reassignptr_int(sclr, name, mem_path, name_target, mem_path_target)
1761  integer(I4B), pointer, intent(inout) :: sclr !< pointer to integer scalar
1762  character(len=*), intent(in) :: name !< variable name
1763  character(len=*), intent(in) :: mem_path !< path where variable is stored
1764  character(len=*), intent(in) :: name_target !< name of target variable
1765  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1766  ! -- local
1767  type(memorytype), pointer :: mt
1768  type(memorytype), pointer :: mt2
1769  logical(LGP) :: found
1770  ! -- code
1771  call get_from_memorystore(name, mem_path, mt, found)
1772  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1773  if (associated(sclr)) then
1775  deallocate (sclr)
1776  end if
1777  sclr => mt2%intsclr
1778  mt%intsclr => sclr
1779  mt%element_size = i4b
1780  mt%isize = 1
1781  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', mt%isize
1782  !
1783  ! -- set master information
1784  mt%master = .false.
1785  mt%mastername = name_target
1786  mt%masterPath = mem_path_target
1787  end subroutine reassignptr_int
1788 
1789  !> @brief Set the pointer for a 1-dimensional integer array to
1790  !< a target array already stored in the memory manager
1791  subroutine reassignptr_int1d(aint, name, mem_path, name_target, mem_path_target)
1792  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< pointer to 1d integer array
1793  character(len=*), intent(in) :: name !< variable name
1794  character(len=*), intent(in) :: mem_path !< path where variable is stored
1795  character(len=*), intent(in) :: name_target !< name of target variable
1796  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1797  ! -- local
1798  type(memorytype), pointer :: mt
1799  type(memorytype), pointer :: mt2
1800  logical(LGP) :: found
1801  ! -- code
1802  call get_from_memorystore(name, mem_path, mt, found)
1803  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1804  if (size(aint) > 0) then
1805  nvalues_aint = nvalues_aint - size(aint)
1806  deallocate (aint)
1807  end if
1808  aint => mt2%aint1d
1809  mt%aint1d => aint
1810  mt%element_size = i4b
1811  mt%isize = size(aint)
1812  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', mt%isize
1813  !
1814  ! -- set master information
1815  mt%master = .false.
1816  mt%mastername = name_target
1817  mt%masterPath = mem_path_target
1818  end subroutine reassignptr_int1d
1819 
1820  !> @brief Set the pointer for a 2-dimensional integer array to
1821  !< a target array already stored in the memory manager
1822  subroutine reassignptr_int2d(aint, name, mem_path, name_target, mem_path_target)
1823  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< pointer to 2d integer array
1824  character(len=*), intent(in) :: name !< variable name
1825  character(len=*), intent(in) :: mem_path !< path where variable is stored
1826  character(len=*), intent(in) :: name_target !< name of target variable
1827  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1828  ! -- local
1829  type(memorytype), pointer :: mt
1830  type(memorytype), pointer :: mt2
1831  logical(LGP) :: found
1832  integer(I4B) :: ncol
1833  integer(I4B) :: nrow
1834  ! -- code
1835  call get_from_memorystore(name, mem_path, mt, found)
1836  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1837  if (size(aint) > 0) then
1838  nvalues_aint = nvalues_aint - size(aint)
1839  deallocate (aint)
1840  end if
1841  aint => mt2%aint2d
1842  mt%aint2d => aint
1843  mt%element_size = i4b
1844  mt%isize = size(aint)
1845  ncol = size(aint, dim=1)
1846  nrow = size(aint, dim=2)
1847  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
1848  !
1849  ! -- set master information
1850  mt%master = .false.
1851  mt%mastername = name_target
1852  mt%masterPath = mem_path_target
1853  end subroutine reassignptr_int2d
1854 
1855  !> @brief Set the pointer for a 1-dimensional real array to
1856  !< a target array already stored in the memory manager
1857  subroutine reassignptr_dbl1d(adbl, name, mem_path, name_target, mem_path_target)
1858  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< pointer to 1d real array
1859  character(len=*), intent(in) :: name !< variable name
1860  character(len=*), intent(in) :: mem_path !< path where variable is stored
1861  character(len=*), intent(in) :: name_target !< name of target variable
1862  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1863  ! -- local
1864  type(memorytype), pointer :: mt
1865  type(memorytype), pointer :: mt2
1866  logical(LGP) :: found
1867  ! -- code
1868  call get_from_memorystore(name, mem_path, mt, found)
1869  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1870  if (size(adbl) > 0) then
1871  nvalues_adbl = nvalues_adbl - size(adbl)
1872  deallocate (adbl)
1873  end if
1874  adbl => mt2%adbl1d
1875  mt%adbl1d => adbl
1876  mt%element_size = dp
1877  mt%isize = size(adbl)
1878  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', mt%isize
1879  !
1880  ! -- set master information
1881  mt%master = .false.
1882  mt%mastername = name_target
1883  mt%masterPath = mem_path_target
1884  end subroutine reassignptr_dbl1d
1885 
1886  !> @brief Set the pointer for a 2-dimensional real array to
1887  !< a target array already stored in the memory manager
1888  subroutine reassignptr_dbl2d(adbl, name, mem_path, name_target, mem_path_target)
1889  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 2d real array
1890  character(len=*), intent(in) :: name !< variable name
1891  character(len=*), intent(in) :: mem_path !< path where variable is stored
1892  character(len=*), intent(in) :: name_target !< name of target variable
1893  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1894  ! -- local
1895  type(memorytype), pointer :: mt
1896  type(memorytype), pointer :: mt2
1897  logical(LGP) :: found
1898  integer(I4B) :: ncol
1899  integer(I4b) :: nrow
1900  ! -- code
1901  call get_from_memorystore(name, mem_path, mt, found)
1902  call get_from_memorystore(name_target, mem_path_target, mt2, found)
1903  if (size(adbl) > 0) then
1904  nvalues_adbl = nvalues_adbl - size(adbl)
1905  deallocate (adbl)
1906  end if
1907  adbl => mt2%adbl2d
1908  mt%adbl2d => adbl
1909  mt%element_size = dp
1910  mt%isize = size(adbl)
1911  ncol = size(adbl, dim=1)
1912  nrow = size(adbl, dim=2)
1913  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
1914  !
1915  ! -- set master information
1916  mt%master = .false.
1917  mt%mastername = name_target
1918  mt%masterPath = mem_path_target
1919  end subroutine reassignptr_dbl2d
1920 
1921  !> @brief Deallocate a variable-length character string
1922  !<
1923  subroutine deallocate_str(sclr, name, mem_path)
1924  character(len=*), pointer, intent(inout) :: sclr !< pointer to string
1925  character(len=*), intent(in), optional :: name !< variable name
1926  character(len=*), intent(in), optional :: mem_path !< path where variable is stored
1927  ! -- local
1928  type(memorytype), pointer :: mt
1929  logical(LGP) :: found
1930  type(memorycontaineriteratortype), allocatable :: itr
1931  ! -- code
1932  found = .false.
1933  if (present(name) .and. present(mem_path)) then
1934  call get_from_memorystore(name, mem_path, mt, found)
1935  nullify (mt%strsclr)
1936  else
1937  itr = memorystore%iterator()
1938  do while (itr%has_next())
1939  call itr%next()
1940  mt => itr%value()
1941  if (associated(mt%strsclr, sclr)) then
1942  nullify (mt%strsclr)
1943  found = .true.
1944  exit
1945  end if
1946  end do
1947  end if
1948  if (.not. found) then
1949  call store_error('Programming error in deallocate_str.', terminate=.true.)
1950  else
1951  if (mt%master) then
1952  deallocate (sclr)
1953  else
1954  nullify (sclr)
1955  end if
1956  end if
1957  end subroutine deallocate_str
1958 
1959  !> @brief Deallocate an array of defined-length character strings
1960  !!
1961  !<
1962  subroutine deallocate_str1d(astr1d, name, mem_path)
1963  character(len=*), dimension(:), pointer, contiguous, intent(inout) :: astr1d !< array of strings
1964  character(len=*), optional, intent(in) :: name !< variable name
1965  character(len=*), optional, intent(in) :: mem_path !< path where variable is stored
1966  ! -- local
1967  type(memorytype), pointer :: mt
1968  logical(LGP) :: found
1969  type(memorycontaineriteratortype), allocatable :: itr
1970  ! -- code
1971  !
1972  ! -- process optional variables
1973  found = .false.
1974  if (present(name) .and. present(mem_path)) then
1975  call get_from_memorystore(name, mem_path, mt, found)
1976  nullify (mt%astr1d)
1977  else
1978  itr = memorystore%iterator()
1979  do while (itr%has_next())
1980  call itr%next()
1981  mt => itr%value()
1982  if (associated(mt%astr1d, astr1d)) then
1983  nullify (mt%astr1d)
1984  found = .true.
1985  exit
1986  end if
1987  end do
1988  end if
1989  if (.not. found .and. size(astr1d) > 0) then
1990  call store_error('programming error in deallocate_str1d', terminate=.true.)
1991  else
1992  if (mt%master) then
1993  deallocate (astr1d)
1994  else
1995  nullify (astr1d)
1996  end if
1997  end if
1998  end subroutine deallocate_str1d
1999 
2000  !> @brief Deallocate an array of deferred-length character strings
2001  !!
2002  !<
2003  subroutine deallocate_charstr1d(astr1d, name, mem_path)
2004  type(characterstringtype), dimension(:), pointer, contiguous, &
2005  intent(inout) :: astr1d !< array of strings
2006  character(len=*), optional, intent(in) :: name !< variable name
2007  character(len=*), optional, intent(in) :: mem_path !< path where variable is stored
2008  ! -- local
2009  type(memorytype), pointer :: mt
2010  logical(LGP) :: found
2011  type(memorycontaineriteratortype), allocatable :: itr
2012  ! -- code
2013  !
2014  ! -- process optional variables
2015  found = .false.
2016  if (present(name) .and. present(mem_path)) then
2017  call get_from_memorystore(name, mem_path, mt, found)
2018  nullify (mt%acharstr1d)
2019  else
2020  itr = memorystore%iterator()
2021  do while (itr%has_next())
2022  call itr%next()
2023  mt => itr%value()
2024  if (associated(mt%acharstr1d, astr1d)) then
2025  nullify (mt%acharstr1d)
2026  found = .true.
2027  exit
2028  end if
2029  end do
2030  end if
2031  if (.not. found .and. size(astr1d) > 0) then
2032  call store_error('programming error in deallocate_charstr1d', &
2033  terminate=.true.)
2034  else
2035  if (mt%master) then
2036  deallocate (astr1d)
2037  else
2038  nullify (astr1d)
2039  end if
2040  end if
2041  end subroutine deallocate_charstr1d
2042 
2043  !> @brief Deallocate a logical scalar
2044  !<
2045  subroutine deallocate_logical(sclr)
2046  logical(LGP), pointer, intent(inout) :: sclr !< logical scalar to deallocate
2047  ! -- local
2048  class(memorytype), pointer :: mt
2049  logical(LGP) :: found
2050  type(memorycontaineriteratortype), allocatable :: itr
2051  ! -- code
2052  found = .false.
2053  itr = memorystore%iterator()
2054  do while (itr%has_next())
2055  call itr%next()
2056  mt => itr%value()
2057  if (associated(mt%logicalsclr, sclr)) then
2058  nullify (mt%logicalsclr)
2059  found = .true.
2060  exit
2061  end if
2062  end do
2063  if (.not. found) then
2064  call store_error('programming error in deallocate_logical', &
2065  terminate=.true.)
2066  else
2067  if (mt%master) then
2068  deallocate (sclr)
2069  else
2070  nullify (sclr)
2071  end if
2072  end if
2073  end subroutine deallocate_logical
2074 
2075  !> @brief Deallocate a integer scalar
2076  !<
2077  subroutine deallocate_int(sclr)
2078  integer(I4B), pointer, intent(inout) :: sclr !< integer variable to deallocate
2079  ! -- local
2080  class(memorytype), pointer :: mt
2081  logical(LGP) :: found
2082  type(memorycontaineriteratortype), allocatable :: itr
2083  ! -- code
2084  found = .false.
2085  itr = memorystore%iterator()
2086  do while (itr%has_next())
2087  call itr%next()
2088  mt => itr%value()
2089  if (associated(mt%intsclr, sclr)) then
2090  nullify (mt%intsclr)
2091  found = .true.
2092  exit
2093  end if
2094  end do
2095  if (.not. found) then
2096  call store_error('Programming error in deallocate_int.', terminate=.true.)
2097  else
2098  if (mt%master) then
2099  deallocate (sclr)
2100  else
2101  nullify (sclr)
2102  end if
2103  end if
2104  end subroutine deallocate_int
2105 
2106  !> @brief Deallocate a real scalar
2107  !<
2108  subroutine deallocate_dbl(sclr)
2109  real(DP), pointer, intent(inout) :: sclr !< real variable to deallocate
2110  ! -- local
2111  class(memorytype), pointer :: mt
2112  logical(LGP) :: found
2113  type(memorycontaineriteratortype), allocatable :: itr
2114  ! -- code
2115  found = .false.
2116  itr = memorystore%iterator()
2117  do while (itr%has_next())
2118  call itr%next()
2119  mt => itr%value()
2120  if (associated(mt%dblsclr, sclr)) then
2121  nullify (mt%dblsclr)
2122  found = .true.
2123  exit
2124  end if
2125  end do
2126  if (.not. found) then
2127  call store_error('Programming error in deallocate_dbl.', terminate=.true.)
2128  else
2129  if (mt%master) then
2130  deallocate (sclr)
2131  else
2132  nullify (sclr)
2133  end if
2134  end if
2135  end subroutine deallocate_dbl
2136 
2137  !> @brief Deallocate a 1-dimensional integer array
2138  !<
2139  subroutine deallocate_int1d(aint, name, mem_path)
2140  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< 1d integer array to deallocate
2141  character(len=*), optional :: name !< variable name
2142  character(len=*), optional :: mem_path !< path where variable is stored
2143  ! -- local
2144  type(memorytype), pointer :: mt
2145  logical(LGP) :: found
2146  type(memorycontaineriteratortype), allocatable :: itr
2147  ! -- code
2148  !
2149  ! -- process optional variables
2150  found = .false.
2151  if (present(name) .and. present(mem_path)) then
2152  call get_from_memorystore(name, mem_path, mt, found)
2153  nullify (mt%aint1d)
2154  else
2155  itr = memorystore%iterator()
2156  do while (itr%has_next())
2157  call itr%next()
2158  mt => itr%value()
2159  if (associated(mt%aint1d, aint)) then
2160  nullify (mt%aint1d)
2161  found = .true.
2162  exit
2163  end if
2164  end do
2165  end if
2166  if (.not. found .and. size(aint) > 0) then
2167  call store_error('programming error in deallocate_int1d', terminate=.true.)
2168  else
2169  if (mt%master) then
2170  deallocate (aint)
2171  else
2172  nullify (aint)
2173  end if
2174  end if
2175  end subroutine deallocate_int1d
2176 
2177  !> @brief Deallocate a 2-dimensional integer array
2178  !<
2179  subroutine deallocate_int2d(aint, name, mem_path)
2180  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< 2d integer array to deallocate
2181  character(len=*), optional :: name !< variable name
2182  character(len=*), optional :: mem_path !< path where variable is stored
2183  ! -- local
2184  type(memorytype), pointer :: mt
2185  logical(LGP) :: found
2186  type(memorycontaineriteratortype), allocatable :: itr
2187  ! -- code
2188  !
2189  ! -- process optional variables
2190  found = .false.
2191  if (present(name) .and. present(mem_path)) then
2192  call get_from_memorystore(name, mem_path, mt, found)
2193  nullify (mt%aint2d)
2194  else
2195  itr = memorystore%iterator()
2196  do while (itr%has_next())
2197  call itr%next()
2198  mt => itr%value()
2199  if (associated(mt%aint2d, aint)) then
2200  nullify (mt%aint2d)
2201  found = .true.
2202  exit
2203  end if
2204  end do
2205  end if
2206  if (.not. found .and. size(aint) > 0) then
2207  call store_error('programming error in deallocate_int2d', terminate=.true.)
2208  else
2209  if (mt%master) then
2210  deallocate (aint)
2211  else
2212  nullify (aint)
2213  end if
2214  end if
2215  end subroutine deallocate_int2d
2216 
2217  !> @brief Deallocate a 3-dimensional integer array
2218  !<
2219  subroutine deallocate_int3d(aint, name, mem_path)
2220  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< 3d integer array to deallocate
2221  character(len=*), optional :: name !< variable name
2222  character(len=*), optional :: mem_path !< path where variable is stored
2223  ! -- local
2224  type(memorytype), pointer :: mt
2225  logical(LGP) :: found
2226  type(memorycontaineriteratortype), allocatable :: itr
2227  ! -- code
2228  !
2229  ! -- process optional variables
2230  found = .false.
2231  if (present(name) .and. present(mem_path)) then
2232  call get_from_memorystore(name, mem_path, mt, found)
2233  nullify (mt%aint3d)
2234  else
2235  itr = memorystore%iterator()
2236  do while (itr%has_next())
2237  call itr%next()
2238  mt => itr%value()
2239  if (associated(mt%aint3d, aint)) then
2240  nullify (mt%aint3d)
2241  found = .true.
2242  exit
2243  end if
2244  end do
2245  end if
2246  if (.not. found .and. size(aint) > 0) then
2247  call store_error('programming error in deallocate_int3d', terminate=.true.)
2248  else
2249  if (mt%master) then
2250  deallocate (aint)
2251  else
2252  nullify (aint)
2253  end if
2254  end if
2255  end subroutine deallocate_int3d
2256 
2257  !> @brief Deallocate a 1-dimensional real array
2258  !<
2259  subroutine deallocate_dbl1d(adbl, name, mem_path)
2260  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< 1d real array to deallocate
2261  character(len=*), optional :: name !< variable name
2262  character(len=*), optional :: mem_path !< path where variable is stored
2263  ! -- local
2264  type(memorytype), pointer :: mt
2265  logical(LGP) :: found
2266  type(memorycontaineriteratortype), allocatable :: itr
2267  ! -- code
2268  !
2269  ! -- process optional variables
2270  found = .false.
2271  if (present(name) .and. present(mem_path)) then
2272  call get_from_memorystore(name, mem_path, mt, found)
2273  nullify (mt%adbl1d)
2274  else
2275  itr = memorystore%iterator()
2276  do while (itr%has_next())
2277  call itr%next()
2278  mt => itr%value()
2279  if (associated(mt%adbl1d, adbl)) then
2280  nullify (mt%adbl1d)
2281  found = .true.
2282  exit
2283  end if
2284  end do
2285  end if
2286  if (.not. found .and. size(adbl) > 0) then
2287  call store_error('programming error in deallocate_dbl1d', terminate=.true.)
2288  else
2289  if (mt%master) then
2290  deallocate (adbl)
2291  else
2292  nullify (adbl)
2293  end if
2294  end if
2295  end subroutine deallocate_dbl1d
2296 
2297  !> @brief Deallocate a 2-dimensional real array
2298  !<
2299  subroutine deallocate_dbl2d(adbl, name, mem_path)
2300  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< 2d real array to deallocate
2301  character(len=*), optional :: name !< variable name
2302  character(len=*), optional :: mem_path !< path where variable is stored
2303  ! -- local
2304  type(memorytype), pointer :: mt
2305  logical(LGP) :: found
2306  type(memorycontaineriteratortype), allocatable :: itr
2307  ! -- code
2308  !
2309  ! -- process optional variables
2310  found = .false.
2311  if (present(name) .and. present(mem_path)) then
2312  call get_from_memorystore(name, mem_path, mt, found)
2313  nullify (mt%adbl2d)
2314  else
2315  itr = memorystore%iterator()
2316  do while (itr%has_next())
2317  call itr%next()
2318  mt => itr%value()
2319  if (associated(mt%adbl2d, adbl)) then
2320  nullify (mt%adbl2d)
2321  found = .true.
2322  exit
2323  end if
2324  end do
2325  end if
2326  if (.not. found .and. size(adbl) > 0) then
2327  call store_error('programming error in deallocate_dbl2d', terminate=.true.)
2328  else
2329  if (mt%master) then
2330  deallocate (adbl)
2331  else
2332  nullify (adbl)
2333  end if
2334  end if
2335  end subroutine deallocate_dbl2d
2336 
2337  !> @brief Deallocate a 3-dimensional real array
2338  !<
2339  subroutine deallocate_dbl3d(adbl, name, mem_path)
2340  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< 3d real array to deallocate
2341  character(len=*), optional :: name !< variable name
2342  character(len=*), optional :: mem_path !< path where variable is stored
2343  ! -- local
2344  type(memorytype), pointer :: mt
2345  logical(LGP) :: found
2346  type(memorycontaineriteratortype), allocatable :: itr
2347  ! -- code
2348  !
2349  ! -- process optional variables
2350  found = .false.
2351  if (present(name) .and. present(mem_path)) then
2352  call get_from_memorystore(name, mem_path, mt, found)
2353  nullify (mt%adbl3d)
2354  else
2355  itr = memorystore%iterator()
2356  do while (itr%has_next())
2357  call itr%next()
2358  mt => itr%value()
2359  if (associated(mt%adbl3d, adbl)) then
2360  nullify (mt%adbl3d)
2361  found = .true.
2362  exit
2363  end if
2364  end do
2365  end if
2366  if (.not. found .and. size(adbl) > 0) then
2367  call store_error('programming error in deallocate_dbl3d', terminate=.true.)
2368  else
2369  if (mt%master) then
2370  deallocate (adbl)
2371  else
2372  nullify (adbl)
2373  end if
2374  end if
2375  end subroutine deallocate_dbl3d
2376 
2377  !> @brief Set the memory print option
2378  !<
2379  subroutine mem_set_print_option(iout, keyword, error_msg)
2380  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2381  character(len=*), intent(in) :: keyword !< memory print option
2382  character(len=*), intent(inout) :: error_msg !< returned error message if keyword is not valid option
2383  ! -- local
2384  ! -- format
2385  ! -- code
2386  select case (keyword)
2387  case ('NONE')
2388  iprmem = 0
2389  write (iout, '(4x, a)') &
2390  'LIMITED MEMORY INFORMATION WILL BE WRITTEN.'
2391  case ('SUMMARY')
2392  iprmem = 1
2393  write (iout, '(4x, a)') &
2394  'A SUMMARY OF SIMULATION MEMORY INFORMATION WILL BE WRITTEN.'
2395  case ('ALL')
2396  iprmem = 2
2397  write (iout, '(4x, a)') &
2398  'ALL SIMULATION MEMORY INFORMATION WILL BE WRITTEN.'
2399  case default
2400  error_msg = "Unknown memory print option '"//trim(keyword)//"."
2401  end select
2402  end subroutine mem_set_print_option
2403 
2404  !> @brief Create a table if memory_print_option is 'SUMMARY'
2405  !<
2406  subroutine mem_summary_table(iout, nrows, cunits)
2407  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2408  integer(I4B), intent(in) :: nrows !< number of table rows
2409  character(len=*), intent(in) :: cunits !< memory units (bytes, kilobytes, megabytes, or gigabytes)
2410  ! -- local
2411  character(len=LINELENGTH) :: title
2412  character(len=LINELENGTH) :: text
2413  integer(I4B) :: nterms
2414  ! -- formats
2415  ! -- code
2416  nterms = 6
2417  !
2418  ! -- set up table title
2419  title = 'SUMMARY INFORMATION ON VARIABLES STORED IN THE MEMORY MANAGER, '// &
2420  'IN '//trim(cunits)
2421  !
2422  ! -- set up stage tableobj
2423  call table_cr(memtab, 'MEM SUM', title)
2424  call memtab%table_df(nrows, nterms, iout)
2425  !
2426  ! -- data type
2427  text = 'COMPONENT'
2428  call memtab%initialize_column(text, 20, alignment=tableft)
2429  !
2430  ! -- memory allocated for characters
2431  text = 'CHARACTER'
2432  call memtab%initialize_column(text, 15, alignment=tabcenter)
2433  !
2434  ! -- memory allocated for logical
2435  text = 'LOGICAL'
2436  call memtab%initialize_column(text, 15, alignment=tabcenter)
2437  !
2438  ! -- memory allocated for integers
2439  text = 'INTEGER'
2440  call memtab%initialize_column(text, 15, alignment=tabcenter)
2441  !
2442  ! -- memory allocated for reals
2443  text = 'REAL'
2444  call memtab%initialize_column(text, 15, alignment=tabcenter)
2445  !
2446  ! -- total memory allocated
2447  text = 'TOTAL'
2448  call memtab%initialize_column(text, 15, alignment=tabcenter)
2449  end subroutine mem_summary_table
2450 
2451  !> @brief Create a table if memory_print_option is 'ALL'
2452  !<
2453  subroutine mem_detailed_table(iout, nrows)
2454  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2455  integer(I4B), intent(in) :: nrows !< number of table rows
2456  ! -- local
2457  character(len=LINELENGTH) :: title
2458  character(len=LINELENGTH) :: text
2459  integer(I4B) :: nterms
2460  ! -- formats
2461  ! -- code
2462  nterms = 5
2463  !
2464  ! -- set up table title
2465  title = 'DETAILED INFORMATION ON VARIABLES STORED IN THE MEMORY MANAGER'
2466  !
2467  ! -- set up stage tableobj
2468  call table_cr(memtab, 'MEM DET', title)
2469  call memtab%table_df(nrows, nterms, iout)
2470  !
2471  ! -- origin
2472  text = 'ORIGIN'
2473  call memtab%initialize_column(text, lenmempath, alignment=tableft)
2474  !
2475  ! -- variable
2476  text = 'VARIABLE NAME'
2477  call memtab%initialize_column(text, lenvarname, alignment=tableft)
2478  !
2479  ! -- data type
2480  text = 'DATA TYPE'
2481  call memtab%initialize_column(text, 16, alignment=tableft)
2482  !
2483  ! -- size
2484  text = 'NUMBER OF ITEMS'
2485  call memtab%initialize_column(text, 20, alignment=tabright)
2486  !
2487  ! -- is it a pointer
2488  text = 'ASSOCIATED VARIABLE'
2489  call memtab%initialize_column(text, lenmemaddress, alignment=tableft)
2490  end subroutine mem_detailed_table
2491 
2492  !> @brief Write a row for the memory_print_option 'SUMMARY' table
2493  !<
2494  subroutine mem_summary_line(component, rchars, rlog, rint, rreal, bytes)
2495  character(len=*), intent(in) :: component !< character defining the program component (e.g. solution)
2496  real(DP), intent(in) :: rchars !< allocated size of characters (in common units)
2497  real(DP), intent(in) :: rlog !< allocated size of logical (in common units)
2498  real(DP), intent(in) :: rint !< allocated size of integer variables (in common units)
2499  real(DP), intent(in) :: rreal !< allocated size of real variables (in common units)
2500  real(DP), intent(in) :: bytes !< total allocated memory in memory manager (in common units)
2501  ! -- formats
2502  ! -- code
2503  !
2504  ! -- write line
2505  call memtab%add_term(component)
2506  call memtab%add_term(rchars)
2507  call memtab%add_term(rlog)
2508  call memtab%add_term(rint)
2509  call memtab%add_term(rreal)
2510  call memtab%add_term(bytes)
2511  end subroutine mem_summary_line
2512 
2513  !> @brief Determine appropriate memory unit and conversion factor
2514  !<
2515  subroutine mem_units(bytes, fact, cunits)
2516  ! -- dummy
2517  real(DP), intent(in) :: bytes !< total nr. of bytes
2518  real(DP), intent(inout) :: fact !< conversion factor
2519  character(len=*), intent(inout) :: cunits !< string with memory unit
2520  ! -- local
2521  ! -- formats
2522  ! -- code
2523  !
2524  ! -- initialize factor and unit string
2525  cunits = 'UNKNOWN'
2526  fact = done
2527  !
2528  ! -- set factor
2529  if (bytes < dep3) then
2530  fact = done
2531  cunits = 'BYTES'
2532  else if (bytes < dep6) then
2533  fact = dem3
2534  cunits = 'KILOBYTES'
2535  else if (bytes < dep9) then
2536  fact = dem6
2537  cunits = 'MEGABYTES'
2538  else
2539  fact = dem9
2540  cunits = 'GIGABYTES'
2541  end if
2542  end subroutine mem_units
2543 
2544  !> @brief Create and fill a table with the total allocated memory
2545  !< in the memory manager
2546  subroutine mem_summary_total(iout, bytes)
2547  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2548  real(DP), intent(in) :: bytes !< total number of bytes allocated in the memory manager
2549  ! -- local
2550  character(len=LINELENGTH) :: title
2551  character(len=LINELENGTH) :: text
2552  character(LEN=10) :: cunits
2553  integer(I4B) :: nterms
2554  integer(I4B) :: nrows
2555  real(DP) :: fact
2556  real(DP) :: smb
2557  ! -- formats
2558  ! -- code
2559  !
2560  ! -- calculate factor and memory units
2561  call mem_units(bytes, fact, cunits)
2562  !
2563  ! -- set table terms
2564  nterms = 2
2565  nrows = 6
2566  !
2567  ! -- set up table title
2568  title = 'MEMORY MANAGER TOTAL STORAGE BY DATA TYPE, IN '//trim(cunits)
2569  !
2570  ! -- set up stage tableobj
2571  call table_cr(memtab, 'MEM TOT', title)
2572  call memtab%table_df(nrows, nterms, iout)
2573  !
2574  ! -- data type
2575  text = 'DATA TYPE'
2576  call memtab%initialize_column(text, 15, alignment=tableft)
2577  !
2578  ! -- number of values
2579  text = 'ALLOCATED MEMORY'
2580  call memtab%initialize_column(text, 15, alignment=tabcenter)
2581  !
2582  ! -- write data
2583  !
2584  ! -- characters
2585  smb = real(nvalues_astr, dp) * fact
2586  call memtab%add_term('Character')
2587  call memtab%add_term(smb)
2588  !
2589  ! -- logicals
2590  smb = real(nvalues_alogical * lgp, dp) * fact
2591  call memtab%add_term('Logical')
2592  call memtab%add_term(smb)
2593  !
2594  ! -- integers
2595  smb = real(nvalues_aint * i4b, dp) * fact
2596  call memtab%add_term('Integer')
2597  call memtab%add_term(smb)
2598  !
2599  ! -- reals
2600  smb = real(nvalues_adbl * dp, dp) * fact
2601  call memtab%add_term('Real')
2602  call memtab%add_term(smb)
2603  !
2604  ! -- total memory usage
2605  call memtab%print_separator()
2606  smb = bytes * fact
2607  call memtab%add_term('Total')
2608  call memtab%add_term(smb)
2609  !
2610  ! -- Virtual memory
2611  smb = calc_virtual_mem() * fact
2612  call memtab%add_term('Virtual')
2613  call memtab%add_term(smb)
2614  !
2615  ! -- deallocate table
2616  call mem_cleanup_table()
2617  end subroutine mem_summary_total
2618 
2619  !> @brief Generic function to clean a memory manager table
2620  !<
2621  subroutine mem_cleanup_table()
2622  ! -- local
2623  ! -- formats
2624  ! -- code
2625  call memtab%table_da()
2626  deallocate (memtab)
2627  nullify (memtab)
2628  end subroutine mem_cleanup_table
2629 
2630  !> @brief Write memory manager memory usage based on the
2631  !! user-specified memory_print_option
2632  !!
2633  !! The total memory usage by data types (int, real, etc.)
2634  !! is written for every simulation.
2635  !<
2636  subroutine mem_write_usage(iout)
2637  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2638  ! -- local
2639  class(memorytype), pointer :: mt
2640  character(len=LENMEMADDRESS), allocatable, dimension(:) :: cunique
2641  ! character(len=LENMEMPATH) :: mem_path
2642  character(len=LENMEMPATH) :: context
2643  character(len=LENCOMPONENTNAME) :: component
2644  character(len=LENCOMPONENTNAME) :: subcomponent
2645  character(len=LENMEMADDRESS) :: context_component
2646  character(LEN=10) :: cunits
2647  type(memorycontaineriteratortype), allocatable :: itr
2648  integer(I4B) :: icomp
2649  integer(I4B) :: ilen
2650  integer(I8B) :: nchars
2651  integer(I8B) :: nlog
2652  integer(I8B) :: nint
2653  integer(I8B) :: nreal
2654  real(dp) :: simbytes
2655  real(dp) :: fact
2656  real(dp) :: rchars
2657  real(dp) :: rlog
2658  real(dp) :: rint
2659  real(dp) :: rreal
2660  real(dp) :: bytes
2661  ! -- formats
2662  ! -- code
2663  !
2664  ! -- Calculate simulation memory allocation
2665  simbytes = (nvalues_astr + &
2666  nvalues_alogical * lgp + &
2667  nvalues_aint * i4b + &
2668  nvalues_adbl * dp)
2669  simbytes = real(simbytes, dp)
2670  !
2671  ! -- calculate factor and memory units
2672  call mem_units(simbytes, fact, cunits)
2673  !
2674  ! -- Write summary table for simulation components
2675  if (iprmem == 1) then
2676  !
2677  ! -- Find unique names of simulation components
2678  call mem_unique_origins(cunique)
2679  call mem_summary_table(iout, size(cunique), cunits)
2680  do icomp = 1, size(cunique)
2681  nchars = 0
2682  nlog = 0
2683  nint = 0
2684  nreal = 0
2685  bytes = dzero
2686  ilen = len_trim(cunique(icomp))
2687  itr = memorystore%iterator()
2688  do while (itr%has_next())
2689  call itr%next()
2690  mt => itr%value()
2691  call split_mem_path(mt%path, component, subcomponent)
2692  context = get_mem_path_context(mt%path)
2693  context_component = trim(context)//component
2694  if (cunique(icomp) /= context_component(1:ilen)) cycle
2695  if (.not. mt%master) cycle
2696  if (mt%memtype(1:6) == 'STRING') then
2697  nchars = nchars + mt%isize * mt%element_size
2698  else if (mt%memtype(1:7) == 'LOGICAL') then
2699  nlog = nlog + mt%isize
2700  else if (mt%memtype(1:7) == 'INTEGER') then
2701  nint = nint + mt%isize
2702  else if (mt%memtype(1:6) == 'DOUBLE') then
2703  nreal = nreal + mt%isize
2704  end if
2705  end do
2706  !
2707  ! -- calculate size of each data type in bytes
2708  rchars = real(nchars, dp) * fact
2709  rlog = real(nlog * lgp, dp) * fact
2710  rint = real(nint * i4b, dp) * fact
2711  rreal = real(nreal * dp, dp) * fact
2712  !
2713  ! -- calculate total storage in bytes
2714  bytes = rchars + rlog + rint + rreal
2715  !
2716  ! -- write data
2717  call mem_summary_line(cunique(icomp), rchars, rlog, rint, rreal, bytes)
2718  end do
2719  call mem_cleanup_table()
2720  end if
2721  !
2722  ! -- Write table with all variables for iprmem == 2
2723  if (iprmem == 2) then
2724  call mem_print_detailed(iout)
2725  end if
2726  !
2727  ! -- Write total memory allocation
2728  call mem_summary_total(iout, simbytes)
2729  end subroutine mem_write_usage
2730 
2731  subroutine mem_print_detailed(iout)
2732  integer(I4B) :: iout
2733  ! local
2734  class(memorytype), pointer :: mt
2735  type(memorycontaineriteratortype), allocatable :: itr
2736 
2737  call mem_detailed_table(iout, memorystore%count())
2738  itr = memorystore%iterator()
2739  do while (itr%has_next())
2740  call itr%next()
2741  mt => itr%value()
2742  call mt%table_entry(memtab)
2743  end do
2744  call mem_cleanup_table()
2745 
2746  end subroutine mem_print_detailed
2747 
2748  !> @brief Sum up virtual memory, i.e. memory
2749  !< that is owned by other processes
2750  function calc_virtual_mem() result(vmem_size)
2751  real(dp) :: vmem_size
2752  ! local
2753  type(memorycontaineriteratortype), allocatable :: itr
2754  type(memorytype), pointer :: mt
2755 
2756  vmem_size = dzero
2757  itr = memorystore%iterator()
2758  do while (itr%has_next())
2759  call itr%next()
2760  mt => itr%value()
2761  if (index(mt%path, "__P") == 1) then
2762  vmem_size = mt%element_size * mt%isize + vmem_size
2763  end if
2764  end do
2765 
2766  end function calc_virtual_mem
2767 
2768  !> @brief Deallocate memory in the memory manager
2769  !<
2770  subroutine mem_da()
2771  ! -- modules
2772  use versionmodule, only: idevelopmode
2773  use inputoutputmodule, only: upcase
2774  ! -- local
2775  class(memorytype), pointer :: mt
2776  character(len=LINELENGTH) :: error_msg
2777  character(len=LENVARNAME) :: ucname
2778  type(memorycontaineriteratortype), allocatable :: itr
2779  ! -- code
2780  itr = memorystore%iterator()
2781  do while (itr%has_next())
2782  call itr%next()
2783  mt => itr%value()
2784  if (idevelopmode == 1) then
2785  !
2786  ! -- check if memory has been deallocated
2787  if (mt%mt_associated() .and. mt%element_size == -1) then
2788  error_msg = trim(adjustl(mt%path))//' '// &
2789  trim(adjustl(mt%name))//' has invalid element size'
2790  call store_error(trim(error_msg))
2791  end if
2792  !
2793  ! -- check if memory has been deallocated
2794  if (mt%mt_associated() .and. mt%isize > 0) then
2795  error_msg = trim(adjustl(mt%path))//' '// &
2796  trim(adjustl(mt%name))//' not deallocated'
2797  call store_error(trim(error_msg))
2798  end if
2799  !
2800  ! -- check case of varname
2801  ucname = mt%name
2802  call upcase(ucname)
2803  if (mt%name /= ucname) then
2804  error_msg = trim(adjustl(mt%path))//' '// &
2805  trim(adjustl(mt%name))//' not upper case'
2806  call store_error(trim(error_msg))
2807  end if
2808  end if
2809  !
2810  ! -- deallocate instance of memory type
2811  deallocate (mt)
2812  end do
2813  call memorystore%clear()
2814  if (count_errors() > 0) then
2815  call store_error('Could not clear memory list.', terminate=.true.)
2816  end if
2817  end subroutine mem_da
2818 
2819  !> @brief Create a array with unique first components from all memory paths.
2820  !! Only the first component of the memory path is evaluated.
2821  !<
2822  subroutine mem_unique_origins(cunique)
2823  ! -- modules
2825  ! -- dummy
2826  character(len=LENMEMADDRESS), allocatable, dimension(:), intent(inout) :: &
2827  cunique !< array with unique first components
2828  ! -- local
2829  class(memorytype), pointer :: mt
2830  character(len=LENMEMPATH) :: context
2831  character(len=LENCOMPONENTNAME) :: component
2832  character(len=LENCOMPONENTNAME) :: subcomponent
2833  character(len=LENMEMADDRESS) :: context_component
2834  type(memorycontaineriteratortype), allocatable :: itr
2835  integer(I4B) :: ipa
2836  ! -- code
2837  !
2838  ! -- initialize cunique
2839  allocate (cunique(0))
2840  !
2841  ! -- find unique origins
2842  itr = memorystore%iterator()
2843  do while (itr%has_next())
2844  call itr%next()
2845  mt => itr%value()
2846  call split_mem_path(mt%path, component, subcomponent)
2847  context = get_mem_path_context(mt%path)
2848  context_component = trim(context)//component
2849  ipa = ifind(cunique, context_component)
2850  if (ipa < 1) then
2851  call expandarray(cunique, 1)
2852  cunique(size(cunique)) = context_component
2853  end if
2854  end do
2855  end subroutine mem_unique_origins
2856 
2857 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.