FE-Project
Loading...
Searching...
No Matches
mod_atmos_vars.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
10!-------------------------------------------------------------------------------
11#include "scaleFElib.h"
13 !-----------------------------------------------------------------------------
14 !
15 !++ Used modules
16 !
17 use scale_precision
18 use scale_io
19 use scale_prc
20 use scale_debug
21 use scale_tracer, only: &
22 qa, tracer_name, tracer_desc, tracer_unit
23 use scale_atmos_hydrometeor, only: &
24 atmos_hydrometeor_dry
25
26 use scale_element_base, only: &
28 use scale_mesh_base, only: meshbase
29 use scale_mesh_base2d, only: &
31 use scale_mesh_base3d, only: &
32 meshbase3d, &
33 dimtype_xyz => meshbase3d_dimtypeid_xyz
34
38 use scale_localmeshfield_base, only: &
40 use scale_meshfield_base, only: &
42
45
47
48 use scale_model_var_manager, only: &
49 modelvarmanager, variableinfo
50
52 prgvar_num, auxvar_num, phytend_num1 => phytend_num, &
59
60 use mod_atmos_mesh, only: atmosmesh
61
62 !-----------------------------------------------------------------------------
63 implicit none
64 private
65
66 !-----------------------------------------------------------------------------
67 !
68 !++ Public type & procedures
69 !
70 type, public :: atmosvars
71 class(atmosmesh), pointer :: mesh
72
73 !- prognostic variables
74 type(meshfield3d), allocatable :: prog_vars(:)
75 type(modelvarmanager) :: progvars_manager
76 integer :: prog_vars_commid
77
78 !- tracer variables
79 type(meshfield3d), allocatable :: qtrc_vars(:)
80 type(modelvarmanager) :: qtrcvars_manager
81 integer :: qtrc_vars_commid
82
83 !- auxiliary variables
84 type(meshfield3d), allocatable :: aux_vars(:)
85 type(modelvarmanager) :: auxvars_manager
86 integer :: aux_vars_commid
87
88 !- auxiliary variables (2D)
89 type(meshfield2d), allocatable :: aux_vars2d(:)
90 type(modelvarmanager) :: auxvars2d_manager
91
92 !-
93 type(modelvarmanager), pointer :: ptr_mp_auxvars2d_manager
94
95 logical :: moist
96 type(meshfield3d), pointer :: qv
97 type(meshfield3d) :: zero
98
99 !- Tendency with physics
100 type(meshfield3d), allocatable :: phy_tend(:)
101 type(modelvarmanager) :: phytends_manager
102 integer :: phytends_commid
103 integer :: phytend_num_tot
104
105 !--
106 integer, allocatable :: diagvars2d_histid(:)
107 integer, allocatable :: diagvars3d_histid(:)
108
110
111 logical :: check_range
112 logical :: check_total
113
114 contains
115 procedure :: init => atmosvars_init
116 procedure :: final => atmosvars_final
117 procedure :: calc_diagnostics => atmosvars_calculatediagnostics
118 procedure :: calc_diagvar => atmosvars_calcdiagvar
119 procedure :: calc_diagvar2d => atmosvars_calcdiagvar2d
120 procedure :: history => atmosvars_history
121 procedure :: check => atmosvars_check
122 procedure :: monitor => atmosvars_monitor
123 procedure :: read_restart_file => atmosvar_read_restart_file
124 procedure :: write_restart_file => atmosvar_write_restart_file
125 procedure :: regist_physvar_manager => atmosvars_regist_physvar_manager
126 end type atmosvars
127
137
138 !-----------------------------------------------------------------------------
139 !
140 !++ Public parameters & variables
141 !
142 real(rp), parameter :: progvars_check_min(prgvar_num) = (/ -1.0_rp, -100.0_rp, -200.0_rp, -200.0_rp, -200.0_rp /)
143 real(rp), parameter :: progvars_check_max(prgvar_num) = (/ 1.0_rp, 100.0_rp, 200.0_rp, 200.0_rp, 200.0_rp /)
144
145
146 ! Surface variables
147
148 integer, public, parameter :: atmos_auxvars2d_prec_id = 1
149 integer, public, parameter :: atmos_auxvars2d_prec_engi_id = 2
150 integer, public, parameter :: atmos_auxvars2d_num = 2
151
152 type(variableinfo), public :: atmos_auxvars2d_vinfo(atmos_auxvars2d_num)
153 DATA atmos_auxvars2d_vinfo / &
154 variableinfo( atmos_auxvars2d_prec_id , 'PREC', 'surface precipitaion flux' , 'kg/m2/s', 2, 'XY', 'precipitation_flux' ), &
155 variableinfo( atmos_auxvars2d_prec_engi_id, 'PREC_ENGI', 'internal energy of precipitation' , 'J/m2', 2, 'XY', '' ) /
156
157
158
159 ! Diagnostic variables
160
161 integer, public, parameter :: atmos_diagvars_dens_id = 1
162 integer, public, parameter :: atmos_diagvars_u_id = 2
163 integer, public, parameter :: atmos_diagvars_v_id = 3
164 integer, public, parameter :: atmos_diagvars_w_id = 4
165 integer, public, parameter :: atmos_diagvars_t_id = 5
166 integer, public, parameter :: atmos_diagvars_umet_id = 6
167 integer, public, parameter :: atmos_diagvars_vmet_id = 7
168 integer, public, parameter :: atmos_diagvars_qdry_id = 8
169 integer, public, parameter :: atmos_diagvars_rh_id = 9
170 integer, public, parameter :: atmos_diagvars_engt_id = 10
171 integer, public, parameter :: atmos_diagvars_engp_id = 11
172 integer, public, parameter :: atmos_diagvars_engk_id = 12
173 integer, public, parameter :: atmos_diagvars_engi_id = 13
174 integer, public, parameter :: atmos_diagvars3d_num = 13
175
176 type(variableinfo), public :: atmos_diagvars3d_vinfo(atmos_diagvars3d_num)
178 variableinfo( atmos_diagvars_dens_id, 'DENS', 'density' , 'kg/m3', 3, 'XYZ', 'air_density' ), &
179 variableinfo( atmos_diagvars_u_id , 'U' , 'velocity u' , 'm/s' , 3, 'XYZ', 'x_wind' ), &
180 variableinfo( atmos_diagvars_v_id , 'V' , 'velocity v' , 'm/s' , 3, 'XYZ', 'y_wind' ), &
181 variableinfo( atmos_diagvars_w_id , 'W' , 'velocity w' , 'm/s' , 3, 'XYZ', 'upward_air_velocity' ), &
182 variableinfo( atmos_diagvars_t_id , 'T' , 'temperature' , 'K' , 3, 'XYZ', 'air_temperature' ), &
183 variableinfo( atmos_diagvars_umet_id, 'Umet', 'eastward velocity' , 'm/s' , 3, 'XYZ', 'x_wind' ), &
184 variableinfo( atmos_diagvars_vmet_id, 'Vmet', 'northward velocity' , 'm/s' , 3, 'XYZ', 'y_wind' ), &
185 variableinfo( atmos_diagvars_qdry_id, 'QDRY', 'dry air' , 'kg/kg', 3, 'XYZ', '' ), &
186 variableinfo( atmos_diagvars_rh_id , 'RH', 'relative humidity(liq)', '%', 3, 'XYZ', 'relative_humidity' ), &
187 variableinfo( atmos_diagvars_engt_id, 'ENGT', 'total energy' , 'J/m3' , 3, 'XYZ', '' ), &
188 variableinfo( atmos_diagvars_engp_id, 'ENGP', 'potential energy' , 'J/m3' , 3, 'XYZ', '' ), &
189 variableinfo( atmos_diagvars_engk_id, 'ENGK', 'kinetic energy' , 'J/m3' , 3, 'XYZ', '' ), &
190 variableinfo( atmos_diagvars_engi_id, 'ENGI', 'internal energy' , 'J/m3' , 3, 'XYZ', '' ) /
191
192
193 integer, public, parameter :: atmos_diagvars_rain_id = 1
194 integer, public, parameter :: atmos_diagvars_snow_id = 2
195 integer, public, parameter :: atmos_diagvars2d_num = 2
196
197 type(variableinfo), public :: atmos_diagvars2d_vinfo(atmos_diagvars2d_num)
199 variableinfo( atmos_diagvars_rain_id, 'RAIN', 'surface rain flux' , 'kg/m2/s', 2, 'XY', 'rainfall_flux' ), &
200 variableinfo( atmos_diagvars_snow_id, 'SNOW', 'surface snow flux' , 'kg/m2/s', 2, 'XY', 'snowfall_flux' ) /
201
202 !-----------------------------------------------------------------------------
203 !
204 !++ Private procedures & variables
205 !
206 !-------------------
207
208 private :: vars_calc_diagnosevar_lc
209 private :: vars_calc_diagnosevar2d_lc
210
211 ! for monitor
212
213 integer, private, parameter :: im_qdry = 1
214 integer, private, parameter :: im_qtot = 2
215 integer, private, parameter :: im_engt = 3
216 integer, private, parameter :: im_engp = 4
217 integer, private, parameter :: im_engk = 5
218 integer, private, parameter :: im_engi = 6
219 integer, private, parameter :: dvm_nmax = 6
220 integer, private :: dv_monit_id(dvm_nmax)
221
222contains
223
224!OCL SERIAL
225 subroutine atmosvars_init( this, atm_mesh )
226 use scale_atmos_hydrometeor, only: &
227 atmos_hydrometeor_dry
229 monitor_reg => file_monitor_meshfield_reg
230
234 implicit none
235
236 class(atmosvars), target, intent(inout) :: this
237 class(atmosmesh), target, intent(inout) :: atm_mesh
238
239 integer :: n
240 integer :: iv
241 integer :: iq
242 logical :: reg_file_hist
243
244 type(modelvarmanager) :: diagvar_manager ! dummy
245 type(meshfield2d) :: diag_vars2d(atmos_diagvars2d_num) ! dummy
246 type(meshfield3d) :: diag_vars3d(atmos_diagvars3d_num) ! dummy
247
248 logical :: check_range = .false.
249 logical :: check_total = .false.
250
251 namelist / param_atmos_vars / &
252 check_range, &
253 check_total
254
255 character(len=H_LONG) :: in_basename = ''
256 logical :: in_postfix_timelabel = .false.
257 character(len=H_LONG) :: out_basename = ''
258 logical :: out_postfix_timelabel = .true.
259 character(len=H_MID) :: out_title = ''
260 character(len=H_SHORT) :: out_dtype = 'DEFAULT'
261
262 namelist / param_atmos_vars_restart / &
263 in_basename, &
264 in_postfix_timelabel, &
265 out_basename, &
266 out_postfix_timelabel, &
267 out_title, &
268 out_dtype
269
270 integer :: ierr
271 logical :: is_specified
272
273 class(meshbase3d), pointer :: mesh3d
274 class(meshbase2d), pointer :: mesh2d
275
276 type(variableinfo) :: prgvar_info(prgvar_num)
277 !--------------------------------------------------
278
279 log_info('AtmosVars_Init',*)
280
281 !- read namelist
282 rewind(io_fid_conf)
283 read(io_fid_conf,nml=param_atmos_vars,iostat=ierr)
284 if( ierr < 0 ) then !--- missing
285 log_info("ATMOS_vars_setup",*) 'Not found namelist. Default used.'
286 elseif( ierr > 0 ) then !--- fatal error
287 log_error("ATMOS_vars_setup",*) 'Invalid names in namelist PARAM_ATMOS_VARS. Check!'
288 call prc_abort
289 endif
290 log_nml(param_atmos_vars)
291
292 !- Set the pointer of mesh
293 this%mesh => atm_mesh
294 mesh3d => atm_mesh%ptr_mesh
295 call mesh3d%GetMesh2D( mesh2d )
296
297 !- Initialize variables associated with dynamical core
298 ! (prognostic variables, tracer variables, 3D auxiliary variables, and tendencies of physical processes)
299
300 call this%PROGVARS_manager%Init()
301 call this%QTRCVARS_manager%Init()
302 call this%AUXVARS_manager%Init()
303 call this%PHYTENDS_manager%Init()
304
305 allocate( this%PROG_VARS(prgvar_num) )
306 allocate( this%QTRC_VARS(0:qa) )
307 allocate( this%AUX_VARS(auxvar_num) )
308
309 this%PHYTEND_NUM_TOT = phytend_num1 + max(1,qa)
310 allocate( this%PHY_TEND(this%PHYTEND_NUM_TOT) )
311
313 this%PROG_VARS, this%QTRC_VARS, this%AUX_VARS, this%PHY_TEND, & ! (inout)
314 this%PROGVARS_manager, this%QTRCVARS_manager, this%AUXVARS_manager, this%PHYTENDS_manager, & ! (inout)
315 this%PHYTEND_NUM_TOT, mesh3d, & ! (in)
316 prgvar_info ) ! (out)
317
318 ! Setup communicator
319
320 call atm_mesh%Create_communicator( &
322 this%PROGVARS_manager, & ! (inout)
323 this%PROG_VARS(:), & ! (in)
324 this%PROG_VARS_commID ) ! (out)
325
326 if ( .not. atmos_hydrometeor_dry ) then
327 call atm_mesh%Create_communicator( &
328 qa, 0, 0, & ! (in)
329 this%QTRCVARS_manager, & ! (inout)
330 this%QTRC_VARS(1:qa), & ! (in)
331 this%QTRC_VARS_commID ) ! (out)
332 end if
333
334 call atm_mesh%Create_communicator( &
335 auxvar_num, 0, 0, & ! (in)
336 this%AUXVARS_manager, & ! (inout)
337 this%AUX_VARS(:), & ! (in)
338 this%AUX_VARS_commID ) ! (out)
339
340 ! Output list of prognostic variables
341
342 log_newline
343 log_info("ATMOS_vars_setup",*) 'List of prognostic variables (ATMOS) '
344 log_info_cont('(1x,A,A24,A,A48,A,A12,A)') &
345 ' |', 'VARNAME ','|', &
346 'DESCRIPTION ', '[', 'UNIT ', ']'
347 do iv = 1, prgvar_num
348 log_info_cont('(1x,A,I3,A,A24,A,A48,A,A12,A)') &
349 'NO.',iv,'|',prgvar_info(iv)%NAME,'|', prgvar_info(iv)%DESC,'[', prgvar_info(iv)%UNIT,']'
350 end do
351 do iv = 1, qa
352 log_info_cont('(1x,A,I3,A,A24,A,A48,A,A12,A)') &
353 'NO.',prgvar_num+iv,'|',tracer_name(iv),'|', tracer_desc(iv),'[', tracer_unit(iv),']'
354 end do
355 log_newline
356
357
358 !- Initialize 2D auxiliary variables
359 call this%AUXVARS2D_manager%Init()
360 allocate( this%AUX_VARS2D(atmos_auxvars2d_num) )
361
362 reg_file_hist = .true.
363 do iv = 1, atmos_auxvars2d_num
364 call this%AUXVARS2D_manager%Regist( &
365 atmos_auxvars2d_vinfo(iv), mesh2d, & ! (in)
366 this%AUX_VARS2D(iv), & ! (inout)
367 reg_file_hist, fill_zero=.true. ) ! (in)
368 end do
369
370
371 !- Initialize diagnostic variables for output
372
373 ! 3D
374 call diagvar_manager%Init()
375 allocate( this%DIAGVARS3D_HISTID(atmos_diagvars3d_num) )
376
377 reg_file_hist = .true.
378 do iv = 1, atmos_diagvars3d_num
379 call diagvar_manager%Regist( &
380 atmos_diagvars3d_vinfo(iv), atm_mesh%ptr_mesh, & ! (in)
381 diag_vars3d(iv), reg_file_hist ) ! (out)
382
383 this%DIAGVARS3D_HISTID(iv) = diag_vars3d(iv)%hist_id
384 end do
385 call diagvar_manager%Final()
386
387 ! 2D
388 call diagvar_manager%Init()
389 allocate( this%DIAGVARS2D_HISTID(atmos_diagvars2d_num) )
390
391 reg_file_hist = .true.
392 do iv = 1, atmos_diagvars2d_num
393 call diagvar_manager%Regist( &
394 atmos_diagvars2d_vinfo(iv), mesh2d, & ! (in)
395 diag_vars2d(iv), reg_file_hist ) ! (out)
396
397 this%DIAGVARS2D_HISTID(iv) = diag_vars2d(iv)%hist_id
398 end do
399 call diagvar_manager%Final()
400
401 !-- Setup information for input/output restart files.
402
403 is_specified = .true.
404 !- read namelist
405 rewind(io_fid_conf)
406 read(io_fid_conf,nml=param_atmos_vars_restart,iostat=ierr)
407 if( ierr < 0 ) then !--- missing
408 log_info("AtmosVars_Init",*) 'Not found namelist PARAM_ATMOS_VARS_RESTART. Default used.'
409 is_specified = .false.
410 elseif( ierr > 0 ) then !--- fatal error
411 log_error("AtmosVars_Init",*) 'Not appropriate names in namelist PARAM_ATMOS_VARS_RESTART. Check!'
412 call prc_abort
413 endif
414 log_nml(param_atmos_vars_restart)
415
416 if (is_specified) then
417 call atm_mesh%Setup_restartfile( this%restart_file, &
418 in_basename, in_postfix_timelabel, &
419 out_basename, out_postfix_timelabel, out_dtype, out_title, &
421 else
422 call atm_mesh%Setup_restartfile( this%restart_file, &
424 end if
425
426 !-----< monitor output setup >-----
427
428 call monitor_reg( 'QTOT', 'water mass', 'kg', & ! (in)
429 dv_monit_id(im_qtot), & ! (out)
430 dim_type='ATM3D', is_tendency=.false. ) ! (in)
431
432 call monitor_reg( 'ENGT', 'total energy', 'J', & ! (in)
433 dv_monit_id(im_engt), & ! (out)
434 dim_type='ATM3D', is_tendency=.false. ) ! (in)
435 call monitor_reg( 'ENGP', 'potential energy', 'J', & ! (in)
436 dv_monit_id(im_engp), & ! (out)
437 dim_type='ATM3D', is_tendency=.false. ) ! (in)
438 call monitor_reg( 'ENGK', 'kinetic energy', 'J', & ! (in)
439 dv_monit_id(im_engk), & ! (out)
440 dim_type='ATM3D', is_tendency=.false. ) ! (in)
441 call monitor_reg( 'ENGI', 'internal energy', 'J', & ! (in)
442 dv_monit_id(im_engi), & ! (out)
443 dim_type='ATM3D', is_tendency=.false. ) ! (in)
444
445
446 !-----< check the range of values >-----
447
448 this%check_range = check_range
449 this%check_total = check_total
450 log_info("ATMOS_vars_setup",*) 'Check value range of variables? : ', check_range
451 log_info("ATMOS_vars_setup",*) 'Check total value of variables? : ', check_total
452
453 return
454 end subroutine atmosvars_init
455
456!OCL SERIAL
457 subroutine atmosvars_final( this )
458 implicit none
459 class(atmosvars), intent(inout) :: this
460
461 !--------------------------------------------------
462
463 log_info('AtmosVars_Final',*)
464
465 call this%restart_file%Final()
466
467 call this%PROGVARS_manager%Final()
468 deallocate( this%PROG_VARS )
469
470 call this%QTRCVARS_manager%Final()
471 deallocate( this%QTRC_VARS )
472
473 call this%AUXVARS_manager%Final()
474 deallocate( this%AUX_VARS )
475
476 call this%AUXVARS2D_manager%Final()
477 deallocate( this%AUX_VARS2D )
478
479 call this%PHYTENDS_manager%Final()
480 deallocate( this%PHY_TEND )
481
482 deallocate( this%DIAGVARS3D_HISTID )
483
484 return
485 end subroutine atmosvars_final
486
487!OCL SERIAL
488 subroutine atmosvars_regist_physvar_manager( this, &
489 mp_AUXVARS2D_manager )
490 implicit none
491
492 class(atmosvars), target, intent(inout) :: this
493 type(modelvarmanager), intent(in), target :: mp_auxvars2d_manager
494 !----------------------------------------------
495
496 this%ptr_MP_AUXVARS2D_manager => mp_auxvars2d_manager
497
498 return
499 end subroutine atmosvars_regist_physvar_manager
500
501!OCL SERIAL
502 subroutine atmosvars_history( this )
504 implicit none
505 class(atmosvars), intent(inout) :: this
506
507 integer :: v
508 integer :: hst_id
509
510 type(meshfield3d) :: tmp_field3d
511 class(meshbase3d), pointer :: mesh3d
512
513 type(meshfield2d) :: tmp_field2d
514 class(meshbase2d), pointer :: mesh2d
515 !-------------------------------------------------------------------------
516
517 mesh3d => this%PROG_VARS(1)%mesh
518 call mesh3d%GetMesh2D(mesh2d)
519
520 do v = 1, prgvar_num
521 hst_id = this%PROG_VARS(v)%hist_id
522 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%PROG_VARS(v) )
523 end do
524
525 do v = 1, qa
526 hst_id = this%QTRC_VARS(v)%hist_id
527 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%QTRC_VARS(v) )
528 end do
529
530 call this%Calc_diagnostics()
531
532 do v = 1, auxvar_num
533 hst_id = this%AUX_VARS(v)%hist_id
534 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%AUX_VARS(v) )
535 end do
536 do v = 1, atmos_auxvars2d_num
537 hst_id = this%AUX_VARS2D(v)%hist_id
538 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%AUX_VARS2D(v) )
539 end do
540
541 do v = 1, this%PHYTEND_NUM_TOT
542 hst_id = this%PHY_TEND(v)%hist_id
543 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%PHY_TEND(v) )
544 end do
545
546 !- Output diagnostic variables
547
548 ! 3D
549 call tmp_field3d%Init( "tmp_field", "", mesh3d)
550 do v = 1, atmos_diagvars3d_num
551 hst_id = this%DIAGVARS3D_HISTID(v)
552 if ( hst_id > 0 ) then
553 call atmosvars_calcdiagvar( this, atmos_diagvars3d_vinfo(v)%NAME, tmp_field3d )
554 call file_history_meshfield_put( hst_id, tmp_field3d )
555 end if
556 end do
557 call tmp_field3d%Final()
558
559 ! 2D
560 call tmp_field2d%Init( "tmp_field", "", mesh2d)
561 do v = 1, atmos_diagvars2d_num
562 hst_id = this%DIAGVARS2D_HISTID(v)
563 if ( hst_id > 0 ) then
564 call atmosvars_calcdiagvar2d( this, atmos_diagvars2d_vinfo(v)%NAME, tmp_field2d )
565 call file_history_meshfield_put( hst_id, tmp_field2d )
566 end if
567 end do
568 call tmp_field2d%Final()
569
570 return
571 end subroutine atmosvars_history
572
573!OCL SERIAL
574 subroutine atmosvar_read_restart_file( this, atmos_mesh, dyncore )
575
579
582 implicit none
583
584 class(atmosvars), intent(inout), target :: this
585 class(atmosmesh), intent(in) :: atmos_mesh
586 class(atmdyndgmdriver_nonhydro3d), intent(inout) :: dyncore
587
588 integer :: iv
589
590 integer :: domid
591 integer :: ke
592 class(meshbase3d), pointer :: mesh3d
593 class(localmesh3d), pointer :: lcmesh3d
594 class(meshfield3d), pointer :: phyd_ref
595 !---------------------------------------
596
597 log_newline
598 log_info("ATMOSVar_read_restart_file",*) 'Open restart file (ATMOS) '
599
600 !- Open restart file
601 call this%restart_file%Open()
602
603 !- Read restart file
604
605 do iv=1, prgvar_num
606 call this%restart_file%Read_var( dimtype_xyz, this%PROG_VARS(iv)%varname, &
607 this%PROG_VARS(iv) )
608 end do
609 do iv=1, auxvar_denshydro_id
610 call this%restart_file%Read_var( dimtype_xyz, this%AUX_VARS(iv)%varname, &
611 this%AUX_VARS(iv) )
612 end do
613 do iv=1, qa
614 call this%restart_file%Read_var( dimtype_xyz, this%QTRC_VARS(iv)%varname, &
615 this%QTRC_VARS(iv) )
616 end do
617
618 !- Close restart file
619 log_info("ATMOSVar_read_restart_file",*) 'Close restart file (ATMOS) '
620 call this%restart_file%Close()
621
622 !-- Diagnostic pressure
623 call vars_calc_specific_heat( this )
624
625 call dyncore%calc_pressure( this%AUX_VARS(auxvar_pres_id), &
626 this%PROGVARS_manager, this%AUXVARS_manager )
627
628 ! Set reference value of hydrostatic pressure
629 phyd_ref => this%AUX_VARS(auxvar_preshydro_ref_id)
630 mesh3d => phyd_ref%mesh
631 do domid=1, mesh3d%LOCAL_MESH_NUM
632 lcmesh3d => mesh3d%lcmesh_list(domid)
633 !$omp parallel do
634 do ke=lcmesh3d%NeS, lcmesh3d%NeE
635 phyd_ref%local(domid)%val(:,ke) = 0.0_rp
636 end do
637 end do
638
639 !-- Check read data
640 call this%Check( force = .true. )
641
642 !-- Calculate diagnostic variables
643 call this%Calc_diagnostics()
644
645 !-- Communicate halo data of hydrostatic & diagnostic variables
646 call this%AUXVARS_manager%MeshFieldComm_Exchange()
647
648 !-- Set horizontal gradient of hydrostatic pressure
649 call dyncore%update_phyd_hgrad( this%AUX_VARS(auxvar_preshydro_id), phyd_ref, &
650 mesh3d, atmos_mesh%element3D_operation )
651
652 return
653 end subroutine atmosvar_read_restart_file
654
655!OCL SERIAL
656 subroutine atmosvar_write_restart_file( this )
659
660 implicit none
661 class(atmosvars), intent(inout) :: this
662
663 type(variableinfo) :: prgvar_info(prgvar_num)
664 type(variableinfo) :: auxvar_info(auxvar_num)
665
666 integer :: iv, rf_vid
667 !---------------------------------------
668
669 log_newline
670 log_info("ATMOSVar_write_restart_file",*) 'Create restart file (ATMOS) '
671
672 !- Check data which will be written to restart file
673 call this%Check( force = .true. )
674
675 !- Create restart file
676 call this%restart_file%Create()
677 call prc_mpibarrier()
678
679 !- Define variables
680
681 call atm_dyn_dgm_nonhydro3d_common_get_varinfo( prgvar_info, auxvar_info )
682
683 do iv=1, prgvar_num
684 rf_vid = iv
685 call this%restart_file%Def_var( this%PROG_VARS(iv), &
686 prgvar_info(iv)%DESC, rf_vid, dimtype_xyz )
687 end do
688 do iv=1, auxvar_denshydro_id
689 rf_vid = prgvar_num + iv
690 call this%restart_file%Def_var( this%AUX_VARS(iv), &
691 auxvar_info(iv)%DESC, rf_vid, dimtype_xyz )
692 end do
693 do iv=1, qa
694 rf_vid = rf_vid + 1
695 call this%restart_file%Def_var( this%QTRC_VARS(iv), &
696 tracer_desc(iv), rf_vid, dimtype_xyz )
697 end do
698
699 call this%restart_file%End_def()
700
701 !- Write restart file
702 do iv=1, prgvar_num
703 rf_vid = iv
704 call this%restart_file%Write_var(rf_vid, this%PROG_VARS(iv) )
705 end do
706 do iv=1, auxvar_denshydro_id
707 rf_vid = prgvar_num + iv
708 call this%restart_file%Write_var(rf_vid, this%AUX_VARS(iv) )
709 end do
710 do iv=1, qa
711 rf_vid = rf_vid + 1
712 call this%restart_file%Write_var(rf_vid, this%QTRC_VARS(iv) )
713 end do
714
715 !- Close restart file
716 log_info("ATMOSVar_write_restart_file",*) 'Close restart file (ATMOS) '
717 call this%restart_file%Close()
718
719 return
720 end subroutine atmosvar_write_restart_file
721
722!OCL SERIAL
723 subroutine atmosvars_check( this, force )
724
725 use scale_meshfield_statistics, only: &
728
729 implicit none
730 class(atmosvars), intent(inout) :: this
731 logical, intent(in), optional :: force
732
733 integer :: iv
734 integer :: iv_diag
735 integer :: n
736 logical :: check
737
738 class(meshbase3d), pointer :: mesh3d
739 class(localmeshbase), pointer :: lcmesh
740 class(localmeshfieldbase), pointer :: lcfield
741 type(elementbase), pointer :: elem
742 character(len=H_MID) :: varname
743
744 type(meshfield3d) :: vel_fields(3)
745 type(meshfield3d) :: work
746 !--------------------------------------------------------------------------
747
748 if ( present(force) ) then
749 check = force
750 else
751 check = this%check_range
752 end if
753
754 if (check) then
755 do iv=1, prgvar_num
756 if ( iv == prgvar_therm_id ) cycle
757
758 mesh3d => this%PROG_VARS(iv)%mesh
759 do n=1, mesh3d%LOCAL_MESH_NUM
760 lcmesh => mesh3d%lcmesh_list(n)
761 elem => lcmesh%refElem
762 call this%PROG_VARS(iv)%GetLocalMeshField(n, lcfield)
763 write(varname,'(a,i3.3,a)') this%PROG_VARS(iv)%varname//'(domID=', n, ')'
764 call valcheck( elem%Np, 1, elem%Np, lcmesh%NeA, lcmesh%NeS, lcmesh%NeE, lcfield%val(:,:), &
765 progvars_check_min(iv), progvars_check_max(iv), trim(varname), __file__, __line__ )
766 end do
767 end do
768
769 mesh3d => this%PROG_VARS(1)%mesh
770 do iv=1, 3
771 iv_diag = atmos_diagvars_u_id + iv - 1
772 call vel_fields(iv)%Init( atmos_diagvars3d_vinfo(iv_diag)%NAME, "", mesh3d )
773 call atmosvars_calcdiagvar( this, vel_fields(iv)%varname, vel_fields(iv) )
774 end do
775 call meshfield_statistics_detail( vel_fields(:) )
776 do iv=1, 3
777 call vel_fields(iv)%Final()
778 end do
779 end if
780
781 if ( present(force) ) then
782 check = force
783 else
784 check = this%check_total
785 end if
786 if (check) then
787 mesh3d => this%PROG_VARS(1)%mesh
788 call work%Init("tmp", "", mesh3d)
789 call work%Final()
790 end if
791
792 return
793 end subroutine atmosvars_check
794
795!OCL SERIAL
796 subroutine atmosvars_monitor( this )
799
800 implicit none
801 class(atmosvars), intent(inout) :: this
802
803 integer :: iv
804 class(meshbase3d), pointer :: mesh3d
805 type(meshfield3d) :: work
806
807 integer :: n
808 integer :: ke
809 class(localmesh3d), pointer :: lcmesh
810 !--------------------------------------------------------------------------
811
812 mesh3d => this%PROG_VARS(1)%mesh
813 call work%Init("tmp", "", mesh3d)
814
815 do iv=1, prgvar_num
816 call file_monitor_meshfield_put( this%PROG_VARS(iv)%monitor_id, this%PROG_VARS(iv) )
817 end do
818
819 do iv=1, qa
820 if ( this%QTRC_VARS(iv)%monitor_id > 0 ) then
821 do n=1, mesh3d%LOCAL_MESH_NUM
822 lcmesh => mesh3d%lcmesh_list(n)
823 !$omp parallel do
824 do ke=lcmesh%NeS, lcmesh%NeE
825 work%local(n)%val(:,ke) = ( this%AUX_VARS(auxvar_denshydro_id)%local(n)%val(:,ke) &
826 + this%PROG_VARS(prgvar_ddens_id )%local(n)%val(:,ke) &
827 ) * this%QTRC_VARS(iv)%local(n)%val(:,ke)
828 end do
829 end do
830 call file_monitor_meshfield_put( this%QTRC_VARS(iv)%monitor_id, work )
831 end if
832 end do
833 if ( dv_monit_id(im_qtot) > 0 ) then
834 do n=1, mesh3d%LOCAL_MESH_NUM
835 lcmesh => mesh3d%lcmesh_list(n)
836 !$omp parallel do
837 do ke=lcmesh%NeS, lcmesh%NeE
838 work%local(n)%val(:,ke) = ( this%AUX_VARS(auxvar_denshydro_id)%local(n)%val(:,ke) &
839 + this%PROG_VARS(prgvar_ddens_id )%local(n)%val(:,ke) &
840 ) * ( 1.0_rp - this%AUX_VARS(auxvar_qdry_id)%local(n)%val(:,ke) )
841 end do
842 end do
843 call file_monitor_meshfield_put( dv_monit_id(im_qtot), work )
844 end if
845
846 !##### Energy Budget #####
847
848 if ( dv_monit_id(im_engt) > 0 ) then
849 call atmosvars_calcdiagvar( this, 'ENGT', work )
850 call file_monitor_meshfield_put( dv_monit_id(im_engt), work )
851 end if
852 if ( dv_monit_id(im_engp) > 0 ) then
853 call atmosvars_calcdiagvar( this, 'ENGP', work )
854 call file_monitor_meshfield_put( dv_monit_id(im_engp), work )
855 end if
856 if ( dv_monit_id(im_engk) > 0 ) then
857 call atmosvars_calcdiagvar( this, 'ENGK', work )
858 call file_monitor_meshfield_put( dv_monit_id(im_engk), work )
859 end if
860 if ( dv_monit_id(im_engi) > 0 ) then
861 call atmosvars_calcdiagvar( this, 'ENGI', work )
862 call file_monitor_meshfield_put( dv_monit_id(im_engi), work )
863 end if
864
865 call work%Final()
866
867 return
868 end subroutine atmosvars_monitor
869
870 !---- Getter ---------------------------------------------------------------------------
871
872!OCL SERIAL
873 subroutine atmosvars_getlocalmeshprgvar( domID, mesh, prgvars_list, auxvars_list, &
874 varid, &
875 var, DENS_hyd, PRES_hyd, lcmesh3D )
876
877 implicit none
878
879 integer, intent(in) :: domid
880 class(meshbase), intent(in) :: mesh
881 class(modelvarmanager), intent(inout) :: prgvars_list
882 class(modelvarmanager), intent(inout) :: auxvars_list
883 integer, intent(in) :: varid
884 class(localmeshfieldbase), pointer, intent(out) :: var
885 class(localmeshfieldbase), pointer, intent(out), optional :: dens_hyd, pres_hyd
886 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
887
888 class(meshfieldbase), pointer :: field
889 class(localmeshbase), pointer :: lcmesh
890 !-------------------------------------------------------
891
892 !--
893 call prgvars_list%Get(varid, field)
894 call field%GetLocalMeshField(domid, var)
895
896 if (present(dens_hyd)) then
897 call auxvars_list%Get(auxvar_denshydro_id, field)
898 call field%GetLocalMeshField(domid, dens_hyd)
899 end if
900 if (present(pres_hyd)) then
901 call auxvars_list%Get(auxvar_preshydro_id, field)
902 call field%GetLocalMeshField(domid, pres_hyd)
903 end if
904
905 if (present(lcmesh3d)) then
906 call mesh%GetLocalMesh( domid, lcmesh )
907 nullify( lcmesh3d )
908
909 select type(lcmesh)
910 type is (localmesh3d)
911 if (present(lcmesh3d)) lcmesh3d => lcmesh
912 end select
913 end if
914
915 return
916 end subroutine atmosvars_getlocalmeshprgvar
917
918!OCL SERIAL
919 subroutine atmosvars_getlocalmeshprgvars( domID, mesh, prgvars_list, auxvars_list, &
920 DDENS, MOMX, MOMY, MOMZ, THERM, &
921 DENS_hyd, PRES_hyd, Rtot, CVtot, CPtot, &
922 lcmesh3D )
923 implicit none
924 integer, intent(in) :: domid
925 class(meshbase), intent(in) :: mesh
926 class(modelvarmanager), intent(inout) :: prgvars_list
927 class(modelvarmanager), intent(inout) :: auxvars_list
928 class(localmeshfieldbase), pointer, intent(out) :: ddens, momx, momy, momz, therm
929 class(localmeshfieldbase), pointer, intent(out) :: dens_hyd, pres_hyd
930 class(localmeshfieldbase), pointer, intent(out) :: rtot, cvtot, cptot
931 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
932
933 class(meshfieldbase), pointer :: field
934 class(localmeshbase), pointer :: lcmesh
935 !-------------------------------------------------------
936
937 !--
938 call prgvars_list%Get(prgvar_ddens_id, field)
939 call field%GetLocalMeshField(domid, ddens)
940
941 call prgvars_list%Get(prgvar_momx_id, field)
942 call field%GetLocalMeshField(domid, momx)
943
944 call prgvars_list%Get(prgvar_momy_id, field)
945 call field%GetLocalMeshField(domid, momy)
946
947 call prgvars_list%Get(prgvar_momz_id, field)
948 call field%GetLocalMeshField(domid, momz)
949
950 call prgvars_list%Get(prgvar_therm_id, field)
951 call field%GetLocalMeshField(domid, therm)
952
953 !--
954 call auxvars_list%Get(auxvar_denshydro_id, field)
955 call field%GetLocalMeshField(domid, dens_hyd)
956
957 call auxvars_list%Get(auxvar_preshydro_id, field)
958 call field%GetLocalMeshField(domid, pres_hyd)
959
960 call auxvars_list%Get(auxvar_rtot_id, field)
961 call field%GetLocalMeshField(domid, rtot)
962
963 call auxvars_list%Get(auxvar_cvtot_id, field)
964 call field%GetLocalMeshField(domid, cvtot)
965
966 call auxvars_list%Get(auxvar_cptot_id, field)
967 call field%GetLocalMeshField(domid, cptot)
968
969 !---
970
971 if ( present(lcmesh3d) ) then
972 call mesh%GetLocalMesh( domid, lcmesh )
973 nullify( lcmesh3d )
974
975 select type(lcmesh)
976 type is (localmesh3d)
977 if (present(lcmesh3d)) lcmesh3d => lcmesh
978 end select
979 end if
980
981 return
982 end subroutine atmosvars_getlocalmeshprgvars
983
984!OCL SERIAL
985 subroutine atmosvars_getlocalmeshsfcvar( domID, mesh, auxvars2D_list, &
986 PREC, PREC_ENGI, lcmesh2D )
987
988 implicit none
989 integer, intent(in) :: domid
990 class(meshbase), intent(in) :: mesh
991 class(modelvarmanager), intent(inout) :: auxvars2d_list
992 class(localmeshfieldbase), pointer, intent(out) :: prec, prec_engi
993 class(localmesh2d), pointer, intent(out), optional :: lcmesh2d
994
995 class(meshfieldbase), pointer :: field
996 class(localmeshbase), pointer :: lcmesh
997 !-------------------------------------------------------
998
999 !--
1000 call auxvars2d_list%Get(atmos_auxvars2d_prec_id, field)
1001 call field%GetLocalMeshField(domid, prec)
1002
1003 call auxvars2d_list%Get(atmos_auxvars2d_prec_engi_id, field)
1004 call field%GetLocalMeshField(domid, prec_engi)
1005
1006 if (present(lcmesh2d)) then
1007 call mesh%GetLocalMesh( domid, lcmesh )
1008 nullify( lcmesh2d )
1009
1010 select type(lcmesh)
1011 type is (localmesh2d)
1012 if (present(lcmesh2d)) lcmesh2d => lcmesh
1013 end select
1014 end if
1015
1016 return
1017 end subroutine atmosvars_getlocalmeshsfcvar
1018
1019!OCL SERIAL
1020 subroutine atmosvars_getlocalmeshqtrcvar( domID, mesh, trcvars_list, &
1021 varid, &
1022 var, lcmesh3D )
1023
1024 implicit none
1025 integer, intent(in) :: domid
1026 class(meshbase), intent(in) :: mesh
1027 class(modelvarmanager), intent(inout) :: trcvars_list
1028 integer, intent(in) :: varid
1029 class(localmeshfieldbase), pointer, intent(out) :: var
1030 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
1031
1032 class(meshfieldbase), pointer :: field
1033 class(localmeshbase), pointer :: lcmesh
1034 !-------------------------------------------------------
1035
1036 !--
1037 call trcvars_list%Get(varid, field)
1038 call field%GetLocalMeshField(domid, var)
1039
1040 if (present(lcmesh3d)) then
1041 call mesh%GetLocalMesh( domid, lcmesh )
1042 nullify( lcmesh3d )
1043
1044 select type(lcmesh)
1045 type is (localmesh3d)
1046 if (present(lcmesh3d)) lcmesh3d => lcmesh
1047 end select
1048 end if
1049
1050 return
1051 end subroutine atmosvars_getlocalmeshqtrcvar
1052
1053!OCL SERIAL
1054 subroutine atmosvars_getlocalmeshqtrc_qv( domID, mesh, trcvars_list, forcing_list, &
1055 var, var_tp, lcmesh3D )
1056
1057 use scale_atmos_hydrometeor, only: &
1058 atmos_hydrometeor_dry, &
1059 i_qv
1060 implicit none
1061 integer, intent(in) :: domid
1062 class(meshbase), intent(in) :: mesh
1063 class(modelvarmanager), intent(inout) :: trcvars_list
1064 class(modelvarmanager), intent(inout) :: forcing_list
1065 class(localmeshfieldbase), pointer, intent(out) :: var
1066 class(localmeshfieldbase), pointer, intent(out) :: var_tp
1067 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
1068
1069 class(meshfieldbase), pointer :: field
1070 class(localmeshbase), pointer :: lcmesh
1071
1072 integer :: iq, tend_iq
1073 !-------------------------------------------------------
1074
1075 !--
1076 if ( atmos_hydrometeor_dry ) then
1077 iq =0; tend_iq = phytend_num1+1
1078 else
1079 iq = i_qv; tend_iq = phytend_num1 + i_qv
1080 end if
1081
1082 call trcvars_list%Get(iq, field)
1083 call field%GetLocalMeshField(domid, var)
1084
1085 call forcing_list%Get(tend_iq, field)
1086 call field%GetLocalMeshField(domid, var_tp)
1087
1088 if (present(lcmesh3d)) then
1089 call mesh%GetLocalMesh( domid, lcmesh )
1090 nullify( lcmesh3d )
1091
1092 select type(lcmesh)
1093 type is (localmesh3d)
1094 if (present(lcmesh3d)) lcmesh3d => lcmesh
1095 end select
1096 end if
1097
1098 return
1099 end subroutine atmosvars_getlocalmeshqtrc_qv
1100
1101!OCL SERIAL
1102 subroutine atmosvars_getlocalmeshqtrcvarlist( domID, mesh, trcvars_list, &
1103 varid_s, &
1104 var_list, lcmesh3D )
1105
1106 implicit none
1107 integer, intent(in) :: domid
1108 class(meshbase), intent(in) :: mesh
1109 class(modelvarmanager), intent(inout) :: trcvars_list
1110 integer, intent(in) :: varid_s
1111 type(localmeshfieldbaselist), intent(out) :: var_list(:)
1112 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
1113
1114 class(meshfieldbase), pointer :: field
1115 class(localmeshbase), pointer :: lcmesh
1116
1117 integer :: iq
1118 !-------------------------------------------------------
1119
1120 !--
1121 do iq = varid_s, varid_s + size(var_list) - 1
1122 call trcvars_list%Get(iq, field)
1123 call field%GetLocalMeshField(domid, var_list(iq-varid_s+1)%ptr)
1124 end do
1125 if (present(lcmesh3d)) then
1126 call mesh%GetLocalMesh( domid, lcmesh )
1127 nullify( lcmesh3d )
1128
1129 select type(lcmesh)
1130 type is (localmesh3d)
1131 if (present(lcmesh3d)) lcmesh3d => lcmesh
1132 end select
1133 end if
1134
1135 return
1137
1138!OCL SERIAL
1139 subroutine atmosvars_getlocalmeshphyauxvars( domID, mesh, phyauxvars_list, &
1140 PRES, PT, &
1141 lcmesh3D )
1142
1143 implicit none
1144 integer, intent(in) :: domid
1145 class(meshbase), intent(in) :: mesh
1146 class(modelvarmanager), intent(inout) :: phyauxvars_list
1147 class(localmeshfieldbase), pointer, intent(out) :: pres, pt
1148 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
1149
1150 class(meshfieldbase), pointer :: field
1151 class(localmeshbase), pointer :: lcmesh
1152 !-------------------------------------------------------
1153
1154 !--
1155 call phyauxvars_list%Get(auxvar_pres_id, field)
1156 call field%GetLocalMeshField(domid, pres)
1157
1158 call phyauxvars_list%Get(auxvar_pt_id, field)
1159 call field%GetLocalMeshField(domid, pt)
1160
1161 !---
1162
1163 if (present(lcmesh3d)) then
1164 call mesh%GetLocalMesh( domid, lcmesh )
1165 nullify( lcmesh3d )
1166
1167 select type(lcmesh)
1168 type is (localmesh3d)
1169 if (present(lcmesh3d)) lcmesh3d => lcmesh
1170 end select
1171 end if
1172
1173 return
1175
1176!OCL SERIAL
1177 subroutine atmosvars_getlocalmeshphytends( domID, mesh, phytends_list, &
1178 DENS_tp, MOMX_tp, MOMY_tp, MOMZ_tp, RHOT_tp, RHOH_p, &
1179 RHOQ_tp, &
1180 lcmesh3D )
1181
1182 implicit none
1183 integer, intent(in) :: domid
1184 class(meshbase), intent(in) :: mesh
1185 class(modelvarmanager), intent(inout) :: phytends_list
1186 class(localmeshfieldbase), pointer, intent(out) :: dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp
1187 class(localmeshfieldbase), pointer, intent(out) :: rhoh_p
1188 type(localmeshfieldbaselist), intent(inout), optional :: rhoq_tp(qa)
1189 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
1190
1191 class(meshfieldbase), pointer :: field
1192 class(localmeshbase), pointer :: lcmesh
1193
1194 integer :: iq
1195 !-------------------------------------------------------
1196
1197 !--
1198 call phytends_list%Get(phytend_dens_id, field)
1199 call field%GetLocalMeshField(domid, dens_tp)
1200
1201 call phytends_list%Get(phytend_momx_id, field)
1202 call field%GetLocalMeshField(domid, momx_tp)
1203
1204 call phytends_list%Get(phytend_momy_id, field)
1205 call field%GetLocalMeshField(domid, momy_tp)
1206
1207 call phytends_list%Get(phytend_momz_id, field)
1208 call field%GetLocalMeshField(domid, momz_tp)
1209
1210 call phytends_list%Get(phytend_rhot_id, field)
1211 call field%GetLocalMeshField(domid, rhot_tp)
1212
1213 call phytends_list%Get(phytend_rhoh_id, field)
1214 call field%GetLocalMeshField(domid, rhoh_p)
1215
1216 if ( present(rhoq_tp) ) then
1217 do iq = 1, qa
1218 call phytends_list%Get(phytend_num1+iq, field)
1219 call field%GetLocalMeshField(domid, rhoq_tp(iq)%ptr)
1220 end do
1221 end if
1222
1223 !---
1224 if ( present(lcmesh3d) ) then
1225 call mesh%GetLocalMesh( domid, lcmesh )
1226 nullify( lcmesh3d )
1227
1228 select type(lcmesh)
1229 type is (localmesh3d)
1230 if ( present(lcmesh3d) ) lcmesh3d => lcmesh
1231 end select
1232 end if
1233
1234 return
1235 end subroutine atmosvars_getlocalmeshphytends
1236
1237!OCL SERIAL
1238 subroutine atmosvars_getlocalmeshqtrcphytend( domID, mesh, phytends_list, &
1239 qtrcid, &
1240 RHOQ_tp )
1241
1242 implicit none
1243 integer, intent(in) :: domid
1244 class(meshbase), intent(in) :: mesh
1245 class(modelvarmanager), intent(inout) :: phytends_list
1246 integer, intent(in) :: qtrcid
1247 class(localmeshfieldbase), pointer, intent(out) :: rhoq_tp
1248
1249 class(meshfieldbase), pointer :: field
1250 class(localmeshbase), pointer :: lcmesh
1251 !-------------------------------------------------------
1252
1253 call phytends_list%Get(phytend_num1 + qtrcid, field)
1254 call field%GetLocalMeshField(domid, rhoq_tp)
1255
1256 return
1257 end subroutine atmosvars_getlocalmeshqtrcphytend
1258
1259 !-----------------------------------------------------------------------------
1261!OCL SERIAL
1262 subroutine atmosvars_calculatediagnostics( this )
1263 use scale_const, only: &
1264 rdry => const_rdry, &
1265 cpdry => const_cpdry, &
1266 cvdry => const_cvdry, &
1267 pres00 => const_pre00
1268 use scale_tracer, only: &
1269 tracer_mass, tracer_r, tracer_cv, tracer_cp
1270 use scale_atmos_thermodyn, only: &
1271 atmos_thermodyn_specific_heat
1272 implicit none
1273 class(atmosvars), intent(inout), target :: this
1274
1275 class(localmesh3d), pointer :: lcmesh3d
1276 integer :: n
1277 integer :: varid
1278 integer :: ke
1279
1280 class(meshfield3d), pointer :: field
1281 class(elementbase3d), pointer :: elem3d
1282
1283 type(localmeshfieldbaselist) :: qtrc(qa)
1284 !-------------------------------------------------------
1285
1286 ! Calculate specific heat
1287 call vars_calc_specific_heat( this )
1288
1289 ! Calculate diagnostic variables
1290 do varid=auxvar_denshydro_id+1, auxvar_pt_id
1291 field => this%AUX_VARS(varid)
1292 do n=1, field%mesh%LOCAL_MESH_NUM
1294 field%mesh, this%QTRCVARS_manager, &
1295 1, qtrc, lcmesh3d )
1296
1297 elem3d => lcmesh3d%refElem3D
1298
1299 call vars_calc_diagnosevar_lc( &
1300 field%varname, field%local(n)%val, &
1301 this%PROG_VARS(prgvar_ddens_id)%local(n)%val, &
1302 this%PROG_VARS(prgvar_momx_id)%local(n)%val, &
1303 this%PROG_VARS(prgvar_momy_id)%local(n)%val, &
1304 this%PROG_VARS(prgvar_momz_id)%local(n)%val, &
1305 this%AUX_VARS(auxvar_pres_id)%local(n)%val, &
1306 this%AUX_VARS(auxvar_qdry_id)%local(n)%val, &
1307 qtrc, &
1308 this%AUX_VARS(auxvar_denshydro_id)%local(n)%val, &
1309 this%AUX_VARS(auxvar_preshydro_id)%local(n)%val, &
1310 this%AUX_VARS(auxvar_rtot_id )%local(n)%val, &
1311 this%AUX_VARS(auxvar_cvtot_id)%local(n)%val, &
1312 this%AUX_VARS(auxvar_cptot_id)%local(n)%val, &
1313 lcmesh3d, lcmesh3d%refElem3D )
1314 end do
1315 end do
1316
1317 return
1318 end subroutine atmosvars_calculatediagnostics
1319
1320!OCL SERIAL
1321 subroutine atmosvars_calcdiagvar( this, field_name, field_work )
1322 use scale_const, only: &
1323 rdry => const_rdry, &
1324 cpdry => const_cpdry, &
1325 cvdry => const_cvdry, &
1326 pres00 => const_pre00
1327 use scale_tracer, only: &
1328 tracer_mass, tracer_r, tracer_cv, tracer_cp
1329 use scale_atmos_thermodyn, only: &
1330 atmos_thermodyn_specific_heat
1331
1332 implicit none
1333 class(atmosvars), intent(inout) :: this
1334 character(*), intent(in) :: field_name
1335 type(meshfield3d), intent(inout) :: field_work
1336
1337 class(localmesh3d), pointer :: lcmesh3d
1338 class(elementbase3d), pointer :: elem3d
1339 integer :: n
1340 integer :: ke
1341 integer :: iq
1342
1343 type(meshfield3d) :: field_work_uvmet(2)
1344 logical :: is_uvmet
1345 integer :: uvmet_i
1346
1347 type(localmeshfieldbaselist) :: qtrc(qa)
1348 !--------------------------------------------------
1349
1350 is_uvmet = .false.
1351 if ( field_name == 'Umet' ) then
1352 is_uvmet = .true.; uvmet_i = 1
1353 else if ( field_name == 'Vmet' ) then
1354 is_uvmet = .true.; uvmet_i = 2
1355 end if
1356
1357 field_work%varname = field_name
1358
1359 do n=1, field_work%mesh%LOCAL_MESH_NUM
1361 field_work%mesh, this%QTRCVARS_manager, &
1362 1, qtrc, lcmesh3d )
1363
1364 if ( .not. is_uvmet ) then
1365 call vars_calc_diagnosevar_lc( field_name, field_work%local(n)%val, &
1366 this%PROG_VARS(prgvar_ddens_id)%local(n)%val, &
1367 this%PROG_VARS(prgvar_momx_id)%local(n)%val, &
1368 this%PROG_VARS(prgvar_momy_id)%local(n)%val, &
1369 this%PROG_VARS(prgvar_momz_id)%local(n)%val, &
1370 this%AUX_VARS(auxvar_pres_id)%local(n)%val, &
1371 this%AUX_VARS(auxvar_qdry_id)%local(n)%val, &
1372 qtrc, &
1373 this%AUX_VARS(auxvar_denshydro_id)%local(n)%val, &
1374 this%AUX_VARS(auxvar_preshydro_id)%local(n)%val, &
1375 this%AUX_VARS(auxvar_rtot_id )%local(n)%val, &
1376 this%AUX_VARS(auxvar_cvtot_id)%local(n)%val, &
1377 this%AUX_VARS(auxvar_cptot_id)%local(n)%val, &
1378 lcmesh3d, lcmesh3d%refElem3D )
1379 else
1380 call field_work_uvmet(1)%Init( 'Umet', '', field_work%mesh )
1381 call field_work_uvmet(2)%Init( 'Vmet', '', field_work%mesh )
1382 call this%mesh%Calc_UVmet( &
1383 this%PROG_VARS(prgvar_momx_id), this%PROG_VARS(prgvar_momy_id), & ! (in)
1384 field_work_uvmet(1), field_work_uvmet(2) ) ! (inout)
1385 !$omp parallel do
1386 do ke=lcmesh3d%NeS, lcmesh3d%NeE
1387 field_work%local(n)%val(:,ke) = field_work_uvmet(uvmet_i)%local(n)%val(:,ke) &
1388 / ( this%AUX_VARS (auxvar_denshydro_id)%local(n)%val(:,ke) &
1389 + this%PROG_VARS(prgvar_ddens_id )%local(n)%val(:,ke) )
1390 end do
1391 call field_work_uvmet(1)%Final()
1392 call field_work_uvmet(2)%Final()
1393 end if
1394 end do
1395
1396 return
1397 end subroutine atmosvars_calcdiagvar
1398
1399!OCL SERIAL
1400 subroutine atmosvars_calcdiagvar2d( this, field_name, field_work )
1401 implicit none
1402 class(atmosvars), intent(inout) :: this
1403 character(*), intent(in) :: field_name
1404 type(meshfield2d), intent(inout), target :: field_work
1405
1406 class(localmesh2d), pointer :: lcmesh2d
1407 class(elementbase2d), pointer :: elem2d
1408
1409 integer :: n
1410 !--------------------------------------------------
1411
1412 field_work%varname = field_name
1413
1414 do n=1, field_work%mesh%LOCAL_MESH_NUM
1415 lcmesh2d => field_work%mesh%lcmesh_list(n)
1416 call vars_calc_diagnosevar2d_lc( field_name, field_work%local(n)%val, &
1417 this%ptr_MP_AUXVARS2D_manager, &
1418 field_work%mesh, lcmesh2d, lcmesh2d%refElem2D )
1419 end do
1420
1421 return
1422 end subroutine atmosvars_calcdiagvar2d
1423
1424!-- private -----------------------------------------------------------------------
1425
1426!OCL SERIAL
1427 subroutine vars_calc_diagnosevar_lc( field_name, var_out, &
1428 DDENS_, MOMX_, MOMY_, MOMZ_, PRES_, QDRY_, QTRC, &
1429 DENS_hyd, PRES_hyd, Rtot, CVtot, CPTot, &
1430 lcmesh, elem )
1431
1432 use scale_const, only: &
1433 grav => const_grav, &
1434 rdry => const_rdry, &
1435 rvap => const_rvap, &
1436 cpdry => const_cpdry, &
1437 cvdry => const_cvdry, &
1438 pres00 => const_pre00
1439 use scale_tracer, only: &
1440 tracer_inq_id, &
1441 tracer_cv, tracer_engi0
1442 use scale_atmos_saturation, only: &
1443 atmos_saturation_psat_liq
1444 implicit none
1445
1446 class(localmesh3d), intent(in) :: lcmesh
1447 class(elementbase3d), intent(in) :: elem
1448 character(*), intent(in) :: field_name
1449 real(rp), intent(out) :: var_out(elem%np,lcmesh%nea)
1450 real(rp), intent(in) :: ddens_(elem%np,lcmesh%nea)
1451 real(rp), intent(in) :: momx_(elem%np,lcmesh%nea)
1452 real(rp), intent(in) :: momy_(elem%np,lcmesh%nea)
1453 real(rp), intent(in) :: momz_(elem%np,lcmesh%nea)
1454 real(rp), intent(in) :: pres_(elem%np,lcmesh%nea)
1455 real(rp), intent(in) :: qdry_(elem%np,lcmesh%nea)
1456 type(localmeshfieldbaselist), intent(in) :: qtrc(qa)
1457 real(rp), intent(in) :: dens_hyd(elem%np,lcmesh%nea)
1458 real(rp), intent(in) :: pres_hyd(elem%np,lcmesh%nea)
1459 real(rp), intent(in) :: rtot (elem%np,lcmesh%nea)
1460 real(rp), intent(in) :: cvtot(elem%np,lcmesh%nea)
1461 real(rp), intent(in) :: cptot(elem%np,lcmesh%nea)
1462
1463 integer :: ke, ke2d
1464 integer :: iq
1465 real(rp) :: dens(elem%np), temp(elem%np)
1466 real(rp) :: mom_u1(elem%np), mom_u2(elem%np), g_11(elem%np), g_12(elem%np), g_22(elem%np)
1467 real(rp) :: psat(elem%np)
1468
1469 integer :: iq_qv
1470 !-------------------------------------------------------------------------
1471
1472 select case(trim(field_name))
1473 case('DENS')
1474 !$omp parallel do
1475 do ke=1, lcmesh%Ne
1476 var_out(:,ke) = ddens_(:,ke) + dens_hyd(:,ke)
1477 end do
1478
1479 case('U')
1480 !$omp parallel do private (DENS)
1481 do ke=1, lcmesh%Ne
1482 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1483 var_out(:,ke) = momx_(:,ke) / dens(:)
1484 end do
1485
1486 case('V')
1487 !$omp parallel do private (DENS)
1488 do ke=1, lcmesh%Ne
1489 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1490 var_out(:,ke) = momy_(:,ke) / dens(:)
1491 end do
1492
1493 case('W')
1494 !$omp parallel do private (DENS)
1495 do ke=1, lcmesh%Ne
1496 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1497 var_out(:,ke) = momz_(:,ke) / dens(:)
1498 end do
1499
1500 case ( 'PRES' )
1501 case('PRES_diff')
1502 !$omp parallel do
1503 do ke=1, lcmesh%Ne
1504 var_out(:,ke) = pres_(:,ke) - pres_hyd(:,ke)
1505 end do
1506
1507 case('T')
1508 !$omp parallel do
1509 do ke=1, lcmesh%Ne
1510 var_out(:,ke) = pres_(:,ke) / (rtot(:,ke) * (ddens_(:,ke) + dens_hyd(:,ke)) )
1511 end do
1512
1513 case('T_diff')
1514 !$omp parallel do
1515 do ke=1, lcmesh%Ne
1516 var_out(:,ke) = pres_(:,ke) / ( rtot(:,ke) * (ddens_(:,ke) + dens_hyd(:,ke)) ) &
1517 - pres_hyd(:,ke) / ( rdry * dens_hyd(:,ke) )
1518 end do
1519
1520 case('PT')
1521 !$omp parallel do private( DENS )
1522 do ke=1, lcmesh%Ne
1523 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1524 var_out(:,ke) = pres_(:,ke) / (rtot(:,ke) * dens(:) ) * ( pres00 / pres_(:,ke) )**( rtot(:,ke) / cptot(:,ke) )
1525 end do
1526
1527 case('PT_diff')
1528 !$omp parallel do private( DENS )
1529 do ke=1, lcmesh%Ne
1530 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1531 var_out(:,ke) = pres_(:,ke) / (rtot(:,ke) * dens(:) ) * ( pres00 / pres_(:,ke) )**( rtot(:,ke) / cptot(:,ke) ) &
1532 - pres00/rdry * (pres_hyd(:,ke)/pres00)**(cvdry/cpdry) / dens_hyd(:,ke)
1533 end do
1534
1535 case( 'RH', 'RHL' )
1536 if ( atmos_hydrometeor_dry ) then
1537 var_out(:,ke) = 0.0_rp
1538 else
1539 call tracer_inq_id( "QV", iq_qv )
1540
1541 !$omp parallel do private (TEMP, PSAT)
1542 do ke=1, lcmesh%Ne
1543 temp(:) = pres_(:,ke) / (rtot(:,ke) * (ddens_(:,ke) + dens_hyd(:,ke)) )
1544
1545 call atmos_saturation_psat_liq( &
1546 elem%Np, 1, elem%Np, temp(:), & ! (in)
1547 psat(:) ) ! (out)
1548
1549 var_out(:,ke) = ( ddens_(:,ke) + dens_hyd(:,ke) ) * qtrc(iq_qv)%ptr%val(:,ke) &
1550 / psat(:) * rvap * temp(:) * 100.0_rp
1551 end do
1552 end if
1553
1554 case('ENGK')
1555 !$omp parallel do private (ke2D, DENS, mom_u1, mom_u2, G_11, G_12, G_22)
1556 do ke=1, lcmesh%Ne
1557 ke2d = lcmesh%EMap3Dto2D(ke)
1558
1559 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1560 g_11(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,1)
1561 g_12(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,2)
1562 g_22(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,2,2)
1563
1564 mom_u1(:) = g_11(:) * momx_(:,ke) + g_12(:) * momy_(:,ke)
1565 mom_u2(:) = g_12(:) * momx_(:,ke) + g_22(:) * momy_(:,ke)
1566
1567 var_out(:,ke) = 0.5_rp * ( momx_(:,ke) * mom_u1(:) + momy_(:,ke) * mom_u2(:) + momz_(:,ke)**2 ) / dens(:)
1568 end do
1569
1570 case('ENGP')
1571 !$omp parallel do private (DENS)
1572 do ke=1, lcmesh%Ne
1573 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1574 var_out(:,ke) = dens(:) * grav * lcmesh%zlev(:,ke)
1575 end do
1576
1577 case('ENGI')
1578 !$omp parallel do private (DENS, iq)
1579 do ke=1, lcmesh%Ne
1580 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1581 var_out(:,ke) = qdry_(:,ke) * pres_(:,ke) / rtot(:,ke) * cvdry
1582 do iq = 1, qa
1583 var_out(:,ke) = var_out(:,ke) &
1584 + qtrc(iq)%ptr%val(:,ke) * ( pres_(:,ke) / rtot(:,ke) * tracer_cv(iq) + dens(:) * tracer_engi0(iq) )
1585 end do
1586 end do
1587
1588 case('ENGT')
1589 !$omp parallel do private (ke2D, DENS, mom_u1, mom_u2, iq, G_11, G_12, G_22)
1590 do ke=1, lcmesh%Ne
1591 ke2d = lcmesh%EMap3Dto2D(ke)
1592
1593 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1594
1595 g_11(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,1)
1596 g_12(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,2)
1597 g_22(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,2,2)
1598 mom_u1(:) = g_11(:) * momx_(:,ke) + g_12(:) * momy_(:,ke)
1599 mom_u2(:) = g_12(:) * momx_(:,ke) + g_22(:) * momy_(:,ke)
1600
1601 ! ENGI
1602 var_out(:,ke) = qdry_(:,ke) * pres_(:,ke) / rtot(:,ke) * cvdry
1603 do iq = 1, qa
1604 var_out(:,ke) = var_out(:,ke) &
1605 + qtrc(iq)%ptr%val(:,ke) * ( pres_(:,ke) / rtot(:,ke) * tracer_cv(iq) + dens(:) * tracer_engi0(iq) )
1606 end do
1607 ! ENGT
1608 var_out(:,ke) = &
1609 0.5_rp * ( momx_(:,ke) * mom_u1(:) + momy_(:,ke) * mom_u2(:) + momz_(:,ke)**2 ) / dens(:) & ! ENGK
1610 + var_out(:,ke) & ! ENGI
1611 + dens(:) * grav * lcmesh%pos_en(:,ke,3) ! ENGP
1612 end do
1613
1614 case default
1615 log_error("AtmosVars_calc_diagnoseVar_lc",*) 'The name of diagnostic variable is not suported. Check!', field_name
1616 call prc_abort
1617
1618 end select
1619
1620 return
1621 end subroutine vars_calc_diagnosevar_lc
1622
1623!OCL SERIAL
1624 subroutine vars_calc_diagnosevar2d_lc( field_name, & ! (in)
1625 var_out, & ! (out)
1626 mp_auxvars2d, mesh2d, lcmesh, elem ) ! (in)
1627
1628 use mod_atmos_phy_mp_vars, only: &
1630
1631 implicit none
1632 class(localmesh2d), intent(in) :: lcmesh
1633 class(elementbase2d), intent(in) :: elem
1634 character(*), intent(in) :: field_name
1635 real(rp), intent(out) :: var_out(elem%np,lcmesh%nea)
1636 class(modelvarmanager), intent(inout) :: mp_auxvars2d
1637 class(meshbase2d), intent(in) :: mesh2d
1638
1639 integer :: ke
1640
1641 class(localmeshfieldbase), pointer :: sflx_rain_mp, sflx_snow_mp, sflx_engi_mp
1642 !-------------------------------------------------------------------------
1643
1644 select case(trim(field_name))
1645 case('RAIN', 'SNOW')
1647 lcmesh%lcdomID, mesh2d, mp_auxvars2d, &
1648 sflx_rain_mp, sflx_snow_mp, sflx_engi_mp )
1649 end select
1650
1651 select case(trim(field_name))
1652 case('RAIN')
1653 !$omp parallel do
1654 do ke=lcmesh%NeS, lcmesh%NeE
1655 var_out(:,ke) = sflx_rain_mp%val(:,ke)
1656 end do
1657 case('SNOW')
1658 !$omp parallel do
1659 do ke=lcmesh%NeS, lcmesh%NeE
1660 var_out(:,ke) = sflx_snow_mp%val(:,ke)
1661 end do
1662 case default
1663 log_error("AtmosVars_calc_diagnoseVar2D_lc",*) 'The name of diagnostic variable is not suported. Check!', field_name
1664 call prc_abort
1665 end select
1666
1667 return
1668 end subroutine vars_calc_diagnosevar2d_lc
1669
1670!OCL SERIAL
1671 subroutine vars_calc_specific_heat( this )
1672 use scale_const, only: &
1673 rdry => const_rdry, &
1674 cpdry => const_cpdry, &
1675 cvdry => const_cvdry, &
1676 pres00 => const_pre00
1677 use scale_tracer, only: &
1678 tracer_mass, tracer_r, tracer_cv, tracer_cp
1679 use scale_atmos_thermodyn, only: &
1680 atmos_thermodyn_specific_heat
1681 implicit none
1682 class(atmosvars), intent(inout), target :: this
1683
1684 class(localmesh3d), pointer :: lcmesh3d
1685 integer :: n
1686 integer :: varid
1687 integer :: ke
1688 integer :: iq
1689
1690 class(elementbase3d), pointer :: elem3d
1691
1692 real(rp), allocatable :: q_tmp(:,:)
1693 !-------------------------------------------------------
1694
1695 ! Calculate specific heat
1696 do n=1, this%AUX_VARS(1)%mesh%LOCAL_MESH_NUM
1697 lcmesh3d => this%AUX_VARS(1)%mesh%lcmesh_list(n)
1698 elem3d => lcmesh3d%refElem3D
1699 allocate( q_tmp(elem3d%Np,qa) )
1700
1701 !$omp parallel do private(ke, iq, q_tmp)
1702 do ke = lcmesh3d%NeS, lcmesh3d%NeE
1703 do iq = 1, qa
1704 q_tmp(:,iq) = this%QTRC_VARS(iq)%local(n)%val(:,ke)
1705 end do
1706 call atmos_thermodyn_specific_heat( &
1707 elem3d%Np, 1, elem3d%Np, qa, & ! (in)
1708 q_tmp(:,:), tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), & ! (in)
1709 this%AUX_VARS(auxvar_qdry_id )%local(n)%val(:,ke), & ! (out)
1710 this%AUX_VARS(auxvar_rtot_id )%local(n)%val(:,ke), & ! (out)
1711 this%AUX_VARS(auxvar_cvtot_id)%local(n)%val(:,ke), & ! (out)
1712 this%AUX_VARS(auxvar_cptot_id)%local(n)%val(:,ke) ) ! (out)
1713 end do
1714 deallocate(q_tmp)
1715 end do
1716 end subroutine vars_calc_specific_heat
1717
1718end module mod_atmos_vars
module Atmosphere / Mesh
module ATMOSPHERE physics / Cloud Microphysics
subroutine, public atmosphympvars_getlocalmeshfields_sfcflx(domid, mesh, sfcflx_list, sflx_rain, sflx_snow, sflx_engi)
module ATMOSPHERE / Variables
integer, parameter, public atmos_diagvars2d_num
real(rp), dimension(prgvar_num), parameter progvars_check_min
integer, parameter, public atmos_diagvars_dens_id
integer, parameter, public atmos_diagvars_w_id
subroutine, public atmosvars_getlocalmeshsfcvar(domid, mesh, auxvars2d_list, prec, prec_engi, lcmesh2d)
type(variableinfo), dimension(atmos_diagvars3d_num), public atmos_diagvars3d_vinfo
subroutine, public atmosvars_getlocalmeshqtrc_qv(domid, mesh, trcvars_list, forcing_list, var, var_tp, lcmesh3d)
subroutine, public atmosvars_getlocalmeshprgvar(domid, mesh, prgvars_list, auxvars_list, varid, var, dens_hyd, pres_hyd, lcmesh3d)
integer, parameter, public atmos_auxvars2d_prec_id
integer, parameter, public atmos_diagvars3d_num
subroutine, public atmosvars_getlocalmeshprgvars(domid, mesh, prgvars_list, auxvars_list, ddens, momx, momy, momz, therm, dens_hyd, pres_hyd, rtot, cvtot, cptot, lcmesh3d)
integer, parameter, public atmos_diagvars_vmet_id
integer, parameter, public atmos_diagvars_rh_id
integer, parameter, public atmos_diagvars_engi_id
integer, parameter, public atmos_diagvars_umet_id
integer, parameter, public atmos_diagvars_qdry_id
subroutine, public atmosvars_getlocalmeshqtrcvarlist(domid, mesh, trcvars_list, varid_s, var_list, lcmesh3d)
integer, parameter, public atmos_auxvars2d_prec_engi_id
integer, parameter, public atmos_diagvars_engk_id
type(variableinfo), dimension(atmos_diagvars2d_num), public atmos_diagvars2d_vinfo
integer, parameter, public atmos_diagvars_v_id
integer, parameter, public atmos_auxvars2d_num
type(variableinfo), dimension(atmos_auxvars2d_num), public atmos_auxvars2d_vinfo
subroutine, public atmosvars_getlocalmeshphytends(domid, mesh, phytends_list, dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp, rhoh_p, rhoq_tp, lcmesh3d)
subroutine, public atmosvars_getlocalmeshqtrcphytend(domid, mesh, phytends_list, qtrcid, rhoq_tp)
integer, parameter, public atmos_diagvars_t_id
subroutine, public atmosvars_getlocalmeshqtrcvar(domid, mesh, trcvars_list, varid, var, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphyauxvars(domid, mesh, phyauxvars_list, pres, pt, lcmesh3d)
integer, parameter, public atmos_diagvars_snow_id
integer, parameter, public atmos_diagvars_engt_id
integer, parameter, public atmos_diagvars_engp_id
integer, parameter, public atmos_diagvars_rain_id
integer, parameter, public atmos_diagvars_u_id
module FElib / Fluid dyn solver / Atmosphere / driver (3D nonhydrostatic model)
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / Common
subroutine, public atm_dyn_dgm_nonhydro3d_common_get_varinfo(prgvar_info, auxvar_info, phytend_info)
subroutine, public atm_dyn_dgm_nonhydro3d_common_setup_variables(prgvars, qtrcvars, auxvars, phytends, prgvar_manager, qtrcvar_manager, auxvar_manager, phytend_manager, phytend_num_tot, mesh3d, prgvar_varinfo)
module FElib / Element / Base
type(file_restart_meshfield), public restart_file
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
integer, public meshbase3d_dimtypeid_xyz
module FElib / Mesh / Base
module FElib / Data / base
module FElib / Data / Statistics
module FElib / Data / Communication base
module FElib / Data / Communication 3D cubic domain
FElib / model framework / variable manager.
Derived type to provide a driver of dynamical core with the atmospheric nonhydrostatic equations.