MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
prtfmimodule Module Reference

Data Types

type  prtfmitype
 

Functions/Subroutines

subroutine, public fmi_cr (fmiobj, name_model, input_mempath, inunit, iout)
 Create a new PrtFmi object. More...
 
subroutine fmi_ad (this)
 Time step advance. More...
 
subroutine prtfmi_df (this, dis, idryinactive)
 Define the flow model interface. More...
 
subroutine accumulate_flows (this)
 Accumulate flows. More...
 
subroutine mark_boundary_face (this, ic, iface)
 Mark a face as a boundary face. More...
 
logical(lgp) function is_boundary_face (this, ic, iface)
 Check if a face is assigned to a boundary package. More...
 
logical(lgp) function is_net_out_boundary_face (this, ic, iface)
 Check if a face is an assigned boundary with net outflow. More...
 
integer(i4b) function iflowface_to_icellface (this, iflowface)
 Convert an iflowface number to a cell face number. Maps bottom (-2) -> max_faces - 1, top (-1) -> max_faces. More...
 

Variables

character(len=lenpackagename) text = ' PRTFMI'
 

Function/Subroutine Documentation

◆ accumulate_flows()

subroutine prtfmimodule::accumulate_flows ( class(prtfmitype this)
private

Definition at line 153 of file prt-fmi.f90.

154  ! dummy
155  class(PrtFmiType) :: this
156  ! local
157  integer(I4B) :: j, i, ip, ib
158  integer(I4B) :: ioffset, iflowface, iauxiflowface, icellface
159  real(DP) :: qbnd
160  character(len=LENAUXNAME) :: auxname
161  integer(I4B) :: naux
162 
163  this%StorageFlows = dzero
164  if (this%igwfstrgss /= 0) &
165  this%StorageFlows = this%StorageFlows + this%gwfstrgss
166  if (this%igwfstrgsy /= 0) &
167  this%StorageFlows = this%StorageFlows + this%gwfstrgsy
168 
169  this%SourceFlows = dzero
170  this%SinkFlows = dzero
171  this%BoundaryFlows = dzero
172  this%BoundaryFaces = 0
173  do ip = 1, this%nflowpack
174  iauxiflowface = 0
175  naux = this%gwfpackages(ip)%naux
176  if (naux > 0) then
177  do j = 1, naux
178  auxname = this%gwfpackages(ip)%auxname(j)
179  if (trim(adjustl(auxname)) == "IFLOWFACE") then
180  iauxiflowface = j
181  exit
182  end if
183  end do
184  end if
185  do ib = 1, this%gwfpackages(ip)%nbound
186  i = this%gwfpackages(ip)%nodelist(ib)
187  if (i <= 0) cycle
188  if (this%ibound(i) <= 0) cycle
189  qbnd = this%gwfpackages(ip)%get_flow(ib)
190  ! todo, after initial release: default iflowface values for different packages
191  iflowface = 0
192  icellface = 0
193  if (iauxiflowface > 0) then
194  iflowface = nint(this%gwfpackages(ip)%auxvar(iauxiflowface, ib))
195  icellface = this%iflowface_to_icellface(iflowface)
196  end if
197  if (icellface > 0) then
198  call this%mark_boundary_face(i, icellface)
199  ioffset = (i - 1) * this%max_faces
200  this%BoundaryFlows(ioffset + icellface) = &
201  this%BoundaryFlows(ioffset + icellface) + qbnd
202  else if (qbnd .gt. dzero) then
203  this%SourceFlows(i) = this%SourceFlows(i) + qbnd
204  else if (qbnd .lt. dzero) then
205  this%SinkFlows(i) = this%SinkFlows(i) + qbnd
206  end if
207  end do
208  end do
209 

◆ fmi_ad()

subroutine prtfmimodule::fmi_ad ( class(prtfmitype this)
private

Definition at line 71 of file prt-fmi.f90.

72  ! modules
73  use constantsmodule, only: dhdry
74  ! dummy
75  class(PrtFmiType) :: this
76  ! local
77  integer(I4B) :: n
78  character(len=15) :: nodestr
79  character(len=*), parameter :: fmtdry = &
80  &"(/1X,'WARNING: DRY CELL ENCOUNTERED AT ',a,'; RESET AS INACTIVE')"
81  character(len=*), parameter :: fmtrewet = &
82  &"(/1X,'DRY CELL REACTIVATED AT ', a)"
83 
84  ! Set flag to indicated that flows are being updated. For the case where
85  ! flows may be reused (only when flows are read from a file) then set
86  ! the flag to zero to indicated that flows were not updated
87  this%iflowsupdated = 1
88 
89  ! If reading flows from a budget file, read the next set of records
90  if (this%iubud /= 0) call this%advance_bfr()
91 
92  ! If reading heads from a head file, read the next set of records
93  if (this%iuhds /= 0) call this%advance_hfr()
94 
95  ! If mover flows are being read from file, read the next set of records
96  if (this%iumvr /= 0) &
97  call this%mvrbudobj%bfr_advance(this%dis, this%iout)
98 
99  ! Accumulate flows
100  call this%accumulate_flows()
101 
102  ! if flow cell is dry, then set this%ibound = 0
103  do n = 1, this%dis%nodes
104  ! Calculate the ibound-like array that has 0 if saturation
105  ! is zero and 1 otherwise
106  if (this%gwfsat(n) > dzero) then
107  this%ibdgwfsat0(n) = 1
108  else
109  this%ibdgwfsat0(n) = 0
110  end if
111 
112  ! Check if active model cell is inactive for flow
113  if (this%ibound(n) > 0) then
114  if (this%gwfhead(n) == dhdry) then
115  ! cell should be made inactive
116  this%ibound(n) = 0
117  call this%dis%noder_to_string(n, nodestr)
118  write (this%iout, fmtdry) trim(nodestr)
119  end if
120  end if
121 
122  ! Convert dry model cell to active if flow has rewet
123  if (this%ibound(n) == 0) then
124  if (this%gwfhead(n) /= dhdry) then
125  ! cell is now wet
126  this%ibound(n) = 1
127  call this%dis%noder_to_string(n, nodestr)
128  write (this%iout, fmtrewet) trim(nodestr)
129  end if
130  end if
131  end do
132 
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dhdry
real dry cell constant
Definition: Constants.f90:94

◆ fmi_cr()

subroutine, public prtfmimodule::fmi_cr ( type(prtfmitype), pointer  fmiobj,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  input_mempath,
integer(i4b), intent(inout)  inunit,
integer(i4b), intent(in)  iout 
)

Definition at line 43 of file prt-fmi.f90.

44  ! dummy
45  type(PrtFmiType), pointer :: fmiobj
46  character(len=*), intent(in) :: name_model
47  character(len=*), intent(in) :: input_mempath
48  integer(I4B), intent(inout) :: inunit
49  integer(I4B), intent(in) :: iout
50 
51  ! Create the object
52  allocate (fmiobj)
53 
54  ! create name and memory path
55  call fmiobj%set_names(1, name_model, 'FMI', 'FMI', input_mempath)
56  fmiobj%text = text
57 
58  ! Allocate scalars
59  call fmiobj%allocate_scalars()
60 
61  ! Set variables
62  fmiobj%inunit = inunit
63  fmiobj%iout = iout
64 
65  ! Assign dependent variable label
66  fmiobj%depvartype = 'TRACKS '
67 
Here is the caller graph for this function:

◆ iflowface_to_icellface()

integer(i4b) function prtfmimodule::iflowface_to_icellface ( class(prtfmitype), intent(inout)  this,
integer(i4b), intent(in)  iflowface 
)
private

Definition at line 285 of file prt-fmi.f90.

286  class(PrtFmiType), intent(inout) :: this
287  integer(I4B), intent(in) :: iflowface
288  integer(I4B) :: iface
289 
290  iface = iflowface
291  if (iface < 0) iface = iface + this%max_faces + 1

◆ is_boundary_face()

logical(lgp) function prtfmimodule::is_boundary_face ( class(prtfmitype this,
integer(i4b), intent(in)  ic,
integer(i4b), intent(in)  iface 
)
private
Parameters
[in]icnode number (reduced)
[in]ifacecell face number

Definition at line 240 of file prt-fmi.f90.

241  class(PrtFmiType) :: this
242  integer(I4B), intent(in) :: ic !< node number (reduced)
243  integer(I4B), intent(in) :: iface !< cell face number
244  logical(LGP) :: is_boundary
245  ! local
246  integer(I4B) :: bit_pos
247 
248  is_boundary = .false.
249  if (ic <= 0 .or. ic > this%dis%nodes) then
250  print *, 'Invalid cell number: ', ic
251  print *, 'Expected a value in range [1, ', this%dis%nodes, ']'
252  call pstop(1)
253  end if
254  if (iface <= 0) then
255  print *, 'Invalid face number: ', iface
256  print *, 'Expected a value in range [1, ', this%max_faces, ']'
257  call pstop(1)
258  end if
259  bit_pos = iface - 1 ! bit position 0-based
260  if (bit_pos < 0 .or. bit_pos > 31) then
261  print *, 'Invalid bitmask position: ', iface
262  print *, 'Expected a value in range [0, 31]'
263  call pstop(1)
264  end if
265  is_boundary = btest(this%BoundaryFaces(ic), bit_pos)
Here is the call graph for this function:

◆ is_net_out_boundary_face()

logical(lgp) function prtfmimodule::is_net_out_boundary_face ( class(prtfmitype this,
integer(i4b), intent(in)  ic,
integer(i4b), intent(in)  iface 
)
private
Parameters
[in]icnode number (reduced)
[in]ifacecell face number

Definition at line 269 of file prt-fmi.f90.

270  class(PrtFmiType) :: this
271  integer(I4B), intent(in) :: ic !< node number (reduced)
272  integer(I4B), intent(in) :: iface !< cell face number
273  logical(LGP) :: is_net_out_boundary
274  ! local
275  integer(I4B) :: ioffset
276 
277  is_net_out_boundary = .false.
278  if (.not. this%is_boundary_face(ic, iface)) return
279  ioffset = (ic - 1) * this%max_faces
280  if (this%BoundaryFlows(ioffset + iface) < dzero) is_net_out_boundary = .true.

◆ mark_boundary_face()

subroutine prtfmimodule::mark_boundary_face ( class(prtfmitype this,
integer(i4b), intent(in)  ic,
integer(i4b), intent(in)  iface 
)
private
Parameters
[in]icnode number (reduced)
[in]ifacecell face number

Definition at line 213 of file prt-fmi.f90.

214  class(PrtFmiType) :: this
215  integer(I4B), intent(in) :: ic !< node number (reduced)
216  integer(I4B), intent(in) :: iface !< cell face number
217  ! local
218  integer(I4B) :: bit_pos
219 
220  if (ic <= 0 .or. ic > this%dis%nodes) then
221  print *, 'Invalid cell number: ', ic
222  print *, 'Expected a value in range [1, ', this%dis%nodes, ']'
223  call pstop(1)
224  end if
225  if (iface <= 0) then
226  print *, 'Invalid face number: ', iface
227  print *, 'Expected a value in range [1, ', this%max_faces, ']'
228  call pstop(1)
229  end if
230  bit_pos = iface - 1 ! bit position 0-based
231  if (bit_pos < 0 .or. bit_pos > 31) then
232  print *, 'Invalid bitmask position: ', iface
233  print *, 'Expected a value in range [0, 31]'
234  call pstop(1)
235  end if
236  this%BoundaryFaces(ic) = ibset(this%BoundaryFaces(ic), bit_pos)
Here is the call graph for this function:

◆ prtfmi_df()

subroutine prtfmimodule::prtfmi_df ( class(prtfmitype this,
class(disbasetype), intent(in), pointer  dis,
integer(i4b), intent(in)  idryinactive 
)

Definition at line 136 of file prt-fmi.f90.

137  class(PrtFmiType) :: this
138  class(DisBaseType), pointer, intent(in) :: dis
139  integer(I4B), intent(in) :: idryinactive
140 
141  call this%FlowModelInterfaceType%fmi_df(dis, idryinactive)
142 
143  this%max_faces = this%dis%get_max_npolyverts() + 2
144  allocate (this%StorageFlows(this%dis%nodes))
145  allocate (this%SourceFlows(this%dis%nodes))
146  allocate (this%SinkFlows(this%dis%nodes))
147  allocate (this%BoundaryFlows(this%dis%nodes * this%max_faces))
148  allocate (this%BoundaryFaces(this%dis%nodes))
149 

Variable Documentation

◆ text

character(len=lenpackagename) prtfmimodule::text = ' PRTFMI'
private

Definition at line 17 of file prt-fmi.f90.

17  character(len=LENPACKAGENAME) :: text = ' PRTFMI'