MODFLOW 6  version 6.7.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 1166 of file gwf.f90.

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

1332  implicit none
1333  class(*), pointer, intent(inout) :: model
1334  class(GwfModelType), pointer :: gwfModel
1335 
1336  gwfmodel => null()
1337  if (.not. associated(model)) return
1338  select type (model)
1339  class is (gwfmodeltype)
1340  gwfmodel => model
1341  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 1346 of file gwf.f90.

1348  ! -- modules
1351  ! -- dummy
1352  class(GwfModelType) :: this
1353  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
1354  type(CharacterStringType), dimension(:), contiguous, &
1355  pointer, intent(inout) :: pkgtypes
1356  type(CharacterStringType), dimension(:), contiguous, &
1357  pointer, intent(inout) :: pkgnames
1358  type(CharacterStringType), dimension(:), contiguous, &
1359  pointer, intent(inout) :: mempaths
1360  integer(I4B), dimension(:), contiguous, &
1361  pointer, intent(inout) :: inunits
1362  ! -- local
1363  integer(I4B) :: ipakid, ipaknum
1364  character(len=LENFTYPE) :: pkgtype, bndptype
1365  character(len=LENPACKAGENAME) :: pkgname
1366  character(len=LENMEMPATH) :: mempath
1367  integer(I4B), pointer :: inunit
1368  integer(I4B) :: n
1369 
1370  if (allocated(bndpkgs)) then
1371  !
1372  ! -- create stress packages
1373  ipakid = 1
1374  bndptype = ''
1375  do n = 1, size(bndpkgs)
1376  !
1377  pkgtype = pkgtypes(bndpkgs(n))
1378  pkgname = pkgnames(bndpkgs(n))
1379  mempath = mempaths(bndpkgs(n))
1380  inunit => inunits(bndpkgs(n))
1381  !
1382  if (bndptype /= pkgtype) then
1383  ipaknum = 1
1384  bndptype = pkgtype
1385  end if
1386  !
1387  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1388  inunit, this%iout)
1389  ipakid = ipakid + 1
1390  ipaknum = ipaknum + 1
1391  end do
1392  !
1393  ! -- cleanup
1394  deallocate (bndpkgs)
1395  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 1400 of file gwf.f90.

