MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
gwf-hfb.f90
Go to the documentation of this file.
1 
3 
4  use kindmodule, only: dp, i4b
5  use simvariablesmodule, only: errmsg
7  use xt3dmodule, only: xt3dtype
8  use gwfvscmodule, only: gwfvsctype
10  use basedismodule, only: disbasetype
12 
13  implicit none
14 
15  private
16  public :: gwfhfbtype
17  public :: hfb_cr
18 
19  type, extends(numericalpackagetype) :: gwfhfbtype
20 
21  type(gwfvsctype), pointer :: vsc => null() !< viscosity object
22  integer(I4B), pointer :: maxhfb => null() !< max number of hfb's
23  integer(I4B), pointer :: nhfb => null() !< number of hfb's
24  integer(I4B), dimension(:), pointer, contiguous :: noden => null() !< first cell
25  integer(I4B), dimension(:), pointer, contiguous :: nodem => null() !< second cell
26  integer(I4B), dimension(:), pointer, contiguous :: idxloc => null() !< position in model ja
27  real(dp), dimension(:), pointer, contiguous :: hydchr => null() !< hydraulic characteristic of the barrier
28  real(dp), dimension(:), pointer, contiguous :: csatsav => null() !< value of condsat prior to hfb modification
29  real(dp), dimension(:), pointer, contiguous :: condsav => null() !< saved conductance of combined npf and hfb
30  type(xt3dtype), pointer :: xt3d => null() !< pointer to xt3d object
31  !
32  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< pointer to model ibound
33  integer(I4B), dimension(:), pointer, contiguous :: icelltype => null() !< pointer to model icelltype
34  integer(I4B), dimension(:), pointer, contiguous :: ihc => null() !< pointer to model ihc
35  integer(I4B), dimension(:), pointer, contiguous :: ia => null() !< pointer to model ia
36  integer(I4B), dimension(:), pointer, contiguous :: ja => null() !< pointer to model ja
37  integer(I4B), dimension(:), pointer, contiguous :: jas => null() !< pointer to model jas
38  integer(I4B), dimension(:), pointer, contiguous :: isym => null() !< pointer to model isym
39  real(dp), dimension(:), pointer, contiguous :: condsat => null() !< pointer to model condsat
40  real(dp), dimension(:), pointer, contiguous :: top => null() !< pointer to model top
41  real(dp), dimension(:), pointer, contiguous :: bot => null() !< pointer to model bot
42  real(dp), dimension(:), pointer, contiguous :: hwva => null() !< pointer to model hwva
43  real(dp), dimension(:), pointer, contiguous :: hnew => null() !< pointer to model xnew
44  !
45  ! -- viscosity flag
46  integer(I4B), pointer :: ivsc => null() !< flag indicating if viscosity is active in the model
47 
48  contains
49 
50  procedure :: hfb_ar
51  procedure :: hfb_rp
52  procedure :: hfb_fc
53  procedure :: hfb_cq
54  procedure :: hfb_da
55  procedure :: allocate_scalars
56  procedure, private :: allocate_arrays
57  procedure, private :: source_options
58  procedure, private :: source_dimensions
59  procedure, private :: source_data
60  procedure, private :: check_data
61  procedure, private :: condsat_reset
62  procedure, private :: condsat_modify
63 
64  end type gwfhfbtype
65 
66 contains
67 
68  !> @brief Create a new hfb object
69  !<
70  subroutine hfb_cr(hfbobj, name_model, input_mempath, inunit, iout)
71  ! -- dummy
72  type(gwfhfbtype), pointer :: hfbobj
73  character(len=*), intent(in) :: name_model
74  character(len=*), intent(in) :: input_mempath
75  integer(I4B), intent(in) :: inunit
76  integer(I4B), intent(in) :: iout
77  !
78  ! -- Create the object
79  allocate (hfbobj)
80  !
81  ! -- create name and memory path
82  call hfbobj%set_names(1, name_model, 'HFB', 'HFB', input_mempath)
83  !
84  ! -- Allocate scalars
85  call hfbobj%allocate_scalars()
86  !
87  ! -- Save unit numbers
88  hfbobj%inunit = inunit
89  hfbobj%iout = iout
90  end subroutine hfb_cr
91 
92  !> @brief Allocate and read
93  !<
94  subroutine hfb_ar(this, ibound, xt3d, dis, invsc, vsc)
95  ! -- modules
98  ! -- dummy
99  class(gwfhfbtype) :: this
100  integer(I4B), dimension(:), pointer, contiguous :: ibound
101  type(xt3dtype), pointer :: xt3d
102  class(disbasetype), pointer, intent(inout) :: dis !< discretization package
103  integer(I4B), pointer :: invsc !< indicates if viscosity package is active
104  type(gwfvsctype), pointer, intent(in) :: vsc !< viscosity package
105  ! -- formats
106  character(len=*), parameter :: fmtheader = &
107  "(1x, /1x, 'HFB -- HORIZONTAL FLOW BARRIER PACKAGE, VERSION 8, ', &
108  &'4/24/2015 INPUT READ FROM MEMPATH: ', a, /)"
109  !
110  ! -- Print a message identifying the node property flow package.
111  write (this%iout, fmtheader) this%input_mempath
112  !
113  ! -- Set pointers
114  this%dis => dis
115  this%ibound => ibound
116  this%xt3d => xt3d
117  !
118  call mem_setptr(this%icelltype, 'ICELLTYPE', &
119  create_mem_path(this%name_model, 'NPF'))
120  call mem_setptr(this%ihc, 'IHC', create_mem_path(this%name_model, 'CON'))
121  call mem_setptr(this%ia, 'IA', create_mem_path(this%name_model, 'CON'))
122  call mem_setptr(this%ja, 'JA', create_mem_path(this%name_model, 'CON'))
123  call mem_setptr(this%jas, 'JAS', create_mem_path(this%name_model, 'CON'))
124  call mem_setptr(this%isym, 'ISYM', create_mem_path(this%name_model, 'CON'))
125  call mem_setptr(this%condsat, 'CONDSAT', create_mem_path(this%name_model, &
126  'NPF'))
127  call mem_setptr(this%top, 'TOP', create_mem_path(this%name_model, 'DIS'))
128  call mem_setptr(this%bot, 'BOT', create_mem_path(this%name_model, 'DIS'))
129  call mem_setptr(this%hwva, 'HWVA', create_mem_path(this%name_model, 'CON'))
130  !
131  call this%source_options()
132  call this%source_dimensions()
133  call this%allocate_arrays()
134  !
135  ! -- If vsc package active, set ivsc
136  if (invsc /= 0) then
137  this%ivsc = 1
138  this%vsc => vsc
139  !
140  ! -- Notify user via listing file viscosity accounted for by HFB package
141  write (this%iout, '(/1x,a,a)') 'Viscosity active in ', &
142  trim(this%filtyp)//' Package calculations: '//trim(adjustl(this%packName))
143  end if
144  end subroutine hfb_ar
145 
146  !> @brief Check for new HFB stress period data
147  !<
148  subroutine hfb_rp(this)
149  ! -- modules
151  use tdismodule, only: kper
152  ! -- dummy
153  class(gwfhfbtype) :: this
154  ! -- local
155  integer(I4B), pointer :: iper
156  ! -- formats
157  character(len=*), parameter :: fmtlsp = &
158  &"(1X,/1X,'REUSING ',A,'S FROM LAST STRESS PERIOD')"
159 
160  call mem_setptr(iper, 'IPER', this%input_mempath)
161  if (iper == kper) then
162  call this%condsat_reset()
163  call this%source_data()
164  call this%condsat_modify()
165  else
166  write (this%iout, fmtlsp) 'HFB'
167  end if
168  end subroutine hfb_rp
169 
170  !> @brief Fill matrix terms
171  !!
172  !! Fill amatsln for the following conditions:
173  !! 1. XT3D
174  !! OR
175  !! 2. Not Newton, and
176  !! 3. Cell type n is convertible or cell type m is convertible
177  !<
178  subroutine hfb_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
179  ! -- modules
180  use constantsmodule, only: dhalf, dzero, done
181  ! -- dummy
182  class(gwfhfbtype) :: this
183  integer(I4B) :: kiter
184  class(matrixbasetype), pointer :: matrix_sln
185  integer(I4B), intent(in), dimension(:) :: idxglo
186  real(DP), intent(inout), dimension(:) :: rhs
187  real(DP), intent(inout), dimension(:) :: hnew
188  ! -- local
189  integer(I4B) :: nodes, nja
190  integer(I4B) :: ihfb, n, m
191  integer(I4B) :: ipos
192  integer(I4B) :: idiag, isymcon
193  integer(I4B) :: ixt3d
194  real(DP) :: cond, condhfb, aterm
195  real(DP) :: fawidth, faheight
196  real(DP) :: topn, topm, botn, botm
197  real(DP) :: viscratio
198  !
199  ! -- initialize variables
200  viscratio = done
201  nodes = this%dis%nodes
202  nja = this%dis%con%nja
203  if (associated(this%xt3d%ixt3d)) then
204  ixt3d = this%xt3d%ixt3d
205  else
206  ixt3d = 0
207  end if
208  !
209  if (ixt3d > 0) then
210  !
211  do ihfb = 1, this%nhfb
212  n = min(this%noden(ihfb), this%nodem(ihfb))
213  m = max(this%noden(ihfb), this%nodem(ihfb))
214  ! -- Skip if either cell is inactive.
215  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
216  !!! if(this%icelltype(n) == 1 .or. this%icelltype(m) == 1) then
217  if (this%ivsc /= 0) then
218  call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
219  end if
220  ! -- Compute scale factor for hfb correction
221  if (this%hydchr(ihfb) > dzero) then
222  if (this%inewton == 0) then
223  ipos = this%idxloc(ihfb)
224  topn = this%top(n)
225  topm = this%top(m)
226  botn = this%bot(n)
227  botm = this%bot(m)
228  if (this%icelltype(n) == 1) then
229  if (hnew(n) < topn) topn = hnew(n)
230  end if
231  if (this%icelltype(m) == 1) then
232  if (hnew(m) < topm) topm = hnew(m)
233  end if
234  if (this%ihc(this%jas(ipos)) == 2) then
235  faheight = min(topn, topm) - max(botn, botm)
236  else
237  faheight = dhalf * ((topn - botn) + (topm - botm))
238  end if
239  fawidth = this%hwva(this%jas(ipos))
240  condhfb = this%hydchr(ihfb) * viscratio * &
241  fawidth * faheight
242  else
243  condhfb = this%hydchr(ihfb) * viscratio
244  end if
245  else
246  condhfb = this%hydchr(ihfb)
247  end if
248  ! -- Make hfb corrections for xt3d
249  call this%xt3d%xt3d_fhfb(kiter, nodes, nja, matrix_sln, idxglo, &
250  rhs, hnew, n, m, condhfb)
251  end do
252  !
253  else
254  !
255  ! -- For Newton, the effect of the barrier is included in condsat.
256  if (this%inewton == 0) then
257  do ihfb = 1, this%nhfb
258  ipos = this%idxloc(ihfb)
259  aterm = matrix_sln%get_value_pos(idxglo(ipos))
260  n = this%noden(ihfb)
261  m = this%nodem(ihfb)
262  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
263  !
264  if (this%ivsc /= 0) then
265  call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
266  end if
267  !
268  if (this%icelltype(n) == 1 .or. this%icelltype(m) == 1 .or. &
269  this%ivsc /= 0) then
270  !
271  ! -- Calculate hfb conductance
272  topn = this%top(n)
273  topm = this%top(m)
274  botn = this%bot(n)
275  botm = this%bot(m)
276  if (this%icelltype(n) == 1) then
277  if (hnew(n) < topn) topn = hnew(n)
278  end if
279  if (this%icelltype(m) == 1) then
280  if (hnew(m) < topm) topm = hnew(m)
281  end if
282  if (this%ihc(this%jas(ipos)) == 2) then
283  faheight = min(topn, topm) - max(botn, botm)
284  else
285  faheight = dhalf * ((topn - botn) + (topm - botm))
286  end if
287  if (this%hydchr(ihfb) > dzero) then
288  fawidth = this%hwva(this%jas(ipos))
289  condhfb = this%hydchr(ihfb) * viscratio * &
290  fawidth * faheight
291  cond = aterm * condhfb / (aterm + condhfb)
292  else
293  cond = -aterm * this%hydchr(ihfb)
294  end if
295  !
296  ! -- Save cond for budget calculation
297  this%condsav(ihfb) = cond
298  !
299  ! -- Fill row n diag and off diag
300  idiag = this%ia(n)
301  call matrix_sln%add_value_pos(idxglo(idiag), aterm - cond)
302  call matrix_sln%set_value_pos(idxglo(ipos), cond)
303  !
304  ! -- Fill row m diag and off diag
305  isymcon = this%isym(ipos)
306  idiag = this%ia(m)
307  call matrix_sln%add_value_pos(idxglo(idiag), aterm - cond)
308  call matrix_sln%set_value_pos(idxglo(isymcon), cond)
309  !
310  end if
311  end do
312  end if
313  !
314  end if
315  end subroutine hfb_fc
316 
317  !> @brief flowja will automatically include the effects of the hfb for
318  !! confined and newton cases when xt3d is not used.
319  !!
320  !! This method recalculates flowja for the other cases.
321  !<
322  subroutine hfb_cq(this, hnew, flowja)
323  ! -- modules
324  use constantsmodule, only: dhalf, dzero, done
325  ! -- dummy
326  class(gwfhfbtype) :: this
327  real(DP), intent(inout), dimension(:) :: hnew
328  real(DP), intent(inout), dimension(:) :: flowja
329  ! -- local
330  integer(I4B) :: ihfb, n, m
331  integer(I4B) :: ipos
332  real(DP) :: qnm
333  real(DP) :: cond
334  integer(I4B) :: ixt3d
335  real(DP) :: condhfb
336  real(DP) :: fawidth, faheight
337  real(DP) :: topn, topm, botn, botm
338  real(DP) :: viscratio
339  !
340  ! -- initialize viscratio
341  viscratio = done
342  !
343  if (associated(this%xt3d%ixt3d)) then
344  ixt3d = this%xt3d%ixt3d
345  else
346  ixt3d = 0
347  end if
348  !
349  if (ixt3d > 0) then
350  !
351  do ihfb = 1, this%nhfb
352  n = min(this%noden(ihfb), this%nodem(ihfb))
353  m = max(this%noden(ihfb), this%nodem(ihfb))
354  ! -- Skip if either cell is inactive.
355  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
356  !!! if(this%icelltype(n) == 1 .or. this%icelltype(m) == 1) then
357  if (this%ivsc /= 0) then
358  call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
359  end if
360  !
361  ! -- Compute scale factor for hfb correction
362  if (this%hydchr(ihfb) > dzero) then
363  if (this%inewton == 0) then
364  ipos = this%idxloc(ihfb)
365  topn = this%top(n)
366  topm = this%top(m)
367  botn = this%bot(n)
368  botm = this%bot(m)
369  if (this%icelltype(n) == 1) then
370  if (hnew(n) < topn) topn = hnew(n)
371  end if
372  if (this%icelltype(m) == 1) then
373  if (hnew(m) < topm) topm = hnew(m)
374  end if
375  if (this%ihc(this%jas(ipos)) == 2) then
376  faheight = min(topn, topm) - max(botn, botm)
377  else
378  faheight = dhalf * ((topn - botn) + (topm - botm))
379  end if
380  fawidth = this%hwva(this%jas(ipos))
381  condhfb = this%hydchr(ihfb) * viscratio * &
382  fawidth * faheight
383  else
384  condhfb = this%hydchr(ihfb)
385  end if
386  else
387  condhfb = this%hydchr(ihfb)
388  end if
389  ! -- Make hfb corrections for xt3d
390  call this%xt3d%xt3d_flowjahfb(n, m, hnew, flowja, condhfb)
391  end do
392  !
393  else
394  !
395  ! -- Recalculate flowja for non-newton unconfined.
396  if (this%inewton == 0) then
397  do ihfb = 1, this%nhfb
398  n = this%noden(ihfb)
399  m = this%nodem(ihfb)
400  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
401  if (this%icelltype(n) == 1 .or. this%icelltype(m) == 1 .or. &
402  this%ivsc /= 0) then
403  ipos = this%dis%con%getjaindex(n, m)
404  !
405  ! -- condsav already accnts for visc adjustment
406  cond = this%condsav(ihfb)
407  qnm = cond * (hnew(m) - hnew(n))
408  flowja(ipos) = qnm
409  ipos = this%dis%con%getjaindex(m, n)
410  flowja(ipos) = -qnm
411  !
412  end if
413  end do
414  end if
415  !
416  end if
417  end subroutine hfb_cq
418 
419  !> @brief Deallocate memory
420  !<
421  subroutine hfb_da(this)
422  ! -- modules
424  ! -- dummy
425  class(gwfhfbtype) :: this
426  !
427  ! -- Scalars
428  call mem_deallocate(this%maxhfb)
429  call mem_deallocate(this%nhfb)
430  call mem_deallocate(this%ivsc)
431  !
432  ! -- Arrays
433  if (this%inunit > 0) then
434  call mem_deallocate(this%noden)
435  call mem_deallocate(this%nodem)
436  call mem_deallocate(this%hydchr)
437  call mem_deallocate(this%idxloc)
438  call mem_deallocate(this%csatsav)
439  call mem_deallocate(this%condsav)
440  end if
441  !
442  ! -- deallocate parent
443  call this%NumericalPackageType%da()
444  !
445  ! -- nullify pointers
446  this%xt3d => null()
447  this%inewton => null()
448  this%ibound => null()
449  this%icelltype => null()
450  this%ihc => null()
451  this%ia => null()
452  this%ja => null()
453  this%jas => null()
454  this%isym => null()
455  this%condsat => null()
456  this%top => null()
457  this%bot => null()
458  this%hwva => null()
459  this%vsc => null()
460  end subroutine hfb_da
461 
462  !> @brief Allocate package scalars
463  !<
464  subroutine allocate_scalars(this)
465  ! -- modules
467  ! -- dummy
468  class(gwfhfbtype) :: this
469  !
470  ! -- allocate scalars in NumericalPackageType
471  call this%NumericalPackageType%allocate_scalars()
472  !
473  ! -- allocate scalars
474  call mem_allocate(this%maxhfb, 'MAXHFB', this%memoryPath)
475  call mem_allocate(this%nhfb, 'NHFB', this%memoryPath)
476  !
477  ! -- allocate flag for determining if vsc active
478  call mem_allocate(this%ivsc, 'IVSC', this%memoryPath)
479  !
480  ! -- initialize
481  this%maxhfb = 0
482  this%nhfb = 0
483  this%ivsc = 0
484  end subroutine allocate_scalars
485 
486  !> @brief Allocate package arrays
487  !<
488  subroutine allocate_arrays(this)
489  ! -- modules
491  ! -- dummy
492  class(gwfhfbtype) :: this
493  ! -- local
494  integer(I4B) :: ihfb
495  !
496  call mem_allocate(this%noden, this%maxhfb, 'NODEN', this%memoryPath)
497  call mem_allocate(this%nodem, this%maxhfb, 'NODEM', this%memoryPath)
498  call mem_allocate(this%hydchr, this%maxhfb, 'HYDCHR', this%memoryPath)
499  call mem_allocate(this%idxloc, this%maxhfb, 'IDXLOC', this%memoryPath)
500  call mem_allocate(this%csatsav, this%maxhfb, 'CSATSAV', this%memoryPath)
501  call mem_allocate(this%condsav, this%maxhfb, 'CONDSAV', this%memoryPath)
502  !
503  ! -- initialize idxloc to 0
504  do ihfb = 1, this%maxhfb
505  this%idxloc(ihfb) = 0
506  end do
507  end subroutine allocate_arrays
508 
509  !> @ brief Source hfb input options
510  !<
511  subroutine source_options(this)
512  ! -- modules
515  ! -- dummy
516  class(gwfhfbtype) :: this
517  ! -- local
518  type(gwfhfbparamfoundtype) :: found
519 
520  ! update options from input context
521  call mem_set_value(this%iprpak, 'PRINT_INPUT', this%input_mempath, &
522  found%print_input)
523 
524  ! log options
525  write (this%iout, '(1x,a)') 'PROCESSING HFB OPTIONS'
526  if (found%print_input) then
527  write (this%iout, '(4x,a)') &
528  'THE LIST OF HFBS WILL BE PRINTED.'
529  end if
530  write (this%iout, '(1x,a)') 'END OF HFB OPTIONS'
531  end subroutine source_options
532 
533  !> @ brief Source hfb input options
534  !<
535  subroutine source_dimensions(this)
536  ! -- modules
539  ! -- dummy
540  class(gwfhfbtype) :: this
541  ! -- local
542  type(gwfhfbparamfoundtype) :: found
543 
544  ! update dimensions from input context
545  call mem_set_value(this%maxhfb, 'MAXBOUND', this%input_mempath, &
546  found%maxbound)
547 
548  ! log dimensions
549  write (this%iout, '(/1x,a)') 'PROCESSING HFB DIMENSIONS'
550  write (this%iout, '(4x,a,i7)') 'MAXHFB = ', this%maxhfb
551  write (this%iout, '(1x,a)') 'END OF HFB DIMENSIONS'
552 
553  ! check dimensions
554  if (this%maxhfb <= 0) then
555  write (errmsg, '(a)') &
556  'MAXHFB must be specified with value greater than zero.'
557  call store_error(errmsg)
558  call store_error_filename(this%input_mempath)
559  end if
560  end subroutine source_dimensions
561 
562  !> @ brief source hfb period data
563  !<
564  subroutine source_data(this)
565  ! -- modules
566  use tdismodule, only: kper
567  use constantsmodule, only: linelength
569  use geomutilmodule, only: get_node
570  ! -- dummy
571  class(gwfhfbtype), intent(inout) :: this
572  ! -- local
573  integer(I4B), dimension(:, :), pointer, contiguous :: cellids1, cellids2
574  integer(I4B), dimension(:), pointer :: cellid1, cellid2
575  real(DP), dimension(:), pointer, contiguous :: hydchr
576  character(len=LINELENGTH) :: nodenstr, nodemstr
577  integer(I4B), pointer :: nbound
578  integer(I4B) :: n, nodeu1, nodeu2, noder1, noder2
579  ! -- formats
580  character(len=*), parameter :: fmthfb = "(i10, 2a10, 1(1pg15.6))"
581 
582  ! set input context pointers
583  call mem_setptr(nbound, 'NBOUND', this%input_mempath)
584  call mem_setptr(cellids1, 'CELLID1', this%input_mempath)
585  call mem_setptr(cellids2, 'CELLID2', this%input_mempath)
586  call mem_setptr(hydchr, 'HYDCHR', this%input_mempath)
587 
588  ! set nhfb
589  this%nhfb = nbound
590 
591  ! log data
592  write (this%iout, '(//,1x,a)') 'READING HFB DATA'
593  if (this%iprpak > 0) then
594  write (this%iout, '(3a10, 1a15)') 'HFB NUM', 'CELL1', 'CELL2', &
595  'HYDCHR'
596  end if
597 
598  ! update state
599  do n = 1, this%nhfb
600 
601  ! set cellid
602  cellid1 => cellids1(:, n)
603  cellid2 => cellids2(:, n)
604 
605  ! set node user
606  if (this%dis%ndim == 1) then
607  nodeu1 = cellid1(1)
608  nodeu2 = cellid2(1)
609  elseif (this%dis%ndim == 2) then
610  nodeu1 = get_node(cellid1(1), 1, cellid1(2), &
611  this%dis%mshape(1), 1, &
612  this%dis%mshape(2))
613  nodeu2 = get_node(cellid2(1), 1, cellid2(2), &
614  this%dis%mshape(1), 1, &
615  this%dis%mshape(2))
616  else
617  nodeu1 = get_node(cellid1(1), cellid1(2), cellid1(3), &
618  this%dis%mshape(1), &
619  this%dis%mshape(2), &
620  this%dis%mshape(3))
621  nodeu2 = get_node(cellid2(1), cellid2(2), cellid2(3), &
622  this%dis%mshape(1), &
623  this%dis%mshape(2), &
624  this%dis%mshape(3))
625  end if
626 
627  ! set nodes
628  noder1 = this%dis%get_nodenumber(nodeu1, 1)
629  noder2 = this%dis%get_nodenumber(nodeu2, 1)
630  if (noder1 <= 0 .or. &
631  noder2 <= 0) then
632  cycle
633  else
634  this%noden(n) = noder1
635  this%nodem(n) = noder2
636  end if
637 
638  this%hydchr(n) = hydchr(n)
639 
640  ! print input if requested
641  if (this%iprpak /= 0) then
642  call this%dis%noder_to_string(this%noden(n), nodenstr)
643  call this%dis%noder_to_string(this%nodem(n), nodemstr)
644  write (this%iout, fmthfb) n, trim(adjustl(nodenstr)), &
645  trim(adjustl(nodemstr)), this%hydchr(n)
646  end if
647  end do
648 
649  ! check errors
650  if (count_errors() > 0) then
651  call store_error('Errors encountered in HFB input file.')
652  call store_error_filename(this%input_fname)
653  end if
654 
655  ! finalize logging
656  write (this%iout, '(3x,i0,a,i0)') this%nhfb, &
657  ' HFBs READ FOR STRESS PERIOD ', kper
658  write (this%iout, '(1x,a)') 'END READING HFB DATA'
659 
660  ! input data check
661  call this%check_data()
662  end subroutine source_data
663 
664  !> @brief Check for hfb's between two unconnected cells and write a warning
665  !!
666  !! Store ipos in idxloc
667  !<
668  subroutine check_data(this)
669  ! -- modules
670  use constantsmodule, only: linelength
671  ! -- dummy
672  class(gwfhfbtype) :: this
673  ! -- local
674  integer(I4B) :: ihfb, n, m
675  integer(I4B) :: ipos
676  character(len=LINELENGTH) :: nodenstr, nodemstr
677  logical :: found
678  ! -- formats
679  character(len=*), parameter :: fmterr = "(1x, 'HFB no. ',i0, &
680  &' is between two unconnected cells: ', a, ' and ', a)"
681  character(len=*), parameter :: fmtverr = "(1x, 'HFB no. ',i0, &
682  &' is between two cells not horizontally connected: ', a, ' and ', a)"
683  !
684  do ihfb = 1, this%nhfb
685  n = this%noden(ihfb)
686  m = this%nodem(ihfb)
687  found = .false.
688  do ipos = this%ia(n) + 1, this%ia(n + 1) - 1
689  if (m == this%ja(ipos)) then
690  found = .true.
691  this%idxloc(ihfb) = ipos
692  exit
693  end if
694  end do
695  !
696  ! -- check to make sure cells are connected
697  if (.not. found) then
698  call this%dis%noder_to_string(n, nodenstr)
699  call this%dis%noder_to_string(m, nodemstr)
700  write (errmsg, fmterr) ihfb, trim(adjustl(nodenstr)), &
701  trim(adjustl(nodemstr))
702  call store_error(errmsg)
703  else
704  !
705  ! -- check to make sure cells are not vertically connected
706  ipos = this%idxloc(ihfb)
707  if (this%ihc(this%jas(ipos)) == 0) then
708  call this%dis%noder_to_string(n, nodenstr)
709  call this%dis%noder_to_string(m, nodemstr)
710  write (errmsg, fmtverr) ihfb, trim(adjustl(nodenstr)), &
711  trim(adjustl(nodemstr))
712  call store_error(errmsg)
713  end if
714  end if
715  end do
716  !
717  ! -- Stop if errors detected
718  if (count_errors() > 0) then
719  call store_error_filename(this%input_fname)
720  end if
721  end subroutine check_data
722 
723  !> @brief Reset condsat to its value prior to being modified by hfb's
724  !<
725  subroutine condsat_reset(this)
726  ! -- dummy
727  class(gwfhfbtype) :: this
728  ! -- local
729  integer(I4B) :: ihfb
730  integer(I4B) :: ipos
731  !
732  do ihfb = 1, this%nhfb
733  ipos = this%idxloc(ihfb)
734  this%condsat(this%jas(ipos)) = this%csatsav(ihfb)
735  end do
736  end subroutine condsat_reset
737 
738  !> @brief Modify condsat
739  !!
740  !! Modify condsat for the following conditions:
741  !! 1. If Newton is active
742  !! 2. If icelltype for n and icelltype for m is 0
743  !<
744  subroutine condsat_modify(this)
745  ! -- modules
746  use constantsmodule, only: dhalf, dzero
747  ! -- dummy
748  class(gwfhfbtype) :: this
749  ! -- local
750  integer(I4B) :: ihfb, n, m
751  integer(I4B) :: ipos
752  real(DP) :: cond, condhfb
753  real(DP) :: fawidth, faheight
754  real(DP) :: topn, topm, botn, botm
755  !
756  do ihfb = 1, this%nhfb
757  ipos = this%idxloc(ihfb)
758  cond = this%condsat(this%jas(ipos))
759  this%csatsav(ihfb) = cond
760  n = this%noden(ihfb)
761  m = this%nodem(ihfb)
762  !
763  if (this%inewton == 1 .or. &
764  (this%icelltype(n) == 0 .and. this%icelltype(m) == 0)) then
765  !
766  ! -- Calculate hfb conductance
767  topn = this%top(n)
768  topm = this%top(m)
769  botn = this%bot(n)
770  botm = this%bot(m)
771  if (this%ihc(this%jas(ipos)) == 2) then
772  faheight = min(topn, topm) - max(botn, botm)
773  else
774  faheight = dhalf * ((topn - botn) + (topm - botm))
775  end if
776  if (this%hydchr(ihfb) > dzero) then
777  fawidth = this%hwva(this%jas(ipos))
778  condhfb = this%hydchr(ihfb) * &
779  fawidth * faheight
780  cond = cond * condhfb / (cond + condhfb)
781  else
782  cond = -cond * this%hydchr(ihfb)
783  end if
784  this%condsat(this%jas(ipos)) = cond
785  end if
786  end do
787  end subroutine condsat_modify
788 
789 end module gwfhfbmodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:83
subroutine hfb_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
Fill matrix terms.
Definition: gwf-hfb.f90:179
subroutine check_data(this)
Check for hfb's between two unconnected cells and write a warning.
Definition: gwf-hfb.f90:669
subroutine hfb_da(this)
Deallocate memory.
Definition: gwf-hfb.f90:422
subroutine hfb_cq(this, hnew, flowja)
flowja will automatically include the effects of the hfb for confined and newton cases when xt3d is n...
Definition: gwf-hfb.f90:323
subroutine, public hfb_cr(hfbobj, name_model, input_mempath, inunit, iout)
Create a new hfb object.
Definition: gwf-hfb.f90:71
subroutine allocate_scalars(this)
Allocate package scalars.
Definition: gwf-hfb.f90:465
subroutine condsat_modify(this)
Modify condsat.
Definition: gwf-hfb.f90:745
subroutine condsat_reset(this)
Reset condsat to its value prior to being modified by hfb's.
Definition: gwf-hfb.f90:726
subroutine source_data(this)
@ brief source hfb period data
Definition: gwf-hfb.f90:565
subroutine source_dimensions(this)
@ brief Source hfb input options
Definition: gwf-hfb.f90:536
subroutine hfb_rp(this)
Check for new HFB stress period data.
Definition: gwf-hfb.f90:149
subroutine allocate_arrays(this)
Allocate package arrays.
Definition: gwf-hfb.f90:489
subroutine hfb_ar(this, ibound, xt3d, dis, invsc, vsc)
Allocate and read.
Definition: gwf-hfb.f90:95
subroutine source_options(this)
@ brief Source hfb input options
Definition: gwf-hfb.f90:512
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
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
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23