MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
gwe-mwe.f90
Go to the documentation of this file.
1 ! -- Multi-Aquifer Well Energy Transport Module
2 ! -- todo: save the mwe temperature into the mwe aux variable?
3 ! -- todo: calculate the maw DENSE aux variable using temperature?
4 !
5 ! MAW flows (flowbudptr) index var MWE term Transport Type
6 !---------------------------------------------------------------------------------
7 
8 ! -- terms from MAW that will be handled by parent APT Package
9 ! FLOW-JA-FACE idxbudfjf FLOW-JA-FACE cv2cv (note that this doesn't exist for MAW)
10 ! GWF (aux FLOW-AREA) idxbudgwf GWF cv2gwf
11 ! STORAGE (aux VOLUME) idxbudsto none used for cv volumes
12 ! FROM-MVR idxbudfmvr FROM-MVR q * text = this%qfrommvr(:)
13 ! TO-MVR idxbudtmvr TO-MVR q * tfeat
14 
15 ! -- MAW terms
16 ! RATE idxbudrate RATE q < 0: q * twell, else q * tuser
17 ! FW-RATE idxbudfwrt FW-RATE q * twell
18 ! RATE-TO-MVR idxbudrtmv RATE-TO-MVR q * twell
19 ! FW-RATE-TO-MVR idxbudfrtm FW-RATE-TO-MVR q * twell
20 ! WELL-AQUIFER CONDUCTION idxbudmwcd MW-CONDUCTION K_t_f * WetArea / thickness
21 
22 ! -- terms from MAW that should be skipped
23 ! CONSTANT-TO-MVR ? CONSTANT-TO-MVR q * twell
24 
25 ! -- terms from a flow file that should be skipped
26 ! CONSTANT none none none
27 ! AUXILIARY none none none
28 
29 ! -- terms that are written to the energy transport budget file
30 ! none none STORAGE (aux ENER) dE/dt
31 ! none none AUXILIARY none
32 ! none none CONSTANT accumulate
33 !
34 !
36 
37  use kindmodule, only: dp, i4b
38  use constantsmodule, only: dzero, linelength
39  use simmodule, only: store_error
40  use bndmodule, only: bndtype, getbndfromlist
41  use tspfmimodule, only: tspfmitype
42  use mawmodule, only: mawtype
43  use observemodule, only: observetype
48 
49  implicit none
50 
51  public mwe_create
52 
53  character(len=*), parameter :: ftype = 'MWE'
54  character(len=*), parameter :: flowtype = 'MAW'
55  character(len=16) :: text = ' MWE'
56 
57  type, extends(tspapttype) :: gwemwetype
58 
59  type(gweinputdatatype), pointer :: gwecommon => null() !< pointer to shared gwe data used by multiple packages but set in mst
60 
61  integer(I4B), pointer :: idxbudrate => null() ! index of well rate terms in flowbudptr
62  integer(I4B), pointer :: idxbudfwrt => null() ! index of flowing well rate terms in flowbudptr
63  integer(I4B), pointer :: idxbudrtmv => null() ! index of rate to mover terms in flowbudptr
64  integer(I4B), pointer :: idxbudfrtm => null() ! index of flowing well rate to mover terms in flowbudptr
65  integer(I4B), pointer :: idxbudmwcd => null() ! index of well bore conduction terms in flowbudptr
66  real(dp), dimension(:), pointer, contiguous :: temprate => null() ! well rate temperature
67 
68  contains
69 
70  procedure :: bnd_da => mwe_da
71  procedure :: allocate_scalars
72  procedure :: apt_allocate_arrays => mwe_allocate_arrays
73  procedure :: find_apt_package => find_mwe_package
74  procedure :: pak_fc_expanded => mwe_fc_expanded
75  procedure :: pak_solve => mwe_solve
76  procedure :: pak_get_nbudterms => mwe_get_nbudterms
77  procedure :: pak_setup_budobj => mwe_setup_budobj
78  procedure :: pak_fill_budobj => mwe_fill_budobj
79  procedure :: mwe_rate_term
80  procedure :: mwe_fwrt_term
81  procedure :: mwe_rtmv_term
82  procedure :: mwe_frtm_term
83  procedure :: pak_df_obs => mwe_df_obs
84  procedure :: pak_rp_obs => mwe_rp_obs
85  procedure :: pak_bd_obs => mwe_bd_obs
86  procedure :: pak_set_stressperiod => mwe_set_stressperiod
87 
88  end type gwemwetype
89 
90 contains
91 
92  !> Create new MWE package
93  !<
94  subroutine mwe_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
95  fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
96  ! -- dummy
97  class(bndtype), pointer :: packobj
98  integer(I4B), intent(in) :: id
99  integer(I4B), intent(in) :: ibcnum
100  integer(I4B), intent(in) :: inunit
101  integer(I4B), intent(in) :: iout
102  character(len=*), intent(in) :: namemodel
103  character(len=*), intent(in) :: pakname
104  type(tspfmitype), pointer :: fmi
105  real(dp), intent(in), pointer :: eqnsclfac !< governing equation scale factor
106  type(gweinputdatatype), intent(in), target :: gwecommon !< shared data container for use by multiple GWE packages
107  character(len=*), intent(in) :: dvt !< For GWE, set to "TEMPERATURE" in TspAptType
108  character(len=*), intent(in) :: dvu !< For GWE, set to "energy" in TspAptType
109  character(len=*), intent(in) :: dvua !< For GWE, set to "E" in TspAptType
110  ! -- local
111  type(gwemwetype), pointer :: mweobj
112  !
113  ! -- Allocate the object and assign values to object variables
114  allocate (mweobj)
115  packobj => mweobj
116  !
117  ! -- Create name and memory path
118  call packobj%set_names(ibcnum, namemodel, pakname, ftype)
119  packobj%text = text
120  !
121  ! -- Allocate scalars
122  call mweobj%allocate_scalars()
123  !
124  ! -- Initialize package
125  call packobj%pack_initialize()
126  !
127  packobj%inunit = inunit
128  packobj%iout = iout
129  packobj%id = id
130  packobj%ibcnum = ibcnum
131  packobj%ncolbnd = 1
132  packobj%iscloc = 1
133  !
134  ! -- Store pointer to flow model interface. When the GwfGwe exchange is
135  ! created, it sets fmi%bndlist so that the GWE model has access to all
136  ! the flow packages
137  mweobj%fmi => fmi
138  !
139  ! -- Store pointer to governing equation scale factor
140  mweobj%eqnsclfac => eqnsclfac
141  !
142  ! -- Store pointer to shared data module for accessing cpw, rhow
143  ! for the budget calculations, and for accessing the latent heat of
144  ! vaporization for evaporative cooling.
145  mweobj%gwecommon => gwecommon
146  !
147  ! -- Set labels that will be used in generalized APT class
148  mweobj%depvartype = dvt
149  mweobj%depvarunit = dvu
150  mweobj%depvarunitabbrev = dvua
151  end subroutine mwe_create
152 
153  !> @brief Find corresponding mwe package
154  !<
155  subroutine find_mwe_package(this)
156  ! -- modules
158  ! -- dummy
159  class(gwemwetype) :: this
160  ! -- local
161  character(len=LINELENGTH) :: errmsg
162  class(bndtype), pointer :: packobj
163  integer(I4B) :: ip, icount
164  integer(I4B) :: nbudterm
165  logical :: found
166  !
167  ! -- Initialize found to false, and error later if flow package cannot
168  ! be found
169  found = .false.
170  !
171  ! -- If user is specifying flows in a binary budget file, then set up
172  ! the budget file reader, otherwise set a pointer to the flow package
173  ! budobj
174  if (this%fmi%flows_from_file) then
175  call this%fmi%set_aptbudobj_pointer(this%flowpackagename, this%flowbudptr)
176  if (associated(this%flowbudptr)) found = .true.
177  !
178  else
179  if (associated(this%fmi%gwfbndlist)) then
180  ! -- Look through gwfbndlist for a flow package with the same name as
181  ! this transport package name
182  do ip = 1, this%fmi%gwfbndlist%Count()
183  packobj => getbndfromlist(this%fmi%gwfbndlist, ip)
184  if (packobj%packName == this%flowpackagename) then
185  found = .true.
186  !
187  ! -- Store BndType pointer to packobj, and then
188  ! use the select type to point to the budobj in flow package
189  this%flowpackagebnd => packobj
190  select type (packobj)
191  type is (mawtype)
192  this%flowbudptr => packobj%budobj
193  end select
194  end if
195  if (found) exit
196  end do
197  end if
198  end if
199  !
200  ! -- Error if flow package not found
201  if (.not. found) then
202  write (errmsg, '(a)') 'Could not find flow package with name '&
203  &//trim(adjustl(this%flowpackagename))//'.'
204  call store_error(errmsg)
205  call this%parser%StoreErrorUnit()
206  end if
207  !
208  ! -- Allocate space for idxbudssm, which indicates whether this is a
209  ! special budget term or one that is a general source and sink
210  nbudterm = this%flowbudptr%nbudterm
211  call mem_allocate(this%idxbudssm, nbudterm, 'IDXBUDSSM', this%memoryPath)
212  !
213  ! -- Process budget terms and identify special budget terms
214  write (this%iout, '(/, a, a)') &
215  'PROCESSING '//ftype//' INFORMATION FOR ', this%packName
216  write (this%iout, '(a)') ' IDENTIFYING FLOW TERMS IN '//flowtype//' PACKAGE'
217  write (this%iout, '(a, i0)') &
218  ' NUMBER OF '//flowtype//' = ', this%flowbudptr%ncv
219  icount = 1
220  do ip = 1, this%flowbudptr%nbudterm
221  select case (trim(adjustl(this%flowbudptr%budterm(ip)%flowtype)))
222  case ('FLOW-JA-FACE')
223  this%idxbudfjf = ip
224  this%idxbudssm(ip) = 0
225  case ('GWF')
226  this%idxbudgwf = ip
227  this%idxbudssm(ip) = 0
228  case ('STORAGE')
229  this%idxbudsto = ip
230  this%idxbudssm(ip) = 0
231  case ('RATE')
232  this%idxbudrate = ip
233  this%idxbudssm(ip) = 0
234  case ('FW-RATE')
235  this%idxbudfwrt = ip
236  this%idxbudssm(ip) = 0
237  case ('RATE-TO-MVR')
238  this%idxbudrtmv = ip
239  this%idxbudssm(ip) = 0
240  case ('FW-RATE-TO-MVR')
241  this%idxbudfrtm = ip
242  this%idxbudssm(ip) = 0
243  case ('TO-MVR')
244  this%idxbudtmvr = ip
245  this%idxbudssm(ip) = 0
246  case ('FROM-MVR')
247  this%idxbudfmvr = ip
248  this%idxbudssm(ip) = 0
249  case ('AUXILIARY')
250  this%idxbudaux = ip
251  this%idxbudssm(ip) = 0
252  case default
253  !
254  ! -- Set idxbudssm equal to a column index for where the temperatures
255  ! are stored in the tempbud(nbudssm, ncv) array
256  this%idxbudssm(ip) = icount
257  icount = icount + 1
258  end select
259  write (this%iout, '(a, i0, " = ", a,/, a, i0)') &
260  ' TERM ', ip, trim(adjustl(this%flowbudptr%budterm(ip)%flowtype)), &
261  ' MAX NO. OF ENTRIES = ', this%flowbudptr%budterm(ip)%maxlist
262  end do
263  write (this%iout, '(a, //)') 'DONE PROCESSING '//ftype//' INFORMATION'
264  !
265  ! -- Streambed conduction term
266  this%idxbudmwcd = this%idxbudgwf
267  end subroutine find_mwe_package
268 
269  !> @brief Add matrix terms related to MWE
270  !!
271  !! This routine is called from TspAptType%apt_fc_expanded() in
272  !! order to add matrix terms specifically for MWE
273  !<
274  subroutine mwe_fc_expanded(this, rhs, ia, idxglo, matrix_sln)
275  ! -- dummy
276  class(gwemwetype) :: this
277  real(DP), dimension(:), intent(inout) :: rhs
278  integer(I4B), dimension(:), intent(in) :: ia
279  integer(I4B), dimension(:), intent(in) :: idxglo
280  class(matrixbasetype), pointer :: matrix_sln
281  ! -- local
282  integer(I4B) :: j, n, n1, n2
283  integer(I4B) :: iloc
284  integer(I4B) :: iposd, iposoffd
285  integer(I4B) :: ipossymd, ipossymoffd
286  integer(I4B) :: auxpos
287  real(DP) :: rrate
288  real(DP) :: rhsval
289  real(DP) :: hcofval
290  real(DP) :: ctherm ! kluge?
291  real(DP) :: wa !< wetted area
292  real(DP) :: ktf !< thermal conductivity of streambed material
293  real(DP) :: s !< thickness of conductive wellbore material
294  !
295  ! -- Add puping rate contribution
296  if (this%idxbudrate /= 0) then
297  do j = 1, this%flowbudptr%budterm(this%idxbudrate)%nlist
298  call this%mwe_rate_term(j, n1, n2, rrate, rhsval, hcofval)
299  iloc = this%idxlocnode(n1)
300  iposd = this%idxpakdiag(n1)
301  call matrix_sln%add_value_pos(iposd, hcofval)
302  rhs(iloc) = rhs(iloc) + rhsval
303  end do
304  end if
305  !
306  ! -- Add flowing well rate contribution
307  if (this%idxbudfwrt /= 0) then
308  do j = 1, this%flowbudptr%budterm(this%idxbudfwrt)%nlist
309  call this%mwe_fwrt_term(j, n1, n2, rrate, rhsval, hcofval)
310  iloc = this%idxlocnode(n1)
311  iposd = this%idxpakdiag(n1)
312  call matrix_sln%add_value_pos(iposd, hcofval)
313  rhs(iloc) = rhs(iloc) + rhsval
314  end do
315  end if
316  !
317  ! -- Add rate to mover contribution
318  if (this%idxbudrtmv /= 0) then
319  do j = 1, this%flowbudptr%budterm(this%idxbudrtmv)%nlist
320  call this%mwe_rtmv_term(j, n1, n2, rrate, rhsval, hcofval)
321  iloc = this%idxlocnode(n1)
322  iposd = this%idxpakdiag(n1)
323  call matrix_sln%add_value_pos(iposd, hcofval)
324  rhs(iloc) = rhs(iloc) + rhsval
325  end do
326  end if
327  !
328  ! -- Add puping rate contribution
329  if (this%idxbudfrtm /= 0) then
330  do j = 1, this%flowbudptr%budterm(this%idxbudfrtm)%nlist
331  call this%mwe_frtm_term(j, n1, n2, rrate, rhsval, hcofval)
332  iloc = this%idxlocnode(n1)
333  iposd = this%idxpakdiag(n1)
334  call matrix_sln%add_value_pos(iposd, hcofval)
335  rhs(iloc) = rhs(iloc) + rhsval
336  end do
337  end if
338  !
339  ! -- Add wellbore conduction contribution
340  do j = 1, this%flowbudptr%budterm(this%idxbudgwf)%nlist
341  !
342  ! -- Set n to feature number and process if active features
343  n = this%flowbudptr%budterm(this%idxbudgwf)%id1(j)
344  if (this%iboundpak(n) /= 0) then
345  !
346  ! -- Set acoef and rhs to negative so they are relative to mwe and not gwe
347  auxpos = this%flowbudptr%budterm(this%idxbudgwf)%naux
348  wa = this%flowbudptr%budterm(this%idxbudgwf)%auxvar(auxpos, j)
349  ktf = this%ktf(n)
350  s = this%rfeatthk(n)
351  ctherm = ktf * wa / s
352  !
353  ! -- Add to mwe row
354  iposd = this%idxdglo(j)
355  iposoffd = this%idxoffdglo(j)
356  call matrix_sln%add_value_pos(iposd, -ctherm) ! kluge note: make sure the signs on ctherm are correct here and below
357  call matrix_sln%add_value_pos(iposoffd, ctherm)
358  !
359  ! -- Add to gwe row for mwe connection
360  ipossymd = this%idxsymdglo(j)
361  ipossymoffd = this%idxsymoffdglo(j)
362  call matrix_sln%add_value_pos(ipossymd, -ctherm)
363  call matrix_sln%add_value_pos(ipossymoffd, ctherm)
364  end if
365  end do
366  end subroutine mwe_fc_expanded
367 
368  !> @brief Add terms specific to multi-aquifer wells to the explicit multi-
369  !! aquifer well energy transport solve
370  !<
371  subroutine mwe_solve(this)
372  ! -- dummy
373  class(gwemwetype) :: this
374  ! -- local
375  integer(I4B) :: j
376  integer(I4B) :: n1, n2
377  real(DP) :: rrate
378  !
379  ! -- Add well pumping contribution
380  if (this%idxbudrate /= 0) then
381  do j = 1, this%flowbudptr%budterm(this%idxbudrate)%nlist
382  call this%mwe_rate_term(j, n1, n2, rrate)
383  this%dbuff(n1) = this%dbuff(n1) + rrate
384  end do
385  end if
386  !
387  ! -- Add flowing well rate contribution
388  if (this%idxbudfwrt /= 0) then
389  do j = 1, this%flowbudptr%budterm(this%idxbudfwrt)%nlist
390  call this%mwe_fwrt_term(j, n1, n2, rrate)
391  this%dbuff(n1) = this%dbuff(n1) + rrate
392  end do
393  end if
394  !
395  ! -- Add well pumping rate to mover contribution
396  if (this%idxbudrtmv /= 0) then
397  do j = 1, this%flowbudptr%budterm(this%idxbudrtmv)%nlist
398  call this%mwe_rtmv_term(j, n1, n2, rrate)
399  this%dbuff(n1) = this%dbuff(n1) + rrate
400  end do
401  end if
402  !
403  ! -- Add flowing well rate to mover contribution
404  if (this%idxbudfrtm /= 0) then
405  do j = 1, this%flowbudptr%budterm(this%idxbudfrtm)%nlist
406  call this%mwe_frtm_term(j, n1, n2, rrate)
407  this%dbuff(n1) = this%dbuff(n1) + rrate
408  end do
409  end if
410  end subroutine mwe_solve
411 
412  !> @brief Function to return the number of budget terms just for this package
413  !!
414  !! This overrides a function in the parent class.
415  !<
416  function mwe_get_nbudterms(this) result(nbudterms)
417  ! -- dummy
418  class(gwemwetype) :: this
419  ! -- return
420  integer(I4B) :: nbudterms
421  !
422  ! -- Number of potential budget terms is 5
423  nbudterms = 1 ! RATE
424  if (this%idxbudfwrt /= 0) nbudterms = nbudterms + 1
425  if (this%idxbudrtmv /= 0) nbudterms = nbudterms + 1
426  if (this%idxbudfrtm /= 0) nbudterms = nbudterms + 1
427  if (this%idxbudmwcd /= 0) nbudterms = nbudterms + 1
428  end function mwe_get_nbudterms
429 
430  !> @brief Set up the budget object that stores all the mwe flows
431  !<
432  subroutine mwe_setup_budobj(this, idx)
433  ! -- modules
434  use constantsmodule, only: lenbudtxt
435  ! -- dummy
436  class(gwemwetype) :: this
437  integer(I4B), intent(inout) :: idx
438  ! -- local
439  integer(I4B) :: n, n1, n2
440  integer(I4B) :: maxlist, naux
441  real(DP) :: q
442  character(len=LENBUDTXT) :: text
443  !
444  ! -- User-specified rate
445  text = ' RATE'
446  idx = idx + 1
447  maxlist = this%flowbudptr%budterm(this%idxbudrate)%maxlist
448  naux = 0
449  call this%budobj%budterm(idx)%initialize(text, &
450  this%name_model, &
451  this%packName, &
452  this%name_model, &
453  this%packName, &
454  maxlist, .false., .false., &
455  naux)
456  !
457  ! -- Flowing well rate
458  if (this%idxbudfwrt /= 0) then
459  text = ' FW-RATE'
460  idx = idx + 1
461  maxlist = this%flowbudptr%budterm(this%idxbudfwrt)%maxlist
462  naux = 0
463  call this%budobj%budterm(idx)%initialize(text, &
464  this%name_model, &
465  this%packName, &
466  this%name_model, &
467  this%packName, &
468  maxlist, .false., .false., &
469  naux)
470  end if
471  !
472  ! -- User-specified flow rate to mover
473  if (this%idxbudrtmv /= 0) then
474  text = ' RATE-TO-MVR'
475  idx = idx + 1
476  maxlist = this%flowbudptr%budterm(this%idxbudrtmv)%maxlist
477  naux = 0
478  call this%budobj%budterm(idx)%initialize(text, &
479  this%name_model, &
480  this%packName, &
481  this%name_model, &
482  this%packName, &
483  maxlist, .false., .false., &
484  naux)
485  end if
486  !
487  ! -- Fflowing well rate to mover
488  if (this%idxbudfrtm /= 0) then
489  text = ' FW-RATE-TO-MVR'
490  idx = idx + 1
491  maxlist = this%flowbudptr%budterm(this%idxbudfrtm)%maxlist
492  naux = 0
493  call this%budobj%budterm(idx)%initialize(text, &
494  this%name_model, &
495  this%packName, &
496  this%name_model, &
497  this%packName, &
498  maxlist, .false., .false., &
499  naux)
500  end if
501  !
502  ! -- Conduction through wellbore (and/or filter pack)
503  text = ' WELLBORE-COND'
504  idx = idx + 1
505  maxlist = this%flowbudptr%budterm(this%idxbudmwcd)%maxlist
506  naux = 0
507  call this%budobj%budterm(idx)%initialize(text, &
508  this%name_model, &
509  this%packName, &
510  this%name_model, &
511  this%packName, &
512  maxlist, .false., .false., &
513  naux)
514  call this%budobj%budterm(idx)%reset(maxlist)
515  q = dzero
516  do n = 1, maxlist
517  n1 = this%flowbudptr%budterm(this%idxbudgwf)%id1(n)
518  n2 = this%flowbudptr%budterm(this%idxbudgwf)%id2(n)
519  call this%budobj%budterm(idx)%update_term(n1, n2, q)
520  end do
521  end subroutine mwe_setup_budobj
522 
523  !> @brief Copy flow terms into this%budobj
524  !<
525  subroutine mwe_fill_budobj(this, idx, x, flowja, ccratin, ccratout)
526  ! -- dummy
527  class(gwemwetype) :: this
528  integer(I4B), intent(inout) :: idx
529  real(DP), dimension(:), intent(in) :: x
530  real(DP), dimension(:), contiguous, intent(inout) :: flowja
531  real(DP), intent(inout) :: ccratin
532  real(DP), intent(inout) :: ccratout
533  ! -- local
534  integer(I4B) :: j, n1, n2
535  integer(I4B) :: nlist
536  integer(I4B) :: igwfnode
537  integer(I4B) :: idiag
538  integer(I4B) :: auxpos
539  real(DP) :: q
540  real(DP) :: ctherm
541  real(DP) :: wa !< wetted area
542  real(DP) :: ktf !< thermal conductivity of streambed material
543  real(DP) :: s !< thickness of conductive streambed materia
544  !
545  ! -- Rate
546  idx = idx + 1
547  nlist = this%flowbudptr%budterm(this%idxbudrate)%nlist
548  call this%budobj%budterm(idx)%reset(nlist)
549  do j = 1, nlist
550  call this%mwe_rate_term(j, n1, n2, q)
551  call this%budobj%budterm(idx)%update_term(n1, n2, q)
552  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
553  end do
554  !
555  ! -- FW-Rate
556  if (this%idxbudfwrt /= 0) then
557  idx = idx + 1
558  nlist = this%flowbudptr%budterm(this%idxbudfwrt)%nlist
559  call this%budobj%budterm(idx)%reset(nlist)
560  do j = 1, nlist
561  call this%mwe_fwrt_term(j, n1, n2, q)
562  call this%budobj%budterm(idx)%update_term(n1, n2, q)
563  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
564  end do
565  end if
566  !
567  ! -- Rate-To-MVR
568  if (this%idxbudrtmv /= 0) then
569  idx = idx + 1
570  nlist = this%flowbudptr%budterm(this%idxbudrtmv)%nlist
571  call this%budobj%budterm(idx)%reset(nlist)
572  do j = 1, nlist
573  call this%mwe_rtmv_term(j, n1, n2, q)
574  call this%budobj%budterm(idx)%update_term(n1, n2, q)
575  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
576  end do
577  end if
578  !
579  ! -- FW-Rate-To-MVR
580  if (this%idxbudfrtm /= 0) then
581  idx = idx + 1
582  nlist = this%flowbudptr%budterm(this%idxbudfrtm)%nlist
583  call this%budobj%budterm(idx)%reset(nlist)
584  do j = 1, nlist
585  call this%mwe_frtm_term(j, n1, n2, q)
586  call this%budobj%budterm(idx)%update_term(n1, n2, q)
587  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
588  end do
589  end if
590  !
591  ! -- Wellbore-Cond
592  idx = idx + 1
593  call this%budobj%budterm(idx)%reset(this%maxbound)
594  do j = 1, this%flowbudptr%budterm(this%idxbudmwcd)%nlist
595  q = dzero
596  n1 = this%flowbudptr%budterm(this%idxbudmwcd)%id1(j)
597  if (this%iboundpak(n1) /= 0) then
598  igwfnode = this%flowbudptr%budterm(this%idxbudmwcd)%id2(j)
599  auxpos = this%flowbudptr%budterm(this%idxbudgwf)%naux
600  wa = this%flowbudptr%budterm(this%idxbudgwf)%auxvar(auxpos, j)
601  ktf = this%ktf(n1)
602  s = this%rfeatthk(n1)
603  ctherm = ktf * wa / s
604  q = ctherm * (x(igwfnode) - this%xnewpak(n1))
605  end if
606  call this%budobj%budterm(idx)%update_term(n1, igwfnode, q)
607  call this%apt_accumulate_ccterm(n1, q, ccratin, ccratout)
608  if (this%iboundpak(n1) /= 0) then
609  ! -- Contribution to gwe cell budget
610  this%simvals(j) = this%simvals(j) - q
611  idiag = this%dis%con%ia(igwfnode)
612  flowja(idiag) = flowja(idiag) - q
613  end if
614  end do
615  end subroutine mwe_fill_budobj
616 
617  !> @brief Allocate scalars specific to the multi-aquifer well energy
618  !! transport (MWE) package.
619  !<
620  subroutine allocate_scalars(this)
621  ! -- modules
623  ! -- dummy
624  class(gwemwetype) :: this
625  !
626  ! -- Allocate scalars in TspAptType
627  call this%TspAptType%allocate_scalars()
628  !
629  ! -- Allocate
630  call mem_allocate(this%idxbudrate, 'IDXBUDRATE', this%memoryPath)
631  call mem_allocate(this%idxbudfwrt, 'IDXBUDFWRT', this%memoryPath)
632  call mem_allocate(this%idxbudrtmv, 'IDXBUDRTMV', this%memoryPath)
633  call mem_allocate(this%idxbudfrtm, 'IDXBUDFRTM', this%memoryPath)
634  call mem_allocate(this%idxbudmwcd, 'IDXBUDMWCD', this%memoryPath)
635  !
636  ! -- Initialize
637  this%idxbudrate = 0
638  this%idxbudfwrt = 0
639  this%idxbudrtmv = 0
640  this%idxbudfrtm = 0
641  this%idxbudmwcd = 0
642  end subroutine allocate_scalars
643 
644  !> @brief Allocate arrays specific to the streamflow mass transport (SFT)
645  !! package
646  !<
647  subroutine mwe_allocate_arrays(this)
648  ! -- modules
650  ! -- dummy
651  class(gwemwetype), intent(inout) :: this
652  ! -- local
653  integer(I4B) :: n
654  !
655  ! -- Time series
656  call mem_allocate(this%temprate, this%ncv, 'TEMPRATE', this%memoryPath)
657  !
658  ! -- Call standard TspAptType allocate arrays
659  call this%TspAptType%apt_allocate_arrays()
660  !
661  ! -- Initialize
662  do n = 1, this%ncv
663  this%temprate(n) = dzero
664  end do
665  end subroutine mwe_allocate_arrays
666 
667  !> @brief Deallocate memory associated with MWE package
668  !<
669  subroutine mwe_da(this)
670  ! -- modules
672  ! -- dummy
673  class(gwemwetype) :: this
674  !
675  ! -- Deallocate scalars
676  call mem_deallocate(this%idxbudrate)
677  call mem_deallocate(this%idxbudfwrt)
678  call mem_deallocate(this%idxbudrtmv)
679  call mem_deallocate(this%idxbudfrtm)
680  call mem_deallocate(this%idxbudmwcd)
681  !
682  ! -- Deallocate time series
683  call mem_deallocate(this%temprate)
684  !
685  ! -- Deallocate scalars in TspAptType
686  call this%TspAptType%bnd_da()
687  end subroutine mwe_da
688 
689  !> @brief Thermal transport matrix term(s) associated with a user-specified
690  !! flow rate (mwe_rate_term)
691  !<
692  subroutine mwe_rate_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
693  ! -- dummy
694  class(gwemwetype) :: this
695  integer(I4B), intent(in) :: ientry
696  integer(I4B), intent(inout) :: n1
697  integer(I4B), intent(inout) :: n2
698  real(DP), intent(inout), optional :: rrate
699  real(DP), intent(inout), optional :: rhsval
700  real(DP), intent(inout), optional :: hcofval
701  ! -- local
702  real(DP) :: qbnd
703  real(DP) :: ctmp
704  real(DP) :: h, r
705  !
706  n1 = this%flowbudptr%budterm(this%idxbudrate)%id1(ientry)
707  n2 = this%flowbudptr%budterm(this%idxbudrate)%id2(ientry)
708  ! -- Note that qbnd is negative for an extracting well
709  qbnd = this%flowbudptr%budterm(this%idxbudrate)%flow(ientry)
710  if (qbnd < dzero) then
711  ctmp = this%xnewpak(n1)
712  h = qbnd
713  r = dzero
714  else
715  ctmp = this%temprate(n1)
716  h = dzero
717  r = -qbnd * ctmp
718  end if
719  if (present(rrate)) rrate = qbnd * ctmp * this%eqnsclfac
720  if (present(rhsval)) rhsval = r * this%eqnsclfac
721  if (present(hcofval)) hcofval = h * this%eqnsclfac
722  end subroutine mwe_rate_term
723 
724  !> @brief Thermal transport matrix term(s) associated with a flowing-
725  !! well rate term associated with pumping (or injection)
726  !<
727  subroutine mwe_fwrt_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
728  ! -- dummy
729  class(gwemwetype) :: this
730  integer(I4B), intent(in) :: ientry
731  integer(I4B), intent(inout) :: n1
732  integer(I4B), intent(inout) :: n2
733  real(DP), intent(inout), optional :: rrate
734  real(DP), intent(inout), optional :: rhsval
735  real(DP), intent(inout), optional :: hcofval
736  ! -- local
737  real(DP) :: qbnd
738  real(DP) :: ctmp
739  !
740  n1 = this%flowbudptr%budterm(this%idxbudfwrt)%id1(ientry)
741  n2 = this%flowbudptr%budterm(this%idxbudfwrt)%id2(ientry)
742  qbnd = this%flowbudptr%budterm(this%idxbudfwrt)%flow(ientry)
743  ctmp = this%xnewpak(n1)
744  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
745  if (present(rhsval)) rhsval = dzero
746  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac
747  end subroutine mwe_fwrt_term
748 
749  !> @brief Thermal transport matrix term(s) associated with pumped-water-
750  !! to-mover term (mwe_rtmv_term)
751  !<
752  subroutine mwe_rtmv_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
753  ! -- dummy
754  class(gwemwetype) :: this
755  integer(I4B), intent(in) :: ientry
756  integer(I4B), intent(inout) :: n1
757  integer(I4B), intent(inout) :: n2
758  real(DP), intent(inout), optional :: rrate
759  real(DP), intent(inout), optional :: rhsval
760  real(DP), intent(inout), optional :: hcofval
761  ! -- local
762  real(DP) :: qbnd
763  real(DP) :: ctmp
764  !
765  n1 = this%flowbudptr%budterm(this%idxbudrtmv)%id1(ientry)
766  n2 = this%flowbudptr%budterm(this%idxbudrtmv)%id2(ientry)
767  qbnd = this%flowbudptr%budterm(this%idxbudrtmv)%flow(ientry)
768  ctmp = this%xnewpak(n1)
769  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
770  if (present(rhsval)) rhsval = dzero
771  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac
772  end subroutine mwe_rtmv_term
773 
774  !> @brief Thermal transport matrix term(s) associated with the flowing-
775  !! well-rate-to-mover term (mwe_frtm_term)
776  !<
777  subroutine mwe_frtm_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
778  ! -- dummy
779  class(gwemwetype) :: this
780  integer(I4B), intent(in) :: ientry
781  integer(I4B), intent(inout) :: n1
782  integer(I4B), intent(inout) :: n2
783  real(DP), intent(inout), optional :: rrate
784  real(DP), intent(inout), optional :: rhsval
785  real(DP), intent(inout), optional :: hcofval
786  ! -- local
787  real(DP) :: qbnd
788  real(DP) :: ctmp
789  !
790  n1 = this%flowbudptr%budterm(this%idxbudfrtm)%id1(ientry)
791  n2 = this%flowbudptr%budterm(this%idxbudfrtm)%id2(ientry)
792  qbnd = this%flowbudptr%budterm(this%idxbudfrtm)%flow(ientry)
793  ctmp = this%xnewpak(n1)
794  if (present(rrate)) rrate = ctmp * qbnd * this%eqnsclfac
795  if (present(rhsval)) rhsval = dzero
796  if (present(hcofval)) hcofval = qbnd * this%eqnsclfac
797  end subroutine mwe_frtm_term
798 
799  !> @brief Observations
800  !!
801  !! Store the observation type supported by the APT package and override
802  !! BndType%bnd_df_obs
803  !<
804  subroutine mwe_df_obs(this)
805  ! -- dummy
806  class(gwemwetype) :: this
807  ! -- local
808  integer(I4B) :: indx
809  !
810  ! -- Store obs type and assign procedure pointer
811  ! for temperature observation type.
812  call this%obs%StoreObsType('temperature', .false., indx)
813  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
814  !
815  ! -- Flow-ja-face not supported for MWE
816  !call this%obs%StoreObsType('flow-ja-face', .true., indx)
817  !this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsID
818  !
819  ! -- Store obs type and assign procedure pointer
820  ! for from-mvr observation type.
821  call this%obs%StoreObsType('from-mvr', .true., indx)
822  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
823  !
824  ! -- To-mvr not supported for mwe
825  !call this%obs%StoreObsType('to-mvr', .true., indx)
826  !this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsID
827  !
828  ! -- Store obs type and assign procedure pointer
829  ! for storage observation type.
830  call this%obs%StoreObsType('storage', .true., indx)
831  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
832  !
833  ! -- Store obs type and assign procedure pointer
834  ! for constant observation type.
835  call this%obs%StoreObsType('constant', .true., indx)
836  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
837  !
838  ! -- Store obs type and assign procedure pointer
839  ! for observation type: mwe
840  call this%obs%StoreObsType('mwe', .true., indx)
841  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid12
842  !
843  ! -- Store obs type and assign procedure pointer
844  ! for rate observation type.
845  call this%obs%StoreObsType('rate', .true., indx)
846  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
847  !
848  ! -- Store obs type and assign procedure pointer
849  ! for observation type.
850  call this%obs%StoreObsType('fw-rate', .true., indx)
851  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
852  !
853  ! -- Store obs type and assign procedure pointer
854  ! for observation type.
855  call this%obs%StoreObsType('rate-to-mvr', .true., indx)
856  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
857  !
858  ! -- Store obs type and assign procedure pointer
859  ! for observation type.
860  call this%obs%StoreObsType('fw-rate-to-mvr', .true., indx)
861  this%obs%obsData(indx)%ProcessIdPtr => apt_process_obsid
862  end subroutine mwe_df_obs
863 
864  !> @brief Process package specific obs
865  !!
866  !! Method to process specific observations for this package.
867  !<
868  subroutine mwe_rp_obs(this, obsrv, found)
869  ! -- dummy
870  class(gwemwetype), intent(inout) :: this !< package class
871  type(observetype), intent(inout) :: obsrv !< observation object
872  logical, intent(inout) :: found !< indicate whether observation was found
873  !
874  found = .true.
875  select case (obsrv%ObsTypeId)
876  case ('RATE')
877  call this%rp_obs_byfeature(obsrv)
878  case ('FW-RATE')
879  call this%rp_obs_byfeature(obsrv)
880  case ('RATE-TO-MVR')
881  call this%rp_obs_byfeature(obsrv)
882  case ('FW-RATE-TO-MVR')
883  call this%rp_obs_byfeature(obsrv)
884  case default
885  found = .false.
886  end select
887  end subroutine mwe_rp_obs
888 
889  !> @brief Calculate observation value and pass it back to APT
890  !<
891  subroutine mwe_bd_obs(this, obstypeid, jj, v, found)
892  ! -- dummy
893  class(gwemwetype), intent(inout) :: this
894  character(len=*), intent(in) :: obstypeid
895  real(DP), intent(inout) :: v
896  integer(I4B), intent(in) :: jj
897  logical, intent(inout) :: found
898  ! -- local
899  integer(I4B) :: n1, n2
900  !
901  found = .true.
902  select case (obstypeid)
903  case ('RATE')
904  if (this%iboundpak(jj) /= 0) then
905  call this%mwe_rate_term(jj, n1, n2, v)
906  end if
907  case ('FW-RATE')
908  if (this%iboundpak(jj) /= 0 .and. this%idxbudfwrt > 0) then
909  call this%mwe_fwrt_term(jj, n1, n2, v)
910  end if
911  case ('RATE-TO-MVR')
912  if (this%iboundpak(jj) /= 0 .and. this%idxbudrtmv > 0) then
913  call this%mwe_rtmv_term(jj, n1, n2, v)
914  end if
915  case ('FW-RATE-TO-MVR')
916  if (this%iboundpak(jj) /= 0 .and. this%idxbudfrtm > 0) then
917  call this%mwe_frtm_term(jj, n1, n2, v)
918  end if
919  case default
920  found = .false.
921  end select
922  end subroutine mwe_bd_obs
923 
924  !> @brief Sets the stress period attributes for keyword use.
925  !<
926  subroutine mwe_set_stressperiod(this, itemno, keyword, found)
927  ! -- modules
929  ! -- dummy
930  class(gwemwetype), intent(inout) :: this
931  integer(I4B), intent(in) :: itemno
932  character(len=*), intent(in) :: keyword
933  logical, intent(inout) :: found
934  ! -- local
935  character(len=LINELENGTH) :: text
936  integer(I4B) :: ierr
937  integer(I4B) :: jj
938  real(DP), pointer :: bndElem => null()
939  !
940  ! RATE <rate>
941  !
942  found = .true.
943  select case (keyword)
944  case ('RATE')
945  ierr = this%apt_check_valid(itemno)
946  if (ierr /= 0) then
947  goto 999
948  end if
949  call this%parser%GetString(text)
950  jj = 1
951  bndelem => this%temprate(itemno)
952  call read_value_or_time_series_adv(text, itemno, jj, bndelem, &
953  this%packName, 'BND', this%tsManager, &
954  this%iprpak, 'RATE')
955  case default
956  !
957  ! -- Keyword not recognized so return to caller with found = .false.
958  found = .false.
959  end select
960  !
961 999 continue
962  end subroutine mwe_set_stressperiod
963 
964 end module gwemwemodule
This module contains the base boundary package.
class(bndtype) function, pointer, public getbndfromlist(list, idx)
Get boundary from package list.
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 dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
subroutine allocate_scalars(this)
Allocate scalars specific to the multi-aquifer well energy transport (MWE) package.
Definition: gwe-mwe.f90:621
character(len= *), parameter flowtype
Definition: gwe-mwe.f90:54
subroutine mwe_rp_obs(this, obsrv, found)
Process package specific obs.
Definition: gwe-mwe.f90:869
subroutine mwe_rtmv_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Thermal transport matrix term(s) associated with pumped-water- to-mover term (mwe_rtmv_term)
Definition: gwe-mwe.f90:753
subroutine mwe_allocate_arrays(this)
Allocate arrays specific to the streamflow mass transport (SFT) package.
Definition: gwe-mwe.f90:648
subroutine, public mwe_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
Create new MWE package.
Definition: gwe-mwe.f90:96
subroutine mwe_bd_obs(this, obstypeid, jj, v, found)
Calculate observation value and pass it back to APT.
Definition: gwe-mwe.f90:892
subroutine find_mwe_package(this)
Find corresponding mwe package.
Definition: gwe-mwe.f90:156
subroutine mwe_fill_budobj(this, idx, x, flowja, ccratin, ccratout)
Copy flow terms into thisbudobj.
Definition: gwe-mwe.f90:526
character(len= *), parameter ftype
Definition: gwe-mwe.f90:53
character(len=16) text
Definition: gwe-mwe.f90:55
subroutine mwe_da(this)
Deallocate memory associated with MWE package.
Definition: gwe-mwe.f90:670
subroutine mwe_setup_budobj(this, idx)
Set up the budget object that stores all the mwe flows.
Definition: gwe-mwe.f90:433
integer(i4b) function mwe_get_nbudterms(this)
Function to return the number of budget terms just for this package.
Definition: gwe-mwe.f90:417
subroutine mwe_frtm_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Thermal transport matrix term(s) associated with the flowing- well-rate-to-mover term (mwe_frtm_term)
Definition: gwe-mwe.f90:778
subroutine mwe_solve(this)
Add terms specific to multi-aquifer wells to the explicit multi- aquifer well energy transport solve.
Definition: gwe-mwe.f90:372
subroutine mwe_fc_expanded(this, rhs, ia, idxglo, matrix_sln)
Add matrix terms related to MWE.
Definition: gwe-mwe.f90:275
subroutine mwe_set_stressperiod(this, itemno, keyword, found)
Sets the stress period attributes for keyword use.
Definition: gwe-mwe.f90:927
subroutine mwe_fwrt_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Thermal transport matrix term(s) associated with a flowing- well rate term associated with pumping (o...
Definition: gwe-mwe.f90:728
subroutine mwe_rate_term(this, ientry, n1, n2, rrate, rhsval, hcofval)
Thermal transport matrix term(s) associated with a user-specified flow rate (mwe_rate_term)
Definition: gwe-mwe.f90:693
subroutine mwe_df_obs(this)
Observations.
Definition: gwe-mwe.f90:805
This module defines variable data types.
Definition: kind.f90:8
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
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, public apt_process_obsid(obsrv, dis, inunitobs, iout)
Process observation IDs for an advanced package.
Definition: tsp-apt.f90:2851
subroutine, public apt_process_obsid12(obsrv, dis, inunitobs, iout)
Process observation IDs for a package.
Definition: tsp-apt.f90:2894
@ brief BndType
Data for sharing among multiple packages. Originally read in from.