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