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