MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
tsp-apt.f90
Go to the documentation of this file.
1 ! -- Advanced Package Transport Module
2 ! -- This module contains most of the routines for simulating transport
3 ! -- through the advanced packages.
4 ! -- Future work:
5 ! * support decay, sorption
6 ! * dispersion in SFT and UZT?
7 !
8 ! AFP flows (flowbudptr) index var ATP term Transport Type
9 !---------------------------------------------------------------------------------
10 
11 ! -- specialized terms in the flow budget
12 ! FLOW-JA-FACE idxbudfjf FLOW-JA-FACE cv2cv
13 ! GWF (aux FLOW-AREA) idxbudgwf GWF cv2gwf
14 ! STORAGE (aux VOLUME) idxbudsto none used for cv volumes
15 ! FROM-MVR idxbudfmvr FROM-MVR q * cext = this%qfrommvr(:) ! rhow*cpw is applied to various terms for heat transport
16 ! TO-MVR idxbudtmvr TO-MVR q * cfeat
17 
18 ! -- generalized source/sink terms (except ET?)
19 ! RAINFALL idxbudrain RAINFALL q * crain
20 ! EVAPORATION idxbudevap EVAPORATION cfeat<cevap: q*cfeat, else: q*cevap ! latent heat may be applied for evaporative cooling
21 ! RUNOFF idxbudroff RUNOFF q * croff
22 ! EXT-INFLOW idxbudiflw EXT-INFLOW q * ciflw
23 ! WITHDRAWAL idxbudwdrl WITHDRAWAL q * cfeat
24 ! EXT-OUTFLOW idxbudoutf EXT-OUTFLOW q * cfeat
25 
26 ! -- terms from a flow file that should be skipped
27 ! CONSTANT none none none
28 ! AUXILIARY none none none
29 
30 ! -- terms that are written to the transport budget file
31 ! none none STORAGE (aux MASS) dM/dt
32 ! none none AUXILIARY none
33 ! none none CONSTANT accumulate
34 !
35 !
37 
38  use kindmodule, only: dp, i4b, lgp
45  use simvariablesmodule, only: errmsg
46  use bndmodule, only: bndtype
47  use tspfmimodule, only: tspfmitype
50  use tablemodule, only: tabletype, table_cr
51  use observemodule, only: observetype
53  use basedismodule, only: disbasetype
55 
56  implicit none
57 
58  public :: tspapttype
59  public :: apt_process_obsid
60  public :: apt_process_obsid12
61 
62  character(len=LENFTYPE) :: ftype = 'APT'
63  character(len=LENVARNAME) :: text = ' APT'
64 
65  type, extends(bndtype) :: tspapttype
66 
67  character(len=LENPACKAGENAME) :: flowpackagename = '' !< name of corresponding flow package
68  character(len=8), &
69  dimension(:), pointer, contiguous :: status => null() !< active, inactive, constant
70  character(len=LENAUXNAME) :: cauxfpconc = '' !< name of aux column in flow package auxvar array for concentration (or temperature)
71  integer(I4B), pointer :: iauxfpconc => null() !< column in flow package bound array to insert concs
72  integer(I4B), pointer :: imatrows => null() !< if active, add new rows to matrix
73  integer(I4B), pointer :: iprconc => null() !< print conc to listing file
74  integer(I4B), pointer :: iconcout => null() !< unit number for conc output file
75  integer(I4B), pointer :: ibudgetout => null() !< unit number for budget output file
76  integer(I4B), pointer :: ibudcsv => null() !< unit number for csv budget output file
77  integer(I4B), pointer :: ncv => null() !< number of control volumes
78  integer(I4B), pointer :: igwfaptpak => null() !< package number of corresponding this package
79  integer(I4B), pointer :: idxprepak => null() !< budget-object index that precedes package-specific budget objects
80  integer(I4B), pointer :: idxlastpak => null() !< budget-object index of last package-specific budget object
81  real(dp), dimension(:), pointer, contiguous :: strt => null() !< starting feature concentration (or temperature)
82  integer(I4B), dimension(:), pointer, contiguous :: idxlocnode => null() !< map position in global rhs and x array of pack entry
83  integer(I4B), dimension(:), pointer, contiguous :: idxpakdiag => null() !< map diag position of feature in global amat
84  integer(I4B), dimension(:), pointer, contiguous :: idxdglo => null() !< map position in global array of package diagonal row entries
85  integer(I4B), dimension(:), pointer, contiguous :: idxoffdglo => null() !< map position in global array of package off diagonal row entries
86  integer(I4B), dimension(:), pointer, contiguous :: idxsymdglo => null() !< map position in global array of package diagonal entries to model rows
87  integer(I4B), dimension(:), pointer, contiguous :: idxsymoffdglo => null() !< map position in global array of package off diagonal entries to model rows
88  integer(I4B), dimension(:), pointer, contiguous :: idxfjfdglo => null() !< map diagonal feature to feature in global amat
89  integer(I4B), dimension(:), pointer, contiguous :: idxfjfoffdglo => null() !< map off diagonal feature to feature in global amat
90  integer(I4B), dimension(:), pointer, contiguous :: iboundpak => null() !< package ibound
91  real(dp), dimension(:), pointer, contiguous :: xnewpak => null() !< feature concentration (or temperature) for current time step
92  real(dp), dimension(:), pointer, contiguous :: xoldpak => null() !< feature concentration (or temperature) from previous time step
93  real(dp), dimension(:), pointer, contiguous :: dbuff => null() !< temporary storage array
94  character(len=LENBOUNDNAME), &
95  dimension(:), pointer, contiguous :: featname => null()
96  real(dp), dimension(:), pointer, contiguous :: concfeat => null() !< concentration (or temperature) of the feature
97  real(dp), dimension(:, :), pointer, contiguous :: lauxvar => null() !< auxiliary variable
98  type(tspfmitype), pointer :: fmi => null() !< pointer to fmi object
99  real(dp), dimension(:), pointer, contiguous :: qsto => null() !< mass (or energy) flux due to storage change
100  real(dp), dimension(:), pointer, contiguous :: ccterm => null() !< mass (or energy) flux required to maintain constant concentration (or temperature)
101  integer(I4B), pointer :: idxbudfjf => null() !< index of flow ja face in flowbudptr
102  integer(I4B), pointer :: idxbudgwf => null() !< index of gwf terms in flowbudptr
103  integer(I4B), pointer :: idxbudsto => null() !< index of storage terms in flowbudptr
104  integer(I4B), pointer :: idxbudtmvr => null() !< index of to mover terms in flowbudptr
105  integer(I4B), pointer :: idxbudfmvr => null() !< index of from mover terms in flowbudptr
106  integer(I4B), pointer :: idxbudaux => null() !< index of auxiliary terms in flowbudptr
107  integer(I4B), dimension(:), pointer, contiguous :: idxbudssm => null() !< flag that flowbudptr%buditem is a general solute source/sink
108  integer(I4B), pointer :: nconcbudssm => null() !< number of concbudssm terms (columns)
109  real(dp), dimension(:, :), pointer, contiguous :: concbudssm => null() !< user specified concentrations (or temperatures) for flow terms
110  real(dp), dimension(:), pointer, contiguous :: qmfrommvr => null() !< a mass or energy flow coming from the mover that needs to be added
111  real(dp), pointer :: eqnsclfac => null() !< governing equation scale factor; =1. for solute; =rhow*cpw for energy
112  character(len=LENVARNAME) :: depvartype = '' !< stores string identifying dependent variable type, depending on model type
113  character(len=LENVARNAME) :: depvarunit = '' !< "mass" or "energy"
114  character(len=LENVARNAME) :: depvarunitabbrev = '' !< "M" or "E"
115  !
116  ! -- pointer to flow package boundary
117  type(bndtype), pointer :: flowpackagebnd => null()
118  !
119  ! -- budget objects
120  type(budgetobjecttype), pointer :: budobj => null() !< apt solute budget object
121  type(budgetobjecttype), pointer :: flowbudptr => null() !< GWF flow budget object
122  !
123  ! -- table objects
124  type(tabletype), pointer :: dvtab => null()
125 
126  contains
127 
128  procedure :: set_pointers => apt_set_pointers
129  procedure :: bnd_ac => apt_ac
130  procedure :: bnd_mc => apt_mc
131  procedure :: bnd_ar => apt_ar
132  procedure :: bnd_rp => apt_rp
133  procedure :: bnd_ad => apt_ad
134  procedure :: bnd_reset => apt_reset
135  procedure :: bnd_fc => apt_fc
136  procedure, public :: apt_fc_expanded ! Made public for uze
137  procedure, public :: apt_ad_chk
138  procedure :: pak_fc_expanded
139  procedure, private :: apt_fc_nonexpanded
140  procedure, public :: apt_cfupdate ! Made public for uze
141  procedure :: apt_check_valid
142  procedure :: apt_set_stressperiod
143  procedure :: pak_set_stressperiod
145  procedure :: bnd_cq => apt_cq
146  procedure :: bnd_ot_package_flows => apt_ot_package_flows
147  procedure :: bnd_ot_dv => apt_ot_dv
148  procedure :: bnd_ot_bdsummary => apt_ot_bdsummary
149  procedure :: bnd_da => apt_da
150  procedure :: allocate_scalars
152  procedure :: apt_allocate_arrays
153  procedure :: find_apt_package
154  procedure :: apt_solve
155  procedure :: pak_solve
156  procedure :: bnd_options => apt_options
157  procedure :: read_dimensions => apt_read_dimensions
158  procedure :: apt_read_cvs
159  procedure :: read_initial_attr => apt_read_initial_attr
160  procedure :: define_listlabel
161  ! -- methods for observations
162  procedure :: bnd_obs_supported => apt_obs_supported
163  procedure :: bnd_df_obs => apt_df_obs
164  procedure :: pak_df_obs
165  procedure :: pak_rp_obs
166  procedure :: bnd_rp_obs => apt_rp_obs
167  procedure :: rp_obs_byfeature
168  procedure :: rp_obs_budterm
169  procedure :: rp_obs_flowjaface
170  procedure :: bnd_bd_obs => apt_bd_obs
171  procedure :: pak_bd_obs
172  procedure :: pak_get_nbudterms
173  procedure :: apt_setup_budobj
174  procedure :: pak_setup_budobj
175  procedure :: apt_fill_budobj
176  procedure :: pak_fill_budobj
177  procedure, public :: apt_get_volumes !< made public for sft/sfe
178  procedure, public :: apt_stor_term
179  procedure, public :: apt_tmvr_term
180  procedure, public :: apt_fmvr_term !< made public for uze
181  procedure, public :: apt_fjf_term !< made public for uze
182  procedure, private :: apt_copy2flowp
183  procedure, private :: apt_setup_tableobj
184  procedure, public :: get_mvr_depvar
185 
186  end type tspapttype
187 
188 contains
189 
190  !> @brief Add package connection to matrix
191  !<
192  subroutine apt_ac(this, moffset, sparse)
194  use sparsemodule, only: sparsematrix
195  ! -- dummy
196  class(tspapttype), intent(inout) :: this
197  integer(I4B), intent(in) :: moffset
198  type(sparsematrix), intent(inout) :: sparse
199  ! -- local
200  integer(I4B) :: i, n
201  integer(I4B) :: jj, jglo
202  integer(I4B) :: nglo
203  ! -- format
204  !
205  ! -- Add package rows to sparse
206  if (this%imatrows /= 0) then
207  !
208  ! -- diagonal
209  do n = 1, this%ncv
210  nglo = moffset + this%dis%nodes + this%ioffset + n
211  call sparse%addconnection(nglo, nglo, 1)
212  end do
213  !
214  ! -- apt-gwf connections
215  do i = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
216  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(i)
217  jj = this%flowbudptr%budterm(this%idxbudgwf)%id2(i)
218  nglo = moffset + this%dis%nodes + this%ioffset + n
219  jglo = jj + moffset
220  call sparse%addconnection(nglo, jglo, 1)
221  call sparse%addconnection(jglo, nglo, 1)
222  end do
223  !
224  ! -- apt-apt connections
225  if (this%idxbudfjf /= 0) then
226  do i = 1, this%flowbudptr%budterm(this%idxbudfjf)%maxlist
227  n = this%flowbudptr%budterm(this%idxbudfjf)%id1(i)
228  jj = this%flowbudptr%budterm(this%idxbudfjf)%id2(i)
229  nglo = moffset + this%dis%nodes + this%ioffset + n
230  jglo = moffset + this%dis%nodes + this%ioffset + jj
231  call sparse%addconnection(nglo, jglo, 1)
232  end do
233  end if
234  end if
235  end subroutine apt_ac
236 
237  !> @brief Advanced package transport map package connections to matrix
238  !<
239  subroutine apt_mc(this, moffset, matrix_sln)
240  use sparsemodule, only: sparsematrix
241  ! -- dummy
242  class(tspapttype), intent(inout) :: this
243  integer(I4B), intent(in) :: moffset
244  class(matrixbasetype), pointer :: matrix_sln
245  ! -- local
246  integer(I4B) :: n, j, iglo, jglo
247  integer(I4B) :: ipos
248  ! -- format
249  !
250  ! -- allocate memory for index arrays
251  call this%apt_allocate_index_arrays()
252  !
253  ! -- store index positions
254  if (this%imatrows /= 0) then
255  !
256  ! -- Find the position of each connection in the global ia, ja structure
257  ! and store them in idxglo. idxglo allows this model to insert or
258  ! retrieve values into or from the global A matrix
259  ! -- apt rows
260  do n = 1, this%ncv
261  this%idxlocnode(n) = this%dis%nodes + this%ioffset + n
262  iglo = moffset + this%dis%nodes + this%ioffset + n
263  this%idxpakdiag(n) = matrix_sln%get_position_diag(iglo)
264  end do
265  do ipos = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
266  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(ipos)
267  j = this%flowbudptr%budterm(this%idxbudgwf)%id2(ipos)
268  iglo = moffset + this%dis%nodes + this%ioffset + n
269  jglo = j + moffset
270  this%idxdglo(ipos) = matrix_sln%get_position_diag(iglo)
271  this%idxoffdglo(ipos) = matrix_sln%get_position(iglo, jglo)
272  end do
273  !
274  ! -- apt contributions to gwf portion of global matrix
275  do ipos = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
276  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(ipos)
277  j = this%flowbudptr%budterm(this%idxbudgwf)%id2(ipos)
278  iglo = j + moffset
279  jglo = moffset + this%dis%nodes + this%ioffset + n
280  this%idxsymdglo(ipos) = matrix_sln%get_position_diag(iglo)
281  this%idxsymoffdglo(ipos) = matrix_sln%get_position(iglo, jglo)
282  end do
283  !
284  ! -- apt-apt contributions to gwf portion of global matrix
285  if (this%idxbudfjf /= 0) then
286  do ipos = 1, this%flowbudptr%budterm(this%idxbudfjf)%nlist
287  n = this%flowbudptr%budterm(this%idxbudfjf)%id1(ipos)
288  j = this%flowbudptr%budterm(this%idxbudfjf)%id2(ipos)
289  iglo = moffset + this%dis%nodes + this%ioffset + n
290  jglo = moffset + this%dis%nodes + this%ioffset + j
291  this%idxfjfdglo(ipos) = matrix_sln%get_position_diag(iglo)
292  this%idxfjfoffdglo(ipos) = matrix_sln%get_position(iglo, jglo)
293  end do
294  end if
295  end if
296  end subroutine apt_mc
297 
298  !> @brief Advanced package transport allocate and read (ar) routine
299  !<
300  subroutine apt_ar(this)
301  ! -- modules
302  ! -- dummy
303  class(tspapttype), intent(inout) :: this
304  ! -- local
305  integer(I4B) :: j
306  logical :: found
307  ! -- formats
308  character(len=*), parameter :: fmtapt = &
309  "(1x,/1x,'APT -- ADVANCED PACKAGE TRANSPORT, VERSION 1, 3/5/2020', &
310  &' INPUT READ FROM UNIT ', i0, //)"
311  !
312  ! -- Get obs setup
313  call this%obs%obs_ar()
314  !
315  ! --print a message identifying the apt package.
316  write (this%iout, fmtapt) this%inunit
317  !
318  ! -- Allocate arrays
319  call this%apt_allocate_arrays()
320  !
321  ! -- read optional initial package parameters
322  call this%read_initial_attr()
323  !
324  ! -- Find the package index in the GWF model or GWF budget file
325  ! for the corresponding apt flow package
326  call this%fmi%get_package_index(this%flowpackagename, this%igwfaptpak)
327  !
328  ! -- Tell fmi that this package is being handled by APT, otherwise
329  ! SSM would handle the flows into GWT from this pack. Then point the
330  ! fmi data for an advanced package to xnewpak and qmfrommvr
331  this%fmi%iatp(this%igwfaptpak) = 1
332  this%fmi%datp(this%igwfaptpak)%concpack => this%get_mvr_depvar()
333  this%fmi%datp(this%igwfaptpak)%qmfrommvr => this%qmfrommvr
334  !
335  ! -- If there is an associated flow package and the user wishes to put
336  ! simulated concentrations (or temperatures) into a aux variable
337  ! column, then find the column number.
338  if (associated(this%flowpackagebnd)) then
339  if (this%cauxfpconc /= '') then
340  found = .false.
341  do j = 1, this%flowpackagebnd%naux
342  if (this%flowpackagebnd%auxname(j) == this%cauxfpconc) then
343  this%iauxfpconc = j
344  found = .true.
345  exit
346  end if
347  end do
348  if (this%iauxfpconc == 0) then
349  errmsg = 'Could not find auxiliary variable '// &
350  trim(adjustl(this%cauxfpconc))//' in flow package '// &
351  trim(adjustl(this%flowpackagename))
352  call store_error(errmsg)
353  call this%parser%StoreErrorUnit()
354  else
355  ! -- tell package not to update this auxiliary variable
356  this%flowpackagebnd%noupdateauxvar(this%iauxfpconc) = 1
357  call this%apt_copy2flowp()
358  end if
359  end if
360  end if
361  end subroutine apt_ar
362 
363  !> @brief Advanced package transport read and prepare (rp) routine
364  !!
365  !! This subroutine calls the attached packages' read and prepare routines.
366  !<
367  subroutine apt_rp(this)
368  use tdismodule, only: kper, nper
369  ! -- dummy
370  class(tspapttype), intent(inout) :: this
371  ! -- local
372  integer(I4B) :: ierr
373  integer(I4B) :: n
374  logical :: isfound, endOfBlock
375  character(len=LINELENGTH) :: title
376  character(len=LINELENGTH) :: line
377  integer(I4B) :: itemno
378  integer(I4B) :: igwfnode
379  ! -- formats
380  character(len=*), parameter :: fmtblkerr = &
381  &"('Error. Looking for BEGIN PERIOD iper. Found ', a, ' instead.')"
382  character(len=*), parameter :: fmtlsp = &
383  &"(1X,/1X,'REUSING ',A,'S FROM LAST STRESS PERIOD')"
384  !
385  ! -- set nbound to maxbound
386  this%nbound = this%maxbound
387  !
388  ! -- Set ionper to the stress period number for which a new block of data
389  ! will be read.
390  if (this%inunit == 0) return
391  !
392  ! -- get stress period data
393  if (this%ionper < kper) then
394  !
395  ! -- get period block
396  call this%parser%GetBlock('PERIOD', isfound, ierr, &
397  supportopenclose=.true., &
398  blockrequired=.false.)
399  if (isfound) then
400  !
401  ! -- read ionper and check for increasing period numbers
402  call this%read_check_ionper()
403  else
404  !
405  ! -- PERIOD block not found
406  if (ierr < 0) then
407  ! -- End of file found; data applies for remainder of simulation.
408  this%ionper = nper + 1
409  else
410  ! -- Found invalid block
411  call this%parser%GetCurrentLine(line)
412  write (errmsg, fmtblkerr) adjustl(trim(line))
413  call store_error(errmsg)
414  call this%parser%StoreErrorUnit()
415  end if
416  end if
417  end if
418  !
419  ! -- Read data if ionper == kper
420  if (this%ionper == kper) then
421  !
422  ! -- setup table for period data
423  if (this%iprpak /= 0) then
424  !
425  ! -- reset the input table object
426  title = trim(adjustl(this%text))//' PACKAGE ('// &
427  trim(adjustl(this%packName))//') DATA FOR PERIOD'
428  write (title, '(a,1x,i6)') trim(adjustl(title)), kper
429  call table_cr(this%inputtab, this%packName, title)
430  call this%inputtab%table_df(1, 4, this%iout, finalize=.false.)
431  text = 'NUMBER'
432  call this%inputtab%initialize_column(text, 10, alignment=tabcenter)
433  text = 'KEYWORD'
434  call this%inputtab%initialize_column(text, 20, alignment=tableft)
435  do n = 1, 2
436  write (text, '(a,1x,i6)') 'VALUE', n
437  call this%inputtab%initialize_column(text, 15, alignment=tabcenter)
438  end do
439  end if
440  !
441  ! -- read data
442  stressperiod: do
443  call this%parser%GetNextLine(endofblock)
444  if (endofblock) exit
445  !
446  ! -- get feature number
447  itemno = this%parser%GetInteger()
448  !
449  ! -- read data from the rest of the line
450  call this%apt_set_stressperiod(itemno)
451  !
452  ! -- write line to table
453  if (this%iprpak /= 0) then
454  call this%parser%GetCurrentLine(line)
455  call this%inputtab%line_to_columns(line)
456  end if
457  end do stressperiod
458 
459  if (this%iprpak /= 0) then
460  call this%inputtab%finalize_table()
461  end if
462  !
463  ! -- using stress period data from the previous stress period
464  else
465  write (this%iout, fmtlsp) trim(this%filtyp)
466  end if
467  !
468  ! -- write summary of stress period error messages
469  ierr = count_errors()
470  if (ierr > 0) then
471  call this%parser%StoreErrorUnit()
472  end if
473  !
474  ! -- fill arrays
475  do n = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
476  igwfnode = this%flowbudptr%budterm(this%idxbudgwf)%id2(n)
477  this%nodelist(n) = igwfnode
478  end do
479  end subroutine apt_rp
480 
481  subroutine apt_ad_chk(this)
482  ! -- dummy
483  class(tspapttype), intent(inout) :: this
484  ! function available for override by packages
485  end subroutine apt_ad_chk
486 
487  !> @brief Advanced package transport set stress period routine.
488  !!
489  !! Set a stress period attribute for an advanced transport package feature
490  !! (itemno) using keywords.
491  !<
492  subroutine apt_set_stressperiod(this, itemno)
493  ! -- module
495  ! -- dummy
496  class(tspapttype), intent(inout) :: this
497  integer(I4B), intent(in) :: itemno
498  ! -- local
499  character(len=LINELENGTH) :: text
500  character(len=LINELENGTH) :: caux
501  character(len=LINELENGTH) :: keyword
502  integer(I4B) :: ierr
503  integer(I4B) :: ii
504  integer(I4B) :: jj
505  real(DP), pointer :: bndElem => null()
506  logical :: found
507  ! -- formats
508  !
509  ! -- Support these general options in LKT, SFT, MWT, UZT
510  ! STATUS <status>
511  ! CONCENTRATION <concentration> or TEMPERATURE <temperature>
512  ! WITHDRAWAL <withdrawal>
513  ! AUXILIARY <auxname> <auxval>
514  !
515  ! -- read line
516  call this%parser%GetStringCaps(keyword)
517  select case (keyword)
518  case ('STATUS')
519  ierr = this%apt_check_valid(itemno)
520  if (ierr /= 0) then
521  goto 999
522  end if
523  call this%parser%GetStringCaps(text)
524  this%status(itemno) = text(1:8)
525  if (text == 'CONSTANT') then
526  this%iboundpak(itemno) = -1
527  else if (text == 'INACTIVE') then
528  this%iboundpak(itemno) = 0
529  else if (text == 'ACTIVE') then
530  this%iboundpak(itemno) = 1
531  else
532  write (errmsg, '(a,a)') &
533  'Unknown '//trim(this%text)//' status keyword: ', text//'.'
534  call store_error(errmsg)
535  end if
536  case ('CONCENTRATION', 'TEMPERATURE')
537  ierr = this%apt_check_valid(itemno)
538  if (ierr /= 0) then
539  goto 999
540  end if
541  call this%parser%GetString(text)
542  jj = 1 ! For feature concentration
543  bndelem => this%concfeat(itemno)
544  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
545  this%packName, 'BND', this%tsManager, &
546  this%iprpak, this%depvartype)
547  case ('AUXILIARY')
548  ierr = this%apt_check_valid(itemno)
549  if (ierr /= 0) then
550  goto 999
551  end if
552  call this%parser%GetStringCaps(caux)
553  do jj = 1, this%naux
554  if (trim(adjustl(caux)) /= trim(adjustl(this%auxname(jj)))) cycle
555  call this%parser%GetString(text)
556  ii = itemno
557  bndelem => this%lauxvar(jj, ii)
558  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
559  this%packName, 'AUX', &
560  this%tsManager, this%iprpak, &
561  this%auxname(jj))
562  exit
563  end do
564  case default
565  !
566  ! -- call the specific package to look for stress period data
567  call this%pak_set_stressperiod(itemno, keyword, found)
568  !
569  ! -- terminate with error if data not valid
570  if (.not. found) then
571  write (errmsg, '(2a)') &
572  'Unknown '//trim(adjustl(this%text))//' data keyword: ', &
573  trim(keyword)//'.'
574  call store_error(errmsg)
575  end if
576  end select
577  !
578  ! -- terminate if any errors were detected
579 999 if (count_errors() > 0) then
580  call this%parser%StoreErrorUnit()
581  end if
582  end subroutine apt_set_stressperiod
583 
584  !> @brief Advanced package transport set stress period routine.
585  !!
586  !! Set a stress period attribute for an individual package. This routine
587  !! must be overridden.
588  !<
589  subroutine pak_set_stressperiod(this, itemno, keyword, found)
590  ! -- dummy
591  class(tspapttype), intent(inout) :: this
592  integer(I4B), intent(in) :: itemno
593  character(len=*), intent(in) :: keyword
594  logical, intent(inout) :: found
595  ! -- local
596 
597  ! -- formats
598  !
599  ! -- this routine should never be called
600  found = .false.
601  call store_error('Program error: pak_set_stressperiod not implemented.', &
602  terminate=.true.)
603  end subroutine pak_set_stressperiod
604 
605  !> @brief Advanced package transport routine
606  !!
607  !! Determine if a valid feature number has been specified.
608  !<
609  function apt_check_valid(this, itemno) result(ierr)
610  ! -- return
611  integer(I4B) :: ierr
612  ! -- dummy
613  class(tspapttype), intent(inout) :: this
614  integer(I4B), intent(in) :: itemno
615  ! -- formats
616  ierr = 0
617  if (itemno < 1 .or. itemno > this%ncv) then
618  write (errmsg, '(a,1x,i6,1x,a,1x,i6)') &
619  'Featureno ', itemno, 'must be > 0 and <= ', this%ncv
620  call store_error(errmsg)
621  ierr = 1
622  end if
623  end function apt_check_valid
624 
625  !> @brief Advanced package transport utility function
626  !!
627  !! Set the concentration (or temperature) to be used by either MVT or MVE
628  !<
629  function get_mvr_depvar(this)
630  ! -- dummy
631  class(tspapttype) :: this
632  ! -- return
633  real(dp), dimension(:), contiguous, pointer :: get_mvr_depvar
634  !
635  get_mvr_depvar => this%xnewpak
636  end function get_mvr_depvar
637 
638  !> @brief Advanced package transport routine
639  !!
640  !! Add package connections to matrix
641  !<
642  subroutine apt_ad(this)
643  ! -- modules
645  ! -- dummy
646  class(tspapttype) :: this
647  ! -- local
648  integer(I4B) :: n
649  integer(I4B) :: j, iaux
650  !
651  ! -- Advance the time series
652  call this%TsManager%ad()
653  !
654  ! -- update auxiliary variables by copying from the derived-type time
655  ! series variable into the bndpackage auxvar variable so that this
656  ! information is properly written to the GWF budget file
657  if (this%naux > 0) then
658  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
659  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
660  do iaux = 1, this%naux
661  this%auxvar(iaux, j) = this%lauxvar(iaux, n)
662  end do
663  end do
664  end if
665  !
666  ! -- copy xnew into xold and set xnewpak to specified concentration (or
667  ! temperature) for constant concentration/temperature features
668  if (ifailedstepretry == 0) then
669  do n = 1, this%ncv
670  this%xoldpak(n) = this%xnewpak(n)
671  if (this%iboundpak(n) < 0) then
672  this%xnewpak(n) = this%concfeat(n)
673  end if
674  end do
675  else
676  do n = 1, this%ncv
677  this%xnewpak(n) = this%xoldpak(n)
678  if (this%iboundpak(n) < 0) then
679  this%xnewpak(n) = this%concfeat(n)
680  end if
681  end do
682  end if
683  !
684  ! -- pakmvrobj ad
685  !if (this%imover == 1) then
686  ! call this%pakmvrobj%ad()
687  !end if
688  !
689  ! -- For each observation, push simulated value and corresponding
690  ! simulation time from "current" to "preceding" and reset
691  ! "current" value.
692  call this%obs%obs_ad()
693  !
694  ! -- run package-specific checks
695  call this%apt_ad_chk()
696  end subroutine apt_ad
697 
698  !> @brief Override bnd reset for custom mover logic
699  subroutine apt_reset(this)
700  class(tspapttype) :: this !< GwtAptType object
701  ! local
702  integer(I4B) :: i
703  !
704  do i = 1, size(this%qmfrommvr)
705  this%qmfrommvr(i) = dzero
706  end do
707  end subroutine apt_reset
708 
709  subroutine apt_fc(this, rhs, ia, idxglo, matrix_sln)
710  ! -- modules
711  ! -- dummy
712  class(tspapttype) :: this
713  real(DP), dimension(:), intent(inout) :: rhs
714  integer(I4B), dimension(:), intent(in) :: ia
715  integer(I4B), dimension(:), intent(in) :: idxglo
716  class(matrixbasetype), pointer :: matrix_sln
717  ! -- local
718  !
719  ! -- Call fc depending on whether or not a matrix is expanded or not
720  if (this%imatrows == 0) then
721  call this%apt_fc_nonexpanded(rhs, ia, idxglo, matrix_sln)
722  else
723  call this%apt_fc_expanded(rhs, ia, idxglo, matrix_sln)
724  end if
725  end subroutine apt_fc
726 
727  !> @brief Advanced package transport fill coefficient (fc) method
728  !!
729  !! Routine to formulate the nonexpanded matrix case in which feature
730  !! concentrations (or temperatures) are solved explicitly
731  !<
732  subroutine apt_fc_nonexpanded(this, rhs, ia, idxglo, matrix_sln)
733  ! -- modules
734  ! -- dummy
735  class(tspapttype) :: this
736  real(DP), dimension(:), intent(inout) :: rhs
737  integer(I4B), dimension(:), intent(in) :: ia
738  integer(I4B), dimension(:), intent(in) :: idxglo
739  class(matrixbasetype), pointer :: matrix_sln
740  ! -- local
741  integer(I4B) :: j, igwfnode, idiag
742  !
743  ! -- solve for concentration (or temperatures) in the features
744  call this%apt_solve()
745  !
746  ! -- add hcof and rhs terms (from apt_solve) to the gwf matrix
747  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
748  igwfnode = this%flowbudptr%budterm(this%idxbudgwf)%id2(j)
749  if (this%ibound(igwfnode) < 1) cycle
750  idiag = idxglo(ia(igwfnode))
751  call matrix_sln%add_value_pos(idiag, this%hcof(j))
752  rhs(igwfnode) = rhs(igwfnode) + this%rhs(j)
753  end do
754  end subroutine apt_fc_nonexpanded
755 
756  !> @brief Advanced package transport fill coefficient (fc) method
757  !!
758  !! Routine to formulate the expanded matrix case in which new rows are added
759  !! to the system of equations for each advanced package transport feature
760  !<
761  subroutine apt_fc_expanded(this, rhs, ia, idxglo, matrix_sln)
762  ! -- modules
763  ! -- dummy
764  class(tspapttype) :: this
765  real(DP), dimension(:), intent(inout) :: rhs
766  integer(I4B), dimension(:), intent(in) :: ia
767  integer(I4B), dimension(:), intent(in) :: idxglo
768  class(matrixbasetype), pointer :: matrix_sln
769  ! -- local
770  integer(I4B) :: j, n, n1, n2
771  integer(I4B) :: iloc
772  integer(I4B) :: iposd, iposoffd
773  integer(I4B) :: ipossymd, ipossymoffd
774  real(DP) :: cold
775  real(DP) :: qbnd, qbnd_scaled
776  real(DP) :: omega
777  real(DP) :: rrate
778  real(DP) :: rhsval
779  real(DP) :: hcofval
780  !
781  ! -- call the specific method for the advanced transport package, such as
782  ! what would be overridden by
783  ! GwtLktType, GwtSftType, GwtMwtType, GwtUztType
784  ! This routine will add terms for rainfall, runoff, or other terms
785  ! specific to the package
786  call this%pak_fc_expanded(rhs, ia, idxglo, matrix_sln)
787  !
788  ! -- mass (or energy) storage in features
789  do n = 1, this%ncv
790  cold = this%xoldpak(n)
791  iloc = this%idxlocnode(n)
792  iposd = this%idxpakdiag(n)
793  call this%apt_stor_term(n, n1, n2, rrate, rhsval, hcofval)
794  call matrix_sln%add_value_pos(iposd, hcofval)
795  rhs(iloc) = rhs(iloc) + rhsval
796  end do
797  !
798  ! -- add to mover contribution
799  if (this%idxbudtmvr /= 0) then
800  do j = 1, this%flowbudptr%budterm(this%idxbudtmvr)%nlist
801  call this%apt_tmvr_term(j, n1, n2, rrate, rhsval, hcofval)
802  iloc = this%idxlocnode(n1)
803  iposd = this%idxpakdiag(n1)
804  call matrix_sln%add_value_pos(iposd, hcofval)
805  rhs(iloc) = rhs(iloc) + rhsval
806  end do
807  end if
808  !
809  ! -- add from mover contribution
810  if (this%idxbudfmvr /= 0) then
811  do n = 1, this%ncv
812  rhsval = this%qmfrommvr(n) ! this will already be in terms of energy for heat transport
813  iloc = this%idxlocnode(n)
814  rhs(iloc) = rhs(iloc) - rhsval
815  end do
816  end if
817  !
818  ! -- go through each apt-gwf connection
819  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
820  !
821  ! -- set n to feature number and process if active feature
822  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
823  if (this%iboundpak(n) /= 0) then
824  !
825  ! -- set acoef and rhs to negative so they are relative to apt and not gwt
826  qbnd = this%flowbudptr%budterm(this%idxbudgwf)%flow(j)
827  omega = dzero
828  if (qbnd < dzero) omega = done
829  qbnd_scaled = qbnd * this%eqnsclfac
830  !
831  ! -- add to apt row
832  iposd = this%idxdglo(j)
833  iposoffd = this%idxoffdglo(j)
834  call matrix_sln%add_value_pos(iposd, omega * qbnd_scaled)
835  call matrix_sln%add_value_pos(iposoffd, (done - omega) * qbnd_scaled)
836  !
837  ! -- add to gwf row for apt connection
838  ipossymd = this%idxsymdglo(j)
839  ipossymoffd = this%idxsymoffdglo(j)
840  call matrix_sln%add_value_pos(ipossymd, -(done - omega) * qbnd_scaled)
841  call matrix_sln%add_value_pos(ipossymoffd, -omega * qbnd_scaled)
842  end if
843  end do
844  !
845  ! -- go through each apt-apt connection
846  if (this%idxbudfjf /= 0) then
847  do j = 1, this%flowbudptr%budterm(this%idxbudfjf)%nlist
848  n1 = this%flowbudptr%budterm(this%idxbudfjf)%id1(j)
849  n2 = this%flowbudptr%budterm(this%idxbudfjf)%id2(j)
850  qbnd = this%flowbudptr%budterm(this%idxbudfjf)%flow(j)
851  if (qbnd <= dzero) then
852  omega = done
853  else
854  omega = dzero
855  end if
856  qbnd_scaled = qbnd * this%eqnsclfac
857  iposd = this%idxfjfdglo(j)
858  iposoffd = this%idxfjfoffdglo(j)
859  call matrix_sln%add_value_pos(iposd, omega * qbnd_scaled)
860  call matrix_sln%add_value_pos(iposoffd, (done - omega) * qbnd_scaled)
861  end do
862  end if
863  end subroutine apt_fc_expanded
864 
865  !> @brief Advanced package transport fill coefficient (fc) method
866  !!
867  !! Routine to allow a subclass advanced transport package to inject
868  !! terms into the matrix assembly. This method must be overridden.
869  !<
870  subroutine pak_fc_expanded(this, rhs, ia, idxglo, matrix_sln)
871  ! -- modules
872  ! -- dummy
873  class(tspapttype) :: this
874  real(DP), dimension(:), intent(inout) :: rhs
875  integer(I4B), dimension(:), intent(in) :: ia
876  integer(I4B), dimension(:), intent(in) :: idxglo
877  class(matrixbasetype), pointer :: matrix_sln
878  ! -- local
879  !
880  ! -- this routine should never be called
881  call store_error('Program error: pak_fc_expanded not implemented.', &
882  terminate=.true.)
883  end subroutine pak_fc_expanded
884 
885  !> @brief Advanced package transport routine
886  !!
887  !! Calculate advanced package transport hcof and rhs so transport budget is
888  !! calculated.
889  !<
890  subroutine apt_cfupdate(this)
891  ! -- modules
892  ! -- dummy
893  class(tspapttype) :: this
894  ! -- local
895  integer(I4B) :: j, n
896  real(DP) :: qbnd
897  real(DP) :: omega
898  !
899  ! -- Calculate hcof and rhs terms so GWF exchanges are calculated correctly
900  ! -- go through each apt-gwf connection and calculate
901  ! rhs and hcof terms for gwt/gwe matrix rows
902  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
903  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
904  this%hcof(j) = dzero
905  this%rhs(j) = dzero
906  if (this%iboundpak(n) /= 0) then
907  qbnd = this%flowbudptr%budterm(this%idxbudgwf)%flow(j)
908  omega = dzero
909  if (qbnd < dzero) omega = done
910  this%hcof(j) = -(done - omega) * qbnd * this%eqnsclfac
911  this%rhs(j) = omega * qbnd * this%xnewpak(n) * this%eqnsclfac
912  end if
913  end do
914  end subroutine apt_cfupdate
915 
916  !> @brief Advanced package transport calculate flows (cq) routine
917  !!
918  !! Calculate flows for the advanced package transport feature
919  !<
920  subroutine apt_cq(this, x, flowja, iadv)
921  ! -- modules
922  ! -- dummy
923  class(tspapttype), intent(inout) :: this
924  real(DP), dimension(:), intent(in) :: x
925  real(DP), dimension(:), contiguous, intent(inout) :: flowja
926  integer(I4B), optional, intent(in) :: iadv
927  ! -- local
928  integer(I4B) :: n, n1, n2
929  real(DP) :: rrate
930  !
931  ! -- Solve the feature concentrations (or temperatures) again or update
932  ! the feature hcof and rhs terms
933  if (this%imatrows == 0) then
934  call this%apt_solve()
935  else
936  call this%apt_cfupdate()
937  end if
938  !
939  ! -- call base functionality in bnd_cq
940  call this%BndType%bnd_cq(x, flowja)
941  !
942  ! -- calculate storage term
943  do n = 1, this%ncv
944  rrate = dzero
945  if (this%iboundpak(n) > 0) then
946  call this%apt_stor_term(n, n1, n2, rrate)
947  end if
948  this%qsto(n) = rrate
949  end do
950  !
951  ! -- Copy concentrations (or temperatures) into the flow package auxiliary variable
952  call this%apt_copy2flowp()
953  !
954  ! -- fill the budget object
955  call this%apt_fill_budobj(x, flowja)
956  end subroutine apt_cq
957 
958  !> @brief Save advanced package flows routine
959  !<
960  subroutine apt_ot_package_flows(this, icbcfl, ibudfl)
961  use tdismodule, only: kstp, kper, delt, pertim, totim
962  class(tspapttype) :: this
963  integer(I4B), intent(in) :: icbcfl
964  integer(I4B), intent(in) :: ibudfl
965  integer(I4B) :: ibinun
966  !
967  ! -- write the flows from the budobj
968  ibinun = 0
969  if (this%ibudgetout /= 0) then
970  ibinun = this%ibudgetout
971  end if
972  if (icbcfl == 0) ibinun = 0
973  if (ibinun > 0) then
974  call this%budobj%save_flows(this%dis, ibinun, kstp, kper, delt, &
975  pertim, totim, this%iout)
976  end if
977  !
978  ! -- Print lake flows table
979  if (ibudfl /= 0 .and. this%iprflow /= 0) then
980  call this%budobj%write_flowtable(this%dis, kstp, kper)
981  end if
982  end subroutine apt_ot_package_flows
983 
984  subroutine apt_ot_dv(this, idvsave, idvprint)
985  ! -- modules
986  use constantsmodule, only: lenbudtxt
987  use tdismodule, only: kstp, kper, pertim, totim
989  use inputoutputmodule, only: ulasav
990  ! -- dummy
991  class(tspapttype) :: this
992  integer(I4B), intent(in) :: idvsave
993  integer(I4B), intent(in) :: idvprint
994  ! -- local
995  integer(I4B) :: ibinun
996  integer(I4B) :: n
997  real(DP) :: c
998  character(len=LENBUDTXT) :: text
999  !
1000  ! -- set unit number for binary dependent variable output
1001  ibinun = 0
1002  if (this%iconcout /= 0) then
1003  ibinun = this%iconcout
1004  end if
1005  if (idvsave == 0) ibinun = 0
1006  !
1007  ! -- write binary output
1008  if (ibinun > 0) then
1009  do n = 1, this%ncv
1010  c = this%xnewpak(n)
1011  if (this%iboundpak(n) == 0) then
1012  c = dhnoflo
1013  end if
1014  this%dbuff(n) = c
1015  end do
1016  write (text, '(a)') str_pad_left(this%depvartype, lenvarname)
1017  call ulasav(this%dbuff, text, kstp, kper, pertim, totim, &
1018  this%ncv, 1, 1, ibinun)
1019  end if
1020  !
1021  ! -- write apt conc table
1022  if (idvprint /= 0 .and. this%iprconc /= 0) then
1023  !
1024  ! -- set table kstp and kper
1025  call this%dvtab%set_kstpkper(kstp, kper)
1026  !
1027  ! -- fill concentration data
1028  do n = 1, this%ncv
1029  if (this%inamedbound == 1) then
1030  call this%dvtab%add_term(this%featname(n))
1031  end if
1032  call this%dvtab%add_term(n)
1033  call this%dvtab%add_term(this%xnewpak(n))
1034  end do
1035  end if
1036  end subroutine apt_ot_dv
1037 
1038  !> @brief Print advanced package transport dependent variables
1039  !<
1040  subroutine apt_ot_bdsummary(this, kstp, kper, iout, ibudfl)
1041  ! -- module
1042  use tdismodule, only: totim, delt
1043  ! -- dummy
1044  class(tspapttype) :: this !< TspAptType object
1045  integer(I4B), intent(in) :: kstp !< time step number
1046  integer(I4B), intent(in) :: kper !< period number
1047  integer(I4B), intent(in) :: iout !< flag and unit number for the model listing file
1048  integer(I4B), intent(in) :: ibudfl !< flag indicating budget should be written
1049  !
1050  call this%budobj%write_budtable(kstp, kper, iout, ibudfl, totim, delt)
1051  end subroutine apt_ot_bdsummary
1052 
1053  !> @ brief Allocate scalars
1054  !!
1055  !! Allocate scalar variables for an advanced package
1056  !<
1057  subroutine allocate_scalars(this)
1058  ! -- modules
1060  ! -- dummy
1061  class(tspapttype) :: this
1062  ! -- local
1063  !
1064  ! -- allocate scalars in NumericalPackageType
1065  call this%BndType%allocate_scalars()
1066  !
1067  ! -- Allocate
1068  call mem_allocate(this%iauxfpconc, 'IAUXFPCONC', this%memoryPath)
1069  call mem_allocate(this%imatrows, 'IMATROWS', this%memoryPath)
1070  call mem_allocate(this%iprconc, 'IPRCONC', this%memoryPath)
1071  call mem_allocate(this%iconcout, 'ICONCOUT', this%memoryPath)
1072  call mem_allocate(this%ibudgetout, 'IBUDGETOUT', this%memoryPath)
1073  call mem_allocate(this%ibudcsv, 'IBUDCSV', this%memoryPath)
1074  call mem_allocate(this%igwfaptpak, 'IGWFAPTPAK', this%memoryPath)
1075  call mem_allocate(this%ncv, 'NCV', this%memoryPath)
1076  call mem_allocate(this%idxbudfjf, 'IDXBUDFJF', this%memoryPath)
1077  call mem_allocate(this%idxbudgwf, 'IDXBUDGWF', this%memoryPath)
1078  call mem_allocate(this%idxbudsto, 'IDXBUDSTO', this%memoryPath)
1079  call mem_allocate(this%idxbudtmvr, 'IDXBUDTMVR', this%memoryPath)
1080  call mem_allocate(this%idxbudfmvr, 'IDXBUDFMVR', this%memoryPath)
1081  call mem_allocate(this%idxbudaux, 'IDXBUDAUX', this%memoryPath)
1082  call mem_allocate(this%nconcbudssm, 'NCONCBUDSSM', this%memoryPath)
1083  call mem_allocate(this%idxprepak, 'IDXPREPAK', this%memoryPath)
1084  call mem_allocate(this%idxlastpak, 'IDXLASTPAK', this%memoryPath)
1085  !
1086  ! -- Initialize
1087  this%iauxfpconc = 0
1088  this%imatrows = 1
1089  this%iprconc = 0
1090  this%iconcout = 0
1091  this%ibudgetout = 0
1092  this%ibudcsv = 0
1093  this%igwfaptpak = 0
1094  this%ncv = 0
1095  this%idxbudfjf = 0
1096  this%idxbudgwf = 0
1097  this%idxbudsto = 0
1098  this%idxbudtmvr = 0
1099  this%idxbudfmvr = 0
1100  this%idxbudaux = 0
1101  this%nconcbudssm = 0
1102  this%idxprepak = 0
1103  this%idxlastpak = 0
1104  !
1105  ! -- set this package as causing asymmetric matrix terms
1106  this%iasym = 1
1107  end subroutine allocate_scalars
1108 
1109  !> @ brief Allocate index arrays
1110  !!
1111  !! Allocate arrays that map to locations in the numerical solution
1112  !<
1113  subroutine apt_allocate_index_arrays(this)
1114  ! -- modules
1116  ! -- dummy
1117  class(tspapttype), intent(inout) :: this
1118  ! -- local
1119  integer(I4B) :: n
1120  !
1121  if (this%imatrows /= 0) then
1122  !
1123  ! -- count number of flow-ja-face connections
1124  n = 0
1125  if (this%idxbudfjf /= 0) then
1126  n = this%flowbudptr%budterm(this%idxbudfjf)%maxlist
1127  end if
1128  !
1129  ! -- allocate pointers to global matrix
1130  call mem_allocate(this%idxlocnode, this%ncv, 'IDXLOCNODE', &
1131  this%memoryPath)
1132  call mem_allocate(this%idxpakdiag, this%ncv, 'IDXPAKDIAG', &
1133  this%memoryPath)
1134  call mem_allocate(this%idxdglo, this%maxbound, 'IDXGLO', &
1135  this%memoryPath)
1136  call mem_allocate(this%idxoffdglo, this%maxbound, 'IDXOFFDGLO', &
1137  this%memoryPath)
1138  call mem_allocate(this%idxsymdglo, this%maxbound, 'IDXSYMDGLO', &
1139  this%memoryPath)
1140  call mem_allocate(this%idxsymoffdglo, this%maxbound, 'IDXSYMOFFDGLO', &
1141  this%memoryPath)
1142  call mem_allocate(this%idxfjfdglo, n, 'IDXFJFDGLO', &
1143  this%memoryPath)
1144  call mem_allocate(this%idxfjfoffdglo, n, 'IDXFJFOFFDGLO', &
1145  this%memoryPath)
1146  else
1147  call mem_allocate(this%idxlocnode, 0, 'IDXLOCNODE', &
1148  this%memoryPath)
1149  call mem_allocate(this%idxpakdiag, 0, 'IDXPAKDIAG', &
1150  this%memoryPath)
1151  call mem_allocate(this%idxdglo, 0, 'IDXGLO', &
1152  this%memoryPath)
1153  call mem_allocate(this%idxoffdglo, 0, 'IDXOFFDGLO', &
1154  this%memoryPath)
1155  call mem_allocate(this%idxsymdglo, 0, 'IDXSYMDGLO', &
1156  this%memoryPath)
1157  call mem_allocate(this%idxsymoffdglo, 0, 'IDXSYMOFFDGLO', &
1158  this%memoryPath)
1159  call mem_allocate(this%idxfjfdglo, 0, 'IDXFJFDGLO', &
1160  this%memoryPath)
1161  call mem_allocate(this%idxfjfoffdglo, 0, 'IDXFJFOFFDGLO', &
1162  this%memoryPath)
1163  end if
1164  end subroutine apt_allocate_index_arrays
1165 
1166  !> @ brief Allocate arrays
1167  !!
1168  !! Allocate advanced package transport arrays
1169  !<
1170  subroutine apt_allocate_arrays(this)
1171  ! -- modules
1173  ! -- dummy
1174  class(tspapttype), intent(inout) :: this
1175  ! -- local
1176  integer(I4B) :: n
1177  !
1178  ! -- call standard BndType allocate scalars
1179  call this%BndType%allocate_arrays()
1180  !
1181  ! -- Allocate
1182  !
1183  ! -- allocate and initialize dbuff
1184  if (this%iconcout > 0) then
1185  call mem_allocate(this%dbuff, this%ncv, 'DBUFF', this%memoryPath)
1186  do n = 1, this%ncv
1187  this%dbuff(n) = dzero
1188  end do
1189  else
1190  call mem_allocate(this%dbuff, 0, 'DBUFF', this%memoryPath)
1191  end if
1192  !
1193  ! -- allocate character array for status
1194  allocate (this%status(this%ncv))
1195  !
1196  ! -- time series
1197  call mem_allocate(this%concfeat, this%ncv, 'CONCFEAT', this%memoryPath)
1198  !
1199  ! -- budget terms
1200  call mem_allocate(this%qsto, this%ncv, 'QSTO', this%memoryPath)
1201  call mem_allocate(this%ccterm, this%ncv, 'CCTERM', this%memoryPath)
1202  !
1203  ! -- concentration for budget terms
1204  call mem_allocate(this%concbudssm, this%nconcbudssm, this%ncv, &
1205  'CONCBUDSSM', this%memoryPath)
1206  !
1207  ! -- mass (or energy) added from the mover transport package
1208  call mem_allocate(this%qmfrommvr, this%ncv, 'QMFROMMVR', this%memoryPath)
1209  !
1210  ! -- initialize arrays
1211  do n = 1, this%ncv
1212  this%status(n) = 'ACTIVE'
1213  this%qsto(n) = dzero
1214  this%ccterm(n) = dzero
1215  this%qmfrommvr(n) = dzero
1216  this%concbudssm(:, n) = dzero
1217  this%concfeat(n) = dzero
1218  end do
1219  end subroutine apt_allocate_arrays
1220 
1221  !> @ brief Deallocate memory
1222  !!
1223  !! Deallocate memory associated with this package
1224  !<
1225  subroutine apt_da(this)
1226  ! -- modules
1228  ! -- dummy
1229  class(tspapttype) :: this
1230  ! -- local
1231  !
1232  ! -- deallocate arrays
1233  call mem_deallocate(this%dbuff)
1234  call mem_deallocate(this%qsto)
1235  call mem_deallocate(this%ccterm)
1236  call mem_deallocate(this%strt)
1237  call mem_deallocate(this%lauxvar)
1238  call mem_deallocate(this%xoldpak)
1239  if (this%imatrows == 0) then
1240  call mem_deallocate(this%iboundpak)
1241  call mem_deallocate(this%xnewpak)
1242  end if
1243  call mem_deallocate(this%concbudssm)
1244  call mem_deallocate(this%concfeat)
1245  call mem_deallocate(this%qmfrommvr)
1246  deallocate (this%status)
1247  deallocate (this%featname)
1248  !
1249  ! -- budobj
1250  call this%budobj%budgetobject_da()
1251  deallocate (this%budobj)
1252  nullify (this%budobj)
1253  !
1254  ! -- conc table
1255  if (this%iprconc > 0) then
1256  call this%dvtab%table_da()
1257  deallocate (this%dvtab)
1258  nullify (this%dvtab)
1259  end if
1260  !
1261  ! -- index pointers
1262  call mem_deallocate(this%idxlocnode)
1263  call mem_deallocate(this%idxpakdiag)
1264  call mem_deallocate(this%idxdglo)
1265  call mem_deallocate(this%idxoffdglo)
1266  call mem_deallocate(this%idxsymdglo)
1267  call mem_deallocate(this%idxsymoffdglo)
1268  call mem_deallocate(this%idxfjfdglo)
1269  call mem_deallocate(this%idxfjfoffdglo)
1270  !
1271  ! -- deallocate scalars
1272  call mem_deallocate(this%iauxfpconc)
1273  call mem_deallocate(this%imatrows)
1274  call mem_deallocate(this%iprconc)
1275  call mem_deallocate(this%iconcout)
1276  call mem_deallocate(this%ibudgetout)
1277  call mem_deallocate(this%ibudcsv)
1278  call mem_deallocate(this%igwfaptpak)
1279  call mem_deallocate(this%ncv)
1280  call mem_deallocate(this%idxbudfjf)
1281  call mem_deallocate(this%idxbudgwf)
1282  call mem_deallocate(this%idxbudsto)
1283  call mem_deallocate(this%idxbudtmvr)
1284  call mem_deallocate(this%idxbudfmvr)
1285  call mem_deallocate(this%idxbudaux)
1286  call mem_deallocate(this%idxbudssm)
1287  call mem_deallocate(this%nconcbudssm)
1288  call mem_deallocate(this%idxprepak)
1289  call mem_deallocate(this%idxlastpak)
1290  !
1291  ! -- deallocate scalars in NumericalPackageType
1292  call this%BndType%bnd_da()
1293  end subroutine apt_da
1294 
1295  !> @brief Find corresponding advanced package transport package
1296  !<
1297  subroutine find_apt_package(this)
1298  ! -- modules
1300  ! -- dummy
1301  class(tspapttype) :: this
1302  ! -- local
1303  !
1304  ! -- this routine should never be called
1305  call store_error('Program error: pak_solve not implemented.', &
1306  terminate=.true.)
1307  end subroutine find_apt_package
1308 
1309  !> @brief Set options specific to the TspAptType
1310  !!
1311  !! This routine overrides BndType%bnd_options
1312  !<
1313  subroutine apt_options(this, option, found)
1314  use constantsmodule, only: maxcharlen, dzero
1315  use openspecmodule, only: access, form
1317  ! -- dummy
1318  class(tspapttype), intent(inout) :: this
1319  character(len=*), intent(inout) :: option
1320  logical, intent(inout) :: found
1321  ! -- local
1322  character(len=MAXCHARLEN) :: fname, keyword
1323  ! -- formats
1324  character(len=*), parameter :: fmtaptbin = &
1325  "(4x, a, 1x, a, 1x, ' WILL BE SAVED TO FILE: ', a, &
1326  &/4x, 'OPENED ON UNIT: ', I0)"
1327  !
1328  found = .true.
1329  select case (option)
1330  case ('FLOW_PACKAGE_NAME')
1331  call this%parser%GetStringCaps(this%flowpackagename)
1332  write (this%iout, '(4x,a)') &
1333  'THIS '//trim(adjustl(this%text))//' PACKAGE CORRESPONDS TO A GWF &
1334  &PACKAGE WITH THE NAME '//trim(adjustl(this%flowpackagename))
1335  case ('FLOW_PACKAGE_AUXILIARY_NAME')
1336  call this%parser%GetStringCaps(this%cauxfpconc)
1337  write (this%iout, '(4x,a)') &
1338  'SIMULATED CONCENTRATIONS WILL BE COPIED INTO THE FLOW PACKAGE &
1339  &AUXILIARY VARIABLE WITH THE NAME '//trim(adjustl(this%cauxfpconc))
1340  case ('DEV_NONEXPANDING_MATRIX')
1341  ! -- use an iterative solution where concentration is not solved
1342  ! as part of the matrix. It is instead solved separately with a
1343  ! general mixing equation and then added to the RHS of the GWT
1344  ! equations
1345  call this%parser%DevOpt()
1346  this%imatrows = 0
1347  write (this%iout, '(4x,a)') &
1348  trim(adjustl(this%text))// &
1349  ' WILL NOT ADD ADDITIONAL ROWS TO THE A MATRIX.'
1350  case ('PRINT_CONCENTRATION', 'PRINT_TEMPERATURE')
1351  this%iprconc = 1
1352  write (this%iout, '(4x,a,1x,a,1x,a)') trim(adjustl(this%text))// &
1353  trim(adjustl(this%depvartype))//'S WILL BE PRINTED TO LISTING &
1354  &FILE.'
1355  case ('CONCENTRATION', 'TEMPERATURE')
1356  call this%parser%GetStringCaps(keyword)
1357  if (keyword == 'FILEOUT') then
1358  call this%parser%GetString(fname)
1359  this%iconcout = getunit()
1360  call openfile(this%iconcout, this%iout, fname, 'DATA(BINARY)', &
1361  form, access, 'REPLACE')
1362  write (this%iout, fmtaptbin) &
1363  trim(adjustl(this%text)), trim(adjustl(this%depvartype)), &
1364  trim(fname), this%iconcout
1365  else
1366  write (errmsg, "('Optional', 1x, a, 1X, 'keyword must &
1367  &be followed by FILEOUT')") this%depvartype
1368  call store_error(errmsg)
1369  end if
1370  case ('BUDGET')
1371  call this%parser%GetStringCaps(keyword)
1372  if (keyword == 'FILEOUT') then
1373  call this%parser%GetString(fname)
1374  call assign_iounit(this%ibudgetout, this%inunit, "BUDGET fileout")
1375  call openfile(this%ibudgetout, this%iout, fname, 'DATA(BINARY)', &
1376  form, access, 'REPLACE')
1377  write (this%iout, fmtaptbin) trim(adjustl(this%text)), 'BUDGET', &
1378  trim(fname), this%ibudgetout
1379  else
1380  call store_error('Optional BUDGET keyword must be followed by FILEOUT')
1381  end if
1382  case ('BUDGETCSV')
1383  call this%parser%GetStringCaps(keyword)
1384  if (keyword == 'FILEOUT') then
1385  call this%parser%GetString(fname)
1386  call assign_iounit(this%ibudcsv, this%inunit, "BUDGETCSV fileout")
1387  call openfile(this%ibudcsv, this%iout, fname, 'CSV', &
1388  filstat_opt='REPLACE')
1389  write (this%iout, fmtaptbin) trim(adjustl(this%text)), 'BUDGET CSV', &
1390  trim(fname), this%ibudcsv
1391  else
1392  call store_error('Optional BUDGETCSV keyword must be followed by &
1393  &FILEOUT')
1394  end if
1395  case default
1396  !
1397  ! -- No options found
1398  found = .false.
1399  end select
1400  end subroutine apt_options
1401 
1402  !> @brief Determine dimensions for this advanced package
1403  !<
1404  subroutine apt_read_dimensions(this)
1405  ! -- dummy
1406  class(tspapttype), intent(inout) :: this
1407  ! -- local
1408  integer(I4B) :: ierr
1409  ! -- format
1410  !
1411  ! -- Set a pointer to the GWF LAK Package budobj
1412  if (this%flowpackagename == '') then
1413  this%flowpackagename = this%packName
1414  write (this%iout, '(4x,a)') &
1415  'THE FLOW PACKAGE NAME FOR '//trim(adjustl(this%text))//' WAS NOT &
1416  &SPECIFIED. SETTING FLOW PACKAGE NAME TO '// &
1417  &trim(adjustl(this%flowpackagename))
1418 
1419  end if
1420  call this%find_apt_package()
1421  !
1422  ! -- Set dimensions from the GWF advanced package
1423  this%ncv = this%flowbudptr%ncv
1424  this%maxbound = this%flowbudptr%budterm(this%idxbudgwf)%maxlist
1425  this%nbound = this%maxbound
1426  write (this%iout, '(a, a)') 'SETTING DIMENSIONS FOR PACKAGE ', this%packName
1427  write (this%iout, '(2x,a,i0)') 'NUMBER OF CONTROL VOLUMES = ', this%ncv
1428  write (this%iout, '(2x,a,i0)') 'MAXBOUND = ', this%maxbound
1429  write (this%iout, '(2x,a,i0)') 'NBOUND = ', this%nbound
1430  if (this%imatrows /= 0) then
1431  this%npakeq = this%ncv
1432  write (this%iout, '(2x,a)') trim(adjustl(this%text))// &
1433  ' SOLVED AS PART OF GWT MATRIX EQUATIONS'
1434  else
1435  write (this%iout, '(2x,a)') trim(adjustl(this%text))// &
1436  ' SOLVED SEPARATELY FROM GWT MATRIX EQUATIONS '
1437  end if
1438  write (this%iout, '(a, //)') 'DONE SETTING DIMENSIONS FOR '// &
1439  trim(adjustl(this%text))
1440  !
1441  ! -- Check for errors
1442  if (this%ncv < 0) then
1443  write (errmsg, '(a)') &
1444  'Number of control volumes could not be determined correctly.'
1445  call store_error(errmsg)
1446  end if
1447  !
1448  ! -- stop if errors were encountered in the DIMENSIONS block
1449  ierr = count_errors()
1450  if (ierr > 0) then
1451  call store_error_unit(this%inunit)
1452  end if
1453  !
1454  ! -- read packagedata block
1455  call this%apt_read_cvs()
1456  !
1457  ! -- Call define_listlabel to construct the list label that is written
1458  ! when PRINT_INPUT option is used.
1459  call this%define_listlabel()
1460  !
1461  ! -- setup the budget object
1462  call this%apt_setup_budobj()
1463  !
1464  ! -- setup the conc table object
1465  call this%apt_setup_tableobj()
1466  end subroutine apt_read_dimensions
1467 
1468  !> @brief Read feature information for this advanced package
1469  !<
1470  subroutine apt_read_cvs(this)
1471  ! -- modules
1474  ! -- dummy
1475  class(tspapttype), intent(inout) :: this
1476  ! -- local
1477  character(len=LINELENGTH) :: text
1478  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
1479  character(len=9) :: cno
1480  character(len=50), dimension(:), allocatable :: caux
1481  integer(I4B) :: ierr
1482  logical :: isfound, endOfBlock
1483  integer(I4B) :: n
1484  integer(I4B) :: ii, jj
1485  integer(I4B) :: iaux
1486  integer(I4B) :: itmp
1487  integer(I4B) :: nlak
1488  integer(I4B) :: nconn
1489  integer(I4B), dimension(:), pointer, contiguous :: nboundchk
1490  real(DP), pointer :: bndElem => null()
1491  !
1492  ! -- initialize itmp
1493  itmp = 0
1494  !
1495  ! -- allocate apt data
1496  call mem_allocate(this%strt, this%ncv, 'STRT', this%memoryPath)
1497  call mem_allocate(this%lauxvar, this%naux, this%ncv, 'LAUXVAR', &
1498  this%memoryPath)
1499  !
1500  ! -- lake boundary and concentrations
1501  if (this%imatrows == 0) then
1502  call mem_allocate(this%iboundpak, this%ncv, 'IBOUND', this%memoryPath)
1503  call mem_allocate(this%xnewpak, this%ncv, 'XNEWPAK', this%memoryPath)
1504  end if
1505  call mem_allocate(this%xoldpak, this%ncv, 'XOLDPAK', this%memoryPath)
1506  !
1507  ! -- allocate character storage not managed by the memory manager
1508  allocate (this%featname(this%ncv)) ! ditch after boundnames allocated??
1509  !allocate(this%status(this%ncv))
1510  !
1511  do n = 1, this%ncv
1512  this%strt(n) = dep20
1513  this%lauxvar(:, n) = dzero
1514  this%xoldpak(n) = dep20
1515  if (this%imatrows == 0) then
1516  this%iboundpak(n) = 1
1517  this%xnewpak(n) = dep20
1518  end if
1519  end do
1520  !
1521  ! -- allocate local storage for aux variables
1522  if (this%naux > 0) then
1523  allocate (caux(this%naux))
1524  end if
1525  !
1526  ! -- allocate and initialize temporary variables
1527  allocate (nboundchk(this%ncv))
1528  do n = 1, this%ncv
1529  nboundchk(n) = 0
1530  end do
1531  !
1532  ! -- get packagedata block
1533  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
1534  supportopenclose=.true.)
1535  !
1536  ! -- parse locations block if detected
1537  if (isfound) then
1538  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text))// &
1539  ' PACKAGEDATA'
1540  nlak = 0
1541  nconn = 0
1542  do
1543  call this%parser%GetNextLine(endofblock)
1544  if (endofblock) exit
1545  n = this%parser%GetInteger()
1546 
1547  if (n < 1 .or. n > this%ncv) then
1548  write (errmsg, '(a,1x,i6)') &
1549  'Itemno must be > 0 and <= ', this%ncv
1550  call store_error(errmsg)
1551  cycle
1552  end if
1553  !
1554  ! -- increment nboundchk
1555  nboundchk(n) = nboundchk(n) + 1
1556  !
1557  ! -- strt
1558  this%strt(n) = this%parser%GetDouble()
1559  !
1560  ! -- get aux data
1561  do iaux = 1, this%naux
1562  call this%parser%GetString(caux(iaux))
1563  end do
1564 
1565  ! -- set default bndName
1566  write (cno, '(i9.9)') n
1567  bndname = 'Feature'//cno
1568 
1569  ! -- featname
1570  if (this%inamedbound /= 0) then
1571  call this%parser%GetStringCaps(bndnametemp)
1572  if (bndnametemp /= '') then
1573  bndname = bndnametemp
1574  end if
1575  end if
1576  this%featname(n) = bndname
1577 
1578  ! -- fill time series aware data
1579  ! -- fill aux data
1580  do jj = 1, this%naux
1581  text = caux(jj)
1582  ii = n
1583  bndelem => this%lauxvar(jj, ii)
1584  call read_value_or_time_series_adv(text, ii, jj, bndelem, &
1585  this%packName, 'AUX', &
1586  this%tsManager, this%iprpak, &
1587  this%auxname(jj))
1588  end do
1589  !
1590  nlak = nlak + 1
1591  end do
1592  !
1593  ! -- check for duplicate or missing lakes
1594  do n = 1, this%ncv
1595  if (nboundchk(n) == 0) then
1596  write (errmsg, '(a,1x,i0)') 'No data specified for feature', n
1597  call store_error(errmsg)
1598  else if (nboundchk(n) > 1) then
1599  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
1600  'Data for feature', n, 'specified', nboundchk(n), 'times'
1601  call store_error(errmsg)
1602  end if
1603  end do
1604  !
1605  write (this%iout, '(1x,a)') &
1606  'END OF '//trim(adjustl(this%text))//' PACKAGEDATA'
1607  else
1608  call store_error('Required packagedata block not found.')
1609  end if
1610  !
1611  ! -- terminate if any errors were detected
1612  if (count_errors() > 0) then
1613  call this%parser%StoreErrorUnit()
1614  end if
1615  !
1616  ! -- deallocate local storage for aux variables
1617  if (this%naux > 0) then
1618  deallocate (caux)
1619  end if
1620  !
1621  ! -- deallocate local storage for nboundchk
1622  deallocate (nboundchk)
1623  end subroutine apt_read_cvs
1624 
1625  !> @brief Read the initial parameters for an advanced package
1626  !<
1627  subroutine apt_read_initial_attr(this)
1628  use constantsmodule, only: linelength
1629  use budgetmodule, only: budget_cr
1630  ! -- dummy
1631  class(tspapttype), intent(inout) :: this
1632  ! -- local
1633  !character(len=LINELENGTH) :: text
1634  integer(I4B) :: j, n
1635 
1636  !
1637  ! -- initialize xnewpak and set feature concentration (or temperature)
1638  ! -- todo: this should be a time series?
1639  do n = 1, this%ncv
1640  this%xnewpak(n) = this%strt(n)
1641  !
1642  ! -- todo: read aux
1643  !
1644  ! -- todo: read boundname
1645  end do
1646  !
1647  ! -- initialize status (iboundpak) of lakes to active
1648  do n = 1, this%ncv
1649  if (this%status(n) == 'CONSTANT') then
1650  this%iboundpak(n) = -1
1651  else if (this%status(n) == 'INACTIVE') then
1652  this%iboundpak(n) = 0
1653  else if (this%status(n) == 'ACTIVE ') then
1654  this%iboundpak(n) = 1
1655  end if
1656  end do
1657  !
1658  ! -- set boundname for each connection
1659  if (this%inamedbound /= 0) then
1660  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
1661  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
1662  this%boundname(j) = this%featname(n)
1663  end do
1664  end if
1665  !
1666  ! -- copy boundname into boundname_cst
1667  call this%copy_boundname()
1668  end subroutine apt_read_initial_attr
1669 
1670  !> @brief Add terms specific to advanced package transport to the explicit
1671  !! solve
1672  !!
1673  !! Explicit solve for concentration (or temperature) in advaced package
1674  !! features, which is an alternative to the iterative implicit solve.
1675  !<
1676  subroutine apt_solve(this)
1677  use constantsmodule, only: linelength
1678  ! -- dummy
1679  class(tspapttype) :: this
1680  ! -- local
1681  integer(I4B) :: n, j, igwfnode
1682  integer(I4B) :: n1, n2
1683  real(DP) :: rrate
1684  real(DP) :: ctmp
1685  real(DP) :: c1, qbnd
1686  real(DP) :: hcofval, rhsval
1687  !
1688  ! -- initialize dbuff
1689  do n = 1, this%ncv
1690  this%dbuff(n) = dzero
1691  end do
1692  !
1693  ! -- call the individual package routines to add terms specific to the
1694  ! advanced transport package
1695  call this%pak_solve()
1696  !
1697  ! -- add to mover contribution
1698  if (this%idxbudtmvr /= 0) then
1699  do j = 1, this%flowbudptr%budterm(this%idxbudtmvr)%nlist
1700  call this%apt_tmvr_term(j, n1, n2, rrate)
1701  this%dbuff(n1) = this%dbuff(n1) + rrate
1702  end do
1703  end if
1704  !
1705  ! -- add from mover contribution
1706  if (this%idxbudfmvr /= 0) then
1707  do n1 = 1, size(this%qmfrommvr)
1708  rrate = this%qmfrommvr(n1) ! Will be in terms of energy for heat transport
1709  this%dbuff(n1) = this%dbuff(n1) + rrate
1710  end do
1711  end if
1712  !
1713  ! -- go through each gwf connection and accumulate
1714  ! total mass (or energy) in dbuff mass
1715  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
1716  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
1717  this%hcof(j) = dzero
1718  this%rhs(j) = dzero
1719  igwfnode = this%flowbudptr%budterm(this%idxbudgwf)%id2(j)
1720  qbnd = this%flowbudptr%budterm(this%idxbudgwf)%flow(j)
1721  if (qbnd <= dzero) then
1722  ctmp = this%xnewpak(n)
1723  this%rhs(j) = qbnd * ctmp * this%eqnsclfac
1724  else
1725  ctmp = this%xnew(igwfnode)
1726  this%hcof(j) = -qbnd * this%eqnsclfac
1727  end if
1728  c1 = qbnd * ctmp * this%eqnsclfac
1729  this%dbuff(n) = this%dbuff(n) + c1
1730  end do
1731  !
1732  ! -- go through each "within apt-apt" connection (e.g., lak-lak) and
1733  ! accumulate total mass (or energy) in dbuff mass
1734  if (this%idxbudfjf /= 0) then
1735  do j = 1, this%flowbudptr%budterm(this%idxbudfjf)%nlist
1736  call this%apt_fjf_term(j, n1, n2, rrate)
1737  c1 = rrate
1738  this%dbuff(n1) = this%dbuff(n1) + c1
1739  end do
1740  end if
1741  !
1742  ! -- calculate the feature concentration/temperature
1743  do n = 1, this%ncv
1744  call this%apt_stor_term(n, n1, n2, rrate, rhsval, hcofval)
1745  !
1746  ! -- at this point, dbuff has q * c for all sources, so now
1747  ! add Vold / dt * Cold
1748  this%dbuff(n) = this%dbuff(n) - rhsval
1749  !
1750  ! -- Now to calculate c, need to divide dbuff by hcofval
1751  c1 = -this%dbuff(n) / hcofval
1752  if (this%iboundpak(n) > 0) then
1753  this%xnewpak(n) = c1
1754  end if
1755  end do
1756  end subroutine apt_solve
1757 
1758  !> @brief Add terms specific to advanced package transport features to the
1759  !! explicit solve routine
1760  !!
1761  !! This routine must be overridden by the specific apt package
1762  !<
1763  subroutine pak_solve(this)
1764  ! -- dummy
1765  class(tspapttype) :: this
1766  ! -- local
1767  !
1768  ! -- this routine should never be called
1769  call store_error('Program error: pak_solve not implemented.', &
1770  terminate=.true.)
1771  end subroutine pak_solve
1772 
1773  !> @brief Accumulate constant concentration (or temperature) terms for budget
1774  !<
1775  subroutine apt_accumulate_ccterm(this, ilak, rrate, ccratin, ccratout)
1776  ! -- dummy
1777  class(tspapttype) :: this
1778  integer(I4B), intent(in) :: ilak
1779  real(DP), intent(in) :: rrate
1780  real(DP), intent(inout) :: ccratin
1781  real(DP), intent(inout) :: ccratout
1782  ! -- locals
1783  real(DP) :: q
1784  ! format
1785  ! code
1786  !
1787  if (this%iboundpak(ilak) < 0) then
1788  q = -rrate
1789  this%ccterm(ilak) = this%ccterm(ilak) + q
1790  !
1791  ! -- See if flow is into lake or out of lake.
1792  if (q < dzero) then
1793  !
1794  ! -- Flow is out of lake subtract rate from ratout.
1795  ccratout = ccratout - q
1796  else
1797  !
1798  ! -- Flow is into lake; add rate to ratin.
1799  ccratin = ccratin + q
1800  end if
1801  end if
1802  end subroutine apt_accumulate_ccterm
1803 
1804  !> @brief Define the list heading that is written to iout when PRINT_INPUT
1805  !! option is used.
1806  !<
1807  subroutine define_listlabel(this)
1808  class(tspapttype), intent(inout) :: this
1809  !
1810  ! -- create the header list label
1811  this%listlabel = trim(this%filtyp)//' NO.'
1812  if (this%dis%ndim == 3) then
1813  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
1814  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
1815  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
1816  elseif (this%dis%ndim == 2) then
1817  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
1818  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
1819  else
1820  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
1821  end if
1822  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
1823  if (this%inamedbound == 1) then
1824  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
1825  end if
1826  end subroutine define_listlabel
1827 
1828  !> @brief Set pointers to model arrays and variables so that a package has
1829  !! access to these items.
1830  !<
1831  subroutine apt_set_pointers(this, neq, ibound, xnew, xold, flowja)
1832  class(tspapttype) :: this
1833  integer(I4B), pointer :: neq
1834  integer(I4B), dimension(:), pointer, contiguous :: ibound
1835  real(DP), dimension(:), pointer, contiguous :: xnew
1836  real(DP), dimension(:), pointer, contiguous :: xold
1837  real(DP), dimension(:), pointer, contiguous :: flowja
1838  ! -- local
1839  integer(I4B) :: istart, iend
1840  !
1841  ! -- call base BndType set_pointers
1842  call this%BndType%set_pointers(neq, ibound, xnew, xold, flowja)
1843  !
1844  ! -- Set the pointers
1845  !
1846  ! -- set package pointers
1847  if (this%imatrows /= 0) then
1848  istart = this%dis%nodes + this%ioffset + 1
1849  iend = istart + this%ncv - 1
1850  this%iboundpak => this%ibound(istart:iend)
1851  this%xnewpak => this%xnew(istart:iend)
1852  end if
1853  end subroutine apt_set_pointers
1854 
1855  !> @brief Return the feature new volume and old volume
1856  !<
1857  subroutine apt_get_volumes(this, icv, vnew, vold, delt)
1858  ! -- modules
1859  ! -- dummy
1860  class(tspapttype) :: this
1861  integer(I4B), intent(in) :: icv
1862  real(DP), intent(inout) :: vnew, vold
1863  real(DP), intent(in) :: delt
1864  ! -- local
1865  real(DP) :: qss
1866  !
1867  ! -- get volumes
1868  vold = dzero
1869  vnew = vold
1870  if (this%idxbudsto /= 0) then
1871  qss = this%flowbudptr%budterm(this%idxbudsto)%flow(icv)
1872  vnew = this%flowbudptr%budterm(this%idxbudsto)%auxvar(1, icv)
1873  vold = vnew + qss * delt
1874  end if
1875  end subroutine apt_get_volumes
1876 
1877  !> @brief Function to return the number of budget terms just for this package
1878  !!
1879  !! This function must be overridden.
1880  !<
1881  function pak_get_nbudterms(this) result(nbudterms)
1882  ! -- modules
1883  ! -- dummy
1884  class(tspapttype) :: this
1885  ! -- return
1886  integer(I4B) :: nbudterms
1887  ! -- local
1888  !
1889  ! -- this routine should never be called
1890  call store_error('Program error: pak_get_nbudterms not implemented.', &
1891  terminate=.true.)
1892  nbudterms = 0
1893  end function pak_get_nbudterms
1894 
1895  !> @brief Set up the budget object that stores advanced package flow terms
1896  !<
1897  subroutine apt_setup_budobj(this)
1898  ! -- modules
1899  use constantsmodule, only: lenbudtxt
1900  ! -- dummy
1901  class(tspapttype) :: this
1902  ! -- local
1903  integer(I4B) :: nbudterm
1904  integer(I4B) :: nlen
1905  integer(I4B) :: n, n1, n2
1906  integer(I4B) :: maxlist, naux
1907  integer(I4B) :: idx
1908  logical :: ordered_id1
1909  real(DP) :: q
1910  character(len=LENBUDTXT) :: bddim_opt
1911  character(len=LENBUDTXT) :: text, textt
1912  character(len=LENBUDTXT), dimension(1) :: auxtxt
1913  !
1914  ! -- initialize nbudterm
1915  nbudterm = 0
1916  !
1917  ! -- Determine if there are flow-ja-face terms
1918  nlen = 0
1919  if (this%idxbudfjf /= 0) then
1920  nlen = this%flowbudptr%budterm(this%idxbudfjf)%maxlist
1921  end if
1922  !
1923  ! -- Determine the number of budget terms associated with apt.
1924  ! These are fixed for the simulation and cannot change
1925  !
1926  ! -- add one if flow-ja-face present
1927  if (this%idxbudfjf /= 0) nbudterm = nbudterm + 1
1928  !
1929  ! -- All the APT packages have GWF, STORAGE, and CONSTANT
1930  nbudterm = nbudterm + 3
1931  !
1932  ! -- add terms for the specific package
1933  nbudterm = nbudterm + this%pak_get_nbudterms()
1934  !
1935  ! -- add for mover terms and auxiliary
1936  if (this%idxbudtmvr /= 0) nbudterm = nbudterm + 1
1937  if (this%idxbudfmvr /= 0) nbudterm = nbudterm + 1
1938  if (this%naux > 0) nbudterm = nbudterm + 1
1939  !
1940  ! -- set up budobj
1941  call budgetobject_cr(this%budobj, this%packName)
1942  !
1943  bddim_opt = this%depvarunitabbrev
1944  call this%budobj%budgetobject_df(this%ncv, nbudterm, 0, 0, &
1945  bddim_opt=bddim_opt, ibudcsv=this%ibudcsv)
1946  idx = 0
1947  !
1948  ! -- Go through and set up each budget term
1949  if (nlen > 0) then
1950  text = ' FLOW-JA-FACE'
1951  idx = idx + 1
1952  maxlist = this%flowbudptr%budterm(this%idxbudfjf)%maxlist
1953  naux = 0
1954  ordered_id1 = this%flowbudptr%budterm(this%idxbudfjf)%ordered_id1
1955  call this%budobj%budterm(idx)%initialize(text, &
1956  this%name_model, &
1957  this%packName, &
1958  this%name_model, &
1959  this%packName, &
1960  maxlist, .false., .false., &
1961  naux, ordered_id1=ordered_id1)
1962  !
1963  ! -- store outlet connectivity
1964  call this%budobj%budterm(idx)%reset(maxlist)
1965  q = dzero
1966  do n = 1, maxlist
1967  n1 = this%flowbudptr%budterm(this%idxbudfjf)%id1(n)
1968  n2 = this%flowbudptr%budterm(this%idxbudfjf)%id2(n)
1969  call this%budobj%budterm(idx)%update_term(n1, n2, q)
1970  end do
1971  end if
1972  !
1973  ! --
1974  text = ' GWF'
1975  idx = idx + 1
1976  maxlist = this%flowbudptr%budterm(this%idxbudgwf)%maxlist
1977  naux = 0
1978  call this%budobj%budterm(idx)%initialize(text, &
1979  this%name_model, &
1980  this%packName, &
1981  this%name_model, &
1982  this%name_model, &
1983  maxlist, .false., .true., &
1984  naux)
1985  call this%budobj%budterm(idx)%reset(maxlist)
1986  q = dzero
1987  do n = 1, maxlist
1988  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(n)
1989  n2 = this%flowbudptr%budterm(this%idxbudgwf)%id2(n)
1990  call this%budobj%budterm(idx)%update_term(n1, n2, q)
1991  end do
1992  !
1993  ! -- Reserve space for the package specific terms
1994  this%idxprepak = idx
1995  call this%pak_setup_budobj(idx)
1996  this%idxlastpak = idx
1997  !
1998  ! --
1999  text = ' STORAGE'
2000  idx = idx + 1
2001  maxlist = this%flowbudptr%budterm(this%idxbudsto)%maxlist
2002  naux = 1
2003  write (textt, '(a)') str_pad_left(this%depvarunit, 16)
2004  auxtxt(1) = textt ! ' MASS' or ' ENERGY'
2005  call this%budobj%budterm(idx)%initialize(text, &
2006  this%name_model, &
2007  this%packName, &
2008  this%name_model, &
2009  this%packName, &
2010  maxlist, .false., .false., &
2011  naux, auxtxt)
2012  if (this%idxbudtmvr /= 0) then
2013  !
2014  ! --
2015  text = ' TO-MVR'
2016  idx = idx + 1
2017  maxlist = this%flowbudptr%budterm(this%idxbudtmvr)%maxlist
2018  naux = 0
2019  ordered_id1 = this%flowbudptr%budterm(this%idxbudtmvr)%ordered_id1
2020  call this%budobj%budterm(idx)%initialize(text, &
2021  this%name_model, &
2022  this%packName, &
2023  this%name_model, &
2024  this%packName, &
2025  maxlist, .false., .false., &
2026  naux, ordered_id1=ordered_id1)
2027  end if
2028  if (this%idxbudfmvr /= 0) then
2029  !
2030  ! --
2031  text = ' FROM-MVR'
2032  idx = idx + 1
2033  maxlist = this%ncv
2034  naux = 0
2035  call this%budobj%budterm(idx)%initialize(text, &
2036  this%name_model, &
2037  this%packName, &
2038  this%name_model, &
2039  this%packName, &
2040  maxlist, .false., .false., &
2041  naux)
2042  end if
2043  !
2044  ! --
2045  text = ' CONSTANT'
2046  idx = idx + 1
2047  maxlist = this%ncv
2048  naux = 0
2049  call this%budobj%budterm(idx)%initialize(text, &
2050  this%name_model, &
2051  this%packName, &
2052  this%name_model, &
2053  this%packName, &
2054  maxlist, .false., .false., &
2055  naux)
2056 
2057  !
2058  ! --
2059  naux = this%naux
2060  if (naux > 0) then
2061  !
2062  ! --
2063  text = ' AUXILIARY'
2064  idx = idx + 1
2065  maxlist = this%ncv
2066  call this%budobj%budterm(idx)%initialize(text, &
2067  this%name_model, &
2068  this%packName, &
2069  this%name_model, &
2070  this%packName, &
2071  maxlist, .false., .false., &
2072  naux, this%auxname)
2073  end if
2074  !
2075  ! -- if flow for each control volume are written to the listing file
2076  if (this%iprflow /= 0) then
2077  call this%budobj%flowtable_df(this%iout)
2078  end if
2079  end subroutine apt_setup_budobj
2080 
2081  !> @brief Set up a budget object that stores an advanced package flows
2082  !!
2083  !! Individual packages set up their budget terms. Must be overridden.
2084  !<
2085  subroutine pak_setup_budobj(this, idx)
2086  ! -- modules
2087  ! -- dummy
2088  class(tspapttype) :: this
2089  integer(I4B), intent(inout) :: idx
2090  ! -- local
2091  !
2092  ! -- this routine should never be called
2093  call store_error('Program error: pak_setup_budobj not implemented.', &
2094  terminate=.true.)
2095  end subroutine pak_setup_budobj
2096 
2097  !> @brief Copy flow terms into this%budobj
2098  !<
2099  subroutine apt_fill_budobj(this, x, flowja)
2100  ! -- modules
2101  use tdismodule, only: delt
2102  ! -- dummy
2103  class(tspapttype) :: this
2104  real(DP), dimension(:), intent(in) :: x
2105  real(DP), dimension(:), contiguous, intent(inout) :: flowja
2106  ! -- local
2107  integer(I4B) :: naux
2108  real(DP), dimension(:), allocatable :: auxvartmp
2109  integer(I4B) :: i, j, n1, n2
2110  integer(I4B) :: idx
2111  integer(I4B) :: nlen
2112  integer(I4B) :: nlist
2113  integer(I4B) :: igwfnode
2114  real(DP) :: q
2115  real(DP) :: v0, v1
2116  real(DP) :: ccratin, ccratout
2117  ! -- formats
2118  !
2119  ! -- initialize counter
2120  idx = 0
2121  !
2122  ! -- initialize ccterm, which is used to sum up all mass (or energy) flows
2123  ! into a constant concentration (or temperature) cell
2124  ccratin = dzero
2125  ccratout = dzero
2126  do n1 = 1, this%ncv
2127  this%ccterm(n1) = dzero
2128  end do
2129  !
2130  ! -- FLOW JA FACE
2131  nlen = 0
2132  if (this%idxbudfjf /= 0) then
2133  nlen = this%flowbudptr%budterm(this%idxbudfjf)%maxlist
2134  end if
2135  if (nlen > 0) then
2136  idx = idx + 1
2137  nlist = this%flowbudptr%budterm(this%idxbudfjf)%maxlist
2138  call this%budobj%budterm(idx)%reset(nlist)
2139  q = dzero
2140  do j = 1, nlist
2141  call this%apt_fjf_term(j, n1, n2, q)
2142  call this%budobj%budterm(idx)%update_term(n1, n2, q)
2143  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
2144  end do
2145  end if
2146  !
2147  ! -- GWF (LEAKAGE)
2148  idx = idx + 1
2149  call this%budobj%budterm(idx)%reset(this%maxbound)
2150  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
2151  q = dzero
2152  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
2153  if (this%iboundpak(n1) /= 0) then
2154  igwfnode = this%flowbudptr%budterm(this%idxbudgwf)%id2(j)
2155  q = this%hcof(j) * x(igwfnode) - this%rhs(j)
2156  q = -q ! flip sign so relative to advanced package feature
2157  end if
2158  call this%budobj%budterm(idx)%update_term(n1, igwfnode, q)
2159  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
2160  end do
2161  !
2162  ! -- skip individual package terms for now and process them last
2163  ! -- in case they depend on the other terms (as for uze)
2164  idx = this%idxlastpak
2165  !
2166  ! -- STORAGE
2167  idx = idx + 1
2168  call this%budobj%budterm(idx)%reset(this%ncv)
2169  allocate (auxvartmp(1))
2170  do n1 = 1, this%ncv
2171  call this%apt_get_volumes(n1, v1, v0, delt)
2172  auxvartmp(1) = v1 * this%xnewpak(n1) ! Note: When GWE is added, check if this needs a factor of eqnsclfac
2173  q = this%qsto(n1)
2174  call this%budobj%budterm(idx)%update_term(n1, n1, q, auxvartmp)
2175  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
2176  end do
2177  deallocate (auxvartmp)
2178  !
2179  ! -- TO MOVER
2180  if (this%idxbudtmvr /= 0) then
2181  idx = idx + 1
2182  nlist = this%flowbudptr%budterm(this%idxbudtmvr)%nlist
2183  call this%budobj%budterm(idx)%reset(nlist)
2184  do j = 1, nlist
2185  call this%apt_tmvr_term(j, n1, n2, q)
2186  call this%budobj%budterm(idx)%update_term(n1, n2, q)
2187  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
2188  end do
2189  end if
2190  !
2191  ! -- FROM MOVER
2192  if (this%idxbudfmvr /= 0) then
2193  idx = idx + 1
2194  nlist = this%ncv
2195  call this%budobj%budterm(idx)%reset(nlist)
2196  do j = 1, nlist
2197  call this%apt_fmvr_term(j, n1, n2, q)
2198  call this%budobj%budterm(idx)%update_term(n1, n1, q)
2199  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
2200  end do
2201  end if
2202  !
2203  ! -- CONSTANT FLOW
2204  idx = idx + 1
2205  call this%budobj%budterm(idx)%reset(this%ncv)
2206  do n1 = 1, this%ncv
2207  q = this%ccterm(n1)
2208  call this%budobj%budterm(idx)%update_term(n1, n1, q)
2209  end do
2210  !
2211  ! -- AUXILIARY VARIABLES
2212  naux = this%naux
2213  if (naux > 0) then
2214  idx = idx + 1
2215  allocate (auxvartmp(naux))
2216  call this%budobj%budterm(idx)%reset(this%ncv)
2217  do n1 = 1, this%ncv
2218  q = dzero
2219  do i = 1, naux
2220  auxvartmp(i) = this%lauxvar(i, n1)
2221  end do
2222  call this%budobj%budterm(idx)%update_term(n1, n1, q, auxvartmp)
2223  end do
2224  deallocate (auxvartmp)
2225  end if
2226  !
2227  ! -- individual package terms processed last
2228  idx = this%idxprepak
2229  call this%pak_fill_budobj(idx, x, flowja, ccratin, ccratout)
2230  !
2231  ! --Terms are filled, now accumulate them for this time step
2232  call this%budobj%accumulate_terms()
2233  end subroutine apt_fill_budobj
2234 
2235  !> @brief Copy flow terms into this%budobj, must be overridden
2236  !<
2237  subroutine pak_fill_budobj(this, idx, x, flowja, ccratin, ccratout)
2238  ! -- modules
2239  ! -- dummy
2240  class(tspapttype) :: this
2241  integer(I4B), intent(inout) :: idx
2242  real(DP), dimension(:), intent(in) :: x
2243  real(DP), dimension(:), contiguous, intent(inout) :: flowja
2244  real(DP), intent(inout) :: ccratin
2245  real(DP), intent(inout) :: ccratout
2246  ! -- local
2247  ! -- formats
2248  !
2249  ! -- this routine should never be called
2250  call store_error('Program error: pak_fill_budobj not implemented.', &
2251  terminate=.true.)
2252  end subroutine pak_fill_budobj
2253 
2254  !> @brief Account for mass or energy storage in advanced package features
2255  !<
2256  subroutine apt_stor_term(this, ientry, n1, n2, rrate, &
2257  rhsval, hcofval)
2258  use tdismodule, only: delt
2259  class(tspapttype) :: this
2260  integer(I4B), intent(in) :: ientry
2261  integer(I4B), intent(inout) :: n1
2262  integer(I4B), intent(inout) :: n2
2263  real(DP), intent(inout), optional :: rrate
2264  real(DP), intent(inout), optional :: rhsval
2265  real(DP), intent(inout), optional :: hcofval
2266  real(DP) :: v0, v1
2267  real(DP) :: c0, c1
2268  !
2269  n1 = ientry
2270  n2 = ientry
2271  call this%apt_get_volumes(n1, v1, v0, delt)
2272  c0 = this%xoldpak(n1)
2273  c1 = this%xnewpak(n1)
2274  if (present(rrate)) then
2275  rrate = (-c1 * v1 / delt + c0 * v0 / delt) * this%eqnsclfac
2276  end if
2277  !
2278  if (present(rhsval)) rhsval = -c0 * v0 * this%eqnsclfac / delt
2279  if (present(hcofval)) hcofval = -v1 * this%eqnsclfac / delt
2280  end subroutine apt_stor_term
2281 
2282  !> @brief Account for mass or energy transferred to the MVR package
2283  !<
2284  subroutine apt_tmvr_term(this, ientry, n1, n2, rrate, &
2285  rhsval, hcofval)
2286  ! -- modules
2287  ! -- dummy
2288  class(tspapttype) :: this
2289  integer(I4B), intent(in) :: ientry
2290  integer(I4B), intent(inout) :: n1
2291  integer(I4B), intent(inout) :: n2
2292  real(DP), intent(inout), optional :: rrate
2293  real(DP), intent(inout), optional :: rhsval
2294  real(DP), intent(inout), optional :: hcofval
2295  ! -- local
2296  real(DP) :: qbnd
2297  real(DP) :: ctmp
2298  !
2299  ! -- Calculate MVR-related terms
2300  n1 = this%flowbudptr%budterm(this%idxbudtmvr)%id1(ientry)
2301  n2 = this%flowbudptr%budterm(this%idxbudtmvr)%id2(ientry)
2302  qbnd = this%flowbudptr%budterm(this%idxbudtmvr)%flow(ientry)
2303  ctmp = this%xnewpak(n1)
2304  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
2305  if (present(rhsval)) rhsval = dzero
2306  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac
2307  end subroutine apt_tmvr_term
2308 
2309  !> @brief Account for mass or energy transferred to this package from the
2310  !! MVR package
2311  !<
2312  subroutine apt_fmvr_term(this, ientry, n1, n2, rrate, &
2313  rhsval, hcofval)
2314  ! -- modules
2315  ! -- dummy
2316  class(tspapttype) :: this
2317  integer(I4B), intent(in) :: ientry
2318  integer(I4B), intent(inout) :: n1
2319  integer(I4B), intent(inout) :: n2
2320  real(DP), intent(inout), optional :: rrate
2321  real(DP), intent(inout), optional :: rhsval
2322  real(DP), intent(inout), optional :: hcofval
2323  !
2324  ! -- Calculate MVR-related terms
2325  n1 = ientry
2326  n2 = n1
2327  if (present(rrate)) rrate = this%qmfrommvr(n1) ! NOTE: When bringing in GWE, ensure this is in terms of energy. Might need to apply eqnsclfac here.
2328  if (present(rhsval)) rhsval = this%qmfrommvr(n1)
2329  if (present(hcofval)) hcofval = dzero
2330  end subroutine apt_fmvr_term
2331 
2332  !> @brief Go through each "within apt-apt" connection (e.g., lkt-lkt, or
2333  !! sft-sft) and accumulate total mass (or energy) in dbuff mass
2334  !<
2335  subroutine apt_fjf_term(this, ientry, n1, n2, rrate, &
2336  rhsval, hcofval)
2337  ! -- modules
2338  ! -- dummy
2339  class(tspapttype) :: this
2340  integer(I4B), intent(in) :: ientry
2341  integer(I4B), intent(inout) :: n1
2342  integer(I4B), intent(inout) :: n2
2343  real(DP), intent(inout), optional :: rrate
2344  real(DP), intent(inout), optional :: rhsval
2345  real(DP), intent(inout), optional :: hcofval
2346  ! -- local
2347  real(DP) :: qbnd
2348  real(DP) :: ctmp
2349  !
2350  n1 = this%flowbudptr%budterm(this%idxbudfjf)%id1(ientry)
2351  n2 = this%flowbudptr%budterm(this%idxbudfjf)%id2(ientry)
2352  qbnd = this%flowbudptr%budterm(this%idxbudfjf)%flow(ientry)
2353  if (qbnd <= 0) then
2354  ctmp = this%xnewpak(n1)
2355  else
2356  ctmp = this%xnewpak(n2)
2357  end if
2358  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
2359  if (present(rhsval)) rhsval = -rrate * this%eqnsclfac
2360  if (present(hcofval)) hcofval = dzero
2361  end subroutine apt_fjf_term
2362 
2363  !> @brief Copy concentrations (or temperatures) into flow package aux
2364  !! variable
2365  !<
2366  subroutine apt_copy2flowp(this)
2367  ! -- modules
2368  ! -- dummy
2369  class(tspapttype) :: this
2370  ! -- local
2371  integer(I4B) :: n, j
2372  !
2373  ! -- copy
2374  if (this%iauxfpconc /= 0) then
2375  !
2376  ! -- go through each apt-gwf connection
2377  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
2378  !
2379  ! -- set n to feature number and process if active feature
2380  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
2381  this%flowpackagebnd%auxvar(this%iauxfpconc, j) = this%xnewpak(n)
2382  end do
2383  end if
2384  end subroutine apt_copy2flowp
2385 
2386  !> @brief Determine whether an obs type is supported
2387  !!
2388  !! This function:
2389  !! - returns true if APT package supports named observation.
2390  !! - overrides BndType%bnd_obs_supported()
2391  !<
2392  logical function apt_obs_supported(this)
2393  ! -- modules
2394  ! -- dummy
2395  class(tspapttype) :: this
2396  !
2397  ! -- Set to true
2398  apt_obs_supported = .true.
2399  end function apt_obs_supported
2400 
2401  !> @brief Define observation type
2402  !!
2403  !! This routine:
2404  !! - stores observation types supported by APT package.
2405  !! - overrides BndType%bnd_df_obs
2406  !<
2407  subroutine apt_df_obs(this)
2408  ! -- modules
2409  ! -- dummy
2410  class(tspapttype) :: this
2411  ! -- local
2412  !
2413  ! -- call additional specific observations for lkt, sft, mwt, and uzt
2414  call this%pak_df_obs()
2415  end subroutine apt_df_obs
2416 
2417  !> @brief Define apt observation type
2418  !!
2419  !! This routine:
2420  !! - stores observations supported by the APT package
2421  !! - must be overridden by child class
2422  subroutine pak_df_obs(this)
2423  ! -- modules
2424  ! -- dummy
2425  class(tspapttype) :: this
2426  ! -- local
2427  !
2428  ! -- this routine should never be called
2429  call store_error('Program error: pak_df_obs not implemented.', &
2430  terminate=.true.)
2431  end subroutine pak_df_obs
2432 
2433  !> @brief Process package specific obs
2434  !!
2435  !! Method to process specific observations for this package.
2436  !<
2437  subroutine pak_rp_obs(this, obsrv, found)
2438  ! -- dummy
2439  class(tspapttype), intent(inout) :: this !< package class
2440  type(observetype), intent(inout) :: obsrv !< observation object
2441  logical, intent(inout) :: found !< indicate whether observation was found
2442  ! -- local
2443  !
2444  ! -- this routine should never be called
2445  call store_error('Program error: pak_rp_obs not implemented.', &
2446  terminate=.true.)
2447  end subroutine pak_rp_obs
2448 
2449  !> @brief Prepare observation
2450  !!
2451  !! Find the indices for this observation assuming they are indexed by
2452  !! feature number
2453  !<
2454  subroutine rp_obs_byfeature(this, obsrv)
2455  class(tspapttype), intent(inout) :: this !< object
2456  type(observetype), intent(inout) :: obsrv !< observation
2457  integer(I4B) :: nn1
2458  integer(I4B) :: j
2459  logical :: jfound
2460  character(len=*), parameter :: fmterr = &
2461  "('Boundary ', a, ' for observation ', a, &
2462  &' is invalid in package ', a)"
2463  nn1 = obsrv%NodeNumber
2464  if (nn1 == namedboundflag) then
2465  jfound = .false.
2466  do j = 1, this%ncv
2467  if (this%featname(j) == obsrv%FeatureName) then
2468  jfound = .true.
2469  call obsrv%AddObsIndex(j)
2470  end if
2471  end do
2472  if (.not. jfound) then
2473  write (errmsg, fmterr) trim(obsrv%FeatureName), trim(obsrv%Name), &
2474  trim(this%packName)
2475  call store_error(errmsg)
2476  end if
2477  else
2478  !
2479  ! -- ensure nn1 is > 0 and < ncv
2480  if (nn1 < 0 .or. nn1 > this%ncv) then
2481  write (errmsg, '(7a, i0, a, i0, a)') &
2482  'Observation ', trim(obsrv%Name), ' of type ', &
2483  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2484  trim(this%packName), ' was assigned ID = ', nn1, &
2485  '. ID must be >= 1 and <= ', this%ncv, '.'
2486  call store_error(errmsg)
2487  end if
2488  call obsrv%AddObsIndex(nn1)
2489  end if
2490  end subroutine rp_obs_byfeature
2491 
2492  !> @brief Prepare observation
2493  !!
2494  !! Find the indices for this observation assuming they are first indexed
2495  !! by feature number and secondly by a connection number
2496  !<
2497  subroutine rp_obs_budterm(this, obsrv, budterm)
2498  class(tspapttype), intent(inout) :: this !< object
2499  type(observetype), intent(inout) :: obsrv !< observation
2500  type(budgettermtype), intent(in) :: budterm !< budget term
2501  integer(I4B) :: nn1
2502  integer(I4B) :: iconn
2503  integer(I4B) :: icv
2504  integer(I4B) :: idx
2505  integer(I4B) :: j
2506  logical :: jfound
2507  character(len=*), parameter :: fmterr = &
2508  "('Boundary ', a, ' for observation ', a, &
2509  &' is invalid in package ', a)"
2510  nn1 = obsrv%NodeNumber
2511  if (nn1 == namedboundflag) then
2512  jfound = .false.
2513  do j = 1, budterm%nlist
2514  icv = budterm%id1(j)
2515  if (this%featname(icv) == obsrv%FeatureName) then
2516  jfound = .true.
2517  call obsrv%AddObsIndex(j)
2518  end if
2519  end do
2520  if (.not. jfound) then
2521  write (errmsg, fmterr) trim(obsrv%FeatureName), trim(obsrv%Name), &
2522  trim(this%packName)
2523  call store_error(errmsg)
2524  end if
2525  else
2526  !
2527  ! -- ensure nn1 is > 0 and < ncv
2528  if (nn1 < 0 .or. nn1 > this%ncv) then
2529  write (errmsg, '(7a, i0, a, i0, a)') &
2530  'Observation ', trim(obsrv%Name), ' of type ', &
2531  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2532  trim(this%packName), ' was assigned ID = ', nn1, &
2533  '. ID must be >= 1 and <= ', this%ncv, '.'
2534  call store_error(errmsg)
2535  end if
2536  iconn = obsrv%NodeNumber2
2537  do j = 1, budterm%nlist
2538  if (budterm%id1(j) == nn1) then
2539  ! -- Look for the first occurrence of nn1, then set indxbnds
2540  ! to the iconn record after that
2541  idx = j + iconn - 1
2542  call obsrv%AddObsIndex(idx)
2543  exit
2544  end if
2545  end do
2546  if (idx < 1 .or. idx > budterm%nlist) then
2547  write (errmsg, '(7a, i0, a, i0, a)') &
2548  'Observation ', trim(obsrv%Name), ' of type ', &
2549  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2550  trim(this%packName), ' specifies iconn = ', iconn, &
2551  ', but this is not a valid connection for ID ', nn1, '.'
2552  call store_error(errmsg)
2553  else if (budterm%id1(idx) /= nn1) then
2554  write (errmsg, '(7a, i0, a, i0, a)') &
2555  'Observation ', trim(obsrv%Name), ' of type ', &
2556  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2557  trim(this%packName), ' specifies iconn = ', iconn, &
2558  ', but this is not a valid connection for ID ', nn1, '.'
2559  call store_error(errmsg)
2560  end if
2561  end if
2562  end subroutine rp_obs_budterm
2563 
2564  !> @brief Prepare observation
2565  !!
2566  !! Find the indices for this observation assuming they are first indexed
2567  !! by a feature number and secondly by a second feature number
2568  !<
2569  subroutine rp_obs_flowjaface(this, obsrv, budterm)
2570  class(tspapttype), intent(inout) :: this !< object
2571  type(observetype), intent(inout) :: obsrv !< observation
2572  type(budgettermtype), intent(in) :: budterm !< budget term
2573  integer(I4B) :: nn1
2574  integer(I4B) :: nn2
2575  integer(I4B) :: icv
2576  integer(I4B) :: j
2577  logical :: jfound
2578  character(len=*), parameter :: fmterr = &
2579  "('Boundary ', a, ' for observation ', a, &
2580  &' is invalid in package ', a)"
2581  nn1 = obsrv%NodeNumber
2582  if (nn1 == namedboundflag) then
2583  jfound = .false.
2584  do j = 1, budterm%nlist
2585  icv = budterm%id1(j)
2586  if (this%featname(icv) == obsrv%FeatureName) then
2587  jfound = .true.
2588  call obsrv%AddObsIndex(j)
2589  end if
2590  end do
2591  if (.not. jfound) then
2592  write (errmsg, fmterr) trim(obsrv%FeatureName), trim(obsrv%Name), &
2593  trim(this%packName)
2594  call store_error(errmsg)
2595  end if
2596  else
2597  !
2598  ! -- ensure nn1 is > 0 and < ncv
2599  if (nn1 < 0 .or. nn1 > this%ncv) then
2600  write (errmsg, '(7a, i0, a, i0, a)') &
2601  'Observation ', trim(obsrv%Name), ' of type ', &
2602  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2603  trim(this%packName), ' was assigned ID = ', nn1, &
2604  '. ID must be >= 1 and <= ', this%ncv, '.'
2605  call store_error(errmsg)
2606  end if
2607  nn2 = obsrv%NodeNumber2
2608  !
2609  ! -- ensure nn2 is > 0 and < ncv
2610  if (nn2 < 0 .or. nn2 > this%ncv) then
2611  write (errmsg, '(7a, i0, a, i0, a)') &
2612  'Observation ', trim(obsrv%Name), ' of type ', &
2613  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2614  trim(this%packName), ' was assigned ID2 = ', nn2, &
2615  '. ID must be >= 1 and <= ', this%ncv, '.'
2616  call store_error(errmsg)
2617  end if
2618  ! -- Look for nn1 and nn2 in id1 and id2
2619  jfound = .false.
2620  do j = 1, budterm%nlist
2621  if (budterm%id1(j) == nn1 .and. budterm%id2(j) == nn2) then
2622  call obsrv%AddObsIndex(j)
2623  jfound = .true.
2624  end if
2625  end do
2626  if (.not. jfound) then
2627  write (errmsg, '(7a, i0, a, i0, a)') &
2628  'Observation ', trim(obsrv%Name), ' of type ', &
2629  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2630  trim(this%packName), &
2631  ' specifies a connection between feature ', nn1, &
2632  ' feature ', nn2, ', but these features are not connected.'
2633  call store_error(errmsg)
2634  end if
2635  end if
2636  end subroutine rp_obs_flowjaface
2637 
2638  !> @brief Read and prepare apt-related observations
2639  !!
2640  !! Method to process specific observations for an apt package
2641  !<
2642  subroutine apt_rp_obs(this)
2643  ! -- modules
2644  use tdismodule, only: kper
2645  ! -- dummy
2646  class(tspapttype), intent(inout) :: this
2647  ! -- local
2648  integer(I4B) :: i
2649  logical :: found
2650  class(observetype), pointer :: obsrv => null()
2651  !
2652  if (kper == 1) then
2653  do i = 1, this%obs%npakobs
2654  obsrv => this%obs%pakobs(i)%obsrv
2655  select case (obsrv%ObsTypeId)
2656  case ('CONCENTRATION', 'TEMPERATURE')
2657  call this%rp_obs_byfeature(obsrv)
2658  !
2659  ! -- catch non-cumulative observation assigned to observation defined
2660  ! by a boundname that is assigned to more than one element
2661  if (obsrv%indxbnds_count > 1) then
2662  write (errmsg, '(a, a, a, a)') &
2663  trim(adjustl(this%depvartype))// &
2664  ' for observation', trim(adjustl(obsrv%Name)), &
2665  ' must be assigned to a feature with a unique boundname.'
2666  call store_error(errmsg)
2667  end if
2668  case ('LKT', 'SFT', 'MWT', 'UZT', 'LKE', 'SFE', 'MWE', 'UZE')
2669  call this%rp_obs_budterm(obsrv, &
2670  this%flowbudptr%budterm(this%idxbudgwf))
2671  case ('FLOW-JA-FACE')
2672  if (this%idxbudfjf > 0) then
2673  call this%rp_obs_flowjaface(obsrv, &
2674  this%flowbudptr%budterm(this%idxbudfjf))
2675  else
2676  write (errmsg, '(7a)') &
2677  'Observation ', trim(obsrv%Name), ' of type ', &
2678  trim(adjustl(obsrv%ObsTypeId)), ' in package ', &
2679  trim(this%packName), &
2680  ' cannot be processed because there are no flow connections.'
2681  call store_error(errmsg)
2682  end if
2683  case ('STORAGE')
2684  call this%rp_obs_byfeature(obsrv)
2685  case ('CONSTANT')
2686  call this%rp_obs_byfeature(obsrv)
2687  case ('FROM-MVR')
2688  call this%rp_obs_byfeature(obsrv)
2689  case default
2690  !
2691  ! -- check the child package for any specific obs
2692  found = .false.
2693  call this%pak_rp_obs(obsrv, found)
2694  !
2695  ! -- if none found then terminate with an error
2696  if (.not. found) then
2697  errmsg = 'Unrecognized observation type "'// &
2698  trim(obsrv%ObsTypeId)//'" for '// &
2699  trim(adjustl(this%text))//' package '// &
2700  trim(this%packName)
2701  call store_error(errmsg, terminate=.true.)
2702  end if
2703  end select
2704 
2705  end do
2706  !
2707  ! -- check for errors
2708  if (count_errors() > 0) then
2709  call store_error_unit(this%obs%inunitobs)
2710  end if
2711  end if
2712  end subroutine apt_rp_obs
2713 
2714  !> @brief Calculate observation values
2715  !!
2716  !! Routine calculates observations common to SFT/LKT/MWT/UZT
2717  !! (or SFE/LKE/MWE/UZE) for as many TspAptType observations that are common
2718  !! among the advanced transport packages
2719  !<
2720  subroutine apt_bd_obs(this)
2721  ! -- modules
2722  ! -- dummy
2723  class(tspapttype) :: this
2724  ! -- local
2725  integer(I4B) :: i
2726  integer(I4B) :: igwfnode
2727  integer(I4B) :: j
2728  integer(I4B) :: jj
2729  integer(I4B) :: n
2730  integer(I4B) :: n1
2731  integer(I4B) :: n2
2732  real(DP) :: v
2733  type(observetype), pointer :: obsrv => null()
2734  logical :: found
2735  !
2736  ! -- Write simulated values for all Advanced Package observations
2737  if (this%obs%npakobs > 0) then
2738  call this%obs%obs_bd_clear()
2739  do i = 1, this%obs%npakobs
2740  obsrv => this%obs%pakobs(i)%obsrv
2741  do j = 1, obsrv%indxbnds_count
2742  v = dnodata
2743  jj = obsrv%indxbnds(j)
2744  select case (obsrv%ObsTypeId)
2745  case ('CONCENTRATION', 'TEMPERATURE')
2746  if (this%iboundpak(jj) /= 0) then
2747  v = this%xnewpak(jj)
2748  end if
2749  case ('LKT', 'SFT', 'MWT', 'UZT', 'LKE', 'SFE', 'MWE', 'UZE')
2750  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(jj)
2751  if (this%iboundpak(n) /= 0) then
2752  igwfnode = this%flowbudptr%budterm(this%idxbudgwf)%id2(jj)
2753  v = this%hcof(jj) * this%xnew(igwfnode) - this%rhs(jj)
2754  v = -v
2755  end if
2756  case ('FLOW-JA-FACE')
2757  n = this%flowbudptr%budterm(this%idxbudfjf)%id1(jj)
2758  if (this%iboundpak(n) /= 0) then
2759  call this%apt_fjf_term(jj, n1, n2, v)
2760  end if
2761  case ('STORAGE')
2762  if (this%iboundpak(jj) /= 0) then
2763  v = this%qsto(jj)
2764  end if
2765  case ('CONSTANT')
2766  if (this%iboundpak(jj) /= 0) then
2767  v = this%ccterm(jj)
2768  end if
2769  case ('FROM-MVR')
2770  if (this%iboundpak(jj) /= 0 .and. this%idxbudfmvr > 0) then
2771  call this%apt_fmvr_term(jj, n1, n2, v)
2772  end if
2773  case ('TO-MVR')
2774  if (this%idxbudtmvr > 0) then
2775  n = this%flowbudptr%budterm(this%idxbudtmvr)%id1(jj)
2776  if (this%iboundpak(n) /= 0) then
2777  call this%apt_tmvr_term(jj, n1, n2, v)
2778  end if
2779  end if
2780  case default
2781  found = .false.
2782  !
2783  ! -- check the child package for any specific obs
2784  call this%pak_bd_obs(obsrv%ObsTypeId, jj, v, found)
2785  !
2786  ! -- if none found then terminate with an error
2787  if (.not. found) then
2788  errmsg = 'Unrecognized observation type "'// &
2789  trim(obsrv%ObsTypeId)//'" for '// &
2790  trim(adjustl(this%text))//' package '// &
2791  trim(this%packName)
2792  call store_error(errmsg, terminate=.true.)
2793  end if
2794  end select
2795  call this%obs%SaveOneSimval(obsrv, v)
2796  end do
2797  end do
2798  !
2799  ! -- write summary of error messages
2800  if (count_errors() > 0) then
2801  call store_error_unit(this%obs%inunitobs)
2802  end if
2803  end if
2804  end subroutine apt_bd_obs
2805 
2806  !> @brief Check if observation exists in an advanced package
2807  !<
2808  subroutine pak_bd_obs(this, obstypeid, jj, v, found)
2809  ! -- dummy
2810  class(tspapttype), intent(inout) :: this
2811  character(len=*), intent(in) :: obstypeid
2812  integer(I4B), intent(in) :: jj
2813  real(DP), intent(inout) :: v
2814  logical, intent(inout) :: found
2815  ! -- local
2816  !
2817  ! -- set found = .false. because obstypeid is not known
2818  found = .false.
2819  end subroutine pak_bd_obs
2820 
2821  !> @brief Process observation IDs for an advanced package
2822  !!
2823  !! Method to process observation ID strings for an APT package.
2824  !! This processor is only for observation types that support ID1
2825  !! and not ID2.
2826  !<
2827  subroutine apt_process_obsid(obsrv, dis, inunitobs, iout)
2828  ! -- dummy variables
2829  type(observetype), intent(inout) :: obsrv !< Observation object
2830  class(disbasetype), intent(in) :: dis !< Discretization object
2831  integer(I4B), intent(in) :: inunitobs !< file unit number for the package observation file
2832  integer(I4B), intent(in) :: iout !< model listing file unit number
2833  ! -- local variables
2834  integer(I4B) :: nn1
2835  integer(I4B) :: icol
2836  integer(I4B) :: istart
2837  integer(I4B) :: istop
2838  character(len=LINELENGTH) :: string
2839  character(len=LENBOUNDNAME) :: bndname
2840  !
2841  ! -- initialize local variables
2842  string = obsrv%IDstring
2843  !
2844  ! -- Extract reach number from string and store it.
2845  ! If 1st item is not an integer(I4B), it should be a
2846  ! boundary name--deal with it.
2847  icol = 1
2848  !
2849  ! -- get reach number or boundary name
2850  call extract_idnum_or_bndname(string, icol, istart, istop, nn1, bndname)
2851  if (nn1 == namedboundflag) then
2852  obsrv%FeatureName = bndname
2853  end if
2854  !
2855  ! -- store reach number (NodeNumber)
2856  obsrv%NodeNumber = nn1
2857  !
2858  ! -- store NodeNumber2 as 1 so that this can be used
2859  ! as the iconn value for SFT. This works for SFT
2860  ! because there is only one reach per GWT connection.
2861  obsrv%NodeNumber2 = 1
2862  end subroutine apt_process_obsid
2863 
2864  !> @brief Process observation IDs for a package
2865  !!
2866  !! Method to process observation ID strings for an APT package. This
2867  !! processor is for the case where if ID1 is an integer then ID2 must be
2868  !! provided.
2869  !<
2870  subroutine apt_process_obsid12(obsrv, dis, inunitobs, iout)
2871  ! -- dummy variables
2872  type(observetype), intent(inout) :: obsrv !< Observation object
2873  class(disbasetype), intent(in) :: dis !< Discretization object
2874  integer(I4B), intent(in) :: inunitobs !< file unit number for the package observation file
2875  integer(I4B), intent(in) :: iout !< model listing file unit number
2876  ! -- local variables
2877  integer(I4B) :: nn1
2878  integer(I4B) :: iconn
2879  integer(I4B) :: icol
2880  integer(I4B) :: istart
2881  integer(I4B) :: istop
2882  character(len=LINELENGTH) :: string
2883  character(len=LENBOUNDNAME) :: bndname
2884  !
2885  ! -- initialize local variables
2886  string = obsrv%IDstring
2887  !
2888  ! -- Extract reach number from string and store it.
2889  ! If 1st item is not an integer(I4B), it should be a
2890  ! boundary name--deal with it.
2891  icol = 1
2892  !
2893  ! -- get reach number or boundary name
2894  call extract_idnum_or_bndname(string, icol, istart, istop, nn1, bndname)
2895  if (nn1 == namedboundflag) then
2896  obsrv%FeatureName = bndname
2897  else
2898  call extract_idnum_or_bndname(string, icol, istart, istop, iconn, bndname)
2899  if (len_trim(bndname) < 1 .and. iconn < 0) then
2900  write (errmsg, '(a,1x,a,a,1x,a,1x,a)') &
2901  'For observation type', trim(adjustl(obsrv%ObsTypeId)), &
2902  ', ID given as an integer and not as boundname,', &
2903  'but ID2 is missing. Either change ID to valid', &
2904  'boundname or supply valid entry for ID2.'
2905  call store_error(errmsg)
2906  end if
2907  obsrv%NodeNumber2 = iconn
2908  end if
2909  !
2910  ! -- store reach number (NodeNumber)
2911  obsrv%NodeNumber = nn1
2912  end subroutine apt_process_obsid12
2913 
2914  !> @brief Setup a table object an advanced package
2915  !!
2916  !! Set up the table object that is used to write the apt concentration
2917  !! (or temperature) data. The terms listed here must correspond in the
2918  !! apt_ot method.
2919  !<
2920  subroutine apt_setup_tableobj(this)
2921  ! -- modules
2923  ! -- dummy
2924  class(tspapttype) :: this
2925  ! -- local
2926  integer(I4B) :: nterms
2927  character(len=LINELENGTH) :: title
2928  character(len=LINELENGTH) :: text_temp
2929  !
2930  ! -- setup well head table
2931  if (this%iprconc > 0) then
2932  !
2933  ! -- Determine the number of head table columns
2934  nterms = 2
2935  if (this%inamedbound == 1) nterms = nterms + 1
2936  !
2937  ! -- set up table title
2938  title = trim(adjustl(this%text))//' PACKAGE ('// &
2939  trim(adjustl(this%packName))// &
2940  ') '//trim(adjustl(this%depvartype))// &
2941  &' FOR EACH CONTROL VOLUME'
2942  !
2943  ! -- set up dv tableobj
2944  call table_cr(this%dvtab, this%packName, title)
2945  call this%dvtab%table_df(this%ncv, nterms, this%iout, &
2946  transient=.true.)
2947  !
2948  ! -- Go through and set up table budget term
2949  if (this%inamedbound == 1) then
2950  text_temp = 'NAME'
2951  call this%dvtab%initialize_column(text_temp, 20, alignment=tableft)
2952  end if
2953  !
2954  ! -- feature number
2955  text_temp = 'NUMBER'
2956  call this%dvtab%initialize_column(text_temp, 10, alignment=tabcenter)
2957  !
2958  ! -- feature conc
2959  text_temp = this%depvartype(1:4)
2960  call this%dvtab%initialize_column(text_temp, 12, alignment=tabcenter)
2961  end if
2962  end subroutine apt_setup_tableobj
2963 
2964 end module tspaptmodule
This module contains the base boundary package.
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public budget_cr(this, name_model)
@ brief Create a new budget object
Definition: Budget.f90:84
subroutine, public budgetobject_cr(this, name)
Create a new budget object.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dhdry
real dry cell constant
Definition: Constants.f90:94
@ tabcenter
centered table column
Definition: Constants.f90:172
@ tabright
right justified table column
Definition: Constants.f90:173
@ tableft
left justified table column
Definition: Constants.f90:171
@ tabucstring
upper case string table data
Definition: Constants.f90:180
@ tabstring
string table data
Definition: Constants.f90:179
@ tabreal
real table data
Definition: Constants.f90:182
@ tabinteger
integer table data
Definition: Constants.f90:181
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter namedboundflag
named bound flag
Definition: Constants.f90:49
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
real(dp), parameter dep20
real constant 1e20
Definition: Constants.f90:91
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
integer(i4b), parameter lenauxname
maximum length of a aux variable
Definition: Constants.f90:35
integer(i4b), parameter lenboundname
maximum length of a bound name
Definition: Constants.f90:36
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine, public assign_iounit(iounit, errunit, description)
@ brief assign io unit number
subroutine, public extract_idnum_or_bndname(line, icol, istart, istop, idnum, bndname)
Starting at position icol, define string as line(istart:istop).
integer(i4b) function, public getunit()
Get a free unit number.
character(len=max(len_trim(str), width)) function, public str_pad_left(str, width)
Function for string manipulation.
subroutine, public ulasav(buf, text, kstp, kper, pertim, totim, ncol, nrow, ilay, ichn)
Save 1 layer array on disk.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module defines variable data types.
Definition: kind.f90:8
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b) ifailedstepretry
current retry for this time step
subroutine, public table_cr(this, name, title)
Definition: Table.f90:87
real(dp), pointer, public pertim
time relative to start of stress period
Definition: tdis.f90:30
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21
subroutine, public read_value_or_time_series_adv(textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, varName)
Call this subroutine from advanced packages to define timeseries link for a variable (varName).
subroutine apt_ot_bdsummary(this, kstp, kper, iout, ibudfl)
Print advanced package transport dependent variables.
Definition: tsp-apt.f90:1041
subroutine apt_cfupdate(this)
Advanced package transport routine.
Definition: tsp-apt.f90:891
subroutine pak_setup_budobj(this, idx)
Set up a budget object that stores an advanced package flows.
Definition: tsp-apt.f90:2086
subroutine apt_ar(this)
Advanced package transport allocate and read (ar) routine.
Definition: tsp-apt.f90:301
subroutine pak_fill_budobj(this, idx, x, flowja, ccratin, ccratout)
Copy flow terms into thisbudobj, must be overridden.
Definition: tsp-apt.f90:2238
subroutine rp_obs_flowjaface(this, obsrv, budterm)
Prepare observation.
Definition: tsp-apt.f90:2570
subroutine pak_set_stressperiod(this, itemno, keyword, found)
Advanced package transport set stress period routine.
Definition: tsp-apt.f90:590
subroutine pak_bd_obs(this, obstypeid, jj, v, found)
Check if observation exists in an advanced package.
Definition: tsp-apt.f90:2809
subroutine apt_fc_expanded(this, rhs, ia, idxglo, matrix_sln)
Advanced package transport fill coefficient (fc) method.
Definition: tsp-apt.f90:762
subroutine pak_fc_expanded(this, rhs, ia, idxglo, matrix_sln)
Advanced package transport fill coefficient (fc) method.
Definition: tsp-apt.f90:871
subroutine apt_fjf_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Go through each "within apt-apt" connection (e.g., lkt-lkt, or sft-sft) and accumulate total mass (or...
Definition: tsp-apt.f90:2337
subroutine rp_obs_budterm(this, obsrv, budterm)
Prepare observation.
Definition: tsp-apt.f90:2498
subroutine allocate_scalars(this)
@ brief Allocate scalars
Definition: tsp-apt.f90:1058
subroutine apt_read_initial_attr(this)
Read the initial parameters for an advanced package.
Definition: tsp-apt.f90:1628
subroutine apt_setup_budobj(this)
Set up the budget object that stores advanced package flow terms.
Definition: tsp-apt.f90:1898
subroutine apt_allocate_index_arrays(this)
@ brief Allocate index arrays
Definition: tsp-apt.f90:1114
subroutine apt_rp_obs(this)
Read and prepare apt-related observations.
Definition: tsp-apt.f90:2643
character(len=lenvarname) text
Definition: tsp-apt.f90:63
subroutine apt_ot_package_flows(this, icbcfl, ibudfl)
Save advanced package flows routine.
Definition: tsp-apt.f90:961
real(dp) function, dimension(:), pointer, contiguous get_mvr_depvar(this)
Advanced package transport utility function.
Definition: tsp-apt.f90:630
subroutine apt_allocate_arrays(this)
@ brief Allocate arrays
Definition: tsp-apt.f90:1171
integer(i4b) function pak_get_nbudterms(this)
Function to return the number of budget terms just for this package.
Definition: tsp-apt.f90:1882
subroutine apt_set_stressperiod(this, itemno)
Advanced package transport set stress period routine.
Definition: tsp-apt.f90:493
subroutine apt_df_obs(this)
Define observation type.
Definition: tsp-apt.f90:2408
character(len=lenftype) ftype
Definition: tsp-apt.f90:62
subroutine apt_fmvr_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Account for mass or energy transferred to this package from the MVR package.
Definition: tsp-apt.f90:2314
subroutine apt_da(this)
@ brief Deallocate memory
Definition: tsp-apt.f90:1226
subroutine apt_read_dimensions(this)
Determine dimensions for this advanced package.
Definition: tsp-apt.f90:1405
subroutine apt_mc(this, moffset, matrix_sln)
Advanced package transport map package connections to matrix.
Definition: tsp-apt.f90:240
subroutine apt_fill_budobj(this, x, flowja)
Copy flow terms into thisbudobj.
Definition: tsp-apt.f90:2100
subroutine pak_rp_obs(this, obsrv, found)
Process package specific obs.
Definition: tsp-apt.f90:2438
logical function apt_obs_supported(this)
Determine whether an obs type is supported.
Definition: tsp-apt.f90:2393
subroutine apt_read_cvs(this)
Read feature information for this advanced package.
Definition: tsp-apt.f90:1471
subroutine apt_bd_obs(this)
Calculate observation values.
Definition: tsp-apt.f90:2721
subroutine apt_solve(this)
Add terms specific to advanced package transport to the explicit solve.
Definition: tsp-apt.f90:1677
subroutine apt_tmvr_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Account for mass or energy transferred to the MVR package.
Definition: tsp-apt.f90:2286
subroutine apt_get_volumes(this, icv, vnew, vold, delt)
Return the feature new volume and old volume.
Definition: tsp-apt.f90:1858
subroutine apt_set_pointers(this, neq, ibound, xnew, xold, flowja)
Set pointers to model arrays and variables so that a package has access to these items.
Definition: tsp-apt.f90:1832
subroutine apt_accumulate_ccterm(this, ilak, rrate, ccratin, ccratout)
Accumulate constant concentration (or temperature) terms for budget.
Definition: tsp-apt.f90:1776
subroutine apt_setup_tableobj(this)
Setup a table object an advanced package.
Definition: tsp-apt.f90:2921
subroutine find_apt_package(this)
Find corresponding advanced package transport package.
Definition: tsp-apt.f90:1298
subroutine apt_rp(this)
Advanced package transport read and prepare (rp) routine.
Definition: tsp-apt.f90:368
subroutine apt_stor_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Account for mass or energy storage in advanced package features.
Definition: tsp-apt.f90:2258
subroutine apt_options(this, option, found)
Set options specific to the TspAptType.
Definition: tsp-apt.f90:1314
subroutine rp_obs_byfeature(this, obsrv)
Prepare observation.
Definition: tsp-apt.f90:2455
subroutine apt_copy2flowp(this)
Copy concentrations (or temperatures) into flow package aux variable.
Definition: tsp-apt.f90:2367
subroutine apt_cq(this, x, flowja, iadv)
Advanced package transport calculate flows (cq) routine.
Definition: tsp-apt.f90:921
subroutine apt_ac(this, moffset, sparse)
Add package connection to matrix.
Definition: tsp-apt.f90:193
subroutine, public apt_process_obsid(obsrv, dis, inunitobs, iout)
Process observation IDs for an advanced package.
Definition: tsp-apt.f90:2828
subroutine apt_ad(this)
Advanced package transport routine.
Definition: tsp-apt.f90:643
integer(i4b) function apt_check_valid(this, itemno)
Advanced package transport routine.
Definition: tsp-apt.f90:610
subroutine apt_ad_chk(this)
Definition: tsp-apt.f90:482
subroutine define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
Definition: tsp-apt.f90:1808
subroutine pak_solve(this)
Add terms specific to advanced package transport features to the explicit solve routine.
Definition: tsp-apt.f90:1764
subroutine apt_fc_nonexpanded(this, rhs, ia, idxglo, matrix_sln)
Advanced package transport fill coefficient (fc) method.
Definition: tsp-apt.f90:733
subroutine apt_reset(this)
Override bnd reset for custom mover logic.
Definition: tsp-apt.f90:700
subroutine apt_fc(this, rhs, ia, idxglo, matrix_sln)
Definition: tsp-apt.f90:710
subroutine apt_ot_dv(this, idvsave, idvprint)
Definition: tsp-apt.f90:985
subroutine pak_df_obs(this)
Define apt observation type.
Definition: tsp-apt.f90:2423
subroutine, public apt_process_obsid12(obsrv, dis, inunitobs, iout)
Process observation IDs for a package.
Definition: tsp-apt.f90:2871
@ brief BndType