MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
VirtualGweModel.f90
Go to the documentation of this file.
2  use kindmodule, only: i4b
8  implicit none
9  private
10 
11  public :: add_virtual_gwe_model
12 
14  ! CND
15  type(virtualinttype), pointer :: cnd_idisp => null()
16  type(virtualdbl1dtype), pointer :: cnd_alh => null()
17  type(virtualdbl1dtype), pointer :: cnd_alv => null()
18  type(virtualdbl1dtype), pointer :: cnd_ath1 => null()
19  type(virtualdbl1dtype), pointer :: cnd_ath2 => null()
20  type(virtualdbl1dtype), pointer :: cnd_atv => null()
21  type(virtualdbl1dtype), pointer :: cnd_ktw => null()
22  type(virtualdbl1dtype), pointer :: cnd_kts => null()
23  ! FMI
24  type(virtualinttype), pointer :: fmi_igwfspdis => null()
25  type(virtualdbl1dtype), pointer :: fmi_gwfhead => null()
26  type(virtualdbl1dtype), pointer :: fmi_gwfsat => null()
27  type(virtualdbl2dtype), pointer :: fmi_gwfspdis => null()
28  type(virtualdbl1dtype), pointer :: fmi_gwfflowja => null()
29  ! EST
30  type(virtualdbl1dtype), pointer :: est_porosity => null()
31  ! GWE Model fields
32  type(virtualinttype), pointer :: incnd => null()
33  type(virtualinttype), pointer :: inest => null()
34  contains
35  ! public
36  procedure :: create => vgwe_create
37  procedure :: prepare_stage => vgwe_prepare_stage
38  procedure :: destroy => vgwe_destroy
39  ! private
40  procedure, private :: init_virtual_data
41  procedure, private :: allocate_data
42  procedure, private :: deallocate_data
43  end type virtualgwemodeltype
44 
45 contains
46 
47  subroutine add_virtual_gwe_model(model_id, model_name, model)
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(virtualgwemodeltype), pointer :: virtual_gwe_model
54  class(*), pointer :: obj
55 
56  allocate (virtual_gwe_model)
57  call virtual_gwe_model%create(model_name, model_id, model)
58 
59  obj => virtual_gwe_model
60  call virtual_model_list%Add(obj)
61 
62  end subroutine add_virtual_gwe_model
63 
64  subroutine vgwe_create(this, name, id, model)
65  class(virtualgwemodeltype) :: 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_gwemodel_type
73 
74  call this%allocate_data()
75  call this%init_virtual_data()
76 
77  end subroutine vgwe_create
78 
79  subroutine init_virtual_data(this)
80  class(virtualgwemodeltype) :: this
81 
82  call this%set(this%cnd_idisp%base(), 'IDISP', 'CND', map_all_type)
83  call this%set(this%cnd_alh%base(), 'ALH', 'CND', map_node_type)
84  call this%set(this%cnd_alv%base(), 'ALV', 'CND', map_node_type)
85  call this%set(this%cnd_ath1%base(), 'ATH1', 'CND', map_node_type)
86  call this%set(this%cnd_ath2%base(), 'ATH2', 'CND', map_node_type)
87  call this%set(this%cnd_atv%base(), 'ATV', 'CND', map_node_type)
88  call this%set(this%cnd_ktw%base(), 'KTW', 'CND', map_node_type)
89  call this%set(this%cnd_kts%base(), 'KTS', 'CND', 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%est_porosity%base(), 'POROSITY', 'EST', map_node_type)
96  call this%set(this%incnd%base(), 'INCND', '', map_all_type)
97  call this%set(this%inest%base(), 'INEST', '', map_all_type)
98 
99  end subroutine init_virtual_data
100 
101  subroutine vgwe_prepare_stage(this, stage)
102  class(virtualgwemodeltype) :: 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%cnd_idisp%base(), (/stg_aft_mdl_df/))
116  call this%map(this%incnd%base(), (/stg_aft_mdl_df/))
117  call this%map(this%inest%base(), (/stg_aft_mdl_df/))
118 
119  else if (stage == stg_bfr_con_ar) then
120 
121  nr_nodes = this%element_maps(map_node_type)%nr_virt_elems
122  nr_conns = this%element_maps(map_conn_type)%nr_virt_elems
123 
124  call this%map(this%x%base(), nr_nodes, &
126  call this%map(this%ibound%base(), nr_nodes, (/stg_bfr_con_ar/))
127 
128  call this%map(this%fmi_igwfspdis%base(), (/stg_bfr_con_ar/))
129 
130  if (this%cnd_idisp%get() > 0) then
131  call this%map(this%cnd_alh%base(), nr_nodes, (/stg_bfr_con_ar/))
132  call this%map(this%cnd_alv%base(), nr_nodes, (/stg_bfr_con_ar/))
133  call this%map(this%cnd_ath1%base(), nr_nodes, (/stg_bfr_con_ar/))
134  call this%map(this%cnd_ath2%base(), nr_nodes, (/stg_bfr_con_ar/))
135  call this%map(this%cnd_atv%base(), nr_nodes, (/stg_bfr_con_ar/))
136  call this%map(this%cnd_ktw%base(), nr_nodes, (/stg_bfr_con_ar/))
137  call this%map(this%cnd_kts%base(), nr_nodes, (/stg_bfr_con_ar/))
138  end if
139 
140  if (this%incnd%get() > 0 .and. this%inest%get() > 0) then
141  call this%map(this%est_porosity%base(), nr_nodes, (/stg_aft_con_ar/))
142  end if
143 
144  else if (stage == stg_aft_con_ar) then
145 
146  nr_nodes = this%element_maps(map_node_type)%nr_virt_elems
147  nr_conns = this%element_maps(map_conn_type)%nr_virt_elems
148 
149  call this%map(this%fmi_gwfhead%base(), nr_nodes, (/stg_bfr_exg_ad/))
150  call this%map(this%fmi_gwfsat%base(), nr_nodes, (/stg_bfr_exg_ad/))
151  if (this%fmi_igwfspdis%get() > 0) then
152  call this%map(this%fmi_gwfspdis%base(), 3, nr_nodes, (/stg_bfr_exg_ad/))
153  else
154  call this%map(this%fmi_gwfspdis%base(), 3, 0, (/stg_never/))
155  end if
156  call this%map(this%fmi_gwfflowja%base(), nr_conns, (/stg_bfr_exg_ad/))
157 
158  end if
159 
160  end subroutine vgwe_prepare_stage
161 
162  subroutine allocate_data(this)
163  class(virtualgwemodeltype) :: this
164 
165  allocate (this%cnd_idisp)
166  allocate (this%cnd_alh)
167  allocate (this%cnd_alv)
168  allocate (this%cnd_ath1)
169  allocate (this%cnd_ath2)
170  allocate (this%cnd_atv)
171  allocate (this%cnd_ktw)
172  allocate (this%cnd_kts)
173  allocate (this%fmi_igwfspdis)
174  allocate (this%fmi_gwfhead)
175  allocate (this%fmi_gwfsat)
176  allocate (this%fmi_gwfspdis)
177  allocate (this%fmi_gwfflowja)
178  allocate (this%est_porosity)
179  allocate (this%incnd)
180  allocate (this%inest)
181 
182  end subroutine allocate_data
183 
184  subroutine deallocate_data(this)
185  class(virtualgwemodeltype) :: this
186 
187  deallocate (this%cnd_idisp)
188  deallocate (this%cnd_alh)
189  deallocate (this%cnd_alv)
190  deallocate (this%cnd_ath1)
191  deallocate (this%cnd_ath2)
192  deallocate (this%cnd_atv)
193  deallocate (this%cnd_ktw)
194  deallocate (this%cnd_kts)
195  deallocate (this%fmi_igwfspdis)
196  deallocate (this%fmi_gwfhead)
197  deallocate (this%fmi_gwfsat)
198  deallocate (this%fmi_gwfspdis)
199  deallocate (this%fmi_gwfflowja)
200  deallocate (this%est_porosity)
201  deallocate (this%incnd)
202  deallocate (this%inest)
203 
204  end subroutine deallocate_data
205 
206  subroutine vgwe_destroy(this)
207  class(virtualgwemodeltype) :: this
208 
209  call this%VirtualModelType%destroy()
210  call this%deallocate_data()
211 
212  end subroutine vgwe_destroy
213 
214 end module virtualgwemodelmodule
This module defines variable data types.
Definition: kind.f90:8
integer(i4b), parameter, public stg_aft_con_ar
afterr connection allocate read
Definition: SimStages.f90:18
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_conn_type
Definition: VirtualBase.f90:15
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_gwemodel_type
type(listtype), public virtual_model_list
subroutine vgwe_prepare_stage(this, stage)
subroutine allocate_data(this)
subroutine vgwe_destroy(this)
subroutine, public add_virtual_gwe_model(model_id, model_name, model)
subroutine init_virtual_data(this)
subroutine vgwe_create(this, name, id, model)
subroutine deallocate_data(this)