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

Data Types

type  gwfnpftype
 

Functions/Subroutines

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

Function/Subroutine Documentation

◆ allocate_arrays()

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

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

1217  ! -- dummy
1218  class(GwfNpftype) :: this
1219  integer(I4B), intent(in) :: ncells
1220  integer(I4B), intent(in) :: njas
1221  ! -- local
1222  integer(I4B) :: n
1223  !
1224  call mem_allocate(this%ithickstartflag, ncells, 'ITHICKSTARTFLAG', &
1225  this%memoryPath)
1226  call mem_allocate(this%icelltype, ncells, 'ICELLTYPE', this%memoryPath)
1227  call mem_allocate(this%k11, ncells, 'K11', this%memoryPath)
1228  call mem_allocate(this%sat, ncells, 'SAT', this%memoryPath)
1229  call mem_allocate(this%condsat, njas, 'CONDSAT', this%memoryPath)
1230  !
1231  ! -- Optional arrays dimensioned to full size initially
1232  call mem_allocate(this%k22, ncells, 'K22', this%memoryPath)
1233  call mem_allocate(this%k33, ncells, 'K33', this%memoryPath)
1234  call mem_allocate(this%wetdry, ncells, 'WETDRY', this%memoryPath)
1235  call mem_allocate(this%angle1, ncells, 'ANGLE1', this%memoryPath)
1236  call mem_allocate(this%angle2, ncells, 'ANGLE2', this%memoryPath)
1237  call mem_allocate(this%angle3, ncells, 'ANGLE3', this%memoryPath)
1238  !
1239  ! -- Optional arrays
1240  call mem_allocate(this%ibotnode, 0, 'IBOTNODE', this%memoryPath)
1241  call mem_allocate(this%nodedge, 0, 'NODEDGE', this%memoryPath)
1242  call mem_allocate(this%ihcedge, 0, 'IHCEDGE', this%memoryPath)
1243  call mem_allocate(this%propsedge, 0, 0, 'PROPSEDGE', this%memoryPath)
1244  call mem_allocate(this%iedge_ptr, 0, 'NREDGESNODE', this%memoryPath)
1245  call mem_allocate(this%edge_idxs, 0, 'EDGEIDXS', this%memoryPath)
1246  !
1247  ! -- Optional arrays only needed when vsc package is active
1248  call mem_allocate(this%k11input, 0, 'K11INPUT', this%memoryPath)
1249  call mem_allocate(this%k22input, 0, 'K22INPUT', this%memoryPath)
1250  call mem_allocate(this%k33input, 0, 'K33INPUT', this%memoryPath)
1251  !
1252  ! -- Specific discharge is (re-)allocated when nedges is known
1253  call mem_allocate(this%spdis, 3, 0, 'SPDIS', this%memoryPath)
1254  !
1255  ! -- Time-varying property flag arrays
1256  call mem_allocate(this%nodekchange, ncells, 'NODEKCHANGE', this%memoryPath)
1257  !
1258  ! -- initialize iangle1, iangle2, iangle3, and wetdry
1259  do n = 1, ncells
1260  this%angle1(n) = dzero
1261  this%angle2(n) = dzero
1262  this%angle3(n) = dzero
1263  this%wetdry(n) = dzero
1264  this%nodekchange(n) = dzero
1265  end do
1266  !
1267  ! -- allocate variable names
1268  allocate (this%aname(this%iname))
1269  this%aname = [' ICELLTYPE', ' K', &
1270  ' K33', ' K22', &
1271  ' WETDRY', ' ANGLE1', &
1272  ' ANGLE2', ' ANGLE3']

◆ allocate_scalars()

subroutine gwfnpfmodule::allocate_scalars ( class(gwfnpftype this)

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

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

1098  ! -- modules
1100  ! -- dummy
1101  class(GwfNpftype) :: this
1102  !
1103  ! -- allocate scalars in NumericalPackageType
1104  call this%NumericalPackageType%allocate_scalars()
1105  !
1106  ! -- Allocate scalars
1107  call mem_allocate(this%iname, 'INAME', this%memoryPath)
1108  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
1109  call mem_allocate(this%ixt3drhs, 'IXT3DRHS', this%memoryPath)
1110  call mem_allocate(this%satomega, 'SATOMEGA', this%memoryPath)
1111  call mem_allocate(this%hnoflo, 'HNOFLO', this%memoryPath)
1112  call mem_allocate(this%hdry, 'HDRY', this%memoryPath)
1113  call mem_allocate(this%icellavg, 'ICELLAVG', this%memoryPath)
1114  call mem_allocate(this%iavgkeff, 'IAVGKEFF', this%memoryPath)
1115  call mem_allocate(this%ik22, 'IK22', this%memoryPath)
1116  call mem_allocate(this%ik33, 'IK33', this%memoryPath)
1117  call mem_allocate(this%ik22overk, 'IK22OVERK', this%memoryPath)
1118  call mem_allocate(this%ik33overk, 'IK33OVERK', this%memoryPath)
1119  call mem_allocate(this%iperched, 'IPERCHED', this%memoryPath)
1120  call mem_allocate(this%ivarcv, 'IVARCV', this%memoryPath)
1121  call mem_allocate(this%idewatcv, 'IDEWATCV', this%memoryPath)
1122  call mem_allocate(this%ithickstrt, 'ITHICKSTRT', this%memoryPath)
1123  call mem_allocate(this%ihighcellsat, 'IHIGHCELLSAT', this%memoryPath)
1124  call mem_allocate(this%icalcspdis, 'ICALCSPDIS', this%memoryPath)
1125  call mem_allocate(this%isavspdis, 'ISAVSPDIS', this%memoryPath)
1126  call mem_allocate(this%isavsat, 'ISAVSAT', this%memoryPath)
1127  call mem_allocate(this%irewet, 'IREWET', this%memoryPath)
1128  call mem_allocate(this%wetfct, 'WETFCT', this%memoryPath)
1129  call mem_allocate(this%iwetit, 'IWETIT', this%memoryPath)
1130  call mem_allocate(this%ihdwet, 'IHDWET', this%memoryPath)
1131  call mem_allocate(this%iangle1, 'IANGLE1', this%memoryPath)
1132  call mem_allocate(this%iangle2, 'IANGLE2', this%memoryPath)
1133  call mem_allocate(this%iangle3, 'IANGLE3', this%memoryPath)
1134  call mem_allocate(this%iwetdry, 'IWETDRY', this%memoryPath)
1135  call mem_allocate(this%nedges, 'NEDGES', this%memoryPath)
1136  call mem_allocate(this%lastedge, 'LASTEDGE', this%memoryPath)
1137  call mem_allocate(this%intvk, 'INTVK', this%memoryPath)
1138  call mem_allocate(this%invsc, 'INVSC', this%memoryPath)
1139  call mem_allocate(this%kchangeper, 'KCHANGEPER', this%memoryPath)
1140  call mem_allocate(this%kchangestp, 'KCHANGESTP', this%memoryPath)
1141  !
1142  ! -- set pointer to inewtonur
1143  call mem_setptr(this%igwfnewtonur, 'INEWTONUR', &
1144  create_mem_path(this%name_model))
1145  !
1146  ! -- Initialize value
1147  this%iname = 8
1148  this%ixt3d = 0
1149  this%ixt3drhs = 0
1150  this%satomega = dzero
1151  this%hnoflo = dhnoflo !1.d30
1152  this%hdry = dhdry !-1.d30
1153  this%icellavg = ccond_hmean
1154  this%iavgkeff = 0
1155  this%ik22 = 0
1156  this%ik33 = 0
1157  this%ik22overk = 0
1158  this%ik33overk = 0
1159  this%iperched = 0
1160  this%ivarcv = 0
1161  this%idewatcv = 0
1162  this%ithickstrt = 0
1163  this%ihighcellsat = 0
1164  this%icalcspdis = 0
1165  this%isavspdis = 0
1166  this%isavsat = 0
1167  this%irewet = 0
1168  this%wetfct = done
1169  this%iwetit = 1
1170  this%ihdwet = 0
1171  this%iangle1 = 0
1172  this%iangle2 = 0
1173  this%iangle3 = 0
1174  this%iwetdry = 0
1175  this%nedges = 0
1176  this%lastedge = 0
1177  this%intvk = 0
1178  this%invsc = 0
1179  this%kchangeper = 0
1180  this%kchangestp = 0
1181  !
1182  ! -- If newton is on, then NPF creates asymmetric matrix
1183  this%iasym = this%inewton
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Here is the call graph for this function:

◆ calc_condsat()

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

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

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

2030  ! -- dummy variables
2031  class(GwfNpfType) :: this
2032  integer(I4B), intent(in) :: node
2033  logical, intent(in) :: upperOnly
2034  ! -- local variables
2035  integer(I4B) :: ii, m, n, ihc, jj
2036  real(DP) :: topm, topn, topnode, botm, botn, botnode, satm, satn, satnode
2037  real(DP) :: hyn, hym, hn, hm, fawidth, csat
2038  !
2039  satnode = this%calc_initial_sat(node)
2040  !
2041  topnode = this%dis%top(node)
2042  botnode = this%dis%bot(node)
2043  !
2044  ! -- Go through the connecting cells
2045  do ii = this%dis%con%ia(node) + 1, this%dis%con%ia(node + 1) - 1
2046  !
2047  ! -- Set the m cell number and cycle if lower triangle connection and
2048  ! -- we're not updating both upper and lower matrix parts for this node
2049  m = this%dis%con%ja(ii)
2050  jj = this%dis%con%jas(ii)
2051  if (m < node) then
2052  if (upperonly) cycle
2053  ! m => node, n => neighbour
2054  n = m
2055  m = node
2056  topm = topnode
2057  botm = botnode
2058  satm = satnode
2059  topn = this%dis%top(n)
2060  botn = this%dis%bot(n)
2061  satn = this%calc_initial_sat(n)
2062  else
2063  ! n => node, m => neighbour
2064  n = node
2065  topn = topnode
2066  botn = botnode
2067  satn = satnode
2068  topm = this%dis%top(m)
2069  botm = this%dis%bot(m)
2070  satm = this%calc_initial_sat(m)
2071  end if
2072  !
2073  ihc = this%dis%con%ihc(jj)
2074  hyn = this%hy_eff(n, m, ihc, ipos=ii)
2075  hym = this%hy_eff(m, n, ihc, ipos=ii)
2076  if (this%ithickstartflag(n) == 0) then
2077  hn = topn
2078  else
2079  hn = this%ic%strt(n)
2080  end if
2081  if (this%ithickstartflag(m) == 0) then
2082  hm = topm
2083  else
2084  hm = this%ic%strt(m)
2085  end if
2086  !
2087  ! -- Calculate conductance depending on whether connection is
2088  ! vertical (0), horizontal (1), or staggered horizontal (2)
2089  if (ihc == c3d_vertical) then
2090  !
2091  ! -- Vertical conductance for fully saturated conditions
2092  csat = vcond(1, 1, 1, 1, 0, 1, 1, done, &
2093  botn, botm, &
2094  hyn, hym, &
2095  satn, satm, &
2096  topn, topm, &
2097  botn, botm, &
2098  this%dis%con%hwva(jj))
2099  else
2100  !
2101  ! -- Horizontal conductance for fully saturated conditions
2102  fawidth = this%dis%con%hwva(jj)
2103  csat = hcond(1, 1, 1, 1, 0, &
2104  ihc, &
2105  this%icellavg, &
2106  done, &
2107  hn, hm, satn, satm, hyn, hym, &
2108  topn, topm, &
2109  botn, botm, &
2110  this%dis%con%cl1(jj), &
2111  this%dis%con%cl2(jj), &
2112  fawidth)
2113  end if
2114  this%condsat(jj) = csat
2115  end do
Here is the call graph for this function:

◆ calc_initial_sat()

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

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

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

2126  ! -- dummy variables
2127  class(GwfNpfType) :: this
2128  integer(I4B), intent(in) :: n
2129  ! -- Return
2130  real(DP) :: satn
2131  !
2132  satn = done
2133  if (this%ibound(n) /= 0 .and. this%ithickstartflag(n) /= 0) then
2134  call this%thksat(n, this%ic%strt(n), satn)
2135  end if

◆ calc_max_conns()

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

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

2760  class(GwfNpfType) :: this
2761  integer(I4B) :: max_conns
2762  ! local
2763  integer(I4B) :: n, m, ic
2764 
2765  max_conns = 0
2766  do n = 1, this%dis%nodes
2767 
2768  ! Count internal model connections
2769  ic = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
2770 
2771  ! Add edge connections
2772  do m = 1, this%nedges
2773  if (this%nodedge(m) == n) then
2774  ic = ic + 1
2775  end if
2776  end do
2777 
2778  ! Set max number of connections for any cell
2779  if (ic > max_conns) max_conns = ic
2780  end do
2781 

◆ calc_spdis()

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

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

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

2861  ! -- dummy
2862  class(GwfNpfType) :: this !< this NPF instance
2863  integer(I4B) :: n !< node n
2864  integer(I4B) :: m !< node m
2865  integer(I4B) :: ihc !< 1 = horizontal connection, 0 for vertical
2866  ! -- return
2867  real(DP) :: satThickness !< saturated thickness
2868  !
2869  satthickness = thksatnm(this%ibound(n), &
2870  this%ibound(m), &
2871  this%icelltype(n), &
2872  this%icelltype(m), &
2873  this%inewton, &
2874  ihc, &
2875  this%hnew(n), &
2876  this%hnew(m), &
2877  this%sat(n), &
2878  this%sat(m), &
2879  this%dis%top(n), &
2880  this%dis%top(m), &
2881  this%dis%bot(n), &
2882  this%dis%bot(m))

◆ check_options()

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

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

1457  ! -- modules
1458  use simmodule, only: store_error, store_warning, &
1460  use constantsmodule, only: linelength
1461  ! -- dummy
1462  class(GwfNpftype) :: this
1463  !
1464  ! -- set omega value used for saturation calculations
1465  if (this%inewton > 0) then
1466  this%satomega = dem6
1467  end if
1468  !
1469  if (this%inewton > 0) then
1470  if (this%iperched > 0) then
1471  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1472  'BE USED WITH PERCHED OPTION.'
1473  call store_error(errmsg)
1474  end if
1475  if (this%ivarcv > 0) then
1476  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1477  'BE USED WITH VARIABLECV OPTION.'
1478  call store_error(errmsg)
1479  end if
1480  if (this%irewet > 0) then
1481  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1482  'BE USED WITH REWET OPTION.'
1483  call store_error(errmsg)
1484  end if
1485  else
1486  if (this%ihighcellsat /= 0) then
1487  write (warnmsg, '(a)') 'HIGHEST_CELL_SATURATION '// &
1488  'option cannot be used when NEWTON option in not specified. '// &
1489  'Resetting HIGHEST_CELL_SATURATION option to off.'
1490  this%ihighcellsat = 0
1491  call store_warning(warnmsg)
1492  end if
1493  end if
1494  !
1495  if (this%ixt3d /= 0) then
1496  if (this%icellavg > 0) then
1497  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. '// &
1498  'ALTERNATIVE_CELL_AVERAGING OPTION '// &
1499  'CANNOT BE USED WITH XT3D OPTION.'
1500  call store_error(errmsg)
1501  end if
1502  if (this%ithickstrt > 0) then
1503  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. THICKSTRT OPTION '// &
1504  'CANNOT BE USED WITH XT3D OPTION.'
1505  call store_error(errmsg)
1506  end if
1507  if (this%iperched > 0) then
1508  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. PERCHED OPTION '// &
1509  'CANNOT BE USED WITH XT3D OPTION.'
1510  call store_error(errmsg)
1511  end if
1512  if (this%ivarcv > 0) then
1513  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. VARIABLECV OPTION '// &
1514  'CANNOT BE USED WITH XT3D OPTION.'
1515  call store_error(errmsg)
1516  end if
1517  end if
1518  !
1519  ! -- Terminate if errors
1520  if (count_errors() > 0) then
1521  call store_error_filename(this%input_fname)
1522  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
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:

◆ highest_cell_saturation()

subroutine gwfnpfmodule::highest_cell_saturation ( class(gwfnpftype this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
real(dp), intent(in)  hn,
real(dp), intent(in)  hm,
real(dp), intent(inout)  satn,
real(dp), intent(inout)  satm 
)

Calculate the saturation based on the maximum cell bottom for two connected cells

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

593  ! dummy
594  class(GwfNpfType) :: this
595  integer(I4B), intent(in) :: n, m
596  real(DP), intent(in) :: hn, hm
597  real(DP), intent(inout) :: satn, satm
598  ! local
599  integer(I4B) :: ihdbot
600  real(DP) :: botn, botm
601  real(DP) :: top, bot
602 
603  botn = this%dis%bot(n)
604  botm = this%dis%bot(m)
605 
606  ihdbot = n
607  if (botm > botn) ihdbot = m
608 
609  ! recalculate saturation if the difference in elevation between
610  ! two cells exceed a threshold value
611  if (abs(botm - botn) >= dem2) then
612  top = this%dis%top(ihdbot)
613  bot = this%dis%bot(ihdbot)
614  satn = squadraticsaturation(top, bot, hn, this%satomega)
615  satm = squadraticsaturation(top, bot, hm, this%satomega)
616  end if
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 2359 of file gwf-npf.f90.

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

◆ increase_edge_count()

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

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

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

2750  ! -- dummy
2751  class(GwfNpfType) :: this
2752  integer(I4B), intent(in) :: nedges
2753  !
2754  this%nedges = this%nedges + nedges

◆ log_griddata()

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

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

1528  ! -- modules
1530  ! -- dummy
1531  class(GwfNpfType) :: this
1532  type(GwfNpfParamFoundType), intent(in) :: found
1533  !
1534  write (this%iout, '(1x,a)') 'Setting NPF Griddata'
1535  !
1536  if (found%icelltype) then
1537  write (this%iout, '(4x,a)') 'ICELLTYPE set from input file'
1538  end if
1539  !
1540  if (found%k) then
1541  write (this%iout, '(4x,a)') 'K set from input file'
1542  end if
1543  !
1544  if (found%k33) then
1545  write (this%iout, '(4x,a)') 'K33 set from input file'
1546  else
1547  write (this%iout, '(4x,a)') 'K33 not provided. Setting K33 = K.'
1548  end if
1549  !
1550  if (found%k22) then
1551  write (this%iout, '(4x,a)') 'K22 set from input file'
1552  else
1553  write (this%iout, '(4x,a)') 'K22 not provided. Setting K22 = K.'
1554  end if
1555  !
1556  if (found%wetdry) then
1557  write (this%iout, '(4x,a)') 'WETDRY set from input file'
1558  end if
1559  !
1560  if (found%angle1) then
1561  write (this%iout, '(4x,a)') 'ANGLE1 set from input file'
1562  end if
1563  !
1564  if (found%angle2) then
1565  write (this%iout, '(4x,a)') 'ANGLE2 set from input file'
1566  end if
1567  !
1568  if (found%angle3) then
1569  write (this%iout, '(4x,a)') 'ANGLE3 set from input file'
1570  end if
1571  !
1572  write (this%iout, '(1x,a,/)') 'End Setting NPF Griddata'

◆ log_options()

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

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

1278  ! -- modules
1279  use kindmodule, only: lgp
1281  ! -- dummy
1282  class(GwfNpftype) :: this
1283  ! -- locals
1284  type(GwfNpfParamFoundType), intent(in) :: found
1285  !
1286  write (this%iout, '(1x,a)') 'Setting NPF Options'
1287  if (found%iprflow) &
1288  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be printed &
1289  &to listing file whenever ICBCFL is not zero.'
1290  if (found%ipakcb) &
1291  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be saved &
1292  &to binary file whenever ICBCFL is not zero.'
1293  if (found%cellavg) &
1294  write (this%iout, '(4x,a,i0)') 'Alternative cell averaging [1=logarithmic, &
1295  &2=AMT-LMK, 3=AMT-HMK] set to: ', &
1296  this%icellavg
1297  if (found%ithickstrt) &
1298  write (this%iout, '(4x,a)') 'THICKSTRT option has been activated.'
1299  if (found%ihighcellsat) &
1300  write (this%iout, '(4x,a)') 'HIGHEST_CELL_SATURATION option &
1301  &has been activated.'
1302  if (found%iperched) &
1303  write (this%iout, '(4x,a)') 'Vertical flow will be adjusted for perched &
1304  &conditions.'
1305  if (found%ivarcv) &
1306  write (this%iout, '(4x,a)') 'Vertical conductance varies with water table.'
1307  if (found%idewatcv) &
1308  write (this%iout, '(4x,a)') 'Vertical conductance is calculated using &
1309  &only the saturated thickness and properties &
1310  &of the overlying cell if the head in the &
1311  &underlying cell is below its top.'
1312  if (found%ixt3d) write (this%iout, '(4x,a)') 'XT3D formulation is selected.'
1313  if (found%ixt3drhs) &
1314  write (this%iout, '(4x,a)') 'XT3D RHS formulation is selected.'
1315  if (found%isavspdis) &
1316  write (this%iout, '(4x,a)') 'Specific discharge will be calculated at cell &
1317  &centers and written to DATA-SPDIS in budget &
1318  &file when requested.'
1319  if (found%isavsat) &
1320  write (this%iout, '(4x,a)') 'Saturation will be written to DATA-SAT in &
1321  &budget file when requested.'
1322  if (found%ik22overk) &
1323  write (this%iout, '(4x,a)') 'Values specified for K22 are anisotropy &
1324  &ratios and will be multiplied by K before &
1325  &being used in calculations.'
1326  if (found%ik33overk) &
1327  write (this%iout, '(4x,a)') 'Values specified for K33 are anisotropy &
1328  &ratios and will be multiplied by K before &
1329  &being used in calculations.'
1330  if (found%inewton) &
1331  write (this%iout, '(4x,a)') 'NEWTON-RAPHSON method disabled for unconfined &
1332  &cells'
1333  if (found%satomega) &
1334  write (this%iout, '(4x,a,1pg15.6)') 'Saturation omega: ', this%satomega
1335  if (found%irewet) &
1336  write (this%iout, '(4x,a)') 'Rewetting is active.'
1337  if (found%wetfct) &
1338  write (this%iout, '(4x,a,1pg15.6)') &
1339  'Wetting factor (WETFCT) has been set to: ', this%wetfct
1340  if (found%iwetit) &
1341  write (this%iout, '(4x,a,i5)') &
1342  'Wetting iteration interval (IWETIT) has been set to: ', this%iwetit
1343  if (found%ihdwet) &
1344  write (this%iout, '(4x,a,i5)') &
1345  'Head rewet equation (IHDWET) has been set to: ', this%ihdwet
1346  write (this%iout, '(1x,a,/)') 'End Setting NPF Options'
This module defines variable data types.
Definition: kind.f90:8

◆ npf_ac()

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

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

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

◆ npf_ad()

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

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

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

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

◆ npf_ar()

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

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

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

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

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

◆ npf_cf()

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

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

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

◆ npf_cq()

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

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

789  ! -- dummy
790  class(GwfNpfType) :: this
791  real(DP), intent(inout), dimension(:) :: hnew
792  real(DP), intent(inout), dimension(:) :: flowja
793  ! -- local
794  integer(I4B) :: n, ipos, m
795  real(DP) :: qnm
796  !
797  ! -- Calculate the flow across each cell face and store in flowja
798  !
799  if (this%ixt3d /= 0) then
800  call this%xt3d%xt3d_flowja(hnew, flowja)
801  else
802  !
803  do n = 1, this%dis%nodes
804  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
805  m = this%dis%con%ja(ipos)
806  if (m < n) cycle
807  call this%qcalc(n, m, hnew(n), hnew(m), ipos, qnm)
808  flowja(ipos) = qnm
809  flowja(this%dis%con%isym(ipos)) = -qnm
810  end do
811  end do
812  !
813  end if

◆ npf_cr()

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

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

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

◆ npf_da()

subroutine gwfnpfmodule::npf_da ( class(gwfnpftype this)

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

996  ! -- modules
999  ! -- dummy
1000  class(GwfNpftype) :: this
1001 
1002  ! free spdis work structure
1003  if (this%icalcspdis == 1 .and. this%spdis_wa%is_created()) &
1004  call this%spdis_wa%destroy()
1005  deallocate (this%spdis_wa)
1006  !
1007  ! -- Deallocate input memory
1008  call memorystore_remove(this%name_model, 'NPF', idm_context)
1009  !
1010  ! -- TVK
1011  if (this%intvk /= 0) then
1012  call this%tvk%da()
1013  deallocate (this%tvk)
1014  end if
1015  !
1016  ! -- VSC
1017  if (this%invsc /= 0) then
1018  nullify (this%vsc)
1019  end if
1020  !
1021  ! -- Strings
1022  !
1023  ! -- Scalars
1024  call mem_deallocate(this%iname)
1025  call mem_deallocate(this%ixt3d)
1026  call mem_deallocate(this%ixt3drhs)
1027  call mem_deallocate(this%satomega)
1028  call mem_deallocate(this%hnoflo)
1029  call mem_deallocate(this%hdry)
1030  call mem_deallocate(this%icellavg)
1031  call mem_deallocate(this%iavgkeff)
1032  call mem_deallocate(this%ik22)
1033  call mem_deallocate(this%ik33)
1034  call mem_deallocate(this%iperched)
1035  call mem_deallocate(this%ivarcv)
1036  call mem_deallocate(this%idewatcv)
1037  call mem_deallocate(this%ithickstrt)
1038  call mem_deallocate(this%ihighcellsat)
1039  call mem_deallocate(this%isavspdis)
1040  call mem_deallocate(this%isavsat)
1041  call mem_deallocate(this%icalcspdis)
1042  call mem_deallocate(this%irewet)
1043  call mem_deallocate(this%wetfct)
1044  call mem_deallocate(this%iwetit)
1045  call mem_deallocate(this%ihdwet)
1046  call mem_deallocate(this%ibotnode)
1047  call mem_deallocate(this%iwetdry)
1048  call mem_deallocate(this%iangle1)
1049  call mem_deallocate(this%iangle2)
1050  call mem_deallocate(this%iangle3)
1051  call mem_deallocate(this%nedges)
1052  call mem_deallocate(this%lastedge)
1053  call mem_deallocate(this%ik22overk)
1054  call mem_deallocate(this%ik33overk)
1055  call mem_deallocate(this%intvk)
1056  call mem_deallocate(this%invsc)
1057  call mem_deallocate(this%kchangeper)
1058  call mem_deallocate(this%kchangestp)
1059  !
1060  ! -- Deallocate arrays
1061  deallocate (this%aname)
1062  call mem_deallocate(this%ithickstartflag)
1063  call mem_deallocate(this%icelltype)
1064  call mem_deallocate(this%k11)
1065  call mem_deallocate(this%k22)
1066  call mem_deallocate(this%k33)
1067  call mem_deallocate(this%k11input)
1068  call mem_deallocate(this%k22input)
1069  call mem_deallocate(this%k33input)
1070  call mem_deallocate(this%sat, 'SAT', this%memoryPath)
1071  call mem_deallocate(this%condsat)
1072  call mem_deallocate(this%wetdry)
1073  call mem_deallocate(this%angle1)
1074  call mem_deallocate(this%angle2)
1075  call mem_deallocate(this%angle3)
1076  call mem_deallocate(this%nodedge)
1077  call mem_deallocate(this%ihcedge)
1078  call mem_deallocate(this%propsedge)
1079  call mem_deallocate(this%iedge_ptr)
1080  call mem_deallocate(this%edge_idxs)
1081  call mem_deallocate(this%spdis, 'SPDIS', this%memoryPath)
1082  call mem_deallocate(this%nodekchange)
1083  !
1084  ! -- deallocate parent
1085  call this%NumericalPackageType%da()
1086 
1087  ! pointers
1088  this%hnew => null()
1089 
subroutine, public memorystore_remove(component, subcomponent, context)
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
Here is the call graph for this function:

◆ npf_df()

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

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

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

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

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

◆ npf_fc()

subroutine gwfnpfmodule::npf_fc ( class(gwfnpftype this,
integer(i4b)  kiter,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs,
real(dp), dimension(:), intent(inout)  hnew 
)
private

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

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

◆ npf_fn()

subroutine gwfnpfmodule::npf_fn ( class(gwfnpftype this,
integer(i4b)  kiter,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs,
real(dp), dimension(:), intent(inout)  hnew 
)
private

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

622  ! -- dummy
623  class(GwfNpfType) :: this
624  integer(I4B) :: kiter
625  class(MatrixBaseType), pointer :: matrix_sln
626  integer(I4B), intent(in), dimension(:) :: idxglo
627  real(DP), intent(inout), dimension(:) :: rhs
628  real(DP), intent(inout), dimension(:) :: hnew
629  ! -- local
630  integer(I4B) :: nodes, nja
631  integer(I4B) :: n, m, ii, idiag
632  integer(I4B) :: isymcon, idiagm
633  integer(I4B) :: iups
634  integer(I4B) :: idn
635  real(DP) :: cond
636  real(DP) :: consterm
637  real(DP) :: filledterm
638  real(DP) :: derv
639  real(DP) :: hds
640  real(DP) :: term
641  real(DP) :: topup
642  real(DP) :: botup
643  !
644  ! -- add newton terms to solution matrix
645  nodes = this%dis%nodes
646  nja = this%dis%con%nja
647  if (this%ixt3d /= 0) then
648  call this%xt3d%xt3d_fn(kiter, nodes, nja, matrix_sln, idxglo, rhs, hnew)
649  else
650  !
651  do n = 1, nodes
652  idiag = this%dis%con%ia(n)
653  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
654  if (this%dis%con%mask(ii) == 0) cycle
655 
656  m = this%dis%con%ja(ii)
657  isymcon = this%dis%con%isym(ii)
658  ! work on upper triangle
659  if (m < n) cycle
660  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0 .and. &
661  this%ivarcv == 0) then
662  !call this%vcond(n,m,hnew(n),hnew(m),ii,cond)
663  ! do nothing
664  else
665  ! determine upstream node
666  iups = m
667  if (hnew(m) < hnew(n)) iups = n
668  idn = n
669  if (iups == n) idn = m
670  !
671  ! -- no newton terms if upstream cell is confined
672  if (this%icelltype(iups) == 0) cycle
673  !
674  ! -- Set the upstream top and bot, and then recalculate for a
675  ! vertically staggered horizontal connection
676  topup = this%dis%top(iups)
677  botup = this%dis%bot(iups)
678  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 2) then
679  topup = min(this%dis%top(n), this%dis%top(m))
680  botup = max(this%dis%bot(n), this%dis%bot(m))
681  end if
682  !
683  ! get saturated conductivity for derivative
684  cond = this%condsat(this%dis%con%jas(ii))
685  !
686  ! compute additional term
687  consterm = -cond * (hnew(iups) - hnew(idn)) !needs to use hwadi instead of hnew(idn)
688  !filledterm = cond
689  filledterm = matrix_sln%get_value_pos(idxglo(ii))
690  derv = squadraticsaturationderivative(topup, botup, hnew(iups), &
691  this%satomega)
692  idiagm = this%dis%con%ia(m)
693  ! fill jacobian for n being the upstream node
694  if (iups == n) then
695  hds = hnew(m)
696  !isymcon = this%dis%con%isym(ii)
697  term = consterm * derv
698  rhs(n) = rhs(n) + term * hnew(n) !+ amat(idxglo(isymcon)) * (dwadi * hds - hds) !need to add dwadi
699  rhs(m) = rhs(m) - term * hnew(n) !- amat(idxglo(isymcon)) * (dwadi * hds - hds) !need to add dwadi
700  ! fill in row of n
701  call matrix_sln%add_value_pos(idxglo(idiag), term)
702  ! fill newton term in off diagonal if active cell
703  if (this%ibound(n) > 0) then
704  filledterm = matrix_sln%get_value_pos(idxglo(ii))
705  call matrix_sln%set_value_pos(idxglo(ii), filledterm) !* dwadi !need to add dwadi
706  end if
707  !fill row of m
708  filledterm = matrix_sln%get_value_pos(idxglo(idiagm))
709  call matrix_sln%set_value_pos(idxglo(idiagm), filledterm) !- filledterm * (dwadi - DONE) !need to add dwadi
710  ! fill newton term in off diagonal if active cell
711  if (this%ibound(m) > 0) then
712  call matrix_sln%add_value_pos(idxglo(isymcon), -term)
713  end if
714  ! fill jacobian for m being the upstream node
715  else
716  hds = hnew(n)
717  term = -consterm * derv
718  rhs(n) = rhs(n) + term * hnew(m) !+ amat(idxglo(ii)) * (dwadi * hds - hds) !need to add dwadi
719  rhs(m) = rhs(m) - term * hnew(m) !- amat(idxglo(ii)) * (dwadi * hds - hds) !need to add dwadi
720  ! fill in row of n
721  filledterm = matrix_sln%get_value_pos(idxglo(idiag))
722  call matrix_sln%set_value_pos(idxglo(idiag), filledterm) !- filledterm * (dwadi - DONE) !need to add dwadi
723  ! fill newton term in off diagonal if active cell
724  if (this%ibound(n) > 0) then
725  call matrix_sln%add_value_pos(idxglo(ii), term)
726  end if
727  !fill row of m
728  call matrix_sln%add_value_pos(idxglo(idiagm), -term)
729  ! fill newton term in off diagonal if active cell
730  if (this%ibound(m) > 0) then
731  filledterm = matrix_sln%get_value_pos(idxglo(isymcon))
732  call matrix_sln%set_value_pos(idxglo(isymcon), filledterm) !* dwadi !need to add dwadi
733  end if
734  end if
735  end if
736 
737  end do
738  end do
739  !
740  end if
Here is the call graph for this function:

◆ npf_mc()

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

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

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

◆ npf_nur()

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

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

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

749  ! -- dummy
750  class(GwfNpfType) :: this
751  integer(I4B), intent(in) :: neqmod
752  real(DP), dimension(neqmod), intent(inout) :: x
753  real(DP), dimension(neqmod), intent(in) :: xtemp
754  real(DP), dimension(neqmod), intent(inout) :: dx
755  integer(I4B), intent(inout) :: inewtonur
756  real(DP), intent(inout) :: dxmax
757  integer(I4B), intent(inout) :: locmax
758  ! -- local
759  integer(I4B) :: n
760  real(DP) :: botm
761  real(DP) :: xx
762  real(DP) :: dxx
763  !
764  ! -- Newton-Raphson under-relaxation
765  do n = 1, this%dis%nodes
766  if (this%ibound(n) < 1) cycle
767  if (this%icelltype(n) > 0) then
768  botm = this%dis%bot(this%ibotnode(n))
769  ! -- only apply Newton-Raphson under-relaxation if
770  ! solution head is below the bottom of the model
771  if (x(n) < botm) then
772  inewtonur = 1
773  xx = xtemp(n) * (done - dp9) + botm * dp9
774  dxx = x(n) - xx
775  if (abs(dxx) > abs(dxmax)) then
776  locmax = n
777  dxmax = dxx
778  end if
779  x(n) = xx
780  dx(n) = dzero
781  end if
782  end if
783  end do

◆ npf_print_model_flows()

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

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

957  ! -- modules
958  use tdismodule, only: kper, kstp
959  use constantsmodule, only: lenbigline
960  ! -- dummy
961  class(GwfNpfType) :: this
962  integer(I4B), intent(in) :: ibudfl
963  real(DP), intent(inout), dimension(:) :: flowja
964  ! -- local
965  character(len=LENBIGLINE) :: line
966  character(len=30) :: tempstr
967  integer(I4B) :: n, ipos, m
968  real(DP) :: qnm
969  ! -- formats
970  character(len=*), parameter :: fmtiprflow = &
971  &"(/,4x,'CALCULATED INTERCELL FLOW FOR PERIOD ', i0, ' STEP ', i0)"
972  !
973  ! -- Write flowja to list file if requested
974  if (ibudfl /= 0 .and. this%iprflow > 0) then
975  write (this%iout, fmtiprflow) kper, kstp
976  do n = 1, this%dis%nodes
977  line = ''
978  call this%dis%noder_to_string(n, tempstr)
979  line = trim(tempstr)//':'
980  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
981  m = this%dis%con%ja(ipos)
982  call this%dis%noder_to_string(m, tempstr)
983  line = trim(line)//' '//trim(tempstr)
984  qnm = flowja(ipos)
985  write (tempstr, '(1pg15.6)') qnm
986  line = trim(line)//' '//trim(adjustl(tempstr))
987  end do
988  write (this%iout, '(a)') trim(line)
989  end do
990  end if
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15

◆ npf_rp()

subroutine gwfnpfmodule::npf_rp ( class(gwfnpftype this)

Read and prepare NPF stress period data.

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

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

◆ npf_save_model_flows()

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

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

920  ! -- dummy
921  class(GwfNpfType) :: this
922  real(DP), dimension(:), intent(in) :: flowja
923  integer(I4B), intent(in) :: icbcfl
924  integer(I4B), intent(in) :: icbcun
925  ! -- local
926  integer(I4B) :: ibinun
927  !
928  ! -- Set unit number for binary output
929  if (this%ipakcb < 0) then
930  ibinun = icbcun
931  elseif (this%ipakcb == 0) then
932  ibinun = 0
933  else
934  ibinun = this%ipakcb
935  end if
936  if (icbcfl == 0) ibinun = 0
937  !
938  ! -- Write the face flows if requested
939  if (ibinun /= 0) then
940  call this%dis%record_connection_array(flowja, ibinun, this%iout)
941  end if
942  !
943  ! -- Calculate specific discharge at cell centers and write, if requested
944  if (this%isavspdis /= 0) then
945  if (ibinun /= 0) call this%sav_spdis(ibinun)
946  end if
947  !
948  ! -- Save saturation, if requested
949  if (this%isavsat /= 0) then
950  if (ibinun /= 0) call this%sav_sat(ibinun)
951  end if

◆ prepare_edge_lookup()

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

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

2816  class(GwfNpfType) :: this
2817  ! local
2818  integer(I4B) :: i, inode, iedge
2819  integer(I4B) :: n, start, end
2820  integer(I4B) :: prev_cnt, strt_idx, ipos
2821 
2822  do i = 1, size(this%iedge_ptr)
2823  this%iedge_ptr(i) = 0
2824  end do
2825  do i = 1, size(this%edge_idxs)
2826  this%edge_idxs(i) = 0
2827  end do
2828 
2829  ! count
2830  do iedge = 1, this%nedges
2831  n = this%nodedge(iedge)
2832  this%iedge_ptr(n) = this%iedge_ptr(n) + 1
2833  end do
2834 
2835  ! determine start indexes
2836  prev_cnt = this%iedge_ptr(1)
2837  this%iedge_ptr(1) = 1
2838  do inode = 2, this%dis%nodes + 1
2839  strt_idx = this%iedge_ptr(inode - 1) + prev_cnt
2840  prev_cnt = this%iedge_ptr(inode)
2841  this%iedge_ptr(inode) = strt_idx
2842  end do
2843 
2844  ! loop over edges to fill lookup table
2845  do iedge = 1, this%nedges
2846  n = this%nodedge(iedge)
2847  start = this%iedge_ptr(n)
2848  end = this%iedge_ptr(n + 1) - 1
2849  do ipos = start, end
2850  if (this%edge_idxs(ipos) > 0) cycle ! go to next
2851  this%edge_idxs(ipos) = iedge
2852  exit
2853  end do
2854  end do
2855 

◆ prepcheck()

subroutine gwfnpfmodule::prepcheck ( class(gwfnpftype this)

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

1667  ! -- modules
1668  use constantsmodule, only: linelength, dpio180
1670  ! -- dummy
1671  class(GwfNpfType) :: this
1672  ! -- local
1673  character(len=24), dimension(:), pointer :: aname
1674  character(len=LINELENGTH) :: cellstr, errmsg
1675  integer(I4B) :: nerr, n
1676  ! -- format
1677  character(len=*), parameter :: fmtkerr = &
1678  &"(1x, 'Hydraulic property ',a,' is <= 0 for cell ',a, ' ', 1pg15.6)"
1679  character(len=*), parameter :: fmtkerr2 = &
1680  &"(1x, '... ', i0,' additional errors not shown for ',a)"
1681  !
1682  ! -- initialize
1683  aname => this%aname
1684  !
1685  ! -- check k11
1686  nerr = 0
1687  do n = 1, size(this%k11)
1688  if (this%k11(n) <= dzero) then
1689  nerr = nerr + 1
1690  if (nerr <= 20) then
1691  call this%dis%noder_to_string(n, cellstr)
1692  write (errmsg, fmtkerr) trim(adjustl(aname(2))), trim(cellstr), &
1693  this%k11(n)
1694  call store_error(errmsg)
1695  end if
1696  end if
1697  end do
1698  if (nerr > 20) then
1699  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(2)))
1700  call store_error(errmsg)
1701  end if
1702  !
1703  ! -- check k33 because it was read
1704  if (this%ik33 /= 0) then
1705  !
1706  ! -- Check to make sure values are greater than or equal to zero
1707  nerr = 0
1708  do n = 1, size(this%k33)
1709  if (this%ik33overk /= 0) this%k33(n) = this%k33(n) * this%k11(n)
1710  if (this%k33(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(3))), trim(cellstr), &
1715  this%k33(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(3)))
1722  call store_error(errmsg)
1723  end if
1724  end if
1725  !
1726  ! -- check k22 because it was read
1727  if (this%ik22 /= 0) then
1728  !
1729  ! -- Check to make sure that angles are available
1730  if (this%dis%con%ianglex == 0) then
1731  write (errmsg, '(a)') 'Error. ANGLDEGX not provided in '// &
1732  'discretization file, but K22 was specified. '
1733  call store_error(errmsg)
1734  end if
1735  !
1736  ! -- Check to make sure values are greater than or equal to zero
1737  nerr = 0
1738  do n = 1, size(this%k22)
1739  if (this%ik22overk /= 0) this%k22(n) = this%k22(n) * this%k11(n)
1740  if (this%k22(n) <= dzero) then
1741  nerr = nerr + 1
1742  if (nerr <= 20) then
1743  call this%dis%noder_to_string(n, cellstr)
1744  write (errmsg, fmtkerr) trim(adjustl(aname(4))), trim(cellstr), &
1745  this%k22(n)
1746  call store_error(errmsg)
1747  end if
1748  end if
1749  end do
1750  if (nerr > 20) then
1751  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(4)))
1752  call store_error(errmsg)
1753  end if
1754  end if
1755  !
1756  ! -- check for wetdry conflicts
1757  if (this%irewet == 1) then
1758  if (this%iwetdry == 0) then
1759  write (errmsg, '(a, a, a)') 'Error in GRIDDATA block: ', &
1760  trim(adjustl(aname(5))), ' not found.'
1761  call store_error(errmsg)
1762  end if
1763  end if
1764  !
1765  ! -- Check for angle conflicts
1766  if (this%iangle1 /= 0) then
1767  do n = 1, size(this%angle1)
1768  this%angle1(n) = this%angle1(n) * dpio180
1769  end do
1770  else
1771  if (this%ixt3d /= 0) then
1772  this%iangle1 = 1
1773  write (this%iout, '(a)') 'XT3D IN USE, BUT ANGLE1 NOT SPECIFIED. '// &
1774  'SETTING ANGLE1 TO ZERO.'
1775  do n = 1, size(this%angle1)
1776  this%angle1(n) = dzero
1777  end do
1778  end if
1779  end if
1780  if (this%iangle2 /= 0) then
1781  if (this%iangle1 == 0) then
1782  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE1. '// &
1783  'ANGLE2 REQUIRES ANGLE1. '
1784  call store_error(errmsg)
1785  end if
1786  if (this%iangle3 == 0) then
1787  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE3. '// &
1788  'SPECIFY BOTH OR NEITHER ONE. '
1789  call store_error(errmsg)
1790  end if
1791  do n = 1, size(this%angle2)
1792  this%angle2(n) = this%angle2(n) * dpio180
1793  end do
1794  end if
1795  if (this%iangle3 /= 0) then
1796  if (this%iangle1 == 0) then
1797  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE1. '// &
1798  'ANGLE3 REQUIRES ANGLE1. '
1799  call store_error(errmsg)
1800  end if
1801  if (this%iangle2 == 0) then
1802  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE2. '// &
1803  'SPECIFY BOTH OR NEITHER ONE. '
1804  call store_error(errmsg)
1805  end if
1806  do n = 1, size(this%angle3)
1807  this%angle3(n) = this%angle3(n) * dpio180
1808  end do
1809  end if
1810  !
1811  ! -- terminate if data errors
1812  if (count_errors() > 0) then
1813  call store_error_filename(this%input_fname)
1814  end if
real(dp), parameter dpio180
real constant
Definition: Constants.f90:130
Here is the call graph for this function:

◆ preprocess_input()

subroutine gwfnpfmodule::preprocess_input ( class(gwfnpftype this)

This routine consists of the following steps:

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

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

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

◆ rewet_check()

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

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

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

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

◆ sav_sat()

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

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

2718  ! -- dummy
2719  class(GwfNpfType) :: this
2720  integer(I4B), intent(in) :: ibinun
2721  ! -- local
2722  character(len=16) :: text
2723  character(len=16), dimension(1) :: auxtxt
2724  real(DP), dimension(1) :: a
2725  integer(I4B) :: n
2726  integer(I4B) :: naux
2727  !
2728  ! -- Write the header
2729  text = ' DATA-SAT'
2730  naux = 1
2731  auxtxt(:) = [' sat']
2732  call this%dis%record_srcdst_list_header(text, this%name_model, &
2733  this%packName, this%name_model, &
2734  this%packName, naux, auxtxt, ibinun, &
2735  this%dis%nodes, this%iout)
2736  !
2737  ! -- Write a zero for Q, and then write saturation as an aux variables
2738  do n = 1, this%dis%nodes
2739  a(1) = this%sat(n)
2740  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, a)
2741  end do

◆ sav_spdis()

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

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

2690  ! -- dummy
2691  class(GwfNpfType) :: this
2692  integer(I4B), intent(in) :: ibinun
2693  ! -- local
2694  character(len=16) :: text
2695  character(len=16), dimension(3) :: auxtxt
2696  integer(I4B) :: n
2697  integer(I4B) :: naux
2698  !
2699  ! -- Write the header
2700  text = ' DATA-SPDIS'
2701  naux = 3
2702  auxtxt(:) = [' qx', ' qy', ' qz']
2703  call this%dis%record_srcdst_list_header(text, this%name_model, &
2704  this%packName, this%name_model, &
2705  this%packName, naux, auxtxt, ibinun, &
2706  this%dis%nodes, this%iout)
2707  !
2708  ! -- Write a zero for Q, and then write qx, qy, qz as aux variables
2709  do n = 1, this%dis%nodes
2710  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, &
2711  this%spdis(:, n))
2712  end do

◆ set_edge_properties()

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

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

2788  ! -- dummy
2789  class(GwfNpfType) :: this
2790  integer(I4B), intent(in) :: nodedge
2791  integer(I4B), intent(in) :: ihcedge
2792  real(DP), intent(in) :: q
2793  real(DP), intent(in) :: area
2794  real(DP), intent(in) :: nx
2795  real(DP), intent(in) :: ny
2796  real(DP), intent(in) :: distance
2797  ! -- local
2798  integer(I4B) :: lastedge
2799  !
2800  this%lastedge = this%lastedge + 1
2801  lastedge = this%lastedge
2802  this%nodedge(lastedge) = nodedge
2803  this%ihcedge(lastedge) = ihcedge
2804  this%propsedge(1, lastedge) = q
2805  this%propsedge(2, lastedge) = area
2806  this%propsedge(3, lastedge) = nx
2807  this%propsedge(4, lastedge) = ny
2808  this%propsedge(5, lastedge) = distance
2809  !
2810  ! -- If this is the last edge, then the next call must be starting a new
2811  ! edge properties assignment loop, so need to reset lastedge to 0
2812  if (this%lastedge == this%nedges) this%lastedge = 0

◆ set_options()

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

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

1439  ! -- dummy
1440  class(GwfNpftype) :: this
1441  type(GwfNpfOptionsType), intent(in) :: options
1442  !
1443  this%ithickstrt = options%ithickstrt
1444  this%ihighcellsat = options%ihighcellsat
1445  this%iperched = options%iperched
1446  this%ivarcv = options%ivarcv
1447  this%idewatcv = options%idewatcv
1448  this%irewet = options%irewet
1449  this%wetfct = options%wetfct
1450  this%iwetit = options%iwetit
1451  this%ihdwet = options%ihdwet