1401  ! -- modules
1404  use arrayhandlersmodule, only: expandarray
1405  use memorymanagermodule, only: mem_setptr
1407  use simvariablesmodule, only: idm_context
1408  use dismodule, only: dis_cr
1409  use disvmodule, only: disv_cr
1410  use disumodule, only: disu_cr
1411  use gwfnpfmodule, only: npf_cr
1412  use xt3dmodule, only: xt3d_cr
1413  use gwfbuymodule, only: buy_cr
1414  use gwfvscmodule, only: vsc_cr
1415  use gwfstomodule, only: sto_cr
1416  use gwfcsubmodule, only: csub_cr
1417  use gwfmvrmodule, only: mvr_cr
1418  use gwfhfbmodule, only: hfb_cr
1419  use gwficmodule, only: ic_cr
1420  use gwfocmodule, only: oc_cr
1421  ! -- dummy
1422  class(GwfModelType) :: this
1423  ! -- local
1424  type(CharacterStringType), dimension(:), contiguous, &
1425  pointer :: pkgtypes => null()
1426  type(CharacterStringType), dimension(:), contiguous, &
1427  pointer :: pkgnames => null()
1428  type(CharacterStringType), dimension(:), contiguous, &
1429  pointer :: mempaths => null()
1430  integer(I4B), dimension(:), contiguous, &
1431  pointer :: inunits => null()
1432  character(len=LENMEMPATH) :: model_mempath
1433  character(len=LENFTYPE) :: pkgtype
1434  character(len=LENPACKAGENAME) :: pkgname
1435  character(len=LENMEMPATH) :: mempath
1436  integer(I4B), pointer :: inunit
1437  integer(I4B), dimension(:), allocatable :: bndpkgs
1438  integer(I4B) :: n
1439  integer(I4B) :: indis = 0 ! DIS enabled flag
1440  character(len=LENMEMPATH) :: mempathnpf = ''
1441  character(len=LENMEMPATH) :: mempathic = ''
1442  character(len=LENMEMPATH) :: mempathsto = ''
1443  !
1444  ! -- set input model memory path
1445  model_mempath = create_mem_path(component=this%name, context=idm_context)
1446  !
1447  ! -- set pointers to model path package info
1448  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1449  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1450  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1451  call mem_setptr(inunits, 'INUNITS', model_mempath)
1452  !
1453  do n = 1, size(pkgtypes)
1454  !
1455  ! attributes for this input package
1456  pkgtype = pkgtypes(n)
1457  pkgname = pkgnames(n)
1458  mempath = mempaths(n)
1459  inunit => inunits(n)
1460  !
1461  ! -- create dis package as it is a prerequisite for other packages
1462  select case (pkgtype)
1463  case ('DIS6')
1464  indis = 1
1465  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1466  case ('DISV6')
1467  indis = 1
1468  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1469  case ('DISU6')
1470  indis = 1
1471  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1472  case ('NPF6')
1473  this%innpf = 1
1474  mempathnpf = mempath
1475  case ('BUY6')
1476  this%inbuy = inunit
1477  case ('VSC6')
1478  this%invsc = inunit
1479  case ('GNC6')
1480  this%ingnc = inunit
1481  case ('HFB6')
1482  this%inhfb = inunit
1483  case ('STO6')
1484  this%insto = 1
1485  mempathsto = mempath
1486  case ('CSUB6')
1487  this%incsub = inunit
1488  case ('IC6')
1489  this%inic = 1
1490  mempathic = mempath
1491  case ('MVR6')
1492  this%inmvr = inunit
1493  case ('OC6')
1494  this%inoc = inunit
1495  case ('OBS6')
1496  this%inobs = inunit
1497  case ('WEL6', 'DRN6', 'RIV6', 'GHB6', 'RCH6', &
1498  'EVT6', 'API6', 'CHD6', 'MAW6', 'SFR6', &
1499  'LAK6', 'UZF6')
1500  call expandarray(bndpkgs)
1501  bndpkgs(size(bndpkgs)) = n
1502  case default
1503  ! TODO
1504  end select
1505  end do
1506  !
1507  ! -- Create packages that are tied directly to model
1508  call npf_cr(this%npf, this%name, mempathnpf, this%innpf, this%iout)
1509  call xt3d_cr(this%xt3d, this%name, this%innpf, this%iout)
1510  call buy_cr(this%buy, this%name, this%inbuy, this%iout)
1511  call vsc_cr(this%vsc, this%name, this%invsc, this%iout)
1512  call gnc_cr(this%gnc, this%name, this%ingnc, this%iout)
1513  call hfb_cr(this%hfb, this%name, this%inhfb, this%iout)
1514  call sto_cr(this%sto, this%name, mempathsto, this%insto, this%iout)
1515  call csub_cr(this%csub, this%name, this%insto, this%sto%packName, &
1516  this%incsub, this%iout)
1517  call ic_cr(this%ic, this%name, mempathic, this%inic, this%iout, this%dis)
1518  call mvr_cr(this%mvr, this%name, this%inmvr, this%iout, this%dis)
1519  call oc_cr(this%oc, this%name, this%inoc, this%iout)
1520  call gwf_obs_cr(this%obs, this%inobs)
1521  !
1522  ! -- Check to make sure that required ftype's have been specified
1523  call this%ftype_check(indis)
1524  !
1525  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: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:103
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:325
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: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: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, mempath, inunit, iout)
@ brief Create a new package object
Definition: gwf-sto.f90:78
subroutine, public vsc_cr(vscobj, name_model, inunit, iout)
@ brief Create a new package object
Definition: gwf-vsc.f90:140
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 1296 of file gwf.f90.

