MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
VirtualGwfModel.f90
Go to the documentation of this file.
2  use kindmodule, only: i4b
8  implicit none
9  private
10 
11  public :: add_virtual_gwf_model
12 
13  type, public, extends(virtualmodeltype) :: virtualgwfmodeltype
14  ! NPF
15  type(virtualinttype), pointer :: npf_iangle1 => null()
16  type(virtualinttype), pointer :: npf_iangle2 => null()
17  type(virtualinttype), pointer :: npf_iangle3 => null()
18  type(virtualinttype), pointer :: npf_iwetdry => null()
19  type(virtualinttype), pointer :: inbuy => null()
20  type(virtualint1dtype), pointer :: npf_icelltype => null()
21  type(virtualdbl1dtype), pointer :: npf_k11 => null()
22  type(virtualdbl1dtype), pointer :: npf_k22 => null()
23  type(virtualdbl1dtype), pointer :: npf_k33 => null()
24  type(virtualdbl1dtype), pointer :: npf_angle1 => null()
25  type(virtualdbl1dtype), pointer :: npf_angle2 => null()
26  type(virtualdbl1dtype), pointer :: npf_angle3 => null()
27  type(virtualdbl1dtype), pointer :: npf_wetdry => null()
28  type(virtualdbl1dtype), pointer :: buy_dense => null()
29  contains
30  ! public
31  procedure :: create => vgwf_create
32  procedure :: destroy => vgwf_destroy
33  procedure :: prepare_stage => vgwf_prepare_stage
34  ! private
35  procedure, private :: init_virtual_data
36  procedure, private :: allocate_data
37  procedure, private :: deallocate_data
38  end type virtualgwfmodeltype
39 
40 contains
41 
42  !> @brief Add virtual GWF model
43  !<
44  subroutine add_virtual_gwf_model(model_id, model_name, model)
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 
59  end subroutine add_virtual_gwf_model
60 
61  subroutine vgwf_create(this, name, id, model)
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 
74  end subroutine vgwf_create
75 
76  subroutine init_virtual_data(this)
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 
94  end subroutine init_virtual_data
95 
96  subroutine vgwf_prepare_stage(this, stage)
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, &
119  call this%map(this%ibound%base(), nr_nodes, &
121  call this%map(this%x_old%base(), nr_nodes, &
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 
158  end subroutine vgwf_prepare_stage
159 
160  subroutine vgwf_destroy(this)
161  class(virtualgwfmodeltype) :: this
162 
163  call this%VirtualModelType%destroy()
164  call this%deallocate_data()
165 
166  end subroutine vgwf_destroy
167 
168  subroutine allocate_data(this)
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 
186  end subroutine allocate_data
187 
188  subroutine deallocate_data(this)
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 
206  end subroutine deallocate_data
207 
208 end module virtualgwfmodelmodule
This module defines variable data types.
Definition: kind.f90:8
integer(i4b), parameter, public stg_aft_mdl_df
after model define
Definition: SimStages.f90:11
integer(i4b), parameter, public stg_bfr_exg_ad
before exchange advance (per solution)
Definition: SimStages.f90:21
integer(i4b), parameter, public stg_bfr_exg_cf
before exchange calculate (per solution)
Definition: SimStages.f90:22
integer(i4b), parameter, public stg_never
never
Definition: SimStages.f90:9
integer(i4b), parameter, public stg_bfr_con_ar
before connection allocate read
Definition: SimStages.f90:17
integer(i4b), parameter, public map_all_type
Definition: VirtualBase.f90:13
integer(i4b), parameter, public map_node_type
Definition: VirtualBase.f90:14
integer(i4b), parameter, public vdc_gwfmodel_type
type(listtype), public virtual_model_list
subroutine vgwf_prepare_stage(this, stage)
subroutine deallocate_data(this)
subroutine, public add_virtual_gwf_model(model_id, model_name, model)
Add virtual GWF model.
subroutine init_virtual_data(this)
subroutine vgwf_create(this, name, id, model)
subroutine allocate_data(this)
subroutine vgwf_destroy(this)