MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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)
 
subroutine gwf_ot_flow (this, icbcfl, ibudfl, icbcun)
 
subroutine gwf_ot_dv (this, idvsave, idvprint, ipflag)
 
subroutine gwf_ot_bdsummary (this, ibudfl, ipflag)
 
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 1217 of file gwf.f90.

1218  ! -- modules
1220  ! -- dummy
1221  class(GwfModelType) :: this
1222  character(len=*), intent(in) :: modelname
1223  !
1224  ! -- allocate members from parent class
1225  call this%NumericalModelType%allocate_scalars(modelname)
1226  !
1227  ! -- allocate members that are part of model class
1228  call mem_allocate(this%inic, 'INIC', this%memoryPath)
1229  call mem_allocate(this%inoc, 'INOC', this%memoryPath)
1230  call mem_allocate(this%innpf, 'INNPF', this%memoryPath)
1231  call mem_allocate(this%inbuy, 'INBUY', this%memoryPath)
1232  call mem_allocate(this%invsc, 'INVSC', this%memoryPath)
1233  call mem_allocate(this%insto, 'INSTO', this%memoryPath)
1234  call mem_allocate(this%incsub, 'INCSUB', this%memoryPath)
1235  call mem_allocate(this%inmvr, 'INMVR', this%memoryPath)
1236  call mem_allocate(this%inhfb, 'INHFB', this%memoryPath)
1237  call mem_allocate(this%ingnc, 'INGNC', this%memoryPath)
1238  call mem_allocate(this%inobs, 'INOBS', this%memoryPath)
1239  call mem_allocate(this%iss, 'ISS', this%memoryPath)
1240  call mem_allocate(this%inewtonur, 'INEWTONUR', this%memoryPath)
1241  !
1242  this%inic = 0
1243  this%inoc = 0
1244  this%innpf = 0
1245  this%inbuy = 0
1246  this%invsc = 0
1247  this%insto = 0
1248  this%incsub = 0
1249  this%inmvr = 0
1250  this%inhfb = 0
1251  this%ingnc = 0
1252  this%inobs = 0
1253  this%iss = 1 !default is steady-state (i.e., no STO package)
1254  this%inewtonur = 0 !default is to not use newton bottom head dampening
1255  !
1256  ! -- return
1257  return

◆ castasgwfmodel()

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

Definition at line 1390 of file gwf.f90.

1391  implicit none
1392  class(*), pointer, intent(inout) :: model
1393  class(GwfModelType), pointer :: gwfModel
1394 
1395  gwfmodel => null()
1396  if (.not. associated(model)) return
1397  select type (model)
1398  class is (gwfmodeltype)
1399  gwfmodel => model
1400  end select
1401  return
1402 
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 1407 of file gwf.f90.

1409  ! -- modules
1412  ! -- dummy
1413  class(GwfModelType) :: this
1414  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
1415  type(CharacterStringType), dimension(:), contiguous, &
1416  pointer, intent(inout) :: pkgtypes
1417  type(CharacterStringType), dimension(:), contiguous, &
1418  pointer, intent(inout) :: pkgnames
1419  type(CharacterStringType), dimension(:), contiguous, &
1420  pointer, intent(inout) :: mempaths
1421  integer(I4B), dimension(:), contiguous, &
1422  pointer, intent(inout) :: inunits
1423  ! -- local
1424  integer(I4B) :: ipakid, ipaknum
1425  character(len=LENFTYPE) :: pkgtype, bndptype
1426  character(len=LENPACKAGENAME) :: pkgname
1427  character(len=LENMEMPATH) :: mempath
1428  integer(I4B), pointer :: inunit
1429  integer(I4B) :: n
1430 
1431  if (allocated(bndpkgs)) then
1432  !
1433  ! -- create stress packages
1434  ipakid = 1
1435  bndptype = ''
1436  do n = 1, size(bndpkgs)
1437  !
1438  pkgtype = pkgtypes(bndpkgs(n))
1439  pkgname = pkgnames(bndpkgs(n))
1440  mempath = mempaths(bndpkgs(n))
1441  inunit => inunits(bndpkgs(n))
1442  !
1443  if (bndptype /= pkgtype) then
1444  ipaknum = 1
1445  bndptype = pkgtype
1446  end if
1447  !
1448  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1449  inunit, this%iout)
1450  ipakid = ipakid + 1
1451  ipaknum = ipaknum + 1
1452  end do
1453  !
1454  ! -- cleanup
1455  deallocate (bndpkgs)
1456  end if
1457  !
1458  ! -- return
1459  return
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:22
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 1464 of file gwf.f90.

1465  ! -- modules
1468  use arrayhandlersmodule, only: expandarray
1469  use memorymanagermodule, only: mem_setptr
1471  use simvariablesmodule, only: idm_context
1472  use dismodule, only: dis_cr
1473  use disvmodule, only: disv_cr
1474  use disumodule, only: disu_cr
1475  use gwfnpfmodule, only: npf_cr
1476  use xt3dmodule, only: xt3d_cr
1477  use gwfbuymodule, only: buy_cr
1478  use gwfvscmodule, only: vsc_cr
1479  use gwfstomodule, only: sto_cr
1480  use gwfcsubmodule, only: csub_cr
1481  use gwfmvrmodule, only: mvr_cr
1482  use gwfhfbmodule, only: hfb_cr
1483  use gwficmodule, only: ic_cr
1484  use gwfocmodule, only: oc_cr
1485  ! -- dummy
1486  class(GwfModelType) :: this
1487  ! -- local
1488  type(CharacterStringType), dimension(:), contiguous, &
1489  pointer :: pkgtypes => null()
1490  type(CharacterStringType), dimension(:), contiguous, &
1491  pointer :: pkgnames => null()
1492  type(CharacterStringType), dimension(:), contiguous, &
1493  pointer :: mempaths => null()
1494  integer(I4B), dimension(:), contiguous, &
1495  pointer :: inunits => null()
1496  character(len=LENMEMPATH) :: model_mempath
1497  character(len=LENFTYPE) :: pkgtype
1498  character(len=LENPACKAGENAME) :: pkgname
1499  character(len=LENMEMPATH) :: mempath
1500  integer(I4B), pointer :: inunit
1501  integer(I4B), dimension(:), allocatable :: bndpkgs
1502  integer(I4B) :: n
1503  integer(I4B) :: indis = 0 ! DIS enabled flag
1504  character(len=LENMEMPATH) :: mempathnpf = ''
1505  character(len=LENMEMPATH) :: mempathic = ''
1506  !
1507  ! -- set input model memory path
1508  model_mempath = create_mem_path(component=this%name, context=idm_context)
1509  !
1510  ! -- set pointers to model path package info
1511  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1512  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1513  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1514  call mem_setptr(inunits, 'INUNITS', model_mempath)
1515  !
1516  do n = 1, size(pkgtypes)
1517  !
1518  ! attributes for this input package
1519  pkgtype = pkgtypes(n)
1520  pkgname = pkgnames(n)
1521  mempath = mempaths(n)
1522  inunit => inunits(n)
1523  !
1524  ! -- create dis package as it is a prerequisite for other packages
1525  select case (pkgtype)
1526  case ('DIS6')
1527  indis = 1
1528  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1529  case ('DISV6')
1530  indis = 1
1531  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1532  case ('DISU6')
1533  indis = 1
1534  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1535  case ('NPF6')
1536  this%innpf = 1
1537  mempathnpf = mempath
1538  case ('BUY6')
1539  this%inbuy = inunit
1540  case ('VSC6')
1541  this%invsc = inunit
1542  case ('GNC6')
1543  this%ingnc = inunit
1544  case ('HFB6')
1545  this%inhfb = inunit
1546  case ('STO6')
1547  this%insto = inunit
1548  case ('CSUB6')
1549  this%incsub = inunit
1550  case ('IC6')
1551  this%inic = 1
1552  mempathic = mempath
1553  case ('MVR6')
1554  this%inmvr = inunit
1555  case ('OC6')
1556  this%inoc = inunit
1557  case ('OBS6')
1558  this%inobs = inunit
1559  case ('WEL6', 'DRN6', 'RIV6', 'GHB6', 'RCH6', &
1560  'EVT6', 'API6', 'CHD6', 'MAW6', 'SFR6', &
1561  'LAK6', 'UZF6')
1562  call expandarray(bndpkgs)
1563  bndpkgs(size(bndpkgs)) = n
1564  case default
1565  ! TODO
1566  end select
1567  end do
1568  !
1569  ! -- Create packages that are tied directly to model
1570  call npf_cr(this%npf, this%name, mempathnpf, this%innpf, this%iout)
1571  call xt3d_cr(this%xt3d, this%name, this%innpf, this%iout)
1572  call buy_cr(this%buy, this%name, this%inbuy, this%iout)
1573  call vsc_cr(this%vsc, this%name, this%invsc, this%iout)
1574  call gnc_cr(this%gnc, this%name, this%ingnc, this%iout)
1575  call hfb_cr(this%hfb, this%name, this%inhfb, this%iout)
1576  call sto_cr(this%sto, this%name, this%insto, this%iout)
1577  call csub_cr(this%csub, this%name, this%insto, this%sto%packName, &
1578  this%incsub, this%iout)
1579  call ic_cr(this%ic, this%name, mempathic, this%inic, this%iout, this%dis)
1580  call mvr_cr(this%mvr, this%name, this%inmvr, this%iout, this%dis)
1581  call oc_cr(this%oc, this%name, this%inoc, this%iout)
1582  call gwf_obs_cr(this%obs, this%inobs)
1583  !
1584  ! -- Check to make sure that required ftype's have been specified
1585  call this%ftype_check(indis)
1586  !
1587  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
1588  !
1589  ! -- return
1590  return
Definition: Dis.f90:1
subroutine, public dis_cr(dis, name_model, input_mempath, inunit, iout)
Create a new structured discretization object.
Definition: Dis.f90:97
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:109
subroutine, public buy_cr(buyobj, name_model, inunit, iout)
Create a new BUY object.
Definition: gwf-buy.f90:106
This module contains the CSUB package methods.
Definition: gwf-csub.f90:9
subroutine, public csub_cr(csubobj, name_model, istounit, stoPckName, inunit, iout)
@ brief Create a new package object
Definition: gwf-csub.f90:322
subroutine, public hfb_cr(hfbobj, name_model, inunit, iout)
Create a new hfb object.
Definition: gwf-hfb.f90:70
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:185
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:153
subroutine, public oc_cr(ocobj, name_model, 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, inunit, iout)
@ brief Create a new package object
Definition: gwf-sto.f90:76
subroutine, public vsc_cr(vscobj, name_model, 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 1352 of file gwf.f90.

1353  ! -- modules
1354  use constantsmodule, only: linelength
1355  use simmodule, only: store_error, count_errors
1356  ! -- dummy
1357  class(GwfModelType) :: this
1358  integer(I4B), intent(in) :: indis
1359  ! -- local
1360  !
1361  ! -- Check for IC8, DIS(u), and NPF. Stop if not present.
1362  if (this%inic == 0) then
1363  write (errmsg, '(a)') &
1364  'Initial Conditions (IC6) package not specified.'
1365  call store_error(errmsg)
1366  end if
1367  if (indis == 0) then
1368  write (errmsg, '(a)') &
1369  'Discretization (DIS6, DISV6, or DISU6) Package not specified.'
1370  call store_error(errmsg)
1371  end if
1372  if (this%innpf == 0) then
1373  write (errmsg, '(a)') &
1374  'Node Property Flow (NPF6) Package not specified.'
1375  call store_error(errmsg)
1376  end if
1377  !
1378  if (count_errors() > 0) then
1379  write (errmsg, '(a)') 'One or more required package(s) not specified.'
1380  call store_error(errmsg)
1381  call store_error_filename(this%filename)
1382  end if
1383  !
1384  ! -- return
1385  return
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 260 of file gwf.f90.

261  ! -- modules
262  use sparsemodule, only: sparsematrix
263  ! -- dummy
264  class(GwfModelType) :: this
265  type(sparsematrix), intent(inout) :: sparse
266  ! -- local
267  class(BndType), pointer :: packobj
268  integer(I4B) :: ip
269  !
270  ! -- Add the primary grid connections of this model to sparse
271  call this%dis%dis_ac(this%moffset, sparse)
272  !
273  ! -- Add any additional connections that NPF may need
274  if (this%innpf > 0) call this%npf%npf_ac(this%moffset, sparse)
275  !
276  ! -- Add any package connections
277  do ip = 1, this%bndlist%Count()
278  packobj => getbndfromlist(this%bndlist, ip)
279  call packobj%bnd_ac(this%moffset, sparse)
280  end do
281  !
282  ! -- If GNC is active, then add the gnc connections to sparse
283  if (this%ingnc > 0) call this%gnc%gnc_ac(sparse)
284  !
285  ! -- return
286  return
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 413 of file gwf.f90.

414  ! -- modules
416  ! -- dummy
417  class(GwfModelType) :: this
418  class(BndType), pointer :: packobj
419  ! -- local
420  integer(I4B) :: irestore
421  integer(I4B) :: ip, n
422  !
423  ! -- Reset state variable
424  irestore = 0
425  if (ifailedstepretry > 0) irestore = 1
426  if (irestore == 0) then
427  !
428  ! -- copy x into xold
429  do n = 1, this%dis%nodes
430  this%xold(n) = this%x(n)
431  end do
432  else
433  !
434  ! -- copy xold into x if this time step is a redo
435  do n = 1, this%dis%nodes
436  this%x(n) = this%xold(n)
437  end do
438  end if
439  !
440  ! -- Advance
441  if (this%invsc > 0) call this%vsc%vsc_ad()
442  if (this%innpf > 0) call this%npf%npf_ad(this%dis%nodes, this%xold, &
443  this%x, irestore)
444  if (this%insto > 0) call this%sto%sto_ad()
445  if (this%incsub > 0) call this%csub%csub_ad(this%dis%nodes, this%x)
446  if (this%inbuy > 0) call this%buy%buy_ad()
447  if (this%inmvr > 0) call this%mvr%mvr_ad()
448  do ip = 1, this%bndlist%Count()
449  packobj => getbndfromlist(this%bndlist, ip)
450  call packobj%bnd_ad()
451  if (this%invsc > 0) call this%vsc%vsc_ad_bnd(packobj, this%x)
452  if (isimcheck > 0) then
453  call packobj%bnd_ck()
454  end if
455  end do
456  !
457  ! -- Push simulated values to preceding time/subtime step
458  call this%obs%obs_ad()
459  !
460  ! -- return
461  return
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 327 of file gwf.f90.

328  ! -- dummy
329  class(GwfModelType) :: this
330  ! -- locals
331  integer(I4B) :: ip
332  class(BndType), pointer :: packobj
333  !
334  ! -- Allocate and read modules attached to model
335  if (this%inic > 0) call this%ic%ic_ar(this%x)
336  if (this%innpf > 0) call this%npf%npf_ar(this%ic, this%vsc, this%ibound, &
337  this%x)
338  if (this%invsc > 0) call this%vsc%vsc_ar(this%ibound)
339  if (this%inbuy > 0) call this%buy%buy_ar(this%npf, this%ibound)
340  if (this%inhfb > 0) call this%hfb%hfb_ar(this%ibound, this%xt3d, this%dis, &
341  this%invsc, this%vsc)
342  if (this%insto > 0) call this%sto%sto_ar(this%dis, this%ibound)
343  if (this%incsub > 0) call this%csub%csub_ar(this%dis, this%ibound)
344  if (this%inmvr > 0) call this%mvr%mvr_ar()
345  if (this%inobs > 0) call this%obs%gwf_obs_ar(this%ic, this%x, this%flowja)
346  !
347  ! -- Call dis_ar to write binary grid file
348  call this%dis%dis_ar(this%npf%icelltype)
349  !
350  ! -- set up output control
351  call this%oc%oc_ar(this%x, this%dis, this%npf%hnoflo)
352  call this%budget%set_ibudcsv(this%oc%ibudcsv)
353  !
354  ! -- Package input files now open, so allocate and read
355  do ip = 1, this%bndlist%Count()
356  packobj => getbndfromlist(this%bndlist, ip)
357  call packobj%set_pointers(this%dis%nodes, this%ibound, this%x, &
358  this%xold, this%flowja)
359  ! -- Read and allocate package
360  call packobj%bnd_ar()
361  if (this%inbuy > 0) call this%buy%buy_ar_bnd(packobj, this%x)
362  if (this%invsc > 0) call this%vsc%vsc_ar_bnd(packobj)
363  end do
364  !
365  ! -- return
366  return
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 824 of file gwf.f90.

825  ! -- modules
826  use sparsemodule, only: csr_diagsum
827  ! -- dummy
828  class(GwfModelType) :: this
829  integer(I4B), intent(in) :: icnvg
830  integer(I4B), intent(in) :: isuppress_output
831  ! -- local
832  integer(I4B) :: ip
833  class(BndType), pointer :: packobj
834  !
835  ! -- Finalize calculation of flowja by adding face flows to the diagonal.
836  ! This results in the flow residual being stored in the diagonal
837  ! position for each cell.
838  call csr_diagsum(this%dis%con%ia, this%flowja)
839  !
840  ! -- Save the solution convergence flag
841  this%icnvg = icnvg
842  !
843  ! -- Budget routines (start by resetting). Sole purpose of this section
844  ! is to add in and outs to model budget. All ins and out for a model
845  ! should be added here to this%budget. In a subsequent exchange call,
846  ! exchange flows might also be added.
847  call this%budget%reset()
848  if (this%insto > 0) call this%sto%sto_bd(isuppress_output, this%budget)
849  if (this%incsub > 0) call this%csub%csub_bd(isuppress_output, this%budget)
850  if (this%inmvr > 0) call this%mvr%mvr_bd()
851  do ip = 1, this%bndlist%Count()
852  packobj => getbndfromlist(this%bndlist, ip)
853  call packobj%bnd_bd(this%budget)
854  end do
855  !
856  ! -- npf velocities have to be calculated here, after gwf-gwf exchanges
857  ! have passed in their contributions from exg_cq()
858  if (this%innpf > 0) then
859  if (this%npf%icalcspdis /= 0) then
860  call this%npf%calc_spdis(this%flowja)
861  end if
862  end if
863  !
864  ! -- Return
865  return
subroutine csr_diagsum(ia, flowja)
Definition: Sparse.f90:281
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 contributions.

(1) adds the entry to the budget object

Definition at line 1166 of file gwf.f90.

1167  ! -- modules
1168  use constantsmodule, only: lenbudtxt
1169  use tdismodule, only: delt
1170  ! -- dummy
1171  class(GwfModelType) :: this
1172  real(DP), dimension(:, :), intent(in) :: budterm
1173  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
1174  character(len=*), intent(in) :: rowlabel
1175  !
1176  call this%budget%addentry(budterm, delt, budtxt, rowlabel=rowlabel)
1177  !
1178  ! -- return
1179  return
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:36
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 586 of file gwf.f90.

587  ! -- dummy
588  class(GwfModelType) :: this
589  integer(I4B), intent(in) :: innertot
590  integer(I4B), intent(in) :: kiter
591  integer(I4B), intent(in) :: iend
592  integer(I4B), intent(in) :: icnvgmod
593  character(len=LENPAKLOC), intent(inout) :: cpak
594  integer(I4B), intent(inout) :: ipak
595  real(DP), intent(inout) :: dpak
596  ! -- local
597  class(BndType), pointer :: packobj
598  integer(I4B) :: ip
599  ! -- formats
600  !
601  ! -- If mover is on, then at least 2 outers required
602  if (this%inmvr > 0) then
603  call this%mvr%mvr_cc(innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
604  end if
605  !
606  ! -- csub convergence check
607  if (this%incsub > 0) then
608  call this%csub%csub_cc(innertot, kiter, iend, icnvgmod, &
609  this%dis%nodes, this%x, this%xold, &
610  cpak, ipak, dpak)
611  end if
612  !
613  ! -- Call package cc routines
614  do ip = 1, this%bndlist%Count()
615  packobj => getbndfromlist(this%bndlist, ip)
616  call packobj%bnd_cc(innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
617  end do
618  !
619  ! -- return
620  return
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 465 of file gwf.f90.

466  ! -- dummy
467  class(GwfModelType) :: this
468  integer(I4B), intent(in) :: kiter
469  ! -- local
470  class(BndType), pointer :: packobj
471  integer(I4B) :: ip
472  !
473  ! -- Call package cf routines
474  if (this%innpf > 0) call this%npf%npf_cf(kiter, this%dis%nodes, this%x)
475  if (this%inbuy > 0) call this%buy%buy_cf(kiter)
476  do ip = 1, this%bndlist%Count()
477  packobj => getbndfromlist(this%bndlist, ip)
478  call packobj%bnd_cf()
479  if (this%inbuy > 0) call this%buy%buy_cf_bnd(packobj, this%x)
480  end do
481  !
482  ! -- return
483  return
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 777 of file gwf.f90.

778  ! -- modules
779  ! -- dummy
780  class(GwfModelType) :: this
781  integer(I4B), intent(in) :: icnvg
782  integer(I4B), intent(in) :: isuppress_output
783  ! -- local
784  integer(I4B) :: i
785  integer(I4B) :: ip
786  class(BndType), pointer :: packobj
787  !
788  ! -- Construct the flowja array. Flowja is calculated each time, even if
789  ! output is suppressed. (flowja is positive into a cell.) The diagonal
790  ! position of the flowja array will contain the flow residual after
791  ! these routines are called, so each package is responsible for adding
792  ! its flow to this diagonal position.
793  do i = 1, this%nja
794  this%flowja(i) = dzero
795  end do
796  if (this%innpf > 0) call this%npf%npf_cq(this%x, this%flowja)
797  if (this%inbuy > 0) call this%buy%buy_cq(this%x, this%flowja)
798  if (this%inhfb > 0) call this%hfb%hfb_cq(this%x, this%flowja)
799  if (this%ingnc > 0) call this%gnc%gnc_cq(this%flowja)
800  if (this%insto > 0) call this%sto%sto_cq(this%flowja, this%x, this%xold)
801  if (this%incsub > 0) call this%csub%csub_cq(this%dis%nodes, this%x, &
802  this%xold, isuppress_output, &
803  this%flowja)
804  !
805  ! -- Go through packages and call cq routines. cf() routines are called
806  ! first to regenerate non-linear terms to be consistent with the final
807  ! head solution.
808  do ip = 1, this%bndlist%Count()
809  packobj => getbndfromlist(this%bndlist, ip)
810  call packobj%bnd_cf()
811  if (this%inbuy > 0) call this%buy%buy_cf_bnd(packobj, this%x)
812  call packobj%bnd_cq(this%x, this%flowja)
813  end do
814  !
815  ! -- Return
816  return
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()
208  !
209  ! -- return
210  return
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:161
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 1082 of file gwf.f90.

1083  ! -- modules
1086  use simvariablesmodule, only: idm_context
1087  ! -- dummy
1088  class(GwfModelType) :: this
1089  ! -- local
1090  integer(I4B) :: ip
1091  class(BndType), pointer :: packobj
1092  !
1093  ! -- Deallocate idm memory
1094  call memorylist_remove(this%name, 'NAM', idm_context)
1095  call memorylist_remove(component=this%name, context=idm_context)
1096  !
1097  ! -- Internal flow packages deallocate
1098  call this%dis%dis_da()
1099  call this%ic%ic_da()
1100  call this%npf%npf_da()
1101  call this%xt3d%xt3d_da()
1102  call this%buy%buy_da()
1103  call this%vsc%vsc_da()
1104  call this%gnc%gnc_da()
1105  call this%sto%sto_da()
1106  call this%csub%csub_da()
1107  call this%budget%budget_da()
1108  call this%hfb%hfb_da()
1109  call this%mvr%mvr_da()
1110  call this%oc%oc_da()
1111  call this%obs%obs_da()
1112  !
1113  ! -- Internal package objects
1114  deallocate (this%dis)
1115  deallocate (this%ic)
1116  deallocate (this%npf)
1117  deallocate (this%xt3d)
1118  deallocate (this%buy)
1119  deallocate (this%vsc)
1120  deallocate (this%gnc)
1121  deallocate (this%sto)
1122  deallocate (this%csub)
1123  deallocate (this%budget)
1124  deallocate (this%hfb)
1125  deallocate (this%mvr)
1126  deallocate (this%obs)
1127  deallocate (this%oc)
1128  !
1129  ! -- Boundary packages
1130  do ip = 1, this%bndlist%Count()
1131  packobj => getbndfromlist(this%bndlist, ip)
1132  call packobj%bnd_da()
1133  deallocate (packobj)
1134  end do
1135  !
1136  ! -- Scalars
1137  call mem_deallocate(this%inic)
1138  call mem_deallocate(this%inoc)
1139  call mem_deallocate(this%inobs)
1140  call mem_deallocate(this%innpf)
1141  call mem_deallocate(this%inbuy)
1142  call mem_deallocate(this%invsc)
1143  call mem_deallocate(this%insto)
1144  call mem_deallocate(this%incsub)
1145  call mem_deallocate(this%inmvr)
1146  call mem_deallocate(this%inhfb)
1147  call mem_deallocate(this%ingnc)
1148  call mem_deallocate(this%iss)
1149  call mem_deallocate(this%inewtonur)
1150  !
1151  ! -- NumericalModelType
1152  call this%NumericalModelType%model_da()
1153  !
1154  ! -- return
1155  return
subroutine, public memorylist_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 219 of file gwf.f90.

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

488  ! -- dummy
489  class(GwfModelType) :: this
490  integer(I4B), intent(in) :: kiter
491  class(MatrixBaseType), pointer :: matrix_sln
492  integer(I4B), intent(in) :: inwtflag
493  ! -- local
494  class(BndType), pointer :: packobj
495  integer(I4B) :: ip
496  integer(I4B) :: inwt, inwtsto, inwtcsub, inwtpak
497  !
498  ! -- newton flags
499  inwt = inwtflag
500  if (inwtflag == 1) inwt = this%npf%inewton
501  inwtsto = inwtflag
502  if (this%insto > 0) then
503  if (inwtflag == 1) inwtsto = this%sto%inewton
504  end if
505  inwtcsub = inwtflag
506  if (this%incsub > 0) then
507  if (inwtflag == 1) inwtcsub = this%csub%inewton
508  end if
509  !
510  ! -- Fill standard conductance terms
511  if (this%innpf > 0) call this%npf%npf_fc(kiter, matrix_sln, this%idxglo, &
512  this%rhs, this%x)
513  if (this%inbuy > 0) call this%buy%buy_fc(kiter, matrix_sln, this%idxglo, &
514  this%rhs, this%x)
515  if (this%inhfb > 0) call this%hfb%hfb_fc(kiter, matrix_sln, this%idxglo, &
516  this%rhs, this%x)
517  if (this%ingnc > 0) call this%gnc%gnc_fc(kiter, matrix_sln)
518  ! -- storage
519  if (this%insto > 0) then
520  call this%sto%sto_fc(kiter, this%xold, this%x, matrix_sln, &
521  this%idxglo, this%rhs)
522  end if
523  ! -- skeletal storage, compaction, and land subsidence
524  if (this%incsub > 0) then
525  call this%csub%csub_fc(kiter, this%xold, this%x, matrix_sln, &
526  this%idxglo, this%rhs)
527  end if
528  if (this%inmvr > 0) call this%mvr%mvr_fc()
529  do ip = 1, this%bndlist%Count()
530  packobj => getbndfromlist(this%bndlist, ip)
531  call packobj%bnd_fc(this%rhs, this%ia, this%idxglo, matrix_sln)
532  end do
533  !
534  !--Fill newton terms
535  if (this%innpf > 0) then
536  if (inwt /= 0) then
537  call this%npf%npf_fn(kiter, matrix_sln, this%idxglo, this%rhs, this%x)
538  end if
539  end if
540  !
541  ! -- Fill newton terms for ghost nodes
542  if (this%ingnc > 0) then
543  if (inwt /= 0) then
544  call this%gnc%gnc_fn(kiter, matrix_sln, this%npf%condsat, &
545  ivarcv_opt=this%npf%ivarcv, &
546  ictm1_opt=this%npf%icelltype, &
547  ictm2_opt=this%npf%icelltype)
548  end if
549  end if
550  !
551  ! -- Fill newton terms for storage
552  if (this%insto > 0) then
553  if (inwtsto /= 0) then
554  call this%sto%sto_fn(kiter, this%xold, this%x, matrix_sln, &
555  this%idxglo, this%rhs)
556  end if
557  end if
558  !
559  ! -- Fill newton terms for skeletal storage, compaction, and land subsidence
560  if (this%incsub > 0) then
561  if (inwtcsub /= 0) then
562  call this%csub%csub_fn(kiter, this%xold, this%x, matrix_sln, &
563  this%idxglo, this%rhs)
564  end if
565  end if
566  !
567  ! -- Fill Newton terms for packages
568  do ip = 1, this%bndlist%Count()
569  packobj => getbndfromlist(this%bndlist, ip)
570  inwtpak = inwtflag
571  if (inwtflag == 1) inwtpak = packobj%inewton
572  if (inwtpak /= 0) then
573  call packobj%bnd_fn(this%rhs, this%ia, this%idxglo, matrix_sln)
574  end if
575  end do
576  !
577  ! -- return
578  return
Here is the call graph for this function:

◆ gwf_fp()

subroutine gwfmodule::gwf_fp ( class(gwfmodeltype this)

Definition at line 1067 of file gwf.f90.

1068  ! -- modules
1069  ! -- dummy
1070  class(GwfModelType) :: this
1071  ! -- local
1072  !
1073  ! -- csub final processing
1074  if (this%incsub > 0) then
1075  call this%csub%csub_fp()
1076  end if
1077  !
1078  return

◆ gwf_get_iasym()

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

Definition at line 1185 of file gwf.f90.

1186  class(GwfModelType) :: this
1187  ! -- local
1188  integer(I4B) :: iasym
1189  integer(I4B) :: ip
1190  class(BndType), pointer :: packobj
1191  !
1192  ! -- Start by setting iasym to zero
1193  iasym = 0
1194  !
1195  ! -- NPF
1196  if (this%innpf > 0) then
1197  if (this%npf%iasym /= 0) iasym = 1
1198  if (this%npf%ixt3d /= 0) iasym = 1
1199  end if
1200  !
1201  ! -- GNC
1202  if (this%ingnc > 0) then
1203  if (this%gnc%iasym /= 0) iasym = 1
1204  end if
1205  !
1206  ! -- Check for any packages that introduce matrix asymmetry
1207  do ip = 1, this%bndlist%Count()
1208  packobj => getbndfromlist(this%bndlist, ip)
1209  if (packobj%iasym /= 0) iasym = 1
1210  end do
1211  !
1212  ! -- return
1213  return
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 292 of file gwf.f90.

293  ! -- dummy
294  class(GwfModelType) :: this
295  class(MatrixBaseType), pointer :: matrix_sln
296  ! -- local
297  class(BndType), pointer :: packobj
298  integer(I4B) :: ip
299  !
300  ! -- Find the position of each connection in the global ia, ja structure
301  ! and store them in idxglo.
302  call this%dis%dis_mc(this%moffset, this%idxglo, matrix_sln)
303  !
304  ! -- Map any additional connections that NPF may need
305  if (this%innpf > 0) call this%npf%npf_mc(this%moffset, matrix_sln)
306  !
307  ! -- Map any package connections
308  do ip = 1, this%bndlist%Count()
309  packobj => getbndfromlist(this%bndlist, ip)
310  call packobj%bnd_mc(this%moffset, matrix_sln)
311  end do
312  !
313  ! -- For implicit gnc, need to store positions of gnc connections
314  ! in solution matrix connection
315  if (this%ingnc > 0) call this%gnc%gnc_mc(matrix_sln)
316  !
317  ! -- return
318  return
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 729 of file gwf.f90.

730  ! modules
731  use constantsmodule, only: done, dp9
732  ! -- dummy
733  class(GwfModelType) :: this
734  integer(I4B), intent(in) :: neqmod
735  real(DP), dimension(neqmod), intent(inout) :: x
736  real(DP), dimension(neqmod), intent(in) :: xtemp
737  real(DP), dimension(neqmod), intent(inout) :: dx
738  integer(I4B), intent(inout) :: inewtonur
739  real(DP), intent(inout) :: dxmax
740  integer(I4B), intent(inout) :: locmax
741  ! -- local
742  integer(I4B) :: i0
743  integer(I4B) :: i1
744  class(BndType), pointer :: packobj
745  integer(I4B) :: ip
746  !
747  ! -- apply Newton-Raphson under-relaxation if model is using
748  ! the Newton-Raphson formulation and this Newton-Raphson
749  ! under-relaxation is turned on.
750  if (this%inewton /= 0 .and. this%inewtonur /= 0) then
751  if (this%innpf > 0) then
752  call this%npf%npf_nur(neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
753  end if
754  !
755  ! -- Call package nur routines
756  i0 = this%dis%nodes + 1
757  do ip = 1, this%bndlist%Count()
758  packobj => getbndfromlist(this%bndlist, ip)
759  if (packobj%npakeq > 0) then
760  i1 = i0 + packobj%npakeq - 1
761  call packobj%bnd_nur(packobj%npakeq, x(i0:i1), xtemp(i0:i1), &
762  dx(i0:i1), inewtonur, dxmax, locmax)
763  i0 = i1 + 1
764  end if
765  end do
766  end if
767  !
768  ! -- return
769  return
real(dp), parameter dp9
real constant 9/10
Definition: Constants.f90:71
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
Here is the call graph for this function:

◆ gwf_ot()

subroutine gwfmodule::gwf_ot ( class(gwfmodeltype this)

Definition at line 869 of file gwf.f90.

870  ! -- modules
871  use tdismodule, only: kstp, kper, tdis_ot, endofperiod
872  ! -- dummy
873  class(GwfModelType) :: this
874  ! -- local
875  integer(I4B) :: idvsave
876  integer(I4B) :: idvprint
877  integer(I4B) :: icbcfl
878  integer(I4B) :: icbcun
879  integer(I4B) :: ibudfl
880  integer(I4B) :: ipflag
881  ! -- formats
882  character(len=*), parameter :: fmtnocnvg = &
883  "(1X,/9X,'****FAILED TO MEET SOLVER CONVERGENCE CRITERIA IN TIME STEP ', &
884  &I0,' OF STRESS PERIOD ',I0,'****')"
885  !
886  ! -- Set write and print flags
887  idvsave = 0
888  idvprint = 0
889  icbcfl = 0
890  ibudfl = 0
891  if (this%oc%oc_save('HEAD')) idvsave = 1
892  if (this%oc%oc_print('HEAD')) idvprint = 1
893  if (this%oc%oc_save('BUDGET')) icbcfl = 1
894  if (this%oc%oc_print('BUDGET')) ibudfl = 1
895  icbcun = this%oc%oc_save_unit('BUDGET')
896  !
897  ! -- Override ibudfl and idvprint flags for nonconvergence
898  ! and end of period
899  ibudfl = this%oc%set_print_flag('BUDGET', this%icnvg, endofperiod)
900  idvprint = this%oc%set_print_flag('HEAD', this%icnvg, endofperiod)
901  !
902  ! Calculate and save observations
903  call this%gwf_ot_obs()
904  !
905  ! Save and print flows
906  call this%gwf_ot_flow(icbcfl, ibudfl, icbcun)
907  !
908  ! Save and print dependent variables
909  call this%gwf_ot_dv(idvsave, idvprint, ipflag)
910  !
911  ! Print budget summaries
912  call this%gwf_ot_bdsummary(ibudfl, ipflag)
913  !
914  ! -- Timing Output; if any dependent variables or budgets
915  ! are printed, then ipflag is set to 1.
916  if (ipflag == 1) call tdis_ot(this%iout)
917  !
918  ! -- Write non-convergence message
919  if (this%icnvg == 0) then
920  write (this%iout, fmtnocnvg) kstp, kper
921  end if
922  !
923  ! -- Return
924  return
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:349
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 1034 of file gwf.f90.

1035  use tdismodule, only: kstp, kper, totim, delt
1036  class(GwfModelType) :: this
1037  integer(I4B), intent(in) :: ibudfl
1038  integer(I4B), intent(inout) :: ipflag
1039  class(BndType), pointer :: packobj
1040  integer(I4B) :: ip
1041 
1042  !
1043  ! -- Package budget summary
1044  do ip = 1, this%bndlist%Count()
1045  packobj => getbndfromlist(this%bndlist, ip)
1046  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
1047  end do
1048 
1049  ! -- mover budget summary
1050  if (this%inmvr > 0) then
1051  call this%mvr%mvr_ot_bdsummary(ibudfl)
1052  end if
1053 
1054  ! -- model budget summary
1055  call this%budget%finalize_step(delt)
1056  if (ibudfl /= 0) then
1057  ipflag = 1
1058  call this%budget%budget_ot(kstp, kper, this%iout)
1059  end if
1060 
1061  ! -- Write to budget csv every time step
1062  call this%budget%writecsv(totim)
1063 
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 1000 of file gwf.f90.

1001  class(GwfModelType) :: this
1002  integer(I4B), intent(in) :: idvsave
1003  integer(I4B), intent(in) :: idvprint
1004  integer(I4B), intent(inout) :: ipflag
1005  class(BndType), pointer :: packobj
1006  integer(I4B) :: ip
1007  !
1008  ! -- Save compaction to binary file
1009  if (this%incsub > 0) call this%csub%csub_ot_dv(idvsave, idvprint)
1010  !
1011  ! -- save density to binary file
1012  if (this%inbuy > 0) then
1013  call this%buy%buy_ot_dv(idvsave)
1014  end if
1015  !
1016  ! -- save viscosity to binary file
1017  if (this%invsc > 0) then
1018  call this%vsc%vsc_ot_dv(idvsave)
1019  end if
1020  !
1021  ! -- Print advanced package dependent variables
1022  do ip = 1, this%bndlist%Count()
1023  packobj => getbndfromlist(this%bndlist, ip)
1024  call packobj%bnd_ot_dv(idvsave, idvprint)
1025  end do
1026  !
1027  ! -- save head and print head
1028  call this%oc%oc_ot(ipflag)
1029  !
1030  ! -- Return
1031  return
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 951 of file gwf.f90.

952  class(GwfModelType) :: this
953  integer(I4B), intent(in) :: icbcfl
954  integer(I4B), intent(in) :: ibudfl
955  integer(I4B), intent(in) :: icbcun
956  class(BndType), pointer :: packobj
957  integer(I4B) :: ip
958 
959  ! -- Save GWF flows
960  if (this%insto > 0) then
961  call this%sto%sto_save_model_flows(icbcfl, icbcun)
962  end if
963  if (this%innpf > 0) then
964  call this%npf%npf_save_model_flows(this%flowja, icbcfl, icbcun)
965  end if
966  if (this%incsub > 0) call this%csub%csub_save_model_flows(icbcfl, icbcun)
967  do ip = 1, this%bndlist%Count()
968  packobj => getbndfromlist(this%bndlist, ip)
969  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=0, icbcun=icbcun)
970  end do
971 
972  ! -- Save advanced package flows
973  do ip = 1, this%bndlist%Count()
974  packobj => getbndfromlist(this%bndlist, ip)
975  call packobj%bnd_ot_package_flows(icbcfl=icbcfl, ibudfl=0)
976  end do
977  if (this%inmvr > 0) then
978  call this%mvr%mvr_ot_saveflow(icbcfl, ibudfl)
979  end if
980 
981  ! -- Print GWF flows
982  if (this%innpf > 0) call this%npf%npf_print_model_flows(ibudfl, this%flowja)
983  if (this%ingnc > 0) call this%gnc%gnc_ot(ibudfl)
984  do ip = 1, this%bndlist%Count()
985  packobj => getbndfromlist(this%bndlist, ip)
986  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=ibudfl, icbcun=0)
987  end do
988 
989  ! -- Print advanced package flows
990  do ip = 1, this%bndlist%Count()
991  packobj => getbndfromlist(this%bndlist, ip)
992  call packobj%bnd_ot_package_flows(icbcfl=0, ibudfl=ibudfl)
993  end do
994  if (this%inmvr > 0) then
995  call this%mvr%mvr_ot_printflow(icbcfl, ibudfl)
996  end if
997 
Here is the call graph for this function:

◆ gwf_ot_obs()

subroutine gwfmodule::gwf_ot_obs ( class(gwfmodeltype this)

Definition at line 927 of file gwf.f90.

928  class(GwfModelType) :: this
929  class(BndType), pointer :: packobj
930  integer(I4B) :: ip
931 
932  ! -- Calculate and save GWF observations
933  call this%obs%obs_bd()
934  call this%obs%obs_ot()
935 
936  ! -- Calculate and save csub observations
937  if (this%incsub > 0) then
938  call this%csub%csub_bd_obs()
939  call this%csub%obs%obs_ot()
940  end if
941 
942  ! -- Calculate and save package observations
943  do ip = 1, this%bndlist%Count()
944  packobj => getbndfromlist(this%bndlist, ip)
945  call packobj%bnd_bd_obs()
946  call packobj%bnd_ot_obs()
947  end do
948 
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 655 of file gwf.f90.

656  ! -- modules
657  use constantsmodule, only: done
658  use tdismodule, only: delt
659  ! -- dummy
660  class(GwfModelType) :: this
661  class(VectorBaseType), pointer :: vec_residual
662  integer(I4B), intent(inout) :: iptc
663  real(DP), intent(inout) :: ptcf
664  ! -- local
665  integer(I4B) :: n
666  integer(I4B) :: iptct
667  real(DP) :: v
668  real(DP) :: resid
669  real(DP) :: ptcdelem1
670  !
671  ! -- set temporary flag indicating if pseudo-transient continuation should
672  ! be used for this model and time step
673  iptct = 0
674  ! -- only apply pseudo-transient continuation to problems using the
675  ! Newton-Raphson formulations for steady-state stress periods
676  if (this%iss > 0) then
677  if (this%inewton > 0) then
678  iptct = this%inewton
679  else
680  iptct = this%npf%inewton
681  end if
682  end if
683  !
684  ! -- calculate pseudo-transient continuation factor for model
685  if (iptct > 0) then
686  !
687  ! -- calculate the pseudo-time step using the residual
688  do n = 1, this%dis%nodes
689  if (this%npf%ibound(n) < 1) cycle
690  !
691  ! -- get the maximum volume of the cell (head at top of cell)
692  v = this%dis%get_cell_volume(n, this%dis%top(n))
693  !
694  ! -- set the residual
695  resid = vec_residual%get_value_local(n)
696  !
697  ! -- calculate the reciprocal of the pseudo-time step
698  ! resid [L3/T] / volume [L3] = [1/T]
699  ptcdelem1 = abs(resid) / v
700  !
701  ! -- set ptcf if the reciprocal of the pseudo-time step
702  ! exceeds the current value (equivalent to using the
703  ! smallest pseudo-time step)
704  if (ptcdelem1 > ptcf) ptcf = ptcdelem1
705  end do
706  !
707  ! -- protection for the case where the residuals are zero
708  if (ptcf == dzero) then
709  ptcf = done / (delt * dten)
710  end if
711  end if
712  !
713  ! -- reset ipc if needed
714  if (iptc == 0) then
715  if (iptct > 0) iptc = 1
716  end if
717  !
718  ! -- return
719  return

◆ 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 628 of file gwf.f90.

629  ! -- dummy
630  class(GwfModelType) :: this
631  integer(I4B), intent(inout) :: iptc
632  !
633  ! -- determine if pseudo-transient continuation should be applied to this
634  ! model - pseudo-transient continuation only applied to problems that
635  ! use the Newton-Raphson formulation during steady-state stress periods
636  iptc = 0
637  if (this%iss > 0) then
638  if (this%inewton > 0) then
639  iptc = this%inewton
640  else
641  iptc = this%npf%inewton
642  end if
643  end if
644  !
645  ! -- return
646  return

◆ gwf_rp()

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

(1) calls package read and prepare routines

Definition at line 374 of file gwf.f90.

375  ! -- modules
376  use tdismodule, only: readnewdata
377  ! -- dummy
378  class(GwfModelType) :: this
379  ! -- local
380  class(BndType), pointer :: packobj
381  integer(I4B) :: ip
382  !
383  ! -- Check with TDIS on whether or not it is time to RP
384  if (.not. readnewdata) return
385  !
386  ! -- Read and prepare
387  if (this%innpf > 0) call this%npf%npf_rp()
388  if (this%inbuy > 0) call this%buy%buy_rp()
389  if (this%invsc > 0) call this%vsc%vsc_rp()
390  if (this%inhfb > 0) call this%hfb%hfb_rp()
391  if (this%inoc > 0) call this%oc%oc_rp()
392  if (this%insto > 0) call this%sto%sto_rp()
393  if (this%incsub > 0) call this%csub%csub_rp()
394  if (this%inmvr > 0) call this%mvr%mvr_rp()
395  do ip = 1, this%bndlist%Count()
396  packobj => getbndfromlist(this%bndlist, ip)
397  call packobj%bnd_rp()
398  call packobj%bnd_rp_obs()
399  end do
400  !
401  ! -- Check for steady state period
402  call this%steady_period_check()
403  !
404  ! -- Return
405  return
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 1595 of file gwf.f90.

1597  class(GwfModelType) :: this
1598  type(GwfNamParamFoundType), intent(in) :: found
1599 
1600  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1601 
1602  if (found%newton) then
1603  write (this%iout, '(4x,a)') &
1604  'NEWTON-RAPHSON method enabled for the model.'
1605  if (found%under_relaxation) then
1606  write (this%iout, '(4x,a,a)') &
1607  'NEWTON-RAPHSON UNDER-RELAXATION based on the bottom ', &
1608  'elevation of the model will be applied to the model.'
1609  end if
1610  end if
1611 
1612  if (found%print_input) then
1613  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1614  'FOR ALL MODEL STRESS PACKAGES'
1615  end if
1616 
1617  if (found%print_flows) then
1618  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1619  'FOR ALL MODEL PACKAGES'
1620  end if
1621 
1622  if (found%save_flows) then
1623  write (this%iout, '(4x,a)') &
1624  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1625  end if
1626 
1627  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 1266 of file gwf.f90.

1268  ! -- modules
1269  use constantsmodule, only: linelength
1270  use simmodule, only: store_error
1271  use chdmodule, only: chd_create
1272  use welmodule, only: wel_create
1273  use drnmodule, only: drn_create
1274  use rivmodule, only: riv_create
1275  use ghbmodule, only: ghb_create
1276  use rchmodule, only: rch_create
1277  use evtmodule, only: evt_create
1278  use mawmodule, only: maw_create
1279  use sfrmodule, only: sfr_create
1280  use lakmodule, only: lak_create
1281  use uzfmodule, only: uzf_create
1282  use apimodule, only: api_create
1283  ! -- dummy
1284  class(GwfModelType) :: this
1285  character(len=*), intent(in) :: filtyp
1286  integer(I4B), intent(in) :: ipakid
1287  integer(I4B), intent(in) :: ipaknum
1288  character(len=*), intent(in) :: pakname
1289  character(len=*), intent(in) :: mempath
1290  integer(I4B), intent(in) :: inunit
1291  integer(I4B), intent(in) :: iout
1292  ! -- local
1293  class(BndType), pointer :: packobj
1294  class(BndType), pointer :: packobj2
1295  integer(I4B) :: ip
1296  !
1297  ! -- This part creates the package object
1298  select case (filtyp)
1299  case ('CHD6')
1300  call chd_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1301  pakname, mempath)
1302  case ('WEL6')
1303  call wel_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1304  pakname, mempath)
1305  case ('DRN6')
1306  call drn_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1307  pakname, mempath)
1308  case ('RIV6')
1309  call riv_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1310  pakname, mempath)
1311  case ('GHB6')
1312  call ghb_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1313  pakname, mempath)
1314  case ('RCH6')
1315  call rch_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1316  pakname, mempath)
1317  case ('EVT6')
1318  call evt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1319  pakname, mempath)
1320  case ('MAW6')
1321  call maw_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1322  case ('SFR6')
1323  call sfr_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1324  case ('LAK6')
1325  call lak_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1326  case ('UZF6')
1327  call uzf_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1328  case ('API6')
1329  call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1330  case default
1331  write (errmsg, *) 'Invalid package type: ', filtyp
1332  call store_error(errmsg, terminate=.true.)
1333  end select
1334  !
1335  ! -- Check to make sure that the package name is unique, then store a
1336  ! pointer to the package in the model bndlist
1337  do ip = 1, this%bndlist%Count()
1338  packobj2 => getbndfromlist(this%bndlist, ip)
1339  if (packobj2%packName == pakname) then
1340  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
1341  'already exists: ', trim(pakname)
1342  call store_error(errmsg, terminate=.true.)
1343  end if
1344  end do
1345  call addbndtolist(this%bndlist, packobj)
1346  !
1347  ! -- return
1348  return
This module contains the API package methods.
Definition: gwf-api.f90:12
subroutine, public api_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname)
@ brief Create a new package object
Definition: gwf-api.f90:49
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:70
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:63
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:236
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:175
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 1637 of file gwf.f90.

1638  ! -- modules
1639  use tdismodule, only: kper
1641  use simvariablesmodule, only: warnmsg
1642  use simmodule, only: store_warning
1643  ! -- dummy
1644  class(GwfModelType) :: this
1645  if (this%iss == 1) then
1646  if (isadaptiveperiod(kper)) then
1647  write (warnmsg, '(a,a,a,i0,a)') &
1648  'GWF Model (', trim(this%name), ') is steady state for period ', &
1649  kper, ' and adaptive time stepping is active. Adaptive time &
1650  &stepping may not work properly for steady-state conditions.'
1651  call store_warning(warnmsg)
1652  end if
1653  end if
1654  return
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