MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
gwf-npf.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
4  use constantsmodule, only: dzero, dem9, dem8, dem7, dem6, dem2, &
5  dhalf, dp9, done, dtwo, &
6  dhnoflo, dhdry, dem10, &
13  use basedismodule, only: disbasetype
14  use gwficmodule, only: gwfictype
15  use gwfvscmodule, only: gwfvsctype
16  use xt3dmodule, only: xt3dtype
19  use tvkmodule, only: tvktype, tvk_cr
24  use hgeoutilmodule, only: hyeff
26  condmean, thksatnm, &
28 
29  implicit none
30 
31  private
32  public :: gwfnpftype
33  public :: npf_cr
34 
35  type, extends(numericalpackagetype) :: gwfnpftype
36 
37  type(gwfictype), pointer :: ic => null() !< initial conditions object
38  type(gwfvsctype), pointer :: vsc => null() !< viscosity object
39  type(xt3dtype), pointer :: xt3d => null() !< xt3d pointer
40  integer(I4B), pointer :: iname => null() !< length of variable names
41  character(len=24), dimension(:), pointer :: aname => null() !< variable names
42  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< pointer to model ibound
43  real(dp), dimension(:), pointer, contiguous :: hnew => null() !< pointer to model xnew
44  integer(I4B), pointer :: ixt3d => null() !< xt3d flag (0 is off, 1 is lhs, 2 is rhs)
45  integer(I4B), pointer :: ixt3drhs => null() !< xt3d rhs flag, xt3d rhs is set active if 1
46  integer(I4B), pointer :: iperched => null() !< vertical flow corrections if 1
47  integer(I4B), pointer :: ivarcv => null() !< CV is function of water table
48  integer(I4B), pointer :: idewatcv => null() !< CV may be a discontinuous function of water table
49  integer(I4B), pointer :: ithickstrt => null() !< thickstrt option flag
50  integer(I4B), pointer :: ihighcellsat => null() !< highest_cell_saturation option flag
51  integer(I4B), pointer :: igwfnewtonur => null() !< newton head dampening using node bottom option flag
52  integer(I4B), pointer :: icalcspdis => null() !< Calculate specific discharge at cell centers
53  integer(I4B), pointer :: isavspdis => null() !< Save specific discharge at cell centers
54  integer(I4B), pointer :: isavsat => null() !< Save sat to budget file
55  real(dp), pointer :: hnoflo => null() !< default is 1.e30
56  real(dp), pointer :: satomega => null() !< newton-raphson saturation omega
57  integer(I4B), pointer :: irewet => null() !< rewetting (0:off, 1:on)
58  integer(I4B), pointer :: iwetit => null() !< wetting interval (default is 1)
59  integer(I4B), pointer :: ihdwet => null() !< (0 or not 0)
60  integer(I4B), pointer :: icellavg => null() !< harmonic(0), logarithmic(1), or arithmetic thick-log K (2)
61  real(dp), pointer :: wetfct => null() !< wetting factor
62  real(dp), pointer :: hdry => null() !< default is -1.d30
63  integer(I4B), dimension(:), pointer, contiguous :: icelltype => null() !< confined (0) or convertible (1)
64  integer(I4B), dimension(:), pointer, contiguous :: ithickstartflag => null() !< array of flags for handling the thickstrt option
65  !
66  ! K properties
67  real(dp), dimension(:), pointer, contiguous :: k11 => null() !< hydraulic conductivity; if anisotropic, then this is Kx prior to rotation
68  real(dp), dimension(:), pointer, contiguous :: k22 => null() !< hydraulic conductivity; if specified then this is Ky prior to rotation
69  real(dp), dimension(:), pointer, contiguous :: k33 => null() !< hydraulic conductivity; if specified then this is Kz prior to rotation
70  real(dp), dimension(:), pointer, contiguous :: k11input => null() !< hydraulic conductivity originally specified by user prior to TVK or VSC modification
71  real(dp), dimension(:), pointer, contiguous :: k22input => null() !< hydraulic conductivity originally specified by user prior to TVK or VSC modification
72  real(dp), dimension(:), pointer, contiguous :: k33input => null() !< hydraulic conductivity originally specified by user prior to TVK or VSC modification
73  integer(I4B), pointer :: iavgkeff => null() !< effective conductivity averaging (0: harmonic, 1: arithmetic)
74  integer(I4B), pointer :: ik22 => null() !< flag that k22 is specified
75  integer(I4B), pointer :: ik33 => null() !< flag that k33 is specified
76  integer(I4B), pointer :: ik22overk => null() !< flag that k22 is specified as anisotropy ratio
77  integer(I4B), pointer :: ik33overk => null() !< flag that k33 is specified as anisotropy ratio
78  integer(I4B), pointer :: iangle1 => null() !< flag to indicate angle1 was read
79  integer(I4B), pointer :: iangle2 => null() !< flag to indicate angle2 was read
80  integer(I4B), pointer :: iangle3 => null() !< flag to indicate angle3 was read
81  real(dp), dimension(:), pointer, contiguous :: angle1 => null() !< k ellipse rotation in xy plane around z axis (yaw)
82  real(dp), dimension(:), pointer, contiguous :: angle2 => null() !< k ellipse rotation up from xy plane around y axis (pitch)
83  real(dp), dimension(:), pointer, contiguous :: angle3 => null() !< k tensor rotation around x axis (roll)
84  !
85  integer(I4B), pointer :: iwetdry => null() !< flag to indicate angle1 was read
86  real(dp), dimension(:), pointer, contiguous :: wetdry => null() !< wetdry array
87  real(dp), dimension(:), pointer, contiguous :: sat => null() !< saturation (0. to 1.) for each cell
88  real(dp), dimension(:), pointer, contiguous :: condsat => null() !< saturated conductance (symmetric array)
89  integer(I4B), dimension(:), pointer, contiguous :: ibotnode => null() !< bottom node used if igwfnewtonur /= 0
90  ! spdis machinery:
91  real(dp), dimension(:, :), pointer, contiguous :: spdis => null() !< specific discharge : qx, qy, qz (nodes, 3)
92  integer(I4B), pointer :: nedges => null() !< number of cell edges
93  integer(I4B), pointer :: lastedge => null() !< last edge number
94  integer(I4B), dimension(:), pointer, contiguous :: nodedge => null() !< array of node numbers that have edges
95  integer(I4B), dimension(:), pointer, contiguous :: ihcedge => null() !< edge type (horizontal or vertical)
96  real(dp), dimension(:, :), pointer, contiguous :: propsedge => null() !< edge properties (Q, area, nx, ny, distance)
97  integer(I4B), dimension(:), pointer, contiguous :: iedge_ptr => null() !< csr pointer into edge index array
98  integer(I4B), dimension(:), pointer, contiguous :: edge_idxs => null() !< sorted edge indexes for faster lookup
99  type(spdisworkarraytype), pointer :: spdis_wa => null() !< work arrays for spdis calculation
100  !
101  integer(I4B), pointer :: intvk => null() ! TVK (time-varying K) unit number (0 if unused)
102  integer(I4B), pointer :: invsc => null() ! VSC (viscosity) unit number (0 if unused); viscosity leads to time-varying K's
103  type(tvktype), pointer :: tvk => null() ! TVK object
104  integer(I4B), pointer :: kchangeper => null() ! last stress period in which any node K (or K22, or K33) values were changed (0 if unchanged from start of simulation)
105  integer(I4B), pointer :: kchangestp => null() ! last time step in which any node K (or K22, or K33) values were changed (0 if unchanged from start of simulation)
106  integer(I4B), dimension(:), pointer, contiguous :: nodekchange => null() ! grid array of flags indicating for each node whether its K (or K22, or K33) value changed (1) at (kchangeper, kchangestp) or not (0)
107 
108  contains
109 
110  procedure :: npf_df
111  procedure :: npf_ac
112  procedure :: npf_mc
113  procedure :: npf_ar
114  procedure :: npf_rp
115  procedure :: npf_ad
116  procedure :: npf_cf
117  procedure :: npf_fc
118  procedure :: npf_fn
119  procedure :: npf_cq
120  procedure :: npf_save_model_flows
121  procedure :: npf_nur
123  procedure :: npf_da
124  procedure, private :: thksat => sgwf_npf_thksat
125  procedure, private :: qcalc => sgwf_npf_qcalc
126  procedure, private :: wd => sgwf_npf_wetdry
127  procedure, private :: wdmsg => sgwf_npf_wdmsg
128  procedure :: allocate_scalars
129  procedure, private :: store_original_k_arrays
130  procedure, private :: allocate_arrays
131  procedure, private :: source_options
132  procedure, private :: source_griddata
133  procedure, private :: log_options
134  procedure, private :: log_griddata
135  procedure, private :: set_options
136  procedure, private :: check_options
137  procedure, private :: prepcheck
138  procedure, private :: preprocess_input
139  procedure, private :: calc_condsat
140  procedure, private :: calc_initial_sat
141  procedure, public :: rewet_check
142  procedure, public :: hy_eff
143  procedure, public :: calc_spdis
144  procedure, public :: sav_spdis
145  procedure, public :: sav_sat
146  procedure, public :: increase_edge_count
147  procedure, public :: set_edge_properties
148  procedure, public :: calcsatthickness
149  procedure, private :: calc_max_conns
150  procedure, private :: prepare_edge_lookup
151  procedure, private :: highest_cell_saturation
152  end type
153 
154 contains
155 
156  !> @brief Create a new NPF object. Pass a inunit value of 0 if npf data will
157  !! initialized from memory
158  !<
159  subroutine npf_cr(npfobj, name_model, input_mempath, inunit, iout)
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 
197  end subroutine npf_cr
198 
199  !> @brief Define the NPF package instance
200  !!
201  !! This is a hybrid routine: it either reads the options for this package
202  !! from the input file, or the optional argument @param npf_options
203  !! should be passed. A consistency check is performed, and finally
204  !! xt3d_df is called, when enabled.
205  !<
206  subroutine npf_df(this, dis, xt3d, ingnc, invsc, npf_options)
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
255  end subroutine npf_df
256 
257  !> @brief Add connections for extended neighbors to the sparse matrix
258  !<
259  subroutine npf_ac(this, moffset, sparse)
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)
269  end subroutine npf_ac
270 
271  !> @brief Map connections and construct iax, jax, and idxglox
272  !<
273  subroutine npf_mc(this, moffset, matrix_sln)
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)
280  end subroutine npf_mc
281 
282  !> @brief Allocate and read this NPF instance
283  !!
284  !! Allocate remaining package arrays, preprocess the input data and
285  !! call *_ar on xt3d, when active.
286  !<
287  subroutine npf_ar(this, ic, vsc, ibound, hnew)
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
361  end subroutine npf_ar
362 
363  !> @brief Read and prepare method for package
364  !!
365  !! Read and prepare NPF stress period data.
366  !<
367  subroutine npf_rp(this)
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
376  end subroutine npf_rp
377 
378  !> @brief Advance
379  !!
380  !! Sets hold (head old) to bot whenever a wettable cell is dry
381  !<
382  subroutine npf_ad(this, nodes, hold, hnew, irestore)
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
442  end subroutine npf_ad
443 
444  !> @brief Routines associated fill coefficients
445  !<
446  subroutine npf_cf(this, kiter, nodes, hnew)
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
472  end subroutine npf_cf
473 
474  !> @brief Formulate coefficients
475  !<
476  subroutine npf_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
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
585  end subroutine npf_fc
586 
587  !> @brief Calculate dry cell saturation
588  !!
589  !! Calculate the saturation based on the maximum cell bottom for
590  !! two connected cells
591  !<
592  subroutine highest_cell_saturation(this, n, m, hn, hm, satn, satm)
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
617  end subroutine highest_cell_saturation
618 
619  !> @brief Fill newton terms
620  !<
621  subroutine npf_fn(this, kiter, matrix_sln, idxglo, rhs, hnew)
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
741  end subroutine npf_fn
742 
743  !> @brief Under-relaxation
744  !!
745  !! Under-relaxation of Groundwater Flow Model Heads for current outer
746  !! iteration using the cell bottoms at the bottom of the model
747  !<
748  subroutine npf_nur(this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
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  integer(I4B) :: ibot
761  real(DP) :: botm
762  real(DP) :: xx
763  real(DP) :: dxx
764  !
765  ! -- Newton-Raphson under-relaxation
766  do n = 1, this%dis%nodes
767  if (this%ibound(n) < 1) cycle
768  ibot = this%ibotnode(n)
769  ! Newton-Raphson under-relaxation is only applied to convertible cells where
770  ! the bottom cell in a stack is convertible
771  if (this%icelltype(n) > 0 .and. this%icelltype(ibot) > 0) then
772  botm = this%dis%bot(ibot)
773  ! Newton-Raphson under-relaxation applied when solution head is
774  ! below the bottom of the model
775  if (x(n) < botm) then
776  inewtonur = 1
777  xx = xtemp(n) * (done - dp9) + botm * dp9
778  dxx = xx - xtemp(n)
779  if (abs(dxx) > abs(dxmax)) then
780  locmax = n
781  dxmax = dxx
782  end if
783  x(n) = xx
784  dx(n) = xtemp(n) - x(n)
785  end if
786  end if
787  end do
788  end subroutine npf_nur
789 
790  !> @brief Calculate flowja
791  !<
792  subroutine npf_cq(this, hnew, flowja)
793  ! -- dummy
794  class(gwfnpftype) :: this
795  real(DP), intent(inout), dimension(:) :: hnew
796  real(DP), intent(inout), dimension(:) :: flowja
797  ! -- local
798  integer(I4B) :: n, ipos, m
799  real(DP) :: qnm
800  !
801  ! -- Calculate the flow across each cell face and store in flowja
802  !
803  if (this%ixt3d /= 0) then
804  call this%xt3d%xt3d_flowja(hnew, flowja)
805  else
806  !
807  do n = 1, this%dis%nodes
808  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
809  m = this%dis%con%ja(ipos)
810  if (m < n) cycle
811  call this%qcalc(n, m, hnew(n), hnew(m), ipos, qnm)
812  flowja(ipos) = qnm
813  flowja(this%dis%con%isym(ipos)) = -qnm
814  end do
815  end do
816  !
817  end if
818  end subroutine npf_cq
819 
820  !> @brief Fractional cell saturation
821  !<
822  subroutine sgwf_npf_thksat(this, n, hn, thksat)
823  ! -- dummy
824  class(gwfnpftype) :: this
825  integer(I4B), intent(in) :: n
826  real(DP), intent(in) :: hn
827  real(DP), intent(inout) :: thksat
828  !
829  ! -- Standard Formulation
830  if (hn >= this%dis%top(n)) then
831  thksat = done
832  else
833  thksat = (hn - this%dis%bot(n)) / (this%dis%top(n) - this%dis%bot(n))
834  end if
835  !
836  ! -- Newton-Raphson Formulation
837  if (this%inewton /= 0) then
838  thksat = squadraticsaturation(this%dis%top(n), this%dis%bot(n), hn, &
839  this%satomega)
840  end if
841  end subroutine sgwf_npf_thksat
842 
843  !> @brief Flow between two cells
844  !<
845  subroutine sgwf_npf_qcalc(this, n, m, hn, hm, icon, qnm)
846  ! -- dummy
847  class(gwfnpftype) :: this
848  integer(I4B), intent(in) :: n
849  integer(I4B), intent(in) :: m
850  real(DP), intent(in) :: hn
851  real(DP), intent(in) :: hm
852  integer(I4B), intent(in) :: icon
853  real(DP), intent(inout) :: qnm
854  ! -- local
855  real(DP) :: hyn, hym
856  real(DP) :: condnm
857  real(DP) :: hntemp, hmtemp
858  real(DP) :: satn, satm
859  integer(I4B) :: ihc
860  !
861  ! -- Initialize
862  ihc = this%dis%con%ihc(this%dis%con%jas(icon))
863  hyn = this%hy_eff(n, m, ihc, ipos=icon)
864  hym = this%hy_eff(m, n, ihc, ipos=icon)
865  !
866  ! -- Calculate conductance
867  if (ihc == c3d_vertical) then
868  condnm = vcond(this%ibound(n), this%ibound(m), &
869  this%icelltype(n), this%icelltype(m), this%inewton, &
870  this%ivarcv, this%idewatcv, &
871  this%condsat(this%dis%con%jas(icon)), hn, hm, &
872  hyn, hym, &
873  this%sat(n), this%sat(m), &
874  this%dis%top(n), this%dis%top(m), &
875  this%dis%bot(n), this%dis%bot(m), &
876  this%dis%con%hwva(this%dis%con%jas(icon)))
877  else
878  satn = this%sat(n)
879  satm = this%sat(m)
880  if (this%ihighcellsat /= 0) then
881  call this%highest_cell_saturation(n, m, hn, hm, satn, satm)
882  end if
883 
884  condnm = hcond(this%ibound(n), this%ibound(m), &
885  this%icelltype(n), this%icelltype(m), &
886  this%inewton, &
887  this%dis%con%ihc(this%dis%con%jas(icon)), &
888  this%icellavg, &
889  this%condsat(this%dis%con%jas(icon)), &
890  hn, hm, satn, satm, hyn, hym, &
891  this%dis%top(n), this%dis%top(m), &
892  this%dis%bot(n), this%dis%bot(m), &
893  this%dis%con%cl1(this%dis%con%jas(icon)), &
894  this%dis%con%cl2(this%dis%con%jas(icon)), &
895  this%dis%con%hwva(this%dis%con%jas(icon)))
896  end if
897  !
898  ! -- Initialize hntemp and hmtemp
899  hntemp = hn
900  hmtemp = hm
901  !
902  ! -- Check and adjust for dewatered conditions
903  if (this%iperched /= 0) then
904  if (this%dis%con%ihc(this%dis%con%jas(icon)) == 0) then
905  if (n > m) then
906  if (this%icelltype(n) /= 0) then
907  if (hn < this%dis%top(n)) hntemp = this%dis%bot(m)
908  end if
909  else
910  if (this%icelltype(m) /= 0) then
911  if (hm < this%dis%top(m)) hmtemp = this%dis%bot(n)
912  end if
913  end if
914  end if
915  end if
916  !
917  ! -- Calculate flow positive into cell n
918  qnm = condnm * (hmtemp - hntemp)
919  end subroutine sgwf_npf_qcalc
920 
921  !> @brief Record flowja and calculate specific discharge if requested
922  !<
923  subroutine npf_save_model_flows(this, flowja, icbcfl, icbcun)
924  ! -- dummy
925  class(gwfnpftype) :: this
926  real(DP), dimension(:), intent(in) :: flowja
927  integer(I4B), intent(in) :: icbcfl
928  integer(I4B), intent(in) :: icbcun
929  ! -- local
930  integer(I4B) :: ibinun
931  !
932  ! -- Set unit number for binary output
933  if (this%ipakcb < 0) then
934  ibinun = icbcun
935  elseif (this%ipakcb == 0) then
936  ibinun = 0
937  else
938  ibinun = this%ipakcb
939  end if
940  if (icbcfl == 0) ibinun = 0
941  !
942  ! -- Write the face flows if requested
943  if (ibinun /= 0) then
944  call this%dis%record_connection_array(flowja, ibinun, this%iout)
945  end if
946  !
947  ! -- Calculate specific discharge at cell centers and write, if requested
948  if (this%isavspdis /= 0) then
949  if (ibinun /= 0) call this%sav_spdis(ibinun)
950  end if
951  !
952  ! -- Save saturation, if requested
953  if (this%isavsat /= 0) then
954  if (ibinun /= 0) call this%sav_sat(ibinun)
955  end if
956  end subroutine npf_save_model_flows
957 
958  !> @brief Print budget
959  !<
960  subroutine npf_print_model_flows(this, ibudfl, flowja)
961  ! -- modules
962  use tdismodule, only: kper, kstp
963  use constantsmodule, only: lenbigline
964  ! -- dummy
965  class(gwfnpftype) :: this
966  integer(I4B), intent(in) :: ibudfl
967  real(DP), intent(inout), dimension(:) :: flowja
968  ! -- local
969  character(len=LENBIGLINE) :: line
970  character(len=30) :: tempstr
971  integer(I4B) :: n, ipos, m
972  real(DP) :: qnm
973  ! -- formats
974  character(len=*), parameter :: fmtiprflow = &
975  &"(/,4x,'CALCULATED INTERCELL FLOW FOR PERIOD ', i0, ' STEP ', i0)"
976  !
977  ! -- Write flowja to list file if requested
978  if (ibudfl /= 0 .and. this%iprflow > 0) then
979  write (this%iout, fmtiprflow) kper, kstp
980  do n = 1, this%dis%nodes
981  line = ''
982  call this%dis%noder_to_string(n, tempstr)
983  line = trim(tempstr)//':'
984  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
985  m = this%dis%con%ja(ipos)
986  call this%dis%noder_to_string(m, tempstr)
987  line = trim(line)//' '//trim(tempstr)
988  qnm = flowja(ipos)
989  write (tempstr, '(1pg15.6)') qnm
990  line = trim(line)//' '//trim(adjustl(tempstr))
991  end do
992  write (this%iout, '(a)') trim(line)
993  end do
994  end if
995  end subroutine npf_print_model_flows
996 
997  !> @brief Deallocate variables
998  !<
999  subroutine npf_da(this)
1000  ! -- modules
1002  use simvariablesmodule, only: idm_context
1003  ! -- dummy
1004  class(gwfnpftype) :: this
1005 
1006  ! free spdis work structure
1007  if (this%icalcspdis == 1 .and. this%spdis_wa%is_created()) &
1008  call this%spdis_wa%destroy()
1009  deallocate (this%spdis_wa)
1010  !
1011  ! -- Deallocate input memory
1012  call memorystore_remove(this%name_model, 'NPF', idm_context)
1013  !
1014  ! -- TVK
1015  if (this%intvk /= 0) then
1016  call this%tvk%da()
1017  deallocate (this%tvk)
1018  end if
1019  !
1020  ! -- VSC
1021  if (this%invsc /= 0) then
1022  nullify (this%vsc)
1023  end if
1024  !
1025  ! -- Strings
1026  !
1027  ! -- Scalars
1028  call mem_deallocate(this%iname)
1029  call mem_deallocate(this%ixt3d)
1030  call mem_deallocate(this%ixt3drhs)
1031  call mem_deallocate(this%satomega)
1032  call mem_deallocate(this%hnoflo)
1033  call mem_deallocate(this%hdry)
1034  call mem_deallocate(this%icellavg)
1035  call mem_deallocate(this%iavgkeff)
1036  call mem_deallocate(this%ik22)
1037  call mem_deallocate(this%ik33)
1038  call mem_deallocate(this%iperched)
1039  call mem_deallocate(this%ivarcv)
1040  call mem_deallocate(this%idewatcv)
1041  call mem_deallocate(this%ithickstrt)
1042  call mem_deallocate(this%ihighcellsat)
1043  call mem_deallocate(this%isavspdis)
1044  call mem_deallocate(this%isavsat)
1045  call mem_deallocate(this%icalcspdis)
1046  call mem_deallocate(this%irewet)
1047  call mem_deallocate(this%wetfct)
1048  call mem_deallocate(this%iwetit)
1049  call mem_deallocate(this%ihdwet)
1050  call mem_deallocate(this%ibotnode)
1051  call mem_deallocate(this%iwetdry)
1052  call mem_deallocate(this%iangle1)
1053  call mem_deallocate(this%iangle2)
1054  call mem_deallocate(this%iangle3)
1055  call mem_deallocate(this%nedges)
1056  call mem_deallocate(this%lastedge)
1057  call mem_deallocate(this%ik22overk)
1058  call mem_deallocate(this%ik33overk)
1059  call mem_deallocate(this%intvk)
1060  call mem_deallocate(this%invsc)
1061  call mem_deallocate(this%kchangeper)
1062  call mem_deallocate(this%kchangestp)
1063  !
1064  ! -- Deallocate arrays
1065  deallocate (this%aname)
1066  call mem_deallocate(this%ithickstartflag)
1067  call mem_deallocate(this%icelltype)
1068  call mem_deallocate(this%k11)
1069  call mem_deallocate(this%k22)
1070  call mem_deallocate(this%k33)
1071  call mem_deallocate(this%k11input)
1072  call mem_deallocate(this%k22input)
1073  call mem_deallocate(this%k33input)
1074  call mem_deallocate(this%sat, 'SAT', this%memoryPath)
1075  call mem_deallocate(this%condsat)
1076  call mem_deallocate(this%wetdry)
1077  call mem_deallocate(this%angle1)
1078  call mem_deallocate(this%angle2)
1079  call mem_deallocate(this%angle3)
1080  call mem_deallocate(this%nodedge)
1081  call mem_deallocate(this%ihcedge)
1082  call mem_deallocate(this%propsedge)
1083  call mem_deallocate(this%iedge_ptr)
1084  call mem_deallocate(this%edge_idxs)
1085  call mem_deallocate(this%spdis, 'SPDIS', this%memoryPath)
1086  call mem_deallocate(this%nodekchange)
1087  !
1088  ! -- deallocate parent
1089  call this%NumericalPackageType%da()
1090 
1091  ! pointers
1092  this%hnew => null()
1093 
1094  end subroutine npf_da
1095 
1096  !> @ brief Allocate scalars
1097  !!
1098  !! Allocate and initialize scalars for the VSC package. The base model
1099  !! allocate scalars method is also called.
1100  !<
1101  subroutine allocate_scalars(this)
1102  ! -- modules
1104  ! -- dummy
1105  class(gwfnpftype) :: this
1106  !
1107  ! -- allocate scalars in NumericalPackageType
1108  call this%NumericalPackageType%allocate_scalars()
1109  !
1110  ! -- Allocate scalars
1111  call mem_allocate(this%iname, 'INAME', this%memoryPath)
1112  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
1113  call mem_allocate(this%ixt3drhs, 'IXT3DRHS', this%memoryPath)
1114  call mem_allocate(this%satomega, 'SATOMEGA', this%memoryPath)
1115  call mem_allocate(this%hnoflo, 'HNOFLO', this%memoryPath)
1116  call mem_allocate(this%hdry, 'HDRY', this%memoryPath)
1117  call mem_allocate(this%icellavg, 'ICELLAVG', this%memoryPath)
1118  call mem_allocate(this%iavgkeff, 'IAVGKEFF', this%memoryPath)
1119  call mem_allocate(this%ik22, 'IK22', this%memoryPath)
1120  call mem_allocate(this%ik33, 'IK33', this%memoryPath)
1121  call mem_allocate(this%ik22overk, 'IK22OVERK', this%memoryPath)
1122  call mem_allocate(this%ik33overk, 'IK33OVERK', this%memoryPath)
1123  call mem_allocate(this%iperched, 'IPERCHED', this%memoryPath)
1124  call mem_allocate(this%ivarcv, 'IVARCV', this%memoryPath)
1125  call mem_allocate(this%idewatcv, 'IDEWATCV', this%memoryPath)
1126  call mem_allocate(this%ithickstrt, 'ITHICKSTRT', this%memoryPath)
1127  call mem_allocate(this%ihighcellsat, 'IHIGHCELLSAT', this%memoryPath)
1128  call mem_allocate(this%icalcspdis, 'ICALCSPDIS', this%memoryPath)
1129  call mem_allocate(this%isavspdis, 'ISAVSPDIS', this%memoryPath)
1130  call mem_allocate(this%isavsat, 'ISAVSAT', this%memoryPath)
1131  call mem_allocate(this%irewet, 'IREWET', this%memoryPath)
1132  call mem_allocate(this%wetfct, 'WETFCT', this%memoryPath)
1133  call mem_allocate(this%iwetit, 'IWETIT', this%memoryPath)
1134  call mem_allocate(this%ihdwet, 'IHDWET', this%memoryPath)
1135  call mem_allocate(this%iangle1, 'IANGLE1', this%memoryPath)
1136  call mem_allocate(this%iangle2, 'IANGLE2', this%memoryPath)
1137  call mem_allocate(this%iangle3, 'IANGLE3', this%memoryPath)
1138  call mem_allocate(this%iwetdry, 'IWETDRY', this%memoryPath)
1139  call mem_allocate(this%nedges, 'NEDGES', this%memoryPath)
1140  call mem_allocate(this%lastedge, 'LASTEDGE', this%memoryPath)
1141  call mem_allocate(this%intvk, 'INTVK', this%memoryPath)
1142  call mem_allocate(this%invsc, 'INVSC', this%memoryPath)
1143  call mem_allocate(this%kchangeper, 'KCHANGEPER', this%memoryPath)
1144  call mem_allocate(this%kchangestp, 'KCHANGESTP', this%memoryPath)
1145  !
1146  ! -- set pointer to inewtonur
1147  call mem_setptr(this%igwfnewtonur, 'INEWTONUR', &
1148  create_mem_path(this%name_model))
1149  !
1150  ! -- Initialize value
1151  this%iname = 8
1152  this%ixt3d = 0
1153  this%ixt3drhs = 0
1154  this%satomega = dzero
1155  this%hnoflo = dhnoflo !1.d30
1156  this%hdry = dhdry !-1.d30
1157  this%icellavg = ccond_hmean
1158  this%iavgkeff = 0
1159  this%ik22 = 0
1160  this%ik33 = 0
1161  this%ik22overk = 0
1162  this%ik33overk = 0
1163  this%iperched = 0
1164  this%ivarcv = 0
1165  this%idewatcv = 0
1166  this%ithickstrt = 0
1167  this%ihighcellsat = 0
1168  this%icalcspdis = 0
1169  this%isavspdis = 0
1170  this%isavsat = 0
1171  this%irewet = 0
1172  this%wetfct = done
1173  this%iwetit = 1
1174  this%ihdwet = 0
1175  this%iangle1 = 0
1176  this%iangle2 = 0
1177  this%iangle3 = 0
1178  this%iwetdry = 0
1179  this%nedges = 0
1180  this%lastedge = 0
1181  this%intvk = 0
1182  this%invsc = 0
1183  this%kchangeper = 0
1184  this%kchangestp = 0
1185  !
1186  ! -- If newton is on, then NPF creates asymmetric matrix
1187  this%iasym = this%inewton
1188  end subroutine allocate_scalars
1189 
1190  !> @ brief Store backup copy of hydraulic conductivity when the VSC
1191  !! package is activate
1192  !!
1193  !! The K arrays (K11, etc.) get multiplied by the viscosity ratio so that
1194  !! subsequent uses of K already take into account the effect of viscosity.
1195  !! Thus the original user-specified K array values are lost unless they are
1196  !! backed up in k11input, for example. In a new stress period/time step,
1197  !! the values in k11input are multiplied by the viscosity ratio, not k11
1198  !! since it contains viscosity-adjusted hydraulic conductivity values.
1199  !<
1200  subroutine store_original_k_arrays(this, ncells, njas)
1201  ! -- modules
1203  ! -- dummy
1204  class(gwfnpftype) :: this
1205  integer(I4B), intent(in) :: ncells
1206  integer(I4B), intent(in) :: njas
1207  ! -- local
1208  integer(I4B) :: n
1209  !
1210  ! -- Retain copy of user-specified K arrays
1211  do n = 1, ncells
1212  this%k11input(n) = this%k11(n)
1213  this%k22input(n) = this%k22(n)
1214  this%k33input(n) = this%k33(n)
1215  end do
1216  end subroutine store_original_k_arrays
1217 
1218  !> @brief Allocate npf arrays
1219  !<
1220  subroutine allocate_arrays(this, ncells, njas)
1221  ! -- dummy
1222  class(gwfnpftype) :: this
1223  integer(I4B), intent(in) :: ncells
1224  integer(I4B), intent(in) :: njas
1225  ! -- local
1226  integer(I4B) :: n
1227  !
1228  call mem_allocate(this%ithickstartflag, ncells, 'ITHICKSTARTFLAG', &
1229  this%memoryPath)
1230  call mem_allocate(this%icelltype, ncells, 'ICELLTYPE', this%memoryPath)
1231  call mem_allocate(this%k11, ncells, 'K11', this%memoryPath)
1232  call mem_allocate(this%sat, ncells, 'SAT', this%memoryPath)
1233  call mem_allocate(this%condsat, njas, 'CONDSAT', this%memoryPath)
1234  !
1235  ! -- Optional arrays dimensioned to full size initially
1236  call mem_allocate(this%k22, ncells, 'K22', this%memoryPath)
1237  call mem_allocate(this%k33, ncells, 'K33', this%memoryPath)
1238  call mem_allocate(this%wetdry, ncells, 'WETDRY', this%memoryPath)
1239  call mem_allocate(this%angle1, ncells, 'ANGLE1', this%memoryPath)
1240  call mem_allocate(this%angle2, ncells, 'ANGLE2', this%memoryPath)
1241  call mem_allocate(this%angle3, ncells, 'ANGLE3', this%memoryPath)
1242  !
1243  ! -- Optional arrays
1244  call mem_allocate(this%ibotnode, 0, 'IBOTNODE', this%memoryPath)
1245  call mem_allocate(this%nodedge, 0, 'NODEDGE', this%memoryPath)
1246  call mem_allocate(this%ihcedge, 0, 'IHCEDGE', this%memoryPath)
1247  call mem_allocate(this%propsedge, 0, 0, 'PROPSEDGE', this%memoryPath)
1248  call mem_allocate(this%iedge_ptr, 0, 'NREDGESNODE', this%memoryPath)
1249  call mem_allocate(this%edge_idxs, 0, 'EDGEIDXS', this%memoryPath)
1250  !
1251  ! -- Optional arrays only needed when vsc package is active
1252  call mem_allocate(this%k11input, 0, 'K11INPUT', this%memoryPath)
1253  call mem_allocate(this%k22input, 0, 'K22INPUT', this%memoryPath)
1254  call mem_allocate(this%k33input, 0, 'K33INPUT', this%memoryPath)
1255  !
1256  ! -- Specific discharge is (re-)allocated when nedges is known
1257  call mem_allocate(this%spdis, 3, 0, 'SPDIS', this%memoryPath)
1258  !
1259  ! -- Time-varying property flag arrays
1260  call mem_allocate(this%nodekchange, ncells, 'NODEKCHANGE', this%memoryPath)
1261  !
1262  ! -- initialize iangle1, iangle2, iangle3, and wetdry
1263  do n = 1, ncells
1264  this%angle1(n) = dzero
1265  this%angle2(n) = dzero
1266  this%angle3(n) = dzero
1267  this%wetdry(n) = dzero
1268  this%nodekchange(n) = dzero
1269  end do
1270  !
1271  ! -- allocate variable names
1272  allocate (this%aname(this%iname))
1273  this%aname = [' ICELLTYPE', ' K', &
1274  ' K33', ' K22', &
1275  ' WETDRY', ' ANGLE1', &
1276  ' ANGLE2', ' ANGLE3']
1277  end subroutine allocate_arrays
1278 
1279  !> @brief Log npf options sourced from the input mempath
1280  !<
1281  subroutine log_options(this, found)
1282  ! -- modules
1283  use kindmodule, only: lgp
1285  ! -- dummy
1286  class(gwfnpftype) :: this
1287  ! -- locals
1288  type(gwfnpfparamfoundtype), intent(in) :: found
1289  !
1290  write (this%iout, '(1x,a)') 'Setting NPF Options'
1291  if (found%iprflow) &
1292  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be printed &
1293  &to listing file whenever ICBCFL is not zero.'
1294  if (found%ipakcb) &
1295  write (this%iout, '(4x,a)') 'Cell-by-cell flow information will be saved &
1296  &to binary file whenever ICBCFL is not zero.'
1297  if (found%cellavg) &
1298  write (this%iout, '(4x,a,i0)') 'Alternative cell averaging [1=logarithmic, &
1299  &2=AMT-LMK, 3=AMT-HMK] set to: ', &
1300  this%icellavg
1301  if (found%ithickstrt) &
1302  write (this%iout, '(4x,a)') 'THICKSTRT option has been activated.'
1303  if (found%ihighcellsat) &
1304  write (this%iout, '(4x,a)') 'HIGHEST_CELL_SATURATION option &
1305  &has been activated.'
1306  if (found%iperched) &
1307  write (this%iout, '(4x,a)') 'Vertical flow will be adjusted for perched &
1308  &conditions.'
1309  if (found%ivarcv) &
1310  write (this%iout, '(4x,a)') 'Vertical conductance varies with water table.'
1311  if (found%idewatcv) &
1312  write (this%iout, '(4x,a)') 'Vertical conductance is calculated using &
1313  &only the saturated thickness and properties &
1314  &of the overlying cell if the head in the &
1315  &underlying cell is below its top.'
1316  if (found%ixt3d) write (this%iout, '(4x,a)') 'XT3D formulation is selected.'
1317  if (found%ixt3drhs) &
1318  write (this%iout, '(4x,a)') 'XT3D RHS formulation is selected.'
1319  if (found%isavspdis) &
1320  write (this%iout, '(4x,a)') 'Specific discharge will be calculated at cell &
1321  &centers and written to DATA-SPDIS in budget &
1322  &file when requested.'
1323  if (found%isavsat) &
1324  write (this%iout, '(4x,a)') 'Saturation will be written to DATA-SAT in &
1325  &budget file when requested.'
1326  if (found%ik22overk) &
1327  write (this%iout, '(4x,a)') 'Values specified for K22 are anisotropy &
1328  &ratios and will be multiplied by K before &
1329  &being used in calculations.'
1330  if (found%ik33overk) &
1331  write (this%iout, '(4x,a)') 'Values specified for K33 are anisotropy &
1332  &ratios and will be multiplied by K before &
1333  &being used in calculations.'
1334  if (found%inewton) &
1335  write (this%iout, '(4x,a)') 'NEWTON-RAPHSON method disabled for unconfined &
1336  &cells'
1337  if (found%satomega) &
1338  write (this%iout, '(4x,a,1pg15.6)') 'Saturation omega: ', this%satomega
1339  if (found%irewet) &
1340  write (this%iout, '(4x,a)') 'Rewetting is active.'
1341  if (found%wetfct) &
1342  write (this%iout, '(4x,a,1pg15.6)') &
1343  'Wetting factor (WETFCT) has been set to: ', this%wetfct
1344  if (found%iwetit) &
1345  write (this%iout, '(4x,a,i5)') &
1346  'Wetting iteration interval (IWETIT) has been set to: ', this%iwetit
1347  if (found%ihdwet) &
1348  write (this%iout, '(4x,a,i5)') &
1349  'Head rewet equation (IHDWET) has been set to: ', this%ihdwet
1350  write (this%iout, '(1x,a,/)') 'End Setting NPF Options'
1351  end subroutine log_options
1352 
1353  !> @brief Update simulation options from input mempath
1354  !<
1355  subroutine source_options(this)
1356  ! -- modules
1362  use sourcecommonmodule, only: filein_fname
1363  ! -- dummy
1364  class(gwfnpftype) :: this
1365  ! -- locals
1366  character(len=LENVARNAME), dimension(3) :: cellavg_method = &
1367  &[character(len=LENVARNAME) :: 'LOGARITHMIC', 'AMT-LMK', 'AMT-HMK']
1368  type(gwfnpfparamfoundtype) :: found
1369  character(len=LINELENGTH) :: tvk6_filename
1370  !
1371  ! -- update defaults with idm sourced values
1372  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
1373  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
1374  call mem_set_value(this%icellavg, 'CELLAVG', this%input_mempath, &
1375  cellavg_method, found%cellavg)
1376  call mem_set_value(this%ithickstrt, 'ITHICKSTRT', this%input_mempath, &
1377  found%ithickstrt)
1378  call mem_set_value(this%ihighcellsat, 'IHIGHCELLSAT', this%input_mempath, &
1379  found%ihighcellsat)
1380  call mem_set_value(this%iperched, 'IPERCHED', this%input_mempath, &
1381  found%iperched)
1382  call mem_set_value(this%ivarcv, 'IVARCV', this%input_mempath, found%ivarcv)
1383  call mem_set_value(this%idewatcv, 'IDEWATCV', this%input_mempath, &
1384  found%idewatcv)
1385  call mem_set_value(this%ixt3d, 'IXT3D', this%input_mempath, found%ixt3d)
1386  call mem_set_value(this%ixt3drhs, 'IXT3DRHS', this%input_mempath, &
1387  found%ixt3drhs)
1388  call mem_set_value(this%isavspdis, 'ISAVSPDIS', this%input_mempath, &
1389  found%isavspdis)
1390  call mem_set_value(this%isavsat, 'ISAVSAT', this%input_mempath, found%isavsat)
1391  call mem_set_value(this%ik22overk, 'IK22OVERK', this%input_mempath, &
1392  found%ik22overk)
1393  call mem_set_value(this%ik33overk, 'IK33OVERK', this%input_mempath, &
1394  found%ik33overk)
1395  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
1396  call mem_set_value(this%satomega, 'SATOMEGA', this%input_mempath, &
1397  found%satomega)
1398  call mem_set_value(this%irewet, 'IREWET', this%input_mempath, found%irewet)
1399  call mem_set_value(this%wetfct, 'WETFCT', this%input_mempath, found%wetfct)
1400  call mem_set_value(this%iwetit, 'IWETIT', this%input_mempath, found%iwetit)
1401  call mem_set_value(this%ihdwet, 'IHDWET', this%input_mempath, found%ihdwet)
1402  !
1403  ! -- save flows option active
1404  if (found%ipakcb) this%ipakcb = -1
1405  !
1406  ! -- xt3d active with rhs
1407  if (found%ixt3d .and. found%ixt3drhs) this%ixt3d = 2
1408  !
1409  ! -- save specific discharge active
1410  if (found%isavspdis) this%icalcspdis = this%isavspdis
1411  !
1412  ! -- no newton specified
1413  if (found%inewton) then
1414  this%inewton = 0
1415  this%iasym = 0
1416  end if
1417  !
1418  ! -- enforce 0 or 1 TVK6_FILENAME entries in option block
1419  if (filein_fname(tvk6_filename, 'TVK6_FILENAME', this%input_mempath, &
1420  this%input_fname)) then
1421  call openfile(this%intvk, this%iout, tvk6_filename, 'TVK')
1422  call tvk_cr(this%tvk, this%name_model, this%intvk, this%iout)
1423  end if
1424  !
1425  ! -- verify ALTERNATIVE_CELL_AVERAGING input value is supported
1426  if (found%cellavg) then
1427  if (this%icellavg == 0) then
1428  errmsg = 'Unrecognized input value for ALTERNATIVE_CELL_AVERAGING option.'
1429  call store_error(errmsg)
1430  call store_error_filename(this%input_fname)
1431  end if
1432  end if
1433  !
1434  ! -- log options
1435  if (this%iout > 0) then
1436  call this%log_options(found)
1437  end if
1438  end subroutine source_options
1439 
1440  !> @brief Set options in the NPF object
1441  !<
1442  subroutine set_options(this, options)
1443  ! -- dummy
1444  class(gwfnpftype) :: this
1445  type(gwfnpfoptionstype), intent(in) :: options
1446  !
1447  this%ithickstrt = options%ithickstrt
1448  this%ihighcellsat = options%ihighcellsat
1449  this%iperched = options%iperched
1450  this%ivarcv = options%ivarcv
1451  this%idewatcv = options%idewatcv
1452  this%irewet = options%irewet
1453  this%wetfct = options%wetfct
1454  this%iwetit = options%iwetit
1455  this%ihdwet = options%ihdwet
1456  end subroutine set_options
1457 
1458  !> @brief Check for conflicting NPF options
1459  !<
1460  subroutine check_options(this)
1461  ! -- modules
1462  use simmodule, only: store_error, store_warning, &
1464  use constantsmodule, only: linelength
1465  ! -- dummy
1466  class(gwfnpftype) :: this
1467  !
1468  ! -- set omega value used for saturation calculations
1469  if (this%inewton > 0) then
1470  this%satomega = dem6
1471  end if
1472  !
1473  if (this%inewton > 0) then
1474  if (this%iperched > 0) then
1475  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1476  'BE USED WITH PERCHED OPTION.'
1477  call store_error(errmsg)
1478  end if
1479  if (this%ivarcv > 0) then
1480  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1481  'BE USED WITH VARIABLECV OPTION.'
1482  call store_error(errmsg)
1483  end if
1484  if (this%irewet > 0) then
1485  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. NEWTON OPTION CANNOT '// &
1486  'BE USED WITH REWET OPTION.'
1487  call store_error(errmsg)
1488  end if
1489  else
1490  if (this%ihighcellsat /= 0) then
1491  write (warnmsg, '(a)') 'HIGHEST_CELL_SATURATION '// &
1492  'option cannot be used when NEWTON option in not specified. '// &
1493  'Resetting HIGHEST_CELL_SATURATION option to off.'
1494  this%ihighcellsat = 0
1495  call store_warning(warnmsg)
1496  end if
1497  end if
1498  !
1499  if (this%ixt3d /= 0) then
1500  if (this%icellavg > 0) then
1501  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. '// &
1502  'ALTERNATIVE_CELL_AVERAGING OPTION '// &
1503  'CANNOT BE USED WITH XT3D OPTION.'
1504  call store_error(errmsg)
1505  end if
1506  if (this%ithickstrt > 0) then
1507  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. THICKSTRT OPTION '// &
1508  'CANNOT BE USED WITH XT3D OPTION.'
1509  call store_error(errmsg)
1510  end if
1511  if (this%iperched > 0) then
1512  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. PERCHED OPTION '// &
1513  'CANNOT BE USED WITH XT3D OPTION.'
1514  call store_error(errmsg)
1515  end if
1516  if (this%ivarcv > 0) then
1517  write (errmsg, '(a)') 'ERROR IN NPF OPTIONS. VARIABLECV OPTION '// &
1518  'CANNOT BE USED WITH XT3D OPTION.'
1519  call store_error(errmsg)
1520  end if
1521  end if
1522  !
1523  ! -- Terminate if errors
1524  if (count_errors() > 0) then
1525  call store_error_filename(this%input_fname)
1526  end if
1527  end subroutine check_options
1528 
1529  !> @brief Write dimensions to list file
1530  !<
1531  subroutine log_griddata(this, found)
1532  ! -- modules
1534  ! -- dummy
1535  class(gwfnpftype) :: this
1536  type(gwfnpfparamfoundtype), intent(in) :: found
1537  !
1538  write (this%iout, '(1x,a)') 'Setting NPF Griddata'
1539  !
1540  if (found%icelltype) then
1541  write (this%iout, '(4x,a)') 'ICELLTYPE set from input file'
1542  end if
1543  !
1544  if (found%k) then
1545  write (this%iout, '(4x,a)') 'K set from input file'
1546  end if
1547  !
1548  if (found%k33) then
1549  write (this%iout, '(4x,a)') 'K33 set from input file'
1550  else
1551  write (this%iout, '(4x,a)') 'K33 not provided. Setting K33 = K.'
1552  end if
1553  !
1554  if (found%k22) then
1555  write (this%iout, '(4x,a)') 'K22 set from input file'
1556  else
1557  write (this%iout, '(4x,a)') 'K22 not provided. Setting K22 = K.'
1558  end if
1559  !
1560  if (found%wetdry) then
1561  write (this%iout, '(4x,a)') 'WETDRY set from input file'
1562  end if
1563  !
1564  if (found%angle1) then
1565  write (this%iout, '(4x,a)') 'ANGLE1 set from input file'
1566  end if
1567  !
1568  if (found%angle2) then
1569  write (this%iout, '(4x,a)') 'ANGLE2 set from input file'
1570  end if
1571  !
1572  if (found%angle3) then
1573  write (this%iout, '(4x,a)') 'ANGLE3 set from input file'
1574  end if
1575  !
1576  write (this%iout, '(1x,a,/)') 'End Setting NPF Griddata'
1577  end subroutine log_griddata
1578 
1579  !> @brief Update simulation griddata from input mempath
1580  !<
1581  subroutine source_griddata(this)
1582  ! -- modules
1583  use simmodule, only: count_errors, store_error
1587  ! -- dummy
1588  class(gwfnpftype) :: this
1589  ! -- locals
1590  character(len=LINELENGTH) :: errmsg
1591  type(gwfnpfparamfoundtype) :: found
1592  logical, dimension(2) :: afound
1593  integer(I4B), dimension(:), pointer, contiguous :: map
1594  !
1595  ! -- set map to convert user input data into reduced data
1596  map => null()
1597  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1598  !
1599  ! -- update defaults with idm sourced values
1600  call mem_set_value(this%icelltype, 'ICELLTYPE', this%input_mempath, map, &
1601  found%icelltype)
1602  call mem_set_value(this%k11, 'K', this%input_mempath, map, found%k)
1603  call mem_set_value(this%k33, 'K33', this%input_mempath, map, found%k33)
1604  call mem_set_value(this%k22, 'K22', this%input_mempath, map, found%k22)
1605  call mem_set_value(this%wetdry, 'WETDRY', this%input_mempath, map, &
1606  found%wetdry)
1607  call mem_set_value(this%angle1, 'ANGLE1', this%input_mempath, map, &
1608  found%angle1)
1609  call mem_set_value(this%angle2, 'ANGLE2', this%input_mempath, map, &
1610  found%angle2)
1611  call mem_set_value(this%angle3, 'ANGLE3', this%input_mempath, map, &
1612  found%angle3)
1613  !
1614  ! -- ensure ICELLTYPE was found
1615  if (.not. found%icelltype) then
1616  write (errmsg, '(a)') 'Error in GRIDDATA block: ICELLTYPE not found.'
1617  call store_error(errmsg)
1618  end if
1619  !
1620  ! -- ensure K was found
1621  if (.not. found%k) then
1622  write (errmsg, '(a)') 'Error in GRIDDATA block: K not found.'
1623  call store_error(errmsg)
1624  end if
1625  !
1626  ! -- set error if ik33overk set with no k33
1627  if (.not. found%k33 .and. this%ik33overk /= 0) then
1628  write (errmsg, '(a)') 'K33OVERK option specified but K33 not specified.'
1629  call store_error(errmsg)
1630  end if
1631  !
1632  ! -- set error if ik22overk set with no k22
1633  if (.not. found%k22 .and. this%ik22overk /= 0) then
1634  write (errmsg, '(a)') 'K22OVERK option specified but K22 not specified.'
1635  call store_error(errmsg)
1636  end if
1637  !
1638  ! -- handle found side effects
1639  if (found%k33) this%ik33 = 1
1640  if (found%k22) this%ik22 = 1
1641  if (found%wetdry) this%iwetdry = 1
1642  if (found%angle1) this%iangle1 = 1
1643  if (found%angle2) this%iangle2 = 1
1644  if (found%angle3) this%iangle3 = 1
1645  !
1646  ! -- handle not found side effects
1647  if (.not. found%k33) then
1648  call mem_set_value(this%k33, 'K', this%input_mempath, map, afound(1))
1649  end if
1650  if (.not. found%k22) then
1651  call mem_set_value(this%k22, 'K', this%input_mempath, map, afound(2))
1652  end if
1653  if (.not. found%wetdry) call mem_reallocate(this%wetdry, 1, 'WETDRY', &
1654  trim(this%memoryPath))
1655  if (.not. found%angle1 .and. this%ixt3d == 0) &
1656  call mem_reallocate(this%angle1, 0, 'ANGLE1', trim(this%memoryPath))
1657  if (.not. found%angle2 .and. this%ixt3d == 0) &
1658  call mem_reallocate(this%angle2, 0, 'ANGLE2', trim(this%memoryPath))
1659  if (.not. found%angle3 .and. this%ixt3d == 0) &
1660  call mem_reallocate(this%angle3, 0, 'ANGLE3', trim(this%memoryPath))
1661  !
1662  ! -- log griddata
1663  if (this%iout > 0) then
1664  call this%log_griddata(found)
1665  end if
1666  end subroutine source_griddata
1667 
1668  !> @brief Initialize and check NPF data
1669  !<
1670  subroutine prepcheck(this)
1671  ! -- modules
1672  use constantsmodule, only: linelength, dpio180
1674  ! -- dummy
1675  class(gwfnpftype) :: this
1676  ! -- local
1677  character(len=24), dimension(:), pointer :: aname
1678  character(len=LINELENGTH) :: cellstr, errmsg
1679  integer(I4B) :: nerr, n
1680  ! -- format
1681  character(len=*), parameter :: fmtkerr = &
1682  &"(1x, 'Hydraulic property ',a,' is <= 0 for cell ',a, ' ', 1pg15.6)"
1683  character(len=*), parameter :: fmtkerr2 = &
1684  &"(1x, '... ', i0,' additional errors not shown for ',a)"
1685  !
1686  ! -- initialize
1687  aname => this%aname
1688  !
1689  ! -- check k11
1690  nerr = 0
1691  do n = 1, size(this%k11)
1692  if (this%k11(n) <= dzero) then
1693  nerr = nerr + 1
1694  if (nerr <= 20) then
1695  call this%dis%noder_to_string(n, cellstr)
1696  write (errmsg, fmtkerr) trim(adjustl(aname(2))), trim(cellstr), &
1697  this%k11(n)
1698  call store_error(errmsg)
1699  end if
1700  end if
1701  end do
1702  if (nerr > 20) then
1703  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(2)))
1704  call store_error(errmsg)
1705  end if
1706  !
1707  ! -- check k33 because it was read
1708  if (this%ik33 /= 0) then
1709  !
1710  ! -- Check to make sure values are greater than or equal to zero
1711  nerr = 0
1712  do n = 1, size(this%k33)
1713  if (this%ik33overk /= 0) this%k33(n) = this%k33(n) * this%k11(n)
1714  if (this%k33(n) <= dzero) then
1715  nerr = nerr + 1
1716  if (nerr <= 20) then
1717  call this%dis%noder_to_string(n, cellstr)
1718  write (errmsg, fmtkerr) trim(adjustl(aname(3))), trim(cellstr), &
1719  this%k33(n)
1720  call store_error(errmsg)
1721  end if
1722  end if
1723  end do
1724  if (nerr > 20) then
1725  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(3)))
1726  call store_error(errmsg)
1727  end if
1728  end if
1729  !
1730  ! -- check k22 because it was read
1731  if (this%ik22 /= 0) then
1732  !
1733  ! -- Check to make sure that angles are available
1734  if (this%dis%con%ianglex == 0) then
1735  write (errmsg, '(a)') 'Error. ANGLDEGX not provided in '// &
1736  'discretization file, but K22 was specified. '
1737  call store_error(errmsg)
1738  end if
1739  !
1740  ! -- Check to make sure values are greater than or equal to zero
1741  nerr = 0
1742  do n = 1, size(this%k22)
1743  if (this%ik22overk /= 0) this%k22(n) = this%k22(n) * this%k11(n)
1744  if (this%k22(n) <= dzero) then
1745  nerr = nerr + 1
1746  if (nerr <= 20) then
1747  call this%dis%noder_to_string(n, cellstr)
1748  write (errmsg, fmtkerr) trim(adjustl(aname(4))), trim(cellstr), &
1749  this%k22(n)
1750  call store_error(errmsg)
1751  end if
1752  end if
1753  end do
1754  if (nerr > 20) then
1755  write (errmsg, fmtkerr2) nerr, trim(adjustl(aname(4)))
1756  call store_error(errmsg)
1757  end if
1758  end if
1759  !
1760  ! -- check for wetdry conflicts
1761  if (this%irewet == 1) then
1762  if (this%iwetdry == 0) then
1763  write (errmsg, '(a, a, a)') 'Error in GRIDDATA block: ', &
1764  trim(adjustl(aname(5))), ' not found.'
1765  call store_error(errmsg)
1766  end if
1767  end if
1768  !
1769  ! -- Check for angle conflicts
1770  if (this%iangle1 /= 0) then
1771  do n = 1, size(this%angle1)
1772  this%angle1(n) = this%angle1(n) * dpio180
1773  end do
1774  else
1775  if (this%ixt3d /= 0) then
1776  this%iangle1 = 1
1777  write (this%iout, '(a)') 'XT3D IN USE, BUT ANGLE1 NOT SPECIFIED. '// &
1778  'SETTING ANGLE1 TO ZERO.'
1779  do n = 1, size(this%angle1)
1780  this%angle1(n) = dzero
1781  end do
1782  end if
1783  end if
1784  if (this%iangle2 /= 0) then
1785  if (this%iangle1 == 0) then
1786  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE1. '// &
1787  'ANGLE2 REQUIRES ANGLE1. '
1788  call store_error(errmsg)
1789  end if
1790  if (this%iangle3 == 0) then
1791  write (errmsg, '(a)') 'ANGLE2 SPECIFIED BUT NOT ANGLE3. '// &
1792  'SPECIFY BOTH OR NEITHER ONE. '
1793  call store_error(errmsg)
1794  end if
1795  do n = 1, size(this%angle2)
1796  this%angle2(n) = this%angle2(n) * dpio180
1797  end do
1798  end if
1799  if (this%iangle3 /= 0) then
1800  if (this%iangle1 == 0) then
1801  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE1. '// &
1802  'ANGLE3 REQUIRES ANGLE1. '
1803  call store_error(errmsg)
1804  end if
1805  if (this%iangle2 == 0) then
1806  write (errmsg, '(a)') 'ANGLE3 SPECIFIED BUT NOT ANGLE2. '// &
1807  'SPECIFY BOTH OR NEITHER ONE. '
1808  call store_error(errmsg)
1809  end if
1810  do n = 1, size(this%angle3)
1811  this%angle3(n) = this%angle3(n) * dpio180
1812  end do
1813  end if
1814  !
1815  ! -- terminate if data errors
1816  if (count_errors() > 0) then
1817  call store_error_filename(this%input_fname)
1818  end if
1819  end subroutine prepcheck
1820 
1821  !> @brief preprocess the NPF input data
1822  !!
1823  !! This routine consists of the following steps:
1824  !!
1825  !! 1. convert cells to noflow when all transmissive parameters equal zero
1826  !! 2. perform initial wetting and drying
1827  !! 3. initialize cell saturation
1828  !! 4. calculate saturated conductance (when not xt3d)
1829  !! 5. If NEWTON under-relaxation, determine lower most node
1830  !<
1831  subroutine preprocess_input(this)
1832  ! -- modules
1833  use constantsmodule, only: linelength
1835  ! -- dummy
1836  class(gwfnpftype) :: this !< the instance of the NPF package
1837  ! -- local
1838  integer(I4B) :: n, m, ii, nn
1839  real(DP) :: hyn, hym
1840  real(DP) :: satn, topn, botn
1841  integer(I4B) :: nextn
1842  real(DP) :: minbot, botm
1843  logical :: finished
1844  character(len=LINELENGTH) :: cellstr, errmsg
1845  ! -- format
1846  character(len=*), parameter :: fmtcnv = &
1847  "(1X,'CELL ', A, &
1848  &' ELIMINATED BECAUSE ALL HYDRAULIC CONDUCTIVITIES TO NODE ARE 0.')"
1849  character(len=*), parameter :: fmtnct = &
1850  &"(1X,'Negative cell thickness at cell ', A)"
1851  character(len=*), parameter :: fmtihbe = &
1852  &"(1X,'Initial head, bottom elevation:',1P,2G13.5)"
1853  character(len=*), parameter :: fmttebe = &
1854  &"(1X,'Top elevation, bottom elevation:',1P,2G13.5)"
1855  !
1856  do n = 1, this%dis%nodes
1857  this%ithickstartflag(n) = 0
1858  end do
1859  !
1860  ! -- Insure that each cell has at least one non-zero transmissive parameter
1861  ! Note that a cell can be deactivated even if it has a valid connection
1862  ! to another model.
1863  nodeloop: do n = 1, this%dis%nodes
1864  !
1865  ! -- Skip if already inactive
1866  if (this%ibound(n) == 0) then
1867  if (this%irewet /= 0) then
1868  if (this%wetdry(n) == dzero) cycle nodeloop
1869  else
1870  cycle nodeloop
1871  end if
1872  end if
1873  !
1874  ! -- Cycle if k11 is not zero
1875  if (this%k11(n) /= dzero) cycle nodeloop
1876  !
1877  ! -- Cycle if at least one vertical connection has non-zero k33
1878  ! for n and m
1879  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
1880  m = this%dis%con%ja(ii)
1881  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0) then
1882  hyn = this%k11(n)
1883  if (this%ik33 /= 0) hyn = this%k33(n)
1884  if (hyn /= dzero) then
1885  hym = this%k11(m)
1886  if (this%ik33 /= 0) hym = this%k33(m)
1887  if (hym /= dzero) cycle
1888  end if
1889  end if
1890  end do
1891  !
1892  ! -- If this part of the loop is reached, then all connections have
1893  ! zero transmissivity, so convert to noflow.
1894  this%ibound(n) = 0
1895  this%hnew(n) = this%hnoflo
1896  if (this%irewet /= 0) this%wetdry(n) = dzero
1897  call this%dis%noder_to_string(n, cellstr)
1898  write (this%iout, fmtcnv) trim(adjustl(cellstr))
1899  !
1900  end do nodeloop
1901  !
1902  ! -- Preprocess cell status and heads based on initial conditions
1903  if (this%inewton == 0) then
1904  !
1905  ! -- For standard formulation (non-Newton) call wetdry routine
1906  call this%wd(0, this%hnew)
1907  else
1908  !
1909  ! -- Newton formulation, so adjust heads to be above bottom
1910  ! (Not used in present formulation because variable cv
1911  ! cannot be used with Newton)
1912  if (this%ivarcv == 1) then
1913  do n = 1, this%dis%nodes
1914  if (this%hnew(n) < this%dis%bot(n)) then
1915  this%hnew(n) = this%dis%bot(n) + dem6
1916  end if
1917  end do
1918  end if
1919  end if
1920  !
1921  ! -- If THCKSTRT is not active, then loop through icelltype and replace
1922  ! any negative values with 1.
1923  if (this%ithickstrt == 0) then
1924  do n = 1, this%dis%nodes
1925  if (this%icelltype(n) < 0) then
1926  this%icelltype(n) = 1
1927  end if
1928  end do
1929  end if
1930  !
1931  ! -- Initialize sat to zero for ibound=0 cells, unless the cell can
1932  ! rewet. Initialize sat to the saturated fraction based on strt
1933  ! if icelltype is negative and the THCKSTRT option is in effect.
1934  ! Initialize sat to 1.0 for all other cells in order to calculate
1935  ! condsat in next section.
1936  do n = 1, this%dis%nodes
1937  if (this%ibound(n) == 0) then
1938  this%sat(n) = done
1939  if (this%icelltype(n) < 0 .and. this%ithickstrt /= 0) then
1940  this%ithickstartflag(n) = 1
1941  this%icelltype(n) = 0
1942  end if
1943  else
1944  topn = this%dis%top(n)
1945  botn = this%dis%bot(n)
1946  if (this%icelltype(n) < 0 .and. this%ithickstrt /= 0) then
1947  call this%thksat(n, this%ic%strt(n), satn)
1948  if (botn > this%ic%strt(n)) then
1949  call this%dis%noder_to_string(n, cellstr)
1950  write (errmsg, fmtnct) trim(adjustl(cellstr))
1951  call store_error(errmsg)
1952  write (errmsg, fmtihbe) this%ic%strt(n), botn
1953  call store_error(errmsg)
1954  end if
1955  this%ithickstartflag(n) = 1
1956  this%icelltype(n) = 0
1957  else
1958  satn = done
1959  if (botn > topn) then
1960  call this%dis%noder_to_string(n, cellstr)
1961  write (errmsg, fmtnct) trim(adjustl(cellstr))
1962  call store_error(errmsg)
1963  write (errmsg, fmttebe) topn, botn
1964  call store_error(errmsg)
1965  end if
1966  end if
1967  this%sat(n) = satn
1968  end if
1969  end do
1970  if (count_errors() > 0) then
1971  call store_error_filename(this%input_fname)
1972  end if
1973  !
1974  ! -- Calculate condsat, but only if xt3d is not active. If xt3d is
1975  ! active, then condsat is allocated to size of zero.
1976  if (this%ixt3d == 0) then
1977  !
1978  ! -- Calculate the saturated conductance for all connections assuming
1979  ! that saturation is 1 (except for case where icelltype was entered
1980  ! as a negative value and THCKSTRT option in effect)
1981  do n = 1, this%dis%nodes
1982  call this%calc_condsat(n, .true.)
1983  end do
1984  !
1985  end if
1986  !
1987  ! -- Determine the lower most node
1988  if (this%igwfnewtonur /= 0) then
1989  call mem_reallocate(this%ibotnode, this%dis%nodes, 'IBOTNODE', &
1990  trim(this%memoryPath))
1991  do n = 1, this%dis%nodes
1992  !
1993  minbot = this%dis%bot(n)
1994  nn = n
1995  finished = .false.
1996  do while (.not. finished)
1997  nextn = 0
1998  !
1999  ! -- Go through the connecting cells
2000  do ii = this%dis%con%ia(nn) + 1, this%dis%con%ia(nn + 1) - 1
2001  !
2002  ! -- Set the m cell number
2003  m = this%dis%con%ja(ii)
2004  botm = this%dis%bot(m)
2005  !
2006  ! -- select vertical connections: ihc == 0
2007  if (this%dis%con%ihc(this%dis%con%jas(ii)) == 0) then
2008  if (m > nn .and. botm < minbot) then
2009  nextn = m
2010  minbot = botm
2011  end if
2012  end if
2013  end do
2014  if (nextn > 0) then
2015  nn = nextn
2016  else
2017  finished = .true.
2018  end if
2019  end do
2020  this%ibotnode(n) = nn
2021  end do
2022  end if
2023  !
2024  ! -- nullify unneeded gwf pointers
2025  this%igwfnewtonur => null()
2026  end subroutine preprocess_input
2027 
2028  !> @brief Calculate CONDSAT array entries for the given node
2029  !!
2030  !! Calculate saturated conductances for all connections of the given node,
2031  !! or optionally for the upper portion of the matrix only.
2032  !<
2033  subroutine calc_condsat(this, node, upperOnly)
2034  ! -- dummy variables
2035  class(gwfnpftype) :: this
2036  integer(I4B), intent(in) :: node
2037  logical, intent(in) :: upperOnly
2038  ! -- local variables
2039  integer(I4B) :: ii, m, n, ihc, jj
2040  real(DP) :: topm, topn, topnode, botm, botn, botnode, satm, satn, satnode
2041  real(DP) :: hyn, hym, hn, hm, fawidth, csat
2042  !
2043  satnode = this%calc_initial_sat(node)
2044  !
2045  topnode = this%dis%top(node)
2046  botnode = this%dis%bot(node)
2047  !
2048  ! -- Go through the connecting cells
2049  do ii = this%dis%con%ia(node) + 1, this%dis%con%ia(node + 1) - 1
2050  !
2051  ! -- Set the m cell number and cycle if lower triangle connection and
2052  ! -- we're not updating both upper and lower matrix parts for this node
2053  m = this%dis%con%ja(ii)
2054  jj = this%dis%con%jas(ii)
2055  if (m < node) then
2056  if (upperonly) cycle
2057  ! m => node, n => neighbour
2058  n = m
2059  m = node
2060  topm = topnode
2061  botm = botnode
2062  satm = satnode
2063  topn = this%dis%top(n)
2064  botn = this%dis%bot(n)
2065  satn = this%calc_initial_sat(n)
2066  else
2067  ! n => node, m => neighbour
2068  n = node
2069  topn = topnode
2070  botn = botnode
2071  satn = satnode
2072  topm = this%dis%top(m)
2073  botm = this%dis%bot(m)
2074  satm = this%calc_initial_sat(m)
2075  end if
2076  !
2077  ihc = this%dis%con%ihc(jj)
2078  hyn = this%hy_eff(n, m, ihc, ipos=ii)
2079  hym = this%hy_eff(m, n, ihc, ipos=ii)
2080  if (this%ithickstartflag(n) == 0) then
2081  hn = topn
2082  else
2083  hn = this%ic%strt(n)
2084  end if
2085  if (this%ithickstartflag(m) == 0) then
2086  hm = topm
2087  else
2088  hm = this%ic%strt(m)
2089  end if
2090  !
2091  ! -- Calculate conductance depending on whether connection is
2092  ! vertical (0), horizontal (1), or staggered horizontal (2)
2093  if (ihc == c3d_vertical) then
2094  !
2095  ! -- Vertical conductance for fully saturated conditions
2096  csat = vcond(1, 1, 1, 1, 0, 1, 1, done, &
2097  botn, botm, &
2098  hyn, hym, &
2099  satn, satm, &
2100  topn, topm, &
2101  botn, botm, &
2102  this%dis%con%hwva(jj))
2103  else
2104  !
2105  ! -- Horizontal conductance for fully saturated conditions
2106  fawidth = this%dis%con%hwva(jj)
2107  csat = hcond(1, 1, 1, 1, 0, &
2108  ihc, &
2109  this%icellavg, &
2110  done, &
2111  hn, hm, satn, satm, hyn, hym, &
2112  topn, topm, &
2113  botn, botm, &
2114  this%dis%con%cl1(jj), &
2115  this%dis%con%cl2(jj), &
2116  fawidth)
2117  end if
2118  this%condsat(jj) = csat
2119  end do
2120  end subroutine calc_condsat
2121 
2122  !> @brief Calculate initial saturation for the given node
2123  !!
2124  !! Calculate saturation as a fraction of thickness for the given node, used
2125  !! for saturated conductance calculations: full thickness by default (1.0) or
2126  !! saturation based on initial conditions if the THICKSTRT option is used.
2127  !!
2128  !<
2129  function calc_initial_sat(this, n) result(satn)
2130  ! -- dummy variables
2131  class(gwfnpftype) :: this
2132  integer(I4B), intent(in) :: n
2133  ! -- Return
2134  real(dp) :: satn
2135  !
2136  satn = done
2137  if (this%ibound(n) /= 0 .and. this%ithickstartflag(n) /= 0) then
2138  call this%thksat(n, this%ic%strt(n), satn)
2139  end if
2140  end function calc_initial_sat
2141 
2142  !> @brief Perform wetting and drying
2143  !<
2144  subroutine sgwf_npf_wetdry(this, kiter, hnew)
2145  ! -- modules
2146  use tdismodule, only: kstp, kper
2148  use constantsmodule, only: linelength
2149  ! -- dummy
2150  class(gwfnpftype) :: this
2151  integer(I4B), intent(in) :: kiter
2152  real(DP), intent(inout), dimension(:) :: hnew
2153  ! -- local
2154  integer(I4B) :: n, m, ii, ihc
2155  real(DP) :: ttop, bbot, thick
2156  integer(I4B) :: ncnvrt, ihdcnv
2157  character(len=30), dimension(5) :: nodcnvrt
2158  character(len=30) :: nodestr
2159  character(len=3), dimension(5) :: acnvrt
2160  character(len=LINELENGTH) :: errmsg
2161  integer(I4B) :: irewet
2162  ! -- formats
2163  character(len=*), parameter :: fmtnct = &
2164  "(1X,/1X,'Negative cell thickness at (layer,row,col)', &
2165  &I4,',',I5,',',I5)"
2166  character(len=*), parameter :: fmttopbot = &
2167  &"(1X,'Top elevation, bottom elevation:',1P,2G13.5)"
2168  character(len=*), parameter :: fmttopbotthk = &
2169  &"(1X,'Top elevation, bottom elevation, thickness:',1P,3G13.5)"
2170  character(len=*), parameter :: fmtdrychd = &
2171  &"(1X,/1X,'CONSTANT-HEAD CELL WENT DRY -- SIMULATION ABORTED')"
2172  character(len=*), parameter :: fmtni = &
2173  &"(1X,'CELLID=',a,' ITERATION=',I0,' TIME STEP=',I0,' STRESS PERIOD=',I0)"
2174  !
2175  ! -- Initialize
2176  ncnvrt = 0
2177  ihdcnv = 0
2178  !
2179  ! -- Convert dry cells to wet
2180  do n = 1, this%dis%nodes
2181  do ii = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
2182  m = this%dis%con%ja(ii)
2183  ihc = this%dis%con%ihc(this%dis%con%jas(ii))
2184  call this%rewet_check(kiter, n, hnew(m), this%ibound(m), ihc, hnew, &
2185  irewet)
2186  if (irewet == 1) then
2187  call this%wdmsg(2, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2188  end if
2189  end do
2190  end do
2191  !
2192  ! -- Perform drying
2193  do n = 1, this%dis%nodes
2194  !
2195  ! -- cycle if inactive or confined
2196  if (this%ibound(n) == 0) cycle
2197  if (this%icelltype(n) == 0) cycle
2198  !
2199  ! -- check for negative cell thickness
2200  bbot = this%dis%bot(n)
2201  ttop = this%dis%top(n)
2202  if (bbot > ttop) then
2203  write (errmsg, fmtnct) n
2204  call store_error(errmsg)
2205  write (errmsg, fmttopbot) ttop, bbot
2206  call store_error(errmsg)
2207  call store_error_filename(this%input_fname)
2208  end if
2209  !
2210  ! -- Calculate saturated thickness
2211  if (this%icelltype(n) /= 0) then
2212  if (hnew(n) < ttop) ttop = hnew(n)
2213  end if
2214  thick = ttop - bbot
2215  !
2216  ! -- If thick<0 print message, set hnew, and ibound
2217  if (thick <= dzero) then
2218  call this%wdmsg(1, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2219  hnew(n) = this%hdry
2220  if (this%ibound(n) < 0) then
2221  errmsg = 'CONSTANT-HEAD CELL WENT DRY -- SIMULATION ABORTED'
2222  call store_error(errmsg)
2223  write (errmsg, fmttopbotthk) ttop, bbot, thick
2224  call store_error(errmsg)
2225  call this%dis%noder_to_string(n, nodestr)
2226  write (errmsg, fmtni) trim(adjustl(nodestr)), kiter, kstp, kper
2227  call store_error(errmsg)
2228  call store_error_filename(this%input_fname)
2229  end if
2230  this%ibound(n) = 0
2231  end if
2232  end do
2233  !
2234  ! -- Print remaining cell conversions
2235  call this%wdmsg(0, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
2236  !
2237  ! -- Change ibound from 30000 to 1
2238  do n = 1, this%dis%nodes
2239  if (this%ibound(n) == 30000) this%ibound(n) = 1
2240  end do
2241  end subroutine sgwf_npf_wetdry
2242 
2243  !> @brief Determine if a cell should rewet
2244  !!
2245  !! This method can be called from any external object that has a head that
2246  !! can be used to rewet the GWF cell node. The ihc value is used to
2247  !! determine if it is a vertical or horizontal connection, which can operate
2248  !! differently depending on user settings.
2249  !<
2250  subroutine rewet_check(this, kiter, node, hm, ibdm, ihc, hnew, irewet)
2251  ! -- dummy
2252  class(gwfnpftype) :: this
2253  integer(I4B), intent(in) :: kiter
2254  integer(I4B), intent(in) :: node
2255  real(DP), intent(in) :: hm
2256  integer(I4B), intent(in) :: ibdm
2257  integer(I4B), intent(in) :: ihc
2258  real(DP), intent(inout), dimension(:) :: hnew
2259  integer(I4B), intent(out) :: irewet
2260  ! -- local
2261  integer(I4B) :: itflg
2262  real(DP) :: wd, awd, turnon, bbot
2263  !
2264  irewet = 0
2265  !
2266  ! -- Convert a dry cell to wet if it meets the criteria
2267  if (this%irewet > 0) then
2268  itflg = mod(kiter, this%iwetit)
2269  if (itflg == 0) then
2270  if (this%ibound(node) == 0 .and. this%wetdry(node) /= dzero) then
2271  !
2272  ! -- Calculate wetting elevation
2273  bbot = this%dis%bot(node)
2274  wd = this%wetdry(node)
2275  awd = wd
2276  if (wd < 0) awd = -wd
2277  turnon = bbot + awd
2278  !
2279  ! -- Check head in adjacent cells to see if wetting elevation has
2280  ! been reached
2281  if (ihc == c3d_vertical) then
2282  !
2283  ! -- check cell below
2284  if (ibdm > 0 .and. hm >= turnon) irewet = 1
2285  else
2286  if (wd > dzero) then
2287  !
2288  ! -- check horizontally adjacent cells
2289  if (ibdm > 0 .and. hm >= turnon) irewet = 1
2290  end if
2291  end if
2292  !
2293  if (irewet == 1) then
2294  ! -- rewet cell; use equation 3a if ihdwet=0; use equation 3b if
2295  ! ihdwet is not 0.
2296  if (this%ihdwet == 0) then
2297  hnew(node) = bbot + this%wetfct * (hm - bbot)
2298  else
2299  hnew(node) = bbot + this%wetfct * awd !(hm - bbot)
2300  end if
2301  this%ibound(node) = 30000
2302  end if
2303  end if
2304  end if
2305  end if
2306  end subroutine rewet_check
2307 
2308  !> @brief Print wet/dry message
2309  !<
2310  subroutine sgwf_npf_wdmsg(this, icode, ncnvrt, nodcnvrt, acnvrt, ihdcnv, &
2311  kiter, n)
2312  ! -- modules
2313  use tdismodule, only: kstp, kper
2314  ! -- dummy
2315  class(gwfnpftype) :: this
2316  integer(I4B), intent(in) :: icode
2317  integer(I4B), intent(inout) :: ncnvrt
2318  character(len=30), dimension(5), intent(inout) :: nodcnvrt
2319  character(len=3), dimension(5), intent(inout) :: acnvrt
2320  integer(I4B), intent(inout) :: ihdcnv
2321  integer(I4B), intent(in) :: kiter
2322  integer(I4B), intent(in) :: n
2323  ! -- local
2324  integer(I4B) :: l
2325  ! -- formats
2326  character(len=*), parameter :: fmtcnvtn = &
2327  "(1X,/1X,'CELL CONVERSIONS FOR ITER.=',I0, &
2328  &' STEP=',I0,' PERIOD=',I0,' (NODE or LRC)')"
2329  character(len=*), parameter :: fmtnode = "(1X,3X,5(A4, A20))"
2330  !
2331  ! -- Keep track of cell conversions
2332  if (icode > 0) then
2333  ncnvrt = ncnvrt + 1
2334  call this%dis%noder_to_string(n, nodcnvrt(ncnvrt))
2335  if (icode == 1) then
2336  acnvrt(ncnvrt) = 'DRY'
2337  else
2338  acnvrt(ncnvrt) = 'WET'
2339  end if
2340  end if
2341  !
2342  ! -- Print a line if 5 conversions have occurred or if icode indicates that a
2343  ! partial line should be printed
2344  if (ncnvrt == 5 .or. (icode == 0 .and. ncnvrt > 0)) then
2345  if (ihdcnv == 0) write (this%iout, fmtcnvtn) kiter, kstp, kper
2346  ihdcnv = 1
2347  write (this%iout, fmtnode) &
2348  (acnvrt(l), trim(adjustl(nodcnvrt(l))), l=1, ncnvrt)
2349  ncnvrt = 0
2350  end if
2351  end subroutine sgwf_npf_wdmsg
2352 
2353  !> @brief Calculate the effective hydraulic conductivity for the n-m connection
2354  !!
2355  !! n is primary node node number
2356  !! m is connected node (not used if vg is provided)
2357  !! ihc is horizontal indicator (0 vertical, 1 horizontal, 2 vertically
2358  !! staggered)
2359  !! ipos_opt is position of connection in ja array
2360  !! vg is the global unit vector that expresses the direction from which to
2361  !! calculate an effective hydraulic conductivity.
2362  !<
2363  function hy_eff(this, n, m, ihc, ipos, vg) result(hy)
2364  ! -- return
2365  real(dp) :: hy
2366  ! -- dummy
2367  class(gwfnpftype) :: this
2368  integer(I4B), intent(in) :: n
2369  integer(I4B), intent(in) :: m
2370  integer(I4B), intent(in) :: ihc
2371  integer(I4B), intent(in), optional :: ipos
2372  real(dp), dimension(3), intent(in), optional :: vg
2373  ! -- local
2374  integer(I4B) :: iipos
2375  real(dp) :: hy11, hy22, hy33
2376  real(dp) :: ang1, ang2, ang3
2377  real(dp) :: vg1, vg2, vg3
2378  !
2379  ! -- Initialize
2380  iipos = 0
2381  if (present(ipos)) iipos = ipos
2382  hy11 = this%k11(n)
2383  hy22 = this%k11(n)
2384  hy33 = this%k11(n)
2385  hy22 = this%k22(n)
2386  hy33 = this%k33(n)
2387  !
2388  ! -- Calculate effective K based on whether connection is vertical
2389  ! or horizontal
2390  if (ihc == c3d_vertical) then
2391  !
2392  ! -- Handle rotated anisotropy case that would affect the effective
2393  ! vertical hydraulic conductivity
2394  hy = hy33
2395  if (this%iangle2 > 0) then
2396  if (present(vg)) then
2397  vg1 = vg(1)
2398  vg2 = vg(2)
2399  vg3 = vg(3)
2400  else
2401  call this%dis%connection_normal(n, m, ihc, vg1, vg2, vg3, iipos)
2402  end if
2403  ang1 = this%angle1(n)
2404  ang2 = this%angle2(n)
2405  ang3 = dzero
2406  if (this%iangle3 > 0) ang3 = this%angle3(n)
2407  hy = hyeff(hy11, hy22, hy33, ang1, ang2, ang3, vg1, vg2, vg3, &
2408  this%iavgkeff)
2409  end if
2410  !
2411  else
2412  !
2413  ! -- Handle horizontal case
2414  hy = hy11
2415  if (this%ik22 > 0) then
2416  if (present(vg)) then
2417  vg1 = vg(1)
2418  vg2 = vg(2)
2419  vg3 = vg(3)
2420  else
2421  call this%dis%connection_normal(n, m, ihc, vg1, vg2, vg3, iipos)
2422  end if
2423  ang1 = dzero
2424  ang2 = dzero
2425  ang3 = dzero
2426  if (this%iangle1 > 0) then
2427  ang1 = this%angle1(n)
2428  if (this%iangle2 > 0) then
2429  ang2 = this%angle2(n)
2430  if (this%iangle3 > 0) ang3 = this%angle3(n)
2431  end if
2432  end if
2433  hy = hyeff(hy11, hy22, hy33, ang1, ang2, ang3, vg1, vg2, vg3, &
2434  this%iavgkeff)
2435  end if
2436  !
2437  end if
2438  end function hy_eff
2439 
2440  !> @brief Calculate the 3 components of specific discharge at the cell center
2441  !<
2442  subroutine calc_spdis(this, flowja)
2443  ! -- modules
2444  use simmodule, only: store_error
2445  ! -- dummy
2446  class(gwfnpftype) :: this
2447  real(DP), intent(in), dimension(:) :: flowja
2448  ! -- local
2449  integer(I4B) :: n
2450  integer(I4B) :: m
2451  integer(I4B) :: ipos
2452  integer(I4B) :: iedge
2453  integer(I4B) :: isympos
2454  integer(I4B) :: ihc
2455  integer(I4B) :: ic
2456  integer(I4B) :: iz
2457  integer(I4B) :: nc
2458  integer(I4B) :: ncz
2459  real(DP) :: qz
2460  real(DP) :: vx
2461  real(DP) :: vy
2462  real(DP) :: vz
2463  real(DP) :: xn
2464  real(DP) :: yn
2465  real(DP) :: zn
2466  real(DP) :: xc
2467  real(DP) :: yc
2468  real(DP) :: zc
2469  real(DP) :: cl1
2470  real(DP) :: cl2
2471  real(DP) :: dltot
2472  real(DP) :: ooclsum
2473  real(DP) :: dsumx
2474  real(DP) :: dsumy
2475  real(DP) :: dsumz
2476  real(DP) :: denom
2477  real(DP) :: area
2478  real(DP) :: dz
2479  real(DP) :: axy
2480  real(DP) :: ayx
2481  logical :: nozee = .true.
2482  type(spdisworkarraytype), pointer :: swa => null() !< pointer to spdis work arrays structure
2483  !
2484  ! -- Ensure dis has necessary information
2485  if (this%icalcspdis /= 0 .and. this%dis%con%ianglex == 0) then
2486  call store_error('Error. ANGLDEGX not provided in '// &
2487  'discretization file. ANGLDEGX required for '// &
2488  'calculation of specific discharge.', terminate=.true.)
2489  end if
2490 
2491  swa => this%spdis_wa
2492  if (.not. swa%is_created()) then
2493  ! prepare work arrays
2494  call this%spdis_wa%create(this%calc_max_conns())
2495 
2496  ! prepare lookup table
2497  if (this%nedges > 0) call this%prepare_edge_lookup()
2498  end if
2499  !
2500  ! -- Go through each cell and calculate specific discharge
2501  do n = 1, this%dis%nodes
2502  !
2503  ! -- first calculate geometric properties for x and y directions and
2504  ! the specific discharge at a face (vi)
2505  ic = 0
2506  iz = 0
2507 
2508  ! reset work arrays
2509  call swa%reset()
2510 
2511  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
2512  m = this%dis%con%ja(ipos)
2513  isympos = this%dis%con%jas(ipos)
2514  ihc = this%dis%con%ihc(isympos)
2515  area = this%dis%con%hwva(isympos)
2516  if (ihc == c3d_vertical) then
2517  !
2518  ! -- vertical connection
2519  iz = iz + 1
2520  !call this%dis%connection_normal(n, m, ihc, xn, yn, zn, ipos)
2521  call this%dis%connection_vector(n, m, nozee, this%sat(n), this%sat(m), &
2522  ihc, xc, yc, zc, dltot)
2523  cl1 = this%dis%con%cl1(isympos)
2524  cl2 = this%dis%con%cl2(isympos)
2525  if (m < n) then
2526  cl1 = this%dis%con%cl2(isympos)
2527  cl2 = this%dis%con%cl1(isympos)
2528  end if
2529  ooclsum = done / (cl1 + cl2)
2530  swa%diz(iz) = dltot * cl1 * ooclsum
2531  qz = flowja(ipos)
2532  if (n > m) qz = -qz
2533  swa%viz(iz) = qz / area
2534  else
2535  !
2536  ! -- horizontal connection
2537  ic = ic + 1
2538  dz = thksatnm(this%ibound(n), this%ibound(m), &
2539  this%icelltype(n), this%icelltype(m), &
2540  this%inewton, ihc, &
2541  this%hnew(n), this%hnew(m), this%sat(n), this%sat(m), &
2542  this%dis%top(n), this%dis%top(m), this%dis%bot(n), &
2543  this%dis%bot(m))
2544  area = area * dz
2545  call this%dis%connection_normal(n, m, ihc, xn, yn, zn, ipos)
2546  call this%dis%connection_vector(n, m, nozee, this%sat(n), this%sat(m), &
2547  ihc, xc, yc, zc, dltot)
2548  cl1 = this%dis%con%cl1(isympos)
2549  cl2 = this%dis%con%cl2(isympos)
2550  if (m < n) then
2551  cl1 = this%dis%con%cl2(isympos)
2552  cl2 = this%dis%con%cl1(isympos)
2553  end if
2554  ooclsum = done / (cl1 + cl2)
2555  swa%nix(ic) = -xn
2556  swa%niy(ic) = -yn
2557  swa%di(ic) = dltot * cl1 * ooclsum
2558  if (area > dzero) then
2559  swa%vi(ic) = flowja(ipos) / area
2560  else
2561  swa%vi(ic) = dzero
2562  end if
2563  end if
2564  end do
2565 
2566  ! add contribution from edge flows (i.e. from exchanges)
2567  if (this%nedges > 0) then
2568  do ipos = this%iedge_ptr(n), this%iedge_ptr(n + 1) - 1
2569  iedge = this%edge_idxs(ipos)
2570 
2571  ! propsedge: (Q, area, nx, ny, distance)
2572  ihc = this%ihcedge(iedge)
2573  area = this%propsedge(2, iedge)
2574  if (ihc == c3d_vertical) then
2575  iz = iz + 1
2576  swa%viz(iz) = this%propsedge(1, iedge) / area
2577  swa%diz(iz) = this%propsedge(5, iedge)
2578  else
2579  ic = ic + 1
2580  swa%nix(ic) = -this%propsedge(3, iedge)
2581  swa%niy(ic) = -this%propsedge(4, iedge)
2582  swa%di(ic) = this%propsedge(5, iedge)
2583  if (area > dzero) then
2584  swa%vi(ic) = this%propsedge(1, iedge) / area
2585  else
2586  swa%vi(ic) = dzero
2587  end if
2588  end if
2589  end do
2590  end if
2591  !
2592  ! -- Assign number of vertical and horizontal connections
2593  ncz = iz
2594  nc = ic
2595  !
2596  ! -- calculate z weight (wiz) and z velocity
2597  if (ncz == 1) then
2598  swa%wiz(1) = done
2599  else
2600  dsumz = dzero
2601  do iz = 1, ncz
2602  dsumz = dsumz + swa%diz(iz)
2603  end do
2604  denom = (ncz - done)
2605  if (denom < dzero) denom = dzero
2606  dsumz = dsumz + dem10 * dsumz
2607  do iz = 1, ncz
2608  if (dsumz > dzero) swa%wiz(iz) = done - swa%diz(iz) / dsumz
2609  if (denom > 0) then
2610  swa%wiz(iz) = swa%wiz(iz) / denom
2611  else
2612  swa%wiz(iz) = dzero
2613  end if
2614  end do
2615  end if
2616  vz = dzero
2617  do iz = 1, ncz
2618  vz = vz + swa%wiz(iz) * swa%viz(iz)
2619  end do
2620  !
2621  ! -- distance-based weighting
2622  nc = ic
2623  dsumx = dzero
2624  dsumy = dzero
2625  dsumz = dzero
2626  do ic = 1, nc
2627  swa%wix(ic) = swa%di(ic) * abs(swa%nix(ic))
2628  swa%wiy(ic) = swa%di(ic) * abs(swa%niy(ic))
2629  dsumx = dsumx + swa%wix(ic)
2630  dsumy = dsumy + swa%wiy(ic)
2631  end do
2632  !
2633  ! -- Finish computing omega weights. Add a tiny bit
2634  ! to dsum so that the normalized omega weight later
2635  ! evaluates to (essentially) 1 in the case of a single
2636  ! relevant connection, avoiding 0/0.
2637  dsumx = dsumx + dem10 * dsumx
2638  dsumy = dsumy + dem10 * dsumy
2639  do ic = 1, nc
2640  swa%wix(ic) = (dsumx - swa%wix(ic)) * abs(swa%nix(ic))
2641  swa%wiy(ic) = (dsumy - swa%wiy(ic)) * abs(swa%niy(ic))
2642  end do
2643  !
2644  ! -- compute B weights
2645  dsumx = dzero
2646  dsumy = dzero
2647  do ic = 1, nc
2648  swa%bix(ic) = swa%wix(ic) * sign(done, swa%nix(ic))
2649  swa%biy(ic) = swa%wiy(ic) * sign(done, swa%niy(ic))
2650  dsumx = dsumx + swa%wix(ic) * abs(swa%nix(ic))
2651  dsumy = dsumy + swa%wiy(ic) * abs(swa%niy(ic))
2652  end do
2653  if (dsumx > dzero) dsumx = done / dsumx
2654  if (dsumy > dzero) dsumy = done / dsumy
2655  axy = dzero
2656  ayx = dzero
2657  do ic = 1, nc
2658  swa%bix(ic) = swa%bix(ic) * dsumx
2659  swa%biy(ic) = swa%biy(ic) * dsumy
2660  axy = axy + swa%bix(ic) * swa%niy(ic)
2661  ayx = ayx + swa%biy(ic) * swa%nix(ic)
2662  end do
2663  !
2664  ! -- Calculate specific discharge. The divide by zero checking below
2665  ! is problematic for cells with only one flow, such as can happen
2666  ! with triangular cells in corners. In this case, the resulting
2667  ! cell velocity will be calculated as zero. The method should be
2668  ! improved so that edge flows of zero are included in these
2669  ! calculations. But this needs to be done with consideration for LGR
2670  ! cases in which flows are submitted from an exchange.
2671  vx = dzero
2672  vy = dzero
2673  do ic = 1, nc
2674  vx = vx + (swa%bix(ic) - axy * swa%biy(ic)) * swa%vi(ic)
2675  vy = vy + (swa%biy(ic) - ayx * swa%bix(ic)) * swa%vi(ic)
2676  end do
2677  denom = done - axy * ayx
2678  if (denom /= dzero) then
2679  vx = vx / denom
2680  vy = vy / denom
2681  end if
2682  !
2683  this%spdis(1, n) = vx
2684  this%spdis(2, n) = vy
2685  this%spdis(3, n) = vz
2686  !
2687  end do
2688 
2689  end subroutine calc_spdis
2690 
2691  !> @brief Save specific discharge in binary format to ibinun
2692  !<
2693  subroutine sav_spdis(this, ibinun)
2694  ! -- dummy
2695  class(gwfnpftype) :: this
2696  integer(I4B), intent(in) :: ibinun
2697  ! -- local
2698  character(len=16) :: text
2699  character(len=16), dimension(3) :: auxtxt
2700  integer(I4B) :: n
2701  integer(I4B) :: naux
2702  !
2703  ! -- Write the header
2704  text = ' DATA-SPDIS'
2705  naux = 3
2706  auxtxt(:) = [' qx', ' qy', ' qz']
2707  call this%dis%record_srcdst_list_header(text, this%name_model, &
2708  this%packName, this%name_model, &
2709  this%packName, naux, auxtxt, ibinun, &
2710  this%dis%nodes, this%iout)
2711  !
2712  ! -- Write a zero for Q, and then write qx, qy, qz as aux variables
2713  do n = 1, this%dis%nodes
2714  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, &
2715  this%spdis(:, n))
2716  end do
2717  end subroutine sav_spdis
2718 
2719  !> @brief Save saturation in binary format to ibinun
2720  !<
2721  subroutine sav_sat(this, ibinun)
2722  ! -- dummy
2723  class(gwfnpftype) :: this
2724  integer(I4B), intent(in) :: ibinun
2725  ! -- local
2726  character(len=16) :: text
2727  character(len=16), dimension(1) :: auxtxt
2728  real(DP), dimension(1) :: a
2729  integer(I4B) :: n
2730  integer(I4B) :: naux
2731  !
2732  ! -- Write the header
2733  text = ' DATA-SAT'
2734  naux = 1
2735  auxtxt(:) = [' sat']
2736  call this%dis%record_srcdst_list_header(text, this%name_model, &
2737  this%packName, this%name_model, &
2738  this%packName, naux, auxtxt, ibinun, &
2739  this%dis%nodes, this%iout)
2740  !
2741  ! -- Write a zero for Q, and then write saturation as an aux variables
2742  do n = 1, this%dis%nodes
2743  a(1) = this%sat(n)
2744  call this%dis%record_mf6_list_entry(ibinun, n, n, dzero, naux, a)
2745  end do
2746  end subroutine sav_sat
2747 
2748  !> @brief Reserve space for nedges cells that have an edge on them.
2749  !!
2750  !! This must be called before the npf%allocate_arrays routine, which is
2751  !! called from npf%ar.
2752  !<
2753  subroutine increase_edge_count(this, nedges)
2754  ! -- dummy
2755  class(gwfnpftype) :: this
2756  integer(I4B), intent(in) :: nedges
2757  !
2758  this%nedges = this%nedges + nedges
2759  end subroutine increase_edge_count
2760 
2761  !> @brief Calculate the maximum number of connections for any cell
2762  !<
2763  function calc_max_conns(this) result(max_conns)
2764  class(gwfnpftype) :: this
2765  integer(I4B) :: max_conns
2766  ! local
2767  integer(I4B) :: n, m, ic
2768 
2769  max_conns = 0
2770  do n = 1, this%dis%nodes
2771 
2772  ! Count internal model connections
2773  ic = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
2774 
2775  ! Add edge connections
2776  do m = 1, this%nedges
2777  if (this%nodedge(m) == n) then
2778  ic = ic + 1
2779  end if
2780  end do
2781 
2782  ! Set max number of connections for any cell
2783  if (ic > max_conns) max_conns = ic
2784  end do
2785 
2786  end function calc_max_conns
2787 
2788  !> @brief Provide the npf package with edge properties
2789  !<
2790  subroutine set_edge_properties(this, nodedge, ihcedge, q, area, nx, ny, &
2791  distance)
2792  ! -- dummy
2793  class(gwfnpftype) :: this
2794  integer(I4B), intent(in) :: nodedge
2795  integer(I4B), intent(in) :: ihcedge
2796  real(DP), intent(in) :: q
2797  real(DP), intent(in) :: area
2798  real(DP), intent(in) :: nx
2799  real(DP), intent(in) :: ny
2800  real(DP), intent(in) :: distance
2801  ! -- local
2802  integer(I4B) :: lastedge
2803  !
2804  this%lastedge = this%lastedge + 1
2805  lastedge = this%lastedge
2806  this%nodedge(lastedge) = nodedge
2807  this%ihcedge(lastedge) = ihcedge
2808  this%propsedge(1, lastedge) = q
2809  this%propsedge(2, lastedge) = area
2810  this%propsedge(3, lastedge) = nx
2811  this%propsedge(4, lastedge) = ny
2812  this%propsedge(5, lastedge) = distance
2813  !
2814  ! -- If this is the last edge, then the next call must be starting a new
2815  ! edge properties assignment loop, so need to reset lastedge to 0
2816  if (this%lastedge == this%nedges) this%lastedge = 0
2817  end subroutine set_edge_properties
2818 
2819  subroutine prepare_edge_lookup(this)
2820  class(gwfnpftype) :: this
2821  ! local
2822  integer(I4B) :: i, inode, iedge
2823  integer(I4B) :: n, start, end
2824  integer(I4B) :: prev_cnt, strt_idx, ipos
2825 
2826  do i = 1, size(this%iedge_ptr)
2827  this%iedge_ptr(i) = 0
2828  end do
2829  do i = 1, size(this%edge_idxs)
2830  this%edge_idxs(i) = 0
2831  end do
2832 
2833  ! count
2834  do iedge = 1, this%nedges
2835  n = this%nodedge(iedge)
2836  this%iedge_ptr(n) = this%iedge_ptr(n) + 1
2837  end do
2838 
2839  ! determine start indexes
2840  prev_cnt = this%iedge_ptr(1)
2841  this%iedge_ptr(1) = 1
2842  do inode = 2, this%dis%nodes + 1
2843  strt_idx = this%iedge_ptr(inode - 1) + prev_cnt
2844  prev_cnt = this%iedge_ptr(inode)
2845  this%iedge_ptr(inode) = strt_idx
2846  end do
2847 
2848  ! loop over edges to fill lookup table
2849  do iedge = 1, this%nedges
2850  n = this%nodedge(iedge)
2851  start = this%iedge_ptr(n)
2852  end = this%iedge_ptr(n + 1) - 1
2853  do ipos = start, end
2854  if (this%edge_idxs(ipos) > 0) cycle ! go to next
2855  this%edge_idxs(ipos) = iedge
2856  exit
2857  end do
2858  end do
2859 
2860  end subroutine prepare_edge_lookup
2861 
2862  !> Calculate saturated thickness between cell n and m
2863  !<
2864  function calcsatthickness(this, n, m, ihc) result(satThickness)
2865  ! -- dummy
2866  class(gwfnpftype) :: this !< this NPF instance
2867  integer(I4B) :: n !< node n
2868  integer(I4B) :: m !< node m
2869  integer(I4B) :: ihc !< 1 = horizontal connection, 0 for vertical
2870  ! -- return
2871  real(dp) :: satthickness !< saturated thickness
2872  !
2873  satthickness = thksatnm(this%ibound(n), &
2874  this%ibound(m), &
2875  this%icelltype(n), &
2876  this%icelltype(m), &
2877  this%inewton, &
2878  ihc, &
2879  this%hnew(n), &
2880  this%hnew(m), &
2881  this%sat(n), &
2882  this%sat(m), &
2883  this%dis%top(n), &
2884  this%dis%top(m), &
2885  this%dis%bot(n), &
2886  this%dis%bot(m))
2887  end function calcsatthickness
2888 
2889 end module gwfnpfmodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
@ c3d_vertical
vertical connection
Definition: Constants.f90:222
real(dp), parameter dhdry
real dry cell constant
Definition: Constants.f90:94
real(dp), parameter dp9
real constant 9/10
Definition: Constants.f90:72
real(dp), parameter dem10
real constant 1e-10
Definition: Constants.f90:113
real(dp), parameter dem7
real constant 1e-7
Definition: Constants.f90:110
real(dp), parameter dem8
real constant 1e-8
Definition: Constants.f90:111
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dpio180
real constant
Definition: Constants.f90:130
real(dp), parameter dem6
real constant 1e-6
Definition: Constants.f90:109
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dem9
real constant 1e-9
Definition: Constants.f90:112
real(dp), parameter dem2
real constant 1e-2
Definition: Constants.f90:105
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module contains stateless conductance functions.
real(dp) function, public thksatnm(ibdn, ibdm, ictn, ictm, iupstream, ihc, hn, hm, satn, satm, topn, topm, botn, botm)
Calculate wetted cell thickness at interface between two cells.
real(dp) function, public condmean(k1, k2, thick1, thick2, cl1, cl2, width, iavgmeth)
Calculate the conductance between two cells.
real(dp) function, public hcond(ibdn, ibdm, ictn, ictm, iupstream, ihc, icellavg, condsat, hn, hm, satn, satm, hkn, hkm, topn, topm, botn, botm, cln, clm, fawidth)
Horizontal conductance between two cells.
@, public ccond_hmean
Harmonic mean.
real(dp) function, public vcond(ibdn, ibdm, ictn, ictm, inewton, ivarcv, idewatcv, condsat, hn, hm, vkn, vkm, satn, satm, topn, topm, botn, botm, flowarea)
Vertical conductance between two cells.
subroutine set_options(this, options)
Set options in the NPF object.
Definition: gwf-npf.f90:1443
subroutine npf_save_model_flows(this, flowja, icbcfl, icbcun)
Record flowja and calculate specific discharge if requested.
Definition: gwf-npf.f90:924
subroutine source_options(this)
Update simulation options from input mempath.
Definition: gwf-npf.f90:1356
real(dp) function calc_initial_sat(this, n)
Calculate initial saturation for the given node.
Definition: gwf-npf.f90:2130
subroutine calc_condsat(this, node, upperOnly)
Calculate CONDSAT array entries for the given node.
Definition: gwf-npf.f90:2034
subroutine rewet_check(this, kiter, node, hm, ibdm, ihc, hnew, irewet)
Determine if a cell should rewet.
Definition: gwf-npf.f90:2251
integer(i4b) function calc_max_conns(this)
Calculate the maximum number of connections for any cell.
Definition: gwf-npf.f90:2764
subroutine npf_mc(this, moffset, matrix_sln)
Map connections and construct iax, jax, and idxglox.
Definition: gwf-npf.f90:274
subroutine npf_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
Formulate coefficients.
Definition: gwf-npf.f90:477
subroutine npf_ac(this, moffset, sparse)
Add connections for extended neighbors to the sparse matrix.
Definition: gwf-npf.f90:260
subroutine sgwf_npf_wetdry(this, kiter, hnew)
Perform wetting and drying.
Definition: gwf-npf.f90:2145
subroutine source_griddata(this)
Update simulation griddata from input mempath.
Definition: gwf-npf.f90:1582
subroutine npf_nur(this, neqmod, x, xtemp, dx, inewtonur, dxmax, locmax)
Under-relaxation.
Definition: gwf-npf.f90:749
subroutine sav_spdis(this, ibinun)
Save specific discharge in binary format to ibinun.
Definition: gwf-npf.f90:2694
subroutine sgwf_npf_thksat(this, n, hn, thksat)
Fractional cell saturation.
Definition: gwf-npf.f90:823
subroutine preprocess_input(this)
preprocess the NPF input data
Definition: gwf-npf.f90:1832
subroutine set_edge_properties(this, nodedge, ihcedge, q, area, nx, ny, distance)
Provide the npf package with edge properties.
Definition: gwf-npf.f90:2792
subroutine prepare_edge_lookup(this)
Definition: gwf-npf.f90:2820
subroutine npf_da(this)
Deallocate variables.
Definition: gwf-npf.f90:1000
subroutine log_griddata(this, found)
Write dimensions to list file.
Definition: gwf-npf.f90:1532
real(dp) function hy_eff(this, n, m, ihc, ipos, vg)
Calculate the effective hydraulic conductivity for the n-m connection.
Definition: gwf-npf.f90:2364
subroutine calc_spdis(this, flowja)
Calculate the 3 components of specific discharge at the cell center.
Definition: gwf-npf.f90:2443
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.
Definition: gwf-npf.f90:160
subroutine increase_edge_count(this, nedges)
Reserve space for nedges cells that have an edge on them.
Definition: gwf-npf.f90:2754
subroutine npf_fn(this, kiter, matrix_sln, idxglo, rhs, hnew)
Fill newton terms.
Definition: gwf-npf.f90:622
subroutine log_options(this, found)
Log npf options sourced from the input mempath.
Definition: gwf-npf.f90:1282
subroutine npf_print_model_flows(this, ibudfl, flowja)
Print budget.
Definition: gwf-npf.f90:961
subroutine allocate_arrays(this, ncells, njas)
Allocate npf arrays.
Definition: gwf-npf.f90:1221
subroutine sav_sat(this, ibinun)
Save saturation in binary format to ibinun.
Definition: gwf-npf.f90:2722
subroutine prepcheck(this)
Initialize and check NPF data.
Definition: gwf-npf.f90:1671
subroutine npf_df(this, dis, xt3d, ingnc, invsc, npf_options)
Define the NPF package instance.
Definition: gwf-npf.f90:207
subroutine highest_cell_saturation(this, n, m, hn, hm, satn, satm)
Calculate dry cell saturation.
Definition: gwf-npf.f90:593
subroutine npf_ad(this, nodes, hold, hnew, irestore)
Advance.
Definition: gwf-npf.f90:383
subroutine sgwf_npf_wdmsg(this, icode, ncnvrt, nodcnvrt, acnvrt, ihdcnv, kiter, n)
Print wet/dry message.
Definition: gwf-npf.f90:2312
subroutine sgwf_npf_qcalc(this, n, m, hn, hm, icon, qnm)
Flow between two cells.
Definition: gwf-npf.f90:846
subroutine npf_rp(this)
Read and prepare method for package.
Definition: gwf-npf.f90:368
real(dp) function calcsatthickness(this, n, m, ihc)
Calculate saturated thickness between cell n and m.
Definition: gwf-npf.f90:2865
subroutine store_original_k_arrays(this, ncells, njas)
@ brief Store backup copy of hydraulic conductivity when the VSC package is activate
Definition: gwf-npf.f90:1201
subroutine allocate_scalars(this)
@ brief Allocate scalars
Definition: gwf-npf.f90:1102
subroutine npf_cq(this, hnew, flowja)
Calculate flowja.
Definition: gwf-npf.f90:793
subroutine npf_ar(this, ic, vsc, ibound, hnew)
Allocate and read this NPF instance.
Definition: gwf-npf.f90:288
subroutine check_options(this)
Check for conflicting NPF options.
Definition: gwf-npf.f90:1461
subroutine npf_cf(this, kiter, nodes, hnew)
Routines associated fill coefficients.
Definition: gwf-npf.f90:447
General-purpose hydrogeologic functions.
Definition: HGeoUtil.f90:2
real(dp) function, public hyeff(k11, k22, k33, ang1, ang2, ang3, vg1, vg2, vg3, iavgmeth)
Calculate the effective horizontal hydraulic conductivity from an ellipse using a specified direction...
Definition: HGeoUtil.f90:31
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
subroutine, public memorystore_remove(component, subcomponent, context)
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the base numerical package type.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
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
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
character(len=linelength) idm_context
character(len=maxcharlen) warnmsg
warning message string
real(dp) function squadraticsaturation(top, bot, x, eps)
@ brief sQuadraticSaturation
real(dp) function squadraticsaturationderivative(top, bot, x, eps)
@ brief Derivative of the quadratic saturation function
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
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
This module contains time-varying conductivity package methods.
Definition: gwf-tvk.f90:8
subroutine, public tvk_cr(tvk, name_model, inunit, iout)
Create a new TvkType object.
Definition: gwf-tvk.f90:53
subroutine, public xt3d_cr(xt3dobj, name_model, inunit, iout, ldispopt)
Create a new xt3d object.
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Data structure and helper methods for passing NPF options into npf_df, as an alternative to reading t...
Helper class with work arrays for the SPDIS calculation in NPF.