MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
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 source_options (this)
 @ brief Source input options for package More...
 
subroutine source_packagedata (this)
 @ brief Source input options for package More...
 
subroutine read_grid (this)
 Read/validate flow model grid. 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 594 of file FlowModelInterface.f90.

595  ! -- modules
596  use tdismodule, only: kstp, kper
597  ! -- dummy
598  class(FlowModelInterfaceType) :: this
599  ! -- local
600  logical :: success
601  integer(I4B) :: n
602  integer(I4B) :: ipos
603  integer(I4B) :: nu, nr
604  integer(I4B) :: ip, i
605  logical :: readnext
606  ! -- format
607  character(len=*), parameter :: fmtkstpkper = &
608  "(1x,/1x,'FMI READING BUDGET TERMS &
609  &FOR KSTP ', i0, ' KPER ', i0)"
610  character(len=*), parameter :: fmtbudkstpkper = &
611  "(1x,/1x, 'FMI SETTING BUDGET TERMS &
612  &FOR KSTP ', i0, ' AND KPER ', &
613  &i0, ' TO BUDGET FILE TERMS FROM &
614  &KSTP ', i0, ' AND KPER ', i0)"
615  !
616  ! -- If the latest record read from the budget file is from a stress
617  ! -- period with only one time step, reuse that record (do not read a
618  ! -- new record) if the running model is still in that same stress period,
619  ! -- or if that record is the last one in the budget file.
620  readnext = .true.
621  if (kstp * kper > 1) then
622  if (this%bfr%header%kstp == 1) then
623  if (this%bfr%endoffile) then
624  readnext = .false.
625  else if (this%bfr%headernext%kper == kper + 1) then
626  readnext = .false.
627  end if
628  else if (this%bfr%endoffile) then
629  write (errmsg, '(4x,a)') 'REACHED END OF GWF BUDGET &
630  &FILE BEFORE READING SUFFICIENT BUDGET INFORMATION FOR THIS &
631  &GWT SIMULATION.'
632  call store_error(errmsg)
633  call store_error_unit(this%iubud)
634  end if
635  end if
636  !
637  ! -- Read the next record
638  if (readnext) then
639  !
640  ! -- Write the current time step and stress period
641  write (this%iout, fmtkstpkper) kstp, kper
642  !
643  ! -- loop through the budget terms for this stress period
644  ! i is the counter for gwf flow packages
645  ip = 1
646  do n = 1, this%bfr%nbudterms
647  call this%bfr%read_record(success, this%iout)
648  if (.not. success) then
649  write (errmsg, '(4x,a)') 'GWF BUDGET READ NOT SUCCESSFUL'
650  call store_error(errmsg)
651  call store_error_unit(this%iubud)
652  end if
653  !
654  ! -- Ensure kper is same between model and budget file
655  if (kper /= this%bfr%header%kper) then
656  write (errmsg, '(4x,a)') 'PERIOD NUMBER IN BUDGET FILE &
657  &DOES NOT MATCH PERIOD NUMBER IN TRANSPORT MODEL. IF THERE &
658  &IS MORE THAN ONE TIME STEP IN THE BUDGET FILE FOR A GIVEN &
659  &STRESS PERIOD, BUDGET FILE TIME STEPS MUST MATCH GWT MODEL &
660  &TIME STEPS ONE-FOR-ONE IN THAT STRESS PERIOD.'
661  call store_error(errmsg)
662  call store_error_unit(this%iubud)
663  end if
664  !
665  ! -- if budget file kstp > 1, then kstp must match
666  if (this%bfr%header%kstp > 1 .and. (kstp /= this%bfr%header%kstp)) then
667  write (errmsg, '(4x,a)') 'TIME STEP NUMBER IN BUDGET FILE &
668  &DOES NOT MATCH TIME STEP NUMBER IN TRANSPORT MODEL. IF THERE &
669  &IS MORE THAN ONE TIME STEP IN THE BUDGET FILE FOR A GIVEN STRESS &
670  &PERIOD, BUDGET FILE TIME STEPS MUST MATCH GWT MODEL TIME STEPS &
671  &ONE-FOR-ONE IN THAT STRESS PERIOD.'
672  call store_error(errmsg)
673  call store_error_unit(this%iubud)
674  end if
675  !
676  ! -- parse based on the type of data, and compress all user node
677  ! numbers into reduced node numbers
678  select type (h => this%bfr%header)
679  type is (budgetfileheadertype)
680  select case (trim(adjustl(h%budtxt)))
681  case ('FLOW-JA-FACE')
682  !
683  ! -- bfr%flowja contains only reduced connections so there is
684  ! a one-to-one match with this%gwfflowja
685  do ipos = 1, size(this%bfr%flowja)
686  this%gwfflowja(ipos) = this%bfr%flowja(ipos)
687  end do
688  case ('DATA-SPDIS')
689  do i = 1, h%nlist
690  nu = this%bfr%nodesrc(i)
691  nr = this%dis%get_nodenumber(nu, 0)
692  if (nr <= 0) cycle
693  this%gwfspdis(1, nr) = this%bfr%auxvar(1, i)
694  this%gwfspdis(2, nr) = this%bfr%auxvar(2, i)
695  this%gwfspdis(3, nr) = this%bfr%auxvar(3, i)
696  end do
697  case ('DATA-SAT')
698  do i = 1, h%nlist
699  nu = this%bfr%nodesrc(i)
700  nr = this%dis%get_nodenumber(nu, 0)
701  if (nr <= 0) cycle
702  this%gwfsat(nr) = this%bfr%auxvar(1, i)
703  end do
704  case ('STO-SS')
705  do nu = 1, this%dis%nodesuser
706  nr = this%dis%get_nodenumber(nu, 0)
707  if (nr <= 0) cycle
708  this%gwfstrgss(nr) = this%bfr%flow(nu)
709  end do
710  case ('STO-SY')
711  do nu = 1, this%dis%nodesuser
712  nr = this%dis%get_nodenumber(nu, 0)
713  if (nr <= 0) cycle
714  this%gwfstrgsy(nr) = this%bfr%flow(nu)
715  end do
716  case default
717  call this%gwfpackages(ip)%copy_values( &
718  h%nlist, &
719  this%bfr%nodesrc, &
720  this%bfr%flow, &
721  this%bfr%auxvar)
722  do i = 1, this%gwfpackages(ip)%nbound
723  nu = this%gwfpackages(ip)%nodelist(i)
724  nr = this%dis%get_nodenumber(nu, 0)
725  this%gwfpackages(ip)%nodelist(i) = nr
726  end do
727  ip = ip + 1
728  end select
729  end select
730  end do
731  else
732  !
733  ! -- write message to indicate that flows are being reused
734  write (this%iout, fmtbudkstpkper) kstp, kper, &
735  this%bfr%header%kstp, this%bfr%header%kper
736  !
737  ! -- set the flag to indicate that flows were not updated
738  this%iflowsupdated = 0
739  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 757 of file FlowModelInterface.f90.

758  ! modules
759  use tdismodule, only: kstp, kper
760  class(FlowModelInterfaceType) :: this
761  integer(I4B) :: nu, nr, i, ilay
762  integer(I4B) :: ncpl
763  real(DP) :: val
764  logical :: readnext
765  logical :: success
766  character(len=*), parameter :: fmtkstpkper = &
767  "(1x,/1x,'FMI READING HEAD FOR &
768  &KSTP ', i0, ' KPER ', i0)"
769  character(len=*), parameter :: fmthdskstpkper = &
770  "(1x,/1x, 'FMI SETTING HEAD FOR KSTP ', i0, ' AND KPER ', &
771  &i0, ' TO BINARY FILE HEADS FROM KSTP ', i0, ' AND KPER ', i0)"
772  !
773  ! -- If the latest record read from the head file is from a stress
774  ! -- period with only one time step, reuse that record (do not read a
775  ! -- new record) if the running model is still in that same stress period,
776  ! -- or if that record is the last one in the head file.
777  readnext = .true.
778  if (kstp * kper > 1) then
779  if (this%hfr%header%kstp == 1) then
780  if (this%hfr%endoffile) then
781  readnext = .false.
782  else if (this%hfr%headernext%kper == kper + 1) then
783  readnext = .false.
784  end if
785  else if (this%hfr%endoffile) then
786  write (errmsg, '(4x,a)') 'REACHED END OF GWF HEAD &
787  &FILE BEFORE READING SUFFICIENT HEAD INFORMATION FOR THIS &
788  &GWT SIMULATION.'
789  call store_error(errmsg)
790  call store_error_unit(this%iuhds)
791  end if
792  end if
793  !
794  ! -- Read the next record
795  if (readnext) then
796  !
797  ! -- write to list file that heads are being read
798  write (this%iout, fmtkstpkper) kstp, kper
799  !
800  ! -- loop through the layered heads for this time step
801  do ilay = 1, this%hfr%nlay
802  !
803  ! -- read next head chunk
804  call this%hfr%read_record(success, this%iout)
805  if (.not. success) then
806  write (errmsg, '(4x,a)') 'GWF HEAD READ NOT SUCCESSFUL'
807  call store_error(errmsg)
808  call store_error_unit(this%iuhds)
809  end if
810  !
811  ! -- Ensure kper is same between model and head file
812  if (kper /= this%hfr%header%kper) then
813  write (errmsg, '(4x,a)') 'PERIOD NUMBER IN HEAD FILE &
814  &DOES NOT MATCH PERIOD NUMBER IN TRANSPORT MODEL. IF THERE &
815  &IS MORE THAN ONE TIME STEP IN THE HEAD FILE FOR A GIVEN STRESS &
816  &PERIOD, HEAD FILE TIME STEPS MUST MATCH GWT MODEL TIME STEPS &
817  &ONE-FOR-ONE IN THAT STRESS PERIOD.'
818  call store_error(errmsg)
819  call store_error_unit(this%iuhds)
820  end if
821  !
822  ! -- if head file kstp > 1, then kstp must match
823  if (this%hfr%header%kstp > 1 .and. (kstp /= this%hfr%header%kstp)) then
824  write (errmsg, '(4x,a)') 'TIME STEP NUMBER IN HEAD FILE &
825  &DOES NOT MATCH TIME STEP NUMBER IN TRANSPORT MODEL. IF THERE &
826  &IS MORE THAN ONE TIME STEP IN THE HEAD FILE FOR A GIVEN STRESS &
827  &PERIOD, HEAD FILE TIME STEPS MUST MATCH GWT MODEL TIME STEPS &
828  &ONE-FOR-ONE IN THAT STRESS PERIOD.'
829  call store_error(errmsg)
830  call store_error_unit(this%iuhds)
831  end if
832  !
833  ! -- fill the head array for this layer and
834  ! compress into reduced form
835  ncpl = size(this%hfr%head)
836  do i = 1, ncpl
837  nu = (ilay - 1) * ncpl + i
838  nr = this%dis%get_nodenumber(nu, 0)
839  val = this%hfr%head(i)
840  if (nr > 0) this%gwfhead(nr) = val
841  end do
842  end do
843  else
844  write (this%iout, fmthdskstpkper) kstp, kper, &
845  this%hfr%header%kstp, this%hfr%header%kper
846  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 247 of file FlowModelInterface.f90.

249  !modules
250  use constantsmodule, only: dzero
251  ! -- dummy
252  class(FlowModelInterfaceType) :: this
253  integer(I4B), intent(in) :: nodes
254  ! -- local
255  integer(I4B) :: n
256  !
257  ! -- Allocate ibdgwfsat0, which is an indicator array marking cells with
258  ! saturation greater than 0.0 with a value of 1
259  call mem_allocate(this%ibdgwfsat0, nodes, 'IBDGWFSAT0', this%memoryPath)
260  do n = 1, nodes
261  this%ibdgwfsat0(n) = 1
262  end do
263  !
264  ! -- Allocate differently depending on whether or not flows are
265  ! being read from a file.
266  if (this%flows_from_file) then
267  call mem_allocate(this%gwfflowja, this%dis%con%nja, &
268  'GWFFLOWJA', this%memoryPath)
269  call mem_allocate(this%gwfsat, nodes, 'GWFSAT', this%memoryPath)
270  call mem_allocate(this%gwfhead, nodes, 'GWFHEAD', this%memoryPath)
271  call mem_allocate(this%gwfspdis, 3, nodes, 'GWFSPDIS', this%memoryPath)
272  do n = 1, nodes
273  this%gwfsat(n) = done
274  this%gwfhead(n) = dzero
275  this%gwfspdis(:, n) = dzero
276  end do
277  do n = 1, size(this%gwfflowja)
278  this%gwfflowja(n) = dzero
279  end do
280  !
281  ! -- allocate and initialize storage arrays
282  if (this%igwfstrgss == 0) then
283  call mem_allocate(this%gwfstrgss, 1, 'GWFSTRGSS', this%memoryPath)
284  else
285  call mem_allocate(this%gwfstrgss, nodes, 'GWFSTRGSS', this%memoryPath)
286  end if
287  if (this%igwfstrgsy == 0) then
288  call mem_allocate(this%gwfstrgsy, 1, 'GWFSTRGSY', this%memoryPath)
289  else
290  call mem_allocate(this%gwfstrgsy, nodes, 'GWFSTRGSY', this%memoryPath)
291  end if
292  do n = 1, size(this%gwfstrgss)
293  this%gwfstrgss(n) = dzero
294  end do
295  do n = 1, size(this%gwfstrgsy)
296  this%gwfstrgsy(n) = dzero
297  end do
298  ! allocate and initialize cell type array. if the FMI is in a separate
299  ! simulation from the GWF model, we expect cell type to have been read
300  ! already if the binary grid file was provided to FMI. otherwise don't
301  ! initialize the cell type array to any default; unless it is received
302  ! from GWF NPF by an EXG it's undefined as indicated by igwfceltyp = 0
303  ! (this is because some coupled models need cell type, but some don't)
304  if (this%igwfceltyp == 0) &
305  call mem_allocate(this%gwfceltyp, nodes, 'GWFCELTYP', this%memoryPath)
306  !
307  ! -- If there is no fmi package, then there are no flows at all or a
308  ! connected GWF model, so allocate gwfpackages to zero
309  if (this%inunit == 0) call this%allocate_gwfpackages(this%nflowpack)
310  end if
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ 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 1021 of file FlowModelInterface.f90.

1022  ! -- modules
1023  use constantsmodule, only: lenmempath
1025  ! -- dummy
1026  class(FlowModelInterfaceType) :: this
1027  integer(I4B), intent(in) :: ngwfterms
1028  ! -- local
1029  integer(I4B) :: n
1030  character(len=LENMEMPATH) :: memPath
1031  !
1032  ! -- direct allocate
1033  allocate (this%gwfpackages(ngwfterms))
1034  allocate (this%flowpacknamearray(ngwfterms))
1035  !
1036  ! -- mem_allocate
1037  call mem_allocate(this%igwfmvrterm, ngwfterms, 'IGWFMVRTERM', this%memoryPath)
1038  !
1039  ! -- initialize
1040  this%nflowpack = ngwfterms
1041  do n = 1, this%nflowpack
1042  this%igwfmvrterm(n) = 0
1043  this%flowpacknamearray(n) = ''
1044  !
1045  ! -- Create a mempath for each individual flow package data set
1046  ! of the form, MODELNAME/FMI-FTn
1047  write (mempath, '(a, i0)') trim(this%memoryPath)//'-FT', n
1048  call this%gwfpackages(n)%initialize(mempath)
1049  end do
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27

◆ allocate_scalars()

subroutine flowmodelinterfacemodule::allocate_scalars ( class(flowmodelinterfacetype this)

Definition at line 204 of file FlowModelInterface.f90.

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

◆ deallocate_gwfpackages()

subroutine flowmodelinterfacemodule::deallocate_gwfpackages ( class(flowmodelinterfacetype this)

Definition at line 1053 of file FlowModelInterface.f90.

1054  class(FlowModelInterfaceType) :: this
1055  integer(I4B) :: n
1056 
1057  do n = 1, this%nflowpack
1058  call this%gwfpackages(n)%da()
1059  end do

◆ finalize_bfr()

subroutine flowmodelinterfacemodule::finalize_bfr ( class(flowmodelinterfacetype this)

Definition at line 743 of file FlowModelInterface.f90.

744  class(FlowModelInterfaceType) :: this
745  call this%bfr%finalize()

◆ finalize_hfr()

subroutine flowmodelinterfacemodule::finalize_hfr ( class(flowmodelinterfacetype this)

Definition at line 850 of file FlowModelInterface.f90.

851  class(FlowModelInterfaceType) :: this
852  close (this%iuhds)

◆ fmi_ar()

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

Definition at line 141 of file FlowModelInterface.f90.

142  ! -- modules
143  ! -- dummy
144  class(FlowModelInterfaceType) :: this
145  integer(I4B), dimension(:), pointer, contiguous :: ibound
146  !
147  ! -- store pointers to arguments that were passed in
148  this%ibound => ibound
149  !
150  ! -- Allocate arrays
151  call this%allocate_arrays(this%dis%nodes)

◆ fmi_da()

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

Definition at line 156 of file FlowModelInterface.f90.

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

◆ fmi_df()

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

Definition at line 85 of file FlowModelInterface.f90.

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

◆ 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 1063 of file FlowModelInterface.f90.

1064  use bndmodule, only: bndtype, getbndfromlist
1065  class(FlowModelInterfaceType) :: this
1066  character(len=*), intent(in) :: name
1067  integer(I4B), intent(inout) :: idx
1068  ! -- local
1069  integer(I4B) :: ip
1070  !
1071  ! -- Look through all the packages and return the index with name
1072  idx = 0
1073  do ip = 1, size(this%flowpacknamearray)
1074  if (this%flowpacknamearray(ip) == name) then
1075  idx = ip
1076  exit
1077  end if
1078  end do
1079  if (idx == 0) then
1080  call store_error('Error in get_package_index. Could not find '//name, &
1081  terminate=.true.)
1082  end if
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 581 of file FlowModelInterface.f90.

582  class(FlowModelInterfaceType) :: this
583  integer(I4B) :: ncrbud
584  call this%bfr%initialize(this%iubud, this%iout, ncrbud)
585  ! todo: need to run through the budget terms
586  ! 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 860 of file FlowModelInterface.f90.

861  ! -- dummy
862  class(FlowModelInterfaceType) :: this
863  ! -- local
864  integer(I4B) :: nflowpack
865  integer(I4B) :: i, ip
866  integer(I4B) :: naux
867  logical :: found_flowja
868  logical :: found_dataspdis
869  logical :: found_datasat
870  logical :: found_stoss
871  logical :: found_stosy
872  integer(I4B), dimension(:), allocatable :: imap
873  !
874  ! -- Calculate the number of gwf flow packages
875  allocate (imap(this%bfr%nbudterms))
876  imap(:) = 0
877  nflowpack = 0
878  found_flowja = .false.
879  found_dataspdis = .false.
880  found_datasat = .false.
881  found_stoss = .false.
882  found_stosy = .false.
883  do i = 1, this%bfr%nbudterms
884  select case (trim(adjustl(this%bfr%budtxtarray(i))))
885  case ('FLOW-JA-FACE')
886  found_flowja = .true.
887  case ('DATA-SPDIS')
888  found_dataspdis = .true.
889  this%igwfspdis = 1
890  case ('DATA-SAT')
891  found_datasat = .true.
892  case ('STO-SS')
893  found_stoss = .true.
894  this%igwfstrgss = 1
895  case ('STO-SY')
896  found_stosy = .true.
897  this%igwfstrgsy = 1
898  case default
899  nflowpack = nflowpack + 1
900  imap(i) = 1
901  end select
902  end do
903  !
904  ! -- allocate gwfpackage arrays
905  call this%allocate_gwfpackages(nflowpack)
906  !
907  ! -- Copy the package name and aux names from budget file reader
908  ! to the gwfpackages derived-type variable
909  ip = 1
910  do i = 1, this%bfr%nbudterms
911  if (imap(i) == 0) cycle
912  call this%gwfpackages(ip)%set_name(this%bfr%dstpackagenamearray(i), &
913  this%bfr%budtxtarray(i))
914  naux = this%bfr%nauxarray(i)
915  call this%gwfpackages(ip)%set_auxname(naux, this%bfr%auxtxtarray(1:naux, i))
916  ip = ip + 1
917  end do
918  !
919  ! -- Copy just the package names for the boundary packages into
920  ! the flowpacknamearray
921  ip = 1
922  do i = 1, size(imap)
923  if (imap(i) == 1) then
924  this%flowpacknamearray(ip) = this%bfr%dstpackagenamearray(i)
925  ip = ip + 1
926  end if
927  end do
928  !
929  ! -- Error if specific discharge, saturation or flowja not found
930  if (.not. found_dataspdis) then
931  write (errmsg, '(4x,a)') 'SPECIFIC DISCHARGE NOT FOUND IN &
932  &BUDGET FILE. SAVE_SPECIFIC_DISCHARGE AND &
933  &SAVE_FLOWS MUST BE ACTIVATED IN THE NPF PACKAGE.'
934  call store_error(errmsg)
935  end if
936  if (.not. found_datasat) then
937  write (errmsg, '(4x,a)') 'SATURATION NOT FOUND IN &
938  &BUDGET FILE. SAVE_SATURATION AND &
939  &SAVE_FLOWS MUST BE ACTIVATED IN THE NPF PACKAGE.'
940  call store_error(errmsg)
941  end if
942  if (.not. found_flowja) then
943  write (errmsg, '(4x,a)') 'FLOWJA NOT FOUND IN &
944  &BUDGET FILE. SAVE_FLOWS MUST &
945  &BE ACTIVATED IN THE NPF PACKAGE.'
946  call store_error(errmsg)
947  end if
948  if (count_errors() > 0) then
949  call store_error_filename(this%input_fname)
950  end if
Here is the call graph for this function:

◆ initialize_gwfterms_from_gwfbndlist()

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

Definition at line 954 of file FlowModelInterface.f90.

955  ! -- modules
956  use bndmodule, only: bndtype, getbndfromlist
957  ! -- dummy
958  class(FlowModelInterfaceType) :: this
959  ! -- local
960  integer(I4B) :: ngwfpack
961  integer(I4B) :: ngwfterms
962  integer(I4B) :: ip
963  integer(I4B) :: imover
964  integer(I4B) :: ntomvr
965  integer(I4B) :: iterm
966  character(len=LENPACKAGENAME) :: budtxt
967  class(BndType), pointer :: packobj => null()
968  !
969  ! -- determine size of gwf terms
970  ngwfpack = this%gwfbndlist%Count()
971  !
972  ! -- Count number of to-mvr terms, but do not include advanced packages
973  ! as those mover terms are not losses from the cell, but rather flows
974  ! within the advanced package
975  ntomvr = 0
976  do ip = 1, ngwfpack
977  packobj => getbndfromlist(this%gwfbndlist, ip)
978  imover = packobj%imover
979  if (packobj%isadvpak /= 0) imover = 0
980  if (imover /= 0) then
981  ntomvr = ntomvr + 1
982  end if
983  end do
984  !
985  ! -- Allocate arrays in fmi of size ngwfterms, which is the number of
986  ! packages plus the number of packages with mover terms.
987  ngwfterms = ngwfpack + ntomvr
988  call this%allocate_gwfpackages(ngwfterms)
989  !
990  ! -- Assign values in the fmi package
991  iterm = 1
992  do ip = 1, ngwfpack
993  !
994  ! -- set and store names
995  packobj => getbndfromlist(this%gwfbndlist, ip)
996  budtxt = adjustl(packobj%text)
997  call this%gwfpackages(iterm)%set_name(packobj%packName, budtxt)
998  this%flowpacknamearray(iterm) = packobj%packName
999  iterm = iterm + 1
1000  !
1001  ! -- if this package has a mover associated with it, then add another
1002  ! term that corresponds to the mover flows
1003  imover = packobj%imover
1004  if (packobj%isadvpak /= 0) imover = 0
1005  if (imover /= 0) then
1006  budtxt = trim(adjustl(packobj%text))//'-TO-MVR'
1007  call this%gwfpackages(iterm)%set_name(packobj%packName, budtxt)
1008  this%flowpacknamearray(iterm) = packobj%packName
1009  this%igwfmvrterm(iterm) = 1
1010  iterm = iterm + 1
1011  end if
1012  end do
Here is the call graph for this function:

◆ initialize_hfr()

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

Definition at line 749 of file FlowModelInterface.f90.

750  class(FlowModelInterfaceType) :: this
751  call this%hfr%initialize(this%iuhds, this%iout)
752  ! todo: need to run through the head terms
753  ! and do some checking

◆ read_grid()

subroutine flowmodelinterfacemodule::read_grid ( class(flowmodelinterfacetype this)

Definition at line 429 of file FlowModelInterface.f90.

430  ! -- modules
431  use dismodule, only: distype
432  use disvmodule, only: disvtype
433  use disumodule, only: disutype
434  use dis2dmodule, only: dis2dtype
435  use disv2dmodule, only: disv2dtype
436  use disv1dmodule, only: disv1dtype
437  ! -- dummy
438  class(FlowModelInterfaceType) :: this
439  ! -- local
440  integer(I4B) :: user_nodes
441  integer(I4B), allocatable :: idomain1d(:), idomain2d(:, :), idomain3d(:, :, :)
442  ! -- formats
443  character(len=*), parameter :: fmticterr = &
444  &"('Error in ',a,': Binary grid file does not contain ICELLTYPE.')"
445  character(len=*), parameter :: fmtdiserr = &
446  "('Error in ',a,': Models do not have the same discretization. &
447  &GWF model has ', i0, ' user nodes, this model has ', i0, '. &
448  &Ensure discretization packages, including IDOMAIN, are identical.')"
449  character(len=*), parameter :: fmtidomerr = &
450  "('Error in ',a,': models do not have the same discretization. &
451  &Models have different IDOMAIN arrays. &
452  &Ensure discretization packages, including IDOMAIN, are identical.')"
453 
454  call this%gfr%initialize(this%iugrb)
455 
456  ! load icelltype array
457  if (.not. this%gfr%has_variable("ICELLTYPE")) then
458  write (errmsg, fmticterr) trim(this%text)
459  call store_error(errmsg, terminate=.true.)
460  end if
461  this%igwfceltyp = 1
462  call mem_allocate(this%gwfceltyp, this%dis%nodesuser, &
463  'GWFCELTYP', this%memoryPath)
464  call this%gfr%read_int_1d_into("ICELLTYPE", this%gwfceltyp)
465 
466  ! check grid equivalence
467  select case (this%gfr%grid_type)
468  case ('DIS')
469  select type (dis => this%dis)
470  type is (distype)
471  user_nodes = this%gfr%read_int("NCELLS")
472  if (user_nodes /= this%dis%nodesuser) then
473  write (errmsg, fmtdiserr) &
474  trim(this%text), user_nodes, this%dis%nodesuser
475  call store_error(errmsg, terminate=.true.)
476  end if
477  idomain1d = this%gfr%read_int_1d("IDOMAIN")
478  idomain3d = reshape(idomain1d, [ &
479  this%gfr%read_int("NCOL"), &
480  this%gfr%read_int("NROW"), &
481  this%gfr%read_int("NLAY") &
482  ])
483  if (.not. all(dis%idomain == idomain3d)) then
484  write (errmsg, fmtidomerr) trim(this%text)
485  call store_error(errmsg, terminate=.true.)
486  end if
487  end select
488  case ('DISV')
489  select type (dis => this%dis)
490  type is (disvtype)
491  user_nodes = this%gfr%read_int("NCELLS")
492  if (user_nodes /= this%dis%nodesuser) then
493  write (errmsg, fmtdiserr) &
494  trim(this%text), user_nodes, this%dis%nodesuser
495  call store_error(errmsg, terminate=.true.)
496  end if
497  idomain1d = this%gfr%read_int_1d("IDOMAIN")
498  idomain2d = reshape(idomain1d, [ &
499  this%gfr%read_int("NCPL"), &
500  this%gfr%read_int("NLAY") &
501  ])
502  if (.not. all(dis%idomain == idomain2d)) then
503  write (errmsg, fmtidomerr) trim(this%text)
504  call store_error(errmsg, terminate=.true.)
505  end if
506  end select
507  case ('DISU')
508  select type (dis => this%dis)
509  type is (disutype)
510  user_nodes = this%gfr%read_int("NODES")
511  if (user_nodes /= this%dis%nodesuser) then
512  write (errmsg, fmtdiserr) &
513  trim(this%text), user_nodes, this%dis%nodesuser
514  call store_error(errmsg, terminate=.true.)
515  end if
516  idomain1d = this%gfr%read_int_1d("IDOMAIN")
517  if (.not. all(dis%idomain == idomain1d)) then
518  write (errmsg, fmtidomerr) trim(this%text)
519  call store_error(errmsg, terminate=.true.)
520  end if
521  end select
522  case ('DIS2D')
523  select type (dis => this%dis)
524  type is (dis2dtype)
525  user_nodes = this%gfr%read_int("NCELLS")
526  if (user_nodes /= this%dis%nodesuser) then
527  write (errmsg, fmtdiserr) &
528  trim(this%text), user_nodes, this%dis%nodesuser
529  call store_error(errmsg, terminate=.true.)
530  end if
531  idomain1d = this%gfr%read_int_1d("IDOMAIN")
532  idomain2d = reshape(idomain1d, [ &
533  this%gfr%read_int("NCOL"), &
534  this%gfr%read_int("NROW") &
535  ])
536  if (.not. all(dis%idomain == idomain2d)) then
537  write (errmsg, fmtidomerr) trim(this%text)
538  call store_error(errmsg, terminate=.true.)
539  end if
540  end select
541  case ('DISV2D')
542  select type (dis => this%dis)
543  type is (disv2dtype)
544  user_nodes = this%gfr%read_int("NODES")
545  if (user_nodes /= this%dis%nodesuser) then
546  write (errmsg, fmtdiserr) &
547  trim(this%text), user_nodes, this%dis%nodesuser
548  call store_error(errmsg, terminate=.true.)
549  end if
550  idomain1d = this%gfr%read_int_1d("IDOMAIN")
551  if (.not. all(dis%idomain == idomain1d)) then
552  write (errmsg, fmtidomerr) trim(this%text)
553  call store_error(errmsg, terminate=.true.)
554  end if
555  end select
556  case ('DISV1D')
557  select type (dis => this%dis)
558  type is (disv1dtype)
559  user_nodes = this%gfr%read_int("NCELLS")
560  if (user_nodes /= this%dis%nodesuser) then
561  write (errmsg, fmtdiserr) &
562  trim(this%text), user_nodes, this%dis%nodesuser
563  call store_error(errmsg, terminate=.true.)
564  end if
565  idomain1d = this%gfr%read_int_1d("IDOMAIN")
566  if (.not. all(dis%idomain == idomain1d)) then
567  write (errmsg, fmtidomerr) trim(this%text)
568  call store_error(errmsg, terminate=.true.)
569  end if
570  end select
571  end select
572 
573  if (allocated(idomain3d)) deallocate (idomain3d)
574  if (allocated(idomain2d)) deallocate (idomain2d)
575  if (allocated(idomain1d)) deallocate (idomain1d)
576 
577  call this%gfr%finalize()
Definition: Dis.f90:1
Structured grid discretization.
Definition: Dis2d.f90:23
Structured grid discretization.
Definition: Dis.f90:23
Unstructured grid discretization.
Definition: Disu.f90:28
Vertex grid discretization.
Definition: Disv2d.f90:24
Vertex grid discretization.
Definition: Disv.f90:24
Here is the call graph for this function:

◆ source_options()

subroutine flowmodelinterfacemodule::source_options ( class(flowmodelinterfacetype this)

Definition at line 315 of file FlowModelInterface.f90.

316  ! -- modules
318  ! -- dummy
319  class(FlowModelInterfaceType) :: this
320  ! -- local
321  logical(LGP) :: found_ipakcb
322  character(len=*), parameter :: fmtisvflow = &
323  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
324  &WHENEVER ICBCFL IS NOT ZERO AND FLOW IMBALANCE CORRECTION ACTIVE.')"
325 
326  ! -- source package input
327  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', this%input_mempath, &
328  found_ipakcb)
329 
330  write (this%iout, '(1x,a)') 'PROCESSING FMI OPTIONS'
331 
332  if (found_ipakcb) then
333  this%ipakcb = -1
334  write (this%iout, fmtisvflow)
335  end if
336 
337  write (this%iout, '(1x,a)') 'END OF FMI OPTIONS'

◆ source_packagedata()

subroutine flowmodelinterfacemodule::source_packagedata ( class(flowmodelinterfacetype this)

Definition at line 342 of file FlowModelInterface.f90.

343  ! -- modules
347  use openspecmodule, only: access, form
350  ! -- dummy
351  class(FlowModelInterfaceType) :: this
352  ! -- local
353  type(CharacterStringType), dimension(:), contiguous, &
354  pointer :: flowtypes
355  type(CharacterStringType), dimension(:), contiguous, &
356  pointer :: fileops
357  type(CharacterStringType), dimension(:), contiguous, &
358  pointer :: fnames
359  character(len=LINELENGTH) :: flowtype, fileop, fname
360  integer(I4B) :: inunit, n
361  logical(LGP) :: exist
362 
363  call mem_setptr(flowtypes, 'FLOWTYPE', this%input_mempath)
364  call mem_setptr(fileops, 'FILEIN', this%input_mempath)
365  call mem_setptr(fnames, 'FNAME', this%input_mempath)
366 
367  write (this%iout, '(1x,a)') 'PROCESSING FMI PACKAGEDATA'
368 
369  do n = 1, size(flowtypes)
370  flowtype = flowtypes(n)
371  fileop = fileops(n)
372  fname = fnames(n)
373 
374  inquire (file=trim(fname), exist=exist)
375  if (.not. exist) then
376  call store_error('Could not find file '//trim(fname))
377  cycle
378  end if
379 
380  if (fileop /= 'FILEIN') then
381  call store_error('Unexpected packagedata input keyword read: "' &
382  //trim(fileop)//'".')
383  cycle
384  end if
385 
386  select case (flowtype)
387  case ('GWFBUDGET')
388  inunit = getunit()
389  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
390  access, 'UNKNOWN')
391  this%iubud = inunit
392  call this%initialize_bfr()
393  case ('GWFHEAD')
394  inunit = getunit()
395  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
396  access, 'UNKNOWN')
397  this%iuhds = inunit
398  call this%initialize_hfr()
399  case ('GWFMOVER')
400  inunit = getunit()
401  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
402  access, 'UNKNOWN')
403  this%iumvr = inunit
404  call budgetobject_cr_bfr(this%mvrbudobj, 'MVT', this%iumvr, &
405  this%iout)
406  call this%mvrbudobj%fill_from_bfr(this%dis, this%iout)
407  case ('GWFGRID')
408  inunit = getunit()
409  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', &
410  form, access, 'UNKNOWN')
411  this%iugrb = inunit
412  call this%read_grid()
413  case default
414  write (errmsg, '(a,3(1x,a))') &
415  'UNKNOWN', trim(adjustl(this%text)), 'PACKAGEDATA:', trim(flowtype)
416  call store_error(errmsg)
417  end select
418  end do
419 
420  write (this%iout, '(1x,a)') 'END OF FMI PACKAGEDATA'
421 
422  if (count_errors() > 0) then
423  call store_error_filename(this%input_fname)
424  end if
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
real(dp), parameter dem6
real constant 1e-6
Definition: Constants.f90:109
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
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function: