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