MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
swf-flw.f90
Go to the documentation of this file.
1 !> @brief This module contains the FLW package methods
2 !!
3 !! This module can be used to represent inflow to streams. It is
4 !! designed similarly to the GWF WEL package.
5 !!
6 !<
8  ! -- modules
9  use kindmodule, only: dp, i4b
10  use constantsmodule, only: dzero, dem1, done, lenftype, dnodata, &
12  use simvariablesmodule, only: errmsg
15  use bndmodule, only: bndtype
16  use bndextmodule, only: bndexttype
19  use observemodule, only: observetype
25  !
26  implicit none
27  !
28  private
29  public :: flw_create
30  !
31  character(len=LENFTYPE) :: ftype = 'FLW' !< package ftype
32  character(len=16) :: text = ' FLW' !< package flow text string
33  !
34  type, extends(bndexttype) :: swfflwtype
35  real(dp), dimension(:), pointer, contiguous :: q => null() !< volumetric rate
36  contains
37  procedure :: allocate_scalars => flw_allocate_scalars
38  procedure :: allocate_arrays => flw_allocate_arrays
39  procedure :: source_options => flw_options
40  procedure :: log_flw_options
41  procedure :: bnd_rp => flw_rp
42  procedure :: bnd_cf => flw_cf
43  procedure :: bnd_fc => flw_fc
44  procedure :: bnd_da => flw_da
45  procedure :: define_listlabel
46  procedure :: bound_value => flw_bound_value
47  procedure :: q_mult
48  ! -- methods for observations
49  procedure, public :: bnd_obs_supported => flw_obs_supported
50  procedure, public :: bnd_df_obs => flw_df_obs
51  procedure, public :: bnd_bd_obs => flw_bd_obs
52  ! -- methods for time series
53  procedure, public :: bnd_rp_ts => flw_rp_ts
54  end type swfflwtype
55 
56 contains
57 
58  !> @ brief Create a new package object
59  !!
60  !! Create a new FLW Package object
61  !!
62  !<
63  subroutine flw_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
64  mempath)
65  ! -- dummy variables
66  class(bndtype), pointer :: packobj !< pointer to default package type
67  integer(I4B), intent(in) :: id !< package id
68  integer(I4B), intent(in) :: ibcnum !< boundary condition number
69  integer(I4B), intent(in) :: inunit !< unit number of FLW package input file
70  integer(I4B), intent(in) :: iout !< unit number of model listing file
71  character(len=*), intent(in) :: namemodel !< model name
72  character(len=*), intent(in) :: pakname !< package name
73  character(len=*), intent(in) :: mempath !< input mempath
74  ! -- local variables
75  type(swfflwtype), pointer :: flwobj
76  !
77  ! -- allocate the object and assign values to object variables
78  allocate (flwobj)
79  packobj => flwobj
80  !
81  ! -- create name and memory path
82  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
83  packobj%text = text
84  !
85  ! -- allocate scalars
86  call flwobj%allocate_scalars()
87  !
88  ! -- initialize package
89  call packobj%pack_initialize()
90 
91  packobj%inunit = inunit
92  packobj%iout = iout
93  packobj%id = id
94  packobj%ibcnum = ibcnum
95  packobj%ncolbnd = 1
96  packobj%iscloc = 1
97  packobj%ictMemPath = ''
98  !
99  ! -- return
100  return
101  end subroutine flw_create
102 
103  !> @ brief Allocate scalars
104  !!
105  !! Allocate and initialize scalars for the FLW package. The base model
106  !! allocate scalars method is also called.
107  !!
108  !<
109  subroutine flw_allocate_scalars(this)
110  ! -- modules
112  ! -- dummy variables
113  class(swfflwtype) :: this !< SwfFlwType object
114  !
115  ! -- call base type allocate scalars
116  call this%BndExtType%allocate_scalars()
117  !
118  ! -- allocate the object and assign values to object variables
119  ! none for this package
120  !
121  ! -- Set values
122  ! none for this package
123  !
124  ! -- return
125  return
126  end subroutine flw_allocate_scalars
127 
128  !> @ brief Allocate arrays
129  !!
130  !! Allocate and initialize arrays for the SWF package
131  !!
132  !<
133  subroutine flw_allocate_arrays(this, nodelist, auxvar)
134  ! -- modules
136  ! -- dummy
137  class(swfflwtype) :: this
138  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
139  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
140  ! -- local
141  !
142  ! -- call BndExtType allocate scalars
143  call this%BndExtType%allocate_arrays(nodelist, auxvar)
144  !
145  ! -- set q array input context pointer
146  call mem_setptr(this%q, 'Q', this%input_mempath)
147  !
148  ! -- checkin q array input context pointer
149  call mem_checkin(this%q, 'Q', this%memoryPath, &
150  'Q', this%input_mempath)
151  !
152  ! -- return
153  return
154  end subroutine flw_allocate_arrays
155 
156  !> @ brief Deallocate package memory
157  !!
158  !! Deallocate SWF package scalars and arrays.
159  !!
160  !<
161  subroutine flw_da(this)
162  ! -- modules
164  ! -- dummy variables
165  class(swfflwtype) :: this !< SwfFlwType object
166  !
167  ! -- Deallocate parent package
168  call this%BndExtType%bnd_da()
169  !
170  ! -- scalars
171  call mem_deallocate(this%q, 'Q', this%memoryPath)
172  !
173  ! -- return
174  return
175  end subroutine flw_da
176 
177  !> @ brief Source additional options for package
178  !!
179  !! Source additional options for SWF package.
180  !!
181  !<
182  subroutine flw_options(this)
183  ! -- modules
184  use inputoutputmodule, only: urword
187  ! -- dummy variables
188  class(swfflwtype), intent(inout) :: this !< SwfFlwType object
189  ! -- local variables
190  type(swfflwparamfoundtype) :: found
191  ! -- formats
192  !
193  ! -- source base BndExtType options
194  call this%BndExtType%source_options()
195  !
196  ! -- source options from input context
197  ! none
198  !
199  ! -- log SWF specific options
200  call this%log_flw_options(found)
201  !
202  ! -- return
203  return
204  end subroutine flw_options
205 
206  !> @ brief Log SWF specific package options
207  !<
208  subroutine log_flw_options(this, found)
209  ! -- modules
211  ! -- dummy variables
212  class(swfflwtype), intent(inout) :: this !< BndExtType object
213  type(swfflwparamfoundtype), intent(in) :: found
214  ! -- local variables
215  ! -- format
216  !
217  ! -- log found options
218  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
219  //' OPTIONS'
220  !
221  ! if (found%mover) then
222  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
223  ! end if
224  !
225  ! -- close logging block
226  write (this%iout, '(1x,a)') &
227  'END OF '//trim(adjustl(this%text))//' OPTIONS'
228  !
229  ! -- return
230  return
231  end subroutine log_flw_options
232 
233  !> @ brief SWF read and prepare
234  !!
235  !<
236  subroutine flw_rp(this)
237  use tdismodule, only: kper
238  ! -- dummy
239  class(swfflwtype), intent(inout) :: this
240  ! -- local
241  !
242  if (this%iper /= kper) return
243  !
244  ! -- Call the parent class read and prepare
245  call this%BndExtType%bnd_rp()
246  !
247  ! -- Write the list to iout if requested
248  if (this%iprpak /= 0) then
249  call this%write_list()
250  end if
251  !
252  ! -- return
253  return
254  end subroutine flw_rp
255 
256  !> @ brief Formulate the package hcof and rhs terms.
257  !!
258  !! Formulate the hcof and rhs terms for the FLW package that will be
259  !! added to the coefficient matrix and right-hand side vector.
260  !!
261  !<
262  subroutine flw_cf(this)
263  ! -- dummy variables
264  class(swfflwtype) :: this !< SwfFlwType object
265  ! -- local variables
266  integer(I4B) :: i, node
267  real(DP) :: q
268  !
269  ! -- Return if no inflows
270  if (this%nbound == 0) return
271  !
272  ! -- Calculate hcof and rhs for each flw entry
273  do i = 1, this%nbound
274  node = this%nodelist(i)
275  this%hcof(i) = dzero
276  if (this%ibound(node) <= 0) then
277  this%rhs(i) = dzero
278  cycle
279  end if
280  q = this%q_mult(i)
281  this%rhs(i) = -q
282  end do
283  !
284  return
285  end subroutine flw_cf
286 
287  !> @ brief Copy hcof and rhs terms into solution.
288  !!
289  !! Add the hcof and rhs terms for the FLW package to the
290  !! coefficient matrix and right-hand side vector.
291  !!
292  !<
293  subroutine flw_fc(this, rhs, ia, idxglo, matrix_sln)
294  ! -- dummy variables
295  class(swfflwtype) :: this !< SwfFlwType object
296  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
297  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
298  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
299  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
300  ! -- local variables
301  integer(I4B) :: i
302  integer(I4B) :: n
303  integer(I4B) :: ipos
304  !
305  ! -- pakmvrobj fc
306  if (this%imover == 1) then
307  call this%pakmvrobj%fc()
308  end if
309  !
310  ! -- Copy package rhs and hcof into solution rhs and amat
311  do i = 1, this%nbound
312  n = this%nodelist(i)
313  rhs(n) = rhs(n) + this%rhs(i)
314  ipos = ia(n)
315  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
316  !
317  ! -- If mover is active and this flw item is discharging,
318  ! store available water (as positive value).
319  if (this%imover == 1 .and. this%rhs(i) > dzero) then
320  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
321  end if
322  end do
323  !
324  ! -- return
325  return
326  end subroutine flw_fc
327 
328  !> @ brief Define the list label for the package
329  !!
330  !! Method defined the list label for the FLW package. The list label is
331  !! the heading that is written to iout when PRINT_INPUT option is used.
332  !!
333  !<
334  subroutine define_listlabel(this)
335  ! -- dummy variables
336  class(swfflwtype), intent(inout) :: this !< SwfFlwType object
337  !
338  ! -- create the header list label
339  this%listlabel = trim(this%filtyp)//' NO.'
340  if (this%dis%ndim == 3) then
341  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
342  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
343  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
344  elseif (this%dis%ndim == 2) then
345  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
346  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
347  else
348  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
349  end if
350  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
351  if (this%inamedbound == 1) then
352  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
353  end if
354  !
355  ! -- return
356  return
357  end subroutine define_listlabel
358 
359  ! -- Procedures related to observations
360 
361  !> @brief Determine if observations are supported.
362  !!
363  !! Function to determine if observations are supported by the FLW package.
364  !! Observations are supported by the FLW package.
365  !!
366  !! @return flw_obs_supported boolean indicating if observations are supported
367  !!
368  !<
369  logical function flw_obs_supported(this)
370  ! -- dummy variables
371  class(swfflwtype) :: this !< SwfFlwType object
372  !
373  ! -- set boolean
374  flw_obs_supported = .true.
375  !
376  ! -- return
377  return
378  end function flw_obs_supported
379 
380  !> @brief Define the observation types available in the package
381  !!
382  !! Method to define the observation types available in the FLW package.
383  !!
384  !<
385  subroutine flw_df_obs(this)
386  ! -- dummy variables
387  class(swfflwtype) :: this !< SwfFlwType object
388  ! -- local variables
389  integer(I4B) :: indx
390  !
391  ! -- initialize observations
392  call this%obs%StoreObsType('flw', .true., indx)
393  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
394  !
395  ! -- Store obs type and assign procedure pointer
396  ! for to-mvr observation type.
397  call this%obs%StoreObsType('to-mvr', .true., indx)
398  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
399  !
400  ! -- return
401  return
402  end subroutine flw_df_obs
403 
404  !> @brief Save observations for the package
405  !!
406  !! Method to save simulated values for the FLW package.
407  !!
408  !<
409  subroutine flw_bd_obs(this)
410  ! -- dummy variables
411  class(swfflwtype) :: this !< SwfFlwType object
412  ! -- local variables
413  integer(I4B) :: i
414  integer(I4B) :: n
415  integer(I4B) :: jj
416  real(DP) :: v
417  type(observetype), pointer :: obsrv => null()
418  !
419  ! -- clear the observations
420  call this%obs%obs_bd_clear()
421  !
422  ! -- Save simulated values for all of package's observations.
423  do i = 1, this%obs%npakobs
424  obsrv => this%obs%pakobs(i)%obsrv
425  if (obsrv%BndFound) then
426  do n = 1, obsrv%indxbnds_count
427  v = dnodata
428  jj = obsrv%indxbnds(n)
429  select case (obsrv%ObsTypeId)
430  case ('TO-MVR')
431  if (this%imover == 1) then
432  v = this%pakmvrobj%get_qtomvr(jj)
433  if (v > dzero) then
434  v = -v
435  end if
436  end if
437  case ('FLW')
438  v = this%simvals(jj)
439  case default
440  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
441  call store_error(errmsg)
442  end select
443  call this%obs%SaveOneSimval(obsrv, v)
444  end do
445  else
446  call this%obs%SaveOneSimval(obsrv, dnodata)
447  end if
448  end do
449  !
450  ! -- return
451  return
452  end subroutine flw_bd_obs
453 
454  ! -- Procedure related to time series
455 
456  !> @brief Assign time series links for the package
457  !!
458  !! Assign the time series links for the FLW package. Only
459  !! the Q variable can be defined with time series.
460  !!
461  !<
462  subroutine flw_rp_ts(this)
463  ! -- dummy variables
464  class(swfflwtype), intent(inout) :: this !< SwfFlwType object
465  ! -- local variables
466  integer(I4B) :: i, nlinks
467  type(timeserieslinktype), pointer :: tslink => null()
468  !
469  ! -- set up the time series links
470  nlinks = this%TsManager%boundtslinks%Count()
471  do i = 1, nlinks
472  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
473  if (associated(tslink)) then
474  if (tslink%JCol == 1) then
475  tslink%Text = 'Q'
476  end if
477  end if
478  end do
479  !
480  ! -- return
481  return
482  end subroutine flw_rp_ts
483 
484  function q_mult(this, row) result(q)
485  ! -- modules
486  use constantsmodule, only: dzero
487  ! -- dummy variables
488  class(swfflwtype), intent(inout) :: this !< BndExtType object
489  integer(I4B), intent(in) :: row
490  ! -- result
491  real(dp) :: q
492  !
493  if (this%iauxmultcol > 0) then
494  q = this%q(row) * this%auxvar(this%iauxmultcol, row)
495  else
496  q = this%q(row)
497  end if
498  !
499  ! -- return
500  return
501  end function q_mult
502 
503  !> @ brief Return a bound value
504  !!
505  !! Return a bound value associated with an ncolbnd index
506  !! and row.
507  !!
508  !<
509  function flw_bound_value(this, col, row) result(bndval)
510  ! -- modules
511  use constantsmodule, only: dzero
512  ! -- dummy variables
513  class(swfflwtype), intent(inout) :: this !< BndExtType object
514  integer(I4B), intent(in) :: col
515  integer(I4B), intent(in) :: row
516  ! -- result
517  real(dp) :: bndval
518  !
519  select case (col)
520  case (1)
521  bndval = this%q_mult(row)
522  case default
523  errmsg = 'Programming error. FLW bound value requested column '&
524  &'outside range of ncolbnd (1).'
525  call store_error(errmsg)
526  call store_error_filename(this%input_fname)
527  end select
528  !
529  ! -- return
530  return
531  end function flw_bound_value
532 
533 end module swfflwmodule
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:94
real(dp), parameter dem1
real constant 1e-1
Definition: Constants.f90:102
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:38
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
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
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:249
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_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
real(dp) function sqsaturationderivative(top, bot, x, c1, c2)
@ brief sQSaturationDerivative
real(dp) function sqsaturation(top, bot, x, c1, c2)
@ brief sQSaturation
This module contains the FLW package methods.
Definition: swf-flw.f90:7
subroutine flw_df_obs(this)
Define the observation types available in the package.
Definition: swf-flw.f90:386
character(len=lenftype) ftype
package ftype
Definition: swf-flw.f90:31
subroutine flw_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: swf-flw.f90:263
subroutine, public flw_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
@ brief Create a new package object
Definition: swf-flw.f90:65
subroutine flw_da(this)
@ brief Deallocate package memory
Definition: swf-flw.f90:162
real(dp) function flw_bound_value(this, col, row)
@ brief Return a bound value
Definition: swf-flw.f90:510
logical function flw_obs_supported(this)
Determine if observations are supported.
Definition: swf-flw.f90:370
subroutine flw_allocate_scalars(this)
@ brief Allocate scalars
Definition: swf-flw.f90:110
character(len=16) text
package flow text string
Definition: swf-flw.f90:32
subroutine flw_rp(this)
@ brief SWF read and prepare
Definition: swf-flw.f90:237
subroutine flw_bd_obs(this)
Save observations for the package.
Definition: swf-flw.f90:410
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: swf-flw.f90:335
subroutine log_flw_options(this, found)
@ brief Log SWF specific package options
Definition: swf-flw.f90:209
subroutine flw_options(this)
@ brief Source additional options for package
Definition: swf-flw.f90:183
real(dp) function q_mult(this, row)
Definition: swf-flw.f90:485
subroutine flw_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: swf-flw.f90:134
subroutine flw_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: swf-flw.f90:294
subroutine flw_rp_ts(this)
Assign time series links for the package.
Definition: swf-flw.f90:463
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType