MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
gwfhfbmodule Module Reference

Data Types

type  gwfhfbtype
 

Functions/Subroutines

subroutine, public hfb_cr (hfbobj, name_model, inunit, iout)
 Create a new hfb object. More...
 
subroutine hfb_ar (this, ibound, xt3d, dis, invsc, vsc)
 Allocate and read. More...
 
subroutine hfb_rp (this)
 Check for new HFB stress period data. More...
 
subroutine hfb_fc (this, kiter, matrix_sln, idxglo, rhs, hnew)
 Fill matrix terms. More...
 
subroutine hfb_cq (this, hnew, flowja)
 flowja will automatically include the effects of the hfb for confined and newton cases when xt3d is not used. More...
 
subroutine hfb_da (this)
 Deallocate memory. More...
 
subroutine allocate_scalars (this)
 Allocate package scalars. More...
 
subroutine allocate_arrays (this)
 Allocate package arrays. More...
 
subroutine read_options (this)
 Read a hfb options block. More...
 
subroutine read_dimensions (this)
 Read the dimensions for this package. More...
 
subroutine read_data (this)
 Read HFB period block. More...
 
subroutine check_data (this)
 Check for hfb's between two unconnected cells and write a warning. More...
 
subroutine condsat_reset (this)
 Reset condsat to its value prior to being modified by hfb's. More...
 
subroutine condsat_modify (this)
 Modify condsat. More...
 

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwfhfbmodule::allocate_arrays ( class(gwfhfbtype this)

Definition at line 521 of file gwf-hfb.f90.

522  ! -- modules
524  ! -- dummy
525  class(GwfHfbType) :: this
526  ! -- local
527  integer(I4B) :: ihfb
528  !
529  call mem_allocate(this%noden, this%maxhfb, 'NODEN', this%memoryPath)
530  call mem_allocate(this%nodem, this%maxhfb, 'NODEM', this%memoryPath)
531  call mem_allocate(this%hydchr, this%maxhfb, 'HYDCHR', this%memoryPath)
532  call mem_allocate(this%idxloc, this%maxhfb, 'IDXLOC', this%memoryPath)
533  call mem_allocate(this%csatsav, this%maxhfb, 'CSATSAV', this%memoryPath)
534  call mem_allocate(this%condsav, this%maxhfb, 'CONDSAV', this%memoryPath)
535  !
536  ! -- initialize idxloc to 0
537  do ihfb = 1, this%maxhfb
538  this%idxloc(ihfb) = 0
539  end do

◆ allocate_scalars()

subroutine gwfhfbmodule::allocate_scalars ( class(gwfhfbtype this)

Definition at line 497 of file gwf-hfb.f90.

498  ! -- modules
500  ! -- dummy
501  class(GwfHfbType) :: this
502  !
503  ! -- allocate scalars in NumericalPackageType
504  call this%NumericalPackageType%allocate_scalars()
505  !
506  ! -- allocate scalars
507  call mem_allocate(this%maxhfb, 'MAXHFB', this%memoryPath)
508  call mem_allocate(this%nhfb, 'NHFB', this%memoryPath)
509  !
510  ! -- allocate flag for determining if vsc active
511  call mem_allocate(this%ivsc, 'IVSC', this%memoryPath)
512  !
513  ! -- initialize
514  this%maxhfb = 0
515  this%nhfb = 0
516  this%ivsc = 0

◆ check_data()

subroutine gwfhfbmodule::check_data ( class(gwfhfbtype this)

Store ipos in idxloc

Definition at line 711 of file gwf-hfb.f90.

712  ! -- modules
713  use constantsmodule, only: linelength
715  ! -- dummy
716  class(GwfHfbType) :: this
717  ! -- local
718  integer(I4B) :: ihfb, n, m
719  integer(I4B) :: ipos
720  character(len=LINELENGTH) :: nodenstr, nodemstr
721  character(len=LINELENGTH) :: errmsg
722  logical :: found
723  ! -- formats
724  character(len=*), parameter :: fmterr = "(1x, 'HFB no. ',i0, &
725  &' is between two unconnected cells: ', a, ' and ', a)"
726  character(len=*), parameter :: fmtverr = "(1x, 'HFB no. ',i0, &
727  &' is between two cells not horizontally connected: ', a, ' and ', a)"
728  !
729  do ihfb = 1, this%nhfb
730  n = this%noden(ihfb)
731  m = this%nodem(ihfb)
732  found = .false.
733  do ipos = this%ia(n) + 1, this%ia(n + 1) - 1
734  if (m == this%ja(ipos)) then
735  found = .true.
736  this%idxloc(ihfb) = ipos
737  exit
738  end if
739  end do
740  !
741  ! -- check to make sure cells are connected
742  if (.not. found) then
743  call this%dis%noder_to_string(n, nodenstr)
744  call this%dis%noder_to_string(m, nodemstr)
745  write (errmsg, fmterr) ihfb, trim(adjustl(nodenstr)), &
746  trim(adjustl(nodemstr))
747  call store_error(errmsg)
748  else
749  !
750  ! -- check to make sure cells are not vertically connected
751  ipos = this%idxloc(ihfb)
752  if (this%ihc(this%jas(ipos)) == 0) then
753  call this%dis%noder_to_string(n, nodenstr)
754  call this%dis%noder_to_string(m, nodemstr)
755  write (errmsg, fmtverr) ihfb, trim(adjustl(nodenstr)), &
756  trim(adjustl(nodemstr))
757  call store_error(errmsg)
758  end if
759  end if
760  end do
761  !
762  ! -- Stop if errors detected
763  if (count_errors() > 0) then
764  call store_error_unit(this%inunit)
765  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
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
Here is the call graph for this function:

◆ condsat_modify()

subroutine gwfhfbmodule::condsat_modify ( class(gwfhfbtype this)
private

Modify condsat for the following conditions:

  1. If Newton is active
  2. If icelltype for n and icelltype for m is 0

Definition at line 789 of file gwf-hfb.f90.

790  ! -- modules
791  use constantsmodule, only: dhalf, dzero
792  ! -- dummy
793  class(GwfHfbType) :: this
794  ! -- local
795  integer(I4B) :: ihfb, n, m
796  integer(I4B) :: ipos
797  real(DP) :: cond, condhfb
798  real(DP) :: fawidth, faheight
799  real(DP) :: topn, topm, botn, botm
800  !
801  do ihfb = 1, this%nhfb
802  ipos = this%idxloc(ihfb)
803  cond = this%condsat(this%jas(ipos))
804  this%csatsav(ihfb) = cond
805  n = this%noden(ihfb)
806  m = this%nodem(ihfb)
807  !
808  if (this%inewton == 1 .or. &
809  (this%icelltype(n) == 0 .and. this%icelltype(m) == 0)) then
810  !
811  ! -- Calculate hfb conductance
812  topn = this%top(n)
813  topm = this%top(m)
814  botn = this%bot(n)
815  botm = this%bot(m)
816  if (this%ihc(this%jas(ipos)) == 2) then
817  faheight = min(topn, topm) - max(botn, botm)
818  else
819  faheight = dhalf * ((topn - botn) + (topm - botm))
820  end if
821  if (this%hydchr(ihfb) > dzero) then
822  fawidth = this%hwva(this%jas(ipos))
823  condhfb = this%hydchr(ihfb) * &
824  fawidth * faheight
825  cond = cond * condhfb / (cond + condhfb)
826  else
827  cond = -cond * this%hydchr(ihfb)
828  end if
829  this%condsat(this%jas(ipos)) = cond
830  end if
831  end do
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ condsat_reset()

subroutine gwfhfbmodule::condsat_reset ( class(gwfhfbtype this)

Definition at line 770 of file gwf-hfb.f90.

771  ! -- dummy
772  class(GwfHfbType) :: this
773  ! -- local
774  integer(I4B) :: ihfb
775  integer(I4B) :: ipos
776  !
777  do ihfb = 1, this%nhfb
778  ipos = this%idxloc(ihfb)
779  this%condsat(this%jas(ipos)) = this%csatsav(ihfb)
780  end do

◆ hfb_ar()

subroutine gwfhfbmodule::hfb_ar ( class(gwfhfbtype this,
integer(i4b), dimension(:), pointer, contiguous  ibound,
type(xt3dtype), pointer  xt3d,
class(disbasetype), intent(inout), pointer  dis,
integer(i4b), pointer  invsc,
type(gwfvsctype), intent(in), pointer  vsc 
)
private
Parameters
[in,out]disdiscretization package
invscindicates if viscosity package is active
[in]vscviscosity package

Definition at line 95 of file gwf-hfb.f90.

96  ! -- modules
99  ! -- dummy
100  class(GwfHfbType) :: this
101  integer(I4B), dimension(:), pointer, contiguous :: ibound
102  type(Xt3dType), pointer :: xt3d
103  class(DisBaseType), pointer, intent(inout) :: dis !< discretization package
104  integer(I4B), pointer :: invsc !< indicates if viscosity package is active
105  type(GwfVscType), pointer, intent(in) :: vsc !< viscosity package
106  ! -- formats
107  character(len=*), parameter :: fmtheader = &
108  "(1x, /1x, 'HFB -- HORIZONTAL FLOW BARRIER PACKAGE, VERSION 8, ', &
109  &'4/24/2015 INPUT READ FROM UNIT ', i4, //)"
110  !
111  ! -- Print a message identifying the node property flow package.
112  write (this%iout, fmtheader) this%inunit
113  !
114  ! -- Set pointers
115  this%dis => dis
116  this%ibound => ibound
117  this%xt3d => xt3d
118  !
119  call mem_setptr(this%icelltype, 'ICELLTYPE', &
120  create_mem_path(this%name_model, 'NPF'))
121  call mem_setptr(this%ihc, 'IHC', create_mem_path(this%name_model, 'CON'))
122  call mem_setptr(this%ia, 'IA', create_mem_path(this%name_model, 'CON'))
123  call mem_setptr(this%ja, 'JA', create_mem_path(this%name_model, 'CON'))
124  call mem_setptr(this%jas, 'JAS', create_mem_path(this%name_model, 'CON'))
125  call mem_setptr(this%isym, 'ISYM', create_mem_path(this%name_model, 'CON'))
126  call mem_setptr(this%condsat, 'CONDSAT', create_mem_path(this%name_model, &
127  'NPF'))
128  call mem_setptr(this%top, 'TOP', create_mem_path(this%name_model, 'DIS'))
129  call mem_setptr(this%bot, 'BOT', create_mem_path(this%name_model, 'DIS'))
130  call mem_setptr(this%hwva, 'HWVA', create_mem_path(this%name_model, 'CON'))
131  !
132  call this%read_options()
133  call this%read_dimensions()
134  call this%allocate_arrays()
135  !
136  ! -- If vsc package active, set ivsc
137  if (invsc /= 0) then
138  this%ivsc = 1
139  this%vsc => vsc
140  !
141  ! -- Notify user via listing file viscosity accounted for by HFB package
142  write (this%iout, '(/1x,a,a)') 'Viscosity active in ', &
143  trim(this%filtyp)//' Package calculations: '//trim(adjustl(this%packName))
144  end if
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Here is the call graph for this function:

◆ hfb_cq()

subroutine gwfhfbmodule::hfb_cq ( class(gwfhfbtype this,
real(dp), dimension(:), intent(inout)  hnew,
real(dp), dimension(:), intent(inout)  flowja 
)

This method recalculates flowja for the other cases.

Definition at line 355 of file gwf-hfb.f90.

356  ! -- modules
357  use constantsmodule, only: dhalf, dzero, done
358  ! -- dummy
359  class(GwfHfbType) :: this
360  real(DP), intent(inout), dimension(:) :: hnew
361  real(DP), intent(inout), dimension(:) :: flowja
362  ! -- local
363  integer(I4B) :: ihfb, n, m
364  integer(I4B) :: ipos
365  real(DP) :: qnm
366  real(DP) :: cond
367  integer(I4B) :: ixt3d
368  real(DP) :: condhfb
369  real(DP) :: fawidth, faheight
370  real(DP) :: topn, topm, botn, botm
371  real(DP) :: viscratio
372  !
373  ! -- initialize viscratio
374  viscratio = done
375  !
376  if (associated(this%xt3d%ixt3d)) then
377  ixt3d = this%xt3d%ixt3d
378  else
379  ixt3d = 0
380  end if
381  !
382  if (ixt3d > 0) then
383  !
384  do ihfb = 1, this%nhfb
385  n = min(this%noden(ihfb), this%nodem(ihfb))
386  m = max(this%noden(ihfb), this%nodem(ihfb))
387  ! -- Skip if either cell is inactive.
388  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
389  !!! if(this%icelltype(n) == 1 .or. this%icelltype(m) == 1) then
390  if (this%ivsc /= 0) then
391  call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
392  end if
393  !
394  ! -- Compute scale factor for hfb correction
395  if (this%hydchr(ihfb) > dzero) then
396  if (this%inewton == 0) then
397  ipos = this%idxloc(ihfb)
398  topn = this%top(n)
399  topm = this%top(m)
400  botn = this%bot(n)
401  botm = this%bot(m)
402  if (this%icelltype(n) == 1) then
403  if (hnew(n) < topn) topn = hnew(n)
404  end if
405  if (this%icelltype(m) == 1) then
406  if (hnew(m) < topm) topm = hnew(m)
407  end if
408  if (this%ihc(this%jas(ipos)) == 2) then
409  faheight = min(topn, topm) - max(botn, botm)
410  else
411  faheight = dhalf * ((topn - botn) + (topm - botm))
412  end if
413  fawidth = this%hwva(this%jas(ipos))
414  condhfb = this%hydchr(ihfb) * viscratio * &
415  fawidth * faheight
416  else
417  condhfb = this%hydchr(ihfb)
418  end if
419  else
420  condhfb = this%hydchr(ihfb)
421  end if
422  ! -- Make hfb corrections for xt3d
423  call this%xt3d%xt3d_flowjahfb(n, m, hnew, flowja, condhfb)
424  end do
425  !
426  else
427  !
428  ! -- Recalculate flowja for non-newton unconfined.
429  if (this%inewton == 0) then
430  do ihfb = 1, this%nhfb
431  n = this%noden(ihfb)
432  m = this%nodem(ihfb)
433  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
434  if (this%icelltype(n) == 1 .or. this%icelltype(m) == 1 .or. &
435  this%ivsc /= 0) then
436  ipos = this%dis%con%getjaindex(n, m)
437  !
438  ! -- condsav already accnts for visc adjustment
439  cond = this%condsav(ihfb)
440  qnm = cond * (hnew(m) - hnew(n))
441  flowja(ipos) = qnm
442  ipos = this%dis%con%getjaindex(m, n)
443  flowja(ipos) = -qnm
444  !
445  end if
446  end do
447  end if
448  !
449  end if
real(dp), parameter done
real constant 1
Definition: Constants.f90:76

◆ hfb_cr()

subroutine, public gwfhfbmodule::hfb_cr ( type(gwfhfbtype), pointer  hfbobj,
character(len=*), intent(in)  name_model,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Definition at line 69 of file gwf-hfb.f90.

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)
Here is the caller graph for this function:

◆ hfb_da()

subroutine gwfhfbmodule::hfb_da ( class(gwfhfbtype this)

Definition at line 454 of file gwf-hfb.f90.

455  ! -- modules
457  ! -- dummy
458  class(GwfHfbType) :: this
459  !
460  ! -- Scalars
461  call mem_deallocate(this%maxhfb)
462  call mem_deallocate(this%nhfb)
463  call mem_deallocate(this%ivsc)
464  !
465  ! -- Arrays
466  if (this%inunit > 0) then
467  call mem_deallocate(this%noden)
468  call mem_deallocate(this%nodem)
469  call mem_deallocate(this%hydchr)
470  call mem_deallocate(this%idxloc)
471  call mem_deallocate(this%csatsav)
472  call mem_deallocate(this%condsav)
473  end if
474  !
475  ! -- deallocate parent
476  call this%NumericalPackageType%da()
477  !
478  ! -- nullify pointers
479  this%xt3d => null()
480  this%inewton => null()
481  this%ibound => null()
482  this%icelltype => null()
483  this%ihc => null()
484  this%ia => null()
485  this%ja => null()
486  this%jas => null()
487  this%isym => null()
488  this%condsat => null()
489  this%top => null()
490  this%bot => null()
491  this%hwva => null()
492  this%vsc => null()

◆ hfb_fc()

subroutine gwfhfbmodule::hfb_fc ( class(gwfhfbtype this,
integer(i4b)  kiter,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs,
real(dp), dimension(:), intent(inout)  hnew 
)

Fill amatsln for the following conditions:

  1. XT3D OR
  2. Not Newton, and
  3. Cell type n is convertible or cell type m is convertible

Definition at line 211 of file gwf-hfb.f90.

212  ! -- modules
213  use constantsmodule, only: dhalf, dzero, done
214  ! -- dummy
215  class(GwfHfbType) :: this
216  integer(I4B) :: kiter
217  class(MatrixBaseType), pointer :: matrix_sln
218  integer(I4B), intent(in), dimension(:) :: idxglo
219  real(DP), intent(inout), dimension(:) :: rhs
220  real(DP), intent(inout), dimension(:) :: hnew
221  ! -- local
222  integer(I4B) :: nodes, nja
223  integer(I4B) :: ihfb, n, m
224  integer(I4B) :: ipos
225  integer(I4B) :: idiag, isymcon
226  integer(I4B) :: ixt3d
227  real(DP) :: cond, condhfb, aterm
228  real(DP) :: fawidth, faheight
229  real(DP) :: topn, topm, botn, botm
230  real(DP) :: viscratio
231  !
232  ! -- initialize variables
233  viscratio = done
234  nodes = this%dis%nodes
235  nja = this%dis%con%nja
236  if (associated(this%xt3d%ixt3d)) then
237  ixt3d = this%xt3d%ixt3d
238  else
239  ixt3d = 0
240  end if
241  !
242  if (ixt3d > 0) then
243  !
244  do ihfb = 1, this%nhfb
245  n = min(this%noden(ihfb), this%nodem(ihfb))
246  m = max(this%noden(ihfb), this%nodem(ihfb))
247  ! -- Skip if either cell is inactive.
248  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
249  !!! if(this%icelltype(n) == 1 .or. this%icelltype(m) == 1) then
250  if (this%ivsc /= 0) then
251  call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
252  end if
253  ! -- Compute scale factor for hfb correction
254  if (this%hydchr(ihfb) > dzero) then
255  if (this%inewton == 0) then
256  ipos = this%idxloc(ihfb)
257  topn = this%top(n)
258  topm = this%top(m)
259  botn = this%bot(n)
260  botm = this%bot(m)
261  if (this%icelltype(n) == 1) then
262  if (hnew(n) < topn) topn = hnew(n)
263  end if
264  if (this%icelltype(m) == 1) then
265  if (hnew(m) < topm) topm = hnew(m)
266  end if
267  if (this%ihc(this%jas(ipos)) == 2) then
268  faheight = min(topn, topm) - max(botn, botm)
269  else
270  faheight = dhalf * ((topn - botn) + (topm - botm))
271  end if
272  fawidth = this%hwva(this%jas(ipos))
273  condhfb = this%hydchr(ihfb) * viscratio * &
274  fawidth * faheight
275  else
276  condhfb = this%hydchr(ihfb) * viscratio
277  end if
278  else
279  condhfb = this%hydchr(ihfb)
280  end if
281  ! -- Make hfb corrections for xt3d
282  call this%xt3d%xt3d_fhfb(kiter, nodes, nja, matrix_sln, idxglo, &
283  rhs, hnew, n, m, condhfb)
284  end do
285  !
286  else
287  !
288  ! -- For Newton, the effect of the barrier is included in condsat.
289  if (this%inewton == 0) then
290  do ihfb = 1, this%nhfb
291  ipos = this%idxloc(ihfb)
292  aterm = matrix_sln%get_value_pos(idxglo(ipos))
293  n = this%noden(ihfb)
294  m = this%nodem(ihfb)
295  if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
296  !
297  if (this%ivsc /= 0) then
298  call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
299  end if
300  !
301  if (this%icelltype(n) == 1 .or. this%icelltype(m) == 1 .or. &
302  this%ivsc /= 0) then
303  !
304  ! -- Calculate hfb conductance
305  topn = this%top(n)
306  topm = this%top(m)
307  botn = this%bot(n)
308  botm = this%bot(m)
309  if (this%icelltype(n) == 1) then
310  if (hnew(n) < topn) topn = hnew(n)
311  end if
312  if (this%icelltype(m) == 1) then
313  if (hnew(m) < topm) topm = hnew(m)
314  end if
315  if (this%ihc(this%jas(ipos)) == 2) then
316  faheight = min(topn, topm) - max(botn, botm)
317  else
318  faheight = dhalf * ((topn - botn) + (topm - botm))
319  end if
320  if (this%hydchr(ihfb) > dzero) then
321  fawidth = this%hwva(this%jas(ipos))
322  condhfb = this%hydchr(ihfb) * viscratio * &
323  fawidth * faheight
324  cond = aterm * condhfb / (aterm + condhfb)
325  else
326  cond = -aterm * this%hydchr(ihfb)
327  end if
328  !
329  ! -- Save cond for budget calculation
330  this%condsav(ihfb) = cond
331  !
332  ! -- Fill row n diag and off diag
333  idiag = this%ia(n)
334  call matrix_sln%add_value_pos(idxglo(idiag), aterm - cond)
335  call matrix_sln%set_value_pos(idxglo(ipos), cond)
336  !
337  ! -- Fill row m diag and off diag
338  isymcon = this%isym(ipos)
339  idiag = this%ia(m)
340  call matrix_sln%add_value_pos(idxglo(idiag), aterm - cond)
341  call matrix_sln%set_value_pos(idxglo(isymcon), cond)
342  !
343  end if
344  end do
345  end if
346  !
347  end if

◆ hfb_rp()

subroutine gwfhfbmodule::hfb_rp ( class(gwfhfbtype this)

Definition at line 149 of file gwf-hfb.f90.

150  ! -- modules
151  use constantsmodule, only: linelength
153  use tdismodule, only: kper, nper
154  ! -- dummy
155  class(GwfHfbType) :: this
156  ! -- local
157  character(len=LINELENGTH) :: line, errmsg
158  integer(I4B) :: ierr
159  logical :: isfound
160  ! -- formats
161  character(len=*), parameter :: fmtblkerr = &
162  &"('Error. Looking for BEGIN PERIOD iper. Found ', a, ' instead.')"
163  character(len=*), parameter :: fmtlsp = &
164  &"(1X,/1X,'REUSING ',A,'S FROM LAST STRESS PERIOD')"
165  !
166  ! -- Set ionper to the stress period number for which a new block of data
167  ! will be read.
168  if (this%ionper < kper) then
169  !
170  ! -- get period block
171  call this%parser%GetBlock('PERIOD', isfound, ierr, &
172  supportopenclose=.true., &
173  blockrequired=.false.)
174  if (isfound) then
175  !
176  ! -- read ionper and check for increasing period numbers
177  call this%read_check_ionper()
178  else
179  !
180  ! -- PERIOD block not found
181  if (ierr < 0) then
182  ! -- End of file found; data applies for remainder of simulation.
183  this%ionper = nper + 1
184  else
185  ! -- Found invalid block
186  call this%parser%GetCurrentLine(line)
187  write (errmsg, fmtblkerr) adjustl(trim(line))
188  call store_error(errmsg)
189  call this%parser%StoreErrorUnit()
190  end if
191  end if
192  end if
193  !
194  if (this%ionper == kper) then
195  call this%condsat_reset()
196  call this%read_data()
197  call this%condsat_modify()
198  else
199  write (this%iout, fmtlsp) 'HFB'
200  end if
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
Here is the call graph for this function:

◆ read_data()

subroutine gwfhfbmodule::read_data ( class(gwfhfbtype this)

Data are in form of: L, IROW1, ICOL1, IROW2, ICOL2, HYDCHR or for unstructured: N1, N2, HYDCHR

Definition at line 639 of file gwf-hfb.f90.

640  ! -- modules
641  use constantsmodule, only: linelength
643  use tdismodule, only: kper
644  ! -- dummy
645  class(GwfHfbType) :: this
646  ! -- local
647  character(len=LINELENGTH) :: nodenstr, nodemstr, cellidm, cellidn
648  integer(I4B) :: ihfb, nerr
649  logical :: endOfBlock
650  ! -- formats
651  character(len=*), parameter :: fmthfb = "(i10, 2a10, 1(1pg15.6))"
652  !
653  write (this%iout, '(//,1x,a)') 'READING HFB DATA'
654  if (this%iprpak > 0) then
655  write (this%iout, '(3a10, 1a15)') 'HFB NUM', 'CELL1', 'CELL2', &
656  'HYDCHR'
657  end if
658  !
659  ihfb = 0
660  this%nhfb = 0
661  readloop: do
662  !
663  ! -- Check for END of block
664  call this%parser%GetNextLine(endofblock)
665  if (endofblock) exit
666  !
667  ! -- Reset lloc and read noden, nodem, and hydchr
668  ihfb = ihfb + 1
669  if (ihfb > this%maxhfb) then
670  call store_error('MAXHFB not large enough.')
671  call this%parser%StoreErrorUnit()
672  end if
673  call this%parser%GetCellid(this%dis%ndim, cellidn)
674  this%noden(ihfb) = this%dis%noder_from_cellid(cellidn, &
675  this%parser%iuactive, &
676  this%iout)
677  call this%parser%GetCellid(this%dis%ndim, cellidm)
678  this%nodem(ihfb) = this%dis%noder_from_cellid(cellidm, &
679  this%parser%iuactive, &
680  this%iout)
681  this%hydchr(ihfb) = this%parser%GetDouble()
682  !
683  ! -- Print input if requested
684  if (this%iprpak /= 0) then
685  call this%dis%noder_to_string(this%noden(ihfb), nodenstr)
686  call this%dis%noder_to_string(this%nodem(ihfb), nodemstr)
687  write (this%iout, fmthfb) ihfb, trim(adjustl(nodenstr)), &
688  trim(adjustl(nodemstr)), this%hydchr(ihfb)
689  end if
690  !
691  this%nhfb = ihfb
692  end do readloop
693  !
694  ! -- Stop if errors
695  nerr = count_errors()
696  if (nerr > 0) then
697  call store_error('Errors encountered in HFB input file.')
698  call this%parser%StoreErrorUnit()
699  end if
700  !
701  write (this%iout, '(3x,i0,a,i0)') this%nhfb, &
702  ' HFBs READ FOR STRESS PERIOD ', kper
703  call this%check_data()
704  write (this%iout, '(1x,a)') 'END READING HFB DATA'
Here is the call graph for this function:

◆ read_dimensions()

subroutine gwfhfbmodule::read_dimensions ( class(gwfhfbtype), intent(inout)  this)

Definition at line 584 of file gwf-hfb.f90.

585  use constantsmodule, only: linelength
587  ! -- dummy
588  class(GwfHfbType), intent(inout) :: this
589  ! -- local
590  character(len=LINELENGTH) :: errmsg, keyword
591  integer(I4B) :: ierr
592  logical :: isfound, endOfBlock
593  !
594  ! -- get dimensions block
595  call this%parser%GetBlock('DIMENSIONS', isfound, ierr, &
596  supportopenclose=.true.)
597  !
598  ! -- parse dimensions block if detected
599  if (isfound) then
600  write (this%iout, '(/1x,a)') 'PROCESSING HFB DIMENSIONS'
601  do
602  call this%parser%GetNextLine(endofblock)
603  if (endofblock) exit
604  call this%parser%GetStringCaps(keyword)
605  select case (keyword)
606  case ('MAXHFB')
607  this%maxhfb = this%parser%GetInteger()
608  write (this%iout, '(4x,a,i7)') 'MAXHFB = ', this%maxhfb
609  case default
610  write (errmsg, '(a,a)') &
611  'Unknown HFB dimension: ', trim(keyword)
612  call store_error(errmsg)
613  call this%parser%StoreErrorUnit()
614  end select
615  end do
616  !
617  write (this%iout, '(1x,a)') 'END OF HFB DIMENSIONS'
618  else
619  call store_error('Required DIMENSIONS block not found.')
620  call this%parser%StoreErrorUnit()
621  end if
622  !
623  ! -- verify dimensions were set
624  if (this%maxhfb <= 0) then
625  write (errmsg, '(a)') &
626  'MAXHFB must be specified with value greater than zero.'
627  call store_error(errmsg)
628  call this%parser%StoreErrorUnit()
629  end if
Here is the call graph for this function:

◆ read_options()

subroutine gwfhfbmodule::read_options ( class(gwfhfbtype this)

Definition at line 544 of file gwf-hfb.f90.

545  ! -- modules
546  use constantsmodule, only: linelength
548  ! -- dummy
549  class(GwfHfbType) :: this
550  ! -- local
551  character(len=LINELENGTH) :: errmsg, keyword
552  integer(I4B) :: ierr
553  logical :: isfound, endOfBlock
554  !
555  ! -- get options block
556  call this%parser%GetBlock('OPTIONS', isfound, ierr, &
557  supportopenclose=.true., blockrequired=.false.)
558  !
559  ! -- parse options block if detected
560  if (isfound) then
561  write (this%iout, '(1x,a)') 'PROCESSING HFB OPTIONS'
562  do
563  call this%parser%GetNextLine(endofblock)
564  if (endofblock) exit
565  call this%parser%GetStringCaps(keyword)
566  select case (keyword)
567  case ('PRINT_INPUT')
568  this%iprpak = 1
569  write (this%iout, '(4x,a)') &
570  'THE LIST OF HFBS WILL BE PRINTED.'
571  case default
572  write (errmsg, '(a,a)') 'Unknown HFB option: ', &
573  trim(keyword)
574  call store_error(errmsg)
575  call this%parser%StoreErrorUnit()
576  end select
577  end do
578  write (this%iout, '(1x,a)') 'END OF HFB OPTIONS'
579  end if
Here is the call graph for this function: