MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
swfgwfexchangemodule Module Reference

This module contains the SwfGwfExchangeModule Module. More...

Data Types

type  swfgwfexchangetype
 

Functions/Subroutines

subroutine, public swfgwf_cr (filename, name, id, m1_id, m2_id, input_mempath)
 @ brief Create SWF GWF exchange More...
 
subroutine swf_gwf_df (this)
 @ brief Define 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...
 
subroutine source_options (this, iout)
 @ brief Source options More...
 
subroutine source_dimensions (this, iout)
 Source dimension from input context. More...
 
integer(i4b) function noder (this, model, cellid, iout)
 
character(len=20) function cellstr (this, model, cellid, iout)
 
subroutine source_data (this, iout)
 Source exchange data from input context. More...
 
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, n1, n2)
 @ brief Calculate flow 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.

Function/Subroutine Documentation

◆ allocate_arrays()

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

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

396  ! -- dummy
397  class(SwfGwfExchangeType) :: this !< instance of exchange object
398  !
399  call mem_allocate(this%nodem1, this%nexg, 'NODEM1', this%memoryPath)
400  call mem_allocate(this%nodem2, this%nexg, 'NODEM2', this%memoryPath)
401  call mem_allocate(this%cond, this%nexg, 'COND', this%memoryPath)
402  call mem_allocate(this%idxglo, this%nexg, 'IDXGLO', this%memoryPath)
403  call mem_allocate(this%idxsymglo, this%nexg, 'IDXSYMGLO', this%memoryPath)
404  call mem_allocate(this%simvals, this%nexg, 'SIMVALS', this%memoryPath)
405  !
406  ! -- Return
407  return

◆ allocate_scalars()

subroutine swfgwfexchangemodule::allocate_scalars ( class(swfgwfexchangetype this)

Allocate scalar variables

Parameters
thisSwfGwfExchangeType

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

371  ! -- modules
372  ! -- dummy
373  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
374  !
375  allocate (this%filename)
376  this%filename = ''
377  !
378  call mem_allocate(this%ipr_input, 'IPR_INPUT', this%memoryPath)
379  call mem_allocate(this%ipr_flow, 'IPR_FLOW', this%memoryPath)
380  call mem_allocate(this%nexg, 'NEXG', this%memoryPath)
381  call mem_allocate(this%inobs, 'INOBS', this%memoryPath)
382  !
383  this%ipr_input = 0
384  this%ipr_flow = 0
385  this%nexg = 0
386  this%inobs = 0
387  !
388  ! -- Return
389  return

◆ cellstr()

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

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

522  ! -- modules
523  ! -- dummy
524  class(SwfGwfExchangeType) :: this !< instance of exchange object
525  class(NumericalModelType), pointer, intent(in) :: model
526  integer(I4B), dimension(:), pointer, intent(in) :: cellid
527  integer(I4B), intent(in) :: iout !< the output file unit
528  character(len=20) :: cellstr
529  character(len=*), parameter :: fmtndim1 = &
530  "('(',i0,')')"
531  character(len=*), parameter :: fmtndim2 = &
532  "('(',i0,',',i0,')')"
533  character(len=*), parameter :: fmtndim3 = &
534  "('(',i0,',',i0,',',i0,')')"
535  !
536  cellstr = ''
537  !
538  select case (model%dis%ndim)
539  case (1)
540  write (cellstr, fmtndim1) cellid(1)
541  case (2)
542  write (cellstr, fmtndim2) cellid(1), cellid(2)
543  case (3)
544  write (cellstr, fmtndim3) cellid(1), cellid(2), cellid(3)
545  case default
546  end select
547  !
548  ! -- return
549  return

◆ noder()

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

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

492  ! -- modules
493  use geomutilmodule, only: get_node
494  ! -- dummy
495  class(SwfGwfExchangeType) :: this !< instance of exchange object
496  class(NumericalModelType), pointer, intent(in) :: model
497  integer(I4B), dimension(:), pointer, intent(in) :: cellid
498  integer(I4B), intent(in) :: iout !< the output file unit
499  integer(I4B) :: noder, node
500  !
501  if (model%dis%ndim == 1) then
502  node = cellid(1)
503  elseif (model%dis%ndim == 2) then
504  node = get_node(cellid(1), 1, cellid(2), &
505  model%dis%mshape(1), 1, &
506  model%dis%mshape(2))
507  else
508  node = get_node(cellid(1), cellid(2), cellid(3), &
509  model%dis%mshape(1), &
510  model%dis%mshape(2), &
511  model%dis%mshape(3))
512  end if
513  noder = model%dis%get_nodenumber(node, 0)
514  !
515  ! -- return
516  return
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:80
Here is the call graph for this function:

◆ qcalc()

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

Calculate the flow for the specified exchange and node numbers

Parameters
thisSwfGwfExchangeType

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

688  ! -- return
689  real(DP) :: qcalc
690  ! -- dummy
691  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
692  integer(I4B), intent(in) :: iexg
693  integer(I4B), intent(in) :: n1
694  integer(I4B), intent(in) :: n2
695  ! -- local
696  !
697  ! -- Calculate flow between nodes in the two models
698  qcalc = this%cond(iexg) * (this%gwfmodel2%x(n2) - this%swfmodel1%x(n1))
699  !
700  ! -- Return
701  return

◆ source_data()

subroutine swfgwfexchangemodule::source_data ( class(swfgwfexchangetype this,
integer(i4b), intent(in)  iout 
)
private
Parameters
thisinstance of exchange object
[in]ioutthe output file unit

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

555  ! -- modules
557  ! -- dummy
558  class(SwfGwfExchangeType) :: this !< instance of exchange object
559  integer(I4B), intent(in) :: iout !< the output file unit
560  ! -- local
561  integer(I4B), dimension(:, :), contiguous, pointer :: cellidm1
562  integer(I4B), dimension(:, :), contiguous, pointer :: cellidm2
563  real(DP), dimension(:), contiguous, pointer :: cond
564  character(len=20) :: cellstr1, cellstr2
565  integer(I4B) :: nerr
566  integer(I4B) :: iexg, nodem1, nodem2
567  ! -- format
568  character(len=*), parameter :: fmtexglabel = "(1x, 3a10, 50(a16))"
569  character(len=*), parameter :: fmtexgdata = &
570  "(5x, a, 1x, a ,50(1pg16.6))"
571  !
572  call mem_setptr(cellidm1, 'CELLIDM1', this%input_mempath)
573  call mem_setptr(cellidm2, 'CELLIDM2', this%input_mempath)
574  call mem_setptr(cond, 'COND', this%input_mempath)
575  !
576  write (iout, '(1x,a)') 'PROCESSING EXCHANGEDATA'
577  !
578  if (this%ipr_input /= 0) then
579  write (iout, fmtexglabel) 'NODEM1', 'NODEM2', 'COND'
580  end if
581  !
582  do iexg = 1, this%nexg
583  !
584  if (associated(this%model1)) then
585  !
586  ! -- Determine user node number
587  nodem1 = this%noder(this%model1, cellidm1(:, iexg), iout)
588  this%nodem1(iexg) = nodem1
589  !
590  else
591  this%nodem1(iexg) = -1
592  end if
593  !
594  if (associated(this%model2)) then
595  !
596  ! -- Determine user node number
597  nodem2 = this%noder(this%model2, cellidm2(:, iexg), iout)
598  this%nodem2(iexg) = nodem2
599  !
600  else
601  this%nodem2(iexg) = -1
602  end if
603  !
604  ! -- Read rest of input line
605  this%cond(iexg) = cond(iexg)
606  !
607  ! -- Write the data to listing file if requested
608  if (this%ipr_input /= 0) then
609  cellstr1 = this%cellstr(this%model1, cellidm1(:, iexg), iout)
610  cellstr2 = this%cellstr(this%model2, cellidm2(:, iexg), iout)
611  write (iout, fmtexgdata) trim(cellstr1), trim(cellstr2), &
612  this%cond(iexg)
613  end if
614  !
615  ! -- Check to see if nodem1 is outside of active domain
616  if (associated(this%model1)) then
617  if (nodem1 <= 0) then
618  cellstr1 = this%cellstr(this%model1, cellidm1(:, iexg), iout)
619  write (errmsg, *) &
620  trim(adjustl(this%model1%name))// &
621  ' Cell is outside active grid domain ('// &
622  trim(adjustl(cellstr1))//').'
623  call store_error(errmsg)
624  end if
625  end if
626  !
627  ! -- Check to see if nodem2 is outside of active domain
628  if (associated(this%model2)) then
629  if (nodem2 <= 0) then
630  cellstr2 = this%cellstr(this%model2, cellidm2(:, iexg), iout)
631  write (errmsg, *) &
632  trim(adjustl(this%model2%name))// &
633  ' Cell is outside active grid domain ('// &
634  trim(adjustl(cellstr2))//').'
635  call store_error(errmsg)
636  end if
637  end if
638  end do
639  !
640  write (iout, '(1x,a)') 'END OF EXCHANGEDATA'
641  !
642  ! -- Stop if errors
643  nerr = count_errors()
644  if (nerr > 0) then
645  call store_error('Errors encountered in exchange input file.')
646  call store_error_filename(this%filename)
647  end if
648  !
649  ! -- Return
650  return
Here is the call graph for this function:

◆ source_dimensions()

subroutine swfgwfexchangemodule::source_dimensions ( class(swfgwfexchangetype this,
integer(i4b), intent(in)  iout 
)
Parameters
thisinstance of exchange object
[in]ioutfor logging

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

465  ! -- modules
468  ! -- dummy
469  class(SwfGwfExchangeType) :: this !< instance of exchange object
470  integer(I4B), intent(in) :: iout !< for logging
471  ! -- local
472  type(ExgSwfgwfParamFoundType) :: found
473  !
474  ! -- update defaults with idm sourced values
475  call mem_set_value(this%nexg, 'NEXG', this%input_mempath, found%nexg)
476  !
477  write (iout, '(1x,a)') 'PROCESSING EXCHANGE DIMENSIONS'
478  !
479  if (found%nexg) then
480  write (iout, '(4x,a,i0)') 'NEXG = ', this%nexg
481  end if
482  !
483  write (iout, '(1x,a)') 'END OF EXCHANGE DIMENSIONS'
484  !
485  ! -- return
486  return

◆ source_options()

subroutine swfgwfexchangemodule::source_options ( class(swfgwfexchangetype this,
integer(i4b), intent(in)  iout 
)
private

Source the options block

Parameters
thisGwfExchangeType

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

415  ! -- modules
416  use constantsmodule, only: lenvarname, dem6
422  ! -- dummy
423  class(SwfGwfExchangeType) :: this !< GwfExchangeType
424  integer(I4B), intent(in) :: iout
425  ! -- local
426  type(ExgSwfgwfParamFoundType) :: found
427  !
428  ! -- update defaults with idm sourced values
429  call mem_set_value(this%ipr_input, 'IPR_INPUT', &
430  this%input_mempath, found%ipr_input)
431  call mem_set_value(this%ipr_flow, 'IPR_FLOW', &
432  this%input_mempath, found%ipr_flow)
433  !
434  write (iout, '(1x,a)') 'PROCESSING SWF-GWF EXCHANGE OPTIONS'
435  !
436  if (found%ipr_input) then
437  write (iout, '(4x,a)') &
438  'THE LIST OF EXCHANGES WILL BE PRINTED.'
439  end if
440  !
441  if (found%ipr_flow) then
442  write (iout, '(4x,a)') &
443  'EXCHANGE FLOWS WILL BE PRINTED TO LIST FILES.'
444  end if
445  !
446  ! -- enforce 0 or 1 OBS6_FILENAME entries in option block
447  ! if (.not. this%is_datacopy) then
448  ! if (filein_fname(this%obs%inputFilename, 'OBS6_FILENAME', &
449  ! this%input_mempath, this%filename)) then
450  ! this%obs%active = .true.
451  ! this%obs%inUnitObs = GetUnit()
452  ! call openfile(this%obs%inUnitObs, iout, this%obs%inputFilename, 'OBS')
453  ! end if
454  ! end if
455  !
456  write (iout, '(1x,a)') 'END OF SWF-GWF EXCHANGE OPTIONS'
457  !
458  ! -- Return
459  return
This module contains simulation constants.
Definition: Constants.f90:9
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:108
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 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:

◆ 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 229 of file exg-swfgwf.f90.

230  ! -- modules
231  use sparsemodule, only: sparsematrix
232  ! -- dummy
233  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
234  type(sparsematrix), intent(inout) :: sparse
235  ! -- local
236  integer(I4B) :: n, iglo, jglo
237  !
238  ! -- add exchange connections
239  do n = 1, this%nexg
240  iglo = this%nodem1(n) + this%swfmodel1%moffset
241  jglo = this%nodem2(n) + this%gwfmodel2%moffset
242  call sparse%addconnection(iglo, jglo, 1)
243  call sparse%addconnection(jglo, iglo, 1)
244  end do
245  !
246  ! -- Return
247  return

◆ swf_gwf_add_to_flowja()

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

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

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

◆ 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 746 of file exg-swfgwf.f90.

747  ! -- modules
749  use budgetmodule, only: rate_accumulator
750  ! -- dummy
751  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
752  integer(I4B), intent(inout) :: icnvg
753  integer(I4B), intent(in) :: isuppress_output
754  integer(I4B), intent(in) :: isolnid
755  ! -- local
756  character(len=LENBUDTXT), dimension(1) :: budtxt
757  real(DP), dimension(2, 1) :: budterm
758  real(DP) :: ratin, ratout
759  !
760  ! -- initialize
761  budtxt(1) = ' FLOW-JA-FACE'
762  !
763  ! -- Calculate ratin/ratout and pass to model budgets
764  call rate_accumulator(this%simvals, ratin, ratout)
765  !
766  ! -- Add the budget terms to model 1
767  if (associated(this%swfmodel1)) then
768  budterm(1, 1) = ratin
769  budterm(2, 1) = ratout
770  call this%swfmodel1%model_bdentry(budterm, budtxt, this%name)
771  end if
772  !
773  ! -- Add the budget terms to model 2
774  if (associated(this%gwfmodel2)) then
775  budterm(1, 1) = ratout
776  budterm(2, 1) = ratin
777  call this%gwfmodel2%model_bdentry(budterm, budtxt, this%name)
778  end if
779  !
780  ! -- Add any flows from one model into a constant head in another model
781  ! as a separate budget term called FLOW-JA-FACE-CHD
782  call this%swf_gwf_chd_bd()
783  !
784  ! -- Return
785  return
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:664
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:22
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:36
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 857 of file exg-swfgwf.f90.

858  ! -- dummy
859  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
860  ! -- local
861  integer(I4B) :: icbcfl, ibudfl
862  ! !
863  ! ! -- budget for model1
864  ! if (associated(this%swfmodel1)) then
865  ! call this%swf_gwf_bdsav_model(this%swfmodel1, this%gwfmodel2%name)
866  ! end if
867  ! !
868  ! ! -- budget for model2
869  ! if (associated(this%gwfmodel2)) then
870  ! call this%swf_gwf_bdsav_model(this%gwfmodel2, this%swfmodel1%name)
871  ! end if
872  !
873  ! -- Set icbcfl, ibudfl to zero so that flows will be printed and
874  ! saved, if the options were set in the MVR package
875  icbcfl = 1
876  ibudfl = 1
877  !
878  ! -- Calculate and write simulated values for observations
879  if (this%inobs /= 0) then
880  call this%swf_gwf_save_simvals()
881  end if
882  !
883  ! -- Return
884  return

◆ swf_gwf_calc_simvals()

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

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

657  ! -- modules
658  use constantsmodule, only: dzero
659  ! -- dummy
660  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
661  ! -- local
662  integer(I4B) :: i
663  integer(I4B) :: n1, n2
664  integer(I4B) :: ibdn1, ibdn2
665  real(DP) :: rrate
666  !
667  do i = 1, this%nexg
668  rrate = dzero
669  n1 = this%nodem1(i)
670  n2 = this%nodem2(i)
671  ibdn1 = this%swfmodel1%ibound(n1)
672  ibdn2 = this%gwfmodel2%ibound(n2)
673  if (ibdn1 /= 0 .and. ibdn2 /= 0) then
674  rrate = this%qcalc(i, n1, n2)
675  end if
676  this%simvals(i) = rrate
677  end do
678  !
679  ! -- Return
680  return

◆ 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 792 of file exg-swfgwf.f90.

793  ! -- modules
795  ! -- dummy
796  class(SwfGwfExchangeType) :: this !< GwfExchangeType
797  ! -- local
798  character(len=LENBUDTXT), dimension(1) :: budtxt
799  integer(I4B) :: n
800  integer(I4B) :: i
801  real(DP), dimension(2, 1) :: budterm
802  real(DP) :: ratin, ratout
803  real(DP) :: q
804  !
805  ! -- initialize
806  budtxt(1) = 'FLOW-JA-FACE-CHD'
807  !
808  ! -- Add the constant-head budget terms for flow from model 2 into model 1
809  if (associated(this%swfmodel1)) then
810  ratin = dzero
811  ratout = dzero
812  do i = 1, this%nexg
813  n = this%nodem1(i)
814  if (this%swfmodel1%ibound(n) < 0) then
815  q = this%simvals(i)
816  if (q > dzero) then
817  ratout = ratout + q
818  else
819  ratin = ratin - q
820  end if
821  end if
822  end do
823  budterm(1, 1) = ratin
824  budterm(2, 1) = ratout
825  call this%swfmodel1%model_bdentry(budterm, budtxt, this%name)
826  end if
827  !
828  ! -- Add the constant-head budget terms for flow from model 1 into model 2
829  if (associated(this%gwfmodel2)) then
830  ratin = dzero
831  ratout = dzero
832  do i = 1, this%nexg
833  n = this%nodem2(i)
834  if (this%gwfmodel2%ibound(n) < 0) then
835  ! -- flip flow sign as flow is relative to model 1
836  q = -this%simvals(i)
837  if (q > dzero) then
838  ratout = ratout + q
839  else
840  ratin = ratin - q
841  end if
842  end if
843  end do
844  budterm(1, 1) = ratin
845  budterm(2, 1) = ratout
846  call this%gwfmodel2%model_bdentry(budterm, budtxt, this%name)
847  end if
848  !
849  ! -- Return
850  return

◆ 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 1142 of file exg-swfgwf.f90.

1143  ! -- dummy
1144  class(SwfGwfExchangeType) :: this !< the instance of the exchange
1145  class(BaseModelType), pointer, intent(in) :: model !< the model to which the exchange might hold a connection
1146  ! -- return
1147  logical(LGP) :: is_connected !< true, when connected
1148  !
1149  is_connected = .false.
1150  select type (model)
1151  class is (gwfmodeltype)
1152  if (associated(this%gwfmodel2, model)) then
1153  is_connected = .true.
1154  end if
1155  class is (swfmodeltype)
1156  if (associated(this%swfmodel1, model)) then
1157  is_connected = .true.
1158  end if
1159  end select
1160  !
1161  ! -- Return
1162  return

◆ 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 
)
private

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

Parameters
thisSwfGwfExchangeType

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

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

◆ swf_gwf_da()

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

Deallocate memory associated with this object

Parameters
thisSwfGwfExchangeType

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

338  ! -- modules
340  ! -- dummy
341  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
342  !
343  ! -- objects
344  call this%obs%obs_da()
345  deallocate (this%obs)
346  !
347  ! -- arrays
348  call mem_deallocate(this%nodem1)
349  call mem_deallocate(this%nodem2)
350  call mem_deallocate(this%cond)
351  call mem_deallocate(this%idxglo)
352  call mem_deallocate(this%idxsymglo)
353  call mem_deallocate(this%simvals)
354  !
355  ! -- scalars
356  deallocate (this%filename)
357  call mem_deallocate(this%ipr_input)
358  call mem_deallocate(this%ipr_flow)
359  call mem_deallocate(this%nexg)
360  call mem_deallocate(this%inobs)
361  !
362  ! -- Return
363  return

◆ swf_gwf_df()

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

Define SWF to GWF exchange object.

Parameters
thisSwfGwfExchangeType

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

180  ! -- modules
181  ! -- dummy
182  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
183  ! -- local
184  !
185  ! -- log the exchange
186  write (iout, '(/a,a)') ' Creating exchange: ', this%name
187  !
188  ! -- Ensure models are in same solution
189  if (associated(this%swfmodel1) .and. associated(this%gwfmodel2)) then
190  if (this%swfmodel1%idsoln /= this%gwfmodel2%idsoln) then
191  call store_error('Two models are connected in a SWF-GWF '// &
192  'exchange but they are in different solutions. '// &
193  'Models must be in same solution: '// &
194  trim(this%swfmodel1%name)//' '// &
195  trim(this%gwfmodel2%name))
196  call this%parser%StoreErrorUnit()
197  end if
198  end if
199  !
200  ! -- source options
201  call this%source_options(iout)
202  !
203  ! -- source dimensions
204  call this%source_dimensions(iout)
205  !
206  ! -- allocate arrays
207  call this%allocate_arrays()
208  !
209  ! -- source exchange data
210  call this%source_data(iout)
211  !
212  ! -- Store obs
213  ! call this%swf_gwf_df_obs()
214  ! if (associated(this%swfmodel1)) then
215  ! call this%obs%obs_df(iout, this%name, 'SWF-GWF', this%swfmodel1%dis)
216  ! end if
217  ! !
218  ! ! -- validate
219  ! call this%validate_exchange()
220  !
221  ! -- Return
222  return
Here is the call graph for this function:

◆ 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 281 of file exg-swfgwf.f90.

282  ! -- modules
283  ! -- dummy
284  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
285  integer(I4B), intent(in) :: kiter
286  class(MatrixBaseType), pointer :: matrix_sln
287  real(DP), dimension(:), intent(inout) :: rhs_sln
288  integer(I4B), optional, intent(in) :: inwtflag
289  ! -- local
290  integer(I4B) :: i, nodem1sln, nodem2sln
291  !
292  ! -- Put this%cond into amatsln
293  do i = 1, this%nexg
294  call matrix_sln%set_value_pos(this%idxglo(i), this%cond(i))
295  call matrix_sln%set_value_pos(this%idxsymglo(i), this%cond(i))
296 
297  nodem1sln = this%nodem1(i) + this%swfmodel1%moffset
298  nodem2sln = this%nodem2(i) + this%gwfmodel2%moffset
299  call matrix_sln%add_diag_value(nodem1sln, -this%cond(i))
300  call matrix_sln%add_diag_value(nodem2sln, -this%cond(i))
301  end do
302  !
303  ! -- Return
304  return

◆ 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 254 of file exg-swfgwf.f90.

255  ! -- modules
256  use sparsemodule, only: sparsematrix
257  ! -- dummy
258  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
259  class(MatrixBaseType), pointer :: matrix_sln !< the system matrix
260  ! -- local
261  integer(I4B) :: n, iglo, jglo
262  !
263  ! -- map exchange connections
264  do n = 1, this%nexg
265  iglo = this%nodem1(n) + this%swfmodel1%moffset
266  jglo = this%nodem2(n) + this%gwfmodel2%moffset
267  this%idxglo(n) = matrix_sln%get_position(iglo, jglo)
268  this%idxsymglo(n) = matrix_sln%get_position(jglo, iglo)
269  end do
270  !
271  ! -- Return
272  return

◆ swf_gwf_ot()

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

Write output

Parameters
thisSwfGwfExchangeType

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

1047  ! -- modules
1048  use simvariablesmodule, only: iout
1049  use constantsmodule, only: dzero, linelength
1050  ! -- dummy
1051  class(SwfGwfExchangeType) :: this !< SwfGwfExchangeType
1052  ! -- local
1053  integer(I4B) :: iexg, n1, n2
1054  real(DP) :: flow
1055  character(len=LINELENGTH) :: node1str, node2str
1056  ! -- format
1057  character(len=*), parameter :: fmtheader2 = &
1058  "(/1x, 'SUMMARY OF EXCHANGE RATES FOR EXCHANGE ', a, ' WITH ID ', i0, /, &
1059  &2a16, 4a16, /, 96('-'))"
1060  character(len=*), parameter :: fmtdata = &
1061  "(2a16, 5(1pg16.6))"
1062  !
1063  ! -- Call bdsave
1064  call this%swf_gwf_bdsav()
1065  !
1066  ! -- Write a table of exchanges
1067  if (this%ipr_flow /= 0) then
1068  write (iout, fmtheader2) trim(adjustl(this%name)), this%id, 'NODEM1', &
1069  'NODEM2', 'COND', 'X_M1', 'X_M2', 'FLOW'
1070  do iexg = 1, this%nexg
1071  n1 = this%nodem1(iexg)
1072  n2 = this%nodem2(iexg)
1073  flow = this%simvals(iexg)
1074  call this%swfmodel1%dis%noder_to_string(n1, node1str)
1075  call this%gwfmodel2%dis%noder_to_string(n2, node2str)
1076  write (iout, fmtdata) trim(adjustl(node1str)), &
1077  trim(adjustl(node2str)), &
1078  this%cond(iexg), this%swfmodel1%x(n1), &
1079  this%gwfmodel2%x(n2), flow
1080  end do
1081  end if
1082  !
1083  ! -- OBS output
1084  call this%obs%obs_ot()
1085  !
1086  ! -- Return
1087  return
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
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 1094 of file exg-swfgwf.f90.

1095  ! -- modules
1097  use simvariablesmodule, only: errmsg
1098  use constantsmodule, only: dzero
1099  use observemodule, only: observetype
1100  ! -- dummy
1101  class(SwfGwfExchangeType), intent(inout) :: this
1102  ! -- local
1103  integer(I4B) :: i
1104  integer(I4B) :: j
1105  integer(I4B) :: n1
1106  integer(I4B) :: n2
1107  integer(I4B) :: iexg
1108  real(DP) :: v
1109  type(ObserveType), pointer :: obsrv => null()
1110  !
1111  ! -- Write simulated values for all gwf-gwf observations
1112  if (this%obs%npakobs > 0) then
1113  call this%obs%obs_bd_clear()
1114  do i = 1, this%obs%npakobs
1115  obsrv => this%obs%pakobs(i)%obsrv
1116  do j = 1, obsrv%indxbnds_count
1117  iexg = obsrv%indxbnds(j)
1118  v = dzero
1119  select case (obsrv%ObsTypeId)
1120  case ('FLOW-JA-FACE')
1121  n1 = this%nodem1(iexg)
1122  n2 = this%nodem2(iexg)
1123  v = this%simvals(iexg)
1124  case default
1125  errmsg = 'Unrecognized observation type: '// &
1126  trim(obsrv%ObsTypeId)
1127  call store_error(errmsg)
1128  call store_error_unit(this%inobs)
1129  end select
1130  call this%obs%SaveOneSimval(obsrv, v)
1131  end do
1132  end do
1133  end if
1134  !
1135  ! -- Return
1136  return
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:

◆ swfgwf_cr()

subroutine, public swfgwfexchangemodule::swfgwf_cr ( 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 SWF 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 96 of file exg-swfgwf.f90.

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