MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
gwe-esl.f90
Go to the documentation of this file.
2  !
3  use kindmodule, only: dp, i4b
5  use bndmodule, only: bndtype
12  !
13  implicit none
14  !
15  private
16  public :: esl_create
17  !
18  character(len=LENFTYPE) :: ftype = 'ESL'
19  character(len=16) :: text = ' ESL'
20  !
21  type, extends(bndtype) :: gweesltype
22 
23  type(gweinputdatatype), pointer :: gwecommon => null() !< pointer to shared gwe data used by multiple packages but set in mst
24 
25  contains
26 
27  procedure :: allocate_scalars => esl_allocate_scalars
28  procedure :: bnd_cf => esl_cf
29  procedure :: bnd_fc => esl_fc
30  procedure :: bnd_da => esl_da
31  procedure :: define_listlabel
32  ! -- methods for observations
33  procedure, public :: bnd_obs_supported => esl_obs_supported
34  procedure, public :: bnd_df_obs => esl_df_obs
35  ! -- methods for time series
36  procedure, public :: bnd_rp_ts => esl_rp_ts
37 
38  end type gweesltype
39 
40 contains
41 
42  !> @brief Create an energy source loading package
43  !!
44  !! This subroutine points bndobj to the newly created package
45  !<
46  subroutine esl_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
47  gwecommon)
48  ! -- dummy
49  class(bndtype), pointer :: packobj
50  integer(I4B), intent(in) :: id
51  integer(I4B), intent(in) :: ibcnum
52  integer(I4B), intent(in) :: inunit
53  integer(I4B), intent(in) :: iout
54  character(len=*), intent(in) :: namemodel
55  character(len=*), intent(in) :: pakname
56  type(gweinputdatatype), intent(in), target :: gwecommon !< shared data container for use by multiple GWE packages
57  ! -- local
58  type(gweesltype), pointer :: eslobj
59  !
60  ! -- Allocate the object and assign values to object variables
61  allocate (eslobj)
62  packobj => eslobj
63  !
64  ! -- Create name and memory path
65  call packobj%set_names(ibcnum, namemodel, pakname, ftype)
66  packobj%text = text
67  !
68  ! -- Allocate scalars
69  call eslobj%allocate_scalars()
70  !
71  ! -- Initialize package
72  call packobj%pack_initialize()
73  !
74  packobj%inunit = inunit
75  packobj%iout = iout
76  packobj%id = id
77  packobj%ibcnum = ibcnum
78  packobj%ncolbnd = 1
79  packobj%iscloc = 1
80  !
81  ! -- Store pointer to shared data module for accessing cpw, rhow
82  ! for the budget calculations, and for accessing the latent heat of
83  ! vaporization for evaporative cooling.
84  eslobj%gwecommon => gwecommon
85  end subroutine esl_create
86 
87  !> @brief Deallocate memory
88  !<
89  subroutine esl_da(this)
90  ! -- modules
92  ! -- dummy
93  class(gweesltype) :: this
94  !
95  ! -- Deallocate parent package
96  call this%BndType%bnd_da()
97  end subroutine esl_da
98 
99  !> @brief Allocate scalars
100  !!
101  !! Allocate scalars specific to this energy source loading package
102  !<
103  subroutine esl_allocate_scalars(this)
104  ! -- modules
106  ! -- dummy
107  class(gweesltype) :: this
108  !
109  ! -- call standard BndType allocate scalars
110  call this%BndType%allocate_scalars()
111  !
112  ! -- allocate the object and assign values to object variables
113  !
114  ! -- Set values
115  end subroutine esl_allocate_scalars
116 
117  !> @brief Formulate the HCOF and RHS terms
118  !!
119  !! This subroutine:
120  !! - calculates hcof and rhs terms
121  !! - skip if no sources
122  !<
123  subroutine esl_cf(this)
124  ! -- dummy
125  class(gweesltype) :: this
126  ! -- local
127  integer(I4B) :: i, node
128  real(DP) :: q
129  !
130  ! -- Return if no sources
131  if (this%nbound == 0) return
132  !
133  ! -- Calculate hcof and rhs for each source entry
134  do i = 1, this%nbound
135  node = this%nodelist(i)
136  this%hcof(i) = dzero
137  if (this%ibound(node) <= 0) then
138  this%rhs(i) = dzero
139  cycle
140  end if
141  q = this%bound(1, i)
142  this%rhs(i) = -q
143  end do
144  end subroutine esl_cf
145 
146  !> @brief Add matrix terms related to specified energy source loading
147  !!
148  !! Copy rhs and hcof into solution rhs and amat
149  !<
150  subroutine esl_fc(this, rhs, ia, idxglo, matrix_sln)
151  ! -- dummy
152  class(gweesltype) :: this
153  real(DP), dimension(:), intent(inout) :: rhs
154  integer(I4B), dimension(:), intent(in) :: ia
155  integer(I4B), dimension(:), intent(in) :: idxglo
156  class(matrixbasetype), pointer :: matrix_sln
157  ! -- local
158  integer(I4B) :: i, n, ipos
159  !
160  ! -- pakmvrobj fc
161  if (this%imover == 1) then
162  call this%pakmvrobj%fc()
163  end if
164  !
165  ! -- Copy package rhs and hcof into solution rhs and amat
166  do i = 1, this%nbound
167  n = this%nodelist(i)
168  rhs(n) = rhs(n) + this%rhs(i)
169  ipos = ia(n)
170  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
171  !
172  ! -- If mover is active and mass is being withdrawn,
173  ! store available mass (as positive value).
174  if (this%imover == 1 .and. this%rhs(i) > dzero) then
175  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
176  end if
177  end do
178  end subroutine esl_fc
179 
180  !> @brief Define list labels
181  !!
182  !! Define the list heading that is written to iout when
183  !! PRINT_INPUT option is used.
184  !<
185  subroutine define_listlabel(this)
186  ! -- dummy
187  class(gweesltype), intent(inout) :: this
188  !
189  ! -- Create the header list label
190  this%listlabel = trim(this%filtyp)//' NO.'
191  if (this%dis%ndim == 3) then
192  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
193  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
194  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
195  elseif (this%dis%ndim == 2) then
196  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
197  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
198  else
199  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
200  end if
201  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
202  if (this%inamedbound == 1) then
203  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
204  end if
205  end subroutine define_listlabel
206 
207  ! -- Procedures related to observations
208 
209  !> @brief Support function for specified energy source loading observations
210  !!
211  !! This function:
212  !! - returns true because ESL package supports observations.
213  !! - overrides BndType%bnd_obs_supported()
214  !<
215  logical function esl_obs_supported(this)
216  implicit none
217  ! -- dummy
218  class(gweesltype) :: this
219  !
220  esl_obs_supported = .true.
221  end function esl_obs_supported
222 
223  !> @brief Define observations
224  !!
225  !! This subroutine:
226  !! - stores observation types supported by ESL package.
227  !! - overrides BndType%bnd_df_obs
228  !<
229  subroutine esl_df_obs(this)
230  implicit none
231  ! -- dummy
232  class(gweesltype) :: this
233  ! -- local
234  integer(I4B) :: indx
235  !
236  call this%obs%StoreObsType('esl', .true., indx)
237  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
238  !
239  ! -- Store obs type and assign procedure pointer
240  ! for to-mvr observation type.
241  call this%obs%StoreObsType('to-mvr', .true., indx)
242  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
243  end subroutine esl_df_obs
244 
245  !> @brief Procedure related to time series
246  !!
247  !! Assign tsLink%Text appropriately for all time series in use by package.
248  !! In the ESL package only the SENERRATE variable can be controlled by time
249  !! series.
250  !<
251  subroutine esl_rp_ts(this)
252  ! -- dummy
253  class(gweesltype), intent(inout) :: this
254  ! -- local
255  integer(I4B) :: i, nlinks
256  type(timeserieslinktype), pointer :: tslink => null()
257  !
258  nlinks = this%TsManager%boundtslinks%Count()
259  do i = 1, nlinks
260  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
261  if (associated(tslink)) then
262  if (tslink%JCol == 1) then
263  tslink%Text = 'SMASSRATE'
264  end if
265  end if
266  end do
267  end subroutine esl_rp_ts
268 
269 end module gweeslmodule
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dem1
real constant 1e-1
Definition: Constants.f90:103
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine, public esl_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, gwecommon)
Create an energy source loading package.
Definition: gwe-esl.f90:48
subroutine esl_da(this)
Deallocate memory.
Definition: gwe-esl.f90:90
subroutine esl_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwe-esl.f90:124
character(len=lenftype) ftype
Definition: gwe-esl.f90:18
subroutine define_listlabel(this)
Define list labels.
Definition: gwe-esl.f90:186
subroutine esl_df_obs(this)
Define observations.
Definition: gwe-esl.f90:230
subroutine esl_fc(this, rhs, ia, idxglo, matrix_sln)
Add matrix terms related to specified energy source loading.
Definition: gwe-esl.f90:151
logical function esl_obs_supported(this)
Support function for specified energy source loading observations.
Definition: gwe-esl.f90:216
character(len=16) text
Definition: gwe-esl.f90:19
subroutine esl_rp_ts(this)
Procedure related to time series.
Definition: gwe-esl.f90:252
subroutine esl_allocate_scalars(this)
Allocate scalars.
Definition: gwe-esl.f90:104
This module defines variable data types.
Definition: kind.f90:8
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:246
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType
Data for sharing among multiple packages. Originally read in from.