40   subroutine chfgwf_cr(filename, name, id, m1_id, m2_id, input_mempath)
 
   43     character(len=*), 
intent(in) :: filename
 
   44     character(len=*) :: name
 
   45     integer(I4B), 
intent(in) :: id
 
   46     integer(I4B), 
intent(in) :: m1_id
 
   47     integer(I4B), 
intent(in) :: m2_id
 
   48     character(len=*), 
intent(in) :: input_mempath
 
   55     baseexchange => exchange
 
   58     call exchange%initialize(filename, name, 
'CHF', id, m1_id, m2_id, &
 
   72     write (
iout, 
'(/a,a)') 
' Creating exchange: ', this%name
 
   75     if (
associated(this%swfmodel) .and. 
associated(this%gwfmodel)) 
then 
   76       if (this%swfmodel%idsoln /= this%gwfmodel%idsoln) 
then 
   77         call store_error(
'Two models are connected in a SWF-GWF '// &
 
   78                          'exchange but they are in different solutions. '// &
 
   79                          'Models must be in same solution: '// &
 
   80                          trim(this%swfmodel%name)//
' '// &
 
   81                          trim(this%gwfmodel%name))
 
   87     call this%source_options(
iout)
 
   90     call this%source_dimensions(
iout)
 
   93     call this%allocate_arrays()
 
   96     call this%source_data(
iout)
 
  121     integer(I4B), 
intent(in) :: iout
 
  127                        this%input_mempath, found%ipr_input)
 
  129                        this%input_mempath, found%ipr_flow)
 
  131                        this%input_mempath, found%ifixedcond)
 
  133     write (iout, 
'(1x,a)') 
'Processing CHF-GWF exchange options' 
  135     if (found%ipr_input) 
then 
  136       write (iout, 
'(4x,a)') &
 
  137         'The list of exchanges will be printed.' 
  140     if (found%ipr_flow) 
then 
  141       write (iout, 
'(4x,a)') &
 
  142         'Exchange flows will be printed to list files.' 
  145     if (found%ifixedcond) 
then 
  146       write (iout, 
'(4x,a)') &
 
  147         'Conductance is fixed as product of BEDLEAK and CFACT.' 
  160     write (iout, 
'(1x,a)') 
'End of CHF-GWF exchange options' 
  172     integer(I4B), 
intent(in) :: iout
 
  177     call mem_set_value(this%nexg, 
'NEXG', this%input_mempath, found%nexg)
 
  179     write (iout, 
'(1x,a)') 
'PROCESSING EXCHANGE DIMENSIONS' 
  182       write (iout, 
'(4x,a,i0)') 
'NEXG = ', this%nexg
 
  185     write (iout, 
'(1x,a)') 
'END OF EXCHANGE DIMENSIONS' 
  196     integer(I4B), 
intent(in) :: iout
 
  198     integer(I4B), 
dimension(:, :), 
contiguous, 
pointer :: cellidm1
 
  199     integer(I4B), 
dimension(:, :), 
contiguous, 
pointer :: cellidm2
 
  200     real(DP), 
dimension(:), 
contiguous, 
pointer :: bedleak
 
  201     real(DP), 
dimension(:), 
contiguous, 
pointer :: cfact
 
  202     character(len=20) :: cellstr1, cellstr2
 
  204     integer(I4B) :: iexg, nodeswf, nodegwf
 
  206     character(len=*), 
parameter :: fmtexglabel = 
"(1x, 3a10, 50(a16))" 
  207     character(len=*), 
parameter :: fmtexgdata = &
 
  208                                    "(5x, a, 1x, a ,50(1pg16.6))" 
  210     call mem_setptr(cellidm1, 
'CELLIDM1', this%input_mempath)
 
  211     call mem_setptr(cellidm2, 
'CELLIDM2', this%input_mempath)
 
  212     call mem_setptr(bedleak, 
'BEDLEAK', this%input_mempath)
 
  213     call mem_setptr(cfact, 
'CFACT', this%input_mempath)
 
  215     write (iout, 
'(1x,a)') 
'PROCESSING EXCHANGEDATA' 
  217     if (this%ipr_input /= 0) 
then 
  218       write (iout, fmtexglabel) 
'NODEM1', 
'NODEM2', 
'BEDLEAK', 
'CFACT' 
  221     do iexg = 1, this%nexg
 
  223       if (
associated(this%model1)) 
then 
  225         nodeswf = this%noder(this%model1, cellidm1(:, iexg), iout)
 
  226         this%nodeswf(iexg) = nodeswf
 
  228         this%nodeswf(iexg) = -1
 
  231       if (
associated(this%model2)) 
then 
  233         nodegwf = this%noder(this%model2, cellidm2(:, iexg), iout)
 
  234         this%nodegwf(iexg) = nodegwf
 
  236         this%nodegwf(iexg) = -1
 
  240       this%bedleak(iexg) = bedleak(iexg)
 
  241       this%cfact(iexg) = cfact(iexg)
 
  244       if (this%ipr_input /= 0) 
then 
  245         cellstr1 = this%cellstr(this%model1, cellidm1(:, iexg), iout)
 
  246         cellstr2 = this%cellstr(this%model2, cellidm2(:, iexg), iout)
 
  247         write (iout, fmtexgdata) trim(cellstr1), trim(cellstr2), &
 
  248           this%bedleak(iexg), this%cfact(iexg)
 
  252       if (
associated(this%model1)) 
then 
  253         if (nodeswf <= 0) 
then 
  254           cellstr1 = this%cellstr(this%model1, cellidm1(:, iexg), iout)
 
  256             trim(adjustl(this%model1%name))// &
 
  257             ' Cell is outside active grid domain ('// &
 
  258             trim(adjustl(cellstr1))//
').' 
  264       if (
associated(this%model2)) 
then 
  265         if (nodegwf <= 0) 
then 
  266           cellstr2 = this%cellstr(this%model2, cellidm2(:, iexg), iout)
 
  268             trim(adjustl(this%model2%name))// &
 
  269             ' Cell is outside active grid domain ('// &
 
  270             trim(adjustl(cellstr2))//
').' 
  276     write (iout, 
'(1x,a)') 
'END OF EXCHANGEDATA' 
  281       call store_error(
'Errors encountered in exchange input file.')
 
subroutine, public addbaseexchangetolist(list, exchange)
Add the exchange object (BaseExchangeType) to a list.
class(basemodeltype) function, pointer, public getbasemodelfromlist(list, idx)
This module contains the ChfGwfExchangeModule Module.
subroutine source_dimensions(this, iout)
Source dimension from input context.
subroutine, public chfgwf_cr(filename, name, id, m1_id, m2_id, input_mempath)
@ brief Create CHF GWF exchange
subroutine source_options(this, iout)
@ brief Source options
subroutine source_data(this, iout)
Source exchange data from input context.
subroutine chf_gwf_df(this)
@ brief Define CHF GWF exchange
This module contains simulation constants.
integer(i4b), parameter lenvarname
maximum length of a variable name
real(dp), parameter dem6
real constant 1e-6
This module defines variable data types.
type(listtype), public basemodellist
type(listtype), public baseexchangelist
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived type ObsType.
subroutine, public obs_cr(obs, inobs)
@ brief Create a new ObsType object
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
integer(i4b) function, public count_errors()
Return number of errors.
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
integer(i4b), dimension(:), allocatable model_loc_idx
equals the local index into the basemodel list (-1 when not available)
integer(i4b) iout
file unit number for simulation output
This module contains the SourceCommonModule.
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
This module contains the SwfGwfExchangeModule Module.
Surface Water Flow (SWF) Module.
Highest level model type. All models extend this parent type.
This class is used to store a single deferred-length character string. It was designed to work in an ...