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

Data Types

type  gwtcnctype
 

Functions/Subroutines

subroutine, public cnc_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, depvartype, mempath)
 Create a new constant concentration or temperature package. More...
 
subroutine cnc_allocate_arrays (this, nodelist, auxvar)
 Allocate arrays specific to the constant concentration/tempeature package. More...
 
subroutine cnc_rp (this)
 Constant concentration/temperature read and prepare (rp) routine. More...
 
subroutine cnc_ad (this)
 Constant concentration/temperature package advance routine. More...
 
subroutine cnc_ck (this)
 Check constant concentration/temperature boundary condition data. More...
 
subroutine cnc_fc (this, rhs, ia, idxglo, matrix_sln)
 Override bnd_fc and do nothing. More...
 
subroutine cnc_cq (this, x, flowja, iadv)
 Calculate flow associated with constant concentration/temperature boundary. More...
 
subroutine cnc_bd (this, model_budget)
 Add package ratin/ratout to model budget. More...
 
subroutine cnc_da (this)
 Deallocate memory. More...
 
subroutine define_listlabel (this)
 Define labels used in list file. More...
 
logical function cnc_obs_supported (this)
 Procedure related to observation processing. More...
 
subroutine cnc_df_obs (this)
 Procedure related to observation processing. More...
 
subroutine cnc_rp_ts (this)
 Procedure related to time series. More...
 
real(dp) function conc_mult (this, row)
 Apply auxiliary multiplier to specified concentration if. More...
 
real(dp) function cnc_bound_value (this, col, row)
 @ brief Return a bound value More...
 

Variables

character(len=lenftype) ftype = 'CNC'
 
character(len=lenpackagename) text = ' CNC'
 

Function/Subroutine Documentation

◆ cnc_ad()

subroutine gwtcncmodule::cnc_ad ( class(gwtcnctype this)

Add package connections to matrix

Definition at line 182 of file gwt-cnc.f90.

183  ! -- modules
184  ! -- dummy
185  class(GwtCncType) :: this
186  ! -- local
187  integer(I4B) :: i, node
188  real(DP) :: cb
189  ! -- formats
190  !
191  ! -- Advance the time series
192  call this%TsManager%ad()
193  !
194  ! -- Process each entry in the constant concentration cell list
195  do i = 1, this%nbound
196  node = this%nodelist(i)
197  cb = this%conc_mult(i)
198  !
199  this%xnew(node) = cb
200  this%xold(node) = this%xnew(node)
201  end do
202  !
203  ! -- For each observation, push simulated value and corresponding
204  ! simulation time from "current" to "preceding" and reset
205  ! "current" value.
206  call this%obs%obs_ad()

◆ cnc_allocate_arrays()

subroutine gwtcncmodule::cnc_allocate_arrays ( class(gwtcnctype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)
private

Definition at line 97 of file gwt-cnc.f90.

98  ! -- modules
100  ! -- dummy
101  class(GwtCncType) :: this
102  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
103  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
104  ! -- local
105  integer(I4B) :: i
106  !
107  ! -- call standard BndType allocate scalars
108  call this%BndExtType%allocate_arrays(nodelist, auxvar)
109  !
110  ! -- allocate ratecncex
111  call mem_allocate(this%ratecncin, this%maxbound, 'RATECNCIN', this%memoryPath)
112  call mem_allocate(this%ratecncout, this%maxbound, 'RATECNCOUT', &
113  this%memoryPath)
114  do i = 1, this%maxbound
115  this%ratecncin(i) = dzero
116  this%ratecncout(i) = dzero
117  end do
118  ! -- set constant head array input context pointer
119  call mem_setptr(this%tspvar, 'TSPVAR', this%input_mempath)
120  !
121  ! -- checkin constant head array input context pointer
122  call mem_checkin(this%tspvar, 'TSPVAR', this%memoryPath, &
123  'TSPVAR', this%input_mempath)
124  !

◆ cnc_bd()

subroutine gwtcncmodule::cnc_bd ( class(gwtcnctype this,
type(budgettype), intent(inout)  model_budget 
)
private

Definition at line 324 of file gwt-cnc.f90.

325  ! -- modules
326  use tdismodule, only: delt
328  ! -- dummy
329  class(GwtCncType) :: this
330  ! -- local
331  type(BudgetType), intent(inout) :: model_budget
332  real(DP) :: ratin
333  real(DP) :: ratout
334  real(DP) :: dum
335  integer(I4B) :: isuppress_output
336  !
337  isuppress_output = 0
338  call rate_accumulator(this%ratecncin(1:this%nbound), ratin, dum)
339  call rate_accumulator(this%ratecncout(1:this%nbound), ratout, dum)
340  call model_budget%addentry(ratin, ratout, delt, this%text, &
341  isuppress_output, this%packName)
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
Derived type for the Budget object.
Definition: Budget.f90:39
Here is the call graph for this function:

◆ cnc_bound_value()

real(dp) function gwtcncmodule::cnc_bound_value ( class(gwtcnctype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)

Return a bound value associated with an ncolbnd index and row.

Parameters
[in,out]thisBndExtType object

Definition at line 471 of file gwt-cnc.f90.

472  ! -- modules
473  use constantsmodule, only: dzero
474  ! -- dummy variables
475  class(GwtCncType), intent(inout) :: this !< BndExtType object
476  integer(I4B), intent(in) :: col
477  integer(I4B), intent(in) :: row
478  ! -- result
479  real(DP) :: bndval
480  !
481  select case (col)
482  case (1)
483  bndval = this%conc_mult(row)
484  case default
485  write (errmsg, '(3a)') 'Programming error. ', &
486  & adjustl(trim(this%filtyp)), ' bound value requested column '&
487  &'outside range of ncolbnd (1).'
488  call store_error(errmsg)
489  call store_error_filename(this%input_fname)
490  end select
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Here is the call graph for this function:

◆ cnc_ck()

subroutine gwtcncmodule::cnc_ck ( class(gwtcnctype), intent(inout)  this)
private

Definition at line 211 of file gwt-cnc.f90.

212  ! -- modules
213  ! -- dummy
214  class(GwtCncType), intent(inout) :: this
215  ! -- local
216  character(len=30) :: nodestr
217  integer(I4B) :: i
218  integer(I4B) :: node
219  ! -- formats
220  character(len=*), parameter :: fmtcncerr = &
221  &"('Specified dependent variable boundary ',i0, &
222  &' conc (',g0,') is less than zero for cell', a)"
223  !
224  ! -- check stress period data
225  do i = 1, this%nbound
226  node = this%nodelist(i)
227  ! -- accumulate errors
228  if (this%conc_mult(i) < dzero) then
229  call this%dis%noder_to_string(node, nodestr)
230  write (errmsg, fmt=fmtcncerr) i, this%tspvar(i), trim(nodestr)
231  call store_error(errmsg)
232  end if
233  end do
234  !
235  ! -- write summary of cnc package error messages
236  if (count_errors() > 0) then
237  call store_error_filename(this%input_fname)
238  end if
Here is the call graph for this function:

◆ cnc_cq()

subroutine gwtcncmodule::cnc_cq ( class(gwtcnctype), intent(inout)  this,
real(dp), dimension(:), intent(in)  x,
real(dp), dimension(:), intent(inout), contiguous  flowja,
integer(i4b), intent(in), optional  iadv 
)
private

This method overrides bnd_cq()

Definition at line 261 of file gwt-cnc.f90.

262  ! -- modules
263  ! -- dummy
264  class(GwtCncType), intent(inout) :: this
265  real(DP), dimension(:), intent(in) :: x
266  real(DP), dimension(:), contiguous, intent(inout) :: flowja
267  integer(I4B), optional, intent(in) :: iadv
268  ! -- local
269  integer(I4B) :: i
270  integer(I4B) :: ipos
271  integer(I4B) :: node
272  integer(I4B) :: n2
273  integer(I4B) :: idiag
274  real(DP) :: rate
275  real(DP) :: ratein, rateout
276  real(DP) :: q
277  !
278  ! -- If no boundaries, skip flow calculations.
279  if (this%nbound > 0) then
280  !
281  ! -- Loop through each boundary calculating flow.
282  do i = 1, this%nbound
283  node = this%nodelist(i)
284  idiag = this%dis%con%ia(node)
285  rate = dzero
286  ratein = dzero
287  rateout = dzero
288  !
289  ! -- Calculate the flow rate into the cell.
290  do ipos = this%dis%con%ia(node) + 1, &
291  this%dis%con%ia(node + 1) - 1
292  q = flowja(ipos)
293  rate = rate - q
294  ! -- only accumulate chin and chout for active
295  ! connected cells
296  n2 = this%dis%con%ja(ipos)
297  if (this%ibound(n2) > 0) then
298  if (q < dzero) then
299  ratein = ratein - q
300  else
301  rateout = rateout + q
302  end if
303  end if
304  end do
305  !
306  ! -- For CNC, store total flow in rhs so it is available for other
307  ! calculations
308  this%rhs(i) = -rate
309  this%hcof(i) = dzero
310  !
311  ! -- Save simulated value to simvals array.
312  this%simvals(i) = rate
313  this%ratecncin(i) = ratein
314  this%ratecncout(i) = rateout
315  flowja(idiag) = flowja(idiag) + rate
316  !
317  end do
318  !
319  end if

◆ cnc_create()

subroutine, public gwtcncmodule::cnc_create ( class(bndtype), pointer  packobj,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  ibcnum,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  namemodel,
character(len=*), intent(in)  pakname,
character(len=lenvarname), intent(in)  depvartype,
character(len=*), intent(in)  mempath 
)

Routine points packobj to the newly created package

Definition at line 55 of file gwt-cnc.f90.

57  ! -- dummy
58  class(BndType), pointer :: packobj
59  integer(I4B), intent(in) :: id
60  integer(I4B), intent(in) :: ibcnum
61  integer(I4B), intent(in) :: inunit
62  integer(I4B), intent(in) :: iout
63  character(len=*), intent(in) :: namemodel
64  character(len=*), intent(in) :: pakname
65  character(len=LENVARNAME), intent(in) :: depvartype
66  character(len=*), intent(in) :: mempath
67  ! -- local
68  type(GwtCncType), pointer :: cncobj
69  !
70  ! -- allocate the object and assign values to object variables
71  allocate (cncobj)
72  packobj => cncobj
73  !
74  ! -- create name and memory path
75  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
76  packobj%text = text
77  !
78  ! -- allocate scalars
79  call cncobj%allocate_scalars()
80  !
81  ! -- initialize package
82  call packobj%pack_initialize()
83  !
84  ! -- store values
85  packobj%inunit = inunit
86  packobj%iout = iout
87  packobj%id = id
88  packobj%ibcnum = ibcnum
89  !
90  ! -- Store the appropriate label based on the dependent variable
91  cncobj%depvartype = depvartype
Here is the caller graph for this function:

◆ cnc_da()

subroutine gwtcncmodule::cnc_da ( class(gwtcnctype this)

Method to deallocate memory for the package.

Definition at line 348 of file gwt-cnc.f90.

349  ! -- modules
351  ! -- dummy
352  class(GwtCncType) :: this
353  !
354  ! -- Deallocate parent package
355  call this%BndExtType%bnd_da()
356  !
357  ! -- arrays
358  call mem_deallocate(this%ratecncin)
359  call mem_deallocate(this%ratecncout)
360  call mem_deallocate(this%tspvar, 'TSPVAR', this%memoryPath)

◆ cnc_df_obs()

subroutine gwtcncmodule::cnc_df_obs ( class(gwtcnctype this)
private

This routine:

  • defines observations
  • stores observation types supported by either of the SDV packages (CNC or CNT),
  • overrides BndExtTypebnd_df_obs

Definition at line 411 of file gwt-cnc.f90.

412  ! -- dummy
413  class(GwtCncType) :: this
414  ! -- local
415  integer(I4B) :: indx
416  !
417  call this%obs%StoreObsType(this%filtyp, .true., indx)
418  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ cnc_fc()

subroutine gwtcncmodule::cnc_fc ( class(gwtcnctype this,
real(dp), dimension(:), intent(inout)  rhs,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), dimension(:), intent(in)  idxglo,
class(matrixbasetype), pointer  matrix_sln 
)
private

For constant concentration/temperature boundary type, the call to bnd_fc needs to be overwritten to prevent logic found in bnd from being executed

Definition at line 246 of file gwt-cnc.f90.

247  ! -- dummy
248  class(GwtCncType) :: this
249  real(DP), dimension(:), intent(inout) :: rhs
250  integer(I4B), dimension(:), intent(in) :: ia
251  integer(I4B), dimension(:), intent(in) :: idxglo
252  class(MatrixBaseType), pointer :: matrix_sln
253  ! -- local

◆ cnc_obs_supported()

logical function gwtcncmodule::cnc_obs_supported ( class(gwtcnctype this)
private

This routine:

  • returns true because the SDV package supports observations,
  • overrides packagetype_obs_supported()

Definition at line 396 of file gwt-cnc.f90.

397  ! -- dummy
398  class(GwtCncType) :: this
399  !
400  cnc_obs_supported = .true.

◆ cnc_rp()

subroutine gwtcncmodule::cnc_rp ( class(gwtcnctype), intent(inout)  this)

Definition at line 129 of file gwt-cnc.f90.

130  ! -- modules
131  use simmodule, only: store_error
132  use inputoutputmodule, only: lowcase
133  implicit none
134  ! -- dummy
135  class(GwtCncType), intent(inout) :: this
136  ! -- local
137  integer(I4B) :: i, node, ibd, ierr
138  character(len=30) :: nodestr
139  character(len=LENVARNAME) :: dvtype
140  !
141  ! -- Reset previous CNCs to active cell
142  do i = 1, this%nbound
143  node = this%nodelist(i)
144  this%ibound(node) = this%ibcnum
145  end do
146  !
147  ! -- Call the parent class read and prepare
148  call this%BndExtType%bnd_rp()
149  !
150  ! -- Set ibound to -(ibcnum + 1) for constant concentration cells
151  ierr = 0
152  do i = 1, this%nbound
153  node = this%nodelist(i)
154  ibd = this%ibound(node)
155  if (ibd < 0) then
156  call this%dis%noder_to_string(node, nodestr)
157  dvtype = trim(this%depvartype)
158  call lowcase(dvtype)
159  call store_error('Cell is already a constant ' &
160  //dvtype//': '//trim(adjustl(nodestr)))
161  ierr = ierr + 1
162  else
163  this%ibound(node) = -this%ibcnum
164  end if
165  end do
166  !
167  ! -- Stop if errors detected
168  if (ierr > 0) then
169  call store_error_filename(this%input_fname)
170  end if
171  !
172  ! -- Write the list to iout if requested
173  if (this%iprpak /= 0) then
174  call this%write_list()
175  end if
subroutine, public lowcase(word)
Convert to lower case.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:

◆ cnc_rp_ts()

subroutine gwtcncmodule::cnc_rp_ts ( class(gwtcnctype), intent(inout)  this)
private

Assign tsLinkText appropriately for all time series in use by package. For any specified dependent variable package, for example either the constant concentration or constant temperature packages, the dependent variable can be controlled by time series.

Definition at line 430 of file gwt-cnc.f90.

431  ! -- dummy
432  class(GwtCncType), intent(inout) :: this
433  ! -- local
434  integer(I4B) :: i, nlinks
435  type(TimeSeriesLinkType), pointer :: tslink => null()
436  !
437  nlinks = this%TsManager%boundtslinks%Count()
438  do i = 1, nlinks
439  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
440  if (associated(tslink)) then
441  select case (tslink%JCol)
442  case (1)
443  tslink%Text = trim(this%depvartype)
444  end select
445  end if
446  end do
Here is the call graph for this function:

◆ conc_mult()

real(dp) function gwtcncmodule::conc_mult ( class(gwtcnctype), intent(inout)  this,
integer(i4b), intent(in)  row 
)
private
Parameters
[in,out]thisBndExtType object

Definition at line 451 of file gwt-cnc.f90.

452  ! -- modules
453  use constantsmodule, only: dzero
454  ! -- dummy variables
455  class(GwtCncType), intent(inout) :: this !< BndExtType object
456  integer(I4B), intent(in) :: row
457  ! -- result
458  real(DP) :: conc
459  !
460  if (this%iauxmultcol > 0) then
461  conc = this%tspvar(row) * this%auxvar(this%iauxmultcol, row)
462  else
463  conc = this%tspvar(row)
464  end if

◆ define_listlabel()

subroutine gwtcncmodule::define_listlabel ( class(gwtcnctype), intent(inout)  this)

Define the list heading that is written to iout when PRINT_INPUT option is used.

Definition at line 368 of file gwt-cnc.f90.

369  ! -- dummy
370  class(GwtCncType), intent(inout) :: this
371  !
372  ! -- create the header list label
373  this%listlabel = trim(this%filtyp)//' NO.'
374  if (this%dis%ndim == 3) then
375  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
376  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
377  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
378  elseif (this%dis%ndim == 2) then
379  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
380  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
381  else
382  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
383  end if
384  write (this%listlabel, '(a, a16)') trim(this%listlabel), &
385  trim(this%depvartype)
386  if (this%inamedbound == 1) then
387  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
388  end if

Variable Documentation

◆ ftype

character(len=lenftype) gwtcncmodule::ftype = 'CNC'
private

Definition at line 21 of file gwt-cnc.f90.

21  character(len=LENFTYPE) :: ftype = 'CNC'

◆ text

character(len=lenpackagename) gwtcncmodule::text = ' CNC'
private

Definition at line 22 of file gwt-cnc.f90.

22  character(len=LENPACKAGENAME) :: text = ' CNC'