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