MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
prt-fmi.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
6  use simvariablesmodule, only: errmsg
8  use basedismodule, only: disbasetype
10  use cellmodule, only: max_poly_cells
11 
12  implicit none
13  private
14  public :: prtfmitype
15  public :: fmi_cr
16 
17  character(len=LENPACKAGENAME) :: text = ' PRTFMI'
18 
20 
21  double precision, allocatable, public :: sourceflows(:) ! cell source flows array
22  double precision, allocatable, public :: sinkflows(:) ! cell sink flows array
23  double precision, allocatable, public :: storageflows(:) ! cell storage flows array
24  double precision, allocatable, public :: boundaryflows(:) ! cell boundary flows array
25 
26  contains
27 
28  procedure :: fmi_ad
29  procedure :: fmi_df => prtfmi_df
30  procedure, private :: accumulate_flows
31 
32  end type prtfmitype
33 
34 contains
35 
36  !> @brief Create a new PrtFmi object
37  subroutine fmi_cr(fmiobj, name_model, inunit, iout)
38  ! -- dummy
39  type(prtfmitype), pointer :: fmiobj
40  character(len=*), intent(in) :: name_model
41  integer(I4B), intent(inout) :: inunit
42  integer(I4B), intent(in) :: iout
43  !
44  ! -- Create the object
45  allocate (fmiobj)
46  !
47  ! -- create name and memory path
48  call fmiobj%set_names(1, name_model, 'FMI', 'FMI')
49  fmiobj%text = text
50  !
51  ! -- Allocate scalars
52  call fmiobj%allocate_scalars()
53  !
54  ! -- Set variables
55  fmiobj%inunit = inunit
56  fmiobj%iout = iout
57  !
58  ! -- Initialize block parser
59  call fmiobj%parser%Initialize(fmiobj%inunit, fmiobj%iout)
60  !
61  ! -- Assign dependent variable label
62  fmiobj%depvartype = 'TRACKS '
63 
64  end subroutine fmi_cr
65 
66  !> @brief Time step advance
67  subroutine fmi_ad(this)
68  ! -- modules
69  use constantsmodule, only: dhdry
70  ! -- dummy
71  class(prtfmitype) :: this
72  ! -- local
73  integer(I4B) :: n
74  character(len=15) :: nodestr
75  character(len=*), parameter :: fmtdry = &
76  &"(/1X,'WARNING: DRY CELL ENCOUNTERED AT ',a,'; RESET AS INACTIVE')"
77  character(len=*), parameter :: fmtrewet = &
78  &"(/1X,'DRY CELL REACTIVATED AT ', a)"
79  !
80  ! -- Set flag to indicated that flows are being updated. For the case where
81  ! flows may be reused (only when flows are read from a file) then set
82  ! the flag to zero to indicated that flows were not updated
83  this%iflowsupdated = 1
84  !
85  ! -- If reading flows from a budget file, read the next set of records
86  if (this%iubud /= 0) then
87  call this%advance_bfr()
88  end if
89  !
90  ! -- If reading heads from a head file, read the next set of records
91  if (this%iuhds /= 0) then
92  call this%advance_hfr()
93  end if
94  !
95  ! -- If mover flows are being read from file, read the next set of records
96  if (this%iumvr /= 0) then
97  call this%mvrbudobj%bfr_advance(this%dis, this%iout)
98  end if
99  !
100  ! -- Accumulate flows
101  call this%accumulate_flows()
102  !
103  ! -- if flow cell is dry, then set this%ibound = 0
104  do n = 1, this%dis%nodes
105  !
106  ! -- Calculate the ibound-like array that has 0 if saturation
107  ! is zero and 1 otherwise
108  if (this%gwfsat(n) > dzero) then
109  this%ibdgwfsat0(n) = 1
110  else
111  this%ibdgwfsat0(n) = 0
112  end if
113  !
114  ! -- Check if active model cell is inactive for flow
115  if (this%ibound(n) > 0) then
116  if (this%gwfhead(n) == dhdry) then
117  ! -- cell should be made inactive
118  this%ibound(n) = 0
119  call this%dis%noder_to_string(n, nodestr)
120  write (this%iout, fmtdry) trim(nodestr)
121  end if
122  end if
123  !
124  ! -- Convert dry model cell to active if flow has rewet
125  if (this%ibound(n) == 0) then
126  if (this%gwfhead(n) /= dhdry) then
127  ! -- cell is now wet
128  this%ibound(n) = 1
129  call this%dis%noder_to_string(n, nodestr)
130  write (this%iout, fmtrewet) trim(nodestr)
131  end if
132  end if
133  end do
134 
135  end subroutine fmi_ad
136 
137  !> @brief Define the flow model interface
138  subroutine prtfmi_df(this, dis, idryinactive)
139  ! -- modules
140  use simmodule, only: store_error
141  ! -- dummy
142  class(prtfmitype) :: this
143  class(disbasetype), pointer, intent(in) :: dis
144  integer(I4B), intent(in) :: idryinactive
145  !
146  ! -- Call parent class define
147  call this%FlowModelInterfaceType%fmi_df(dis, idryinactive)
148  !
149  ! -- Allocate arrays
150  allocate (this%StorageFlows(this%dis%nodes))
151  allocate (this%SourceFlows(this%dis%nodes))
152  allocate (this%SinkFlows(this%dis%nodes))
153  allocate (this%BoundaryFlows(this%dis%nodes * max_poly_cells))
154 
155  end subroutine prtfmi_df
156 
157  !> @brief Accumulate flows
158  subroutine accumulate_flows(this)
159  implicit none
160  ! -- dummy
161  class(prtfmitype) :: this
162  ! -- local
163  integer(I4B) :: j, i, ip, ib
164  integer(I4B) :: ioffset, iflowface, iauxiflowface
165  real(DP) :: qbnd
166  character(len=LENAUXNAME) :: auxname
167  integer(I4B) :: naux
168  !
169  this%StorageFlows = dzero
170  if (this%igwfstrgss /= 0) &
171  this%StorageFlows = this%StorageFlows + &
172  this%gwfstrgss
173  if (this%igwfstrgsy /= 0) &
174  this%StorageFlows = this%StorageFlows + &
175  this%gwfstrgsy
176 
177  this%SourceFlows = dzero
178  this%SinkFlows = dzero
179  this%BoundaryFlows = dzero
180  do ip = 1, this%nflowpack
181  iauxiflowface = 0
182  naux = this%gwfpackages(ip)%naux
183  if (naux > 0) then
184  do j = 1, naux
185  auxname = this%gwfpackages(ip)%auxname(j)
186  if (trim(adjustl(auxname)) == "IFLOWFACE") then
187  iauxiflowface = j
188  exit
189  end if
190  end do
191  end if
192  do ib = 1, this%gwfpackages(ip)%nbound
193  i = this%gwfpackages(ip)%nodelist(ib)
194  if (i <= 0) cycle
195  if (this%ibound(i) <= 0) cycle
196  qbnd = this%gwfpackages(ip)%get_flow(ib)
197  ! todo, after initial release: default iflowface values for different packages
198  iflowface = 0
199  if (iauxiflowface > 0) then
200  iflowface = nint(this%gwfpackages(ip)%auxvar(iauxiflowface, ib))
201  ! this maps bot -2 -> 9, top -1 -> 10; see note re: max faces below
202  if (iflowface < 0) iflowface = iflowface + max_poly_cells + 1
203  end if
204  if (iflowface .gt. 0) then
205  ioffset = (i - 1) * max_poly_cells
206  this%BoundaryFlows(ioffset + iflowface) = &
207  this%BoundaryFlows(ioffset + iflowface) + qbnd
208  else if (qbnd .gt. dzero) then
209  this%SourceFlows(i) = this%SourceFlows(i) + qbnd
210  else if (qbnd .lt. dzero) then
211  this%SinkFlows(i) = this%SinkFlows(i) + qbnd
212  end if
213  end do
214  end do
215 
216  end subroutine accumulate_flows
217 
218 end module prtfmimodule
integer(i4b), parameter, public max_poly_cells
Definition: Cell.f90:9
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dhdry
real dry cell constant
Definition: Constants.f90:94
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter lenauxname
maximum length of a aux variable
Definition: Constants.f90:35
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
This module defines variable data types.
Definition: kind.f90:8
subroutine, public fmi_cr(fmiobj, name_model, inunit, iout)
Create a new PrtFmi object.
Definition: prt-fmi.f90:38
subroutine accumulate_flows(this)
Accumulate flows.
Definition: prt-fmi.f90:159
subroutine fmi_ad(this)
Time step advance.
Definition: prt-fmi.f90:68
subroutine prtfmi_df(this, dis, idryinactive)
Define the flow model interface.
Definition: prt-fmi.f90:139
character(len=lenpackagename) text
Definition: prt-fmi.f90:17
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string