◆ sgwf_npf_qcalc()

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

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

842  ! -- dummy
843  class(GwfNpfType) :: this
844  integer(I4B), intent(in) :: n
845  integer(I4B), intent(in) :: m
846  real(DP), intent(in) :: hn
847  real(DP), intent(in) :: hm
848  integer(I4B), intent(in) :: icon
849  real(DP), intent(inout) :: qnm
850  ! -- local
851  real(DP) :: hyn, hym
852  real(DP) :: condnm
853  real(DP) :: hntemp, hmtemp
854  real(DP) :: satn, satm
855  integer(I4B) :: ihc
856  !
857  ! -- Initialize
858  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
859  hyn = this%hy_eff(n, m, ihc, ipos=icon)
860  hym = this%hy_eff(m, n, ihc, ipos=icon)
861  !
862  ! -- Calculate conductance
863  if (ihc == c3d_vertical) then
864  condnm = vcond(this%ibound(n), this%ibound(m), &
865  this%icelltype(n), this%icelltype(m), this%inewton, &
866  this%ivarcv, this%idewatcv, &
867  this%condsat(this%dis%con%jas(icon)), hn, hm, &
868  hyn, hym, &
869  this%sat(n), this%sat(m), &
870  this%dis%top(n), this%dis%top(m), &
871  this%dis%bot(n), this%dis%bot(m), &
872  this%dis%con%hwva(this%dis%con%jas(icon)))
873  else
874  satn = this%sat(n)
875  satm = this%sat(m)
876  if (this%ihighcellsat /= 0) then
877  call this%highest_cell_saturation(n, m, hn, hm, satn, satm)
878  end if
879 
880  condnm = hcond(this%ibound(n), this%ibound(m), &
881  this%icelltype(n), this%icelltype(m), &
882  this%inewton, &
883  this%dis%con%ihc(this%dis%con%jas(icon)), &
884  this%icellavg, &
885  this%condsat(this%dis%con%jas(icon)), &
886  hn, hm, satn, satm, hyn, hym, &
887  this%dis%top(n), this%dis%top(m), &
888  this%dis%bot(n), this%dis%bot(m), &
889  this%dis%con%cl1(this%dis%con%jas(icon)), &
890  this%dis%con%cl2(this%dis%con%jas(icon)), &
891  this%dis%con%hwva(this%dis%con%jas(icon)))
892  end if
893  !
894  ! -- Initialize hntemp and hmtemp
895  hntemp = hn
896  hmtemp = hm
897  !
898  ! -- Check and adjust for dewatered conditions
899  if (this%iperched /= 0) then
900  if (this%dis%con%ihc(this%dis%con%jas(icon)) == 0) then
901  if (n > m) then
902  if (this%icelltype(n) /= 0) then
903  if (hn < this%dis%top(n)) hntemp = this%dis%bot(m)
904  end if
905  else
906  if (this%icelltype(m) /= 0) then
907  if (hm < this%dis%top(m)) hmtemp = this%dis%bot(n)
908  end if
909  end if
910  end if
911  end if
912  !
913  ! -- Calculate flow positive into cell n
914  qnm = condnm * (hmtemp - hntemp)
Here is the call graph for this function:

◆ sgwf_npf_thksat()

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

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

819  ! -- dummy
820  class(GwfNpfType) :: this
821  integer(I4B), intent(in) :: n
822  real(DP), intent(in) :: hn
823  real(DP), intent(inout) :: thksat
824  !
825  ! -- Standard Formulation
826  if (hn >= this%dis%top(n)) then
827  thksat = done
828  else
829  thksat = (hn - this%dis%bot(n)) / (this%dis%top(n) - this%dis%bot(n))
830  end if
831  !
832  ! -- Newton-Raphson Formulation
833  if (this%inewton /= 0) then
834  thksat = squadraticsaturation(this%dis%top(n), this%dis%bot(n), hn, &
835  this%satomega)
836  end if
Here is the call graph for this function:

◆ sgwf_npf_wdmsg()

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

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

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

◆ sgwf_npf_wetdry()

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

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

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

◆ source_griddata()

subroutine gwfnpfmodule::source_griddata ( class(gwfnpftype this)

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

1578  ! -- modules
1579  use simmodule, only: count_errors, store_error
1583  ! -- dummy
1584  class(GwfNpftype) :: this
1585  ! -- locals
1586  character(len=LINELENGTH) :: errmsg
1587  type(GwfNpfParamFoundType) :: found
1588  logical, dimension(2) :: afound
1589  integer(I4B), dimension(:), pointer, contiguous :: map
1590  !
1591  ! -- set map to convert user input data into reduced data
1592  map => null()
1593  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1594  !
1595  ! -- update defaults with idm sourced values
1596  call mem_set_value(this%icelltype, 'ICELLTYPE', this%input_mempath, map, &
1597  found%icelltype)
1598  call mem_set_value(this%k11, 'K', this%input_mempath, map, found%k)
1599  call mem_set_value(this%k33, 'K33', this%input_mempath, map, found%k33)
1600  call mem_set_value(this%k22, 'K22', this%input_mempath, map, found%k22)
1601  call mem_set_value(this%wetdry, 'WETDRY', this%input_mempath, map, &
1602  found%wetdry)
1603  call mem_set_value(this%angle1, 'ANGLE1', this%input_mempath, map, &
1604  found%angle1)
1605  call mem_set_value(this%angle2, 'ANGLE2', this%input_mempath, map, &
1606  found%angle2)
1607  call mem_set_value(this%angle3, 'ANGLE3', this%input_mempath, map, &
1608  found%angle3)
1609  !
1610  ! -- ensure ICELLTYPE was found
1611  if (.not. found%icelltype) then
1612  write (errmsg, '(a)') 'Error in GRIDDATA block: ICELLTYPE not found.'
1613  call store_error(errmsg)
1614  end if
1615  !
1616  ! -- ensure K was found
1617  if (.not. found%k) then
1618  write (errmsg, '(a)') 'Error in GRIDDATA block: K not found.'
1619  call store_error(errmsg)
1620  end if
1621  !
1622  ! -- set error if ik33overk set with no k33
1623  if (.not. found%k33 .and. this%ik33overk /= 0) then
1624  write (errmsg, '(a)') 'K33OVERK option specified but K33 not specified.'
1625  call store_error(errmsg)
1626  end if
1627  !
1628  ! -- set error if ik22overk set with no k22
1629  if (.not. found%k22 .and. this%ik22overk /= 0) then
1630  write (errmsg, '(a)') 'K22OVERK option specified but K22 not specified.'
1631  call store_error(errmsg)
1632  end if
1633  !
1634  ! -- handle found side effects
1635  if (found%k33) this%ik33 = 1
1636  if (found%k22) this%ik22 = 1
1637  if (found%wetdry) this%iwetdry = 1
1638  if (found%angle1) this%iangle1 = 1
1639  if (found%angle2) this%iangle2 = 1
1640  if (found%angle3) this%iangle3 = 1
1641  !
1642  ! -- handle not found side effects
1643  if (.not. found%k33) then
1644  call mem_set_value(this%k33, 'K', this%input_mempath, map, afound(1))
1645  end if
1646  if (.not. found%k22) then
1647  call mem_set_value(this%k22, 'K', this%input_mempath, map, afound(2))
1648  end if
1649  if (.not. found%wetdry) call mem_reallocate(this%wetdry, 1, 'WETDRY', &
1650  trim(this%memoryPath))
1651  if (.not. found%angle1 .and. this%ixt3d == 0) &
1652  call mem_reallocate(this%angle1, 0, 'ANGLE1', trim(this%memoryPath))
1653  if (.not. found%angle2 .and. this%ixt3d == 0) &
1654  call mem_reallocate(this%angle2, 0, 'ANGLE2', trim(this%memoryPath))
1655  if (.not. found%angle3 .and. this%ixt3d == 0) &
1656  call mem_reallocate(this%angle3, 0, 'ANGLE3', trim(this%memoryPath))
1657  !
1658  ! -- log griddata
1659  if (this%iout > 0) then
1660  call this%log_griddata(found)
1661  end if
Here is the call graph for this function:

◆ source_options()

subroutine gwfnpfmodule::source_options ( class(gwfnpftype this)

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

1352  ! -- modules
1358  use sourcecommonmodule, only: filein_fname
1359  ! -- dummy
1360  class(GwfNpftype) :: this
1361  ! -- locals
1362  character(len=LENVARNAME), dimension(3) :: cellavg_method = &
1363  &[character(len=LENVARNAME) :: 'LOGARITHMIC', 'AMT-LMK', 'AMT-HMK']
1364  type(gwfnpfparamfoundtype) :: found
1365  character(len=LINELENGTH) :: tvk6_filename
1366  !
1367  ! -- update defaults with idm sourced values
1368  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
1369  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
1370  call mem_set_value(this%icellavg, 'CELLAVG', this%input_mempath, &
1371  cellavg_method, found%cellavg)
1372  call mem_set_value(this%ithickstrt, 'ITHICKSTRT', this%input_mempath, &
1373  found%ithickstrt)
1374  call mem_set_value(this%ihighcellsat, 'IHIGHCELLSAT', this%input_mempath, &
1375  found%ihighcellsat)
1376  call mem_set_value(this%iperched, 'IPERCHED', this%input_mempath, &
1377  found%iperched)
1378  call mem_set_value(this%ivarcv, 'IVARCV', this%input_mempath, found%ivarcv)
1379  call mem_set_value(this%idewatcv, 'IDEWATCV', this%input_mempath, &
1380  found%idewatcv)
1381  call mem_set_value(this%ixt3d, 'IXT3D', this%input_mempath, found%ixt3d)
1382  call mem_set_value(this%ixt3drhs, 'IXT3DRHS', this%input_mempath, &
1383  found%ixt3drhs)
1384  call mem_set_value(this%isavspdis, 'ISAVSPDIS', this%input_mempath, &
1385  found%isavspdis)
1386  call mem_set_value(this%isavsat, 'ISAVSAT', this%input_mempath, found%isavsat)
1387  call mem_set_value(this%ik22overk, 'IK22OVERK', this%input_mempath, &
1388  found%ik22overk)
1389  call mem_set_value(this%ik33overk, 'IK33OVERK', this%input_mempath, &
1390  found%ik33overk)
1391  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
1392  call mem_set_value(this%satomega, 'SATOMEGA', this%input_mempath, &
1393  found%satomega)
1394  call mem_set_value(this%irewet, 'IREWET', this%input_mempath, found%irewet)
1395  call mem_set_value(this%wetfct, 'WETFCT', this%input_mempath, found%wetfct)
1396  call mem_set_value(this%iwetit, 'IWETIT', this%input_mempath, found%iwetit)
1397  call mem_set_value(this%ihdwet, 'IHDWET', this%input_mempath, found%ihdwet)
1398  !
1399  ! -- save flows option active
1400  if (found%ipakcb) this%ipakcb = -1
1401  !
1402  ! -- xt3d active with rhs
1403  if (found%ixt3d .and. found%ixt3drhs) this%ixt3d = 2
1404  !
1405  ! -- save specific discharge active
1406  if (found%isavspdis) this%icalcspdis = this%isavspdis
1407  !
1408  ! -- no newton specified
1409  if (found%inewton) then
1410  this%inewton = 0
1411  this%iasym = 0
1412  end if
1413  !
1414  ! -- enforce 0 or 1 TVK6_FILENAME entries in option block
1415  if (filein_fname(tvk6_filename, 'TVK6_FILENAME', this%input_mempath, &
1416  this%input_fname)) then
1417  call openfile(this%intvk, this%iout, tvk6_filename, 'TVK')
1418  call tvk_cr(this%tvk, this%name_model, this%intvk, this%iout)
1419  end if
1420  !
1421  ! -- verify ALTERNATIVE_CELL_AVERAGING input value is supported
1422  if (found%cellavg) then
1423  if (this%icellavg == 0) then
1424  errmsg = 'Unrecognized input value for ALTERNATIVE_CELL_AVERAGING option.'
1425  call store_error(errmsg)
1426  call store_error_filename(this%input_fname)
1427  end if
1428  end if
1429  !
1430  ! -- log options
1431  if (this%iout > 0) then
1432  call this%log_options(found)
1433  end if
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function:

◆ store_original_k_arrays()

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

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

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

1197  ! -- modules
1199  ! -- dummy
1200  class(GwfNpftype) :: this
1201  integer(I4B), intent(in) :: ncells
1202  integer(I4B), intent(in) :: njas
1203  ! -- local
1204  integer(I4B) :: n
1205  !
1206  ! -- Retain copy of user-specified K arrays
1207  do n = 1, ncells
1208  this%k11input(n) = this%k11(n)
1209  this%k22input(n) = this%k22(n)
1210  this%k33input(n) = this%k33(n)
1211  end do