MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
virtualgwtmodelmodule Module Reference

Data Types

type  virtualgwtmodeltype
 

Functions/Subroutines

subroutine, public add_virtual_gwt_model (model_id, model_name, model)
 
subroutine vgwt_create (this, name, id, model)
 
subroutine init_virtual_data (this)
 
subroutine vgwt_prepare_stage (this, stage)
 
subroutine allocate_data (this)
 
subroutine deallocate_data (this)
 
subroutine vgwt_destroy (this)
 

Function/Subroutine Documentation

◆ add_virtual_gwt_model()

subroutine, public virtualgwtmodelmodule::add_virtual_gwt_model ( integer(i4b)  model_id,
character(len=*)  model_name,
class(numericalmodeltype), pointer  model 
)
Parameters
model_idglobal model id
model_namemodel name
modelthe actual model (can be null() when remote)

Definition at line 47 of file VirtualGwtModel.f90.

49  integer(I4B) :: model_id !< global model id
50  character(len=*) :: model_name !< model name
51  class(NumericalModelType), pointer :: model !< the actual model (can be null() when remote)
52  ! local
53  class(VirtualGwtModelType), pointer :: virtual_gwt_model
54  class(*), pointer :: obj
55 
56  allocate (virtual_gwt_model)
57  call virtual_gwt_model%create(model_name, model_id, model)
58 
59  obj => virtual_gwt_model
60  call virtual_model_list%Add(obj)
61 
type(listtype), public virtual_model_list
Here is the caller graph for this function:

◆ allocate_data()

subroutine virtualgwtmodelmodule::allocate_data ( class(virtualgwtmodeltype this)
private

Definition at line 164 of file VirtualGwtModel.f90.

165  class(VirtualGwtModelType) :: this
166 
167  allocate (this%dsp_idiffc)
168  allocate (this%dsp_idisp)
169  allocate (this%dsp_diffc)
170  allocate (this%dsp_alh)
171  allocate (this%dsp_alv)
172  allocate (this%dsp_ath1)
173  allocate (this%dsp_ath2)
174  allocate (this%dsp_atv)
175  allocate (this%fmi_igwfspdis)
176  allocate (this%fmi_gwfhead)
177  allocate (this%fmi_gwfsat)
178  allocate (this%fmi_gwfspdis)
179  allocate (this%fmi_gwfflowja)
180  allocate (this%mst_thetam)
181  allocate (this%indsp)
182  allocate (this%inmst)
183 

◆ deallocate_data()

subroutine virtualgwtmodelmodule::deallocate_data ( class(virtualgwtmodeltype this)
private

Definition at line 186 of file VirtualGwtModel.f90.

187  class(VirtualGwtModelType) :: this
188 
189  deallocate (this%dsp_idiffc)
190  deallocate (this%dsp_idisp)
191  deallocate (this%dsp_diffc)
192  deallocate (this%dsp_alh)
193  deallocate (this%dsp_alv)
194  deallocate (this%dsp_ath1)
195  deallocate (this%dsp_ath2)
196  deallocate (this%dsp_atv)
197  deallocate (this%fmi_igwfspdis)
198  deallocate (this%fmi_gwfhead)
199  deallocate (this%fmi_gwfsat)
200  deallocate (this%fmi_gwfspdis)
201  deallocate (this%fmi_gwfflowja)
202  deallocate (this%mst_thetam)
203  deallocate (this%indsp)
204  deallocate (this%inmst)
205 

◆ init_virtual_data()

subroutine virtualgwtmodelmodule::init_virtual_data ( class(virtualgwtmodeltype this)
private

Definition at line 79 of file VirtualGwtModel.f90.

80  class(VirtualGwtModelType) :: this
81 
82  call this%set(this%dsp_idiffc%base(), 'IDIFFC', 'DSP', map_all_type)
83  call this%set(this%dsp_idisp%base(), 'IDISP', 'DSP', map_all_type)
84  call this%set(this%dsp_diffc%base(), 'DIFFC', 'DSP', map_node_type)
85  call this%set(this%dsp_alh%base(), 'ALH', 'DSP', map_node_type)
86  call this%set(this%dsp_alv%base(), 'ALV', 'DSP', map_node_type)
87  call this%set(this%dsp_ath1%base(), 'ATH1', 'DSP', map_node_type)
88  call this%set(this%dsp_ath2%base(), 'ATH2', 'DSP', map_node_type)
89  call this%set(this%dsp_atv%base(), 'ATV', 'DSP', map_node_type)
90  call this%set(this%fmi_igwfspdis%base(), 'IGWFSPDIS', 'FMI', map_all_type)
91  call this%set(this%fmi_gwfhead%base(), 'GWFHEAD', 'FMI', map_node_type)
92  call this%set(this%fmi_gwfsat%base(), 'GWFSAT', 'FMI', map_node_type)
93  call this%set(this%fmi_gwfspdis%base(), 'GWFSPDIS', 'FMI', map_node_type)
94  call this%set(this%fmi_gwfflowja%base(), 'GWFFLOWJA', 'FMI', map_conn_type)
95  call this%set(this%mst_thetam%base(), 'THETAM', 'MST', map_node_type)
96  call this%set(this%indsp%base(), 'INDSP', '', map_all_type)
97  call this%set(this%inmst%base(), 'INMST', '', map_all_type)
98 

◆ vgwt_create()

subroutine virtualgwtmodelmodule::vgwt_create ( class(virtualgwtmodeltype this,
character(len=*)  name,
integer(i4b)  id,
class(numericalmodeltype), pointer  model 
)

Definition at line 64 of file VirtualGwtModel.f90.

65  class(VirtualGwtModelType) :: this
66  character(len=*) :: name
67  integer(I4B) :: id
68  class(NumericalModelType), pointer :: model
69 
70  ! create base
71  call this%VirtualModelType%create(name, id, model)
72  this%container_type = vdc_gwtmodel_type
73 
74  call this%allocate_data()
75  call this%init_virtual_data()
76 

◆ vgwt_destroy()

subroutine virtualgwtmodelmodule::vgwt_destroy ( class(virtualgwtmodeltype this)
private

Definition at line 208 of file VirtualGwtModel.f90.

209  class(VirtualGwtModelType) :: this
210 
211  call this%VirtualModelType%destroy()
212  call this%deallocate_data()
213 

◆ vgwt_prepare_stage()

subroutine virtualgwtmodelmodule::vgwt_prepare_stage ( class(virtualgwtmodeltype this,
integer(i4b)  stage 
)
private

Definition at line 101 of file VirtualGwtModel.f90.

102  class(VirtualGwtModelType) :: this
103  integer(I4B) :: stage
104  ! local
105  integer(I4B) :: nr_nodes, nr_conns
106 
107  ! prepare base (=numerical) model data items
108  call this%VirtualModelType%prepare_stage(stage)
109 
110  nr_nodes = 0
111  nr_conns = 0
112 
113  if (stage == stg_aft_mdl_df) then
114 
115  call this%map(this%dsp_idiffc%base(), (/stg_aft_mdl_df/))
116  call this%map(this%dsp_idisp%base(), (/stg_aft_mdl_df/))
117  call this%map(this%indsp%base(), (/stg_aft_mdl_df/))
118  call this%map(this%inmst%base(), (/stg_aft_mdl_df/))
119 
120  else if (stage == stg_bfr_con_ar) then
121 
122  nr_nodes = this%element_maps(map_node_type)%nr_virt_elems
123  nr_conns = this%element_maps(map_conn_type)%nr_virt_elems
124 
125  call this%map(this%fmi_igwfspdis%base(), (/stg_bfr_con_ar/))
126  call this%map(this%x%base(), nr_nodes, &
127  (/stg_bfr_con_ar, stg_bfr_exg_ad, stg_bfr_exg_cf/))
128  call this%map(this%ibound%base(), nr_nodes, (/stg_bfr_con_ar/))
129 
130  if (this%dsp_idiffc%get() > 0) then
131  call this%map(this%dsp_diffc%base(), nr_nodes, (/stg_bfr_con_ar/))
132  end if
133 
134  if (this%dsp_idisp%get() > 0) then
135  call this%map(this%dsp_alh%base(), nr_nodes, (/stg_bfr_con_ar/))
136  call this%map(this%dsp_alv%base(), nr_nodes, (/stg_bfr_con_ar/))
137  call this%map(this%dsp_ath1%base(), nr_nodes, (/stg_bfr_con_ar/))
138  call this%map(this%dsp_ath2%base(), nr_nodes, (/stg_bfr_con_ar/))
139  call this%map(this%dsp_atv%base(), nr_nodes, (/stg_bfr_con_ar/))
140  end if
141 
142  if (this%indsp%get() > 0 .and. this%inmst%get() > 0) then
143  call this%map(this%mst_thetam%base(), nr_nodes, (/stg_aft_con_ar/))
144  end if
145 
146  else if (stage == stg_aft_con_ar) then
147 
148  nr_nodes = this%element_maps(map_node_type)%nr_virt_elems
149  nr_conns = this%element_maps(map_conn_type)%nr_virt_elems
150 
151  call this%map(this%fmi_gwfhead%base(), nr_nodes, (/stg_bfr_exg_ad/))
152  call this%map(this%fmi_gwfsat%base(), nr_nodes, (/stg_bfr_exg_ad/))
153  if (this%fmi_igwfspdis%get() > 0) then
154  call this%map(this%fmi_gwfspdis%base(), 3, nr_nodes, (/stg_bfr_exg_ad/))
155  else
156  call this%map(this%fmi_gwfspdis%base(), 3, 0, (/stg_never/))
157  end if
158  call this%map(this%fmi_gwfflowja%base(), nr_conns, (/stg_bfr_exg_ad/))
159 
160  end if
161