MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
prtfmimodule Module Reference

Data Types

type  prtfmitype
 

Enumerations

enum  
 IFLOWFACE numbers for top and bottom faces. More...
 

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, icellface)
 Mark a face as a boundary face. More...
 
logical(lgp) function is_boundary_face (this, ic, icellface)
 Check if a face is assigned to a boundary package. More...
 
logical(lgp) function is_net_out_boundary_face (this, ic, icellface)
 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'
 
@, public iflowface_top = -1
 
@, public iflowface_bottom = -2
 

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
private

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

Function/Subroutine Documentation

◆ accumulate_flows()

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

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

170  ! dummy
171  class(PrtFmiType) :: this
172  ! local
173  integer(I4B) :: j, i, ip, ib
174  integer(I4B) :: iflowface, iauxiflowface, icellface
175  real(DP) :: qbnd
176  character(len=LENAUXNAME) :: auxname
177  integer(I4B) :: naux
178 
179  this%StorageFlows = dzero
180  if (this%igwfstrgss /= 0) &
181  this%StorageFlows = this%StorageFlows + this%gwfstrgss
182  if (this%igwfstrgsy /= 0) &
183  this%StorageFlows = this%StorageFlows + this%gwfstrgsy
184 
185  this%SourceFlows = dzero
186  this%SinkFlows = dzero
187  this%BoundaryFlows = dzero
188  this%BoundaryFaces = 0
189  do ip = 1, this%nflowpack
190  iauxiflowface = 0
191  naux = this%gwfpackages(ip)%naux
192  if (naux > 0) then
193  do j = 1, naux
194  auxname = this%gwfpackages(ip)%auxname(j)
195  if (trim(adjustl(auxname)) == "IFLOWFACE") then
196  iauxiflowface = j
197  exit
198  end if
199  end do
200  end if
201  do ib = 1, this%gwfpackages(ip)%nbound
202  i = this%gwfpackages(ip)%nodelist(ib)
203  if (i <= 0) cycle
204  if (this%ibound(i) <= 0) cycle
205  qbnd = this%gwfpackages(ip)%get_flow(ib)
206  ! todo, after initial release: default iflowface values for different packages
207  iflowface = 0
208  icellface = 0
209  if (iauxiflowface > 0) then
210  iflowface = nint(this%gwfpackages(ip)%auxvar(iauxiflowface, ib))
211  icellface = this%iflowface_to_icellface(iflowface)
212  end if
213  if (icellface > 0) then
214  call this%mark_boundary_face(i, icellface)
215  this%BoundaryFlows(i, icellface) = &
216  this%BoundaryFlows(i, icellface) + qbnd
217  else if (qbnd .gt. dzero) then
218  this%SourceFlows(i) = this%SourceFlows(i) + qbnd
219  else if (qbnd .lt. dzero) then
220  this%SinkFlows(i) = this%SinkFlows(i) + qbnd
221  end if
222  end do
223  end do
224 

◆ fmi_ad()

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

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

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

51  ! dummy
52  type(PrtFmiType), pointer :: fmiobj
53  character(len=*), intent(in) :: name_model
54  character(len=*), intent(in) :: input_mempath
55  integer(I4B), intent(inout) :: inunit
56  integer(I4B), intent(in) :: iout
57 
58  ! Create the object
59  allocate (fmiobj)
60 
61  ! create name and memory path
62  call fmiobj%set_names(1, name_model, 'FMI', 'FMI', input_mempath)
63  fmiobj%text = text
64 
65  ! Allocate scalars
66  call fmiobj%allocate_scalars()
67 
68  ! Set variables
69  fmiobj%inunit = inunit
70  fmiobj%iout = iout
71 
72  ! Assign dependent variable label
73  fmiobj%depvartype = 'TRACKS '
74 
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 299 of file prt-fmi.f90.

300  class(PrtFmiType), intent(inout) :: this
301  integer(I4B), intent(in) :: iflowface
302  integer(I4B) :: icellface
303 
304  icellface = iflowface
305  if (icellface < 0) icellface = icellface + this%max_faces - iflowface_top

◆ is_boundary_face()

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

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

256  class(PrtFmiType) :: this
257  integer(I4B), intent(in) :: ic !< node number (reduced)
258  integer(I4B), intent(in) :: icellface !< cell face number
259  logical(LGP) :: is_boundary
260  ! local
261  integer(I4B) :: bit_pos
262 
263  is_boundary = .false.
264  if (ic <= 0 .or. ic > this%dis%nodes) then
265  print *, 'Invalid cell number: ', ic
266  print *, 'Expected a value in range [1, ', this%dis%nodes, ']'
267  call pstop(1)
268  end if
269  if (icellface <= 0) then
270  print *, 'Invalid face number: ', icellface
271  print *, 'Expected a value in range [1, ', this%max_faces, ']'
272  call pstop(1)
273  end if
274  bit_pos = icellface - 1 ! bit position 0-based
275  if (bit_pos < 0 .or. bit_pos > 31) then
276  print *, 'Invalid bitmask position: ', bit_pos
277  print *, 'Expected a value in range [0, 31]'
278  call pstop(1)
279  end if
280  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)  icellface 
)
private
Parameters
[in]icnode number (reduced)
[in]icellfacecell face number

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

286  class(PrtFmiType) :: this
287  integer(I4B), intent(in) :: ic !< node number (reduced)
288  integer(I4B), intent(in) :: icellface !< cell face number
289  logical(LGP) :: is_net_out_boundary
290 
291  is_net_out_boundary = .false.
292  if (.not. this%is_boundary_face(ic, icellface)) return
293  if (this%BoundaryFlows(ic, icellface) < dzero) &
294  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)  icellface 
)
private
Parameters
[in]icnode number (reduced)
[in]icellfacecell face number

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

229  class(PrtFmiType) :: this
230  integer(I4B), intent(in) :: ic !< node number (reduced)
231  integer(I4B), intent(in) :: icellface !< cell face number
232  ! local
233  integer(I4B) :: bit_pos
234 
235  if (ic <= 0 .or. ic > this%dis%nodes) then
236  print *, 'Invalid cell number: ', ic
237  print *, 'Expected a value in range [1, ', this%dis%nodes, ']'
238  call pstop(1)
239  end if
240  if (icellface <= 0) then
241  print *, 'Invalid face number: ', icellface
242  print *, 'Expected a value in range [1, ', this%max_faces, ']'
243  call pstop(1)
244  end if
245  bit_pos = icellface - 1 ! bit position 0-based
246  if (bit_pos < 0 .or. bit_pos > 31) then
247  print *, 'Invalid bitmask position: ', bit_pos
248  print *, 'Expected a value in range [0, 31]'
249  call pstop(1)
250  end if
251  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 143 of file prt-fmi.f90.

144  class(PrtFmiType) :: this
145  class(DisBaseType), pointer, intent(in) :: dis
146  integer(I4B), intent(in) :: idryinactive
147 
148  call this%FlowModelInterfaceType%fmi_df(dis, idryinactive)
149 
150  this%max_faces = this%dis%get_max_npolyverts() + 2
151  if (this%max_faces > 32) then
152  write (errmsg, '(a,i0,a,i0,a)') &
153  'DISV grid contains a cell with ', this%max_faces - 2, &
154  ' lateral faces. Cells may have at most 30 lateral faces.'
155  call store_error(errmsg)
156  call this%parser%StoreErrorUnit()
157  return
158  end if
159 
160  allocate (this%StorageFlows(this%dis%nodes))
161  allocate (this%SourceFlows(this%dis%nodes))
162  allocate (this%SinkFlows(this%dis%nodes))
163  allocate (this%BoundaryFlows(this%dis%nodes, this%max_faces))
164  allocate (this%BoundaryFaces(this%dis%nodes))
165 
Here is the call graph for this function:

Variable Documentation

◆ iflowface_bottom

@, public prtfmimodule::iflowface_bottom = -2

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

23  enumerator :: IFLOWFACE_BOTTOM = -2

◆ iflowface_top

@, public prtfmimodule::iflowface_top = -1

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

22  enumerator :: IFLOWFACE_TOP = -1

◆ text

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

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

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