MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
gwfmodule Module Reference

Data Types

type  gwfmodeltype
 

Functions/Subroutines

subroutine, public gwf_cr (filename, id, modelname)
 Create a new groundwater flow model object. More...
 
subroutine gwf_df (this)
 Define packages of the model. More...
 
subroutine gwf_ac (this, sparse)
 Add the internal connections of this model to the sparse matrix. More...
 
subroutine gwf_mc (this, matrix_sln)
 Map the positions of this models connections in the numerical solution coefficient matrix. More...
 
subroutine gwf_ar (this)
 GroundWater Flow Model Allocate and Read. More...
 
subroutine gwf_rp (this)
 GroundWater Flow Model Read and Prepare. More...
 
subroutine gwf_ad (this)
 GroundWater Flow Model Time Step Advance. More...
 
subroutine gwf_cf (this, kiter)
 GroundWater Flow Model calculate coefficients. More...
 
subroutine gwf_fc (this, kiter, matrix_sln, inwtflag)
 GroundWater Flow Model fill coefficients. More...
 
subroutine gwf_cc (this, innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
 GroundWater Flow Model Final Convergence Check for Boundary Packages. More...
 
subroutine gwf_ptcchk (this, iptc)
 check if pseudo-transient continuation factor should be used More...
 
subroutine gwf_ptc (this, vec_residual, iptc, ptcf)
 calculate maximum pseudo-transient continuation factor More...
 
subroutine gwf_nur (this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
 under-relaxation More...
 
subroutine gwf_cq (this, icnvg, isuppress_output)
 Groundwater flow model calculate flow. More...
 
subroutine gwf_bd (this, icnvg, isuppress_output)
 GroundWater Flow Model Budget. More...
 
subroutine gwf_ot (this)
 GroundWater Flow Model Output. More...
 
subroutine gwf_ot_obs (this)
 GroundWater Flow Model output observations. More...
 
subroutine gwf_ot_flow (this, icbcfl, ibudfl, icbcun)
 Groundwater Flow Model output flows. More...
 
subroutine gwf_ot_dv (this, idvsave, idvprint, ipflag)
 Groundwater Flow Model output dependent variable. More...
 
subroutine gwf_ot_bdsummary (this, ibudfl, ipflag)
 Groundwater Flow Model output budget summary. More...
 
subroutine gwf_fp (this)
 Final processing. More...
 
subroutine gwf_da (this)
 Deallocate. More...
 
subroutine gwf_bdentry (this, budterm, budtxt, rowlabel)
 GroundWater Flow Model Budget Entry. More...
 
integer(i4b) function gwf_get_iasym (this)
 return 1 if any package causes the matrix to be asymmetric. Otherwise return 0. More...
 
subroutine allocate_scalars (this, modelname)
 Allocate memory for non-allocatable members. More...
 
subroutine package_create (this, filtyp, ipakid, ipaknum, pakname, mempath, inunit, iout)
 Create boundary condition packages for this model. More...
 
subroutine ftype_check (this, indis)
 Check to make sure required input files have been specified. More...
 
class(gwfmodeltype) function, pointer, public castasgwfmodel (model)
 Cast to GWF model. More...
 
subroutine create_bndpkgs (this, bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
 Source package info and begin to process. More...
 
subroutine create_packages (this)
 Source package info and begin to process. More...
 
subroutine log_namfile_options (this, found)
 Write model namfile options to list file. More...
 
subroutine steady_period_check (this)
 Check for steady state period. More...
 

Variables

integer(i4b), parameter, public gwf_nbasepkg = 50
 GWF base package array descriptors. More...
 
character(len=lenpackagetype), dimension(gwf_nbasepkg), public gwf_basepkg
 
integer(i4b), parameter, public gwf_nmultipkg = 50
 GWF multi package array descriptors. More...
 
character(len=lenpackagetype), dimension(gwf_nmultipkg), public gwf_multipkg
 
integer(i4b), parameter niunit_gwf = GWF_NBASEPKG + GWF_NMULTIPKG
 

Function/Subroutine Documentation

◆ allocate_scalars()

subroutine gwfmodule::allocate_scalars ( class(gwfmodeltype this,
character(len=*), intent(in)  modelname 
)
private

Definition at line 1167 of file gwf.f90.

1168  ! -- modules
1170  ! -- dummy
1171  class(GwfModelType) :: this
1172  character(len=*), intent(in) :: modelname
1173  !
1174  ! -- allocate members from parent class
1175  call this%NumericalModelType%allocate_scalars(modelname)
1176  !
1177  ! -- allocate members that are part of model class
1178  call mem_allocate(this%inic, 'INIC', this%memoryPath)
1179  call mem_allocate(this%inoc, 'INOC', this%memoryPath)
1180  call mem_allocate(this%innpf, 'INNPF', this%memoryPath)
1181  call mem_allocate(this%inbuy, 'INBUY', this%memoryPath)
1182  call mem_allocate(this%invsc, 'INVSC', this%memoryPath)
1183  call mem_allocate(this%insto, 'INSTO', this%memoryPath)
1184  call mem_allocate(this%incsub, 'INCSUB', this%memoryPath)
1185  call mem_allocate(this%inmvr, 'INMVR', this%memoryPath)
1186  call mem_allocate(this%inhfb, 'INHFB', this%memoryPath)
1187  call mem_allocate(this%ingnc, 'INGNC', this%memoryPath)
1188  call mem_allocate(this%inobs, 'INOBS', this%memoryPath)
1189  call mem_allocate(this%iss, 'ISS', this%memoryPath)
1190  call mem_allocate(this%inewtonur, 'INEWTONUR', this%memoryPath)
1191  !
1192  this%inic = 0
1193  this%inoc = 0
1194  this%innpf = 0
1195  this%inbuy = 0
1196  this%invsc = 0
1197  this%insto = 0
1198  this%incsub = 0
1199  this%inmvr = 0
1200  this%inhfb = 0
1201  this%ingnc = 0
1202  this%inobs = 0
1203  this%iss = 1 !default is steady-state (i.e., no STO package)
1204  this%inewtonur = 0 !default is to not use newton bottom head dampening

◆ castasgwfmodel()

class(gwfmodeltype) function, pointer, public gwfmodule::castasgwfmodel ( class(*), intent(inout), pointer  model)

Definition at line 1333 of file gwf.f90.

1334  implicit none
1335  class(*), pointer, intent(inout) :: model
1336  class(GwfModelType), pointer :: gwfModel
1337 
1338  gwfmodel => null()
1339  if (.not. associated(model)) return
1340  select type (model)
1341  class is (gwfmodeltype)
1342  gwfmodel => model
1343  end select
Here is the caller graph for this function:

◆ create_bndpkgs()

subroutine gwfmodule::create_bndpkgs ( class(gwfmodeltype this,
integer(i4b), dimension(:), intent(inout), allocatable  bndpkgs,
type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  pkgtypes,
type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  pkgnames,
type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  mempaths,
integer(i4b), dimension(:), intent(inout), pointer, contiguous  inunits 
)
private

Definition at line 1348 of file gwf.f90.

1350  ! -- modules
1353  ! -- dummy
1354  class(GwfModelType) :: this
1355  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
1356  type(CharacterStringType), dimension(:), contiguous, &
1357  pointer, intent(inout) :: pkgtypes
1358  type(CharacterStringType), dimension(:), contiguous, &
1359  pointer, intent(inout) :: pkgnames
1360  type(CharacterStringType), dimension(:), contiguous, &
1361  pointer, intent(inout) :: mempaths
1362  integer(I4B), dimension(:), contiguous, &
1363  pointer, intent(inout) :: inunits
1364  ! -- local
1365  integer(I4B) :: ipakid, ipaknum
1366  character(len=LENFTYPE) :: pkgtype, bndptype
1367  character(len=LENPACKAGENAME) :: pkgname
1368  character(len=LENMEMPATH) :: mempath
1369  integer(I4B), pointer :: inunit
1370  integer(I4B) :: n
1371 
1372  if (allocated(bndpkgs)) then
1373  !
1374  ! -- create stress packages
1375  ipakid = 1
1376  bndptype = ''
1377  do n = 1, size(bndpkgs)
1378  !
1379  pkgtype = pkgtypes(bndpkgs(n))
1380  pkgname = pkgnames(bndpkgs(n))
1381  mempath = mempaths(bndpkgs(n))
1382  inunit => inunits(bndpkgs(n))
1383  !
1384  if (bndptype /= pkgtype) then
1385  ipaknum = 1
1386  bndptype = pkgtype
1387  end if
1388  !
1389  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1390  inunit, this%iout)
1391  ipakid = ipakid + 1
1392  ipaknum = ipaknum + 1
1393  end do
1394  !
1395  ! -- cleanup
1396  deallocate (bndpkgs)
1397  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23

◆ create_packages()

subroutine gwfmodule::create_packages ( class(gwfmodeltype this)

Definition at line 1402 of file gwf.f90.

1403  ! -- modules
1406  use arrayhandlersmodule, only: expandarray
1407  use memorymanagermodule, only: mem_setptr
1409  use simvariablesmodule, only: idm_context
1410  use dismodule, only: dis_cr
1411  use disvmodule, only: disv_cr
1412  use disumodule, only: disu_cr
1413  use gwfnpfmodule, only: npf_cr
1414  use xt3dmodule, only: xt3d_cr
1415  use gwfbuymodule, only: buy_cr
1416  use gwfvscmodule, only: vsc_cr
1417  use gwfstomodule, only: sto_cr
1418  use gwfcsubmodule, only: csub_cr
1419  use gwfmvrmodule, only: mvr_cr
1420  use gwfhfbmodule, only: hfb_cr
1421  use gwficmodule, only: ic_cr
1422  use gwfocmodule, only: oc_cr
1423  ! -- dummy
1424  class(GwfModelType) :: this
1425  ! -- local
1426  type(CharacterStringType), dimension(:), contiguous, &
1427  pointer :: pkgtypes => null()
1428  type(CharacterStringType), dimension(:), contiguous, &
1429  pointer :: pkgnames => null()
1430  type(CharacterStringType), dimension(:), contiguous, &
1431  pointer :: mempaths => null()
1432  integer(I4B), dimension(:), contiguous, &
1433  pointer :: inunits => null()
1434  character(len=LENMEMPATH) :: model_mempath
1435  character(len=LENFTYPE) :: pkgtype
1436  character(len=LENPACKAGENAME) :: pkgname
1437  character(len=LENMEMPATH) :: mempath
1438  integer(I4B), pointer :: inunit
1439  integer(I4B), dimension(:), allocatable :: bndpkgs
1440  integer(I4B) :: n
1441  integer(I4B) :: indis = 0 ! DIS enabled flag
1442  character(len=LENMEMPATH) :: mempathbuy = ''
1443  character(len=LENMEMPATH) :: mempathcsub = ''
1444  character(len=LENMEMPATH) :: mempathhfb = ''
1445  character(len=LENMEMPATH) :: mempathic = ''
1446  character(len=LENMEMPATH) :: mempathnpf = ''
1447  character(len=LENMEMPATH) :: mempathoc = ''
1448  character(len=LENMEMPATH) :: mempathsto = ''
1449  character(len=LENMEMPATH) :: mempathvsc = ''
1450  !
1451  ! -- set input model memory path
1452  model_mempath = create_mem_path(component=this%name, context=idm_context)
1453  !
1454  ! -- set pointers to model path package info
1455  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1456  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1457  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1458  call mem_setptr(inunits, 'INUNITS', model_mempath)
1459  !
1460  do n = 1, size(pkgtypes)
1461  !
1462  ! attributes for this input package
1463  pkgtype = pkgtypes(n)
1464  pkgname = pkgnames(n)
1465  mempath = mempaths(n)
1466  inunit => inunits(n)
1467  !
1468  ! -- create dis package as it is a prerequisite for other packages
1469  select case (pkgtype)
1470  case ('DIS6')
1471  indis = 1
1472  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1473  case ('DISV6')
1474  indis = 1
1475  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1476  case ('DISU6')
1477  indis = 1
1478  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1479  case ('NPF6')
1480  this%innpf = 1
1481  mempathnpf = mempath
1482  case ('BUY6')
1483  this%inbuy = 1
1484  mempathbuy = mempath
1485  case ('VSC6')
1486  this%invsc = 1
1487  mempathvsc = mempath
1488  case ('GNC6')
1489  this%ingnc = inunit
1490  case ('HFB6')
1491  this%inhfb = 1
1492  mempathhfb = mempath
1493  case ('STO6')
1494  this%insto = 1
1495  mempathsto = mempath
1496  case ('CSUB6')
1497  this%incsub = 1
1498  mempathcsub = mempath
1499  case ('IC6')
1500  this%inic = 1
1501  mempathic = mempath
1502  case ('MVR6')
1503  this%inmvr = inunit
1504  case ('OC6')
1505  this%inoc = 1
1506  mempathoc = mempath
1507  case ('OBS6')
1508  this%inobs = inunit
1509  case ('WEL6', 'DRN6', 'RIV6', 'GHB6', 'RCH6', &
1510  'EVT6', 'API6', 'CHD6', 'MAW6', 'SFR6', &
1511  'LAK6', 'UZF6')
1512  call expandarray(bndpkgs)
1513  bndpkgs(size(bndpkgs)) = n
1514  case default
1515  ! TODO
1516  end select
1517  end do
1518  !
1519  ! -- Create packages that are tied directly to model
1520  call npf_cr(this%npf, this%name, mempathnpf, this%innpf, this%iout)
1521  call xt3d_cr(this%xt3d, this%name, this%innpf, this%iout)
1522  call buy_cr(this%buy, this%name, mempathbuy, this%inbuy, this%iout)
1523  call vsc_cr(this%vsc, this%name, mempathvsc, this%invsc, this%iout)
1524  call gnc_cr(this%gnc, this%name, this%ingnc, this%iout)
1525  call hfb_cr(this%hfb, this%name, mempathhfb, this%inhfb, this%iout)
1526  call sto_cr(this%sto, this%name, mempathsto, this%insto, this%iout)
1527  call csub_cr(this%csub, this%name, mempathcsub, this%insto, &
1528  this%sto%packName, this%incsub, this%iout)
1529  call ic_cr(this%ic, this%name, mempathic, this%inic, this%iout, this%dis)
1530  call mvr_cr(this%mvr, this%name, this%inmvr, this%iout, this%dis)
1531  call oc_cr(this%oc, this%name, mempathoc, this%inoc, this%iout)
1532  call gwf_obs_cr(this%obs, this%inobs)
1533  !
1534  ! -- Check to make sure that required ftype's have been specified
1535  call this%ftype_check(indis)
1536  !
1537  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
Definition: Dis.f90:1
subroutine, public dis_cr(dis, name_model, input_mempath, inunit, iout)
Create a new structured discretization object.
Definition: Dis.f90:99
subroutine, public disu_cr(dis, name_model, input_mempath, inunit, iout)
Create a new unstructured discretization object.
Definition: Disu.f90:127
subroutine, public disv_cr(dis, name_model, input_mempath, inunit, iout)
Create a new discretization by vertices object.
Definition: Disv.f90:111
subroutine, public buy_cr(buyobj, name_model, input_mempath, inunit, iout)
Create a new BUY object.
Definition: gwf-buy.f90:104
This module contains the CSUB package methods.
Definition: gwf-csub.f90:9
subroutine, public csub_cr(csubobj, name_model, mempath, istounit, stoPckName, inunit, iout)
@ brief Create a new package object
Definition: gwf-csub.f90:320
subroutine, public hfb_cr(hfbobj, name_model, input_mempath, inunit, iout)
Create a new hfb object.
Definition: gwf-hfb.f90:71
subroutine, public ic_cr(ic, name_model, input_mempath, inunit, iout, dis)
Create a new initial conditions object.
Definition: gwf-ic.f90:33
subroutine, public mvr_cr(mvrobj, name_parent, inunit, iout, dis, iexgmvr)
Create a new mvr object.
Definition: gwf-mvr.f90:186
subroutine, public npf_cr(npfobj, name_model, input_mempath, inunit, iout)
Create a new NPF object. Pass a inunit value of 0 if npf data will initialized from memory.
Definition: gwf-npf.f90:160
subroutine, public oc_cr(ocobj, name_model, input_mempath, inunit, iout)
@ brief Create GwfOcType
Definition: gwf-oc.f90:32
This module contains the storage package methods.
Definition: gwf-sto.f90:8
subroutine, public sto_cr(stoobj, name_model, mempath, inunit, iout)
@ brief Create a new package object
Definition: gwf-sto.f90:78
subroutine, public vsc_cr(vscobj, name_model, input_mempath, inunit, iout)
@ brief Create a new package object
Definition: gwf-vsc.f90:143
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
subroutine, public xt3d_cr(xt3dobj, name_model, inunit, iout, ldispopt)
Create a new xt3d object.
Here is the call graph for this function:

◆ ftype_check()

subroutine gwfmodule::ftype_check ( class(gwfmodeltype this,
integer(i4b), intent(in)  indis 
)

Definition at line 1298 of file gwf.f90.

1299  ! -- modules
1300  use constantsmodule, only: linelength
1301  use simmodule, only: store_error, count_errors
1302  ! -- dummy
1303  class(GwfModelType) :: this
1304  integer(I4B), intent(in) :: indis
1305  ! -- local
1306  !
1307  ! -- Check for IC8, DIS(u), and NPF. Stop if not present.
1308  if (this%inic == 0) then
1309  write (errmsg, '(a)') &
1310  'Initial Conditions (IC6) package not specified.'
1311  call store_error(errmsg)
1312  end if
1313  if (indis == 0) then
1314  write (errmsg, '(a)') &
1315  'Discretization (DIS6, DISV6, or DISU6) Package not specified.'
1316  call store_error(errmsg)
1317  end if
1318  if (this%innpf == 0) then
1319  write (errmsg, '(a)') &
1320  'Node Property Flow (NPF6) Package not specified.'
1321  call store_error(errmsg)
1322  end if
1323  !
1324  if (count_errors() > 0) then
1325  write (errmsg, '(a)') 'One or more required package(s) not specified.'
1326  call store_error(errmsg)
1327  call store_error_filename(this%filename)
1328  end if
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
Here is the call graph for this function:

◆ gwf_ac()

subroutine gwfmodule::gwf_ac ( class(gwfmodeltype this,
type(sparsematrix), intent(inout)  sparse 
)
private

Definition at line 255 of file gwf.f90.

256  ! -- modules
257  use sparsemodule, only: sparsematrix
258  ! -- dummy
259  class(GwfModelType) :: this
260  type(sparsematrix), intent(inout) :: sparse
261  ! -- local
262  class(BndType), pointer :: packobj
263  integer(I4B) :: ip
264  !
265  ! -- Add the primary grid connections of this model to sparse
266  call this%dis%dis_ac(this%moffset, sparse)
267  !
268  ! -- Add any additional connections that NPF may need
269  if (this%innpf > 0) call this%npf%npf_ac(this%moffset, sparse)
270  !
271  ! -- Add any package connections
272  do ip = 1, this%bndlist%Count()
273  packobj => getbndfromlist(this%bndlist, ip)
274  call packobj%bnd_ac(this%moffset, sparse)
275  end do
276  !
277  ! -- If GNC is active, then add the gnc connections to sparse
278  if (this%ingnc > 0) call this%gnc%gnc_ac(sparse)
Here is the call graph for this function:

◆ gwf_ad()

subroutine gwfmodule::gwf_ad ( class(gwfmodeltype this)

(1) calls package advance subroutines

Definition at line 397 of file gwf.f90.

398  ! -- modules
400  ! -- dummy
401  class(GwfModelType) :: this
402  class(BndType), pointer :: packobj
403  ! -- local
404  integer(I4B) :: irestore
405  integer(I4B) :: ip, n
406  !
407  ! -- Reset state variable
408  irestore = 0
409  if (ifailedstepretry > 0) irestore = 1
410  if (irestore == 0) then
411  !
412  ! -- copy x into xold
413  do n = 1, this%dis%nodes
414  this%xold(n) = this%x(n)
415  end do
416  else
417  !
418  ! -- copy xold into x if this time step is a redo
419  do n = 1, this%dis%nodes
420  this%x(n) = this%xold(n)
421  end do
422  end if
423  !
424  ! -- Advance
425  if (this%invsc > 0) call this%vsc%vsc_ad()
426  if (this%innpf > 0) call this%npf%npf_ad(this%dis%nodes, this%xold, &
427  this%x, irestore)
428  if (this%insto > 0) call this%sto%sto_ad()
429  if (this%incsub > 0) call this%csub%csub_ad(this%dis%nodes, this%x)
430  if (this%inbuy > 0) call this%buy%buy_ad()
431  if (this%inmvr > 0) call this%mvr%mvr_ad()
432  do ip = 1, this%bndlist%Count()
433  packobj => getbndfromlist(this%bndlist, ip)
434  call packobj%bnd_ad()
435  if (this%invsc > 0) call this%vsc%vsc_ad_bnd(packobj, this%x)
436  if (isimcheck > 0) then
437  call packobj%bnd_ck()
438  end if
439  end do
440  !
441  ! -- Push simulated values to preceding time/subtime step
442  call this%obs%obs_ad()
integer(i4b) isimcheck
simulation input check flag (1) to check input, (0) to ignore checks
integer(i4b) ifailedstepretry
current retry for this time step
Here is the call graph for this function:

◆ gwf_ar()

subroutine gwfmodule::gwf_ar ( class(gwfmodeltype this)
private

(1) allocates and reads packages part of this model, (2) allocates memory for arrays part of this model object

Definition at line 316 of file gwf.f90.

317  ! -- dummy
318  class(GwfModelType) :: this
319  ! -- locals
320  integer(I4B) :: ip
321  class(BndType), pointer :: packobj
322  !
323  ! -- Allocate and read modules attached to model
324  if (this%inic > 0) call this%ic%ic_ar(this%x)
325  if (this%innpf > 0) call this%npf%npf_ar(this%ic, this%vsc, this%ibound, &
326  this%x)
327  if (this%invsc > 0) call this%vsc%vsc_ar(this%ibound)
328  if (this%inbuy > 0) call this%buy%buy_ar(this%npf, this%ibound)
329  if (this%inhfb > 0) call this%hfb%hfb_ar(this%ibound, this%xt3d, this%dis, &
330  this%invsc, this%vsc)
331  if (this%insto > 0) call this%sto%sto_ar(this%dis, this%ibound)
332  if (this%incsub > 0) call this%csub%csub_ar(this%dis, this%ibound)
333  if (this%inmvr > 0) call this%mvr%mvr_ar()
334  if (this%inobs > 0) call this%obs%gwf_obs_ar(this%ic, this%x, this%flowja)
335  !
336  ! -- Call dis_ar to write binary grid file
337  call this%dis%dis_ar(this%npf%icelltype)
338  !
339  ! -- set up output control
340  call this%oc%oc_ar(this%x, this%dis, this%npf%hnoflo)
341  call this%budget%set_ibudcsv(this%oc%ibudcsv)
342  !
343  ! -- Package input files now open, so allocate and read
344  do ip = 1, this%bndlist%Count()
345  packobj => getbndfromlist(this%bndlist, ip)
346  call packobj%set_pointers(this%dis%nodes, this%ibound, this%x, &
347  this%xold, this%flowja)
348  ! -- Read and allocate package
349  call packobj%bnd_ar()
350  if (this%inbuy > 0) call this%buy%buy_ar_bnd(packobj, this%x)
351  if (this%invsc > 0) call this%vsc%vsc_ar_bnd(packobj)
352  end do
Here is the call graph for this function:

◆ gwf_bd()

subroutine gwfmodule::gwf_bd ( class(gwfmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)
private

(1) Calculate stress package contributions to model budget

Definition at line 786 of file gwf.f90.

787  ! -- modules
788  use sparsemodule, only: csr_diagsum
789  ! -- dummy
790  class(GwfModelType) :: this
791  integer(I4B), intent(in) :: icnvg
792  integer(I4B), intent(in) :: isuppress_output
793  ! -- local
794  integer(I4B) :: ip
795  class(BndType), pointer :: packobj
796  !
797  ! -- Finalize calculation of flowja by adding face flows to the diagonal.
798  ! This results in the flow residual being stored in the diagonal
799  ! position for each cell.
800  call csr_diagsum(this%dis%con%ia, this%flowja)
801  !
802  ! -- Save the solution convergence flag
803  this%icnvg = icnvg
804  !
805  ! -- Budget routines (start by resetting). Sole purpose of this section
806  ! is to add in and outs to model budget. All ins and out for a model
807  ! should be added here to this%budget. In a subsequent exchange call,
808  ! exchange flows might also be added.
809  call this%budget%reset()
810  if (this%insto > 0) call this%sto%sto_bd(isuppress_output, this%budget)
811  if (this%incsub > 0) call this%csub%csub_bd(isuppress_output, this%budget)
812  if (this%inmvr > 0) call this%mvr%mvr_bd()
813  do ip = 1, this%bndlist%Count()
814  packobj => getbndfromlist(this%bndlist, ip)
815  call packobj%bnd_bd(this%budget)
816  end do
817  !
818  ! -- npf velocities have to be calculated here, after gwf-gwf exchanges
819  ! have passed in their contributions from exg_cq()
820  if (this%innpf > 0) then
821  if (this%npf%icalcspdis /= 0) then
822  call this%npf%calc_spdis(this%flowja)
823  end if
824  end if
subroutine csr_diagsum(ia, flowja)
Definition: Sparse.f90:263
Here is the call graph for this function:

◆ gwf_bdentry()

subroutine gwfmodule::gwf_bdentry ( class(gwfmodeltype this,
real(dp), dimension(:, :), intent(in)  budterm,
character(len=lenbudtxt), dimension(:), intent(in)  budtxt,
character(len=*), intent(in)  rowlabel 
)

This subroutine adds a budget entry to the flow budget. It was added as a method for the gwf model object so that the exchange object could add its

Definition at line 1121 of file gwf.f90.

1122  ! -- modules
1123  use constantsmodule, only: lenbudtxt
1124  use tdismodule, only: delt
1125  ! -- dummy
1126  class(GwfModelType) :: this
1127  real(DP), dimension(:, :), intent(in) :: budterm
1128  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
1129  character(len=*), intent(in) :: rowlabel
1130  !
1131  call this%budget%addentry(budterm, delt, budtxt, rowlabel=rowlabel)
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29

◆ gwf_cc()

subroutine gwfmodule::gwf_cc ( class(gwfmodeltype this,
integer(i4b), intent(in)  innertot,
integer(i4b), intent(in)  kiter,
integer(i4b), intent(in)  iend,
integer(i4b), intent(in)  icnvgmod,
character(len=lenpakloc), intent(inout)  cpak,
integer(i4b), intent(inout)  ipak,
real(dp), intent(inout)  dpak 
)
private

(1) calls package cc routines

Definition at line 563 of file gwf.f90.

564  ! -- dummy
565  class(GwfModelType) :: this
566  integer(I4B), intent(in) :: innertot
567  integer(I4B), intent(in) :: kiter
568  integer(I4B), intent(in) :: iend
569  integer(I4B), intent(in) :: icnvgmod
570  character(len=LENPAKLOC), intent(inout) :: cpak
571  integer(I4B), intent(inout) :: ipak
572  real(DP), intent(inout) :: dpak
573  ! -- local
574  class(BndType), pointer :: packobj
575  integer(I4B) :: ip
576  ! -- formats
577  !
578  ! -- If mover is on, then at least 2 outers required
579  if (this%inmvr > 0) then
580  call this%mvr%mvr_cc(innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
581  end if
582  !
583  ! -- csub convergence check
584  if (this%incsub > 0) then
585  call this%csub%csub_cc(innertot, kiter, iend, icnvgmod, &
586  this%dis%nodes, this%x, this%xold, &
587  cpak, ipak, dpak)
588  end if
589  !
590  ! -- Call package cc routines
591  do ip = 1, this%bndlist%Count()
592  packobj => getbndfromlist(this%bndlist, ip)
593  call packobj%bnd_cc(innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
594  end do
Here is the call graph for this function:

◆ gwf_cf()

subroutine gwfmodule::gwf_cf ( class(gwfmodeltype this,
integer(i4b), intent(in)  kiter 
)

Definition at line 447 of file gwf.f90.

448  ! -- dummy
449  class(GwfModelType) :: this
450  integer(I4B), intent(in) :: kiter
451  ! -- local
452  class(BndType), pointer :: packobj
453  integer(I4B) :: ip
454  !
455  ! -- Call package cf routines
456  if (this%innpf > 0) call this%npf%npf_cf(kiter, this%dis%nodes, this%x)
457  if (this%inbuy > 0) call this%buy%buy_cf(kiter)
458  do ip = 1, this%bndlist%Count()
459  packobj => getbndfromlist(this%bndlist, ip)
460  call packobj%bnd_cf()
461  if (this%inbuy > 0) call this%buy%buy_cf_bnd(packobj, this%x)
462  end do
Here is the call graph for this function:

◆ gwf_cq()

subroutine gwfmodule::gwf_cq ( class(gwfmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)

(1) Calculate intercell flows (flowja)

Definition at line 742 of file gwf.f90.

743  ! -- modules
744  ! -- dummy
745  class(GwfModelType) :: this
746  integer(I4B), intent(in) :: icnvg
747  integer(I4B), intent(in) :: isuppress_output
748  ! -- local
749  integer(I4B) :: i
750  integer(I4B) :: ip
751  class(BndType), pointer :: packobj
752  !
753  ! -- Construct the flowja array. Flowja is calculated each time, even if
754  ! output is suppressed. (flowja is positive into a cell.) The diagonal
755  ! position of the flowja array will contain the flow residual after
756  ! these routines are called, so each package is responsible for adding
757  ! its flow to this diagonal position.
758  do i = 1, this%nja
759  this%flowja(i) = dzero
760  end do
761  if (this%innpf > 0) call this%npf%npf_cq(this%x, this%flowja)
762  if (this%inbuy > 0) call this%buy%buy_cq(this%x, this%flowja)
763  if (this%inhfb > 0) call this%hfb%hfb_cq(this%x, this%flowja)
764  if (this%ingnc > 0) call this%gnc%gnc_cq(this%flowja)
765  if (this%insto > 0) call this%sto%sto_cq(this%flowja, this%x, this%xold)
766  if (this%incsub > 0) call this%csub%csub_cq(this%dis%nodes, this%x, &
767  this%xold, isuppress_output, &
768  this%flowja)
769  !
770  ! -- Go through packages and call cq routines. cf() routines are called
771  ! first to regenerate non-linear terms to be consistent with the final
772  ! head solution.
773  do ip = 1, this%bndlist%Count()
774  packobj => getbndfromlist(this%bndlist, ip)
775  call packobj%bnd_cf()
776  if (this%inbuy > 0) call this%buy%buy_cf_bnd(packobj, this%x)
777  call packobj%bnd_cq(this%x, this%flowja)
778  end do
Here is the call graph for this function:

◆ gwf_cr()

subroutine, public gwfmodule::gwf_cr ( character(len=*), intent(in)  filename,
integer(i4b), intent(in)  id,
character(len=*), intent(in)  modelname 
)

(1) creates model object and add to modellist (2) assign values

Parameters
[in]filenameinput file
[in]idconsecutive model number listed in mfsim.nam
[in]modelnamename of the model

Definition at line 137 of file gwf.f90.

138  ! -- modules
139  use listsmodule, only: basemodellist
141  use constantsmodule, only: linelength
146  use budgetmodule, only: budget_cr
147  ! -- dummy
148  character(len=*), intent(in) :: filename !< input file
149  integer(I4B), intent(in) :: id !< consecutive model number listed in mfsim.nam
150  character(len=*), intent(in) :: modelname !< name of the model
151  ! -- local
152  type(GwfModelType), pointer :: this
153  class(BaseModelType), pointer :: model
154  character(len=LENMEMPATH) :: input_mempath
155  character(len=LINELENGTH) :: lst_fname
156  type(GwfNamParamFoundType) :: found
157  ! -- format
158  !
159  ! -- Allocate a new GWF Model (this) and add it to basemodellist
160  allocate (this)
161  !
162  ! -- Set memory path before allocation in memory manager can be done
163  this%memoryPath = create_mem_path(modelname)
164  !
165  call this%allocate_scalars(modelname)
166  model => this
167  call addbasemodeltolist(basemodellist, model)
168  !
169  ! -- Assign values
170  this%filename = filename
171  this%name = modelname
172  this%macronym = 'GWF'
173  this%id = id
174  !
175  ! -- set input model namfile memory path
176  input_mempath = create_mem_path(modelname, 'NAM', idm_context)
177  !
178  ! -- copy option params from input context
179  call mem_set_value(lst_fname, 'LIST', input_mempath, found%list)
180  call mem_set_value(this%inewton, 'NEWTON', input_mempath, found%newton)
181  call mem_set_value(this%inewtonur, 'UNDER_RELAXATION', input_mempath, &
182  found%under_relaxation)
183  call mem_set_value(this%iprpak, 'PRINT_INPUT', input_mempath, &
184  found%print_input)
185  call mem_set_value(this%iprflow, 'PRINT_FLOWS', input_mempath, &
186  found%print_flows)
187  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', input_mempath, found%save_flows)
188  !
189  ! -- create the list file
190  call this%create_lstfile(lst_fname, filename, found%list, &
191  'GROUNDWATER FLOW MODEL (GWF)')
192  !
193  ! -- activate save_flows if found
194  if (found%save_flows) then
195  this%ipakcb = -1
196  end if
197  !
198  ! -- log set options
199  if (this%iout > 0) then
200  call this%log_namfile_options(found)
201  end if
202  !
203  ! -- Create utility objects
204  call budget_cr(this%budget, this%name)
205  !
206  ! -- create model packages
207  call this%create_packages()
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:160
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public budget_cr(this, name_model)
@ brief Create a new budget object
Definition: Budget.f90:84
type(listtype), public basemodellist
Definition: mf6lists.f90:16
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gwf_da()

subroutine gwfmodule::gwf_da ( class(gwfmodeltype this)
private

Definition at line 1043 of file gwf.f90.

1044  ! -- modules
1047  use simvariablesmodule, only: idm_context
1048  ! -- dummy
1049  class(GwfModelType) :: this
1050  ! -- local
1051  integer(I4B) :: ip
1052  class(BndType), pointer :: packobj
1053  !
1054  ! -- Deallocate idm memory
1055  call memorystore_remove(this%name, 'NAM', idm_context)
1056  call memorystore_remove(component=this%name, context=idm_context)
1057  !
1058  ! -- Internal flow packages deallocate
1059  call this%dis%dis_da()
1060  call this%ic%ic_da()
1061  call this%npf%npf_da()
1062  call this%xt3d%xt3d_da()
1063  call this%buy%buy_da()
1064  call this%vsc%vsc_da()
1065  call this%gnc%gnc_da()
1066  call this%sto%sto_da()
1067  call this%csub%csub_da()
1068  call this%budget%budget_da()
1069  call this%hfb%hfb_da()
1070  call this%mvr%mvr_da()
1071  call this%oc%oc_da()
1072  call this%obs%obs_da()
1073  !
1074  ! -- Internal package objects
1075  deallocate (this%dis)
1076  deallocate (this%ic)
1077  deallocate (this%npf)
1078  deallocate (this%xt3d)
1079  deallocate (this%buy)
1080  deallocate (this%vsc)
1081  deallocate (this%gnc)
1082  deallocate (this%sto)
1083  deallocate (this%csub)
1084  deallocate (this%budget)
1085  deallocate (this%hfb)
1086  deallocate (this%mvr)
1087  deallocate (this%obs)
1088  deallocate (this%oc)
1089  !
1090  ! -- Boundary packages
1091  do ip = 1, this%bndlist%Count()
1092  packobj => getbndfromlist(this%bndlist, ip)
1093  call packobj%bnd_da()
1094  deallocate (packobj)
1095  end do
1096  !
1097  ! -- Scalars
1098  call mem_deallocate(this%inic)
1099  call mem_deallocate(this%inoc)
1100  call mem_deallocate(this%inobs)
1101  call mem_deallocate(this%innpf)
1102  call mem_deallocate(this%inbuy)
1103  call mem_deallocate(this%invsc)
1104  call mem_deallocate(this%insto)
1105  call mem_deallocate(this%incsub)
1106  call mem_deallocate(this%inmvr)
1107  call mem_deallocate(this%inhfb)
1108  call mem_deallocate(this%ingnc)
1109  call mem_deallocate(this%iss)
1110  call mem_deallocate(this%inewtonur)
1111  !
1112  ! -- NumericalModelType
1113  call this%NumericalModelType%model_da()
subroutine, public memorystore_remove(component, subcomponent, context)
Here is the call graph for this function:

◆ gwf_df()

subroutine gwfmodule::gwf_df ( class(gwfmodeltype this)

(1) call df routines for each package (2) set gwf variables and pointers

Definition at line 216 of file gwf.f90.

217  ! -- modules
218  ! -- dummy
219  class(GwfModelType) :: this
220  ! -- local
221  integer(I4B) :: ip
222  class(BndType), pointer :: packobj
223  !
224  ! -- Define packages and utility objects
225  call this%dis%dis_df()
226  call this%npf%npf_df(this%dis, this%xt3d, this%ingnc, this%invsc)
227  call this%oc%oc_df()
228  call this%budget%budget_df(niunit_gwf, 'VOLUME', 'L**3')
229  if (this%inbuy > 0) call this%buy%buy_df(this%dis)
230  if (this%invsc > 0) call this%vsc%vsc_df(this%dis)
231  if (this%ingnc > 0) call this%gnc%gnc_df(this)
232  !
233  ! -- Assign or point model members to dis members
234  ! this%neq will be incremented if packages add additional unknowns
235  this%neq = this%dis%nodes
236  this%nja = this%dis%nja
237  this%ia => this%dis%con%ia
238  this%ja => this%dis%con%ja
239  !
240  ! -- Allocate model arrays, now that neq and nja are known
241  call this%allocate_arrays()
242  !
243  ! -- Define packages and assign iout for time series managers
244  do ip = 1, this%bndlist%Count()
245  packobj => getbndfromlist(this%bndlist, ip)
246  call packobj%bnd_df(this%neq, this%dis)
247  end do
248  !
249  ! -- Store information needed for observations
250  call this%obs%obs_df(this%iout, this%name, 'GWF', this%dis)
Here is the call graph for this function:

◆ gwf_fc()

subroutine gwfmodule::gwf_fc ( class(gwfmodeltype this,
integer(i4b), intent(in)  kiter,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), intent(in)  inwtflag 
)
private

Definition at line 467 of file gwf.f90.

468  ! -- dummy
469  class(GwfModelType) :: this
470  integer(I4B), intent(in) :: kiter
471  class(MatrixBaseType), pointer :: matrix_sln
472  integer(I4B), intent(in) :: inwtflag
473  ! -- local
474  class(BndType), pointer :: packobj
475  integer(I4B) :: ip
476  integer(I4B) :: inwt, inwtsto, inwtcsub, inwtpak
477  !
478  ! -- newton flags
479  inwt = inwtflag
480  if (inwtflag == 1) inwt = this%npf%inewton
481  inwtsto = inwtflag
482  if (this%insto > 0) then
483  if (inwtflag == 1) inwtsto = this%sto%inewton
484  end if
485  inwtcsub = inwtflag
486  if (this%incsub > 0) then
487  if (inwtflag == 1) inwtcsub = this%csub%inewton
488  end if
489  !
490  ! -- Fill standard conductance terms
491  if (this%innpf > 0) call this%npf%npf_fc(kiter, matrix_sln, this%idxglo, &
492  this%rhs, this%x)
493  if (this%inbuy > 0) call this%buy%buy_fc(kiter, matrix_sln, this%idxglo, &
494  this%rhs, this%x)
495  if (this%inhfb > 0) call this%hfb%hfb_fc(kiter, matrix_sln, this%idxglo, &
496  this%rhs, this%x)
497  if (this%ingnc > 0) call this%gnc%gnc_fc(kiter, matrix_sln)
498  ! -- storage
499  if (this%insto > 0) then
500  call this%sto%sto_fc(kiter, this%xold, this%x, matrix_sln, &
501  this%idxglo, this%rhs)
502  end if
503  ! -- skeletal storage, compaction, and land subsidence
504  if (this%incsub > 0) then
505  call this%csub%csub_fc(kiter, this%xold, this%x, matrix_sln, &
506  this%idxglo, this%rhs)
507  end if
508  if (this%inmvr > 0) call this%mvr%mvr_fc()
509  do ip = 1, this%bndlist%Count()
510  packobj => getbndfromlist(this%bndlist, ip)
511  call packobj%bnd_fc(this%rhs, this%ia, this%idxglo, matrix_sln)
512  end do
513  !
514  !--Fill newton terms
515  if (this%innpf > 0) then
516  if (inwt /= 0) then
517  call this%npf%npf_fn(kiter, matrix_sln, this%idxglo, this%rhs, this%x)
518  end if
519  end if
520  !
521  ! -- Fill newton terms for ghost nodes
522  if (this%ingnc > 0) then
523  if (inwt /= 0) then
524  call this%gnc%gnc_fn(kiter, matrix_sln, this%npf%condsat, &
525  ivarcv_opt=this%npf%ivarcv, &
526  ictm1_opt=this%npf%icelltype, &
527  ictm2_opt=this%npf%icelltype)
528  end if
529  end if
530  !
531  ! -- Fill newton terms for storage
532  if (this%insto > 0) then
533  if (inwtsto /= 0) then
534  call this%sto%sto_fn(kiter, this%xold, this%x, matrix_sln, &
535  this%idxglo, this%rhs)
536  end if
537  end if
538  !
539  ! -- Fill newton terms for skeletal storage, compaction, and land subsidence
540  if (this%incsub > 0) then
541  if (inwtcsub /= 0) then
542  call this%csub%csub_fn(kiter, this%xold, this%x, matrix_sln, &
543  this%idxglo, this%rhs)
544  end if
545  end if
546  !
547  ! -- Fill Newton terms for packages
548  do ip = 1, this%bndlist%Count()
549  packobj => getbndfromlist(this%bndlist, ip)
550  inwtpak = inwtflag
551  if (inwtflag == 1) inwtpak = packobj%inewton
552  if (inwtpak /= 0) then
553  call packobj%bnd_fn(this%rhs, this%ia, this%idxglo, matrix_sln)
554  end if
555  end do
Here is the call graph for this function:

◆ gwf_fp()

subroutine gwfmodule::gwf_fp ( class(gwfmodeltype this)

Definition at line 1029 of file gwf.f90.

1030  ! -- modules
1031  ! -- dummy
1032  class(GwfModelType) :: this
1033  ! -- local
1034  !
1035  ! -- csub final processing
1036  if (this%incsub > 0) then
1037  call this%csub%csub_fp()
1038  end if

◆ gwf_get_iasym()

integer(i4b) function gwfmodule::gwf_get_iasym ( class(gwfmodeltype this)

Definition at line 1137 of file gwf.f90.

1138  class(GwfModelType) :: this
1139  ! -- local
1140  integer(I4B) :: iasym
1141  integer(I4B) :: ip
1142  class(BndType), pointer :: packobj
1143  !
1144  ! -- Start by setting iasym to zero
1145  iasym = 0
1146  !
1147  ! -- NPF
1148  if (this%innpf > 0) then
1149  if (this%npf%iasym /= 0) iasym = 1
1150  if (this%npf%ixt3d /= 0) iasym = 1
1151  end if
1152  !
1153  ! -- GNC
1154  if (this%ingnc > 0) then
1155  if (this%gnc%iasym /= 0) iasym = 1
1156  end if
1157  !
1158  ! -- Check for any packages that introduce matrix asymmetry
1159  do ip = 1, this%bndlist%Count()
1160  packobj => getbndfromlist(this%bndlist, ip)
1161  if (packobj%iasym /= 0) iasym = 1
1162  end do
Here is the call graph for this function:

◆ gwf_mc()

subroutine gwfmodule::gwf_mc ( class(gwfmodeltype this,
class(matrixbasetype), pointer  matrix_sln 
)

Definition at line 284 of file gwf.f90.

285  ! -- dummy
286  class(GwfModelType) :: this
287  class(MatrixBaseType), pointer :: matrix_sln
288  ! -- local
289  class(BndType), pointer :: packobj
290  integer(I4B) :: ip
291  !
292  ! -- Find the position of each connection in the global ia, ja structure
293  ! and store them in idxglo.
294  call this%dis%dis_mc(this%moffset, this%idxglo, matrix_sln)
295  !
296  ! -- Map any additional connections that NPF may need
297  if (this%innpf > 0) call this%npf%npf_mc(this%moffset, matrix_sln)
298  !
299  ! -- Map any package connections
300  do ip = 1, this%bndlist%Count()
301  packobj => getbndfromlist(this%bndlist, ip)
302  call packobj%bnd_mc(this%moffset, matrix_sln)
303  end do
304  !
305  ! -- For implicit gnc, need to store positions of gnc connections
306  ! in solution matrix connection
307  if (this%ingnc > 0) call this%gnc%gnc_mc(matrix_sln)
Here is the call graph for this function:

◆ gwf_nur()

subroutine gwfmodule::gwf_nur ( class(gwfmodeltype this,
integer(i4b), intent(in)  neqmod,
real(dp), dimension(neqmod), intent(inout)  x,
real(dp), dimension(neqmod), intent(in)  xtemp,
real(dp), dimension(neqmod), intent(inout)  dx,
integer(i4b), intent(inout)  inewtonur,
real(dp), intent(inout)  dxmax,
integer(i4b), intent(inout)  locmax 
)

(1) Under-relaxation of Groundwater Flow Model Heads for current outer iteration using the cell bottoms at the bottom of the model

Definition at line 697 of file gwf.f90.

698  ! modules
699  use constantsmodule, only: done, dp9
700  ! -- dummy
701  class(GwfModelType) :: this
702  integer(I4B), intent(in) :: neqmod
703  real(DP), dimension(neqmod), intent(inout) :: x
704  real(DP), dimension(neqmod), intent(in) :: xtemp
705  real(DP), dimension(neqmod), intent(inout) :: dx
706  integer(I4B), intent(inout) :: inewtonur
707  real(DP), intent(inout) :: dxmax
708  integer(I4B), intent(inout) :: locmax
709  ! -- local
710  integer(I4B) :: i0
711  integer(I4B) :: i1
712  class(BndType), pointer :: packobj
713  integer(I4B) :: ip
714  !
715  ! -- apply Newton-Raphson under-relaxation if model is using
716  ! the Newton-Raphson formulation and this Newton-Raphson
717  ! under-relaxation is turned on.
718  if (this%inewton /= 0 .and. this%inewtonur /= 0) then
719  if (this%innpf > 0) then
720  call this%npf%npf_nur(neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
721  end if
722  !
723  ! -- Call package nur routines
724  i0 = this%dis%nodes + 1
725  do ip = 1, this%bndlist%Count()
726  packobj => getbndfromlist(this%bndlist, ip)
727  if (packobj%npakeq > 0) then
728  i1 = i0 + packobj%npakeq - 1
729  call packobj%bnd_nur(packobj%npakeq, x(i0:i1), xtemp(i0:i1), &
730  dx(i0:i1), inewtonur, dxmax, locmax)
731  i0 = i1 + 1
732  end if
733  end do
734  end if
real(dp), parameter dp9
real constant 9/10
Definition: Constants.f90:72
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
Here is the call graph for this function:

◆ gwf_ot()

subroutine gwfmodule::gwf_ot ( class(gwfmodeltype this)

Definition at line 829 of file gwf.f90.

830  ! -- modules
831  use tdismodule, only: kstp, kper, tdis_ot, endofperiod
832  ! -- dummy
833  class(GwfModelType) :: this
834  ! -- local
835  integer(I4B) :: idvsave
836  integer(I4B) :: idvprint
837  integer(I4B) :: icbcfl
838  integer(I4B) :: icbcun
839  integer(I4B) :: ibudfl
840  integer(I4B) :: ipflag
841  ! -- formats
842  character(len=*), parameter :: fmtnocnvg = &
843  "(1X,/9X,'****FAILED TO MEET SOLVER CONVERGENCE CRITERIA IN TIME STEP ', &
844  &I0,' OF STRESS PERIOD ',I0,'****')"
845  !
846  ! -- Set write and print flags
847  idvsave = 0
848  idvprint = 0
849  icbcfl = 0
850  ibudfl = 0
851  if (this%oc%oc_save('HEAD')) idvsave = 1
852  if (this%oc%oc_print('HEAD')) idvprint = 1
853  if (this%oc%oc_save('BUDGET')) icbcfl = 1
854  if (this%oc%oc_print('BUDGET')) ibudfl = 1
855  icbcun = this%oc%oc_save_unit('BUDGET')
856  !
857  ! -- Override ibudfl and idvprint flags for nonconvergence
858  ! and end of period
859  ibudfl = this%oc%set_print_flag('BUDGET', this%icnvg, endofperiod)
860  idvprint = this%oc%set_print_flag('HEAD', this%icnvg, endofperiod)
861  !
862  ! Calculate and save observations
863  call this%gwf_ot_obs()
864  !
865  ! Save and print flows
866  call this%gwf_ot_flow(icbcfl, ibudfl, icbcun)
867  !
868  ! Save and print dependent variables
869  call this%gwf_ot_dv(idvsave, idvprint, ipflag)
870  !
871  ! Print budget summaries
872  call this%gwf_ot_bdsummary(ibudfl, ipflag)
873  !
874  ! -- Timing Output; if any dependent variables or budgets
875  ! are printed, then ipflag is set to 1.
876  if (ipflag == 1) call tdis_ot(this%iout)
877  !
878  ! -- Write non-convergence message
879  if (this%icnvg == 0) then
880  write (this%iout, fmtnocnvg) kstp, kper
881  end if
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:27
subroutine, public tdis_ot(iout)
Print simulation time.
Definition: tdis.f90:274
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ gwf_ot_bdsummary()

subroutine gwfmodule::gwf_ot_bdsummary ( class(gwfmodeltype this,
integer(i4b), intent(in)  ibudfl,
integer(i4b), intent(inout)  ipflag 
)
private

Definition at line 996 of file gwf.f90.

997  use tdismodule, only: kstp, kper, totim, delt
998  class(GwfModelType) :: this
999  integer(I4B), intent(in) :: ibudfl
1000  integer(I4B), intent(inout) :: ipflag
1001  class(BndType), pointer :: packobj
1002  integer(I4B) :: ip
1003 
1004  ! -- Package budget summary
1005  do ip = 1, this%bndlist%Count()
1006  packobj => getbndfromlist(this%bndlist, ip)
1007  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
1008  end do
1009 
1010  ! -- mover budget summary
1011  if (this%inmvr > 0) then
1012  call this%mvr%mvr_ot_bdsummary(ibudfl)
1013  end if
1014 
1015  ! -- model budget summary
1016  call this%budget%finalize_step(delt)
1017  if (ibudfl /= 0) then
1018  ipflag = 1
1019  call this%budget%budget_ot(kstp, kper, this%iout)
1020  end if
1021 
1022  ! -- Write to budget csv every time step
1023  call this%budget%writecsv(totim)
1024 
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
Here is the call graph for this function:

◆ gwf_ot_dv()

subroutine gwfmodule::gwf_ot_dv ( class(gwfmodeltype this,
integer(i4b), intent(in)  idvsave,
integer(i4b), intent(in)  idvprint,
integer(i4b), intent(inout)  ipflag 
)
private

Definition at line 963 of file gwf.f90.

964  class(GwfModelType) :: this
965  integer(I4B), intent(in) :: idvsave
966  integer(I4B), intent(in) :: idvprint
967  integer(I4B), intent(inout) :: ipflag
968  class(BndType), pointer :: packobj
969  integer(I4B) :: ip
970  !
971  ! -- Save compaction to binary file
972  if (this%incsub > 0) call this%csub%csub_ot_dv(idvsave, idvprint)
973  !
974  ! -- save density to binary file
975  if (this%inbuy > 0) then
976  call this%buy%buy_ot_dv(idvsave)
977  end if
978  !
979  ! -- save viscosity to binary file
980  if (this%invsc > 0) then
981  call this%vsc%vsc_ot_dv(idvsave)
982  end if
983  !
984  ! -- Print advanced package dependent variables
985  do ip = 1, this%bndlist%Count()
986  packobj => getbndfromlist(this%bndlist, ip)
987  call packobj%bnd_ot_dv(idvsave, idvprint)
988  end do
989  !
990  ! -- save head and print head
991  call this%oc%oc_ot(ipflag)
Here is the call graph for this function:

◆ gwf_ot_flow()

subroutine gwfmodule::gwf_ot_flow ( class(gwfmodeltype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  ibudfl,
integer(i4b), intent(in)  icbcun 
)
private

Definition at line 912 of file gwf.f90.

913  class(GwfModelType) :: this
914  integer(I4B), intent(in) :: icbcfl
915  integer(I4B), intent(in) :: ibudfl
916  integer(I4B), intent(in) :: icbcun
917  class(BndType), pointer :: packobj
918  integer(I4B) :: ip
919 
920  ! -- Save GWF flows
921  if (this%insto > 0) then
922  call this%sto%sto_save_model_flows(icbcfl, icbcun)
923  end if
924  if (this%innpf > 0) then
925  call this%npf%npf_save_model_flows(this%flowja, icbcfl, icbcun)
926  end if
927  if (this%incsub > 0) call this%csub%csub_save_model_flows(icbcfl, icbcun)
928  do ip = 1, this%bndlist%Count()
929  packobj => getbndfromlist(this%bndlist, ip)
930  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=0, icbcun=icbcun)
931  end do
932 
933  ! -- Save advanced package flows
934  do ip = 1, this%bndlist%Count()
935  packobj => getbndfromlist(this%bndlist, ip)
936  call packobj%bnd_ot_package_flows(icbcfl=icbcfl, ibudfl=0)
937  end do
938  if (this%inmvr > 0) then
939  call this%mvr%mvr_ot_saveflow(icbcfl, ibudfl)
940  end if
941 
942  ! -- Print GWF flows
943  if (this%innpf > 0) call this%npf%npf_print_model_flows(ibudfl, this%flowja)
944  if (this%ingnc > 0) call this%gnc%gnc_ot(ibudfl)
945  do ip = 1, this%bndlist%Count()
946  packobj => getbndfromlist(this%bndlist, ip)
947  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=ibudfl, icbcun=0)
948  end do
949 
950  ! -- Print advanced package flows
951  do ip = 1, this%bndlist%Count()
952  packobj => getbndfromlist(this%bndlist, ip)
953  call packobj%bnd_ot_package_flows(icbcfl=0, ibudfl=ibudfl)
954  end do
955  if (this%inmvr > 0) then
956  call this%mvr%mvr_ot_printflow(icbcfl, ibudfl)
957  end if
958 
Here is the call graph for this function:

◆ gwf_ot_obs()

subroutine gwfmodule::gwf_ot_obs ( class(gwfmodeltype this)

Definition at line 886 of file gwf.f90.

887  class(GwfModelType) :: this
888  class(BndType), pointer :: packobj
889  integer(I4B) :: ip
890 
891  ! -- Calculate and save GWF observations
892  call this%obs%obs_bd()
893  call this%obs%obs_ot()
894 
895  ! -- Calculate and save csub observations
896  if (this%incsub > 0) then
897  call this%csub%csub_bd_obs()
898  call this%csub%obs%obs_ot()
899  end if
900 
901  ! -- Calculate and save package observations
902  do ip = 1, this%bndlist%Count()
903  packobj => getbndfromlist(this%bndlist, ip)
904  call packobj%bnd_bd_obs()
905  call packobj%bnd_ot_obs()
906  end do
907 
Here is the call graph for this function:

◆ gwf_ptc()

subroutine gwfmodule::gwf_ptc ( class(gwfmodeltype this,
class(vectorbasetype), pointer  vec_residual,
integer(i4b), intent(inout)  iptc,
real(dp), intent(inout)  ptcf 
)
private

(1) Calculate maximum pseudo-transient continuation factor for the current outer iteration

Definition at line 626 of file gwf.f90.

627  ! -- modules
628  use constantsmodule, only: done
629  use tdismodule, only: delt
630  ! -- dummy
631  class(GwfModelType) :: this
632  class(VectorBaseType), pointer :: vec_residual
633  integer(I4B), intent(inout) :: iptc
634  real(DP), intent(inout) :: ptcf
635  ! -- local
636  integer(I4B) :: n
637  integer(I4B) :: iptct
638  real(DP) :: v
639  real(DP) :: resid
640  real(DP) :: ptcdelem1
641  !
642  ! -- set temporary flag indicating if pseudo-transient continuation should
643  ! be used for this model and time step
644  iptct = 0
645  ! -- only apply pseudo-transient continuation to problems using the
646  ! Newton-Raphson formulations for steady-state stress periods
647  if (this%iss > 0) then
648  if (this%inewton > 0) then
649  iptct = this%inewton
650  else
651  iptct = this%npf%inewton
652  end if
653  end if
654  !
655  ! -- calculate pseudo-transient continuation factor for model
656  if (iptct > 0) then
657  !
658  ! -- calculate the pseudo-time step using the residual
659  do n = 1, this%dis%nodes
660  if (this%npf%ibound(n) < 1) cycle
661  !
662  ! -- get the maximum volume of the cell (head at top of cell)
663  v = this%dis%get_cell_volume(n, this%dis%top(n))
664  !
665  ! -- set the residual
666  resid = vec_residual%get_value_local(n)
667  !
668  ! -- calculate the reciprocal of the pseudo-time step
669  ! resid [L3/T] / volume [L3] = [1/T]
670  ptcdelem1 = abs(resid) / v
671  !
672  ! -- set ptcf if the reciprocal of the pseudo-time step
673  ! exceeds the current value (equivalent to using the
674  ! smallest pseudo-time step)
675  if (ptcdelem1 > ptcf) ptcf = ptcdelem1
676  end do
677  !
678  ! -- protection for the case where the residuals are zero
679  if (ptcf == dzero) then
680  ptcf = done / (delt * dten)
681  end if
682  end if
683  !
684  ! -- reset ipc if needed
685  if (iptc == 0) then
686  if (iptct > 0) iptc = 1
687  end if

◆ gwf_ptcchk()

subroutine gwfmodule::gwf_ptcchk ( class(gwfmodeltype this,
integer(i4b), intent(inout)  iptc 
)
private

(1) Check if pseudo-transient continuation factor should be used

Definition at line 602 of file gwf.f90.

603  ! -- dummy
604  class(GwfModelType) :: this
605  integer(I4B), intent(inout) :: iptc
606  !
607  ! -- determine if pseudo-transient continuation should be applied to this
608  ! model - pseudo-transient continuation only applied to problems that
609  ! use the Newton-Raphson formulation during steady-state stress periods
610  iptc = 0
611  if (this%iss > 0) then
612  if (this%inewton > 0) then
613  iptc = this%inewton
614  else
615  iptc = this%npf%inewton
616  end if
617  end if

◆ gwf_rp()

subroutine gwfmodule::gwf_rp ( class(gwfmodeltype this)
private

(1) calls package read and prepare routines

Definition at line 360 of file gwf.f90.

361  ! -- modules
362  use tdismodule, only: readnewdata
363  ! -- dummy
364  class(GwfModelType) :: this
365  ! -- local
366  class(BndType), pointer :: packobj
367  integer(I4B) :: ip
368  !
369  ! -- Check with TDIS on whether or not it is time to RP
370  if (.not. readnewdata) return
371  !
372  ! -- Read and prepare
373  if (this%innpf > 0) call this%npf%npf_rp()
374  if (this%inbuy > 0) call this%buy%buy_rp()
375  if (this%invsc > 0) call this%vsc%vsc_rp()
376  if (this%inhfb > 0) call this%hfb%hfb_rp()
377  if (this%inoc > 0) call this%oc%oc_rp()
378  if (this%insto > 0) call this%sto%sto_rp()
379  if (this%incsub > 0) call this%csub%csub_rp()
380  if (this%inmvr > 0) call this%mvr%mvr_rp()
381  do ip = 1, this%bndlist%Count()
382  packobj => getbndfromlist(this%bndlist, ip)
383  call packobj%bnd_rp()
384  call packobj%bnd_rp_log()
385  call packobj%bnd_rp_obs()
386  end do
387  !
388  ! -- Check for steady state period
389  call this%steady_period_check()
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:26
Here is the call graph for this function:

◆ log_namfile_options()

subroutine gwfmodule::log_namfile_options ( class(gwfmodeltype this,
type(gwfnamparamfoundtype), intent(in)  found 
)

Definition at line 1542 of file gwf.f90.

1544  class(GwfModelType) :: this
1545  type(GwfNamParamFoundType), intent(in) :: found
1546 
1547  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1548 
1549  if (found%newton) then
1550  write (this%iout, '(4x,a)') &
1551  'NEWTON-RAPHSON method enabled for the model.'
1552  if (found%under_relaxation) then
1553  write (this%iout, '(4x,a,a)') &
1554  'NEWTON-RAPHSON UNDER-RELAXATION based on the bottom ', &
1555  'elevation of the model will be applied to the model.'
1556  end if
1557  end if
1558 
1559  if (found%print_input) then
1560  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1561  'FOR ALL MODEL STRESS PACKAGES'
1562  end if
1563 
1564  if (found%print_flows) then
1565  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1566  'FOR ALL MODEL PACKAGES'
1567  end if
1568 
1569  if (found%save_flows) then
1570  write (this%iout, '(4x,a)') &
1571  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1572  end if
1573 
1574  write (this%iout, '(1x,a)') 'END NAMEFILE OPTIONS:'

◆ package_create()

subroutine gwfmodule::package_create ( class(gwfmodeltype this,
character(len=*), intent(in)  filtyp,
integer(i4b), intent(in)  ipakid,
integer(i4b), intent(in)  ipaknum,
character(len=*), intent(in)  pakname,
character(len=*), intent(in)  mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

(1) create new-style package (2) add a pointer to the package

Definition at line 1213 of file gwf.f90.

1215  ! -- modules
1216  use constantsmodule, only: linelength
1217  use simmodule, only: store_error
1218  use chdmodule, only: chd_create
1219  use welmodule, only: wel_create
1220  use drnmodule, only: drn_create
1221  use rivmodule, only: riv_create
1222  use ghbmodule, only: ghb_create
1223  use rchmodule, only: rch_create
1224  use evtmodule, only: evt_create
1225  use mawmodule, only: maw_create
1226  use sfrmodule, only: sfr_create
1227  use lakmodule, only: lak_create
1228  use uzfmodule, only: uzf_create
1229  use apimodule, only: api_create
1230  ! -- dummy
1231  class(GwfModelType) :: this
1232  character(len=*), intent(in) :: filtyp
1233  integer(I4B), intent(in) :: ipakid
1234  integer(I4B), intent(in) :: ipaknum
1235  character(len=*), intent(in) :: pakname
1236  character(len=*), intent(in) :: mempath
1237  integer(I4B), intent(in) :: inunit
1238  integer(I4B), intent(in) :: iout
1239  ! -- local
1240  class(BndType), pointer :: packobj
1241  class(BndType), pointer :: packobj2
1242  integer(I4B) :: ip
1243  !
1244  ! -- This part creates the package object
1245  select case (filtyp)
1246  case ('CHD6')
1247  call chd_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1248  pakname, mempath)
1249  case ('WEL6')
1250  call wel_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1251  pakname, mempath)
1252  case ('DRN6')
1253  call drn_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1254  pakname, mempath)
1255  case ('RIV6')
1256  call riv_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1257  pakname, mempath)
1258  case ('GHB6')
1259  call ghb_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1260  pakname, mempath)
1261  case ('RCH6')
1262  call rch_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1263  pakname, mempath)
1264  case ('EVT6')
1265  call evt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1266  pakname, mempath)
1267  case ('MAW6')
1268  call maw_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1269  case ('SFR6')
1270  call sfr_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1271  case ('LAK6')
1272  call lak_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1273  case ('UZF6')
1274  call uzf_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1275  case ('API6')
1276  call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1277  pakname, mempath)
1278  case default
1279  write (errmsg, *) 'Invalid package type: ', filtyp
1280  call store_error(errmsg, terminate=.true.)
1281  end select
1282  !
1283  ! -- Check to make sure that the package name is unique, then store a
1284  ! pointer to the package in the model bndlist
1285  do ip = 1, this%bndlist%Count()
1286  packobj2 => getbndfromlist(this%bndlist, ip)
1287  if (packobj2%packName == pakname) then
1288  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
1289  'already exists: ', trim(pakname)
1290  call store_error(errmsg, terminate=.true.)
1291  end if
1292  end do
1293  call addbndtolist(this%bndlist, packobj)
This module contains the API package methods.
Definition: gwf-api.f90:12
subroutine, public api_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
@ brief Create a new package object
Definition: gwf-api.f90:51
subroutine, public chd_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a new constant head package.
Definition: gwf-chd.f90:56
subroutine, public drn_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Drn Package and point packobj to the new package.
Definition: gwf-drn.f90:60
subroutine, public evt_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a new Evapotranspiration Segments Package and point pakobj to the new package.
Definition: gwf-evt.f90:65
subroutine, public ghb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Ghb Package and point bndobj to the new package.
Definition: gwf-ghb.f90:48
subroutine, public lak_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname)
Create a new LAK Package and point bndobj to the new package.
Definition: gwf-lak.f90:291
subroutine, public maw_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname)
Create a New Multi-Aquifer Well (MAW) Package.
Definition: gwf-maw.f90:234
subroutine, public rch_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Recharge Package.
Definition: gwf-rch.f90:56
subroutine, public riv_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Riv Package and point packobj to the new package.
Definition: gwf-riv.f90:51
This module contains the SFR package methods.
Definition: gwf-sfr.f90:7
subroutine, public sfr_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname)
@ brief Create a new package object
Definition: gwf-sfr.f90:295
subroutine, public uzf_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname)
Create a New UZF Package and point packobj to the new package.
Definition: gwf-uzf.f90:176
This module contains the WEL package methods.
Definition: gwf-wel.f90:15
subroutine, public wel_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
@ brief Create a new package object
Definition: gwf-wel.f90:75
Here is the call graph for this function:

◆ steady_period_check()

subroutine gwfmodule::steady_period_check ( class(gwfmodeltype this)

Write warning message if steady state period and adaptive time stepping is active for the period

Definition at line 1584 of file gwf.f90.

1585  ! -- modules
1586  use tdismodule, only: kper
1588  use simvariablesmodule, only: warnmsg
1589  use simmodule, only: store_warning
1590  ! -- dummy
1591  class(GwfModelType) :: this
1592  if (this%iss == 1) then
1593  if (isadaptiveperiod(kper)) then
1594  write (warnmsg, '(a,a,a,i0,a)') &
1595  'GWF Model (', trim(this%name), ') is steady state for period ', &
1596  kper, ' and adaptive time stepping is active. Adaptive time &
1597  &stepping may not work properly for steady-state conditions.'
1598  call store_warning(warnmsg)
1599  end if
1600  end if
logical(lgp) function, public isadaptiveperiod(kper)
@ brief Determine if period is adaptive
Definition: ats.f90:45
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
character(len=maxcharlen) warnmsg
warning message string
Here is the call graph for this function:

Variable Documentation

◆ gwf_basepkg

character(len=lenpackagetype), dimension(gwf_nbasepkg), public gwfmodule::gwf_basepkg

Definition at line 107 of file gwf.f90.

107  character(len=LENPACKAGETYPE), dimension(GWF_NBASEPKG) :: GWF_BASEPKG

◆ gwf_multipkg

character(len=lenpackagetype), dimension(gwf_nmultipkg), public gwfmodule::gwf_multipkg

Definition at line 120 of file gwf.f90.

120  character(len=LENPACKAGETYPE), dimension(GWF_NMULTIPKG) :: GWF_MULTIPKG

◆ gwf_nbasepkg

integer(i4b), parameter, public gwfmodule::gwf_nbasepkg = 50

GWF6 model base package types. Only listed packages are candidates for input and these will be loaded in the order specified.

Definition at line 106 of file gwf.f90.

106  integer(I4B), parameter :: GWF_NBASEPKG = 50

◆ gwf_nmultipkg

integer(i4b), parameter, public gwfmodule::gwf_nmultipkg = 50

GWF6 model multi-instance package types. Only listed packages are candidates for input and these will be loaded in the order specified.

Definition at line 119 of file gwf.f90.

119  integer(I4B), parameter :: GWF_NMULTIPKG = 50

◆ niunit_gwf

integer(i4b), parameter gwfmodule::niunit_gwf = GWF_NBASEPKG + GWF_NMULTIPKG
private

Definition at line 127 of file gwf.f90.

127  integer(I4B), parameter :: NIUNIT_GWF = gwf_nbasepkg + gwf_nmultipkg