MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
virtualgwfmodelmodule Module Reference

Data Types

type  virtualgwfmodeltype
 

Functions/Subroutines

subroutine, public add_virtual_gwf_model (model_id, model_name, model)
 Add virtual GWF model. More...
 
subroutine vgwf_create (this, name, id, model)
 
subroutine init_virtual_data (this)
 
subroutine vgwf_prepare_stage (this, stage)
 
subroutine vgwf_destroy (this)
 
subroutine allocate_data (this)
 
subroutine deallocate_data (this)
 

Function/Subroutine Documentation

◆ add_virtual_gwf_model()

subroutine, public virtualgwfmodelmodule::add_virtual_gwf_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 44 of file VirtualGwfModel.f90.

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

◆ allocate_data()

subroutine virtualgwfmodelmodule::allocate_data ( class(virtualgwfmodeltype this)
private

Definition at line 168 of file VirtualGwfModel.f90.

169  class(VirtualGwfModelType) :: this
170 
171  allocate (this%npf_iangle1)
172  allocate (this%npf_iangle2)
173  allocate (this%npf_iangle3)
174  allocate (this%npf_iwetdry)
175  allocate (this%inbuy)
176  allocate (this%npf_icelltype)
177  allocate (this%npf_k11)
178  allocate (this%npf_k22)
179  allocate (this%npf_k33)
180  allocate (this%npf_angle1)
181  allocate (this%npf_angle2)
182  allocate (this%npf_angle3)
183  allocate (this%npf_wetdry)
184  allocate (this%buy_dense)
185 

◆ deallocate_data()

subroutine virtualgwfmodelmodule::deallocate_data ( class(virtualgwfmodeltype this)
private

Definition at line 188 of file VirtualGwfModel.f90.

189  class(VirtualGwfModelType) :: this
190 
191  deallocate (this%npf_iangle1)
192  deallocate (this%npf_iangle2)
193  deallocate (this%npf_iangle3)
194  deallocate (this%npf_iwetdry)
195  deallocate (this%inbuy)
196  deallocate (this%npf_icelltype)
197  deallocate (this%npf_k11)
198  deallocate (this%npf_k22)
199  deallocate (this%npf_k33)
200  deallocate (this%npf_angle1)
201  deallocate (this%npf_angle2)
202  deallocate (this%npf_angle3)
203  deallocate (this%npf_wetdry)
204  deallocate (this%buy_dense)
205 

◆ init_virtual_data()

subroutine virtualgwfmodelmodule::init_virtual_data ( class(virtualgwfmodeltype this)
private

Definition at line 76 of file VirtualGwfModel.f90.

77  class(VirtualGwfModelType) :: this
78 
79  call this%set(this%npf_iangle1%base(), 'IANGLE1', 'NPF', map_all_type)
80  call this%set(this%npf_iangle2%base(), 'IANGLE2', 'NPF', map_all_type)
81  call this%set(this%npf_iangle3%base(), 'IANGLE3', 'NPF', map_all_type)
82  call this%set(this%npf_iwetdry%base(), 'IWETDRY', 'NPF', map_all_type)
83  call this%set(this%inbuy%base(), 'INBUY', '', map_all_type)
84  call this%set(this%npf_icelltype%base(), 'ICELLTYPE', 'NPF', map_node_type)
85  call this%set(this%npf_k11%base(), 'K11', 'NPF', map_node_type)
86  call this%set(this%npf_k22%base(), 'K22', 'NPF', map_node_type)
87  call this%set(this%npf_k33%base(), 'K33', 'NPF', map_node_type)
88  call this%set(this%npf_angle1%base(), 'ANGLE1', 'NPF', map_node_type)
89  call this%set(this%npf_angle2%base(), 'ANGLE2', 'NPF', map_node_type)
90  call this%set(this%npf_angle3%base(), 'ANGLE3', 'NPF', map_node_type)
91  call this%set(this%npf_wetdry%base(), 'WETDRY', 'NPF', map_node_type)
92  call this%set(this%buy_dense%base(), 'DENSE', 'BUY', map_node_type)
93 

◆ vgwf_create()

subroutine virtualgwfmodelmodule::vgwf_create ( class(virtualgwfmodeltype this,
character(len=*)  name,
integer(i4b)  id,
class(numericalmodeltype), pointer  model 
)

Definition at line 61 of file VirtualGwfModel.f90.

62  class(VirtualGwfModelType) :: this
63  character(len=*) :: name
64  integer(I4B) :: id
65  class(NumericalModelType), pointer :: model
66 
67  ! create base
68  call this%VirtualModelType%create(name, id, model)
69  this%container_type = vdc_gwfmodel_type
70 
71  call this%allocate_data()
72  call this%init_virtual_data()
73 

◆ vgwf_destroy()

subroutine virtualgwfmodelmodule::vgwf_destroy ( class(virtualgwfmodeltype this)
private

Definition at line 160 of file VirtualGwfModel.f90.

161  class(VirtualGwfModelType) :: this
162 
163  call this%VirtualModelType%destroy()
164  call this%deallocate_data()
165 

◆ vgwf_prepare_stage()

subroutine virtualgwfmodelmodule::vgwf_prepare_stage ( class(virtualgwfmodeltype this,
integer(i4b)  stage 
)
private

Definition at line 96 of file VirtualGwfModel.f90.

97  class(VirtualGwfModelType) :: this
98  integer(I4B) :: stage
99  ! local
100  integer(I4B) :: nr_nodes
101 
102  ! prepare base (=numerical) model data items
103  call this%VirtualModelType%prepare_stage(stage)
104 
105  if (stage == stg_aft_mdl_df) then
106 
107  call this%map(this%npf_iangle1%base(), (/stg_aft_mdl_df/))
108  call this%map(this%npf_iangle2%base(), (/stg_aft_mdl_df/))
109  call this%map(this%npf_iangle3%base(), (/stg_aft_mdl_df/))
110  call this%map(this%npf_iwetdry%base(), (/stg_aft_mdl_df/))
111  call this%map(this%inbuy%base(), (/stg_aft_mdl_df/))
112 
113  else if (stage == stg_bfr_con_ar) then
114 
115  nr_nodes = this%element_maps(map_node_type)%nr_virt_elems
116  ! Num. model data
117  call this%map(this%x%base(), nr_nodes, &
118  (/stg_bfr_con_ar, stg_bfr_exg_ad, stg_bfr_exg_cf/))
119  call this%map(this%ibound%base(), nr_nodes, &
120  (/stg_bfr_con_ar, stg_bfr_exg_ad, stg_bfr_exg_cf/))
121  call this%map(this%x_old%base(), nr_nodes, &
122  (/stg_bfr_exg_ad, stg_bfr_exg_cf/))
123  ! NPF
124  call this%map(this%npf_icelltype%base(), nr_nodes, (/stg_bfr_con_ar/))
125  call this%map(this%npf_k11%base(), nr_nodes, (/stg_bfr_con_ar/))
126  call this%map(this%npf_k22%base(), nr_nodes, (/stg_bfr_con_ar/))
127  call this%map(this%npf_k33%base(), nr_nodes, (/stg_bfr_con_ar/))
128 
129  if (this%npf_iangle1%get() > 0) then
130  call this%map(this%npf_angle1%base(), nr_nodes, (/stg_bfr_con_ar/))
131  else
132  call this%map(this%npf_angle1%base(), 0, (/stg_never/))
133  end if
134  if (this%npf_iangle2%get() > 0) then
135  call this%map(this%npf_angle2%base(), nr_nodes, (/stg_bfr_con_ar/))
136  else
137  call this%map(this%npf_angle2%base(), 0, (/stg_never/))
138  end if
139  if (this%npf_iangle3%get() > 0) then
140  call this%map(this%npf_angle3%base(), nr_nodes, (/stg_bfr_con_ar/))
141  else
142  call this%map(this%npf_angle3%base(), 0, (/stg_never/))
143  end if
144  if (this%npf_iwetdry%get() > 0) then
145  call this%map(this%npf_wetdry%base(), nr_nodes, (/stg_bfr_con_ar/))
146  else
147  call this%map(this%npf_wetdry%base(), 0, (/stg_never/))
148  end if
149 
150  if (this%inbuy%get() > 0) then
151  call this%map(this%buy_dense%base(), nr_nodes, (/stg_bfr_exg_cf/))
152  else
153  call this%map(this%buy_dense%base(), 0, (/stg_never/))
154  end if
155 
156  end if
157