MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
gwf-buy.f90
Go to the documentation of this file.
1 ! Buoyancy Package for representing variable-density groundwater flow
2 ! The BUY Package does not work yet with the NPF XT3D option
3 
5 
6  use kindmodule, only: dp, i4b
13  use basedismodule, only: disbasetype
14  use gwfnpfmodule, only: gwfnpftype
17 
18  implicit none
19 
20  private
21  public :: gwfbuytype
22  public :: buy_cr
23 
25  real(dp), dimension(:), pointer :: conc => null() !< pointer to concentration array
26  integer(I4B), dimension(:), pointer :: icbund => null() !< store pointer to gwt ibound array
27  end type concentrationpointer
28 
29  type, extends(numericalpackagetype) :: gwfbuytype
30  type(gwfnpftype), pointer :: npf => null() !< npf object
31  integer(I4B), pointer :: ioutdense => null() !< unit number for saving density
32  integer(I4B), pointer :: iform => null() !< formulation: 0 freshwater head, 1 hh rhs, 2 hydraulic head
33  integer(I4B), pointer :: ireadelev => null() !< if 1 then elev has been allocated and filled
34  integer(I4B), pointer :: ireadconcbuy => null() !< if 1 then dense has been read from this buy input file
35  integer(I4B), pointer :: iconcset => null() !< if 1 then conc is pointed to a gwt model%x
36  real(dp), pointer :: denseref => null() !< reference fluid density
37  real(dp), dimension(:), pointer, contiguous :: dense => null() !< density
38  real(dp), dimension(:), pointer, contiguous :: concbuy => null() !< concentration array if specified in buy package
39  real(dp), dimension(:), pointer, contiguous :: elev => null() !< cell center elevation (optional; if not specified, then use (top+bot)/2)
40  integer(I4B), dimension(:), pointer :: ibound => null() !< store pointer to ibound
41 
42  integer(I4B), pointer :: nrhospecies => null() !< number of species used in equation of state to calculate density
43  real(dp), dimension(:), pointer, contiguous :: drhodc => null() !< change in density with change in concentration
44  real(dp), dimension(:), pointer, contiguous :: crhoref => null() !< reference concentration used in equation of state
45  real(dp), dimension(:), pointer, contiguous :: ctemp => null() !< temporary array of size (nrhospec) to pass to calcdens
46  character(len=LENMODELNAME), dimension(:), allocatable :: cmodelname !< names of gwt models used in equation of state
47  character(len=LENAUXNAME), dimension(:), allocatable :: cauxspeciesname !< names of gwt models used in equation of state
48 
49  type(concentrationpointer), allocatable, dimension(:) :: modelconc !< concentration pointer for each transport model
50 
51  contains
52  procedure :: buy_df
53  procedure :: buy_ar
54  procedure :: buy_ar_bnd
55  procedure :: buy_rp
56  procedure :: buy_ad
57  procedure :: buy_cf
58  procedure :: buy_cf_bnd
59  procedure :: buy_fc
60  procedure :: buy_ot_dv
61  procedure :: buy_cq
62  procedure :: buy_da
63  procedure, private :: calcbuy
64  procedure, private :: calchhterms
65  procedure, private :: buy_calcdens
66  procedure, private :: buy_calcelev
67  procedure :: allocate_scalars
68  procedure, private :: allocate_arrays
69  procedure, private :: source_options
70  procedure, private :: set_options
71  procedure, private :: log_options
72  procedure, private :: source_dimensions
73  procedure, private :: source_packagedata
74  procedure, private :: set_packagedata
76  end type gwfbuytype
77 
78 contains
79 
80  !> @brief Generic function to calculate fluid density from concentration
81  !<
82  function calcdens(denseref, drhodc, crhoref, conc) result(dense)
83  ! -- dummy
84  real(dp), intent(in) :: denseref
85  real(dp), dimension(:), intent(in) :: drhodc
86  real(dp), dimension(:), intent(in) :: crhoref
87  real(dp), dimension(:), intent(in) :: conc
88  ! -- return
89  real(dp) :: dense
90  ! -- local
91  integer(I4B) :: nrhospec
92  integer(I4B) :: i
93  !
94  nrhospec = size(drhodc)
95  dense = denseref
96  do i = 1, nrhospec
97  dense = dense + drhodc(i) * (conc(i) - crhoref(i))
98  end do
99  end function calcdens
100 
101  !> @brief Create a new BUY object
102  !<
103  subroutine buy_cr(buyobj, name_model, input_mempath, inunit, iout)
104  ! -- dummy
105  type(gwfbuytype), pointer :: buyobj
106  character(len=*), intent(in) :: name_model
107  character(len=*), intent(in) :: input_mempath
108  integer(I4B), intent(in) :: inunit
109  integer(I4B), intent(in) :: iout
110  !
111  ! -- Create the object
112  allocate (buyobj)
113  !
114  ! -- create name and memory path
115  call buyobj%set_names(1, name_model, 'BUY', 'BUY', input_mempath)
116  !
117  ! -- Allocate scalars
118  call buyobj%allocate_scalars()
119  !
120  ! -- Set variables
121  buyobj%inunit = inunit
122  buyobj%iout = iout
123  end subroutine buy_cr
124 
125  !> @brief Read options and package data, or set from argument
126  !<
127  subroutine buy_df(this, dis, buy_input)
128  ! -- dummy
129  class(gwfbuytype) :: this !< this buoyancy package
130  class(disbasetype), pointer, intent(in) :: dis !< pointer to discretization
131  type(gwfbuyinputdatatype), optional, intent(in) :: buy_input !< optional buy input data, otherwise read from file
132  ! -- local
133  ! -- formats
134  character(len=*), parameter :: fmtbuy = &
135  "(1x,/1x,'BUY -- Buoyancy Package, Version 1, 5/16/2018', &
136  &' input read from mempath: ', a, /)"
137  !
138  ! --print a message identifying the buoyancy package.
139  if (.not. present(buy_input)) then
140  write (this%iout, fmtbuy) this%input_mempath
141  end if
142  !
143  ! -- store pointers to arguments that were passed in
144  this%dis => dis
145 
146  if (.not. present(buy_input)) then
147  !
148  ! -- Read buoyancy options
149  call this%source_options()
150  !
151  ! -- Read buoyancy dimensions
152  call this%source_dimensions()
153  else
154  ! set from input data instead
155  call this%set_options(buy_input)
156  this%nrhospecies = buy_input%nrhospecies
157  end if
158  !
159  ! -- Allocate arrays
160  call this%allocate_arrays(dis%nodes)
161 
162  if (.not. present(buy_input)) then
163  !
164  ! -- Read buoyancy packagedata
165  call this%source_packagedata()
166  else
167  ! set from input data instead
168  call this%set_packagedata(buy_input)
169  end if
170  end subroutine buy_df
171 
172  !> @brief Allocate and Read
173  !<
174  subroutine buy_ar(this, npf, ibound)
175  ! -- dummy
176  class(gwfbuytype) :: this
177  type(gwfnpftype), pointer, intent(in) :: npf
178  integer(I4B), dimension(:), pointer :: ibound
179  !
180  ! -- store pointers to arguments that were passed in
181  this%npf => npf
182  this%ibound => ibound
183  !
184  ! -- Ensure NPF XT3D is not on
185  if (this%npf%ixt3d /= 0) then
186  call store_error('Error in model '//trim(this%name_model)// &
187  '. The XT3D option cannot be used with the BUY Package.')
188  call store_error_filename(this%input_fname)
189  end if
190  !
191  ! -- Calculate cell elevations
192  call this%buy_calcelev()
193  end subroutine buy_ar
194 
195  !> @brief Buoyancy ar_bnd routine to activate density in packages
196  !!
197  !! This routine is called from gwf_ar() as it goes through each package
198  !<
199  subroutine buy_ar_bnd(this, packobj, hnew)
200  ! -- modules
201  use bndmodule, only: bndtype
202  use lakmodule, only: laktype
203  use sfrmodule, only: sfrtype
204  use mawmodule, only: mawtype
205  ! -- dummy
206  class(gwfbuytype) :: this
207  class(bndtype), pointer :: packobj
208  real(DP), intent(in), dimension(:) :: hnew
209  !
210  ! -- Add density terms based on boundary package type
211  select case (packobj%filtyp)
212  case ('LAK')
213  !
214  ! -- activate density for lake package
215  select type (packobj)
216  type is (laktype)
217  call packobj%lak_activate_density()
218  end select
219  !
220  case ('SFR')
221  !
222  ! -- activate density for sfr package
223  select type (packobj)
224  type is (sfrtype)
225  call packobj%sfr_activate_density()
226  end select
227  !
228  case ('MAW')
229  !
230  ! -- activate density for maw package
231  select type (packobj)
232  type is (mawtype)
233  call packobj%maw_activate_density()
234  end select
235  !
236  case default
237  !
238  ! -- nothing
239  end select
240  end subroutine buy_ar_bnd
241 
242  !> @brief Check for new buy period data
243  !<
244  subroutine buy_rp(this)
245  ! -- modules
246  use tdismodule, only: kstp, kper
247  ! -- dummy
248  class(gwfbuytype) :: this
249  ! -- local
250  character(len=LINELENGTH) :: errmsg
251  integer(I4B) :: i
252  ! -- formats
253  character(len=*), parameter :: fmtc = &
254  "('Buoyancy Package does not have a concentration set &
255  &for species ',i0,'. One or more model names may be specified &
256  &incorrectly in the PACKAGEDATA block or a gwf-gwt exchange may need &
257  &to be activated.')"
258  !
259  ! -- Check to make sure all concentration pointers have been set
260  if (kstp * kper == 1) then
261  do i = 1, this%nrhospecies
262  if (.not. associated(this%modelconc(i)%conc)) then
263  write (errmsg, fmtc) i
264  call store_error(errmsg)
265  end if
266  end do
267  if (count_errors() > 0) then
268  call store_error_filename(this%input_fname)
269  end if
270  end if
271  end subroutine buy_rp
272 
273  !> @brief Advance
274  !<
275  subroutine buy_ad(this)
276  ! -- dummy
277  class(gwfbuytype) :: this
278  !
279  ! -- update density using the last concentration
280  call this%buy_calcdens()
281  end subroutine buy_ad
282 
283  !> @brief Fill coefficients
284  !<
285  subroutine buy_cf(this, kiter)
286  ! -- dummy
287  class(gwfbuytype) :: this
288  integer(I4B) :: kiter
289  !
290  ! -- Recalculate the elev array for this iteration
291  if (this%ireadelev == 0) then
292  if (this%iform == 1 .or. this%iform == 2) then
293  call this%buy_calcelev()
294  end if
295  end if
296  end subroutine buy_cf
297 
298  !> @brief Fill coefficients
299  !<
300  subroutine buy_cf_bnd(this, packobj, hnew) !, hcof, rhs, auxnam, auxvar)
301  ! -- modules
302  use bndmodule, only: bndtype
303  ! -- dummy
304  class(gwfbuytype) :: this
305  class(bndtype), pointer :: packobj
306  real(DP), intent(in), dimension(:) :: hnew
307  ! -- local
308  integer(I4B) :: i, j
309  integer(I4B) :: n, locdense, locelev
310  integer(I4B), dimension(:), allocatable :: locconc
311  !
312  ! -- Return if freshwater head formulation; all boundary heads must be
313  ! entered as freshwater equivalents
314  if (this%iform == 0) return
315  !
316  ! -- initialize
317  locdense = 0
318  locelev = 0
319  allocate (locconc(this%nrhospecies))
320  locconc(:) = 0
321  !
322  ! -- Add buoyancy terms for head-dependent boundaries
323  do n = 1, packobj%naux
324  if (packobj%auxname(n) == 'DENSITY') then
325  locdense = n
326  else if (packobj%auxname(n) == 'ELEVATION') then
327  locelev = n
328  end if
329  end do
330  !
331  ! -- find aux columns for concentrations that affect density
332  do i = 1, this%nrhospecies
333  locconc(i) = 0
334  do j = 1, packobj%naux
335  if (this%cauxspeciesname(i) == packobj%auxname(j)) then
336  locconc(i) = j
337  exit
338  end if
339  end do
340  if (locconc(i) == 0) then
341  ! -- one not found, so don't use and mark all as 0
342  locconc(:) = 0
343  exit
344  end if
345  end do
346  !
347  ! -- Add density terms based on boundary package type
348  select case (packobj%filtyp)
349  case ('GHB')
350  !
351  ! -- general head boundary
352  call buy_cf_ghb(packobj, hnew, this%dense, this%elev, this%denseref, &
353  locelev, locdense, locconc, this%drhodc, this%crhoref, &
354  this%ctemp, this%iform)
355  case ('RIV')
356  !
357  ! -- river
358  call buy_cf_riv(packobj, hnew, this%dense, this%elev, this%denseref, &
359  locelev, locdense, locconc, this%drhodc, this%crhoref, &
360  this%ctemp, this%iform)
361  case ('DRN')
362  !
363  ! -- drain
364  call buy_cf_drn(packobj, hnew, this%dense, this%denseref)
365  case ('LAK')
366  !
367  ! -- lake
368  call buy_cf_lak(packobj, hnew, this%dense, this%elev, this%denseref, &
369  locdense, locconc, this%drhodc, this%crhoref, &
370  this%ctemp, this%iform)
371  case ('SFR')
372  !
373  ! -- sfr
374  call buy_cf_sfr(packobj, hnew, this%dense, this%elev, this%denseref, &
375  locdense, locconc, this%drhodc, this%crhoref, &
376  this%ctemp, this%iform)
377  case ('MAW')
378  !
379  ! -- maw
380  call buy_cf_maw(packobj, hnew, this%dense, this%elev, this%denseref, &
381  locdense, locconc, this%drhodc, this%crhoref, &
382  this%ctemp, this%iform)
383  case default
384  !
385  ! -- nothing
386  end select
387  !
388  ! -- deallocate
389  deallocate (locconc)
390  end subroutine buy_cf_bnd
391 
392  !> @brief Return the density of the boundary package using one of several
393  !! different options in the following order of priority:
394  !! 1. Assign as aux variable in column with name 'DENSITY'
395  !! 2. Calculate using equation of state and nrhospecies aux columns
396  !! 3. If neither of those, then assign as denseref
397  !<
398  function get_bnd_density(n, locdense, locconc, denseref, drhodc, crhoref, &
399  ctemp, auxvar) result(densebnd)
400  ! -- dummy
401  integer(I4B), intent(in) :: n
402  integer(I4B), intent(in) :: locdense
403  integer(I4B), dimension(:), intent(in) :: locconc
404  real(dp), intent(in) :: denseref
405  real(dp), dimension(:), intent(in) :: drhodc
406  real(dp), dimension(:), intent(in) :: crhoref
407  real(dp), dimension(:), intent(inout) :: ctemp
408  real(dp), dimension(:, :), intent(in) :: auxvar
409  ! -- return
410  real(dp) :: densebnd
411  ! -- local
412  integer(I4B) :: i
413  !
414  ! -- assign boundary density based on one of three options
415  if (locdense > 0) then
416  ! -- assign density to an aux column named 'DENSITY'
417  densebnd = auxvar(locdense, n)
418  else if (locconc(1) > 0) then
419  ! -- calculate density using one or more concentration auxcolumns
420  do i = 1, size(locconc)
421  ctemp(i) = dzero
422  if (locconc(i) > 0) then
423  ctemp(i) = auxvar(locconc(i), n)
424  end if
425  end do
426  densebnd = calcdens(denseref, drhodc, crhoref, ctemp)
427  else
428  ! -- neither of the above, so assign as denseref
429  densebnd = denseref
430  end if
431  end function get_bnd_density
432 
433  !> @brief Fill ghb coefficients
434  !<
435  subroutine buy_cf_ghb(packobj, hnew, dense, elev, denseref, locelev, &
436  locdense, locconc, drhodc, crhoref, ctemp, &
437  iform)
438  ! -- modules
439  use bndmodule, only: bndtype
440  use ghbmodule, only: ghbtype
441  class(bndtype), pointer :: packobj
442  ! -- dummy
443  real(DP), intent(in), dimension(:) :: hnew
444  real(DP), intent(in), dimension(:) :: dense
445  real(DP), intent(in), dimension(:) :: elev
446  real(DP), intent(in) :: denseref
447  integer(I4B), intent(in) :: locelev
448  integer(I4B), intent(in) :: locdense
449  integer(I4B), dimension(:), intent(in) :: locconc
450  real(DP), dimension(:), intent(in) :: drhodc
451  real(DP), dimension(:), intent(in) :: crhoref
452  real(DP), dimension(:), intent(inout) :: ctemp
453  integer(I4B), intent(in) :: iform
454  ! -- local
455  integer(I4B) :: n
456  integer(I4B) :: node
457  real(DP) :: denseghb
458  real(DP) :: elevghb
459  real(DP) :: hghb
460  real(DP) :: cond
461  real(DP) :: hcofterm, rhsterm
462  !
463  ! -- Process density terms for each GHB
464  select type (packobj)
465  type is (ghbtype)
466  do n = 1, packobj%nbound
467  node = packobj%nodelist(n)
468  if (packobj%ibound(node) <= 0) cycle
469  !
470  ! -- density
471  denseghb = get_bnd_density(n, locdense, locconc, denseref, &
472  drhodc, crhoref, ctemp, packobj%auxvar)
473  !
474  ! -- elevation
475  elevghb = elev(node)
476  if (locelev > 0) elevghb = packobj%auxvar(locelev, n)
477  !
478  ! -- boundary head and conductance
479  hghb = packobj%bhead(n)
480  cond = packobj%cond(n)
481  !
482  ! -- calculate HCOF and RHS terms
483  call calc_ghb_hcof_rhs_terms(denseref, denseghb, dense(node), &
484  elevghb, elev(node), hghb, hnew(node), &
485  cond, iform, rhsterm, hcofterm)
486  packobj%hcof(n) = packobj%hcof(n) + hcofterm
487  packobj%rhs(n) = packobj%rhs(n) - rhsterm
488  !
489  end do
490  end select
491  end subroutine buy_cf_ghb
492 
493  !> @brief Calculate density hcof and rhs terms for ghb conditions
494  !<
495  subroutine calc_ghb_hcof_rhs_terms(denseref, denseghb, densenode, &
496  elevghb, elevnode, hghb, hnode, &
497  cond, iform, rhsterm, hcofterm)
498  ! -- dummy
499  real(DP), intent(in) :: denseref
500  real(DP), intent(in) :: denseghb
501  real(DP), intent(in) :: densenode
502  real(DP), intent(in) :: elevghb
503  real(DP), intent(in) :: elevnode
504  real(DP), intent(in) :: hghb
505  real(DP), intent(in) :: hnode
506  real(DP), intent(in) :: cond
507  integer(I4B), intent(in) :: iform
508  real(DP), intent(inout) :: rhsterm
509  real(DP), intent(inout) :: hcofterm
510  ! -- local
511  real(DP) :: t1, t2
512  real(DP) :: avgdense, avgelev
513  !
514  ! -- Calculate common terms
515  avgdense = dhalf * denseghb + dhalf * densenode
516  avgelev = dhalf * elevghb + dhalf * elevnode
517  t1 = avgdense / denseref - done
518  t2 = (denseghb - densenode) / denseref
519  !
520  ! -- Add hcof terms
521  hcofterm = -cond * t1
522  if (iform == 2) then
523  !
524  ! -- this term goes on RHS for iform == 1
525  hcofterm = hcofterm + dhalf * cond * t2
526  end if
527  !
528  ! -- Add rhs terms
529  rhsterm = cond * t1 * hghb
530  rhsterm = rhsterm - cond * t2 * avgelev
531  rhsterm = rhsterm + dhalf * cond * t2 * hghb
532  if (iform == 1) then
533  !
534  ! -- this term goes on LHS for iform == 2
535  rhsterm = rhsterm + dhalf * cond * t2 * hnode
536  end if
537  end subroutine calc_ghb_hcof_rhs_terms
538 
539  !> @brief Fill riv coefficients
540  !<
541  subroutine buy_cf_riv(packobj, hnew, dense, elev, denseref, locelev, &
542  locdense, locconc, drhodc, crhoref, ctemp, &
543  iform)
544  ! -- modules
545  use bndmodule, only: bndtype
546  use rivmodule, only: rivtype
547  class(bndtype), pointer :: packobj
548  ! -- dummy
549  real(DP), intent(in), dimension(:) :: hnew
550  real(DP), intent(in), dimension(:) :: dense
551  real(DP), intent(in), dimension(:) :: elev
552  real(DP), intent(in) :: denseref
553  integer(I4B), intent(in) :: locelev
554  integer(I4B), intent(in) :: locdense
555  integer(I4B), dimension(:), intent(in) :: locconc
556  real(DP), dimension(:), intent(in) :: drhodc
557  real(DP), dimension(:), intent(in) :: crhoref
558  real(DP), dimension(:), intent(inout) :: ctemp
559  integer(I4B), intent(in) :: iform
560  ! -- local
561  integer(I4B) :: n
562  integer(I4B) :: node
563  real(DP) :: denseriv
564  real(DP) :: elevriv
565  real(DP) :: hriv
566  real(DP) :: rbot
567  real(DP) :: cond
568  real(DP) :: hcofterm
569  real(DP) :: rhsterm
570  !
571  ! -- Process density terms for each RIV
572  select type (packobj)
573  type is (rivtype)
574  do n = 1, packobj%nbound
575  node = packobj%nodelist(n)
576  if (packobj%ibound(node) <= 0) cycle
577  !
578  ! -- density
579  denseriv = get_bnd_density(n, locdense, locconc, denseref, &
580  drhodc, crhoref, ctemp, packobj%auxvar)
581  !
582  ! -- elevation
583  elevriv = elev(node)
584  if (locelev > 0) elevriv = packobj%auxvar(locelev, n)
585  !
586  ! -- boundary head and conductance
587  hriv = packobj%stage(n)
588  cond = packobj%cond(n)
589  rbot = packobj%rbot(n)
590  !
591  ! -- calculate and add terms depending on whether head is above rbot
592  if (hnew(node) > rbot) then
593  !
594  ! --calculate HCOF and RHS terms, similar to GHB in this case
595  call calc_ghb_hcof_rhs_terms(denseref, denseriv, dense(node), &
596  elevriv, elev(node), hriv, hnew(node), &
597  cond, iform, rhsterm, hcofterm)
598  else
599  hcofterm = dzero
600  rhsterm = cond * (denseriv / denseref - done) * (hriv - rbot)
601  end if
602  !
603  ! -- Add terms to package hcof and rhs accumulators
604  packobj%hcof(n) = packobj%hcof(n) + hcofterm
605  packobj%rhs(n) = packobj%rhs(n) - rhsterm
606  end do
607  end select
608  end subroutine buy_cf_riv
609 
610  !> @brief Fill drn coefficients
611  !<
612  subroutine buy_cf_drn(packobj, hnew, dense, denseref)
613  ! -- modules
614  use bndmodule, only: bndtype
615  use drnmodule, only: drntype
616  class(bndtype), pointer :: packobj
617  ! -- dummy
618  real(DP), intent(in), dimension(:) :: hnew
619  real(DP), intent(in), dimension(:) :: dense
620  real(DP), intent(in) :: denseref
621  ! -- local
622  integer(I4B) :: n
623  integer(I4B) :: node
624  real(DP) :: rho
625  real(DP) :: hbnd
626  real(DP) :: cond
627  real(DP) :: hcofterm
628  real(DP) :: rhsterm
629  !
630  ! -- Process density terms for each DRN
631  select type (packobj)
632  type is (drntype)
633  do n = 1, packobj%nbound
634  node = packobj%nodelist(n)
635  if (packobj%ibound(node) <= 0) cycle
636  rho = dense(node)
637  hbnd = packobj%elev(n)
638  cond = packobj%cond(n)
639  if (hnew(node) > hbnd) then
640  hcofterm = -cond * (rho / denseref - done)
641  rhsterm = hcofterm * hbnd
642  packobj%hcof(n) = packobj%hcof(n) + hcofterm
643  packobj%rhs(n) = packobj%rhs(n) + rhsterm
644  end if
645  end do
646  end select
647  end subroutine buy_cf_drn
648 
649  !> @brief Pass density information into lak package; density terms are
650  !! calculated in the lake package as part of lak_calculate_density_exchange
651  !! method
652  !<
653  subroutine buy_cf_lak(packobj, hnew, dense, elev, denseref, locdense, &
654  locconc, drhodc, crhoref, ctemp, iform)
655  ! -- modules
656  use bndmodule, only: bndtype
657  use lakmodule, only: laktype
658  class(bndtype), pointer :: packobj
659  ! -- dummy
660  real(DP), intent(in), dimension(:) :: hnew
661  real(DP), intent(in), dimension(:) :: dense
662  real(DP), intent(in), dimension(:) :: elev
663  real(DP), intent(in) :: denseref
664  integer(I4B), intent(in) :: locdense
665  integer(I4B), dimension(:), intent(in) :: locconc
666  real(DP), dimension(:), intent(in) :: drhodc
667  real(DP), dimension(:), intent(in) :: crhoref
668  real(DP), dimension(:), intent(inout) :: ctemp
669  integer(I4B), intent(in) :: iform
670  ! -- local
671  integer(I4B) :: n
672  integer(I4B) :: node
673  real(DP) :: denselak
674  !
675  ! -- Insert the lake and gwf relative densities into col 1 and 2 and the
676  ! gwf elevation into col 3 of the lake package denseterms array
677  select type (packobj)
678  type is (laktype)
679  do n = 1, packobj%nbound
680  !
681  ! -- get gwf node number
682  node = packobj%nodelist(n)
683  if (packobj%ibound(node) <= 0) cycle
684  !
685  ! -- Determine lak density
686  denselak = get_bnd_density(n, locdense, locconc, denseref, &
687  drhodc, crhoref, ctemp, packobj%auxvar)
688  !
689  ! -- fill lak relative density into column 1 of denseterms
690  packobj%denseterms(1, n) = denselak / denseref
691  !
692  ! -- fill gwf relative density into column 2 of denseterms
693  packobj%denseterms(2, n) = dense(node) / denseref
694  !
695  ! -- fill gwf elevation into column 3 of denseterms
696  packobj%denseterms(3, n) = elev(node)
697  !
698  end do
699  end select
700  end subroutine buy_cf_lak
701 
702  !> @brief Pass density information into sfr package; density terms are
703  !! calculated in the sfr package as part of sfr_calculate_density_exchange
704  !! method
705  !<
706  subroutine buy_cf_sfr(packobj, hnew, dense, elev, denseref, locdense, &
707  locconc, drhodc, crhoref, ctemp, iform)
708  ! -- modules
709  use bndmodule, only: bndtype
710  use sfrmodule, only: sfrtype
711  class(bndtype), pointer :: packobj
712  ! -- dummy
713  real(DP), intent(in), dimension(:) :: hnew
714  real(DP), intent(in), dimension(:) :: dense
715  real(DP), intent(in), dimension(:) :: elev
716  real(DP), intent(in) :: denseref
717  integer(I4B), intent(in) :: locdense
718  integer(I4B), dimension(:), intent(in) :: locconc
719  real(DP), dimension(:), intent(in) :: drhodc
720  real(DP), dimension(:), intent(in) :: crhoref
721  real(DP), dimension(:), intent(inout) :: ctemp
722  integer(I4B), intent(in) :: iform
723  ! -- local
724  integer(I4B) :: n
725  integer(I4B) :: node
726  real(DP) :: densesfr
727  !
728  ! -- Insert the sfr and gwf relative densities into col 1 and 2 and the
729  ! gwf elevation into col 3 of the sfr package denseterms array
730  select type (packobj)
731  type is (sfrtype)
732  do n = 1, packobj%nbound
733  !
734  ! -- get gwf node number
735  node = packobj%nodelist(n)
736  if (packobj%ibound(node) <= 0) cycle
737  !
738  ! -- Determine sfr density
739  densesfr = get_bnd_density(n, locdense, locconc, denseref, &
740  drhodc, crhoref, ctemp, packobj%auxvar)
741  !
742  ! -- fill sfr relative density into column 1 of denseterms
743  packobj%denseterms(1, n) = densesfr / denseref
744  !
745  ! -- fill gwf relative density into column 2 of denseterms
746  packobj%denseterms(2, n) = dense(node) / denseref
747  !
748  ! -- fill gwf elevation into column 3 of denseterms
749  packobj%denseterms(3, n) = elev(node)
750  !
751  end do
752  end select
753  end subroutine buy_cf_sfr
754 
755  !> @brief Pass density information into maw package; density terms are
756  !! calculated in the maw package as part of maw_calculate_density_exchange
757  !! method
758  !<
759  subroutine buy_cf_maw(packobj, hnew, dense, elev, denseref, locdense, &
760  locconc, drhodc, crhoref, ctemp, iform)
761  ! -- modules
762  use bndmodule, only: bndtype
763  use mawmodule, only: mawtype
764  class(bndtype), pointer :: packobj
765  ! -- dummy
766  real(DP), intent(in), dimension(:) :: hnew
767  real(DP), intent(in), dimension(:) :: dense
768  real(DP), intent(in), dimension(:) :: elev
769  real(DP), intent(in) :: denseref
770  integer(I4B), intent(in) :: locdense
771  integer(I4B), dimension(:), intent(in) :: locconc
772  real(DP), dimension(:), intent(in) :: drhodc
773  real(DP), dimension(:), intent(in) :: crhoref
774  real(DP), dimension(:), intent(inout) :: ctemp
775  integer(I4B), intent(in) :: iform
776  ! -- local
777  integer(I4B) :: n
778  integer(I4B) :: node
779  real(DP) :: densemaw
780  !
781  ! -- Insert the maw and gwf relative densities into col 1 and 2 and the
782  ! gwf elevation into col 3 of the maw package denseterms array
783  select type (packobj)
784  type is (mawtype)
785  do n = 1, packobj%nbound
786  !
787  ! -- get gwf node number
788  node = packobj%nodelist(n)
789  if (packobj%ibound(node) <= 0) cycle
790  !
791  ! -- Determine maw density
792  densemaw = get_bnd_density(n, locdense, locconc, denseref, &
793  drhodc, crhoref, ctemp, packobj%auxvar)
794  !
795  ! -- fill maw relative density into column 1 of denseterms
796  packobj%denseterms(1, n) = densemaw / denseref
797  !
798  ! -- fill gwf relative density into column 2 of denseterms
799  packobj%denseterms(2, n) = dense(node) / denseref
800  !
801  ! -- fill gwf elevation into column 3 of denseterms
802  packobj%denseterms(3, n) = elev(node)
803  !
804  end do
805  end select
806  end subroutine buy_cf_maw
807 
808  !> @brief Fill coefficients
809  !<
810  subroutine buy_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
811  ! -- dummy
812  class(gwfbuytype) :: this
813  integer(I4B) :: kiter
814  class(matrixbasetype), pointer :: matrix_sln
815  integer(I4B), intent(in), dimension(:) :: idxglo
816  real(DP), dimension(:), intent(inout) :: rhs
817  real(DP), intent(inout), dimension(:) :: hnew
818  ! -- local
819  integer(I4B) :: n, m, ipos, idiag
820  real(DP) :: rhsterm, amatnn, amatnm
821  !
822  ! -- initialize
823  amatnn = dzero
824  amatnm = dzero
825  !
826  ! -- fill buoyancy flow term
827  do n = 1, this%dis%nodes
828  if (this%ibound(n) == 0) cycle
829  idiag = this%dis%con%ia(n)
830  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
831  m = this%dis%con%ja(ipos)
832  if (this%ibound(m) == 0) cycle
833  if (this%iform == 0) then
834  call this%calcbuy(n, m, ipos, hnew(n), hnew(m), rhsterm)
835  else
836  call this%calchhterms(n, m, ipos, hnew(n), hnew(m), rhsterm, &
837  amatnn, amatnm)
838  end if
839  !
840  ! -- Add terms to rhs, diagonal, and off diagonal
841  rhs(n) = rhs(n) - rhsterm
842  call matrix_sln%add_value_pos(idxglo(idiag), -amatnn)
843  call matrix_sln%add_value_pos(idxglo(ipos), amatnm)
844  end do
845  end do
846  end subroutine buy_fc
847 
848  !> @brief Save density array to binary file
849  !<
850  subroutine buy_ot_dv(this, idvfl)
851  ! -- dummy
852  class(gwfbuytype) :: this
853  integer(I4B), intent(in) :: idvfl
854  ! -- local
855  character(len=1) :: cdatafmp = ' ', editdesc = ' '
856  integer(I4B) :: ibinun
857  integer(I4B) :: iprint
858  integer(I4B) :: nvaluesp
859  integer(I4B) :: nwidthp
860  real(DP) :: dinact
861  !
862  ! -- Set unit number for density output
863  if (this%ioutdense /= 0) then
864  ibinun = 1
865  else
866  ibinun = 0
867  end if
868  if (idvfl == 0) ibinun = 0
869  !
870  ! -- save density array
871  if (ibinun /= 0) then
872  iprint = 0
873  dinact = dhnoflo
874  !
875  ! -- write density to binary file
876  if (this%ioutdense /= 0) then
877  ibinun = this%ioutdense
878  call this%dis%record_array(this%dense, this%iout, iprint, ibinun, &
879  ' DENSITY', cdatafmp, nvaluesp, &
880  nwidthp, editdesc, dinact)
881  end if
882  end if
883  end subroutine buy_ot_dv
884 
885  !> @brief Add buy term to flowja
886  !<
887  subroutine buy_cq(this, hnew, flowja)
888  implicit none
889  class(gwfbuytype) :: this
890  real(DP), intent(in), dimension(:) :: hnew
891  real(DP), intent(inout), dimension(:) :: flowja
892  integer(I4B) :: n, m, ipos
893  real(DP) :: deltaQ
894  real(DP) :: rhsterm, amatnn, amatnm
895  !
896  ! -- Calculate the flow across each cell face and store in flowja
897  do n = 1, this%dis%nodes
898  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
899  m = this%dis%con%ja(ipos)
900  if (m < n) cycle
901  if (this%iform == 0) then
902  ! -- equivalent freshwater head formulation
903  call this%calcbuy(n, m, ipos, hnew(n), hnew(m), deltaq)
904  else
905  ! -- hydraulic head formulation
906  call this%calchhterms(n, m, ipos, hnew(n), hnew(m), rhsterm, &
907  amatnn, amatnm)
908  deltaq = amatnm * hnew(m) - amatnn * hnew(n) + rhsterm
909  end if
910  flowja(ipos) = flowja(ipos) + deltaq
911  flowja(this%dis%con%isym(ipos)) = flowja(this%dis%con%isym(ipos)) - &
912  deltaq
913  end do
914  end do
915  end subroutine buy_cq
916 
917  !> @brief Deallocate
918  !<
919  subroutine buy_da(this)
920  ! -- dummy
921  class(gwfbuytype) :: this
922  !
923  ! -- Deallocate arrays if package was active
924  if (this%inunit > 0) then
925  call mem_deallocate(this%elev)
926  call mem_deallocate(this%dense)
927  call mem_deallocate(this%concbuy)
928  call mem_deallocate(this%drhodc)
929  call mem_deallocate(this%crhoref)
930  call mem_deallocate(this%ctemp)
931  deallocate (this%cmodelname)
932  deallocate (this%cauxspeciesname)
933  deallocate (this%modelconc)
934  end if
935  !
936  ! -- Scalars
937  call mem_deallocate(this%ioutdense)
938  call mem_deallocate(this%iform)
939  call mem_deallocate(this%ireadelev)
940  call mem_deallocate(this%ireadconcbuy)
941  call mem_deallocate(this%iconcset)
942  call mem_deallocate(this%denseref)
943  !
944  call mem_deallocate(this%nrhospecies)
945  !
946  ! -- deallocate parent
947  call this%NumericalPackageType%da()
948  end subroutine buy_da
949 
950  !> @ brief Source dimensions for package
951  !<
952  subroutine source_dimensions(this)
953  ! -- modules
956  ! -- dummy
957  class(gwfbuytype), intent(inout) :: this
958  ! -- local
959  type(gwfbuyparamfoundtype) :: found
960 
961  ! update defaults from input context
962  call mem_set_value(this%nrhospecies, 'NRHOSPECIES', this%input_mempath, &
963  found%nrhospecies)
964 
965  write (this%iout, '(/1x,a)') 'Processing BUY DIMENSIONS block'
966  write (this%iout, '(4x,a,i0)') 'NRHOSPECIES = ', this%nrhospecies
967  write (this%iout, '(1x,a)') 'End of BUY DIMENSIONS block'
968 
969  ! check dimension
970  if (this%nrhospecies < 1) then
971  call store_error('NRHOSPECIES must be greater than zero.')
972  call store_error_filename(this%input_fname)
973  end if
974  end subroutine source_dimensions
975 
976  !> @ brief source packagedata for package
977  !<
978  subroutine source_packagedata(this)
979  ! -- modules
983  ! -- dummy
984  class(gwfbuytype), intent(inout) :: this
985  ! -- local
986  integer(I4B), dimension(:), pointer, contiguous :: irhospec
987  type(characterstringtype), dimension(:), pointer, &
988  contiguous :: modelnames, auxspeciesnames
989  real(DP), dimension(:), pointer, contiguous :: drhodc, crhoref
990  integer(I4B), dimension(:), allocatable :: itemp
991  character(len=LINELENGTH) :: modelname, auxspeciesname, line, errmsg
992  character(len=10) :: c10
993  character(len=16) :: c16
994  integer(I4B) :: n
995  ! -- format
996  character(len=*), parameter :: fmterr = &
997  "('Invalid value for IRHOSPEC (',i0,') detected in BUY Package. &
998  &IRHOSPEC must be > 0 and <= NRHOSPECIES, and duplicate values &
999  &are not allowed.')"
1000 
1001  ! initialize
1002  allocate (itemp(this%nrhospecies))
1003  itemp(:) = 0
1004 
1005  ! set input context pointers
1006  call mem_setptr(irhospec, 'IRHOSPEC', this%input_mempath)
1007  call mem_setptr(drhodc, 'DRHODC', this%input_mempath)
1008  call mem_setptr(crhoref, 'CRHOREF', this%input_mempath)
1009  call mem_setptr(modelnames, 'MODELNAME', this%input_mempath)
1010  call mem_setptr(auxspeciesnames, 'AUXSPECIESNAME', this%input_mempath)
1011 
1012  ! process package data
1013  do n = 1, size(irhospec)
1014  modelname = modelnames(n)
1015  auxspeciesname = auxspeciesnames(n)
1016 
1017  if (irhospec(n) < 1 .or. irhospec(n) > this%nrhospecies) then
1018  write (errmsg, fmterr) irhospec(n)
1019  call store_error(errmsg)
1020  end if
1021  if (itemp(irhospec(n)) /= 0) then
1022  write (errmsg, fmterr) irhospec(n)
1023  call store_error(errmsg)
1024  end if
1025  itemp(irhospec(n)) = 1
1026 
1027  this%drhodc(irhospec(n)) = drhodc(n)
1028  this%crhoref(irhospec(n)) = crhoref(n)
1029  this%cmodelname(irhospec(n)) = trim(modelname)
1030  this%cauxspeciesname(irhospec(n)) = trim(auxspeciesname)
1031  end do
1032 
1033  ! Check for errors.
1034  if (count_errors() > 0) then
1035  call store_error_filename(this%input_fname)
1036  end if
1037 
1038  ! log package data
1039  write (this%iout, '(/,1x,a)') 'Processing BUY PACKAGEDATA block'
1040 
1041  ! write packagedata information
1042  write (this%iout, '(1x,a)') 'Summary of species information in BUY Package'
1043  write (this%iout, '(1a11, 4a17)') &
1044  'SPECIES', 'DRHODC', 'CRHOREF', 'MODEL', &
1045  'AUXSPECIESNAME'
1046  do n = 1, this%nrhospecies
1047  write (c10, '(i0)') n
1048  line = ' '//adjustr(c10)
1049  write (c16, '(g15.6)') this%drhodc(n)
1050  line = trim(line)//' '//adjustr(c16)
1051  write (c16, '(g15.6)') this%crhoref(n)
1052  line = trim(line)//' '//adjustr(c16)
1053  write (c16, '(a)') this%cmodelname(n)
1054  line = trim(line)//' '//adjustr(c16)
1055  write (c16, '(a)') this%cauxspeciesname(n)
1056  line = trim(line)//' '//adjustr(c16)
1057  write (this%iout, '(a)') trim(line)
1058  end do
1059 
1060  write (this%iout, '(1x,a)') 'End of BUY PACKAGEDATA block'
1061 
1062  ! cleanup
1063  deallocate (itemp)
1064  call memorystore_release('IRHOSPEC', this%input_mempath)
1065  call memorystore_release('DRHODC', this%input_mempath)
1066  call memorystore_release('CRHOREF', this%input_mempath)
1067  call memorystore_release('MODELNAME', this%input_mempath)
1068  call memorystore_release('AUXSPECIESNAME', this%input_mempath)
1069  end subroutine source_packagedata
1070 
1071  !> @brief Sets package data instead of reading from file
1072  !<
1073  subroutine set_packagedata(this, input_data)
1074  ! -- dummy
1075  class(gwfbuytype) :: this !< this buyoancy pkg
1076  type(gwfbuyinputdatatype), intent(in) :: input_data !< the input data to be set
1077  ! -- local
1078  integer(I4B) :: ispec
1079 
1080  do ispec = 1, this%nrhospecies
1081  this%drhodc(ispec) = input_data%drhodc(ispec)
1082  this%crhoref(ispec) = input_data%crhoref(ispec)
1083  this%cmodelname(ispec) = input_data%cmodelname(ispec)
1084  this%cauxspeciesname(ispec) = input_data%cauxspeciesname(ispec)
1085  end do
1086  end subroutine set_packagedata
1087 
1088  !> @brief Calculate buyancy term for this connection
1089  !<
1090  subroutine calcbuy(this, n, m, icon, hn, hm, buy)
1091  ! -- modules
1093  ! -- dummy
1094  class(gwfbuytype) :: this
1095  integer(I4B), intent(in) :: n
1096  integer(I4B), intent(in) :: m
1097  integer(I4B), intent(in) :: icon
1098  real(DP), intent(in) :: hn
1099  real(DP), intent(in) :: hm
1100  real(DP), intent(inout) :: buy
1101  ! -- local
1102  integer(I4B) :: ihc
1103  real(DP) :: densen, densem, cl1, cl2, avgdense, wt, elevn, elevm, &
1104  cond, tp, bt
1105  real(DP) :: hyn
1106  real(DP) :: hym
1107  !
1108  ! -- Average density
1109  densen = this%dense(n)
1110  densem = this%dense(m)
1111  if (m > n) then
1112  cl1 = this%dis%con%cl1(this%dis%con%jas(icon))
1113  cl2 = this%dis%con%cl2(this%dis%con%jas(icon))
1114  else
1115  cl1 = this%dis%con%cl2(this%dis%con%jas(icon))
1116  cl2 = this%dis%con%cl1(this%dis%con%jas(icon))
1117  end if
1118  wt = cl1 / (cl1 + cl2)
1119  avgdense = wt * densen + (done - wt) * densem
1120  !
1121  ! -- Elevations
1122  if (this%ireadelev == 0) then
1123  tp = this%dis%top(n)
1124  bt = this%dis%bot(n)
1125  elevn = bt + dhalf * this%npf%sat(n) * (tp - bt)
1126  tp = this%dis%top(m)
1127  bt = this%dis%bot(m)
1128  elevm = bt + dhalf * this%npf%sat(m) * (tp - bt)
1129  else
1130  elevn = this%elev(n)
1131  elevm = this%elev(m)
1132  end if
1133  !
1134  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
1135  hyn = this%npf%hy_eff(n, m, ihc, ipos=icon)
1136  hym = this%npf%hy_eff(m, n, ihc, ipos=icon)
1137  !
1138  ! -- Conductance
1139  if (this%dis%con%ihc(this%dis%con%jas(icon)) == 0) then
1140  cond = vcond(this%ibound(n), this%ibound(m), &
1141  this%npf%icelltype(n), this%npf%icelltype(m), &
1142  this%npf%inewton, &
1143  this%npf%ivarcv, this%npf%idewatcv, &
1144  this%npf%condsat(this%dis%con%jas(icon)), hn, hm, &
1145  hyn, hym, &
1146  this%npf%sat(n), this%npf%sat(m), &
1147  this%dis%top(n), this%dis%top(m), &
1148  this%dis%bot(n), this%dis%bot(m), &
1149  this%dis%con%hwva(this%dis%con%jas(icon)))
1150  else
1151  cond = hcond(this%ibound(n), this%ibound(m), &
1152  this%npf%icelltype(n), this%npf%icelltype(m), &
1153  this%npf%inewton, &
1154  this%dis%con%ihc(this%dis%con%jas(icon)), &
1155  this%npf%icellavg, &
1156  this%npf%condsat(this%dis%con%jas(icon)), &
1157  hn, hm, this%npf%sat(n), this%npf%sat(m), &
1158  hyn, hym, &
1159  this%dis%top(n), this%dis%top(m), &
1160  this%dis%bot(n), this%dis%bot(m), &
1161  this%dis%con%cl1(this%dis%con%jas(icon)), &
1162  this%dis%con%cl2(this%dis%con%jas(icon)), &
1163  this%dis%con%hwva(this%dis%con%jas(icon)))
1164  end if
1165  !
1166  ! -- Calculate buoyancy term
1167  buy = cond * (avgdense - this%denseref) / this%denseref * (elevm - elevn)
1168  end subroutine calcbuy
1169 
1170  !> @brief Calculate hydraulic head term for this connection
1171  !<
1172  subroutine calchhterms(this, n, m, icon, hn, hm, rhsterm, amatnn, amatnm)
1173  ! -- modules
1175  ! -- dummy
1176  class(gwfbuytype) :: this
1177  integer(I4B), intent(in) :: n
1178  integer(I4B), intent(in) :: m
1179  integer(I4B), intent(in) :: icon
1180  real(DP), intent(in) :: hn
1181  real(DP), intent(in) :: hm
1182  real(DP), intent(inout) :: rhsterm
1183  real(DP), intent(inout) :: amatnn
1184  real(DP), intent(inout) :: amatnm
1185  ! -- local
1186  integer(I4B) :: ihc
1187  real(DP) :: densen, densem, cl1, cl2, avgdense, wt, elevn, elevm, cond
1188  real(DP) :: rhonormn, rhonormm
1189  real(DP) :: rhoterm
1190  real(DP) :: elevnm
1191  real(DP) :: hphi
1192  real(DP) :: hyn
1193  real(DP) :: hym
1194  !
1195  ! -- Average density
1196  densen = this%dense(n)
1197  densem = this%dense(m)
1198  if (m > n) then
1199  cl1 = this%dis%con%cl1(this%dis%con%jas(icon))
1200  cl2 = this%dis%con%cl2(this%dis%con%jas(icon))
1201  else
1202  cl1 = this%dis%con%cl2(this%dis%con%jas(icon))
1203  cl2 = this%dis%con%cl1(this%dis%con%jas(icon))
1204  end if
1205  wt = cl1 / (cl1 + cl2)
1206  avgdense = wt * densen + (1.0 - wt) * densem
1207  !
1208  ! -- Elevations
1209  elevn = this%elev(n)
1210  elevm = this%elev(m)
1211  elevnm = (done - wt) * elevn + wt * elevm
1212  !
1213  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
1214  hyn = this%npf%hy_eff(n, m, ihc, ipos=icon)
1215  hym = this%npf%hy_eff(m, n, ihc, ipos=icon)
1216  !
1217  ! -- Conductance
1218  if (ihc == 0) then
1219  cond = vcond(this%ibound(n), this%ibound(m), &
1220  this%npf%icelltype(n), this%npf%icelltype(m), &
1221  this%npf%inewton, &
1222  this%npf%ivarcv, this%npf%idewatcv, &
1223  this%npf%condsat(this%dis%con%jas(icon)), hn, hm, &
1224  hyn, hym, &
1225  this%npf%sat(n), this%npf%sat(m), &
1226  this%dis%top(n), this%dis%top(m), &
1227  this%dis%bot(n), this%dis%bot(m), &
1228  this%dis%con%hwva(this%dis%con%jas(icon)))
1229  else
1230  cond = hcond(this%ibound(n), this%ibound(m), &
1231  this%npf%icelltype(n), this%npf%icelltype(m), &
1232  this%npf%inewton, &
1233  this%dis%con%ihc(this%dis%con%jas(icon)), &
1234  this%npf%icellavg, &
1235  this%npf%condsat(this%dis%con%jas(icon)), &
1236  hn, hm, this%npf%sat(n), this%npf%sat(m), &
1237  hyn, hym, &
1238  this%dis%top(n), this%dis%top(m), &
1239  this%dis%bot(n), this%dis%bot(m), &
1240  this%dis%con%cl1(this%dis%con%jas(icon)), &
1241  this%dis%con%cl2(this%dis%con%jas(icon)), &
1242  this%dis%con%hwva(this%dis%con%jas(icon)))
1243  end if
1244  !
1245  ! -- Calculate terms
1246  rhonormn = densen / this%denseref
1247  rhonormm = densem / this%denseref
1248  rhoterm = wt * rhonormn + (done - wt) * rhonormm
1249  amatnn = cond * (rhoterm - done)
1250  amatnm = amatnn
1251  rhsterm = -cond * (rhonormm - rhonormn) * elevnm
1252  if (this%iform == 1) then
1253  ! -- rhs (lag the h terms and keep matrix symmetric)
1254  hphi = (done - wt) * hn + wt * hm
1255  rhsterm = rhsterm + cond * hphi * (rhonormm - rhonormn)
1256  else
1257  ! -- lhs, results in asymmetric matrix due to weight term
1258  amatnn = amatnn - cond * (done - wt) * (rhonormm - rhonormn)
1259  amatnm = amatnm + cond * wt * (rhonormm - rhonormn)
1260  end if
1261  end subroutine calchhterms
1262 
1263  !> @brief calculate fluid density from concentration
1264  !<
1265  subroutine buy_calcdens(this)
1266  ! -- dummy
1267  class(gwfbuytype) :: this
1268 
1269  ! -- local
1270  integer(I4B) :: n
1271  integer(I4B) :: i
1272  !
1273  ! -- Calculate the density using the specified concentration array
1274  do n = 1, this%dis%nodes
1275  do i = 1, this%nrhospecies
1276  if (this%modelconc(i)%icbund(n) == 0) then
1277  this%ctemp(i) = dzero
1278  else
1279  this%ctemp(i) = this%modelconc(i)%conc(n)
1280  end if
1281  end do
1282  this%dense(n) = calcdens(this%denseref, this%drhodc, this%crhoref, &
1283  this%ctemp)
1284  end do
1285  end subroutine buy_calcdens
1286 
1287  !> @brief Calculate cell elevations to use in density flow equations
1288  !<
1289  subroutine buy_calcelev(this)
1290  ! -- dummy
1291  class(gwfbuytype) :: this
1292  ! -- local
1293  integer(I4B) :: n
1294  real(DP) :: tp, bt, frac
1295  !
1296  ! -- Calculate the elev array
1297  do n = 1, this%dis%nodes
1298  tp = this%dis%top(n)
1299  bt = this%dis%bot(n)
1300  frac = this%npf%sat(n)
1301  this%elev(n) = bt + dhalf * frac * (tp - bt)
1302  end do
1303  end subroutine buy_calcelev
1304 
1305  !> @brief Allocate scalars used by the package
1306  !<
1307  subroutine allocate_scalars(this)
1308  ! -- modules
1309  use constantsmodule, only: dzero
1310  ! -- dummy
1311  class(gwfbuytype) :: this
1312  ! -- local
1313  !
1314  ! -- allocate scalars in NumericalPackageType
1315  call this%NumericalPackageType%allocate_scalars()
1316  !
1317  ! -- Allocate
1318  call mem_allocate(this%ioutdense, 'IOUTDENSE', this%memoryPath)
1319  call mem_allocate(this%iform, 'IFORM', this%memoryPath)
1320  call mem_allocate(this%ireadelev, 'IREADELEV', this%memoryPath)
1321  call mem_allocate(this%ireadconcbuy, 'IREADCONCBUY', this%memoryPath)
1322  call mem_allocate(this%iconcset, 'ICONCSET', this%memoryPath)
1323  call mem_allocate(this%denseref, 'DENSEREF', this%memoryPath)
1324  !
1325  call mem_allocate(this%nrhospecies, 'NRHOSPECIES', this%memoryPath)
1326  !
1327  ! -- Initialize
1328  this%ioutdense = 0
1329  this%ireadelev = 0
1330  this%iconcset = 0
1331  this%ireadconcbuy = 0
1332  this%denseref = 1000.d0
1333  !
1334  this%nrhospecies = 0
1335  !
1336  ! -- Initialize default to LHS implementation of hydraulic head formulation
1337  this%iform = 2
1338  this%iasym = 1
1339  end subroutine allocate_scalars
1340 
1341  !> @brief Allocate arrays used by the package
1342  !<
1343  subroutine allocate_arrays(this, nodes)
1344  ! -- dummy
1345  class(gwfbuytype) :: this
1346  integer(I4B), intent(in) :: nodes
1347  ! -- local
1348  integer(I4B) :: i
1349  !
1350  ! -- Allocate
1351  call mem_allocate(this%dense, nodes, 'DENSE', this%memoryPath)
1352  call mem_allocate(this%concbuy, 0, 'CONCBUY', this%memoryPath)
1353  call mem_allocate(this%elev, nodes, 'ELEV', this%memoryPath)
1354  call mem_allocate(this%drhodc, this%nrhospecies, 'DRHODC', this%memoryPath)
1355  call mem_allocate(this%crhoref, this%nrhospecies, 'CRHOREF', this%memoryPath)
1356  call mem_allocate(this%ctemp, this%nrhospecies, 'CTEMP', this%memoryPath)
1357  allocate (this%cmodelname(this%nrhospecies))
1358  allocate (this%cauxspeciesname(this%nrhospecies))
1359  allocate (this%modelconc(this%nrhospecies))
1360  !
1361  ! -- Initialize
1362  do i = 1, nodes
1363  this%dense(i) = this%denseref
1364  this%elev(i) = dzero
1365  end do
1366  !
1367  ! -- Initialize nrhospecies arrays
1368  do i = 1, this%nrhospecies
1369  this%drhodc(i) = dzero
1370  this%crhoref(i) = dzero
1371  this%ctemp(i) = dzero
1372  this%cmodelname(i) = ''
1373  this%cauxspeciesname(i) = ''
1374  end do
1375  end subroutine allocate_arrays
1376 
1377  !> @ brief Source input options
1378  !<
1379  subroutine source_options(this)
1380  ! -- modules
1382  use openspecmodule, only: access, form
1383  use inputoutputmodule, only: getunit, openfile
1385  ! -- dummy
1386  class(gwfbuytype), intent(inout) :: this
1387  ! -- local
1388  character(len=LINELENGTH) :: densityfile
1389  type(gwfbuyparamfoundtype) :: found
1390 
1391  ! initialize variables
1392  densityfile = ''
1393 
1394  ! update defaults from input context
1395  call mem_set_value(this%iform, 'HHFORM_RHS', this%input_mempath, &
1396  found%hhform_rhs)
1397  call mem_set_value(this%denseref, 'DENSEREF', this%input_mempath, &
1398  found%denseref)
1399  call mem_set_value(this%iform, 'DEV_EFH_FORM', this%input_mempath, &
1400  found%dev_efh_form)
1401  call mem_set_value(densityfile, 'DENSITYFILE', this%input_mempath, &
1402  found%densityfile)
1403 
1404  ! update input dependent internal state
1405  if (found%hhform_rhs) this%iasym = 0
1406  if (found%dev_efh_form) then
1407  this%iform = 0
1408  this%iasym = 0
1409  end if
1410 
1411  ! fileout options
1412  if (found%densityfile) then
1413  this%ioutdense = getunit()
1414  call openfile(this%ioutdense, this%iout, densityfile, 'DATA(BINARY)', &
1415  form, access, 'REPLACE')
1416  end if
1417 
1418  ! log options
1419  call this%log_options(found, densityfile)
1420  end subroutine source_options
1421 
1422  !> @ brief log input options
1423  !<
1424  subroutine log_options(this, found, densityfile)
1425  ! -- modules
1427  ! -- dummy variables
1428  class(gwfbuytype), intent(inout) :: this
1429  type(gwfbuyparamfoundtype), intent(in) :: found
1430  character(len=*), intent(in) :: densityfile
1431  ! -- local variables
1432  ! -- formats
1433  character(len=*), parameter :: fmtfileout = &
1434  "(4x, 'BUY ', 1x, a, 1x, ' will be saved to file: ', &
1435  &a, ' opened on unit: ', I7)"
1436  !
1437  write (this%iout, '(1x,a)') 'Processing BUY OPTIONS block'
1438 
1439  if (found%hhform_rhs) then
1440  write (this%iout, '(4x,a)') &
1441  'Hydraulic head formulation set to right-hand side'
1442  end if
1443  if (found%denseref) then
1444  write (this%iout, '(4x,a,1pg15.6)') &
1445  'Reference density has been set to: ', this%denseref
1446  end if
1447  if (found%dev_efh_form) then
1448  write (this%iout, '(4x,a)') &
1449  'Formulation set to equivalent freshwater head'
1450  end if
1451  if (found%densityfile) then
1452  write (this%iout, fmtfileout) &
1453  'DENSITY', trim(densityfile), this%ioutdense
1454  end if
1455 
1456  write (this%iout, '(1x,a)') 'End of BUY OPTIONS block'
1457  end subroutine log_options
1458 
1459  !> @brief Sets options as opposed to reading them from a file
1460  !<
1461  subroutine set_options(this, input_data)
1462  ! -- dummy
1463  class(gwfbuytype) :: this
1464  type(gwfbuyinputdatatype), intent(in) :: input_data !< the input data to be set
1465  !
1466  this%iform = input_data%iform
1467  this%denseref = input_data%denseref
1468  !
1469  ! derived option:
1470  ! if not iform==2, there is no asymmetry
1471  if (this%iform == 0 .or. this%iform == 1) then
1472  this%iasym = 0
1473  end if
1474  end subroutine set_options
1475 
1476  !> @brief Pass in a gwt model name, concentration array and ibound, and store
1477  !! a pointer to these in the BUY package so that density can be calculated
1478  !! from them
1479  !!
1480  !! This routine is called from the gwfgwt exchange in the exg_ar() method
1481  !<
1482  subroutine set_concentration_pointer(this, modelname, conc, icbund)
1483  ! -- dummy
1484  class(gwfbuytype) :: this
1485  character(len=LENMODELNAME), intent(in) :: modelname
1486  real(DP), dimension(:), pointer :: conc
1487  integer(I4B), dimension(:), pointer :: icbund
1488  ! -- local
1489  integer(I4B) :: i
1490  logical :: found
1491  !
1492  this%iconcset = 1
1493  found = .false.
1494  do i = 1, this%nrhospecies
1495  if (this%cmodelname(i) == modelname) then
1496  this%modelconc(i)%conc => conc
1497  this%modelconc(i)%icbund => icbund
1498  found = .true.
1499  exit
1500  end if
1501  end do
1502  end subroutine set_concentration_pointer
1503 
1504 end module gwfbuymodule
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenmodelname
maximum length of the model name
Definition: Constants.f90:22
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
integer(i4b), parameter lenauxname
maximum length of a aux variable
Definition: Constants.f90:35
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
real(dp) function calcdens(denseref, drhodc, crhoref, conc)
Generic function to calculate fluid density from concentration.
Definition: gwf-buy.f90:83
subroutine buy_cf(this, kiter)
Fill coefficients.
Definition: gwf-buy.f90:286
subroutine buy_cf_bnd(this, packobj, hnew)
Fill coefficients.
Definition: gwf-buy.f90:301
subroutine buy_cf_lak(packobj, hnew, dense, elev, denseref, locdense, locconc, drhodc, crhoref, ctemp, iform)
Pass density information into lak package; density terms are calculated in the lake package as part o...
Definition: gwf-buy.f90:655
subroutine set_options(this, input_data)
Sets options as opposed to reading them from a file.
Definition: gwf-buy.f90:1462
subroutine calchhterms(this, n, m, icon, hn, hm, rhsterm, amatnn, amatnm)
Calculate hydraulic head term for this connection.
Definition: gwf-buy.f90:1173
subroutine buy_rp(this)
Check for new buy period data.
Definition: gwf-buy.f90:245
subroutine buy_cf_maw(packobj, hnew, dense, elev, denseref, locdense, locconc, drhodc, crhoref, ctemp, iform)
Pass density information into maw package; density terms are calculated in the maw package as part of...
Definition: gwf-buy.f90:761
subroutine buy_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
Fill coefficients.
Definition: gwf-buy.f90:811
subroutine, public buy_cr(buyobj, name_model, input_mempath, inunit, iout)
Create a new BUY object.
Definition: gwf-buy.f90:104
subroutine calc_ghb_hcof_rhs_terms(denseref, denseghb, densenode, elevghb, elevnode, hghb, hnode, cond, iform, rhsterm, hcofterm)
Calculate density hcof and rhs terms for ghb conditions.
Definition: gwf-buy.f90:498
subroutine allocate_scalars(this)
Allocate scalars used by the package.
Definition: gwf-buy.f90:1308
subroutine buy_cq(this, hnew, flowja)
Add buy term to flowja.
Definition: gwf-buy.f90:888
subroutine buy_calcelev(this)
Calculate cell elevations to use in density flow equations.
Definition: gwf-buy.f90:1290
subroutine set_concentration_pointer(this, modelname, conc, icbund)
Pass in a gwt model name, concentration array and ibound, and store a pointer to these in the BUY pac...
Definition: gwf-buy.f90:1483
subroutine buy_cf_riv(packobj, hnew, dense, elev, denseref, locelev, locdense, locconc, drhodc, crhoref, ctemp, iform)
Fill riv coefficients.
Definition: gwf-buy.f90:544
subroutine set_packagedata(this, input_data)
Sets package data instead of reading from file.
Definition: gwf-buy.f90:1074
subroutine buy_df(this, dis, buy_input)
Read options and package data, or set from argument.
Definition: gwf-buy.f90:128
subroutine source_dimensions(this)
@ brief Source dimensions for package
Definition: gwf-buy.f90:953
subroutine buy_cf_sfr(packobj, hnew, dense, elev, denseref, locdense, locconc, drhodc, crhoref, ctemp, iform)
Pass density information into sfr package; density terms are calculated in the sfr package as part of...
Definition: gwf-buy.f90:708
subroutine buy_ot_dv(this, idvfl)
Save density array to binary file.
Definition: gwf-buy.f90:851
real(dp) function get_bnd_density(n, locdense, locconc, denseref, drhodc, crhoref, ctemp, auxvar)
Return the density of the boundary package using one of several different options in the following or...
Definition: gwf-buy.f90:400
subroutine buy_cf_ghb(packobj, hnew, dense, elev, denseref, locelev, locdense, locconc, drhodc, crhoref, ctemp, iform)
Fill ghb coefficients.
Definition: gwf-buy.f90:438
subroutine buy_ar(this, npf, ibound)
Allocate and Read.
Definition: gwf-buy.f90:175
subroutine buy_da(this)
Deallocate.
Definition: gwf-buy.f90:920
subroutine buy_cf_drn(packobj, hnew, dense, denseref)
Fill drn coefficients.
Definition: gwf-buy.f90:613
subroutine source_options(this)
@ brief Source input options
Definition: gwf-buy.f90:1380
subroutine buy_ad(this)
Advance.
Definition: gwf-buy.f90:276
subroutine log_options(this, found, densityfile)
@ brief log input options
Definition: gwf-buy.f90:1425
subroutine source_packagedata(this)
@ brief source packagedata for package
Definition: gwf-buy.f90:979
subroutine allocate_arrays(this, nodes)
Allocate arrays used by the package.
Definition: gwf-buy.f90:1344
subroutine buy_ar_bnd(this, packobj, hnew)
Buoyancy ar_bnd routine to activate density in packages.
Definition: gwf-buy.f90:200
subroutine calcbuy(this, n, m, icon, hn, hm, buy)
Calculate buyancy term for this connection.
Definition: gwf-buy.f90:1091
subroutine buy_calcdens(this)
calculate fluid density from concentration
Definition: gwf-buy.f90:1266
This module contains stateless conductance functions.
real(dp) function, public hcond(ibdn, ibdm, ictn, ictm, iupstream, ihc, icellavg, condsat, hn, hm, satn, satm, hkn, hkm, topn, topm, botn, botm, cln, clm, fawidth)
Horizontal conductance between two cells.
real(dp) function, public vcond(ibdn, ibdm, ictn, ictm, inewton, ivarcv, idewatcv, condsat, hn, hm, vkn, vkm, satn, satm, topn, topm, botn, botm, flowarea)
Vertical conductance between two cells.
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module defines variable data types.
Definition: kind.f90:8
subroutine, public memorystore_release(varname, memory_path)
Release a single variable from the memory store.
This module contains the base numerical package type.
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
This module contains the SFR package methods.
Definition: gwf-sfr.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_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:27
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:26
@ brief BndType
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Data structure to transfer input configuration to the.