MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
flowmodelinterfacemodule Module Reference

Data Types

type  flowmodelinterfacetype
 

Functions/Subroutines

subroutine fmi_df (this, dis, idryinactive)
 Define the flow model interface. More...
 
subroutine fmi_ar (this, ibound)
 Allocate the package. More...
 
subroutine fmi_da (this)
 Deallocate variables. More...
 
subroutine allocate_scalars (this)
 Allocate scalars. More...
 
subroutine allocate_arrays (this, nodes)
 Allocate arrays. More...
 
subroutine read_options (this)
 Read options from input file. More...
 
subroutine read_packagedata (this)
 Read packagedata block from input file. More...
 
subroutine initialize_bfr (this)
 Initialize the budget file reader. More...
 
subroutine advance_bfr (this)
 Advance the budget file reader. More...
 
subroutine finalize_bfr (this)
 Finalize the budget file reader. More...
 
subroutine initialize_hfr (this)
 Initialize the head file reader. More...
 
subroutine advance_hfr (this)
 Advance the head file reader. More...
 
subroutine finalize_hfr (this)
 Finalize the head file reader. More...
 
subroutine initialize_gwfterms_from_bfr (this)
 Initialize gwf terms from budget file. More...
 
subroutine initialize_gwfterms_from_gwfbndlist (this)
 Initialize gwf terms from a GWF exchange. More...
 
subroutine allocate_gwfpackages (this, ngwfterms)
 Allocate budget packages. More...
 
subroutine deallocate_gwfpackages (this)
 Deallocate memory in the gwfpackages array. More...
 
subroutine get_package_index (this, name, idx)
 Find the package index for the package with the given name. More...
 

Function/Subroutine Documentation

◆ advance_bfr()

