MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 192 of file gwt-cnc.f90.

193  ! -- modules
194  ! -- dummy
195  class(GwtCncType) :: this
196  ! -- local
197  integer(I4B) :: i, node
198  real(DP) :: cb
199  ! -- formats
200  !
201  ! -- Advance the time series
202  call this%TsManager%ad()
203  !
204  ! -- Process each entry in the constant concentration cell list
205  do i = 1, this%nbound
206  node = this%nodelist(i)
207  cb = this%conc_mult(i)
208  !
209  this%xnew(node) = cb
210  this%xold(node) = this%xnew(node)
211  end do
212  !
213  ! -- For each observation, push simulated value and corresponding
214  ! simulation time from "current" to "preceding" and reset
215  ! "current" value.
216  call this%obs%obs_ad()
217  !
218  ! -- Return
219  return

◆ 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 101 of file gwt-cnc.f90.

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

◆ cnc_bd()

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

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

347  ! -- modules
348  use tdismodule, only: delt
350  ! -- dummy
351  class(GwtCncType) :: this
352  ! -- local
353  type(BudgetType), intent(inout) :: model_budget
354  real(DP) :: ratin
355  real(DP) :: ratout
356  real(DP) :: dum
357  integer(I4B) :: isuppress_output
358  !
359  isuppress_output = 0
360  call rate_accumulator(this%ratecncin(1:this%nbound), ratin, dum)
361  call rate_accumulator(this%ratecncout(1:this%nbound), ratout, dum)
362  call model_budget%addentry(ratin, ratout, delt, this%text, &
363  isuppress_output, this%packName)
364  !
365  ! -- Return
366  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
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 514 of file gwt-cnc.f90.

515  ! -- modules
516  use constantsmodule, only: dzero
517  ! -- dummy variables
518  class(GwtCncType), intent(inout) :: this !< BndExtType object
519  integer(I4B), intent(in) :: col
520  integer(I4B), intent(in) :: row
521  ! -- result
522  real(DP) :: bndval
523  !
524  select case (col)
525  case (1)
526  bndval = this%conc_mult(row)
527  case default
528  write (errmsg, '(3a)') 'Programming error. ', &
529  & adjustl(trim(this%filtyp)), ' bound value requested column '&
530  &'outside range of ncolbnd (1).'
531  call store_error(errmsg)
532  call store_error_filename(this%input_fname)
533  end select
534  !
535  ! -- Return
536  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
Here is the call graph for this function:

◆ cnc_ck()

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

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

225  ! -- modules
226  ! -- dummy
227  class(GwtCncType), intent(inout) :: this
228  ! -- local
229  character(len=30) :: nodestr
230  integer(I4B) :: i
231  integer(I4B) :: node
232  ! -- formats
233  character(len=*), parameter :: fmtcncerr = &
234  &"('Specified dependent variable boundary ',i0, &
235  &' conc (',g0,') is less than zero for cell', a)"
236  !
237  ! -- check stress period data
238  do i = 1, this%nbound
239  node = this%nodelist(i)
240  ! -- accumulate errors
241  if (this%conc_mult(i) < dzero) then
242  call this%dis%noder_to_string(node, nodestr)
243  write (errmsg, fmt=fmtcncerr) i, this%tspvar(i), trim(nodestr)
244  call store_error(errmsg)
245  end if
246  end do
247  !
248  ! -- write summary of cnc package error messages
249  if (count_errors() > 0) then
250  call store_error_filename(this%input_fname)
251  end if
252  !
253  ! -- Return
254  return
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 280 of file gwt-cnc.f90.

281  ! -- modules
282  ! -- dummy
283  class(GwtCncType), intent(inout) :: this
284  real(DP), dimension(:), intent(in) :: x
285  real(DP), dimension(:), contiguous, intent(inout) :: flowja
286  integer(I4B), optional, intent(in) :: iadv
287  ! -- local
288  integer(I4B) :: i
289  integer(I4B) :: ipos
290  integer(I4B) :: node
291  integer(I4B) :: n2
292  integer(I4B) :: idiag
293  real(DP) :: rate
294  real(DP) :: ratein, rateout
295  real(DP) :: q
296  !
297  ! -- If no boundaries, skip flow calculations.
298  if (this%nbound > 0) then
299  !
300  ! -- Loop through each boundary calculating flow.
301  do i = 1, this%nbound
302  node = this%nodelist(i)
303  idiag = this%dis%con%ia(node)
304  rate = dzero
305  ratein = dzero
306  rateout = dzero
307  !
308  ! -- Calculate the flow rate into the cell.
309  do ipos = this%dis%con%ia(node) + 1, &
310  this%dis%con%ia(node + 1) - 1
311  q = flowja(ipos)
312  rate = rate - q
313  ! -- only accumulate chin and chout for active
314  ! connected cells
315  n2 = this%dis%con%ja(ipos)
316  if (this%ibound(n2) > 0) then
317  if (q < dzero) then
318  ratein = ratein - q
319  else
320  rateout = rateout + q
321  end if
322  end if
323  end do
324  !
325  ! -- For CNC, store total flow in rhs so it is available for other
326  ! calculations
327  this%rhs(i) = -rate
328  this%hcof(i) = dzero
329  !
330  ! -- Save simulated value to simvals array.
331  this%simvals(i) = rate
332  this%ratecncin(i) = ratein
333  this%ratecncout(i) = rateout
334  flowja(idiag) = flowja(idiag) + rate
335  !
336  end do
337  !
338  end if
339  !
340  ! -- Return
341  return

◆ 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 56 of file gwt-cnc.f90.

58  ! -- dummy
59  class(BndType), pointer :: packobj
60  integer(I4B), intent(in) :: id
61  integer(I4B), intent(in) :: ibcnum
62  integer(I4B), intent(in) :: inunit
63  integer(I4B), intent(in) :: iout
64  character(len=*), intent(in) :: namemodel
65  character(len=*), intent(in) :: pakname
66  character(len=LENVARNAME), intent(in) :: depvartype
67  character(len=*), intent(in) :: mempath
68  ! -- local
69  type(GwtCncType), pointer :: cncobj
70  !
71  ! -- allocate the object and assign values to object variables
72  allocate (cncobj)
73  packobj => cncobj
74  !
75  ! -- create name and memory path
76  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
77  packobj%text = text
78  !
79  ! -- allocate scalars
80  call cncobj%allocate_scalars()
81  !
82  ! -- initialize package
83  call packobj%pack_initialize()
84  !
85  ! -- store values
86  packobj%inunit = inunit
87  packobj%iout = iout
88  packobj%id = id
89  packobj%ibcnum = ibcnum
90  !
91  ! -- Store the appropriate label based on the dependent variable
92  cncobj%depvartype = depvartype
93  !
94  ! -- Return
95  return
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 373 of file gwt-cnc.f90.

374  ! -- modules
376  ! -- dummy
377  class(GwtCncType) :: this
378  !
379  ! -- Deallocate parent package
380  call this%BndExtType%bnd_da()
381  !
382  ! -- arrays
383  call mem_deallocate(this%ratecncin)
384  call mem_deallocate(this%ratecncout)
385  call mem_deallocate(this%tspvar, 'TSPVAR', this%memoryPath)
386  !
387  ! -- Return
388  return

◆ 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 445 of file gwt-cnc.f90.

446  ! -- dummy
447  class(GwtCncType) :: this
448  ! -- local
449  integer(I4B) :: indx
450  !
451  call this%obs%StoreObsType(this%filtyp, .true., indx)
452  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
453  !
454  ! -- Return
455  return
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 262 of file gwt-cnc.f90.

263  ! -- dummy
264  class(GwtCncType) :: this
265  real(DP), dimension(:), intent(inout) :: rhs
266  integer(I4B), dimension(:), intent(in) :: ia
267  integer(I4B), dimension(:), intent(in) :: idxglo
268  class(MatrixBaseType), pointer :: matrix_sln
269  ! -- local
270  !
271  ! -- Return
272  return

◆ 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 427 of file gwt-cnc.f90.

428  ! -- dummy
429  class(GwtCncType) :: this
430  !
431  cnc_obs_supported = .true.
432  !
433  ! -- Return
434  return

◆ cnc_rp()

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

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

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

468  ! -- dummy
469  class(GwtCncType), intent(inout) :: this
470  ! -- local
471  integer(I4B) :: i, nlinks
472  type(TimeSeriesLinkType), pointer :: tslink => null()
473  !
474  nlinks = this%TsManager%boundtslinks%Count()
475  do i = 1, nlinks
476  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
477  if (associated(tslink)) then
478  select case (tslink%JCol)
479  case (1)
480  tslink%Text = trim(this%depvartype)
481  end select
482  end if
483  end do
484  !
485  ! -- Return
486  return
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 491 of file gwt-cnc.f90.

492  ! -- modules
493  use constantsmodule, only: dzero
494  ! -- dummy variables
495  class(GwtCncType), intent(inout) :: this !< BndExtType object
496  integer(I4B), intent(in) :: row
497  ! -- result
498  real(DP) :: conc
499  !
500  if (this%iauxmultcol > 0) then
501  conc = this%tspvar(row) * this%auxvar(this%iauxmultcol, row)
502  else
503  conc = this%tspvar(row)
504  end if
505  !
506  ! -- Return
507  return

◆ 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 396 of file gwt-cnc.f90.

397  ! -- dummy
398  class(GwtCncType), intent(inout) :: this
399  !
400  ! -- create the header list label
401  this%listlabel = trim(this%filtyp)//' NO.'
402  if (this%dis%ndim == 3) then
403  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
404  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
405  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
406  elseif (this%dis%ndim == 2) then
407  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
408  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
409  else
410  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
411  end if
412  write (this%listlabel, '(a, a16)') trim(this%listlabel), &
413  trim(this%depvartype)
414  if (this%inamedbound == 1) then
415  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
416  end if
417  !
418  ! -- Return
419  return

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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