1297  ! -- modules
1298  use constantsmodule, only: linelength
1299  use simmodule, only: store_error, count_errors
1300  ! -- dummy
1301  class(GwfModelType) :: this
1302  integer(I4B), intent(in) :: indis
1303  ! -- local
1304  !
1305  ! -- Check for IC8, DIS(u), and NPF. Stop if not present.
1306  if (this%inic == 0) then
1307  write (errmsg, '(a)') &
1308  'Initial Conditions (IC6) package not specified.'
1309  call store_error(errmsg)
1310  end if
1311  if (indis == 0) then
1312  write (errmsg, '(a)') &
1313  'Discretization (DIS6, DISV6, or DISU6) Package not specified.'
1314  call store_error(errmsg)
1315  end if
1316  if (this%innpf == 0) then
1317  write (errmsg, '(a)') &
1318  'Node Property Flow (NPF6) Package not specified.'
1319  call store_error(errmsg)
1320  end if
1321  !
1322  if (count_errors() > 0) then
1323  write (errmsg, '(a)') 'One or more required package(s) not specified.'
1324  call store_error(errmsg)
1325  call store_error_filename(this%filename)
1326  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 396 of file gwf.f90.

397  ! -- modules
399  ! -- dummy
400  class(GwfModelType) :: this
401  class(BndType), pointer :: packobj
402  ! -- local
403  integer(I4B) :: irestore
404  integer(I4B) :: ip, n
405  !
406  ! -- Reset state variable
407  irestore = 0
408  if (ifailedstepretry > 0) irestore = 1
409  if (irestore == 0) then
410  !
411  ! -- copy x into xold
412  do n = 1, this%dis%nodes
413  this%xold(n) = this%x(n)
414  end do
415  else
416  !
417  ! -- copy xold into x if this time step is a redo
418  do n = 1, this%dis%nodes
419  this%x(n) = this%xold(n)
420  end do
421  end if
422  !
423  ! -- Advance
424  if (this%invsc > 0) call this%vsc%vsc_ad()
425  if (this%innpf > 0) call this%npf%npf_ad(this%dis%nodes, this%xold, &
426  this%x, irestore)
427  if (this%insto > 0) call this%sto%sto_ad()
428  if (this%incsub > 0) call this%csub%csub_ad(this%dis%nodes, this%x)
429  if (this%inbuy > 0) call this%buy%buy_ad()
430  if (this%inmvr > 0) call this%mvr%mvr_ad()
431  do ip = 1, this%bndlist%Count()
432  packobj => getbndfromlist(this%bndlist, ip)
433  call packobj%bnd_ad()
434  if (this%invsc > 0) call this%vsc%vsc_ad_bnd(packobj, this%x)
435  if (isimcheck > 0) then
436  call packobj%bnd_ck()
437  end if
438  end do
439  !
440  ! -- Push simulated values to preceding time/subtime step
441  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 785 of file gwf.f90.

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

1121  ! -- modules
1122  use constantsmodule, only: lenbudtxt
1123  use tdismodule, only: delt
1124  ! -- dummy
1125  class(GwfModelType) :: this
1126  real(DP), dimension(:, :), intent(in) :: budterm
1127  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
1128  character(len=*), intent(in) :: rowlabel
1129  !
1130  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 562 of file gwf.f90.

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

447  ! -- dummy
448  class(GwfModelType) :: this
449  integer(I4B), intent(in) :: kiter
450  ! -- local
451  class(BndType), pointer :: packobj
452  integer(I4B) :: ip
453  !
454  ! -- Call package cf routines
455  if (this%innpf > 0) call this%npf%npf_cf(kiter, this%dis%nodes, this%x)
456  if (this%inbuy > 0) call this%buy%buy_cf(kiter)
457  do ip = 1, this%bndlist%Count()
458  packobj => getbndfromlist(this%bndlist, ip)
459  call packobj%bnd_cf()
460  if (this%inbuy > 0) call this%buy%buy_cf_bnd(packobj, this%x)
461  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 741 of file gwf.f90.

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

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

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

◆ gwf_fp()

subroutine gwfmodule::gwf_fp ( class(gwfmodeltype this)

Definition at line 1028 of file gwf.f90.

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

◆ gwf_get_iasym()

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

Definition at line 1136 of file gwf.f90.

1137  class(GwfModelType) :: this
1138  ! -- local
1139  integer(I4B) :: iasym
1140  integer(I4B) :: ip
1141  class(BndType), pointer :: packobj
1142  !
1143  ! -- Start by setting iasym to zero
1144  iasym = 0
1145  !
1146  ! -- NPF
1147  if (this%innpf > 0) then
1148  if (this%npf%iasym /= 0) iasym = 1
1149  if (this%npf%ixt3d /= 0) iasym = 1
1150  end if
1151  !
1152  ! -- GNC
1153  if (this%ingnc > 0) then
1154  if (this%gnc%iasym /= 0) iasym = 1
1155  end if
1156  !
1157  ! -- Check for any packages that introduce matrix asymmetry
1158  do ip = 1, this%bndlist%Count()
1159  packobj => getbndfromlist(this%bndlist, ip)
1160  if (packobj%iasym /= 0) iasym = 1
1161  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 696 of file gwf.f90.

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

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

996  use tdismodule, only: kstp, kper, totim, delt
997  class(GwfModelType) :: this
998  integer(I4B), intent(in) :: ibudfl
999  integer(I4B), intent(inout) :: ipflag
1000  class(BndType), pointer :: packobj
1001  integer(I4B) :: ip
1002 
1003  ! -- Package budget summary
1004  do ip = 1, this%bndlist%Count()
1005  packobj => getbndfromlist(this%bndlist, ip)
1006  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
1007  end do
1008 
1009  ! -- mover budget summary
1010  if (this%inmvr > 0) then
1011  call this%mvr%mvr_ot_bdsummary(ibudfl)
1012  end if
1013 
1014  ! -- model budget summary
1015  call this%budget%finalize_step(delt)
1016  if (ibudfl /= 0) then
1017  ipflag = 1
1018  call this%budget%budget_ot(kstp, kper, this%iout)
1019  end if
1020 
1021  ! -- Write to budget csv every time step
1022  call this%budget%writecsv(totim)
1023 
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 962 of file gwf.f90.

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

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

◆ gwf_ot_obs()

subroutine gwfmodule::gwf_ot_obs ( class(gwfmodeltype this)

Definition at line 885 of file gwf.f90.

886  class(GwfModelType) :: this
887  class(BndType), pointer :: packobj
888  integer(I4B) :: ip
889 
890  ! -- Calculate and save GWF observations
891  call this%obs%obs_bd()
892  call this%obs%obs_ot()
893 
894  ! -- Calculate and save csub observations
895  if (this%incsub > 0) then
896  call this%csub%csub_bd_obs()
897  call this%csub%obs%obs_ot()
898  end if
899 
900  ! -- Calculate and save package observations
901  do ip = 1, this%bndlist%Count()
902  packobj => getbndfromlist(this%bndlist, ip)
903  call packobj%bnd_bd_obs()
904  call packobj%bnd_ot_obs()
905  end do
906 
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 625 of file gwf.f90.

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

602  ! -- dummy
603  class(GwfModelType) :: this
604  integer(I4B), intent(inout) :: iptc
605  !
606  ! -- determine if pseudo-transient continuation should be applied to this
607  ! model - pseudo-transient continuation only applied to problems that
608  ! use the Newton-Raphson formulation during steady-state stress periods
609  iptc = 0
610  if (this%iss > 0) then
611  if (this%inewton > 0) then
612  iptc = this%inewton
613  else
614  iptc = this%npf%inewton
615  end if
616  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_obs()
385  end do
386  !
387  ! -- Check for steady state period
388  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 1530 of file gwf.f90.

1532  class(GwfModelType) :: this
1533  type(GwfNamParamFoundType), intent(in) :: found
1534 
1535  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1536 
1537  if (found%newton) then
1538  write (this%iout, '(4x,a)') &
1539  'NEWTON-RAPHSON method enabled for the model.'
1540  if (found%under_relaxation) then
1541  write (this%iout, '(4x,a,a)') &
1542  'NEWTON-RAPHSON UNDER-RELAXATION based on the bottom ', &
1543  'elevation of the model will be applied to the model.'
1544  end if
1545  end if
1546 
1547  if (found%print_input) then
1548  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1549  'FOR ALL MODEL STRESS PACKAGES'
1550  end if
1551 
1552  if (found%print_flows) then
1553  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1554  'FOR ALL MODEL PACKAGES'
1555  end if
1556 
1557  if (found%save_flows) then
1558  write (this%iout, '(4x,a)') &
1559  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1560  end if
1561 
1562  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 1212 of file gwf.f90.

1214  ! -- modules
1215  use constantsmodule, only: linelength
1216  use simmodule, only: store_error
1217  use chdmodule, only: chd_create
1218  use welmodule, only: wel_create
1219  use drnmodule, only: drn_create
1220  use rivmodule, only: riv_create
1221  use ghbmodule, only: ghb_create
1222  use rchmodule, only: rch_create
1223  use evtmodule, only: evt_create
1224  use mawmodule, only: maw_create
1225  use sfrmodule, only: sfr_create
1226  use lakmodule, only: lak_create
1227  use uzfmodule, only: uzf_create
1228  use apimodule, only: api_create
1229  ! -- dummy
1230  class(GwfModelType) :: this
1231  character(len=*), intent(in) :: filtyp
1232  integer(I4B), intent(in) :: ipakid
1233  integer(I4B), intent(in) :: ipaknum
1234  character(len=*), intent(in) :: pakname
1235  character(len=*), intent(in) :: mempath
1236  integer(I4B), intent(in) :: inunit
1237  integer(I4B), intent(in) :: iout
1238  ! -- local
1239  class(BndType), pointer :: packobj
1240  class(BndType), pointer :: packobj2
1241  integer(I4B) :: ip
1242  !
1243  ! -- This part creates the package object
1244  select case (filtyp)
1245  case ('CHD6')
1246  call chd_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1247  pakname, mempath)
1248  case ('WEL6')
1249  call wel_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1250  pakname, mempath)
1251  case ('DRN6')
1252  call drn_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1253  pakname, mempath)
1254  case ('RIV6')
1255  call riv_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1256  pakname, mempath)
1257  case ('GHB6')
1258  call ghb_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1259  pakname, mempath)
1260  case ('RCH6')
1261  call rch_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1262  pakname, mempath)
1263  case ('EVT6')
1264  call evt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
1265  pakname, mempath)
1266  case ('MAW6')
1267  call maw_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1268  case ('SFR6')
1269  call sfr_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1270  case ('LAK6')
1271  call lak_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1272  case ('UZF6')
1273  call uzf_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1274  case ('API6')
1275  call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
1276  case default
1277  write (errmsg, *) 'Invalid package type: ', filtyp
1278  call store_error(errmsg, terminate=.true.)
1279  end select
1280  !
1281  ! -- Check to make sure that the package name is unique, then store a
1282  ! pointer to the package in the model bndlist
1283  do ip = 1, this%bndlist%Count()
1284  packobj2 => getbndfromlist(this%bndlist, ip)
1285  if (packobj2%packName == pakname) then
1286  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
1287  'already exists: ', trim(pakname)
1288  call store_error(errmsg, terminate=.true.)
1289  end if
1290  end do
1291  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)
@ 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:69
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:296
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 1572 of file gwf.f90.

1573  ! -- modules
1574  use tdismodule, only: kper
1576  use simvariablesmodule, only: warnmsg
1577  use simmodule, only: store_warning
1578  ! -- dummy
1579  class(GwfModelType) :: this
1580  if (this%iss == 1) then
1581  if (isadaptiveperiod(kper)) then
1582  write (warnmsg, '(a,a,a,i0,a)') &
1583  'GWF Model (', trim(this%name), ') is steady state for period ', &
1584  kper, ' and adaptive time stepping is active. Adaptive time &
1585  &stepping may not work properly for steady-state conditions.'
1586  call store_warning(warnmsg)
1587  end if
1588  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