FE-Project
Loading...
Searching...
No Matches
mod_experiment.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2#include "scaleFElib.h"
4 !-----------------------------------------------------------------------------
5 !
6 !++ Used modules
7 !
8 use scale_precision
9 use scale_io
10 use scale_prc
11
12 use scale_const, only: &
13 pi => const_pi, &
14 grav => const_grav, &
15 rdry => const_rdry, &
16 cpdry => const_cpdry, &
17 cvdry => const_cvdry, &
18 pres00 => const_pre00
19
31
32 !-----------------------------------------------------------------------------
33 implicit none
34 private
35 !-----------------------------------------------------------------------------
36 !
37 !++ Public type & procedures
38 !
39
40 type, public :: experiment
41 character(len=H_SHORT) :: label
42 procedure(exp_setinitcond_lc), pointer :: setinitcond_lc => null()
43 procedure(exp_geostrophic_balance_correction_lc), pointer :: geostrophic_balance_correction_lc => null()
44 contains
45 procedure, public :: init_base => experiment_init
46 generic :: init => init_base
47 procedure, public :: final_base => experiment_final
48 generic :: final => final_base
49 procedure, public :: setinitcond => experiment_setinitcond
50 procedure, public :: regist_setinitcond => experiment_regist_set_initcond
51 procedure, public :: regist_geostrophic_balance_correction => experiment_regist_geostrophic_balance_correction
52 end type experiment
53
55 class(localmeshfieldbase), pointer :: ptr
56 end type
57
58 interface
59 subroutine exp_setinitcond_lc( &
60 this, DENS_hyd, PRES_hyd, DDENS, MOMX, MOMY, MOMZ, DRHOT, &
61 tracer_field_list, &
62 x, y, z, dom_xmin, dom_xmax, dom_ymin, dom_ymax, dom_zmin, &
63 dom_zmax, lcmesh, elem )
64
65 import experiment
66 import localmesh3d
67 import elementbase3d
69 import rp
70
71 class(experiment), intent(inout) :: this
72 type(localmesh3d), intent(in) :: lcmesh
73 class(elementbase3d), intent(in) :: elem
74 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh%NeA)
75 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh%NeA)
76 real(RP), intent(out) :: DDENS(elem%Np,lcmesh%NeA)
77 real(RP), intent(out) :: MOMX(elem%Np,lcmesh%NeA)
78 real(RP), intent(out) :: MOMY(elem%Np,lcmesh%NeA)
79 real(RP), intent(out) :: MOMZ(elem%Np,lcmesh%NeA)
80 real(RP), intent(out) :: DRHOT(elem%Np,lcmesh%NeA)
81 type(tracerlocalmeshfield_ptr), intent(inout) :: tracer_field_list(:)
82 real(RP), intent(in) :: x(elem%Np,lcmesh%Ne)
83 real(RP), intent(in) :: y(elem%Np,lcmesh%Ne)
84 real(RP), intent(in) :: z(elem%Np,lcmesh%Ne)
85 real(RP), intent(in) :: dom_xmin, dom_xmax
86 real(RP), intent(in) :: dom_ymin, dom_ymax
87 real(RP), intent(in) :: dom_zmin, dom_zmax
88 end subroutine exp_setinitcond_lc
89
90 subroutine exp_geostrophic_balance_correction_lc( this, &
91 DENS_hyd, PRES_hyd, DDENS, MOMX, MOMY, MOMZ, DRHOT, &
92 lcmesh, elem )
93
94 import experiment
95 import localmesh3d
96 import elementbase3d
97 import rp
98
99 class(experiment), intent(inout) :: this
100 type(localmesh3d), intent(in) :: lcmesh
101 class(elementbase3d), intent(in) :: elem
102 real(RP), intent(inout) :: DENS_hyd(elem%Np,lcmesh%NeA)
103 real(RP), intent(in) :: PRES_hyd(elem%Np,lcmesh%NeA)
104 real(RP), intent(inout) :: DDENS(elem%Np,lcmesh%NeA)
105 real(RP), intent(inout) :: MOMX(elem%Np,lcmesh%NeA)
106 real(RP), intent(inout) :: MOMY(elem%Np,lcmesh%NeA)
107 real(RP), intent(inout) :: MOMZ(elem%Np,lcmesh%NeA)
108 real(RP), intent(inout) :: DRHOT(elem%Np,lcmesh%NeA)
109 end subroutine exp_geostrophic_balance_correction_lc
110 end interface
111
112 !-----------------------------------------------------------------------------
113 !
114 !++ Public parameters & variables
115 !
116
117 !-----------------------------------------------------------------------------
118 !
119 !++ Private procedures
120 !
121 !-------------------
122
123 !-----------------------------------------------------------------------------
124 !
125 !++ Private parameters & variables
126 !
127
128contains
129 subroutine experiment_init( this, exp_name )
130 implicit none
131 class(experiment), intent(inout) :: this
132 character(len=*), intent(in) :: exp_name
133 !----------------------------------------------------------------------
134
135 this%label = exp_name
136
137 this%setInitCond_lc => experiment_setinitcond_lc_dummy
138 this%geostrophic_balance_correction_lc => experiment_geostrophic_balance_correction_lc_dummy
139
140 return
141 end subroutine experiment_init
142
143 subroutine experiment_final( this )
144 implicit none
145 class(experiment), intent(inout) :: this
146 !----------------------------------------------------------------------
147
148 return
149 end subroutine experiment_final
150
151 subroutine experiment_regist_set_initcond( this, exp_SetInitCond_lc )
152 implicit none
153 class(experiment), intent(inout) :: this
154 interface
155 subroutine exp_setinitcond_lc( &
156 this, &
157 DENS_hyd, PRES_hyd, DDENS, MOMX, MOMY, MOMZ, DRHOT, &
158 tracer_field_list, &
159 x, y, z, dom_xmin, dom_xmax, dom_ymin, dom_ymax, dom_zmin, &
160 dom_zmax, lcmesh, elem )
161
162 import experiment
163 import localmesh3d
164 import elementbase3d
166 import rp
167
168 class(experiment), intent(inout) :: this
169 type(localmesh3d), intent(in) :: lcmesh
170 class(elementbase3d), intent(in) :: elem
171 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh%NeA)
172 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh%NeA)
173 real(RP), intent(out) :: DDENS(elem%Np,lcmesh%NeA)
174 real(RP), intent(out) :: MOMX(elem%Np,lcmesh%NeA)
175 real(RP), intent(out) :: MOMY(elem%Np,lcmesh%NeA)
176 real(RP), intent(out) :: MOMZ(elem%Np,lcmesh%NeA)
177 real(RP), intent(out) :: DRHOT(elem%Np,lcmesh%NeA)
178 type(tracerlocalmeshfield_ptr), intent(inout) :: tracer_field_list(:)
179 real(RP), intent(in) :: x(elem%Np,lcmesh%Ne)
180 real(RP), intent(in) :: y(elem%Np,lcmesh%Ne)
181 real(RP), intent(in) :: z(elem%Np,lcmesh%Ne)
182 real(RP), intent(in) :: dom_xmin, dom_xmax
183 real(RP), intent(in) :: dom_ymin, dom_ymax
184 real(RP), intent(in) :: dom_zmin, dom_zmax
185 end subroutine exp_setinitcond_lc
186 end interface
187 !----------------------------------------------------------------------
188
189 this%setInitCond_lc => exp_setinitcond_lc
190 return
191 end subroutine experiment_regist_set_initcond
192
193 subroutine experiment_regist_geostrophic_balance_correction( this, exp_geostrophic_balance_correction_lc )
194 implicit none
195 class(experiment), intent(inout) :: this
196 interface
197 subroutine exp_geostrophic_balance_correction_lc( this, &
198 DENS_hyd, PRES_hyd, DDENS, MOMX, MOMY, MOMZ, DRHOT, &
199 lcmesh, elem )
200
201 import experiment
202 import localmesh3d
203 import elementbase3d
204 import rp
205
206 class(experiment), intent(inout) :: this
207 type(localmesh3d), intent(in) :: lcmesh
208 class(elementbase3d), intent(in) :: elem
209 real(RP), intent(inout) :: DENS_hyd(elem%Np,lcmesh%NeA)
210 real(RP), intent(in) :: PRES_hyd(elem%Np,lcmesh%NeA)
211 real(RP), intent(inout) :: DDENS(elem%Np,lcmesh%NeA)
212 real(RP), intent(inout) :: MOMX(elem%Np,lcmesh%NeA)
213 real(RP), intent(inout) :: MOMY(elem%Np,lcmesh%NeA)
214 real(RP), intent(inout) :: MOMZ(elem%Np,lcmesh%NeA)
215 real(RP), intent(inout) :: DRHOT(elem%Np,lcmesh%NeA)
216 end subroutine exp_geostrophic_balance_correction_lc
217 end interface
218 !----------------------------------------------------------------------
219
220 this%geostrophic_balance_correction_lc => exp_geostrophic_balance_correction_lc
221 return
222 end subroutine experiment_regist_geostrophic_balance_correction
223
224 subroutine experiment_setinitcond( this, &
225 model_mesh, atm_prgvars_manager, atm_auxvars_manager, atm_trcvars_manager )
226
227 use scale_tracer, only: qa
228
232
236
237 use mod_atmos_vars, only: &
240 use mod_atmos_mesh, only: atmosmesh
241
242 implicit none
243
244 class(experiment), intent(inout) :: this
245 class(atmosmesh), target, intent(in) :: model_mesh
246 class(modelvarmanager), intent(inout) :: atm_prgvars_manager
247 class(modelvarmanager), intent(inout) :: atm_auxvars_manager
248 class(modelvarmanager), intent(inout) :: atm_trcvars_manager
249
250 class(localmeshfieldbase), pointer :: DDENS, MOMX, MOMY, MOMZ, DRHOT
251 class(localmeshfieldbase), pointer :: DENS_hyd, PRES_hyd
252 class(localmeshfieldbase), pointer :: Rtot, CVtot, CPtot
253
254 integer :: n
255 class(localmesh3d), pointer :: lcmesh3D
256 class(meshbase3d), pointer :: mesh
257
258 type(meshfieldcommcubedom3d), target :: hydvars_comm_rm
259 type(meshfieldcommcubedspheredom3d), target :: hydvars_comm_gm
260 class(meshfieldcommbase), pointer :: hydvars_comm
261 type(meshfieldcontainer) :: hydvars_comm_list(1)
262 class(meshfieldbase), pointer :: field_ptr
263
264 type(tracerlocalmeshfield_ptr) :: tracer_field_list(max(1,QA))
265 integer :: iq
266 !----------------------------------------------------------------------
267
268 mesh => model_mesh%ptr_mesh
269
270 do n=1, mesh%LOCAL_MESH_NUM
272 mesh, atm_prgvars_manager, atm_auxvars_manager, &
273 ddens, momx, momy, momz, drhot, &
274 dens_hyd, pres_hyd, rtot, cvtot, cptot, &
275 lcmesh3d )
276
277 do iq=1, qa
278 call atmosvars_getlocalmeshqtrcvar( n, mesh, atm_trcvars_manager, &
279 iq, tracer_field_list(iq)%ptr )
280 end do
281
282 select type (mesh)
283 type is (meshcubedom3d)
284 call this%setInitCond_lc( &
285 dens_hyd%val, pres_hyd%val, & ! (out)
286 ddens%val, momx%val, momy%val, momz%val, drhot%val, & ! (out)
287 tracer_field_list, & ! (inout)
288 lcmesh3d%pos_en(:,:,1), lcmesh3d%pos_en(:,:,2), lcmesh3d%pos_en(:,:,3), & ! (in)
289 mesh%xmin_gl, mesh%xmax_gl, mesh%ymin_gl, mesh%ymax_gl, mesh%zmin_gl, mesh%zmax_gl, & ! (in)
290 lcmesh3d, lcmesh3d%refElem3D ) ! (in)
291 type is (meshcubedspheredom3d)
292 call this%setInitCond_lc( &
293 dens_hyd%val, pres_hyd%val, & ! (out)
294 ddens%val, momx%val, momy%val, momz%val, drhot%val, & ! (out)
295 tracer_field_list, & ! (inout)
296 lcmesh3d%pos_en(:,:,1), lcmesh3d%pos_en(:,:,2), lcmesh3d%pos_en(:,:,3), & ! (in)
297 mesh%xmin_gl, mesh%xmax_gl, mesh%ymin_gl, mesh%ymax_gl, mesh%zmin_gl, mesh%zmax_gl, & ! (in)
298 lcmesh3d, lcmesh3d%refElem3D ) ! (in)
299 end select
300 end do
301
302 !------------------------------------------------------
303
304 call atm_auxvars_manager%Get(auxvar_preshydro_id, field_ptr)
305 select type(field_ptr)
306 type is (meshfield3d)
307 hydvars_comm_list(1)%field3d => field_ptr
308 end select
309 select type (mesh)
310 type is (meshcubedom3d)
311 call hydvars_comm_rm%Init(1, 0, 0, mesh)
312 hydvars_comm => hydvars_comm_rm
313 type is (meshcubedspheredom3d)
314 call hydvars_comm_gm%Init(1, 0, 0, mesh)
315 hydvars_comm => hydvars_comm_gm
316 end select
317
318 call hydvars_comm%Put(hydvars_comm_list, 1)
319 call hydvars_comm%Exchange()
320 call hydvars_comm%Get(hydvars_comm_list, 1)
321
322 !--------------------------------
323 do n=1, mesh%LOCAL_MESH_NUM
325 mesh, atm_prgvars_manager, atm_auxvars_manager, &
326 ddens, momx, momy, momz, drhot, &
327 dens_hyd, pres_hyd, rtot, cvtot, cptot, &
328 lcmesh3d )
329
330 call this%geostrophic_balance_correction_lc( &
331 dens_hyd%val, pres_hyd%val, & ! (out)
332 ddens%val, momx%val, momy%val, momz%val, drhot%val, & ! (out)
333 lcmesh3d, lcmesh3d%refElem3D ) ! (in)
334 end do
335
336 call atm_auxvars_manager%Get(auxvar_denshydro_id, field_ptr)
337 select type(field_ptr)
338 type is (meshfield3d)
339 hydvars_comm_list(1)%field3d => field_ptr
340 end select
341 call hydvars_comm%Put(hydvars_comm_list, 1)
342 call hydvars_comm%Exchange()
343 call hydvars_comm%Get(hydvars_comm_list, 1)
344
345 select type (mesh)
346 type is (meshcubedom3d)
347 call hydvars_comm_rm%Final()
348 type is (meshcubedspheredom3d)
349 call hydvars_comm_gm%Final()
350 end select
351
352 return
353 end subroutine experiment_setinitcond
354
355 !------
356!OCL SERIAL
357 subroutine experiment_setinitcond_lc_dummy( this, &
358 DENS_hyd, PRES_hyd, DDENS, MOMX, MOMY, MOMZ, DRHOT, tracer_field_list, &
359 x, y, z, dom_xmin, dom_xmax, dom_ymin, dom_ymax, dom_zmin, dom_zmax, &
360 lcmesh, elem )
361 implicit none
362
363 class(experiment), intent(inout) :: this
364 type(localmesh3d), intent(in) :: lcmesh
365 class(elementbase3d), intent(in) :: elem
366 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh%NeA)
367 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh%NeA)
368 real(RP), intent(out) :: DDENS(elem%Np,lcmesh%NeA)
369 real(RP), intent(out) :: MOMX(elem%Np,lcmesh%NeA)
370 real(RP), intent(out) :: MOMY(elem%Np,lcmesh%NeA)
371 real(RP), intent(out) :: MOMZ(elem%Np,lcmesh%NeA)
372 real(RP), intent(out) :: DRHOT(elem%Np,lcmesh%NeA)
373 type(tracerlocalmeshfield_ptr), intent(inout) :: tracer_field_list(:)
374 real(RP), intent(in) :: x(elem%Np,lcmesh%Ne)
375 real(RP), intent(in) :: y(elem%Np,lcmesh%Ne)
376 real(RP), intent(in) :: z(elem%Np,lcmesh%Ne)
377 real(RP), intent(in) :: dom_xmin, dom_xmax
378 real(RP), intent(in) :: dom_ymin, dom_ymax
379 real(RP), intent(in) :: dom_zmin, dom_zmax
380 !---------------------------------------------------
381
382 return
383 end subroutine experiment_setinitcond_lc_dummy
384
385 subroutine experiment_geostrophic_balance_correction_lc_dummy( this, &
386 DENS_hyd, PRES_hyd, DDENS, MOMX, MOMY, MOMZ, DRHOT, &
387 lcmesh, elem )
388
389 implicit none
390 class(experiment), intent(inout) :: this
391 type(localmesh3d), intent(in) :: lcmesh
392 class(elementbase3d), intent(in) :: elem
393 real(RP), intent(inout) :: DENS_hyd(elem%Np,lcmesh%NeA)
394 real(RP), intent(in) :: PRES_hyd(elem%Np,lcmesh%NeA)
395 real(RP), intent(inout) :: DDENS(elem%Np,lcmesh%NeA)
396 real(RP), intent(inout) :: MOMX(elem%Np,lcmesh%NeA)
397 real(RP), intent(inout) :: MOMY(elem%Np,lcmesh%NeA)
398 real(RP), intent(inout) :: MOMZ(elem%Np,lcmesh%NeA)
399 real(RP), intent(inout) :: DRHOT(elem%Np,lcmesh%NeA)
400 !---------------------------------------------------
401 return
402 end subroutine experiment_geostrophic_balance_correction_lc_dummy
403
404end module mod_experiment
module Atmosphere / Mesh
module ATMOSPHERE / Variables
subroutine, public atmosvars_getlocalmeshprgvars(domid, mesh, prgvars_list, auxvars_list, ddens, momx, momy, momz, therm, dens_hyd, pres_hyd, rtot, cvtot, cptot, lcmesh3d)
subroutine, public atmosvars_getlocalmeshqtrcvar(domid, mesh, trcvars_list, varid, var, lcmesh3d)
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / Common
module FElib / Element / Base
module FElib / Element / hexahedron
module FElib / Mesh / Local 3D
module FElib / Mesh / Base 3D
module FElib / Mesh / Cubic 3D domain
module FElib / Mesh / Cubed-sphere 3D domain
module FElib / Data / base
module FElib / Data / Communication base
module FElib / Data / Communication 3D cubic domain
module FElib / Data / Communication in 3D cubed-sphere domain
FElib / model framework / variable manager.
Base derived type to manage data communication.