MODFLOW 6  version 6.7.0.dev1
USGS Modular Hydrologic Model
tspfmimodule Module Reference

Data Types

type  dataadvancedpackagetype
 
type  budobjptrarray
 
type  tspfmitype
 

Functions/Subroutines

subroutine, public fmi_cr (fmiobj, name_model, inunit, iout, eqnsclfac, depvartype)
 Create a new FMI object. More...
 
subroutine fmi_rp (this, inmvr)
 Read and prepare. More...
 
subroutine fmi_ad (this, cnew)
 Advance routine for FMI object. More...
 
subroutine fmi_fc (this, nodes, cold, nja, matrix_sln, idxglo, rhs)
 Calculate coefficients and fill matrix and rhs terms associated with FMI object. More...
 
subroutine fmi_cq (this, cnew, flowja)
 Calculate flow correction. More...
 
subroutine fmi_bd (this, isuppress_output, model_budget)
 Calculate budget terms associated with FMI object. More...
 
subroutine fmi_ot_flow (this, icbcfl, icbcun)
 Save budget terms associated with FMI object. More...
 
subroutine gwtfmi_da (this)
 Deallocate variables. More...
 
subroutine gwtfmi_allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine gwtfmi_allocate_arrays (this, nodes)
 @ brief Allocate arrays for FMI object More...
 
subroutine set_active_status (this, cnew)
 Set gwt transport cell status. More...
 
real(dp) function gwfsatold (this, n, delt)
 Calculate the previous saturation level. More...
 
subroutine gwtfmi_read_options (this)
 Read options from input file. More...
 
subroutine gwtfmi_read_packagedata (this)
 Read PACKAGEDATA block. More...
 
subroutine set_aptbudobj_pointer (this, name, budobjptr)
 Set the pointer to a budget object. More...
 
subroutine initialize_gwfterms_from_bfr (this)
 Initialize the groundwater flow terms based on the budget file reader. More...
 
subroutine initialize_gwfterms_from_gwfbndlist (this)
 Initialize groundwater flow terms from the groundwater budget. More...
 
subroutine gwtfmi_allocate_gwfpackages (this, ngwfterms)
 Initialize an array for storing PackageBudget objects. More...
 
subroutine gwtfmi_deallocate_gwfpackages (this)
 Deallocate memory. More...
 

Variables

character(len=lenpackagename) text = ' GWTFMI'
 
integer(i4b), parameter nbditems = 2
 
character(len=lenbudtxt), dimension(nbditemsbudtxt
 

Function/Subroutine Documentation

◆ fmi_ad()

subroutine tspfmimodule::fmi_ad ( class(tspfmitype this,
real(dp), dimension(:), intent(inout)  cnew 
)

Definition at line 138 of file tsp-fmi.f90.

139  ! -- modules
140  use constantsmodule, only: dhdry
141  ! -- dummy
142  class(TspFmiType) :: this
143  real(DP), intent(inout), dimension(:) :: cnew
144  ! -- local
145  integer(I4B) :: n
146  character(len=*), parameter :: fmtdry = &
147  &"(/1X,'WARNING: DRY CELL ENCOUNTERED AT ',a,'; RESET AS INACTIVE &
148  &WITH DRY CONCENTRATION = ', G13.5)"
149  character(len=*), parameter :: fmtrewet = &
150  &"(/1X,'DRY CELL REACTIVATED AT ', a,&
151  &' WITH STARTING CONCENTRATION =',G13.5)"
152  !
153  ! -- Set flag to indicated that flows are being updated. For the case where
154  ! flows may be reused (only when flows are read from a file) then set
155  ! the flag to zero to indicate that flows were not updated
156  this%iflowsupdated = 1
157  !
158  ! -- If reading flows from a budget file, read the next set of records
159  if (this%iubud /= 0) then
160  call this%advance_bfr()
161  end if
162  !
163  ! -- If reading heads from a head file, read the next set of records
164  if (this%iuhds /= 0) then
165  call this%advance_hfr()
166  end if
167  !
168  ! -- If mover flows are being read from file, read the next set of records
169  if (this%iumvr /= 0) then
170  call this%mvrbudobj%bfr_advance(this%dis, this%iout)
171  end if
172  !
173  ! -- If advanced package flows are being read from file, read the next set of records
174  if (this%flows_from_file .and. this%inunit /= 0) then
175  do n = 1, size(this%aptbudobj)
176  call this%aptbudobj(n)%ptr%bfr_advance(this%dis, this%iout)
177  end do
178  end if
179  !
180  ! -- set inactive transport cell status
181  if (this%idryinactive /= 0) then
182  call this%set_active_status(cnew)
183  end if
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dhdry
real dry cell constant
Definition: Constants.f90:94

◆ fmi_bd()

subroutine tspfmimodule::fmi_bd ( class(tspfmitype this,
integer(i4b), intent(in)  isuppress_output,
type(budgettype), intent(inout)  model_budget 
)
private

Definition at line 251 of file tsp-fmi.f90.

252  ! -- modules
253  use tdismodule, only: delt
255  ! -- dummy
256  class(TspFmiType) :: this
257  integer(I4B), intent(in) :: isuppress_output
258  type(BudgetType), intent(inout) :: model_budget
259  ! -- local
260  real(DP) :: rin
261  real(DP) :: rout
262  !
263  ! -- flow correction
264  if (this%iflowerr /= 0) then
265  call rate_accumulator(this%flowcorrect, rin, rout)
266  call model_budget%addentry(rin, rout, delt, budtxt(2), isuppress_output)
267  end if
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
Derived type for the Budget object.
Definition: Budget.f90:39
Here is the call graph for this function:

◆ fmi_cq()

subroutine tspfmimodule::fmi_cq ( class(tspfmitype this,
real(dp), dimension(:), intent(in)  cnew,
real(dp), dimension(:), intent(inout), contiguous  flowja 
)
private

Where there is a flow imbalance for a given cell, a correction may be applied if selected

Definition at line 222 of file tsp-fmi.f90.

223  ! -- modules
224  ! -- dummy
225  class(TspFmiType) :: this
226  real(DP), intent(in), dimension(:) :: cnew
227  real(DP), dimension(:), contiguous, intent(inout) :: flowja
228  ! -- local
229  integer(I4B) :: n
230  integer(I4B) :: idiag
231  real(DP) :: rate
232  !
233  ! -- If not adding flow error correction, return
234  if (this%iflowerr /= 0) then
235  !
236  ! -- Accumulate the flow correction term
237  do n = 1, this%dis%nodes
238  rate = dzero
239  idiag = this%dis%con%ia(n)
240  if (this%ibound(n) > 0) then
241  rate = -this%gwfflowja(idiag) * cnew(n) * this%eqnsclfac
242  end if
243  this%flowcorrect(n) = rate
244  flowja(idiag) = flowja(idiag) + rate
245  end do
246  end if

◆ fmi_cr()

subroutine, public tspfmimodule::fmi_cr ( type(tspfmitype), pointer  fmiobj,
character(len=*), intent(in)  name_model,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
real(dp), intent(in), pointer  eqnsclfac,
character(len=lenvarname), intent(in)  depvartype 
)
Parameters
[in]eqnsclfacgoverning equation scale factor

Definition at line 74 of file tsp-fmi.f90.

75  ! -- dummy
76  type(TspFmiType), pointer :: fmiobj
77  character(len=*), intent(in) :: name_model
78  integer(I4B), intent(in) :: inunit
79  integer(I4B), intent(in) :: iout
80  real(DP), intent(in), pointer :: eqnsclfac !< governing equation scale factor
81  character(len=LENVARNAME), intent(in) :: depvartype
82  !
83  ! -- Create the object
84  allocate (fmiobj)
85  !
86  ! -- create name and memory path
87  call fmiobj%set_names(1, name_model, 'FMI', 'FMI')
88  fmiobj%text = text
89  !
90  ! -- Allocate scalars
91  call fmiobj%allocate_scalars()
92  !
93  ! -- Set variables
94  fmiobj%inunit = inunit
95  fmiobj%iout = iout
96  !
97  ! -- Initialize block parser
98  call fmiobj%parser%Initialize(fmiobj%inunit, fmiobj%iout)
99  !
100  ! -- Assign label based on dependent variable
101  fmiobj%depvartype = depvartype
102  !
103  ! -- Store pointer to governing equation scale factor
104  fmiobj%eqnsclfac => eqnsclfac
Here is the caller graph for this function:

◆ fmi_fc()

subroutine tspfmimodule::fmi_fc ( class(tspfmitype this,
integer, intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cold,
integer(i4b), intent(in)  nja,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(nja), intent(in)  idxglo,
real(dp), dimension(nodes), intent(inout)  rhs 
)

Definition at line 189 of file tsp-fmi.f90.

190  ! -- dummy
191  class(TspFmiType) :: this
192  integer, intent(in) :: nodes
193  real(DP), intent(in), dimension(nodes) :: cold
194  integer(I4B), intent(in) :: nja
195  class(MatrixBaseType), pointer :: matrix_sln
196  integer(I4B), intent(in), dimension(nja) :: idxglo
197  real(DP), intent(inout), dimension(nodes) :: rhs
198  ! -- local
199  integer(I4B) :: n, idiag, idiag_sln
200  real(DP) :: qcorr
201  !
202  ! -- Calculate the flow imbalance error and make a correction for it
203  if (this%iflowerr /= 0) then
204  !
205  ! -- Correct the transport solution for the flow imbalance by adding
206  ! the flow residual to the diagonal
207  do n = 1, nodes
208  idiag = this%dis%con%ia(n)
209  idiag_sln = idxglo(idiag)
210  !call matrix_sln%add_value_pos(idiag_sln, -this%gwfflowja(idiag))
211  qcorr = -this%gwfflowja(idiag) * this%eqnsclfac
212  call matrix_sln%add_value_pos(idiag_sln, qcorr)
213  end do
214  end if

◆ fmi_ot_flow()

subroutine tspfmimodule::fmi_ot_flow ( class(tspfmitype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  icbcun 
)

Definition at line 272 of file tsp-fmi.f90.

273  ! -- dummy
274  class(TspFmiType) :: this
275  integer(I4B), intent(in) :: icbcfl
276  integer(I4B), intent(in) :: icbcun
277  ! -- local
278  integer(I4B) :: ibinun
279  integer(I4B) :: iprint, nvaluesp, nwidthp
280  character(len=1) :: cdatafmp = ' ', editdesc = ' '
281  real(DP) :: dinact
282  !
283  ! -- Set unit number for binary output
284  if (this%ipakcb < 0) then
285  ibinun = icbcun
286  elseif (this%ipakcb == 0) then
287  ibinun = 0
288  else
289  ibinun = this%ipakcb
290  end if
291  if (icbcfl == 0) ibinun = 0
292  !
293  ! -- Do not save flow corrections if not active
294  if (this%iflowerr == 0) ibinun = 0
295  !
296  ! -- Record the storage rates if requested
297  if (ibinun /= 0) then
298  iprint = 0
299  dinact = dzero
300  !
301  ! -- flow correction
302  call this%dis%record_array(this%flowcorrect, this%iout, iprint, -ibinun, &
303  budtxt(2), cdatafmp, nvaluesp, &
304  nwidthp, editdesc, dinact)
305  end if

◆ fmi_rp()

subroutine tspfmimodule::fmi_rp ( class(tspfmitype this,
integer(i4b), intent(in)  inmvr 
)
private

Definition at line 109 of file tsp-fmi.f90.

110  ! -- modules
111  use tdismodule, only: kper, kstp
112  ! -- dummy
113  class(TspFmiType) :: this
114  integer(I4B), intent(in) :: inmvr
115  ! -- local
116  ! -- formats
117  !
118  ! --Check to make sure MVT Package is active if mvr flows are available.
119  ! This cannot be checked until RP because exchange doesn't set a pointer
120  ! to mvrbudobj until exg_ar().
121  if (kper * kstp == 1) then
122  if (associated(this%mvrbudobj) .and. inmvr == 0) then
123  write (errmsg, '(a)') 'GWF water mover is active but the GWT MVT &
124  &package has not been specified. activate GWT MVT package.'
125  call store_error(errmsg, terminate=.true.)
126  end if
127  if (.not. associated(this%mvrbudobj) .and. inmvr > 0) then
128  write (errmsg, '(a)') 'GWF water mover terms are not available &
129  &but the GWT MVT package has been activated. Activate GWF-GWT &
130  &exchange or specify GWFMOVER in FMI PACKAGEDATA.'
131  call store_error(errmsg, terminate=.true.)
132  end if
133  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:

◆ gwfsatold()

real(dp) function tspfmimodule::gwfsatold ( class(tspfmitype this,
integer(i4b), intent(in)  n,
real(dp), intent(in)  delt 
)

Calculate the groundwater cell head saturation for the end of the last time step

Definition at line 504 of file tsp-fmi.f90.

505  ! -- modules
506  ! -- dummy
507  class(TspFmiType) :: this
508  integer(I4B), intent(in) :: n
509  real(DP), intent(in) :: delt
510  ! -- result
511  real(DP) :: satold
512  ! -- local
513  real(DP) :: vcell
514  real(DP) :: vnew
515  real(DP) :: vold
516  !
517  ! -- calculate the value
518  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
519  vnew = vcell * this%gwfsat(n)
520  vold = vnew
521  if (this%igwfstrgss /= 0) vold = vold + this%gwfstrgss(n) * delt
522  if (this%igwfstrgsy /= 0) vold = vold + this%gwfstrgsy(n) * delt
523  satold = vold / vcell

◆ gwtfmi_allocate_arrays()

subroutine tspfmimodule::gwtfmi_allocate_arrays ( class(tspfmitype this,
integer(i4b), intent(in)  nodes 
)

Method to allocate arrays for the FMI package.

Definition at line 391 of file tsp-fmi.f90.

393  ! -- modules
394  use constantsmodule, only: dzero
395  ! -- dummy
396  class(TspFmiType) :: this
397  integer(I4B), intent(in) :: nodes
398  ! -- local
399  integer(I4B) :: n
400  !
401  ! -- allocate parent arrays
402  call this%FlowModelInterfaceType%allocate_arrays(nodes)
403  !
404  ! -- Allocate variables needed for all cases
405  if (this%iflowerr == 0) then
406  call mem_allocate(this%flowcorrect, 1, 'FLOWCORRECT', this%memoryPath)
407  else
408  call mem_allocate(this%flowcorrect, nodes, 'FLOWCORRECT', this%memoryPath)
409  end if
410  do n = 1, size(this%flowcorrect)
411  this%flowcorrect(n) = dzero
412  end do
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ gwtfmi_allocate_gwfpackages()

subroutine tspfmimodule::gwtfmi_allocate_gwfpackages ( class(tspfmitype this,
integer(i4b), intent(in)  ngwfterms 
)

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

Definition at line 905 of file tsp-fmi.f90.

906  ! -- modules
907  use constantsmodule, only: lenmempath
909  ! -- dummy
910  class(TspFmiType) :: this
911  integer(I4B), intent(in) :: ngwfterms
912  ! -- local
913  integer(I4B) :: n
914  character(len=LENMEMPATH) :: memPath
915  !
916  ! -- direct allocate
917  allocate (this%gwfpackages(ngwfterms))
918  allocate (this%flowpacknamearray(ngwfterms))
919  allocate (this%datp(ngwfterms))
920  !
921  ! -- mem_allocate
922  call mem_allocate(this%iatp, ngwfterms, 'IATP', this%memoryPath)
923  call mem_allocate(this%igwfmvrterm, ngwfterms, 'IGWFMVRTERM', this%memoryPath)
924  !
925  ! -- initialize
926  this%nflowpack = ngwfterms
927  do n = 1, this%nflowpack
928  this%iatp(n) = 0
929  this%igwfmvrterm(n) = 0
930  this%flowpacknamearray(n) = ''
931  !
932  ! -- Create a mempath for each individual flow package data set
933  ! of the form, MODELNAME/FMI-FTn
934  write (mempath, '(a, i0)') trim(this%memoryPath)//'-FT', n
935  call this%gwfpackages(n)%initialize(mempath)
936  end do
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27

◆ gwtfmi_allocate_scalars()

subroutine tspfmimodule::gwtfmi_allocate_scalars ( class(tspfmitype this)

Allocate scalar variables for an FMI object

Definition at line 366 of file tsp-fmi.f90.

367  ! -- modules
369  ! -- dummy
370  class(TspFmiType) :: this
371  ! -- local
372  !
373  ! -- allocate scalars in parent
374  call this%FlowModelInterfaceType%allocate_scalars()
375  !
376  ! -- Allocate
377  call mem_allocate(this%iflowerr, 'IFLOWERR', this%memoryPath)
378  !
379  ! -- Although not a scalar, allocate the advanced package transport
380  ! budget object to zero so that it can be dynamically resized later
381  allocate (this%aptbudobj(0))
382  !
383  ! -- Initialize
384  this%iflowerr = 0

◆ gwtfmi_da()

subroutine tspfmimodule::gwtfmi_da ( class(tspfmitype this)
private

Deallocate memory associated with FMI object

Definition at line 312 of file tsp-fmi.f90.

313  ! -- modules
315  ! -- dummy
316  class(TspFmiType) :: this
317  ! -- todo: finalize hfr and bfr either here or in a finalize routine
318  !
319  ! -- deallocate any memory stored with gwfpackages
320  call this%deallocate_gwfpackages()
321  !
322  ! -- deallocate fmi arrays
323  if (associated(this%datp)) then
324  deallocate (this%datp)
325  deallocate (this%gwfpackages)
326  deallocate (this%flowpacknamearray)
327  call mem_deallocate(this%iatp)
328  call mem_deallocate(this%igwfmvrterm)
329  end if
330 
331  deallocate (this%aptbudobj)
332  call mem_deallocate(this%flowcorrect)
333  call mem_deallocate(this%ibdgwfsat0)
334  if (this%flows_from_file) then
335  call mem_deallocate(this%gwfstrgss)
336  call mem_deallocate(this%gwfstrgsy)
337  end if
338  !
339  ! -- special treatment, these could be from mem_checkin
340  call mem_deallocate(this%gwfhead, 'GWFHEAD', this%memoryPath)
341  call mem_deallocate(this%gwfsat, 'GWFSAT', this%memoryPath)
342  call mem_deallocate(this%gwfspdis, 'GWFSPDIS', this%memoryPath)
343  call mem_deallocate(this%gwfflowja, 'GWFFLOWJA', this%memoryPath)
344  !
345  ! -- deallocate scalars
346  call mem_deallocate(this%flows_from_file)
347  call mem_deallocate(this%iflowsupdated)
348  call mem_deallocate(this%iflowerr)
349  call mem_deallocate(this%igwfstrgss)
350  call mem_deallocate(this%igwfstrgsy)
351  call mem_deallocate(this%iubud)
352  call mem_deallocate(this%iuhds)
353  call mem_deallocate(this%iumvr)
354  call mem_deallocate(this%iugrb)
355  call mem_deallocate(this%nflowpack)
356  call mem_deallocate(this%idryinactive)
357  !
358  ! -- deallocate parent
359  call this%NumericalPackageType%da()

◆ gwtfmi_deallocate_gwfpackages()

subroutine tspfmimodule::gwtfmi_deallocate_gwfpackages ( class(tspfmitype this)

Deallocate memory that stores the gwfpackages array

Definition at line 943 of file tsp-fmi.f90.

944  ! -- modules
945  ! -- dummy
946  class(TspFmiType) :: this
947  ! -- local
948  integer(I4B) :: n
949  !
950  ! -- initialize
951  do n = 1, this%nflowpack
952  call this%gwfpackages(n)%da()
953  end do

◆ gwtfmi_read_options()

subroutine tspfmimodule::gwtfmi_read_options ( class(tspfmitype this)
private

Definition at line 528 of file tsp-fmi.f90.

529  ! -- modules
530  use constantsmodule, only: linelength, dem6
533  ! -- dummy
534  class(TspFmiType) :: this
535  ! -- local
536  character(len=LINELENGTH) :: keyword
537  integer(I4B) :: ierr
538  logical :: isfound, endOfBlock
539  character(len=*), parameter :: fmtisvflow = &
540  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
541  &WHENEVER ICBCFL IS NOT ZERO AND FLOW IMBALANCE CORRECTION ACTIVE.')"
542  character(len=*), parameter :: fmtifc = &
543  &"(4x,'MASS WILL BE ADDED OR REMOVED TO COMPENSATE FOR FLOW IMBALANCE.')"
544  !
545  ! -- get options block
546  call this%parser%GetBlock('OPTIONS', isfound, ierr, blockrequired=.false., &
547  supportopenclose=.true.)
548  !
549  ! -- parse options block if detected
550  if (isfound) then
551  write (this%iout, '(1x,a)') 'PROCESSING FMI OPTIONS'
552  do
553  call this%parser%GetNextLine(endofblock)
554  if (endofblock) exit
555  call this%parser%GetStringCaps(keyword)
556  select case (keyword)
557  case ('SAVE_FLOWS')
558  this%ipakcb = -1
559  write (this%iout, fmtisvflow)
560  case ('FLOW_IMBALANCE_CORRECTION')
561  write (this%iout, fmtifc)
562  this%iflowerr = 1
563  case default
564  write (errmsg, '(a,a)') 'Unknown FMI option: ', &
565  trim(keyword)
566  call store_error(errmsg)
567  call this%parser%StoreErrorUnit()
568  end select
569  end do
570  write (this%iout, '(1x,a)') 'END OF FMI OPTIONS'
571  end if
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
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
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
Here is the call graph for this function:

◆ gwtfmi_read_packagedata()

subroutine tspfmimodule::gwtfmi_read_packagedata ( class(tspfmitype this)

Read packagedata block from input file

Definition at line 578 of file tsp-fmi.f90.

579  ! -- modules
580  use openspecmodule, only: access, form
584  ! -- dummy
585  class(TspFmiType) :: this
586  ! -- local
587  type(BudgetObjectType), pointer :: budobjptr
588  character(len=LINELENGTH) :: keyword, fname
589  character(len=LENPACKAGENAME) :: pname
590  integer(I4B) :: i
591  integer(I4B) :: ierr
592  integer(I4B) :: inunit
593  integer(I4B) :: iapt
594  logical :: isfound, endOfBlock
595  logical :: blockrequired
596  logical :: exist
597  type(BudObjPtrArray), dimension(:), allocatable :: tmpbudobj
598  !
599  ! -- initialize
600  iapt = 0
601  blockrequired = .true.
602  !
603  ! -- get options block
604  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
605  blockrequired=blockrequired, &
606  supportopenclose=.true.)
607  !
608  ! -- parse options block if detected
609  if (isfound) then
610  write (this%iout, '(1x,a)') 'PROCESSING FMI PACKAGEDATA'
611  do
612  call this%parser%GetNextLine(endofblock)
613  if (endofblock) exit
614  call this%parser%GetStringCaps(keyword)
615  select case (keyword)
616  case ('GWFBUDGET')
617  call this%parser%GetStringCaps(keyword)
618  if (keyword /= 'FILEIN') then
619  call store_error('GWFBUDGET keyword must be followed by '// &
620  '"FILEIN" then by filename.')
621  call this%parser%StoreErrorUnit()
622  end if
623  call this%parser%GetString(fname)
624  inunit = getunit()
625  inquire (file=trim(fname), exist=exist)
626  if (.not. exist) then
627  call store_error('Could not find file '//trim(fname))
628  call this%parser%StoreErrorUnit()
629  end if
630  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
631  access, 'UNKNOWN')
632  this%iubud = inunit
633  call this%initialize_bfr()
634  case ('GWFHEAD')
635  call this%parser%GetStringCaps(keyword)
636  if (keyword /= 'FILEIN') then
637  call store_error('GWFHEAD keyword must be followed by '// &
638  '"FILEIN" then by filename.')
639  call this%parser%StoreErrorUnit()
640  end if
641  call this%parser%GetString(fname)
642  inquire (file=trim(fname), exist=exist)
643  if (.not. exist) then
644  call store_error('Could not find file '//trim(fname))
645  call this%parser%StoreErrorUnit()
646  end if
647  inunit = getunit()
648  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
649  access, 'UNKNOWN')
650  this%iuhds = inunit
651  call this%initialize_hfr()
652  case ('GWFMOVER')
653  call this%parser%GetStringCaps(keyword)
654  if (keyword /= 'FILEIN') then
655  call store_error('GWFMOVER keyword must be followed by '// &
656  '"FILEIN" then by filename.')
657  call this%parser%StoreErrorUnit()
658  end if
659  call this%parser%GetString(fname)
660  inunit = getunit()
661  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
662  access, 'UNKNOWN')
663  this%iumvr = inunit
664  call budgetobject_cr_bfr(this%mvrbudobj, 'MVT', this%iumvr, &
665  this%iout)
666  call this%mvrbudobj%fill_from_bfr(this%dis, this%iout)
667  case default
668  !
669  ! --expand the size of aptbudobj, which stores a pointer to the budobj
670  allocate (tmpbudobj(iapt))
671  do i = 1, size(this%aptbudobj)
672  tmpbudobj(i)%ptr => this%aptbudobj(i)%ptr
673  end do
674  deallocate (this%aptbudobj)
675  allocate (this%aptbudobj(iapt + 1))
676  do i = 1, size(tmpbudobj)
677  this%aptbudobj(i)%ptr => tmpbudobj(i)%ptr
678  end do
679  deallocate (tmpbudobj)
680  !
681  ! -- Open the budget file and start filling it
682  iapt = iapt + 1
683  pname = keyword(1:lenpackagename)
684  call this%parser%GetStringCaps(keyword)
685  if (keyword /= 'FILEIN') then
686  call store_error('Package name must be followed by '// &
687  '"FILEIN" then by filename.')
688  call this%parser%StoreErrorUnit()
689  end if
690  call this%parser%GetString(fname)
691  inunit = getunit()
692  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
693  access, 'UNKNOWN')
694  call budgetobject_cr_bfr(budobjptr, pname, inunit, &
695  this%iout, colconv2=['GWF '])
696  call budobjptr%fill_from_bfr(this%dis, this%iout)
697  this%aptbudobj(iapt)%ptr => budobjptr
698  end select
699  end do
700  write (this%iout, '(1x,a)') 'END OF FMI PACKAGEDATA'
701  end if
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Here is the call graph for this function:

◆ initialize_gwfterms_from_bfr()

subroutine tspfmimodule::initialize_gwfterms_from_bfr ( class(tspfmitype this)
private

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

Definition at line 735 of file tsp-fmi.f90.

736  ! -- modules
739  ! -- dummy
740  class(TspFmiType) :: this
741  ! -- local
742  integer(I4B) :: nflowpack
743  integer(I4B) :: i, ip
744  integer(I4B) :: naux
745  logical :: found_flowja
746  logical :: found_dataspdis
747  logical :: found_datasat
748  logical :: found_stoss
749  logical :: found_stosy
750  integer(I4B), dimension(:), allocatable :: imap
751  !
752  ! -- Calculate the number of gwf flow packages
753  allocate (imap(this%bfr%nbudterms))
754  imap(:) = 0
755  nflowpack = 0
756  found_flowja = .false.
757  found_dataspdis = .false.
758  found_datasat = .false.
759  found_stoss = .false.
760  found_stosy = .false.
761  do i = 1, this%bfr%nbudterms
762  select case (trim(adjustl(this%bfr%budtxtarray(i))))
763  case ('FLOW-JA-FACE')
764  found_flowja = .true.
765  case ('DATA-SPDIS')
766  found_dataspdis = .true.
767  case ('DATA-SAT')
768  found_datasat = .true.
769  case ('STO-SS')
770  found_stoss = .true.
771  this%igwfstrgss = 1
772  case ('STO-SY')
773  found_stosy = .true.
774  this%igwfstrgsy = 1
775  case default
776  nflowpack = nflowpack + 1
777  imap(i) = 1
778  end select
779  end do
780  !
781  ! -- allocate gwfpackage arrays (gwfpackages, iatp, datp, ...)
782  call this%allocate_gwfpackages(nflowpack)
783  !
784  ! -- Copy the package name and aux names from budget file reader
785  ! to the gwfpackages derived-type variable
786  ip = 1
787  do i = 1, this%bfr%nbudterms
788  if (imap(i) == 0) cycle
789  call this%gwfpackages(ip)%set_name(this%bfr%dstpackagenamearray(i), &
790  this%bfr%budtxtarray(i))
791  naux = this%bfr%nauxarray(i)
792  call this%gwfpackages(ip)%set_auxname(naux, &
793  this%bfr%auxtxtarray(1:naux, i))
794  ip = ip + 1
795  end do
796  !
797  ! -- Copy just the package names for the boundary packages into
798  ! the flowpacknamearray
799  ip = 1
800  do i = 1, size(imap)
801  if (imap(i) == 1) then
802  this%flowpacknamearray(ip) = this%bfr%dstpackagenamearray(i)
803  ip = ip + 1
804  end if
805  end do
806  !
807  ! -- Error if specific discharge, saturation or flowja not found
808  if (.not. found_dataspdis) then
809  write (errmsg, '(a)') 'Specific discharge not found in &
810  &budget file. SAVE_SPECIFIC_DISCHARGE and &
811  &SAVE_FLOWS must be activated in the NPF package.'
812  call store_error(errmsg)
813  end if
814  if (.not. found_datasat) then
815  write (errmsg, '(a)') 'Saturation not found in &
816  &budget file. SAVE_SATURATION and &
817  &SAVE_FLOWS must be activated in the NPF package.'
818  call store_error(errmsg)
819  end if
820  if (.not. found_flowja) then
821  write (errmsg, '(a)') 'FLOWJA not found in &
822  &budget file. SAVE_FLOWS must &
823  &be activated in the NPF package.'
824  call store_error(errmsg)
825  end if
826  if (count_errors() > 0) then
827  call this%parser%StoreErrorUnit()
828  end if
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
Here is the call graph for this function:

◆ initialize_gwfterms_from_gwfbndlist()

subroutine tspfmimodule::initialize_gwfterms_from_gwfbndlist ( class(tspfmitype this)

Flows are coming from a gwf-gwt exchange object

Definition at line 835 of file tsp-fmi.f90.

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

◆ set_active_status()

subroutine tspfmimodule::set_active_status ( class(tspfmitype this,
real(dp), dimension(:), intent(inout)  cnew 
)

Dry GWF cells are treated differently by GWT and GWE. Transport does not occur in deactivated GWF cells; however, GWE still simulates conduction through dry cells.

Definition at line 421 of file tsp-fmi.f90.

422  ! -- modules
423  use constantsmodule, only: dhdry
424  ! -- dummy
425  class(TspFmiType) :: this
426  real(DP), intent(inout), dimension(:) :: cnew
427  ! -- local
428  integer(I4B) :: n
429  integer(I4B) :: m
430  integer(I4B) :: ipos
431  real(DP) :: crewet, tflow, flownm
432  character(len=15) :: nodestr
433  ! -- formats
434  character(len=*), parameter :: fmtoutmsg1 = &
435  "(1x,'WARNING: DRY CELL ENCOUNTERED AT ', a,'; RESET AS INACTIVE WITH &
436  &DRY ', a, '=', G13.5)"
437  character(len=*), parameter :: fmtoutmsg2 = &
438  &"(1x,'DRY CELL REACTIVATED AT', a, 'WITH STARTING', a, '=', G13.5)"
439  !
440  do n = 1, this%dis%nodes
441  ! -- Calculate the ibound-like array that has 0 if saturation
442  ! is zero and 1 otherwise
443  if (this%gwfsat(n) > dzero) then
444  this%ibdgwfsat0(n) = 1
445  else
446  this%ibdgwfsat0(n) = 0
447  end if
448  !
449  ! -- Check if active transport cell is inactive for flow
450  if (this%ibound(n) > 0) then
451  if (this%gwfhead(n) == dhdry) then
452  ! -- transport cell should be made inactive
453  this%ibound(n) = 0
454  cnew(n) = dhdry
455  call this%dis%noder_to_string(n, nodestr)
456  write (this%iout, fmtoutmsg1) &
457  trim(nodestr), trim(adjustl(this%depvartype)), dhdry
458  end if
459  end if
460  end do
461  !
462  ! -- if flow cell is dry, then set gwt%ibound = 0 and conc to dry
463  do n = 1, this%dis%nodes
464  !
465  ! -- Convert dry transport cell to active if flow has rewet
466  if (cnew(n) == dhdry) then
467  if (this%gwfhead(n) /= dhdry) then
468  !
469  ! -- obtain weighted concentration/temperature
470  crewet = dzero
471  tflow = dzero
472  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
473  m = this%dis%con%ja(ipos)
474  flownm = this%gwfflowja(ipos)
475  if (flownm > 0) then
476  if (this%ibound(m) /= 0) then
477  crewet = crewet + cnew(m) * flownm ! kluge note: apparently no need to multiply flows by eqnsclfac
478  tflow = tflow + this%gwfflowja(ipos) ! since it will divide out below anyway
479  end if
480  end if
481  end do
482  if (tflow > dzero) then
483  crewet = crewet / tflow
484  else
485  crewet = dzero
486  end if
487  !
488  ! -- cell is now wet
489  this%ibound(n) = 1
490  cnew(n) = crewet
491  call this%dis%noder_to_string(n, nodestr)
492  write (this%iout, fmtoutmsg2) &
493  trim(nodestr), trim(adjustl(this%depvartype)), crewet
494  end if
495  end if
496  end do

◆ set_aptbudobj_pointer()

subroutine tspfmimodule::set_aptbudobj_pointer ( class(tspfmitype this,
character(len=*), intent(in)  name,
type(budgetobjecttype), pointer  budobjptr 
)

An advanced transport can pass in a name and a pointer budget object, and this routine will look through the budget objects managed by FMI and point to the one with the same name, such as LAK-1, SFR-1, etc.

Definition at line 711 of file tsp-fmi.f90.

712  ! -- modules
713  class(TspFmiType) :: this
714  ! -- dumm
715  character(len=*), intent(in) :: name
716  type(BudgetObjectType), pointer :: budobjptr
717  ! -- local
718  integer(I4B) :: i
719  !
720  ! -- find and set the pointer
721  do i = 1, size(this%aptbudobj)
722  if (this%aptbudobj(i)%ptr%name == name) then
723  budobjptr => this%aptbudobj(i)%ptr
724  exit
725  end if
726  end do

Variable Documentation

◆ budtxt

character(len=lenbudtxt), dimension(nbditems) tspfmimodule::budtxt
private

Definition at line 25 of file tsp-fmi.f90.

25  character(len=LENBUDTXT), dimension(NBDITEMS) :: budtxt

◆ nbditems

integer(i4b), parameter tspfmimodule::nbditems = 2
private

Definition at line 24 of file tsp-fmi.f90.

24  integer(I4B), parameter :: NBDITEMS = 2

◆ text

character(len=lenpackagename) tspfmimodule::text = ' GWTFMI'
private

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

22  character(len=LENPACKAGENAME) :: text = ' GWTFMI'