MODFLOW 6
version 6.7.0.dev0
USGS Modular Hydrologic Model
|
This module contains stateless conductance functions. More...
Enumerations | |
enum | { ccond_lmean = 1 , ccond_amtlmk = 2 , ccond_amthmk = 3 } |
enumerator that defines the conductance options More... | |
Functions/Subroutines | |
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. More... | |
real(dp) function | convertible_upstream (hn, hm, satn, satm, condsat) |
Convertible cell(s) with upstream weighted horizontal conductance. More... | |
real(dp) function | convertible_standard (ihc, icellavg, satn, satm, hkn, hkm, topn, topm, botn, botm, cln, clm, fawidth) |
Convertible cell(s) with standard weighted horizontal conductance. More... | |
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. More... | |
real(dp) function, public | condmean (k1, k2, thick1, thick2, cl1, cl2, width, iavgmeth) |
Calculate the conductance between two cells. More... | |
real(dp) function | logmean (d1, d2) |
Calculate the the logarithmic mean of two double precision numbers. More... | |
real(dp) function, public | thksatnm (ibdn, ibdm, ictn, ictm, iupstream, ihc, hn, hm, satn, satm, topn, topm, botn, botm) |
Calculate wetted cell thickness at interface between two cells. More... | |
real(dp) function, public | staggered_thkfrac (top, bot, sat, topc, botc) |
Calculate the thickness fraction for staggered grids. More... | |
Variables | |
@, public | ccond_hmean = 0 |
Harmonic mean. More... | |
This module contains the functions to calculate the horizontal and vertical conductance between two cells that are used in the the node property flow (NPF) package. It also contains functions to calculate the wetted cell fraction. This module does not depend on the NPF package.
|
private |
Definition at line 27 of file GwfConductanceUtils.f90.
real(dp) function, public gwfconductanceutilsmodule::condmean | ( | real(dp), intent(in) | k1, |
real(dp), intent(in) | k2, | ||
real(dp), intent(in) | thick1, | ||
real(dp), intent(in) | thick2, | ||
real(dp), intent(in) | cl1, | ||
real(dp), intent(in) | cl2, | ||
real(dp), intent(in) | width, | ||
integer(i4b), intent(in) | iavgmeth | ||
) |
[in] | k1 | hydraulic conductivity for cell n (in the direction of cell m) |
[in] | k2 | hydraulic conductivity for cell m (in the direction of celln) |
[in] | thick1 | saturated thickness for cell 1 |
[in] | thick2 | saturated thickness for cell 2 |
[in] | cl1 | distance from the center of cell n to the shared face with cell m |
[in] | cl2 | distance from the center of cell m to the shared face with cell n |
[in] | width | width of cell perpendicular to flow |
[in] | iavgmeth | averaging method |
Definition at line 226 of file GwfConductanceUtils.f90.
|
private |
[in] | ihc | connection type |
[in] | icellavg | cell averaging option |
[in] | satn | cell n wetted cell fraction |
[in] | satm | cell m wetted cell fraction |
[in] | hkn | horizontal hydraulic conductivity for cell n (in the direction of cell m) |
[in] | hkm | horizontal hydraulic conductivity for cell m (in the direction of cell n) |
[in] | topn | top of cell n |
[in] | topm | top of cell m |
[in] | botn | bottom of cell n |
[in] | botm | bottom of cell m |
[in] | cln | distance from the center of cell n to the shared face with cell m |
[in] | clm | distance from the center of cell m to the shared face with cell n |
[in] | fawidth | width of cell perpendicular to flow |
Definition at line 113 of file GwfConductanceUtils.f90.
|
private |
[in] | condsat | saturated conductance |
[in] | hn | cell n head |
[in] | hm | cell m head |
[in] | satn | cell n wetted cell fraction |
[in] | satm | cell m wetted cell fraction |
Definition at line 90 of file GwfConductanceUtils.f90.
real(dp) function, public gwfconductanceutilsmodule::hcond | ( | integer(i4b), intent(in) | ibdn, |
integer(i4b), intent(in) | ibdm, | ||
integer(i4b), intent(in) | ictn, | ||
integer(i4b), intent(in) | ictm, | ||
integer(i4b), intent(in) | iupstream, | ||
integer(i4b), intent(in) | ihc, | ||
integer(i4b), intent(in) | icellavg, | ||
real(dp), intent(in) | condsat, | ||
real(dp), intent(in) | hn, | ||
real(dp), intent(in) | hm, | ||
real(dp), intent(in) | satn, | ||
real(dp), intent(in) | satm, | ||
real(dp), intent(in) | hkn, | ||
real(dp), intent(in) | hkm, | ||
real(dp), intent(in) | topn, | ||
real(dp), intent(in) | topm, | ||
real(dp), intent(in) | botn, | ||
real(dp), intent(in) | botm, | ||
real(dp), intent(in) | cln, | ||
real(dp), intent(in) | clm, | ||
real(dp), intent(in) | fawidth | ||
) |
This function uses a weighted transmissivity in the harmonic mean conductance calculations. This differs from the MODFLOW-NWT and MODFLOW-USG conductance calculations for the Newton-Raphson formulation which use a weighted hydraulic conductivity.
[in] | ibdn | cell n active flag |
[in] | ibdm | cell m active flag |
[in] | ictn | cell n convertible cell flag |
[in] | ictm | cell m convertible cell flag |
[in] | iupstream | flag for upstream weighting |
[in] | ihc | connection type |
[in] | icellavg | cell averaging option |
[in] | condsat | saturated conductance |
[in] | hn | cell n head |
[in] | hm | cell m head |
[in] | satn | cell n wetted cell fraction |
[in] | satm | cell m wetted cell fraction |
[in] | hkn | horizontal hydraulic conductivity for cell n (in the direction of cell m) |
[in] | hkm | horizontal hydraulic conductivity for cell m (in the direction of cell n) |
[in] | topn | top of cell n |
[in] | topm | top of cell m |
[in] | botn | bottom of cell n |
[in] | botm | bottom of cell m |
[in] | cln | distance from the center of cell n to the shared face with cell m |
[in] | clm | distance from the center of cell m to the shared face with cell n |
[in] | fawidth | width of cell perpendicular to flow |
Definition at line 43 of file GwfConductanceUtils.f90.
|
private |
Use an approximation if the ratio is near 1
[in] | d1 | first number |
[in] | d2 | second number |
Definition at line 290 of file GwfConductanceUtils.f90.
real(dp) function, public gwfconductanceutilsmodule::staggered_thkfrac | ( | real(dp) | top, |
real(dp) | bot, | ||
real(dp) | sat, | ||
real(dp) | topc, | ||
real(dp) | botc | ||
) |
top | top of cell |
bot | bottom of cell |
sat | cell saturation |
topc | top of connected cell |
botc | bottom of connected cells |
Definition at line 375 of file GwfConductanceUtils.f90.
real(dp) function, public gwfconductanceutilsmodule::thksatnm | ( | integer(i4b), intent(in) | ibdn, |
integer(i4b), intent(in) | ibdm, | ||
integer(i4b), intent(in) | ictn, | ||
integer(i4b), intent(in) | ictm, | ||
integer(i4b), intent(in) | iupstream, | ||
integer(i4b), intent(in) | ihc, | ||
real(dp), intent(in) | hn, | ||
real(dp), intent(in) | hm, | ||
real(dp), intent(in) | satn, | ||
real(dp), intent(in) | satm, | ||
real(dp), intent(in) | topn, | ||
real(dp), intent(in) | topm, | ||
real(dp), intent(in) | botn, | ||
real(dp), intent(in) | botm | ||
) |
[in] | ibdn | cell n active flag |
[in] | ibdm | cell m active flag |
[in] | ictn | cell n convertible cell flag |
[in] | ictm | cell m convertible cell flag |
[in] | iupstream | flag for upstream weighting |
[in] | ihc | connection type |
[in] | hn | cell n head |
[in] | hm | cell m head |
[in] | satn | cell n wetted cell fraction |
[in] | satm | cell m wetted cell fraction |
[in] | topn | top of cell n |
[in] | topm | top of cell m |
[in] | botn | bottom of cell n |
[in] | botm | bottom of cell m |
Definition at line 309 of file GwfConductanceUtils.f90.
real(dp) function, public gwfconductanceutilsmodule::vcond | ( | integer(i4b), intent(in) | ibdn, |
integer(i4b), intent(in) | ibdm, | ||
integer(i4b), intent(in) | ictn, | ||
integer(i4b), intent(in) | ictm, | ||
integer(i4b), intent(in) | inewton, | ||
integer(i4b), intent(in) | ivarcv, | ||
integer(i4b), intent(in) | idewatcv, | ||
real(dp), intent(in) | condsat, | ||
real(dp), intent(in) | hn, | ||
real(dp), intent(in) | hm, | ||
real(dp), intent(in) | vkn, | ||
real(dp), intent(in) | vkm, | ||
real(dp), intent(in) | satn, | ||
real(dp), intent(in) | satm, | ||
real(dp), intent(in) | topn, | ||
real(dp), intent(in) | topm, | ||
real(dp), intent(in) | botn, | ||
real(dp), intent(in) | botm, | ||
real(dp), intent(in) | flowarea | ||
) |
[in] | ibdn | cell n active flag |
[in] | ibdm | cell m active flag |
[in] | ictn | cell n convertible cell flag |
[in] | ictm | cell m convertible cell flag |
[in] | inewton | flag for Newton-Raphson formulation |
[in] | ivarcv | variable vertical conductance flag |
[in] | idewatcv | dewatered vertical conductance flag |
[in] | condsat | saturated conductance |
[in] | hn | cell n head |
[in] | hm | cell m head |
[in] | vkn | vertical hydraulic conductivity for cell n (in the direction of cell m) |
[in] | vkm | vertical hydraulic conductivity for cell m (in the direction of cell n) |
[in] | satn | cell n wetted cell fraction |
[in] | satm | cell m wetted cell fraction |
[in] | topn | top of cell n |
[in] | topm | top of cell m |
[in] | botn | bottom of cell n |
[in] | botm | bottom of cell m |
[in] | flowarea | flow area between cell n and m |
Definition at line 149 of file GwfConductanceUtils.f90.
@, public gwfconductanceutilsmodule::ccond_hmean = 0 |
Definition at line 28 of file GwfConductanceUtils.f90.