MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
memorymanagermodule::mem_reallocate Interface Reference
Collaboration diagram for memorymanagermodule::mem_reallocate:
Collaboration graph

Private Member Functions

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 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...
 

Detailed Description

Definition at line 78 of file MemoryManager.f90.

Member Function/Subroutine Documentation

◆ reallocate_charstr1d()

subroutine memorymanagermodule::mem_reallocate::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
Here is the call graph for this function:

◆ reallocate_dbl1d()

subroutine memorymanagermodule::mem_reallocate::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
Here is the call graph for this function:

◆ reallocate_dbl2d()

subroutine memorymanagermodule::mem_reallocate::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
Here is the call graph for this function:

◆ reallocate_int1d()

subroutine memorymanagermodule::mem_reallocate::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
Here is the call graph for this function:

◆ reallocate_int2d()

subroutine memorymanagermodule::mem_reallocate::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
Here is the call graph for this function:

◆ reallocate_str1d()

subroutine memorymanagermodule::mem_reallocate::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
Here is the call graph for this function:

The documentation for this interface was generated from the following file: