MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
swfstomodule Module Reference

This module contains the storage package methods. More...

Data Types

type  swfstotype
 

Functions/Subroutines

subroutine, public sto_cr (stoobj, name_model, mempath, inunit, iout, cxs)
 @ brief Create a new package object More...
 
subroutine sto_ar (this, dis, ibound)
 @ brief Allocate and read method for package More...
 
subroutine sto_rp (this)
 @ brief Read and prepare method for package More...
 
subroutine sto_ad (this)
 @ brief Advance the package More...
 
subroutine sto_fc (this, kiter, stage_old, stage_new, matrix_sln, idxglo, rhs)
 @ brief Fill A and right-hand side for the package More...
 
subroutine sto_fc_dis1d (this, kiter, stage_old, stage_new, matrix_sln, idxglo, rhs)
 @ brief Fill A and right-hand side for the package More...
 
subroutine sto_fc_dis2d (this, kiter, stage_old, stage_new, matrix_sln, idxglo, rhs)
 @ brief Fill A and right-hand side for the package More...
 
subroutine sto_cq (this, flowja, stage_new, stage_old)
 @ brief Calculate flows for package More...
 
subroutine calc_storage_dis1d (this, n, stage_new, stage_old, dx, qsto, derv)
 
subroutine calc_storage_dis2d (this, n, stage_new, stage_old, qsto, derv)
 
subroutine sto_bd (this, isuppress_output, model_budget)
 @ brief Model budget calculation for package More...
 
subroutine sto_save_model_flows (this, icbcfl, icbcun)
 @ brief Save model flows for package More...
 
subroutine sto_da (this)
 @ brief Deallocate package memory More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine allocate_arrays (this, nodes)
 @ brief Allocate package arrays More...
 
subroutine source_options (this)
 @ brief Source input options for package More...
 
subroutine log_options (this, found)
 @ brief Log found options for package More...
 
subroutine source_data (this)
 @ brief Source input data for package More...
 
subroutine set_dfw_pointers (this)
 Set pointers to channel properties in DFW Package. More...
 
real(dp) function, dimension(:), pointer reach_length_pointer (this)
 

Variables

character(len=lenbudtxt), dimension(1) budtxt = [' STORAGE']
 

Detailed Description

This module contains the methods used to add the effects of storage on the surface water flow equation.

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine swfstomodule::allocate_arrays ( class(swfstotype), target  this,
integer(i4b), intent(in)  nodes 
)

Allocate and initialize STO package arrays.

Parameters
thisSwfStoType object
[in]nodesactive model nodes

Definition at line 545 of file swf-sto.f90.

546  ! -- modules
548  ! -- dummy variables
549  class(SwfStoType), target :: this !< SwfStoType object
550  integer(I4B), intent(in) :: nodes !< active model nodes
551  ! -- local variables
552  integer(I4B) :: n
553  !
554  ! -- Allocate arrays
555  call mem_allocate(this%qsto, nodes, 'STRGSS', this%memoryPath)
556  !
557  ! -- set input context pointers
558  if (this%inunit > 0) then
559  call mem_setptr(this%iper, 'IPER', this%input_mempath)
560  call mem_setptr(this%storage, 'STORAGE', this%input_mempath)
561  end if
562  !
563  ! -- Initialize arrays
564  this%iss = 0
565  do n = 1, nodes
566  this%qsto(n) = dzero
567  end do

◆ allocate_scalars()

subroutine swfstomodule::allocate_scalars ( class(swfstotype this)

Allocate and initialize scalars for the STO package. The base numerical package allocate scalars method is also called.

Parameters
thisSwfStoType object

Definition at line 524 of file swf-sto.f90.

525  ! -- modules
527  ! -- dummy variables
528  class(SwfStoType) :: this !< SwfStoType object
529  !
530  ! -- allocate scalars in NumericalPackageType
531  call this%NumericalPackageType%allocate_scalars()
532  !
533  ! -- allocate scalars
534  !call mem_allocate(this%xxx, 'XXX', this%memoryPath)
535  !
536  ! -- initialize scalars
537  !this%xxx = 0

◆ calc_storage_dis1d()

subroutine swfstomodule::calc_storage_dis1d ( class(swfstotype this,
integer(i4b), intent(in)  n,
real(dp), intent(in)  stage_new,
real(dp), intent(in)  stage_old,
real(dp), intent(in)  dx,
real(dp), intent(inout)  qsto,
real(dp), intent(inout), optional  derv 
)

Definition at line 360 of file swf-sto.f90.

361  ! module
362  use tdismodule, only: delt
364  ! dummy
365  class(SwfStoType) :: this
366  integer(I4B), intent(in) :: n
367  real(DP), intent(in) :: stage_new
368  real(DP), intent(in) :: stage_old
369  real(DP), intent(in) :: dx
370  real(DP), intent(inout) :: qsto
371  real(DP), intent(inout), optional :: derv
372  ! local
373  real(DP) :: depth_new
374  real(DP) :: depth_old
375  real(DP) :: width_n
376  real(DP) :: width_m
377  real(DP) :: cxs_area_new
378  real(DP) :: cxs_area_old
379  real(DP) :: cxs_area_eps
380  real(DP) :: eps
381 
382  call this%dis%get_flow_width(n, n, 0, width_n, width_m)
383  depth_new = stage_new - this%dis%bot(n)
384  depth_old = stage_old - this%dis%bot(n)
385  cxs_area_new = this%cxs%get_area(this%idcxs(n), width_n, depth_new)
386  cxs_area_old = this%cxs%get_area(this%idcxs(n), width_n, depth_old)
387  qsto = (cxs_area_new - cxs_area_old) * dx / delt
388  if (present(derv)) then
389  eps = get_perturbation(depth_new)
390  cxs_area_eps = this%cxs%get_area(this%idcxs(n), width_n, depth_new + eps)
391  derv = (cxs_area_eps - cxs_area_new) * dx / delt / eps
392  end if
393 
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:372
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
Here is the call graph for this function:

◆ calc_storage_dis2d()

subroutine swfstomodule::calc_storage_dis2d ( class(swfstotype this,
integer(i4b), intent(in)  n,
real(dp), intent(in)  stage_new,
real(dp), intent(in)  stage_old,
real(dp), intent(inout)  qsto,
real(dp), intent(inout), optional  derv 
)

Definition at line 396 of file swf-sto.f90.

397  ! module
398  use tdismodule, only: delt
400  ! dummy
401  class(SwfStoType) :: this
402  integer(I4B), intent(in) :: n
403  real(DP), intent(in) :: stage_new
404  real(DP), intent(in) :: stage_old
405  real(DP), intent(inout) :: qsto
406  real(DP), intent(inout), optional :: derv
407  ! local
408  real(DP) :: area
409  real(DP) :: depth_new
410  real(DP) :: depth_old
411  real(DP) :: depth_eps
412  real(DP) :: volume_new
413  real(DP) :: volume_old
414  real(DP) :: eps
415 
416  area = this%dis%get_area(n)
417  depth_new = stage_new - this%dis%bot(n)
418  depth_old = stage_old - this%dis%bot(n)
419  volume_new = area * depth_new
420  volume_old = area * depth_old
421  qsto = (volume_new - volume_old) / delt
422 
423  if (present(derv)) then
424  eps = get_perturbation(depth_new)
425  depth_eps = depth_new + eps
426  derv = (depth_eps - depth_new) * area / delt / eps
427  end if
428 
Here is the call graph for this function:

◆ log_options()

subroutine swfstomodule::log_options ( class(swfstotype this,
type(swfstoparamfoundtype), intent(in)  found 
)

Log options block for STO package.

Parameters
thisSwfStoType object

Definition at line 601 of file swf-sto.f90.

602  ! -- modules
605  ! -- dummy variables
606  class(SwfStoType) :: this !< SwfStoType object
607  type(SwfStoParamFoundType), intent(in) :: found
608  ! -- local variables
609  ! -- formats
610  character(len=*), parameter :: fmtisvflow = &
611  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
612  &WHENEVER ICBCFL IS NOT ZERO.')"
613  !
614  write (this%iout, '(1x,a)') 'PROCESSING STORAGE OPTIONS'
615  !
616  if (found%ipakcb) then
617  write (this%iout, fmtisvflow)
618  end if
619  !
620  write (this%iout, '(1x,a)') 'END OF STORAGE OPTIONS'

◆ reach_length_pointer()

real(dp) function, dimension(:), pointer swfstomodule::reach_length_pointer ( class(swfstotype this)
Parameters
thisthis instance

Definition at line 666 of file swf-sto.f90.

667  ! dummy
668  class(SwfStoType) :: this !< this instance
669  ! return
670  real(DP), dimension(:), pointer :: ptr
671  ! local
672  class(DisBaseType), pointer :: dis
673 
674  ptr => null()
675  dis => this%dis
676  select type (dis)
677  type is (disv1dtype)
678  ptr => dis%length
679  end select
680 

◆ set_dfw_pointers()

subroutine swfstomodule::set_dfw_pointers ( class(swfstotype this)
Parameters
thisthis instance

Definition at line 653 of file swf-sto.f90.

654  ! -- modules
656  ! -- dummy
657  class(SwfStoType) :: this !< this instance
658  ! -- local
659  character(len=LENMEMPATH) :: dfw_mem_path
660 
661  dfw_mem_path = create_mem_path(this%name_model, 'DFW')
662  call mem_setptr(this%idcxs, 'IDCXS', dfw_mem_path)
663 
Here is the call graph for this function:

◆ source_data()

subroutine swfstomodule::source_data ( class(swfstotype this)

Source griddata block parameters for STO package.

Parameters
thisSwfStoType object

Definition at line 628 of file swf-sto.f90.

629  ! -- modules
632  ! -- dummy variables
633  class(SwfStotype) :: this !< SwfStoType object
634  ! -- local variables
635  character(len=24), dimension(1) :: aname
636  integer(I4B), dimension(:), pointer, contiguous :: map
637  !type(SwfStoParamFoundType) :: found
638  !
639  ! -- initialize data
640  data aname(1)/' XXX'/
641  !
642  ! -- set map to reduce data input arrays
643  map => null()
644  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
645  !
646  ! -- log griddata
647  write (this%iout, '(1x,a)') 'PROCESSING GRIDDATA'
648  write (this%iout, '(1x,a)') 'END PROCESSING GRIDDATA'

◆ source_options()

subroutine swfstomodule::source_options ( class(swfstotype this)

Source options block parameters for STO package.

Parameters
thisSwfStoType object

Definition at line 575 of file swf-sto.f90.

576  ! -- modules
580  ! -- dummy variables
581  class(SwfStoType) :: this !< SwfStoType object
582  ! -- local variables
583  type(SwfStoParamFoundType) :: found
584  !
585  ! -- source package input
586  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
587  !
588  if (found%ipakcb) then
589  this%ipakcb = -1
590  end if
591  !
592  ! -- log found options
593  call this%log_options(found)
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
Here is the call graph for this function:

◆ sto_ad()

subroutine swfstomodule::sto_ad ( class(swfstotype this)

Advance data in the STO package.

Parameters
thisSwfStoType object

Definition at line 189 of file swf-sto.f90.

190  ! -- modules
191  ! -- dummy variables
192  class(SwfStoType) :: this !< SwfStoType object

◆ sto_ar()

subroutine swfstomodule::sto_ar ( class(swfstotype this,
class(disbasetype), intent(in), pointer  dis,
integer(i4b), dimension(:), pointer, contiguous  ibound 
)

Method to allocate and read static data for the STO package.

Parameters
thisSwfStoType object
[in]dismodel discretization object
iboundmodel ibound array

Definition at line 102 of file swf-sto.f90.

103  ! -- modules
106  ! -- dummy variables
107  class(SwfStoType) :: this !< SwfStoType object
108  class(DisBaseType), pointer, intent(in) :: dis !< model discretization object
109  integer(I4B), dimension(:), pointer, contiguous :: ibound !< model ibound array
110  ! -- local variables
111  ! -- formats
112  character(len=*), parameter :: fmtsto = &
113  "(1x,/1x,'STO -- STORAGE PACKAGE, VERSION 1, 10/27/2023', &
114  &' INPUT READ FROM UNIT ', i0, //)"
115  !
116  ! --print a message identifying the storage package.
117  write (this%iout, fmtsto) this%inunit
118 
119  ! -- set pointers to data in dfw package
120  call this%set_dfw_pointers()
121 
122  !
123  ! -- store pointers to arguments that were passed in
124  this%dis => dis
125  this%ibound => ibound
126  !
127  ! -- set pointer to model iss
128  call mem_setptr(this%iss, 'ISS', create_mem_path(this%name_model))
129  !
130  ! -- Allocate arrays
131  call this%allocate_arrays(dis%nodes)
132  !
133  ! -- Read storage options
134  call this%source_options()
135  !
136  ! -- read the data block
137  ! no griddata at the moment for SWF Storage Package
138  ! call this%source_data()
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Here is the call graph for this function:

◆ sto_bd()

subroutine swfstomodule::sto_bd ( class(swfstotype this,
integer(i4b), intent(in)  isuppress_output,
type(budgettype), intent(inout)  model_budget 
)

Budget calculation for the STO package components. Components include specific storage and specific yield storage.

Parameters
thisSwfStoType object
[in]isuppress_outputflag to suppress model output
[in,out]model_budgetmodel budget object

Definition at line 437 of file swf-sto.f90.

438  ! -- modules
439  use tdismodule, only: delt
441  ! -- dummy variables
442  class(SwfStoType) :: this !< SwfStoType object
443  integer(I4B), intent(in) :: isuppress_output !< flag to suppress model output
444  type(BudgetType), intent(inout) :: model_budget !< model budget object
445  ! -- local variables
446  real(DP) :: rin
447  real(DP) :: rout
448  !
449  ! -- Add storage rates to model budget
450  call rate_accumulator(this%qsto, rin, rout)
451  call model_budget%addentry(rin, rout, delt, ' STO', &
452  isuppress_output, ' STORAGE')
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
Derived type for the Budget object.
Definition: Budget.f90:39
Here is the call graph for this function:

◆ sto_cq()

subroutine swfstomodule::sto_cq ( class(swfstotype this,
real(dp), dimension(:), intent(inout)  flowja,
real(dp), dimension(:), intent(inout)  stage_new,
real(dp), dimension(:), intent(inout)  stage_old 
)

Definition at line 320 of file swf-sto.f90.

321  ! -- dummy
322  class(SwfStoType) :: this
323  real(DP), intent(inout), dimension(:) :: flowja
324  real(DP), intent(inout), dimension(:) :: stage_new
325  real(DP), intent(inout), dimension(:) :: stage_old
326  ! -- local
327  real(DP), dimension(:), pointer :: reach_length
328  integer(I4B) :: n
329  integer(I4B) :: idiag
330  real(DP) :: dx
331  real(DP) :: q
332 
333  ! -- test if steady-state stress period
334  if (this%iss /= 0) return
335 
336  ! Set pointer to reach_length for 1d
337  reach_length => this%reach_length_pointer()
338 
339  ! -- Calculate storage term
340  do n = 1, this%dis%nodes
341  !
342  ! -- skip if constant stage
343  if (this%ibound(n) < 0) cycle
344 
345  ! Calculate storage for either the DIS1D or DIS2D cases and
346  ! add to flowja
347  if (associated(reach_length)) then
348  dx = reach_length(n)
349  call this%calc_storage_dis1d(n, stage_new(n), stage_old(n), dx, q)
350  else
351  call this%calc_storage_dis2d(n, stage_new(n), stage_old(n), q)
352  end if
353  this%qsto(n) = -q
354  idiag = this%dis%con%ia(n)
355  flowja(idiag) = flowja(idiag) + this%qsto(n)
356 
357  end do

◆ sto_cr()

subroutine, public swfstomodule::sto_cr ( type(swfstotype), pointer  stoobj,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
type(swfcxstype), intent(in), pointer  cxs 
)

Create a new storage (STO) object

Parameters
stoobjSwfStoType object
[in]name_modelname of model
[in]mempathinput context mem path
[in]inunitpackage input file unit
[in]ioutmodel listing file unit
[in]cxsthe pointer to the cxs package

Definition at line 70 of file swf-sto.f90.

71  ! -- dummy variables
72  type(SwfStoType), pointer :: stoobj !< SwfStoType object
73  character(len=*), intent(in) :: name_model !< name of model
74  character(len=*), intent(in) :: mempath !< input context mem path
75  integer(I4B), intent(in) :: inunit !< package input file unit
76  integer(I4B), intent(in) :: iout !< model listing file unit
77  type(SwfCxsType), pointer, intent(in) :: cxs !< the pointer to the cxs package
78  !
79  ! -- Create the object
80  allocate (stoobj)
81  !
82  ! -- create name and memory path
83  call stoobj%set_names(1, name_model, 'STO', 'STO', mempath)
84  !
85  ! -- Allocate scalars
86  call stoobj%allocate_scalars()
87  !
88  ! -- Set variables
89  stoobj%inunit = inunit
90  stoobj%iout = iout
91 
92  ! -- store pointers
93  stoobj%cxs => cxs
94 
Here is the caller graph for this function:

◆ sto_da()

subroutine swfstomodule::sto_da ( class(swfstotype this)

Deallocate STO package scalars and arrays.

Parameters
thisSwfStoType object

Definition at line 498 of file swf-sto.f90.

499  ! -- modules
501  ! -- dummy variables
502  class(SwfStoType) :: this !< SwfStoType object
503  !
504  ! -- Deallocate arrays if package is active
505  if (this%inunit > 0) then
506  call mem_deallocate(this%qsto)
507  nullify (this%idcxs)
508  nullify (this%iper)
509  nullify (this%storage)
510  end if
511  !
512  ! -- Deallocate scalars
513  !
514  ! -- deallocate parent
515  call this%NumericalPackageType%da()

◆ sto_fc()

subroutine swfstomodule::sto_fc ( class(swfstotype this,
integer(i4b)  kiter,
real(dp), dimension(:), intent(inout)  stage_old,
real(dp), dimension(:), intent(inout)  stage_new,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs 
)

Fill the coefficient matrix and right-hand side with the STO package terms.

Definition at line 200 of file swf-sto.f90.

201  ! -- modules
202  use tdismodule, only: delt
203  ! -- dummy
204  class(SwfStoType) :: this
205  integer(I4B) :: kiter
206  real(DP), intent(inout), dimension(:) :: stage_old
207  real(DP), intent(inout), dimension(:) :: stage_new
208  class(MatrixBaseType), pointer :: matrix_sln
209  integer(I4B), intent(in), dimension(:) :: idxglo
210  real(DP), intent(inout), dimension(:) :: rhs
211  ! -- local
212  character(len=LINELENGTH) :: distype = ''
213  ! -- formats
214  character(len=*), parameter :: fmtsperror = &
215  &"('Detected time step length of zero. SWF Storage Package cannot be ', &
216  &'used unless delt is non-zero.')"
217  !
218  ! -- test if steady-state stress period
219  if (this%iss /= 0) return
220  !
221  ! -- Ensure time step length is not zero
222  if (delt == dzero) then
223  write (errmsg, fmtsperror)
224  call store_error(errmsg, terminate=.true.)
225  end if
226 
227  call this%dis%get_dis_type(distype)
228  if (distype == 'DISV1D') then
229  call this%sto_fc_dis1d(kiter, stage_old, stage_new, matrix_sln, idxglo, rhs)
230  else
231  call this%sto_fc_dis2d(kiter, stage_old, stage_new, matrix_sln, idxglo, rhs)
232  end if
233 
Here is the call graph for this function:

◆ sto_fc_dis1d()

subroutine swfstomodule::sto_fc_dis1d ( class(swfstotype this,
integer(i4b)  kiter,
real(dp), dimension(:), intent(inout)  stage_old,
real(dp), dimension(:), intent(inout)  stage_new,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs 
)

Fill the coefficient matrix and right-hand side with the STO package terms.

Definition at line 241 of file swf-sto.f90.

243  ! -- modules
244  ! -- dummy
245  class(SwfStoType) :: this
246  integer(I4B) :: kiter
247  real(DP), intent(inout), dimension(:) :: stage_old
248  real(DP), intent(inout), dimension(:) :: stage_new
249  class(MatrixBaseType), pointer :: matrix_sln
250  integer(I4B), intent(in), dimension(:) :: idxglo
251  real(DP), intent(inout), dimension(:) :: rhs
252  ! -- local
253  integer(I4B) :: n, idiag
254  real(DP) :: derv
255  real(DP) :: qsto
256  real(DP), dimension(:), pointer :: reach_length
257 
258  ! Set pointer to reach_length for 1d
259  reach_length => this%reach_length_pointer()
260 
261  ! -- Calculate coefficients and put into amat
262  do n = 1, this%dis%nodes
263 
264  ! -- skip if constant stage
265  if (this%ibound(n) < 0) cycle
266 
267  call this%calc_storage_dis1d(n, stage_new(n), stage_old(n), &
268  reach_length(n), qsto, derv)
269 
270  ! -- Fill amat and rhs
271  idiag = this%dis%con%ia(n)
272  call matrix_sln%add_value_pos(idxglo(idiag), -derv)
273  rhs(n) = rhs(n) + qsto - derv * stage_new(n)
274 
275  end do

◆ sto_fc_dis2d()

subroutine swfstomodule::sto_fc_dis2d ( class(swfstotype this,
integer(i4b)  kiter,
real(dp), dimension(:), intent(inout)  stage_old,
real(dp), dimension(:), intent(inout)  stage_new,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs 
)

Fill the coefficient matrix and right-hand side with the STO package terms.

Definition at line 283 of file swf-sto.f90.

285  ! -- modules
286  ! -- dummy
287  class(SwfStoType) :: this
288  integer(I4B) :: kiter
289  real(DP), intent(inout), dimension(:) :: stage_old
290  real(DP), intent(inout), dimension(:) :: stage_new
291  class(MatrixBaseType), pointer :: matrix_sln
292  integer(I4B), intent(in), dimension(:) :: idxglo
293  real(DP), intent(inout), dimension(:) :: rhs
294  ! -- local
295  integer(I4B) :: n, idiag
296  real(DP) :: derv
297  real(DP) :: qsto
298 
299  ! -- Calculate coefficients and put into amat
300  do n = 1, this%dis%nodes
301  !
302  ! -- skip if constant stage
303  if (this%ibound(n) < 0) cycle
304 
305  ! Calculate storage and derivative term
306  call this%calc_storage_dis2d(n, stage_new(n), stage_old(n), &
307  qsto, derv)
308 
309  ! -- Fill amat and rhs
310  idiag = this%dis%con%ia(n)
311  call matrix_sln%add_value_pos(idxglo(idiag), -derv)
312  rhs(n) = rhs(n) + qsto - derv * stage_new(n)
313 
314  end do
315 

◆ sto_rp()

subroutine swfstomodule::sto_rp ( class(swfstotype this)

Method to read and prepare stress period data for the STO package.

Parameters
thisSwfStoType object

Definition at line 146 of file swf-sto.f90.

147  ! -- modules
148  use tdismodule, only: kper
149  implicit none
150  ! -- dummy variables
151  class(SwfStoType) :: this !< SwfStoType object
152  ! -- local variables
153  character(len=16) :: css(0:1)
154  ! -- data
155  data css(0)/' TRANSIENT'/
156  data css(1)/' STEADY-STATE'/
157  !
158  ! -- confirm package is active
159  if (this%inunit <= 0) return
160  !
161  ! -- confirm loaded iper
162  if (this%iper /= kper) return
163  !
164  write (this%iout, '(//,1x,a)') 'PROCESSING STORAGE PERIOD DATA'
165  !
166  ! -- set period iss
167  if (this%storage == 'STEADY-STATE') then
168  this%iss = 1
169  else if (this%storage == 'TRANSIENT') then
170  this%iss = 0
171  else
172  write (errmsg, '(a,a)') 'Unknown STORAGE data tag: ', &
173  trim(this%storage)
174  call store_error(errmsg)
175  call store_error_filename(this%input_fname)
176  end if
177  !
178  write (this%iout, '(1x,a)') 'END PROCESSING STORAGE PERIOD DATA'
179  !
180  write (this%iout, '(//1X,A,I0,A,A,/)') &
181  'STRESS PERIOD ', kper, ' IS ', trim(adjustl(css(this%iss)))
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ sto_save_model_flows()

subroutine swfstomodule::sto_save_model_flows ( class(swfstotype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  icbcun 
)

Save cell-by-cell budget terms for the STO package.

Parameters
thisSwfStoType object
[in]icbcflflag to output budget data
[in]icbcuncell-by-cell file unit number

Definition at line 460 of file swf-sto.f90.

461  ! -- dummy variables
462  class(SwfStoType) :: this !< SwfStoType object
463  integer(I4B), intent(in) :: icbcfl !< flag to output budget data
464  integer(I4B), intent(in) :: icbcun !< cell-by-cell file unit number
465  ! -- local variables
466  integer(I4B) :: ibinun
467  integer(I4B) :: iprint, nvaluesp, nwidthp
468  character(len=1) :: cdatafmp = ' ', editdesc = ' '
469  real(DP) :: dinact
470  !
471  ! -- Set unit number for binary output
472  if (this%ipakcb < 0) then
473  ibinun = icbcun
474  elseif (this%ipakcb == 0) then
475  ibinun = 0
476  else
477  ibinun = this%ipakcb
478  end if
479  if (icbcfl == 0) ibinun = 0
480  !
481  ! -- Record the storage rates if requested
482  if (ibinun /= 0) then
483  iprint = 0
484  dinact = dzero
485  !
486  ! -- qsto
487  call this%dis%record_array(this%qsto, this%iout, iprint, -ibinun, &
488  budtxt(1), cdatafmp, nvaluesp, &
489  nwidthp, editdesc, dinact)
490  end if

Variable Documentation

◆ budtxt

character(len=lenbudtxt), dimension(1) swfstomodule::budtxt = [' STORAGE']

Definition at line 23 of file swf-sto.f90.

23  character(len=LENBUDTXT), dimension(1) :: budtxt = & !< text labels for budget terms
24  &[' STORAGE']