MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
GwtInterfaceModel.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 gwtdspmodule, only: dsp_cr, gwtdsptype
14  use gwtmstmodule, only: mst_cr
15  use tspobsmodule, only: tsp_obs_cr
17 
18  implicit none
19  private
20 
21  !> The GWT Interface Model is a utility to calculate the solution's exchange
22  !! coefficients from the interface between a GWT model and its GWT neighbors.
23  !! The interface model itself will not be part of the solution, it is not
24  !! being solved.
25  !<
26  type, public, extends(gwtmodeltype) :: gwtinterfacemodeltype
27 
28  integer(i4B), pointer :: iadvscheme => null() !< the advection scheme: 0 = up, 1 = central, 2 = tvd
29  integer(i4B), pointer :: ixt3d => null() !< xt3d setting: 0 = off, 1 = lhs, 2 = rhs
30  real(dp), pointer :: ieqnsclfac => null() !< governing eqn scaling factor: 1: GWT, >1: GWE
31 
32  class(gridconnectiontype), pointer :: gridconnection => null() !< The grid connection class will provide the interface grid
33  class(gwtmodeltype), private, pointer :: owner => null() !< the real GWT model for which the exchange coefficients
34  !! are calculated with this interface model
35 
36  contains
37  procedure, pass(this) :: gwtifmod_cr
38  procedure :: model_df => gwtifmod_df
39  procedure :: model_ar => gwtifmod_ar
40  procedure :: model_da => gwtifmod_da
41  procedure, public :: allocate_fmi
42  procedure :: allocate_scalars
43  end type gwtinterfacemodeltype
44 
45 contains
46 
47  !> @brief Create the interface model, analogously to what
48  !< happens in gwt_cr
49  subroutine gwtifmod_cr(this, name, iout, gridConn)
50  ! -- dummy
51  class(gwtinterfacemodeltype) :: this !< the GWT interface model
52  character(len=*), intent(in) :: name !< the interface model's name
53  integer(I4B), intent(in) :: iout !< the output unit
54  class(gridconnectiontype), pointer, intent(in) :: gridConn !< the grid connection data for creating a DISU
55  ! local
56  class(*), pointer :: modelPtr
57  integer(I4B), target :: inobs
58  integer(I4B) :: adv_unit, dsp_unit
59  !
60  this%memoryPath = create_mem_path(name)
61  call this%allocate_scalars(name)
62  !
63  ! defaults
64  this%iAdvScheme = 0
65  this%ixt3d = 0
66  this%ieqnsclfac = done
67  !
68  this%iout = iout
69  this%gridConnection => gridconn
70  modelptr => gridconn%model
71  this%owner => castasgwtmodel(modelptr)
72  !
73  inobs = 0
74  adv_unit = 0
75  dsp_unit = 0
76  if (this%owner%inadv > 0) then
77  this%inadv = huge(1_i4b)
78  adv_unit = huge(1_i4b)
79  end if
80  if (this%owner%indsp > 0) then
81  this%indsp = huge(1_i4b)
82  dsp_unit = huge(1_i4b)
83  end if
84  !
85  ! create dis and packages
86  call disu_cr(this%dis, this%name, '', -1, this%iout)
87  call fmi_cr(this%fmi, this%name, 0, this%iout, this%ieqnsclfac, &
88  this%depvartype)
89  call adv_cr(this%adv, this%name, adv_unit, this%iout, this%fmi, &
90  this%ieqnsclfac)
91  call dsp_cr(this%dsp, this%name, '', -dsp_unit, this%iout, this%fmi)
92  call tsp_obs_cr(this%obs, inobs, this%depvartype)
93  end subroutine gwtifmod_cr
94 
95  !> @brief Allocate scalars associated with the interface model object
96  !<
97  subroutine allocate_scalars(this, modelname)
98  ! -- dummy
99  class(gwtinterfacemodeltype) :: this !< the GWT interface model
100  character(len=*), intent(in) :: modelname !< the model name
101  !
102  call this%GwtModelType%allocate_scalars(modelname)
103  !
104  call mem_allocate(this%iAdvScheme, 'ADVSCHEME', this%memoryPath)
105  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
106  call mem_allocate(this%ieqnsclfac, 'IEQNSCLFAC', this%memoryPath)
107  end subroutine allocate_scalars
108 
109  !> @brief Allocate a Flow Model Interface (FMI) object for the interface model
110  !<
111  subroutine allocate_fmi(this)
112  ! -- dummy
113  class(gwtinterfacemodeltype) :: this !< the GWT interface model
114  !
115  call mem_allocate(this%fmi%gwfflowja, this%nja, 'GWFFLOWJA', &
116  this%fmi%memoryPath)
117  call mem_allocate(this%fmi%gwfhead, this%neq, 'GWFHEAD', &
118  this%fmi%memoryPath)
119  call mem_allocate(this%fmi%gwfsat, this%neq, 'GWFSAT', &
120  this%fmi%memoryPath)
121  call mem_allocate(this%fmi%gwfspdis, 3, this%neq, 'GWFSPDIS', &
122  this%fmi%memoryPath)
123  end subroutine allocate_fmi
124 
125  !> @brief Define the GWT interface model
126  !<
127  subroutine gwtifmod_df(this)
128  ! -- dummy
129  class(gwtinterfacemodeltype) :: this !< the GWT interface model
130  ! -- local
131  class(*), pointer :: disPtr
132  type(tspadvoptionstype) :: adv_options
133  type(gwtdspoptionstype) :: dsp_options
134  !
135  this%moffset = 0
136  adv_options%iAdvScheme = this%iAdvScheme
137  dsp_options%ixt3d = this%ixt3d
138  !
139  ! define DISU
140  disptr => this%dis
141  call this%gridConnection%getDiscretization(castasdisutype(disptr))
142  call this%fmi%fmi_df(this%dis, 1)
143  !
144  if (this%inadv > 0) then
145  call this%adv%adv_df(adv_options)
146  end if
147  if (this%indsp > 0) then
148  this%dsp%idiffc = this%owner%dsp%idiffc
149  this%dsp%idisp = this%owner%dsp%idisp
150  call this%dsp%dsp_df(this%dis, dsp_options)
151  if (this%dsp%idiffc > 0) then
152  call mem_reallocate(this%dsp%diffc, this%dis%nodes, 'DIFFC', &
153  trim(this%dsp%memoryPath))
154  end if
155  if (this%dsp%idisp > 0) then
156  call mem_reallocate(this%dsp%alh, this%dis%nodes, 'ALH', &
157  trim(this%dsp%memoryPath))
158  call mem_reallocate(this%dsp%alv, this%dis%nodes, 'ALV', &
159  trim(this%dsp%memoryPath))
160  call mem_reallocate(this%dsp%ath1, this%dis%nodes, 'ATH1', &
161  trim(this%dsp%memoryPath))
162  call mem_reallocate(this%dsp%ath2, this%dis%nodes, 'ATH2', &
163  trim(this%dsp%memoryPath))
164  call mem_reallocate(this%dsp%atv, this%dis%nodes, 'ATV', &
165  trim(this%dsp%memoryPath))
166  end if
167  allocate (this%mst)
168  call mem_allocate(this%mst%thetam, this%dis%nodes, &
169  'THETAM', create_mem_path(this%name, 'MST'))
170  end if
171  !
172  ! assign or point model members to dis members
173  this%neq = this%dis%nodes
174  this%nja = this%dis%nja
175  this%ia => this%dis%con%ia
176  this%ja => this%dis%con%ja
177  !
178  ! allocate model arrays, now that neq and nja are assigned
179  call this%allocate_arrays()
180  end subroutine gwtifmod_df
181 
182  !> @brief Override allocate and read the GWT interface model and its packages
183  !! so that we can create stuff from memory instead of input files
184  !<
185  subroutine gwtifmod_ar(this)
186  ! -- dummy
187  class(gwtinterfacemodeltype) :: this !< the GWT interface model
188  !
189  call this%fmi%fmi_ar(this%ibound)
190  if (this%inadv > 0) then
191  call this%adv%adv_ar(this%dis, this%ibound)
192  end if
193  if (this%indsp > 0) then
194  call this%dsp%dsp_ar(this%ibound, this%mst%thetam)
195  end if
196  end subroutine gwtifmod_ar
197 
198  !> @brief Clean up resources
199  !<
200  subroutine gwtifmod_da(this)
201  ! -- dummy
202  class(gwtinterfacemodeltype) :: this !< the GWT interface model
203  !
204 
205  ! this
206  call mem_deallocate(this%iAdvScheme)
207  call mem_deallocate(this%ixt3d)
208  call mem_deallocate(this%ieqnsclfac)
209  !
210  ! gwt packages
211  call this%dis%dis_da()
212  call this%fmi%fmi_da()
213  call this%adv%adv_da()
214  call this%dsp%dsp_da()
215  !
216  deallocate (this%dis)
217  deallocate (this%fmi)
218  deallocate (this%adv)
219  deallocate (this%dsp)
220  !
221  if (associated(this%mst)) then
222  call mem_deallocate(this%mst%thetam)
223  deallocate (this%mst)
224  end if
225  !
226  ! gwt scalars
227  call mem_deallocate(this%inic)
228  call mem_deallocate(this%infmi)
229  call mem_deallocate(this%inadv)
230  call mem_deallocate(this%indsp)
231  call mem_deallocate(this%inssm)
232  call mem_deallocate(this%inmst)
233  call mem_deallocate(this%inmvt)
234  call mem_deallocate(this%inoc)
235  call mem_deallocate(this%inobs)
236  call mem_deallocate(this%eqnsclfac)
237  !
238  ! base
239  call this%NumericalModelType%model_da()
240  end subroutine gwtifmod_da
241 
242 end module gwtinterfacemodelmodule
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 dsp_cr(dspobj, name_model, input_mempath, inunit, iout, fmi)
Create a DSP object.
Definition: gwt-dsp.f90:78
subroutine allocate_scalars(this, modelname)
Allocate scalars associated with the interface model object.
subroutine gwtifmod_ar(this)
Override allocate and read the GWT interface model and its packages so that we can create stuff from ...
subroutine gwtifmod_da(this)
Clean up resources.
subroutine gwtifmod_df(this)
Define the GWT interface model.
subroutine gwtifmod_cr(this, name, iout, gridConn)
Create the interface model, analogously to what.
subroutine allocate_fmi(this)
Allocate a Flow Model Interface (FMI) object for the interface model.
Definition: gwt.f90:8
class(gwtmodeltype) function, pointer, public castasgwtmodel(model)
Cast to GwtModelType.
Definition: gwt.f90:824
– @ brief Mobile Storage and Transfer (MST) Module
Definition: gwt-mst.f90:10
subroutine, public mst_cr(mstobj, name_model, inunit, iout, fmi)
@ brief Create a new package object
Definition: gwt-mst.f90:103
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 dsp option data
The GWT Interface Model is a utility to calculate the solution's exchange coefficients from the inter...