MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwtspcmodule Module Reference

This module contains the GwtSpc Module. More...

Data Types

type  gwtspctype
 Derived type for managing SPC input. More...
 

Functions/Subroutines

subroutine initialize (this, dis, id, inunit, iout, name_model, packNameFlow)
 @ brief Initialize the SPC type More...
 
subroutine allocate_scalars (this)
 @ brief Allocate package scalars More...
 
subroutine read_options (this)
 @ brief Read options for package More...
 
subroutine read_dimensions (this)
 @ brief Read dimensions for package More...
 
subroutine allocate_arrays (this)
 @ brief Allocate package arrays More...
 
real(dp) function get_value (this, ientry, nbound_flow)
 @ brief Get the data value from this package More...
 
subroutine spc_rp (this)
 @ brief Read and prepare More...
 
subroutine spc_rp_list (this)
 @ brief spc_rp_list More...
 
subroutine spc_rp_array (this, line)
 @ brief spc_rp_array More...
 
subroutine spc_ad (this, nbound_flowpack, budtxt)
 @ brief Advance More...
 
subroutine spc_da (this)
 @ brief Deallocate variables More...
 
subroutine read_check_ionper (this)
 @ brief Check ionper More...
 
subroutine set_value (this, ival)
 @ brief Set the data value from the input file More...
 
subroutine check_flow_package (this, nbound_flowpack, budtxt)
 @ brief check_flow_package More...
 

Variables

character(len=lenftype) ftype = 'SPC'
 
character(len=lenpackagename) text = 'STRESS PACK CONC'
 

Detailed Description

This module contains the code for reading and storing a generic input file of source and sink concentrations.

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwtspcmodule::allocate_arrays ( class(gwtspctype this)
private

Allocate and initialize package arrays.

Parameters
thisGwtSpcType object

Definition at line 321 of file GwtSpc.f90.

322  ! -- modules
324  ! -- dummy variables
325  class(GwtSpcType) :: this !< GwtSpcType object
326  ! -- local
327  integer(I4B) :: i
328  !
329  ! -- allocate array
330  call mem_allocate(this%dblvec, this%maxbound, 'DBLVEC', this%memoryPath)
331  !
332  ! -- initialize dblvec to zero
333  do i = 1, this%maxbound
334  this%dblvec(i) = dzero
335  end do
336  !
337  ! -- return
338  return

◆ allocate_scalars()

subroutine gwtspcmodule::allocate_scalars ( class(gwtspctype this)
private

Allocate and initialize package scalars.

Parameters
thisGwtSpcType object

Definition at line 148 of file GwtSpc.f90.

149  ! -- modules
151  ! -- dummy variables
152  class(GwtSpcType) :: this !< GwtSpcType object
153  !
154  ! -- allocate scalars in memory manager
155  call mem_allocate(this%id, 'ID', this%memoryPath)
156  call mem_allocate(this%inunit, 'INUNIT', this%memoryPath)
157  call mem_allocate(this%iout, 'IOUT', this%memoryPath)
158  call mem_allocate(this%maxbound, 'MAXBOUND', this%memoryPath)
159  call mem_allocate(this%ionper, 'IONPER', this%memoryPath)
160  call mem_allocate(this%lastonper, 'LASTONPER', this%memoryPath)
161  call mem_allocate(this%iprpak, 'IPRPAK', this%memoryPath)
162  call mem_allocate(this%readasarrays, 'READASARRAYS', this%memoryPath)
163  !
164  ! -- allocate special derived types
165  allocate (this%TsManager)
166  allocate (this%TasManager)
167  !
168  ! -- initialize
169  this%id = 0
170  this%inunit = 0
171  this%iout = 0
172  this%maxbound = 0
173  this%ionper = 0
174  this%lastonper = 0
175  this%iprpak = 0
176  this%readasarrays = .false.
177  !
178  ! -- return
179  return

◆ check_flow_package()

subroutine gwtspcmodule::check_flow_package ( class(gwtspctype), intent(inout)  this,
integer(i4b), intent(in)  nbound_flowpack,
character(len=*), intent(in)  budtxt 
)

Check to make sure that flow package information is consistent with this SPC information.

Parameters
[in,out]thisGwtSpcType object

Definition at line 745 of file GwtSpc.f90.

746  ! -- modules
747  ! -- dummy
748  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
749  integer(I4B), intent(in) :: nbound_flowpack
750  character(len=*), intent(in) :: budtxt
751  ! -- local
752  !
753  ! -- Check and make sure MAXBOUND is not less than nbound_flowpack
754  if (this%maxbound < nbound_flowpack) then
755  write (errmsg, '(a, a, a, i0, a, i0, a)') &
756  'The SPC Package corresponding to flow package ', &
757  trim(this%packNameFlow), &
758  ' has MAXBOUND set less than the number of boundaries &
759  &active in this package. Found MAXBOUND equal ', &
760  this%maxbound, &
761  ' and number of flow boundaries (NBOUND) equal ', &
762  nbound_flowpack, &
763  '. Increase MAXBOUND in the SPC input file for this package.'
764  call store_error(errmsg)
765  call this%parser%StoreErrorUnit()
766  end if
767  !
768  ! -- If budtxt is RCHA or EVTA, then readasarrays must be used, otherwise
769  ! readasarrays cannot be used
770  select case (trim(adjustl(budtxt)))
771  case ('RCHA')
772  if (.not. this%readasarrays) then
773  write (errmsg, '(a, a, a)') &
774  'Array-based recharge must be used with array-based stress package &
775  &concentrations. GWF Package ', trim(this%packNameFlow), ' is being &
776  &used with list-based SPC6 input. Use array-based SPC6 input instead.'
777  call store_error(errmsg)
778  call this%parser%StoreErrorUnit()
779  end if
780  case ('EVTA')
781  if (.not. this%readasarrays) then
782  write (errmsg, '(a, a, a)') &
783  'Array-based evapotranspiration must be used with array-based stress &
784  &package concentrations. GWF Package ', trim(this%packNameFlow), &
785  &' is being used with list-based SPC6 input. Use array-based SPC6 &
786  &input instead.'
787  call store_error(errmsg)
788  call this%parser%StoreErrorUnit()
789  end if
790  case default
791  if (this%readasarrays) then
792  write (errmsg, '(a, a, a)') &
793  'List-based packages must be used with list-based stress &
794  &package concentrations. GWF Package ', trim(this%packNameFlow), &
795  &' is being used with array-based SPC6 input. Use list-based SPC6 &
796  &input instead.'
797  call store_error(errmsg)
798  call this%parser%StoreErrorUnit()
799  end if
800  end select
801  !
802  ! -- return
803  return
Here is the call graph for this function:

◆ get_value()

real(dp) function gwtspcmodule::get_value ( class(gwtspctype this,
integer(i4b), intent(in)  ientry,
integer(i4b), intent(in)  nbound_flow 
)

Get the floating point value from the dblvec array.

Parameters
thisGwtSpcType object
[in]ientryindex of the data to return
[in]nbound_flowsize of bound list in flow package

Definition at line 346 of file GwtSpc.f90.

347  class(GwtSpcType) :: this !< GwtSpcType object
348  integer(I4B), intent(in) :: ientry !< index of the data to return
349  integer(I4B), intent(in) :: nbound_flow !< size of bound list in flow package
350  real(DP) :: value
351  integer(I4B) :: nu
352  if (this%readasarrays) then
353  ! Special handling for reduced grids and readasarrays
354  ! if flow and transport are in the same simulation, then
355  ! ientry is a user node number and it corresponds to the
356  ! correct position in the dblvec array. But if flow and
357  ! transport are not in the same simulation, then ientry is
358  ! a reduced node number, because the list of flows in the
359  ! budget file do not include idomain < 1 entries. In this
360  ! case, ientry must be converted to a user node number so
361  ! that it corresponds to a user array, which includes
362  ! idomain < 1 values.
363  if (nbound_flow == this%maxbound) then
364  ! flow and transport are in the same simulation or there
365  ! are no idomain < 1 cells.
366  value = this%dblvec(ientry)
367  else
368  ! This identifies case where flow and transport must be
369  ! in a separate simulation, because nbound_flow is not
370  ! the same as this%maxbound. Under these conditions, we
371  ! must assume that ientry corresponds to a flow list that
372  ! would be of size ncpl if flow and transport were in the
373  ! same simulation, but because boundary cells with
374  ! idomain < 1 are not written to binary budget file, the
375  ! list size is smaller.
376  nu = this%dis%get_nodeuser(ientry)
377  value = this%dblvec(nu)
378  end if
379  else
380  value = this%dblvec(ientry)
381  end if
382  return

◆ initialize()

subroutine gwtspcmodule::initialize ( class(gwtspctype this,
class(disbasetype), intent(in), pointer  dis,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  packNameFlow 
)
private

Initialize the SPC object by setting up the parser, and time series manager, reading options and dimensions, and allocating memory.

Parameters
thisGwtSpcType
[in]disdiscretization package
[in]idid number for this spc package
[in]inunitunit number for input
[in]ioutunit number for output
[in]name_modelcharacter string containing model name
[in]packnameflowcharacter string containing name of corresponding flow package

Definition at line 87 of file GwtSpc.f90.

88  ! -- dummy variables
89  class(GwtSpcType) :: this !< GwtSpcType
90  class(DisBaseType), pointer, intent(in) :: dis !< discretization package
91  integer(I4B), intent(in) :: id !< id number for this spc package
92  integer(I4B), intent(in) :: inunit !< unit number for input
93  integer(I4B), intent(in) :: iout !< unit number for output
94  character(len=*), intent(in) :: name_model !< character string containing model name
95  character(len=*), intent(in) :: packNameflow !< character string containing name of corresponding flow package
96  ! -- local
97  !
98  ! -- construct the memory path
99  write (this%packName, '(a, i0)') 'SPC'//'-', id
100  this%name_model = name_model
101  this%memoryPath = create_mem_path(this%name_model, this%packName)
102  !
103  ! -- allocate scalar variables
104  call this%allocate_scalars()
105  !
106  ! -- assign member values
107  this%id = id
108  this%inunit = inunit
109  this%iout = iout
110  this%packNameFlow = packnameflow
111  !
112  ! -- set pointers
113  this%dis => dis
114  !
115  ! -- Setup the parser
116  call this%parser%Initialize(this%inunit, this%iout)
117  !
118  ! -- Setup the time series manager
119  call tsmanager_cr(this%TsManager, this%iout)
120  call tasmanager_cr(this%TasManager, dis, name_model, this%iout)
121  !
122  ! -- read options
123  call this%read_options()
124  !
125  ! -- read dimensions
126  if (this%readasarrays) then
127  this%maxbound = this%dis%get_ncpl()
128  else
129  call this%read_dimensions()
130  end if
131  !
132  ! -- allocate arrays
133  call this%allocate_arrays()
134  !
135  ! -- Now that time series are read, call define
136  call this%tsmanager%tsmanager_df()
137  call this%tasmanager%tasmanager_df()
138  !
139  ! -- return
140  return
Here is the call graph for this function:

◆ read_check_ionper()

subroutine gwtspcmodule::read_check_ionper ( class(gwtspctype), intent(inout)  this)

Generic method to read and check ionperiod, which is used to determine if new period data should be read from the input file. The check of ionperiod also makes sure periods are increasing in subsequent period data blocks. Copied from NumericalPackage

Parameters
[in,out]thisGwtSpcType object

Definition at line 683 of file GwtSpc.f90.

684  ! -- modules
685  use tdismodule, only: kper
686  ! -- dummy variables
687  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
688  !
689  ! -- save last value and read period number
690  this%lastonper = this%ionper
691  this%ionper = this%parser%GetInteger()
692  !
693  ! -- make check
694  if (this%ionper <= this%lastonper) then
695  write (errmsg, '(a, i0, a, i0, a, i0, a)') &
696  'Error in stress period ', kper, &
697  '. Period numbers not increasing. Found ', this%ionper, &
698  ' but last period block was assigned ', this%lastonper, '.'
699  call store_error(errmsg)
700  call this%parser%StoreErrorUnit()
701  end if
702  !
703  ! -- return
704  return
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ read_dimensions()

subroutine gwtspcmodule::read_dimensions ( class(gwtspctype), intent(inout)  this)
private

Read dimensions for this package.

Parameters
[in,out]thisGwtSpcType object

Definition at line 263 of file GwtSpc.f90.

264  ! -- dummy variables
265  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
266  ! -- local variables
267  character(len=LINELENGTH) :: keyword
268  logical(LGP) :: isfound
269  logical(LGP) :: endOfBlock
270  integer(I4B) :: ierr
271  !
272  ! -- get dimensions block
273  call this%parser%GetBlock('DIMENSIONS', isfound, ierr, &
274  supportopenclose=.true.)
275  !
276  ! -- parse dimensions block if detected
277  if (isfound) then
278  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(text))// &
279  ' DIMENSIONS'
280  do
281  call this%parser%GetNextLine(endofblock)
282  if (endofblock) exit
283  call this%parser%GetStringCaps(keyword)
284  select case (keyword)
285  case ('MAXBOUND')
286  this%maxbound = this%parser%GetInteger()
287  write (this%iout, '(4x,a,i7)') 'MAXBOUND = ', this%maxbound
288  case default
289  write (errmsg, '(a,3(1x,a))') &
290  'Unknown', trim(text), 'dimension:', trim(keyword)
291  call store_error(errmsg)
292  end select
293  end do
294  !
295  write (this%iout, '(1x,a)') 'END OF '//trim(adjustl(text))//' DIMENSIONS'
296  else
297  call store_error('Required DIMENSIONS block not found.')
298  call this%parser%StoreErrorUnit()
299  end if
300  !
301  ! -- verify dimensions were set
302  if (this%maxbound <= 0) then
303  write (errmsg, '(a)') 'MAXBOUND must be an integer greater than zero.'
304  call store_error(errmsg)
305  end if
306  !
307  ! -- terminate if there are errors
308  if (count_errors() > 0) then
309  call this%parser%StoreErrorUnit()
310  end if
311  !
312  ! -- return
313  return
Here is the call graph for this function:

◆ read_options()

subroutine gwtspcmodule::read_options ( class(gwtspctype this)

Read options for this package.

Definition at line 187 of file GwtSpc.f90.

188  ! -- modules
189  ! -- dummy
190  class(GwtSpcType) :: this
191  ! -- local
192  character(len=LINELENGTH) :: keyword, fname
193  integer(I4B) :: ierr
194  logical :: isfound, endOfBlock
195  ! -- formats
196  character(len=*), parameter :: fmtiprpak = &
197  &"(4x,'SPC INFORMATION WILL BE PRINTED TO LISTING FILE.')"
198  character(len=*), parameter :: fmtreadasarrays = &
199  "(4x,'SPC INFORMATION WILL BE READ AS ARRAYS RATHER THAN IN LIST &
200  &FORMAT.')"
201  character(len=*), parameter :: fmtts = &
202  &"(4x, 'TIME-SERIES DATA WILL BE READ FROM FILE: ', a)"
203  character(len=*), parameter :: fmttas = &
204  &"(4x, 'TIME-ARRAY SERIES DATA WILL BE READ FROM FILE: ', a)"
205  !
206  ! -- get options block
207  call this%parser%GetBlock('OPTIONS', isfound, ierr, blockrequired=.false., &
208  supportopenclose=.true.)
209  !
210  ! -- parse options block if detected
211  if (isfound) then
212  write (this%iout, '(1x,a)') 'PROCESSING SPC OPTIONS'
213  do
214  call this%parser%GetNextLine(endofblock)
215  if (endofblock) exit
216  call this%parser%GetStringCaps(keyword)
217  select case (keyword)
218  case ('PRINT_INPUT')
219  this%iprpak = 1
220  write (this%iout, fmtiprpak)
221  case ('READASARRAYS')
222  this%readasarrays = .true.
223  write (this%iout, fmtreadasarrays)
224  case ('TS6')
225  call this%parser%GetStringCaps(keyword)
226  if (trim(adjustl(keyword)) /= 'FILEIN') then
227  errmsg = 'TS6 keyword must be followed by "FILEIN" '// &
228  'then by filename.'
229  call store_error(errmsg)
230  end if
231  call this%parser%GetString(fname)
232  write (this%iout, fmtts) trim(fname)
233  call this%TsManager%add_tsfile(fname, this%inunit)
234  case ('TAS6')
235  call this%parser%GetStringCaps(keyword)
236  if (trim(adjustl(keyword)) /= 'FILEIN') then
237  errmsg = 'TAS6 keyword must be followed by "FILEIN" '// &
238  'then by filename.'
239  call store_error(errmsg)
240  call this%parser%StoreErrorUnit()
241  end if
242  call this%parser%GetString(fname)
243  write (this%iout, fmttas) trim(fname)
244  call this%TasManager%add_tasfile(fname)
245  case default
246  write (errmsg, '(a,a)') 'Unknown SPC option: ', trim(keyword)
247  call store_error(errmsg)
248  call this%parser%StoreErrorUnit()
249  end select
250  end do
251  write (this%iout, '(1x,a)') 'END OF SPC OPTIONS'
252  end if
253  !
254  ! -- Return
255  return
Here is the call graph for this function:

◆ set_value()

subroutine gwtspcmodule::set_value ( class(gwtspctype), intent(inout)  this,
integer(i4b), intent(in)  ival 
)

Set the floating point value in the dblvec array using strings returned from the parser. Allow for time series names.

Parameters
[in,out]thisGwtSpcType object

Definition at line 713 of file GwtSpc.f90.

714  ! -- modules
716  ! -- dummy
717  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
718  integer(I4B), intent(in) :: ival
719  ! -- local
720  character(len=LINELENGTH) :: keyword
721  integer(I4B) :: jj
722  real(DP), pointer :: bndElem => null()
723  !
724  ! -- read remainder of variables on the line
725  call this%parser%GetStringCaps(keyword)
726  select case (keyword)
727  case ('CONCENTRATION')
728  call this%parser%GetString(text)
729  jj = 1 ! For CONCENTRATION
730  bndelem => this%dblvec(ival)
731  call read_value_or_time_series_adv(text, ival, jj, bndelem, this%packName, &
732  'BND', this%tsManager, this%iprpak, &
733  'CONCENTRATION')
734 
735  end select
736  return
subroutine, public read_value_or_time_series_adv(textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, varName)
Call this subroutine from advanced packages to define timeseries link for a variable (varName).
Here is the call graph for this function:

◆ spc_ad()

subroutine gwtspcmodule::spc_ad ( class(gwtspctype), intent(inout)  this,
integer(i4b), intent(in)  nbound_flowpack,
character(len=*), intent(in)  budtxt 
)

Call the advance method on the time series so that new values are interpolated and entered into dblvec

Parameters
[in,out]thisGwtSpcType object

Definition at line 623 of file GwtSpc.f90.

624  ! -- modules
625  ! -- dummy
626  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
627  integer(I4B), intent(in) :: nbound_flowpack
628  character(len=*), intent(in) :: budtxt
629  ! -- local
630  !
631  ! -- Advance the time series
632  call this%TsManager%ad()
633  call this%TasManager%ad()
634  !
635  ! -- Check flow package consistency
636  call this%check_flow_package(nbound_flowpack, budtxt)
637  !
638  ! -- return
639  return

◆ spc_da()

subroutine gwtspcmodule::spc_da ( class(gwtspctype this)
private

Deallocate and nullify package variables.

Parameters
thisGwtSpcType object

Definition at line 647 of file GwtSpc.f90.

648  ! -- modules
650  ! -- dummy variables
651  class(GwtSpcType) :: this !< GwtSpcType object
652  !
653  ! -- deallocate arrays in memory manager
654  call mem_deallocate(this%dblvec)
655  !
656  ! -- deallocate scalars in memory manager
657  call mem_deallocate(this%id)
658  call mem_deallocate(this%inunit)
659  call mem_deallocate(this%iout)
660  call mem_deallocate(this%maxbound)
661  call mem_deallocate(this%ionper)
662  call mem_deallocate(this%lastonper)
663  call mem_deallocate(this%iprpak)
664  call mem_deallocate(this%readasarrays)
665  !
666  ! -- deallocate derived types
667  call this%TsManager%da()
668  deallocate (this%TsManager)
669  nullify (this%TsManager)
670  !
671  ! -- return
672  return

◆ spc_rp()

subroutine gwtspcmodule::spc_rp ( class(gwtspctype), intent(inout)  this)
private

Read and prepare the period data block and fill dblvec if the next period block corresponds to this time step.

Parameters
[in,out]thisGwtSpcType object

Definition at line 391 of file GwtSpc.f90.

392  ! -- modules
393  use tdismodule, only: kper, nper
394  ! -- dummy
395  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
396  ! -- local
397  character(len=LINELENGTH) :: line
398  logical :: isfound
399  integer(I4B) :: ierr
400  ! -- formats
401  character(len=*), parameter :: fmtblkerr = &
402  &"('Looking for BEGIN PERIOD iper. Found ', a, ' instead.')"
403  character(len=*), parameter :: fmtlsp = &
404  &"(1X,/1X,'REUSING ',A,'S FROM LAST STRESS PERIOD')"
405  !
406  ! -- Set ionper to the stress period number for which a new block of data
407  ! will be read.
408  if (this%inunit == 0) return
409  !
410  ! -- get stress period data
411  if (this%ionper < kper) then
412  !
413  ! -- get period block
414  call this%parser%GetBlock('PERIOD', isfound, ierr, &
415  supportopenclose=.true., &
416  blockrequired=.false.)
417  if (isfound) then
418  !
419  ! -- read ionper and check for increasing period numbers
420  call this%read_check_ionper()
421  else
422  !
423  ! -- PERIOD block not found
424  if (ierr < 0) then
425  ! -- End of file found; data applies for remainder of simulation.
426  this%ionper = nper + 1
427  else
428  ! -- Found invalid block
429  call this%parser%GetCurrentLine(line)
430  write (errmsg, fmtblkerr) adjustl(trim(line))
431  call store_error(errmsg, terminate=.true.)
432  end if
433  end if
434  end if
435  !
436  ! -- Read data if ionper == kper
437  if (this%ionper == kper) then
438  !
439  ! -- Remove all time-series and time-array-series links associated with
440  ! this package.
441  ! Do not reset as we are using a "settings" approach here in which the
442  ! settings remain the same until the user changes them.
443  ! call this%TsManager%Reset(this%packName)
444  call this%TasManager%Reset(this%packName)
445  if (this%readasarrays) then
446  call this%spc_rp_array(line)
447  else
448  call this%spc_rp_list()
449  end if
450  !
451  ! -- using data from the last stress period
452  else
453  write (this%iout, fmtlsp) trim(ftype)
454  end if
455  !
456  ! -- write summary of maw well stress period error messages
457  if (count_errors() > 0) then
458  call this%parser%StoreErrorUnit()
459  end if
460  !
461  ! -- return
462  return
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21
Here is the call graph for this function:

◆ spc_rp_array()

subroutine gwtspcmodule::spc_rp_array ( class(gwtspctype), intent(inout)  this,
character(len=linelength), intent(inout)  line 
)

Read the stress period data in array format

Parameters
[in,out]thisGwtSpcType object

Definition at line 538 of file GwtSpc.f90.

540  use simmodule, only: store_error
541  use arrayhandlersmodule, only: ifind
542  ! -- dummy
543  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
544  character(len=LINELENGTH), intent(inout) :: line
545  ! -- local
546  integer(I4B) :: n
547  integer(I4B) :: ncolbnd
548  integer(I4B) :: jauxcol, ivarsread
549  integer(I4B), dimension(:), allocatable, target :: nodelist
550  character(len=LENTIMESERIESNAME) :: tasName
551  character(len=24), dimension(1) :: aname
552  character(len=LINELENGTH) :: keyword
553  logical :: endOfBlock
554  logical :: convertFlux
555  !
556  ! -- these time array series pointers need to be non-contiguous
557  ! because a slice of bound is passed
558  real(DP), dimension(:), pointer :: bndArrayPtr => null()
559  ! -- formats
560  ! -- data
561  data aname(1)/' CONCENTRATION'/
562  !
563 ! ------------------------------------------------------------------------------
564  !
565  ! -- Initialize
566  jauxcol = 0
567  ivarsread = 0
568  ncolbnd = 1
569  allocate (nodelist(this%maxbound))
570  do n = 1, size(nodelist)
571  nodelist(n) = n
572  end do
573  !
574  ! -- Read CONCENTRATION variables as arrays
575  do
576  call this%parser%GetNextLine(endofblock)
577  if (endofblock) exit
578  call this%parser%GetStringCaps(keyword)
579  !
580  ! -- Parse the keywords
581  select case (keyword)
582  case ('CONCENTRATION')
583  !
584  ! -- Look for keyword TIMEARRAYSERIES and time-array series
585  ! name on line, following RECHARGE
586  call this%parser%GetStringCaps(keyword)
587  if (keyword == 'TIMEARRAYSERIES') then
588  ! -- Get time-array series name
589  call this%parser%GetStringCaps(tasname)
590  bndarrayptr => this%dblvec(:)
591  ! Make a time-array-series link and add it to the list of links
592  ! contained in the TimeArraySeriesManagerType object.
593  convertflux = .false.
594  call this%TasManager%MakeTasLink(this%packName, bndarrayptr, &
595  this%iprpak, tasname, &
596  'CONCENTRATION', &
597  convertflux, nodelist, &
598  this%parser%iuactive)
599  else
600  !
601  ! -- Read the concentration array
602  call this%dis%read_layer_array(nodelist, this%dblvec, ncolbnd, &
603  this%maxbound, 1, aname(1), &
604  this%parser%iuactive, this%iout)
605  end if
606  !
607  case default
608  call store_error('Looking for CONCENTRATION. Found: '//trim(line))
609  call this%parser%StoreErrorUnit()
610  end select
611 
612  end do
613  !
614  return
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lentimeseriesname
maximum length of a time series name
Definition: Constants.f90:41
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:

◆ spc_rp_list()

subroutine gwtspcmodule::spc_rp_list ( class(gwtspctype), intent(inout)  this)

Read the stress period data in list format

Parameters
[in,out]thisGwtSpcType object

Definition at line 470 of file GwtSpc.f90.

471  ! -- modules
472  use tdismodule, only: kper
473  ! -- dummy
474  class(GwtSpcType), intent(inout) :: this !< GwtSpcType object
475  ! -- local
476  character(len=LINELENGTH) :: line
477  character(len=LINELENGTH) :: title
478  character(len=LINELENGTH) :: tabletext
479  logical :: endOfBlock
480  integer(I4B) :: ival
481  !
482  !
483  ! -- setup table for period data
484  if (this%iprpak /= 0) then
485  !
486  ! -- reset the input table object
487  title = trim(adjustl(text))//' PACKAGE ('// &
488  'SPC'//') DATA FOR PERIOD'
489  write (title, '(a,1x,i6)') trim(adjustl(title)), kper
490  call table_cr(this%inputtab, ftype, title)
491  call this%inputtab%table_df(1, 3, this%iout, finalize=.false.)
492  tabletext = 'NUMBER'
493  call this%inputtab%initialize_column(tabletext, 10, alignment=tabcenter)
494  tabletext = 'DATA TYPE'
495  call this%inputtab%initialize_column(tabletext, 20, alignment=tableft)
496  write (tabletext, '(a,1x,i6)') 'VALUE'
497  call this%inputtab%initialize_column(tabletext, 15, alignment=tabcenter)
498  end if
499  !
500  ! -- read data
501  do
502  call this%parser%GetNextLine(endofblock)
503  if (endofblock) exit
504 
505  ival = this%parser%GetInteger()
506  if (ival < 1 .or. ival > this%maxbound) then
507  write (errmsg, '(2(a,1x),i0,a)') &
508  'IVAL must be greater than 0 and', &
509  'less than or equal to ', this%maxbound, '.'
510  call store_error(errmsg)
511  cycle
512  end if
513  !
514  ! -- set stress period data
515  call this%set_value(ival)
516  !
517  ! -- write line to table
518  if (this%iprpak /= 0) then
519  call this%parser%GetCurrentLine(line)
520  call this%inputtab%line_to_columns(line)
521  end if
522  end do
523  !
524  ! -- finalize the table
525  if (this%iprpak /= 0) then
526  call this%inputtab%finalize_table()
527  end if
528  !
529  ! -- return
530  return
Here is the call graph for this function:

Variable Documentation

◆ ftype

character(len=lenftype) gwtspcmodule::ftype = 'SPC'
private

Definition at line 27 of file GwtSpc.f90.

27  character(len=LENFTYPE) :: ftype = 'SPC'

◆ text

character(len=lenpackagename) gwtspcmodule::text = 'STRESS PACK CONC'
private

Definition at line 28 of file GwtSpc.f90.

28  character(len=LENPACKAGENAME) :: text = 'STRESS PACK CONC'