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

Data Types

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

Functions/Subroutines

subroutine gwtifmod_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 gwtifmod_df (this)
 Define the GWT interface model. More...
 
subroutine gwtifmod_ar (this)
 Override allocate and read the GWT interface model and its packages so that we can create stuff from memory instead of input files. More...
 
subroutine gwtifmod_da (this)
 Clean up resources. More...
 

Function/Subroutine Documentation

◆ allocate_fmi()

subroutine gwtinterfacemodelmodule::allocate_fmi ( class(gwtinterfacemodeltype this)
private
Parameters
thisthe GWT interface model

Definition at line 117 of file GwtInterfaceModel.f90.

118  ! -- dummy
119  class(GwtInterfaceModelType) :: this !< the GWT interface model
120  !
121  call mem_allocate(this%fmi%gwfflowja, this%nja, 'GWFFLOWJA', &
122  this%fmi%memoryPath)
123  call mem_allocate(this%fmi%gwfhead, this%neq, 'GWFHEAD', &
124  this%fmi%memoryPath)
125  call mem_allocate(this%fmi%gwfsat, this%neq, 'GWFSAT', &
126  this%fmi%memoryPath)
127  call mem_allocate(this%fmi%gwfspdis, 3, this%neq, 'GWFSPDIS', &
128  this%fmi%memoryPath)
129  !
130  ! -- Return
131  return

◆ allocate_scalars()

subroutine gwtinterfacemodelmodule::allocate_scalars ( class(gwtinterfacemodeltype this,
character(len=*), intent(in)  modelname 
)
private
Parameters
thisthe GWT interface model
[in]modelnamethe model name

Definition at line 100 of file GwtInterfaceModel.f90.

101  ! -- dummy
102  class(GwtInterfaceModelType) :: this !< the GWT interface model
103  character(len=*), intent(in) :: modelname !< the model name
104  !
105  call this%GwtModelType%allocate_scalars(modelname)
106  !
107  call mem_allocate(this%iAdvScheme, 'ADVSCHEME', this%memoryPath)
108  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
109  call mem_allocate(this%ieqnsclfac, 'IEQNSCLFAC', this%memoryPath)
110  !
111  ! -- Return
112  return

◆ gwtifmod_ar()

subroutine gwtinterfacemodelmodule::gwtifmod_ar ( class(gwtinterfacemodeltype this)
private
Parameters
thisthe GWT interface model

Definition at line 197 of file GwtInterfaceModel.f90.

198  ! -- dummy
199  class(GwtInterfaceModelType) :: this !< the GWT interface model
200  !
201  call this%fmi%fmi_ar(this%ibound)
202  if (this%inadv > 0) then
203  call this%adv%adv_ar(this%dis, this%ibound)
204  end if
205  if (this%indsp > 0) then
206  call this%dsp%dsp_ar(this%ibound, this%mst%thetam)
207  end if
208  !
209  ! -- Return
210  return

◆ gwtifmod_cr()

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

Definition at line 49 of file GwtInterfaceModel.f90.

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)
93  !
94  ! -- Return
95  return
Here is the call graph for this function:

◆ gwtifmod_da()

subroutine gwtinterfacemodelmodule::gwtifmod_da ( class(gwtinterfacemodeltype this)
private
Parameters
thisthe GWT interface model

Definition at line 215 of file GwtInterfaceModel.f90.

216  ! -- dummy
217  class(GwtInterfaceModelType) :: this !< the GWT interface model
218  !
219 
220  ! this
221  call mem_deallocate(this%iAdvScheme)
222  call mem_deallocate(this%ixt3d)
223  call mem_deallocate(this%ieqnsclfac)
224  !
225  ! gwt packages
226  call this%dis%dis_da()
227  call this%fmi%fmi_da()
228  call this%adv%adv_da()
229  call this%dsp%dsp_da()
230  !
231  deallocate (this%dis)
232  deallocate (this%fmi)
233  deallocate (this%adv)
234  deallocate (this%dsp)
235  !
236  if (associated(this%mst)) then
237  call mem_deallocate(this%mst%thetam)
238  deallocate (this%mst)
239  end if
240  !
241  ! gwt scalars
242  call mem_deallocate(this%inic)
243  call mem_deallocate(this%infmi)
244  call mem_deallocate(this%inadv)
245  call mem_deallocate(this%indsp)
246  call mem_deallocate(this%inssm)
247  call mem_deallocate(this%inmst)
248  call mem_deallocate(this%inmvt)
249  call mem_deallocate(this%inoc)
250  call mem_deallocate(this%inobs)
251  call mem_deallocate(this%eqnsclfac)
252  !
253  ! base
254  call this%NumericalModelType%model_da()
255  !
256  ! -- Return
257  return

◆ gwtifmod_df()

subroutine gwtinterfacemodelmodule::gwtifmod_df ( class(gwtinterfacemodeltype this)
private
Parameters
thisthe GWT interface model

Definition at line 136 of file GwtInterfaceModel.f90.

137  ! -- dummy
138  class(GwtInterfaceModelType) :: this !< the GWT interface model
139  ! -- local
140  class(*), pointer :: disPtr
141  type(TspAdvOptionsType) :: adv_options
142  type(GwtDspOptionsType) :: dsp_options
143  !
144  this%moffset = 0
145  adv_options%iAdvScheme = this%iAdvScheme
146  dsp_options%ixt3d = this%ixt3d
147  !
148  ! define DISU
149  disptr => this%dis
150  call this%gridConnection%getDiscretization(castasdisutype(disptr))
151  call this%fmi%fmi_df(this%dis, 1)
152  !
153  if (this%inadv > 0) then
154  call this%adv%adv_df(adv_options)
155  end if
156  if (this%indsp > 0) then
157  this%dsp%idiffc = this%owner%dsp%idiffc
158  this%dsp%idisp = this%owner%dsp%idisp
159  call this%dsp%dsp_df(this%dis, dsp_options)
160  if (this%dsp%idiffc > 0) then
161  call mem_reallocate(this%dsp%diffc, this%dis%nodes, 'DIFFC', &
162  trim(this%dsp%memoryPath))
163  end if
164  if (this%dsp%idisp > 0) then
165  call mem_reallocate(this%dsp%alh, this%dis%nodes, 'ALH', &
166  trim(this%dsp%memoryPath))
167  call mem_reallocate(this%dsp%alv, this%dis%nodes, 'ALV', &
168  trim(this%dsp%memoryPath))
169  call mem_reallocate(this%dsp%ath1, this%dis%nodes, 'ATH1', &
170  trim(this%dsp%memoryPath))
171  call mem_reallocate(this%dsp%ath2, this%dis%nodes, 'ATH2', &
172  trim(this%dsp%memoryPath))
173  call mem_reallocate(this%dsp%atv, this%dis%nodes, 'ATV', &
174  trim(this%dsp%memoryPath))
175  end if
176  allocate (this%mst)
177  call mem_allocate(this%mst%thetam, this%dis%nodes, &
178  'THETAM', create_mem_path(this%name, 'MST'))
179  end if
180  !
181  ! assign or point model members to dis members
182  this%neq = this%dis%nodes
183  this%nja = this%dis%nja
184  this%ia => this%dis%con%ia
185  this%ja => this%dis%con%ja
186  !
187  ! allocate model arrays, now that neq and nja are assigned
188  call this%allocate_arrays()
189  !
190  ! -- Return
191  return
Here is the call graph for this function: