MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
gwfnpfmodule Module Reference

Data Types

type  gwfnpftype
 

Functions/Subroutines

subroutine, public npf_cr (npfobj, name_model, input_mempath, inunit, iout)
 Create a new NPF object. Pass a inunit value of 0 if npf data will initialized from memory. More...
 
subroutine npf_df (this, dis, xt3d, ingnc, invsc, npf_options)
 Define the NPF package instance. More...
 
subroutine npf_ac (this, moffset, sparse)
 Add connections for extended neighbors to the sparse matrix. More...
 
subroutine npf_mc (this, moffset, matrix_sln)
 Map connections and construct iax, jax, and idxglox. More...
 
subroutine npf_ar (this, ic, vsc, ibound, hnew)
 Allocate and read this NPF instance. More...
 
subroutine npf_rp (this)
 Read and prepare method for package. More...
 
subroutine npf_ad (this, nodes, hold, hnew, irestore)
 Advance. More...
 
subroutine npf_cf (this, kiter, nodes, hnew)
 Routines associated fill coefficients. More...
 
subroutine npf_fc (this, kiter, matrix_sln, idxglo, rhs, hnew)
 Formulate coefficients. More...
 
subroutine npf_fn (this, kiter, matrix_sln, idxglo, rhs, hnew)
 Fill newton terms. More...
 
subroutine npf_nur (this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
 Under-relaxation. More...
 
subroutine npf_cq (this, hnew, flowja)
 Calculate flowja. More...
 
subroutine sgwf_npf_thksat (this, n, hn, thksat)
 Fractional cell saturation. More...
 
subroutine sgwf_npf_qcalc (this, n, m, hn, hm, icon, qnm)
 Flow between two cells. More...
 
subroutine npf_save_model_flows (this, flowja, icbcfl, icbcun)
 Record flowja and calculate specific discharge if requested. More...
 
subroutine npf_print_model_flows (this, ibudfl, flowja)
 Print budget. More...
 
subroutine npf_da (this)
 Deallocate variables. More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine store_original_k_arrays (this, ncells, njas)
 @ brief Store backup copy of hydraulic conductivity when the VSC package is activate More...
 
subroutine allocate_arrays (this, ncells, njas)
 Allocate npf arrays. More...
 
subroutine log_options (this, found)
 Log npf options sourced from the input mempath. More...
 
subroutine source_options (this)
 Update simulation options from input mempath. More...
 
subroutine set_options (this, options)
 Set options in the NPF object. More...
 
subroutine check_options (this)
 Check for conflicting NPF options. More...
 
subroutine log_griddata (this, found)
 Write dimensions to list file. More...
 
subroutine source_griddata (this)
 Update simulation griddata from input mempath. More...
 
subroutine prepcheck (this)
 Initialize and check NPF data. More...
 
subroutine preprocess_input (this)
 preprocess the NPF input data More...
 
subroutine calc_condsat (this, node, upperOnly)
 Calculate CONDSAT array entries for the given node. More...
 
real(dp) function calc_initial_sat (this, n)
 Calculate initial saturation for the given node. More...
 
subroutine sgwf_npf_wetdry (this, kiter, hnew)
 Perform wetting and drying. More...
 
subroutine rewet_check (this, kiter, node, hm, ibdm, ihc, hnew, irewet)
 Determine if a cell should rewet. More...
 
subroutine sgwf_npf_wdmsg (this, icode, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
 Print wet/dry message. More...
 
real(dp) function hy_eff (this, n, m, ihc, ipos, vg)
 Calculate the effective hydraulic conductivity for the n-m connection. More...
 
subroutine calc_spdis (this, flowja)
 Calculate the 3 components of specific discharge at the cell center. More...
 
subroutine sav_spdis (this, ibinun)
 Save specific discharge in binary format to ibinun. More...
 
subroutine sav_sat (this, ibinun)
 Save saturation in binary format to ibinun. More...
 
subroutine increase_edge_count (this, nedges)
 Reserve space for nedges cells that have an edge on them. More...
 
integer(i4b) function calc_max_conns (this)
 Calculate the maximum number of connections for any cell. More...
 
subroutine set_edge_properties (this, nodedge, ihcedge, q, area, nx, ny, distance)
 Provide the npf package with edge properties. More...
 
subroutine prepare_edge_lookup (this)
 
real(dp) function calcsatthickness (this, n, m, ihc)
 Calculate saturated thickness between cell n and m. More...
 

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwfnpfmodule::allocate_arrays ( class(gwfnpftype this,
integer(i4b), intent(in)  ncells,
integer(i4b), intent(in)  njas 
)

Definition at line 1163 of file gwf-npf.f90.

1164  ! -- dummy
1165  class(GwfNpftype) :: this
1166  integer(I4B), intent(in) :: ncells
1167  integer(I4B), intent(in) :: njas
1168  ! -- local
1169  integer(I4B) :: n
1170  !
1171  call mem_allocate(this%ithickstartflag, ncells, 'ITHICKSTARTFLAG', &
1172  this%memoryPath)
1173  call mem_allocate(this%icelltype, ncells, 'ICELLTYPE', this%memoryPath)
1174  call mem_allocate(this%k11, ncells, 'K11', this%memoryPath)
1175  call mem_allocate(this%sat, ncells, 'SAT', this%memoryPath)
1176  call mem_allocate(this%condsat, njas, 'CONDSAT', this%memoryPath)
1177  !
1178  ! -- Optional arrays dimensioned to full size initially
1179  call mem_allocate(this%k22, ncells, 'K22', this%memoryPath)
1180  call mem_allocate(this%k33, ncells, 'K33', this%memoryPath)
1181  call mem_allocate(this%wetdry, ncells, 'WETDRY', this%memoryPath)
1182  call mem_allocate(this%angle1, ncells, 'ANGLE1', this%memoryPath)
1183  call mem_allocate(this%angle2, ncells, 'ANGLE2', this%memoryPath)
1184  call mem_allocate(this%angle3, ncells, 'ANGLE3', this%memoryPath)
1185  !
1186  ! -- Optional arrays
1187  call mem_allocate(this%ibotnode, 0, 'IBOTNODE', this%memoryPath)
1188  call mem_allocate(this%nodedge, 0, 'NODEDGE', this%memoryPath)
1189  call mem_allocate(this%ihcedge, 0, 'IHCEDGE', this%memoryPath)
1190  call mem_allocate(this%propsedge, 0, 0, 'PROPSEDGE', this%memoryPath)
1191  call mem_allocate(this%iedge_ptr, 0, 'NREDGESNODE', this%memoryPath)
1192  call mem_allocate(this%edge_idxs, 0, 'EDGEIDXS', this%memoryPath)
1193  !
1194  ! -- Optional arrays only needed when vsc package is active
1195  call mem_allocate(this%k11input, 0, 'K11INPUT', this%memoryPath)
1196  call mem_allocate(this%k22input, 0, 'K22INPUT', this%memoryPath)
1197  call mem_allocate(this%k33input, 0, 'K33INPUT', this%memoryPath)
1198  !
1199  ! -- Specific discharge is (re-)allocated when nedges is known
1200  call mem_allocate(this%spdis, 3, 0, 'SPDIS', this%memoryPath)
1201  !
1202  ! -- Time-varying property flag arrays
1203  call mem_allocate(this%nodekchange, ncells, 'NODEKCHANGE', this%memoryPath)
1204  !
1205  ! -- initialize iangle1, iangle2, iangle3, and wetdry
1206  do n = 1, ncells
1207  this%angle1(n) = dzero
1208  this%angle2(n) = dzero
1209  this%angle3(n) = dzero
1210  this%wetdry(n) = dzero
1211  this%nodekchange(n) = dzero
1212  end do
1213  !
1214  ! -- allocate variable names
1215  allocate (this%aname(this%iname))
1216  this%aname = [' ICELLTYPE', ' K', &
1217  ' K33', ' K22', &
1218  ' WETDRY', ' ANGLE1', &
1219  ' ANGLE2', ' ANGLE3']

◆ allocate_scalars()

subroutine gwfnpfmodule::allocate_scalars ( class(gwfnpftype this)

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

Definition at line 1046 of file gwf-npf.f90.

1047  ! -- modules
1049  ! -- dummy
1050  class(GwfNpftype) :: this
1051  !
1052  ! -- allocate scalars in NumericalPackageType
1053  call this%NumericalPackageType%allocate_scalars()
1054  !
1055  ! -- Allocate scalars
1056  call mem_allocate(this%iname, 'INAME', this%memoryPath)
1057  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
1058  call mem_allocate(this%ixt3drhs, 'IXT3DRHS', this%memoryPath)
1059  call mem_allocate(this%satomega, 'SATOMEGA', this%memoryPath)
1060  call mem_allocate(this%hnoflo, 'HNOFLO', this%memoryPath)
1061  call mem_allocate(this%hdry, 'HDRY', this%memoryPath)
1062  call mem_allocate(this%icellavg, 'ICELLAVG', this%memoryPath)
1063  call mem_allocate(this%iavgkeff, 'IAVGKEFF', this%memoryPath)
1064  call mem_allocate(this%ik22, 'IK22', this%memoryPath)
1065  call mem_allocate(this%ik33, 'IK33', this%memoryPath)
1066  call mem_allocate(this%ik22overk, 'IK22OVERK', this%memoryPath)
1067  call mem_allocate(this%ik33overk, 'IK33OVERK', this%memoryPath)
1068  call mem_allocate(this%iperched, 'IPERCHED', this%memoryPath)
1069  call mem_allocate(this%ivarcv, 'IVARCV', this%memoryPath)
1070  call mem_allocate(this%idewatcv, 'IDEWATCV', this%memoryPath)
1071  call mem_allocate(this%ithickstrt, 'ITHICKSTRT', this%memoryPath)
1072  call mem_allocate(this%icalcspdis, 'ICALCSPDIS', this%memoryPath)
1073  call mem_allocate(this%isavspdis, 'ISAVSPDIS', this%memoryPath)
1074  call mem_allocate(this%isavsat, 'ISAVSAT', this%memoryPath)
1075  call mem_allocate(this%irewet, 'IREWET', this%memoryPath)
1076  call mem_allocate(this%wetfct, 'WETFCT', this%memoryPath)
1077  call mem_allocate(this%iwetit, 'IWETIT', this%memoryPath)
1078  call mem_allocate(this%ihdwet, 'IHDWET', this%memoryPath)
1079  call mem_allocate(this%iangle1, 'IANGLE1', this%memoryPath)
1080  call mem_allocate(this%iangle2, 'IANGLE2', this%memoryPath)
1081  call mem_allocate(this%iangle3, 'IANGLE3', this%memoryPath)
1082  call mem_allocate(this%iwetdry, 'IWETDRY', this%memoryPath)
1083  call mem_allocate(this%nedges, 'NEDGES', this%memoryPath)
1084  call mem_allocate(this%lastedge, 'LASTEDGE', this%memoryPath)
1085  call mem_allocate(this%intvk, 'INTVK', this%memoryPath)
1086  call mem_allocate(this%invsc, 'INVSC', this%memoryPath)
1087  call mem_allocate(this%kchangeper, 'KCHANGEPER', this%memoryPath)
1088  call mem_allocate(this%kchangestp, 'KCHANGESTP', this%memoryPath)
1089  !
1090  ! -- set pointer to inewtonur
1091  call mem_setptr(this%igwfnewtonur, 'INEWTONUR', &
1092  create_mem_path(this%name_model))
1093  !
1094  ! -- Initialize value
1095  this%iname = 8
1096  this%ixt3d = 0
1097  this%ixt3drhs = 0
1098  this%satomega = dzero
1099  this%hnoflo = dhnoflo !1.d30
1100  this%hdry = dhdry !-1.d30
1101  this%icellavg = ccond_hmean
1102  this%iavgkeff = 0
1103  this%ik22 = 0
1104  this%ik33 = 0
1105  this%ik22overk = 0
1106  this%ik33overk = 0
1107  this%iperched = 0
1108  this%ivarcv = 0
1109  this%idewatcv = 0
1110  this%ithickstrt = 0
1111  this%icalcspdis = 0
1112  this%isavspdis = 0
1113  this%isavsat = 0
1114  this%irewet = 0
1115  this%wetfct = done
1116  this%iwetit = 1
1117  this%ihdwet = 0
1118  this%iangle1 = 0
1119  this%iangle2 = 0
1120  this%iangle3 = 0
1121  this%iwetdry = 0
1122  this%nedges = 0
1123  this%lastedge = 0
1124  this%intvk = 0
1125  this%invsc = 0
1126  this%kchangeper = 0
1127  this%kchangestp = 0
1128  !
1129  ! -- If newton is on, then NPF creates asymmetric matrix
1130  this%iasym = this%inewton
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Here is the call graph for this function:

◆ calc_condsat()

subroutine gwfnpfmodule::calc_condsat ( class(gwfnpftype this,
integer(i4b), intent(in)  node,
logical, intent(in)  upperOnly 
)

Calculate saturated conductances for all connections of the given node, or optionally for the upper portion of the matrix only.

Definition at line 1953 of file gwf-npf.f90.

1954  ! -- dummy variables
1955  class(GwfNpfType) :: this
1956  integer(I4B), intent(in) :: node
1957  logical, intent(in) :: upperOnly
1958  ! -- local variables
1959  integer(I4B) :: ii, m, n, ihc, jj
1960  real(DP) :: topm, topn, topnode, botm, botn, botnode, satm, satn, satnode
1961  real(DP) :: hyn, hym, hn, hm, fawidth, csat
1962  !
1963  satnode = this%calc_initial_sat(node)
1964  !
1965  topnode = this%dis%top(node)
1966  botnode = this%dis%bot(node)
1967  !
1968  ! -- Go through the connecting cells
1969  do ii = this%dis%con%ia(node) + 1, this%dis%con%ia(node + 1) - 1
1970  !
1971  ! -- Set the m cell number and cycle if lower triangle connection and
1972  ! -- we're not updating both upper and lower matrix parts for this node
1973  m = this%dis%con%ja(ii)
1974  jj = this%dis%con%jas(ii)
1975  if (m < node) then
1976  if (upperonly) cycle
1977  ! m => node, n => neighbour
1978  n = m
1979  m = node
1980  topm = topnode
1981  botm = botnode
1982  satm = satnode
1983  topn = this%dis%top(n)
1984  botn = this%dis%bot(n)
1985  satn = this%calc_initial_sat(n)
1986  else
1987  ! n => node, m => neighbour
1988  n = node
1989  topn = topnode
1990  botn = botnode
1991  satn = satnode
1992  topm = this%dis%top(m)
1993  botm = this%dis%bot(m)
1994  satm = this%calc_initial_sat(m)
1995  end if
1996  !
1997  ihc = this%dis%con%ihc(jj)
1998  hyn = this%hy_eff(n, m, ihc, ipos=ii)
1999  hym = this%hy_eff(m, n, ihc, ipos=ii)
2000  if (this%ithickstartflag(n) == 0) then
2001  hn = topn
2002  else
2003  hn = this%ic%strt(n)
2004  end if
2005  if (this%ithickstartflag(m) == 0) then
2006  hm = topm
2007  else
2008  hm = this%ic%strt(m)
2009  end if
2010  !
2011  ! -- Calculate conductance depending on whether connection is
2012  ! vertical (0), horizontal (1), or staggered horizontal (2)
2013  if (ihc == c3d_vertical) then
2014  !
2015  ! -- Vertical conductance for fully saturated conditions
2016  csat = vcond(1, 1, 1, 1, 0, 1, 1, done, &
2017  botn, botm, &
2018  hyn, hym, &
2019  satn, satm, &
2020  topn, topm, &
2021  botn, botm, &
2022  this%dis%con%hwva(jj))
2023  else
2024  !
2025  ! -- Horizontal conductance for fully saturated conditions
2026  fawidth = this%dis%con%hwva(jj)
2027  csat = hcond(1, 1, 1, 1, 0, &
2028  ihc, &
2029  this%icellavg, &
2030  done, &
2031  hn, hm, satn, satm, hyn, hym, &
2032  topn, topm, &
2033  botn, botm, &
2034  this%dis%con%cl1(jj), &
2035  this%dis%con%cl2(jj), &
2036  fawidth)
2037  end if
2038  this%condsat(jj) = csat
2039  end do
Here is the call graph for this function:

◆ calc_initial_sat()

real(dp) function gwfnpfmodule::calc_initial_sat ( class(gwfnpftype this,
integer(i4b), intent(in)  n 
)
private

Calculate saturation as a fraction of thickness for the given node, used for saturated conductance calculations: full thickness by default (1.0) or saturation based on initial conditions if the THICKSTRT option is used.

Definition at line 2049 of file gwf-npf.f90.

2050  ! -- dummy variables
2051  class(GwfNpfType) :: this
2052  integer(I4B), intent(in) :: n
2053  ! -- Return
2054  real(DP) :: satn
2055  !
2056  satn = done
2057  if (this%ibound(n) /= 0 .and. this%ithickstartflag(n) /= 0) then
2058  call this%thksat(n, this%ic%strt(n), satn)
2059  end if

◆ calc_max_conns()

integer(i4b) function gwfnpfmodule::calc_max_conns ( class(gwfnpftype this)
private

Definition at line 2683 of file gwf-npf.f90.

2684  class(GwfNpfType) :: this
2685  integer(I4B) :: max_conns
2686  ! local
2687  integer(I4B) :: n, m, ic
2688 
2689  max_conns = 0
2690  do n = 1, this%dis%nodes
2691 
2692  ! Count internal model connections
2693  ic = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
2694 
2695  ! Add edge connections
2696  do m = 1, this%nedges
2697  if (this%nodedge(m) == n) then
2698  ic = ic + 1
2699  end if
2700  end do
2701 
2702  ! Set max number of connections for any cell
2703  if (ic > max_conns) max_conns = ic
2704  end do
2705 

◆ calc_spdis()

subroutine gwfnpfmodule::calc_spdis ( class(gwfnpftype this,
real(dp), dimension(:), intent(in)  flowja 
)
private

Definition at line 2362 of file gwf-npf.f90.

2363  ! -- modules
2364  use simmodule, only: store_error
2365  ! -- dummy
2366  class(GwfNpfType) :: this
2367  real(DP), intent(in), dimension(:) :: flowja
2368  ! -- local
2369  integer(I4B) :: n
2370  integer(I4B) :: m
2371  integer(I4B) :: ipos
2372  integer(I4B) :: iedge
2373  integer(I4B) :: isympos
2374  integer(I4B) :: ihc
2375  integer(I4B) :: ic
2376  integer(I4B) :: iz
2377  integer(I4B) :: nc
2378  integer(I4B) :: ncz
2379  real(DP) :: qz
2380  real(DP) :: vx
2381  real(DP) :: vy
2382  real(DP) :: vz
2383  real(DP) :: xn
2384  real(DP) :: yn
2385  real(DP) :: zn
2386  real(DP) :: xc
2387  real(DP) :: yc
2388  real(DP) :: zc
2389  real(DP) :: cl1
2390  real(DP) :: cl2
2391  real(DP) :: dltot
2392  real(DP) :: ooclsum
2393  real(DP) :: dsumx
2394  real(DP) :: dsumy
2395  real(DP) :: dsumz
2396  real(DP) :: denom
2397  real(DP) :: area
2398  real(DP) :: dz
2399  real(DP) :: axy
2400  real(DP) :: ayx
2401  logical :: nozee = .true.
2402  type(SpdisWorkArrayType), pointer :: swa => null() !< pointer to spdis work arrays structure
2403  !
2404  ! -- Ensure dis has necessary information
2405  if (this%icalcspdis /= 0 .and. this%dis%con%ianglex == 0) then
2406  call store_error('Error. ANGLDEGX not provided in '// &
2407  'discretization file. ANGLDEGX required for '// &
2408  'calculation of specific discharge.', terminate=.true.)
2409  end if
2410 
2411  swa => this%spdis_wa
2412  if (.not. swa%is_created()) then
2413  ! prepare work arrays
2414  call this%spdis_wa%create(this%calc_max_conns())
2415 
2416  ! prepare lookup table
2417  if (this%nedges > 0) call this%prepare_edge_lookup()
2418  end if
2419  !
2420  ! -- Go through each cell and calculate specific discharge
2421  do n = 1, this%dis%nodes
2422  !
2423  ! -- first calculate geometric properties for x and y directions and
2424  ! the specific discharge at a face (vi)
2425  ic = 0
2426  iz = 0
2427 
2428  ! reset work arrays
2429  call swa%reset()
2430 
2431  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
2432  m = this%dis%con%ja(ipos)
2433  isympos = this%dis%con%jas(ipos)
2434  ihc = this%dis%con%ihc(isympos)
2435  area = this%dis%con%hwva(isympos)
2436  if (ihc == c3d_vertical) then
2437  !
2438  ! -- vertical connection
2439  iz = iz + 1
2440  !call this%dis%connection_normal(n, m, ihc, xn, yn, zn, ipos)
2441  call this%dis%connection_vector(n, m, nozee, this%sat(n), this%sat(m), &
2442  ihc, xc, yc, zc, dltot)
2443  cl1 = this%dis%con%cl1(isympos)
2444  cl2 = this%dis%con%cl2(isympos)
2445  if (m < n) then
2446  cl1 = this%dis%con%cl2(isympos)
2447  cl2 = this%dis%con%cl1(isympos)
2448  end if
2449  ooclsum = done / (cl1 + cl2)
2450  swa%diz(iz) = dltot * cl1 * ooclsum
2451  qz = flowja(ipos)
2452  if (n > m) qz = -qz
2453  swa%viz(iz) = qz / area
2454  else
2455  !
2456  ! -- horizontal connection
2457  ic = ic + 1
2458  dz = thksatnm(this%ibound(n), this%ibound(m), &
2459  this%icelltype(n), this%icelltype(m), &
2460  this%inewton, ihc, &
2461  this%hnew(n), this%hnew(m), this%sat(n), this%sat(m), &
2462  this%dis%top(n), this%dis%top(m), this%dis%bot(n), &
2463  this%dis%bot(m))
2464  area = area * dz
2465  call this%dis%connection_normal(n, m, ihc, xn, yn, zn, ipos)
2466  call this%dis%connection_vector(n, m, nozee, this%sat(n), this%sat(m), &
2467  ihc, xc, yc, zc, dltot)
2468  cl1 = this%dis%con%cl1(isympos)
2469  cl2 = this%dis%con%cl2(isympos)
2470  if (m < n) then
2471  cl1 = this%dis%con%cl2(isympos)
2472  cl2 = this%dis%con%cl1(isympos)
2473  end if
2474  ooclsum = done / (cl1 + cl2)
2475  swa%nix(ic) = -xn
2476  swa%niy(ic) = -yn
2477  swa%di(ic) = dltot * cl1 * ooclsum
2478  if (area > dzero) then
2479  swa%vi(ic) = flowja(ipos) / area
2480  else
2481  swa%vi(ic) = dzero
2482  end if
2483  end if
2484  end do
2485 
2486  ! add contribution from edge flows (i.e. from exchanges)
2487  if (this%nedges > 0) then
2488  do ipos = this%iedge_ptr(n), this%iedge_ptr(n + 1) - 1
2489  iedge = this%edge_idxs(ipos)
2490 
2491  ! propsedge: (Q, area, nx, ny, distance)
2492  ihc = this%ihcedge(iedge)
2493  area = this%propsedge(2, iedge)
2494  if (ihc == c3d_vertical) then
2495  iz = iz + 1
2496  swa%viz(iz) = this%propsedge(1, iedge) / area
2497  swa%diz(iz) = this%propsedge(5, iedge)
2498  else
2499  ic = ic + 1
2500  swa%nix(ic) = -this%propsedge(3, iedge)
2501  swa%niy(ic) = -this%propsedge(4, iedge)
2502  swa%di(ic) = this%propsedge(5, iedge)
2503  if (area > dzero) then
2504  swa%vi(ic) = this%propsedge(1, iedge) / area
2505  else
2506  swa%vi(ic) = dzero
2507  end if
2508  end if
2509  end do
2510  end if
2511  !
2512  ! -- Assign number of vertical and horizontal connections
2513  ncz = iz
2514  nc = ic
2515  !
2516  ! -- calculate z weight (wiz) and z velocity
2517  if (ncz == 1) then
2518  swa%wiz(1) = done
2519  else
2520  dsumz = dzero
2521  do iz = 1, ncz
2522  dsumz = dsumz + swa%diz(iz)
2523  end do
2524  denom = (ncz - done)
2525  if (denom < dzero) denom = dzero
2526  dsumz = dsumz + dem10 * dsumz
2527  do iz = 1, ncz
2528  if (dsumz > dzero) swa%wiz(iz) = done - swa%diz(iz) / dsumz
2529  if (denom > 0) then
2530  swa%wiz(iz) = swa%wiz(iz) / denom
2531  else
2532  swa%wiz(iz) = dzero
2533  end if
2534  end do
2535  end if
2536  vz = dzero
2537  do iz = 1, ncz
2538  vz = vz + swa%wiz(iz) * swa%viz(iz)
2539  end do
2540  !
2541  ! -- distance-based weighting
2542  nc = ic
2543  dsumx = dzero
2544  dsumy = dzero
2545  dsumz = dzero
2546  do ic = 1, nc
2547  swa%wix(ic) = swa%di(ic) * abs(swa%nix(ic))
2548  swa%wiy(ic) = swa%di(ic) * abs(swa%niy(ic))
2549  dsumx = dsumx + swa%wix(ic)
2550  dsumy = dsumy + swa%wiy(ic)
2551  end do
2552  !
2553  ! -- Finish computing omega weights. Add a tiny bit
2554  ! to dsum so that the normalized omega weight later
2555  ! evaluates to (essentially) 1 in the case of a single
2556  ! relevant connection, avoiding 0/0.
2557  dsumx = dsumx + dem10 * dsumx
2558  dsumy = dsumy + dem10 * dsumy
2559  do ic = 1, nc
2560  swa%wix(ic) = (dsumx - swa%wix(ic)) * abs(swa%nix(ic))
2561  swa%wiy(ic) = (dsumy - swa%wiy(ic)) * abs(swa%niy(ic))
2562  end do
2563  !
2564  ! -- compute B weights
2565  dsumx = dzero
2566  dsumy = dzero
2567  do ic = 1, nc
2568  swa%bix(ic) = swa%wix(ic) * sign(done, swa%nix(ic))
2569  swa%biy(ic) = swa%wiy(ic) * sign(done, swa%niy(ic))
2570  dsumx = dsumx + swa%wix(ic) * abs(swa%nix(ic))
2571  dsumy = dsumy + swa%wiy(ic) * abs(swa%niy(ic))
2572  end do
2573  if (dsumx > dzero) dsumx = done / dsumx
2574  if (dsumy > dzero) dsumy = done / dsumy
2575  axy = dzero
2576  ayx = dzero
2577  do ic = 1, nc
2578  swa%bix(ic) = swa%bix(ic) * dsumx
2579  swa%biy(ic) = swa%biy(ic) * dsumy
2580  axy = axy + swa%bix(ic) * swa%niy(ic)
2581  ayx = ayx + swa%biy(ic) * swa%nix(ic)
2582  end do
2583  !
2584  ! -- Calculate specific discharge. The divide by zero checking below
2585  ! is problematic for cells with only one flow, such as can happen
2586  ! with triangular cells in corners. In this case, the resulting
2587  ! cell velocity will be calculated as zero. The method should be
2588  ! improved so that edge flows of zero are included in these
2589  ! calculations. But this needs to be done with consideration for LGR
2590  ! cases in which flows are submitted from an exchange.
2591  vx = dzero
2592  vy = dzero
2593  do ic = 1, nc
2594  vx = vx + (swa%bix(ic) - axy * swa%biy(ic)) * swa%vi(ic)
2595  vy = vy + (swa%biy(ic) - ayx * swa%bix(ic)) * swa%vi(ic)
2596  end do
2597  denom = done - axy * ayx
2598  if (denom /= dzero) then
2599  vx = vx / denom
2600  vy = vy / denom
2601  end if
2602  !
2603  this%spdis(1, n) = vx
2604  this%spdis(2, n) = vy
2605  this%spdis(3, n) = vz
2606  !
2607  end do
2608 
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:

◆ calcsatthickness()

real(dp) function gwfnpfmodule::calcsatthickness ( class(gwfnpftype this,
integer(i4b)  n,
integer(i4b)  m,
integer(i4b)  ihc 
)
private
Parameters
thisthis NPF instance
nnode n
mnode m
ihc1 = horizontal connection, 0 for vertical
Returns
saturated thickness

Definition at line 2784 of file gwf-npf.f90.

2785  ! -- dummy
2786  class(GwfNpfType) :: this !< this NPF instance
2787  integer(I4B) :: n !< node n
2788  integer(I4B) :: m !< node m
2789  integer(I4B) :: ihc !< 1 = horizontal connection, 0 for vertical
2790  ! -- return
2791  real(DP) :: satThickness !< saturated thickness
2792  !
2793  satthickness = thksatnm(this%ibound(n), &
2794  this%ibound(m), &
2795  this%icelltype(n), &
2796  this%icelltype(m), &
2797  this%inewton, &
2798  ihc, &
2799  this%hnew(n), &
2800  this%hnew(m), &
2801  this%sat(n), &
2802  this%sat(m), &
2803  this%dis%top(n), &
2804  this%dis%top(m), &
2805  this%dis%bot(n), &
2806  this%dis%bot(m))

◆ check_options()

subroutine gwfnpfmodule::check_options ( class(gwfnpftype this)
private

Definition at line 1389 of file gwf-npf.f90.

1390  ! -- modules
1392  use constantsmodule, only: linelength
1393  ! -- dummy
1394  class(GwfNpftype) :: this
1395  !
1396  ! -- set omega value used for saturation calculations
1397  if (this%inewton > 0) then
1398  this%satomega = dem6
1399  end if
1400  !
1401  if (this%inewton > 0) then
1402  if (this%iperched > 0) then
1403  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1404  'BE USED WITH PERCHED OPTION.'
1405  call store_error(errmsg)
1406  end if
1407  if (this%ivarcv > 0) then
1408  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1409  'BE USED WITH VARIABLECV OPTION.'
1410  call store_error(errmsg)
1411  end if
1412  if (this%irewet > 0) then
1413  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1414  'BE USED WITH REWET OPTION.'
1415  call store_error(errmsg)
1416  end if
1417  end if
1418  !
1419  if (this%ixt3d /= 0) then
1420  if (this%icellavg > 0) then
1421  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. '// &
1422  'ALTERNATIVE_CELL_AVERAGING OPTION '// &
1423  'CANNOT BE USED WITH XT3D OPTION.'
1424  call store_error(errmsg)
1425  end if
1426  if (this%ithickstrt > 0) then
1427  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. THICKSTRT OPTION '// &
1428  'CANNOT BE USED WITH XT3D OPTION.'
1429  call store_error(errmsg)
1430  end if
1431  if (this%iperched > 0) then
1432  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. PERCHED OPTION '// &
1433  'CANNOT BE USED WITH XT3D OPTION.'
1434  call store_error(errmsg)
1435  end if
1436  if (this%ivarcv > 0) then
1437  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. VARIABLECV OPTION '// &
1438  'CANNOT BE USED WITH XT3D OPTION.'
1439  call store_error(errmsg)
1440  end if
1441  end if
1442  !
1443  ! -- Terminate if errors
1444  if (count_errors() > 0) then
1445  call store_error_filename(this%input_fname)
1446  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
Here is the call graph for this function:

◆ hy_eff()

real(dp) function gwfnpfmodule::hy_eff ( class(gwfnpftype this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
integer(i4b), intent(in)  ihc,
integer(i4b), intent(in), optional  ipos,
real(dp), dimension(3), intent(in), optional  vg 
)

n is primary node node number m is connected node (not used if vg is provided) ihc is horizontal indicator (0 vertical, 1 horizontal, 2 vertically staggered) ipos_opt is position of connection in ja array vg is the global unit vector that expresses the direction from which to calculate an effective hydraulic conductivity.

Definition at line 2283 of file gwf-npf.f90.

2284  ! -- return
2285  real(DP) :: hy
2286  ! -- dummy
2287  class(GwfNpfType) :: this
2288  integer(I4B), intent(in) :: n
2289  integer(I4B), intent(in) :: m
2290  integer(I4B), intent(in) :: ihc
2291  integer(I4B), intent(in), optional :: ipos
2292  real(DP), dimension(3), intent(in), optional :: vg
2293  ! -- local
2294  integer(I4B) :: iipos
2295  real(DP) :: hy11, hy22, hy33
2296  real(DP) :: ang1, ang2, ang3
2297  real(DP) :: vg1, vg2, vg3
2298  !
2299  ! -- Initialize
2300  iipos = 0
2301  if (present(ipos)) iipos = ipos
2302  hy11 = this%k11(n)
2303  hy22 = this%k11(n)
2304  hy33 = this%k11(n)
2305  hy22 = this%k22(n)
2306  hy33 = this%k33(n)
2307  !
2308  ! -- Calculate effective K based on whether connection is vertical
2309  ! or horizontal
2310  if (ihc == c3d_vertical) then
2311  !
2312  ! -- Handle rotated anisotropy case that would affect the effective
2313  ! vertical hydraulic conductivity
2314  hy = hy33
2315  if (this%iangle2 > 0) then
2316  if (present(vg)) then
2317  vg1 = vg(1)
2318  vg2 = vg(2)
2319  vg3 = vg(3)
2320  else
2321  call this%dis%connection_normal(n, m, ihc, vg1, vg2, vg3, iipos)
2322  end if
2323  ang1 = this%angle1(n)
2324  ang2 = this%angle2(n)
2325  ang3 = dzero
2326  if (this%iangle3 > 0) ang3 = this%angle3(n)
2327  hy = hyeff(hy11, hy22, hy33, ang1, ang2, ang3, vg1, vg2, vg3, &
2328  this%iavgkeff)
2329  end if
2330  !
2331  else
2332  !
2333  ! -- Handle horizontal case
2334  hy = hy11
2335  if (this%ik22 > 0) then
2336  if (present(vg)) then
2337  vg1 = vg(1)
2338  vg2 = vg(2)
2339  vg3 = vg(3)
2340  else
2341  call this%dis%connection_normal(n, m, ihc, vg1, vg2, vg3, iipos)
2342  end if
2343  ang1 = dzero
2344  ang2 = dzero
2345  ang3 = dzero
2346  if (this%iangle1 > 0) then
2347  ang1 = this%angle1(n)
2348  if (this%iangle2 > 0) then
2349  ang2 = this%angle2(n)
2350  if (this%iangle3 > 0) ang3 = this%angle3(n)
2351  end if
2352  end if
2353  hy = hyeff(hy11, hy22, hy33, ang1, ang2, ang3, vg1, vg2, vg3, &
2354  this%iavgkeff)
2355  end if
2356  !
2357  end if
Here is the call graph for this function:

◆ increase_edge_count()

subroutine gwfnpfmodule::increase_edge_count ( class(gwfnpftype this,
integer(i4b), intent(in)  nedges 
)
private

This must be called before the npfallocate_arrays routine, which is called from npfar.

Definition at line 2673 of file gwf-npf.f90.

2674  ! -- dummy
2675  class(GwfNpfType) :: this
2676  integer(I4B), intent(in) :: nedges
2677  !
2678  this%nedges = this%nedges + nedges

◆ log_griddata()

subroutine gwfnpfmodule::log_griddata ( class(gwfnpftype this,
type(gwfnpfparamfoundtype), intent(in)  found 
)

Definition at line 1451 of file gwf-npf.f90.

1452  ! -- modules
1454  ! -- dummy
1455  class(GwfNpfType) :: this
1456  type(GwfNpfParamFoundType), intent(in) :: found
1457  !
1458  write (this%iout, '(1x,a)') 'Setting NPF Griddata'
1459  !
1460  if (found%icelltype) then
1461  write (this%iout, '(4x,a)') 'ICELLTYPE set from input file'
1462  end if
1463  !
1464  if (found%k) then
1465  write (this%iout, '(4x,a)') 'K set from input file'
1466  end if
1467  !
1468  if (found%k33) then
1469  write (this%iout, '(4x,a)') 'K33 set from input file'
1470  else
1471  write (this%iout, '(4x,a)') 'K33 not provided. Setting K33 = K.'
1472  end if
1473  !
1474  if (found%k22) then
1475  write (this%iout, '(4x,a)') 'K22 set from input file'
1476  else
1477  write (this%iout, '(4x,a)') 'K22 not provided. Setting K22 = K.'
1478  end if
1479  !
1480  if (found%wetdry) then
1481  write (this%iout, '(4x,a)') 'WETDRY set from input file'
1482  end if
1483  !
1484  if (found%angle1) then
1485  write (this%iout, '(4x,a)') 'ANGLE1 set from input file'
1486  end if
1487  !
1488  if (found%angle2) then
1489  write (this%iout, '(4x,a)') 'ANGLE2 set from input file'
1490  end if
1491  !
1492  if (found%angle3) then
1493  write (this%iout, '(4x,a)') 'ANGLE3 set from input file'
1494  end if
1495  !
1496  write (this%iout, '(1x,a,/)') 'End Setting NPF Griddata'

◆ log_options()

subroutine gwfnpfmodule::log_options ( class(gwfnpftype this,
type(gwfnpfparamfoundtype), intent(in)  found 
)
private

Definition at line 1224 of file gwf-npf.f90.

1225  ! -- modules
1226  use kindmodule, only: lgp
1228  ! -- dummy
1229  class(GwfNpftype) :: this
1230  ! -- locals
1231  type(GwfNpfParamFoundType), intent(in) :: found
1232  !
1233  write (this%iout, '(1x,a)') 'Setting NPF Options'
1234  if (found%iprflow) &
1235  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be printed &
1236  &to listing file whenever ICBCFL is not zero.'
1237  if (found%ipakcb) &
1238  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be saved &
1239  &to binary file whenever ICBCFL is not zero.'
1240  if (found%cellavg) &
1241  write (this%iout, '(4x,a,i0)') 'Alternative cell averaging [1=logarithmic, &
1242  &2=AMT-LMK, 3=AMT-HMK] set to: ', &
1243  this%icellavg
1244  if (found%ithickstrt) &
1245  write (this%iout, '(4x,a)') 'THICKSTRT option has been activated.'
1246  if (found%iperched) &
1247  write (this%iout, '(4x,a)') 'Vertical flow will be adjusted for perched &
1248  &conditions.'
1249  if (found%ivarcv) &
1250  write (this%iout, '(4x,a)') 'Vertical conductance varies with water table.'
1251  if (found%idewatcv) &
1252  write (this%iout, '(4x,a)') 'Vertical conductance is calculated using &
1253  &only the saturated thickness and properties &
1254  &of the overlying cell if the head in the &
1255  &underlying cell is below its top.'
1256  if (found%ixt3d) write (this%iout, '(4x,a)') 'XT3D formulation is selected.'
1257  if (found%ixt3drhs) &
1258  write (this%iout, '(4x,a)') 'XT3D RHS formulation is selected.'
1259  if (found%isavspdis) &
1260  write (this%iout, '(4x,a)') 'Specific discharge will be calculated at cell &
1261  &centers and written to DATA-SPDIS in budget &
1262  &file when requested.'
1263  if (found%isavsat) &
1264  write (this%iout, '(4x,a)') 'Saturation will be written to DATA-SAT in &
1265  &budget file when requested.'
1266  if (found%ik22overk) &
1267  write (this%iout, '(4x,a)') 'Values specified for K22 are anisotropy &
1268  &ratios and will be multiplied by K before &
1269  &being used in calculations.'
1270  if (found%ik33overk) &
1271  write (this%iout, '(4x,a)') 'Values specified for K33 are anisotropy &
1272  &ratios and will be multiplied by K before &
1273  &being used in calculations.'
1274  if (found%inewton) &
1275  write (this%iout, '(4x,a)') 'NEWTON-RAPHSON method disabled for unconfined &
1276  &cells'
1277  if (found%satomega) &
1278  write (this%iout, '(4x,a,1pg15.6)') 'Saturation omega: ', this%satomega
1279  if (found%irewet) &
1280  write (this%iout, '(4x,a)') 'Rewetting is active.'
1281  if (found%wetfct) &
1282  write (this%iout, '(4x,a,1pg15.6)') &
1283  'Wetting factor (WETFCT) has been set to: ', this%wetfct
1284  if (found%iwetit) &
1285  write (this%iout, '(4x,a,i5)') &
1286  'Wetting iteration interval (IWETIT) has been set to: ', this%iwetit
1287  if (found%ihdwet) &
1288  write (this%iout, '(4x,a,i5)') &
1289  'Head rewet equation (IHDWET) has been set to: ', this%ihdwet
1290  write (this%iout, '(1x,a,/)') 'End Setting NPF Options'
This module defines variable data types.
Definition: kind.f90:8

◆ npf_ac()

subroutine gwfnpfmodule::npf_ac ( class(gwfnpftype this,
integer(i4b), intent(in)  moffset,
type(sparsematrix), intent(inout)  sparse 
)

Definition at line 257 of file gwf-npf.f90.

258  ! -- modules
259  use sparsemodule, only: sparsematrix
260  ! -- dummy
261  class(GwfNpftype) :: this
262  integer(I4B), intent(in) :: moffset
263  type(sparsematrix), intent(inout) :: sparse
264  !
265  ! -- Add extended neighbors (neighbors of neighbors)
266  if (this%ixt3d /= 0) call this%xt3d%xt3d_ac(moffset, sparse)

◆ npf_ad()

subroutine gwfnpfmodule::npf_ad ( class(gwfnpftype this,
integer(i4b), intent(in)  nodes,
real(dp), dimension(nodes), intent(inout)  hold,
real(dp), dimension(nodes), intent(inout)  hnew,
integer(i4b), intent(in)  irestore 
)
private

Sets hold (head old) to bot whenever a wettable cell is dry

Definition at line 380 of file gwf-npf.f90.

381  ! -- modules
382  use tdismodule, only: kper, kstp
383  !
384  implicit none
385  ! -- dummy
386  class(GwfNpfType) :: this
387  integer(I4B), intent(in) :: nodes
388  real(DP), dimension(nodes), intent(inout) :: hold
389  real(DP), dimension(nodes), intent(inout) :: hnew
390  integer(I4B), intent(in) :: irestore
391  ! -- local
392  integer(I4B) :: n
393  !
394  ! -- loop through all cells and set hold=bot if wettable cell is dry
395  if (this%irewet > 0) then
396  do n = 1, this%dis%nodes
397  if (this%wetdry(n) == dzero) cycle
398  if (this%ibound(n) /= 0) cycle
399  hold(n) = this%dis%bot(n)
400  end do
401  !
402  ! -- if restore state, then set hnew to DRY if it is a dry wettable cell
403  do n = 1, this%dis%nodes
404  if (this%wetdry(n) == dzero) cycle
405  if (this%ibound(n) /= 0) cycle
406  hnew(n) = dhdry
407  end do
408  end if
409  !
410  ! -- TVK
411  if (this%intvk /= 0) then
412  call this%tvk%ad()
413  end if
414  !
415  ! -- VSC
416  ! -- Hit the TVK-updated K's with VSC correction before calling/updating condsat
417  if (this%invsc /= 0) then
418  call this%vsc%update_k_with_vsc()
419  end if
420  !
421  ! -- If any K values have changed, we need to update CONDSAT or XT3D arrays
422  if (this%kchangeper == kper .and. this%kchangestp == kstp) then
423  if (this%ixt3d == 0) then
424  !
425  ! -- Update the saturated conductance for all connections
426  ! -- of the affected nodes
427  do n = 1, this%dis%nodes
428  if (this%nodekchange(n) == 1) then
429  call this%calc_condsat(n, .false.)
430  end if
431  end do
432  else
433  !
434  ! -- Recompute XT3D coefficients for permanently confined connections
435  if (this%xt3d%lamatsaved .and. .not. this%xt3d%ldispersion) then
436  call this%xt3d%xt3d_fcpc(this%dis%nodes, .true.)
437  end if
438  end if
439  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

◆ npf_ar()

subroutine gwfnpfmodule::npf_ar ( class(gwfnpftype this,
type(gwfictype), intent(in), pointer  ic,
type(gwfvsctype), intent(in), pointer  vsc,
integer(i4b), dimension(:), intent(inout), pointer, contiguous  ibound,
real(dp), dimension(:), intent(inout), pointer, contiguous  hnew 
)
private

Allocate remaining package arrays, preprocess the input data and call *_ar on xt3d, when active.

Parameters
thisinstance of the NPF package
[in]icinitial conditions
[in]vscviscosity package
[in,out]iboundmodel ibound array
[in,out]hnewpointer to model head array

Definition at line 285 of file gwf-npf.f90.

286  ! -- modules
288  ! -- dummy
289  class(GwfNpftype) :: this !< instance of the NPF package
290  type(GwfIcType), pointer, intent(in) :: ic !< initial conditions
291  type(GwfVscType), pointer, intent(in) :: vsc !< viscosity package
292  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: ibound !< model ibound array
293  real(DP), dimension(:), pointer, contiguous, intent(inout) :: hnew !< pointer to model head array
294  ! -- local
295  integer(I4B) :: n
296  !
297  ! -- Store pointers to arguments that were passed in
298  this%ic => ic
299  this%ibound => ibound
300  this%hnew => hnew
301  !
302  if (this%icalcspdis == 1) then
303  call mem_reallocate(this%spdis, 3, this%dis%nodes, 'SPDIS', this%memoryPath)
304  call mem_reallocate(this%nodedge, this%nedges, 'NODEDGE', this%memoryPath)
305  call mem_reallocate(this%ihcedge, this%nedges, 'IHCEDGE', this%memoryPath)
306  call mem_reallocate(this%propsedge, 5, this%nedges, 'PROPSEDGE', &
307  this%memoryPath)
308  call mem_reallocate(this%iedge_ptr, this%dis%nodes + 1, &
309  'NREDGESNODE', this%memoryPath)
310  call mem_reallocate(this%edge_idxs, this%nedges, &
311  'EDGEIDXS', this%memoryPath)
312 
313  do n = 1, this%nedges
314  this%edge_idxs(n) = 0
315  end do
316  do n = 1, this%dis%nodes
317  this%iedge_ptr(n) = 0
318  this%spdis(:, n) = dzero
319  end do
320  end if
321  !
322  ! -- Store pointer to VSC if active
323  if (this%invsc /= 0) then
324  this%vsc => vsc
325  end if
326  !
327  ! -- allocate arrays to store original user input in case TVK/VSC modify them
328  if (this%invsc > 0) then
329  !
330  ! -- Reallocate arrays that store user-input values.
331  call mem_reallocate(this%k11input, this%dis%nodes, 'K11INPUT', &
332  this%memoryPath)
333  call mem_reallocate(this%k22input, this%dis%nodes, 'K22INPUT', &
334  this%memoryPath)
335  call mem_reallocate(this%k33input, this%dis%nodes, 'K33INPUT', &
336  this%memoryPath)
337  ! Allocate arrays that will store the original K values. When VSC active,
338  ! the current Kxx arrays carry the viscosity-adjusted K values.
339  ! This approach leverages existing functionality that makes use of K.
340  call this%store_original_k_arrays(this%dis%nodes, this%dis%njas)
341  end if
342  !
343  ! -- preprocess data
344  call this%preprocess_input()
345  !
346  ! -- xt3d
347  if (this%ixt3d /= 0) then
348  call this%xt3d%xt3d_ar(ibound, this%k11, this%ik33, this%k33, &
349  this%sat, this%ik22, this%k22, &
350  this%iangle1, this%iangle2, this%iangle3, &
351  this%angle1, this%angle2, this%angle3, &
352  this%inewton, this%icelltype)
353  end if
354  !
355  ! -- TVK
356  if (this%intvk /= 0) then
357  call this%tvk%ar(this%dis)
358  end if

◆ npf_cf()

subroutine gwfnpfmodule::npf_cf ( class(gwfnpftype this,
integer(i4b)  kiter,
integer(i4b), intent(in)  nodes,
real(dp), dimension(nodes), intent(inout)  hnew 
)

Definition at line 444 of file gwf-npf.f90.

445  ! -- dummy
446  class(GwfNpfType) :: this
447  integer(I4B) :: kiter
448  integer(I4B), intent(in) :: nodes
449  real(DP), intent(inout), dimension(nodes) :: hnew
450  ! -- local
451  integer(I4B) :: n
452  real(DP) :: satn
453  !
454  ! -- Perform wetting and drying
455  if (this%inewton /= 1) then
456  call this%wd(kiter, hnew)
457  end if
458  !
459  ! -- Calculate saturation for convertible cells
460  do n = 1, this%dis%nodes
461  if (this%icelltype(n) /= 0) then
462  if (this%ibound(n) == 0) then
463  satn = dzero
464  else
465  call this%thksat(n, hnew(n), satn)
466  end if
467  this%sat(n) = satn
468  end if
469  end do

◆ npf_cq()

subroutine gwfnpfmodule::npf_cq ( class(gwfnpftype this,
real(dp), dimension(:), intent(inout)  hnew,
real(dp), dimension(:), intent(inout)  flowja 
)
private

Definition at line 745 of file gwf-npf.f90.

746  ! -- dummy
747  class(GwfNpfType) :: this
748  real(DP), intent(inout), dimension(:) :: hnew
749  real(DP), intent(inout), dimension(:) :: flowja
750  ! -- local
751  integer(I4B) :: n, ipos, m
752  real(DP) :: qnm
753  !
754  ! -- Calculate the flow across each cell face and store in flowja
755  !
756  if (this%ixt3d /= 0) then
757  call this%xt3d%xt3d_flowja(hnew, flowja)
758  else
759  !
760  do n = 1, this%dis%nodes
761  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
762  m = this%dis%con%ja(ipos)
763  if (m < n) cycle
764  call this%qcalc(n, m, hnew(n), hnew(m), ipos, qnm)
765  flowja(ipos) = qnm
766  flowja(this%dis%con%isym(ipos)) = -qnm
767  end do
768  end do
769  !
770  end if

◆ npf_cr()

subroutine, public gwfnpfmodule::npf_cr ( type(gwfnpftype), pointer  npfobj,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  input_mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Definition at line 157 of file gwf-npf.f90.

158  ! -- modules
159  use kindmodule, only: lgp
161  ! -- dummy
162  type(GwfNpfType), pointer :: npfobj
163  character(len=*), intent(in) :: name_model
164  character(len=*), intent(in) :: input_mempath
165  integer(I4B), intent(in) :: inunit
166  integer(I4B), intent(in) :: iout
167  ! -- formats
168  character(len=*), parameter :: fmtheader = &
169  "(1x, /1x, 'NPF -- NODE PROPERTY FLOW PACKAGE, VERSION 1, 3/30/2015', &
170  &' INPUT READ FROM MEMPATH: ', A, /)"
171  !
172  ! -- Create the object
173  allocate (npfobj)
174  !
175  ! -- create name and memory path
176  call npfobj%set_names(1, name_model, 'NPF', 'NPF', input_mempath)
177  !
178  ! -- Allocate scalars
179  call npfobj%allocate_scalars()
180  !
181  ! -- Set variables
182  npfobj%inunit = inunit
183  npfobj%iout = iout
184  !
185  ! -- check if npf is enabled
186  if (inunit > 0) then
187  !
188  ! -- Print a message identifying the node property flow package.
189  write (iout, fmtheader) input_mempath
190  end if
191 
192  ! allocate spdis structure
193  allocate (npfobj%spdis_wa)
194 
Here is the caller graph for this function:

◆ npf_da()

subroutine gwfnpfmodule::npf_da ( class(gwfnpftype this)

Definition at line 945 of file gwf-npf.f90.

946  ! -- modules
949  ! -- dummy
950  class(GwfNpftype) :: this
951 
952  ! free spdis work structure
953  if (this%icalcspdis == 1 .and. this%spdis_wa%is_created()) &
954  call this%spdis_wa%destroy()
955  deallocate (this%spdis_wa)
956  !
957  ! -- Deallocate input memory
958  call memorystore_remove(this%name_model, 'NPF', idm_context)
959  !
960  ! -- TVK
961  if (this%intvk /= 0) then
962  call this%tvk%da()
963  deallocate (this%tvk)
964  end if
965  !
966  ! -- VSC
967  if (this%invsc /= 0) then
968  nullify (this%vsc)
969  end if
970  !
971  ! -- Strings
972  !
973  ! -- Scalars
974  call mem_deallocate(this%iname)
975  call mem_deallocate(this%ixt3d)
976  call mem_deallocate(this%ixt3drhs)
977  call mem_deallocate(this%satomega)
978  call mem_deallocate(this%hnoflo)
979  call mem_deallocate(this%hdry)
980  call mem_deallocate(this%icellavg)
981  call mem_deallocate(this%iavgkeff)
982  call mem_deallocate(this%ik22)
983  call mem_deallocate(this%ik33)
984  call mem_deallocate(this%iperched)
985  call mem_deallocate(this%ivarcv)
986  call mem_deallocate(this%idewatcv)
987  call mem_deallocate(this%ithickstrt)
988  call mem_deallocate(this%isavspdis)
989  call mem_deallocate(this%isavsat)
990  call mem_deallocate(this%icalcspdis)
991  call mem_deallocate(this%irewet)
992  call mem_deallocate(this%wetfct)
993  call mem_deallocate(this%iwetit)
994  call mem_deallocate(this%ihdwet)
995  call mem_deallocate(this%ibotnode)
996  call mem_deallocate(this%iwetdry)
997  call mem_deallocate(this%iangle1)
998  call mem_deallocate(this%iangle2)
999  call mem_deallocate(this%iangle3)
1000  call mem_deallocate(this%nedges)
1001  call mem_deallocate(this%lastedge)
1002  call mem_deallocate(this%ik22overk)
1003  call mem_deallocate(this%ik33overk)
1004  call mem_deallocate(this%intvk)
1005  call mem_deallocate(this%invsc)
1006  call mem_deallocate(this%kchangeper)
1007  call mem_deallocate(this%kchangestp)
1008  !
1009  ! -- Deallocate arrays
1010  deallocate (this%aname)
1011  call mem_deallocate(this%ithickstartflag)
1012  call mem_deallocate(this%icelltype)
1013  call mem_deallocate(this%k11)
1014  call mem_deallocate(this%k22)
1015  call mem_deallocate(this%k33)
1016  call mem_deallocate(this%k11input)
1017  call mem_deallocate(this%k22input)
1018  call mem_deallocate(this%k33input)
1019  call mem_deallocate(this%sat, 'SAT', this%memoryPath)
1020  call mem_deallocate(this%condsat)
1021  call mem_deallocate(this%wetdry)
1022  call mem_deallocate(this%angle1)
1023  call mem_deallocate(this%angle2)
1024  call mem_deallocate(this%angle3)
1025  call mem_deallocate(this%nodedge)
1026  call mem_deallocate(this%ihcedge)
1027  call mem_deallocate(this%propsedge)
1028  call mem_deallocate(this%iedge_ptr)
1029  call mem_deallocate(this%edge_idxs)
1030  call mem_deallocate(this%spdis, 'SPDIS', this%memoryPath)
1031  call mem_deallocate(this%nodekchange)
1032  !
1033  ! -- deallocate parent
1034  call this%NumericalPackageType%da()
1035 
1036  ! pointers
1037  this%hnew => null()
1038 
subroutine, public memorystore_remove(component, subcomponent, context)
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
Here is the call graph for this function:

◆ npf_df()

subroutine gwfnpfmodule::npf_df ( class(gwfnpftype this,
class(disbasetype), intent(inout), pointer  dis,
type(xt3dtype), pointer  xt3d,
integer(i4b), intent(in)  ingnc,
integer(i4b), intent(in)  invsc,
type(gwfnpfoptionstype), intent(in), optional  npf_options 
)

This is a hybrid routine: it either reads the options for this package from the input file, or the optional argument

Parameters
npf_optionsshould be passed. A consistency check is performed, and finally xt3d_df is called, when enabled.
thisinstance of the NPF package
[in,out]disthe pointer to the discretization
xt3dthe pointer to the XT3D 'package'
[in]ingncghostnodes enabled? (>0 means yes)
[in]invscviscosity enabled? (>0 means yes)
[in]npf_optionsthe optional options, for when not constructing from file

Definition at line 204 of file gwf-npf.f90.

205  ! -- modules
206  use simmodule, only: store_error
207  use xt3dmodule, only: xt3d_cr
208  ! -- dummy
209  class(GwfNpftype) :: this !< instance of the NPF package
210  class(DisBaseType), pointer, intent(inout) :: dis !< the pointer to the discretization
211  type(Xt3dType), pointer :: xt3d !< the pointer to the XT3D 'package'
212  integer(I4B), intent(in) :: ingnc !< ghostnodes enabled? (>0 means yes)
213  integer(I4B), intent(in) :: invsc !< viscosity enabled? (>0 means yes)
214  type(GwfNpfOptionsType), optional, intent(in) :: npf_options !< the optional options, for when not constructing from file
215  !
216  ! -- Set a pointer to dis
217  this%dis => dis
218  !
219  ! -- Set flag signifying whether vsc is active
220  if (invsc > 0) this%invsc = invsc
221  !
222  if (.not. present(npf_options)) then
223  !
224  ! -- source options
225  call this%source_options()
226  !
227  ! -- allocate arrays
228  call this%allocate_arrays(this%dis%nodes, this%dis%njas)
229  !
230  ! -- source griddata, set, and convert/check the input
231  call this%source_griddata()
232  call this%prepcheck()
233  else
234  call this%set_options(npf_options)
235  !
236  ! -- allocate arrays
237  call this%allocate_arrays(this%dis%nodes, this%dis%njas)
238  end if
239  !
240  call this%check_options()
241  !
242  ! -- Save pointer to xt3d object
243  this%xt3d => xt3d
244  if (this%ixt3d /= 0) xt3d%ixt3d = this%ixt3d
245  call this%xt3d%xt3d_df(dis)
246  !
247  ! -- Ensure GNC and XT3D are not both on at the same time
248  if (this%ixt3d /= 0 .and. ingnc > 0) then
249  call store_error('Error in model '//trim(this%name_model)// &
250  '. The XT3D option cannot be used with the GNC &
251  &Package.', terminate=.true.)
252  end if
subroutine, public xt3d_cr(xt3dobj, name_model, inunit, iout, ldispopt)
Create a new xt3d object.
Here is the call graph for this function:

◆ npf_fc()

subroutine gwfnpfmodule::npf_fc ( class(gwfnpftype 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 
)
private

Definition at line 474 of file gwf-npf.f90.

475  ! -- modules
476  use constantsmodule, only: done
477  ! -- dummy
478  class(GwfNpfType) :: this
479  integer(I4B) :: kiter
480  class(MatrixBaseType), pointer :: matrix_sln
481  integer(I4B), intent(in), dimension(:) :: idxglo
482  real(DP), intent(inout), dimension(:) :: rhs
483  real(DP), intent(inout), dimension(:) :: hnew
484  ! -- local
485  integer(I4B) :: n, m, ii, idiag, ihc
486  integer(I4B) :: isymcon, idiagm
487  real(DP) :: hyn, hym
488  real(DP) :: cond
489  !
490  ! -- Calculate conductance and put into amat
491  !
492  if (this%ixt3d /= 0) then
493  call this%xt3d%xt3d_fc(kiter, matrix_sln, idxglo, rhs, hnew)
494  else
495  do n = 1, this%dis%nodes
496  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
497  if (this%dis%con%mask(ii) == 0) cycle
498 
499  m = this%dis%con%ja(ii)
500  !
501  ! -- Calculate conductance only for upper triangle but insert into
502  ! upper and lower parts of amat.
503  if (m < n) cycle
504  ihc = this%dis%con%ihc(this%dis%con%jas(ii))
505  hyn = this%hy_eff(n, m, ihc, ipos=ii)
506  hym = this%hy_eff(m, n, ihc, ipos=ii)
507  !
508  ! -- Vertical connection
509  if (ihc == c3d_vertical) then
510  !
511  ! -- Calculate vertical conductance
512  cond = vcond(this%ibound(n), this%ibound(m), &
513  this%icelltype(n), this%icelltype(m), this%inewton, &
514  this%ivarcv, this%idewatcv, &
515  this%condsat(this%dis%con%jas(ii)), hnew(n), hnew(m), &
516  hyn, hym, &
517  this%sat(n), this%sat(m), &
518  this%dis%top(n), this%dis%top(m), &
519  this%dis%bot(n), this%dis%bot(m), &
520  this%dis%con%hwva(this%dis%con%jas(ii)))
521  !
522  ! -- Vertical flow for perched conditions
523  if (this%iperched /= 0) then
524  if (this%icelltype(m) /= 0) then
525  if (hnew(m) < this%dis%top(m)) then
526  !
527  ! -- Fill row n
528  idiag = this%dis%con%ia(n)
529  rhs(n) = rhs(n) - cond * this%dis%bot(n)
530  call matrix_sln%add_value_pos(idxglo(idiag), -cond)
531  !
532  ! -- Fill row m
533  isymcon = this%dis%con%isym(ii)
534  call matrix_sln%add_value_pos(idxglo(isymcon), cond)
535  rhs(m) = rhs(m) + cond * this%dis%bot(n)
536  !
537  ! -- cycle the connection loop
538  cycle
539  end if
540  end if
541  end if
542  !
543  else
544  !
545  ! -- Horizontal conductance
546  cond = hcond(this%ibound(n), this%ibound(m), &
547  this%icelltype(n), this%icelltype(m), &
548  this%inewton, &
549  this%dis%con%ihc(this%dis%con%jas(ii)), &
550  this%icellavg, &
551  this%condsat(this%dis%con%jas(ii)), &
552  hnew(n), hnew(m), this%sat(n), this%sat(m), hyn, hym, &
553  this%dis%top(n), this%dis%top(m), &
554  this%dis%bot(n), this%dis%bot(m), &
555  this%dis%con%cl1(this%dis%con%jas(ii)), &
556  this%dis%con%cl2(this%dis%con%jas(ii)), &
557  this%dis%con%hwva(this%dis%con%jas(ii)))
558  end if
559  !
560  ! -- Fill row n
561  idiag = this%dis%con%ia(n)
562  call matrix_sln%add_value_pos(idxglo(ii), cond)
563  call matrix_sln%add_value_pos(idxglo(idiag), -cond)
564  !
565  ! -- Fill row m
566  isymcon = this%dis%con%isym(ii)
567  idiagm = this%dis%con%ia(m)
568  call matrix_sln%add_value_pos(idxglo(isymcon), cond)
569  call matrix_sln%add_value_pos(idxglo(idiagm), -cond)
570  end do
571  end do
572  !
573  end if
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
Here is the call graph for this function:

◆ npf_fn()

subroutine gwfnpfmodule::npf_fn ( class(gwfnpftype 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 578 of file gwf-npf.f90.

579  ! -- dummy
580  class(GwfNpfType) :: this
581  integer(I4B) :: kiter
582  class(MatrixBaseType), pointer :: matrix_sln
583  integer(I4B), intent(in), dimension(:) :: idxglo
584  real(DP), intent(inout), dimension(:) :: rhs
585  real(DP), intent(inout), dimension(:) :: hnew
586  ! -- local
587  integer(I4B) :: nodes, nja
588  integer(I4B) :: n, m, ii, idiag
589  integer(I4B) :: isymcon, idiagm
590  integer(I4B) :: iups
591  integer(I4B) :: idn
592  real(DP) :: cond
593  real(DP) :: consterm
594  real(DP) :: filledterm
595  real(DP) :: derv
596  real(DP) :: hds
597  real(DP) :: term
598  real(DP) :: topup
599  real(DP) :: botup
600  !
601  ! -- add newton terms to solution matrix
602  nodes = this%dis%nodes
603  nja = this%dis%con%nja
604  if (this%ixt3d /= 0) then
605  call this%xt3d%xt3d_fn(kiter, nodes, nja, matrix_sln, idxglo, rhs, hnew)
606  else
607  !
608  do n = 1, nodes
609  idiag = this%dis%con%ia(n)
610  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
611  if (this%dis%con%mask(ii) == 0) cycle
612 
613  m = this%dis%con%ja(ii)
614  isymcon = this%dis%con%isym(ii)
615  ! work on upper triangle
616  if (m < n) cycle
617  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0 .and. &
618  this%ivarcv == 0) then
619  !call this%vcond(n,m,hnew(n),hnew(m),ii,cond)
620  ! do nothing
621  else
622  ! determine upstream node
623  iups = m
624  if (hnew(m) < hnew(n)) iups = n
625  idn = n
626  if (iups == n) idn = m
627  !
628  ! -- no newton terms if upstream cell is confined
629  if (this%icelltype(iups) == 0) cycle
630  !
631  ! -- Set the upstream top and bot, and then recalculate for a
632  ! vertically staggered horizontal connection
633  topup = this%dis%top(iups)
634  botup = this%dis%bot(iups)
635  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 2) then
636  topup = min(this%dis%top(n), this%dis%top(m))
637  botup = max(this%dis%bot(n), this%dis%bot(m))
638  end if
639  !
640  ! get saturated conductivity for derivative
641  cond = this%condsat(this%dis%con%jas(ii))
642  !
643  ! compute additional term
644  consterm = -cond * (hnew(iups) - hnew(idn)) !needs to use hwadi instead of hnew(idn)
645  !filledterm = cond
646  filledterm = matrix_sln%get_value_pos(idxglo(ii))
647  derv = squadraticsaturationderivative(topup, botup, hnew(iups), &
648  this%satomega)
649  idiagm = this%dis%con%ia(m)
650  ! fill jacobian for n being the upstream node
651  if (iups == n) then
652  hds = hnew(m)
653  !isymcon = this%dis%con%isym(ii)
654  term = consterm * derv
655  rhs(n) = rhs(n) + term * hnew(n) !+ amat(idxglo(isymcon)) * (dwadi * hds - hds) !need to add dwadi
656  rhs(m) = rhs(m) - term * hnew(n) !- amat(idxglo(isymcon)) * (dwadi * hds - hds) !need to add dwadi
657  ! fill in row of n
658  call matrix_sln%add_value_pos(idxglo(idiag), term)
659  ! fill newton term in off diagonal if active cell
660  if (this%ibound(n) > 0) then
661  filledterm = matrix_sln%get_value_pos(idxglo(ii))
662  call matrix_sln%set_value_pos(idxglo(ii), filledterm) !* dwadi !need to add dwadi
663  end if
664  !fill row of m
665  filledterm = matrix_sln%get_value_pos(idxglo(idiagm))
666  call matrix_sln%set_value_pos(idxglo(idiagm), filledterm) !- filledterm * (dwadi - DONE) !need to add dwadi
667  ! fill newton term in off diagonal if active cell
668  if (this%ibound(m) > 0) then
669  call matrix_sln%add_value_pos(idxglo(isymcon), -term)
670  end if
671  ! fill jacobian for m being the upstream node
672  else
673  hds = hnew(n)
674  term = -consterm * derv
675  rhs(n) = rhs(n) + term * hnew(m) !+ amat(idxglo(ii)) * (dwadi * hds - hds) !need to add dwadi
676  rhs(m) = rhs(m) - term * hnew(m) !- amat(idxglo(ii)) * (dwadi * hds - hds) !need to add dwadi
677  ! fill in row of n
678  filledterm = matrix_sln%get_value_pos(idxglo(idiag))
679  call matrix_sln%set_value_pos(idxglo(idiag), filledterm) !- filledterm * (dwadi - DONE) !need to add dwadi
680  ! fill newton term in off diagonal if active cell
681  if (this%ibound(n) > 0) then
682  call matrix_sln%add_value_pos(idxglo(ii), term)
683  end if
684  !fill row of m
685  call matrix_sln%add_value_pos(idxglo(idiagm), -term)
686  ! fill newton term in off diagonal if active cell
687  if (this%ibound(m) > 0) then
688  filledterm = matrix_sln%get_value_pos(idxglo(isymcon))
689  call matrix_sln%set_value_pos(idxglo(isymcon), filledterm) !* dwadi !need to add dwadi
690  end if
691  end if
692  end if
693 
694  end do
695  end do
696  !
697  end if
Here is the call graph for this function:

◆ npf_mc()

subroutine gwfnpfmodule::npf_mc ( class(gwfnpftype this,
integer(i4b), intent(in)  moffset,
class(matrixbasetype), pointer  matrix_sln 
)

Definition at line 271 of file gwf-npf.f90.

272  ! -- dummy
273  class(GwfNpftype) :: this
274  integer(I4B), intent(in) :: moffset
275  class(MatrixBaseType), pointer :: matrix_sln
276  !
277  if (this%ixt3d /= 0) call this%xt3d%xt3d_mc(moffset, matrix_sln)

◆ npf_nur()

subroutine gwfnpfmodule::npf_nur ( class(gwfnpftype this,
integer(i4b), intent(in)  neqmod,
real(dp), dimension(neqmod), intent(inout)  x,
real(dp), dimension(neqmod), intent(in)  xtemp,
real(dp), dimension(neqmod), intent(inout)  dx,
integer(i4b), intent(inout)  inewtonur,
real(dp), intent(inout)  dxmax,
integer(i4b), intent(inout)  locmax 
)
private

Under-relaxation of Groundwater Flow Model Heads for current outer iteration using the cell bottoms at the bottom of the model

Definition at line 705 of file gwf-npf.f90.

706  ! -- dummy
707  class(GwfNpfType) :: this
708  integer(I4B), intent(in) :: neqmod
709  real(DP), dimension(neqmod), intent(inout) :: x
710  real(DP), dimension(neqmod), intent(in) :: xtemp
711  real(DP), dimension(neqmod), intent(inout) :: dx
712  integer(I4B), intent(inout) :: inewtonur
713  real(DP), intent(inout) :: dxmax
714  integer(I4B), intent(inout) :: locmax
715  ! -- local
716  integer(I4B) :: n
717  real(DP) :: botm
718  real(DP) :: xx
719  real(DP) :: dxx
720  !
721  ! -- Newton-Raphson under-relaxation
722  do n = 1, this%dis%nodes
723  if (this%ibound(n) < 1) cycle
724  if (this%icelltype(n) > 0) then
725  botm = this%dis%bot(this%ibotnode(n))
726  ! -- only apply Newton-Raphson under-relaxation if
727  ! solution head is below the bottom of the model
728  if (x(n) < botm) then
729  inewtonur = 1
730  xx = xtemp(n) * (done - dp9) + botm * dp9
731  dxx = x(n) - xx
732  if (abs(dxx) > abs(dxmax)) then
733  locmax = n
734  dxmax = dxx
735  end if
736  x(n) = xx
737  dx(n) = dzero
738  end if
739  end if
740  end do

◆ npf_print_model_flows()

subroutine gwfnpfmodule::npf_print_model_flows ( class(gwfnpftype this,
integer(i4b), intent(in)  ibudfl,
real(dp), dimension(:), intent(inout)  flowja 
)
private

Definition at line 906 of file gwf-npf.f90.

907  ! -- modules
908  use tdismodule, only: kper, kstp
909  use constantsmodule, only: lenbigline
910  ! -- dummy
911  class(GwfNpfType) :: this
912  integer(I4B), intent(in) :: ibudfl
913  real(DP), intent(inout), dimension(:) :: flowja
914  ! -- local
915  character(len=LENBIGLINE) :: line
916  character(len=30) :: tempstr
917  integer(I4B) :: n, ipos, m
918  real(DP) :: qnm
919  ! -- formats
920  character(len=*), parameter :: fmtiprflow = &
921  &"(/,4x,'CALCULATED INTERCELL FLOW FOR PERIOD ', i0, ' STEP ', i0)"
922  !
923  ! -- Write flowja to list file if requested
924  if (ibudfl /= 0 .and. this%iprflow > 0) then
925  write (this%iout, fmtiprflow) kper, kstp
926  do n = 1, this%dis%nodes
927  line = ''
928  call this%dis%noder_to_string(n, tempstr)
929  line = trim(tempstr)//':'
930  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
931  m = this%dis%con%ja(ipos)
932  call this%dis%noder_to_string(m, tempstr)
933  line = trim(line)//' '//trim(tempstr)
934  qnm = flowja(ipos)
935  write (tempstr, '(1pg15.6)') qnm
936  line = trim(line)//' '//trim(adjustl(tempstr))
937  end do
938  write (this%iout, '(a)') trim(line)
939  end do
940  end if
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15

◆ npf_rp()

subroutine gwfnpfmodule::npf_rp ( class(gwfnpftype this)

Read and prepare NPF stress period data.

Definition at line 365 of file gwf-npf.f90.

366  implicit none
367  ! -- dummy
368  class(GwfNpfType) :: this
369  !
370  ! -- TVK
371  if (this%intvk /= 0) then
372  call this%tvk%rp()
373  end if

◆ npf_save_model_flows()

subroutine gwfnpfmodule::npf_save_model_flows ( class(gwfnpftype this,
real(dp), dimension(:), intent(in)  flowja,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  icbcun 
)
private

Definition at line 869 of file gwf-npf.f90.

870  ! -- dummy
871  class(GwfNpfType) :: this
872  real(DP), dimension(:), intent(in) :: flowja
873  integer(I4B), intent(in) :: icbcfl
874  integer(I4B), intent(in) :: icbcun
875  ! -- local
876  integer(I4B) :: ibinun
877  !
878  ! -- Set unit number for binary output
879  if (this%ipakcb < 0) then
880  ibinun = icbcun
881  elseif (this%ipakcb == 0) then
882  ibinun = 0
883  else
884  ibinun = this%ipakcb
885  end if
886  if (icbcfl == 0) ibinun = 0
887  !
888  ! -- Write the face flows if requested
889  if (ibinun /= 0) then
890  call this%dis%record_connection_array(flowja, ibinun, this%iout)
891  end if
892  !
893  ! -- Calculate specific discharge at cell centers and write, if requested
894  if (this%isavspdis /= 0) then
895  if (ibinun /= 0) call this%sav_spdis(ibinun)
896  end if
897  !
898  ! -- Save saturation, if requested
899  if (this%isavsat /= 0) then
900  if (ibinun /= 0) call this%sav_sat(ibinun)
901  end if

◆ prepare_edge_lookup()

subroutine gwfnpfmodule::prepare_edge_lookup ( class(gwfnpftype this)
private

Definition at line 2739 of file gwf-npf.f90.

2740  class(GwfNpfType) :: this
2741  ! local
2742  integer(I4B) :: i, inode, iedge
2743  integer(I4B) :: n, start, end
2744  integer(I4B) :: prev_cnt, strt_idx, ipos
2745 
2746  do i = 1, size(this%iedge_ptr)
2747  this%iedge_ptr(i) = 0
2748  end do
2749  do i = 1, size(this%edge_idxs)
2750  this%edge_idxs(i) = 0
2751  end do
2752 
2753  ! count
2754  do iedge = 1, this%nedges
2755  n = this%nodedge(iedge)
2756  this%iedge_ptr(n) = this%iedge_ptr(n) + 1
2757  end do
2758 
2759  ! determine start indexes
2760  prev_cnt = this%iedge_ptr(1)
2761  this%iedge_ptr(1) = 1
2762  do inode = 2, this%dis%nodes + 1
2763  strt_idx = this%iedge_ptr(inode - 1) + prev_cnt
2764  prev_cnt = this%iedge_ptr(inode)
2765  this%iedge_ptr(inode) = strt_idx
2766  end do
2767 
2768  ! loop over edges to fill lookup table
2769  do iedge = 1, this%nedges
2770  n = this%nodedge(iedge)
2771  start = this%iedge_ptr(n)
2772  end = this%iedge_ptr(n + 1) - 1
2773  do ipos = start, end
2774  if (this%edge_idxs(ipos) > 0) cycle ! go to next
2775  this%edge_idxs(ipos) = iedge
2776  exit
2777  end do
2778  end do
2779 

◆ prepcheck()

subroutine gwfnpfmodule::prepcheck ( class(gwfnpftype this)

Definition at line 1590 of file gwf-npf.f90.

1591  ! -- modules
1592  use constantsmodule, only: linelength, dpio180
1594  ! -- dummy
1595  class(GwfNpfType) :: this
1596  ! -- local
1597  character(len=24), dimension(:), pointer :: aname
1598  character(len=LINELENGTH) :: cellstr, errmsg
1599  integer(I4B) :: nerr, n
1600  ! -- format
1601  character(len=*), parameter :: fmtkerr = &
1602  &"(1x, 'Hydraulic property ',a,' is <= 0 for cell ',a, ' ', 1pg15.6)"
1603  character(len=*), parameter :: fmtkerr2 = &
1604  &"(1x, '... ', i0,' additional errors not shown for ',a)"
1605  !
1606  ! -- initialize
1607  aname => this%aname
1608  !
1609  ! -- check k11
1610  nerr = 0
1611  do n = 1, size(this%k11)
1612  if (this%k11(n) <= dzero) then
1613  nerr = nerr + 1
1614  if (nerr <= 20) then
1615  call this%dis%noder_to_string(n, cellstr)
1616  write (errmsg, fmtkerr) trim(adjustl(aname(2))), trim(cellstr), &
1617  this%k11(n)
1618  call store_error(errmsg)
1619  end if
1620  end if
1621  end do
1622  if (nerr > 20) then
1623  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(2)))
1624  call store_error(errmsg)
1625  end if
1626  !
1627  ! -- check k33 because it was read
1628  if (this%ik33 /= 0) then
1629  !
1630  ! -- Check to make sure values are greater than or equal to zero
1631  nerr = 0
1632  do n = 1, size(this%k33)
1633  if (this%ik33overk /= 0) this%k33(n) = this%k33(n) * this%k11(n)
1634  if (this%k33(n) <= dzero) then
1635  nerr = nerr + 1
1636  if (nerr <= 20) then
1637  call this%dis%noder_to_string(n, cellstr)
1638  write (errmsg, fmtkerr) trim(adjustl(aname(3))), trim(cellstr), &
1639  this%k33(n)
1640  call store_error(errmsg)
1641  end if
1642  end if
1643  end do
1644  if (nerr > 20) then
1645  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(3)))
1646  call store_error(errmsg)
1647  end if
1648  end if
1649  !
1650  ! -- check k22 because it was read
1651  if (this%ik22 /= 0) then
1652  !
1653  ! -- Check to make sure that angles are available
1654  if (this%dis%con%ianglex == 0) then
1655  write (errmsg, '(a)') 'Error. ANGLDEGX not provided in '// &
1656  'discretization file, but K22 was specified. '
1657  call store_error(errmsg)
1658  end if
1659  !
1660  ! -- Check to make sure values are greater than or equal to zero
1661  nerr = 0
1662  do n = 1, size(this%k22)
1663  if (this%ik22overk /= 0) this%k22(n) = this%k22(n) * this%k11(n)
1664  if (this%k22(n) <= dzero) then
1665  nerr = nerr + 1
1666  if (nerr <= 20) then
1667  call this%dis%noder_to_string(n, cellstr)
1668  write (errmsg, fmtkerr) trim(adjustl(aname(4))), trim(cellstr), &
1669  this%k22(n)
1670  call store_error(errmsg)
1671  end if
1672  end if
1673  end do
1674  if (nerr > 20) then
1675  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(4)))
1676  call store_error(errmsg)
1677  end if
1678  end if
1679  !
1680  ! -- check for wetdry conflicts
1681  if (this%irewet == 1) then
1682  if (this%iwetdry == 0) then
1683  write (errmsg, '(a, a, a)') 'Error in GRIDDATA block: ', &
1684  trim(adjustl(aname(5))), ' not found.'
1685  call store_error(errmsg)
1686  end if
1687  end if
1688  !
1689  ! -- Check for angle conflicts
1690  if (this%iangle1 /= 0) then
1691  do n = 1, size(this%angle1)
1692  this%angle1(n) = this%angle1(n) * dpio180
1693  end do
1694  else
1695  if (this%ixt3d /= 0) then
1696  this%iangle1 = 1
1697  write (this%iout, '(a)') 'XT3D IN USE, BUT ANGLE1 NOT SPECIFIED. '// &
1698  'SETTING ANGLE1 TO ZERO.'
1699  do n = 1, size(this%angle1)
1700  this%angle1(n) = dzero
1701  end do
1702  end if
1703  end if
1704  if (this%iangle2 /= 0) then
1705  if (this%iangle1 == 0) then
1706  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE1. '// &
1707  'ANGLE2 REQUIRES ANGLE1. '
1708  call store_error(errmsg)
1709  end if
1710  if (this%iangle3 == 0) then
1711  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE3. '// &
1712  'SPECIFY BOTH OR NEITHER ONE. '
1713  call store_error(errmsg)
1714  end if
1715  do n = 1, size(this%angle2)
1716  this%angle2(n) = this%angle2(n) * dpio180
1717  end do
1718  end if
1719  if (this%iangle3 /= 0) then
1720  if (this%iangle1 == 0) then
1721  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE1. '// &
1722  'ANGLE3 REQUIRES ANGLE1. '
1723  call store_error(errmsg)
1724  end if
1725  if (this%iangle2 == 0) then
1726  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE2. '// &
1727  'SPECIFY BOTH OR NEITHER ONE. '
1728  call store_error(errmsg)
1729  end if
1730  do n = 1, size(this%angle3)
1731  this%angle3(n) = this%angle3(n) * dpio180
1732  end do
1733  end if
1734  !
1735  ! -- terminate if data errors
1736  if (count_errors() > 0) then
1737  call store_error_filename(this%input_fname)
1738  end if
real(dp), parameter dpio180
real constant
Definition: Constants.f90:130
Here is the call graph for this function:

◆ preprocess_input()

subroutine gwfnpfmodule::preprocess_input ( class(gwfnpftype this)

This routine consists of the following steps:

  1. convert cells to noflow when all transmissive parameters equal zero
  2. perform initial wetting and drying
  3. initialize cell saturation
  4. calculate saturated conductance (when not xt3d)
  5. If NEWTON under-relaxation, determine lower most node
    Parameters
    thisthe instance of the NPF package

Definition at line 1751 of file gwf-npf.f90.

1752  ! -- modules
1753  use constantsmodule, only: linelength
1755  ! -- dummy
1756  class(GwfNpfType) :: this !< the instance of the NPF package
1757  ! -- local
1758  integer(I4B) :: n, m, ii, nn
1759  real(DP) :: hyn, hym
1760  real(DP) :: satn, topn, botn
1761  integer(I4B) :: nextn
1762  real(DP) :: minbot, botm
1763  logical :: finished
1764  character(len=LINELENGTH) :: cellstr, errmsg
1765  ! -- format
1766  character(len=*), parameter :: fmtcnv = &
1767  "(1X,'CELL ', A, &
1768  &' ELIMINATED BECAUSE ALL HYDRAULIC CONDUCTIVITIES TO NODE ARE 0.')"
1769  character(len=*), parameter :: fmtnct = &
1770  &"(1X,'Negative cell thickness at cell ', A)"
1771  character(len=*), parameter :: fmtihbe = &
1772  &"(1X,'Initial head, bottom elevation:',1P,2G13.5)"
1773  character(len=*), parameter :: fmttebe = &
1774  &"(1X,'Top elevation, bottom elevation:',1P,2G13.5)"
1775  !
1776  do n = 1, this%dis%nodes
1777  this%ithickstartflag(n) = 0
1778  end do
1779  !
1780  ! -- Insure that each cell has at least one non-zero transmissive parameter
1781  ! Note that a cell can be deactivated even if it has a valid connection
1782  ! to another model.
1783  nodeloop: do n = 1, this%dis%nodes
1784  !
1785  ! -- Skip if already inactive
1786  if (this%ibound(n) == 0) then
1787  if (this%irewet /= 0) then
1788  if (this%wetdry(n) == dzero) cycle nodeloop
1789  else
1790  cycle nodeloop
1791  end if
1792  end if
1793  !
1794  ! -- Cycle if k11 is not zero
1795  if (this%k11(n) /= dzero) cycle nodeloop
1796  !
1797  ! -- Cycle if at least one vertical connection has non-zero k33
1798  ! for n and m
1799  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
1800  m = this%dis%con%ja(ii)
1801  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0) then
1802  hyn = this%k11(n)
1803  if (this%ik33 /= 0) hyn = this%k33(n)
1804  if (hyn /= dzero) then
1805  hym = this%k11(m)
1806  if (this%ik33 /= 0) hym = this%k33(m)
1807  if (hym /= dzero) cycle
1808  end if
1809  end if
1810  end do
1811  !
1812  ! -- If this part of the loop is reached, then all connections have
1813  ! zero transmissivity, so convert to noflow.
1814  this%ibound(n) = 0
1815  this%hnew(n) = this%hnoflo
1816  if (this%irewet /= 0) this%wetdry(n) = dzero
1817  call this%dis%noder_to_string(n, cellstr)
1818  write (this%iout, fmtcnv) trim(adjustl(cellstr))
1819  !
1820  end do nodeloop
1821  !
1822  ! -- Preprocess cell status and heads based on initial conditions
1823  if (this%inewton == 0) then
1824  !
1825  ! -- For standard formulation (non-Newton) call wetdry routine
1826  call this%wd(0, this%hnew)
1827  else
1828  !
1829  ! -- Newton formulation, so adjust heads to be above bottom
1830  ! (Not used in present formulation because variable cv
1831  ! cannot be used with Newton)
1832  if (this%ivarcv == 1) then
1833  do n = 1, this%dis%nodes
1834  if (this%hnew(n) < this%dis%bot(n)) then
1835  this%hnew(n) = this%dis%bot(n) + dem6
1836  end if
1837  end do
1838  end if
1839  end if
1840  !
1841  ! -- If THCKSTRT is not active, then loop through icelltype and replace
1842  ! any negative values with 1.
1843  if (this%ithickstrt == 0) then
1844  do n = 1, this%dis%nodes
1845  if (this%icelltype(n) < 0) then
1846  this%icelltype(n) = 1
1847  end if
1848  end do
1849  end if
1850  !
1851  ! -- Initialize sat to zero for ibound=0 cells, unless the cell can
1852  ! rewet. Initialize sat to the saturated fraction based on strt
1853  ! if icelltype is negative and the THCKSTRT option is in effect.
1854  ! Initialize sat to 1.0 for all other cells in order to calculate
1855  ! condsat in next section.
1856  do n = 1, this%dis%nodes
1857  if (this%ibound(n) == 0) then
1858  this%sat(n) = done
1859  if (this%icelltype(n) < 0 .and. this%ithickstrt /= 0) then
1860  this%ithickstartflag(n) = 1
1861  this%icelltype(n) = 0
1862  end if
1863  else
1864  topn = this%dis%top(n)
1865  botn = this%dis%bot(n)
1866  if (this%icelltype(n) < 0 .and. this%ithickstrt /= 0) then
1867  call this%thksat(n, this%ic%strt(n), satn)
1868  if (botn > this%ic%strt(n)) then
1869  call this%dis%noder_to_string(n, cellstr)
1870  write (errmsg, fmtnct) trim(adjustl(cellstr))
1871  call store_error(errmsg)
1872  write (errmsg, fmtihbe) this%ic%strt(n), botn
1873  call store_error(errmsg)
1874  end if
1875  this%ithickstartflag(n) = 1
1876  this%icelltype(n) = 0
1877  else
1878  satn = done
1879  if (botn > topn) then
1880  call this%dis%noder_to_string(n, cellstr)
1881  write (errmsg, fmtnct) trim(adjustl(cellstr))
1882  call store_error(errmsg)
1883  write (errmsg, fmttebe) topn, botn
1884  call store_error(errmsg)
1885  end if
1886  end if
1887  this%sat(n) = satn
1888  end if
1889  end do
1890  if (count_errors() > 0) then
1891  call store_error_filename(this%input_fname)
1892  end if
1893  !
1894  ! -- Calculate condsat, but only if xt3d is not active. If xt3d is
1895  ! active, then condsat is allocated to size of zero.
1896  if (this%ixt3d == 0) then
1897  !
1898  ! -- Calculate the saturated conductance for all connections assuming
1899  ! that saturation is 1 (except for case where icelltype was entered
1900  ! as a negative value and THCKSTRT option in effect)
1901  do n = 1, this%dis%nodes
1902  call this%calc_condsat(n, .true.)
1903  end do
1904  !
1905  end if
1906  !
1907  ! -- Determine the lower most node
1908  if (this%igwfnewtonur /= 0) then
1909  call mem_reallocate(this%ibotnode, this%dis%nodes, 'IBOTNODE', &
1910  trim(this%memoryPath))
1911  do n = 1, this%dis%nodes
1912  !
1913  minbot = this%dis%bot(n)
1914  nn = n
1915  finished = .false.
1916  do while (.not. finished)
1917  nextn = 0
1918  !
1919  ! -- Go through the connecting cells
1920  do ii = this%dis%con%ia(nn) + 1, this%dis%con%ia(nn + 1) - 1
1921  !
1922  ! -- Set the m cell number
1923  m = this%dis%con%ja(ii)
1924  botm = this%dis%bot(m)
1925  !
1926  ! -- select vertical connections: ihc == 0
1927  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0) then
1928  if (m > nn .and. botm < minbot) then
1929  nextn = m
1930  minbot = botm
1931  end if
1932  end if
1933  end do
1934  if (nextn > 0) then
1935  nn = nextn
1936  else
1937  finished = .true.
1938  end if
1939  end do
1940  this%ibotnode(n) = nn
1941  end do
1942  end if
1943  !
1944  ! -- nullify unneeded gwf pointers
1945  this%igwfnewtonur => null()
Here is the call graph for this function:

◆ rewet_check()

subroutine gwfnpfmodule::rewet_check ( class(gwfnpftype this,
integer(i4b), intent(in)  kiter,
integer(i4b), intent(in)  node,
real(dp), intent(in)  hm,
integer(i4b), intent(in)  ibdm,
integer(i4b), intent(in)  ihc,
real(dp), dimension(:), intent(inout)  hnew,
integer(i4b), intent(out)  irewet 
)

This method can be called from any external object that has a head that can be used to rewet the GWF cell node. The ihc value is used to determine if it is a vertical or horizontal connection, which can operate differently depending on user settings.

Definition at line 2170 of file gwf-npf.f90.

2171  ! -- dummy
2172  class(GwfNpfType) :: this
2173  integer(I4B), intent(in) :: kiter
2174  integer(I4B), intent(in) :: node
2175  real(DP), intent(in) :: hm
2176  integer(I4B), intent(in) :: ibdm
2177  integer(I4B), intent(in) :: ihc
2178  real(DP), intent(inout), dimension(:) :: hnew
2179  integer(I4B), intent(out) :: irewet
2180  ! -- local
2181  integer(I4B) :: itflg
2182  real(DP) :: wd, awd, turnon, bbot
2183  !
2184  irewet = 0
2185  !
2186  ! -- Convert a dry cell to wet if it meets the criteria
2187  if (this%irewet > 0) then
2188  itflg = mod(kiter, this%iwetit)
2189  if (itflg == 0) then
2190  if (this%ibound(node) == 0 .and. this%wetdry(node) /= dzero) then
2191  !
2192  ! -- Calculate wetting elevation
2193  bbot = this%dis%bot(node)
2194  wd = this%wetdry(node)
2195  awd = wd
2196  if (wd < 0) awd = -wd
2197  turnon = bbot + awd
2198  !
2199  ! -- Check head in adjacent cells to see if wetting elevation has
2200  ! been reached
2201  if (ihc == c3d_vertical) then
2202  !
2203  ! -- check cell below
2204  if (ibdm > 0 .and. hm >= turnon) irewet = 1
2205  else
2206  if (wd > dzero) then
2207  !
2208  ! -- check horizontally adjacent cells
2209  if (ibdm > 0 .and. hm >= turnon) irewet = 1
2210  end if
2211  end if
2212  !
2213  if (irewet == 1) then
2214  ! -- rewet cell; use equation 3a if ihdwet=0; use equation 3b if
2215  ! ihdwet is not 0.
2216  if (this%ihdwet == 0) then
2217  hnew(node) = bbot + this%wetfct * (hm - bbot)
2218  else
2219  hnew(node) = bbot + this%wetfct * awd !(hm - bbot)
2220  end if
2221  this%ibound(node) = 30000
2222  end if
2223  end if
2224  end if
2225  end if

◆ sav_sat()

subroutine gwfnpfmodule::sav_sat ( class(gwfnpftype this,
integer(i4b), intent(in)  ibinun 
)
private

Definition at line 2641 of file gwf-npf.f90.

2642  ! -- dummy
2643  class(GwfNpfType) :: this
2644  integer(I4B), intent(in) :: ibinun
2645  ! -- local
2646  character(len=16) :: text
2647  character(len=16), dimension(1) :: auxtxt
2648  real(DP), dimension(1) :: a
2649  integer(I4B) :: n
2650  integer(I4B) :: naux
2651  !
2652  ! -- Write the header
2653  text = ' DATA-SAT'
2654  naux = 1
2655  auxtxt(:) = [' sat']
2656  call this%dis%record_srcdst_list_header(text, this%name_model, &
2657  this%packName, this%name_model, &
2658  this%packName, naux, auxtxt, ibinun, &
2659  this%dis%nodes, this%iout)
2660  !
2661  ! -- Write a zero for Q, and then write saturation as an aux variables
2662  do n = 1, this%dis%nodes
2663  a(1) = this%sat(n)
2664  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, a)
2665  end do

◆ sav_spdis()

subroutine gwfnpfmodule::sav_spdis ( class(gwfnpftype this,
integer(i4b), intent(in)  ibinun 
)

Definition at line 2613 of file gwf-npf.f90.

2614  ! -- dummy
2615  class(GwfNpfType) :: this
2616  integer(I4B), intent(in) :: ibinun
2617  ! -- local
2618  character(len=16) :: text
2619  character(len=16), dimension(3) :: auxtxt
2620  integer(I4B) :: n
2621  integer(I4B) :: naux
2622  !
2623  ! -- Write the header
2624  text = ' DATA-SPDIS'
2625  naux = 3
2626  auxtxt(:) = [' qx', ' qy', ' qz']
2627  call this%dis%record_srcdst_list_header(text, this%name_model, &
2628  this%packName, this%name_model, &
2629  this%packName, naux, auxtxt, ibinun, &
2630  this%dis%nodes, this%iout)
2631  !
2632  ! -- Write a zero for Q, and then write qx, qy, qz as aux variables
2633  do n = 1, this%dis%nodes
2634  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, &
2635  this%spdis(:, n))
2636  end do

◆ set_edge_properties()

subroutine gwfnpfmodule::set_edge_properties ( class(gwfnpftype this,
integer(i4b), intent(in)  nodedge,
integer(i4b), intent(in)  ihcedge,
real(dp), intent(in)  q,
real(dp), intent(in)  area,
real(dp), intent(in)  nx,
real(dp), intent(in)  ny,
real(dp), intent(in)  distance 
)
private

Definition at line 2710 of file gwf-npf.f90.

2712  ! -- dummy
2713  class(GwfNpfType) :: this
2714  integer(I4B), intent(in) :: nodedge
2715  integer(I4B), intent(in) :: ihcedge
2716  real(DP), intent(in) :: q
2717  real(DP), intent(in) :: area
2718  real(DP), intent(in) :: nx
2719  real(DP), intent(in) :: ny
2720  real(DP), intent(in) :: distance
2721  ! -- local
2722  integer(I4B) :: lastedge
2723  !
2724  this%lastedge = this%lastedge + 1
2725  lastedge = this%lastedge
2726  this%nodedge(lastedge) = nodedge
2727  this%ihcedge(lastedge) = ihcedge
2728  this%propsedge(1, lastedge) = q
2729  this%propsedge(2, lastedge) = area
2730  this%propsedge(3, lastedge) = nx
2731  this%propsedge(4, lastedge) = ny
2732  this%propsedge(5, lastedge) = distance
2733  !
2734  ! -- If this is the last edge, then the next call must be starting a new
2735  ! edge properties assignment loop, so need to reset lastedge to 0
2736  if (this%lastedge == this%nedges) this%lastedge = 0

◆ set_options()

subroutine gwfnpfmodule::set_options ( class(gwfnpftype this,
type(gwfnpfoptionstype), intent(in)  options 
)

Definition at line 1371 of file gwf-npf.f90.

1372  ! -- dummy
1373  class(GwfNpftype) :: this
1374  type(GwfNpfOptionsType), intent(in) :: options
1375  !
1376  this%icellavg = options%icellavg
1377  this%ithickstrt = options%ithickstrt
1378  this%iperched = options%iperched
1379  this%ivarcv = options%ivarcv
1380  this%idewatcv = options%idewatcv
1381  this%irewet = options%irewet
1382  this%wetfct = options%wetfct
1383  this%iwetit = options%iwetit
1384  this%ihdwet = options%ihdwet

◆ sgwf_npf_qcalc()

subroutine gwfnpfmodule::sgwf_npf_qcalc ( class(gwfnpftype this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
integer(i4b), intent(in)  icon,
real(dp), intent(inout)  qnm 
)
private

Definition at line 798 of file gwf-npf.f90.

799  ! -- dummy
800  class(GwfNpfType) :: this
801  integer(I4B), intent(in) :: n
802  integer(I4B), intent(in) :: m
803  real(DP), intent(in) :: hn
804  real(DP), intent(in) :: hm
805  integer(I4B), intent(in) :: icon
806  real(DP), intent(inout) :: qnm
807  ! -- local
808  real(DP) :: hyn, hym
809  real(DP) :: condnm
810  real(DP) :: hntemp, hmtemp
811  integer(I4B) :: ihc
812  !
813  ! -- Initialize
814  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
815  hyn = this%hy_eff(n, m, ihc, ipos=icon)
816  hym = this%hy_eff(m, n, ihc, ipos=icon)
817  !
818  ! -- Calculate conductance
819  if (ihc == c3d_vertical) then
820  condnm = vcond(this%ibound(n), this%ibound(m), &
821  this%icelltype(n), this%icelltype(m), this%inewton, &
822  this%ivarcv, this%idewatcv, &
823  this%condsat(this%dis%con%jas(icon)), hn, hm, &
824  hyn, hym, &
825  this%sat(n), this%sat(m), &
826  this%dis%top(n), this%dis%top(m), &
827  this%dis%bot(n), this%dis%bot(m), &
828  this%dis%con%hwva(this%dis%con%jas(icon)))
829  else
830  condnm = hcond(this%ibound(n), this%ibound(m), &
831  this%icelltype(n), this%icelltype(m), &
832  this%inewton, &
833  this%dis%con%ihc(this%dis%con%jas(icon)), &
834  this%icellavg, &
835  this%condsat(this%dis%con%jas(icon)), &
836  hn, hm, this%sat(n), this%sat(m), hyn, hym, &
837  this%dis%top(n), this%dis%top(m), &
838  this%dis%bot(n), this%dis%bot(m), &
839  this%dis%con%cl1(this%dis%con%jas(icon)), &
840  this%dis%con%cl2(this%dis%con%jas(icon)), &
841  this%dis%con%hwva(this%dis%con%jas(icon)))
842  end if
843  !
844  ! -- Initialize hntemp and hmtemp
845  hntemp = hn
846  hmtemp = hm
847  !
848  ! -- Check and adjust for dewatered conditions
849  if (this%iperched /= 0) then
850  if (this%dis%con%ihc(this%dis%con%jas(icon)) == 0) then
851  if (n > m) then
852  if (this%icelltype(n) /= 0) then
853  if (hn < this%dis%top(n)) hntemp = this%dis%bot(m)
854  end if
855  else
856  if (this%icelltype(m) /= 0) then
857  if (hm < this%dis%top(m)) hmtemp = this%dis%bot(n)
858  end if
859  end if
860  end if
861  end if
862  !
863  ! -- Calculate flow positive into cell n
864  qnm = condnm * (hmtemp - hntemp)
Here is the call graph for this function:

◆ sgwf_npf_thksat()

subroutine gwfnpfmodule::sgwf_npf_thksat ( class(gwfnpftype this,
integer(i4b), intent(in)  n,
real(dp), intent(in)  hn,
real(dp), intent(inout)  thksat 
)
private

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

776  ! -- dummy
777  class(GwfNpfType) :: this
778  integer(I4B), intent(in) :: n
779  real(DP), intent(in) :: hn
780  real(DP), intent(inout) :: thksat
781  !
782  ! -- Standard Formulation
783  if (hn >= this%dis%top(n)) then
784  thksat = done
785  else
786  thksat = (hn - this%dis%bot(n)) / (this%dis%top(n) - this%dis%bot(n))
787  end if
788  !
789  ! -- Newton-Raphson Formulation
790  if (this%inewton /= 0) then
791  thksat = squadraticsaturation(this%dis%top(n), this%dis%bot(n), hn, &
792  this%satomega)
793  end if
Here is the call graph for this function:

◆ sgwf_npf_wdmsg()

subroutine gwfnpfmodule::sgwf_npf_wdmsg ( class(gwfnpftype this,
integer(i4b), intent(in)  icode,
integer(i4b), intent(inout)  ncnvrt,
character(len=30), dimension(5), intent(inout)  nodcnvrt,
character(len=3), dimension(5), intent(inout)  acnvrt,
integer(i4b), intent(inout)  ihdcnv,
integer(i4b), intent(in)  kiter,
integer(i4b), intent(in)  n 
)
private

Definition at line 2230 of file gwf-npf.f90.

2232  ! -- modules
2233  use tdismodule, only: kstp, kper
2234  ! -- dummy
2235  class(GwfNpfType) :: this
2236  integer(I4B), intent(in) :: icode
2237  integer(I4B), intent(inout) :: ncnvrt
2238  character(len=30), dimension(5), intent(inout) :: nodcnvrt
2239  character(len=3), dimension(5), intent(inout) :: acnvrt
2240  integer(I4B), intent(inout) :: ihdcnv
2241  integer(I4B), intent(in) :: kiter
2242  integer(I4B), intent(in) :: n
2243  ! -- local
2244  integer(I4B) :: l
2245  ! -- formats
2246  character(len=*), parameter :: fmtcnvtn = &
2247  "(1X,/1X,'CELL CONVERSIONS FOR ITER.=',I0, &
2248  &' STEP=',I0,' PERIOD=',I0,' (NODE or LRC)')"
2249  character(len=*), parameter :: fmtnode = "(1X,3X,5(A4, A20))"
2250  !
2251  ! -- Keep track of cell conversions
2252  if (icode > 0) then
2253  ncnvrt = ncnvrt + 1
2254  call this%dis%noder_to_string(n, nodcnvrt(ncnvrt))
2255  if (icode == 1) then
2256  acnvrt(ncnvrt) = 'DRY'
2257  else
2258  acnvrt(ncnvrt) = 'WET'
2259  end if
2260  end if
2261  !
2262  ! -- Print a line if 5 conversions have occurred or if icode indicates that a
2263  ! partial line should be printed
2264  if (ncnvrt == 5 .or. (icode == 0 .and. ncnvrt > 0)) then
2265  if (ihdcnv == 0) write (this%iout, fmtcnvtn) kiter, kstp, kper
2266  ihdcnv = 1
2267  write (this%iout, fmtnode) &
2268  (acnvrt(l), trim(adjustl(nodcnvrt(l))), l=1, ncnvrt)
2269  ncnvrt = 0
2270  end if

◆ sgwf_npf_wetdry()

subroutine gwfnpfmodule::sgwf_npf_wetdry ( class(gwfnpftype this,
integer(i4b), intent(in)  kiter,
real(dp), dimension(:), intent(inout)  hnew 
)
private

Definition at line 2064 of file gwf-npf.f90.

2065  ! -- modules
2066  use tdismodule, only: kstp, kper
2068  use constantsmodule, only: linelength
2069  ! -- dummy
2070  class(GwfNpfType) :: this
2071  integer(I4B), intent(in) :: kiter
2072  real(DP), intent(inout), dimension(:) :: hnew
2073  ! -- local
2074  integer(I4B) :: n, m, ii, ihc
2075  real(DP) :: ttop, bbot, thick
2076  integer(I4B) :: ncnvrt, ihdcnv
2077  character(len=30), dimension(5) :: nodcnvrt
2078  character(len=30) :: nodestr
2079  character(len=3), dimension(5) :: acnvrt
2080  character(len=LINELENGTH) :: errmsg
2081  integer(I4B) :: irewet
2082  ! -- formats
2083  character(len=*), parameter :: fmtnct = &
2084  "(1X,/1X,'Negative cell thickness at (layer,row,col)', &
2085  &I4,',',I5,',',I5)"
2086  character(len=*), parameter :: fmttopbot = &
2087  &"(1X,'Top elevation, bottom elevation:',1P,2G13.5)"
2088  character(len=*), parameter :: fmttopbotthk = &
2089  &"(1X,'Top elevation, bottom elevation, thickness:',1P,3G13.5)"
2090  character(len=*), parameter :: fmtdrychd = &
2091  &"(1X,/1X,'CONSTANT-HEAD CELL WENT DRY -- SIMULATION ABORTED')"
2092  character(len=*), parameter :: fmtni = &
2093  &"(1X,'CELLID=',a,' ITERATION=',I0,' TIME STEP=',I0,' STRESS PERIOD=',I0)"
2094  !
2095  ! -- Initialize
2096  ncnvrt = 0
2097  ihdcnv = 0
2098  !
2099  ! -- Convert dry cells to wet
2100  do n = 1, this%dis%nodes
2101  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
2102  m = this%dis%con%ja(ii)
2103  ihc = this%dis%con%ihc(this%dis%con%jas(ii))
2104  call this%rewet_check(kiter, n, hnew(m), this%ibound(m), ihc, hnew, &
2105  irewet)
2106  if (irewet == 1) then
2107  call this%wdmsg(2, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2108  end if
2109  end do
2110  end do
2111  !
2112  ! -- Perform drying
2113  do n = 1, this%dis%nodes
2114  !
2115  ! -- cycle if inactive or confined
2116  if (this%ibound(n) == 0) cycle
2117  if (this%icelltype(n) == 0) cycle
2118  !
2119  ! -- check for negative cell thickness
2120  bbot = this%dis%bot(n)
2121  ttop = this%dis%top(n)
2122  if (bbot > ttop) then
2123  write (errmsg, fmtnct) n
2124  call store_error(errmsg)
2125  write (errmsg, fmttopbot) ttop, bbot
2126  call store_error(errmsg)
2127  call store_error_filename(this%input_fname)
2128  end if
2129  !
2130  ! -- Calculate saturated thickness
2131  if (this%icelltype(n) /= 0) then
2132  if (hnew(n) < ttop) ttop = hnew(n)
2133  end if
2134  thick = ttop - bbot
2135  !
2136  ! -- If thick<0 print message, set hnew, and ibound
2137  if (thick <= dzero) then
2138  call this%wdmsg(1, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2139  hnew(n) = this%hdry
2140  if (this%ibound(n) < 0) then
2141  errmsg = 'CONSTANT-HEAD CELL WENT DRY -- SIMULATION ABORTED'
2142  call store_error(errmsg)
2143  write (errmsg, fmttopbotthk) ttop, bbot, thick
2144  call store_error(errmsg)
2145  call this%dis%noder_to_string(n, nodestr)
2146  write (errmsg, fmtni) trim(adjustl(nodestr)), kiter, kstp, kper
2147  call store_error(errmsg)
2148  call store_error_filename(this%input_fname)
2149  end if
2150  this%ibound(n) = 0
2151  end if
2152  end do
2153  !
2154  ! -- Print remaining cell conversions
2155  call this%wdmsg(0, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2156  !
2157  ! -- Change ibound from 30000 to 1
2158  do n = 1, this%dis%nodes
2159  if (this%ibound(n) == 30000) this%ibound(n) = 1
2160  end do
Here is the call graph for this function:

◆ source_griddata()

subroutine gwfnpfmodule::source_griddata ( class(gwfnpftype this)

Definition at line 1501 of file gwf-npf.f90.

1502  ! -- modules
1503  use simmodule, only: count_errors, store_error
1507  ! -- dummy
1508  class(GwfNpftype) :: this
1509  ! -- locals
1510  character(len=LINELENGTH) :: errmsg
1511  type(GwfNpfParamFoundType) :: found
1512  logical, dimension(2) :: afound
1513  integer(I4B), dimension(:), pointer, contiguous :: map
1514  !
1515  ! -- set map to convert user input data into reduced data
1516  map => null()
1517  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1518  !
1519  ! -- update defaults with idm sourced values
1520  call mem_set_value(this%icelltype, 'ICELLTYPE', this%input_mempath, map, &
1521  found%icelltype)
1522  call mem_set_value(this%k11, 'K', this%input_mempath, map, found%k)
1523  call mem_set_value(this%k33, 'K33', this%input_mempath, map, found%k33)
1524  call mem_set_value(this%k22, 'K22', this%input_mempath, map, found%k22)
1525  call mem_set_value(this%wetdry, 'WETDRY', this%input_mempath, map, &
1526  found%wetdry)
1527  call mem_set_value(this%angle1, 'ANGLE1', this%input_mempath, map, &
1528  found%angle1)
1529  call mem_set_value(this%angle2, 'ANGLE2', this%input_mempath, map, &
1530  found%angle2)
1531  call mem_set_value(this%angle3, 'ANGLE3', this%input_mempath, map, &
1532  found%angle3)
1533  !
1534  ! -- ensure ICELLTYPE was found
1535  if (.not. found%icelltype) then
1536  write (errmsg, '(a)') 'Error in GRIDDATA block: ICELLTYPE not found.'
1537  call store_error(errmsg)
1538  end if
1539  !
1540  ! -- ensure K was found
1541  if (.not. found%k) then
1542  write (errmsg, '(a)') 'Error in GRIDDATA block: K not found.'
1543  call store_error(errmsg)
1544  end if
1545  !
1546  ! -- set error if ik33overk set with no k33
1547  if (.not. found%k33 .and. this%ik33overk /= 0) then
1548  write (errmsg, '(a)') 'K33OVERK option specified but K33 not specified.'
1549  call store_error(errmsg)
1550  end if
1551  !
1552  ! -- set error if ik22overk set with no k22
1553  if (.not. found%k22 .and. this%ik22overk /= 0) then
1554  write (errmsg, '(a)') 'K22OVERK option specified but K22 not specified.'
1555  call store_error(errmsg)
1556  end if
1557  !
1558  ! -- handle found side effects
1559  if (found%k33) this%ik33 = 1
1560  if (found%k22) this%ik22 = 1
1561  if (found%wetdry) this%iwetdry = 1
1562  if (found%angle1) this%iangle1 = 1
1563  if (found%angle2) this%iangle2 = 1
1564  if (found%angle3) this%iangle3 = 1
1565  !
1566  ! -- handle not found side effects
1567  if (.not. found%k33) then
1568  call mem_set_value(this%k33, 'K', this%input_mempath, map, afound(1))
1569  end if
1570  if (.not. found%k22) then
1571  call mem_set_value(this%k22, 'K', this%input_mempath, map, afound(2))
1572  end if
1573  if (.not. found%wetdry) call mem_reallocate(this%wetdry, 1, 'WETDRY', &
1574  trim(this%memoryPath))
1575  if (.not. found%angle1 .and. this%ixt3d == 0) &
1576  call mem_reallocate(this%angle1, 0, 'ANGLE1', trim(this%memoryPath))
1577  if (.not. found%angle2 .and. this%ixt3d == 0) &
1578  call mem_reallocate(this%angle2, 0, 'ANGLE2', trim(this%memoryPath))
1579  if (.not. found%angle3 .and. this%ixt3d == 0) &
1580  call mem_reallocate(this%angle3, 0, 'ANGLE3', trim(this%memoryPath))
1581  !
1582  ! -- log griddata
1583  if (this%iout > 0) then
1584  call this%log_griddata(found)
1585  end if
Here is the call graph for this function:

◆ source_options()

subroutine gwfnpfmodule::source_options ( class(gwfnpftype this)

Definition at line 1295 of file gwf-npf.f90.

1296  ! -- modules
1302  use sourcecommonmodule, only: filein_fname
1303  ! -- dummy
1304  class(GwfNpftype) :: this
1305  ! -- locals
1306  character(len=LENVARNAME), dimension(3) :: cellavg_method = &
1307  &[character(len=LENVARNAME) :: 'LOGARITHMIC', 'AMT-LMK', 'AMT-HMK']
1308  type(gwfnpfparamfoundtype) :: found
1309  character(len=LINELENGTH) :: tvk6_filename
1310  !
1311  ! -- update defaults with idm sourced values
1312  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
1313  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
1314  call mem_set_value(this%icellavg, 'CELLAVG', this%input_mempath, &
1315  cellavg_method, found%cellavg)
1316  call mem_set_value(this%ithickstrt, 'ITHICKSTRT', this%input_mempath, &
1317  found%ithickstrt)
1318  call mem_set_value(this%iperched, 'IPERCHED', this%input_mempath, &
1319  found%iperched)
1320  call mem_set_value(this%ivarcv, 'IVARCV', this%input_mempath, found%ivarcv)
1321  call mem_set_value(this%idewatcv, 'IDEWATCV', this%input_mempath, &
1322  found%idewatcv)
1323  call mem_set_value(this%ixt3d, 'IXT3D', this%input_mempath, found%ixt3d)
1324  call mem_set_value(this%ixt3drhs, 'IXT3DRHS', this%input_mempath, &
1325  found%ixt3drhs)
1326  call mem_set_value(this%isavspdis, 'ISAVSPDIS', this%input_mempath, &
1327  found%isavspdis)
1328  call mem_set_value(this%isavsat, 'ISAVSAT', this%input_mempath, found%isavsat)
1329  call mem_set_value(this%ik22overk, 'IK22OVERK', this%input_mempath, &
1330  found%ik22overk)
1331  call mem_set_value(this%ik33overk, 'IK33OVERK', this%input_mempath, &
1332  found%ik33overk)
1333  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
1334  call mem_set_value(this%satomega, 'SATOMEGA', this%input_mempath, &
1335  found%satomega)
1336  call mem_set_value(this%irewet, 'IREWET', this%input_mempath, found%irewet)
1337  call mem_set_value(this%wetfct, 'WETFCT', this%input_mempath, found%wetfct)
1338  call mem_set_value(this%iwetit, 'IWETIT', this%input_mempath, found%iwetit)
1339  call mem_set_value(this%ihdwet, 'IHDWET', this%input_mempath, found%ihdwet)
1340  !
1341  ! -- save flows option active
1342  if (found%ipakcb) this%ipakcb = -1
1343  !
1344  ! -- xt3d active with rhs
1345  if (found%ixt3d .and. found%ixt3drhs) this%ixt3d = 2
1346  !
1347  ! -- save specific discharge active
1348  if (found%isavspdis) this%icalcspdis = this%isavspdis
1349  !
1350  ! -- no newton specified
1351  if (found%inewton) then
1352  this%inewton = 0
1353  this%iasym = 0
1354  end if
1355  !
1356  ! -- enforce 0 or 1 TVK6_FILENAME entries in option block
1357  if (filein_fname(tvk6_filename, 'TVK6_FILENAME', this%input_mempath, &
1358  this%input_fname)) then
1359  call openfile(this%intvk, this%iout, tvk6_filename, 'TVK')
1360  call tvk_cr(this%tvk, this%name_model, this%intvk, this%iout)
1361  end if
1362  !
1363  ! -- log options
1364  if (this%iout > 0) then
1365  call this%log_options(found)
1366  end if
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
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:

◆ store_original_k_arrays()

subroutine gwfnpfmodule::store_original_k_arrays ( class(gwfnpftype this,
integer(i4b), intent(in)  ncells,
integer(i4b), intent(in)  njas 
)

The K arrays (K11, etc.) get multiplied by the viscosity ratio so that subsequent uses of K already take into account the effect of viscosity. Thus the original user-specified K array values are lost unless they are backed up in k11input, for example. In a new stress period/time step, the values in k11input are multiplied by the viscosity ratio, not k11 since it contains viscosity-adjusted hydraulic conductivity values.

Definition at line 1143 of file gwf-npf.f90.

1144  ! -- modules
1146  ! -- dummy
1147  class(GwfNpftype) :: this
1148  integer(I4B), intent(in) :: ncells
1149  integer(I4B), intent(in) :: njas
1150  ! -- local
1151  integer(I4B) :: n
1152  !
1153  ! -- Retain copy of user-specified K arrays
1154  do n = 1, ncells
1155  this%k11input(n) = this%k11(n)
1156  this%k22input(n) = this%k22(n)
1157  this%k33input(n) = this%k33(n)
1158  end do