MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
GweInterfaceModel.f90
Go to the documentation of this file.
2  use kindmodule, only: i4b, dp
3  use constantsmodule, only: done
9  use tspfmimodule, only: fmi_cr, tspfmitype
10  use tspadvmodule, only: adv_cr, tspadvtype
12  use gwecndmodule, only: cnd_cr, gwecndtype
14  use gweestmodule, only: est_cr
15  use tspobsmodule, only: tsp_obs_cr
18 
19  implicit none
20  private
21 
22  !> The GWE Interface Model is a utility to calculate the solution's exchange
23  !! coefficients from the interface between a GWE model and its GWE neighbors.
24  !! The interface model itself will not be part of the solution, it is not
25  !! being solved.
26  !<
27  type, public, extends(gwemodeltype) :: gweinterfacemodeltype
28 
29  integer(i4B), pointer :: iadvscheme => null() !< the advection scheme: 0 = up, 1 = central, 2 = tvd
30  integer(i4B), pointer :: ixt3d => null() !< xt3d setting: 0 = off, 1 = lhs, 2 = rhs
31  real(dp), pointer :: ieqnsclfac => null() !< governing eqn scaling factor: 1: GWT, >1: GWE
32 
33  class(gridconnectiontype), pointer :: gridconnection => null() !< The grid connection class will provide the interface grid
34  class(gwemodeltype), private, pointer :: owner => null() !< the real GWE model for which the exchange coefficients
35  !! are calculated with this interface model
36 
37  real(dp), dimension(:), pointer, contiguous :: porosity => null() !< to be filled with EST porosity
38 
39  contains
40 
41  procedure, pass(this) :: gweifmod_cr
42  procedure :: model_df => gweifmod_df
43  procedure :: model_ar => gweifmod_ar
44  procedure :: model_da => gweifmod_da
45  procedure, public :: allocate_fmi
46  procedure :: allocate_scalars
47  end type gweinterfacemodeltype
48 
49 contains
50 
51  !> @brief Create the interface model, analogously to what
52  !< happens in gwe_cr
53  subroutine gweifmod_cr(this, name, iout, gridConn)
54  ! -- modules
56  ! -- dummy
57  class(gweinterfacemodeltype) :: this !< the GWE interface model
58  character(len=*), intent(in) :: name !< the interface model's name
59  integer(I4B), intent(in) :: iout !< the output unit
60  class(gridconnectiontype), pointer, intent(in) :: gridConn !< the grid connection data for creating a DISU
61  ! -- local
62  class(*), pointer :: modelPtr
63  integer(I4B), target :: inobs
64  integer(I4B) :: adv_unit, cnd_unit
65  !
66  this%memoryPath = create_mem_path(name)
67  call this%allocate_scalars(name)
68  !
69  ! -- Instantiate shared data container
70  call gweshared_dat_cr(this%gwecommon)
71  !
72  ! -- Defaults
73  this%iAdvScheme = 0
74  this%ixt3d = 0
75  this%ieqnsclfac = done
76  !
77  this%iout = iout
78  this%gridConnection => gridconn
79  modelptr => gridconn%model
80  this%owner => castasgwemodel(modelptr)
81  !
82  inobs = 0
83  adv_unit = 0
84  cnd_unit = 0
85  if (this%owner%inadv > 0) then
86  this%inadv = huge(1_i4b)
87  adv_unit = huge(1_i4b)
88  end if
89  if (this%owner%incnd > 0) then
90  this%incnd = huge(1_i4b)
91  cnd_unit = huge(1_i4b)
92  end if
93  !
94  ! -- Create dis and packages
95  call disu_cr(this%dis, this%name, '', -1, this%iout)
96  call fmi_cr(this%fmi, this%name, 0, this%iout, this%ieqnsclfac, &
97  this%depvartype)
98  call adv_cr(this%adv, this%name, adv_unit, this%iout, this%fmi, &
99  this%ieqnsclfac)
100  call cnd_cr(this%cnd, this%name, '', -cnd_unit, this%iout, this%fmi, &
101  this%ieqnsclfac, this%gwecommon)
102  call tsp_obs_cr(this%obs, inobs, this%depvartype)
103  end subroutine gweifmod_cr
104 
105  !> @brief Allocate scalars associated with the interface model object
106  !<
107  subroutine allocate_scalars(this, modelname)
108  ! -- dummy
109  class(gweinterfacemodeltype) :: this !< the GWE interface model
110  character(len=*), intent(in) :: modelname !< the model name
111  !
112  call this%GweModelType%allocate_scalars(modelname)
113  !
114  call mem_allocate(this%iAdvScheme, 'ADVSCHEME', this%memoryPath)
115  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
116  call mem_allocate(this%ieqnsclfac, 'IEQNSCLFAC', this%memoryPath)
117  end subroutine allocate_scalars
118 
119  !> @brief Allocate a Flow Model Interface (FMI) object for the interface model
120  !<
121  subroutine allocate_fmi(this)
122  ! -- dummy
123  class(gweinterfacemodeltype) :: this !< the GWT interface model
124  !
125  call mem_allocate(this%fmi%gwfflowja, this%nja, 'GWFFLOWJA', &
126  this%fmi%memoryPath)
127  call mem_allocate(this%fmi%gwfhead, this%neq, 'GWFHEAD', &
128  this%fmi%memoryPath)
129  call mem_allocate(this%fmi%gwfsat, this%neq, 'GWFSAT', &
130  this%fmi%memoryPath)
131  call mem_allocate(this%fmi%gwfspdis, 3, this%neq, 'GWFSPDIS', &
132  this%fmi%memoryPath)
133  end subroutine allocate_fmi
134 
135  !> @brief Define the GWE interface model
136  !<
137  subroutine gweifmod_df(this)
138  ! -- dummy
139  class(gweinterfacemodeltype) :: this !< the GWE interface model
140  ! -- local
141  class(*), pointer :: disPtr
142  type(tspadvoptionstype) :: adv_options
143  type(gwecndoptionstype) :: cnd_options
144  !
145  this%moffset = 0
146  adv_options%iAdvScheme = this%iAdvScheme
147  cnd_options%ixt3d = this%ixt3d
148  !
149  ! -- Define DISU
150  disptr => this%dis
151  call this%gridConnection%getDiscretization(castasdisutype(disptr))
152  call this%fmi%fmi_df(this%dis, 0)
153  !
154  if (this%inadv > 0) then
155  call this%adv%adv_df(adv_options)
156  end if
157  !
158  if (this%incnd > 0) then
159  this%cnd%idisp = this%owner%cnd%idisp
160  call this%cnd%cnd_df(this%dis, cnd_options)
161  !
162  if (this%cnd%idisp > 0) then
163  call mem_reallocate(this%cnd%alh, this%dis%nodes, 'ALH', &
164  trim(this%cnd%memoryPath))
165  call mem_reallocate(this%cnd%alv, this%dis%nodes, 'ALV', &
166  trim(this%cnd%memoryPath))
167  call mem_reallocate(this%cnd%ath1, this%dis%nodes, 'ATH1', &
168  trim(this%cnd%memoryPath))
169  call mem_reallocate(this%cnd%ath2, this%dis%nodes, 'ATH2', &
170  trim(this%cnd%memoryPath))
171  call mem_reallocate(this%cnd%atv, this%dis%nodes, 'ATV', &
172  trim(this%cnd%memoryPath))
173  call mem_reallocate(this%cnd%ktw, this%dis%nodes, 'KTW', &
174  trim(this%cnd%memoryPath))
175  call mem_reallocate(this%cnd%kts, this%dis%nodes, 'KTS', &
176  trim(this%cnd%memoryPath))
177  end if
178  allocate (this%est)
179  call mem_allocate(this%est%porosity, this%dis%nodes, &
180  'POROSITY', create_mem_path(this%name, 'EST'))
181  end if
182  !
183  ! -- Assign or point model members to dis members
184  this%neq = this%dis%nodes
185  this%nja = this%dis%nja
186  this%ia => this%dis%con%ia
187  this%ja => this%dis%con%ja
188  !
189  ! -- Allocate model arrays, now that neq and nja are assigned
190  call this%allocate_arrays()
191  end subroutine gweifmod_df
192 
193  !> @brief Override allocate and read the GWE interface model and its packages
194  !! so that we can create stuff from memory instead of input files
195  !<
196  subroutine gweifmod_ar(this)
197  ! -- dummy
198  class(gweinterfacemodeltype) :: this !< the GWE interface model
199  !
200  call this%fmi%fmi_ar(this%ibound)
201  if (this%inadv > 0) then
202  call this%adv%adv_ar(this%dis, this%ibound)
203  end if
204  if (this%incnd > 0) then
205  call this%cnd%cnd_ar(this%ibound, this%est%porosity)
206  end if
207  end subroutine gweifmod_ar
208 
209  !> @brief Clean up resources
210  !<
211  subroutine gweifmod_da(this)
212  ! -- dummy
213  class(gweinterfacemodeltype) :: this !< the GWE interface model
214  !
215  ! -- members specified in the interface model input for use by the
216  ! interface model
217  call mem_deallocate(this%iAdvScheme)
218  call mem_deallocate(this%ixt3d)
219  call mem_deallocate(this%ieqnsclfac)
220  !
221  ! -- gwe packages
222  call this%dis%dis_da()
223  call this%fmi%fmi_da()
224  call this%adv%adv_da()
225  call this%cnd%cnd_da()
226  !
227  deallocate (this%dis)
228  deallocate (this%fmi)
229  deallocate (this%adv)
230  deallocate (this%cnd)
231  !
232  if (associated(this%est)) then
233  call mem_deallocate(this%est%porosity)
234  deallocate (this%est)
235  end if
236  !
237  ! -- gwe scalars
238  call mem_deallocate(this%inic)
239  call mem_deallocate(this%infmi)
240  call mem_deallocate(this%inadv)
241  call mem_deallocate(this%incnd)
242  call mem_deallocate(this%inssm)
243  call mem_deallocate(this%inest)
244  call mem_deallocate(this%inmvt)
245  call mem_deallocate(this%inoc)
246  call mem_deallocate(this%inobs)
247  call mem_deallocate(this%eqnsclfac)
248  !
249  ! -- base
250  call this%NumericalModelType%model_da()
251  end subroutine gweifmod_da
252 
253 end module gweinterfacemodelmodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine, public disu_cr(dis, name_model, input_mempath, inunit, iout)
Create a new unstructured discretization object.
Definition: Disu.f90:127
class(disutype) function, pointer, public castasdisutype(dis)
Cast base to DISU.
Definition: Disu.f90:1592
Refactoring issues towards parallel:
subroutine, public cnd_cr(cndobj, name_model, input_mempath, inunit, iout, fmi, eqnsclfac, gwecommon)
Create a new CND object.
Definition: gwe-cnd.f90:86
– @ brief Energy Storage and Transfer (EST) Module
Definition: gwe-est.f90:13
subroutine, public est_cr(estobj, name_model, inunit, iout, fmi, eqnsclfac, gwecommon)
@ brief Create a new EST package object
Definition: gwe-est.f90:87
subroutine, public gweshared_dat_cr(gwe_dat)
Allocate the shared data.
subroutine gweifmod_df(this)
Define the GWE interface model.
subroutine gweifmod_ar(this)
Override allocate and read the GWE interface model and its packages so that we can create stuff from ...
subroutine allocate_fmi(this)
Allocate a Flow Model Interface (FMI) object for the interface model.
subroutine gweifmod_cr(this, name, iout, gridConn)
Create the interface model, analogously to what.
subroutine gweifmod_da(this)
Clean up resources.
subroutine allocate_scalars(this, modelname)
Allocate scalars associated with the interface model object.
Definition: gwe.f90:3
class(gwemodeltype) function, pointer, public castasgwemodel(model)
Cast to GweModelType.
Definition: gwe.f90:771
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
subroutine, public adv_cr(advobj, name_model, inunit, iout, fmi, eqnsclfac)
@ brief Create a new ADV object
Definition: tsp-adv.f90:50
subroutine, public fmi_cr(fmiobj, name_model, inunit, iout, eqnsclfac, depvartype)
Create a new FMI object.
Definition: tsp-fmi.f90:75
subroutine, public tsp_obs_cr(obs, inobs, dvt)
Create a new TspObsType object.
Definition: tsp-obs.f90:44
This class is used to construct the connections object for the interface model's spatial discretizati...
data structure (and helpers) for passing cnd option data
Data for sharing among multiple packages. Originally read in from.
The GWE Interface Model is a utility to calculate the solution's exchange coefficients from the inter...