MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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...
 
subroutine set_edge_properties (this, nodedge, ihcedge, q, area, nx, ny, distance)
 Provide the npf package with edge properties. More...
 
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 1192 of file gwf-npf.f90.

1193  ! -- dummy
1194  class(GwfNpftype) :: this
1195  integer(I4B), intent(in) :: ncells
1196  integer(I4B), intent(in) :: njas
1197  ! -- local
1198  integer(I4B) :: n
1199  !
1200  call mem_allocate(this%ithickstartflag, ncells, 'ITHICKSTARTFLAG', &
1201  this%memoryPath)
1202  call mem_allocate(this%icelltype, ncells, 'ICELLTYPE', this%memoryPath)
1203  call mem_allocate(this%k11, ncells, 'K11', this%memoryPath)
1204  call mem_allocate(this%sat, ncells, 'SAT', this%memoryPath)
1205  call mem_allocate(this%condsat, njas, 'CONDSAT', this%memoryPath)
1206  !
1207  ! -- Optional arrays dimensioned to full size initially
1208  call mem_allocate(this%k22, ncells, 'K22', this%memoryPath)
1209  call mem_allocate(this%k33, ncells, 'K33', this%memoryPath)
1210  call mem_allocate(this%wetdry, ncells, 'WETDRY', this%memoryPath)
1211  call mem_allocate(this%angle1, ncells, 'ANGLE1', this%memoryPath)
1212  call mem_allocate(this%angle2, ncells, 'ANGLE2', this%memoryPath)
1213  call mem_allocate(this%angle3, ncells, 'ANGLE3', this%memoryPath)
1214  !
1215  ! -- Optional arrays
1216  call mem_allocate(this%ibotnode, 0, 'IBOTNODE', this%memoryPath)
1217  call mem_allocate(this%nodedge, 0, 'NODEDGE', this%memoryPath)
1218  call mem_allocate(this%ihcedge, 0, 'IHCEDGE', this%memoryPath)
1219  call mem_allocate(this%propsedge, 0, 0, 'PROPSEDGE', this%memoryPath)
1220  !
1221  ! -- Optional arrays only needed when vsc package is active
1222  call mem_allocate(this%k11input, 0, 'K11INPUT', this%memoryPath)
1223  call mem_allocate(this%k22input, 0, 'K22INPUT', this%memoryPath)
1224  call mem_allocate(this%k33input, 0, 'K33INPUT', this%memoryPath)
1225  !
1226  ! -- Specific discharge is (re-)allocated when nedges is known
1227  call mem_allocate(this%spdis, 3, 0, 'SPDIS', this%memoryPath)
1228  !
1229  ! -- Time-varying property flag arrays
1230  call mem_allocate(this%nodekchange, ncells, 'NODEKCHANGE', this%memoryPath)
1231  !
1232  ! -- initialize iangle1, iangle2, iangle3, and wetdry
1233  do n = 1, ncells
1234  this%angle1(n) = dzero
1235  this%angle2(n) = dzero
1236  this%angle3(n) = dzero
1237  this%wetdry(n) = dzero
1238  this%nodekchange(n) = dzero
1239  end do
1240  !
1241  ! -- allocate variable names
1242  allocate (this%aname(this%iname))
1243  this%aname = [' ICELLTYPE', ' K', &
1244  ' K33', ' K22', &
1245  ' WETDRY', ' ANGLE1', &
1246  ' ANGLE2', ' ANGLE3']
1247  !
1248  ! -- Return
1249  return

◆ 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 1069 of file gwf-npf.f90.

1070  ! -- modules
1072  ! -- dummy
1073  class(GwfNpftype) :: this
1074  !
1075  ! -- allocate scalars in NumericalPackageType
1076  call this%NumericalPackageType%allocate_scalars()
1077  !
1078  ! -- Allocate scalars
1079  call mem_allocate(this%iname, 'INAME', this%memoryPath)
1080  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
1081  call mem_allocate(this%ixt3drhs, 'IXT3DRHS', this%memoryPath)
1082  call mem_allocate(this%satomega, 'SATOMEGA', this%memoryPath)
1083  call mem_allocate(this%hnoflo, 'HNOFLO', this%memoryPath)
1084  call mem_allocate(this%hdry, 'HDRY', this%memoryPath)
1085  call mem_allocate(this%icellavg, 'ICELLAVG', this%memoryPath)
1086  call mem_allocate(this%iavgkeff, 'IAVGKEFF', this%memoryPath)
1087  call mem_allocate(this%ik22, 'IK22', this%memoryPath)
1088  call mem_allocate(this%ik33, 'IK33', this%memoryPath)
1089  call mem_allocate(this%ik22overk, 'IK22OVERK', this%memoryPath)
1090  call mem_allocate(this%ik33overk, 'IK33OVERK', this%memoryPath)
1091  call mem_allocate(this%iperched, 'IPERCHED', this%memoryPath)
1092  call mem_allocate(this%ivarcv, 'IVARCV', this%memoryPath)
1093  call mem_allocate(this%idewatcv, 'IDEWATCV', this%memoryPath)
1094  call mem_allocate(this%ithickstrt, 'ITHICKSTRT', this%memoryPath)
1095  call mem_allocate(this%icalcspdis, 'ICALCSPDIS', this%memoryPath)
1096  call mem_allocate(this%isavspdis, 'ISAVSPDIS', this%memoryPath)
1097  call mem_allocate(this%isavsat, 'ISAVSAT', this%memoryPath)
1098  call mem_allocate(this%irewet, 'IREWET', this%memoryPath)
1099  call mem_allocate(this%wetfct, 'WETFCT', this%memoryPath)
1100  call mem_allocate(this%iwetit, 'IWETIT', this%memoryPath)
1101  call mem_allocate(this%ihdwet, 'IHDWET', this%memoryPath)
1102  call mem_allocate(this%iangle1, 'IANGLE1', this%memoryPath)
1103  call mem_allocate(this%iangle2, 'IANGLE2', this%memoryPath)
1104  call mem_allocate(this%iangle3, 'IANGLE3', this%memoryPath)
1105  call mem_allocate(this%iwetdry, 'IWETDRY', this%memoryPath)
1106  call mem_allocate(this%nedges, 'NEDGES', this%memoryPath)
1107  call mem_allocate(this%lastedge, 'LASTEDGE', this%memoryPath)
1108  call mem_allocate(this%intvk, 'INTVK', this%memoryPath)
1109  call mem_allocate(this%invsc, 'INVSC', this%memoryPath)
1110  call mem_allocate(this%kchangeper, 'KCHANGEPER', this%memoryPath)
1111  call mem_allocate(this%kchangestp, 'KCHANGESTP', this%memoryPath)
1112  !
1113  ! -- set pointer to inewtonur
1114  call mem_setptr(this%igwfnewtonur, 'INEWTONUR', &
1115  create_mem_path(this%name_model))
1116  !
1117  ! -- Initialize value
1118  this%iname = 8
1119  this%ixt3d = 0
1120  this%ixt3drhs = 0
1121  this%satomega = dzero
1122  this%hnoflo = dhnoflo !1.d30
1123  this%hdry = dhdry !-1.d30
1124  this%icellavg = ccond_hmean
1125  this%iavgkeff = 0
1126  this%ik22 = 0
1127  this%ik33 = 0
1128  this%ik22overk = 0
1129  this%ik33overk = 0
1130  this%iperched = 0
1131  this%ivarcv = 0
1132  this%idewatcv = 0
1133  this%ithickstrt = 0
1134  this%icalcspdis = 0
1135  this%isavspdis = 0
1136  this%isavsat = 0
1137  this%irewet = 0
1138  this%wetfct = done
1139  this%iwetit = 1
1140  this%ihdwet = 0
1141  this%iangle1 = 0
1142  this%iangle2 = 0
1143  this%iangle3 = 0
1144  this%iwetdry = 0
1145  this%nedges = 0
1146  this%lastedge = 0
1147  this%intvk = 0
1148  this%invsc = 0
1149  this%kchangeper = 0
1150  this%kchangestp = 0
1151  !
1152  ! -- If newton is on, then NPF creates asymmetric matrix
1153  this%iasym = this%inewton
1154  !
1155  ! -- Return
1156  return
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 2005 of file gwf-npf.f90.

2006  ! -- dummy variables
2007  class(GwfNpfType) :: this
2008  integer(I4B), intent(in) :: node
2009  logical, intent(in) :: upperOnly
2010  ! -- local variables
2011  integer(I4B) :: ii, m, n, ihc, jj
2012  real(DP) :: topm, topn, topnode, botm, botn, botnode, satm, satn, satnode
2013  real(DP) :: hyn, hym, hn, hm, fawidth, csat
2014  !
2015  satnode = this%calc_initial_sat(node)
2016  !
2017  topnode = this%dis%top(node)
2018  botnode = this%dis%bot(node)
2019  !
2020  ! -- Go through the connecting cells
2021  do ii = this%dis%con%ia(node) + 1, this%dis%con%ia(node + 1) - 1
2022  !
2023  ! -- Set the m cell number and cycle if lower triangle connection and
2024  ! -- we're not updating both upper and lower matrix parts for this node
2025  m = this%dis%con%ja(ii)
2026  jj = this%dis%con%jas(ii)
2027  if (m < node) then
2028  if (upperonly) cycle
2029  ! m => node, n => neighbour
2030  n = m
2031  m = node
2032  topm = topnode
2033  botm = botnode
2034  satm = satnode
2035  topn = this%dis%top(n)
2036  botn = this%dis%bot(n)
2037  satn = this%calc_initial_sat(n)
2038  else
2039  ! n => node, m => neighbour
2040  n = node
2041  topn = topnode
2042  botn = botnode
2043  satn = satnode
2044  topm = this%dis%top(m)
2045  botm = this%dis%bot(m)
2046  satm = this%calc_initial_sat(m)
2047  end if
2048  !
2049  ihc = this%dis%con%ihc(jj)
2050  hyn = this%hy_eff(n, m, ihc, ipos=ii)
2051  hym = this%hy_eff(m, n, ihc, ipos=ii)
2052  if (this%ithickstartflag(n) == 0) then
2053  hn = topn
2054  else
2055  hn = this%ic%strt(n)
2056  end if
2057  if (this%ithickstartflag(m) == 0) then
2058  hm = topm
2059  else
2060  hm = this%ic%strt(m)
2061  end if
2062  !
2063  ! -- Calculate conductance depending on whether connection is
2064  ! vertical (0), horizontal (1), or staggered horizontal (2)
2065  if (ihc == c3d_vertical) then
2066  !
2067  ! -- Vertical conductance for fully saturated conditions
2068  csat = vcond(1, 1, 1, 1, 0, 1, 1, done, &
2069  botn, botm, &
2070  hyn, hym, &
2071  satn, satm, &
2072  topn, topm, &
2073  botn, botm, &
2074  this%dis%con%hwva(jj))
2075  else
2076  !
2077  ! -- Horizontal conductance for fully saturated conditions
2078  fawidth = this%dis%con%hwva(jj)
2079  csat = hcond(1, 1, 1, 1, 0, &
2080  ihc, &
2081  this%icellavg, &
2082  done, &
2083  hn, hm, satn, satm, hyn, hym, &
2084  topn, topm, &
2085  botn, botm, &
2086  this%dis%con%cl1(jj), &
2087  this%dis%con%cl2(jj), &
2088  fawidth)
2089  end if
2090  this%condsat(jj) = csat
2091  end do
2092  !
2093  ! -- Return
2094  return
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 2104 of file gwf-npf.f90.

2105  ! -- dummy variables
2106  class(GwfNpfType) :: this
2107  integer(I4B), intent(in) :: n
2108  ! -- Return
2109  real(DP) :: satn
2110  !
2111  satn = done
2112  if (this%ibound(n) /= 0 .and. this%ithickstartflag(n) /= 0) then
2113  call this%thksat(n, this%ic%strt(n), satn)
2114  end if
2115  !
2116  return

◆ calc_spdis()

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

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

2432  ! -- modules
2433  use simmodule, only: store_error
2434  ! -- dummy
2435  class(GwfNpfType) :: this
2436  real(DP), intent(in), dimension(:) :: flowja
2437  ! -- local
2438  integer(I4B) :: n
2439  integer(I4B) :: m
2440  integer(I4B) :: ipos
2441  integer(I4B) :: isympos
2442  integer(I4B) :: ihc
2443  integer(I4B) :: ic
2444  integer(I4B) :: iz
2445  integer(I4B) :: nc
2446  integer(I4B) :: ncz
2447  real(DP) :: qz
2448  real(DP) :: vx
2449  real(DP) :: vy
2450  real(DP) :: vz
2451  real(DP) :: xn
2452  real(DP) :: yn
2453  real(DP) :: zn
2454  real(DP) :: xc
2455  real(DP) :: yc
2456  real(DP) :: zc
2457  real(DP) :: cl1
2458  real(DP) :: cl2
2459  real(DP) :: dltot
2460  real(DP) :: ooclsum
2461  real(DP) :: dsumx
2462  real(DP) :: dsumy
2463  real(DP) :: dsumz
2464  real(DP) :: denom
2465  real(DP) :: area
2466  real(DP) :: dz
2467  real(DP) :: axy
2468  real(DP) :: ayx
2469  real(DP), allocatable, dimension(:) :: vi
2470  real(DP), allocatable, dimension(:) :: di
2471  real(DP), allocatable, dimension(:) :: viz
2472  real(DP), allocatable, dimension(:) :: diz
2473  real(DP), allocatable, dimension(:) :: nix
2474  real(DP), allocatable, dimension(:) :: niy
2475  real(DP), allocatable, dimension(:) :: wix
2476  real(DP), allocatable, dimension(:) :: wiy
2477  real(DP), allocatable, dimension(:) :: wiz
2478  real(DP), allocatable, dimension(:) :: bix
2479  real(DP), allocatable, dimension(:) :: biy
2480  logical :: nozee = .true.
2481  !
2482  ! -- Ensure dis has necessary information
2483  if (this%icalcspdis /= 0 .and. this%dis%con%ianglex == 0) then
2484  call store_error('Error. ANGLDEGX not provided in '// &
2485  'discretization file. ANGLDEGX required for '// &
2486  'calculation of specific discharge.', terminate=.true.)
2487  end if
2488  !
2489  ! -- Find max number of connections and allocate weight arrays
2490  nc = 0
2491  do n = 1, this%dis%nodes
2492  !
2493  ! -- Count internal model connections
2494  ic = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
2495  !
2496  ! -- Count edge connections
2497  do m = 1, this%nedges
2498  if (this%nodedge(m) == n) then
2499  ic = ic + 1
2500  end if
2501  end do
2502  !
2503  ! -- Set max number of connections for any cell
2504  if (ic > nc) nc = ic
2505  end do
2506  !
2507  ! -- Allocate storage arrays needed for cell-centered spdis calculation
2508  allocate (vi(nc))
2509  allocate (di(nc))
2510  allocate (viz(nc))
2511  allocate (diz(nc))
2512  allocate (nix(nc))
2513  allocate (niy(nc))
2514  allocate (wix(nc))
2515  allocate (wiy(nc))
2516  allocate (wiz(nc))
2517  allocate (bix(nc))
2518  allocate (biy(nc))
2519  !
2520  ! -- Go through each cell and calculate specific discharge
2521  do n = 1, this%dis%nodes
2522  !
2523  ! -- first calculate geometric properties for x and y directions and
2524  ! the specific discharge at a face (vi)
2525  ic = 0
2526  iz = 0
2527  vi(:) = dzero
2528  di(:) = dzero
2529  viz(:) = dzero
2530  diz(:) = dzero
2531  nix(:) = dzero
2532  niy(:) = dzero
2533  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
2534  m = this%dis%con%ja(ipos)
2535  isympos = this%dis%con%jas(ipos)
2536  ihc = this%dis%con%ihc(isympos)
2537  area = this%dis%con%hwva(isympos)
2538  if (ihc == c3d_vertical) then
2539  !
2540  ! -- vertical connection
2541  iz = iz + 1
2542  !call this%dis%connection_normal(n, m, ihc, xn, yn, zn, ipos)
2543  call this%dis%connection_vector(n, m, nozee, this%sat(n), this%sat(m), &
2544  ihc, xc, yc, zc, dltot)
2545  cl1 = this%dis%con%cl1(isympos)
2546  cl2 = this%dis%con%cl2(isympos)
2547  if (m < n) then
2548  cl1 = this%dis%con%cl2(isympos)
2549  cl2 = this%dis%con%cl1(isympos)
2550  end if
2551  ooclsum = done / (cl1 + cl2)
2552  diz(iz) = dltot * cl1 * ooclsum
2553  qz = flowja(ipos)
2554  if (n > m) qz = -qz
2555  viz(iz) = qz / area
2556  else
2557  !
2558  ! -- horizontal connection
2559  ic = ic + 1
2560  dz = thksatnm(this%ibound(n), this%ibound(m), &
2561  this%icelltype(n), this%icelltype(m), &
2562  this%inewton, ihc, &
2563  this%hnew(n), this%hnew(m), this%sat(n), this%sat(m), &
2564  this%dis%top(n), this%dis%top(m), this%dis%bot(n), &
2565  this%dis%bot(m))
2566  area = area * dz
2567  call this%dis%connection_normal(n, m, ihc, xn, yn, zn, ipos)
2568  call this%dis%connection_vector(n, m, nozee, this%sat(n), this%sat(m), &
2569  ihc, xc, yc, zc, dltot)
2570  cl1 = this%dis%con%cl1(isympos)
2571  cl2 = this%dis%con%cl2(isympos)
2572  if (m < n) then
2573  cl1 = this%dis%con%cl2(isympos)
2574  cl2 = this%dis%con%cl1(isympos)
2575  end if
2576  ooclsum = done / (cl1 + cl2)
2577  nix(ic) = -xn
2578  niy(ic) = -yn
2579  di(ic) = dltot * cl1 * ooclsum
2580  if (area > dzero) then
2581  vi(ic) = flowja(ipos) / area
2582  else
2583  vi(ic) = dzero
2584  end if
2585  end if
2586  end do
2587  !
2588  ! -- Look through edge flows that may have been provided by an exchange
2589  ! and incorporate them into the averaging arrays
2590  do m = 1, this%nedges
2591  if (this%nodedge(m) == n) then
2592  !
2593  ! -- propsedge: (Q, area, nx, ny, distance)
2594  ihc = this%ihcedge(m)
2595  area = this%propsedge(2, m)
2596  if (ihc == c3d_vertical) then
2597  iz = iz + 1
2598  viz(iz) = this%propsedge(1, m) / area
2599  diz(iz) = this%propsedge(5, m)
2600  else
2601  ic = ic + 1
2602  nix(ic) = -this%propsedge(3, m)
2603  niy(ic) = -this%propsedge(4, m)
2604  di(ic) = this%propsedge(5, m)
2605  if (area > dzero) then
2606  vi(ic) = this%propsedge(1, m) / area
2607  else
2608  vi(ic) = dzero
2609  end if
2610  end if
2611  end if
2612  end do
2613  !
2614  ! -- Assign number of vertical and horizontal connections
2615  ncz = iz
2616  nc = ic
2617  !
2618  ! -- calculate z weight (wiz) and z velocity
2619  if (ncz == 1) then
2620  wiz(1) = done
2621  else
2622  dsumz = dzero
2623  do iz = 1, ncz
2624  dsumz = dsumz + diz(iz)
2625  end do
2626  denom = (ncz - done)
2627  if (denom < dzero) denom = dzero
2628  dsumz = dsumz + dem10 * dsumz
2629  do iz = 1, ncz
2630  if (dsumz > dzero) wiz(iz) = done - diz(iz) / dsumz
2631  if (denom > 0) then
2632  wiz(iz) = wiz(iz) / denom
2633  else
2634  wiz(iz) = dzero
2635  end if
2636  end do
2637  end if
2638  vz = dzero
2639  do iz = 1, ncz
2640  vz = vz + wiz(iz) * viz(iz)
2641  end do
2642  !
2643  ! -- distance-based weighting
2644  nc = ic
2645  dsumx = dzero
2646  dsumy = dzero
2647  dsumz = dzero
2648  do ic = 1, nc
2649  wix(ic) = di(ic) * abs(nix(ic))
2650  wiy(ic) = di(ic) * abs(niy(ic))
2651  dsumx = dsumx + wix(ic)
2652  dsumy = dsumy + wiy(ic)
2653  end do
2654  !
2655  ! -- Finish computing omega weights. Add a tiny bit
2656  ! to dsum so that the normalized omega weight later
2657  ! evaluates to (essentially) 1 in the case of a single
2658  ! relevant connection, avoiding 0/0.
2659  dsumx = dsumx + dem10 * dsumx
2660  dsumy = dsumy + dem10 * dsumy
2661  do ic = 1, nc
2662  wix(ic) = (dsumx - wix(ic)) * abs(nix(ic))
2663  wiy(ic) = (dsumy - wiy(ic)) * abs(niy(ic))
2664  end do
2665  !
2666  ! -- compute B weights
2667  dsumx = dzero
2668  dsumy = dzero
2669  do ic = 1, nc
2670  bix(ic) = wix(ic) * sign(done, nix(ic))
2671  biy(ic) = wiy(ic) * sign(done, niy(ic))
2672  dsumx = dsumx + wix(ic) * abs(nix(ic))
2673  dsumy = dsumy + wiy(ic) * abs(niy(ic))
2674  end do
2675  if (dsumx > dzero) dsumx = done / dsumx
2676  if (dsumy > dzero) dsumy = done / dsumy
2677  axy = dzero
2678  ayx = dzero
2679  do ic = 1, nc
2680  bix(ic) = bix(ic) * dsumx
2681  biy(ic) = biy(ic) * dsumy
2682  axy = axy + bix(ic) * niy(ic)
2683  ayx = ayx + biy(ic) * nix(ic)
2684  end do
2685  !
2686  ! -- Calculate specific discharge. The divide by zero checking below
2687  ! is problematic for cells with only one flow, such as can happen
2688  ! with triangular cells in corners. In this case, the resulting
2689  ! cell velocity will be calculated as zero. The method should be
2690  ! improved so that edge flows of zero are included in these
2691  ! calculations. But this needs to be done with consideration for LGR
2692  ! cases in which flows are submitted from an exchange.
2693  vx = dzero
2694  vy = dzero
2695  do ic = 1, nc
2696  vx = vx + (bix(ic) - axy * biy(ic)) * vi(ic)
2697  vy = vy + (biy(ic) - ayx * bix(ic)) * vi(ic)
2698  end do
2699  denom = done - axy * ayx
2700  if (denom /= dzero) then
2701  vx = vx / denom
2702  vy = vy / denom
2703  end if
2704  !
2705  this%spdis(1, n) = vx
2706  this%spdis(2, n) = vy
2707  this%spdis(3, n) = vz
2708  !
2709  end do
2710  !
2711  ! -- cleanup
2712  deallocate (vi)
2713  deallocate (di)
2714  deallocate (nix)
2715  deallocate (niy)
2716  deallocate (wix)
2717  deallocate (wiy)
2718  deallocate (wiz)
2719  deallocate (bix)
2720  deallocate (biy)
2721  !
2722  ! -- Return
2723  return
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 2841 of file gwf-npf.f90.

2842  ! -- dummy
2843  class(GwfNpfType) :: this !< this NPF instance
2844  integer(I4B) :: n !< node n
2845  integer(I4B) :: m !< node m
2846  integer(I4B) :: ihc !< 1 = horizontal connection, 0 for vertical
2847  ! -- return
2848  real(DP) :: satThickness !< saturated thickness
2849  !
2850  satthickness = thksatnm(this%ibound(n), &
2851  this%ibound(m), &
2852  this%icelltype(n), &
2853  this%icelltype(m), &
2854  this%inewton, &
2855  ihc, &
2856  this%hnew(n), &
2857  this%hnew(m), &
2858  this%sat(n), &
2859  this%sat(m), &
2860  this%dis%top(n), &
2861  this%dis%top(m), &
2862  this%dis%bot(n), &
2863  this%dis%bot(m))
2864  !
2865  ! -- Return
2866  return
Here is the call graph for this function:

◆ check_options()

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

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

1427  ! -- modules
1429  use constantsmodule, only: linelength
1430  ! -- dummy
1431  class(GwfNpftype) :: this
1432  !
1433  ! -- set omega value used for saturation calculations
1434  if (this%inewton > 0) then
1435  this%satomega = dem6
1436  end if
1437  !
1438  if (this%inewton > 0) then
1439  if (this%iperched > 0) then
1440  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1441  'BE USED WITH PERCHED OPTION.'
1442  call store_error(errmsg)
1443  end if
1444  if (this%ivarcv > 0) then
1445  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1446  'BE USED WITH VARIABLECV OPTION.'
1447  call store_error(errmsg)
1448  end if
1449  if (this%irewet > 0) then
1450  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1451  'BE USED WITH REWET OPTION.'
1452  call store_error(errmsg)
1453  end if
1454  end if
1455  !
1456  if (this%ixt3d /= 0) then
1457  if (this%icellavg > 0) then
1458  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. '// &
1459  'ALTERNATIVE_CELL_AVERAGING OPTION '// &
1460  'CANNOT BE USED WITH XT3D OPTION.'
1461  call store_error(errmsg)
1462  end if
1463  if (this%ithickstrt > 0) then
1464  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. THICKSTRT OPTION '// &
1465  'CANNOT BE USED WITH XT3D OPTION.'
1466  call store_error(errmsg)
1467  end if
1468  if (this%iperched > 0) then
1469  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. PERCHED OPTION '// &
1470  'CANNOT BE USED WITH XT3D OPTION.'
1471  call store_error(errmsg)
1472  end if
1473  if (this%ivarcv > 0) then
1474  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. VARIABLECV OPTION '// &
1475  'CANNOT BE USED WITH XT3D OPTION.'
1476  call store_error(errmsg)
1477  end if
1478  end if
1479  !
1480  ! -- Terminate if errors
1481  if (count_errors() > 0) then
1482  call store_error_filename(this%input_fname)
1483  end if
1484  !
1485  ! -- Return
1486  return
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
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 2349 of file gwf-npf.f90.

2350  ! -- return
2351  real(DP) :: hy
2352  ! -- dummy
2353  class(GwfNpfType) :: this
2354  integer(I4B), intent(in) :: n
2355  integer(I4B), intent(in) :: m
2356  integer(I4B), intent(in) :: ihc
2357  integer(I4B), intent(in), optional :: ipos
2358  real(DP), dimension(3), intent(in), optional :: vg
2359  ! -- local
2360  integer(I4B) :: iipos
2361  real(DP) :: hy11, hy22, hy33
2362  real(DP) :: ang1, ang2, ang3
2363  real(DP) :: vg1, vg2, vg3
2364  !
2365  ! -- Initialize
2366  iipos = 0
2367  if (present(ipos)) iipos = ipos
2368  hy11 = this%k11(n)
2369  hy22 = this%k11(n)
2370  hy33 = this%k11(n)
2371  hy22 = this%k22(n)
2372  hy33 = this%k33(n)
2373  !
2374  ! -- Calculate effective K based on whether connection is vertical
2375  ! or horizontal
2376  if (ihc == c3d_vertical) then
2377  !
2378  ! -- Handle rotated anisotropy case that would affect the effective
2379  ! vertical hydraulic conductivity
2380  hy = hy33
2381  if (this%iangle2 > 0) then
2382  if (present(vg)) then
2383  vg1 = vg(1)
2384  vg2 = vg(2)
2385  vg3 = vg(3)
2386  else
2387  call this%dis%connection_normal(n, m, ihc, vg1, vg2, vg3, iipos)
2388  end if
2389  ang1 = this%angle1(n)
2390  ang2 = this%angle2(n)
2391  ang3 = dzero
2392  if (this%iangle3 > 0) ang3 = this%angle3(n)
2393  hy = hyeff(hy11, hy22, hy33, ang1, ang2, ang3, vg1, vg2, vg3, &
2394  this%iavgkeff)
2395  end if
2396  !
2397  else
2398  !
2399  ! -- Handle horizontal case
2400  hy = hy11
2401  if (this%ik22 > 0) then
2402  if (present(vg)) then
2403  vg1 = vg(1)
2404  vg2 = vg(2)
2405  vg3 = vg(3)
2406  else
2407  call this%dis%connection_normal(n, m, ihc, vg1, vg2, vg3, iipos)
2408  end if
2409  ang1 = dzero
2410  ang2 = dzero
2411  ang3 = dzero
2412  if (this%iangle1 > 0) then
2413  ang1 = this%angle1(n)
2414  if (this%iangle2 > 0) then
2415  ang2 = this%angle2(n)
2416  if (this%iangle3 > 0) ang3 = this%angle3(n)
2417  end if
2418  end if
2419  hy = hyeff(hy11, hy22, hy33, ang1, ang2, ang3, vg1, vg2, vg3, &
2420  this%iavgkeff)
2421  end if
2422  !
2423  end if
2424  !
2425  ! -- Return
2426  return
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 2794 of file gwf-npf.f90.

2795  ! -- dummy
2796  class(GwfNpfType) :: this
2797  integer(I4B), intent(in) :: nedges
2798  !
2799  this%nedges = this%nedges + nedges
2800  !
2801  ! -- Return
2802  return

◆ log_griddata()

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

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

1492  ! -- modules
1494  ! -- dummy
1495  class(GwfNpfType) :: this
1496  type(GwfNpfParamFoundType), intent(in) :: found
1497  !
1498  write (this%iout, '(1x,a)') 'Setting NPF Griddata'
1499  !
1500  if (found%icelltype) then
1501  write (this%iout, '(4x,a)') 'ICELLTYPE set from input file'
1502  end if
1503  !
1504  if (found%k) then
1505  write (this%iout, '(4x,a)') 'K set from input file'
1506  end if
1507  !
1508  if (found%k33) then
1509  write (this%iout, '(4x,a)') 'K33 set from input file'
1510  else
1511  write (this%iout, '(4x,a)') 'K33 not provided. Setting K33 = K.'
1512  end if
1513  !
1514  if (found%k22) then
1515  write (this%iout, '(4x,a)') 'K22 set from input file'
1516  else
1517  write (this%iout, '(4x,a)') 'K22 not provided. Setting K22 = K.'
1518  end if
1519  !
1520  if (found%wetdry) then
1521  write (this%iout, '(4x,a)') 'WETDRY set from input file'
1522  end if
1523  !
1524  if (found%angle1) then
1525  write (this%iout, '(4x,a)') 'ANGLE1 set from input file'
1526  end if
1527  !
1528  if (found%angle2) then
1529  write (this%iout, '(4x,a)') 'ANGLE2 set from input file'
1530  end if
1531  !
1532  if (found%angle3) then
1533  write (this%iout, '(4x,a)') 'ANGLE3 set from input file'
1534  end if
1535  !
1536  write (this%iout, '(1x,a,/)') 'End Setting NPF Griddata'
1537  !
1538  ! -- Return
1539  return

◆ log_options()

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

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

1255  ! -- modules
1256  use kindmodule, only: lgp
1258  ! -- dummy
1259  class(GwfNpftype) :: this
1260  ! -- locals
1261  type(GwfNpfParamFoundType), intent(in) :: found
1262  !
1263  write (this%iout, '(1x,a)') 'Setting NPF Options'
1264  if (found%iprflow) &
1265  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be printed &
1266  &to listing file whenever ICBCFL is not zero.'
1267  if (found%ipakcb) &
1268  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be saved &
1269  &to binary file whenever ICBCFL is not zero.'
1270  if (found%cellavg) &
1271  write (this%iout, '(4x,a,i0)') 'Alternative cell averaging [1=logarithmic, &
1272  &2=AMT-LMK, 3=AMT-HMK] set to: ', &
1273  this%icellavg
1274  if (found%ithickstrt) &
1275  write (this%iout, '(4x,a)') 'THICKSTRT option has been activated.'
1276  if (found%iperched) &
1277  write (this%iout, '(4x,a)') 'Vertical flow will be adjusted for perched &
1278  &conditions.'
1279  if (found%ivarcv) &
1280  write (this%iout, '(4x,a)') 'Vertical conductance varies with water table.'
1281  if (found%idewatcv) &
1282  write (this%iout, '(4x,a)') 'Vertical conductance accounts for dewatered &
1283  &portion of an underlying cell.'
1284  if (found%ixt3d) write (this%iout, '(4x,a)') 'XT3D formulation is selected.'
1285  if (found%ixt3drhs) &
1286  write (this%iout, '(4x,a)') 'XT3D RHS formulation is selected.'
1287  if (found%isavspdis) &
1288  write (this%iout, '(4x,a)') 'Specific discharge will be calculated at cell &
1289  &centers and written to DATA-SPDIS in budget &
1290  &file when requested.'
1291  if (found%isavsat) &
1292  write (this%iout, '(4x,a)') 'Saturation will be written to DATA-SAT in &
1293  &budget file when requested.'
1294  if (found%ik22overk) &
1295  write (this%iout, '(4x,a)') 'Values specified for K22 are anisotropy &
1296  &ratios and will be multiplied by K before &
1297  &being used in calculations.'
1298  if (found%ik33overk) &
1299  write (this%iout, '(4x,a)') 'Values specified for K33 are anisotropy &
1300  &ratios and will be multiplied by K before &
1301  &being used in calculations.'
1302  if (found%inewton) &
1303  write (this%iout, '(4x,a)') 'NEWTON-RAPHSON method disabled for unconfined &
1304  &cells'
1305  if (found%satomega) &
1306  write (this%iout, '(4x,a,1pg15.6)') 'Saturation omega: ', this%satomega
1307  if (found%irewet) &
1308  write (this%iout, '(4x,a)') 'Rewetting is active.'
1309  if (found%wetfct) &
1310  write (this%iout, '(4x,a,1pg15.6)') &
1311  'Wetting factor (WETFCT) has been set to: ', this%wetfct
1312  if (found%iwetit) &
1313  write (this%iout, '(4x,a,i5)') &
1314  'Wetting iteration interval (IWETIT) has been set to: ', this%iwetit
1315  if (found%ihdwet) &
1316  write (this%iout, '(4x,a,i5)') &
1317  'Head rewet equation (IHDWET) has been set to: ', this%ihdwet
1318  write (this%iout, '(1x,a,/)') 'End Setting NPF Options'
1319  !
1320  ! -- Return
1321  return
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 254 of file gwf-npf.f90.

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

◆ 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
440  !
441  ! -- Return
442  return
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ 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 288 of file gwf-npf.f90.

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

◆ 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 447 of file gwf-npf.f90.

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

◆ 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 760 of file gwf-npf.f90.

761  ! -- dummy
762  class(GwfNpfType) :: this
763  real(DP), intent(inout), dimension(:) :: hnew
764  real(DP), intent(inout), dimension(:) :: flowja
765  ! -- local
766  integer(I4B) :: n, ipos, m
767  real(DP) :: qnm
768  !
769  ! -- Calculate the flow across each cell face and store in flowja
770  !
771  if (this%ixt3d /= 0) then
772  call this%xt3d%xt3d_flowja(hnew, flowja)
773  else
774  !
775  do n = 1, this%dis%nodes
776  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
777  m = this%dis%con%ja(ipos)
778  if (m < n) cycle
779  call this%qcalc(n, m, hnew(n), hnew(m), ipos, qnm)
780  flowja(ipos) = qnm
781  flowja(this%dis%con%isym(ipos)) = -qnm
782  end do
783  end do
784  !
785  end if
786  !
787  ! -- Return
788  return

◆ 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 152 of file gwf-npf.f90.

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

◆ npf_da()

subroutine gwfnpfmodule::npf_da ( class(gwfnpftype this)

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

976  ! -- modules
979  ! -- dummy
980  class(GwfNpftype) :: this
981  !
982  ! -- Deallocate input memory
983  call memorylist_remove(this%name_model, 'NPF', idm_context)
984  !
985  ! -- TVK
986  if (this%intvk /= 0) then
987  call this%tvk%da()
988  deallocate (this%tvk)
989  end if
990  !
991  ! -- VSC
992  if (this%invsc /= 0) then
993  nullify (this%vsc)
994  end if
995  !
996  ! -- Strings
997  !
998  ! -- Scalars
999  call mem_deallocate(this%iname)
1000  call mem_deallocate(this%ixt3d)
1001  call mem_deallocate(this%ixt3drhs)
1002  call mem_deallocate(this%satomega)
1003  call mem_deallocate(this%hnoflo)
1004  call mem_deallocate(this%hdry)
1005  call mem_deallocate(this%icellavg)
1006  call mem_deallocate(this%iavgkeff)
1007  call mem_deallocate(this%ik22)
1008  call mem_deallocate(this%ik33)
1009  call mem_deallocate(this%iperched)
1010  call mem_deallocate(this%ivarcv)
1011  call mem_deallocate(this%idewatcv)
1012  call mem_deallocate(this%ithickstrt)
1013  call mem_deallocate(this%isavspdis)
1014  call mem_deallocate(this%isavsat)
1015  call mem_deallocate(this%icalcspdis)
1016  call mem_deallocate(this%irewet)
1017  call mem_deallocate(this%wetfct)
1018  call mem_deallocate(this%iwetit)
1019  call mem_deallocate(this%ihdwet)
1020  call mem_deallocate(this%ibotnode)
1021  call mem_deallocate(this%iwetdry)
1022  call mem_deallocate(this%iangle1)
1023  call mem_deallocate(this%iangle2)
1024  call mem_deallocate(this%iangle3)
1025  call mem_deallocate(this%nedges)
1026  call mem_deallocate(this%lastedge)
1027  call mem_deallocate(this%ik22overk)
1028  call mem_deallocate(this%ik33overk)
1029  call mem_deallocate(this%intvk)
1030  call mem_deallocate(this%invsc)
1031  call mem_deallocate(this%kchangeper)
1032  call mem_deallocate(this%kchangestp)
1033  !
1034  ! -- Deallocate arrays
1035  deallocate (this%aname)
1036  call mem_deallocate(this%ithickstartflag)
1037  call mem_deallocate(this%icelltype)
1038  call mem_deallocate(this%k11)
1039  call mem_deallocate(this%k22)
1040  call mem_deallocate(this%k33)
1041  call mem_deallocate(this%k11input)
1042  call mem_deallocate(this%k22input)
1043  call mem_deallocate(this%k33input)
1044  call mem_deallocate(this%sat)
1045  call mem_deallocate(this%condsat)
1046  call mem_deallocate(this%wetdry)
1047  call mem_deallocate(this%angle1)
1048  call mem_deallocate(this%angle2)
1049  call mem_deallocate(this%angle3)
1050  call mem_deallocate(this%nodedge)
1051  call mem_deallocate(this%ihcedge)
1052  call mem_deallocate(this%propsedge)
1053  call mem_deallocate(this%spdis)
1054  call mem_deallocate(this%nodekchange)
1055  !
1056  ! -- deallocate parent
1057  call this%NumericalPackageType%da()
1058 
1059  ! pointers
1060  this%hnew => null()
1061 
subroutine, public memorylist_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 198 of file gwf-npf.f90.

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

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

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

◆ 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 717 of file gwf-npf.f90.

718  ! -- dummy
719  class(GwfNpfType) :: this
720  integer(I4B), intent(in) :: neqmod
721  real(DP), dimension(neqmod), intent(inout) :: x
722  real(DP), dimension(neqmod), intent(in) :: xtemp
723  real(DP), dimension(neqmod), intent(inout) :: dx
724  integer(I4B), intent(inout) :: inewtonur
725  real(DP), intent(inout) :: dxmax
726  integer(I4B), intent(inout) :: locmax
727  ! -- local
728  integer(I4B) :: n
729  real(DP) :: botm
730  real(DP) :: xx
731  real(DP) :: dxx
732  !
733  ! -- Newton-Raphson under-relaxation
734  do n = 1, this%dis%nodes
735  if (this%ibound(n) < 1) cycle
736  if (this%icelltype(n) > 0) then
737  botm = this%dis%bot(this%ibotnode(n))
738  ! -- only apply Newton-Raphson under-relaxation if
739  ! solution head is below the bottom of the model
740  if (x(n) < botm) then
741  inewtonur = 1
742  xx = xtemp(n) * (done - dp9) + botm * dp9
743  dxx = x(n) - xx
744  if (abs(dxx) > abs(dxmax)) then
745  locmax = n
746  dxmax = dxx
747  end if
748  x(n) = xx
749  dx(n) = dzero
750  end if
751  end if
752  end do
753  !
754  ! -- Return
755  return

◆ 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 933 of file gwf-npf.f90.

934  ! -- modules
935  use tdismodule, only: kper, kstp
936  use constantsmodule, only: lenbigline
937  ! -- dummy
938  class(GwfNpfType) :: this
939  integer(I4B), intent(in) :: ibudfl
940  real(DP), intent(inout), dimension(:) :: flowja
941  ! -- local
942  character(len=LENBIGLINE) :: line
943  character(len=30) :: tempstr
944  integer(I4B) :: n, ipos, m
945  real(DP) :: qnm
946  ! -- formats
947  character(len=*), parameter :: fmtiprflow = &
948  &"(/,4x,'CALCULATED INTERCELL FLOW FOR PERIOD ', i0, ' STEP ', i0)"
949  !
950  ! -- Write flowja to list file if requested
951  if (ibudfl /= 0 .and. this%iprflow > 0) then
952  write (this%iout, fmtiprflow) kper, kstp
953  do n = 1, this%dis%nodes
954  line = ''
955  call this%dis%noder_to_string(n, tempstr)
956  line = trim(tempstr)//':'
957  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
958  m = this%dis%con%ja(ipos)
959  call this%dis%noder_to_string(m, tempstr)
960  line = trim(line)//' '//trim(tempstr)
961  qnm = flowja(ipos)
962  write (tempstr, '(1pg15.6)') qnm
963  line = trim(line)//' '//trim(adjustl(tempstr))
964  end do
965  write (this%iout, '(a)') trim(line)
966  end do
967  end if
968  !
969  ! -- Return
970  return
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 362 of file gwf-npf.f90.

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

◆ 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 893 of file gwf-npf.f90.

894  ! -- dummy
895  class(GwfNpfType) :: this
896  real(DP), dimension(:), intent(in) :: flowja
897  integer(I4B), intent(in) :: icbcfl
898  integer(I4B), intent(in) :: icbcun
899  ! -- local
900  integer(I4B) :: ibinun
901  !
902  ! -- Set unit number for binary output
903  if (this%ipakcb < 0) then
904  ibinun = icbcun
905  elseif (this%ipakcb == 0) then
906  ibinun = 0
907  else
908  ibinun = this%ipakcb
909  end if
910  if (icbcfl == 0) ibinun = 0
911  !
912  ! -- Write the face flows if requested
913  if (ibinun /= 0) then
914  call this%dis%record_connection_array(flowja, ibinun, this%iout)
915  end if
916  !
917  ! -- Calculate specific discharge at cell centers and write, if requested
918  if (this%isavspdis /= 0) then
919  if (ibinun /= 0) call this%sav_spdis(ibinun)
920  end if
921  !
922  ! -- Save saturation, if requested
923  if (this%isavsat /= 0) then
924  if (ibinun /= 0) call this%sav_sat(ibinun)
925  end if
926  !
927  ! -- Return
928  return

◆ prepcheck()

subroutine gwfnpfmodule::prepcheck ( class(gwfnpftype this)

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

1637  ! -- modules
1638  use constantsmodule, only: linelength, dpio180
1640  ! -- dummy
1641  class(GwfNpfType) :: this
1642  ! -- local
1643  character(len=24), dimension(:), pointer :: aname
1644  character(len=LINELENGTH) :: cellstr, errmsg
1645  integer(I4B) :: nerr, n
1646  ! -- format
1647  character(len=*), parameter :: fmtkerr = &
1648  &"(1x, 'Hydraulic property ',a,' is <= 0 for cell ',a, ' ', 1pg15.6)"
1649  character(len=*), parameter :: fmtkerr2 = &
1650  &"(1x, '... ', i0,' additional errors not shown for ',a)"
1651  !
1652  ! -- initialize
1653  aname => this%aname
1654  !
1655  ! -- check k11
1656  nerr = 0
1657  do n = 1, size(this%k11)
1658  if (this%k11(n) <= dzero) then
1659  nerr = nerr + 1
1660  if (nerr <= 20) then
1661  call this%dis%noder_to_string(n, cellstr)
1662  write (errmsg, fmtkerr) trim(adjustl(aname(2))), trim(cellstr), &
1663  this%k11(n)
1664  call store_error(errmsg)
1665  end if
1666  end if
1667  end do
1668  if (nerr > 20) then
1669  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(2)))
1670  call store_error(errmsg)
1671  end if
1672  !
1673  ! -- check k33 because it was read
1674  if (this%ik33 /= 0) then
1675  !
1676  ! -- Check to make sure values are greater than or equal to zero
1677  nerr = 0
1678  do n = 1, size(this%k33)
1679  if (this%ik33overk /= 0) this%k33(n) = this%k33(n) * this%k11(n)
1680  if (this%k33(n) <= dzero) then
1681  nerr = nerr + 1
1682  if (nerr <= 20) then
1683  call this%dis%noder_to_string(n, cellstr)
1684  write (errmsg, fmtkerr) trim(adjustl(aname(3))), trim(cellstr), &
1685  this%k33(n)
1686  call store_error(errmsg)
1687  end if
1688  end if
1689  end do
1690  if (nerr > 20) then
1691  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(3)))
1692  call store_error(errmsg)
1693  end if
1694  end if
1695  !
1696  ! -- check k22 because it was read
1697  if (this%ik22 /= 0) then
1698  !
1699  ! -- Check to make sure that angles are available
1700  if (this%dis%con%ianglex == 0) then
1701  write (errmsg, '(a)') 'Error. ANGLDEGX not provided in '// &
1702  'discretization file, but K22 was specified. '
1703  call store_error(errmsg)
1704  end if
1705  !
1706  ! -- Check to make sure values are greater than or equal to zero
1707  nerr = 0
1708  do n = 1, size(this%k22)
1709  if (this%ik22overk /= 0) this%k22(n) = this%k22(n) * this%k11(n)
1710  if (this%k22(n) <= dzero) then
1711  nerr = nerr + 1
1712  if (nerr <= 20) then
1713  call this%dis%noder_to_string(n, cellstr)
1714  write (errmsg, fmtkerr) trim(adjustl(aname(4))), trim(cellstr), &
1715  this%k22(n)
1716  call store_error(errmsg)
1717  end if
1718  end if
1719  end do
1720  if (nerr > 20) then
1721  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(4)))
1722  call store_error(errmsg)
1723  end if
1724  end if
1725  !
1726  ! -- check for wetdry conflicts
1727  if (this%irewet == 1) then
1728  if (this%iwetdry == 0) then
1729  write (errmsg, '(a, a, a)') 'Error in GRIDDATA block: ', &
1730  trim(adjustl(aname(5))), ' not found.'
1731  call store_error(errmsg)
1732  end if
1733  end if
1734  !
1735  ! -- Check for angle conflicts
1736  if (this%iangle1 /= 0) then
1737  do n = 1, size(this%angle1)
1738  this%angle1(n) = this%angle1(n) * dpio180
1739  end do
1740  else
1741  if (this%ixt3d /= 0) then
1742  this%iangle1 = 1
1743  write (this%iout, '(a)') 'XT3D IN USE, BUT ANGLE1 NOT SPECIFIED. '// &
1744  'SETTING ANGLE1 TO ZERO.'
1745  do n = 1, size(this%angle1)
1746  this%angle1(n) = dzero
1747  end do
1748  end if
1749  end if
1750  if (this%iangle2 /= 0) then
1751  if (this%iangle1 == 0) then
1752  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE1. '// &
1753  'ANGLE2 REQUIRES ANGLE1. '
1754  call store_error(errmsg)
1755  end if
1756  if (this%iangle3 == 0) then
1757  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE3. '// &
1758  'SPECIFY BOTH OR NEITHER ONE. '
1759  call store_error(errmsg)
1760  end if
1761  do n = 1, size(this%angle2)
1762  this%angle2(n) = this%angle2(n) * dpio180
1763  end do
1764  end if
1765  if (this%iangle3 /= 0) then
1766  if (this%iangle1 == 0) then
1767  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE1. '// &
1768  'ANGLE3 REQUIRES ANGLE1. '
1769  call store_error(errmsg)
1770  end if
1771  if (this%iangle2 == 0) then
1772  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE2. '// &
1773  'SPECIFY BOTH OR NEITHER ONE. '
1774  call store_error(errmsg)
1775  end if
1776  do n = 1, size(this%angle3)
1777  this%angle3(n) = this%angle3(n) * dpio180
1778  end do
1779  end if
1780  !
1781  ! -- terminate if data errors
1782  if (count_errors() > 0) then
1783  call store_error_filename(this%input_fname)
1784  end if
1785  !
1786  ! -- Return
1787  return
real(dp), parameter dpio180
real constant
Definition: Constants.f90:129
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 1800 of file gwf-npf.f90.

1801  ! -- modules
1802  use constantsmodule, only: linelength
1804  ! -- dummy
1805  class(GwfNpfType) :: this !< the instance of the NPF package
1806  ! -- local
1807  integer(I4B) :: n, m, ii, nn
1808  real(DP) :: hyn, hym
1809  real(DP) :: satn, topn, botn
1810  integer(I4B) :: nextn
1811  real(DP) :: minbot, botm
1812  logical :: finished
1813  character(len=LINELENGTH) :: cellstr, errmsg
1814  ! -- format
1815  character(len=*), parameter :: fmtcnv = &
1816  "(1X,'CELL ', A, &
1817  &' ELIMINATED BECAUSE ALL HYDRAULIC CONDUCTIVITIES TO NODE ARE 0.')"
1818  character(len=*), parameter :: fmtnct = &
1819  &"(1X,'Negative cell thickness at cell ', A)"
1820  character(len=*), parameter :: fmtihbe = &
1821  &"(1X,'Initial head, bottom elevation:',1P,2G13.5)"
1822  character(len=*), parameter :: fmttebe = &
1823  &"(1X,'Top elevation, bottom elevation:',1P,2G13.5)"
1824  !
1825  do n = 1, this%dis%nodes
1826  this%ithickstartflag(n) = 0
1827  end do
1828  !
1829  ! -- Insure that each cell has at least one non-zero transmissive parameter
1830  ! Note that a cell can be deactivated even if it has a valid connection
1831  ! to another model.
1832  nodeloop: do n = 1, this%dis%nodes
1833  !
1834  ! -- Skip if already inactive
1835  if (this%ibound(n) == 0) then
1836  if (this%irewet /= 0) then
1837  if (this%wetdry(n) == dzero) cycle nodeloop
1838  else
1839  cycle nodeloop
1840  end if
1841  end if
1842  !
1843  ! -- Cycle if k11 is not zero
1844  if (this%k11(n) /= dzero) cycle nodeloop
1845  !
1846  ! -- Cycle if at least one vertical connection has non-zero k33
1847  ! for n and m
1848  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
1849  m = this%dis%con%ja(ii)
1850  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0) then
1851  hyn = this%k11(n)
1852  if (this%ik33 /= 0) hyn = this%k33(n)
1853  if (hyn /= dzero) then
1854  hym = this%k11(m)
1855  if (this%ik33 /= 0) hym = this%k33(m)
1856  if (hym /= dzero) cycle
1857  end if
1858  end if
1859  end do
1860  !
1861  ! -- If this part of the loop is reached, then all connections have
1862  ! zero transmissivity, so convert to noflow.
1863  this%ibound(n) = 0
1864  this%hnew(n) = this%hnoflo
1865  if (this%irewet /= 0) this%wetdry(n) = dzero
1866  call this%dis%noder_to_string(n, cellstr)
1867  write (this%iout, fmtcnv) trim(adjustl(cellstr))
1868  !
1869  end do nodeloop
1870  !
1871  ! -- Preprocess cell status and heads based on initial conditions
1872  if (this%inewton == 0) then
1873  !
1874  ! -- For standard formulation (non-Newton) call wetdry routine
1875  call this%wd(0, this%hnew)
1876  else
1877  !
1878  ! -- Newton formulation, so adjust heads to be above bottom
1879  ! (Not used in present formulation because variable cv
1880  ! cannot be used with Newton)
1881  if (this%ivarcv == 1) then
1882  do n = 1, this%dis%nodes
1883  if (this%hnew(n) < this%dis%bot(n)) then
1884  this%hnew(n) = this%dis%bot(n) + dem6
1885  end if
1886  end do
1887  end if
1888  end if
1889  !
1890  ! -- If THCKSTRT is not active, then loop through icelltype and replace
1891  ! any negative values with 1.
1892  if (this%ithickstrt == 0) then
1893  do n = 1, this%dis%nodes
1894  if (this%icelltype(n) < 0) then
1895  this%icelltype(n) = 1
1896  end if
1897  end do
1898  end if
1899  !
1900  ! -- Initialize sat to zero for ibound=0 cells, unless the cell can
1901  ! rewet. Initialize sat to the saturated fraction based on strt
1902  ! if icelltype is negative and the THCKSTRT option is in effect.
1903  ! Initialize sat to 1.0 for all other cells in order to calculate
1904  ! condsat in next section.
1905  do n = 1, this%dis%nodes
1906  if (this%ibound(n) == 0) then
1907  this%sat(n) = done
1908  if (this%icelltype(n) < 0 .and. this%ithickstrt /= 0) then
1909  this%ithickstartflag(n) = 1
1910  this%icelltype(n) = 0
1911  end if
1912  else
1913  topn = this%dis%top(n)
1914  botn = this%dis%bot(n)
1915  if (this%icelltype(n) < 0 .and. this%ithickstrt /= 0) then
1916  call this%thksat(n, this%ic%strt(n), satn)
1917  if (botn > this%ic%strt(n)) then
1918  call this%dis%noder_to_string(n, cellstr)
1919  write (errmsg, fmtnct) trim(adjustl(cellstr))
1920  call store_error(errmsg)
1921  write (errmsg, fmtihbe) this%ic%strt(n), botn
1922  call store_error(errmsg)
1923  end if
1924  this%ithickstartflag(n) = 1
1925  this%icelltype(n) = 0
1926  else
1927  satn = done
1928  if (botn > topn) then
1929  call this%dis%noder_to_string(n, cellstr)
1930  write (errmsg, fmtnct) trim(adjustl(cellstr))
1931  call store_error(errmsg)
1932  write (errmsg, fmttebe) topn, botn
1933  call store_error(errmsg)
1934  end if
1935  end if
1936  this%sat(n) = satn
1937  end if
1938  end do
1939  if (count_errors() > 0) then
1940  call store_error_filename(this%input_fname)
1941  end if
1942  !
1943  ! -- Calculate condsat, but only if xt3d is not active. If xt3d is
1944  ! active, then condsat is allocated to size of zero.
1945  if (this%ixt3d == 0) then
1946  !
1947  ! -- Calculate the saturated conductance for all connections assuming
1948  ! that saturation is 1 (except for case where icelltype was entered
1949  ! as a negative value and THCKSTRT option in effect)
1950  do n = 1, this%dis%nodes
1951  call this%calc_condsat(n, .true.)
1952  end do
1953  !
1954  end if
1955  !
1956  ! -- Determine the lower most node
1957  if (this%igwfnewtonur /= 0) then
1958  call mem_reallocate(this%ibotnode, this%dis%nodes, 'IBOTNODE', &
1959  trim(this%memoryPath))
1960  do n = 1, this%dis%nodes
1961  !
1962  minbot = this%dis%bot(n)
1963  nn = n
1964  finished = .false.
1965  do while (.not. finished)
1966  nextn = 0
1967  !
1968  ! -- Go through the connecting cells
1969  do ii = this%dis%con%ia(nn) + 1, this%dis%con%ia(nn + 1) - 1
1970  !
1971  ! -- Set the m cell number
1972  m = this%dis%con%ja(ii)
1973  botm = this%dis%bot(m)
1974  !
1975  ! -- select vertical connections: ihc == 0
1976  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0) then
1977  if (m > nn .and. botm < minbot) then
1978  nextn = m
1979  minbot = botm
1980  end if
1981  end if
1982  end do
1983  if (nextn > 0) then
1984  nn = nextn
1985  else
1986  finished = .true.
1987  end if
1988  end do
1989  this%ibotnode(n) = nn
1990  end do
1991  end if
1992  !
1993  ! -- nullify unneeded gwf pointers
1994  this%igwfnewtonur => null()
1995  !
1996  ! -- Return
1997  return
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 2230 of file gwf-npf.f90.

2231  ! -- dummy
2232  class(GwfNpfType) :: this
2233  integer(I4B), intent(in) :: kiter
2234  integer(I4B), intent(in) :: node
2235  real(DP), intent(in) :: hm
2236  integer(I4B), intent(in) :: ibdm
2237  integer(I4B), intent(in) :: ihc
2238  real(DP), intent(inout), dimension(:) :: hnew
2239  integer(I4B), intent(out) :: irewet
2240  ! -- local
2241  integer(I4B) :: itflg
2242  real(DP) :: wd, awd, turnon, bbot
2243  !
2244  irewet = 0
2245  !
2246  ! -- Convert a dry cell to wet if it meets the criteria
2247  if (this%irewet > 0) then
2248  itflg = mod(kiter, this%iwetit)
2249  if (itflg == 0) then
2250  if (this%ibound(node) == 0 .and. this%wetdry(node) /= dzero) then
2251  !
2252  ! -- Calculate wetting elevation
2253  bbot = this%dis%bot(node)
2254  wd = this%wetdry(node)
2255  awd = wd
2256  if (wd < 0) awd = -wd
2257  turnon = bbot + awd
2258  !
2259  ! -- Check head in adjacent cells to see if wetting elevation has
2260  ! been reached
2261  if (ihc == c3d_vertical) then
2262  !
2263  ! -- check cell below
2264  if (ibdm > 0 .and. hm >= turnon) irewet = 1
2265  else
2266  if (wd > dzero) then
2267  !
2268  ! -- check horizontally adjacent cells
2269  if (ibdm > 0 .and. hm >= turnon) irewet = 1
2270  end if
2271  end if
2272  !
2273  if (irewet == 1) then
2274  ! -- rewet cell; use equation 3a if ihdwet=0; use equation 3b if
2275  ! ihdwet is not 0.
2276  if (this%ihdwet == 0) then
2277  hnew(node) = bbot + this%wetfct * (hm - bbot)
2278  else
2279  hnew(node) = bbot + this%wetfct * awd !(hm - bbot)
2280  end if
2281  this%ibound(node) = 30000
2282  end if
2283  end if
2284  end if
2285  end if
2286  !
2287  ! -- Return
2288  return

◆ sav_sat()

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

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

2760  ! -- dummy
2761  class(GwfNpfType) :: this
2762  integer(I4B), intent(in) :: ibinun
2763  ! -- local
2764  character(len=16) :: text
2765  character(len=16), dimension(1) :: auxtxt
2766  real(DP), dimension(1) :: a
2767  integer(I4B) :: n
2768  integer(I4B) :: naux
2769  !
2770  ! -- Write the header
2771  text = ' DATA-SAT'
2772  naux = 1
2773  auxtxt(:) = [' sat']
2774  call this%dis%record_srcdst_list_header(text, this%name_model, &
2775  this%packName, this%name_model, &
2776  this%packName, naux, auxtxt, ibinun, &
2777  this%dis%nodes, this%iout)
2778  !
2779  ! -- Write a zero for Q, and then write saturation as an aux variables
2780  do n = 1, this%dis%nodes
2781  a(1) = this%sat(n)
2782  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, a)
2783  end do
2784  !
2785  ! -- Return
2786  return

◆ sav_spdis()

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

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

2729  ! -- dummy
2730  class(GwfNpfType) :: this
2731  integer(I4B), intent(in) :: ibinun
2732  ! -- local
2733  character(len=16) :: text
2734  character(len=16), dimension(3) :: auxtxt
2735  integer(I4B) :: n
2736  integer(I4B) :: naux
2737  !
2738  ! -- Write the header
2739  text = ' DATA-SPDIS'
2740  naux = 3
2741  auxtxt(:) = [' qx', ' qy', ' qz']
2742  call this%dis%record_srcdst_list_header(text, this%name_model, &
2743  this%packName, this%name_model, &
2744  this%packName, naux, auxtxt, ibinun, &
2745  this%dis%nodes, this%iout)
2746  !
2747  ! -- Write a zero for Q, and then write qx, qy, qz as aux variables
2748  do n = 1, this%dis%nodes
2749  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, &
2750  this%spdis(:, n))
2751  end do
2752  !
2753  ! -- Return
2754  return

◆ 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 2807 of file gwf-npf.f90.

2809  ! -- dummy
2810  class(GwfNpfType) :: this
2811  integer(I4B), intent(in) :: nodedge
2812  integer(I4B), intent(in) :: ihcedge
2813  real(DP), intent(in) :: q
2814  real(DP), intent(in) :: area
2815  real(DP), intent(in) :: nx
2816  real(DP), intent(in) :: ny
2817  real(DP), intent(in) :: distance
2818  ! -- local
2819  integer(I4B) :: lastedge
2820  !
2821  this%lastedge = this%lastedge + 1
2822  lastedge = this%lastedge
2823  this%nodedge(lastedge) = nodedge
2824  this%ihcedge(lastedge) = ihcedge
2825  this%propsedge(1, lastedge) = q
2826  this%propsedge(2, lastedge) = area
2827  this%propsedge(3, lastedge) = nx
2828  this%propsedge(4, lastedge) = ny
2829  this%propsedge(5, lastedge) = distance
2830  !
2831  ! -- If this is the last edge, then the next call must be starting a new
2832  ! edge properties assignment loop, so need to reset lastedge to 0
2833  if (this%lastedge == this%nedges) this%lastedge = 0
2834  !
2835  ! -- Return
2836  return

◆ set_options()

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

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

1406  ! -- dummy
1407  class(GwfNpftype) :: this
1408  type(GwfNpfOptionsType), intent(in) :: options
1409  !
1410  this%icellavg = options%icellavg
1411  this%ithickstrt = options%ithickstrt
1412  this%iperched = options%iperched
1413  this%ivarcv = options%ivarcv
1414  this%idewatcv = options%idewatcv
1415  this%irewet = options%irewet
1416  this%wetfct = options%wetfct
1417  this%iwetit = options%iwetit
1418  this%ihdwet = options%ihdwet
1419  !
1420  ! -- Return
1421  return

◆ 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 819 of file gwf-npf.f90.

820  ! -- dummy
821  class(GwfNpfType) :: this
822  integer(I4B), intent(in) :: n
823  integer(I4B), intent(in) :: m
824  real(DP), intent(in) :: hn
825  real(DP), intent(in) :: hm
826  integer(I4B), intent(in) :: icon
827  real(DP), intent(inout) :: qnm
828  ! -- local
829  real(DP) :: hyn, hym
830  real(DP) :: condnm
831  real(DP) :: hntemp, hmtemp
832  integer(I4B) :: ihc
833  !
834  ! -- Initialize
835  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
836  hyn = this%hy_eff(n, m, ihc, ipos=icon)
837  hym = this%hy_eff(m, n, ihc, ipos=icon)
838  !
839  ! -- Calculate conductance
840  if (ihc == c3d_vertical) then
841  condnm = vcond(this%ibound(n), this%ibound(m), &
842  this%icelltype(n), this%icelltype(m), this%inewton, &
843  this%ivarcv, this%idewatcv, &
844  this%condsat(this%dis%con%jas(icon)), hn, hm, &
845  hyn, hym, &
846  this%sat(n), this%sat(m), &
847  this%dis%top(n), this%dis%top(m), &
848  this%dis%bot(n), this%dis%bot(m), &
849  this%dis%con%hwva(this%dis%con%jas(icon)))
850  else
851  condnm = hcond(this%ibound(n), this%ibound(m), &
852  this%icelltype(n), this%icelltype(m), &
853  this%inewton, &
854  this%dis%con%ihc(this%dis%con%jas(icon)), &
855  this%icellavg, &
856  this%condsat(this%dis%con%jas(icon)), &
857  hn, hm, this%sat(n), this%sat(m), hyn, hym, &
858  this%dis%top(n), this%dis%top(m), &
859  this%dis%bot(n), this%dis%bot(m), &
860  this%dis%con%cl1(this%dis%con%jas(icon)), &
861  this%dis%con%cl2(this%dis%con%jas(icon)), &
862  this%dis%con%hwva(this%dis%con%jas(icon)))
863  end if
864  !
865  ! -- Initialize hntemp and hmtemp
866  hntemp = hn
867  hmtemp = hm
868  !
869  ! -- Check and adjust for dewatered conditions
870  if (this%iperched /= 0) then
871  if (this%dis%con%ihc(this%dis%con%jas(icon)) == 0) then
872  if (n > m) then
873  if (this%icelltype(n) /= 0) then
874  if (hn < this%dis%top(n)) hntemp = this%dis%bot(m)
875  end if
876  else
877  if (this%icelltype(m) /= 0) then
878  if (hm < this%dis%top(m)) hmtemp = this%dis%bot(n)
879  end if
880  end if
881  end if
882  end if
883  !
884  ! -- Calculate flow positive into cell n
885  qnm = condnm * (hmtemp - hntemp)
886  !
887  ! -- Return
888  return
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 793 of file gwf-npf.f90.

794  ! -- dummy
795  class(GwfNpfType) :: this
796  integer(I4B), intent(in) :: n
797  real(DP), intent(in) :: hn
798  real(DP), intent(inout) :: thksat
799  !
800  ! -- Standard Formulation
801  if (hn >= this%dis%top(n)) then
802  thksat = done
803  else
804  thksat = (hn - this%dis%bot(n)) / (this%dis%top(n) - this%dis%bot(n))
805  end if
806  !
807  ! -- Newton-Raphson Formulation
808  if (this%inewton /= 0) then
809  thksat = squadraticsaturation(this%dis%top(n), this%dis%bot(n), hn, &
810  this%satomega)
811  end if
812  !
813  ! -- Return
814  return
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 2293 of file gwf-npf.f90.

2295  ! -- modules
2296  use tdismodule, only: kstp, kper
2297  ! -- dummy
2298  class(GwfNpfType) :: this
2299  integer(I4B), intent(in) :: icode
2300  integer(I4B), intent(inout) :: ncnvrt
2301  character(len=30), dimension(5), intent(inout) :: nodcnvrt
2302  character(len=3), dimension(5), intent(inout) :: acnvrt
2303  integer(I4B), intent(inout) :: ihdcnv
2304  integer(I4B), intent(in) :: kiter
2305  integer(I4B), intent(in) :: n
2306  ! -- local
2307  integer(I4B) :: l
2308  ! -- formats
2309  character(len=*), parameter :: fmtcnvtn = &
2310  "(1X,/1X,'CELL CONVERSIONS FOR ITER.=',I0, &
2311  &' STEP=',I0,' PERIOD=',I0,' (NODE or LRC)')"
2312  character(len=*), parameter :: fmtnode = "(1X,3X,5(A4, A20))"
2313  !
2314  ! -- Keep track of cell conversions
2315  if (icode > 0) then
2316  ncnvrt = ncnvrt + 1
2317  call this%dis%noder_to_string(n, nodcnvrt(ncnvrt))
2318  if (icode == 1) then
2319  acnvrt(ncnvrt) = 'DRY'
2320  else
2321  acnvrt(ncnvrt) = 'WET'
2322  end if
2323  end if
2324  !
2325  ! -- Print a line if 5 conversions have occurred or if icode indicates that a
2326  ! partial line should be printed
2327  if (ncnvrt == 5 .or. (icode == 0 .and. ncnvrt > 0)) then
2328  if (ihdcnv == 0) write (this%iout, fmtcnvtn) kiter, kstp, kper
2329  ihdcnv = 1
2330  write (this%iout, fmtnode) &
2331  (acnvrt(l), trim(adjustl(nodcnvrt(l))), l=1, ncnvrt)
2332  ncnvrt = 0
2333  end if
2334  !
2335  ! -- Return
2336  return

◆ 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 2121 of file gwf-npf.f90.

2122  ! -- modules
2123  use tdismodule, only: kstp, kper
2125  use constantsmodule, only: linelength
2126  ! -- dummy
2127  class(GwfNpfType) :: this
2128  integer(I4B), intent(in) :: kiter
2129  real(DP), intent(inout), dimension(:) :: hnew
2130  ! -- local
2131  integer(I4B) :: n, m, ii, ihc
2132  real(DP) :: ttop, bbot, thick
2133  integer(I4B) :: ncnvrt, ihdcnv
2134  character(len=30), dimension(5) :: nodcnvrt
2135  character(len=30) :: nodestr
2136  character(len=3), dimension(5) :: acnvrt
2137  character(len=LINELENGTH) :: errmsg
2138  integer(I4B) :: irewet
2139  ! -- formats
2140  character(len=*), parameter :: fmtnct = &
2141  "(1X,/1X,'Negative cell thickness at (layer,row,col)', &
2142  &I4,',',I5,',',I5)"
2143  character(len=*), parameter :: fmttopbot = &
2144  &"(1X,'Top elevation, bottom elevation:',1P,2G13.5)"
2145  character(len=*), parameter :: fmttopbotthk = &
2146  &"(1X,'Top elevation, bottom elevation, thickness:',1P,3G13.5)"
2147  character(len=*), parameter :: fmtdrychd = &
2148  &"(1X,/1X,'CONSTANT-HEAD CELL WENT DRY -- SIMULATION ABORTED')"
2149  character(len=*), parameter :: fmtni = &
2150  &"(1X,'CELLID=',a,' ITERATION=',I0,' TIME STEP=',I0,' STRESS PERIOD=',I0)"
2151  !
2152  ! -- Initialize
2153  ncnvrt = 0
2154  ihdcnv = 0
2155  !
2156  ! -- Convert dry cells to wet
2157  do n = 1, this%dis%nodes
2158  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
2159  m = this%dis%con%ja(ii)
2160  ihc = this%dis%con%ihc(this%dis%con%jas(ii))
2161  call this%rewet_check(kiter, n, hnew(m), this%ibound(m), ihc, hnew, &
2162  irewet)
2163  if (irewet == 1) then
2164  call this%wdmsg(2, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2165  end if
2166  end do
2167  end do
2168  !
2169  ! -- Perform drying
2170  do n = 1, this%dis%nodes
2171  !
2172  ! -- cycle if inactive or confined
2173  if (this%ibound(n) == 0) cycle
2174  if (this%icelltype(n) == 0) cycle
2175  !
2176  ! -- check for negative cell thickness
2177  bbot = this%dis%bot(n)
2178  ttop = this%dis%top(n)
2179  if (bbot > ttop) then
2180  write (errmsg, fmtnct) n
2181  call store_error(errmsg)
2182  write (errmsg, fmttopbot) ttop, bbot
2183  call store_error(errmsg)
2184  call store_error_filename(this%input_fname)
2185  end if
2186  !
2187  ! -- Calculate saturated thickness
2188  if (this%icelltype(n) /= 0) then
2189  if (hnew(n) < ttop) ttop = hnew(n)
2190  end if
2191  thick = ttop - bbot
2192  !
2193  ! -- If thick<0 print message, set hnew, and ibound
2194  if (thick <= dzero) then
2195  call this%wdmsg(1, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2196  hnew(n) = this%hdry
2197  if (this%ibound(n) < 0) then
2198  errmsg = 'CONSTANT-HEAD CELL WENT DRY -- SIMULATION ABORTED'
2199  call store_error(errmsg)
2200  write (errmsg, fmttopbotthk) ttop, bbot, thick
2201  call store_error(errmsg)
2202  call this%dis%noder_to_string(n, nodestr)
2203  write (errmsg, fmtni) trim(adjustl(nodestr)), kiter, kstp, kper
2204  call store_error(errmsg)
2205  call store_error_filename(this%input_fname)
2206  end if
2207  this%ibound(n) = 0
2208  end if
2209  end do
2210  !
2211  ! -- Print remaining cell conversions
2212  call this%wdmsg(0, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2213  !
2214  ! -- Change ibound from 30000 to 1
2215  do n = 1, this%dis%nodes
2216  if (this%ibound(n) == 30000) this%ibound(n) = 1
2217  end do
2218  !
2219  ! -- Return
2220  return
Here is the call graph for this function:

◆ source_griddata()

subroutine gwfnpfmodule::source_griddata ( class(gwfnpftype this)

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

1545  ! -- modules
1546  use simmodule, only: count_errors, store_error
1550  ! -- dummy
1551  class(GwfNpftype) :: this
1552  ! -- locals
1553  character(len=LINELENGTH) :: errmsg
1554  type(GwfNpfParamFoundType) :: found
1555  logical, dimension(2) :: afound
1556  integer(I4B), dimension(:), pointer, contiguous :: map
1557  !
1558  ! -- set map to convert user input data into reduced data
1559  map => null()
1560  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1561  !
1562  ! -- update defaults with idm sourced values
1563  call mem_set_value(this%icelltype, 'ICELLTYPE', this%input_mempath, map, &
1564  found%icelltype)
1565  call mem_set_value(this%k11, 'K', this%input_mempath, map, found%k)
1566  call mem_set_value(this%k33, 'K33', this%input_mempath, map, found%k33)
1567  call mem_set_value(this%k22, 'K22', this%input_mempath, map, found%k22)
1568  call mem_set_value(this%wetdry, 'WETDRY', this%input_mempath, map, &
1569  found%wetdry)
1570  call mem_set_value(this%angle1, 'ANGLE1', this%input_mempath, map, &
1571  found%angle1)
1572  call mem_set_value(this%angle2, 'ANGLE2', this%input_mempath, map, &
1573  found%angle2)
1574  call mem_set_value(this%angle3, 'ANGLE3', this%input_mempath, map, &
1575  found%angle3)
1576  !
1577  ! -- ensure ICELLTYPE was found
1578  if (.not. found%icelltype) then
1579  write (errmsg, '(a)') 'Error in GRIDDATA block: ICELLTYPE not found.'
1580  call store_error(errmsg)
1581  end if
1582  !
1583  ! -- ensure K was found
1584  if (.not. found%k) then
1585  write (errmsg, '(a)') 'Error in GRIDDATA block: K not found.'
1586  call store_error(errmsg)
1587  end if
1588  !
1589  ! -- set error if ik33overk set with no k33
1590  if (.not. found%k33 .and. this%ik33overk /= 0) then
1591  write (errmsg, '(a)') 'K33OVERK option specified but K33 not specified.'
1592  call store_error(errmsg)
1593  end if
1594  !
1595  ! -- set error if ik22overk set with no k22
1596  if (.not. found%k22 .and. this%ik22overk /= 0) then
1597  write (errmsg, '(a)') 'K22OVERK option specified but K22 not specified.'
1598  call store_error(errmsg)
1599  end if
1600  !
1601  ! -- handle found side effects
1602  if (found%k33) this%ik33 = 1
1603  if (found%k22) this%ik22 = 1
1604  if (found%wetdry) this%iwetdry = 1
1605  if (found%angle1) this%iangle1 = 1
1606  if (found%angle2) this%iangle2 = 1
1607  if (found%angle3) this%iangle3 = 1
1608  !
1609  ! -- handle not found side effects
1610  if (.not. found%k33) then
1611  call mem_set_value(this%k33, 'K', this%input_mempath, map, afound(1))
1612  end if
1613  if (.not. found%k22) then
1614  call mem_set_value(this%k22, 'K', this%input_mempath, map, afound(2))
1615  end if
1616  if (.not. found%wetdry) call mem_reallocate(this%wetdry, 1, 'WETDRY', &
1617  trim(this%memoryPath))
1618  if (.not. found%angle1 .and. this%ixt3d == 0) &
1619  call mem_reallocate(this%angle1, 0, 'ANGLE1', trim(this%memoryPath))
1620  if (.not. found%angle2 .and. this%ixt3d == 0) &
1621  call mem_reallocate(this%angle2, 0, 'ANGLE2', trim(this%memoryPath))
1622  if (.not. found%angle3 .and. this%ixt3d == 0) &
1623  call mem_reallocate(this%angle3, 0, 'ANGLE3', trim(this%memoryPath))
1624  !
1625  ! -- log griddata
1626  if (this%iout > 0) then
1627  call this%log_griddata(found)
1628  end if
1629  !
1630  ! -- Return
1631  return
Here is the call graph for this function:

◆ source_options()

subroutine gwfnpfmodule::source_options ( class(gwfnpftype this)

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

1327  ! -- modules
1333  use sourcecommonmodule, only: filein_fname
1334  ! -- dummy
1335  class(GwfNpftype) :: this
1336  ! -- locals
1337  character(len=LENVARNAME), dimension(3) :: cellavg_method = &
1338  &[character(len=LENVARNAME) :: 'LOGARITHMIC', 'AMT-LMK', 'AMT-HMK']
1339  type(gwfnpfparamfoundtype) :: found
1340  character(len=LINELENGTH) :: tvk6_filename
1341  !
1342  ! -- update defaults with idm sourced values
1343  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
1344  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
1345  call mem_set_value(this%icellavg, 'CELLAVG', this%input_mempath, &
1346  cellavg_method, found%cellavg)
1347  call mem_set_value(this%ithickstrt, 'ITHICKSTRT', this%input_mempath, &
1348  found%ithickstrt)
1349  call mem_set_value(this%iperched, 'IPERCHED', this%input_mempath, &
1350  found%iperched)
1351  call mem_set_value(this%ivarcv, 'IVARCV', this%input_mempath, found%ivarcv)
1352  call mem_set_value(this%idewatcv, 'IDEWATCV', this%input_mempath, &
1353  found%idewatcv)
1354  call mem_set_value(this%ixt3d, 'IXT3D', this%input_mempath, found%ixt3d)
1355  call mem_set_value(this%ixt3drhs, 'IXT3DRHS', this%input_mempath, &
1356  found%ixt3drhs)
1357  call mem_set_value(this%isavspdis, 'ISAVSPDIS', this%input_mempath, &
1358  found%isavspdis)
1359  call mem_set_value(this%isavsat, 'ISAVSAT', this%input_mempath, found%isavsat)
1360  call mem_set_value(this%ik22overk, 'IK22OVERK', this%input_mempath, &
1361  found%ik22overk)
1362  call mem_set_value(this%ik33overk, 'IK33OVERK', this%input_mempath, &
1363  found%ik33overk)
1364  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
1365  call mem_set_value(this%satomega, 'SATOMEGA', this%input_mempath, &
1366  found%satomega)
1367  call mem_set_value(this%irewet, 'IREWET', this%input_mempath, found%irewet)
1368  call mem_set_value(this%wetfct, 'WETFCT', this%input_mempath, found%wetfct)
1369  call mem_set_value(this%iwetit, 'IWETIT', this%input_mempath, found%iwetit)
1370  call mem_set_value(this%ihdwet, 'IHDWET', this%input_mempath, found%ihdwet)
1371  !
1372  ! -- save flows option active
1373  if (found%ipakcb) this%ipakcb = -1
1374  !
1375  ! -- xt3d active with rhs
1376  if (found%ixt3d .and. found%ixt3drhs) this%ixt3d = 2
1377  !
1378  ! -- save specific discharge active
1379  if (found%isavspdis) this%icalcspdis = this%isavspdis
1380  !
1381  ! -- no newton specified
1382  if (found%inewton) then
1383  this%inewton = 0
1384  this%iasym = 0
1385  end if
1386  !
1387  ! -- enforce 0 or 1 TVK6_FILENAME entries in option block
1388  if (filein_fname(tvk6_filename, 'TVK6_FILENAME', this%input_mempath, &
1389  this%input_fname)) then
1390  call openfile(this%intvk, this%iout, tvk6_filename, 'TVK')
1391  call tvk_cr(this%tvk, this%name_model, this%intvk, this%iout)
1392  end if
1393  !
1394  ! -- log options
1395  if (this%iout > 0) then
1396  call this%log_options(found)
1397  end if
1398  !
1399  ! -- Return
1400  return
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 1169 of file gwf-npf.f90.

1170  ! -- modules
1172  ! -- dummy
1173  class(GwfNpftype) :: this
1174  integer(I4B), intent(in) :: ncells
1175  integer(I4B), intent(in) :: njas
1176  ! -- local
1177  integer(I4B) :: n
1178  !
1179  ! -- Retain copy of user-specified K arrays
1180  do n = 1, ncells
1181  this%k11input(n) = this%k11(n)
1182  this%k22input(n) = this%k22(n)
1183  this%k33input(n) = this%k33(n)
1184  end do
1185  !
1186  ! -- Return
1187  return