MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwfbuymodule Module Reference

Data Types

type  concentrationpointer
 
type  gwfbuytype
 

Functions/Subroutines

real(dp) function calcdens (denseref, drhodc, crhoref, conc)
 Generic function to calculate fluid density from concentration. More...
 
subroutine, public buy_cr (buyobj, name_model, inunit, iout)
 Create a new BUY object. More...
 
subroutine buy_df (this, dis, buy_input)
 Read options and package data, or set from argument. More...
 
subroutine buy_ar (this, npf, ibound)
 Allocate and Read. More...
 
subroutine buy_ar_bnd (this, packobj, hnew)
 Buoyancy ar_bnd routine to activate density in packages. More...
 
subroutine buy_rp (this)
 Check for new buy period data. More...
 
subroutine buy_ad (this)
 Advance. More...
 
subroutine buy_cf (this, kiter)
 Fill coefficients. More...
 
subroutine buy_cf_bnd (this, packobj, hnew)
 Fill coefficients. More...
 
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 order of priority: More...
 
subroutine buy_cf_ghb (packobj, hnew, dense, elev, denseref, locelev, locdense, locconc, drhodc, crhoref, ctemp, iform)
 Fill ghb coefficients. More...
 
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. More...
 
subroutine buy_cf_riv (packobj, hnew, dense, elev, denseref, locelev, locdense, locconc, drhodc, crhoref, ctemp, iform)
 Fill riv coefficients. More...
 
subroutine buy_cf_drn (packobj, hnew, dense, denseref)
 Fill drn coefficients. More...
 
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 of lak_calculate_density_exchange method. More...
 
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 sfr_calculate_density_exchange method. More...
 
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 maw_calculate_density_exchange method. More...
 
subroutine buy_fc (this, kiter, matrix_sln, idxglo, rhs, hnew)
 Fill coefficients. More...
 
subroutine buy_ot_dv (this, idvfl)
 Save density array to binary file. More...
 
subroutine buy_cq (this, hnew, flowja)
 Add buy term to flowja. More...
 
subroutine buy_da (this)
 Deallocate. More...
 
subroutine read_dimensions (this)
 Read the dimensions for this package. More...
 
subroutine read_packagedata (this)
 Read PACKAGEDATA block. More...
 
subroutine set_packagedata (this, input_data)
 Sets package data instead of reading from file. More...
 
subroutine calcbuy (this, n, m, icon, hn, hm, buy)
 Calculate buyancy term for this connection. More...
 
subroutine calchhterms (this, n, m, icon, hn, hm, rhsterm, amatnn, amatnm)
 Calculate hydraulic head term for this connection. More...
 
subroutine buy_calcdens (this)
 calculate fluid density from concentration More...
 
subroutine buy_calcelev (this)
 Calculate cell elevations to use in density flow equations. More...
 
subroutine allocate_scalars (this)
 Allocate scalars used by the package. More...
 
subroutine allocate_arrays (this, nodes)
 Allocate arrays used by the package. More...
 
subroutine read_options (this)
 Read package options. More...
 
subroutine set_options (this, input_data)
 Sets options as opposed to reading them from a file. More...
 
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 package so that density can be calculated from them. More...
 

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwfbuymodule::allocate_arrays ( class(gwfbuytype this,
integer(i4b), intent(in)  nodes 
)

Definition at line 1440 of file gwf-buy.f90.

1441  ! -- dummy
1442  class(GwfBuyType) :: this
1443  integer(I4B), intent(in) :: nodes
1444  ! -- local
1445  integer(I4B) :: i
1446  !
1447  ! -- Allocate
1448  call mem_allocate(this%dense, nodes, 'DENSE', this%memoryPath)
1449  call mem_allocate(this%concbuy, 0, 'CONCBUY', this%memoryPath)
1450  call mem_allocate(this%elev, nodes, 'ELEV', this%memoryPath)
1451  call mem_allocate(this%drhodc, this%nrhospecies, 'DRHODC', this%memoryPath)
1452  call mem_allocate(this%crhoref, this%nrhospecies, 'CRHOREF', this%memoryPath)
1453  call mem_allocate(this%ctemp, this%nrhospecies, 'CTEMP', this%memoryPath)
1454  allocate (this%cmodelname(this%nrhospecies))
1455  allocate (this%cauxspeciesname(this%nrhospecies))
1456  allocate (this%modelconc(this%nrhospecies))
1457  !
1458  ! -- Initialize
1459  do i = 1, nodes
1460  this%dense(i) = this%denseref
1461  this%elev(i) = dzero
1462  end do
1463  !
1464  ! -- Initialize nrhospecies arrays
1465  do i = 1, this%nrhospecies
1466  this%drhodc(i) = dzero
1467  this%crhoref(i) = dzero
1468  this%ctemp(i) = dzero
1469  this%cmodelname(i) = ''
1470  this%cauxspeciesname(i) = ''
1471  end do
1472  !
1473  ! -- Return
1474  return

◆ allocate_scalars()

subroutine gwfbuymodule::allocate_scalars ( class(gwfbuytype this)
private

Definition at line 1401 of file gwf-buy.f90.

1402  ! -- modules
1403  use constantsmodule, only: dzero
1404  ! -- dummy
1405  class(GwfBuyType) :: this
1406  ! -- local
1407  !
1408  ! -- allocate scalars in NumericalPackageType
1409  call this%NumericalPackageType%allocate_scalars()
1410  !
1411  ! -- Allocate
1412  call mem_allocate(this%ioutdense, 'IOUTDENSE', this%memoryPath)
1413  call mem_allocate(this%iform, 'IFORM', this%memoryPath)
1414  call mem_allocate(this%ireadelev, 'IREADELEV', this%memoryPath)
1415  call mem_allocate(this%ireadconcbuy, 'IREADCONCBUY', this%memoryPath)
1416  call mem_allocate(this%iconcset, 'ICONCSET', this%memoryPath)
1417  call mem_allocate(this%denseref, 'DENSEREF', this%memoryPath)
1418  !
1419  call mem_allocate(this%nrhospecies, 'NRHOSPECIES', this%memoryPath)
1420  !
1421  ! -- Initialize
1422  this%ioutdense = 0
1423  this%ireadelev = 0
1424  this%iconcset = 0
1425  this%ireadconcbuy = 0
1426  this%denseref = 1000.d0
1427  !
1428  this%nrhospecies = 0
1429  !
1430  ! -- Initialize default to LHS implementation of hydraulic head formulation
1431  this%iform = 2
1432  this%iasym = 1
1433  !
1434  ! -- Return
1435  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64

◆ buy_ad()

subroutine gwfbuymodule::buy_ad ( class(gwfbuytype this)

Definition at line 292 of file gwf-buy.f90.

293  ! -- dummy
294  class(GwfBuyType) :: this
295  !
296  ! -- update density using the last concentration
297  call this%buy_calcdens()
298  !
299  ! -- Return
300  return

◆ buy_ar()

subroutine gwfbuymodule::buy_ar ( class(gwfbuytype this,
type(gwfnpftype), intent(in), pointer  npf,
integer(i4b), dimension(:), pointer  ibound 
)
private

Definition at line 182 of file gwf-buy.f90.

183  ! -- dummy
184  class(GwfBuyType) :: this
185  type(GwfNpfType), pointer, intent(in) :: npf
186  integer(I4B), dimension(:), pointer :: ibound
187  !
188  ! -- store pointers to arguments that were passed in
189  this%npf => npf
190  this%ibound => ibound
191  !
192  ! -- Ensure NPF XT3D is not on
193  if (this%npf%ixt3d /= 0) then
194  call store_error('Error in model '//trim(this%name_model)// &
195  '. The XT3D option cannot be used with the BUY Package.')
196  call this%parser%StoreErrorUnit()
197  end if
198  !
199  ! -- Calculate cell elevations
200  call this%buy_calcelev()
201  !
202  ! -- Return
203  return
Here is the call graph for this function:

◆ buy_ar_bnd()

subroutine gwfbuymodule::buy_ar_bnd ( class(gwfbuytype this,
class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew 
)
private

This routine is called from gwf_ar() as it goes through each package

Definition at line 210 of file gwf-buy.f90.

211  ! -- modules
212  use bndmodule, only: bndtype
213  use lakmodule, only: laktype
214  use sfrmodule, only: sfrtype
215  use mawmodule, only: mawtype
216  ! -- dummy
217  class(GwfBuyType) :: this
218  class(BndType), pointer :: packobj
219  real(DP), intent(in), dimension(:) :: hnew
220  !
221  ! -- Add density terms based on boundary package type
222  select case (packobj%filtyp)
223  case ('LAK')
224  !
225  ! -- activate density for lake package
226  select type (packobj)
227  type is (laktype)
228  call packobj%lak_activate_density()
229  end select
230  !
231  case ('SFR')
232  !
233  ! -- activate density for sfr package
234  select type (packobj)
235  type is (sfrtype)
236  call packobj%sfr_activate_density()
237  end select
238  !
239  case ('MAW')
240  !
241  ! -- activate density for maw package
242  select type (packobj)
243  type is (mawtype)
244  call packobj%maw_activate_density()
245  end select
246  !
247  case default
248  !
249  ! -- nothing
250  end select
251  !
252  ! -- Return
253  return
This module contains the base boundary package.
This module contains the SFR package methods.
Definition: gwf-sfr.f90:7
@ brief BndType

◆ buy_calcdens()

subroutine gwfbuymodule::buy_calcdens ( class(gwfbuytype this)

Definition at line 1353 of file gwf-buy.f90.

1354  ! -- dummy
1355  class(GwfBuyType) :: this
1356 
1357  ! -- local
1358  integer(I4B) :: n
1359  integer(I4B) :: i
1360  !
1361  ! -- Calculate the density using the specified concentration array
1362  do n = 1, this%dis%nodes
1363  do i = 1, this%nrhospecies
1364  if (this%modelconc(i)%icbund(n) == 0) then
1365  this%ctemp = dzero
1366  else
1367  this%ctemp(i) = this%modelconc(i)%conc(n)
1368  end if
1369  end do
1370  this%dense(n) = calcdens(this%denseref, this%drhodc, this%crhoref, &
1371  this%ctemp)
1372  end do
1373  !
1374  ! -- Return
1375  return
Here is the call graph for this function:

◆ buy_calcelev()

subroutine gwfbuymodule::buy_calcelev ( class(gwfbuytype this)
private

Definition at line 1380 of file gwf-buy.f90.

1381  ! -- dummy
1382  class(GwfBuyType) :: this
1383  ! -- local
1384  integer(I4B) :: n
1385  real(DP) :: tp, bt, frac
1386  !
1387  ! -- Calculate the elev array
1388  do n = 1, this%dis%nodes
1389  tp = this%dis%top(n)
1390  bt = this%dis%bot(n)
1391  frac = this%npf%sat(n)
1392  this%elev(n) = bt + dhalf * frac * (tp - bt)
1393  end do
1394  !
1395  ! -- Return
1396  return

◆ buy_cf()

subroutine gwfbuymodule::buy_cf ( class(gwfbuytype this,
integer(i4b)  kiter 
)
private

Definition at line 305 of file gwf-buy.f90.

306  ! -- dummy
307  class(GwfBuyType) :: this
308  integer(I4B) :: kiter
309  !
310  ! -- Recalculate the elev array for this iteration
311  if (this%ireadelev == 0) then
312  if (this%iform == 1 .or. this%iform == 2) then
313  call this%buy_calcelev()
314  end if
315  end if
316  !
317  ! -- Return
318  return

◆ buy_cf_bnd()

subroutine gwfbuymodule::buy_cf_bnd ( class(gwfbuytype this,
class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew 
)
private

Definition at line 323 of file gwf-buy.f90.

324  ! -- modules
325  use bndmodule, only: bndtype
326  ! -- dummy
327  class(GwfBuyType) :: this
328  class(BndType), pointer :: packobj
329  real(DP), intent(in), dimension(:) :: hnew
330  ! -- local
331  integer(I4B) :: i, j
332  integer(I4B) :: n, locdense, locelev
333  integer(I4B), dimension(:), allocatable :: locconc
334  !
335  ! -- Return if freshwater head formulation; all boundary heads must be
336  ! entered as freshwater equivalents
337  if (this%iform == 0) return
338  !
339  ! -- initialize
340  locdense = 0
341  locelev = 0
342  allocate (locconc(this%nrhospecies))
343  locconc(:) = 0
344  !
345  ! -- Add buoyancy terms for head-dependent boundaries
346  do n = 1, packobj%naux
347  if (packobj%auxname(n) == 'DENSITY') then
348  locdense = n
349  else if (packobj%auxname(n) == 'ELEVATION') then
350  locelev = n
351  end if
352  end do
353  !
354  ! -- find aux columns for concentrations that affect density
355  do i = 1, this%nrhospecies
356  locconc(i) = 0
357  do j = 1, packobj%naux
358  if (this%cauxspeciesname(i) == packobj%auxname(j)) then
359  locconc(i) = j
360  exit
361  end if
362  end do
363  if (locconc(i) == 0) then
364  ! -- one not found, so don't use and mark all as 0
365  locconc(:) = 0
366  exit
367  end if
368  end do
369  !
370  ! -- Add density terms based on boundary package type
371  select case (packobj%filtyp)
372  case ('GHB')
373  !
374  ! -- general head boundary
375  call buy_cf_ghb(packobj, hnew, this%dense, this%elev, this%denseref, &
376  locelev, locdense, locconc, this%drhodc, this%crhoref, &
377  this%ctemp, this%iform)
378  case ('RIV')
379  !
380  ! -- river
381  call buy_cf_riv(packobj, hnew, this%dense, this%elev, this%denseref, &
382  locelev, locdense, locconc, this%drhodc, this%crhoref, &
383  this%ctemp, this%iform)
384  case ('DRN')
385  !
386  ! -- drain
387  call buy_cf_drn(packobj, hnew, this%dense, this%denseref)
388  case ('LAK')
389  !
390  ! -- lake
391  call buy_cf_lak(packobj, hnew, this%dense, this%elev, this%denseref, &
392  locdense, locconc, this%drhodc, this%crhoref, &
393  this%ctemp, this%iform)
394  case ('SFR')
395  !
396  ! -- sfr
397  call buy_cf_sfr(packobj, hnew, this%dense, this%elev, this%denseref, &
398  locdense, locconc, this%drhodc, this%crhoref, &
399  this%ctemp, this%iform)
400  case ('MAW')
401  !
402  ! -- maw
403  call buy_cf_maw(packobj, hnew, this%dense, this%elev, this%denseref, &
404  locdense, locconc, this%drhodc, this%crhoref, &
405  this%ctemp, this%iform)
406  case default
407  !
408  ! -- nothing
409  end select
410  !
411  ! -- deallocate
412  deallocate (locconc)
413  !
414  ! -- Return
415  return
Here is the call graph for this function:

◆ buy_cf_drn()

subroutine gwfbuymodule::buy_cf_drn ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(in)  dense,
real(dp), intent(in)  denseref 
)

Definition at line 650 of file gwf-buy.f90.

651  ! -- modules
652  use bndmodule, only: bndtype
653  use drnmodule, only: drntype
654  class(BndType), pointer :: packobj
655  ! -- dummy
656  real(DP), intent(in), dimension(:) :: hnew
657  real(DP), intent(in), dimension(:) :: dense
658  real(DP), intent(in) :: denseref
659  ! -- local
660  integer(I4B) :: n
661  integer(I4B) :: node
662  real(DP) :: rho
663  real(DP) :: hbnd
664  real(DP) :: cond
665  real(DP) :: hcofterm
666  real(DP) :: rhsterm
667  !
668  ! -- Process density terms for each DRN
669  select type (packobj)
670  type is (drntype)
671  do n = 1, packobj%nbound
672  node = packobj%nodelist(n)
673  if (packobj%ibound(node) <= 0) cycle
674  rho = dense(node)
675  hbnd = packobj%elev(n)
676  cond = packobj%cond(n)
677  if (hnew(node) > hbnd) then
678  hcofterm = -cond * (rho / denseref - done)
679  rhsterm = hcofterm * hbnd
680  packobj%hcof(n) = packobj%hcof(n) + hcofterm
681  packobj%rhs(n) = packobj%rhs(n) + rhsterm
682  end if
683  end do
684  end select
685  !
686  ! -- Return
687  return
Here is the caller graph for this function:

◆ buy_cf_ghb()

subroutine gwfbuymodule::buy_cf_ghb ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(in)  dense,
real(dp), dimension(:), intent(in)  elev,
real(dp), intent(in)  denseref,
integer(i4b), intent(in)  locelev,
integer(i4b), intent(in)  locdense,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  drhodc,
real(dp), dimension(:), intent(in)  crhoref,
real(dp), dimension(:), intent(inout)  ctemp,
integer(i4b), intent(in)  iform 
)
private

Definition at line 464 of file gwf-buy.f90.

467  ! -- modules
468  use bndmodule, only: bndtype
469  use ghbmodule, only: ghbtype
470  class(BndType), pointer :: packobj
471  ! -- dummy
472  real(DP), intent(in), dimension(:) :: hnew
473  real(DP), intent(in), dimension(:) :: dense
474  real(DP), intent(in), dimension(:) :: elev
475  real(DP), intent(in) :: denseref
476  integer(I4B), intent(in) :: locelev
477  integer(I4B), intent(in) :: locdense
478  integer(I4B), dimension(:), intent(in) :: locconc
479  real(DP), dimension(:), intent(in) :: drhodc
480  real(DP), dimension(:), intent(in) :: crhoref
481  real(DP), dimension(:), intent(inout) :: ctemp
482  integer(I4B), intent(in) :: iform
483  ! -- local
484  integer(I4B) :: n
485  integer(I4B) :: node
486  real(DP) :: denseghb
487  real(DP) :: elevghb
488  real(DP) :: hghb
489  real(DP) :: cond
490  real(DP) :: hcofterm, rhsterm
491  !
492  ! -- Process density terms for each GHB
493  select type (packobj)
494  type is (ghbtype)
495  do n = 1, packobj%nbound
496  node = packobj%nodelist(n)
497  if (packobj%ibound(node) <= 0) cycle
498  !
499  ! -- density
500  denseghb = get_bnd_density(n, locdense, locconc, denseref, &
501  drhodc, crhoref, ctemp, packobj%auxvar)
502  !
503  ! -- elevation
504  elevghb = elev(node)
505  if (locelev > 0) elevghb = packobj%auxvar(locelev, n)
506  !
507  ! -- boundary head and conductance
508  hghb = packobj%bhead(n)
509  cond = packobj%cond(n)
510  !
511  ! -- calculate HCOF and RHS terms
512  call calc_ghb_hcof_rhs_terms(denseref, denseghb, dense(node), &
513  elevghb, elev(node), hghb, hnew(node), &
514  cond, iform, rhsterm, hcofterm)
515  packobj%hcof(n) = packobj%hcof(n) + hcofterm
516  packobj%rhs(n) = packobj%rhs(n) - rhsterm
517  !
518  end do
519  end select
520  !
521  ! -- Return
522  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ buy_cf_lak()

subroutine gwfbuymodule::buy_cf_lak ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(in)  dense,
real(dp), dimension(:), intent(in)  elev,
real(dp), intent(in)  denseref,
integer(i4b), intent(in)  locdense,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  drhodc,
real(dp), dimension(:), intent(in)  crhoref,
real(dp), dimension(:), intent(inout)  ctemp,
integer(i4b), intent(in)  iform 
)

Definition at line 694 of file gwf-buy.f90.

696  ! -- modules
697  use bndmodule, only: bndtype
698  use lakmodule, only: laktype
699  class(BndType), pointer :: packobj
700  ! -- dummy
701  real(DP), intent(in), dimension(:) :: hnew
702  real(DP), intent(in), dimension(:) :: dense
703  real(DP), intent(in), dimension(:) :: elev
704  real(DP), intent(in) :: denseref
705  integer(I4B), intent(in) :: locdense
706  integer(I4B), dimension(:), intent(in) :: locconc
707  real(DP), dimension(:), intent(in) :: drhodc
708  real(DP), dimension(:), intent(in) :: crhoref
709  real(DP), dimension(:), intent(inout) :: ctemp
710  integer(I4B), intent(in) :: iform
711  ! -- local
712  integer(I4B) :: n
713  integer(I4B) :: node
714  real(DP) :: denselak
715  !
716  ! -- Insert the lake and gwf relative densities into col 1 and 2 and the
717  ! gwf elevation into col 3 of the lake package denseterms array
718  select type (packobj)
719  type is (laktype)
720  do n = 1, packobj%nbound
721  !
722  ! -- get gwf node number
723  node = packobj%nodelist(n)
724  if (packobj%ibound(node) <= 0) cycle
725  !
726  ! -- Determine lak density
727  denselak = get_bnd_density(n, locdense, locconc, denseref, &
728  drhodc, crhoref, ctemp, packobj%auxvar)
729  !
730  ! -- fill lak relative density into column 1 of denseterms
731  packobj%denseterms(1, n) = denselak / denseref
732  !
733  ! -- fill gwf relative density into column 2 of denseterms
734  packobj%denseterms(2, n) = dense(node) / denseref
735  !
736  ! -- fill gwf elevation into column 3 of denseterms
737  packobj%denseterms(3, n) = elev(node)
738  !
739  end do
740  end select
741  !
742  ! -- Return
743  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ buy_cf_maw()

subroutine gwfbuymodule::buy_cf_maw ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(in)  dense,
real(dp), dimension(:), intent(in)  elev,
real(dp), intent(in)  denseref,
integer(i4b), intent(in)  locdense,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  drhodc,
real(dp), dimension(:), intent(in)  crhoref,
real(dp), dimension(:), intent(inout)  ctemp,
integer(i4b), intent(in)  iform 
)

Definition at line 806 of file gwf-buy.f90.

808  ! -- modules
809  use bndmodule, only: bndtype
810  use mawmodule, only: mawtype
811  class(BndType), pointer :: packobj
812  ! -- dummy
813  real(DP), intent(in), dimension(:) :: hnew
814  real(DP), intent(in), dimension(:) :: dense
815  real(DP), intent(in), dimension(:) :: elev
816  real(DP), intent(in) :: denseref
817  integer(I4B), intent(in) :: locdense
818  integer(I4B), dimension(:), intent(in) :: locconc
819  real(DP), dimension(:), intent(in) :: drhodc
820  real(DP), dimension(:), intent(in) :: crhoref
821  real(DP), dimension(:), intent(inout) :: ctemp
822  integer(I4B), intent(in) :: iform
823  ! -- local
824  integer(I4B) :: n
825  integer(I4B) :: node
826  real(DP) :: densemaw
827  !
828  ! -- Insert the maw and gwf relative densities into col 1 and 2 and the
829  ! gwf elevation into col 3 of the maw package denseterms array
830  select type (packobj)
831  type is (mawtype)
832  do n = 1, packobj%nbound
833  !
834  ! -- get gwf node number
835  node = packobj%nodelist(n)
836  if (packobj%ibound(node) <= 0) cycle
837  !
838  ! -- Determine maw density
839  densemaw = get_bnd_density(n, locdense, locconc, denseref, &
840  drhodc, crhoref, ctemp, packobj%auxvar)
841  !
842  ! -- fill maw relative density into column 1 of denseterms
843  packobj%denseterms(1, n) = densemaw / denseref
844  !
845  ! -- fill gwf relative density into column 2 of denseterms
846  packobj%denseterms(2, n) = dense(node) / denseref
847  !
848  ! -- fill gwf elevation into column 3 of denseterms
849  packobj%denseterms(3, n) = elev(node)
850  !
851  end do
852  end select
853  !
854  ! -- Return
855  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ buy_cf_riv()

subroutine gwfbuymodule::buy_cf_riv ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(in)  dense,
real(dp), dimension(:), intent(in)  elev,
real(dp), intent(in)  denseref,
integer(i4b), intent(in)  locelev,
integer(i4b), intent(in)  locdense,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  drhodc,
real(dp), dimension(:), intent(in)  crhoref,
real(dp), dimension(:), intent(inout)  ctemp,
integer(i4b), intent(in)  iform 
)
private

Definition at line 576 of file gwf-buy.f90.

579  ! -- modules
580  use bndmodule, only: bndtype
581  use rivmodule, only: rivtype
582  class(BndType), pointer :: packobj
583  ! -- dummy
584  real(DP), intent(in), dimension(:) :: hnew
585  real(DP), intent(in), dimension(:) :: dense
586  real(DP), intent(in), dimension(:) :: elev
587  real(DP), intent(in) :: denseref
588  integer(I4B), intent(in) :: locelev
589  integer(I4B), intent(in) :: locdense
590  integer(I4B), dimension(:), intent(in) :: locconc
591  real(DP), dimension(:), intent(in) :: drhodc
592  real(DP), dimension(:), intent(in) :: crhoref
593  real(DP), dimension(:), intent(inout) :: ctemp
594  integer(I4B), intent(in) :: iform
595  ! -- local
596  integer(I4B) :: n
597  integer(I4B) :: node
598  real(DP) :: denseriv
599  real(DP) :: elevriv
600  real(DP) :: hriv
601  real(DP) :: rbot
602  real(DP) :: cond
603  real(DP) :: hcofterm
604  real(DP) :: rhsterm
605  !
606  ! -- Process density terms for each RIV
607  select type (packobj)
608  type is (rivtype)
609  do n = 1, packobj%nbound
610  node = packobj%nodelist(n)
611  if (packobj%ibound(node) <= 0) cycle
612  !
613  ! -- density
614  denseriv = get_bnd_density(n, locdense, locconc, denseref, &
615  drhodc, crhoref, ctemp, packobj%auxvar)
616  !
617  ! -- elevation
618  elevriv = elev(node)
619  if (locelev > 0) elevriv = packobj%auxvar(locelev, n)
620  !
621  ! -- boundary head and conductance
622  hriv = packobj%stage(n)
623  cond = packobj%cond(n)
624  rbot = packobj%rbot(n)
625  !
626  ! -- calculate and add terms depending on whether head is above rbot
627  if (hnew(node) > rbot) then
628  !
629  ! --calculate HCOF and RHS terms, similar to GHB in this case
630  call calc_ghb_hcof_rhs_terms(denseref, denseriv, dense(node), &
631  elevriv, elev(node), hriv, hnew(node), &
632  cond, iform, rhsterm, hcofterm)
633  else
634  hcofterm = dzero
635  rhsterm = cond * (denseriv / denseref - done) * (hriv - rbot)
636  end if
637  !
638  ! -- Add terms to package hcof and rhs accumulators
639  packobj%hcof(n) = packobj%hcof(n) + hcofterm
640  packobj%rhs(n) = packobj%rhs(n) - rhsterm
641  end do
642  end select
643  !
644  ! -- Return
645  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ buy_cf_sfr()

subroutine gwfbuymodule::buy_cf_sfr ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(in)  dense,
real(dp), dimension(:), intent(in)  elev,
real(dp), intent(in)  denseref,
integer(i4b), intent(in)  locdense,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  drhodc,
real(dp), dimension(:), intent(in)  crhoref,
real(dp), dimension(:), intent(inout)  ctemp,
integer(i4b), intent(in)  iform 
)

Definition at line 750 of file gwf-buy.f90.

752  ! -- modules
753  use bndmodule, only: bndtype
754  use sfrmodule, only: sfrtype
755  class(BndType), pointer :: packobj
756  ! -- dummy
757  real(DP), intent(in), dimension(:) :: hnew
758  real(DP), intent(in), dimension(:) :: dense
759  real(DP), intent(in), dimension(:) :: elev
760  real(DP), intent(in) :: denseref
761  integer(I4B), intent(in) :: locdense
762  integer(I4B), dimension(:), intent(in) :: locconc
763  real(DP), dimension(:), intent(in) :: drhodc
764  real(DP), dimension(:), intent(in) :: crhoref
765  real(DP), dimension(:), intent(inout) :: ctemp
766  integer(I4B), intent(in) :: iform
767  ! -- local
768  integer(I4B) :: n
769  integer(I4B) :: node
770  real(DP) :: densesfr
771  !
772  ! -- Insert the sfr and gwf relative densities into col 1 and 2 and the
773  ! gwf elevation into col 3 of the sfr package denseterms array
774  select type (packobj)
775  type is (sfrtype)
776  do n = 1, packobj%nbound
777  !
778  ! -- get gwf node number
779  node = packobj%nodelist(n)
780  if (packobj%ibound(node) <= 0) cycle
781  !
782  ! -- Determine sfr density
783  densesfr = get_bnd_density(n, locdense, locconc, denseref, &
784  drhodc, crhoref, ctemp, packobj%auxvar)
785  !
786  ! -- fill sfr relative density into column 1 of denseterms
787  packobj%denseterms(1, n) = densesfr / denseref
788  !
789  ! -- fill gwf relative density into column 2 of denseterms
790  packobj%denseterms(2, n) = dense(node) / denseref
791  !
792  ! -- fill gwf elevation into column 3 of denseterms
793  packobj%denseterms(3, n) = elev(node)
794  !
795  end do
796  end select
797  !
798  ! -- Return
799  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ buy_cq()

subroutine gwfbuymodule::buy_cq ( class(gwfbuytype this,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(inout)  flowja 
)
private

Definition at line 943 of file gwf-buy.f90.

944  implicit none
945  class(GwfBuyType) :: this
946  real(DP), intent(in), dimension(:) :: hnew
947  real(DP), intent(inout), dimension(:) :: flowja
948  integer(I4B) :: n, m, ipos
949  real(DP) :: deltaQ
950  real(DP) :: rhsterm, amatnn, amatnm
951  !
952  ! -- Calculate the flow across each cell face and store in flowja
953  do n = 1, this%dis%nodes
954  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
955  m = this%dis%con%ja(ipos)
956  if (m < n) cycle
957  if (this%iform == 0) then
958  ! -- equivalent freshwater head formulation
959  call this%calcbuy(n, m, ipos, hnew(n), hnew(m), deltaq)
960  else
961  ! -- hydraulic head formulation
962  call this%calchhterms(n, m, ipos, hnew(n), hnew(m), rhsterm, &
963  amatnn, amatnm)
964  deltaq = amatnm * hnew(m) - amatnn * hnew(n) + rhsterm
965  end if
966  flowja(ipos) = flowja(ipos) + deltaq
967  flowja(this%dis%con%isym(ipos)) = flowja(this%dis%con%isym(ipos)) - &
968  deltaq
969  end do
970  end do
971  !
972  ! -- Return
973  return

◆ buy_cr()

subroutine, public gwfbuymodule::buy_cr ( type(gwfbuytype), pointer  buyobj,
character(len=*), intent(in)  name_model,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Definition at line 105 of file gwf-buy.f90.

106  ! -- dummy
107  type(GwfBuyType), pointer :: buyobj
108  character(len=*), intent(in) :: name_model
109  integer(I4B), intent(in) :: inunit
110  integer(I4B), intent(in) :: iout
111  !
112  ! -- Create the object
113  allocate (buyobj)
114  !
115  ! -- create name and memory path
116  call buyobj%set_names(1, name_model, 'BUY', 'BUY')
117  !
118  ! -- Allocate scalars
119  call buyobj%allocate_scalars()
120  !
121  ! -- Set variables
122  buyobj%inunit = inunit
123  buyobj%iout = iout
124  !
125  ! -- Initialize block parser
126  call buyobj%parser%Initialize(buyobj%inunit, buyobj%iout)
127  !
128  ! -- Return
129  return
Here is the caller graph for this function:

◆ buy_da()

subroutine gwfbuymodule::buy_da ( class(gwfbuytype this)
private

Definition at line 978 of file gwf-buy.f90.

979  ! -- dummy
980  class(GwfBuyType) :: this
981  !
982  ! -- Deallocate arrays if package was active
983  if (this%inunit > 0) then
984  call mem_deallocate(this%elev)
985  call mem_deallocate(this%dense)
986  call mem_deallocate(this%concbuy)
987  call mem_deallocate(this%drhodc)
988  call mem_deallocate(this%crhoref)
989  call mem_deallocate(this%ctemp)
990  deallocate (this%cmodelname)
991  deallocate (this%cauxspeciesname)
992  deallocate (this%modelconc)
993  end if
994  !
995  ! -- Scalars
996  call mem_deallocate(this%ioutdense)
997  call mem_deallocate(this%iform)
998  call mem_deallocate(this%ireadelev)
999  call mem_deallocate(this%ireadconcbuy)
1000  call mem_deallocate(this%iconcset)
1001  call mem_deallocate(this%denseref)
1002  !
1003  call mem_deallocate(this%nrhospecies)
1004  !
1005  ! -- deallocate parent
1006  call this%NumericalPackageType%da()
1007  !
1008  ! -- Return
1009  return

◆ buy_df()

subroutine gwfbuymodule::buy_df ( class(gwfbuytype this,
class(disbasetype), intent(in), pointer  dis,
type(gwfbuyinputdatatype), intent(in), optional  buy_input 
)
private
Parameters
thisthis buoyancy package
[in]dispointer to discretization
[in]buy_inputoptional buy input data, otherwise read from file

Definition at line 134 of file gwf-buy.f90.

135  ! -- dummy
136  class(GwfBuyType) :: this !< this buoyancy package
137  class(DisBaseType), pointer, intent(in) :: dis !< pointer to discretization
138  type(GwfBuyInputDataType), optional, intent(in) :: buy_input !< optional buy input data, otherwise read from file
139  ! -- local
140  ! -- formats
141  character(len=*), parameter :: fmtbuy = &
142  "(1x,/1x,'BUY -- Buoyancy Package, Version 1, 5/16/2018', &
143  &' input read from unit ', i0, //)"
144  !
145  ! --print a message identifying the buoyancy package.
146  write (this%iout, fmtbuy) this%inunit
147  !
148  ! -- store pointers to arguments that were passed in
149  this%dis => dis
150 
151  if (.not. present(buy_input)) then
152  !
153  ! -- Read buoyancy options
154  call this%read_options()
155  !
156  ! -- Read buoyancy dimensions
157  call this%read_dimensions()
158  else
159  ! set from input data instead
160  call this%set_options(buy_input)
161  this%nrhospecies = buy_input%nrhospecies
162  end if
163  !
164  ! -- Allocate arrays
165  call this%allocate_arrays(dis%nodes)
166 
167  if (.not. present(buy_input)) then
168  !
169  ! -- Read buoyancy packagedata
170  call this%read_packagedata()
171  else
172  ! set from input data instead
173  call this%set_packagedata(buy_input)
174  end if
175  !
176  ! -- Return
177  return

◆ buy_fc()

subroutine gwfbuymodule::buy_fc ( class(gwfbuytype this,
integer(i4b)  kiter,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs,
real(dp), dimension(:), intent(inout)  hnew 
)

Definition at line 860 of file gwf-buy.f90.

861  ! -- dummy
862  class(GwfBuyType) :: this
863  integer(I4B) :: kiter
864  class(MatrixBaseType), pointer :: matrix_sln
865  integer(I4B), intent(in), dimension(:) :: idxglo
866  real(DP), dimension(:), intent(inout) :: rhs
867  real(DP), intent(inout), dimension(:) :: hnew
868  ! -- local
869  integer(I4B) :: n, m, ipos, idiag
870  real(DP) :: rhsterm, amatnn, amatnm
871  !
872  ! -- initialize
873  amatnn = dzero
874  amatnm = dzero
875  !
876  ! -- fill buoyancy flow term
877  do n = 1, this%dis%nodes
878  if (this%ibound(n) == 0) cycle
879  idiag = this%dis%con%ia(n)
880  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
881  m = this%dis%con%ja(ipos)
882  if (this%ibound(m) == 0) cycle
883  if (this%iform == 0) then
884  call this%calcbuy(n, m, ipos, hnew(n), hnew(m), rhsterm)
885  else
886  call this%calchhterms(n, m, ipos, hnew(n), hnew(m), rhsterm, &
887  amatnn, amatnm)
888  end if
889  !
890  ! -- Add terms to rhs, diagonal, and off diagonal
891  rhs(n) = rhs(n) - rhsterm
892  call matrix_sln%add_value_pos(idxglo(idiag), -amatnn)
893  call matrix_sln%add_value_pos(idxglo(ipos), amatnm)
894  end do
895  end do
896  !
897  ! -- Return
898  return

◆ buy_ot_dv()

subroutine gwfbuymodule::buy_ot_dv ( class(gwfbuytype this,
integer(i4b), intent(in)  idvfl 
)
private

Definition at line 903 of file gwf-buy.f90.

904  ! -- dummy
905  class(GwfBuyType) :: this
906  integer(I4B), intent(in) :: idvfl
907  ! -- local
908  character(len=1) :: cdatafmp = ' ', editdesc = ' '
909  integer(I4B) :: ibinun
910  integer(I4B) :: iprint
911  integer(I4B) :: nvaluesp
912  integer(I4B) :: nwidthp
913  real(DP) :: dinact
914  !
915  ! -- Set unit number for density output
916  if (this%ioutdense /= 0) then
917  ibinun = 1
918  else
919  ibinun = 0
920  end if
921  if (idvfl == 0) ibinun = 0
922  !
923  ! -- save density array
924  if (ibinun /= 0) then
925  iprint = 0
926  dinact = dhnoflo
927  !
928  ! -- write density to binary file
929  if (this%ioutdense /= 0) then
930  ibinun = this%ioutdense
931  call this%dis%record_array(this%dense, this%iout, iprint, ibinun, &
932  ' DENSITY', cdatafmp, nvaluesp, &
933  nwidthp, editdesc, dinact)
934  end if
935  end if
936  !
937  ! -- Return
938  return

◆ buy_rp()

subroutine gwfbuymodule::buy_rp ( class(gwfbuytype this)

Definition at line 258 of file gwf-buy.f90.

259  ! -- modules
260  use tdismodule, only: kstp, kper
261  ! -- dummy
262  class(GwfBuyType) :: this
263  ! -- local
264  character(len=LINELENGTH) :: errmsg
265  integer(I4B) :: i
266  ! -- formats
267  character(len=*), parameter :: fmtc = &
268  "('Buoyancy Package does not have a concentration set &
269  &for species ',i0,'. One or more model names may be specified &
270  &incorrectly in the PACKAGEDATA block or a gwf-gwt exchange may need &
271  &to be activated.')"
272  !
273  ! -- Check to make sure all concentration pointers have been set
274  if (kstp * kper == 1) then
275  do i = 1, this%nrhospecies
276  if (.not. associated(this%modelconc(i)%conc)) then
277  write (errmsg, fmtc) i
278  call store_error(errmsg)
279  end if
280  end do
281  if (count_errors() > 0) then
282  call this%parser%StoreErrorUnit()
283  end if
284  end if
285  !
286  ! -- Return
287  return
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ calc_ghb_hcof_rhs_terms()

subroutine gwfbuymodule::calc_ghb_hcof_rhs_terms ( real(dp), intent(in)  denseref,
real(dp), intent(in)  denseghb,
real(dp), intent(in)  densenode,
real(dp), intent(in)  elevghb,
real(dp), intent(in)  elevnode,
real(dp), intent(in)  hghb,
real(dp), intent(in)  hnode,
real(dp), intent(in)  cond,
integer(i4b), intent(in)  iform,
real(dp), intent(inout)  rhsterm,
real(dp), intent(inout)  hcofterm 
)

Definition at line 527 of file gwf-buy.f90.

530  ! -- dummy
531  real(DP), intent(in) :: denseref
532  real(DP), intent(in) :: denseghb
533  real(DP), intent(in) :: densenode
534  real(DP), intent(in) :: elevghb
535  real(DP), intent(in) :: elevnode
536  real(DP), intent(in) :: hghb
537  real(DP), intent(in) :: hnode
538  real(DP), intent(in) :: cond
539  integer(I4B), intent(in) :: iform
540  real(DP), intent(inout) :: rhsterm
541  real(DP), intent(inout) :: hcofterm
542  ! -- local
543  real(DP) :: t1, t2
544  real(DP) :: avgdense, avgelev
545  !
546  ! -- Calculate common terms
547  avgdense = dhalf * denseghb + dhalf * densenode
548  avgelev = dhalf * elevghb + dhalf * elevnode
549  t1 = avgdense / denseref - done
550  t2 = (denseghb - densenode) / denseref
551  !
552  ! -- Add hcof terms
553  hcofterm = -cond * t1
554  if (iform == 2) then
555  !
556  ! -- this term goes on RHS for iform == 1
557  hcofterm = hcofterm + dhalf * cond * t2
558  end if
559  !
560  ! -- Add rhs terms
561  rhsterm = cond * t1 * hghb
562  rhsterm = rhsterm - cond * t2 * avgelev
563  rhsterm = rhsterm + dhalf * cond * t2 * hghb
564  if (iform == 1) then
565  !
566  ! -- this term goes on LHS for iform == 2
567  rhsterm = rhsterm + dhalf * cond * t2 * hnode
568  end if
569  !
570  ! -- Return
571  return
Here is the caller graph for this function:

◆ calcbuy()

subroutine gwfbuymodule::calcbuy ( class(gwfbuytype this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
integer(i4b), intent(in)  icon,
real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
real(dp), intent(inout)  buy 
)
private

Definition at line 1172 of file gwf-buy.f90.

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) :: buy
1183  ! -- local
1184  integer(I4B) :: ihc
1185  real(DP) :: densen, densem, cl1, cl2, avgdense, wt, elevn, elevm, &
1186  cond, tp, bt
1187  real(DP) :: hyn
1188  real(DP) :: hym
1189  !
1190  ! -- Average density
1191  densen = this%dense(n)
1192  densem = this%dense(m)
1193  if (m > n) then
1194  cl1 = this%dis%con%cl1(this%dis%con%jas(icon))
1195  cl2 = this%dis%con%cl2(this%dis%con%jas(icon))
1196  else
1197  cl1 = this%dis%con%cl2(this%dis%con%jas(icon))
1198  cl2 = this%dis%con%cl1(this%dis%con%jas(icon))
1199  end if
1200  wt = cl1 / (cl1 + cl2)
1201  avgdense = wt * densen + (done - wt) * densem
1202  !
1203  ! -- Elevations
1204  if (this%ireadelev == 0) then
1205  tp = this%dis%top(n)
1206  bt = this%dis%bot(n)
1207  elevn = bt + dhalf * this%npf%sat(n) * (tp - bt)
1208  tp = this%dis%top(m)
1209  bt = this%dis%bot(m)
1210  elevm = bt + dhalf * this%npf%sat(m) * (tp - bt)
1211  else
1212  elevn = this%elev(n)
1213  elevm = this%elev(m)
1214  end if
1215  !
1216  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
1217  hyn = this%npf%hy_eff(n, m, ihc, ipos=icon)
1218  hym = this%npf%hy_eff(m, n, ihc, ipos=icon)
1219  !
1220  ! -- Conductance
1221  if (this%dis%con%ihc(this%dis%con%jas(icon)) == 0) then
1222  cond = vcond(this%ibound(n), this%ibound(m), &
1223  this%npf%icelltype(n), this%npf%icelltype(m), &
1224  this%npf%inewton, &
1225  this%npf%ivarcv, this%npf%idewatcv, &
1226  this%npf%condsat(this%dis%con%jas(icon)), hn, hm, &
1227  hyn, hym, &
1228  this%npf%sat(n), this%npf%sat(m), &
1229  this%dis%top(n), this%dis%top(m), &
1230  this%dis%bot(n), this%dis%bot(m), &
1231  this%dis%con%hwva(this%dis%con%jas(icon)))
1232  else
1233  cond = hcond(this%ibound(n), this%ibound(m), &
1234  this%npf%icelltype(n), this%npf%icelltype(m), &
1235  this%npf%inewton, &
1236  this%dis%con%ihc(this%dis%con%jas(icon)), &
1237  this%npf%icellavg, &
1238  this%npf%condsat(this%dis%con%jas(icon)), &
1239  hn, hm, this%npf%sat(n), this%npf%sat(m), &
1240  hyn, hym, &
1241  this%dis%top(n), this%dis%top(m), &
1242  this%dis%bot(n), this%dis%bot(m), &
1243  this%dis%con%cl1(this%dis%con%jas(icon)), &
1244  this%dis%con%cl2(this%dis%con%jas(icon)), &
1245  this%dis%con%hwva(this%dis%con%jas(icon)))
1246  end if
1247  !
1248  ! -- Calculate buoyancy term
1249  buy = cond * (avgdense - this%denseref) / this%denseref * (elevm - elevn)
1250  !
1251  ! -- Return
1252  return
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.
Here is the call graph for this function:

◆ calcdens()

real(dp) function gwfbuymodule::calcdens ( real(dp), intent(in)  denseref,
real(dp), dimension(:), intent(in)  drhodc,
real(dp), dimension(:), intent(in)  crhoref,
real(dp), dimension(:), intent(in)  conc 
)
private

Definition at line 81 of file gwf-buy.f90.

82  ! -- dummy
83  real(DP), intent(in) :: denseref
84  real(DP), dimension(:), intent(in) :: drhodc
85  real(DP), dimension(:), intent(in) :: crhoref
86  real(DP), dimension(:), intent(in) :: conc
87  ! -- return
88  real(DP) :: dense
89  ! -- local
90  integer(I4B) :: nrhospec
91  integer(I4B) :: i
92  !
93  nrhospec = size(drhodc)
94  dense = denseref
95  do i = 1, nrhospec
96  dense = dense + drhodc(i) * (conc(i) - crhoref(i))
97  end do
98  !
99  ! -- Return
100  return
Here is the caller graph for this function:

◆ calchhterms()

subroutine gwfbuymodule::calchhterms ( class(gwfbuytype this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
integer(i4b), intent(in)  icon,
real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
real(dp), intent(inout)  rhsterm,
real(dp), intent(inout)  amatnn,
real(dp), intent(inout)  amatnm 
)

Definition at line 1257 of file gwf-buy.f90.

1258  ! -- modules
1260  ! -- dummy
1261  class(GwfBuyType) :: this
1262  integer(I4B), intent(in) :: n
1263  integer(I4B), intent(in) :: m
1264  integer(I4B), intent(in) :: icon
1265  real(DP), intent(in) :: hn
1266  real(DP), intent(in) :: hm
1267  real(DP), intent(inout) :: rhsterm
1268  real(DP), intent(inout) :: amatnn
1269  real(DP), intent(inout) :: amatnm
1270  ! -- local
1271  integer(I4B) :: ihc
1272  real(DP) :: densen, densem, cl1, cl2, avgdense, wt, elevn, elevm, cond
1273  real(DP) :: rhonormn, rhonormm
1274  real(DP) :: rhoterm
1275  real(DP) :: elevnm
1276  real(DP) :: hphi
1277  real(DP) :: hyn
1278  real(DP) :: hym
1279  !
1280  ! -- Average density
1281  densen = this%dense(n)
1282  densem = this%dense(m)
1283  if (m > n) then
1284  cl1 = this%dis%con%cl1(this%dis%con%jas(icon))
1285  cl2 = this%dis%con%cl2(this%dis%con%jas(icon))
1286  else
1287  cl1 = this%dis%con%cl2(this%dis%con%jas(icon))
1288  cl2 = this%dis%con%cl1(this%dis%con%jas(icon))
1289  end if
1290  wt = cl1 / (cl1 + cl2)
1291  avgdense = wt * densen + (1.0 - wt) * densem
1292  !
1293  ! -- Elevations
1294  elevn = this%elev(n)
1295  elevm = this%elev(m)
1296  elevnm = (done - wt) * elevn + wt * elevm
1297  !
1298  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
1299  hyn = this%npf%hy_eff(n, m, ihc, ipos=icon)
1300  hym = this%npf%hy_eff(m, n, ihc, ipos=icon)
1301  !
1302  ! -- Conductance
1303  if (ihc == 0) then
1304  cond = vcond(this%ibound(n), this%ibound(m), &
1305  this%npf%icelltype(n), this%npf%icelltype(m), &
1306  this%npf%inewton, &
1307  this%npf%ivarcv, this%npf%idewatcv, &
1308  this%npf%condsat(this%dis%con%jas(icon)), hn, hm, &
1309  hyn, hym, &
1310  this%npf%sat(n), this%npf%sat(m), &
1311  this%dis%top(n), this%dis%top(m), &
1312  this%dis%bot(n), this%dis%bot(m), &
1313  this%dis%con%hwva(this%dis%con%jas(icon)))
1314  else
1315  cond = hcond(this%ibound(n), this%ibound(m), &
1316  this%npf%icelltype(n), this%npf%icelltype(m), &
1317  this%npf%inewton, &
1318  this%dis%con%ihc(this%dis%con%jas(icon)), &
1319  this%npf%icellavg, &
1320  this%npf%condsat(this%dis%con%jas(icon)), &
1321  hn, hm, this%npf%sat(n), this%npf%sat(m), &
1322  hyn, hym, &
1323  this%dis%top(n), this%dis%top(m), &
1324  this%dis%bot(n), this%dis%bot(m), &
1325  this%dis%con%cl1(this%dis%con%jas(icon)), &
1326  this%dis%con%cl2(this%dis%con%jas(icon)), &
1327  this%dis%con%hwva(this%dis%con%jas(icon)))
1328  end if
1329  !
1330  ! -- Calculate terms
1331  rhonormn = densen / this%denseref
1332  rhonormm = densem / this%denseref
1333  rhoterm = wt * rhonormn + (done - wt) * rhonormm
1334  amatnn = cond * (rhoterm - done)
1335  amatnm = amatnn
1336  rhsterm = -cond * (rhonormm - rhonormn) * elevnm
1337  if (this%iform == 1) then
1338  ! -- rhs (lag the h terms and keep matrix symmetric)
1339  hphi = (done - wt) * hn + wt * hm
1340  rhsterm = rhsterm + cond * hphi * (rhonormm - rhonormn)
1341  else
1342  ! -- lhs, results in asymmetric matrix due to weight term
1343  amatnn = amatnn - cond * (done - wt) * (rhonormm - rhonormn)
1344  amatnm = amatnm + cond * wt * (rhonormm - rhonormn)
1345  end if
1346  !
1347  ! -- Return
1348  return
Here is the call graph for this function:

◆ get_bnd_density()

real(dp) function gwfbuymodule::get_bnd_density ( integer(i4b), intent(in)  n,
integer(i4b), intent(in)  locdense,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), intent(in)  denseref,
real(dp), dimension(:), intent(in)  drhodc,
real(dp), dimension(:), intent(in)  crhoref,
real(dp), dimension(:), intent(inout)  ctemp,
real(dp), dimension(:, :), intent(in)  auxvar 
)
  1. Assign as aux variable in column with name 'DENSITY'
  2. Calculate using equation of state and nrhospecies aux columns
  3. If neither of those, then assign as denseref

Definition at line 424 of file gwf-buy.f90.

426  ! -- dummy
427  integer(I4B), intent(in) :: n
428  integer(I4B), intent(in) :: locdense
429  integer(I4B), dimension(:), intent(in) :: locconc
430  real(DP), intent(in) :: denseref
431  real(DP), dimension(:), intent(in) :: drhodc
432  real(DP), dimension(:), intent(in) :: crhoref
433  real(DP), dimension(:), intent(inout) :: ctemp
434  real(DP), dimension(:, :), intent(in) :: auxvar
435  ! -- return
436  real(DP) :: densebnd
437  ! -- local
438  integer(I4B) :: i
439  !
440  ! -- assign boundary density based on one of three options
441  if (locdense > 0) then
442  ! -- assign density to an aux column named 'DENSITY'
443  densebnd = auxvar(locdense, n)
444  else if (locconc(1) > 0) then
445  ! -- calculate density using one or more concentration auxcolumns
446  do i = 1, size(locconc)
447  ctemp(i) = dzero
448  if (locconc(i) > 0) then
449  ctemp(i) = auxvar(locconc(i), n)
450  end if
451  end do
452  densebnd = calcdens(denseref, drhodc, crhoref, ctemp)
453  else
454  ! -- neither of the above, so assign as denseref
455  densebnd = denseref
456  end if
457  !
458  ! -- Return
459  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_dimensions()

subroutine gwfbuymodule::read_dimensions ( class(gwfbuytype), intent(inout)  this)
private

Definition at line 1014 of file gwf-buy.f90.

1015  ! -- dummy
1016  class(GwfBuyType), intent(inout) :: this
1017  ! -- local
1018  character(len=LINELENGTH) :: errmsg, keyword
1019  integer(I4B) :: ierr
1020  logical :: isfound, endOfBlock
1021  ! -- format
1022  !
1023  ! -- get dimensions block
1024  call this%parser%GetBlock('DIMENSIONS', isfound, ierr, &
1025  supportopenclose=.true.)
1026  !
1027  ! -- parse dimensions block if detected
1028  if (isfound) then
1029  write (this%iout, '(/1x,a)') 'Processing BUY DIMENSIONS block'
1030  do
1031  call this%parser%GetNextLine(endofblock)
1032  if (endofblock) exit
1033  call this%parser%GetStringCaps(keyword)
1034  select case (keyword)
1035  case ('NRHOSPECIES')
1036  this%nrhospecies = this%parser%GetInteger()
1037  write (this%iout, '(4x,a,i0)') 'NRHOSPECIES = ', this%nrhospecies
1038  case default
1039  write (errmsg, '(a,a)') &
1040  'Unknown BUY dimension: ', trim(keyword)
1041  call store_error(errmsg)
1042  call this%parser%StoreErrorUnit()
1043  end select
1044  end do
1045  write (this%iout, '(1x,a)') 'End of BUY DIMENSIONS block'
1046  else
1047  call store_error('Required BUY DIMENSIONS block not found.')
1048  call this%parser%StoreErrorUnit()
1049  end if
1050  !
1051  ! -- check dimension
1052  if (this%nrhospecies < 1) then
1053  call store_error('NRHOSPECIES must be greater than zero.')
1054  call this%parser%StoreErrorUnit()
1055  end if
1056  !
1057  ! -- Return
1058  return
Here is the call graph for this function:

◆ read_options()

subroutine gwfbuymodule::read_options ( class(gwfbuytype this)
private

Definition at line 1479 of file gwf-buy.f90.

1480  ! -- modules
1481  use openspecmodule, only: access, form
1483  ! -- dummy
1484  class(GwfBuyType) :: this
1485  ! -- local
1486  character(len=LINELENGTH) :: errmsg, keyword
1487  character(len=MAXCHARLEN) :: fname
1488  integer(I4B) :: ierr
1489  logical :: isfound, endOfBlock
1490  ! -- formats
1491  character(len=*), parameter :: fmtfileout = &
1492  "(4x, 'BUY ', 1x, a, 1x, ' will be saved to file: ', &
1493  &a, /4x, 'opened on unit: ', I7)"
1494  !
1495  ! -- get options block
1496  call this%parser%GetBlock('OPTIONS', isfound, ierr, &
1497  supportopenclose=.true., blockrequired=.false.)
1498  !
1499  ! -- parse options block if detected
1500  if (isfound) then
1501  write (this%iout, '(1x,a)') 'Processing BUY OPTIONS block'
1502  do
1503  call this%parser%GetNextLine(endofblock)
1504  if (endofblock) exit
1505  call this%parser%GetStringCaps(keyword)
1506  select case (keyword)
1507  case ('HHFORMULATION_RHS')
1508  this%iform = 1
1509  this%iasym = 0
1510  write (this%iout, '(4x,a)') &
1511  'Hydraulic head formulation set to right-hand side'
1512  case ('DENSEREF')
1513  this%denseref = this%parser%GetDouble()
1514  write (this%iout, '(4x,a,1pg15.6)') &
1515  'Reference density has been set to: ', &
1516  this%denseref
1517  case ('DEV_EFH_FORMULATION')
1518  call this%parser%DevOpt()
1519  this%iform = 0
1520  this%iasym = 0
1521  write (this%iout, '(4x,a)') &
1522  'Formulation set to equivalent freshwater head'
1523  case ('DENSITY')
1524  call this%parser%GetStringCaps(keyword)
1525  if (keyword == 'FILEOUT') then
1526  call this%parser%GetString(fname)
1527  this%ioutdense = getunit()
1528  call openfile(this%ioutdense, this%iout, fname, 'DATA(BINARY)', &
1529  form, access, 'REPLACE')
1530  write (this%iout, fmtfileout) &
1531  'DENSITY', fname, this%ioutdense
1532  else
1533  errmsg = 'Optional density keyword must be '// &
1534  'followed by FILEOUT'
1535  call store_error(errmsg)
1536  end if
1537  case default
1538  write (errmsg, '(a,a)') 'Unknown BUY option: ', &
1539  trim(keyword)
1540  call store_error(errmsg)
1541  call this%parser%StoreErrorUnit()
1542  end select
1543  end do
1544  write (this%iout, '(1x,a)') 'End of BUY OPTIONS block'
1545  end if
1546  !
1547  ! -- Return
1548  return
subroutine, public urdaux(naux, inunit, iout, lloc, istart, istop, auxname, line, text)
Read auxiliary variables from an input line.
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
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Here is the call graph for this function:

◆ read_packagedata()

subroutine gwfbuymodule::read_packagedata ( class(gwfbuytype this)
private

Definition at line 1063 of file gwf-buy.f90.

1064  ! -- dummy
1065  class(GwfBuyType) :: this
1066  ! -- local
1067  character(len=LINELENGTH) :: errmsg
1068  character(len=LINELENGTH) :: line
1069  integer(I4B) :: ierr
1070  integer(I4B) :: irhospec
1071  logical :: isfound, endOfBlock
1072  logical :: blockrequired
1073  integer(I4B), dimension(:), allocatable :: itemp
1074  character(len=10) :: c10
1075  character(len=16) :: c16
1076  ! -- format
1077  character(len=*), parameter :: fmterr = &
1078  "('Invalid value for IRHOSPEC (',i0,') detected in BUY Package. &
1079  &IRHOSPEC must be > 0 and <= NRHOSPECIES, and duplicate values &
1080  &are not allowed.')"
1081  !
1082  ! -- initialize
1083  allocate (itemp(this%nrhospecies))
1084  itemp(:) = 0
1085  !
1086  ! -- get packagedata block
1087  blockrequired = .true.
1088  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
1089  blockrequired=blockrequired, &
1090  supportopenclose=.true.)
1091  !
1092  ! -- parse packagedata block
1093  if (isfound) then
1094  write (this%iout, '(1x,a)') 'Processing BUY PACKAGEDATA block'
1095  do
1096  call this%parser%GetNextLine(endofblock)
1097  if (endofblock) exit
1098  irhospec = this%parser%GetInteger()
1099  if (irhospec < 1 .or. irhospec > this%nrhospecies) then
1100  write (errmsg, fmterr) irhospec
1101  call store_error(errmsg)
1102  end if
1103  if (itemp(irhospec) /= 0) then
1104  write (errmsg, fmterr) irhospec
1105  call store_error(errmsg)
1106  end if
1107  itemp(irhospec) = 1
1108  this%drhodc(irhospec) = this%parser%GetDouble()
1109  this%crhoref(irhospec) = this%parser%GetDouble()
1110  call this%parser%GetStringCaps(this%cmodelname(irhospec))
1111  call this%parser%GetStringCaps(this%cauxspeciesname(irhospec))
1112  end do
1113  write (this%iout, '(1x,a)') 'End of BUY PACKAGEDATA block'
1114  else
1115  call store_error('Required BUY PACKAGEDATA block not found.')
1116  call this%parser%StoreErrorUnit()
1117  end if
1118  !
1119  ! -- Check for errors.
1120  if (count_errors() > 0) then
1121  call this%parser%StoreErrorUnit()
1122  end if
1123  !
1124  ! -- write packagedata information
1125  write (this%iout, '(/,a)') 'Summary of species information in BUY Package'
1126  write (this%iout, '(1a11, 4a17)') &
1127  'SPECIES', 'DRHODC', 'CRHOREF', 'MODEL', &
1128  'AUXSPECIESNAME'
1129  do irhospec = 1, this%nrhospecies
1130  write (c10, '(i0)') irhospec
1131  line = ' '//adjustr(c10)
1132  write (c16, '(g15.6)') this%drhodc(irhospec)
1133  line = trim(line)//' '//adjustr(c16)
1134  write (c16, '(g15.6)') this%crhoref(irhospec)
1135  line = trim(line)//' '//adjustr(c16)
1136  write (c16, '(a)') this%cmodelname(irhospec)
1137  line = trim(line)//' '//adjustr(c16)
1138  write (c16, '(a)') this%cauxspeciesname(irhospec)
1139  line = trim(line)//' '//adjustr(c16)
1140  write (this%iout, '(a)') trim(line)
1141  end do
1142  !
1143  ! -- deallocate
1144  deallocate (itemp)
1145  !
1146  ! -- Return
1147  return
Here is the call graph for this function:

◆ set_concentration_pointer()

subroutine gwfbuymodule::set_concentration_pointer ( class(gwfbuytype this,
character(len=lenmodelname), intent(in)  modelname,
real(dp), dimension(:), pointer  conc,
integer(i4b), dimension(:), pointer  icbund 
)
private

This routine is called from the gwfgwt exchange in the exg_ar() method

Definition at line 1577 of file gwf-buy.f90.

1578  ! -- dummy
1579  class(GwfBuyType) :: this
1580  character(len=LENMODELNAME), intent(in) :: modelname
1581  real(DP), dimension(:), pointer :: conc
1582  integer(I4B), dimension(:), pointer :: icbund
1583  ! -- local
1584  integer(I4B) :: i
1585  logical :: found
1586  !
1587  this%iconcset = 1
1588  found = .false.
1589  do i = 1, this%nrhospecies
1590  if (this%cmodelname(i) == modelname) then
1591  this%modelconc(i)%conc => conc
1592  this%modelconc(i)%icbund => icbund
1593  found = .true.
1594  exit
1595  end if
1596  end do
1597  !
1598  ! -- Return
1599  return

◆ set_options()

subroutine gwfbuymodule::set_options ( class(gwfbuytype this,
type(gwfbuyinputdatatype), intent(in)  input_data 
)
Parameters
[in]input_datathe input data to be set

Definition at line 1553 of file gwf-buy.f90.

1554  ! -- dummy
1555  class(GwfBuyType) :: this
1556  type(GwfBuyInputDataType), intent(in) :: input_data !< the input data to be set
1557  !
1558  this%iform = input_data%iform
1559  this%denseref = input_data%denseref
1560  !
1561  ! derived option:
1562  ! if not iform==2, there is no asymmetry
1563  if (this%iform == 0 .or. this%iform == 1) then
1564  this%iasym = 0
1565  end if
1566  !
1567  ! -- Return
1568  return

◆ set_packagedata()

subroutine gwfbuymodule::set_packagedata ( class(gwfbuytype this,
type(gwfbuyinputdatatype), intent(in)  input_data 
)
private
Parameters
thisthis buyoancy pkg
[in]input_datathe input data to be set

Definition at line 1152 of file gwf-buy.f90.

1153  ! -- dummy
1154  class(GwfBuyType) :: this !< this buyoancy pkg
1155  type(GwfBuyInputDataType), intent(in) :: input_data !< the input data to be set
1156  ! -- local
1157  integer(I4B) :: ispec
1158 
1159  do ispec = 1, this%nrhospecies
1160  this%drhodc(ispec) = input_data%drhodc(ispec)
1161  this%crhoref(ispec) = input_data%crhoref(ispec)
1162  this%cmodelname(ispec) = input_data%cmodelname(ispec)
1163  this%cauxspeciesname(ispec) = input_data%cauxspeciesname(ispec)
1164  end do
1165  !
1166  ! -- Return
1167  return