MODFLOW 6  version 6.8.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 178 of file gwt-cnc.f90.

179  ! -- modules
180  ! -- dummy
181  class(GwtCncType) :: this
182  ! -- local
183  integer(I4B) :: i, node
184  real(DP) :: cb
185  ! -- formats
186  !
187  ! -- Advance the time series
188  call this%TsManager%ad()
189  !
190  ! -- Process each entry in the constant concentration cell list
191  do i = 1, this%nbound
192  node = this%nodelist(i)
193  cb = this%conc_mult(i)
194  !
195  this%xnew(node) = cb
196  this%xold(node) = this%xnew(node)
197  end do
198  !
199  ! -- For each observation, push simulated value and corresponding
200  ! simulation time from "current" to "preceding" and reset
201  ! "current" value.
202  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 98 of file gwt-cnc.f90.

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

◆ cnc_bd()

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

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

321  ! -- modules
322  use tdismodule, only: delt
324  ! -- dummy
325  class(GwtCncType) :: this
326  ! -- local
327  type(BudgetType), intent(inout) :: model_budget
328  real(DP) :: ratin
329  real(DP) :: ratout
330  real(DP) :: dum
331  integer(I4B) :: isuppress_output
332  !
333  isuppress_output = 0
334  call rate_accumulator(this%ratecncin(1:this%nbound), ratin, dum)
335  call rate_accumulator(this%ratecncout(1:this%nbound), ratout, dum)
336  call model_budget%addentry(ratin, ratout, delt, this%text, &
337  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.

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

468  ! -- modules
469  use constantsmodule, only: dzero
470  ! -- dummy variables
471  class(GwtCncType), intent(inout) :: this
472  integer(I4B), intent(in) :: col
473  integer(I4B), intent(in) :: row
474  ! -- result
475  real(DP) :: bndval
476  !
477  select case (col)
478  case (1)
479  bndval = this%conc_mult(row)
480  case default
481  write (errmsg, '(3a)') 'Programming error. ', &
482  & adjustl(trim(this%filtyp)), ' bound value requested column '&
483  &'outside range of ncolbnd (1).'
484  call store_error(errmsg)
485  call store_error_filename(this%input_fname)
486  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 207 of file gwt-cnc.f90.

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

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

345  ! -- modules
347  ! -- dummy
348  class(GwtCncType) :: this
349  !
350  ! -- Deallocate parent package
351  call this%BndExtType%bnd_da()
352  !
353  ! -- arrays
354  call mem_deallocate(this%ratecncin)
355  call mem_deallocate(this%ratecncout)
356  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 407 of file gwt-cnc.f90.

408  ! -- dummy
409  class(GwtCncType) :: this
410  ! -- local
411  integer(I4B) :: indx
412  !
413  call this%obs%StoreObsType(this%filtyp, .true., indx)
414  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 242 of file gwt-cnc.f90.

243  ! -- dummy
244  class(GwtCncType) :: this
245  real(DP), dimension(:), intent(inout) :: rhs
246  integer(I4B), dimension(:), intent(in) :: ia
247  integer(I4B), dimension(:), intent(in) :: idxglo
248  class(MatrixBaseType), pointer :: matrix_sln
249  ! -- 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 392 of file gwt-cnc.f90.

393  ! -- dummy
394  class(GwtCncType) :: this
395  !
396  cnc_obs_supported = .true.

◆ cnc_rp()

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

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

131  ! -- modules
132  use simmodule, only: store_error
133  use inputoutputmodule, only: lowcase
134  implicit none
135  ! -- dummy
136  class(GwtCncType), intent(inout) :: this
137  ! -- local
138  integer(I4B) :: i, node, ibd, ierr
139  character(len=30) :: nodestr
140  character(len=LENVARNAME) :: dvtype
141  !
142  ! -- Reset previous CNCs to active cell
143  do i = 1, this%nbound
144  node = this%nodelist(i)
145  this%ibound(node) = this%ibcnum
146  end do
147  !
148  ! -- Call the parent class read and prepare
149  call this%BndExtType%bnd_rp()
150  !
151  ! -- Set ibound to -(ibcnum + 1) for constant concentration cells
152  ierr = 0
153  do i = 1, this%nbound
154  node = this%nodelist(i)
155  ibd = this%ibound(node)
156  if (ibd < 0) then
157  call this%dis%noder_to_string(node, nodestr)
158  dvtype = trim(this%depvartype)
159  call lowcase(dvtype)
160  call store_error('Cell is already a constant ' &
161  //dvtype//': '//trim(adjustl(nodestr)))
162  ierr = ierr + 1
163  else
164  this%ibound(node) = -this%ibcnum
165  end if
166  end do
167  !
168  ! -- Stop if errors detected
169  if (ierr > 0) then
170  call store_error_filename(this%input_fname)
171  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 426 of file gwt-cnc.f90.

427  ! -- dummy
428  class(GwtCncType), intent(inout) :: this
429  ! -- local
430  integer(I4B) :: i, nlinks
431  type(TimeSeriesLinkType), pointer :: tslink => null()
432  !
433  nlinks = this%TsManager%boundtslinks%Count()
434  do i = 1, nlinks
435  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
436  if (associated(tslink)) then
437  select case (tslink%JCol)
438  case (1)
439  tslink%Text = trim(this%depvartype)
440  end select
441  end if
442  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

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

448  ! -- modules
449  use constantsmodule, only: dzero
450  ! -- dummy variables
451  class(GwtCncType), intent(inout) :: this
452  integer(I4B), intent(in) :: row
453  ! -- result
454  real(DP) :: conc
455  !
456  if (this%iauxmultcol > 0) then
457  conc = this%tspvar(row) * this%auxvar(this%iauxmultcol, row)
458  else
459  conc = this%tspvar(row)
460  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 364 of file gwt-cnc.f90.

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