MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
gwfvscmodule Module Reference

Data Types

type  concentrationpointer
 
type  gwfvsctype
 

Functions/Subroutines

real(dp) function calc_visc (ivisc, viscref, dviscdc, cviscref, conc, a2, a3, a4)
 Generic function to calculate changes in fluid viscosity using a linear formulation. More...
 
subroutine, public vsc_cr (vscobj, name_model, input_mempath, inunit, iout)
 @ brief Create a new package object More...
 
subroutine vsc_df (this, dis, vsc_input)
 @ brief Define viscosity package options and dimensions More...
 
subroutine vsc_ar (this, ibound)
 @ brief Allocate and read method for viscosity package More...
 
subroutine vsc_ar_bnd (this, packobj)
 Activate viscosity in advanced packages. More...
 
subroutine set_npf_pointers (this)
 Set pointers to NPF variables. More...
 
subroutine vsc_rp (this)
 @ brief Read new period data in viscosity package More...
 
subroutine vsc_ad (this)
 @ brief Advance the viscosity package More...
 
subroutine vsc_ad_bnd (this, packobj, hnew)
 Advance the boundary packages when viscosity is active. More...
 
subroutine vsc_ad_standard_bnd (packobj, hnew, visc, viscref, locelev, locvisc, locconc, dviscdc, cviscref, ivisc, a2, a3, a4, ctemp)
 advance ghb while accounting for viscosity More...
 
subroutine vsc_ad_sfr (packobj, visc, viscref, elev, locvisc, locconc, dviscdc, cviscref, ivisc, a2, a3, a4, ctemp)
 Update sfr-related viscosity ratios. More...
 
subroutine vsc_ad_lak (packobj, visc, viscref, elev, locvisc, locconc, dviscdc, cviscref, ivisc, a2, a3, a4, ctemp)
 Update lak-related viscosity ratios. More...
 
subroutine vsc_ad_maw (packobj, visc, viscref, elev, locvisc, locconc, dviscdc, cviscref, ivisc, a2, a3, a4, ctemp)
 Update maw-related viscosity ratios. More...
 
real(dp) function update_bnd_cond (bndvisc, viscref, spcfdcond)
 Apply viscosity to the conductance term. More...
 
real(dp) function calc_vsc_ratio (viscref, bndvisc)
 calculate and return the viscosity ratio More...
 
real(dp) function calc_bnd_viscosity (n, locvisc, locconc, viscref, dviscdc, cviscref, ctemp, ivisc, a2, a3, a4, auxvar)
 @ brief Calculate the boundary viscosity More...
 
subroutine get_visc_ratio (this, n, m, gwhdn, gwhdm, viscratio)
 Calculate the viscosity ratio. More...
 
subroutine calc_q_visc (this, cellid, viscratio)
 Account for viscosity in the aquiferhorizontal flow barriers. More...
 
subroutine update_k_with_vsc (this)
 Appled the viscosity ratio (mu_o/mu) to the hydraulic conductivity. More...
 
subroutine vsc_set_changed_at (this, kper, kstp)
 Mark K changes as having occurred at (kper, kstp) More...
 
subroutine vsc_ot_dv (this, idvfl)
 @ brief Output viscosity package dependent-variable terms. More...
 
subroutine vsc_da (this)
 @ brief Deallocate viscosity package memory More...
 
subroutine source_dimensions (this)
 @ brief Source dimensions for package More...
 
subroutine source_packagedata (this)
 @ brief source packagedata for package More...
 
subroutine set_packagedata (this, input_data)
 Sets package data instead of reading from file. More...
 
subroutine vsc_calcvisc (this)
 Calculate fluid viscosity. More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine allocate_arrays (this, nodes)
 @ brief Allocate arrays More...
 
subroutine source_options (this)
 @ brief Source VSC options More...
 
subroutine log_options (this, found, viscosityfile)
 @ brief log VSC 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, istmpr)
 @ brief Set pointers to concentration(s) More...
 

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwfvscmodule::allocate_arrays ( class(gwfvsctype this,
integer(i4b), intent(in)  nodes 
)

Allocate and initialize arrays for the VSC package.

Definition at line 1167 of file gwf-vsc.f90.

1168  ! -- dummy
1169  class(GwfVscType) :: this
1170  integer(I4B), intent(in) :: nodes
1171  ! -- local
1172  integer(I4B) :: i
1173  !
1174  ! -- Allocate
1175  call mem_allocate(this%visc, nodes, 'VISC', this%memoryPath)
1176  call mem_allocate(this%ivisc, this%nviscspecies, 'IVISC', this%memoryPath)
1177  call mem_allocate(this%dviscdc, this%nviscspecies, 'DRHODC', &
1178  this%memoryPath)
1179  call mem_allocate(this%cviscref, this%nviscspecies, 'CRHOREF', &
1180  this%memoryPath)
1181  call mem_allocate(this%ctemp, this%nviscspecies, 'CTEMP', this%memoryPath)
1182  allocate (this%cmodelname(this%nviscspecies))
1183  allocate (this%cauxspeciesname(this%nviscspecies))
1184  allocate (this%modelconc(this%nviscspecies))
1185  !
1186  ! -- Initialize
1187  do i = 1, nodes
1188  this%visc(i) = this%viscref
1189  end do
1190  !
1191  ! -- Initialize nviscspecies arrays
1192  do i = 1, this%nviscspecies
1193  this%ivisc(i) = 1
1194  this%dviscdc(i) = dzero
1195  this%cviscref(i) = dzero
1196  this%ctemp(i) = dzero
1197  this%cmodelname(i) = ''
1198  this%cauxspeciesname(i) = ''
1199  end do

◆ allocate_scalars()

subroutine gwfvscmodule::allocate_scalars ( class(gwfvsctype this)
private

Allocate and initialize scalars for the VSC package. The base model allocate scalars method is also called.

Definition at line 1127 of file gwf-vsc.f90.

1128  ! -- modules
1129  use constantsmodule, only: dzero, dten, dep3
1130  ! -- dummy
1131  class(GwfVscType) :: this
1132  !
1133  ! -- allocate scalars in NumericalPackageType
1134  call this%NumericalPackageType%allocate_scalars()
1135  !
1136  ! -- Allocate
1137  call mem_allocate(this%thermivisc, 'THERMIVISC', this%memoryPath)
1138  call mem_allocate(this%idxtmpr, 'IDXTMPR', this%memoryPath)
1139  call mem_allocate(this%ioutvisc, 'IOUTVISC', this%memoryPath)
1140  call mem_allocate(this%ireadelev, 'IREADELEV', this%memoryPath)
1141  call mem_allocate(this%iconcset, 'ICONCSET', this%memoryPath)
1142  call mem_allocate(this%viscref, 'VISCREF', this%memoryPath)
1143  call mem_allocate(this%a2, 'A2', this%memoryPath)
1144  call mem_allocate(this%a3, 'A3', this%memoryPath)
1145  call mem_allocate(this%a4, 'A4', this%memoryPath)
1146  !
1147  call mem_allocate(this%nviscspecies, 'NVISCSPECIES', this%memoryPath)
1148  !
1149  ! -- Initialize
1150  this%thermivisc = 0
1151  this%idxtmpr = 0
1152  this%ioutvisc = 0
1153  this%ireadelev = 0
1154  this%iconcset = 0
1155  this%viscref = dep3
1156  this%A2 = dten
1157  this%A3 = 248.37_dp
1158  this%A4 = 133.15_dp
1159  !
1160  this%nviscspecies = 0
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dep3
real constant 1000
Definition: Constants.f90:88
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dten
real constant 10
Definition: Constants.f90:84

◆ calc_bnd_viscosity()

real(dp) function gwfvscmodule::calc_bnd_viscosity ( integer(i4b), intent(in)  n,
integer(i4b), intent(in)  locvisc,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), intent(in)  viscref,
real(dp), dimension(:), intent(in)  dviscdc,
real(dp), dimension(:), intent(in)  cviscref,
real(dp), dimension(:), intent(inout)  ctemp,
integer(i4b), dimension(:), intent(in)  ivisc,
real(dp), intent(in)  a2,
real(dp), intent(in)  a3,
real(dp), intent(in)  a4,
real(dp), dimension(:, :), intent(in)  auxvar 
)
private

Return the viscosity of the boundary package using one of the options in the following order of priority:

  1. Assign as aux variable in column with name 'VISCOSITY'
  2. Calculate using viscosity equation and nviscspecies aux columns
  3. If neither of those, then assign as viscref !!

Definition at line 730 of file gwf-vsc.f90.

732  ! -- modules
733  ! -- dummy
734  integer(I4B), intent(in) :: n
735  integer(I4B), intent(in) :: locvisc
736  real(DP), intent(in) :: a2, a3, a4
737  integer(I4B), dimension(:), intent(in) :: ivisc
738  integer(I4B), dimension(:), intent(in) :: locconc
739  real(DP), intent(in) :: viscref
740  real(DP), dimension(:), intent(in) :: dviscdc
741  real(DP), dimension(:), intent(in) :: cviscref
742  real(DP), dimension(:), intent(inout) :: ctemp
743  real(DP), dimension(:, :), intent(in) :: auxvar
744  ! -- Return
745  real(DP) :: viscbnd
746  ! -- local
747  integer(I4B) :: i
748  !
749  ! -- assign boundary viscosity based on one of three options
750  if (locvisc > 0) then
751  ! -- assign viscosity to an aux column named 'VISCOSITY'
752  viscbnd = auxvar(locvisc, n)
753  else if (locconc(1) > 0) then
754  ! -- calculate viscosity using one or more concentration auxcolumns
755  do i = 1, size(locconc)
756  ctemp(i) = dzero
757  if (locconc(i) > 0) then
758  ctemp(i) = auxvar(locconc(i), n)
759  end if
760  end do
761  viscbnd = calc_visc(ivisc, viscref, dviscdc, cviscref, ctemp, a2, a3, a4)
762  else
763  ! -- neither of the above, so assign as viscref
764  viscbnd = viscref
765  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_q_visc()

subroutine gwfvscmodule::calc_q_visc ( class(gwfvsctype this,
integer(i4b), intent(in)  cellid,
real(dp), intent(inout)  viscratio 
)

Will return the viscosity associated with the upgradient node (cell) to the HFB package for adjusting the hydraulic characteristic (hydchr) of the barrier

Definition at line 801 of file gwf-vsc.f90.

802  ! -- dummy variables
803  class(GwfVscType) :: this
804  integer(I4B), intent(in) :: cellid
805  ! -- Return
806  real(DP), intent(inout) :: viscratio
807  ! -- local
808  real(DP) :: visc
809  !
810  ! -- Retrieve viscosity for the passed node number
811  visc = this%visc(cellid)
812  !
813  ! -- Calculate the viscosity ratio for the
814  viscratio = calc_vsc_ratio(this%viscref, visc)
Here is the call graph for this function:

◆ calc_visc()

real(dp) function gwfvscmodule::calc_visc ( integer(i4b), dimension(:), intent(in)  ivisc,
real(dp), intent(in)  viscref,
real(dp), dimension(:), intent(in)  dviscdc,
real(dp), dimension(:), intent(in)  cviscref,
real(dp), dimension(:), intent(in)  conc,
real(dp), intent(in)  a2,
real(dp), intent(in)  a3,
real(dp), intent(in)  a4 
)
private

Definition at line 100 of file gwf-vsc.f90.

102  ! -- dummy
103  integer(I4B), dimension(:), intent(in) :: ivisc
104  real(DP), intent(in) :: viscref
105  real(DP), dimension(:), intent(in) :: dviscdc
106  real(DP), dimension(:), intent(in) :: cviscref
107  real(DP), dimension(:), intent(in) :: conc
108  real(DP), intent(in) :: a2, a3, a4
109  ! -- return
110  real(DP) :: visc
111  ! -- local
112  integer(I4B) :: nviscspec
113  integer(I4B) :: i
114  real(DP) :: mu_t
115  real(DP) :: expon
116  !
117  nviscspec = size(dviscdc)
118  visc = viscref
119  !
120  do i = 1, nviscspec
121  if (ivisc(i) == 1) then
122  visc = visc + dviscdc(i) * (conc(i) - cviscref(i))
123  else
124  expon = -1 * a3 * ((conc(i) - cviscref(i)) / &
125  ((conc(i) + a4) * (cviscref(i) + a4)))
126  mu_t = viscref * a2**expon
127  ! Order matters!! (This assumes we apply the temperature correction after
128  ! accounting for solute concentrations)
129  ! If a nonlinear correction is applied, then b/c it takes into account
130  ! viscref, need to subtract it in this case
131  ! At most, there will only ever be 1 nonlinear correction
132  visc = (visc - viscref) + mu_t
133  end if
134  ! end if
135  end do
Here is the caller graph for this function:

◆ calc_vsc_ratio()

real(dp) function gwfvscmodule::calc_vsc_ratio ( real(dp), intent(in)  viscref,
real(dp), intent(in)  bndvisc 
)
private

Definition at line 712 of file gwf-vsc.f90.

713  ! -- dummy
714  real(DP), intent(in) :: viscref
715  real(DP), intent(in) :: bndvisc
716  ! -- local
717  real(DP) :: viscratio
718  !
719  viscratio = viscref / bndvisc
Here is the caller graph for this function:

◆ get_visc_ratio()

subroutine gwfvscmodule::get_visc_ratio ( class(gwfvsctype this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
real(dp), intent(in)  gwhdn,
real(dp), intent(in)  gwhdm,
real(dp), intent(inout)  viscratio 
)
private

Calculate the viscosity ratio applied to the hydraulic characteristic provided by the user. The viscosity ratio is applicable only when the hydraulic characteristic is specified as positive and will not be applied when the hydchr is negative

Definition at line 775 of file gwf-vsc.f90.

776  ! -- modules
777  use constantsmodule, only: done
778  ! -- dummy
779  class(GwfVscType) :: this
780  integer(I4B), intent(in) :: n, m
781  real(DP), intent(in) :: gwhdn, gwhdm
782  real(DP), intent(inout) :: viscratio
783  ! -- loca
784  integer(I4B) :: cellid
785  !
786  viscratio = done
787  if (gwhdm > gwhdn) then
788  cellid = m
789  else if (gwhdn >= gwhdm) then
790  cellid = n
791  end if
792  call this%calc_q_visc(cellid, viscratio)
real(dp), parameter done
real constant 1
Definition: Constants.f90:76

◆ log_options()

subroutine gwfvscmodule::log_options ( class(gwfvsctype), intent(inout)  this,
type(gwfvscparamfoundtype), intent(in)  found,
character(len=*), intent(in)  viscosityfile 
)

Definition at line 1282 of file gwf-vsc.f90.

1283  ! -- modules
1285  ! -- dummy variables
1286  class(GwfVscType), intent(inout) :: this
1287  type(GwfVscParamFoundType), intent(in) :: found
1288  character(len=*), intent(in) :: viscosityfile
1289  ! -- local variables
1290  ! -- formats
1291  character(len=*), parameter :: fmtfileout = &
1292  "(4x, 'VSC', 1x, a, 1x, 'Will be saved to file: ', &
1293  &a, ' opened on unit: ', I7)"
1294  character(len=*), parameter :: fmtlinear = &
1295  "(4x,'Viscosity will vary linearly with temperature &
1296  &change ')"
1297  character(len=*), parameter :: fmtnonlinear = &
1298  "(4x,'Viscosity will vary non-linearly with temperature &
1299  &change ')"
1300 
1301  write (this%iout, '(1x,a)') 'Processing VSC OPTIONS block'
1302 
1303  if (found%viscref) then
1304  write (this%iout, '(4x,a,1pg15.6)') &
1305  'Reference viscosity has been set to: ', this%viscref
1306  end if
1307  if (found%viscosityfile) then
1308  write (this%iout, fmtfileout) &
1309  'VISCOSITY', trim(viscosityfile), this%ioutvisc
1310  end if
1311  if (found%temp_specname) then
1312  write (this%iout, '(4x, a)') 'Temperature species name set to: '// &
1313  trim(this%name_temp_spec)
1314  end if
1315  if (found%thermal_form) then
1316  select case (this%thermivisc)
1317  case (1)
1318  write (this%iout, fmtlinear)
1319  case (2)
1320  write (this%iout, fmtnonlinear)
1321  end select
1322  end if
1323  if (found%thermal_a2) then
1324  if (this%thermivisc == 2) then
1325  write (this%iout, '(4x,a,1pg15.6)') &
1326  'A2 in nonlinear viscosity formulation has been set to: ', &
1327  this%a2
1328  else
1329  write (this%iout, '(4x,a,/,4x,a,/,4x,a)') 'THERMAL_A2 specified by user &
1330  &in VSC Package input file. LINEAR viscosity ', 'formulation also &
1331  &specified. THERMAL_A2 will not affect ', 'viscosity calculations.'
1332  end if
1333  end if
1334  if (found%thermal_a3) then
1335  if (this%thermivisc == 2) then
1336  write (this%iout, '(4x,a,1pg15.6)') &
1337  'A3 in nonlinear viscosity formulation has been set to: ', &
1338  this%a3
1339  else
1340  write (this%iout, '(4x,a,/,4x,a,/,4x,a)') 'THERMAL_A3 specified by user &
1341  &in VSC Package input file. LINEAR viscosity ', 'formulation also &
1342  &specified. THERMAL_A3 will not affect ', 'viscosity calculations.'
1343  end if
1344  end if
1345  if (found%thermal_a4) then
1346  if (this%thermivisc == 2) then
1347  write (this%iout, '(4x,a,1pg15.6)') &
1348  'A4 in nonlinear viscosity formulation has been set to: ', &
1349  this%a4
1350  else
1351  write (this%iout, '(4x,a,/,4x,a,/,4x,a)') 'THERMAL_A4 specified by user &
1352  &in VSC Package input file. LINEAR viscosity ', 'formulation also &
1353  &specified. THERMAL_A4 will not affect ', 'viscosity calculations.'
1354  end if
1355  end if
1356 
1357  write (this%iout, '(1x,a)') 'end of VSC options block'

◆ set_concentration_pointer()

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

Pass in a gwt model name, concentration array, and ibound, and store a pointer to these in the VSC package so that viscosity can be calculated from them. This routine is called from the gwfgwt exchange in the exg_ar() method.

Definition at line 1377 of file gwf-vsc.f90.

1378  ! -- dummy
1379  class(GwfVscType) :: this
1380  character(len=LENMODELNAME), intent(in) :: modelname
1381  real(DP), dimension(:), pointer :: conc
1382  integer(I4B), dimension(:), pointer :: icbund
1383  integer(I4B), optional, intent(in) :: istmpr
1384  ! -- local
1385  integer(I4B) :: i
1386  logical :: found
1387  !
1388  this%iconcset = 1
1389  found = .false.
1390  do i = 1, this%nviscspecies
1391  if (this%cmodelname(i) == modelname) then
1392  this%modelconc(i)%conc => conc
1393  this%modelconc(i)%icbund => icbund
1394  found = .true.
1395  exit
1396  end if
1397  end do

◆ set_npf_pointers()

subroutine gwfvscmodule::set_npf_pointers ( class(gwfvsctype this)

Set array and variable pointers from the NPF package for access by VSC.

Definition at line 301 of file gwf-vsc.f90.

302  ! -- dummy
303  class(GwfVscType) :: this
304  ! -- local
305  character(len=LENMEMPATH) :: npfMemoryPath
306  !
307  ! -- Set pointers to other package variables
308  ! -- NPF
309  npfmemorypath = create_mem_path(this%name_model, 'NPF')
310  call mem_setptr(this%k11, 'K11', npfmemorypath)
311  call mem_setptr(this%k22, 'K22', npfmemorypath)
312  call mem_setptr(this%k33, 'K33', npfmemorypath)
313  call mem_setptr(this%k11input, 'K11INPUT', npfmemorypath)
314  call mem_setptr(this%k22input, 'K22INPUT', npfmemorypath)
315  call mem_setptr(this%k33input, 'K33INPUT', npfmemorypath)
316  call mem_setptr(this%kchangeper, 'KCHANGEPER', npfmemorypath)
317  call mem_setptr(this%kchangestp, 'KCHANGESTP', npfmemorypath)
318  call mem_setptr(this%nodekchange, 'NODEKCHANGE', npfmemorypath)
Here is the call graph for this function:

◆ set_options()

subroutine gwfvscmodule::set_options ( class(gwfvsctype this,
type(gwfvscinputdatatype), intent(in)  input_data 
)
Parameters
[in]input_datathe input data to be set

Definition at line 1362 of file gwf-vsc.f90.

1363  ! -- dummy
1364  class(GwfVscType) :: this
1365  type(GwfVscInputDataType), intent(in) :: input_data !< the input data to be set
1366  !
1367  this%viscref = input_data%viscref

◆ set_packagedata()

subroutine gwfvscmodule::set_packagedata ( class(gwfvsctype this,
type(gwfvscinputdatatype), intent(in)  input_data 
)
Parameters
thisthis vscoancy pkg
[in]input_datathe input data to be set

Definition at line 1078 of file gwf-vsc.f90.

1079  ! -- dummy
1080  class(GwfVscType) :: this !< this vscoancy pkg
1081  type(GwfVscInputDataType), intent(in) :: input_data !< the input data to be set
1082  ! -- local
1083  integer(I4B) :: ispec
1084  !
1085  do ispec = 1, this%nviscspecies
1086  this%dviscdc(ispec) = input_data%dviscdc(ispec)
1087  this%cviscref(ispec) = input_data%cviscref(ispec)
1088  this%cmodelname(ispec) = input_data%cmodelname(ispec)
1089  this%cauxspeciesname(ispec) = input_data%cauxspeciesname(ispec)
1090  end do

◆ source_dimensions()

subroutine gwfvscmodule::source_dimensions ( class(gwfvsctype), intent(inout)  this)
private

Definition at line 949 of file gwf-vsc.f90.

950  ! -- modules
953  ! -- dummy
954  class(GwfVscType), intent(inout) :: this
955  ! -- local
956  type(GwfVscParamFoundType) :: found
957 
958  ! update defaults from input context
959  call mem_set_value(this%nviscspecies, 'NVISCSPECIES', this%input_mempath, &
960  found%nviscspecies)
961 
962  ! log dimensions
963  write (this%iout, '(/1x,a)') 'Processing VSC DIMENSIONS block'
964  write (this%iout, '(4x,a,i0)') 'NVISCSPECIES = ', this%nviscspecies
965  write (this%iout, '(1x,a)') 'End of VSC DIMENSIONS block'
966 
967  ! check dimension
968  if (this%nviscspecies < 1) then
969  call store_error('NVISCSPECIES must be greater than zero.')
970  call store_error_filename(this%input_fname)
971  end if
Here is the call graph for this function:

◆ source_options()

subroutine gwfvscmodule::source_options ( class(gwfvsctype), intent(inout)  this)
private

Definition at line 1204 of file gwf-vsc.f90.

1205  ! -- modules
1206  use constantsmodule, only: lenvarname
1208  use openspecmodule, only: access, form
1209  use inputoutputmodule, only: getunit, openfile
1211  ! -- dummy
1212  class(GwfVscType), intent(inout) :: this
1213  ! -- local
1214  character(len=LENVARNAME), dimension(2) :: thermal_form = &
1215  &[character(len=LENVARNAME) :: 'LINEAR', 'NONLINEAR']
1216  character(len=linelength) :: viscosityfile
1217  type(GwfVscParamFoundType) :: found
1218 
1219  ! initialize variables
1220  viscosityfile = ''
1221 
1222  ! update defaults from input context
1223  call mem_set_value(this%viscref, 'VISCREF', this%input_mempath, &
1224  found%viscref)
1225  call mem_set_value(viscosityfile, 'VISCOSITYFILE', this%input_mempath, &
1226  found%viscosityfile)
1227  call mem_set_value(this%name_temp_spec, 'TEMP_SPECNAME', this%input_mempath, &
1228  found%temp_specname)
1229  call mem_set_value(this%thermivisc, 'THERMAL_FORM', this%input_mempath, &
1230  thermal_form, found%thermal_form)
1231  call mem_set_value(this%a2, 'THERMAL_A2', this%input_mempath, &
1232  found%thermal_a2)
1233  call mem_set_value(this%a3, 'THERMAL_A3', this%input_mempath, &
1234  found%thermal_a3)
1235  call mem_set_value(this%a4, 'THERMAL_A4', this%input_mempath, &
1236  found%thermal_a4)
1237 
1238  ! fileout options
1239  if (found%viscosityfile) then
1240  this%ioutvisc = getunit()
1241  call openfile(this%ioutvisc, this%iout, viscosityfile, 'DATA(BINARY)', &
1242  form, access, 'REPLACE')
1243  end if
1244 
1245  ! set warnings and errors
1246  if (this%thermivisc == 1) then
1247  if (this%a2 == 0.0) then
1248  write (errmsg, '(a)') 'LINEAR option selected for varying &
1249  &viscosity with temperature, but A1, a surrogate for &
1250  &dVISC/dT, set equal to 0.0'
1251  call store_error(errmsg)
1252  call store_error_filename(this%input_fname)
1253  end if
1254  end if
1255  if (this%thermivisc > 1) then
1256  if (this%a2 == 0) then
1257  write (warnmsg, '(a)') 'NONLINEAR option selected for &
1258  &varying viscosity with temperature, but A2 set equal to &
1259  &zero which may lead to unintended values for viscosity'
1260  call store_warning(warnmsg)
1261  end if
1262  if (this%a3 == 0) then
1263  write (warnmsg, '(a)') 'NONLINEAR option selected for &
1264  &varying viscosity with temperature,, but A3 set equal to &
1265  &zero which may lead to unintended values for viscosity'
1266  call store_warning(warnmsg)
1267  end if
1268  if (this%a4 == 0) then
1269  write (warnmsg, '(a)') 'NONLINEAR option selected for &
1270  &varying viscosity with temperature, BUT A4 SET EQUAL TO &
1271  &zero which may lead to unintended values for viscosity'
1272  call store_warning(warnmsg)
1273  end if
1274  end if
1275 
1276  ! log options
1277  call this%log_options(found, viscosityfile)
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
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
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Here is the call graph for this function:

◆ source_packagedata()

subroutine gwfvscmodule::source_packagedata ( class(gwfvsctype), intent(inout)  this)

Definition at line 976 of file gwf-vsc.f90.

977  ! -- modules
980  ! -- dummy
981  class(GwfVscType), intent(inout) :: this
982  integer(I4B), dimension(:), pointer, contiguous :: iviscspec
983  type(CharacterStringType), dimension(:), pointer, &
984  contiguous :: modelnames, auxspeciesnames
985  real(DP), dimension(:), pointer, contiguous :: dviscdc, cviscref
986  integer(I4B), dimension(:), allocatable :: itemp
987  character(len=LINELENGTH) :: modelname, auxspeciesname, line
988  character(len=10) :: c10
989  character(len=16) :: c16
990  integer(I4B) :: n
991  ! -- format
992  character(len=*), parameter :: fmterr = &
993  "('Invalid value for IRHOSPEC (',i0,') detected in VSC Package. &
994  &IRHOSPEC must be > 0 and <= NVISCSPECIES, and duplicate values &
995  &are not allowed.')"
996 
997  ! initialize
998  allocate (itemp(this%nviscspecies))
999  itemp(:) = 0
1000 
1001  ! set input context pointers
1002  call mem_setptr(iviscspec, 'IVISCSPEC', this%input_mempath)
1003  call mem_setptr(dviscdc, 'DVISCDC', this%input_mempath)
1004  call mem_setptr(cviscref, 'CVISCREF', this%input_mempath)
1005  call mem_setptr(modelnames, 'MODELNAME', this%input_mempath)
1006  call mem_setptr(auxspeciesnames, 'AUXSPECIESNAME', this%input_mempath)
1007 
1008  ! process package data
1009  do n = 1, size(iviscspec)
1010  modelname = modelnames(n)
1011  auxspeciesname = auxspeciesnames(n)
1012 
1013  if (iviscspec(n) < 1 .or. iviscspec(n) > this%nviscspecies) then
1014  write (errmsg, fmterr) iviscspec(n)
1015  call store_error(errmsg)
1016  end if
1017  if (itemp(iviscspec(n)) /= 0) then
1018  write (errmsg, fmterr) iviscspec(n)
1019  call store_error(errmsg)
1020  end if
1021  itemp(iviscspec(n)) = 1
1022  !
1023  this%dviscdc(iviscspec(n)) = dviscdc(n)
1024  this%cviscref(iviscspec(n)) = cviscref(n)
1025  this%cmodelname(iviscspec(n)) = trim(modelname)
1026  this%cauxspeciesname(iviscspec(n)) = trim(auxspeciesname)
1027  !
1028  if (auxspeciesname == this%name_temp_spec) then
1029  if (this%idxtmpr > 0) then
1030  write (errmsg, '(a)') 'More than one species in VSC input identified &
1031  &as '//trim(this%name_temp_spec)//'. Only one species may be &
1032  &designated to represent temperature.'
1033  call store_error(errmsg)
1034  else
1035  this%idxtmpr = iviscspec(n)
1036  if (this%thermivisc == 2) then
1037  this%ivisc(iviscspec(n)) = 2
1038  end if
1039  end if
1040  end if
1041  end do
1042 
1043  ! Check for errors.
1044  if (count_errors() > 0) then
1045  call store_error_filename(this%input_fname)
1046  end if
1047 
1048  ! log package data
1049  write (this%iout, '(/,1x,a)') 'Processing VSC PACKAGEDATA block'
1050 
1051  ! write packagedata information
1052  write (this%iout, '(1x,a)') 'Summary of species information in VSC Package'
1053  write (this%iout, '(1a11,5a17)') &
1054  'Species', 'DVISCDC', 'CVISCREF', 'Model', 'AUXSPECIESNAME'
1055  do n = 1, this%nviscspecies
1056  write (c10, '(i0)') n
1057  line = ' '//adjustr(c10)
1058 
1059  write (c16, '(g15.6)') this%dviscdc(n)
1060  line = trim(line)//' '//adjustr(c16)
1061  write (c16, '(g15.6)') this%cviscref(n)
1062  line = trim(line)//' '//adjustr(c16)
1063  write (c16, '(a)') this%cmodelname(n)
1064  line = trim(line)//' '//adjustr(c16)
1065  write (c16, '(a)') this%cauxspeciesname(n)
1066  line = trim(line)//' '//adjustr(c16)
1067  write (this%iout, '(a)') trim(line)
1068  end do
1069 
1070  write (this%iout, '(1x,a)') 'End of VSC PACKAGEDATA block'
1071 
1072  ! deallocate
1073  deallocate (itemp)
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function:

◆ update_bnd_cond()

real(dp) function gwfvscmodule::update_bnd_cond ( real(dp), intent(in)  bndvisc,
real(dp), intent(in)  viscref,
real(dp), intent(in)  spcfdcond 
)

When the viscosity package is active apply the viscosity ratio to the active boundary package's conductance term.

Definition at line 695 of file gwf-vsc.f90.

696  ! -- dummy
697  real(DP), intent(in) :: viscref
698  real(DP), intent(in) :: bndvisc
699  real(DP), intent(in) :: spcfdcond
700  ! -- local
701  real(DP) :: vscratio
702  real(DP) :: updatedcond
703  !
704  vscratio = calc_vsc_ratio(viscref, bndvisc)
705  !
706  ! -- calculate new conductance here
707  updatedcond = vscratio * spcfdcond
Here is the call graph for this function:
Here is the caller graph for this function:

◆ update_k_with_vsc()

subroutine gwfvscmodule::update_k_with_vsc ( class(gwfvsctype this)
private

This routine called after updating the viscosity values using the latest concentration and/or temperature values. The ratio mu_o/mu, reference viscosity divided by the updated viscosity value, is multiplied by K for each cell.

Definition at line 824 of file gwf-vsc.f90.

825  ! -- dummy
826  class(GwfVscType) :: this
827  ! -- local
828  integer(I4B) :: n
829  real(DP) :: viscratio
830  !
831  ! -- For viscosity-based K's, apply change of K to K11 by starting with
832  ! user-specified K values and not the K's leftover from the last viscosity
833  ! update.
834  do n = 1, this%dis%nodes
835  call this%calc_q_visc(n, viscratio)
836  this%k11(n) = this%k11input(n) * viscratio
837  this%k22(n) = this%k22input(n) * viscratio
838  this%k33(n) = this%k33input(n) * viscratio
839  this%nodekchange(n) = 1
840  end do
841  !
842  ! -- Flag kchange
843  call this%vsc_set_changed_at(kper, kstp)

◆ vsc_ad()

subroutine gwfvscmodule::vsc_ad ( class(gwfvsctype this)

Advance data in the VSC package. The method sets or advances time series, time array series, and observation data.

Definition at line 358 of file gwf-vsc.f90.

359  ! -- dummy
360  class(GwfVscType) :: this
361  !
362  ! -- update viscosity using the latest concentration/temperature
363  call this%vsc_calcvisc()

◆ vsc_ad_bnd()

subroutine gwfvscmodule::vsc_ad_bnd ( class(gwfvsctype this,
class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew 
)
private

Update the conductance values associate with inflow from a boundary when VSC package is active.

Definition at line 371 of file gwf-vsc.f90.

372  ! -- modules
373  use bndmodule, only: bndtype
374  ! -- dummy
375  class(GwfVscType) :: this
376  class(BndType), pointer :: packobj
377  real(DP), intent(in), dimension(:) :: hnew
378  ! -- local
379  integer(I4B) :: i, j
380  integer(I4B) :: n, locvisc, locelev
381  integer(I4B), dimension(:), allocatable :: locconc
382  !
383  ! -- initialize
384  locvisc = 0
385  locelev = 0
386  allocate (locconc(this%nviscspecies))
387  locconc(:) = 0
388  !
389  ! -- Add viscosity terms for conductance-dependent boundaries
390  do n = 1, packobj%naux
391  if (packobj%auxname(n) == 'VISCOSITY') then
392  locvisc = n
393  else if (packobj%auxname(n) == 'ELEVATION') then
394  locelev = n
395  end if
396  end do
397  !
398  ! -- find aux columns for conc (or temp.) that affect viscosity
399  do i = 1, this%nviscspecies
400  locconc(i) = 0
401  do j = 1, packobj%naux
402  if (this%cauxspeciesname(i) == packobj%auxname(j)) then
403  locconc(i) = j
404  exit
405  end if
406  end do
407  if (locconc(i) == 0) then
408  ! -- one not found, so don't use and mark all as 0
409  locconc(:) = 0
410  exit
411  end if
412  end do
413  !
414  ! -- apply viscosity terms to inflow from boundary based on package type
415  select case (packobj%filtyp)
416  case ('GHB', 'DRN', 'RIV')
417  !
418  ! -- general head, drain, and river boundary
419  call vsc_ad_standard_bnd(packobj, hnew, this%visc, this%viscref, &
420  locelev, locvisc, locconc, this%dviscdc, &
421  this%cviscref, this%ivisc, this%a2, this%a3, &
422  this%a4, this%ctemp)
423  case ('LAK')
424  !
425  ! -- lake
426  ! Update 'viscratios' internal to lak such that they are
427  ! automatically applied in the LAK calc_cond() routine
428  call vsc_ad_lak(packobj, this%visc, this%viscref, this%elev, locvisc, &
429  locconc, this%dviscdc, this%cviscref, this%ivisc, &
430  this%a2, this%a3, this%a4, this%ctemp)
431  case ('SFR')
432  !
433  ! -- streamflow routing
434  ! Update 'viscratios' internal to sfr such that they are
435  ! automatically applied in the SFR calc_cond() routine
436  call vsc_ad_sfr(packobj, this%visc, this%viscref, this%elev, locvisc, &
437  locconc, this%dviscdc, this%cviscref, this%ivisc, &
438  this%a2, this%a3, this%a4, this%ctemp)
439  case ('MAW')
440  !
441  ! -- multi-aquifer well
442  call vsc_ad_maw(packobj, this%visc, this%viscref, this%elev, locvisc, &
443  locconc, this%dviscdc, this%cviscref, this%ivisc, &
444  this%a2, this%a3, this%a4, this%ctemp)
445  case ('UZF')
446  !
447  ! -- unsaturated-zone flow
448  case default
449  !
450  ! -- nothing
451  end select
452  !
453  ! -- deallocate
454  deallocate (locconc)
This module contains the base boundary package.
@ brief BndType
Here is the call graph for this function:

◆ vsc_ad_lak()

subroutine gwfvscmodule::vsc_ad_lak ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  visc,
real(dp), intent(in)  viscref,
real(dp), dimension(:), intent(in)  elev,
integer(i4b), intent(in)  locvisc,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  dviscdc,
real(dp), dimension(:), intent(in)  cviscref,
integer(i4b), dimension(:), intent(in)  ivisc,
real(dp), intent(in)  a2,
real(dp), intent(in)  a3,
real(dp), intent(in)  a4,
real(dp), dimension(:), intent(inout)  ctemp 
)

When the viscosity package is active, update the viscosity ratio that is applied to the lakebed conductance calculated in the LAK package

Definition at line 587 of file gwf-vsc.f90.

589  ! -- modules
590  use bndmodule, only: bndtype
591  use lakmodule, only: laktype
592  class(BndType), pointer :: packobj
593  ! -- dummy
594  real(DP), intent(in) :: viscref
595  real(DP), intent(in) :: a2, a3, a4
596  integer(I4B), intent(in) :: locvisc
597  integer(I4B), dimension(:), intent(in) :: locconc
598  integer(I4B), dimension(:), intent(in) :: ivisc
599  real(DP), dimension(:), intent(in) :: visc
600  real(DP), dimension(:), intent(in) :: elev
601  real(DP), dimension(:), intent(in) :: dviscdc
602  real(DP), dimension(:), intent(in) :: cviscref
603  real(DP), dimension(:), intent(inout) :: ctemp
604  ! -- local
605  integer(I4B) :: n
606  integer(I4B) :: node
607  real(DP) :: visclak
608  !
609  ! -- update viscosity ratios for updating hyd. cond (and conductance)
610  select type (packobj)
611  type is (laktype)
612  do n = 1, packobj%nbound
613  !
614  ! -- get gwf node number
615  node = packobj%nodelist(n)
616  !
617  ! -- Check if boundary cell is active, cycle if not
618  if (packobj%ibound(node) <= 0) cycle
619  !
620  ! --
621  !
622  ! -- calculate the viscosity associated with the boundary
623  visclak = calc_bnd_viscosity(n, locvisc, locconc, viscref, dviscdc, &
624  cviscref, ctemp, ivisc, a2, a3, a4, &
625  packobj%auxvar)
626  !
627  ! -- fill lak relative viscosity into column 1 of viscratios
628  packobj%viscratios(1, n) = calc_vsc_ratio(viscref, visclak)
629  !
630  ! -- fill gwf relative viscosity into column 2 of viscratios
631  packobj%viscratios(2, n) = calc_vsc_ratio(viscref, visc(node))
632  end do
633  end select
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vsc_ad_maw()

subroutine gwfvscmodule::vsc_ad_maw ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  visc,
real(dp), intent(in)  viscref,
real(dp), dimension(:), intent(in)  elev,
integer(i4b), intent(in)  locvisc,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  dviscdc,
real(dp), dimension(:), intent(in)  cviscref,
integer(i4b), dimension(:), intent(in)  ivisc,
real(dp), intent(in)  a2,
real(dp), intent(in)  a3,
real(dp), intent(in)  a4,
real(dp), dimension(:), intent(inout)  ctemp 
)

When the viscosity package is active, update the viscosity ratio that is applied to the conductance calculated in the MAW package

Definition at line 641 of file gwf-vsc.f90.

643  ! -- modules
644  use bndmodule, only: bndtype
645  use mawmodule, only: mawtype
646  class(BndType), pointer :: packobj
647  ! -- dummy
648  real(DP), intent(in) :: viscref
649  real(DP), intent(in) :: a2, a3, a4
650  integer(I4B), intent(in) :: locvisc
651  integer(I4B), dimension(:), intent(in) :: locconc
652  integer(I4B), dimension(:), intent(in) :: ivisc
653  real(DP), dimension(:), intent(in) :: visc
654  real(DP), dimension(:), intent(in) :: elev
655  real(DP), dimension(:), intent(in) :: dviscdc
656  real(DP), dimension(:), intent(in) :: cviscref
657  real(DP), dimension(:), intent(inout) :: ctemp
658  ! -- local
659  integer(I4B) :: n
660  integer(I4B) :: node
661  real(DP) :: viscmaw
662  !
663  ! -- update viscosity ratios for updating hyd. cond (and conductance)
664  select type (packobj)
665  type is (mawtype)
666  do n = 1, packobj%nbound
667  !
668  ! -- get gwf node number
669  node = packobj%nodelist(n)
670  !
671  ! -- Check if boundary cell is active, cycle if not
672  if (packobj%ibound(node) <= 0) cycle
673  !
674  ! --
675  !
676  ! -- calculate the viscosity associated with the boundary
677  viscmaw = calc_bnd_viscosity(n, locvisc, locconc, viscref, dviscdc, &
678  cviscref, ctemp, ivisc, a2, a3, a4, &
679  packobj%auxvar)
680  !
681  ! -- fill lak relative viscosity into column 1 of viscratios
682  packobj%viscratios(1, n) = calc_vsc_ratio(viscref, viscmaw)
683  !
684  ! -- fill gwf relative viscosity into column 2 of viscratios
685  packobj%viscratios(2, n) = calc_vsc_ratio(viscref, visc(node))
686  end do
687  end select
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vsc_ad_sfr()

subroutine gwfvscmodule::vsc_ad_sfr ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  visc,
real(dp), intent(in)  viscref,
real(dp), dimension(:), intent(in)  elev,
integer(i4b), intent(in)  locvisc,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  dviscdc,
real(dp), dimension(:), intent(in)  cviscref,
integer(i4b), dimension(:), intent(in)  ivisc,
real(dp), intent(in)  a2,
real(dp), intent(in)  a3,
real(dp), intent(in)  a4,
real(dp), dimension(:), intent(inout)  ctemp 
)

When the viscosity package is active, update the viscosity ratio that is applied to the hydraulic conductivity specified in the SFR package

Definition at line 533 of file gwf-vsc.f90.

535  ! -- modules
536  use bndmodule, only: bndtype
537  use sfrmodule, only: sfrtype
538  class(BndType), pointer :: packobj
539  ! -- dummy
540  real(DP), intent(in) :: viscref
541  real(DP), intent(in) :: a2, a3, a4
542  integer(I4B), intent(in) :: locvisc
543  integer(I4B), dimension(:), intent(in) :: locconc
544  integer(I4B), dimension(:), intent(in) :: ivisc
545  real(DP), dimension(:), intent(in) :: visc
546  real(DP), dimension(:), intent(in) :: elev
547  real(DP), dimension(:), intent(in) :: dviscdc
548  real(DP), dimension(:), intent(in) :: cviscref
549  real(DP), dimension(:), intent(inout) :: ctemp
550  ! -- local
551  integer(I4B) :: n
552  integer(I4B) :: node
553  real(DP) :: viscsfr
554  !
555  ! -- update viscosity ratios for updating hyd. cond (and conductance)
556  select type (packobj)
557  type is (sfrtype)
558  do n = 1, packobj%nbound
559  !
560  ! -- get gwf node number
561  node = packobj%nodelist(n)
562  !
563  ! -- Check if boundary cell is active, cycle if not
564  if (packobj%ibound(node) <= 0) cycle
565  !
566  ! --
567  !
568  ! -- calculate the viscosity associated with the boundary
569  viscsfr = calc_bnd_viscosity(n, locvisc, locconc, viscref, dviscdc, &
570  cviscref, ctemp, ivisc, a2, a3, a4, &
571  packobj%auxvar)
572  !
573  ! -- fill sfr relative viscosity into column 1 of viscratios
574  packobj%viscratios(1, n) = calc_vsc_ratio(viscref, viscsfr)
575  !
576  ! -- fill gwf relative viscosity into column 2 of viscratios
577  packobj%viscratios(2, n) = calc_vsc_ratio(viscref, visc(node))
578  end do
579  end select
This module contains the SFR package methods.
Definition: gwf-sfr.f90:7
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vsc_ad_standard_bnd()

subroutine gwfvscmodule::vsc_ad_standard_bnd ( class(bndtype), pointer  packobj,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(in)  visc,
real(dp), intent(in)  viscref,
integer(i4b), intent(in)  locelev,
integer(i4b), intent(in)  locvisc,
integer(i4b), dimension(:), intent(in)  locconc,
real(dp), dimension(:), intent(in)  dviscdc,
real(dp), dimension(:), intent(in)  cviscref,
integer(i4b), dimension(:), intent(in)  ivisc,
real(dp), intent(in)  a2,
real(dp), intent(in)  a3,
real(dp), intent(in)  a4,
real(dp), dimension(:), intent(inout)  ctemp 
)

When flow enters from ghb boundary type, take into account the effects of viscosity on the user-specified conductance terms

Definition at line 462 of file gwf-vsc.f90.

465  ! -- modules
466  use bndmodule, only: bndtype
467  use drnmodule, only: drntype
468  use rivmodule, only: rivtype
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(:) :: visc
474  real(DP), intent(in) :: a2, a3, a4
475  real(DP), intent(in) :: viscref
476  integer(I4B), intent(in) :: locelev
477  integer(I4B), intent(in) :: locvisc
478  integer(I4B), dimension(:), intent(in) :: locconc
479  integer(I4B), dimension(:), intent(in) :: ivisc
480  real(DP), dimension(:), intent(in) :: dviscdc
481  real(DP), dimension(:), intent(in) :: cviscref
482  real(DP), dimension(:), intent(inout) :: ctemp
483  ! -- local
484  integer(I4B) :: n
485  integer(I4B) :: node
486  real(DP) :: viscbnd
487  !
488  ! -- Process density terms for each GHB
489  do n = 1, packobj%nbound
490  node = packobj%nodelist(n)
491  !
492  ! -- Check if boundary cell is active, cycle if not
493  if (packobj%ibound(node) <= 0) cycle
494  !
495  ! -- calculate the viscosity associated with the boundary
496  viscbnd = calc_bnd_viscosity(n, locvisc, locconc, viscref, dviscdc, &
497  cviscref, ctemp, ivisc, a2, a3, a4, &
498  packobj%auxvar)
499  !
500  ! -- update boundary conductance based on viscosity effects
501  select case (packobj%filtyp)
502  case ('DRN')
503  select type (packobj)
504  type is (drntype)
505  packobj%cond(n) = update_bnd_cond(viscbnd, viscref, &
506  packobj%condinput(n))
507  end select
508  case ('GHB')
509  select type (packobj)
510  type is (ghbtype)
511  packobj%cond(n) = update_bnd_cond(viscbnd, viscref, &
512  packobj%condinput(n))
513  end select
514  case ('RIV')
515  select type (packobj)
516  type is (rivtype)
517  packobj%cond(n) = update_bnd_cond(viscbnd, viscref, &
518  packobj%condinput(n))
519  end select
520  case default
521  packobj%bound(2, n) = update_bnd_cond(viscbnd, viscref, &
522  packobj%condinput(n))
523  end select
524  !
525  end do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vsc_ar()

subroutine gwfvscmodule::vsc_ar ( class(gwfvsctype this,
integer(i4b), dimension(:), pointer  ibound 
)
private

Generic method to allocate and read static data for the viscosity package available within the GWF model type.

Definition at line 215 of file gwf-vsc.f90.

216  ! -- dummy
217  class(GwfVscType) :: this
218  integer(I4B), dimension(:), pointer :: ibound
219  !
220  ! -- store pointers to arguments that were passed in
221  this%ibound => ibound
222  !
223  ! -- Set pointers to npf variables
224  call this%set_npf_pointers()

◆ vsc_ar_bnd()

subroutine gwfvscmodule::vsc_ar_bnd ( class(gwfvsctype this,
class(bndtype), pointer  packobj 
)
private

Viscosity ar_bnd routine to activate viscosity in the advanced packages. This routine is called from gwf_ar() as it moves through each package

Definition at line 233 of file gwf-vsc.f90.

234  ! -- modules
235  use bndmodule, only: bndtype
236  use drnmodule, only: drntype
237  use ghbmodule, only: ghbtype
238  use rivmodule, only: rivtype
239  use lakmodule, only: laktype
240  use sfrmodule, only: sfrtype
241  use mawmodule, only: mawtype
242  ! -- dummy
243  class(GwfVscType) :: this
244  class(BndType), pointer :: packobj
245  !
246  ! -- Add density terms based on boundary package type
247  select case (packobj%filtyp)
248  case ('DRN')
249  !
250  ! -- activate viscosity for the drain package
251  select type (packobj)
252  type is (drntype)
253  call packobj%bnd_activate_viscosity()
254  end select
255  case ('GHB')
256  !
257  ! -- activate viscosity for the drain package
258  select type (packobj)
259  type is (ghbtype)
260  call packobj%bnd_activate_viscosity()
261  end select
262  case ('RIV')
263  !
264  ! -- activate viscosity for the drain package
265  select type (packobj)
266  type is (rivtype)
267  call packobj%bnd_activate_viscosity()
268  end select
269  case ('LAK')
270  !
271  ! -- activate viscosity for lake package
272  select type (packobj)
273  type is (laktype)
274  call packobj%lak_activate_viscosity()
275  end select
276  case ('SFR')
277  !
278  ! -- activate viscosity for sfr package
279  select type (packobj)
280  type is (sfrtype)
281  call packobj%sfr_activate_viscosity()
282  end select
283  case ('MAW')
284  !
285  ! -- activate viscosity for maw package
286  select type (packobj)
287  type is (mawtype)
288  call packobj%maw_activate_viscosity()
289  end select
290  case default
291  !
292  ! -- nothing
293  end select

◆ vsc_calcvisc()

subroutine gwfvscmodule::vsc_calcvisc ( class(gwfvsctype this)
private

Calculates fluid viscosity based on concentration or temperature

Definition at line 1098 of file gwf-vsc.f90.

1099  ! -- dummy
1100  class(GwfVscType) :: this
1101  ! -- local
1102  integer(I4B) :: n
1103  integer(I4B) :: i
1104  !
1105  ! -- Calculate the viscosity using the specified concentration and/or
1106  ! temperature arrays
1107  do n = 1, this%dis%nodes
1108  do i = 1, this%nviscspecies
1109  if (this%modelconc(i)%icbund(n) == 0) then
1110  this%ctemp(i) = dzero
1111  else
1112  this%ctemp(i) = this%modelconc(i)%conc(n)
1113  end if
1114  end do
1115  !
1116  this%visc(n) = calc_visc(this%ivisc, this%viscref, this%dviscdc, &
1117  this%cviscref, this%ctemp, this%a2, &
1118  this%a3, this%a4)
1119  end do
Here is the call graph for this function:

◆ vsc_cr()

subroutine, public gwfvscmodule::vsc_cr ( type(gwfvsctype), pointer  vscobj,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  input_mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Create a new VSC Package object.

Definition at line 142 of file gwf-vsc.f90.

143  ! -- dummy
144  type(GwfVscType), pointer :: vscobj
145  character(len=*), intent(in) :: name_model
146  character(len=*), intent(in) :: input_mempath
147  integer(I4B), intent(in) :: inunit
148  integer(I4B), intent(in) :: iout
149  !
150  ! -- Create the object
151  allocate (vscobj)
152  !
153  ! -- create name and memory path
154  call vscobj%set_names(1, name_model, 'VSC', 'VSC', input_mempath)
155  !
156  ! -- Allocate scalars
157  call vscobj%allocate_scalars()
158  !
159  ! -- Set variables
160  vscobj%inunit = inunit
161  vscobj%iout = iout
Here is the caller graph for this function:

◆ vsc_da()

subroutine gwfvscmodule::vsc_da ( class(gwfvsctype this)
private

Deallocate viscosity package scalars and arrays.

Definition at line 904 of file gwf-vsc.f90.

905  ! -- dummy
906  class(GwfVscType) :: this
907  !
908  ! -- Deallocate arrays if package was active
909  if (this%inunit > 0) then
910  call mem_deallocate(this%visc)
911  call mem_deallocate(this%ivisc)
912  call mem_deallocate(this%dviscdc)
913  call mem_deallocate(this%cviscref)
914  call mem_deallocate(this%ctemp)
915  deallocate (this%cmodelname)
916  deallocate (this%cauxspeciesname)
917  deallocate (this%modelconc)
918  end if
919  !
920  ! -- Scalars
921  call mem_deallocate(this%thermivisc)
922  call mem_deallocate(this%idxtmpr)
923  call mem_deallocate(this%ioutvisc)
924  call mem_deallocate(this%ireadelev)
925  call mem_deallocate(this%iconcset)
926  call mem_deallocate(this%viscref)
927  call mem_deallocate(this%nviscspecies)
928  call mem_deallocate(this%a2)
929  call mem_deallocate(this%a3)
930  call mem_deallocate(this%a4)
931  !
932  ! -- Nullify pointers to other package variables
933  nullify (this%k11)
934  nullify (this%k22)
935  nullify (this%k33)
936  nullify (this%k11input)
937  nullify (this%k22input)
938  nullify (this%k33input)
939  nullify (this%kchangeper)
940  nullify (this%kchangestp)
941  nullify (this%nodekchange)
942  !
943  ! -- deallocate parent
944  call this%NumericalPackageType%da()

◆ vsc_df()

subroutine gwfvscmodule::vsc_df ( class(gwfvsctype this,
class(disbasetype), intent(in), pointer  dis,
type(gwfvscinputdatatype), intent(in), optional  vsc_input 
)
private

Define viscosity package options and dimensions

Parameters
thisthis viscosity package
[in]dispointer to discretization
[in]vsc_inputoptional vsc input data, otherwise read from file

Definition at line 168 of file gwf-vsc.f90.

169  ! -- dummy
170  class(GwfVscType) :: this !< this viscosity package
171  class(DisBaseType), pointer, intent(in) :: dis !< pointer to discretization
172  type(GwfVscInputDataType), optional, intent(in) :: vsc_input !< optional vsc input data, otherwise read from file
173  ! -- formats
174  character(len=*), parameter :: fmtvsc = &
175  "(1x,/1x,'VSC -- Viscosity Package, version 1, 11/15/2022', &
176  &' input read from mempath: ', a, /)"
177  !
178  ! --print a message identifying the viscosity package
179  write (this%iout, fmtvsc) this%input_mempath
180  !
181  ! -- store pointers to arguments that were passed in
182  this%dis => dis
183  !
184  if (.not. present(vsc_input)) then
185  !
186  ! -- Read viscosity options
187  call this%source_options()
188  !
189  ! -- Read viscosity dimensions
190  call this%source_dimensions()
191  else
192  ! set from input data instead
193  call this%set_options(vsc_input)
194  this%nviscspecies = vsc_input%nviscspecies
195  end if
196  !
197  ! -- Allocate arrays
198  call this%allocate_arrays(dis%nodes)
199  !
200  if (.not. present(vsc_input)) then
201  !
202  ! -- Read viscosity packagedata
203  call this%source_packagedata()
204  else
205  ! set from input data instead
206  call this%set_packagedata(vsc_input)
207  end if

◆ vsc_ot_dv()

subroutine gwfvscmodule::vsc_ot_dv ( class(gwfvsctype this,
integer(i4b), intent(in)  idvfl 
)
private

Save calculated viscosity array to binary file

Definition at line 865 of file gwf-vsc.f90.

866  ! -- dummy
867  class(GwfVscType) :: this
868  integer(I4B), intent(in) :: idvfl
869  ! -- local
870  character(len=1) :: cdatafmp = ' ', editdesc = ' '
871  integer(I4B) :: ibinun
872  integer(I4B) :: iprint
873  integer(I4B) :: nvaluesp
874  integer(I4B) :: nwidthp
875  real(DP) :: dinact
876  !
877  ! -- Set unit number for viscosity output
878  if (this%ioutvisc /= 0) then
879  ibinun = 1
880  else
881  ibinun = 0
882  end if
883  if (idvfl == 0) ibinun = 0
884  !
885  ! -- save viscosity array
886  if (ibinun /= 0) then
887  iprint = 0
888  dinact = dhnoflo
889  !
890  ! -- write viscosity to binary file
891  if (this%ioutvisc /= 0) then
892  ibinun = this%ioutvisc
893  call this%dis%record_array(this%visc, this%iout, iprint, ibinun, &
894  ' VISCOSITY', cdatafmp, nvaluesp, &
895  nwidthp, editdesc, dinact)
896  end if
897  end if

◆ vsc_rp()

subroutine gwfvscmodule::vsc_rp ( class(gwfvsctype this)
private

Method to read and prepare period data for the VSC package.

Definition at line 325 of file gwf-vsc.f90.

326  ! -- modules
327  use tdismodule, only: kstp, kper
328  ! -- dummy
329  class(GwfVscType) :: this
330  ! -- local
331  integer(I4B) :: i
332  ! -- formats
333  character(len=*), parameter :: fmtc = &
334  "('Viscosity Package does not have a concentration set &
335  &for species ',i0,'. One or more model names may be specified &
336  &incorrectly in the PACKAGEDATA block or a GWF-GWT exchange may need &
337  &to be activated.')"
338  !
339  ! -- Check to make sure all concentration pointers have been set
340  if (kstp * kper == 1) then
341  do i = 1, this%nviscspecies
342  if (.not. associated(this%modelconc(i)%conc)) then
343  write (errmsg, fmtc) i
344  call store_error(errmsg)
345  end if
346  end do
347  if (count_errors() > 0) then
348  call store_error_filename(this%input_fname)
349  end if
350  end if
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:

◆ vsc_set_changed_at()

subroutine gwfvscmodule::vsc_set_changed_at ( class(gwfvsctype this,
integer(i4b), intent(in)  kper,
integer(i4b), intent(in)  kstp 
)
private

Procedure called by VSC code when K updated due to viscosity changes. K values changed at (kper, kstp).

Definition at line 851 of file gwf-vsc.f90.

852  ! -- dummy variables
853  class(GwfVscType) :: this
854  integer(I4B), intent(in) :: kper
855  integer(I4B), intent(in) :: kstp
856  !
857  this%kchangeper = kper
858  this%kchangestp = kstp