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

Data Types

type  gweinterfacemodeltype
 The GWE Interface Model is a utility to calculate the solution's exchange coefficients from the interface between a GWE model and its GWE neighbors. The interface model itself will not be part of the solution, it is not being solved. More...
 

Functions/Subroutines

subroutine gweifmod_cr (this, name, iout, gridConn)
 Create the interface model, analogously to what. More...
 
subroutine allocate_scalars (this, modelname)
 Allocate scalars associated with the interface model object. More...
 
subroutine allocate_fmi (this)
 Allocate a Flow Model Interface (FMI) object for the interface model. More...
 
subroutine gweifmod_df (this)
 Define the GWE interface model. More...
 
subroutine gweifmod_ar (this)
 Override allocate and read the GWE interface model and its packages so that we can create stuff from memory instead of input files. More...
 
subroutine gweifmod_da (this)
 Clean up resources. More...
 

Function/Subroutine Documentation

◆ allocate_fmi()

subroutine gweinterfacemodelmodule::allocate_fmi ( class(gweinterfacemodeltype this)
private
Parameters
thisthe GWT interface model

Definition at line 127 of file GweInterfaceModel.f90.

128  ! -- dummy
129  class(GweInterfaceModelType) :: this !< the GWT interface model
130  !
131  call mem_allocate(this%fmi%gwfflowja, this%nja, 'GWFFLOWJA', &
132  this%fmi%memoryPath)
133  call mem_allocate(this%fmi%gwfhead, this%neq, 'GWFHEAD', &
134  this%fmi%memoryPath)
135  call mem_allocate(this%fmi%gwfsat, this%neq, 'GWFSAT', &
136  this%fmi%memoryPath)
137  call mem_allocate(this%fmi%gwfspdis, 3, this%neq, 'GWFSPDIS', &
138  this%fmi%memoryPath)
139  !
140  ! -- Return
141  return

◆ allocate_scalars()

subroutine gweinterfacemodelmodule::allocate_scalars ( class(gweinterfacemodeltype this,
character(len=*), intent(in)  modelname 
)
Parameters
thisthe GWE interface model
[in]modelnamethe model name

Definition at line 110 of file GweInterfaceModel.f90.

111  ! -- dummy
112  class(GweInterfaceModelType) :: this !< the GWE interface model
113  character(len=*), intent(in) :: modelname !< the model name
114  !
115  call this%GweModelType%allocate_scalars(modelname)
116  !
117  call mem_allocate(this%iAdvScheme, 'ADVSCHEME', this%memoryPath)
118  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
119  call mem_allocate(this%ieqnsclfac, 'IEQNSCLFAC', this%memoryPath)
120  !
121  ! -- Return
122  return

◆ gweifmod_ar()

subroutine gweinterfacemodelmodule::gweifmod_ar ( class(gweinterfacemodeltype this)
private
Parameters
thisthe GWE interface model

Definition at line 211 of file GweInterfaceModel.f90.

212  ! -- dummy
213  class(GweInterfaceModelType) :: this !< the GWE interface model
214  !
215  call this%fmi%fmi_ar(this%ibound)
216  if (this%inadv > 0) then
217  call this%adv%adv_ar(this%dis, this%ibound)
218  end if
219  if (this%incnd > 0) then
220  call this%cnd%cnd_ar(this%ibound, this%est%porosity)
221  end if
222  !
223  ! -- Return
224  return

◆ gweifmod_cr()

subroutine gweinterfacemodelmodule::gweifmod_cr ( class(gweinterfacemodeltype this,
character(len=*), intent(in)  name,
integer(i4b), intent(in)  iout,
class(gridconnectiontype), intent(in), pointer  gridConn 
)
Parameters
thisthe GWE interface model
[in]namethe interface model's name
[in]ioutthe output unit
[in]gridconnthe grid connection data for creating a DISU

Definition at line 53 of file GweInterfaceModel.f90.

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)
103  !
104  ! -- Return
105  return
subroutine, public gweshared_dat_cr(gwe_dat)
Allocate the shared data.
Here is the call graph for this function:

◆ gweifmod_da()

subroutine gweinterfacemodelmodule::gweifmod_da ( class(gweinterfacemodeltype this)
private
Parameters
thisthe GWE interface model

Definition at line 229 of file GweInterfaceModel.f90.

230  ! -- dummy
231  class(GweInterfaceModelType) :: this !< the GWE interface model
232  !
233  ! -- members specified in the interface model input for use by the
234  ! interface model
235  call mem_deallocate(this%iAdvScheme)
236  call mem_deallocate(this%ixt3d)
237  call mem_deallocate(this%ieqnsclfac)
238  !
239  ! -- gwe packages
240  call this%dis%dis_da()
241  call this%fmi%fmi_da()
242  call this%adv%adv_da()
243  call this%cnd%cnd_da()
244  !
245  deallocate (this%dis)
246  deallocate (this%fmi)
247  deallocate (this%adv)
248  deallocate (this%cnd)
249  !
250  if (associated(this%est)) then
251  call mem_deallocate(this%est%porosity)
252  deallocate (this%est)
253  end if
254  !
255  ! -- gwe scalars
256  call mem_deallocate(this%inic)
257  call mem_deallocate(this%infmi)
258  call mem_deallocate(this%inadv)
259  call mem_deallocate(this%incnd)
260  call mem_deallocate(this%inssm)
261  call mem_deallocate(this%inest)
262  call mem_deallocate(this%inmvt)
263  call mem_deallocate(this%inoc)
264  call mem_deallocate(this%inobs)
265  call mem_deallocate(this%eqnsclfac)
266  !
267  ! -- base
268  call this%NumericalModelType%model_da()
269  !
270  ! -- Return
271  return

◆ gweifmod_df()

subroutine gweinterfacemodelmodule::gweifmod_df ( class(gweinterfacemodeltype this)
private
Parameters
thisthe GWE interface model

Definition at line 146 of file GweInterfaceModel.f90.

147  ! -- dummy
148  class(GweInterfaceModelType) :: this !< the GWE interface model
149  ! -- local
150  class(*), pointer :: disPtr
151  type(TspAdvOptionsType) :: adv_options
152  type(GweCndOptionsType) :: cnd_options
153  !
154  this%moffset = 0
155  adv_options%iAdvScheme = this%iAdvScheme
156  cnd_options%ixt3d = this%ixt3d
157  !
158  ! -- Define DISU
159  disptr => this%dis
160  call this%gridConnection%getDiscretization(castasdisutype(disptr))
161  call this%fmi%fmi_df(this%dis, 0)
162  !
163  if (this%inadv > 0) then
164  call this%adv%adv_df(adv_options)
165  end if
166  !
167  if (this%incnd > 0) then
168  this%cnd%idisp = this%owner%cnd%idisp
169  call this%cnd%cnd_df(this%dis, cnd_options)
170  !
171  if (this%cnd%idisp > 0) then
172  call mem_reallocate(this%cnd%alh, this%dis%nodes, 'ALH', &
173  trim(this%cnd%memoryPath))
174  call mem_reallocate(this%cnd%alv, this%dis%nodes, 'ALV', &
175  trim(this%cnd%memoryPath))
176  call mem_reallocate(this%cnd%ath1, this%dis%nodes, 'ATH1', &
177  trim(this%cnd%memoryPath))
178  call mem_reallocate(this%cnd%ath2, this%dis%nodes, 'ATH2', &
179  trim(this%cnd%memoryPath))
180  call mem_reallocate(this%cnd%atv, this%dis%nodes, 'ATV', &
181  trim(this%cnd%memoryPath))
182  call mem_reallocate(this%cnd%ktw, this%dis%nodes, 'KTW', &
183  trim(this%cnd%memoryPath))
184  call mem_reallocate(this%cnd%kts, this%dis%nodes, 'KTS', &
185  trim(this%cnd%memoryPath))
186  end if
187  allocate (this%est)
188  call mem_allocate(this%est%porosity, this%dis%nodes, &
189  'POROSITY', create_mem_path(this%name, 'EST'))
190  end if
191  !
192  ! -- Assign or point model members to dis members
193  this%neq = this%dis%nodes
194  this%nja = this%dis%nja
195  this%ia => this%dis%con%ia
196  this%ja => this%dis%con%ja
197  !
198  ! -- Define shared data (cpw, rhow, latent heat of vaporization)
199  call this%gwecommon%gweshared_dat_df(this%neq)
200  !
201  ! -- Allocate model arrays, now that neq and nja are assigned
202  call this%allocate_arrays()
203  !
204  ! -- Return
205  return
Here is the call graph for this function: