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

This module contains the SwfGwfExchangeModule Module. More...

Data Types

type  swfgwfexchangetype
 

Functions/Subroutines

subroutine initialize (this, filename, name, swf_ftype, id, m1_id, m2_id, input_mempath)
 @ brief Initialize SWF GWF exchange More...
 
subroutine swf_gwf_ac (this, sparse)
 @ brief Add connections More...
 
subroutine swf_gwf_mc (this, matrix_sln)
 @ brief Map connections More...
 
subroutine swf_gwf_fc (this, kiter, matrix_sln, rhs_sln, inwtflag)
 @ brief Fill coefficients More...
 
subroutine swf_gwf_cq (this, icnvg, isuppress_output, isolnid)
 @ brief Calculate flow More...
 
subroutine swf_gwf_da (this)
 @ brief Deallocate More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine allocate_arrays (this)
 Allocate array data, using the number of connected nodes. More...
 
integer(i4b) function noder (this, model, cellid, iout)
 
character(len=20) function cellstr (this, model, cellid, iout)
 
subroutine swf_gwf_calc_simvals (this)
 Calculate flow rates for the exchanges and store them in a member array. More...
 
real(dp) function qcalc (this, iexg, hswf, hgwf)
 @ brief Calculate flow More...
 
real(dp) function get_cond (this, iexg, hswf, hgwf)
 @ brief Calculate conductance More...
 
real(dp) function get_wetted_perimeter (this, nodeswf, depth)
 @ brief Get wetted perimeter for swf channel model More...
 
subroutine swf_gwf_add_to_flowja (this)
 Add exchange flow to each model flowja diagonal position so that residual is calculated correctly. More...
 
subroutine swf_gwf_bd (this, icnvg, isuppress_output, isolnid)
 @ brief Budget More...
 
subroutine swf_gwf_chd_bd (this)
 @ brief swf-gwf-chd-bd More...
 
subroutine swf_gwf_bdsav (this)
 @ brief Budget save More...
 
subroutine swf_gwf_ot (this)
 @ brief Output More...
 
subroutine swf_gwf_save_simvals (this)
 @ brief Save simulated flow observations More...
 
logical(lgp) function swf_gwf_connects_model (this, model)
 Should return true when the exchange should be added to the solution where the model resides. More...
 

Detailed Description

This module contains the code for connecting a SWF model with a GWF model. The SwfGwfExchangeType class is a parent to the CHF and OLF models and should not be used directly.

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine swfgwfexchangemodule::allocate_arrays ( class(swfgwfexchangetype this)
private
Parameters
nexg
thisinstance of exchange object

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

391  ! -- dummy
392  class(SwfGwfExchangeType) :: this !< instance of exchange object
393  !
394  call mem_allocate(this%nodeswf, this%nexg, 'NODEM1', this%memoryPath)
395  call mem_allocate(this%nodegwf, this%nexg, 'NODEM2', this%memoryPath)
396  call mem_allocate(this%bedleak, this%nexg, 'BEDLEAK', this%memoryPath)
397  call mem_allocate(this%cfact, this%nexg, 'CFACT', this%memoryPath)
398  call mem_allocate(this%idxglo, this%nexg, 'IDXGLO', this%memoryPath)
399  call mem_allocate(this%idxsymglo, this%nexg, 'IDXSYMGLO', this%memoryPath)
400  call mem_allocate(this%simvals, this%nexg, 'SIMVALS', this%memoryPath)

◆ allocate_scalars()

subroutine swfgwfexchangemodule::allocate_scalars ( class(swfgwfexchangetype this)

Allocate scalar variables

Parameters
thisSwfGwfExchangeType

Definition at line 364 of file exg-swfgwf.f90.

365  ! -- modules
366  ! -- dummy
367  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
368  !
369  allocate (this%swf_ftype)
370  this%swf_ftype = ''
371  allocate (this%filename)
372  this%filename = ''
373  !
374  call mem_allocate(this%ipr_input, 'IPR_INPUT', this%memoryPath)
375  call mem_allocate(this%ipr_flow, 'IPR_FLOW', this%memoryPath)
376  call mem_allocate(this%ifixedcond, 'IFIXEDCOND', this%memoryPath)
377  call mem_allocate(this%nexg, 'NEXG', this%memoryPath)
378  call mem_allocate(this%inobs, 'INOBS', this%memoryPath)
379  !
380  this%ipr_input = 0
381  this%ipr_flow = 0
382  this%ifixedcond = 0
383  this%nexg = 0
384  this%inobs = 0

◆ cellstr()

character(len=20) function swfgwfexchangemodule::cellstr ( class(swfgwfexchangetype this,
class(numericalmodeltype), intent(in), pointer  model,
integer(i4b), dimension(:), intent(in)  cellid,
integer(i4b), intent(in)  iout 
)
Parameters
thisinstance of exchange object
[in]ioutthe output file unit

Definition at line 432 of file exg-swfgwf.f90.

433  ! -- modules
434  ! -- dummy
435  class(SwfGwfExchangeType) :: this !< instance of exchange object
436  class(NumericalModelType), pointer, intent(in) :: model
437  integer(I4B), dimension(:), intent(in) :: cellid
438  integer(I4B), intent(in) :: iout !< the output file unit
439  character(len=20) :: cellstr
440  character(len=*), parameter :: fmtndim1 = &
441  "('(',i0,')')"
442  character(len=*), parameter :: fmtndim2 = &
443  "('(',i0,',',i0,')')"
444  character(len=*), parameter :: fmtndim3 = &
445  "('(',i0,',',i0,',',i0,')')"
446  !
447  cellstr = ''
448  !
449  select case (model%dis%ndim)
450  case (1)
451  write (cellstr, fmtndim1) cellid(1)
452  case (2)
453  write (cellstr, fmtndim2) cellid(1), cellid(2)
454  case (3)
455  write (cellstr, fmtndim3) cellid(1), cellid(2), cellid(3)
456  case default
457  end select

◆ get_cond()

real(dp) function swfgwfexchangemodule::get_cond ( class(swfgwfexchangetype this,
integer(i4b), intent(in)  iexg,
real(dp), intent(in)  hswf,
real(dp), intent(in)  hgwf 
)
private

Calculate the conductance between the surface water cell and the underlying groundwater cell.

Parameters
thisSwfGwfExchangeType
[in]iexgexchange number
[in]hswfsurface water model head
[in]hgwfgroundwater model head

Definition at line 518 of file exg-swfgwf.f90.

519  ! module
520  use smoothingmodule, only: squadratic
521  ! return
522  real(DP) :: get_cond
523  ! dummy
524  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
525  integer(I4B), intent(in) :: iexg !< exchange number
526  real(DP), intent(in) :: hswf !< surface water model head
527  real(DP), intent(in) :: hgwf !< groundwater model head
528  ! local
529  integer(I4B) :: nodeswf
530  real(DP) :: range = 1.d-6
531  real(DP) :: depth_ups
532  real(DP) :: dydx
533  real(DP) :: smooth_factor
534  real(DP) :: area
535  real(DP) :: perimeter
536 
537  ! -- Calculate or return conductance
538  area = this%cfact(iexg)
539  if (this%ifixedcond == 1) then
540  get_cond = this%bedleak(iexg) * area
541  return
542  end if
543 
544  ! Calculate smooth factor between zero, when the upstream-weighted
545  ! depth is zero, and 1.0, when the upstream weighted depth is
546  ! greater than or equal to the smoothening depth
547  nodeswf = this%nodeswf(iexg)
548  depth_ups = max(hswf, hgwf) - this%swfmodel%dis%bot(nodeswf)
549  call squadratic(depth_ups, range, dydx, smooth_factor)
550 
551  ! For channel model calculate the interaction area as product
552  ! of cfact and upstream-wetted perimeter
553  if (this%swfmodel%dfw%is2d == 0) then
554  perimeter = this%get_wetted_perimeter(nodeswf, depth_ups)
555  area = area * perimeter
556  end if
557 
558  ! Calculate conductance
559  get_cond = smooth_factor * this%bedleak(iexg) * area
subroutine squadratic(x, range, dydx, y)
@ brief sQuadratic
Here is the call graph for this function:

◆ get_wetted_perimeter()

real(dp) function swfgwfexchangemodule::get_wetted_perimeter ( class(swfgwfexchangetype this,
integer(i4b), intent(in)  nodeswf,
real(dp), intent(in)  depth 
)
Parameters
thisSwfGwfExchangeType
[in]nodeswfnode number for surface water model cell
[in]depthwater depth in surface water model cell

Definition at line 564 of file exg-swfgwf.f90.

565  ! return
566  real(DP) :: wp
567  ! dummy
568  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
569  integer(I4B), intent(in) :: nodeswf !< node number for surface water model cell
570  real(DP), intent(in) :: depth !< water depth in surface water model cell
571  ! local
572  integer(I4B) :: idcxs
573  real(DP) :: width
574  real(DP) :: dummy
575 
576  idcxs = this%swfmodel%dfw%idcxs(nodeswf)
577  call this%swfmodel%dis%get_flow_width(nodeswf, nodeswf, 0, width, dummy)
578  wp = this%swfmodel%cxs%get_wetted_perimeter(idcxs, width, depth)
579 

◆ initialize()

subroutine swfgwfexchangemodule::initialize ( class(swfgwfexchangetype this,
character(len=*), intent(in)  filename,
character(len=*)  name,
character(len=*)  swf_ftype,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  m1_id,
integer(i4b), intent(in)  m2_id,
character(len=*), intent(in)  input_mempath 
)
private
Parameters
thisSwfGwfExchangeType
[in]filenamefilename for reading
nameexchange name
swf_ftypetype of swf model, CHF or OLF
[in]idid for the exchange
[in]m1_idid for model 1
[in]m2_idid for model 2

Definition at line 96 of file exg-swfgwf.f90.

98  ! dummy
99  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
100  character(len=*), intent(in) :: filename !< filename for reading
101  character(len=*) :: name !< exchange name
102  character(len=*) :: swf_ftype !< type of swf model, CHF or OLF
103  integer(I4B), intent(in) :: id !< id for the exchange
104  integer(I4B), intent(in) :: m1_id !< id for model 1
105  integer(I4B), intent(in) :: m2_id !< id for model 2
106  character(len=*), intent(in) :: input_mempath
107  ! -- local
108  class(BaseModelType), pointer :: mb
109  integer(I4B) :: m1_index, m2_index
110 
111  ! Assign id and name
112  this%id = id
113  this%name = name
114  this%memoryPath = create_mem_path(this%name)
115  this%input_mempath = input_mempath
116 
117  ! allocate scalars and set defaults
118  call this%allocate_scalars()
119  this%filename = filename
120  this%swf_ftype = swf_ftype
121  this%typename = trim(swf_ftype)//'-GWF'
122 
123  ! set swfmodel
124  m1_index = model_loc_idx(m1_id)
125  if (m1_index > 0) then
126  mb => getbasemodelfromlist(basemodellist, m1_index)
127  select type (mb)
128  class is (swfmodeltype)
129  this%model1 => mb
130  this%swfmodel => mb
131  end select
132  end if
133  ! this%v_model1 => get_virtual_model(m1_id)
134  ! this%is_datacopy = .not. this%v_model1%is_local
135 
136  ! set gwfmodel
137  m2_index = model_loc_idx(m2_id)
138  if (m2_index > 0) then
139  mb => getbasemodelfromlist(basemodellist, m2_index)
140  select type (mb)
141  type is (gwfmodeltype)
142  this%model2 => mb
143  this%gwfmodel => mb
144  end select
145  end if
146  ! this%v_model2 => get_virtual_model(m2_id)
147 
148  ! Verify that the surface water model is of the correct type
149  if (.not. associated(this%swfmodel) .and. m1_index > 0) then
150  write (errmsg, '(7a)') &
151  'Problem with ', &
152  trim(this%typename), &
153  ' exchange ', &
154  trim(this%name), &
155  '. Specified ', &
156  trim(this%swf_ftype), &
157  ' model does not appear to be of the correct type.'
158  call store_error(errmsg, terminate=.true.)
159  end if
160 
161  ! Verify that gwf model is of the correct type
162  if (.not. associated(this%gwfmodel) .and. m2_index > 0) then
163  write (errmsg, '(3a)') 'Problem with SWF-GWF exchange ', &
164  trim(this%name), &
165  '. Specified GWF model does not appear to be of the correct type.'
166  call store_error(errmsg, terminate=.true.)
167  end if
168 
169  ! Create the obs package
170  call obs_cr(this%obs, this%inobs)
171 
Here is the call graph for this function:

◆ noder()

integer(i4b) function swfgwfexchangemodule::noder ( class(swfgwfexchangetype this,
class(numericalmodeltype), intent(in), pointer  model,
integer(i4b), dimension(:), intent(in)  cellid,
integer(i4b), intent(in)  iout 
)
private
Parameters
thisinstance of exchange object
[in]ioutthe output file unit

Definition at line 405 of file exg-swfgwf.f90.

406  ! -- modules
407  use geomutilmodule, only: get_node
408  ! -- dummy
409  class(SwfGwfExchangeType) :: this !< instance of exchange object
410  class(NumericalModelType), pointer, intent(in) :: model
411  integer(I4B), dimension(:), intent(in) :: cellid
412  integer(I4B), intent(in) :: iout !< the output file unit
413  integer(I4B) :: noder, node
414  !
415  if (model%dis%ndim == 1) then
416  node = cellid(1)
417  elseif (model%dis%ndim == 2) then
418  node = get_node(cellid(1), 1, cellid(2), &
419  model%dis%mshape(1), 1, &
420  model%dis%mshape(2))
421  else
422  node = get_node(cellid(1), cellid(2), cellid(3), &
423  model%dis%mshape(1), &
424  model%dis%mshape(2), &
425  model%dis%mshape(3))
426  end if
427  noder = model%dis%get_nodenumber(node, 0)
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:83
Here is the call graph for this function:

◆ qcalc()

real(dp) function swfgwfexchangemodule::qcalc ( class(swfgwfexchangetype this,
integer(i4b), intent(in)  iexg,
real(dp), intent(in)  hswf,
real(dp), intent(in)  hgwf 
)

Calculate the flow for the specified exchange and node numbers. Flow is positive into the surface water model

Parameters
thisSwfGwfExchangeType

Definition at line 496 of file exg-swfgwf.f90.

497  ! -- return
498  real(DP) :: qcalc
499  ! -- dummy
500  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
501  integer(I4B), intent(in) :: iexg
502  real(DP), intent(in) :: hswf
503  real(DP), intent(in) :: hgwf
504  ! -- local
505  real(DP) :: cond
506 
507  ! Calculate flow between swf and gwf models; positive into swf
508  cond = this%get_cond(iexg, hswf, hgwf)
509  qcalc = cond * (hgwf - hswf)
510 

◆ swf_gwf_ac()

subroutine swfgwfexchangemodule::swf_gwf_ac ( class(swfgwfexchangetype this,
type(sparsematrix), intent(inout)  sparse 
)
private

Override parent exg_ac so that gnc can add connections here.

Parameters
thisSwfGwfExchangeType

Definition at line 178 of file exg-swfgwf.f90.

179  ! -- modules
180  use sparsemodule, only: sparsematrix
181  ! -- dummy
182  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
183  type(sparsematrix), intent(inout) :: sparse
184  ! -- local
185  integer(I4B) :: n, iglo, jglo
186  !
187  ! -- add exchange connections
188  do n = 1, this%nexg
189  iglo = this%nodeswf(n) + this%swfmodel%moffset
190  jglo = this%nodegwf(n) + this%gwfmodel%moffset
191  call sparse%addconnection(iglo, jglo, 1)
192  call sparse%addconnection(jglo, iglo, 1)
193  end do

◆ swf_gwf_add_to_flowja()

subroutine swfgwfexchangemodule::swf_gwf_add_to_flowja ( class(swfgwfexchangetype this)
private
Parameters
thisSwfGwfExchangeType

Definition at line 585 of file exg-swfgwf.f90.

586  ! -- modules
587  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
588  ! -- local
589  integer(I4B) :: i
590  integer(I4B) :: n
591  integer(I4B) :: idiag
592  real(DP) :: flow
593  !
594  do i = 1, this%nexg
595  !
596  if (associated(this%swfmodel)) then
597  n = this%nodeswf(i)
598  if (this%swfmodel%ibound(n) > 0) then
599  flow = this%simvals(i)
600  idiag = this%swfmodel%ia(n)
601  this%swfmodel%flowja(idiag) = this%swfmodel%flowja(idiag) + flow
602  end if
603  end if
604  !
605  if (associated(this%gwfmodel)) then
606  n = this%nodegwf(i)
607  if (this%gwfmodel%ibound(n) > 0) then
608  flow = -this%simvals(i)
609  idiag = this%gwfmodel%ia(n)
610  this%gwfmodel%flowja(idiag) = this%gwfmodel%flowja(idiag) + flow
611  end if
612  end if
613  !
614  end do

◆ swf_gwf_bd()

subroutine swfgwfexchangemodule::swf_gwf_bd ( class(swfgwfexchangetype this,
integer(i4b), intent(inout)  icnvg,
integer(i4b), intent(in)  isuppress_output,
integer(i4b), intent(in)  isolnid 
)
private

Accumulate budget terms

Parameters
thisSwfGwfExchangeType

Definition at line 621 of file exg-swfgwf.f90.

622  ! -- modules
624  use budgetmodule, only: rate_accumulator
625  ! -- dummy
626  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
627  integer(I4B), intent(inout) :: icnvg
628  integer(I4B), intent(in) :: isuppress_output
629  integer(I4B), intent(in) :: isolnid
630  ! -- local
631  character(len=LENBUDTXT), dimension(1) :: budtxt
632  real(DP), dimension(2, 1) :: budterm
633  real(DP) :: ratin, ratout
634  !
635  ! -- initialize
636  budtxt(1) = ' FLOW-JA-FACE'
637  !
638  ! -- Calculate ratin/ratout and pass to model budgets
639  call rate_accumulator(this%simvals, ratin, ratout)
640  !
641  ! -- Add the budget terms to model 1
642  if (associated(this%swfmodel)) then
643  budterm(1, 1) = ratin
644  budterm(2, 1) = ratout
645  call this%swfmodel%model_bdentry(budterm, budtxt, this%name)
646  end if
647  !
648  ! -- Add the budget terms to model 2
649  if (associated(this%gwfmodel)) then
650  budterm(1, 1) = ratout
651  budterm(2, 1) = ratin
652  call this%gwfmodel%model_bdentry(budterm, budtxt, this%name)
653  end if
654  !
655  ! -- Add any flows from one model into a constant head in another model
656  ! as a separate budget term called FLOW-JA-FACE-CHD
657  call this%swf_gwf_chd_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
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
Here is the call graph for this function:

◆ swf_gwf_bdsav()

subroutine swfgwfexchangemodule::swf_gwf_bdsav ( class(swfgwfexchangetype this)

Output individual flows to listing file and binary budget files

Parameters
thisSwfGwfExchangeType

Definition at line 726 of file exg-swfgwf.f90.

727  ! -- dummy
728  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
729  ! -- local
730  integer(I4B) :: icbcfl, ibudfl
731  ! !
732  ! ! -- budget for model1
733  ! if (associated(this%swfmodel1)) then
734  ! call this%swf_gwf_bdsav_model(this%swfmodel1, this%gwfmodel2%name)
735  ! end if
736  ! !
737  ! ! -- budget for model2
738  ! if (associated(this%gwfmodel2)) then
739  ! call this%swf_gwf_bdsav_model(this%gwfmodel2, this%swfmodel1%name)
740  ! end if
741  !
742  ! -- Set icbcfl, ibudfl to zero so that flows will be printed and
743  ! saved, if the options were set in the MVR package
744  icbcfl = 1
745  ibudfl = 1
746  !
747  ! -- Calculate and write simulated values for observations
748  if (this%inobs /= 0) then
749  call this%swf_gwf_save_simvals()
750  end if

◆ swf_gwf_calc_simvals()

subroutine swfgwfexchangemodule::swf_gwf_calc_simvals ( class(swfgwfexchangetype this)
private
Parameters
thisSwfGwfExchangeType

Definition at line 463 of file exg-swfgwf.f90.

464  ! -- modules
465  use constantsmodule, only: dzero
466  ! -- dummy
467  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
468  ! -- local
469  integer(I4B) :: iexg
470  integer(I4B) :: nodeswf, nodegwf
471  integer(I4B) :: ibdn1, ibdn2
472  real(DP) :: hswf
473  real(DP) :: hgwf
474  real(DP) :: rrate
475  !
476  do iexg = 1, this%nexg
477  rrate = dzero
478  nodeswf = this%nodeswf(iexg)
479  nodegwf = this%nodegwf(iexg)
480  ibdn1 = this%swfmodel%ibound(nodeswf)
481  ibdn2 = this%gwfmodel%ibound(nodegwf)
482  hswf = this%swfmodel%x(nodeswf)
483  hgwf = this%gwfmodel%x(nodegwf)
484  if (ibdn1 /= 0 .and. ibdn2 /= 0) then
485  rrate = this%qcalc(iexg, hswf, hgwf)
486  end if
487  this%simvals(iexg) = rrate
488  end do

◆ swf_gwf_chd_bd()

subroutine swfgwfexchangemodule::swf_gwf_chd_bd ( class(swfgwfexchangetype this)

Account for flow from an external model into a chd cell

Parameters
thisGwfExchangeType

Definition at line 664 of file exg-swfgwf.f90.

665  ! -- modules
667  ! -- dummy
668  class(SwfGwfExchangeType) :: this !< GwfExchangeType
669  ! -- local
670  character(len=LENBUDTXT), dimension(1) :: budtxt
671  integer(I4B) :: n
672  integer(I4B) :: i
673  real(DP), dimension(2, 1) :: budterm
674  real(DP) :: ratin, ratout
675  real(DP) :: q
676  !
677  ! -- initialize
678  budtxt(1) = 'FLOW-JA-FACE-CHD'
679  !
680  ! -- Add the constant-head budget terms for flow from model 2 into model 1
681  if (associated(this%swfmodel)) then
682  ratin = dzero
683  ratout = dzero
684  do i = 1, this%nexg
685  n = this%nodeswf(i)
686  if (this%swfmodel%ibound(n) < 0) then
687  q = this%simvals(i)
688  if (q > dzero) then
689  ratout = ratout + q
690  else
691  ratin = ratin - q
692  end if
693  end if
694  end do
695  budterm(1, 1) = ratin
696  budterm(2, 1) = ratout
697  call this%swfmodel%model_bdentry(budterm, budtxt, this%name)
698  end if
699  !
700  ! -- Add the constant-head budget terms for flow from model 1 into model 2
701  if (associated(this%gwfmodel)) then
702  ratin = dzero
703  ratout = dzero
704  do i = 1, this%nexg
705  n = this%nodegwf(i)
706  if (this%gwfmodel%ibound(n) < 0) then
707  ! -- flip flow sign as flow is relative to model 1
708  q = -this%simvals(i)
709  if (q > dzero) then
710  ratout = ratout + q
711  else
712  ratin = ratin - q
713  end if
714  end if
715  end do
716  budterm(1, 1) = ratin
717  budterm(2, 1) = ratout
718  call this%gwfmodel%model_bdentry(budterm, budtxt, this%name)
719  end if

◆ swf_gwf_connects_model()

logical(lgp) function swfgwfexchangemodule::swf_gwf_connects_model ( class(swfgwfexchangetype this,
class(basemodeltype), intent(in), pointer  model 
)
Parameters
thisthe instance of the exchange
[in]modelthe model to which the exchange might hold a connection
Returns
true, when connected

Definition at line 1006 of file exg-swfgwf.f90.

1007  ! -- dummy
1008  class(SwfGwfExchangeType) :: this !< the instance of the exchange
1009  class(BaseModelType), pointer, intent(in) :: model !< the model to which the exchange might hold a connection
1010  ! -- return
1011  logical(LGP) :: is_connected !< true, when connected
1012  !
1013  is_connected = .false.
1014  select type (model)
1015  class is (gwfmodeltype)
1016  if (associated(this%gwfmodel, model)) then
1017  is_connected = .true.
1018  end if
1019  class is (swfmodeltype)
1020  if (associated(this%swfmodel, model)) then
1021  is_connected = .true.
1022  end if
1023  end select

◆ swf_gwf_cq()

subroutine swfgwfexchangemodule::swf_gwf_cq ( class(swfgwfexchangetype 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
thisSwfGwfExchangeType

Definition at line 309 of file exg-swfgwf.f90.

310  ! -- dummy
311  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
312  integer(I4B), intent(inout) :: icnvg
313  integer(I4B), intent(in) :: isuppress_output
314  integer(I4B), intent(in) :: isolnid
315  !
316  ! -- calculate flow and store in simvals
317  call this%swf_gwf_calc_simvals()
318  !
319  ! -- set flows to model edges in NPF
320  ! todo: do we add these flows for specific discharge calculation?
321  !call this%swf_gwf_set_flow_to_npf()
322  !
323  ! -- add exchange flows to model's flowja diagonal
324  call this%swf_gwf_add_to_flowja()

◆ swf_gwf_da()

subroutine swfgwfexchangemodule::swf_gwf_da ( class(swfgwfexchangetype this)
private

Deallocate memory associated with this object

Parameters
thisSwfGwfExchangeType

Definition at line 331 of file exg-swfgwf.f90.

332  ! -- modules
334  ! -- dummy
335  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
336  !
337  ! -- objects
338  call this%obs%obs_da()
339  deallocate (this%obs)
340  !
341  ! -- arrays
342  call mem_deallocate(this%nodeswf)
343  call mem_deallocate(this%nodegwf)
344  call mem_deallocate(this%bedleak)
345  call mem_deallocate(this%cfact)
346  call mem_deallocate(this%idxglo)
347  call mem_deallocate(this%idxsymglo)
348  call mem_deallocate(this%simvals)
349  !
350  ! -- scalars
351  deallocate (this%swf_ftype)
352  deallocate (this%filename)
353  call mem_deallocate(this%ipr_input)
354  call mem_deallocate(this%ipr_flow)
355  call mem_deallocate(this%ifixedcond)
356  call mem_deallocate(this%nexg)
357  call mem_deallocate(this%inobs)

◆ swf_gwf_fc()

subroutine swfgwfexchangemodule::swf_gwf_fc ( class(swfgwfexchangetype this,
integer(i4b), intent(in)  kiter,
class(matrixbasetype), pointer  matrix_sln,
real(dp), dimension(:), intent(inout)  rhs_sln,
integer(i4b), intent(in), optional  inwtflag 
)

Fill conductance into coefficient matrix. For now assume all connections are vertical and no newton correction is needed.

Parameters
thisSwfGwfExchangeType

Definition at line 224 of file exg-swfgwf.f90.

225  ! modules
227  ! dummy
228  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
229  integer(I4B), intent(in) :: kiter
230  class(MatrixBaseType), pointer :: matrix_sln
231  real(DP), dimension(:), intent(inout) :: rhs_sln
232  integer(I4B), optional, intent(in) :: inwtflag
233  ! -- local
234  integer(I4B) :: iexg
235  integer(I4B) :: nodeswf
236  integer(I4B) :: nodegwf
237  integer(I4B) :: nodeswf_sln
238  integer(I4B) :: nodegwf_sln
239  integer(I4B) :: ibdn1
240  integer(I4B) :: ibdn2
241  real(DP) :: hswf
242  real(DP) :: hgwf
243  real(DP) :: qnm
244  real(DP) :: qeps
245  real(DP) :: eps
246  real(DP) :: derv
247 
248  ! Fill terms into solution matrix and rhs vector
249  do iexg = 1, this%nexg
250 
251  nodeswf = this%nodeswf(iexg)
252  nodegwf = this%nodegwf(iexg)
253  nodeswf_sln = this%nodeswf(iexg) + this%swfmodel%moffset
254  nodegwf_sln = this%nodegwf(iexg) + this%gwfmodel%moffset
255  ibdn1 = this%swfmodel%ibound(nodeswf)
256  ibdn2 = this%gwfmodel%ibound(nodegwf)
257  hswf = this%swfmodel%x(nodeswf)
258  hgwf = this%gwfmodel%x(nodegwf)
259 
260  ! First add these terms to the row for the surface water model
261 
262  ! Fill the qnm term on the right-hand side
263  qnm = this%qcalc(iexg, hswf, hgwf)
264  rhs_sln(nodeswf_sln) = rhs_sln(nodeswf_sln) - qnm
265 
266  ! Derivative calculation and fill of n terms
267  eps = get_perturbation(hswf)
268  qeps = this%qcalc(iexg, hswf + eps, hgwf)
269  derv = (qeps - qnm) / eps
270  call matrix_sln%add_diag_value(nodeswf_sln, derv)
271  rhs_sln(nodeswf_sln) = rhs_sln(nodeswf_sln) + derv * hswf
272 
273  ! Derivative calculation and fill of m terms
274  eps = get_perturbation(hgwf)
275  qeps = this%qcalc(iexg, hswf, hgwf + eps)
276  derv = (qeps - qnm) / eps
277  call matrix_sln%add_value_pos(this%idxglo(iexg), derv)
278  rhs_sln(nodeswf_sln) = rhs_sln(nodeswf_sln) + derv * hgwf
279 
280  ! now add these terms to the row for the groundwater model
281 
282  ! Fill the qnm term on the right-hand side
283  qnm = -this%qcalc(iexg, hswf, hgwf)
284  rhs_sln(nodegwf_sln) = rhs_sln(nodegwf_sln) - qnm
285 
286  ! Derivative calculation and fill of n terms
287  eps = get_perturbation(hgwf)
288  qeps = -this%qcalc(iexg, hswf, hgwf + eps)
289  derv = (qeps - qnm) / eps
290  call matrix_sln%add_diag_value(nodegwf_sln, derv)
291  rhs_sln(nodegwf_sln) = rhs_sln(nodegwf_sln) + derv * hgwf
292 
293  ! Derivative calculation and fill of m terms
294  eps = get_perturbation(hswf)
295  qeps = -this%qcalc(iexg, hswf + eps, hgwf)
296  derv = (qeps - qnm) / eps
297  call matrix_sln%add_value_pos(this%idxsymglo(iexg), derv)
298  rhs_sln(nodegwf_sln) = rhs_sln(nodegwf_sln) + derv * hswf
299 
300  end do
301 
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:372
Here is the call graph for this function:

◆ swf_gwf_mc()

subroutine swfgwfexchangemodule::swf_gwf_mc ( class(swfgwfexchangetype this,
class(matrixbasetype), pointer  matrix_sln 
)

Map the connections in the global matrix

Parameters
thisSwfGwfExchangeType
matrix_slnthe system matrix

Definition at line 200 of file exg-swfgwf.f90.

201  ! -- modules
202  use sparsemodule, only: sparsematrix
203  ! -- dummy
204  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
205  class(MatrixBaseType), pointer :: matrix_sln !< the system matrix
206  ! -- local
207  integer(I4B) :: n, iglo, jglo
208  !
209  ! -- map exchange connections
210  do n = 1, this%nexg
211  iglo = this%nodeswf(n) + this%swfmodel%moffset
212  jglo = this%nodegwf(n) + this%gwfmodel%moffset
213  this%idxglo(n) = matrix_sln%get_position(iglo, jglo)
214  this%idxsymglo(n) = matrix_sln%get_position(jglo, iglo)
215  end do

◆ swf_gwf_ot()

subroutine swfgwfexchangemodule::swf_gwf_ot ( class(swfgwfexchangetype this)
private

Write output

Parameters
thisSwfGwfExchangeType

Definition at line 916 of file exg-swfgwf.f90.

917  ! -- modules
918  use simvariablesmodule, only: iout
919  use constantsmodule, only: dzero, linelength
920  ! -- dummy
921  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
922  ! -- local
923  integer(I4B) :: iexg, n1, n2
924  real(DP) :: flow
925  character(len=LINELENGTH) :: node1str, node2str
926  ! -- format
927  character(len=*), parameter :: fmtheader2 = &
928  "(/1x, 'SUMMARY OF EXCHANGE RATES FOR EXCHANGE ', a, ' WITH ID ', i0, /, &
929  &2a16, 4a16, /, 96('-'))"
930  character(len=*), parameter :: fmtdata = &
931  "(2a16, 5(1pg16.6))"
932  !
933  ! -- Call bdsave
934  call this%swf_gwf_bdsav()
935  !
936  ! -- Write a table of exchanges
937  if (this%ipr_flow /= 0) then
938  write (iout, fmtheader2) trim(adjustl(this%name)), this%id, 'NODEM1', &
939  'NODEM2', 'COND', 'X_M1', 'X_M2', 'FLOW'
940  do iexg = 1, this%nexg
941  n1 = this%nodeswf(iexg)
942  n2 = this%nodegwf(iexg)
943  flow = this%simvals(iexg)
944  call this%swfmodel%dis%noder_to_string(n1, node1str)
945  call this%gwfmodel%dis%noder_to_string(n2, node2str)
946  write (iout, fmtdata) trim(adjustl(node1str)), &
947  trim(adjustl(node2str)), &
948  this%bedleak(iexg), this%swfmodel%x(n1), &
949  this%gwfmodel%x(n2), flow
950  end do
951  end if
952  !
953  ! -- OBS output
954  call this%obs%obs_ot()
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) iout
file unit number for simulation output

◆ swf_gwf_save_simvals()

subroutine swfgwfexchangemodule::swf_gwf_save_simvals ( class(swfgwfexchangetype), intent(inout)  this)

Save the simulated flows for each exchange

Definition at line 961 of file exg-swfgwf.f90.

962  ! -- modules
964  use simvariablesmodule, only: errmsg
965  use constantsmodule, only: dzero
966  use observemodule, only: observetype
967  ! -- dummy
968  class(SwfGwfExchangeType), intent(inout) :: this
969  ! -- local
970  integer(I4B) :: i
971  integer(I4B) :: j
972  integer(I4B) :: n1
973  integer(I4B) :: n2
974  integer(I4B) :: iexg
975  real(DP) :: v
976  type(ObserveType), pointer :: obsrv => null()
977  !
978  ! -- Write simulated values for all gwf-gwf observations
979  if (this%obs%npakobs > 0) then
980  call this%obs%obs_bd_clear()
981  do i = 1, this%obs%npakobs
982  obsrv => this%obs%pakobs(i)%obsrv
983  do j = 1, obsrv%indxbnds_count
984  iexg = obsrv%indxbnds(j)
985  v = dzero
986  select case (obsrv%ObsTypeId)
987  case ('FLOW-JA-FACE')
988  n1 = this%nodeswf(iexg)
989  n2 = this%nodegwf(iexg)
990  v = this%simvals(iexg)
991  case default
992  errmsg = 'Unrecognized observation type: '// &
993  trim(obsrv%ObsTypeId)
994  call store_error(errmsg)
995  call store_error_unit(this%inobs)
996  end select
997  call this%obs%SaveOneSimval(obsrv, v)
998  end do
999  end do
1000  end if
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
character(len=maxcharlen) errmsg
error message string
Here is the call graph for this function: