MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
tsp-fmi.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b, lgp
7  use simvariablesmodule, only: errmsg
9  use basedismodule, only: disbasetype
10  use listmodule, only: listtype
16 
17  implicit none
18  private
19  public :: tspfmitype
20  public :: fmi_cr
21 
22  character(len=LENPACKAGENAME) :: text = ' GWTFMI'
23 
24  integer(I4B), parameter :: nbditems = 2
25  character(len=LENBUDTXT), dimension(NBDITEMS) :: budtxt
26  data budtxt/' FLOW-ERROR', ' FLOW-CORRECTION'/
27 
29  real(dp), dimension(:), contiguous, pointer :: concpack => null()
30  real(dp), dimension(:), contiguous, pointer :: qmfrommvr => null()
31  end type
32 
34  type(budgetobjecttype), pointer :: ptr
35  end type budobjptrarray
36 
38 
39  integer(I4B), dimension(:), pointer, contiguous :: iatp => null() !< advanced transport package applied to gwfpackages
40  integer(I4B), pointer :: iflowerr => null() !< add the flow error correction
41  real(dp), dimension(:), pointer, contiguous :: flowcorrect => null() !< mass flow correction
42  real(dp), pointer :: eqnsclfac => null() !< governing equation scale factor; =1. for solute; =rhow*cpw for energy
44  dimension(:), pointer, contiguous :: datp => null()
45  type(budobjptrarray), dimension(:), allocatable :: aptbudobj !< flow budget objects for the advanced packages
46 
47  contains
48 
49  procedure :: allocate_arrays => gwtfmi_allocate_arrays
50  procedure :: allocate_gwfpackages => gwtfmi_allocate_gwfpackages
51  procedure :: allocate_scalars => gwtfmi_allocate_scalars
52  procedure :: deallocate_gwfpackages => gwtfmi_deallocate_gwfpackages
53  procedure :: fmi_rp
54  procedure :: fmi_ad
55  procedure :: fmi_fc
56  procedure :: fmi_cq
57  procedure :: fmi_bd
58  procedure :: fmi_ot_flow
59  procedure :: fmi_da => gwtfmi_da
60  procedure :: gwfsatold
63  procedure :: source_options => gwtfmi_source_options
64  procedure :: set_aptbudobj_pointer
65  procedure :: source_packagedata => gwtfmi_source_packagedata
66  procedure :: set_active_status
67 
68  end type tspfmitype
69 
70 contains
71 
72  !> @brief Create a new FMI object
73  !<
74  subroutine fmi_cr(fmiobj, name_model, input_mempath, inunit, iout, eqnsclfac, &
75  depvartype)
76  ! -- dummy
77  type(tspfmitype), pointer :: fmiobj
78  character(len=*), intent(in) :: name_model
79  character(len=*), intent(in) :: input_mempath
80  integer(I4B), intent(in) :: inunit
81  integer(I4B), intent(in) :: iout
82  real(dp), intent(in), pointer :: eqnsclfac !< governing equation scale factor
83  character(len=LENVARNAME), intent(in) :: depvartype
84  !
85  ! -- Create the object
86  allocate (fmiobj)
87  !
88  ! -- create name and memory path
89  call fmiobj%set_names(1, name_model, 'FMI', 'FMI', input_mempath)
90  fmiobj%text = text
91  !
92  ! -- Allocate scalars
93  call fmiobj%allocate_scalars()
94  !
95  ! -- Set variables
96  fmiobj%inunit = inunit
97  fmiobj%iout = iout
98  !
99  ! -- Assign label based on dependent variable
100  fmiobj%depvartype = depvartype
101  !
102  ! -- Store pointer to governing equation scale factor
103  fmiobj%eqnsclfac => eqnsclfac
104  end subroutine fmi_cr
105 
106  !> @brief Read and prepare
107  !<
108  subroutine fmi_rp(this, inmvr)
109  ! -- modules
110  use tdismodule, only: kper, kstp
111  ! -- dummy
112  class(tspfmitype) :: this
113  integer(I4B), intent(in) :: inmvr
114  ! -- local
115  ! -- formats
116  !
117  ! --Check to make sure MVT Package is active if mvr flows are available.
118  ! This cannot be checked until RP because exchange doesn't set a pointer
119  ! to mvrbudobj until exg_ar().
120  if (kper * kstp == 1) then
121  if (associated(this%mvrbudobj) .and. inmvr == 0) then
122  write (errmsg, '(a)') 'GWF water mover is active but the GWT MVT &
123  &package has not been specified. activate GWT MVT package.'
124  call store_error(errmsg, terminate=.true.)
125  end if
126  if (.not. associated(this%mvrbudobj) .and. inmvr > 0) then
127  write (errmsg, '(a)') 'GWF water mover terms are not available &
128  &but the GWT MVT package has been activated. Activate GWF-GWT &
129  &exchange or specify GWFMOVER in FMI PACKAGEDATA.'
130  call store_error(errmsg, terminate=.true.)
131  end if
132  end if
133  end subroutine fmi_rp
134 
135  !> @brief Advance routine for FMI object
136  !<
137  subroutine fmi_ad(this, cnew)
138  ! -- modules
139  use constantsmodule, only: dhdry
140  ! -- dummy
141  class(tspfmitype) :: this
142  real(DP), intent(inout), dimension(:) :: cnew
143  ! -- local
144  integer(I4B) :: n
145  character(len=*), parameter :: fmtdry = &
146  &"(/1X,'WARNING: DRY CELL ENCOUNTERED AT ',a,'; RESET AS INACTIVE &
147  &WITH DRY CONCENTRATION = ', G13.5)"
148  character(len=*), parameter :: fmtrewet = &
149  &"(/1X,'DRY CELL REACTIVATED AT ', a,&
150  &' WITH STARTING CONCENTRATION =',G13.5)"
151  !
152  ! -- Set flag to indicated that flows are being updated. For the case where
153  ! flows may be reused (only when flows are read from a file) then set
154  ! the flag to zero to indicate that flows were not updated
155  this%iflowsupdated = 1
156  !
157  ! -- If reading flows from a budget file, read the next set of records
158  if (this%iubud /= 0) then
159  call this%advance_bfr()
160  end if
161  !
162  ! -- If reading heads from a head file, read the next set of records
163  if (this%iuhds /= 0) then
164  call this%advance_hfr()
165  end if
166  !
167  ! -- If mover flows are being read from file, read the next set of records
168  if (this%iumvr /= 0) then
169  call this%mvrbudobj%bfr_advance(this%dis, this%iout)
170  end if
171  !
172  ! -- If advanced package flows are being read from file, read the next set of records
173  if (this%flows_from_file .and. this%inunit /= 0) then
174  do n = 1, size(this%aptbudobj)
175  call this%aptbudobj(n)%ptr%bfr_advance(this%dis, this%iout)
176  end do
177  end if
178  !
179  ! -- set inactive transport cell status
180  if (this%idryinactive /= 0) then
181  call this%set_active_status(cnew)
182  end if
183  end subroutine fmi_ad
184 
185  !> @brief Calculate coefficients and fill matrix and rhs terms associated
186  !! with FMI object
187  !<
188  subroutine fmi_fc(this, nodes, cold, nja, matrix_sln, idxglo, rhs)
189  ! -- dummy
190  class(tspfmitype) :: this
191  integer, intent(in) :: nodes
192  real(DP), intent(in), dimension(nodes) :: cold
193  integer(I4B), intent(in) :: nja
194  class(matrixbasetype), pointer :: matrix_sln
195  integer(I4B), intent(in), dimension(nja) :: idxglo
196  real(DP), intent(inout), dimension(nodes) :: rhs
197  ! -- local
198  integer(I4B) :: n, idiag, idiag_sln
199  real(DP) :: qcorr
200  !
201  ! -- Calculate the flow imbalance error and make a correction for it
202  if (this%iflowerr /= 0) then
203  !
204  ! -- Correct the transport solution for the flow imbalance by adding
205  ! the flow residual to the diagonal
206  do n = 1, nodes
207  idiag = this%dis%con%ia(n)
208  idiag_sln = idxglo(idiag)
209  !call matrix_sln%add_value_pos(idiag_sln, -this%gwfflowja(idiag))
210  qcorr = -this%gwfflowja(idiag) * this%eqnsclfac
211  call matrix_sln%add_value_pos(idiag_sln, qcorr)
212  end do
213  end if
214  end subroutine fmi_fc
215 
216  !> @brief Calculate flow correction
217  !!
218  !! Where there is a flow imbalance for a given cell, a correction may be
219  !! applied if selected
220  !<
221  subroutine fmi_cq(this, cnew, flowja)
222  ! -- modules
223  ! -- dummy
224  class(tspfmitype) :: this
225  real(DP), intent(in), dimension(:) :: cnew
226  real(DP), dimension(:), contiguous, intent(inout) :: flowja
227  ! -- local
228  integer(I4B) :: n
229  integer(I4B) :: idiag
230  real(DP) :: rate
231  !
232  ! -- If not adding flow error correction, return
233  if (this%iflowerr /= 0) then
234  !
235  ! -- Accumulate the flow correction term
236  do n = 1, this%dis%nodes
237  rate = dzero
238  idiag = this%dis%con%ia(n)
239  if (this%ibound(n) > 0) then
240  rate = -this%gwfflowja(idiag) * cnew(n) * this%eqnsclfac
241  end if
242  this%flowcorrect(n) = rate
243  flowja(idiag) = flowja(idiag) + rate
244  end do
245  end if
246  end subroutine fmi_cq
247 
248  !> @brief Calculate budget terms associated with FMI object
249  !<
250  subroutine fmi_bd(this, isuppress_output, model_budget)
251  ! -- modules
252  use tdismodule, only: delt
254  ! -- dummy
255  class(tspfmitype) :: this
256  integer(I4B), intent(in) :: isuppress_output
257  type(budgettype), intent(inout) :: model_budget
258  ! -- local
259  real(DP) :: rin
260  real(DP) :: rout
261  !
262  ! -- flow correction
263  if (this%iflowerr /= 0) then
264  call rate_accumulator(this%flowcorrect, rin, rout)
265  call model_budget%addentry(rin, rout, delt, budtxt(2), isuppress_output)
266  end if
267  end subroutine fmi_bd
268 
269  !> @brief Save budget terms associated with FMI object
270  !<
271  subroutine fmi_ot_flow(this, icbcfl, icbcun)
272  ! -- dummy
273  class(tspfmitype) :: this
274  integer(I4B), intent(in) :: icbcfl
275  integer(I4B), intent(in) :: icbcun
276  ! -- local
277  integer(I4B) :: ibinun
278  integer(I4B) :: iprint, nvaluesp, nwidthp
279  character(len=1) :: cdatafmp = ' ', editdesc = ' '
280  real(DP) :: dinact
281  !
282  ! -- Set unit number for binary output
283  if (this%ipakcb < 0) then
284  ibinun = icbcun
285  elseif (this%ipakcb == 0) then
286  ibinun = 0
287  else
288  ibinun = this%ipakcb
289  end if
290  if (icbcfl == 0) ibinun = 0
291  !
292  ! -- Do not save flow corrections if not active
293  if (this%iflowerr == 0) ibinun = 0
294  !
295  ! -- Record the storage rates if requested
296  if (ibinun /= 0) then
297  iprint = 0
298  dinact = dzero
299  !
300  ! -- flow correction
301  call this%dis%record_array(this%flowcorrect, this%iout, iprint, -ibinun, &
302  budtxt(2), cdatafmp, nvaluesp, &
303  nwidthp, editdesc, dinact)
304  end if
305  end subroutine fmi_ot_flow
306 
307  !> @brief Deallocate variables
308  !!
309  !! Deallocate memory associated with FMI object
310  !<
311  subroutine gwtfmi_da(this)
312  ! -- modules
314  ! -- dummy
315  class(tspfmitype) :: this
316  ! -- todo: finalize hfr and bfr either here or in a finalize routine
317  !
318  ! -- deallocate any memory stored with gwfpackages
319  call this%deallocate_gwfpackages()
320  !
321  ! -- deallocate fmi arrays
322  if (associated(this%datp)) then
323  deallocate (this%datp)
324  deallocate (this%gwfpackages)
325  deallocate (this%flowpacknamearray)
326  call mem_deallocate(this%iatp)
327  call mem_deallocate(this%igwfmvrterm)
328  end if
329 
330  deallocate (this%aptbudobj)
331  call mem_deallocate(this%flowcorrect)
332  call mem_deallocate(this%ibdgwfsat0)
333  if (this%flows_from_file) then
334  call mem_deallocate(this%gwfstrgss)
335  call mem_deallocate(this%gwfstrgsy)
336  end if
337  !
338  ! -- special treatment, these could be from mem_checkin
339  call mem_deallocate(this%gwfhead, 'GWFHEAD', this%memoryPath)
340  call mem_deallocate(this%gwfsat, 'GWFSAT', this%memoryPath)
341  call mem_deallocate(this%gwfspdis, 'GWFSPDIS', this%memoryPath)
342  call mem_deallocate(this%gwfflowja, 'GWFFLOWJA', this%memoryPath)
343  !
344  ! -- deallocate scalars
345  call mem_deallocate(this%flows_from_file)
346  call mem_deallocate(this%iflowsupdated)
347  call mem_deallocate(this%iflowerr)
348  call mem_deallocate(this%igwfstrgss)
349  call mem_deallocate(this%igwfstrgsy)
350  call mem_deallocate(this%iubud)
351  call mem_deallocate(this%iuhds)
352  call mem_deallocate(this%iumvr)
353  call mem_deallocate(this%iugrb)
354  call mem_deallocate(this%nflowpack)
355  call mem_deallocate(this%idryinactive)
356  !
357  ! -- deallocate parent
358  call this%NumericalPackageType%da()
359  end subroutine gwtfmi_da
360 
361  !> @ brief Allocate scalars
362  !!
363  !! Allocate scalar variables for an FMI object
364  !<
365  subroutine gwtfmi_allocate_scalars(this)
366  ! -- modules
368  ! -- dummy
369  class(tspfmitype) :: this
370  ! -- local
371  !
372  ! -- allocate scalars in parent
373  call this%FlowModelInterfaceType%allocate_scalars()
374  !
375  ! -- Allocate
376  call mem_allocate(this%iflowerr, 'IFLOWERR', this%memoryPath)
377  !
378  ! -- Although not a scalar, allocate the advanced package transport
379  ! budget object to zero so that it can be dynamically resized later
380  allocate (this%aptbudobj(0))
381  !
382  ! -- Initialize
383  this%iflowerr = 0
384  end subroutine gwtfmi_allocate_scalars
385 
386  !> @ brief Allocate arrays for FMI object
387  !!
388  !! Method to allocate arrays for the FMI package.
389  !<
390  subroutine gwtfmi_allocate_arrays(this, nodes)
392  ! -- modules
393  use constantsmodule, only: dzero
394  ! -- dummy
395  class(tspfmitype) :: this
396  integer(I4B), intent(in) :: nodes
397  ! -- local
398  integer(I4B) :: n
399  !
400  ! -- allocate parent arrays
401  call this%FlowModelInterfaceType%allocate_arrays(nodes)
402  !
403  ! -- Allocate variables needed for all cases
404  if (this%iflowerr == 0) then
405  call mem_allocate(this%flowcorrect, 1, 'FLOWCORRECT', this%memoryPath)
406  else
407  call mem_allocate(this%flowcorrect, nodes, 'FLOWCORRECT', this%memoryPath)
408  end if
409  do n = 1, size(this%flowcorrect)
410  this%flowcorrect(n) = dzero
411  end do
412  end subroutine gwtfmi_allocate_arrays
413 
414  !> @brief Set gwt transport cell status
415  !!
416  !! Dry GWF cells are treated differently by GWT and GWE. Transport does not
417  !! occur in deactivated GWF cells; however, GWE still simulates conduction
418  !! through dry cells.
419  !<
420  subroutine set_active_status(this, cnew)
421  ! -- modules
422  use constantsmodule, only: dhdry
423  ! -- dummy
424  class(tspfmitype) :: this
425  real(DP), intent(inout), dimension(:) :: cnew
426  ! -- local
427  integer(I4B) :: n
428  integer(I4B) :: m
429  integer(I4B) :: ipos
430  real(DP) :: crewet, tflow, flownm
431  character(len=15) :: nodestr
432  ! -- formats
433  character(len=*), parameter :: fmtoutmsg1 = &
434  "(1x,'WARNING: DRY CELL ENCOUNTERED AT ', a,'; RESET AS INACTIVE WITH &
435  &DRY ', a, '=', G13.5)"
436  character(len=*), parameter :: fmtoutmsg2 = &
437  &"(1x,'DRY CELL REACTIVATED AT', a, 'WITH STARTING', a, '=', G13.5)"
438  !
439  do n = 1, this%dis%nodes
440  ! -- Calculate the ibound-like array that has 0 if saturation
441  ! is zero and 1 otherwise
442  if (this%gwfsat(n) > dzero) then
443  this%ibdgwfsat0(n) = 1
444  else
445  this%ibdgwfsat0(n) = 0
446  end if
447  !
448  ! -- Check if active transport cell is inactive for flow
449  if (this%ibound(n) > 0) then
450  if (this%gwfhead(n) == dhdry) then
451  ! -- transport cell should be made inactive
452  this%ibound(n) = 0
453  cnew(n) = dhdry
454  call this%dis%noder_to_string(n, nodestr)
455  write (this%iout, fmtoutmsg1) &
456  trim(nodestr), trim(adjustl(this%depvartype)), dhdry
457  end if
458  end if
459  end do
460  !
461  ! -- if flow cell is dry, then set gwt%ibound = 0 and conc to dry
462  do n = 1, this%dis%nodes
463  !
464  ! -- Convert dry transport cell to active if flow has rewet
465  if (cnew(n) == dhdry) then
466  if (this%gwfhead(n) /= dhdry) then
467  !
468  ! -- obtain weighted concentration/temperature
469  crewet = dzero
470  tflow = dzero
471  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
472  m = this%dis%con%ja(ipos)
473  flownm = this%gwfflowja(ipos)
474  if (flownm > 0) then
475  if (this%ibound(m) /= 0) then
476  crewet = crewet + cnew(m) * flownm ! kluge note: apparently no need to multiply flows by eqnsclfac
477  tflow = tflow + this%gwfflowja(ipos) ! since it will divide out below anyway
478  end if
479  end if
480  end do
481  if (tflow > dzero) then
482  crewet = crewet / tflow
483  else
484  crewet = dzero
485  end if
486  !
487  ! -- cell is now wet
488  this%ibound(n) = 1
489  cnew(n) = crewet
490  call this%dis%noder_to_string(n, nodestr)
491  write (this%iout, fmtoutmsg2) &
492  trim(nodestr), trim(adjustl(this%depvartype)), crewet
493  end if
494  end if
495  end do
496  end subroutine set_active_status
497 
498  !> @brief Calculate the previous saturation level
499  !!
500  !! Calculate the groundwater cell head saturation for the end of
501  !! the last time step
502  !<
503  function gwfsatold(this, n, delt) result(satold)
504  ! -- modules
505  ! -- dummy
506  class(tspfmitype) :: this
507  integer(I4B), intent(in) :: n
508  real(dp), intent(in) :: delt
509  ! -- result
510  real(dp) :: satold
511  ! -- local
512  real(dp) :: vcell
513  real(dp) :: vnew
514  real(dp) :: vold
515  !
516  ! -- calculate the value
517  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
518  vnew = vcell * this%gwfsat(n)
519  vold = vnew
520  if (this%igwfstrgss /= 0) vold = vold + this%gwfstrgss(n) * delt
521  if (this%igwfstrgsy /= 0) vold = vold + this%gwfstrgsy(n) * delt
522  satold = vold / vcell
523  end function gwfsatold
524 
525  !> @ brief Source input options for package
526  !<
527  subroutine gwtfmi_source_options(this)
528  ! -- modules
530  ! -- dummy
531  class(tspfmitype) :: this
532  ! -- local
533  logical(LGP) :: found_ipakcb, found_flowerr
534  character(len=*), parameter :: fmtisvflow = &
535  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
536  &WHENEVER ICBCFL IS NOT ZERO AND FLOW IMBALANCE CORRECTION ACTIVE.')"
537  character(len=*), parameter :: fmtifc = &
538  &"(4x,'MASS WILL BE ADDED OR REMOVED TO COMPENSATE FOR FLOW IMBALANCE.')"
539 
540  write (this%iout, '(1x,a)') 'PROCESSING FMI OPTIONS'
541 
542  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', this%input_mempath, &
543  found_ipakcb)
544  call mem_set_value(this%iflowerr, 'IMBALANCECORRECT', this%input_mempath, &
545  found_flowerr)
546 
547  if (found_ipakcb) then
548  this%ipakcb = -1
549  write (this%iout, fmtisvflow)
550  end if
551  if (found_flowerr) write (this%iout, fmtifc)
552 
553  write (this%iout, '(1x,a)') 'END OF FMI OPTIONS'
554  end subroutine gwtfmi_source_options
555 
556  !> @ brief Source input options for package
557  !<
558  subroutine gwtfmi_source_packagedata(this)
559  ! -- modules
563  use openspecmodule, only: access, form
566  ! -- dummy
567  class(tspfmitype) :: this
568  ! -- local
569  type(characterstringtype), dimension(:), contiguous, &
570  pointer :: flowtypes
571  type(characterstringtype), dimension(:), contiguous, &
572  pointer :: fileops
573  type(characterstringtype), dimension(:), contiguous, &
574  pointer :: fnames
575  type(budgetobjecttype), pointer :: budobjptr
576  type(budobjptrarray), dimension(:), allocatable :: tmpbudobj
577  character(len=LINELENGTH) :: flowtype, fileop, fname
578  integer(I4B) :: iapt, inunit, n, i
579  logical(LGP) :: exist
580 
581  iapt = 0
582 
583  call mem_setptr(flowtypes, 'FLOWTYPE', this%input_mempath)
584  call mem_setptr(fileops, 'FILEIN', this%input_mempath)
585  call mem_setptr(fnames, 'FNAME', this%input_mempath)
586 
587  do n = 1, size(flowtypes)
588  flowtype = flowtypes(n)
589  fileop = fileops(n)
590  fname = fnames(n)
591 
592  inquire (file=trim(fname), exist=exist)
593  if (.not. exist) then
594  call store_error('Could not find file '//trim(fname))
595  cycle
596  end if
597 
598  if (fileop /= 'FILEIN') then
599  call store_error('Unexpected packagedata input keyword read: "' &
600  //trim(fileop)//'".')
601  cycle
602  end if
603 
604  select case (flowtype)
605  case ('GWFBUDGET')
606  inunit = getunit()
607  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
608  access, 'UNKNOWN')
609  this%iubud = inunit
610  call this%initialize_bfr()
611  case ('GWFHEAD')
612  inunit = getunit()
613  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
614  access, 'UNKNOWN')
615  this%iuhds = inunit
616  call this%initialize_hfr()
617  case ('GWFMOVER')
618  inunit = getunit()
619  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
620  access, 'UNKNOWN')
621  this%iumvr = inunit
622  call budgetobject_cr_bfr(this%mvrbudobj, 'MVT', this%iumvr, &
623  this%iout)
624  call this%mvrbudobj%fill_from_bfr(this%dis, this%iout)
625  case default
626  !
627  ! --expand the size of aptbudobj, which stores a pointer to the budobj
628  allocate (tmpbudobj(iapt))
629  do i = 1, size(this%aptbudobj)
630  tmpbudobj(i)%ptr => this%aptbudobj(i)%ptr
631  end do
632  deallocate (this%aptbudobj)
633  allocate (this%aptbudobj(iapt + 1))
634  do i = 1, size(tmpbudobj)
635  this%aptbudobj(i)%ptr => tmpbudobj(i)%ptr
636  end do
637  deallocate (tmpbudobj)
638  !
639  ! -- Open the budget file and start filling it
640  iapt = iapt + 1
641  inunit = getunit()
642  call openfile(inunit, this%iout, fname, 'DATA(BINARY)', form, &
643  access, 'UNKNOWN')
644  call budgetobject_cr_bfr(budobjptr, flowtype, inunit, &
645  this%iout, colconv2=['GWF '])
646  call budobjptr%fill_from_bfr(this%dis, this%iout)
647  this%aptbudobj(iapt)%ptr => budobjptr
648  end select
649  end do
650 
651  if (count_errors() > 0) then
652  call store_error_filename(this%input_fname)
653  end if
654 
655  call memorystore_release('FLOWTYPE', this%input_mempath)
656  call memorystore_release('FILEIN', this%input_mempath)
657  call memorystore_release('FNAME', this%input_mempath)
658  end subroutine gwtfmi_source_packagedata
659 
660  !> @brief Set the pointer to a budget object
661  !!
662  !! An advanced transport can pass in a name and a
663  !! pointer budget object, and this routine will look through the budget
664  !! objects managed by FMI and point to the one with the same name, such as
665  !! LAK-1, SFR-1, etc.
666  !<
667  subroutine set_aptbudobj_pointer(this, name, budobjptr)
668  ! -- modules
669  class(tspfmitype) :: this
670  ! -- dumm
671  character(len=*), intent(in) :: name
672  type(budgetobjecttype), pointer :: budobjptr
673  ! -- local
674  integer(I4B) :: i
675  !
676  ! -- find and set the pointer
677  do i = 1, size(this%aptbudobj)
678  if (this%aptbudobj(i)%ptr%name == name) then
679  budobjptr => this%aptbudobj(i)%ptr
680  exit
681  end if
682  end do
683  end subroutine set_aptbudobj_pointer
684 
685  !> @brief Initialize the groundwater flow terms based on the budget file
686  !! reader
687  !!
688  !! Initialize terms and figure out how many different terms and packages
689  !! are contained within the file
690  !<
692  ! -- modules
694  ! -- dummy
695  class(tspfmitype) :: this
696  ! -- local
697  integer(I4B) :: nflowpack
698  integer(I4B) :: i, ip
699  integer(I4B) :: naux
700  logical :: found_flowja
701  logical :: found_dataspdis
702  logical :: found_datasat
703  logical :: found_stoss
704  logical :: found_stosy
705  integer(I4B), dimension(:), allocatable :: imap
706  !
707  ! -- Calculate the number of gwf flow packages
708  allocate (imap(this%bfr%nbudterms))
709  imap(:) = 0
710  nflowpack = 0
711  found_flowja = .false.
712  found_dataspdis = .false.
713  found_datasat = .false.
714  found_stoss = .false.
715  found_stosy = .false.
716  do i = 1, this%bfr%nbudterms
717  select case (trim(adjustl(this%bfr%budtxtarray(i))))
718  case ('FLOW-JA-FACE')
719  found_flowja = .true.
720  case ('DATA-SPDIS')
721  found_dataspdis = .true.
722  case ('DATA-SAT')
723  found_datasat = .true.
724  case ('STO-SS')
725  found_stoss = .true.
726  this%igwfstrgss = 1
727  case ('STO-SY')
728  found_stosy = .true.
729  this%igwfstrgsy = 1
730  case default
731  nflowpack = nflowpack + 1
732  imap(i) = 1
733  end select
734  end do
735  !
736  ! -- allocate gwfpackage arrays (gwfpackages, iatp, datp, ...)
737  call this%allocate_gwfpackages(nflowpack)
738  !
739  ! -- Copy the package name and aux names from budget file reader
740  ! to the gwfpackages derived-type variable
741  ip = 1
742  do i = 1, this%bfr%nbudterms
743  if (imap(i) == 0) cycle
744  call this%gwfpackages(ip)%set_name(this%bfr%dstpackagenamearray(i), &
745  this%bfr%budtxtarray(i))
746  naux = this%bfr%nauxarray(i)
747  call this%gwfpackages(ip)%set_auxname(naux, &
748  this%bfr%auxtxtarray(1:naux, i))
749  ip = ip + 1
750  end do
751  !
752  ! -- Copy just the package names for the boundary packages into
753  ! the flowpacknamearray
754  ip = 1
755  do i = 1, size(imap)
756  if (imap(i) == 1) then
757  this%flowpacknamearray(ip) = this%bfr%dstpackagenamearray(i)
758  ip = ip + 1
759  end if
760  end do
761  !
762  ! -- Error if specific discharge, saturation or flowja not found
763  if (.not. found_dataspdis) then
764  write (errmsg, '(a)') 'Specific discharge not found in &
765  &budget file. SAVE_SPECIFIC_DISCHARGE and &
766  &SAVE_FLOWS must be activated in the NPF package.'
767  call store_error(errmsg)
768  end if
769  if (.not. found_datasat) then
770  write (errmsg, '(a)') 'Saturation not found in &
771  &budget file. SAVE_SATURATION and &
772  &SAVE_FLOWS must be activated in the NPF package.'
773  call store_error(errmsg)
774  end if
775  if (.not. found_flowja) then
776  write (errmsg, '(a)') 'FLOWJA not found in &
777  &budget file. SAVE_FLOWS must &
778  &be activated in the NPF package.'
779  call store_error(errmsg)
780  end if
781  if (count_errors() > 0) then
782  call store_error_filename(this%input_fname)
783  end if
784  end subroutine initialize_gwfterms_from_bfr
785 
786  !> @brief Initialize groundwater flow terms from the groundwater budget
787  !!
788  !! Flows are coming from a gwf-gwt exchange object
789  !<
791  ! -- modules
792  use bndmodule, only: bndtype, getbndfromlist
793  ! -- dummy
794  class(tspfmitype) :: this
795  ! -- local
796  integer(I4B) :: ngwfpack
797  integer(I4B) :: ngwfterms
798  integer(I4B) :: ip
799  integer(I4B) :: imover
800  integer(I4B) :: ntomvr
801  integer(I4B) :: iterm
802  character(len=LENPACKAGENAME) :: budtxt
803  class(bndtype), pointer :: packobj => null()
804  !
805  ! -- determine size of gwf terms
806  ngwfpack = this%gwfbndlist%Count()
807  !
808  ! -- Count number of to-mvr terms, but do not include advanced packages
809  ! as those mover terms are not losses from the cell, but rather flows
810  ! within the advanced package
811  ntomvr = 0
812  do ip = 1, ngwfpack
813  packobj => getbndfromlist(this%gwfbndlist, ip)
814  imover = packobj%imover
815  if (packobj%isadvpak /= 0) imover = 0
816  if (imover /= 0) then
817  ntomvr = ntomvr + 1
818  end if
819  end do
820  !
821  ! -- Allocate arrays in fmi of size ngwfterms, which is the number of
822  ! packages plus the number of packages with mover terms.
823  ngwfterms = ngwfpack + ntomvr
824  call this%allocate_gwfpackages(ngwfterms)
825  !
826  ! -- Assign values in the fmi package
827  iterm = 1
828  do ip = 1, ngwfpack
829  !
830  ! -- set and store names
831  packobj => getbndfromlist(this%gwfbndlist, ip)
832  budtxt = adjustl(packobj%text)
833  call this%gwfpackages(iterm)%set_name(packobj%packName, budtxt)
834  this%flowpacknamearray(iterm) = packobj%packName
835  call this%gwfpackages(iterm)%set_auxname(packobj%naux, &
836  packobj%auxname)
837  iterm = iterm + 1
838  !
839  ! -- if this package has a mover associated with it, then add another
840  ! term that corresponds to the mover flows
841  imover = packobj%imover
842  if (packobj%isadvpak /= 0) imover = 0
843  if (imover /= 0) then
844  budtxt = trim(adjustl(packobj%text))//'-TO-MVR'
845  call this%gwfpackages(iterm)%set_name(packobj%packName, budtxt)
846  this%flowpacknamearray(iterm) = packobj%packName
847  call this%gwfpackages(iterm)%set_auxname(packobj%naux, &
848  packobj%auxname)
849  this%igwfmvrterm(iterm) = 1
850  iterm = iterm + 1
851  end if
852  end do
854 
855  !> @brief Initialize an array for storing PackageBudget objects.
856  !!
857  !! This routine allocates gwfpackages (an array of PackageBudget
858  !! objects) to the proper size and initializes member variables.
859  !<
860  subroutine gwtfmi_allocate_gwfpackages(this, ngwfterms)
861  ! -- modules
862  use constantsmodule, only: lenmempath
864  ! -- dummy
865  class(tspfmitype) :: this
866  integer(I4B), intent(in) :: ngwfterms
867  ! -- local
868  integer(I4B) :: n
869  character(len=LENMEMPATH) :: memPath
870  !
871  ! -- direct allocate
872  allocate (this%gwfpackages(ngwfterms))
873  allocate (this%flowpacknamearray(ngwfterms))
874  allocate (this%datp(ngwfterms))
875  !
876  ! -- mem_allocate
877  call mem_allocate(this%iatp, ngwfterms, 'IATP', this%memoryPath)
878  call mem_allocate(this%igwfmvrterm, ngwfterms, 'IGWFMVRTERM', this%memoryPath)
879  !
880  ! -- initialize
881  this%nflowpack = ngwfterms
882  do n = 1, this%nflowpack
883  this%iatp(n) = 0
884  this%igwfmvrterm(n) = 0
885  this%flowpacknamearray(n) = ''
886  !
887  ! -- Create a mempath for each individual flow package data set
888  ! of the form, MODELNAME/FMI-FTn
889  write (mempath, '(a, i0)') trim(this%memoryPath)//'-FT', n
890  call this%gwfpackages(n)%initialize(mempath)
891  end do
892  end subroutine gwtfmi_allocate_gwfpackages
893 
894  !> @brief Deallocate memory
895  !!
896  !! Deallocate memory that stores the gwfpackages array
897  !<
899  ! -- modules
900  ! -- dummy
901  class(tspfmitype) :: this
902  ! -- local
903  integer(I4B) :: n
904  !
905  ! -- initialize
906  do n = 1, this%nflowpack
907  call this%gwfpackages(n)%da()
908  end do
909  end subroutine gwtfmi_deallocate_gwfpackages
910 
911 end module tspfmimodule
This module contains the base boundary package.
class(bndtype) function, pointer, public getbndfromlist(list, idx)
Get boundary from package list.
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
subroutine, public budgetobject_cr_bfr(this, name, ibinun, iout, colconv1, colconv2)
Create a new budget object from a binary flow file.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
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 lenvarname
maximum length of a variable name
Definition: Constants.f90:17
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dem6
real constant 1e-6
Definition: Constants.f90:109
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
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 defines variable data types.
Definition: kind.f90:8
subroutine, public memorystore_release(varname, memory_path)
Release a single variable from the memory store.
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
This module contains the PackageBudgetModule Module.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:27
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:26
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:32
real(dp) function gwfsatold(this, n, delt)
Calculate the previous saturation level.
Definition: tsp-fmi.f90:504
integer(i4b), parameter nbditems
Definition: tsp-fmi.f90:24
subroutine fmi_bd(this, isuppress_output, model_budget)
Calculate budget terms associated with FMI object.
Definition: tsp-fmi.f90:251
subroutine gwtfmi_deallocate_gwfpackages(this)
Deallocate memory.
Definition: tsp-fmi.f90:899
character(len=lenbudtxt), dimension(nbditems) budtxt
Definition: tsp-fmi.f90:25
subroutine gwtfmi_allocate_scalars(this)
@ brief Allocate scalars
Definition: tsp-fmi.f90:366
subroutine gwtfmi_allocate_arrays(this, nodes)
@ brief Allocate arrays for FMI object
Definition: tsp-fmi.f90:391
subroutine, public fmi_cr(fmiobj, name_model, input_mempath, inunit, iout, eqnsclfac, depvartype)
Create a new FMI object.
Definition: tsp-fmi.f90:76
subroutine gwtfmi_source_options(this)
@ brief Source input options for package
Definition: tsp-fmi.f90:528
subroutine gwtfmi_source_packagedata(this)
@ brief Source input options for package
Definition: tsp-fmi.f90:559
subroutine initialize_gwfterms_from_gwfbndlist(this)
Initialize groundwater flow terms from the groundwater budget.
Definition: tsp-fmi.f90:791
subroutine gwtfmi_allocate_gwfpackages(this, ngwfterms)
Initialize an array for storing PackageBudget objects.
Definition: tsp-fmi.f90:861
subroutine fmi_fc(this, nodes, cold, nja, matrix_sln, idxglo, rhs)
Calculate coefficients and fill matrix and rhs terms associated with FMI object.
Definition: tsp-fmi.f90:189
subroutine initialize_gwfterms_from_bfr(this)
Initialize the groundwater flow terms based on the budget file reader.
Definition: tsp-fmi.f90:692
subroutine fmi_rp(this, inmvr)
Read and prepare.
Definition: tsp-fmi.f90:109
subroutine set_aptbudobj_pointer(this, name, budobjptr)
Set the pointer to a budget object.
Definition: tsp-fmi.f90:668
subroutine set_active_status(this, cnew)
Set gwt transport cell status.
Definition: tsp-fmi.f90:421
subroutine fmi_ot_flow(this, icbcfl, icbcun)
Save budget terms associated with FMI object.
Definition: tsp-fmi.f90:272
character(len=lenpackagename) text
Definition: tsp-fmi.f90:22
subroutine fmi_ad(this, cnew)
Advance routine for FMI object.
Definition: tsp-fmi.f90:138
subroutine fmi_cq(this, cnew, flowja)
Calculate flow correction.
Definition: tsp-fmi.f90:222
subroutine gwtfmi_da(this)
Deallocate variables.
Definition: tsp-fmi.f90:312
@ brief BndType
Derived type for the Budget object.
Definition: Budget.f90:39
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
A generic heterogeneous doubly-linked list.
Definition: List.f90:14
Derived type for storing flows.