MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
memorymanagermodule Module Reference

Data Types

interface  mem_allocate
 
interface  mem_checkin
 
interface  mem_reallocate
 
interface  mem_setptr
 
interface  mem_copyptr
 
interface  mem_reassignptr
 
interface  mem_deallocate
 

Functions/Subroutines

subroutine, public get_mem_type (name, mem_path, var_type)
 @ brief Get the variable memory type More...
 
subroutine, public get_mem_rank (name, mem_path, rank)
 @ brief Get the variable rank More...
 
subroutine, public get_mem_elem_size (name, mem_path, size)
 @ brief Get the memory size of a single element of the stored variable More...
 
subroutine, public get_mem_shape (name, mem_path, mem_shape)
 @ brief Get the variable memory shape More...
 
subroutine, public get_isize (name, mem_path, isize)
 @ brief Get the number of elements for this variable More...
 
subroutine, public get_from_memorylist (name, mem_path, mt, found, check)
 @ brief Get a memory type entry from the memory list More...
 
subroutine allocate_error (varname, mem_path, istat, isize)
 Issue allocation error message and stop program execution. More...
 
subroutine allocate_logical (sclr, name, mem_path)
 Allocate a logical scalar. More...
 
subroutine allocate_str (sclr, ilen, name, mem_path)
 Allocate a character string. More...
 
subroutine allocate_str1d (astr1d, ilen, nrow, name, mem_path)
 Allocate a 1-dimensional defined length string array. More...
 
subroutine allocate_charstr1d (acharstr1d, ilen, nrow, name, mem_path)
 Allocate a 1-dimensional array of deferred-length CharacterStringType. More...
 
subroutine allocate_int (sclr, name, mem_path)
 Allocate a integer scalar. More...
 
subroutine allocate_int1d (aint, nrow, name, mem_path)
 Allocate a 1-dimensional integer array. More...
 
subroutine allocate_int2d (aint, ncol, nrow, name, mem_path)
 Allocate a 2-dimensional integer array. More...
 
subroutine allocate_int3d (aint, ncol, nrow, nlay, name, mem_path)
 Allocate a 3-dimensional integer array. More...
 
subroutine allocate_dbl (sclr, name, mem_path)
 Allocate a real scalar. More...
 
subroutine allocate_dbl1d (adbl, nrow, name, mem_path)
 Allocate a 1-dimensional real array. More...
 
subroutine allocate_dbl2d (adbl, ncol, nrow, name, mem_path)
 Allocate a 2-dimensional real array. More...
 
subroutine allocate_dbl3d (adbl, ncol, nrow, nlay, name, mem_path)
 Allocate a 3-dimensional real array. More...
 
subroutine checkin_int1d (aint, name, mem_path, name2, mem_path2)
 Check in an existing 1d integer array with a new address (name + path) More...
 
subroutine checkin_int2d (aint2d, name, mem_path, name2, mem_path2)
 Check in an existing 2d integer array with a new address (name + path) More...
 
subroutine checkin_dbl1d (adbl, name, mem_path, name2, mem_path2)
 Check in an existing 1d double precision array with a new address (name + path) More...
 
subroutine checkin_dbl2d (adbl2d, name, mem_path, name2, mem_path2)
 Check in an existing 2d double precision array with a new address (name + path) More...
 
subroutine checkin_charstr1d (acharstr1d, ilen, name, mem_path, name2, mem_path2)
 Check in an existing 1d CharacterStringType array with a new address (name + path) More...
 
subroutine reallocate_str1d (astr, ilen, nrow, name, mem_path)
 Reallocate a 1-dimensional defined length string array. More...
 
subroutine reallocate_charstr1d (acharstr1d, ilen, nrow, name, mem_path)
 Reallocate a 1-dimensional deferred length string array. More...
 
subroutine reallocate_int1d (aint, nrow, name, mem_path)
 Reallocate a 1-dimensional integer array. More...
 
subroutine reallocate_int2d (aint, ncol, nrow, name, mem_path)
 Reallocate a 2-dimensional integer array. More...
 
subroutine reallocate_dbl1d (adbl, nrow, name, mem_path)
 Reallocate a 1-dimensional real array. More...
 
subroutine reallocate_dbl2d (adbl, ncol, nrow, name, mem_path)
 Reallocate a 2-dimensional real array. More...
 
subroutine setptr_logical (sclr, name, mem_path)
 Set pointer to a logical scalar. More...
 
subroutine setptr_int (sclr, name, mem_path)
 Set pointer to integer scalar. More...
 
subroutine setptr_int1d (aint, name, mem_path)
 Set pointer to 1d integer array. More...
 
subroutine setptr_int2d (aint, name, mem_path)
 Set pointer to 2d integer array. More...
 
subroutine setptr_int3d (aint, name, mem_path)
 Set pointer to 3d integer array. More...
 
subroutine setptr_dbl (sclr, name, mem_path)
 Set pointer to a real scalar. More...
 
subroutine setptr_dbl1d (adbl, name, mem_path)
 Set pointer to a 1d real array. More...
 
subroutine setptr_dbl2d (adbl, name, mem_path)
 Set pointer to a 2d real array. More...
 
subroutine setptr_dbl3d (adbl, name, mem_path)
 Set pointer to a 3d real array. More...
 
subroutine setptr_str (asrt, name, mem_path)
 Set pointer to a string (scalar) More...
 
subroutine setptr_str1d (astr1d, name, mem_path)
 Set pointer to a fixed-length string array. More...
 
subroutine setptr_charstr1d (acharstr1d, name, mem_path)
 Set pointer to an array of CharacterStringType. More...
 
subroutine copyptr_int1d (aint, name, mem_path, mem_path_copy)
 Make a copy of a 1-dimensional integer array. More...
 
subroutine copyptr_int2d (aint, name, mem_path, mem_path_copy)
 Make a copy of a 2-dimensional integer array. More...
 
subroutine copyptr_dbl1d (adbl, name, mem_path, mem_path_copy)
 Make a copy of a 1-dimensional real array. More...
 
subroutine copyptr_dbl2d (adbl, name, mem_path, mem_path_copy)
 Make a copy of a 2-dimensional real array. More...
 
subroutine, public copy_dbl1d (adbl, name, mem_path)
 Copy values from a 1-dimensional real array in the memory. More...
 
subroutine reassignptr_int (sclr, name, mem_path, name_target, mem_path_target)
 Set the pointer for an integer scalar to. More...
 
subroutine reassignptr_int1d (aint, name, mem_path, name_target, mem_path_target)
 Set the pointer for a 1-dimensional integer array to. More...
 
subroutine reassignptr_int2d (aint, name, mem_path, name_target, mem_path_target)
 Set the pointer for a 2-dimensional integer array to. More...
 
subroutine reassignptr_dbl1d (adbl, name, mem_path, name_target, mem_path_target)
 Set the pointer for a 1-dimensional real array to. More...
 
subroutine reassignptr_dbl2d (adbl, name, mem_path, name_target, mem_path_target)
 Set the pointer for a 2-dimensional real array to. More...
 
subroutine deallocate_str (sclr, name, mem_path)
 Deallocate a variable-length character string. More...
 
subroutine deallocate_str1d (astr1d, name, mem_path)
 Deallocate an array of defined-length character strings. More...
 
subroutine deallocate_charstr1d (astr1d, name, mem_path)
 Deallocate an array of deferred-length character strings. More...
 
subroutine deallocate_logical (sclr)
 Deallocate a logical scalar. More...
 
subroutine deallocate_int (sclr)
 Deallocate a integer scalar. More...
 
subroutine deallocate_dbl (sclr)
 Deallocate a real scalar. More...
 
subroutine deallocate_int1d (aint, name, mem_path)
 Deallocate a 1-dimensional integer array. More...
 
subroutine deallocate_int2d (aint, name, mem_path)
 Deallocate a 2-dimensional integer array. More...
 
subroutine deallocate_int3d (aint, name, mem_path)
 Deallocate a 3-dimensional integer array. More...
 
subroutine deallocate_dbl1d (adbl, name, mem_path)
 Deallocate a 1-dimensional real array. More...
 
subroutine deallocate_dbl2d (adbl, name, mem_path)
 Deallocate a 2-dimensional real array. More...
 
subroutine deallocate_dbl3d (adbl, name, mem_path)
 Deallocate a 3-dimensional real array. More...
 
subroutine, public mem_set_print_option (iout, keyword, error_msg)
 Set the memory print option. More...
 
subroutine mem_summary_table (iout, nrows, cunits)
 Create a table if memory_print_option is 'SUMMARY'. More...
 
subroutine mem_detailed_table (iout, nrows)
 Create a table if memory_print_option is 'ALL'. More...
 
subroutine mem_summary_line (component, rchars, rlog, rint, rreal, bytes)
 Write a row for the memory_print_option 'SUMMARY' table. More...
 
subroutine mem_units (bytes, fact, cunits)
 Determine appropriate memory unit and conversion factor. More...
 
subroutine mem_summary_total (iout, bytes)
 Create and fill a table with the total allocated memory. More...
 
subroutine mem_cleanup_table ()
 Generic function to clean a memory manager table. More...
 
subroutine, public mem_write_usage (iout)
 Write memory manager memory usage based on the user-specified memory_print_option. More...
 
subroutine, public mem_print_detailed (iout)
 
real(dp) function calc_virtual_mem ()
 Sum up virtual memory, i.e. memory. More...
 
subroutine, public mem_da ()
 Deallocate memory in the memory manager. More...
 
subroutine mem_unique_origins (cunique)
 Create a array with unique first components from all memory paths. Only the first component of the memory path is evaluated. More...
 

Variables

type(memorylisttype), public memorylist
 
type(tabletype), pointer memtab => null()
 
integer(i8b) nvalues_alogical = 0
 
integer(i8b) nvalues_astr = 0
 
integer(i8b) nvalues_aint = 0
 
integer(i8b) nvalues_adbl = 0
 
integer(i4b) iprmem = 0
 

Function/Subroutine Documentation

◆ allocate_charstr1d()

subroutine memorymanagermodule::allocate_charstr1d ( type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  acharstr1d,
integer(i4b), intent(in)  ilen,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]acharstr1dvariable for allocation
[in]ilenstring length
[in]nrownumber of strings in array
[in]namevariable name
[in]mem_pathpath where the variable is stored

Definition at line 539 of file MemoryManager.f90.

540  type(CharacterStringType), dimension(:), &
541  pointer, contiguous, intent(inout) :: acharstr1d !< variable for allocation
542  integer(I4B), intent(in) :: ilen !< string length
543  integer(I4B), intent(in) :: nrow !< number of strings in array
544  character(len=*), intent(in) :: name !< variable name
545  character(len=*), intent(in) :: mem_path !< path where the variable is stored
546  ! -- local variables
547  character(len=ilen) :: string
548  type(MemoryType), pointer :: mt
549  integer(I4B) :: n
550  integer(I4B) :: istat
551  integer(I4B) :: isize
552  ! -- code
553  !
554  ! -- initialize string
555  string = ''
556  !
557  ! -- check variable name length
558  call mem_check_length(name, lenvarname, "variable")
559  !
560  ! -- calculate isize
561  isize = nrow
562  !
563  ! -- allocate deferred length string array
564  allocate (acharstr1d(nrow), stat=istat, errmsg=errmsg)
565  !
566  ! -- check for error condition
567  if (istat /= 0) then
568  call allocate_error(name, mem_path, istat, isize)
569  end if
570  !
571  ! -- fill deferred length string with empty string
572  do n = 1, nrow
573  acharstr1d(n) = string
574  end do
575  !
576  ! -- update counter
577  nvalues_astr = nvalues_astr + isize
578  !
579  ! -- allocate memory type
580  allocate (mt)
581  !
582  ! -- set memory type
583  mt%acharstr1d => acharstr1d
584  mt%element_size = ilen
585  mt%isize = isize
586  mt%name = name
587  mt%path = mem_path
588  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
589  !
590  ! -- add deferred length character array to the memory manager list
591  call memorylist%add(mt)
592  !
593  ! -- return
594  return

◆ allocate_dbl()

subroutine memorymanagermodule::allocate_dbl ( real(dp), intent(inout), pointer  sclr,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]sclrvariable for allocation
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 781 of file MemoryManager.f90.

782  real(DP), pointer, intent(inout) :: sclr !< variable for allocation
783  character(len=*), intent(in) :: name !< variable name
784  character(len=*), intent(in) :: mem_path !< path where variable is stored
785  ! -- local
786  type(MemoryType), pointer :: mt
787  integer(I4B) :: istat
788  ! -- code
789  !
790  ! -- check variable name length
791  call mem_check_length(name, lenvarname, "variable")
792  !
793  ! -- allocate real scalar
794  allocate (sclr, stat=istat, errmsg=errmsg)
795  if (istat /= 0) then
796  call allocate_error(name, mem_path, istat, 1)
797  end if
798  !
799  ! -- update counter
800  nvalues_aint = nvalues_aint + 1
801  !
802  ! -- allocate memory type
803  allocate (mt)
804  !
805  ! -- set memory type
806  mt%dblsclr => sclr
807  mt%element_size = dp
808  mt%isize = 1
809  mt%name = name
810  mt%path = mem_path
811  write (mt%memtype, "(a)") 'DOUBLE'
812  !
813  ! -- add memory type to the memory list
814  call memorylist%add(mt)
815  !
816  ! -- return
817  return

◆ allocate_dbl1d()

subroutine memorymanagermodule::allocate_dbl1d ( real(dp), dimension(:), intent(inout), pointer, contiguous  adbl,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblvariable for allocation
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 822 of file MemoryManager.f90.

823  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< variable for allocation
824  integer(I4B), intent(in) :: nrow !< number of rows
825  character(len=*), intent(in) :: name !< variable name
826  character(len=*), intent(in) :: mem_path !< path where variable is stored
827  ! -- local
828  type(MemoryType), pointer :: mt
829  integer(I4B) :: istat
830  integer(I4B) :: isize
831  ! -- code
832  !
833  ! -- check the variable name length
834  call mem_check_length(name, lenvarname, "variable")
835  !
836  ! -- set isize
837  isize = nrow
838  !
839  ! -- allocate the real array
840  allocate (adbl(nrow), stat=istat, errmsg=errmsg)
841  if (istat /= 0) then
842  call allocate_error(name, mem_path, istat, isize)
843  end if
844  !
845  ! -- update counter
846  nvalues_adbl = nvalues_adbl + isize
847  !
848  ! -- allocate memory type
849  allocate (mt)
850  !
851  ! -- set memory type
852  mt%adbl1d => adbl
853  mt%element_size = dp
854  mt%isize = isize
855  mt%name = name
856  mt%path = mem_path
857  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', isize
858  !
859  ! -- add memory type to the memory list
860  call memorylist%add(mt)
861  !
862  ! -- return
863  return
Here is the caller graph for this function:

◆ allocate_dbl2d()

subroutine memorymanagermodule::allocate_dbl2d ( real(dp), dimension(:, :), intent(inout), pointer, contiguous  adbl,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblvariable for allocation
[in]ncolnumber of columns
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 868 of file MemoryManager.f90.

869  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< variable for allocation
870  integer(I4B), intent(in) :: ncol !< number of columns
871  integer(I4B), intent(in) :: nrow !< number of rows
872  character(len=*), intent(in) :: name !< variable name
873  character(len=*), intent(in) :: mem_path !< path where variable is stored
874  ! -- local
875  type(MemoryType), pointer :: mt
876  integer(I4B) :: istat
877  integer(I4B) :: isize
878  ! -- code
879  !
880  ! -- check the variable name length
881  call mem_check_length(name, lenvarname, "variable")
882  !
883  ! -- set isize
884  isize = ncol * nrow
885  !
886  ! -- allocate the real array
887  allocate (adbl(ncol, nrow), stat=istat, errmsg=errmsg)
888  if (istat /= 0) then
889  call allocate_error(name, mem_path, istat, isize)
890  end if
891  !
892  ! -- update counter
893  nvalues_adbl = nvalues_adbl + isize
894  !
895  ! -- allocate memory type
896  allocate (mt)
897  !
898  ! -- set memory type
899  mt%adbl2d => adbl
900  mt%element_size = dp
901  mt%isize = isize
902  mt%name = name
903  mt%path = mem_path
904  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
905  !
906  ! -- add memory type to the memory list
907  call memorylist%add(mt)
908  !
909  ! -- return
910  return
Here is the caller graph for this function:

◆ allocate_dbl3d()

subroutine memorymanagermodule::allocate_dbl3d ( real(dp), dimension(:, :, :), intent(inout), pointer, contiguous  adbl,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  nrow,
integer(i4b), intent(in)  nlay,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblvariable for allocation
[in]ncolnumber of columns
[in]nrownumber of rows
[in]nlaynumber of layers
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 915 of file MemoryManager.f90.

916  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< variable for allocation
917  integer(I4B), intent(in) :: ncol !< number of columns
918  integer(I4B), intent(in) :: nrow !< number of rows
919  integer(I4B), intent(in) :: nlay !< number of layers
920  character(len=*), intent(in) :: name !< variable name
921  character(len=*), intent(in) :: mem_path !< path where variable is stored
922  ! -- local
923  type(MemoryType), pointer :: mt
924  integer(I4B) :: istat
925  integer(I4B) :: isize
926  ! -- code
927  !
928  ! -- check the variable name length
929  call mem_check_length(name, lenvarname, "variable")
930  !
931  ! -- set isize
932  isize = ncol * nrow * nlay
933  !
934  ! -- allocate the real array
935  allocate (adbl(ncol, nrow, nlay), stat=istat, errmsg=errmsg)
936  if (istat /= 0) then
937  call allocate_error(name, mem_path, istat, isize)
938  end if
939  !
940  ! -- update the counter
941  nvalues_adbl = nvalues_adbl + isize
942  !
943  ! -- allocate memory type
944  allocate (mt)
945  !
946  ! -- set memory type
947  mt%adbl3d => adbl
948  mt%element_size = dp
949  mt%isize = isize
950  mt%name = name
951  mt%path = mem_path
952  write (mt%memtype, "(a,' (',i0,',',i0,',',i0,')')") 'DOUBLE', ncol, &
953  nrow, nlay
954  !
955  ! -- add memory type to the memory list
956  call memorylist%add(mt)
957  !
958  ! -- return
959  return

◆ allocate_error()

subroutine memorymanagermodule::allocate_error ( character(len=*), intent(in)  varname,
character(len=*), intent(in)  mem_path,
integer(i4b), intent(in)  istat,
integer(i4b), intent(in)  isize 
)
private
Parameters
[in]varnamevariable name
[in]mem_pathpath where the variable is stored
[in]istatstatus code
[in]isizesize of allocation

Definition at line 351 of file MemoryManager.f90.

352  character(len=*), intent(in) :: varname !< variable name
353  character(len=*), intent(in) :: mem_path !< path where the variable is stored
354  integer(I4B), intent(in) :: istat !< status code
355  integer(I4B), intent(in) :: isize !< size of allocation
356  ! -- local
357  character(len=20) :: csize
358  character(len=20) :: cstat
359  ! -- code
360  !
361  ! -- initialize character variables
362  write (csize, '(i0)') isize
363  write (cstat, '(i0)') istat
364  !
365  ! -- create error message
366  errmsg = "Error trying to allocate memory. Path '"//trim(mem_path)// &
367  "' variable name '"//trim(varname)//"' size '"//trim(csize)// &
368  "'. Error message is '"//trim(adjustl(errmsg))// &
369  "'. Status code is "//trim(cstat)//'.'
370  !
371  ! -- store error and stop program execution
372  call store_error(errmsg, terminate=.true.)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ allocate_int()

subroutine memorymanagermodule::allocate_int ( integer(i4b), intent(inout), pointer  sclr,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]sclrvariable for allocation
[in]namevariable name
[in]mem_pathpath where the variable is stored

Definition at line 599 of file MemoryManager.f90.

600  integer(I4B), pointer, intent(inout) :: sclr !< variable for allocation
601  character(len=*), intent(in) :: name !< variable name
602  character(len=*), intent(in) :: mem_path !< path where the variable is stored
603  ! -- local
604  type(MemoryType), pointer :: mt
605  integer(I4B) :: istat
606  ! -- code
607  !
608  ! -- check variable name length
609  call mem_check_length(name, lenvarname, "variable")
610  !
611  ! -- allocate integer scalar
612  allocate (sclr, stat=istat, errmsg=errmsg)
613  if (istat /= 0) then
614  call allocate_error(name, mem_path, istat, 1)
615  end if
616  !
617  ! -- update counter
618  nvalues_aint = nvalues_aint + 1
619  !
620  ! -- allocate memory type
621  allocate (mt)
622  !
623  ! -- set memory type
624  mt%intsclr => sclr
625  mt%element_size = i4b
626  mt%isize = 1
627  mt%name = name
628  mt%path = mem_path
629  write (mt%memtype, "(a)") 'INTEGER'
630  !
631  ! -- add memory type to the memory list
632  call memorylist%add(mt)
633  !
634  ! -- return
635  return

◆ allocate_int1d()

subroutine memorymanagermodule::allocate_int1d ( integer(i4b), dimension(:), intent(inout), pointer, contiguous  aint,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintvariable for allocation
[in]nrowinteger array number of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 640 of file MemoryManager.f90.

641  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< variable for allocation
642  integer(I4B), intent(in) :: nrow !< integer array number of rows
643  character(len=*), intent(in) :: name !< variable name
644  character(len=*), intent(in) :: mem_path !< path where variable is stored
645  ! --local
646  type(MemoryType), pointer :: mt
647  integer(I4B) :: istat
648  integer(I4B) :: isize
649  ! -- code
650  !
651  ! -- check variable name length
652  call mem_check_length(name, lenvarname, "variable")
653  !
654  ! -- set isize
655  isize = nrow
656  !
657  ! -- allocate integer array
658  allocate (aint(nrow), stat=istat, errmsg=errmsg)
659  if (istat /= 0) then
660  call allocate_error(name, mem_path, istat, isize)
661  end if
662  !
663  ! -- update counter
664  nvalues_aint = nvalues_aint + isize
665  !
666  ! -- allocate memory type
667  allocate (mt)
668  !
669  ! -- set memory type
670  mt%aint1d => aint
671  mt%element_size = i4b
672  mt%isize = isize
673  mt%name = name
674  mt%path = mem_path
675  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', isize
676  !
677  ! -- add memory type to the memory list
678  call memorylist%add(mt)
679  !
680  ! -- return
681  return
Here is the caller graph for this function:

◆ allocate_int2d()

subroutine memorymanagermodule::allocate_int2d ( integer(i4b), dimension(:, :), intent(inout), pointer, contiguous  aint,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintvariable for allocation
[in]ncolnumber of columns
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 686 of file MemoryManager.f90.

687  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< variable for allocation
688  integer(I4B), intent(in) :: ncol !< number of columns
689  integer(I4B), intent(in) :: nrow !< number of rows
690  character(len=*), intent(in) :: name !< variable name
691  character(len=*), intent(in) :: mem_path !< path where variable is stored
692  ! -- local
693  type(MemoryType), pointer :: mt
694  integer(I4B) :: istat
695  integer(I4B) :: isize
696  ! -- code
697  !
698  ! -- check the variable name length
699  call mem_check_length(name, lenvarname, "variable")
700  !
701  ! -- set isize
702  isize = ncol * nrow
703  !
704  ! -- allocate the integer array
705  allocate (aint(ncol, nrow), stat=istat, errmsg=errmsg)
706  if (istat /= 0) then
707  call allocate_error(name, mem_path, istat, isize)
708  end if
709  !
710  ! -- update the counter
711  nvalues_aint = nvalues_aint + isize
712  !
713  ! -- allocate memory type
714  allocate (mt)
715  !
716  ! -- set memory type
717  mt%aint2d => aint
718  mt%element_size = i4b
719  mt%isize = isize
720  mt%name = name
721  mt%path = mem_path
722  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
723  !
724  ! -- add memory type to the memory list
725  call memorylist%add(mt)
726  !
727  ! -- return
Here is the caller graph for this function:

◆ allocate_int3d()

subroutine memorymanagermodule::allocate_int3d ( integer(i4b), dimension(:, :, :), intent(inout), pointer, contiguous  aint,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  nrow,
integer(i4b), intent(in)  nlay,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintvariable for allocation
[in]ncolnumber of columns
[in]nrownumber of rows
[in]nlaynumber of layers
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 732 of file MemoryManager.f90.

733  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< variable for allocation
734  integer(I4B), intent(in) :: ncol !< number of columns
735  integer(I4B), intent(in) :: nrow !< number of rows
736  integer(I4B), intent(in) :: nlay !< number of layers
737  character(len=*), intent(in) :: name !< variable name
738  character(len=*), intent(in) :: mem_path !< path where variable is stored
739  ! -- local
740  type(MemoryType), pointer :: mt
741  integer(I4B) :: istat
742  integer(I4B) :: isize
743  ! -- code
744  !
745  ! -- check variable name length
746  call mem_check_length(name, lenvarname, "variable")
747  !
748  ! -- set isize
749  isize = ncol * nrow * nlay
750  !
751  ! -- allocate integer array
752  allocate (aint(ncol, nrow, nlay), stat=istat, errmsg=errmsg)
753  if (istat /= 0) then
754  call allocate_error(name, mem_path, istat, isize)
755  end if
756  !
757  ! -- update counter
758  nvalues_aint = nvalues_aint + isize
759  !
760  ! -- allocate memory type
761  allocate (mt)
762  !
763  ! -- set memory type
764  mt%aint3d => aint
765  mt%element_size = i4b
766  mt%isize = isize
767  mt%name = name
768  mt%path = mem_path
769  write (mt%memtype, "(a,' (',i0,',',i0,',',i0,')')") 'INTEGER', ncol, &
770  nrow, nlay
771  !
772  ! -- add memory type to the memory list
773  call memorylist%add(mt)
774  !
775  ! -- return
776  return

◆ allocate_logical()

subroutine memorymanagermodule::allocate_logical ( logical(lgp), intent(inout), pointer  sclr,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]sclrvariable for allocation
[in]namevariable name
[in]mem_pathpath where the variable is stored

Definition at line 377 of file MemoryManager.f90.

378  logical(LGP), pointer, intent(inout) :: sclr !< variable for allocation
379  character(len=*), intent(in) :: name !< variable name
380  character(len=*), intent(in) :: mem_path !< path where the variable is stored
381  ! -- local
382  integer(I4B) :: istat
383  type(MemoryType), pointer :: mt
384  ! -- code
385  !
386  ! -- check variable name length
387  call mem_check_length(name, lenvarname, "variable")
388  !
389  ! -- allocate the logical scalar
390  allocate (sclr, stat=istat, errmsg=errmsg)
391  if (istat /= 0) then
392  call allocate_error(name, mem_path, istat, 1)
393  end if
394  !
395  ! -- update counter
396  nvalues_alogical = nvalues_alogical + 1
397  !
398  ! -- allocate memory type
399  allocate (mt)
400  !
401  ! -- set memory type
402  mt%logicalsclr => sclr
403  mt%element_size = lgp
404  mt%isize = 1
405  mt%name = name
406  mt%path = mem_path
407  write (mt%memtype, "(a)") 'LOGICAL'
408  !
409  ! -- add memory type to the memory list
410  call memorylist%add(mt)
411  !
412  ! -- return
413  return

◆ allocate_str()

subroutine memorymanagermodule::allocate_str ( character(len=ilen), intent(inout), pointer  sclr,
integer(i4b), intent(in)  ilen,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in]ilenstring length
[in,out]sclrvariable for allocation
[in]namevariable name
[in]mem_pathpath where the variable is stored

Definition at line 418 of file MemoryManager.f90.

419  integer(I4B), intent(in) :: ilen !< string length
420  character(len=ilen), pointer, intent(inout) :: sclr !< variable for allocation
421  character(len=*), intent(in) :: name !< variable name
422  character(len=*), intent(in) :: mem_path !< path where the variable is stored
423  ! -- local
424  integer(I4B) :: istat
425  type(MemoryType), pointer :: mt
426  ! -- format
427  ! -- code
428  !
429  ! -- make sure ilen is greater than 0
430  if (ilen < 1) then
431  errmsg = 'Programming error in allocate_str. ILEN must be greater than 0.'
432  call store_error(errmsg, terminate=.true.)
433  end if
434  !
435  ! -- check variable name length
436  call mem_check_length(name, lenvarname, "variable")
437  !
438  ! -- allocate string
439  allocate (character(len=ilen) :: sclr, stat=istat, errmsg=errmsg)
440  if (istat /= 0) then
441  call allocate_error(name, mem_path, istat, 1)
442  end if
443  !
444  ! -- set sclr to a empty string
445  sclr = ' '
446  !
447  ! -- update counter
448  nvalues_astr = nvalues_astr + ilen
449  !
450  ! -- allocate memory type
451  allocate (mt)
452  !
453  ! -- set memory type
454  mt%strsclr => sclr
455  mt%element_size = ilen
456  mt%isize = 1
457  mt%name = name
458  mt%path = mem_path
459  write (mt%memtype, "(a,' LEN=',i0)") 'STRING', ilen
460  !
461  ! -- add defined length string to the memory manager list
462  call memorylist%add(mt)
463  !
464  ! -- return
465  return

◆ allocate_str1d()

subroutine memorymanagermodule::allocate_str1d ( character(len=ilen), dimension(:), intent(inout), pointer, contiguous  astr1d,
integer(i4b), intent(in)  ilen,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in]ilenstring length
[in,out]astr1dvariable for allocation
[in]nrownumber of strings in array
[in]namevariable name
[in]mem_pathpath where the variable is stored

Definition at line 470 of file MemoryManager.f90.

471  integer(I4B), intent(in) :: ilen !< string length
472  character(len=ilen), dimension(:), &
473  pointer, contiguous, intent(inout) :: astr1d !< variable for allocation
474  integer(I4B), intent(in) :: nrow !< number of strings in array
475  character(len=*), intent(in) :: name !< variable name
476  character(len=*), intent(in) :: mem_path !< path where the variable is stored
477  ! -- local variables
478  type(MemoryType), pointer :: mt
479  character(len=ilen) :: string
480  integer(I4B) :: n
481  integer(I4B) :: istat
482  integer(I4B) :: isize
483  ! -- code
484  !
485  ! -- initialize string
486  string = ''
487  !
488  ! -- make sure ilen is greater than 0
489  if (ilen < 1) then
490  errmsg = 'Programming error in allocate_str1d. '// &
491  'ILEN must be greater than 0.'
492  call store_error(errmsg, terminate=.true.)
493  end if
494  !
495  ! -- check variable name length
496  call mem_check_length(name, lenvarname, "variable")
497  !
498  ! -- calculate isize
499  isize = nrow
500  !
501  ! -- allocate defined length string array
502  allocate (character(len=ilen) :: astr1d(nrow), stat=istat, errmsg=errmsg)
503  !
504  ! -- check for error condition
505  if (istat /= 0) then
506  call allocate_error(name, mem_path, istat, isize)
507  end if
508  !
509  ! -- fill deferred length string with empty string
510  do n = 1, nrow
511  astr1d(n) = string
512  end do
513  !
514  ! -- update counter
515  nvalues_astr = nvalues_astr + isize
516  !
517  ! -- allocate memory type
518  allocate (mt)
519  !
520  ! -- set memory type
521  ! this does not work with gfortran 11.3 and 12.1
522  ! so we have to disable the pointing to astr1d
523  ! mt%astr1d => astr1d
524  mt%element_size = ilen
525  mt%isize = isize
526  mt%name = name
527  mt%path = mem_path
528  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
529  !
530  ! -- add deferred length character array to the memory manager list
531  call memorylist%add(mt)
532  !
533  ! -- return
534  return

◆ calc_virtual_mem()

real(dp) function memorymanagermodule::calc_virtual_mem
private

Definition at line 2942 of file MemoryManager.f90.

2943  real(DP) :: vmem_size
2944  ! local
2945  integer(I4B) :: i
2946  type(MemoryType), pointer :: mt
2947 
2948  vmem_size = dzero
2949  do i = 1, memorylist%count()
2950  mt => memorylist%Get(i)
2951  if (index(mt%path, "__P") == 1) then
2952  vmem_size = mt%element_size * mt%isize + vmem_size
2953  end if
2954  end do
2955 
Here is the caller graph for this function:

◆ checkin_charstr1d()

subroutine memorymanagermodule::checkin_charstr1d ( type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  acharstr1d,
integer(i4b), intent(in)  ilen,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name2,
character(len=*), intent(in)  mem_path2 
)
private
Parameters
[in,out]acharstr1dthe existing array
[in]namenew variable name
[in]mem_pathnew path where variable is stored
[in]name2existing variable name
[in]mem_path2existing path where variable is stored

Definition at line 1134 of file MemoryManager.f90.

1135  type(CharacterStringType), dimension(:), &
1136  pointer, contiguous, intent(inout) :: acharstr1d !< the existing array
1137  integer(I4B), intent(in) :: ilen
1138  character(len=*), intent(in) :: name !< new variable name
1139  character(len=*), intent(in) :: mem_path !< new path where variable is stored
1140  character(len=*), intent(in) :: name2 !< existing variable name
1141  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
1142  ! --local
1143  type(MemoryType), pointer :: mt
1144  integer(I4B) :: isize
1145  ! -- code
1146  !
1147  ! -- check variable name length
1148  call mem_check_length(name, lenvarname, "variable")
1149  !
1150  ! -- set isize
1151  isize = size(acharstr1d)
1152  !
1153  ! -- allocate memory type
1154  allocate (mt)
1155  !
1156  ! -- set memory type
1157  mt%acharstr1d => acharstr1d
1158  mt%element_size = ilen
1159  mt%isize = isize
1160  mt%name = name
1161  mt%path = mem_path
1162  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, isize
1163  !
1164  ! -- set master information
1165  mt%master = .false.
1166  mt%mastername = name2
1167  mt%masterPath = mem_path2
1168  !
1169  ! -- add memory type to the memory list
1170  call memorylist%add(mt)
1171  !
1172  ! -- return
1173  return

◆ checkin_dbl1d()

subroutine memorymanagermodule::checkin_dbl1d ( real(dp), dimension(:), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name2,
character(len=*), intent(in)  mem_path2 
)
private
Parameters
[in,out]adblthe existing array
[in]namenew variable name
[in]mem_pathnew path where variable is stored
[in]name2existing variable name
[in]mem_path2existing path where variable is stored

Definition at line 1049 of file MemoryManager.f90.

1050  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< the existing array
1051  character(len=*), intent(in) :: name !< new variable name
1052  character(len=*), intent(in) :: mem_path !< new path where variable is stored
1053  character(len=*), intent(in) :: name2 !< existing variable name
1054  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
1055  ! -- local
1056  type(MemoryType), pointer :: mt
1057  integer(I4B) :: isize
1058  ! -- code
1059  !
1060  ! -- check the variable name length
1061  call mem_check_length(name, lenvarname, "variable")
1062  !
1063  ! -- set isize
1064  isize = size(adbl)
1065  !
1066  ! -- allocate memory type
1067  allocate (mt)
1068  !
1069  ! -- set memory type
1070  mt%adbl1d => adbl
1071  mt%element_size = dp
1072  mt%isize = isize
1073  mt%name = name
1074  mt%path = mem_path
1075  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', isize
1076  !
1077  ! -- set master information
1078  mt%master = .false.
1079  mt%mastername = name2
1080  mt%masterPath = mem_path2
1081  !
1082  ! -- add memory type to the memory list
1083  call memorylist%add(mt)
1084  !
1085  ! -- return
1086  return

◆ checkin_dbl2d()

subroutine memorymanagermodule::checkin_dbl2d ( real(dp), dimension(:, :), intent(inout), pointer, contiguous  adbl2d,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name2,
character(len=*), intent(in)  mem_path2 
)
private
Parameters
[in,out]adbl2dthe existing 2d array
[in]namenew variable name
[in]mem_pathnew path where variable is stored
[in]name2existing variable name
[in]mem_path2existing path where variable is stored

Definition at line 1091 of file MemoryManager.f90.

1092  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl2d !< the existing 2d array
1093  character(len=*), intent(in) :: name !< new variable name
1094  character(len=*), intent(in) :: mem_path !< new path where variable is stored
1095  character(len=*), intent(in) :: name2 !< existing variable name
1096  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
1097  ! -- local
1098  type(MemoryType), pointer :: mt
1099  integer(I4B) :: ncol, nrow, isize
1100  ! -- code
1101  !
1102  ! -- check the variable name length
1103  call mem_check_length(name, lenvarname, "variable")
1104  !
1105  ! -- set isize
1106  ncol = size(adbl2d, dim=1)
1107  nrow = size(adbl2d, dim=2)
1108  isize = ncol * nrow
1109  !
1110  ! -- allocate memory type
1111  allocate (mt)
1112  !
1113  ! -- set memory type
1114  mt%adbl2d => adbl2d
1115  mt%isize = isize
1116  mt%name = name
1117  mt%path = mem_path
1118  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
1119  !
1120  ! -- set master information
1121  mt%master = .false.
1122  mt%mastername = name2
1123  mt%masterPath = mem_path2
1124  !
1125  ! -- add memory type to the memory list
1126  call memorylist%add(mt)
1127  !
1128  ! -- return
1129  return

◆ checkin_int1d()

subroutine memorymanagermodule::checkin_int1d ( integer(i4b), dimension(:), intent(in), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name2,
character(len=*), intent(in)  mem_path2 
)
private
Parameters
[in]aintthe existing array
[in]namenew variable name
[in]mem_pathnew path where variable is stored
[in]name2existing variable name
[in]mem_path2existing path where variable is stored

Definition at line 964 of file MemoryManager.f90.

965  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: aint !< the existing array
966  character(len=*), intent(in) :: name !< new variable name
967  character(len=*), intent(in) :: mem_path !< new path where variable is stored
968  character(len=*), intent(in) :: name2 !< existing variable name
969  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
970  ! --local
971  type(MemoryType), pointer :: mt
972  integer(I4B) :: isize
973  ! -- code
974  !
975  ! -- check variable name length
976  call mem_check_length(name, lenvarname, "variable")
977  !
978  ! -- set isize
979  isize = size(aint)
980  !
981  ! -- allocate memory type
982  allocate (mt)
983  !
984  ! -- set memory type
985  mt%aint1d => aint
986  mt%element_size = i4b
987  mt%isize = isize
988  mt%name = name
989  mt%path = mem_path
990  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', isize
991  !
992  ! -- set master information
993  mt%master = .false.
994  mt%mastername = name2
995  mt%masterPath = mem_path2
996  !
997  ! -- add memory type to the memory list
998  call memorylist%add(mt)
999  !
1000  ! -- return
1001  return

◆ checkin_int2d()

subroutine memorymanagermodule::checkin_int2d ( integer(i4b), dimension(:, :), intent(inout), pointer, contiguous  aint2d,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name2,
character(len=*), intent(in)  mem_path2 
)
private
Parameters
[in,out]aint2dthe existing 2d array
[in]namenew variable name
[in]mem_pathnew path where variable is stored
[in]name2existing variable name
[in]mem_path2existing path where variable is stored

Definition at line 1006 of file MemoryManager.f90.

1007  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint2d !< the existing 2d array
1008  character(len=*), intent(in) :: name !< new variable name
1009  character(len=*), intent(in) :: mem_path !< new path where variable is stored
1010  character(len=*), intent(in) :: name2 !< existing variable name
1011  character(len=*), intent(in) :: mem_path2 !< existing path where variable is stored
1012  ! -- local
1013  type(MemoryType), pointer :: mt
1014  integer(I4B) :: ncol, nrow, isize
1015  ! -- code
1016  !
1017  ! -- check the variable name length
1018  call mem_check_length(name, lenvarname, "variable")
1019  !
1020  ! -- set isize
1021  ncol = size(aint2d, dim=1)
1022  nrow = size(aint2d, dim=2)
1023  isize = ncol * nrow
1024  !
1025  ! -- allocate memory type
1026  allocate (mt)
1027  !
1028  ! -- set memory type
1029  mt%aint2d => aint2d
1030  mt%isize = isize
1031  mt%name = name
1032  mt%path = mem_path
1033  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
1034  !
1035  ! -- set master information
1036  mt%master = .false.
1037  mt%mastername = name2
1038  mt%masterPath = mem_path2
1039  !
1040  ! -- add memory type to the memory list
1041  call memorylist%add(mt)
1042  !
1043  ! -- return
1044  return

◆ copy_dbl1d()

subroutine, public memorymanagermodule::copy_dbl1d ( real(dp), dimension(:), intent(inout)  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
Parameters
[in,out]adbltarget array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1887 of file MemoryManager.f90.

1888  real(DP), dimension(:), intent(inout) :: adbl !< target array
1889  character(len=*), intent(in) :: name !< variable name
1890  character(len=*), intent(in) :: mem_path !< path where variable is stored
1891  ! -- local
1892  type(MemoryType), pointer :: mt
1893  logical(LGP) :: found
1894  integer(I4B) :: n
1895  ! -- code
1896  call get_from_memorylist(name, mem_path, mt, found)
1897  do n = 1, size(mt%adbl1d)
1898  adbl(n) = mt%adbl1d(n)
1899  end do
1900  !
1901  ! -- return
1902  return
Here is the call graph for this function:

◆ copyptr_dbl1d()

subroutine memorymanagermodule::copyptr_dbl1d ( real(dp), dimension(:), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in), optional  mem_path_copy 
)
private
Parameters
[in,out]adblreturned copy of 1d real array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]mem_path_copyoptional path where the copy will be stored, if passed then the copy is added to the memory manager

Definition at line 1818 of file MemoryManager.f90.

1819  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< returned copy of 1d real array
1820  character(len=*), intent(in) :: name !< variable name
1821  character(len=*), intent(in) :: mem_path !< path where variable is stored
1822  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1823  !! if passed then the copy is added to the
1824  !! memory manager
1825  ! -- local
1826  type(MemoryType), pointer :: mt
1827  logical(LGP) :: found
1828  integer(I4B) :: n
1829  ! -- code
1830  call get_from_memorylist(name, mem_path, mt, found)
1831  adbl => null()
1832  ! -- check the copy into the memory manager
1833  if (present(mem_path_copy)) then
1834  call allocate_dbl1d(adbl, size(mt%adbl1d), mt%name, mem_path_copy)
1835  ! -- create a local copy
1836  else
1837  allocate (adbl(size(mt%adbl1d)))
1838  end if
1839  do n = 1, size(mt%adbl1d)
1840  adbl(n) = mt%adbl1d(n)
1841  end do
1842  !
1843  ! -- return
1844  return

◆ copyptr_dbl2d()

subroutine memorymanagermodule::copyptr_dbl2d ( real(dp), dimension(:, :), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in), optional  mem_path_copy 
)
private
Parameters
[in,out]adblreturned copy of 2d real array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]mem_path_copyoptional path where the copy will be stored, if passed then the copy is added to the memory manager

Definition at line 1849 of file MemoryManager.f90.

1850  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< returned copy of 2d real array
1851  character(len=*), intent(in) :: name !< variable name
1852  character(len=*), intent(in) :: mem_path !< path where variable is stored
1853  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1854  !! if passed then the copy is added to the
1855  !! memory manager
1856  ! -- local
1857  type(MemoryType), pointer :: mt
1858  logical(LGP) :: found
1859  integer(I4B) :: i
1860  integer(I4B) :: j
1861  integer(I4B) :: ncol
1862  integer(I4B) :: nrow
1863  ! -- code
1864  call get_from_memorylist(name, mem_path, mt, found)
1865  adbl => null()
1866  ncol = size(mt%adbl2d, dim=1)
1867  nrow = size(mt%adbl2d, dim=2)
1868  ! -- check the copy into the memory manager
1869  if (present(mem_path_copy)) then
1870  call allocate_dbl2d(adbl, ncol, nrow, mt%name, mem_path_copy)
1871  ! -- create a local copy
1872  else
1873  allocate (adbl(ncol, nrow))
1874  end if
1875  do i = 1, nrow
1876  do j = 1, ncol
1877  adbl(j, i) = mt%adbl2d(j, i)
1878  end do
1879  end do
1880  !
1881  ! -- return
1882  return

◆ copyptr_int1d()

subroutine memorymanagermodule::copyptr_int1d ( integer(i4b), dimension(:), intent(inout), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in), optional  mem_path_copy 
)
private
Parameters
[in,out]aintreturned copy of 1d integer array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]mem_path_copyoptional path where the copy will be stored, if passed then the copy is added to the memory manager

Definition at line 1749 of file MemoryManager.f90.

1750  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< returned copy of 1d integer array
1751  character(len=*), intent(in) :: name !< variable name
1752  character(len=*), intent(in) :: mem_path !< path where variable is stored
1753  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1754  !! if passed then the copy is added to the
1755  !! memory manager
1756  ! -- local
1757  type(MemoryType), pointer :: mt
1758  logical(LGP) :: found
1759  integer(I4B) :: n
1760  ! -- code
1761  call get_from_memorylist(name, mem_path, mt, found)
1762  aint => null()
1763  ! -- check the copy into the memory manager
1764  if (present(mem_path_copy)) then
1765  call allocate_int1d(aint, size(mt%aint1d), mt%name, mem_path_copy)
1766  ! -- create a local copy
1767  else
1768  allocate (aint(size(mt%aint1d)))
1769  end if
1770  do n = 1, size(mt%aint1d)
1771  aint(n) = mt%aint1d(n)
1772  end do
1773  !
1774  ! -- return
1775  return

◆ copyptr_int2d()

subroutine memorymanagermodule::copyptr_int2d ( integer(i4b), dimension(:, :), intent(inout), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in), optional  mem_path_copy 
)
private
Parameters
[in,out]aintreturned copy of 2d integer array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]mem_path_copyoptional path where the copy will be stored, if passed then the copy is added to the memory manager

Definition at line 1780 of file MemoryManager.f90.

1781  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< returned copy of 2d integer array
1782  character(len=*), intent(in) :: name !< variable name
1783  character(len=*), intent(in) :: mem_path !< path where variable is stored
1784  character(len=*), intent(in), optional :: mem_path_copy !< optional path where the copy will be stored,
1785  !! if passed then the copy is added to the
1786  !! memory manager
1787  ! -- local
1788  type(MemoryType), pointer :: mt
1789  logical(LGP) :: found
1790  integer(I4B) :: i
1791  integer(I4B) :: j
1792  integer(I4B) :: ncol
1793  integer(I4B) :: nrow
1794  ! -- code
1795  call get_from_memorylist(name, mem_path, mt, found)
1796  aint => null()
1797  ncol = size(mt%aint2d, dim=1)
1798  nrow = size(mt%aint2d, dim=2)
1799  ! -- check the copy into the memory manager
1800  if (present(mem_path_copy)) then
1801  call allocate_int2d(aint, ncol, nrow, mt%name, mem_path_copy)
1802  ! -- create a local copy
1803  else
1804  allocate (aint(ncol, nrow))
1805  end if
1806  do i = 1, nrow
1807  do j = 1, ncol
1808  aint(j, i) = mt%aint2d(j, i)
1809  end do
1810  end do
1811  !
1812  ! -- return
1813  return

◆ deallocate_charstr1d()

subroutine memorymanagermodule::deallocate_charstr1d ( type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  astr1d,
character(len=*), intent(in), optional  name,
character(len=*), intent(in), optional  mem_path 
)
private
Parameters
[in,out]astr1darray of strings
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 2167 of file MemoryManager.f90.

2168  type(CharacterStringType), dimension(:), pointer, contiguous, &
2169  intent(inout) :: astr1d !< array of strings
2170  character(len=*), optional, intent(in) :: name !< variable name
2171  character(len=*), optional, intent(in) :: mem_path !< path where variable is stored
2172  ! -- local
2173  type(MemoryType), pointer :: mt
2174  logical(LGP) :: found
2175  integer(I4B) :: ipos
2176  ! -- code
2177  !
2178  ! -- process optional variables
2179  found = .false.
2180  if (present(name) .and. present(mem_path)) then
2181  call get_from_memorylist(name, mem_path, mt, found)
2182  nullify (mt%acharstr1d)
2183  else
2184  do ipos = 1, memorylist%count()
2185  mt => memorylist%Get(ipos)
2186  if (associated(mt%acharstr1d, astr1d)) then
2187  nullify (mt%acharstr1d)
2188  found = .true.
2189  exit
2190  end if
2191  end do
2192  end if
2193  if (.not. found .and. size(astr1d) > 0) then
2194  call store_error('programming error in deallocate_charstr1d', &
2195  terminate=.true.)
2196  else
2197  if (mt%master) then
2198  deallocate (astr1d)
2199  else
2200  nullify (astr1d)
2201  end if
2202  end if
2203  !
2204  ! -- return
2205  return

◆ deallocate_dbl()

subroutine memorymanagermodule::deallocate_dbl ( real(dp), intent(inout), pointer  sclr)
private
Parameters
[in,out]sclrreal variable to deallocate

Definition at line 2275 of file MemoryManager.f90.

2276  real(DP), pointer, intent(inout) :: sclr !< real variable to deallocate
2277  ! -- local
2278  class(MemoryType), pointer :: mt
2279  logical(LGP) :: found
2280  integer(I4B) :: ipos
2281  ! -- code
2282  found = .false.
2283  do ipos = 1, memorylist%count()
2284  mt => memorylist%Get(ipos)
2285  if (associated(mt%dblsclr, sclr)) then
2286  nullify (mt%dblsclr)
2287  found = .true.
2288  exit
2289  end if
2290  end do
2291  if (.not. found) then
2292  call store_error('Programming error in deallocate_dbl.', terminate=.true.)
2293  else
2294  if (mt%master) then
2295  deallocate (sclr)
2296  else
2297  nullify (sclr)
2298  end if
2299  end if
2300  !
2301  ! -- return
2302  return

◆ deallocate_dbl1d()

subroutine memorymanagermodule::deallocate_dbl1d ( real(dp), dimension(:), intent(inout), pointer, contiguous  adbl,
character(len=*), optional  name,
character(len=*), optional  mem_path 
)
private
Parameters
[in,out]adbl1d real array to deallocate
namevariable name
mem_pathpath where variable is stored

Definition at line 2430 of file MemoryManager.f90.

2431  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< 1d real array to deallocate
2432  character(len=*), optional :: name !< variable name
2433  character(len=*), optional :: mem_path !< path where variable is stored
2434  ! -- local
2435  type(MemoryType), pointer :: mt
2436  logical(LGP) :: found
2437  integer(I4B) :: ipos
2438  ! -- code
2439  !
2440  ! -- process optional variables
2441  found = .false.
2442  if (present(name) .and. present(mem_path)) then
2443  call get_from_memorylist(name, mem_path, mt, found)
2444  nullify (mt%adbl1d)
2445  else
2446  do ipos = 1, memorylist%count()
2447  mt => memorylist%Get(ipos)
2448  if (associated(mt%adbl1d, adbl)) then
2449  nullify (mt%adbl1d)
2450  found = .true.
2451  exit
2452  end if
2453  end do
2454  end if
2455  if (.not. found .and. size(adbl) > 0) then
2456  call store_error('programming error in deallocate_dbl1d', terminate=.true.)
2457  else
2458  if (mt%master) then
2459  deallocate (adbl)
2460  else
2461  nullify (adbl)
2462  end if
2463  end if
2464  !
2465  ! -- return
2466  return

◆ deallocate_dbl2d()

subroutine memorymanagermodule::deallocate_dbl2d ( real(dp), dimension(:, :), intent(inout), pointer, contiguous  adbl,
character(len=*), optional  name,
character(len=*), optional  mem_path 
)
private
Parameters
[in,out]adbl2d real array to deallocate
namevariable name
mem_pathpath where variable is stored

Definition at line 2471 of file MemoryManager.f90.

2472  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< 2d real array to deallocate
2473  character(len=*), optional :: name !< variable name
2474  character(len=*), optional :: mem_path !< path where variable is stored
2475  ! -- local
2476  type(MemoryType), pointer :: mt
2477  logical(LGP) :: found
2478  integer(I4B) :: ipos
2479  ! -- code
2480  !
2481  ! -- process optional variables
2482  found = .false.
2483  if (present(name) .and. present(mem_path)) then
2484  call get_from_memorylist(name, mem_path, mt, found)
2485  nullify (mt%adbl2d)
2486  else
2487  do ipos = 1, memorylist%count()
2488  mt => memorylist%Get(ipos)
2489  if (associated(mt%adbl2d, adbl)) then
2490  nullify (mt%adbl2d)
2491  found = .true.
2492  exit
2493  end if
2494  end do
2495  end if
2496  if (.not. found .and. size(adbl) > 0) then
2497  call store_error('programming error in deallocate_dbl2d', terminate=.true.)
2498  else
2499  if (mt%master) then
2500  deallocate (adbl)
2501  else
2502  nullify (adbl)
2503  end if
2504  end if
2505  !
2506  ! -- return
2507  return

◆ deallocate_dbl3d()

subroutine memorymanagermodule::deallocate_dbl3d ( real(dp), dimension(:, :, :), intent(inout), pointer, contiguous  adbl,
character(len=*), optional  name,
character(len=*), optional  mem_path 
)
private
Parameters
[in,out]adbl3d real array to deallocate
namevariable name
mem_pathpath where variable is stored

Definition at line 2512 of file MemoryManager.f90.

2513  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< 3d real array to deallocate
2514  character(len=*), optional :: name !< variable name
2515  character(len=*), optional :: mem_path !< path where variable is stored
2516  ! -- local
2517  type(MemoryType), pointer :: mt
2518  logical(LGP) :: found
2519  integer(I4B) :: ipos
2520  ! -- code
2521  !
2522  ! -- process optional variables
2523  found = .false.
2524  if (present(name) .and. present(mem_path)) then
2525  call get_from_memorylist(name, mem_path, mt, found)
2526  nullify (mt%adbl3d)
2527  else
2528  do ipos = 1, memorylist%count()
2529  mt => memorylist%Get(ipos)
2530  if (associated(mt%adbl3d, adbl)) then
2531  nullify (mt%adbl3d)
2532  found = .true.
2533  exit
2534  end if
2535  end do
2536  end if
2537  if (.not. found .and. size(adbl) > 0) then
2538  call store_error('programming error in deallocate_dbl3d', terminate=.true.)
2539  else
2540  if (mt%master) then
2541  deallocate (adbl)
2542  else
2543  nullify (adbl)
2544  end if
2545  end if
2546  !
2547  ! -- return
2548  return

◆ deallocate_int()

subroutine memorymanagermodule::deallocate_int ( integer(i4b), intent(inout), pointer  sclr)
private
Parameters
[in,out]sclrinteger variable to deallocate

Definition at line 2243 of file MemoryManager.f90.

2244  integer(I4B), pointer, intent(inout) :: sclr !< integer variable to deallocate
2245  ! -- local
2246  class(MemoryType), pointer :: mt
2247  logical(LGP) :: found
2248  integer(I4B) :: ipos
2249  ! -- code
2250  found = .false.
2251  do ipos = 1, memorylist%count()
2252  mt => memorylist%Get(ipos)
2253  if (associated(mt%intsclr, sclr)) then
2254  nullify (mt%intsclr)
2255  found = .true.
2256  exit
2257  end if
2258  end do
2259  if (.not. found) then
2260  call store_error('Programming error in deallocate_int.', terminate=.true.)
2261  else
2262  if (mt%master) then
2263  deallocate (sclr)
2264  else
2265  nullify (sclr)
2266  end if
2267  end if
2268  !
2269  ! -- return
2270  return

◆ deallocate_int1d()

subroutine memorymanagermodule::deallocate_int1d ( integer(i4b), dimension(:), intent(inout), pointer, contiguous  aint,
character(len=*), optional  name,
character(len=*), optional  mem_path 
)
private
Parameters
[in,out]aint1d integer array to deallocate
namevariable name
mem_pathpath where variable is stored

Definition at line 2307 of file MemoryManager.f90.

2308  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< 1d integer array to deallocate
2309  character(len=*), optional :: name !< variable name
2310  character(len=*), optional :: mem_path !< path where variable is stored
2311  ! -- local
2312  type(MemoryType), pointer :: mt
2313  logical(LGP) :: found
2314  integer(I4B) :: ipos
2315  ! -- code
2316  !
2317  ! -- process optional variables
2318  found = .false.
2319  if (present(name) .and. present(mem_path)) then
2320  call get_from_memorylist(name, mem_path, mt, found)
2321  nullify (mt%aint1d)
2322  else
2323  do ipos = 1, memorylist%count()
2324  mt => memorylist%Get(ipos)
2325  if (associated(mt%aint1d, aint)) then
2326  nullify (mt%aint1d)
2327  found = .true.
2328  exit
2329  end if
2330  end do
2331  end if
2332  if (.not. found .and. size(aint) > 0) then
2333  call store_error('programming error in deallocate_int1d', terminate=.true.)
2334  else
2335  if (mt%master) then
2336  deallocate (aint)
2337  else
2338  nullify (aint)
2339  end if
2340  end if
2341  !
2342  ! -- return
2343  return

◆ deallocate_int2d()

subroutine memorymanagermodule::deallocate_int2d ( integer(i4b), dimension(:, :), intent(inout), pointer, contiguous  aint,
character(len=*), optional  name,
character(len=*), optional  mem_path 
)
private
Parameters
[in,out]aint2d integer array to deallocate
namevariable name
mem_pathpath where variable is stored

Definition at line 2348 of file MemoryManager.f90.

2349  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< 2d integer array to deallocate
2350  character(len=*), optional :: name !< variable name
2351  character(len=*), optional :: mem_path !< path where variable is stored
2352  ! -- local
2353  type(MemoryType), pointer :: mt
2354  logical(LGP) :: found
2355  integer(I4B) :: ipos
2356  ! -- code
2357  !
2358  ! -- process optional variables
2359  found = .false.
2360  if (present(name) .and. present(mem_path)) then
2361  call get_from_memorylist(name, mem_path, mt, found)
2362  nullify (mt%aint2d)
2363  else
2364  do ipos = 1, memorylist%count()
2365  mt => memorylist%Get(ipos)
2366  if (associated(mt%aint2d, aint)) then
2367  nullify (mt%aint2d)
2368  found = .true.
2369  exit
2370  end if
2371  end do
2372  end if
2373  if (.not. found .and. size(aint) > 0) then
2374  call store_error('programming error in deallocate_int2d', terminate=.true.)
2375  else
2376  if (mt%master) then
2377  deallocate (aint)
2378  else
2379  nullify (aint)
2380  end if
2381  end if
2382  !
2383  ! -- return
2384  return

◆ deallocate_int3d()

subroutine memorymanagermodule::deallocate_int3d ( integer(i4b), dimension(:, :, :), intent(inout), pointer, contiguous  aint,
character(len=*), optional  name,
character(len=*), optional  mem_path 
)
private
Parameters
[in,out]aint3d integer array to deallocate
namevariable name
mem_pathpath where variable is stored

Definition at line 2389 of file MemoryManager.f90.

2390  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< 3d integer array to deallocate
2391  character(len=*), optional :: name !< variable name
2392  character(len=*), optional :: mem_path !< path where variable is stored
2393  ! -- local
2394  type(MemoryType), pointer :: mt
2395  logical(LGP) :: found
2396  integer(I4B) :: ipos
2397  ! -- code
2398  !
2399  ! -- process optional variables
2400  found = .false.
2401  if (present(name) .and. present(mem_path)) then
2402  call get_from_memorylist(name, mem_path, mt, found)
2403  nullify (mt%aint3d)
2404  else
2405  do ipos = 1, memorylist%count()
2406  mt => memorylist%Get(ipos)
2407  if (associated(mt%aint3d, aint)) then
2408  nullify (mt%aint3d)
2409  found = .true.
2410  exit
2411  end if
2412  end do
2413  end if
2414  if (.not. found .and. size(aint) > 0) then
2415  call store_error('programming error in deallocate_int3d', terminate=.true.)
2416  else
2417  if (mt%master) then
2418  deallocate (aint)
2419  else
2420  nullify (aint)
2421  end if
2422  end if
2423  !
2424  ! -- return
2425  return

◆ deallocate_logical()

subroutine memorymanagermodule::deallocate_logical ( logical(lgp), intent(inout), pointer  sclr)
private
Parameters
[in,out]sclrlogical scalar to deallocate

Definition at line 2210 of file MemoryManager.f90.

2211  logical(LGP), pointer, intent(inout) :: sclr !< logical scalar to deallocate
2212  ! -- local
2213  class(MemoryType), pointer :: mt
2214  logical(LGP) :: found
2215  integer(I4B) :: ipos
2216  ! -- code
2217  found = .false.
2218  do ipos = 1, memorylist%count()
2219  mt => memorylist%Get(ipos)
2220  if (associated(mt%logicalsclr, sclr)) then
2221  nullify (mt%logicalsclr)
2222  found = .true.
2223  exit
2224  end if
2225  end do
2226  if (.not. found) then
2227  call store_error('programming error in deallocate_logical', &
2228  terminate=.true.)
2229  else
2230  if (mt%master) then
2231  deallocate (sclr)
2232  else
2233  nullify (sclr)
2234  end if
2235  end if
2236  !
2237  ! -- return
2238  return

◆ deallocate_str()

subroutine memorymanagermodule::deallocate_str ( character(len=*), intent(inout), pointer  sclr,
character(len=*), intent(in), optional  name,
character(len=*), intent(in), optional  mem_path 
)
private
Parameters
[in,out]sclrpointer to string
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 2085 of file MemoryManager.f90.

2086  character(len=*), pointer, intent(inout) :: sclr !< pointer to string
2087  character(len=*), intent(in), optional :: name !< variable name
2088  character(len=*), intent(in), optional :: mem_path !< path where variable is stored
2089  ! -- local
2090  type(MemoryType), pointer :: mt
2091  logical(LGP) :: found
2092  integer(I4B) :: ipos
2093  ! -- code
2094  found = .false.
2095  if (present(name) .and. present(mem_path)) then
2096  call get_from_memorylist(name, mem_path, mt, found)
2097  nullify (mt%strsclr)
2098  else
2099  do ipos = 1, memorylist%count()
2100  mt => memorylist%Get(ipos)
2101  if (associated(mt%strsclr, sclr)) then
2102  nullify (mt%strsclr)
2103  found = .true.
2104  exit
2105  end if
2106  end do
2107  end if
2108  if (.not. found) then
2109  call store_error('Programming error in deallocate_str.', terminate=.true.)
2110  else
2111  if (mt%master) then
2112  deallocate (sclr)
2113  else
2114  nullify (sclr)
2115  end if
2116  end if
2117  !
2118  ! -- return
2119  return

◆ deallocate_str1d()

subroutine memorymanagermodule::deallocate_str1d ( character(len=*), dimension(:), intent(inout), pointer, contiguous  astr1d,
character(len=*), intent(in), optional  name,
character(len=*), intent(in), optional  mem_path 
)
private
Parameters
[in,out]astr1darray of strings
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 2125 of file MemoryManager.f90.

2126  character(len=*), dimension(:), pointer, contiguous, intent(inout) :: astr1d !< array of strings
2127  character(len=*), optional, intent(in) :: name !< variable name
2128  character(len=*), optional, intent(in) :: mem_path !< path where variable is stored
2129  ! -- local
2130  type(MemoryType), pointer :: mt
2131  logical(LGP) :: found
2132  integer(I4B) :: ipos
2133  ! -- code
2134  !
2135  ! -- process optional variables
2136  found = .false.
2137  if (present(name) .and. present(mem_path)) then
2138  call get_from_memorylist(name, mem_path, mt, found)
2139  nullify (mt%astr1d)
2140  else
2141  do ipos = 1, memorylist%count()
2142  mt => memorylist%Get(ipos)
2143  if (associated(mt%astr1d, astr1d)) then
2144  nullify (mt%astr1d)
2145  found = .true.
2146  exit
2147  end if
2148  end do
2149  end if
2150  if (.not. found .and. size(astr1d) > 0) then
2151  call store_error('programming error in deallocate_str1d', terminate=.true.)
2152  else
2153  if (mt%master) then
2154  deallocate (astr1d)
2155  else
2156  nullify (astr1d)
2157  end if
2158  end if
2159  !
2160  ! -- return
2161  return

◆ get_from_memorylist()

subroutine, public memorymanagermodule::get_from_memorylist ( character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
type(memorytype), intent(inout), pointer  mt,
logical(lgp), intent(out)  found,
logical(lgp), intent(in), optional  check 
)

Default value for

check is .true. which means that this
routine will kill the program when the memory entry cannot be found.
Parameters
[in]namevariable name
[in]mem_pathpath where the variable is stored
[in,out]mtmemory type entry
[out]foundset to .true. when found
[in]checkto suppress aborting the program when not found, set check = .false.

Definition at line 308 of file MemoryManager.f90.

309  character(len=*), intent(in) :: name !< variable name
310  character(len=*), intent(in) :: mem_path !< path where the variable is stored
311  type(MemoryType), pointer, intent(inout) :: mt !< memory type entry
312  logical(LGP), intent(out) :: found !< set to .true. when found
313  logical(LGP), intent(in), optional :: check !< to suppress aborting the program when not found,
314  !! set check = .false.
315  ! -- local
316  integer(I4B) :: ipos
317  logical(LGP) check_opt
318  ! -- code
319  !
320  ! -- initialize
321  mt => null()
322  found = .false.
323  !
324  ! -- iterate over the memory list
325  do ipos = 1, memorylist%count()
326  mt => memorylist%Get(ipos)
327  if (mt%name == name .and. mt%path == mem_path) then
328  found = .true.
329  exit
330  end if
331  end do
332  check_opt = .true.
333  if (present(check)) then
334  check_opt = check
335  end if
336  if (check_opt) then
337  if (.not. found) then
338  errmsg = "Programming error in memory manager. Variable '"// &
339  trim(name)//"' in '"//trim(mem_path)//"' cannot be "// &
340  "assigned because it does not exist in memory manager."
341  call store_error(errmsg, terminate=.true.)
342  end if
343  end if
344  !
345  ! -- return
346  return
Here is the call graph for this function:

◆ get_isize()

subroutine, public memorymanagermodule::get_isize ( character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
integer(i4b), intent(out)  isize 
)

Returns with isize = -1 when not found. Return 1 for scalars.

Parameters
[in]namevariable name
[in]mem_pathpath where the variable is stored
[out]isizenumber of elements (flattened)

Definition at line 275 of file MemoryManager.f90.

276  character(len=*), intent(in) :: name !< variable name
277  character(len=*), intent(in) :: mem_path !< path where the variable is stored
278  integer(I4B), intent(out) :: isize !< number of elements (flattened)
279  ! -- local
280  type(MemoryType), pointer :: mt => null()
281  logical(LGP) :: found
282  logical(LGP) :: terminate
283  ! -- code
284  !
285  ! -- initialize isize to a value to communicate failure
286  isize = -1
287  !
288  ! -- don't exit program if variable not found
289  terminate = .false.
290  !
291  ! -- get the entry from the memory manager
292  call get_from_memorylist(name, mem_path, mt, found, terminate)
293  !
294  ! -- set isize
295  if (found) then
296  isize = mt%isize
297  end if
298  !
299  ! -- return
300  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_mem_elem_size()

subroutine, public memorymanagermodule::get_mem_elem_size ( character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
integer(i4b), intent(out)  size 
)

Memory size in bytes, returns size = -1 when not found. This is

Parameters
[in]namevariable name
[in]mem_pathpath where the variable is stored
[out]sizesize of the variable in bytes

Definition at line 206 of file MemoryManager.f90.

207  character(len=*), intent(in) :: name !< variable name
208  character(len=*), intent(in) :: mem_path !< path where the variable is stored
209  integer(I4B), intent(out) :: size !< size of the variable in bytes
210  ! -- local
211  type(MemoryType), pointer :: mt => null()
212  logical(LGP) :: found
213  ! -- code
214  !
215  ! -- initialize size to a value to communicate failure
216  size = -1
217  !
218  ! -- get the entry from the memory manager
219  call get_from_memorylist(name, mem_path, mt, found)
220  !
221  ! -- set memory size
222  if (found) then
223  size = mt%element_size
224  end if
225  !
226  ! -- return
227  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_mem_rank()

subroutine, public memorymanagermodule::get_mem_rank ( character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
integer(i4b), intent(out)  rank 
)

Returns rank = -1 when not found.

Parameters
[in]namevariable name

Definition at line 167 of file MemoryManager.f90.

168  character(len=*), intent(in) :: name !< variable name
169  character(len=*), intent(in) :: mem_path !< mem_path
170  integer(I4B), intent(out) :: rank !< rank
171  ! -- local
172  type(MemoryType), pointer :: mt => null()
173  logical(LGP) :: found
174  ! -- code
175  !
176  ! -- initialize rank to a value to communicate failure
177  rank = -1
178  !
179  ! -- get the entry from the memory manager
180  call get_from_memorylist(name, mem_path, mt, found)
181  !
182  ! -- set rank
183  if (found) then
184  if (associated(mt%logicalsclr)) rank = 0
185  if (associated(mt%intsclr)) rank = 0
186  if (associated(mt%dblsclr)) rank = 0
187  if (associated(mt%aint1d)) rank = 1
188  if (associated(mt%aint2d)) rank = 2
189  if (associated(mt%aint3d)) rank = 3
190  if (associated(mt%adbl1d)) rank = 1
191  if (associated(mt%adbl2d)) rank = 2
192  if (associated(mt%adbl3d)) rank = 3
193  if (associated(mt%strsclr)) rank = 0
194  if (associated(mt%astr1d)) rank = 1
195  if (associated(mt%acharstr1d)) rank = 1
196  end if
197  !
198  ! -- return
199  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_mem_shape()

subroutine, public memorymanagermodule::get_mem_shape ( character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
integer(i4b), dimension(:), intent(out)  mem_shape 
)

Returns an integer array with the shape (Fortran ordering), and set shape(1) = -1 when not found.

Parameters
[in]namevariable name
[in]mem_pathpath where the variable is stored
[out]mem_shapeshape of the variable

Definition at line 235 of file MemoryManager.f90.

236  character(len=*), intent(in) :: name !< variable name
237  character(len=*), intent(in) :: mem_path !< path where the variable is stored
238  integer(I4B), dimension(:), intent(out) :: mem_shape !< shape of the variable
239  ! -- local
240  type(MemoryType), pointer :: mt => null()
241  logical(LGP) :: found
242  ! -- code
243  !
244  ! -- get the entry from the memory manager
245  call get_from_memorylist(name, mem_path, mt, found)
246  !
247  ! -- set shape
248  if (found) then
249  if (associated(mt%logicalsclr)) mem_shape = shape(mt%logicalsclr)
250  if (associated(mt%intsclr)) mem_shape = shape(mt%logicalsclr)
251  if (associated(mt%dblsclr)) mem_shape = shape(mt%dblsclr)
252  if (associated(mt%aint1d)) mem_shape = shape(mt%aint1d)
253  if (associated(mt%aint2d)) mem_shape = shape(mt%aint2d)
254  if (associated(mt%aint3d)) mem_shape = shape(mt%aint3d)
255  if (associated(mt%adbl1d)) mem_shape = shape(mt%adbl1d)
256  if (associated(mt%adbl2d)) mem_shape = shape(mt%adbl2d)
257  if (associated(mt%adbl3d)) mem_shape = shape(mt%adbl3d)
258  if (associated(mt%strsclr)) mem_shape = shape(mt%strsclr)
259  if (associated(mt%astr1d)) mem_shape = shape(mt%astr1d)
260  if (associated(mt%acharstr1d)) mem_shape = shape(mt%acharstr1d)
261  ! -- to communicate failure
262  else
263  mem_shape(1) = -1
264  end if
265  !
266  ! -- return
267  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_mem_type()

subroutine, public memorymanagermodule::get_mem_type ( character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=lenmemtype), intent(out)  var_type 
)

Returns any of 'LOGICAL', 'INTEGER', 'DOUBLE', 'STRING'. returns 'UNKNOWN' when the variable is not found.

Parameters
[in]namevariable name
[in]mem_pathpath where the variable is stored
[out]var_typememory type

Definition at line 144 of file MemoryManager.f90.

145  character(len=*), intent(in) :: name !< variable name
146  character(len=*), intent(in) :: mem_path !< path where the variable is stored
147  character(len=LENMEMTYPE), intent(out) :: var_type !< memory type
148  ! -- local
149  type(MemoryType), pointer :: mt
150  logical(LGP) :: found
151  ! -- code
152  mt => null()
153  var_type = 'UNKNOWN'
154  call get_from_memorylist(name, mem_path, mt, found)
155  if (found) then
156  var_type = mt%memtype
157  end if
158  !
159  ! -- return
160  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mem_cleanup_table()

subroutine memorymanagermodule::mem_cleanup_table
private

Definition at line 2811 of file MemoryManager.f90.

2812  ! -- local
2813  ! -- formats
2814  ! -- code
2815  call memtab%table_da()
2816  deallocate (memtab)
2817  nullify (memtab)
2818  !
2819  ! -- return
2820  return
Here is the caller graph for this function:

◆ mem_da()

subroutine, public memorymanagermodule::mem_da

Definition at line 2960 of file MemoryManager.f90.

2961  ! -- modules
2962  use versionmodule, only: idevelopmode
2963  use inputoutputmodule, only: upcase
2964  ! -- local
2965  class(MemoryType), pointer :: mt
2966  character(len=LINELENGTH) :: error_msg
2967  character(len=LENVARNAME) :: ucname
2968  integer(I4B) :: ipos
2969  ! -- code
2970  do ipos = 1, memorylist%count()
2971  mt => memorylist%Get(ipos)
2972  if (idevelopmode == 1) then
2973  !
2974  ! -- check if memory has been deallocated
2975  if (mt%mt_associated() .and. mt%element_size == -1) then
2976  error_msg = trim(adjustl(mt%path))//' '// &
2977  trim(adjustl(mt%name))//' has invalid element size'
2978  call store_error(trim(error_msg))
2979  end if
2980  !
2981  ! -- check if memory has been deallocated
2982  if (mt%mt_associated() .and. mt%isize > 0) then
2983  error_msg = trim(adjustl(mt%path))//' '// &
2984  trim(adjustl(mt%name))//' not deallocated'
2985  call store_error(trim(error_msg))
2986  end if
2987  !
2988  ! -- check case of varname
2989  ucname = mt%name
2990  call upcase(ucname)
2991  if (mt%name /= ucname) then
2992  error_msg = trim(adjustl(mt%path))//' '// &
2993  trim(adjustl(mt%name))//' not upper case'
2994  call store_error(trim(error_msg))
2995  end if
2996  end if
2997  !
2998  ! -- deallocate instance of memory type
2999  deallocate (mt)
3000  end do
3001  call memorylist%clear()
3002  if (count_errors() > 0) then
3003  call store_error('Could not clear memory list.', terminate=.true.)
3004  end if
3005  !
3006  ! -- return
3007  return
subroutine, public upcase(word)
Convert to upper case.
This module contains version information.
Definition: version.f90:7
integer(i4b), parameter idevelopmode
Definition: version.f90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mem_detailed_table()

subroutine memorymanagermodule::mem_detailed_table ( integer(i4b), intent(in)  iout,
integer(i4b), intent(in)  nrows 
)
private
Parameters
[in]ioutunit number for mfsim.lst
[in]nrowsnumber of table rows

Definition at line 2631 of file MemoryManager.f90.

2632  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2633  integer(I4B), intent(in) :: nrows !< number of table rows
2634  ! -- local
2635  character(len=LINELENGTH) :: title
2636  character(len=LINELENGTH) :: text
2637  integer(I4B) :: nterms
2638  ! -- formats
2639  ! -- code
2640  nterms = 5
2641  !
2642  ! -- set up table title
2643  title = 'DETAILED INFORMATION ON VARIABLES STORED IN THE MEMORY MANAGER'
2644  !
2645  ! -- set up stage tableobj
2646  call table_cr(memtab, 'MEM DET', title)
2647  call memtab%table_df(nrows, nterms, iout)
2648  !
2649  ! -- origin
2650  text = 'ORIGIN'
2651  call memtab%initialize_column(text, lenmempath, alignment=tableft)
2652  !
2653  ! -- variable
2654  text = 'VARIABLE NAME'
2655  call memtab%initialize_column(text, lenvarname, alignment=tableft)
2656  !
2657  ! -- data type
2658  text = 'DATA TYPE'
2659  call memtab%initialize_column(text, 16, alignment=tableft)
2660  !
2661  ! -- size
2662  text = 'NUMBER OF ITEMS'
2663  call memtab%initialize_column(text, 20, alignment=tabright)
2664  !
2665  ! -- is it a pointer
2666  text = 'ASSOCIATED VARIABLE'
2667  call memtab%initialize_column(text, lenmemaddress, alignment=tableft)
2668  !
2669  ! -- return
2670  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mem_print_detailed()

subroutine, public memorymanagermodule::mem_print_detailed ( integer(i4b)  iout)

Definition at line 2925 of file MemoryManager.f90.

2926  integer(I4B) :: iout
2927  ! local
2928  class(MemoryType), pointer :: mt
2929  integer(I4B) :: ipos
2930 
2931  call mem_detailed_table(iout, memorylist%count())
2932  do ipos = 1, memorylist%count()
2933  mt => memorylist%Get(ipos)
2934  call mt%table_entry(memtab)
2935  end do
2936  call mem_cleanup_table()
2937 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mem_set_print_option()

subroutine, public memorymanagermodule::mem_set_print_option ( integer(i4b), intent(in)  iout,
character(len=*), intent(in)  keyword,
character(len=*), intent(inout)  error_msg 
)
Parameters
[in]ioutunit number for mfsim.lst
[in]keywordmemory print option
[in,out]error_msgreturned error message if keyword is not valid option

Definition at line 2553 of file MemoryManager.f90.

2554  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2555  character(len=*), intent(in) :: keyword !< memory print option
2556  character(len=*), intent(inout) :: error_msg !< returned error message if keyword is not valid option
2557  ! -- local
2558  ! -- format
2559  ! -- code
2560  select case (keyword)
2561  case ('NONE')
2562  iprmem = 0
2563  write (iout, '(4x, a)') &
2564  'LIMITED MEMORY INFORMATION WILL BE WRITTEN.'
2565  case ('SUMMARY')
2566  iprmem = 1
2567  write (iout, '(4x, a)') &
2568  'A SUMMARY OF SIMULATION MEMORY INFORMATION WILL BE WRITTEN.'
2569  case ('ALL')
2570  iprmem = 2
2571  write (iout, '(4x, a)') &
2572  'ALL SIMULATION MEMORY INFORMATION WILL BE WRITTEN.'
2573  case default
2574  error_msg = "Unknown memory print option '"//trim(keyword)//"."
2575  end select
2576  return
Here is the caller graph for this function:

◆ mem_summary_line()

subroutine memorymanagermodule::mem_summary_line ( character(len=*), intent(in)  component,
real(dp), intent(in)  rchars,
real(dp), intent(in)  rlog,
real(dp), intent(in)  rint,
real(dp), intent(in)  rreal,
real(dp), intent(in)  bytes 
)
private
Parameters
[in]componentcharacter defining the program component (e.g. solution)
[in]rcharsallocated size of characters (in common units)
[in]rlogallocated size of logical (in common units)
[in]rintallocated size of integer variables (in common units)
[in]rrealallocated size of real variables (in common units)
[in]bytestotal allocated memory in memory manager (in common units)

Definition at line 2675 of file MemoryManager.f90.

2676  character(len=*), intent(in) :: component !< character defining the program component (e.g. solution)
2677  real(DP), intent(in) :: rchars !< allocated size of characters (in common units)
2678  real(DP), intent(in) :: rlog !< allocated size of logical (in common units)
2679  real(DP), intent(in) :: rint !< allocated size of integer variables (in common units)
2680  real(DP), intent(in) :: rreal !< allocated size of real variables (in common units)
2681  real(DP), intent(in) :: bytes !< total allocated memory in memory manager (in common units)
2682  ! -- formats
2683  ! -- code
2684  !
2685  ! -- write line
2686  call memtab%add_term(component)
2687  call memtab%add_term(rchars)
2688  call memtab%add_term(rlog)
2689  call memtab%add_term(rint)
2690  call memtab%add_term(rreal)
2691  call memtab%add_term(bytes)
2692  !
2693  ! -- return
2694  return
Here is the caller graph for this function:

◆ mem_summary_table()

subroutine memorymanagermodule::mem_summary_table ( integer(i4b), intent(in)  iout,
integer(i4b), intent(in)  nrows,
character(len=*), intent(in)  cunits 
)
private
Parameters
[in]ioutunit number for mfsim.lst
[in]nrowsnumber of table rows
[in]cunitsmemory units (bytes, kilobytes, megabytes, or gigabytes)

Definition at line 2581 of file MemoryManager.f90.

2582  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2583  integer(I4B), intent(in) :: nrows !< number of table rows
2584  character(len=*), intent(in) :: cunits !< memory units (bytes, kilobytes, megabytes, or gigabytes)
2585  ! -- local
2586  character(len=LINELENGTH) :: title
2587  character(len=LINELENGTH) :: text
2588  integer(I4B) :: nterms
2589  ! -- formats
2590  ! -- code
2591  nterms = 6
2592  !
2593  ! -- set up table title
2594  title = 'SUMMARY INFORMATION ON VARIABLES STORED IN THE MEMORY MANAGER, '// &
2595  'IN '//trim(cunits)
2596  !
2597  ! -- set up stage tableobj
2598  call table_cr(memtab, 'MEM SUM', title)
2599  call memtab%table_df(nrows, nterms, iout)
2600  !
2601  ! -- data type
2602  text = 'COMPONENT'
2603  call memtab%initialize_column(text, 20, alignment=tableft)
2604  !
2605  ! -- memory allocated for characters
2606  text = 'CHARACTER'
2607  call memtab%initialize_column(text, 15, alignment=tabcenter)
2608  !
2609  ! -- memory allocated for logical
2610  text = 'LOGICAL'
2611  call memtab%initialize_column(text, 15, alignment=tabcenter)
2612  !
2613  ! -- memory allocated for integers
2614  text = 'INTEGER'
2615  call memtab%initialize_column(text, 15, alignment=tabcenter)
2616  !
2617  ! -- memory allocated for reals
2618  text = 'REAL'
2619  call memtab%initialize_column(text, 15, alignment=tabcenter)
2620  !
2621  ! -- total memory allocated
2622  text = 'TOTAL'
2623  call memtab%initialize_column(text, 15, alignment=tabcenter)
2624  !
2625  ! -- return
2626  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mem_summary_total()

subroutine memorymanagermodule::mem_summary_total ( integer(i4b), intent(in)  iout,
real(dp), intent(in)  bytes 
)
private
Parameters
[in]ioutunit number for mfsim.lst
[in]bytestotal number of bytes allocated in the memory manager

Definition at line 2733 of file MemoryManager.f90.

2734  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2735  real(DP), intent(in) :: bytes !< total number of bytes allocated in the memory manager
2736  ! -- local
2737  character(len=LINELENGTH) :: title
2738  character(len=LINELENGTH) :: text
2739  character(LEN=10) :: cunits
2740  integer(I4B) :: nterms
2741  integer(I4B) :: nrows
2742  real(DP) :: fact
2743  real(DP) :: smb
2744  ! -- formats
2745  ! -- code
2746  !
2747  ! -- calculate factor and memory units
2748  call mem_units(bytes, fact, cunits)
2749  !
2750  ! -- set table terms
2751  nterms = 2
2752  nrows = 6
2753  !
2754  ! -- set up table title
2755  title = 'MEMORY MANAGER TOTAL STORAGE BY DATA TYPE, IN '//trim(cunits)
2756  !
2757  ! -- set up stage tableobj
2758  call table_cr(memtab, 'MEM TOT', title)
2759  call memtab%table_df(nrows, nterms, iout)
2760  !
2761  ! -- data type
2762  text = 'DATA TYPE'
2763  call memtab%initialize_column(text, 15, alignment=tableft)
2764  !
2765  ! -- number of values
2766  text = 'ALLOCATED MEMORY'
2767  call memtab%initialize_column(text, 15, alignment=tabcenter)
2768  !
2769  ! -- write data
2770  !
2771  ! -- characters
2772  smb = real(nvalues_astr, dp) * fact
2773  call memtab%add_term('Character')
2774  call memtab%add_term(smb)
2775  !
2776  ! -- logicals
2777  smb = real(nvalues_alogical * lgp, dp) * fact
2778  call memtab%add_term('Logical')
2779  call memtab%add_term(smb)
2780  !
2781  ! -- integers
2782  smb = real(nvalues_aint * i4b, dp) * fact
2783  call memtab%add_term('Integer')
2784  call memtab%add_term(smb)
2785  !
2786  ! -- reals
2787  smb = real(nvalues_adbl * dp, dp) * fact
2788  call memtab%add_term('Real')
2789  call memtab%add_term(smb)
2790  !
2791  ! -- total memory usage
2792  call memtab%print_separator()
2793  smb = bytes * fact
2794  call memtab%add_term('Total')
2795  call memtab%add_term(smb)
2796  !
2797  ! -- Virtual memory
2798  smb = calc_virtual_mem() * fact
2799  call memtab%add_term('Virtual')
2800  call memtab%add_term(smb)
2801  !
2802  ! -- deallocate table
2803  call mem_cleanup_table()
2804  !
2805  ! -- return
2806  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mem_unique_origins()

subroutine memorymanagermodule::mem_unique_origins ( character(len=lenmemaddress), dimension(:), intent(inout), allocatable  cunique)
Parameters
[in,out]cuniquearray with unique first components

Definition at line 3013 of file MemoryManager.f90.

3014  ! -- modules
3016  ! -- dummy
3017  character(len=LENMEMADDRESS), allocatable, dimension(:), intent(inout) :: &
3018  cunique !< array with unique first components
3019  ! -- local
3020  class(MemoryType), pointer :: mt
3021  character(len=LENMEMPATH) :: context
3022  character(len=LENCOMPONENTNAME) :: component
3023  character(len=LENCOMPONENTNAME) :: subcomponent
3024  character(len=LENMEMADDRESS) :: context_component
3025  integer(I4B) :: ipos
3026  integer(I4B) :: ipa
3027  ! -- code
3028  !
3029  ! -- initialize cunique
3030  allocate (cunique(0))
3031  !
3032  ! -- find unique origins
3033  do ipos = 1, memorylist%count()
3034  mt => memorylist%Get(ipos)
3035  call split_mem_path(mt%path, component, subcomponent)
3036  context = get_mem_path_context(mt%path)
3037  context_component = trim(context)//component
3038  ipa = ifind(cunique, context_component)
3039  if (ipa < 1) then
3040  call expandarray(cunique, 1)
3041  cunique(size(cunique)) = context_component
3042  end if
3043  end do
3044  !
3045  ! -- return
3046  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mem_units()

subroutine memorymanagermodule::mem_units ( real(dp), intent(in)  bytes,
real(dp), intent(inout)  fact,
character(len=*), intent(inout)  cunits 
)
private
Parameters
[in]bytestotal nr. of bytes
[in,out]factconversion factor
[in,out]cunitsstring with memory unit

Definition at line 2699 of file MemoryManager.f90.

2700  ! -- dummy
2701  real(DP), intent(in) :: bytes !< total nr. of bytes
2702  real(DP), intent(inout) :: fact !< conversion factor
2703  character(len=*), intent(inout) :: cunits !< string with memory unit
2704  ! -- local
2705  ! -- formats
2706  ! -- code
2707  !
2708  ! -- initialize factor and unit string
2709  cunits = 'UNKNOWN'
2710  fact = done
2711  !
2712  ! -- set factor
2713  if (bytes < dep3) then
2714  fact = done
2715  cunits = 'BYTES'
2716  else if (bytes < dep6) then
2717  fact = dem3
2718  cunits = 'KILOBYTES'
2719  else if (bytes < dep9) then
2720  fact = dem6
2721  cunits = 'MEGABYTES'
2722  else
2723  fact = dem9
2724  cunits = 'GIGABYTES'
2725  end if
2726  !
2727  ! -- return
2728  return
Here is the caller graph for this function:

◆ mem_write_usage()

subroutine, public memorymanagermodule::mem_write_usage ( integer(i4b), intent(in)  iout)

The total memory usage by data types (int, real, etc.) is written for every simulation.

Parameters
[in]ioutunit number for mfsim.lst

Definition at line 2829 of file MemoryManager.f90.

2830  integer(I4B), intent(in) :: iout !< unit number for mfsim.lst
2831  ! -- local
2832  class(MemoryType), pointer :: mt
2833  character(len=LENMEMADDRESS), allocatable, dimension(:) :: cunique
2834  ! character(len=LENMEMPATH) :: mem_path
2835  character(len=LENMEMPATH) :: context
2836  character(len=LENCOMPONENTNAME) :: component
2837  character(len=LENCOMPONENTNAME) :: subcomponent
2838  character(len=LENMEMADDRESS) :: context_component
2839  character(LEN=10) :: cunits
2840  integer(I4B) :: ipos
2841  integer(I4B) :: icomp
2842  integer(I4B) :: ilen
2843  integer(I8B) :: nchars
2844  integer(I8B) :: nlog
2845  integer(I8B) :: nint
2846  integer(I8B) :: nreal
2847  real(DP) :: simbytes
2848  real(DP) :: fact
2849  real(DP) :: rchars
2850  real(DP) :: rlog
2851  real(DP) :: rint
2852  real(DP) :: rreal
2853  real(DP) :: bytes
2854  ! -- formats
2855  ! -- code
2856  !
2857  ! -- Calculate simulation memory allocation
2858  simbytes = (nvalues_astr + &
2859  nvalues_alogical * lgp + &
2860  nvalues_aint * i4b + &
2861  nvalues_adbl * dp)
2862  simbytes = real(simbytes, dp)
2863  !
2864  ! -- calculate factor and memory units
2865  call mem_units(simbytes, fact, cunits)
2866  !
2867  ! -- Write summary table for simulation components
2868  if (iprmem == 1) then
2869  !
2870  ! -- Find unique names of simulation components
2871  call mem_unique_origins(cunique)
2872  call mem_summary_table(iout, size(cunique), cunits)
2873  do icomp = 1, size(cunique)
2874  nchars = 0
2875  nlog = 0
2876  nint = 0
2877  nreal = 0
2878  bytes = dzero
2879  ilen = len_trim(cunique(icomp))
2880  do ipos = 1, memorylist%count()
2881  mt => memorylist%Get(ipos)
2882  call split_mem_path(mt%path, component, subcomponent)
2883  context = get_mem_path_context(mt%path)
2884  context_component = trim(context)//component
2885  if (cunique(icomp) /= context_component(1:ilen)) cycle
2886  if (.not. mt%master) cycle
2887  if (mt%memtype(1:6) == 'STRING') then
2888  nchars = nchars + mt%isize * mt%element_size
2889  else if (mt%memtype(1:7) == 'LOGICAL') then
2890  nlog = nlog + mt%isize
2891  else if (mt%memtype(1:7) == 'INTEGER') then
2892  nint = nint + mt%isize
2893  else if (mt%memtype(1:6) == 'DOUBLE') then
2894  nreal = nreal + mt%isize
2895  end if
2896  end do
2897  !
2898  ! -- calculate size of each data type in bytes
2899  rchars = real(nchars, dp) * fact
2900  rlog = real(nlog * lgp, dp) * fact
2901  rint = real(nint * i4b, dp) * fact
2902  rreal = real(nreal * dp, dp) * fact
2903  !
2904  ! -- calculate total storage in bytes
2905  bytes = rchars + rlog + rint + rreal
2906  !
2907  ! -- write data
2908  call mem_summary_line(cunique(icomp), rchars, rlog, rint, rreal, bytes)
2909  end do
2910  call mem_cleanup_table()
2911  end if
2912  !
2913  ! -- Write table with all variables for iprmem == 2
2914  if (iprmem == 2) then
2915  call mem_print_detailed(iout)
2916  end if
2917  !
2918  ! -- Write total memory allocation
2919  call mem_summary_total(iout, simbytes)
2920  !
2921  ! -- return
2922  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ reallocate_charstr1d()

subroutine memorymanagermodule::reallocate_charstr1d ( type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  acharstr1d,
integer(i4b), intent(in)  ilen,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]acharstr1dthe reallocated charstring array
[in]ilenstring length
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1262 of file MemoryManager.f90.

1263  type(CharacterStringType), dimension(:), pointer, contiguous, &
1264  intent(inout) :: acharstr1d !< the reallocated charstring array
1265  integer(I4B), intent(in) :: ilen !< string length
1266  integer(I4B), intent(in) :: nrow !< number of rows
1267  character(len=*), intent(in) :: name !< variable name
1268  character(len=*), intent(in) :: mem_path !< path where variable is stored
1269  ! -- local
1270  type(MemoryType), pointer :: mt
1271  logical(LGP) :: found
1272  type(CharacterStringType), dimension(:), allocatable :: astrtemp
1273  character(len=ilen) :: string
1274  integer(I4B) :: istat
1275  integer(I4B) :: isize
1276  integer(I4B) :: isize_old
1277  integer(I4B) :: nrow_old
1278  integer(I4B) :: n
1279  !
1280  ! -- Initialize string
1281  string = ''
1282  !
1283  ! -- Find and assign mt
1284  call get_from_memorylist(name, mem_path, mt, found)
1285  !
1286  ! -- reallocate astr1d
1287  if (found) then
1288  isize_old = mt%isize
1289  if (isize_old > 0) then
1290  nrow_old = size(acharstr1d)
1291  else
1292  nrow_old = 0
1293  end if
1294  !
1295  ! -- calculate isize
1296  isize = nrow
1297  !
1298  ! -- allocate astrtemp
1299  allocate (astrtemp(nrow), stat=istat, errmsg=errmsg)
1300  if (istat /= 0) then
1301  call allocate_error(name, mem_path, istat, isize)
1302  end if
1303  !
1304  ! -- copy existing values
1305  do n = 1, nrow_old
1306  astrtemp(n) = acharstr1d(n)
1307  end do
1308  !
1309  ! -- fill new values with missing values
1310  do n = nrow_old + 1, nrow
1311  astrtemp(n) = string
1312  end do
1313  !
1314  ! -- deallocate mt pointer, repoint, recalculate isize
1315  deallocate (acharstr1d)
1316  !
1317  ! -- allocate astr1d
1318  allocate (acharstr1d(nrow), stat=istat, errmsg=errmsg)
1319  if (istat /= 0) then
1320  call allocate_error(name, mem_path, istat, isize)
1321  end if
1322  !
1323  ! -- fill the reallocated character array
1324  do n = 1, nrow
1325  acharstr1d(n) = astrtemp(n)
1326  end do
1327  !
1328  ! -- deallocate temporary storage
1329  deallocate (astrtemp)
1330  !
1331  ! -- reset memory manager values
1332  mt%acharstr1d => acharstr1d
1333  mt%element_size = ilen
1334  mt%isize = isize
1335  mt%nrealloc = mt%nrealloc + 1
1336  mt%master = .true.
1337  nvalues_astr = nvalues_astr + isize - isize_old
1338  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
1339  else
1340  errmsg = "Programming error, variable '"//trim(name)//"' from '"// &
1341  trim(mem_path)//"' is not defined in the memory manager. Use "// &
1342  "mem_allocate instead."
1343  call store_error(errmsg, terminate=.true.)
1344  end if
1345  !
1346  ! -- return
1347  return

◆ reallocate_dbl1d()

subroutine memorymanagermodule::reallocate_dbl1d ( real(dp), dimension(:), intent(inout), pointer, contiguous  adbl,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblthe reallocated 1d real array
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1447 of file MemoryManager.f90.

1448  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< the reallocated 1d real array
1449  integer(I4B), intent(in) :: nrow !< number of rows
1450  character(len=*), intent(in) :: name !< variable name
1451  character(len=*), intent(in) :: mem_path !< path where variable is stored
1452  ! -- local
1453  type(MemoryType), pointer :: mt
1454  integer(I4B) :: istat
1455  integer(I4B) :: isize
1456  integer(I4B) :: i
1457  integer(I4B) :: isizeold
1458  integer(I4B) :: ifill
1459  logical(LGP) :: found
1460  ! -- code
1461  !
1462  ! -- Find and assign mt
1463  call get_from_memorylist(name, mem_path, mt, found)
1464  !
1465  ! -- Allocate adbl and then refill
1466  isize = nrow
1467  isizeold = size(mt%adbl1d)
1468  ifill = min(isizeold, isize)
1469  allocate (adbl(nrow), stat=istat, errmsg=errmsg)
1470  if (istat /= 0) then
1471  call allocate_error(name, mem_path, istat, isize)
1472  end if
1473  do i = 1, ifill
1474  adbl(i) = mt%adbl1d(i)
1475  end do
1476  !
1477  ! -- deallocate mt pointer, repoint, recalculate isize
1478  deallocate (mt%adbl1d)
1479  mt%adbl1d => adbl
1480  mt%element_size = dp
1481  mt%isize = isize
1482  mt%nrealloc = mt%nrealloc + 1
1483  mt%master = .true.
1484  nvalues_adbl = nvalues_adbl + isize - isizeold
1485  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', isize
1486  !
1487  ! -- return
1488  return

◆ reallocate_dbl2d()

subroutine memorymanagermodule::reallocate_dbl2d ( real(dp), dimension(:, :), intent(inout), pointer, contiguous  adbl,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblthe reallocated 2d real array
[in]ncolnumber of columns
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1493 of file MemoryManager.f90.

1494  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< the reallocated 2d real array
1495  integer(I4B), intent(in) :: ncol !< number of columns
1496  integer(I4B), intent(in) :: nrow !< number of rows
1497  character(len=*), intent(in) :: name !< variable name
1498  character(len=*), intent(in) :: mem_path !< path where variable is stored
1499  ! -- local
1500  type(MemoryType), pointer :: mt
1501  logical(LGP) :: found
1502  integer(I4B) :: istat
1503  integer(I4B), dimension(2) :: ishape
1504  integer(I4B) :: i
1505  integer(I4B) :: j
1506  integer(I4B) :: isize
1507  integer(I4B) :: isizeold
1508  ! -- code
1509  !
1510  ! -- Find and assign mt
1511  call get_from_memorylist(name, mem_path, mt, found)
1512  !
1513  ! -- Allocate adbl and then refill
1514  ishape = shape(mt%adbl2d)
1515  isize = nrow * ncol
1516  isizeold = ishape(1) * ishape(2)
1517  allocate (adbl(ncol, nrow), stat=istat, errmsg=errmsg)
1518  if (istat /= 0) then
1519  call allocate_error(name, mem_path, istat, isize)
1520  end if
1521  do i = 1, ishape(2)
1522  do j = 1, ishape(1)
1523  adbl(j, i) = mt%adbl2d(j, i)
1524  end do
1525  end do
1526  !
1527  ! -- deallocate mt pointer, repoint, recalculate isize
1528  deallocate (mt%adbl2d)
1529  mt%adbl2d => adbl
1530  mt%element_size = dp
1531  mt%isize = isize
1532  mt%nrealloc = mt%nrealloc + 1
1533  mt%master = .true.
1534  nvalues_adbl = nvalues_adbl + isize - isizeold
1535  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
1536  !
1537  ! -- return
1538  return

◆ reallocate_int1d()

subroutine memorymanagermodule::reallocate_int1d ( integer(i4b), dimension(:), intent(inout), pointer, contiguous  aint,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintthe reallocated integer array
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1352 of file MemoryManager.f90.

1353  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< the reallocated integer array
1354  integer(I4B), intent(in) :: nrow !< number of rows
1355  character(len=*), intent(in) :: name !< variable name
1356  character(len=*), intent(in) :: mem_path !< path where variable is stored
1357  ! -- local
1358  type(MemoryType), pointer :: mt
1359  logical(LGP) :: found
1360  integer(I4B) :: istat
1361  integer(I4B) :: isize
1362  integer(I4B) :: i
1363  integer(I4B) :: isizeold
1364  integer(I4B) :: ifill
1365  ! -- code
1366  !
1367  ! -- Find and assign mt
1368  call get_from_memorylist(name, mem_path, mt, found)
1369  !
1370  ! -- Allocate aint and then refill
1371  isize = nrow
1372  isizeold = size(mt%aint1d)
1373  ifill = min(isizeold, isize)
1374  allocate (aint(nrow), stat=istat, errmsg=errmsg)
1375  if (istat /= 0) then
1376  call allocate_error(name, mem_path, istat, isize)
1377  end if
1378  do i = 1, ifill
1379  aint(i) = mt%aint1d(i)
1380  end do
1381  !
1382  ! -- deallocate mt pointer, repoint, recalculate isize
1383  deallocate (mt%aint1d)
1384  mt%aint1d => aint
1385  mt%element_size = i4b
1386  mt%isize = isize
1387  mt%nrealloc = mt%nrealloc + 1
1388  mt%master = .true.
1389  nvalues_aint = nvalues_aint + isize - isizeold
1390  !
1391  ! -- return
1392  return

◆ reallocate_int2d()

subroutine memorymanagermodule::reallocate_int2d ( integer(i4b), dimension(:, :), intent(inout), pointer, contiguous  aint,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintthe reallocated 2d integer array
[in]ncolnumber of columns
[in]nrownumber of rows
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1397 of file MemoryManager.f90.

1398  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< the reallocated 2d integer array
1399  integer(I4B), intent(in) :: ncol !< number of columns
1400  integer(I4B), intent(in) :: nrow !< number of rows
1401  character(len=*), intent(in) :: name !< variable name
1402  character(len=*), intent(in) :: mem_path !< path where variable is stored
1403  ! -- local
1404  type(MemoryType), pointer :: mt
1405  logical(LGP) :: found
1406  integer(I4B) :: istat
1407  integer(I4B), dimension(2) :: ishape
1408  integer(I4B) :: i
1409  integer(I4B) :: j
1410  integer(I4B) :: isize
1411  integer(I4B) :: isizeold
1412  ! -- code
1413  !
1414  ! -- Find and assign mt
1415  call get_from_memorylist(name, mem_path, mt, found)
1416  !
1417  ! -- Allocate aint and then refill
1418  ishape = shape(mt%aint2d)
1419  isize = nrow * ncol
1420  isizeold = ishape(1) * ishape(2)
1421  allocate (aint(ncol, nrow), stat=istat, errmsg=errmsg)
1422  if (istat /= 0) then
1423  call allocate_error(name, mem_path, istat, isize)
1424  end if
1425  do i = 1, ishape(2)
1426  do j = 1, ishape(1)
1427  aint(j, i) = mt%aint2d(j, i)
1428  end do
1429  end do
1430  !
1431  ! -- deallocate mt pointer, repoint, recalculate isize
1432  deallocate (mt%aint2d)
1433  mt%aint2d => aint
1434  mt%element_size = i4b
1435  mt%isize = isize
1436  mt%nrealloc = mt%nrealloc + 1
1437  mt%master = .true.
1438  nvalues_aint = nvalues_aint + isize - isizeold
1439  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
1440  !
1441  ! -- return
1442  return

◆ reallocate_str1d()

subroutine memorymanagermodule::reallocate_str1d ( character(len=ilen), dimension(:), intent(inout), pointer, contiguous  astr,
integer(i4b), intent(in)  ilen,
integer(i4b), intent(in)  nrow,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in]ilenstring length
[in]nrownumber of rows
[in,out]astrthe reallocated string array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1178 of file MemoryManager.f90.

1179  integer(I4B), intent(in) :: ilen !< string length
1180  integer(I4B), intent(in) :: nrow !< number of rows
1181  character(len=ilen), dimension(:), pointer, contiguous, intent(inout) :: astr !< the reallocated string array
1182  character(len=*), intent(in) :: name !< variable name
1183  character(len=*), intent(in) :: mem_path !< path where variable is stored
1184  ! -- local
1185  type(MemoryType), pointer :: mt
1186  logical(LGP) :: found
1187  character(len=ilen), dimension(:), allocatable :: astrtemp
1188  integer(I4B) :: istat
1189  integer(I4B) :: isize
1190  integer(I4B) :: isize_old
1191  integer(I4B) :: nrow_old
1192  integer(I4B) :: n
1193  !
1194  ! -- Find and assign mt
1195  call get_from_memorylist(name, mem_path, mt, found)
1196  !
1197  ! -- reallocate astr1d
1198  if (found) then
1199  isize_old = mt%isize
1200  if (isize_old > 0) then
1201  nrow_old = size(astr)
1202  else
1203  nrow_old = 0
1204  end if
1205  !
1206  ! -- calculate isize
1207  isize = nrow
1208  !
1209  ! -- allocate astrtemp
1210  allocate (astrtemp(nrow), stat=istat, errmsg=errmsg)
1211  if (istat /= 0) then
1212  call allocate_error(name, mem_path, istat, isize)
1213  end if
1214  !
1215  ! -- copy existing values
1216  do n = 1, nrow_old
1217  astrtemp(n) = astr(n)
1218  end do
1219  !
1220  ! -- fill new values with missing values
1221  do n = nrow_old + 1, nrow
1222  astrtemp(n) = ''
1223  end do
1224  !
1225  ! -- deallocate mt pointer, repoint, recalculate isize
1226  deallocate (astr)
1227  !
1228  ! -- allocate astr1d
1229  allocate (astr(nrow), stat=istat, errmsg=errmsg)
1230  if (istat /= 0) then
1231  call allocate_error(name, mem_path, istat, isize)
1232  end if
1233  !
1234  ! -- fill the reallocate character array
1235  do n = 1, nrow
1236  astr(n) = astrtemp(n)
1237  end do
1238  !
1239  ! -- deallocate temporary storage
1240  deallocate (astrtemp)
1241  !
1242  ! -- reset memory manager values
1243  mt%element_size = ilen
1244  mt%isize = isize
1245  mt%nrealloc = mt%nrealloc + 1
1246  mt%master = .true.
1247  nvalues_astr = nvalues_astr + isize - isize_old
1248  write (mt%memtype, "(a,' LEN=',i0,' (',i0,')')") 'STRING', ilen, nrow
1249  else
1250  errmsg = "Programming error, variable '"//trim(name)//"' from '"// &
1251  trim(mem_path)//"' is not defined in the memory manager. Use "// &
1252  "mem_allocate instead."
1253  call store_error(errmsg, terminate=.true.)
1254  end if
1255  !
1256  ! -- return
1257  return

◆ reassignptr_dbl1d()

subroutine memorymanagermodule::reassignptr_dbl1d ( real(dp), dimension(:), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name_target,
character(len=*), intent(in)  mem_path_target 
)
private
Parameters
[in,out]adblpointer to 1d real array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]name_targetname of target variable
[in]mem_path_targetpath where target variable is stored

Definition at line 2013 of file MemoryManager.f90.

2014  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< pointer to 1d real array
2015  character(len=*), intent(in) :: name !< variable name
2016  character(len=*), intent(in) :: mem_path !< path where variable is stored
2017  character(len=*), intent(in) :: name_target !< name of target variable
2018  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
2019  ! -- local
2020  type(MemoryType), pointer :: mt
2021  type(MemoryType), pointer :: mt2
2022  logical(LGP) :: found
2023  ! -- code
2024  call get_from_memorylist(name, mem_path, mt, found)
2025  call get_from_memorylist(name_target, mem_path_target, mt2, found)
2026  if (size(adbl) > 0) then
2027  nvalues_adbl = nvalues_adbl - size(adbl)
2028  deallocate (adbl)
2029  end if
2030  adbl => mt2%adbl1d
2031  mt%adbl1d => adbl
2032  mt%element_size = dp
2033  mt%isize = size(adbl)
2034  write (mt%memtype, "(a,' (',i0,')')") 'DOUBLE', mt%isize
2035  !
2036  ! -- set master information
2037  mt%master = .false.
2038  mt%mastername = name_target
2039  mt%masterPath = mem_path_target
2040  !
2041  ! -- return
2042  return

◆ reassignptr_dbl2d()

subroutine memorymanagermodule::reassignptr_dbl2d ( real(dp), dimension(:, :), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name_target,
character(len=*), intent(in)  mem_path_target 
)
private
Parameters
[in,out]adblpointer to 2d real array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]name_targetname of target variable
[in]mem_path_targetpath where target variable is stored

Definition at line 2047 of file MemoryManager.f90.

2048  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 2d real array
2049  character(len=*), intent(in) :: name !< variable name
2050  character(len=*), intent(in) :: mem_path !< path where variable is stored
2051  character(len=*), intent(in) :: name_target !< name of target variable
2052  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
2053  ! -- local
2054  type(MemoryType), pointer :: mt
2055  type(MemoryType), pointer :: mt2
2056  logical(LGP) :: found
2057  integer(I4B) :: ncol
2058  integer(I4b) :: nrow
2059  ! -- code
2060  call get_from_memorylist(name, mem_path, mt, found)
2061  call get_from_memorylist(name_target, mem_path_target, mt2, found)
2062  if (size(adbl) > 0) then
2063  nvalues_adbl = nvalues_adbl - size(adbl)
2064  deallocate (adbl)
2065  end if
2066  adbl => mt2%adbl2d
2067  mt%adbl2d => adbl
2068  mt%element_size = dp
2069  mt%isize = size(adbl)
2070  ncol = size(adbl, dim=1)
2071  nrow = size(adbl, dim=2)
2072  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'DOUBLE', ncol, nrow
2073  !
2074  ! -- set master information
2075  mt%master = .false.
2076  mt%mastername = name_target
2077  mt%masterPath = mem_path_target
2078  !
2079  ! -- return
2080  return

◆ reassignptr_int()

subroutine memorymanagermodule::reassignptr_int ( integer(i4b), intent(inout), pointer  sclr,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name_target,
character(len=*), intent(in)  mem_path_target 
)
private
Parameters
[in,out]sclrpointer to integer scalar
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]name_targetname of target variable
[in]mem_path_targetpath where target variable is stored

Definition at line 1907 of file MemoryManager.f90.

1908  integer(I4B), pointer, intent(inout) :: sclr !< pointer to integer scalar
1909  character(len=*), intent(in) :: name !< variable name
1910  character(len=*), intent(in) :: mem_path !< path where variable is stored
1911  character(len=*), intent(in) :: name_target !< name of target variable
1912  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1913  ! -- local
1914  type(MemoryType), pointer :: mt
1915  type(MemoryType), pointer :: mt2
1916  logical(LGP) :: found
1917  ! -- code
1918  call get_from_memorylist(name, mem_path, mt, found)
1919  call get_from_memorylist(name_target, mem_path_target, mt2, found)
1920  if (associated(sclr)) then
1921  nvalues_aint = nvalues_aint - 1
1922  deallocate (sclr)
1923  end if
1924  sclr => mt2%intsclr
1925  mt%intsclr => sclr
1926  mt%element_size = i4b
1927  mt%isize = 1
1928  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', mt%isize
1929  !
1930  ! -- set master information
1931  mt%master = .false.
1932  mt%mastername = name_target
1933  mt%masterPath = mem_path_target
1934  !
1935  ! -- return
1936  return

◆ reassignptr_int1d()

subroutine memorymanagermodule::reassignptr_int1d ( integer(i4b), dimension(:), intent(inout), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name_target,
character(len=*), intent(in)  mem_path_target 
)
private
Parameters
[in,out]aintpointer to 1d integer array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]name_targetname of target variable
[in]mem_path_targetpath where target variable is stored

Definition at line 1941 of file MemoryManager.f90.

1942  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< pointer to 1d integer array
1943  character(len=*), intent(in) :: name !< variable name
1944  character(len=*), intent(in) :: mem_path !< path where variable is stored
1945  character(len=*), intent(in) :: name_target !< name of target variable
1946  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1947  ! -- local
1948  type(MemoryType), pointer :: mt
1949  type(MemoryType), pointer :: mt2
1950  logical(LGP) :: found
1951  ! -- code
1952  call get_from_memorylist(name, mem_path, mt, found)
1953  call get_from_memorylist(name_target, mem_path_target, mt2, found)
1954  if (size(aint) > 0) then
1955  nvalues_aint = nvalues_aint - size(aint)
1956  deallocate (aint)
1957  end if
1958  aint => mt2%aint1d
1959  mt%aint1d => aint
1960  mt%element_size = i4b
1961  mt%isize = size(aint)
1962  write (mt%memtype, "(a,' (',i0,')')") 'INTEGER', mt%isize
1963  !
1964  ! -- set master information
1965  mt%master = .false.
1966  mt%mastername = name_target
1967  mt%masterPath = mem_path_target
1968  !
1969  ! -- return
1970  return

◆ reassignptr_int2d()

subroutine memorymanagermodule::reassignptr_int2d ( integer(i4b), dimension(:, :), intent(inout), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path,
character(len=*), intent(in)  name_target,
character(len=*), intent(in)  mem_path_target 
)
private
Parameters
[in,out]aintpointer to 2d integer array
[in]namevariable name
[in]mem_pathpath where variable is stored
[in]name_targetname of target variable
[in]mem_path_targetpath where target variable is stored

Definition at line 1975 of file MemoryManager.f90.

1976  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< pointer to 2d integer array
1977  character(len=*), intent(in) :: name !< variable name
1978  character(len=*), intent(in) :: mem_path !< path where variable is stored
1979  character(len=*), intent(in) :: name_target !< name of target variable
1980  character(len=*), intent(in) :: mem_path_target !< path where target variable is stored
1981  ! -- local
1982  type(MemoryType), pointer :: mt
1983  type(MemoryType), pointer :: mt2
1984  logical(LGP) :: found
1985  integer(I4B) :: ncol
1986  integer(I4B) :: nrow
1987  ! -- code
1988  call get_from_memorylist(name, mem_path, mt, found)
1989  call get_from_memorylist(name_target, mem_path_target, mt2, found)
1990  if (size(aint) > 0) then
1991  nvalues_aint = nvalues_aint - size(aint)
1992  deallocate (aint)
1993  end if
1994  aint => mt2%aint2d
1995  mt%aint2d => aint
1996  mt%element_size = i4b
1997  mt%isize = size(aint)
1998  ncol = size(aint, dim=1)
1999  nrow = size(aint, dim=2)
2000  write (mt%memtype, "(a,' (',i0,',',i0,')')") 'INTEGER', ncol, nrow
2001  !
2002  ! -- set master information
2003  mt%master = .false.
2004  mt%mastername = name_target
2005  mt%masterPath = mem_path_target
2006  !
2007  ! -- return
2008  return

◆ setptr_charstr1d()

subroutine memorymanagermodule::setptr_charstr1d ( type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  acharstr1d,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]acharstr1dthe reallocated charstring array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1731 of file MemoryManager.f90.

1732  type(CharacterStringType), dimension(:), pointer, contiguous, &
1733  intent(inout) :: acharstr1d !< the reallocated charstring array
1734  character(len=*), intent(in) :: name !< variable name
1735  character(len=*), intent(in) :: mem_path !< path where variable is stored
1736  ! -- local
1737  type(MemoryType), pointer :: mt
1738  logical(LGP) :: found
1739  ! -- code
1740  call get_from_memorylist(name, mem_path, mt, found)
1741  acharstr1d => mt%acharstr1d
1742  !
1743  ! -- return
1744  return

◆ setptr_dbl()

subroutine memorymanagermodule::setptr_dbl ( real(dp), intent(inout), pointer  sclr,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]sclrpointer to a real scalar
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1628 of file MemoryManager.f90.

1629  real(DP), pointer, intent(inout) :: sclr !< pointer to a real scalar
1630  character(len=*), intent(in) :: name !< variable name
1631  character(len=*), intent(in) :: mem_path !< path where variable is stored
1632  ! -- local
1633  type(MemoryType), pointer :: mt
1634  logical(LGP) :: found
1635  ! -- code
1636  call get_from_memorylist(name, mem_path, mt, found)
1637  sclr => mt%dblsclr
1638  !
1639  ! -- return
1640  return

◆ setptr_dbl1d()

subroutine memorymanagermodule::setptr_dbl1d ( real(dp), dimension(:), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblpointer to 1d real array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1645 of file MemoryManager.f90.

1646  real(DP), dimension(:), pointer, contiguous, intent(inout) :: adbl !< pointer to 1d real array
1647  character(len=*), intent(in) :: name !< variable name
1648  character(len=*), intent(in) :: mem_path !< path where variable is stored
1649  ! -- local
1650  type(MemoryType), pointer :: mt
1651  logical(LGP) :: found
1652  ! -- code
1653  call get_from_memorylist(name, mem_path, mt, found)
1654  adbl => mt%adbl1d
1655  !
1656  ! -- return
1657  return

◆ setptr_dbl2d()

subroutine memorymanagermodule::setptr_dbl2d ( real(dp), dimension(:, :), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblpointer to 2d real array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1662 of file MemoryManager.f90.

1663  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 2d real array
1664  character(len=*), intent(in) :: name !< variable name
1665  character(len=*), intent(in) :: mem_path !< path where variable is stored
1666  ! -- local
1667  type(MemoryType), pointer :: mt
1668  logical(LGP) :: found
1669  ! -- code
1670  call get_from_memorylist(name, mem_path, mt, found)
1671  adbl => mt%adbl2d
1672  !
1673  ! -- return
1674  return

◆ setptr_dbl3d()

subroutine memorymanagermodule::setptr_dbl3d ( real(dp), dimension(:, :, :), intent(inout), pointer, contiguous  adbl,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]adblpointer to 3d real array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1679 of file MemoryManager.f90.

1680  real(DP), dimension(:, :, :), pointer, contiguous, intent(inout) :: adbl !< pointer to 3d real array
1681  character(len=*), intent(in) :: name !< variable name
1682  character(len=*), intent(in) :: mem_path !< path where variable is stored
1683  ! -- local
1684  type(MemoryType), pointer :: mt
1685  logical(LGP) :: found
1686  ! -- code
1687  call get_from_memorylist(name, mem_path, mt, found)
1688  adbl => mt%adbl3d
1689  !
1690  ! -- return
1691  return

◆ setptr_int()

subroutine memorymanagermodule::setptr_int ( integer(i4b), intent(inout), pointer  sclr,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]sclrpointer to integer scalar
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1560 of file MemoryManager.f90.

1561  integer(I4B), pointer, intent(inout) :: sclr !< pointer to integer scalar
1562  character(len=*), intent(in) :: name !< variable name
1563  character(len=*), intent(in) :: mem_path !< path where variable is stored
1564  ! -- local
1565  type(MemoryType), pointer :: mt
1566  logical(LGP) :: found
1567  ! -- code
1568  call get_from_memorylist(name, mem_path, mt, found)
1569  sclr => mt%intsclr
1570  !
1571  ! -- return
1572  return

◆ setptr_int1d()

subroutine memorymanagermodule::setptr_int1d ( integer(i4b), dimension(:), intent(inout), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintpointer to 1d integer array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1577 of file MemoryManager.f90.

1578  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: aint !< pointer to 1d integer array
1579  character(len=*), intent(in) :: name !< variable name
1580  character(len=*), intent(in) :: mem_path !< path where variable is stored
1581  ! -- local
1582  type(MemoryType), pointer :: mt
1583  logical(LGP) :: found
1584  ! -- code
1585  call get_from_memorylist(name, mem_path, mt, found)
1586  aint => mt%aint1d
1587  !
1588  ! -- return
1589  return

◆ setptr_int2d()

subroutine memorymanagermodule::setptr_int2d ( integer(i4b), dimension(:, :), intent(inout), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintpointer to 2d integer array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1594 of file MemoryManager.f90.

1595  integer(I4B), dimension(:, :), pointer, contiguous, intent(inout) :: aint !< pointer to 2d integer array
1596  character(len=*), intent(in) :: name !< variable name
1597  character(len=*), intent(in) :: mem_path !< path where variable is stored
1598  ! -- local
1599  type(MemoryType), pointer :: mt
1600  logical(LGP) :: found
1601  ! -- code
1602  call get_from_memorylist(name, mem_path, mt, found)
1603  aint => mt%aint2d
1604  !
1605  ! -- return
1606  return

◆ setptr_int3d()

subroutine memorymanagermodule::setptr_int3d ( integer(i4b), dimension(:, :, :), intent(inout), pointer, contiguous  aint,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]aintpointer to 3d integer array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1611 of file MemoryManager.f90.

1612  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(inout) :: aint !< pointer to 3d integer array
1613  character(len=*), intent(in) :: name !< variable name
1614  character(len=*), intent(in) :: mem_path !< path where variable is stored
1615  ! -- local
1616  type(MemoryType), pointer :: mt
1617  logical(LGP) :: found
1618  ! -- code
1619  call get_from_memorylist(name, mem_path, mt, found)
1620  aint => mt%aint3d
1621  !
1622  ! -- return
1623  return

◆ setptr_logical()

subroutine memorymanagermodule::setptr_logical ( logical(lgp), intent(inout), pointer  sclr,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]sclrpointer to logical scalar
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1543 of file MemoryManager.f90.

1544  logical(LGP), pointer, intent(inout) :: sclr !< pointer to logical scalar
1545  character(len=*), intent(in) :: name !< variable name
1546  character(len=*), intent(in) :: mem_path !< path where variable is stored
1547  ! -- local
1548  type(MemoryType), pointer :: mt
1549  logical(LGP) :: found
1550  ! -- code
1551  call get_from_memorylist(name, mem_path, mt, found)
1552  sclr => mt%logicalsclr
1553  !
1554  ! -- return
1555  return

◆ setptr_str()

subroutine memorymanagermodule::setptr_str ( character(len=:), pointer  asrt,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
asrtpointer to the character string
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1696 of file MemoryManager.f90.

1697  character(len=:), pointer :: asrt !< pointer to the character string
1698  character(len=*), intent(in) :: name !< variable name
1699  character(len=*), intent(in) :: mem_path !< path where variable is stored
1700  ! -- local
1701  type(MemoryType), pointer :: mt
1702  logical(LGP) :: found
1703  ! -- code
1704  call get_from_memorylist(name, mem_path, mt, found)
1705  asrt => mt%strsclr
1706  !
1707  ! -- return
1708  return

◆ setptr_str1d()

subroutine memorymanagermodule::setptr_str1d ( character(len=:), dimension(:), intent(inout), pointer, contiguous  astr1d,
character(len=*), intent(in)  name,
character(len=*), intent(in)  mem_path 
)
private
Parameters
[in,out]astr1dpointer to the string array
[in]namevariable name
[in]mem_pathpath where variable is stored

Definition at line 1713 of file MemoryManager.f90.

1714  character(len=:), dimension(:), &
1715  pointer, contiguous, intent(inout) :: astr1d !< pointer to the string array
1716  character(len=*), intent(in) :: name !< variable name
1717  character(len=*), intent(in) :: mem_path !< path where variable is stored
1718  ! -- local
1719  type(MemoryType), pointer :: mt
1720  logical(LGP) :: found
1721  ! -- code
1722  call get_from_memorylist(name, mem_path, mt, found)
1723  astr1d => mt%astr1d
1724  !
1725  ! -- return
1726  return

Variable Documentation

◆ iprmem

integer(i4b) memorymanagermodule::iprmem = 0
private

Definition at line 51 of file MemoryManager.f90.

51  integer(I4B) :: iprmem = 0

◆ memorylist

type(memorylisttype), public memorymanagermodule::memorylist

Definition at line 45 of file MemoryManager.f90.

45  type(MemoryListType) :: memorylist

◆ memtab

type(tabletype), pointer memorymanagermodule::memtab => null()
private

Definition at line 46 of file MemoryManager.f90.

46  type(TableType), pointer :: memtab => null()

◆ nvalues_adbl

integer(i8b) memorymanagermodule::nvalues_adbl = 0
private

Definition at line 50 of file MemoryManager.f90.

50  integer(I8B) :: nvalues_adbl = 0

◆ nvalues_aint

integer(i8b) memorymanagermodule::nvalues_aint = 0
private

Definition at line 49 of file MemoryManager.f90.

49  integer(I8B) :: nvalues_aint = 0

◆ nvalues_alogical

integer(i8b) memorymanagermodule::nvalues_alogical = 0
private

Definition at line 47 of file MemoryManager.f90.

47  integer(I8B) :: nvalues_alogical = 0

◆ nvalues_astr

integer(i8b) memorymanagermodule::nvalues_astr = 0
private

Definition at line 48 of file MemoryManager.f90.

48  integer(I8B) :: nvalues_astr = 0