MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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  !
86  ! -- Return
87  return
88  end subroutine esl_create
89 
90  !> @brief Deallocate memory
91  !<
92  subroutine esl_da(this)
93  ! -- modules
95  ! -- dummy
96  class(gweesltype) :: this
97  !
98  ! -- Deallocate parent package
99  call this%BndType%bnd_da()
100  !
101  ! -- Return
102  return
103  end subroutine esl_da
104 
105  !> @brief Allocate scalars
106  !!
107  !! Allocate scalars specific to this energy source loading package
108  !<
109  subroutine esl_allocate_scalars(this)
110  ! -- modules
112  ! -- dummy
113  class(gweesltype) :: this
114  !
115  ! -- call standard BndType allocate scalars
116  call this%BndType%allocate_scalars()
117  !
118  ! -- allocate the object and assign values to object variables
119  !
120  ! -- Set values
121  !
122  ! -- Return
123  return
124  end subroutine esl_allocate_scalars
125 
126  !> @brief Formulate the HCOF and RHS terms
127  !!
128  !! This subroutine:
129  !! - calculates hcof and rhs terms
130  !! - skip if no sources
131  !<
132  subroutine esl_cf(this)
133  ! -- dummy
134  class(gweesltype) :: this
135  ! -- local
136  integer(I4B) :: i, node
137  real(DP) :: q
138  !
139  ! -- Return if no sources
140  if (this%nbound == 0) return
141  !
142  ! -- Calculate hcof and rhs for each source entry
143  do i = 1, this%nbound
144  node = this%nodelist(i)
145  this%hcof(i) = dzero
146  if (this%ibound(node) <= 0) then
147  this%rhs(i) = dzero
148  cycle
149  end if
150  q = this%bound(1, i)
151  this%rhs(i) = -q
152  end do
153  !
154  ! -- Return
155  return
156  end subroutine esl_cf
157 
158  !> @brief Add matrix terms related to specified energy source loading
159  !!
160  !! Copy rhs and hcof into solution rhs and amat
161  !<
162  subroutine esl_fc(this, rhs, ia, idxglo, matrix_sln)
163  ! -- dummy
164  class(gweesltype) :: this
165  real(DP), dimension(:), intent(inout) :: rhs
166  integer(I4B), dimension(:), intent(in) :: ia
167  integer(I4B), dimension(:), intent(in) :: idxglo
168  class(matrixbasetype), pointer :: matrix_sln
169  ! -- local
170  integer(I4B) :: i, n, ipos
171  !
172  ! -- pakmvrobj fc
173  if (this%imover == 1) then
174  call this%pakmvrobj%fc()
175  end if
176  !
177  ! -- Copy package rhs and hcof into solution rhs and amat
178  do i = 1, this%nbound
179  n = this%nodelist(i)
180  rhs(n) = rhs(n) + this%rhs(i)
181  ipos = ia(n)
182  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
183  !
184  ! -- If mover is active and mass is being withdrawn,
185  ! store available mass (as positive value).
186  if (this%imover == 1 .and. this%rhs(i) > dzero) then
187  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
188  end if
189  end do
190  !
191  ! -- Return
192  return
193  end subroutine esl_fc
194 
195  !> @brief Define list labels
196  !!
197  !! Define the list heading that is written to iout when
198  !! PRINT_INPUT option is used.
199  !<
200  subroutine define_listlabel(this)
201  ! -- dummy
202  class(gweesltype), intent(inout) :: this
203  !
204  ! -- Create the header list label
205  this%listlabel = trim(this%filtyp)//' NO.'
206  if (this%dis%ndim == 3) then
207  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
208  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
209  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
210  elseif (this%dis%ndim == 2) then
211  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
212  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
213  else
214  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
215  end if
216  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
217  if (this%inamedbound == 1) then
218  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
219  end if
220  !
221  ! -- Return
222  return
223  end subroutine define_listlabel
224 
225  ! -- Procedures related to observations
226 
227  !> @brief Support function for specified energy source loading observations
228  !!
229  !! This function:
230  !! - returns true because ESL package supports observations.
231  !! - overrides BndType%bnd_obs_supported()
232  !<
233  logical function esl_obs_supported(this)
234  implicit none
235  ! -- dummy
236  class(gweesltype) :: this
237  !
238  esl_obs_supported = .true.
239  !
240  ! -- Return
241  return
242  end function esl_obs_supported
243 
244  !> @brief Define observations
245  !!
246  !! This subroutine:
247  !! - stores observation types supported by ESL package.
248  !! - overrides BndType%bnd_df_obs
249  !<
250  subroutine esl_df_obs(this)
251  implicit none
252  ! -- dummy
253  class(gweesltype) :: this
254  ! -- local
255  integer(I4B) :: indx
256  !
257  call this%obs%StoreObsType('esl', .true., indx)
258  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
259  !
260  ! -- Store obs type and assign procedure pointer
261  ! for to-mvr observation type.
262  call this%obs%StoreObsType('to-mvr', .true., indx)
263  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
264  !
265  ! -- Return
266  return
267  end subroutine esl_df_obs
268 
269  !> @brief Procedure related to time series
270  !!
271  !! Assign tsLink%Text appropriately for all time series in use by package.
272  !! In the ESL package only the SENERRATE variable can be controlled by time
273  !! series.
274  !<
275  subroutine esl_rp_ts(this)
276  ! -- dummy
277  class(gweesltype), intent(inout) :: this
278  ! -- local
279  integer(I4B) :: i, nlinks
280  type(timeserieslinktype), pointer :: tslink => null()
281  !
282  nlinks = this%TsManager%boundtslinks%Count()
283  do i = 1, nlinks
284  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
285  if (associated(tslink)) then
286  if (tslink%JCol == 1) then
287  tslink%Text = 'SMASSRATE'
288  end if
289  end if
290  end do
291  !
292  ! -- Return
293  return
294  end subroutine esl_rp_ts
295 
296 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: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
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:93
subroutine esl_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwe-esl.f90:133
character(len=lenftype) ftype
Definition: gwe-esl.f90:18
subroutine define_listlabel(this)
Define list labels.
Definition: gwe-esl.f90:201
subroutine esl_df_obs(this)
Define observations.
Definition: gwe-esl.f90:251
subroutine esl_fc(this, rhs, ia, idxglo, matrix_sln)
Add matrix terms related to specified energy source loading.
Definition: gwe-esl.f90:163
logical function esl_obs_supported(this)
Support function for specified energy source loading observations.
Definition: gwe-esl.f90:234
character(len=16) text
Definition: gwe-esl.f90:19
subroutine esl_rp_ts(this)
Procedure related to time series.
Definition: gwe-esl.f90:276
subroutine esl_allocate_scalars(this)
Allocate scalars.
Definition: gwe-esl.f90:110
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:249
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.