subroutine flowmodelinterfacemodule::advance_bfr ( class(flowmodelinterfacetype this)
private

Advance the budget file reader by reading the next chunk of information for the current time step and stress period.

Definition at line 465 of file FlowModelInterface.f90.

466  ! -- modules
467  use tdismodule, only: kstp, kper
468  ! -- dummy
469  class(FlowModelInterfaceType) :: this
470  ! -- local
471  logical :: success
472  integer(I4B) :: n
473  integer(I4B) :: ipos
474  integer(I4B) :: nu, nr
475  integer(I4B) :: ip, i
476  logical :: readnext
477  ! -- format
478  character(len=*), parameter :: fmtkstpkper = &
479  "(1x,/1x,'FMI READING BUDGET TERMS &
480  &FOR KSTP ', i0, ' KPER ', i0)"
481  character(len=*), parameter :: fmtbudkstpkper = &
482  "(1x,/1x, 'FMI SETTING BUDGET TERMS &
483  &FOR KSTP ', i0, ' AND KPER ', &
484  &i0, ' TO BUDGET FILE TERMS FROM &
485  &KSTP ', i0, ' AND KPER ', i0)"
486  !
487  ! -- If the latest record read from the budget file is from a stress
488  ! -- period with only one time step, reuse that record (do not read a
489  ! -- new record) if the running model is still in that same stress period,
490  ! -- or if that record is the last one in the budget file.
491  readnext = .true.
492  if (kstp * kper > 1) then
493  if (this%bfr%kstp == 1) then
494  if (this%bfr%kpernext == kper + 1) then
495  readnext = .false.
496  else if (this%bfr%endoffile) then
497  readnext = .false.
498  end if
499  else if (this%bfr%endoffile) then
500  write (errmsg, '(4x,a)') 'REACHED END OF GWF BUDGET &
501  &FILE BEFORE READING SUFFICIENT BUDGET INFORMATION FOR THIS &
502  &GWT SIMULATION.'
503  call store_error(errmsg)
504  call store_error_unit(this%iubud)
505  end if
506  end if
507  !
508  ! -- Read the next record
509  if (readnext) then
510  !
511  ! -- Write the current time step and stress period
512  write (this%iout, fmtkstpkper) kstp, kper
513  !
514  ! -- loop through the budget terms for this stress period
515  ! i is the counter for gwf flow packages
516  ip = 1
517  do n = 1, this%bfr%nbudterms
518  call this%bfr%read_record(success, this%iout)
519  if (.not. success) then
520  write (errmsg, '(4x,a)') 'GWF BUDGET READ NOT SUCCESSFUL'
521  call store_error(errmsg)
522  call store_error_unit(this%iubud)
523  end if
524  !
525  ! -- Ensure kper is same between model and budget file
526  if (kper /= this%bfr%kper) then
527  write (errmsg, '(4x,a)') 'PERIOD NUMBER IN BUDGET FILE &
528  &DOES NOT MATCH PERIOD NUMBER IN TRANSPORT MODEL. IF THERE &
529  &IS MORE THAN ONE TIME STEP IN THE BUDGET FILE FOR A GIVEN &
530  &STRESS PERIOD, BUDGET FILE TIME STEPS MUST MATCH GWT MODEL &
531  &TIME STEPS ONE-FOR-ONE IN THAT STRESS PERIOD.'
532  call store_error(errmsg)
533  call store_error_unit(this%iubud)
534  end if
535  !
536  ! -- if budget file kstp > 1, then kstp must match
537  if (this%bfr%kstp > 1 .and. (kstp /= this%bfr%kstp)) then
538  write (errmsg, '(4x,a)') 'TIME STEP NUMBER IN BUDGET FILE &
539  &DOES NOT MATCH TIME STEP NUMBER IN TRANSPORT MODEL. IF THERE &
540  &IS MORE THAN ONE TIME STEP IN THE BUDGET FILE FOR A GIVEN STRESS &
541  &PERIOD, BUDGET FILE TIME STEPS MUST MATCH GWT MODEL TIME STEPS &
542  &ONE-FOR-ONE IN THAT STRESS PERIOD.'
543  call store_error(errmsg)
544  call store_error_unit(this%iubud)
545  end if
546  !
547  ! -- parse based on the type of data, and compress all user node
548  ! numbers into reduced node numbers
549  select case (trim(adjustl(this%bfr%budtxt)))
550  case ('FLOW-JA-FACE')
551  !
552  ! -- bfr%flowja contains only reduced connections so there is
553  ! a one-to-one match with this%gwfflowja
554  do ipos = 1, size(this%bfr%flowja)
555  this%gwfflowja(ipos) = this%bfr%flowja(ipos)
556  end do
557  case ('DATA-SPDIS')
558  do i = 1, this%bfr%nlist
559  nu = this%bfr%nodesrc(i)
560  nr = this%dis%get_nodenumber(nu, 0)
561  if (nr <= 0) cycle
562  this%gwfspdis(1, nr) = this%bfr%auxvar(1, i)
563  this%gwfspdis(2, nr) = this%bfr%auxvar(2, i)
564  this%gwfspdis(3, nr) = this%bfr%auxvar(3, i)
565  end do
566  case ('DATA-SAT')
567  do i = 1, this%bfr%nlist
568  nu = this%bfr%nodesrc(i)
569  nr = this%dis%get_nodenumber(nu, 0)
570  if (nr <= 0) cycle
571  this%gwfsat(nr) = this%bfr%auxvar(1, i)
572  end do
573  case ('STO-SS')
574  do nu = 1, this%dis%nodesuser
575  nr = this%dis%get_nodenumber(nu, 0)
576  if (nr <= 0) cycle
577  this%gwfstrgss(nr) = this%bfr%flow(nu)
578  end do
579  case ('STO-SY')
580  do nu = 1, this%dis%nodesuser
581  nr = this%dis%get_nodenumber(nu, 0)
582  if (nr <= 0) cycle
583  this%gwfstrgsy(nr) = this%bfr%flow(nu)
584  end do
585  case default
586  call this%gwfpackages(ip)%copy_values( &
587  this%bfr%nlist, &
588  this%bfr%nodesrc, &
589  this%bfr%flow, &
590  this%bfr%auxvar)
591  do i = 1, this%gwfpackages(ip)%nbound
592  nu = this%gwfpackages(ip)%nodelist(i)
593  nr = this%dis%get_nodenumber(nu, 0)
594  this%gwfpackages(ip)%nodelist(i) = nr
595  end do
596  ip = ip + 1
597  end select
598  end do
599  else
600  !
601  ! -- write message to indicate that flows are being reused
602  write (this%iout, fmtbudkstpkper) kstp, kper, this%bfr%kstp, this%bfr%kper
603  !
604  ! -- set the flag to indicate that flows were not updated
605  this%iflowsupdated = 0
606  end if
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
Here is the call graph for this function:

◆ advance_hfr()

subroutine flowmodelinterfacemodule::advance_hfr ( class(flowmodelinterfacetype this)
private

Definition at line 637 of file FlowModelInterface.f90.

638  ! -- modules
639  use tdismodule, only: kstp, kper
640  class(FlowModelInterfaceType) :: this
641  integer(I4B) :: nu, nr, i, ilay
642  integer(I4B) :: ncpl
643  real(DP) :: val
644  logical :: readnext
645  logical :: success
646  character(len=*), parameter :: fmtkstpkper = &
647  "(1x,/1x,'FMI READING HEAD FOR &
648  &KSTP ', i0, ' KPER ', i0)"
649  character(len=*), parameter :: fmthdskstpkper = &
650  "(1x,/1x, 'FMI SETTING HEAD FOR KSTP ', i0, ' AND KPER ', &
651  &i0, ' TO BINARY FILE HEADS FROM KSTP ', i0, ' AND KPER ', i0)"
652  !
653  ! -- If the latest record read from the head file is from a stress
654  ! -- period with only one time step, reuse that record (do not read a
655  ! -- new record) if the running model is still in that same stress period,
656  ! -- or if that record is the last one in the head file.
657  readnext = .true.
658  if (kstp * kper > 1) then
659  if (this%hfr%kstp == 1) then
660  if (this%hfr%kpernext == kper + 1) then
661  readnext = .false.
662  else if (this%hfr%endoffile) then
663  readnext = .false.
664  end if
665  else if (this%hfr%endoffile) then
666  write (errmsg, '(4x,a)') 'REACHED END OF GWF HEAD &
667  &FILE BEFORE READING SUFFICIENT HEAD INFORMATION FOR THIS &
668  &GWT SIMULATION.'
669  call store_error(errmsg)
670  call store_error_unit(this%iuhds)
671  end if
672  end if
673  !
674  ! -- Read the next record
675  if (readnext) then
676  !
677  ! -- write to list file that heads are being read
678  write (this%iout, fmtkstpkper) kstp, kper
679  !
680  ! -- loop through the layered heads for this time step
681  do ilay = 1, this%hfr%nlay
682  !
683  ! -- read next head chunk
684  call this%hfr%read_record(success, this%iout)
685  if (.not. success) then
686  write (errmsg, '(4x,a)') 'GWF HEAD READ NOT SUCCESSFUL'
687  call store_error(errmsg)
688  call store_error_unit(this%iuhds)
689  end if
690  !
691  ! -- Ensure kper is same between model and head file
692  if (kper /= this%hfr%kper) then
693  write (errmsg, '(4x,a)') 'PERIOD NUMBER IN HEAD FILE &
694  &DOES NOT MATCH PERIOD NUMBER IN TRANSPORT MODEL. IF THERE &
695  &IS MORE THAN ONE TIME STEP IN THE HEAD FILE FOR A GIVEN STRESS &
696  &PERIOD, HEAD FILE TIME STEPS MUST MATCH GWT MODEL TIME STEPS &
697  &ONE-FOR-ONE IN THAT STRESS PERIOD.'
698  call store_error(errmsg)
699  call store_error_unit(this%iuhds)
700  end if
701  !
702  ! -- if head file kstp > 1, then kstp must match
703  if (this%hfr%kstp > 1 .and. (kstp /= this%hfr%kstp)) then
704  write (errmsg, '(4x,a)') 'TIME STEP NUMBER IN HEAD FILE &
705  &DOES NOT MATCH TIME STEP NUMBER IN TRANSPORT MODEL. IF THERE &
706  &IS MORE THAN ONE TIME STEP IN THE HEAD FILE FOR A GIVEN STRESS &
707  &PERIOD, HEAD FILE TIME STEPS MUST MATCH GWT MODEL TIME STEPS &
708  &ONE-FOR-ONE IN THAT STRESS PERIOD.'
709  call store_error(errmsg)
710  call store_error_unit(this%iuhds)
711  end if
712  !
713  ! -- fill the head array for this layer and
714  ! compress into reduced form
715  ncpl = size(this%hfr%head)
716  do i = 1, ncpl
717  nu = (ilay - 1) * ncpl + i
718  nr = this%dis%get_nodenumber(nu, 0)
719  val = this%hfr%head(i)
720  if (nr > 0) this%gwfhead(nr) = val
721  end do
722  end do
723  else
724  write (this%iout, fmthdskstpkper) kstp, kper, this%hfr%kstp, this%hfr%kper
725  end if
Here is the call graph for this function:

◆ allocate_arrays()

subroutine flowmodelinterfacemodule::allocate_arrays ( class(flowmodelinterfacetype this,
integer(i4b), intent(in)  nodes 
)

Definition at line 241 of file FlowModelInterface.f90.

243  !modules
244  use constantsmodule, only: dzero
245  ! -- dummy
246  class(FlowModelInterfaceType) :: this
247  integer(I4B), intent(in) :: nodes
248  ! -- local
249  integer(I4B) :: n
250  !
251  ! -- Allocate ibdgwfsat0, which is an indicator array marking cells with
252  ! saturation greater than 0.0 with a value of 1
253  call mem_allocate(this%ibdgwfsat0, nodes, 'IBDGWFSAT0', this%memoryPath)
254  do n = 1, nodes
255  this%ibdgwfsat0(n) = 1
256  end do
257  !
258  ! -- Allocate differently depending on whether or not flows are
259  ! being read from a file.
260  if (this%flows_from_file) then
261  call mem_allocate(this%gwfflowja, this%dis%con%nja, &
262  'GWFFLOWJA', this%memoryPath)
263  call mem_allocate(this%gwfsat, nodes, 'GWFSAT', this%memoryPath)
264  call mem_allocate(this%gwfhead, nodes, 'GWFHEAD', this%memoryPath)
265  call mem_allocate(this%gwfspdis, 3, nodes, 'GWFSPDIS', this%memoryPath)
266  do n = 1, nodes
267  this%gwfsat(n) = done
268  this%gwfhead(n) = dzero
269  this%gwfspdis(:, n) = dzero
270  end do
271  do n = 1, size(this%gwfflowja)
272  this%gwfflowja(n) = dzero
273  end do
274  !
275  ! -- allocate and initialize storage arrays
276  if (this%igwfstrgss == 0) then
277  call mem_allocate(this%gwfstrgss, 1, 'GWFSTRGSS', this%memoryPath)
278  else
279  call mem_allocate(this%gwfstrgss, nodes, 'GWFSTRGSS', this%memoryPath)
280  end if
281  if (this%igwfstrgsy == 0) then
282  call mem_allocate(this%gwfstrgsy, 1, 'GWFSTRGSY', this%memoryPath)
283  else
284  call mem_allocate(this%gwfstrgsy, nodes, 'GWFSTRGSY', this%memoryPath)
285  end if
286  do n = 1, size(this%gwfstrgss)
287  this%gwfstrgss(n) = dzero
288  end do
289  do n = 1, size(this%gwfstrgsy)
290  this%gwfstrgsy(n) = dzero
291  end do
292  !
293  ! -- If there is no fmi package, then there are no flows at all or a
294  ! connected GWF model, so allocate gwfpackages to zero
295  if (this%inunit == 0) call this%allocate_gwfpackages(this%nflowpack)
296  end if
297  !
298  ! -- Return
299  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64

◆ allocate_gwfpackages()

subroutine flowmodelinterfacemodule::allocate_gwfpackages ( class(flowmodelinterfacetype this,
integer(i4b), intent(in)  ngwfterms 
)

gwfpackages is an array of PackageBudget objects. This routine allocates gwfpackages to the proper size and initializes some member variables.

Definition at line 914 of file FlowModelInterface.f90.

915  ! -- modules
916  use constantsmodule, only: lenmempath
918  ! -- dummy
919  class(FlowModelInterfaceType) :: this
920  integer(I4B), intent(in) :: ngwfterms
921  ! -- local
922  integer(I4B) :: n
923  character(len=LENMEMPATH) :: memPath
924  !
925  ! -- direct allocate
926  allocate (this%gwfpackages(ngwfterms))
927  allocate (this%flowpacknamearray(ngwfterms))
928  !
929  ! -- mem_allocate
930  call mem_allocate(this%igwfmvrterm, ngwfterms, 'IGWFMVRTERM', this%memoryPath)
931  !
932  ! -- initialize
933  this%nflowpack = ngwfterms
934  do n = 1, this%nflowpack
935  this%igwfmvrterm(n) = 0
936  this%flowpacknamearray(n) = ''
937  !
938  ! -- Create a mempath for each individual flow package data set
939  ! of the form, MODELNAME/FMI-FTn
940  write (mempath, '(a, i0)') trim(this%memoryPath)//'-FT', n
941  call this%gwfpackages(n)%initialize(mempath)
942  end do
943  !
944  ! -- return
945  return
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26

◆ allocate_scalars()

subroutine flowmodelinterfacemodule::allocate_scalars ( class(flowmodelinterfacetype this)

Definition at line 202 of file FlowModelInterface.f90.

203  ! -- modules
205  ! -- dummy
206  class(FlowModelInterfaceType) :: this
207  ! -- local
208  !
209  ! -- allocate scalars in NumericalPackageType
210  call this%NumericalPackageType%allocate_scalars()
211  !
212  ! -- Allocate
213  call mem_allocate(this%flows_from_file, 'FLOWS_FROM_FILE', this%memoryPath)
214  call mem_allocate(this%iflowsupdated, 'IFLOWSUPDATED', this%memoryPath)
215  call mem_allocate(this%igwfstrgss, 'IGWFSTRGSS', this%memoryPath)
216  call mem_allocate(this%igwfstrgsy, 'IGWFSTRGSY', this%memoryPath)
217  call mem_allocate(this%iubud, 'IUBUD', this%memoryPath)
218  call mem_allocate(this%iuhds, 'IUHDS', this%memoryPath)
219  call mem_allocate(this%iumvr, 'IUMVR', this%memoryPath)
220  call mem_allocate(this%nflowpack, 'NFLOWPACK', this%memoryPath)
221  call mem_allocate(this%idryinactive, "IDRYINACTIVE", this%memoryPath)
222  !
223  ! !
224  ! -- Initialize
225  this%flows_from_file = .true.
226  this%iflowsupdated = 1
227  this%igwfstrgss = 0
228  this%igwfstrgsy = 0
229  this%iubud = 0
230  this%iuhds = 0
231  this%iumvr = 0
232  this%nflowpack = 0
233  this%idryinactive = 1
234  !
235  ! -- Return
236  return

◆ deallocate_gwfpackages()

subroutine flowmodelinterfacemodule::deallocate_gwfpackages ( class(flowmodelinterfacetype this)

Definition at line 950 of file FlowModelInterface.f90.

951  ! -- modules
952  ! -- dummy
953  class(FlowModelInterfaceType) :: this
954  ! -- local
955  integer(I4B) :: n
956  !
957  ! -- initialize
958  do n = 1, this%nflowpack
959  call this%gwfpackages(n)%da()
960  end do
961  !
962  ! -- return
963  return

◆ finalize_bfr()

subroutine flowmodelinterfacemodule::finalize_bfr ( class(flowmodelinterfacetype this)

Definition at line 611 of file FlowModelInterface.f90.

612  ! -- modules
613  class(FlowModelInterfaceType) :: this
614  ! -- dummy
615  !
616  ! -- Finalize the budget file reader
617  call this%bfr%finalize()
618  !

◆ finalize_hfr()

subroutine flowmodelinterfacemodule::finalize_hfr ( class(flowmodelinterfacetype this)

Definition at line 730 of file FlowModelInterface.f90.

731  ! -- modules
732  class(FlowModelInterfaceType) :: this
733  ! -- dummy
734  !
735  ! -- Finalize the head file reader
736  close (this%iuhds)
737  !

◆ fmi_ar()

subroutine flowmodelinterfacemodule::fmi_ar ( class(flowmodelinterfacetype this,
integer(i4b), dimension(:), pointer, contiguous  ibound 
)

Definition at line 136 of file FlowModelInterface.f90.

137  ! -- modules
138  use simmodule, only: store_error
139  ! -- dummy
140  class(FlowModelInterfaceType) :: this
141  integer(I4B), dimension(:), pointer, contiguous :: ibound
142  !
143  ! -- store pointers to arguments that were passed in
144  this%ibound => ibound
145  !
146  ! -- Allocate arrays
147  call this%allocate_arrays(this%dis%nodes)
148  !
149  ! -- Return
150  return
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:

◆ fmi_da()

subroutine flowmodelinterfacemodule::fmi_da ( class(flowmodelinterfacetype this)

Definition at line 155 of file FlowModelInterface.f90.

156  ! -- modules
158  ! -- dummy
159  class(FlowModelInterfaceType) :: this
160  ! -- todo: finalize hfr and bfr either here or in a finalize routine
161  !
162  ! -- deallocate any memory stored with gwfpackages
163  call this%deallocate_gwfpackages()
164  !
165  ! -- deallocate fmi arrays
166  deallocate (this%gwfpackages)
167  deallocate (this%flowpacknamearray)
168  call mem_deallocate(this%igwfmvrterm)
169  call mem_deallocate(this%ibdgwfsat0)
170  !
171  if (this%flows_from_file) then
172  call mem_deallocate(this%gwfstrgss)
173  call mem_deallocate(this%gwfstrgsy)
174  end if
175  !
176  ! -- special treatment, these could be from mem_checkin
177  call mem_deallocate(this%gwfhead, 'GWFHEAD', this%memoryPath)
178  call mem_deallocate(this%gwfsat, 'GWFSAT', this%memoryPath)
179  call mem_deallocate(this%gwfspdis, 'GWFSPDIS', this%memoryPath)
180  call mem_deallocate(this%gwfflowja, 'GWFFLOWJA', this%memoryPath)
181  !
182  ! -- deallocate scalars
183  call mem_deallocate(this%flows_from_file)
184  call mem_deallocate(this%iflowsupdated)
185  call mem_deallocate(this%igwfstrgss)
186  call mem_deallocate(this%igwfstrgsy)
187  call mem_deallocate(this%iubud)
188  call mem_deallocate(this%iuhds)
189  call mem_deallocate(this%iumvr)
190  call mem_deallocate(this%nflowpack)
191  call mem_deallocate(this%idryinactive)
192  !
193  ! -- deallocate parent
194  call this%NumericalPackageType%da()
195  !
196  ! -- Return
197  return

◆ fmi_df()

subroutine flowmodelinterfacemodule::fmi_df ( class(flowmodelinterfacetype this,
class(disbasetype), intent(in), pointer  dis,
integer(i4b), intent(in)  idryinactive 
)
private

Definition at line 76 of file FlowModelInterface.f90.

77  ! -- modules
78  use simmodule, only: store_error
79  ! -- dummy
80  class(FlowModelInterfaceType) :: this
81  class(DisBaseType), pointer, intent(in) :: dis
82  integer(I4B), intent(in) :: idryinactive
83  ! -- formats
84  character(len=*), parameter :: fmtfmi = &
85  "(1x,/1x,'FMI -- FLOW MODEL INTERFACE, VERSION 2, 8/17/2023', &
86  &' INPUT READ FROM UNIT ', i0, //)"
87  character(len=*), parameter :: fmtfmi0 = &
88  "(1x,/1x,'FMI -- FLOW MODEL INTERFACE,'&
89  &' VERSION 2, 8/17/2023')"
90  !
91  ! --print a message identifying the FMI package.
92  if (this%iout > 0) then
93  if (this%inunit /= 0) then
94  write (this%iout, fmtfmi) this%inunit
95  else
96  write (this%iout, fmtfmi0)
97  if (this%flows_from_file) then
98  write (this%iout, '(a)') ' FLOWS ARE ASSUMED TO BE ZERO.'
99  else
100  write (this%iout, '(a)') ' FLOWS PROVIDED BY A GWF MODEL IN THIS &
101  &SIMULATION'
102  end if
103  end if
104  end if
105  !
106  ! -- Store pointers
107  this%dis => dis
108  !
109  ! -- Read fmi options
110  if (this%inunit /= 0) then
111  call this%read_options()
112  end if
113  !
114  ! -- Read packagedata options
115  if (this%inunit /= 0 .and. this%flows_from_file) then
116  call this%read_packagedata()
117  call this%initialize_gwfterms_from_bfr()
118  end if
119  !
120  ! -- If GWF-Model exchange is active, setup flow terms
121  if (.not. this%flows_from_file) then
122  call this%initialize_gwfterms_from_gwfbndlist()
123  end if
124  !
125  ! -- Set flag that stops dry flows from being deactivated in a GWE
126  ! transport model since conduction will still be simulated.
127  ! 0: GWE (skip deactivation step); 1: GWT (default: use existing code)
128  this%idryinactive = idryinactive
129  !
130  ! -- Return
131  return
Here is the call graph for this function:

◆ get_package_index()

subroutine flowmodelinterfacemodule::get_package_index ( class(flowmodelinterfacetype this,
character(len=*), intent(in)  name,
integer(i4b), intent(inout)  idx 
)
private

Definition at line 968 of file FlowModelInterface.f90.

969  use bndmodule, only: bndtype, getbndfromlist
970  class(FlowModelInterfaceType) :: this
971  character(len=*), intent(in) :: name
972  integer(I4B), intent(inout) :: idx
973  ! -- local
974  integer(I4B) :: ip
975  !
976  ! -- Look through all the packages and return the index with name
977  idx = 0
978  do ip = 1, size(this%flowpacknamearray)
979  if (this%flowpacknamearray(ip) == name) then
980  idx = ip
981  exit
982  end if
983  end do
984  if (idx == 0) then
985  call store_error('Error in get_package_index. Could not find '//name, &
986  terminate=.true.)
987  end if
988  !
989  ! -- return
990  return
This module contains the base boundary package.
class(bndtype) function, pointer, public getbndfromlist(list, idx)
Get boundary from package list.
@ brief BndType
Here is the call graph for this function:

◆ initialize_bfr()

subroutine flowmodelinterfacemodule::initialize_bfr ( class(flowmodelinterfacetype this)

Definition at line 446 of file FlowModelInterface.f90.

447  ! -- modules
448  class(FlowModelInterfaceType) :: this
449  ! -- dummy
450  integer(I4B) :: ncrbud
451  !
452  ! -- Initialize the budget file reader
453  call this%bfr%initialize(this%iubud, this%iout, ncrbud)
454  !
455  ! -- todo: need to run through the budget terms
456  ! and do some checking

◆ initialize_gwfterms_from_bfr()

subroutine flowmodelinterfacemodule::initialize_gwfterms_from_bfr ( class(flowmodelinterfacetype this)
private

initialize terms and figure out how many different terms and packages are contained within the file

Definition at line 746 of file FlowModelInterface.f90.

747  ! -- modules
750  ! -- dummy
751  class(FlowModelInterfaceType) :: this
752  ! -- local
753  integer(I4B) :: nflowpack
754  integer(I4B) :: i, ip
755  integer(I4B) :: naux
756  logical :: found_flowja
757  logical :: found_dataspdis
758  logical :: found_datasat
759  logical :: found_stoss
760  logical :: found_stosy
761  integer(I4B), dimension(:), allocatable :: imap
762  !
763  ! -- Calculate the number of gwf flow packages
764  allocate (imap(this%bfr%nbudterms))
765  imap(:) = 0
766  nflowpack = 0
767  found_flowja = .false.
768  found_dataspdis = .false.
769  found_datasat = .false.
770  found_stoss = .false.
771  found_stosy = .false.
772  do i = 1, this%bfr%nbudterms
773  select case (trim(adjustl(this%bfr%budtxtarray(i))))
774  case ('FLOW-JA-FACE')
775  found_flowja = .true.
776  case ('DATA-SPDIS')
777  found_dataspdis = .true.
778  case ('DATA-SAT')
779  found_datasat = .true.
780  case ('STO-SS')
781  found_stoss = .true.
782  this%igwfstrgss = 1
783  case ('STO-SY')
784  found_stosy = .true.
785  this%igwfstrgsy = 1
786  case default
787  nflowpack = nflowpack + 1
788  imap(i) = 1
789  end select
790  end do
791  !
792  ! -- allocate gwfpackage arrays
793  call this%allocate_gwfpackages(nflowpack)
794  !
795  ! -- Copy the package name and aux names from budget file reader
796  ! to the gwfpackages derived-type variable
797  ip = 1
798  do i = 1, this%bfr%nbudterms
799  if (imap(i) == 0) cycle
800  call this%gwfpackages(ip)%set_name(this%bfr%dstpackagenamearray(i), &
801  this%bfr%budtxtarray(i))
802  naux = this%bfr%nauxarray(i)
803  call this%gwfpackages(ip)%set_auxname(naux, this%bfr%auxtxtarray(1:naux, i))
804  ip = ip + 1
805  end do
806  !
807  ! -- Copy just the package names for the boundary packages into
808  ! the flowpacknamearray
809  ip = 1
810  do i = 1, size(imap)
811  if (imap(i) == 1) then
812  this%flowpacknamearray(ip) = this%bfr%dstpackagenamearray(i)
813  ip = ip + 1
814  end if
815  end do
816  !
817  ! -- Error if specific discharge, saturation or flowja not found
818  if (.not. found_dataspdis) then
819  write (errmsg, '(4x,a)') 'SPECIFIC DISCHARGE NOT FOUND IN &
820  &BUDGET FILE. SAVE_SPECIFIC_DISCHARGE AND &
821  &SAVE_FLOWS MUST BE ACTIVATED IN THE NPF PACKAGE.'
822  call store_error(errmsg)
823  end if
824  if (.not. found_datasat) then
825  write (errmsg, '(4x,a)') 'SATURATION NOT FOUND IN &
826  &BUDGET FILE. SAVE_SATURATION AND &
827  &SAVE_FLOWS MUST BE ACTIVATED IN THE NPF PACKAGE.'
828  call store_error(errmsg)
829  end if
830  if (.not. found_flowja) then
831  write (errmsg, '(4x,a)') 'FLOWJA NOT FOUND IN &
832  &BUDGET FILE. SAVE_FLOWS MUST &
833  &BE ACTIVATED IN THE NPF PACKAGE.'
834  call store_error(errmsg)
835  end if
836  if (count_errors() > 0) then
837  call this%parser%StoreErrorUnit()
838  end if
839  !
840  ! -- return
841  return
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:

◆ initialize_gwfterms_from_gwfbndlist()

subroutine flowmodelinterfacemodule::initialize_gwfterms_from_gwfbndlist ( class(flowmodelinterfacetype this)

Definition at line 846 of file FlowModelInterface.f90.

847  ! -- modules
848  use bndmodule, only: bndtype, getbndfromlist
849  ! -- dummy
850  class(FlowModelInterfaceType) :: this
851  ! -- local
852  integer(I4B) :: ngwfpack
853  integer(I4B) :: ngwfterms
854  integer(I4B) :: ip
855  integer(I4B) :: imover
856  integer(I4B) :: ntomvr
857  integer(I4B) :: iterm
858  character(len=LENPACKAGENAME) :: budtxt
859  class(BndType), pointer :: packobj => null()
860  !
861  ! -- determine size of gwf terms
862  ngwfpack = this%gwfbndlist%Count()
863  !
864  ! -- Count number of to-mvr terms, but do not include advanced packages
865  ! as those mover terms are not losses from the cell, but rather flows
866  ! within the advanced package
867  ntomvr = 0
868  do ip = 1, ngwfpack
869  packobj => getbndfromlist(this%gwfbndlist, ip)
870  imover = packobj%imover
871  if (packobj%isadvpak /= 0) imover = 0
872  if (imover /= 0) then
873  ntomvr = ntomvr + 1
874  end if
875  end do
876  !
877  ! -- Allocate arrays in fmi of size ngwfterms, which is the number of
878  ! packages plus the number of packages with mover terms.
879  ngwfterms = ngwfpack + ntomvr
880  call this%allocate_gwfpackages(ngwfterms)
881  !
882  ! -- Assign values in the fmi package
883  iterm = 1
884  do ip = 1, ngwfpack
885  !
886  ! -- set and store names
887  packobj => getbndfromlist(this%gwfbndlist, ip)
888  budtxt = adjustl(packobj%text)
889  call this%gwfpackages(iterm)%set_name(packobj%packName, budtxt)
890  this%flowpacknamearray(iterm) = packobj%packName
891  iterm = iterm + 1
892  !
893  ! -- if this package has a mover associated with it, then add another
894  ! term that corresponds to the mover flows
895  imover = packobj%imover
896  if (packobj%isadvpak /= 0) imover = 0
897  if (imover /= 0) then
898  budtxt = trim(adjustl(packobj%text))//'-TO-MVR'
899  call this%gwfpackages(iterm)%set_name(packobj%packName, budtxt)
900  this%flowpacknamearray(iterm) = packobj%packName
901  this%igwfmvrterm(iterm) = 1
902  iterm = iterm + 1
903  end if
904  end do
905  return
Here is the call graph for this function:

◆ initialize_hfr()

subroutine flowmodelinterfacemodule::initialize_hfr ( class(flowmodelinterfacetype this)
private

Definition at line 623 of file FlowModelInterface.f90.

624  ! -- modules
625  class(FlowModelInterfaceType) :: this
626  ! -- dummy
627  !
628  ! -- Initialize the budget file reader
629  call this%hfr%initialize(this%iuhds, this%iout)
630  !
631  ! -- todo: need to run through the head terms
632  ! and do some checking

◆ read_options()

subroutine flowmodelinterfacemodule::read_options ( class(flowmodelinterfacetype this)

Definition at line 304 of file FlowModelInterface.f90.

305  ! -- modules
306  use constantsmodule, only: linelength, dem6
309  ! -- dummy
310  class(FlowModelInterfaceType) :: this
311  ! -- local
312  character(len=LINELENGTH) :: keyword
313  integer(I4B) :: ierr
314  logical :: isfound, endOfBlock
315  !
316  ! -- get options block
317  call this%parser%GetBlock('OPTIONS', isfound, ierr, blockrequired=.false., &
318  supportopenclose=.true.)
319  !
320  ! -- parse options block if detected
321  if (isfound) then
322  write (this%iout, '(1x,a)') 'PROCESSING FMI OPTIONS'
323  do
324  call this%parser%GetNextLine(endofblock)
325  if (endofblock) exit
326  call this%parser%GetStringCaps(keyword)
327  select case (keyword)
328  case ('SAVE_FLOWS')
329  this%ipakcb = -1
330  case default
331  write (errmsg, '(a,3(1x,a))') &
332  'UNKNOWN', trim(adjustl(this%text)), 'OPTION:', trim(keyword)
333  call store_error(errmsg)
334  call this%parser%StoreErrorUnit()
335  end select
336  end do
337  write (this%iout, '(1x,a)') 'END OF FMI OPTIONS'
338  end if
339  !
340  ! -- return
341  return
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dem6
real constant 1e-6
Definition: Constants.f90:108
subroutine, public urdaux(naux, inunit, iout, lloc, istart, istop, auxname, line, text)
Read auxiliary variables from an input line.
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
Here is the call graph for this function:

◆ read_packagedata()

subroutine flowmodelinterfacemodule::read_packagedata ( class(flowmodelinterfacetype this)

Definition at line 346 of file FlowModelInterface.f90.

347  ! -- modules
348  use openspecmodule, only: access, form
352  ! -- dummy
353  class(FlowModelInterfaceType) :: this
354  ! -- local
355  character(len=LINELENGTH) :: keyword, fname
356  integer(I4B) :: ierr
357  integer(I4B) :: inunit
358  integer(I4B) :: iapt
359  logical :: isfound, endOfBlock
360  logical :: blockrequired
361  logical :: exist
362  !
363  ! -- initialize
364  iapt = 0
365  blockrequired = .true.
366  !
367  ! -- get options block
368  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
369  blockrequired=blockrequired, &
370  supportopenclose=.true.)
371  !
372  ! -- parse options block if detected
373  if (isfound) then
374  write (this%iout, '(1x,a)') 'PROCESSING FMI PACKAGEDATA'
375  do
376  call this%parser%GetNextLine(endofblock)
377  if (endofblock) exit
378  call this%parser%GetStringCaps(keyword)
379  select case (keyword)
380  case ('GWFBUDGET')
381  call this%parser%GetStringCaps(keyword)
382  if (keyword /= 'FILEIN') then
383  call store_error('GWFBUDGET KEYWORD MUST BE FOLLOWED BY '// &
384  '"FILEIN" then by filename.')
385  call this%parser%StoreErrorUnit()
386  end if
387  call this%parser%GetString(fname)
388  inunit = getunit()
389  inquire (file=trim(fname), exist=exist)
390  if (.not. exist) then
391  call store_error('Could not find file '//trim(fname))
392  call this%parser%StoreErrorUnit()
393  end if
394  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
395  access, 'UNKNOWN')
396  this%iubud = inunit
397  call this%initialize_bfr()
398  case ('GWFHEAD')
399  call this%parser%GetStringCaps(keyword)
400  if (keyword /= 'FILEIN') then
401  call store_error('GWFHEAD KEYWORD MUST BE FOLLOWED BY '// &
402  '"FILEIN" then by filename.')
403  call this%parser%StoreErrorUnit()
404  end if
405  call this%parser%GetString(fname)
406  inquire (file=trim(fname), exist=exist)
407  if (.not. exist) then
408  call store_error('Could not find file '//trim(fname))
409  call this%parser%StoreErrorUnit()
410  end if
411  inunit = getunit()
412  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
413  access, 'UNKNOWN')
414  this%iuhds = inunit
415  call this%initialize_hfr()
416  case ('GWFMOVER')
417  call this%parser%GetStringCaps(keyword)
418  if (keyword /= 'FILEIN') then
419  call store_error('GWFMOVER KEYWORD MUST BE FOLLOWED BY '// &
420  '"FILEIN" then by filename.')
421  call this%parser%StoreErrorUnit()
422  end if
423  call this%parser%GetString(fname)
424  inunit = getunit()
425  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
426  access, 'UNKNOWN')
427  this%iumvr = inunit
428  call budgetobject_cr_bfr(this%mvrbudobj, 'MVT', this%iumvr, &
429  this%iout)
430  call this%mvrbudobj%fill_from_bfr(this%dis, this%iout)
431  case default
432  write (errmsg, '(a,3(1x,a))') &
433  'UNKNOWN', trim(adjustl(this%text)), 'PACKAGEDATA:', trim(keyword)
434  call store_error(errmsg)
435  end select
436  end do
437  write (this%iout, '(1x,a)') 'END OF FMI PACKAGEDATA'
438  end if
439  !
440  ! -- return
441  return
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:22
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Here is the call graph for this function: