MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
gwfgwfexchangemodule Module Reference

This module contains the GwfGwfExchangeModule Module. More...

Data Types

type  gwfexchangetype
 Derived type for GwfExchangeType. More...
 

Functions/Subroutines

subroutine, public gwfexchange_create (filename, name, id, m1_id, m2_id, input_mempath)
 @ brief Create GWF GWF exchange More...
 
subroutine gwf_gwf_df (this)
 @ brief Define GWF GWF exchange More...
 
subroutine validate_exchange (this)
 validate exchange data after reading More...
 
subroutine gwf_gwf_ac (this, sparse)
 @ brief Add connections More...
 
subroutine gwf_gwf_mc (this, matrix_sln)
 @ brief Map connections More...
 
subroutine gwf_gwf_ar (this)
 @ brief Allocate and read More...
 
subroutine gwf_gwf_rp (this)
 @ brief Read and prepare More...
 
subroutine gwf_gwf_ad (this)
 @ brief Advance More...
 
subroutine gwf_gwf_cf (this, kiter)
 @ brief Calculate coefficients More...
 
subroutine gwf_gwf_fc (this, kiter, matrix_sln, rhs_sln, inwtflag)
 @ brief Fill coefficients More...
 
subroutine gwf_gwf_fn (this, kiter, matrix_sln)
 @ brief Fill Newton More...
 
subroutine gwf_gwf_cq (this, icnvg, isuppress_output, isolnid)
 @ brief Calculate flow More...
 
subroutine gwf_gwf_calc_simvals (this)
 Calculate flow rates for the exchanges and store them in a member array. More...
 
subroutine gwf_gwf_add_to_flowja (this)
 Add exchange flow to each model flowja diagonal position so that residual is calculated correctly. More...
 
subroutine gwf_gwf_set_flow_to_npf (this)
 Set flow rates to the edges in the models. More...
 
subroutine gwf_gwf_bd (this, icnvg, isuppress_output, isolnid)
 @ brief Budget More...
 
subroutine gwf_gwf_chd_bd (this)
 @ brief gwf-gwf-chd-bd More...
 
subroutine gwf_gwf_bdsav (this)
 @ brief Budget save More...
 
subroutine gwf_gwf_bdsav_model (this, model)
 
subroutine gwf_gwf_ot (this)
 @ brief Output More...
 
subroutine source_options (this, iout)
 @ brief Source options More...
 
subroutine read_gnc (this)
 @ brief Read ghost nodes More...
 
subroutine read_mvr (this, iout)
 @ brief Read mover More...
 
subroutine rewet (this, kiter)
 @ brief Rewet More...
 
subroutine calc_cond_sat (this)
 
subroutine condcalc (this)
 @ brief Calculate the conductance More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine gwf_gwf_da (this)
 @ brief Deallocate More...
 
subroutine allocate_arrays (this)
 @ brief Allocate arrays More...
 
subroutine gwf_gwf_df_obs (this)
 @ brief Define observations More...
 
subroutine gwf_gwf_rp_obs (this)
 @ brief Read and prepare observations More...
 
subroutine gwf_gwf_fp (this)
 @ brief Final processing More...
 
real(dp) function qcalc (this, iexg, n1, n2)
 @ brief Calculate flow More...
 
integer(i4b) function gwf_gwf_get_iasym (this)
 @ brief Set symmetric flag More...
 
logical(lgp) function gwf_gwf_connects_model (this, model)
 Return true when this exchange provides matrix coefficients for solving. More...
 
logical(lgp) function use_interface_model (this)
 Should interface model be used for this exchange. More...
 
subroutine gwf_gwf_save_simvals (this)
 @ brief Save simulated flow observations More...
 
subroutine gwf_gwf_process_obsid (obsrv, dis, inunitobs, iout)
 @ brief Obs ID processor More...
 
class(gwfexchangetype) function, pointer, public castasgwfexchange (obj)
 @ brief Cast polymorphic object as exchange More...
 
class(gwfexchangetype) function, pointer, public getgwfexchangefromlist (list, idx)
 @ brief Get exchange from list More...
 

Detailed Description

This module contains the code for connecting two GWF Models. The methods are based on the simple two point flux approximation with the option to use ghost nodes to improve accuracy. This exchange is used by GwfGwfConnection with the more sophisticated interface model coupling approach when XT3D is needed.

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwfgwfexchangemodule::allocate_arrays ( class(gwfexchangetype this)

Allocate arrays

Parameters
thisGwfExchangeType

Definition at line 1732 of file exg-gwfgwf.f90.

1733  ! -- modules
1735  ! -- dummy
1736  class(GwfExchangeType) :: this !< GwfExchangeType
1737  ! -- local
1738  character(len=LINELENGTH) :: text
1739  integer(I4B) :: ntabcol, i
1740  !
1741  call this%DisConnExchangeType%allocate_arrays()
1742  !
1743  call mem_allocate(this%cond, this%nexg, 'COND', this%memoryPath)
1744  call mem_allocate(this%idxglo, this%nexg, 'IDXGLO', this%memoryPath)
1745  call mem_allocate(this%idxsymglo, this%nexg, 'IDXSYMGLO', this%memoryPath) !
1746  call mem_allocate(this%condsat, this%nexg, 'CONDSAT', this%memoryPath)
1747  call mem_allocate(this%simvals, this%nexg, 'SIMVALS', this%memoryPath)
1748  !
1749  ! -- Initialize
1750  do i = 1, this%nexg
1751  this%cond(i) = dnodata
1752  end do
1753  !
1754  ! -- allocate and initialize the output table
1755  if (this%iprflow /= 0) then
1756  !
1757  ! -- dimension table
1758  ntabcol = 3
1759  if (this%inamedbound > 0) then
1760  ntabcol = ntabcol + 1
1761  end if
1762  !
1763  ! -- initialize the output table objects
1764  ! outouttab1
1765  if (this%v_model1%is_local) then
1766  call table_cr(this%outputtab1, this%name, ' ')
1767  call this%outputtab1%table_df(this%nexg, ntabcol, this%gwfmodel1%iout, &
1768  transient=.true.)
1769  text = 'NUMBER'
1770  call this%outputtab1%initialize_column(text, 10, alignment=tabcenter)
1771  text = 'CELLID'
1772  call this%outputtab1%initialize_column(text, 20, alignment=tableft)
1773  text = 'RATE'
1774  call this%outputtab1%initialize_column(text, 15, alignment=tabcenter)
1775  if (this%inamedbound > 0) then
1776  text = 'NAME'
1777  call this%outputtab1%initialize_column(text, 20, alignment=tableft)
1778  end if
1779  end if
1780  ! outouttab2
1781  if (this%v_model2%is_local) then
1782  call table_cr(this%outputtab2, this%name, ' ')
1783  call this%outputtab2%table_df(this%nexg, ntabcol, this%gwfmodel2%iout, &
1784  transient=.true.)
1785  text = 'NUMBER'
1786  call this%outputtab2%initialize_column(text, 10, alignment=tabcenter)
1787  text = 'CELLID'
1788  call this%outputtab2%initialize_column(text, 20, alignment=tableft)
1789  text = 'RATE'
1790  call this%outputtab2%initialize_column(text, 15, alignment=tabcenter)
1791  if (this%inamedbound > 0) then
1792  text = 'NAME'
1793  call this%outputtab2%initialize_column(text, 20, alignment=tableft)
1794  end if
1795  end if
1796  end if
Here is the call graph for this function:

◆ allocate_scalars()

subroutine gwfgwfexchangemodule::allocate_scalars ( class(gwfexchangetype this)

Allocate scalar variables

Parameters
thisGwfExchangeType

Definition at line 1644 of file exg-gwfgwf.f90.

1645  ! -- modules
1647  use constantsmodule, only: dzero
1648  ! -- dummy
1649  class(GwfExchangeType) :: this !< GwfExchangeType
1650  !
1651  call this%DisConnExchangeType%allocate_scalars()
1652  !
1653  call mem_allocate(this%icellavg, 'ICELLAVG', this%memoryPath)
1654  call mem_allocate(this%ivarcv, 'IVARCV', this%memoryPath)
1655  call mem_allocate(this%idewatcv, 'IDEWATCV', this%memoryPath)
1656  call mem_allocate(this%inewton, 'INEWTON', this%memoryPath)
1657  call mem_allocate(this%ingnc, 'INGNC', this%memoryPath)
1658  call mem_allocate(this%inmvr, 'INMVR', this%memoryPath)
1659  call mem_allocate(this%inobs, 'INOBS', this%memoryPath)
1660  call mem_allocate(this%satomega, 'SATOMEGA', this%memoryPath)
1661  this%icellavg = 0
1662  this%ivarcv = 0
1663  this%idewatcv = 0
1664  this%inewton = 0
1665  this%ingnc = 0
1666  this%inmvr = 0
1667  this%inobs = 0
1668  this%satomega = dzero
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ calc_cond_sat()

subroutine gwfgwfexchangemodule::calc_cond_sat ( class(gwfexchangetype this)
Parameters
thisGwfExchangeType

Definition at line 1462 of file exg-gwfgwf.f90.

1463  ! -- modules
1466  ! -- dummy
1467  class(GwfExchangeType) :: this !< GwfExchangeType
1468  ! -- local
1469  integer(I4B) :: iexg
1470  integer(I4B) :: n, m, ihc
1471  real(DP) :: topn, topm
1472  real(DP) :: botn, botm
1473  real(DP) :: satn, satm
1474  real(DP) :: thickn, thickm
1475  real(DP) :: angle, hyn, hym
1476  real(DP) :: csat
1477  real(DP) :: fawidth
1478  real(DP), dimension(3) :: vg
1479  !
1480  do iexg = 1, this%nexg
1481  !
1482  ihc = this%ihc(iexg)
1483  n = this%nodem1(iexg)
1484  m = this%nodem2(iexg)
1485  topn = this%gwfmodel1%dis%top(n)
1486  topm = this%gwfmodel2%dis%top(m)
1487  botn = this%gwfmodel1%dis%bot(n)
1488  botm = this%gwfmodel2%dis%bot(m)
1489  satn = this%gwfmodel1%npf%sat(n)
1490  satm = this%gwfmodel2%npf%sat(m)
1491  thickn = (topn - botn) * satn
1492  thickm = (topm - botm) * satm
1493  !
1494  ! -- Calculate conductance depending on connection orientation
1495  if (ihc == 0) then
1496  !
1497  ! -- Vertical conductance for fully saturated conditions
1498  vg(1) = dzero
1499  vg(2) = dzero
1500  vg(3) = done
1501  hyn = this%gwfmodel1%npf%hy_eff(n, 0, ihc, vg=vg)
1502  hym = this%gwfmodel2%npf%hy_eff(m, 0, ihc, vg=vg)
1503  csat = vcond(1, 1, 1, 1, 0, 1, 1, done, &
1504  botn, botm, &
1505  hyn, hym, &
1506  satn, satm, &
1507  topn, topm, &
1508  botn, botm, &
1509  this%hwva(iexg))
1510  else
1511  !
1512  ! -- Calculate horizontal conductance
1513  hyn = this%gwfmodel1%npf%k11(n)
1514  hym = this%gwfmodel2%npf%k11(m)
1515  !
1516  ! -- Check for anisotropy in models, and recalculate hyn and hym
1517  if (this%ianglex > 0) then
1518  angle = this%auxvar(this%ianglex, iexg) * dpio180
1519  vg(1) = abs(cos(angle))
1520  vg(2) = abs(sin(angle))
1521  vg(3) = dzero
1522  !
1523  ! -- anisotropy in model 1
1524  if (this%gwfmodel1%npf%ik22 /= 0) then
1525  hyn = this%gwfmodel1%npf%hy_eff(n, 0, ihc, vg=vg)
1526  end if
1527  !
1528  ! -- anisotropy in model 2
1529  if (this%gwfmodel2%npf%ik22 /= 0) then
1530  hym = this%gwfmodel2%npf%hy_eff(m, 0, ihc, vg=vg)
1531  end if
1532  end if
1533  !
1534  fawidth = this%hwva(iexg)
1535  csat = hcond(1, 1, 1, 1, 0, ihc, &
1536  this%icellavg, done, &
1537  topn, topm, satn, satm, hyn, hym, &
1538  topn, topm, &
1539  botn, botm, &
1540  this%cl1(iexg), this%cl2(iexg), &
1541  fawidth)
1542  end if
1543  !
1544  ! -- store csat in condsat
1545  this%condsat(iexg) = csat
1546  end do
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dpio180
real constant
Definition: Constants.f90:130
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module contains stateless conductance functions.
real(dp) function, public condmean(k1, k2, thick1, thick2, cl1, cl2, width, iavgmeth)
Calculate the conductance between two cells.
real(dp) function, public hcond(ibdn, ibdm, ictn, ictm, iupstream, ihc, icellavg, condsat, hn, hm, satn, satm, hkn, hkm, topn, topm, botn, botm, cln, clm, fawidth)
Horizontal conductance between two cells.
real(dp) function, public vcond(ibdn, ibdm, ictn, ictm, inewton, ivarcv, idewatcv, condsat, hn, hm, vkn, vkm, satn, satm, topn, topm, botn, botm, flowarea)
Vertical conductance between two cells.
Here is the call graph for this function:

◆ castasgwfexchange()

class(gwfexchangetype) function, pointer, public gwfgwfexchangemodule::castasgwfexchange ( class(*), intent(inout), pointer  obj)

Cast polymorphic object as exchange

Definition at line 2069 of file exg-gwfgwf.f90.

2070  implicit none
2071  ! -- dummy
2072  class(*), pointer, intent(inout) :: obj
2073  ! -- return
2074  class(GwfExchangeType), pointer :: res
2075  !
2076  res => null()
2077  if (.not. associated(obj)) return
2078  !
2079  select type (obj)
2080  class is (gwfexchangetype)
2081  res => obj
2082  end select
Here is the caller graph for this function:

◆ condcalc()

subroutine gwfgwfexchangemodule::condcalc ( class(gwfexchangetype this)

Calculate the conductance based on state

Parameters
thisGwfExchangeType

Definition at line 1553 of file exg-gwfgwf.f90.

1554  ! -- modules
1555  use constantsmodule, only: dhalf, dzero, done
1557  ! -- dummy
1558  class(GwfExchangeType) :: this !< GwfExchangeType
1559  ! -- local
1560  integer(I4B) :: iexg
1561  integer(I4B) :: n, m, ihc
1562  integer(I4B) :: ibdn, ibdm
1563  integer(I4B) :: ictn, ictm
1564  real(DP) :: topn, topm
1565  real(DP) :: botn, botm
1566  real(DP) :: satn, satm
1567  real(DP) :: hyn, hym
1568  real(DP) :: angle
1569  real(DP) :: hn, hm
1570  real(DP) :: cond
1571  real(DP) :: fawidth
1572  real(DP), dimension(3) :: vg
1573  !
1574  ! -- Calculate conductance and put into amat
1575  do iexg = 1, this%nexg
1576  ihc = this%ihc(iexg)
1577  n = this%nodem1(iexg)
1578  m = this%nodem2(iexg)
1579  ibdn = this%gwfmodel1%ibound(n)
1580  ibdm = this%gwfmodel2%ibound(m)
1581  ictn = this%gwfmodel1%npf%icelltype(n)
1582  ictm = this%gwfmodel2%npf%icelltype(m)
1583  topn = this%gwfmodel1%dis%top(n)
1584  topm = this%gwfmodel2%dis%top(m)
1585  botn = this%gwfmodel1%dis%bot(n)
1586  botm = this%gwfmodel2%dis%bot(m)
1587  satn = this%gwfmodel1%npf%sat(n)
1588  satm = this%gwfmodel2%npf%sat(m)
1589  hn = this%gwfmodel1%x(n)
1590  hm = this%gwfmodel2%x(m)
1591  !
1592  ! -- Calculate conductance depending on connection orientation
1593  if (ihc == 0) then
1594  !
1595  ! -- Vertical connection
1596  vg(1) = dzero
1597  vg(2) = dzero
1598  vg(3) = done
1599  hyn = this%gwfmodel1%npf%hy_eff(n, 0, ihc, vg=vg)
1600  hym = this%gwfmodel2%npf%hy_eff(m, 0, ihc, vg=vg)
1601  cond = vcond(ibdn, ibdm, ictn, ictm, this%inewton, this%ivarcv, &
1602  this%idewatcv, this%condsat(iexg), hn, hm, hyn, hym, &
1603  satn, satm, topn, topm, botn, botm, this%hwva(iexg))
1604  else
1605  !
1606  ! -- Horizontal Connection
1607  hyn = this%gwfmodel1%npf%k11(n)
1608  hym = this%gwfmodel2%npf%k11(m)
1609  !
1610  ! -- Check for anisotropy in models, and recalculate hyn and hym
1611  if (this%ianglex > 0) then
1612  angle = this%auxvar(this%ianglex, iexg)
1613  vg(1) = abs(cos(angle))
1614  vg(2) = abs(sin(angle))
1615  vg(3) = dzero
1616  !
1617  ! -- anisotropy in model 1
1618  if (this%gwfmodel1%npf%ik22 /= 0) then
1619  hyn = this%gwfmodel1%npf%hy_eff(n, 0, ihc, vg=vg)
1620  end if
1621  !
1622  ! -- anisotropy in model 2
1623  if (this%gwfmodel2%npf%ik22 /= 0) then
1624  hym = this%gwfmodel2%npf%hy_eff(m, 0, ihc, vg=vg)
1625  end if
1626  end if
1627  !
1628  fawidth = this%hwva(iexg)
1629  cond = hcond(ibdn, ibdm, ictn, ictm, this%inewton, &
1630  this%ihc(iexg), this%icellavg, this%condsat(iexg), &
1631  hn, hm, satn, satm, hyn, hym, topn, topm, botn, botm, &
1632  this%cl1(iexg), this%cl2(iexg), fawidth)
1633  end if
1634  !
1635  this%cond(iexg) = cond
1636  !
1637  end do
Here is the call graph for this function:

◆ getgwfexchangefromlist()

class(gwfexchangetype) function, pointer, public gwfgwfexchangemodule::getgwfexchangefromlist ( type(listtype), intent(inout)  list,
integer(i4b), intent(in)  idx 
)

Return an exchange from the list for specified index

Definition at line 2089 of file exg-gwfgwf.f90.

2090  implicit none
2091  ! -- dummy
2092  type(ListType), intent(inout) :: list
2093  integer(I4B), intent(in) :: idx
2094  ! -- return
2095  class(GwfExchangeType), pointer :: res
2096  ! -- local
2097  class(*), pointer :: obj
2098  !
2099  obj => list%GetItem(idx)
2100  res => castasgwfexchange(obj)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gwf_gwf_ac()

subroutine gwfgwfexchangemodule::gwf_gwf_ac ( class(gwfexchangetype this,
type(sparsematrix), intent(inout)  sparse 
)
private

Override parent exg_ac so that gnc can add connections here.

Parameters
thisGwfExchangeType

Definition at line 363 of file exg-gwfgwf.f90.

364  ! -- modules
365  use sparsemodule, only: sparsematrix
366  ! -- dummy
367  class(GwfExchangeType) :: this !< GwfExchangeType
368  type(sparsematrix), intent(inout) :: sparse
369  ! -- local
370  integer(I4B) :: n, iglo, jglo
371  !
372  ! -- add exchange connections
373  do n = 1, this%nexg
374  iglo = this%nodem1(n) + this%gwfmodel1%moffset
375  jglo = this%nodem2(n) + this%gwfmodel2%moffset
376  call sparse%addconnection(iglo, jglo, 1)
377  call sparse%addconnection(jglo, iglo, 1)
378  end do
379  !
380  ! -- add gnc connections
381  if (this%ingnc > 0) then
382  call this%gnc%gnc_ac(sparse)
383  end if

◆ gwf_gwf_ad()

subroutine gwfgwfexchangemodule::gwf_gwf_ad ( class(gwfexchangetype this)

Advance mover and obs

Parameters
thisGwfExchangeType

Definition at line 455 of file exg-gwfgwf.f90.

456  ! -- dummy
457  class(GwfExchangeType) :: this !< GwfExchangeType
458  !
459  ! -- Advance mover
460  if (this%inmvr > 0) call this%mvr%mvr_ad()
461  !
462  ! -- Push simulated values to preceding time step
463  call this%obs%obs_ad()

◆ gwf_gwf_add_to_flowja()

subroutine gwfgwfexchangemodule::gwf_gwf_add_to_flowja ( class(gwfexchangetype this)
Parameters
thisGwfExchangeType

Definition at line 713 of file exg-gwfgwf.f90.

714  ! -- modules
715  class(GwfExchangeType) :: this !< GwfExchangeType
716  ! -- local
717  integer(I4B) :: i
718  integer(I4B) :: n
719  integer(I4B) :: idiag
720  real(DP) :: flow
721  !
722  do i = 1, this%nexg
723  !
724  if (associated(this%gwfmodel1)) then
725  n = this%nodem1(i)
726  if (this%gwfmodel1%ibound(n) > 0) then
727  flow = this%simvals(i)
728  idiag = this%gwfmodel1%ia(n)
729  this%gwfmodel1%flowja(idiag) = this%gwfmodel1%flowja(idiag) + flow
730  end if
731  end if
732  !
733  if (associated(this%gwfmodel2)) then
734  n = this%nodem2(i)
735  if (this%gwfmodel2%ibound(n) > 0) then
736  flow = -this%simvals(i)
737  idiag = this%gwfmodel2%ia(n)
738  this%gwfmodel2%flowja(idiag) = this%gwfmodel2%flowja(idiag) + flow
739  end if
740  end if
741  !
742  end do

◆ gwf_gwf_ar()

subroutine gwfgwfexchangemodule::gwf_gwf_ar ( class(gwfexchangetype this)

Allocated and read and calculate saturated conductance

Parameters
thisGwfExchangeType

Definition at line 417 of file exg-gwfgwf.f90.

418  ! -- dummy
419  class(GwfExchangeType) :: this !< GwfExchangeType
420  !
421  ! -- If mover is active, then call ar routine
422  if (this%inmvr > 0) call this%mvr%mvr_ar()
423  !
424  ! -- Calculate the saturated conductance for all connections
425  if (.not. this%use_interface_model()) call this%calc_cond_sat()
426  !
427  ! -- Observation AR
428  call this%obs%obs_ar()

◆ gwf_gwf_bd()

subroutine gwfgwfexchangemodule::gwf_gwf_bd ( class(gwfexchangetype this,
integer(i4b), intent(inout)  icnvg,
integer(i4b), intent(in)  isuppress_output,
integer(i4b), intent(in)  isolnid 
)

Accumulate budget terms

Parameters
thisGwfExchangeType

Definition at line 858 of file exg-gwfgwf.f90.

859  ! -- modules
861  use budgetmodule, only: rate_accumulator
862  ! -- dummy
863  class(GwfExchangeType) :: this !< GwfExchangeType
864  integer(I4B), intent(inout) :: icnvg
865  integer(I4B), intent(in) :: isuppress_output
866  integer(I4B), intent(in) :: isolnid
867  ! -- local
868  character(len=LENBUDTXT), dimension(1) :: budtxt
869  real(DP), dimension(2, 1) :: budterm
870  real(DP) :: ratin, ratout
871  !
872  ! -- initialize
873  budtxt(1) = ' FLOW-JA-FACE'
874  !
875  ! -- Calculate ratin/ratout and pass to model budgets
876  call rate_accumulator(this%simvals, ratin, ratout)
877  !
878  ! -- Add the budget terms to model 1
879  if (associated(this%gwfmodel1)) then
880  budterm(1, 1) = ratin
881  budterm(2, 1) = ratout
882  call this%gwfmodel1%model_bdentry(budterm, budtxt, this%name)
883  end if
884  !
885  ! -- Add the budget terms to model 2
886  if (associated(this%gwfmodel2)) then
887  budterm(1, 1) = ratout
888  budterm(2, 1) = ratin
889  call this%gwfmodel2%model_bdentry(budterm, budtxt, this%name)
890  end if
891  !
892  ! -- Add any flows from one model into a constant head in another model
893  ! as a separate budget term called FLOW-JA-FACE-CHD
894  call this%gwf_gwf_chd_bd()
895  !
896  ! -- Call mvr bd routine
897  if (this%inmvr > 0) call this%mvr%mvr_bd()
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
Here is the call graph for this function:

◆ gwf_gwf_bdsav()

subroutine gwfgwfexchangemodule::gwf_gwf_bdsav ( class(gwfexchangetype this)

Output individual flows to listing file and binary budget files

Parameters
thisGwfExchangeType

Definition at line 967 of file exg-gwfgwf.f90.

968  ! -- dummy
969  class(GwfExchangeType) :: this !< GwfExchangeType
970  ! -- local
971  integer(I4B) :: icbcfl, ibudfl
972  !
973  ! -- budget for model1
974  if (associated(this%gwfmodel1)) then
975  call this%gwf_gwf_bdsav_model(this%gwfmodel1)
976  end if
977  !
978  ! -- budget for model2
979  if (associated(this%gwfmodel2)) then
980  call this%gwf_gwf_bdsav_model(this%gwfmodel2)
981  end if
982  !
983  ! -- Set icbcfl, ibudfl to zero so that flows will be printed and
984  ! saved, if the options were set in the MVR package
985  icbcfl = 1
986  ibudfl = 1
987  !
988  ! -- Call mvr bd routine
989  if (this%inmvr > 0) call this%mvr%mvr_bdsav(icbcfl, ibudfl, 0)
990  !
991  ! -- Calculate and write simulated values for observations
992  if (this%inobs /= 0) then
993  call this%gwf_gwf_save_simvals()
994  end if

◆ gwf_gwf_bdsav_model()

subroutine gwfgwfexchangemodule::gwf_gwf_bdsav_model ( class(gwfexchangetype this,
class(gwfmodeltype), pointer  model 
)
private
Parameters
thisthis exchange
modelthe model to save budget for

Definition at line 997 of file exg-gwfgwf.f90.

998  ! -- modules
1000  use tdismodule, only: kstp, kper
1001  ! -- dummy
1002  class(GwfExchangeType) :: this !< this exchange
1003  class(GwfModelType), pointer :: model !< the model to save budget for
1004  ! -- local
1005  character(len=LENPACKAGENAME + 4) :: packname
1006  character(len=LENBUDTXT), dimension(1) :: budtxt
1007  type(TableType), pointer :: output_tab
1008  class(VirtualModelType), pointer :: nbr_model
1009  character(len=20) :: nodestr
1010  character(len=LENBOUNDNAME) :: bname
1011  integer(I4B) :: ntabrows
1012  integer(I4B) :: nodeu
1013  integer(I4B) :: i, n1, n2, n1u, n2u
1014  integer(I4B) :: ibinun
1015  real(DP) :: ratin, ratout, rrate
1016  logical(LGP) :: is_for_model1
1017  real(DP), dimension(this%naux) :: auxrow
1018  !
1019  budtxt(1) = ' FLOW-JA-FACE'
1020  packname = 'EXG '//this%name
1021  packname = adjustr(packname)
1022  if (associated(model, this%gwfmodel1)) then
1023  output_tab => this%outputtab1
1024  nbr_model => this%v_model2
1025  is_for_model1 = .true.
1026  else
1027  output_tab => this%outputtab2
1028  nbr_model => this%v_model1
1029  is_for_model1 = .false.
1030  end if
1031  !
1032  ! -- update output tables
1033  if (this%iprflow /= 0) then
1034  !
1035  ! -- update titles
1036  if (model%oc%oc_save('BUDGET')) then
1037  call output_tab%set_title(packname)
1038  end if
1039  !
1040  ! -- set table kstp and kper
1041  call output_tab%set_kstpkper(kstp, kper)
1042  !
1043  ! -- update maxbound of tables
1044  ntabrows = 0
1045  do i = 1, this%nexg
1046  n1 = this%nodem1(i)
1047  n2 = this%nodem2(i)
1048  !
1049  ! -- If both cells are active then calculate flow rate
1050  if (this%v_model1%ibound%get(n1) /= 0 .and. &
1051  this%v_model2%ibound%get(n2) /= 0) then
1052  ntabrows = ntabrows + 1
1053  end if
1054  end do
1055  if (ntabrows > 0) then
1056  call output_tab%set_maxbound(ntabrows)
1057  end if
1058  end if
1059  !
1060  ! -- Print and write budget terms
1061  !
1062  ! -- Set binary unit numbers for saving flows
1063  if (this%ipakcb /= 0) then
1064  ibinun = model%oc%oc_save_unit('BUDGET')
1065  else
1066  ibinun = 0
1067  end if
1068  !
1069  ! -- If save budget flag is zero for this stress period, then
1070  ! shut off saving
1071  if (.not. model%oc%oc_save('BUDGET')) ibinun = 0
1072  !
1073  ! -- If cell-by-cell flows will be saved as a list, write header.
1074  if (ibinun /= 0) then
1075  call model%dis%record_srcdst_list_header(budtxt(1), &
1076  model%name, &
1077  this%name, &
1078  nbr_model%name, &
1079  this%name, &
1080  this%naux, this%auxname, &
1081  ibinun, this%nexg, &
1082  model%iout)
1083  end if
1084  !
1085  ! Initialize accumulators
1086  ratin = dzero
1087  ratout = dzero
1088  !
1089  ! -- Loop through all exchanges
1090  do i = 1, this%nexg
1091  !
1092  ! -- Assign boundary name
1093  if (this%inamedbound > 0) then
1094  bname = this%boundname(i)
1095  else
1096  bname = ''
1097  end if
1098  !
1099  ! -- Calculate the flow rate between n1 and n2
1100  rrate = dzero
1101  n1 = this%nodem1(i)
1102  n2 = this%nodem2(i)
1103  !
1104  ! -- If both cells are active then calculate flow rate
1105  if (this%v_model1%ibound%get(n1) /= 0 .and. &
1106  this%v_model2%ibound%get(n2) /= 0) then
1107  rrate = this%simvals(i)
1108  !
1109  ! -- Print the individual rates to model list files if requested
1110  if (this%iprflow /= 0) then
1111  if (model%oc%oc_save('BUDGET')) then
1112  !
1113  ! -- set nodestr and write outputtab table
1114  if (is_for_model1) then
1115  nodeu = model%dis%get_nodeuser(n1)
1116  call model%dis%nodeu_to_string(nodeu, nodestr)
1117  call output_tab%print_list_entry(i, trim(adjustl(nodestr)), &
1118  rrate, bname)
1119  else
1120  nodeu = model%dis%get_nodeuser(n2)
1121  call model%dis%nodeu_to_string(nodeu, nodestr)
1122  call output_tab%print_list_entry(i, trim(adjustl(nodestr)), &
1123  -rrate, bname)
1124  end if
1125  end if
1126  end if
1127  if (rrate < dzero) then
1128  ratout = ratout - rrate
1129  else
1130  ratin = ratin + rrate
1131  end if
1132  end if
1133  !
1134  ! -- If saving cell-by-cell flows in list, write flow
1135  n1u = this%v_model1%dis_get_nodeuser(n1)
1136  n2u = this%v_model2%dis_get_nodeuser(n2)
1137  if (ibinun /= 0) then
1138  if (this%naux > 0) then
1139  auxrow(:) = this%auxvar(:, i)
1140  end if
1141  if (is_for_model1) then
1142  call model%dis%record_mf6_list_entry(ibinun, n1u, n2u, rrate, &
1143  this%naux, auxrow, &
1144  .false., .false.)
1145  else
1146  call model%dis%record_mf6_list_entry(ibinun, n2u, n1u, -rrate, &
1147  this%naux, auxrow, &
1148  .false., .false.)
1149  end if
1150  end if
1151  !
1152  end do
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

◆ gwf_gwf_calc_simvals()

subroutine gwfgwfexchangemodule::gwf_gwf_calc_simvals ( class(gwfexchangetype this)
private
Parameters
thisGwfExchangeType

Definition at line 683 of file exg-gwfgwf.f90.

684  ! -- modules
685  use constantsmodule, only: dzero
686  ! -- dummy
687  class(GwfExchangeType) :: this !< GwfExchangeType
688  ! -- local
689  integer(I4B) :: i
690  integer(I4B) :: n1, n2
691  integer(I4B) :: ibdn1, ibdn2
692  real(DP) :: rrate
693  !
694  do i = 1, this%nexg
695  rrate = dzero
696  n1 = this%nodem1(i)
697  n2 = this%nodem2(i)
698  ibdn1 = this%gwfmodel1%ibound(n1)
699  ibdn2 = this%gwfmodel2%ibound(n2)
700  if (ibdn1 /= 0 .and. ibdn2 /= 0) then
701  rrate = this%qcalc(i, n1, n2)
702  if (this%ingnc > 0) then
703  rrate = rrate + this%gnc%deltaqgnc(i)
704  end if
705  end if
706  this%simvals(i) = rrate
707  end do

◆ gwf_gwf_cf()

subroutine gwfgwfexchangemodule::gwf_gwf_cf ( class(gwfexchangetype this,
integer(i4b), intent(in)  kiter 
)
private

Rewet as necessary

Parameters
thisGwfExchangeType

Definition at line 470 of file exg-gwfgwf.f90.

471  ! -- dummy
472  class(GwfExchangeType) :: this !< GwfExchangeType
473  integer(I4B), intent(in) :: kiter
474  ! -- local
475  !
476  ! -- Call mvr fc routine
477  if (this%inmvr > 0) call this%mvr%xmvr_cf()
478  !
479  ! -- Rewet cells across models using the wetdry parameters in each model's
480  ! npf package, and the head in the connected model.
481  call this%rewet(kiter)

◆ gwf_gwf_chd_bd()

subroutine gwfgwfexchangemodule::gwf_gwf_chd_bd ( class(gwfexchangetype this)

Account for flow from an external model into a chd cell

Parameters
thisGwfExchangeType

Definition at line 904 of file exg-gwfgwf.f90.

905  ! -- modules
907  use budgetmodule, only: rate_accumulator
908  ! -- dummy
909  class(GwfExchangeType) :: this !< GwfExchangeType
910  ! -- local
911  character(len=LENBUDTXT), dimension(1) :: budtxt
912  integer(I4B) :: n
913  integer(I4B) :: i
914  real(DP), dimension(2, 1) :: budterm
915  real(DP) :: ratin, ratout
916  real(DP) :: q
917  !
918  ! -- initialize
919  budtxt(1) = 'FLOW-JA-FACE-CHD'
920  !
921  ! -- Add the constant-head budget terms for flow from model 2 into model 1
922  if (associated(this%gwfmodel1)) then
923  ratin = dzero
924  ratout = dzero
925  do i = 1, this%nexg
926  n = this%nodem1(i)
927  if (this%gwfmodel1%ibound(n) < 0) then
928  q = this%simvals(i)
929  if (q > dzero) then
930  ratout = ratout + q
931  else
932  ratin = ratin - q
933  end if
934  end if
935  end do
936  budterm(1, 1) = ratin
937  budterm(2, 1) = ratout
938  call this%gwfmodel1%model_bdentry(budterm, budtxt, this%name)
939  end if
940  !
941  ! -- Add the constant-head budget terms for flow from model 1 into model 2
942  if (associated(this%gwfmodel2)) then
943  ratin = dzero
944  ratout = dzero
945  do i = 1, this%nexg
946  n = this%nodem2(i)
947  if (this%gwfmodel2%ibound(n) < 0) then
948  ! -- flip flow sign as flow is relative to model 1
949  q = -this%simvals(i)
950  if (q > dzero) then
951  ratout = ratout + q
952  else
953  ratin = ratin - q
954  end if
955  end if
956  end do
957  budterm(1, 1) = ratin
958  budterm(2, 1) = ratout
959  call this%gwfmodel2%model_bdentry(budterm, budtxt, this%name)
960  end if
Here is the call graph for this function:

◆ gwf_gwf_connects_model()

logical(lgp) function gwfgwfexchangemodule::gwf_gwf_connects_model ( class(gwfexchangetype this,
class(basemodeltype), intent(in), pointer  model 
)
private
Parameters
model
thisGwfExchangeType
[in]modelthe model to which the exchange might hold a connection
Returns
true, when connected

Definition at line 1939 of file exg-gwfgwf.f90.

1940  ! -- dummy
1941  class(GwfExchangeType) :: this !< GwfExchangeType
1942  class(BaseModelType), pointer, intent(in) :: model !< the model to which the exchange might hold a connection
1943  ! -- return
1944  logical(LGP) :: is_connected !< true, when connected
1945  !
1946  is_connected = .false.
1947  !
1948  ! only connected when model is GwfModelType of course
1949  select type (model)
1950  class is (gwfmodeltype)
1951  if (associated(this%gwfmodel1, model)) then
1952  is_connected = .true.
1953  else if (associated(this%gwfmodel2, model)) then
1954  is_connected = .true.
1955  end if
1956  end select

◆ gwf_gwf_cq()

subroutine gwfgwfexchangemodule::gwf_gwf_cq ( class(gwfexchangetype this,
integer(i4b), intent(inout)  icnvg,
integer(i4b), intent(in)  isuppress_output,
integer(i4b), intent(in)  isolnid 
)

Calculate flow between two cells and store in simvals, also set information needed for specific discharge calculation

Parameters
thisGwfExchangeType

Definition at line 663 of file exg-gwfgwf.f90.

664  ! -- dummy
665  class(GwfExchangeType) :: this !< GwfExchangeType
666  integer(I4B), intent(inout) :: icnvg
667  integer(I4B), intent(in) :: isuppress_output
668  integer(I4B), intent(in) :: isolnid
669  !
670  ! -- calculate flow and store in simvals
671  call this%gwf_gwf_calc_simvals()
672  !
673  ! -- set flows to model edges in NPF
674  call this%gwf_gwf_set_flow_to_npf()
675  !
676  ! -- add exchange flows to model's flowja diagonal
677  call this%gwf_gwf_add_to_flowja()

◆ gwf_gwf_da()

subroutine gwfgwfexchangemodule::gwf_gwf_da ( class(gwfexchangetype this)

Deallocate memory associated with this object

Parameters
thisGwfExchangeType

Definition at line 1675 of file exg-gwfgwf.f90.

1676  ! -- modules
1678  ! -- dummy
1679  class(GwfExchangeType) :: this !< GwfExchangeType
1680  !
1681  ! -- objects
1682  if (this%ingnc > 0) then
1683  call this%gnc%gnc_da()
1684  deallocate (this%gnc)
1685  end if
1686  if (this%inmvr > 0) then
1687  call this%mvr%mvr_da()
1688  deallocate (this%mvr)
1689  end if
1690  call this%obs%obs_da()
1691  deallocate (this%obs)
1692  !
1693  ! -- arrays
1694  call mem_deallocate(this%cond)
1695  call mem_deallocate(this%condsat)
1696  call mem_deallocate(this%idxglo)
1697  call mem_deallocate(this%idxsymglo)
1698  call mem_deallocate(this%simvals, 'SIMVALS', this%memoryPath)
1699  !
1700  ! -- output table objects
1701  if (associated(this%outputtab1)) then
1702  call this%outputtab1%table_da()
1703  deallocate (this%outputtab1)
1704  nullify (this%outputtab1)
1705  end if
1706  if (associated(this%outputtab2)) then
1707  call this%outputtab2%table_da()
1708  deallocate (this%outputtab2)
1709  nullify (this%outputtab2)
1710  end if
1711  !
1712  ! -- scalars
1713  deallocate (this%filename)
1714  !
1715  call mem_deallocate(this%icellavg)
1716  call mem_deallocate(this%ivarcv)
1717  call mem_deallocate(this%idewatcv)
1718  call mem_deallocate(this%inewton)
1719  call mem_deallocate(this%ingnc)
1720  call mem_deallocate(this%inmvr)
1721  call mem_deallocate(this%inobs)
1722  call mem_deallocate(this%satomega)
1723  !
1724  ! -- deallocate base
1725  call this%DisConnExchangeType%disconnex_da()

◆ gwf_gwf_df()

subroutine gwfgwfexchangemodule::gwf_gwf_df ( class(gwfexchangetype this)

Define GWF to GWF exchange object.

Parameters
thisGwfExchangeType

Definition at line 206 of file exg-gwfgwf.f90.

207  ! -- modules
208  use simvariablesmodule, only: iout
210  use ghostnodemodule, only: gnc_cr
211  ! -- dummy
212  class(GwfExchangeType) :: this !< GwfExchangeType
213  ! -- local
214  !
215  ! -- log the exchange
216  write (iout, '(/a,a)') ' Creating exchange: ', this%name
217  !
218  ! -- Ensure models are in same solution
219  if (this%v_model1%idsoln%get() /= this%v_model2%idsoln%get()) then
220  call store_error('Two models are connected in a GWF '// &
221  'exchange but they are in different solutions. '// &
222  'GWF models must be in same solution: '// &
223  trim(this%v_model1%name)//' '// &
224  trim(this%v_model2%name))
225  call store_error_filename(this%filename)
226  end if
227  !
228  ! -- source options
229  call this%source_options(iout)
230  !
231  ! -- source dimensions
232  call this%source_dimensions(iout)
233  !
234  ! -- allocate arrays
235  call this%allocate_arrays()
236  !
237  ! -- source exchange data
238  call this%source_data(iout)
239  !
240  ! -- call each model and increase the edge count
241  if (associated(this%gwfmodel1)) then
242  call this%gwfmodel1%npf%increase_edge_count(this%nexg)
243  end if
244  if (associated(this%gwfmodel2)) then
245  call this%gwfmodel2%npf%increase_edge_count(this%nexg)
246  end if
247  !
248  ! -- Create and read ghost node information
249  if (this%ingnc > 0) then
250  if (associated(this%gwfmodel1) .and. associated(this%gwfmodel2)) then
251  call gnc_cr(this%gnc, this%name, this%ingnc, iout)
252  call this%read_gnc()
253  end if
254  end if
255  !
256  ! -- Read mover information
257  if (this%inmvr > 0) then
258  call this%read_mvr(iout)
259  end if
260  !
261  ! -- Store obs
262  call this%gwf_gwf_df_obs()
263  if (associated(this%gwfmodel1)) then
264  call this%obs%obs_df(iout, this%name, 'GWF-GWF', this%gwfmodel1%dis)
265  end if
266  !
267  ! -- validate
268  call this%validate_exchange()
subroutine, public gnc_cr(gncobj, name_parent, inunit, iout)
Create new GNC exchange object.
Definition: GhostNode.f90:61
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) iout
file unit number for simulation output
Here is the call graph for this function:

◆ gwf_gwf_df_obs()

subroutine gwfgwfexchangemodule::gwf_gwf_df_obs ( class(gwfexchangetype this)

Define the observations associated with this object

Parameters
thisGwfExchangeType

Definition at line 1803 of file exg-gwfgwf.f90.

1804  ! -- dummy
1805  class(GwfExchangeType) :: this !< GwfExchangeType
1806  ! -- local
1807  integer(I4B) :: indx
1808  !
1809  ! -- Store obs type and assign procedure pointer
1810  ! for gwf-gwf observation type.
1811  call this%obs%StoreObsType('flow-ja-face', .true., indx)
1812  this%obs%obsData(indx)%ProcessIdPtr => gwf_gwf_process_obsid
Here is the call graph for this function:

◆ gwf_gwf_fc()

subroutine gwfgwfexchangemodule::gwf_gwf_fc ( class(gwfexchangetype this,
integer(i4b), intent(in)  kiter,
class(matrixbasetype), pointer  matrix_sln,
real(dp), dimension(:), intent(inout)  rhs_sln,
integer(i4b), intent(in), optional  inwtflag 
)
private

Calculate conductance and fill coefficient matrix

Parameters
thisGwfExchangeType

Definition at line 488 of file exg-gwfgwf.f90.

489  ! -- modules
490  use constantsmodule, only: dhalf
492  ! -- dummy
493  class(GwfExchangeType) :: this !< GwfExchangeType
494  integer(I4B), intent(in) :: kiter
495  class(MatrixBaseType), pointer :: matrix_sln
496  real(DP), dimension(:), intent(inout) :: rhs_sln
497  integer(I4B), optional, intent(in) :: inwtflag
498  ! -- local
499  integer(I4B) :: inwt, iexg
500  integer(I4B) :: i, nodem1sln, nodem2sln
501  !
502  ! -- calculate the conductance for each exchange connection
503  call this%condcalc()
504  !
505  ! -- if gnc is active, then copy cond into gnc cond (might consider a
506  ! pointer here in the future)
507  if (this%ingnc > 0) then
508  do iexg = 1, this%nexg
509  this%gnc%cond(iexg) = this%cond(iexg)
510  end do
511  end if
512  !
513  ! -- Put this%cond into amatsln
514  do i = 1, this%nexg
515  call matrix_sln%set_value_pos(this%idxglo(i), this%cond(i))
516  call matrix_sln%set_value_pos(this%idxsymglo(i), this%cond(i))
517 
518  nodem1sln = this%nodem1(i) + this%gwfmodel1%moffset
519  nodem2sln = this%nodem2(i) + this%gwfmodel2%moffset
520  call matrix_sln%add_diag_value(nodem1sln, -this%cond(i))
521  call matrix_sln%add_diag_value(nodem2sln, -this%cond(i))
522  end do
523  !
524  ! -- Fill the gnc terms in the solution matrix
525  if (this%ingnc > 0) then
526  call this%gnc%gnc_fc(kiter, matrix_sln)
527  end if
528  !
529  ! -- Call mvr fc routine
530  if (this%inmvr > 0) call this%mvr%mvr_fc()
531  !
532  ! -- Set inwt to exchange newton, but shut off if requested by caller
533  inwt = this%inewton
534  if (present(inwtflag)) then
535  if (inwtflag == 0) inwt = 0
536  end if
537  if (inwt /= 0) then
538  call this%exg_fn(kiter, matrix_sln)
539  end if
540  !
541  ! -- Ghost node Newton-Raphson
542  if (this%ingnc > 0) then
543  if (inwt /= 0) then
544  call this%gnc%gnc_fn(kiter, matrix_sln, this%condsat, &
545  ihc_opt=this%ihc, ivarcv_opt=this%ivarcv, &
546  ictm1_opt=this%gwfmodel1%npf%icelltype, &
547  ictm2_opt=this%gwfmodel2%npf%icelltype)
548  end if
549  end if
Here is the call graph for this function:

◆ gwf_gwf_fn()

subroutine gwfgwfexchangemodule::gwf_gwf_fn ( class(gwfexchangetype this,
integer(i4b), intent(in)  kiter,
class(matrixbasetype), pointer  matrix_sln 
)

Fill amatsln with Newton terms

Parameters
thisGwfExchangeType

Definition at line 556 of file exg-gwfgwf.f90.

557  ! -- modules
559  ! -- dummy
560  class(GwfExchangeType) :: this !< GwfExchangeType
561  integer(I4B), intent(in) :: kiter
562  class(MatrixBaseType), pointer :: matrix_sln
563  ! -- local
564  logical :: nisup
565  integer(I4B) :: iexg
566  integer(I4B) :: n, m
567  integer(I4B) :: nodensln, nodemsln
568  integer(I4B) :: ibdn, ibdm
569  real(DP) :: topn, topm
570  real(DP) :: botn, botm
571  real(DP) :: topup, botup
572  real(DP) :: hn, hm
573  real(DP) :: hup, hdn
574  real(DP) :: cond
575  real(DP) :: term
576  real(DP) :: consterm
577  real(DP) :: derv
578  !
579  do iexg = 1, this%nexg
580  n = this%nodem1(iexg)
581  m = this%nodem2(iexg)
582  nodensln = this%nodem1(iexg) + this%gwfmodel1%moffset
583  nodemsln = this%nodem2(iexg) + this%gwfmodel2%moffset
584  ibdn = this%gwfmodel1%ibound(n)
585  ibdm = this%gwfmodel2%ibound(m)
586  topn = this%gwfmodel1%dis%top(n)
587  topm = this%gwfmodel2%dis%top(m)
588  botn = this%gwfmodel1%dis%bot(n)
589  botm = this%gwfmodel2%dis%bot(m)
590  hn = this%gwfmodel1%x(n)
591  hm = this%gwfmodel2%x(m)
592  if (this%ihc(iexg) == 0) then
593  ! -- vertical connection, newton not supported
594  else
595  ! -- determine upstream node
596  nisup = .false.
597  if (hm < hn) nisup = .true.
598  !
599  ! -- set upstream top and bot
600  if (nisup) then
601  topup = topn
602  botup = botn
603  hup = hn
604  hdn = hm
605  else
606  topup = topm
607  botup = botm
608  hup = hm
609  hdn = hn
610  end if
611  !
612  ! -- no newton terms if upstream cell is confined
613  if (nisup) then
614  if (this%gwfmodel1%npf%icelltype(n) == 0) cycle
615  else
616  if (this%gwfmodel2%npf%icelltype(m) == 0) cycle
617  end if
618  !
619  ! -- set topup and botup
620  if (this%ihc(iexg) == 2) then
621  topup = min(topn, topm)
622  botup = max(botn, botm)
623  end if
624  !
625  ! get saturated conductivity for derivative
626  cond = this%condsat(iexg)
627  !
628  ! -- TO DO deal with MODFLOW-NWT upstream weighting option
629  !
630  ! -- compute terms
631  consterm = -cond * (hup - hdn)
632  derv = squadraticsaturationderivative(topup, botup, hup)
633  if (nisup) then
634  !
635  ! -- fill jacobian with n being upstream
636  term = consterm * derv
637  this%gwfmodel1%rhs(n) = this%gwfmodel1%rhs(n) + term * hn
638  this%gwfmodel2%rhs(m) = this%gwfmodel2%rhs(m) - term * hn
639  call matrix_sln%add_diag_value(nodensln, term)
640  if (ibdm > 0) then
641  call matrix_sln%add_value_pos(this%idxsymglo(iexg), -term)
642  end if
643  else
644  !
645  ! -- fill jacobian with m being upstream
646  term = -consterm * derv
647  this%gwfmodel1%rhs(n) = this%gwfmodel1%rhs(n) + term * hm
648  this%gwfmodel2%rhs(m) = this%gwfmodel2%rhs(m) - term * hm
649  call matrix_sln%add_diag_value(nodemsln, -term)
650  if (ibdn > 0) then
651  call matrix_sln%add_value_pos(this%idxglo(iexg), term)
652  end if
653  end if
654  end if
655  end do
real(dp) function squadraticsaturationderivative(top, bot, x, eps)
@ brief Derivative of the quadratic saturation function
Here is the call graph for this function:

◆ gwf_gwf_fp()

subroutine gwfgwfexchangemodule::gwf_gwf_fp ( class(gwfexchangetype this)

Conduct any final processing

Parameters
thisGwfExchangeType

Definition at line 1890 of file exg-gwfgwf.f90.

1891  ! -- dummy
1892  class(GwfExchangeType) :: this !< GwfExchangeType

◆ gwf_gwf_get_iasym()

integer(i4b) function gwfgwfexchangemodule::gwf_gwf_get_iasym ( class(gwfexchangetype this)
private

Return flag indicating whether or not this exchange will cause the coefficient matrix to be asymmetric.

Parameters
thisGwfExchangeType

Definition at line 1918 of file exg-gwfgwf.f90.

1919  ! -- dummy
1920  class(GwfExchangeType) :: this !< GwfExchangeType
1921  ! -- local
1922  integer(I4B) :: iasym
1923  !
1924  ! -- Start by setting iasym to zero
1925  iasym = 0
1926  !
1927  ! -- Groundwater flow
1928  if (this%inewton /= 0) iasym = 1
1929  !
1930  ! -- GNC
1931  if (this%ingnc > 0) then
1932  if (this%gnc%iasym /= 0) iasym = 1
1933  end if

◆ gwf_gwf_mc()

subroutine gwfgwfexchangemodule::gwf_gwf_mc ( class(gwfexchangetype this,
class(matrixbasetype), pointer  matrix_sln 
)

Map the connections in the global matrix

Parameters
thisGwfExchangeType
matrix_slnthe system matrix

Definition at line 390 of file exg-gwfgwf.f90.

391  ! -- modules
392  use sparsemodule, only: sparsematrix
393  ! -- dummy
394  class(GwfExchangeType) :: this !< GwfExchangeType
395  class(MatrixBaseType), pointer :: matrix_sln !< the system matrix
396  ! -- local
397  integer(I4B) :: n, iglo, jglo
398  !
399  ! -- map exchange connections
400  do n = 1, this%nexg
401  iglo = this%nodem1(n) + this%gwfmodel1%moffset
402  jglo = this%nodem2(n) + this%gwfmodel2%moffset
403  this%idxglo(n) = matrix_sln%get_position(iglo, jglo)
404  this%idxsymglo(n) = matrix_sln%get_position(jglo, iglo)
405  end do
406  !
407  ! -- map gnc connections
408  if (this%ingnc > 0) then
409  call this%gnc%gnc_mc(matrix_sln)
410  end if

◆ gwf_gwf_ot()

subroutine gwfgwfexchangemodule::gwf_gwf_ot ( class(gwfexchangetype this)

Write output

Parameters
thisGwfExchangeType

Definition at line 1159 of file exg-gwfgwf.f90.

1160  ! -- modules
1161  use simvariablesmodule, only: iout
1162  use constantsmodule, only: dzero
1163  ! -- dummy
1164  class(GwfExchangeType) :: this !< GwfExchangeType
1165  ! -- local
1166  integer(I4B) :: iexg, n1, n2
1167  integer(I4B) :: ibudfl
1168  real(DP) :: flow, deltaqgnc
1169  character(len=LINELENGTH) :: node1str, node2str
1170  ! -- format
1171  character(len=*), parameter :: fmtheader = &
1172  "(/1x, 'SUMMARY OF EXCHANGE RATES FOR EXCHANGE ', a, ' WITH ID ', i0, /, &
1173  &2a16, 5a16, /, 112('-'))"
1174  character(len=*), parameter :: fmtheader2 = &
1175  "(/1x, 'SUMMARY OF EXCHANGE RATES FOR EXCHANGE ', a, ' WITH ID ', i0, /, &
1176  &2a16, 4a16, /, 96('-'))"
1177  character(len=*), parameter :: fmtdata = &
1178  "(2a16, 5(1pg16.6))"
1179  !
1180  ! -- Call bdsave
1181  call this%gwf_gwf_bdsav()
1182  !
1183  ! -- Initialize
1184  deltaqgnc = dzero
1185  !
1186  ! -- Write a table of exchanges
1187  if (this%iprflow /= 0) then
1188  if (this%ingnc > 0) then
1189  write (iout, fmtheader) trim(adjustl(this%name)), this%id, 'NODEM1', &
1190  'NODEM2', 'COND', 'X_M1', 'X_M2', 'DELTAQGNC', &
1191  'FLOW'
1192  else
1193  write (iout, fmtheader2) trim(adjustl(this%name)), this%id, 'NODEM1', &
1194  'NODEM2', 'COND', 'X_M1', 'X_M2', 'FLOW'
1195  end if
1196  do iexg = 1, this%nexg
1197  n1 = this%nodem1(iexg)
1198  n2 = this%nodem2(iexg)
1199  flow = this%simvals(iexg)
1200  call this%v_model1%dis_noder_to_string(n1, node1str)
1201  call this%v_model2%dis_noder_to_string(n2, node2str)
1202 
1203  if (this%ingnc > 0) then
1204  deltaqgnc = this%gnc%deltaqgnc(iexg)
1205  write (iout, fmtdata) trim(adjustl(node1str)), &
1206  trim(adjustl(node2str)), &
1207  this%cond(iexg), this%v_model1%x%get(n1), &
1208  this%v_model2%x%get(n2), deltaqgnc, flow
1209  else
1210  write (iout, fmtdata) trim(adjustl(node1str)), &
1211  trim(adjustl(node2str)), &
1212  this%cond(iexg), this%v_model1%x%get(n1), &
1213  this%v_model2%x%get(n2), flow
1214  end if
1215  end do
1216  end if
1217  !
1218  ! -- Mover budget output
1219  ibudfl = 1
1220  if (this%inmvr > 0) call this%mvr%mvr_ot_bdsummary(ibudfl)
1221  !
1222  ! -- OBS output
1223  call this%obs%obs_ot()

◆ gwf_gwf_process_obsid()

subroutine gwfgwfexchangemodule::gwf_gwf_process_obsid ( type(observetype), intent(inout)  obsrv,
class(disbasetype), intent(in)  dis,
integer(i4b), intent(in)  inunitobs,
integer(i4b), intent(in)  iout 
)

Process observations for this exchange

Definition at line 2030 of file exg-gwfgwf.f90.

2031  ! -- modules
2032  use constantsmodule, only: linelength
2033  use inputoutputmodule, only: urword
2034  use observemodule, only: observetype
2035  use basedismodule, only: disbasetype
2036  ! -- dummy
2037  type(ObserveType), intent(inout) :: obsrv
2038  class(DisBaseType), intent(in) :: dis
2039  integer(I4B), intent(in) :: inunitobs
2040  integer(I4B), intent(in) :: iout
2041  ! -- local
2042  integer(I4B) :: n, iexg, istat
2043  integer(I4B) :: icol, istart, istop
2044  real(DP) :: r
2045  character(len=LINELENGTH) :: string
2046  !
2047  string = obsrv%IDstring
2048  icol = 1
2049  ! -- get exchange index
2050  call urword(string, icol, istart, istop, 1, n, r, iout, inunitobs)
2051  read (string(istart:istop), '(i10)', iostat=istat) iexg
2052  if (istat == 0) then
2053  obsrv%intPak1 = iexg
2054  else
2055  ! Integer can't be read from string; it's presumed to be an exchange
2056  ! boundary name (already converted to uppercase)
2057  obsrv%FeatureName = trim(adjustl(string))
2058  ! -- Observation may require summing rates from multiple exchange
2059  ! boundaries, so assign intPak1 as a value that indicates observation
2060  ! is for a named exchange boundary or group of exchange boundaries.
2061  obsrv%intPak1 = namedboundflag
2062  end if
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gwf_gwf_rp()

subroutine gwfgwfexchangemodule::gwf_gwf_rp ( class(gwfexchangetype this)
private

Read new data for mover and obs

Parameters
thisGwfExchangeType

Definition at line 435 of file exg-gwfgwf.f90.

436  ! -- modules
437  use tdismodule, only: readnewdata
438  ! -- dummy
439  class(GwfExchangeType) :: this !< GwfExchangeType
440  !
441  ! -- Check with TDIS on whether or not it is time to RP
442  if (.not. readnewdata) return
443  !
444  ! -- Read and prepare for mover
445  if (this%inmvr > 0) call this%mvr%mvr_rp()
446  !
447  ! -- Read and prepare for observations
448  call this%gwf_gwf_rp_obs()
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:26

◆ gwf_gwf_rp_obs()

subroutine gwfgwfexchangemodule::gwf_gwf_rp_obs ( class(gwfexchangetype this)
private

Handle observation exchanges exchange-boundary names.

Parameters
thisGwfExchangeType

Definition at line 1819 of file exg-gwfgwf.f90.

1820  ! -- modules
1821  use constantsmodule, only: dzero
1822  ! -- dummy
1823  class(GwfExchangeType) :: this !< GwfExchangeType
1824  ! -- local
1825  integer(I4B) :: i
1826  integer(I4B) :: j
1827  class(ObserveType), pointer :: obsrv => null()
1828  character(len=LENBOUNDNAME) :: bname
1829  logical :: jfound
1830  ! -- formats
1831 10 format('Exchange "', a, '" for observation "', a, &
1832  '" is invalid in package "', a, '"')
1833 20 format('Exchange id "', i0, '" for observation "', a, &
1834  '" is invalid in package "', a, '"')
1835  !
1836  do i = 1, this%obs%npakobs
1837  obsrv => this%obs%pakobs(i)%obsrv
1838  !
1839  ! -- indxbnds needs to be reset each stress period because
1840  ! list of boundaries can change each stress period.
1841  ! -- Not true for exchanges, but leave this in for now anyway.
1842  call obsrv%ResetObsIndex()
1843  obsrv%BndFound = .false.
1844  !
1845  bname = obsrv%FeatureName
1846  if (bname /= '') then
1847  ! -- Observation location(s) is(are) based on a boundary name.
1848  ! Iterate through all boundaries to identify and store
1849  ! corresponding index(indices) in bound array.
1850  jfound = .false.
1851  do j = 1, this%nexg
1852  if (this%boundname(j) == bname) then
1853  jfound = .true.
1854  obsrv%BndFound = .true.
1855  obsrv%CurrentTimeStepEndValue = dzero
1856  call obsrv%AddObsIndex(j)
1857  end if
1858  end do
1859  if (.not. jfound) then
1860  write (errmsg, 10) trim(bname), trim(obsrv%ObsTypeId), trim(this%name)
1861  call store_error(errmsg)
1862  end if
1863  else
1864  ! -- Observation location is a single exchange number
1865  if (obsrv%intPak1 <= this%nexg .and. obsrv%intPak1 > 0) then
1866  jfound = .true.
1867  obsrv%BndFound = .true.
1868  obsrv%CurrentTimeStepEndValue = dzero
1869  call obsrv%AddObsIndex(obsrv%intPak1)
1870  else
1871  jfound = .false.
1872  end if
1873  if (.not. jfound) then
1874  write (errmsg, 20) obsrv%intPak1, trim(obsrv%ObsTypeId), trim(this%name)
1875  call store_error(errmsg)
1876  end if
1877  end if
1878  end do
1879  !
1880  ! -- write summary of error messages
1881  if (count_errors() > 0) then
1882  call store_error_filename(this%obs%inputFilename)
1883  end if
Here is the call graph for this function:

◆ gwf_gwf_save_simvals()

subroutine gwfgwfexchangemodule::gwf_gwf_save_simvals ( class(gwfexchangetype), intent(inout)  this)

Save the simulated flows for each exchange

Definition at line 1985 of file exg-gwfgwf.f90.

1986  ! -- modules
1987  use simvariablesmodule, only: errmsg
1988  use constantsmodule, only: dzero
1989  use observemodule, only: observetype
1990  ! -- dummy
1991  class(GwfExchangeType), intent(inout) :: this
1992  ! -- local
1993  integer(I4B) :: i
1994  integer(I4B) :: j
1995  integer(I4B) :: n1
1996  integer(I4B) :: n2
1997  integer(I4B) :: iexg
1998  real(DP) :: v
1999  type(ObserveType), pointer :: obsrv => null()
2000  !
2001  ! -- Write simulated values for all gwf-gwf observations
2002  if (this%obs%npakobs > 0) then
2003  call this%obs%obs_bd_clear()
2004  do i = 1, this%obs%npakobs
2005  obsrv => this%obs%pakobs(i)%obsrv
2006  do j = 1, obsrv%indxbnds_count
2007  iexg = obsrv%indxbnds(j)
2008  v = dzero
2009  select case (obsrv%ObsTypeId)
2010  case ('FLOW-JA-FACE')
2011  n1 = this%nodem1(iexg)
2012  n2 = this%nodem2(iexg)
2013  v = this%simvals(iexg)
2014  case default
2015  errmsg = 'Unrecognized observation type: '// &
2016  trim(obsrv%ObsTypeId)
2017  call store_error(errmsg)
2018  call store_error_filename(this%obs%inputFilename)
2019  end select
2020  call this%obs%SaveOneSimval(obsrv, v)
2021  end do
2022  end do
2023  end if
character(len=maxcharlen) errmsg
error message string
Here is the call graph for this function:

◆ gwf_gwf_set_flow_to_npf()

subroutine gwfgwfexchangemodule::gwf_gwf_set_flow_to_npf ( class(gwfexchangetype this)
private
Parameters
thisGwfExchangeType

Definition at line 747 of file exg-gwfgwf.f90.

748  ! -- modules
749  use constantsmodule, only: dzero, dpio180
751  ! -- dummy
752  class(GwfExchangeType) :: this !< GwfExchangeType
753  ! -- local
754  integer(I4B) :: i
755  integer(I4B) :: n1, n2
756  integer(I4B) :: ibdn1, ibdn2
757  integer(I4B) :: ictn1, ictn2
758  integer(I4B) :: ihc
759  real(DP) :: rrate
760  real(DP) :: topn1, topn2
761  real(DP) :: botn1, botn2
762  real(DP) :: satn1, satn2
763  real(DP) :: hn1, hn2
764  real(DP) :: nx, ny
765  real(DP) :: distance
766  real(DP) :: dltot
767  real(DP) :: hwva
768  real(DP) :: area
769  real(DP) :: thksat
770  real(DP) :: angle
771  !
772  ! -- Return if there neither model needs to calculate specific discharge
773  if (this%gwfmodel1%npf%icalcspdis == 0 .and. &
774  this%gwfmodel2%npf%icalcspdis == 0) return
775  !
776  ! -- Loop through all exchanges using the flow rate
777  ! stored in simvals
778  do i = 1, this%nexg
779  rrate = this%simvals(i)
780  n1 = this%nodem1(i)
781  n2 = this%nodem2(i)
782  ihc = this%ihc(i)
783  hwva = this%hwva(i)
784  ibdn1 = this%gwfmodel1%ibound(n1)
785  ibdn2 = this%gwfmodel2%ibound(n2)
786  ictn1 = this%gwfmodel1%npf%icelltype(n1)
787  ictn2 = this%gwfmodel2%npf%icelltype(n2)
788  topn1 = this%gwfmodel1%dis%top(n1)
789  topn2 = this%gwfmodel2%dis%top(n2)
790  botn1 = this%gwfmodel1%dis%bot(n1)
791  botn2 = this%gwfmodel2%dis%bot(n2)
792  satn1 = this%gwfmodel1%npf%sat(n1)
793  satn2 = this%gwfmodel2%npf%sat(n2)
794  hn1 = this%gwfmodel1%x(n1)
795  hn2 = this%gwfmodel2%x(n2)
796  !
797  ! -- Calculate face normal components
798  if (ihc == 0) then
799  nx = dzero
800  ny = dzero
801  area = hwva
802  if (botn1 < botn2) then
803  ! -- n1 is beneath n2, so rate is positive downward. Flip rate
804  ! upward so that points in positive z direction
805  rrate = -rrate
806  end if
807  else
808  if (this%ianglex > 0) then
809  angle = this%auxvar(this%ianglex, i) * dpio180
810  nx = cos(angle)
811  ny = sin(angle)
812  else
813  ! error?
814  call store_error('error in gwf_gwf_cq', terminate=.true.)
815  end if
816  !
817  ! -- Calculate the saturated thickness at interface between n1 and n2
818  thksat = thksatnm(ibdn1, ibdn2, ictn1, ictn2, this%inewton, ihc, &
819  hn1, hn2, satn1, satn2, &
820  topn1, topn2, botn1, botn2)
821  area = hwva * thksat
822  end if
823  !
824  ! -- Submit this connection and flow information to the npf
825  ! package of gwfmodel1
826  if (this%icdist > 0) then
827  dltot = this%auxvar(this%icdist, i)
828  else
829  call store_error('error in gwf_gwf_cq', terminate=.true.)
830  end if
831  distance = dltot * this%cl1(i) / (this%cl1(i) + this%cl2(i))
832  if (this%gwfmodel1%npf%icalcspdis == 1) then
833  call this%gwfmodel1%npf%set_edge_properties(n1, ihc, rrate, area, &
834  nx, ny, distance)
835  end if
836  !
837  ! -- Submit this connection and flow information to the npf
838  ! package of gwfmodel2
839  if (this%icdist > 0) then
840  dltot = this%auxvar(this%icdist, i)
841  else
842  call store_error('error in gwf_gwf_cq', terminate=.true.)
843  end if
844  if (this%gwfmodel2%npf%icalcspdis == 1) then
845  distance = dltot * this%cl2(i) / (this%cl1(i) + this%cl2(i))
846  if (ihc /= 0) rrate = -rrate
847  call this%gwfmodel2%npf%set_edge_properties(n2, ihc, rrate, area, &
848  -nx, -ny, distance)
849  end if
850  !
851  end do
real(dp) function, public thksatnm(ibdn, ibdm, ictn, ictm, iupstream, ihc, hn, hm, satn, satm, topn, topm, botn, botm)
Calculate wetted cell thickness at interface between two cells.
Here is the call graph for this function:

◆ gwfexchange_create()

subroutine, public gwfgwfexchangemodule::gwfexchange_create ( character(len=*), intent(in)  filename,
character(len=*)  name,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  m1_id,
integer(i4b), intent(in)  m2_id,
character(len=*), intent(in)  input_mempath 
)

Create a new GWF to GWF exchange object.

Parameters
[in]filenamefilename for reading
nameexchange name
[in]idid for the exchange
[in]m1_idid for model 1
[in]m2_idid for model 2

Definition at line 121 of file exg-gwfgwf.f90.

122  ! -- modules
123  use basemodelmodule, only: basemodeltype
125  use listsmodule, only: baseexchangelist
126  use obsmodule, only: obs_cr
128  ! -- dummy
129  character(len=*), intent(in) :: filename !< filename for reading
130  character(len=*) :: name !< exchange name
131  integer(I4B), intent(in) :: id !< id for the exchange
132  integer(I4B), intent(in) :: m1_id !< id for model 1
133  integer(I4B), intent(in) :: m2_id !< id for model 2
134  character(len=*), intent(in) :: input_mempath
135  ! -- local
136  type(GwfExchangeType), pointer :: exchange
137  class(BaseModelType), pointer :: mb
138  class(BaseExchangeType), pointer :: baseexchange
139  integer(I4B) :: m1_index, m2_index
140  !
141  ! -- Create a new exchange and add it to the baseexchangelist container
142  allocate (exchange)
143  baseexchange => exchange
144  call addbaseexchangetolist(baseexchangelist, baseexchange)
145  !
146  ! -- Assign id and name
147  exchange%id = id
148  exchange%name = name
149  exchange%memoryPath = create_mem_path(exchange%name)
150  exchange%input_mempath = input_mempath
151  !
152  ! -- allocate scalars and set defaults
153  call exchange%allocate_scalars()
154  exchange%filename = filename
155  exchange%typename = 'GWF-GWF'
156  !
157  ! -- set gwfmodel1
158  m1_index = model_loc_idx(m1_id)
159  if (m1_index > 0) then
160  mb => getbasemodelfromlist(basemodellist, m1_index)
161  select type (mb)
162  type is (gwfmodeltype)
163  exchange%model1 => mb
164  exchange%gwfmodel1 => mb
165  end select
166  end if
167  exchange%v_model1 => get_virtual_model(m1_id)
168  exchange%is_datacopy = .not. exchange%v_model1%is_local
169  !
170  ! -- set gwfmodel2
171  m2_index = model_loc_idx(m2_id)
172  if (m2_index > 0) then
173  mb => getbasemodelfromlist(basemodellist, m2_index)
174  select type (mb)
175  type is (gwfmodeltype)
176  exchange%model2 => mb
177  exchange%gwfmodel2 => mb
178  end select
179  end if
180  exchange%v_model2 => get_virtual_model(m2_id)
181  !
182  ! -- Verify that gwf model1 is of the correct type
183  if (.not. associated(exchange%gwfmodel1) .and. m1_index > 0) then
184  write (errmsg, '(3a)') 'Problem with GWF-GWF exchange ', &
185  trim(exchange%name), &
186  '. First specified GWF Model does not appear to be of the correct type.'
187  call store_error(errmsg, terminate=.true.)
188  end if
189  !
190  ! -- Verify that gwf model2 is of the correct type
191  if (.not. associated(exchange%gwfmodel2) .and. m2_index > 0) then
192  write (errmsg, '(3a)') 'Problem with GWF-GWF exchange ', &
193  trim(exchange%name), &
194  '. Second specified GWF Model does not appear to be of the correct type.'
195  call store_error(errmsg, terminate=.true.)
196  end if
197  !
198  ! -- Create the obs package
199  call obs_cr(exchange%obs, exchange%inobs)
type(listtype), public baseexchangelist
Definition: mf6lists.f90:25
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public obs_cr(obs, inobs)
@ brief Create a new ObsType object
Definition: Obs.f90:225
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
Here is the call graph for this function:
Here is the caller graph for this function:

◆ qcalc()

real(dp) function gwfgwfexchangemodule::qcalc ( class(gwfexchangetype this,
integer(i4b), intent(in)  iexg,
integer(i4b), intent(in)  n1,
integer(i4b), intent(in)  n2 
)
private

Calculate the flow for the specified exchange and node numbers

Parameters
thisGwfExchangeType

Definition at line 1899 of file exg-gwfgwf.f90.

1900  ! -- return
1901  real(DP) :: qcalc
1902  ! -- dummy
1903  class(GwfExchangeType) :: this !< GwfExchangeType
1904  integer(I4B), intent(in) :: iexg
1905  integer(I4B), intent(in) :: n1
1906  integer(I4B), intent(in) :: n2
1907  ! -- local
1908  !
1909  ! -- Calculate flow between nodes in the two models
1910  qcalc = this%cond(iexg) * (this%gwfmodel2%x(n2) - this%gwfmodel1%x(n1))

◆ read_gnc()

subroutine gwfgwfexchangemodule::read_gnc ( class(gwfexchangetype this)

Read and process ghost nodes

Parameters
thisGwfExchangeType

Definition at line 1330 of file exg-gwfgwf.f90.

1331  ! -- modules
1332  use constantsmodule, only: linelength
1333  ! -- dummy
1334  class(GwfExchangeType) :: this !< GwfExchangeType
1335  ! -- local
1336  integer(I4B) :: i, nm1, nm2, nmgnc1, nmgnc2
1337  character(len=*), parameter :: fmterr = &
1338  "('EXCHANGE NODES ', i0, ' AND ', i0,"// &
1339  "' NOT CONSISTENT WITH GNC NODES ', "// &
1340  "i0, ' AND ', i0)"
1341  !
1342  ! -- If exchange has ghost nodes, then initialize ghost node object
1343  ! This will read the ghost node blocks from the gnc input file.
1344  call this%gnc%gnc_df(this%gwfmodel1, m2=this%gwfmodel2)
1345  !
1346  ! -- Verify gnc is implicit if exchange has Newton Terms
1347  if (.not. this%gnc%implicit .and. this%inewton /= 0) then
1348  call store_error('GNC is explicit, but GWF exchange has active newton.')
1349  call store_error('Add implicit option to GNC or remove NEWTON from '// &
1350  'GWF exchange.')
1351  call store_error_unit(this%ingnc)
1352  end if
1353  !
1354  ! -- Perform checks to ensure GNCs match with GWF-GWF nodes
1355  if (this%nexg /= this%gnc%nexg) then
1356  call store_error('Number of exchanges does not match number of GNCs')
1357  call store_error_unit(this%ingnc)
1358  end if
1359  !
1360  ! -- Go through each entry and confirm
1361  do i = 1, this%nexg
1362  if (this%nodem1(i) /= this%gnc%nodem1(i) .or. &
1363  this%nodem2(i) /= this%gnc%nodem2(i)) then
1364  nm1 = this%gwfmodel1%dis%get_nodeuser(this%nodem1(i))
1365  nm2 = this%gwfmodel2%dis%get_nodeuser(this%nodem2(i))
1366  nmgnc1 = this%gwfmodel1%dis%get_nodeuser(this%gnc%nodem1(i))
1367  nmgnc2 = this%gwfmodel2%dis%get_nodeuser(this%gnc%nodem2(i))
1368  write (errmsg, fmterr) nm1, nm2, nmgnc1, nmgnc2
1369  call store_error(errmsg)
1370  end if
1371  end do
1372  if (count_errors() > 0) then
1373  call store_error_unit(this%ingnc)
1374  end if
1375  !
1376  ! -- close the file
1377  close (this%ingnc)
Here is the call graph for this function:

◆ read_mvr()

subroutine gwfgwfexchangemodule::read_mvr ( class(gwfexchangetype this,
integer(i4b), intent(in)  iout 
)

Read and process movers

Parameters
thisGwfExchangeType

Definition at line 1384 of file exg-gwfgwf.f90.

1385  ! -- modules
1386  use gwfexgmovermodule, only: exg_mvr_cr
1387  ! -- dummy
1388  class(GwfExchangeType) :: this !< GwfExchangeType
1389  integer(I4B), intent(in) :: iout
1390  class(DisBaseType), pointer :: dis
1391  ! -- local
1392  !
1393  ! -- Create and initialize the mover object Here, dis is set to the one
1394  ! for gwfmodel1 so that a call to save flows has an associated dis
1395  ! object. Because the conversion flags for the mover are both false,
1396  ! the dis object does not convert from reduced to user node numbers.
1397  ! So in this case, the dis object is just writing unconverted package
1398  ! numbers to the binary budget file.
1399  dis => null()
1400  if (this%v_model1%is_local) then
1401  dis => this%gwfmodel1%dis
1402  else if (this%v_model2%is_local) then
1403  dis => this%gwfmodel2%dis
1404  end if
1405  call exg_mvr_cr(this%mvr, this%name, this%inmvr, iout, dis)
1406  this%mvr%model1 => this%v_model1
1407  this%mvr%model2 => this%v_model2
1408  this%mvr%suppress_fileout = this%is_datacopy
subroutine, public exg_mvr_cr(exg_mvr, name_parent, inunit, iout, dis)
Here is the call graph for this function:

◆ rewet()

subroutine gwfgwfexchangemodule::rewet ( class(gwfexchangetype this,
integer(i4b), intent(in)  kiter 
)

Check if rewetting should propagate from one model to another

Parameters
thisGwfExchangeType

Definition at line 1415 of file exg-gwfgwf.f90.

1416  ! -- modules
1417  use tdismodule, only: kper, kstp
1418  ! -- dummy
1419  class(GwfExchangeType) :: this !< GwfExchangeType
1420  integer(I4B), intent(in) :: kiter
1421  ! -- local
1422  integer(I4B) :: iexg
1423  integer(I4B) :: n, m
1424  integer(I4B) :: ibdn, ibdm
1425  integer(I4B) :: ihc
1426  real(DP) :: hn, hm
1427  integer(I4B) :: irewet
1428  character(len=30) :: nodestrn, nodestrm
1429  character(len=*), parameter :: fmtrwt = &
1430  "(1x, 'CELL ',A,' REWET FROM GWF MODEL ',A,' CELL ',A, &
1431  &' FOR ITER. ',I0, ' STEP ',I0, ' PERIOD ', I0)"
1432  !
1433  ! -- Use model 1 to rewet model 2 and vice versa
1434  do iexg = 1, this%nexg
1435  n = this%nodem1(iexg)
1436  m = this%nodem2(iexg)
1437  hn = this%gwfmodel1%x(n)
1438  hm = this%gwfmodel2%x(m)
1439  ibdn = this%gwfmodel1%ibound(n)
1440  ibdm = this%gwfmodel2%ibound(m)
1441  ihc = this%ihc(iexg)
1442  call this%gwfmodel1%npf%rewet_check(kiter, n, hm, ibdm, ihc, &
1443  this%gwfmodel1%x, irewet)
1444  if (irewet == 1) then
1445  call this%gwfmodel1%dis%noder_to_string(n, nodestrn)
1446  call this%gwfmodel2%dis%noder_to_string(m, nodestrm)
1447  write (this%gwfmodel1%iout, fmtrwt) trim(nodestrn), &
1448  trim(this%gwfmodel2%name), trim(nodestrm), kiter, kstp, kper
1449  end if
1450  call this%gwfmodel2%npf%rewet_check(kiter, m, hn, ibdn, ihc, &
1451  this%gwfmodel2%x, irewet)
1452  if (irewet == 1) then
1453  call this%gwfmodel1%dis%noder_to_string(n, nodestrm)
1454  call this%gwfmodel2%dis%noder_to_string(m, nodestrn)
1455  write (this%gwfmodel2%iout, fmtrwt) trim(nodestrn), &
1456  trim(this%gwfmodel1%name), trim(nodestrm), kiter, kstp, kper
1457  end if
1458  !
1459  end do

◆ source_options()

subroutine gwfgwfexchangemodule::source_options ( class(gwfexchangetype this,
integer(i4b), intent(in)  iout 
)

Source the options block

Parameters
thisGwfExchangeType

Definition at line 1230 of file exg-gwfgwf.f90.

1231  ! -- modules
1232  use constantsmodule, only: lenvarname, dem6
1233  use inputoutputmodule, only: getunit, openfile
1237  use sourcecommonmodule, only: filein_fname
1238  ! -- dummy
1239  class(GwfExchangeType) :: this !< GwfExchangeType
1240  integer(I4B), intent(in) :: iout
1241  ! -- local
1242  type(ExgGwfgwfParamFoundType) :: found
1243  character(len=LENVARNAME), dimension(3) :: cellavg_method = &
1244  &[character(len=LENVARNAME) :: 'HARMONIC', 'LOGARITHMIC', 'AMT-LMK']
1245  character(len=linelength) :: gnc_fname, mvr_fname
1246  !
1247  ! -- update defaults with idm sourced values
1248  call mem_set_value(this%icellavg, 'CELL_AVERAGING', this%input_mempath, &
1249  cellavg_method, found%cell_averaging)
1250  call mem_set_value(this%inewton, 'NEWTON', this%input_mempath, found%newton)
1251  call mem_set_value(this%ixt3d, 'XT3D', this%input_mempath, found%xt3d)
1252  call mem_set_value(this%ivarcv, 'VARIABLECV', this%input_mempath, &
1253  found%variablecv)
1254  call mem_set_value(this%idewatcv, 'DEWATERED', this%input_mempath, &
1255  found%dewatered)
1256  !
1257  write (iout, '(1x,a)') 'PROCESSING GWF-GWF EXCHANGE OPTIONS'
1258  !
1259  ! -- source base class options
1260  call this%DisConnExchangeType%source_options(iout)
1261  !
1262  if (found%cell_averaging) then
1263  ! -- count from 0
1264  this%icellavg = this%icellavg - 1
1265  write (iout, '(4x,a,a)') &
1266  'CELL AVERAGING METHOD HAS BEEN SET TO: ', &
1267  trim(cellavg_method(this%icellavg + 1))
1268  end if
1269  !
1270  if (found%newton) then
1271  write (iout, '(4x,a)') &
1272  'NEWTON-RAPHSON method used for unconfined cells'
1273  end if
1274  !
1275  if (found%xt3d) then
1276  write (iout, '(4x,a)') 'XT3D WILL BE APPLIED ON THE INTERFACE'
1277  end if
1278  !
1279  if (found%variablecv) then
1280  write (iout, '(4x,a)') &
1281  'VERTICAL CONDUCTANCE VARIES WITH WATER TABLE.'
1282  end if
1283  !
1284  if (found%dewatered) then
1285  write (iout, '(4x,a)') &
1286  'VERTICAL CONDUCTANCE ACCOUNTS FOR DEWATERED PORTION OF '// &
1287  'AN UNDERLYING CELL.'
1288  end if
1289  !
1290  ! -- enforce 0 or 1 GNC6_FILENAME entries in option block
1291  if (filein_fname(gnc_fname, 'GNC6_FILENAME', this%input_mempath, &
1292  this%filename)) then
1293  this%ingnc = getunit()
1294  call openfile(this%ingnc, iout, gnc_fname, 'GNC')
1295  write (iout, '(4x,a)') &
1296  'GHOST NODES WILL BE READ FROM ', trim(gnc_fname)
1297  end if
1298  !
1299  ! -- enforce 0 or 1 MVR6_FILENAME entries in option block
1300  if (filein_fname(mvr_fname, 'MVR6_FILENAME', this%input_mempath, &
1301  this%filename)) then
1302  this%inmvr = getunit()
1303  call openfile(this%inmvr, iout, mvr_fname, 'MVR')
1304  write (iout, '(4x,a)') &
1305  'WATER MOVER INFORMATION WILL BE READ FROM ', trim(mvr_fname)
1306  end if
1307  !
1308  ! -- enforce 0 or 1 OBS6_FILENAME entries in option block
1309  if (.not. this%is_datacopy) then
1310  if (filein_fname(this%obs%inputFilename, 'OBS6_FILENAME', &
1311  this%input_mempath, this%filename)) then
1312  this%obs%active = .true.
1313  this%obs%inUnitObs = getunit()
1314  call openfile(this%obs%inUnitObs, iout, this%obs%inputFilename, 'OBS')
1315  end if
1316  end if
1317  !
1318  write (iout, '(1x,a)') 'END OF GWF-GWF EXCHANGE OPTIONS'
1319  !
1320  ! -- set omega value used for saturation calculations
1321  if (this%inewton > 0) then
1322  this%satomega = dem6
1323  end if
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
real(dp), parameter dem6
real constant 1e-6
Definition: Constants.f90:109
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function:

◆ use_interface_model()

logical(lgp) function gwfgwfexchangemodule::use_interface_model ( class(gwfexchangetype this)
private
Parameters
thisGwfExchangeType
Returns
true when interface model should be used

Definition at line 1961 of file exg-gwfgwf.f90.

1963  ! -- dummy
1964  class(GwfExchangeType) :: this !< GwfExchangeType
1965  ! -- return
1966  logical(LGP) :: use_im !< true when interface model should be used
1967  ! -- local
1968  integer(I4B) :: inbuy_m1
1969 
1970  use_im = this%DisConnExchangeType%use_interface_model()
1971  use_im = use_im .or. (this%ixt3d > 0)
1972 
1973  inbuy_m1 = 0
1974  select type (m => this%v_model1)
1975  class is (virtualgwfmodeltype)
1976  inbuy_m1 = m%inbuy%get()
1977  end select
1978  use_im = use_im .or. (inbuy_m1 > 0)

◆ validate_exchange()

subroutine gwfgwfexchangemodule::validate_exchange ( class(gwfexchangetype this)
Parameters
thisGwfExchangeType

Definition at line 273 of file exg-gwfgwf.f90.

274  ! -- modules
275  class(GwfExchangeType) :: this !< GwfExchangeType
276  ! -- local
277  logical(LGP) :: has_k22, has_spdis, has_vsc
278  !
279  ! Periodic boundary condition in exchange don't allow XT3D (=interface model)
280  if (this%v_model1 == this%v_model2) then
281  if (this%ixt3d > 0) then
282  write (errmsg, '(3a)') 'GWF-GWF exchange ', trim(this%name), &
283  ' is a periodic boundary condition which cannot'// &
284  ' be configured with XT3D'
285  call store_error(errmsg, terminate=.true.)
286  end if
287  end if
288  !
289  ! XT3D needs angle information
290  if (this%ixt3d > 0 .and. this%ianglex == 0) then
291  write (errmsg, '(3a)') 'GWF-GWF exchange ', trim(this%name), &
292  ' requires that ANGLDEGX be specified as an'// &
293  ' auxiliary variable because XT3D is enabled'
294  call store_error(errmsg, terminate=.true.)
295  end if
296  !
297  ! determine if specific functionality is demanded,
298  ! in model 1 or model 2 (in parallel, only one of
299  ! the models is checked, but the exchange is duplicated)
300  has_k22 = .false.
301  has_spdis = .false.
302  has_vsc = .false.
303  if (associated(this%gwfmodel1)) then
304  has_k22 = (this%gwfmodel1%npf%ik22 /= 0)
305  has_spdis = (this%gwfmodel1%npf%icalcspdis /= 0)
306  has_vsc = (this%gwfmodel1%npf%invsc /= 0)
307  end if
308  if (associated(this%gwfmodel2)) then
309  has_k22 = has_k22 .or. (this%gwfmodel2%npf%ik22 /= 0)
310  has_spdis = has_spdis .or. (this%gwfmodel2%npf%icalcspdis /= 0)
311  has_vsc = has_vsc .or. (this%gwfmodel2%npf%invsc /= 0)
312  end if
313  !
314  ! If horizontal anisotropy is in either model1 or model2,
315  ! ANGLDEGX must be provided as an auxiliary variable for this
316  ! GWF-GWF exchange (this%ianglex > 0).
317  if (has_k22) then
318  if (this%ianglex == 0) then
319  write (errmsg, '(3a)') 'GWF-GWF exchange ', trim(this%name), &
320  ' requires that ANGLDEGX be specified as an'// &
321  ' auxiliary variable because K22 was specified'// &
322  ' in one or both groundwater models.'
323  call store_error(errmsg, terminate=.true.)
324  end if
325  end if
326  !
327  ! If specific discharge is needed for model1 or model2,
328  ! ANGLDEGX must be provided as an auxiliary variable for this
329  ! GWF-GWF exchange (this%ianglex > 0).
330  if (has_spdis) then
331  if (this%ianglex == 0) then
332  write (errmsg, '(3a)') 'GWF-GWF exchange ', trim(this%name), &
333  ' requires that ANGLDEGX be specified as an'// &
334  ' auxiliary variable because specific discharge'// &
335  ' is being calculated in one or both'// &
336  ' groundwater models.'
337  call store_error(errmsg, terminate=.true.)
338  end if
339  if (this%icdist == 0) then
340  write (errmsg, '(3a)') 'GWF-GWF exchange ', trim(this%name), &
341  ' requires that CDIST be specified as an'// &
342  ' auxiliary variable because specific discharge'// &
343  ' is being calculated in one or both'// &
344  ' groundwater models.'
345  call store_error(errmsg, terminate=.true.)
346  end if
347  end if
348  !
349  ! If viscosity is on in either model, then terminate with an
350  ! error as viscosity package doesn't work yet with exchanges.
351  if (has_vsc) then
352  write (errmsg, '(3a)') 'GWF-GWF exchange ', trim(this%name), &
353  ' requires that the Viscosity Package is inactive'// &
354  ' in both of the connected models.'
355  call store_error(errmsg, terminate=.true.)
356  end if
Here is the call graph for this function: