SCALE-RM
mod_atmos_vars.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
20  use scale_debug
22  use scale_index
23  use scale_tracer
24  !-----------------------------------------------------------------------------
25  implicit none
26  private
27  !-----------------------------------------------------------------------------
28  !
29  !++ Public procedure
30  !
31  public :: atmos_vars_setup
32  public :: atmos_vars_fillhalo
33  public :: atmos_vars_restart_read
34  public :: atmos_vars_restart_write
35  public :: atmos_vars_restart_check
37  public :: atmos_vars_history
38  public :: atmos_vars_check
40  public :: atmos_vars_get_diagnostic
41  public :: atmos_vars_monitor
42  public :: atmos_vars_finalize
43 
45  public :: atmos_vars_restart_open
48  public :: atmos_vars_restart_close
49 
50  interface atmos_vars_get_diagnostic
51  procedure ATMOS_vars_get_diagnostic_3D
52  procedure ATMOS_vars_get_diagnostic_2D
53  procedure ATMOS_vars_get_diagnostic_1D
54  end interface atmos_vars_get_diagnostic
55 
56  !-----------------------------------------------------------------------------
57  !
58  !++ Public parameters & variables
59  !
60  logical, public :: atmos_restart_output = .false.
61 
62  character(len=H_LONG), public :: atmos_restart_in_basename = ''
63  logical, public :: atmos_restart_in_aggregate
64  logical, public :: atmos_restart_in_postfix_timelabel = .false.
65  character(len=H_LONG), public :: atmos_restart_out_basename = ''
66  logical, public :: atmos_restart_out_aggregate
67  logical, public :: atmos_restart_out_postfix_timelabel = .true.
68  character(len=H_MID), public :: atmos_restart_out_title = 'ATMOS restart'
69  character(len=H_SHORT), public :: atmos_restart_out_dtype = 'DEFAULT'
70 
71  logical, public :: atmos_restart_check = .false.
72  character(len=H_LONG), public :: atmos_restart_check_basename = 'restart_check'
73  real(rp), public :: atmos_restart_check_criterion = 1.e-6_rp
74 
75  ! prognostic variables
76  real(rp), public, target, allocatable :: dens(:,:,:) ! Density [kg/m3]
77  real(rp), public, target, allocatable :: momz(:,:,:) ! momentum z [kg/m2/s]
78  real(rp), public, target, allocatable :: momx(:,:,:) ! momentum x [kg/m2/s]
79  real(rp), public, target, allocatable :: momy(:,:,:) ! momentum y [kg/m2/s]
80  real(rp), public, target, allocatable :: rhot(:,:,:) ! DENS * POTT [K*kg/m3]
81  real(rp), public, target, allocatable :: qtrc(:,:,:,:) ! ratio of mass of tracer to total mass[kg/kg]
82 
83  real(rp), public, target, allocatable :: dens_avw(:,:,:)
84  real(rp), public, target, allocatable :: momz_avw(:,:,:)
85  real(rp), public, target, allocatable :: momx_avw(:,:,:)
86  real(rp), public, target, allocatable :: momy_avw(:,:,:)
87  real(rp), public, target, allocatable :: rhot_avw(:,:,:)
88  real(rp), public, target, allocatable :: qtrc_avw(:,:,:,:)
89 
90  real(rp), public, pointer :: dens_av(:,:,:)
91  real(rp), public, pointer :: momz_av(:,:,:)
92  real(rp), public, pointer :: momx_av(:,:,:)
93  real(rp), public, pointer :: momy_av(:,:,:)
94  real(rp), public, pointer :: rhot_av(:,:,:)
95  real(rp), public, pointer :: qtrc_av(:,:,:,:)
96 
97  real(rp), public, pointer :: qv(:,:,:)
98  real(rp), public, pointer :: qc(:,:,:)
99  real(rp), public, pointer :: qr(:,:,:)
100  real(rp), public, pointer :: qi(:,:,:)
101  real(rp), public, pointer :: qs(:,:,:)
102  real(rp), public, pointer :: qg(:,:,:)
103  real(rp), public, pointer :: qh(:,:,:)
104 
105  real(rp), public, target, allocatable :: qe(:,:,:,:)
106 
107  ! reference state
108  real(rp), public, allocatable :: dens_ref(:,:,:)
109  real(rp), public, allocatable :: pott_ref(:,:,:)
110  real(rp), public, allocatable :: temp_ref(:,:,:)
111  real(rp), public, allocatable :: pres_ref(:,:,:)
112  real(rp), public, allocatable :: qv_ref(:,:,:)
113 
114  ! tendency by physical processes
115  real(rp), public, allocatable :: dens_tp(:,:,:)
116  real(rp), public, allocatable :: momz_tp(:,:,:)
117  real(rp), public, allocatable :: rhou_tp(:,:,:)
118  real(rp), public, allocatable :: rhov_tp(:,:,:)
119  real(rp), public, allocatable :: rhot_tp(:,:,:)
120  real(rp), public, allocatable :: rhoh_p (:,:,:)
121  real(rp), public, allocatable :: rhoq_tp(:,:,:,:)
122 
123  ! (obsolute)
124  real(rp), public, allocatable :: momx_tp(:,:,:)
125  real(rp), public, allocatable :: momy_tp(:,:,:)
126 
127 
128  ! public diagnostic variables
129  real(rp), public, allocatable, target :: w (:,:,:)
130  real(rp), public, allocatable, target :: u (:,:,:)
131  real(rp), public, allocatable, target :: v (:,:,:)
132 
133  real(rp), public, allocatable, target :: pott (:,:,:)
134  real(rp), public, allocatable, target :: temp (:,:,:)
135  real(rp), public, allocatable, target :: pres (:,:,:)
136  real(rp), public, allocatable, target :: exner(:,:,:)
137  real(rp), public, allocatable, target :: phyd (:,:,:)
138  real(rp), public, allocatable, target :: phydh(:,:,:)
139 
140  real(rp), public, allocatable, target :: qdry (:,:,:)
141  real(rp), public, allocatable, target :: rtot (:,:,:)
142  real(rp), public, allocatable, target :: cvtot(:,:,:)
143  real(rp), public, allocatable, target :: cptot(:,:,:)
144 
145  real(rp), public, allocatable, target :: prec (:,:)
146  real(rp), public, allocatable :: prec_engi(:,:)
147  !-----------------------------------------------------------------------------
148  !
149  !++ Private procedure
150  !
151  !-----------------------------------------------------------------------------
152  !
153  !++ Private parameters & variables
154  !
155  logical, private :: atmos_vars_checkrange = .false.
156  real(rp), private :: atmos_vars_checkcfl_soft = 1.0_rp
157  real(rp), private :: atmos_vars_checkcfl_hard = 2.0_rp
158 
159  type vinfo
160  character(len=H_SHORT) :: name
161  character(len=H_MID) :: desc
162  character(len=H_SHORT) :: unit
163  integer :: ndims
164  character(len=H_SHORT) :: dim_type
165  character(len=H_MID) :: stdname
166  end type vinfo
167 
168  ! prognostic variables
169  integer, private, parameter :: pv_nmax = 5
170  type(vinfo), private :: pv_info(pv_nmax)
171  integer, private, allocatable :: pv_id(:)
172 
173  data pv_info / &
174  vinfo( 'DENS', 'density', 'kg/m3', 3, 'ZXY', 'air_density' ), &
175  vinfo( 'MOMZ', 'momentum z', 'kg/m2/s', 3, 'ZHXY', 'upward_mass_flux_of_air' ), &
176  vinfo( 'MOMX', 'momentum x', 'kg/m2/s', 3, 'ZXHY', 'eastward_mass_flux_of_air' ), &
177  vinfo( 'MOMY', 'momentum y', 'kg/m2/s', 3, 'ZXYH', 'northward_mass_flux_of_air' ), &
178  vinfo( 'RHOT', 'rho * theta', 'kg/m3*K', 3, 'ZXY', '' ) /
179 
180 
181  ! private diagnostic variables
182  real(rp), allocatable, target :: lhv (:,:,:)
183  real(rp), allocatable, target :: lhs (:,:,:)
184  real(rp), allocatable, target :: lhf (:,:,:)
185 
186  real(rp), allocatable, target :: potv (:,:,:)
187  real(rp), allocatable, target :: teml (:,:,:)
188  real(rp), allocatable, target :: potl (:,:,:)
189  real(rp), allocatable, target :: pote (:,:,:)
190 
191  real(rp), allocatable, target :: qtot (:,:,:)
192  real(rp), allocatable, target :: qhyd (:,:,:)
193  real(rp), allocatable, target :: qliq (:,:,:)
194  real(rp), allocatable, target :: qice (:,:,:)
195 
196  real(rp), allocatable, target :: lwp (:,:)
197  real(rp), allocatable, target :: iwp (:,:)
198  real(rp), allocatable, target :: pw (:,:)
199 
200  real(rp), allocatable, target :: rain (:,:)
201  real(rp), allocatable, target :: snow (:,:)
202 
203  real(rp), allocatable, target :: qsat (:,:,:)
204  real(rp), allocatable, target :: rha (:,:,:)
205  real(rp), allocatable, target :: rhl (:,:,:)
206  real(rp), allocatable, target :: rhi (:,:,:)
207 
208  real(rp), allocatable, target :: vor (:,:,:)
209  real(rp), allocatable, target :: div (:,:,:)
210  real(rp), allocatable, target :: hdiv (:,:,:)
211  real(rp), allocatable, target :: uabs (:,:,:)
212 
213  real(rp), allocatable, target :: n2 (:,:,:)
214  real(rp), allocatable, target :: pblh (:,:)
215 
216  real(rp), allocatable, target :: mse (:,:,:)
217  real(rp), allocatable, target :: tdew (:,:,:)
218 
219  real(rp), allocatable, target :: cape (:,:)
220  real(rp), allocatable, target :: cin (:,:)
221  real(rp), allocatable, target :: lcl (:,:)
222  real(rp), allocatable, target :: lfc (:,:)
223  real(rp), allocatable, target :: lnb (:,:)
224 
225  real(rp), allocatable, target :: engt (:,:,:)
226  real(rp), allocatable, target :: engp (:,:,:)
227  real(rp), allocatable, target :: engk (:,:,:)
228  real(rp), allocatable, target :: engi (:,:,:)
229 
230  real(rp), allocatable, target :: dens_mean(:)
231  real(rp), allocatable, target :: w_mean (:)
232  real(rp), allocatable, target :: u_mean (:)
233  real(rp), allocatable, target :: v_mean (:)
234  real(rp), allocatable, target :: pt_mean (:)
235  real(rp), allocatable, target :: t_mean (:)
236  real(rp), allocatable, target :: qv_mean (:)
237  real(rp), allocatable, target :: qhyd_mean(:)
238  real(rp), allocatable, target :: qliq_mean(:)
239  real(rp), allocatable, target :: qice_mean(:)
240 
241  real(rp), allocatable, target :: dens_prim(:,:,:)
242  real(rp), allocatable, target :: w_prim (:,:,:)
243  real(rp), allocatable, target :: u_prim (:,:,:)
244  real(rp), allocatable, target :: v_prim (:,:,:)
245  real(rp), allocatable, target :: pt_prim (:,:,:)
246  real(rp), allocatable, target :: w_prim2 (:,:,:)
247  real(rp), allocatable, target :: pt_w_prim(:,:,:)
248  real(rp), allocatable, target :: w_prim3 (:,:,:)
249  real(rp), allocatable, target :: tke_rs (:,:,:)
250 
251  real(rp), allocatable, target :: velz (:,:,:)
252  real(rp), allocatable, target :: velx (:,:,:)
253  real(rp), allocatable, target :: vely (:,:,:)
254  real(rp), allocatable, target :: umet (:,:,:)
255  real(rp), allocatable, target :: vmet (:,:,:)
256 
257  ! id of diagnostic variables
258  !! public
259  integer, private, parameter :: i_w = 1
260  integer, private, parameter :: i_u = 2
261  integer, private, parameter :: i_v = 3
262  integer, private, parameter :: i_pott = 4
263  integer, private, parameter :: i_temp = 5
264  integer, private, parameter :: i_pres = 6
265  integer, private, parameter :: i_exner = 7
266  integer, private, parameter :: i_phyd = 8
267  integer, private, parameter :: i_qdry = 9
268  integer, private, parameter :: i_rtot = 10
269  integer, private, parameter :: i_cvtot = 11
270  integer, private, parameter :: i_cptot = 12
271  !! private
272  integer, private, parameter :: i_lhv = 13
273  integer, private, parameter :: i_lhs = 14
274  integer, private, parameter :: i_lhf = 15
275  integer, private, parameter :: i_potv = 16
276  integer, private, parameter :: i_teml = 17
277  integer, private, parameter :: i_potl = 18
278  integer, private, parameter :: i_pote = 19
279  integer, private, parameter :: i_qtot = 20
280  integer, private, parameter :: i_qhyd = 21
281  integer, private, parameter :: i_qliq = 22
282  integer, private, parameter :: i_qice = 23
283  integer, private, parameter :: i_lwp = 24
284  integer, private, parameter :: i_iwp = 25
285  integer, private, parameter :: i_pw = 26
286  integer, private, parameter :: i_prec = 27
287  integer, private, parameter :: i_rain = 28
288  integer, private, parameter :: i_snow = 29
289  integer, private, parameter :: i_qsat = 30
290  integer, private, parameter :: i_rha = 31
291  integer, private, parameter :: i_rhl = 32
292  integer, private, parameter :: i_rhi = 33
293  integer, private, parameter :: i_vor = 34
294  integer, private, parameter :: i_div = 35
295  integer, private, parameter :: i_hdiv = 36
296  integer, private, parameter :: i_uabs = 37
297  integer, private, parameter :: i_n2 = 38
298  integer, private, parameter :: i_pblh = 39
299  integer, private, parameter :: i_mse = 40
300  integer, private, parameter :: i_tdew = 41
301  integer, private, parameter :: i_cape = 42
302  integer, private, parameter :: i_cin = 43
303  integer, private, parameter :: i_lcl = 44
304  integer, private, parameter :: i_lfc = 45
305  integer, private, parameter :: i_lnb = 46
306  integer, private, parameter :: i_engt = 47
307  integer, private, parameter :: i_engp = 48
308  integer, private, parameter :: i_engk = 49
309  integer, private, parameter :: i_engi = 50
310  integer, private, parameter :: i_dens_mean = 51
311  integer, private, parameter :: i_w_mean = 52
312  integer, private, parameter :: i_u_mean = 53
313  integer, private, parameter :: i_v_mean = 54
314  integer, private, parameter :: i_pt_mean = 55
315  integer, private, parameter :: i_t_mean = 56
316  integer, private, parameter :: i_qv_mean = 57
317  integer, private, parameter :: i_qhyd_mean = 58
318  integer, private, parameter :: i_qliq_mean = 59
319  integer, private, parameter :: i_qice_mean = 60
320  integer, private, parameter :: i_dens_prim = 61
321  integer, private, parameter :: i_w_prim = 62
322  integer, private, parameter :: i_u_prim = 63
323  integer, private, parameter :: i_v_prim = 64
324  integer, private, parameter :: i_pt_prim = 65
325  integer, private, parameter :: i_w_prim2 = 66
326  integer, private, parameter :: i_pt_w_prim = 67
327  integer, private, parameter :: i_w_prim3 = 68
328  integer, private, parameter :: i_tke_rs = 69
329  integer, private, parameter :: i_velz = 70
330  integer, private, parameter :: i_velx = 71
331  integer, private, parameter :: i_vely = 72
332  integer, private, parameter :: i_umet = 73
333  integer, private, parameter :: i_vmet = 74
334 
335  integer, private, parameter :: dv_nmax = 74
336  type(vinfo), private :: dv_info(dv_nmax)
337  logical, private :: dv_calculated(dv_nmax)
338 
339  data dv_info / &
340  vinfo( 'W', 'velocity w', 'm/s', 3, 'ZXY', 'upward_air_velocity' ), &
341  vinfo( 'U', 'velocity u', 'm/s', 3, 'ZXY', 'x_wind' ), &
342  vinfo( 'V', 'velocity v', 'm/s', 3, 'ZXY', 'y_wind' ), &
343  vinfo( 'PT', 'potential temp.', 'K', 3, 'ZXY', 'air_potential_temperature' ), &
344  vinfo( 'T', 'temperature', 'K', 3, 'ZXY', 'air_temperature' ), &
345  vinfo( 'PRES', 'pressure', 'Pa', 3, 'ZXY', 'air_pressure' ), &
346  vinfo( 'EXNER', 'Exner function', '1', 3, 'ZXY', 'dimensionless_exner_function' ), &
347  vinfo( 'PHYD', 'hydrostatic pressure', 'Pa', 3, 'ZXY', '' ), &
348  vinfo( 'QDRY', 'dry air', 'kg/kg', 3, 'ZXY', '' ), &
349  vinfo( 'RTOT', 'Total gas constant', 'J/kg/K', 3, 'ZXY', '' ), &
350  vinfo( 'CVTOT', 'Total heat capacity', 'J/kg/K', 3, 'ZXY', '' ), &
351  vinfo( 'CPTOT', 'Total heat capacity', 'J/kg/K', 3, 'ZXY', '' ), &
352  vinfo( 'LHV', 'latent heat for vaporization', 'J/kg', 3, 'ZXY', '' ), &
353  vinfo( 'LHS', 'latent heat for sublimation', 'J/kg', 3, 'ZXY', '' ), &
354  vinfo( 'LHF', 'latent heat for fusion', 'J/kg', 3, 'ZXY', '' ), &
355  vinfo( 'POTV', 'virtual potential temp.', 'K', 3, 'ZXY', '' ), &
356  vinfo( 'TEML', 'liquid water temperature', 'K', 3, 'ZXY', '' ), &
357  vinfo( 'POTL', 'liquid water potential temp.', 'K', 3, 'ZXY', '' ), &
358  vinfo( 'POTE', 'equivalent potential temp.', 'K', 3, 'ZXY', 'pseudo_equivalent_potential_temperature' ), &
359  vinfo( 'QTOT', 'total water', 'kg/kg', 3, 'ZXY', 'mass_fraction_of_water_in_air' ), &
360  vinfo( 'QHYD', 'total hydrometeors', 'kg/kg', 3, 'ZXY', 'mass_fraction_of_cloud_condensed_water_in_air' ), &
361  vinfo( 'QLIQ', 'total liquid water', 'kg/kg', 3, 'ZXY', '' ), &
362  vinfo( 'QICE', 'total ice water', 'kg/kg', 3, 'ZXY', '' ), &
363  vinfo( 'LWP', 'liquid water path', 'g/m2', 2, 'XY', 'atmosphere_mass_content_of_cloud_liquid_water' ), &
364  vinfo( 'IWP', 'ice water path', 'g/m2', 2, 'XY', '' ), &
365  vinfo( 'PW', 'precipitable water', 'g/m2', 2, 'XY', 'atmosphere_mass_content_of_vapor' ), &
366  vinfo( 'PREC', 'surface precipitation flux', 'kg/m2/s', 2, 'XY', 'precipitation_flux' ), &
367  vinfo( 'RAIN', 'surface rain flux', 'kg/m2/s', 2, 'XY', 'rainfall_flux' ), &
368  vinfo( 'SNOW', 'surface snow flux', 'kg/m2/s', 2, 'XY', 'snowfall_flux' ), &
369  vinfo( 'QSAT', 'saturation specific humidity', 'kg/kg', 3, 'ZXY', '' ), &
370  vinfo( 'RHA', 'relative humidity(liq+ice)', '%', 3, 'ZXY', '' ), &
371  vinfo( 'RH', 'relative humidity(liq)', '%', 3, 'ZXY', 'relative_humidity' ), &
372  vinfo( 'RHI', 'relative humidity(ice)', '%', 3, 'ZXY', '' ), &
373  vinfo( 'VOR', 'vertical vorticity', '1/s', 3, 'ZXY', 'atmosphere_relative_vorticity' ), &
374  vinfo( 'DIV', 'divergence', '1/s', 3, 'ZXY', 'divergence_of_wind' ), &
375  vinfo( 'HDIV', 'horizontal divergence', '1/s', 3, 'ZXY', '' ), &
376  vinfo( 'Uabs', 'absolute velocity', 'm/s', 3, 'ZXY', 'wind_speed' ), &
377  vinfo( 'N2', 'squared Brunt-Vaisala frequency', '1/s2', 3, 'ZXY', 'square_of_brunt_vaisala_frequency_in_air' ), &
378  vinfo( 'PBLH', 'PBL height', 'm', 2, 'XY', 'atmosphere_boundary_layer_thickness' ), &
379  vinfo( 'MSE', 'moist static energy', 'm2/s2', 3, 'ZXY', '' ), &
380  vinfo( 'TDEW', 'dew point', 'K', 3, 'ZXY', 'dew_point_temperature' ), &
381  vinfo( 'CAPE', 'convective avail. pot. energy', 'm2/s2', 2, 'XY', 'atmosphere_specific_convective_available_potential_energy' ), &
382  vinfo( 'CIN', 'convection inhibition', 'm2/s2', 2, 'XY', '' ), &
383  vinfo( 'LCL', 'lifted condensation level', 'm', 2, 'XY', 'atmosphere_lifting_condensation_level' ), &
384  vinfo( 'LFC', 'level of free convection', 'm', 2, 'XY', 'atmosphere_level_of_free_convection' ), &
385  vinfo( 'LNB', 'level of neutral buoyancy', 'm', 2, 'XY', '' ), &
386  vinfo( 'ENGT', 'total energy', 'J/m3', 3, 'ZXY', '' ), &
387  vinfo( 'ENGP', 'potential energy', 'J/m3', 3, 'ZXY', '' ), &
388  vinfo( 'ENGK', 'kinetic energy', 'J/m3', 3, 'ZXY', '' ), &
389  vinfo( 'ENGI', 'internal energy', 'J/m3', 3, 'ZXY', '' ), &
390  vinfo( 'DENS_MEAN', 'horiz. mean of density', 'kg/m3', 1, 'Z', '' ), &
391  vinfo( 'W_MEAN', 'horiz. mean of w', 'm/s', 1, 'Z', '' ), &
392  vinfo( 'U_MEAN', 'horiz. mean of u', 'm/s', 1, 'Z', '' ), &
393  vinfo( 'V_MEAN', 'horiz. mean of v', 'm/s', 1, 'Z', '' ), &
394  vinfo( 'PT_MEAN', 'horiz. mean of pot.', 'K', 1, 'Z', '' ), &
395  vinfo( 'T_MEAN', 'horiz. mean of t', 'K', 1, 'Z', '' ), &
396  vinfo( 'QV_MEAN', 'horiz. mean of QV', '1', 1, 'Z', '' ), &
397  vinfo( 'QHYD_MEAN', 'horiz. mean of QHYD', '1', 1, 'Z', '' ), &
398  vinfo( 'QLIQ_MEAN', 'horiz. mean of QLIQ', '1', 1, 'Z', '' ), &
399  vinfo( 'QICE_MEAN', 'horiz. mean of QICE', '1', 1, 'Z', '' ), &
400  vinfo( 'DENS_PRIM', 'horiz. deviation of density', 'kg/m3', 3, 'ZXY', '' ), &
401  vinfo( 'W_PRIM', 'horiz. deviation of w', 'm/s', 3, 'ZXY', '' ), &
402  vinfo( 'U_PRIM', 'horiz. deviation of u', 'm/s', 3, 'ZXY', '' ), &
403  vinfo( 'V_PRIM', 'horiz. deviation of v', 'm/s', 3, 'ZXY', '' ), &
404  vinfo( 'PT_PRIM', 'horiz. deviation of pot. temp.', 'K', 3, 'ZXY', '' ), &
405  vinfo( 'W_PRIM2', 'variance of w', 'm2/s2', 3, 'ZXY', '' ), &
406  vinfo( 'PT_W_PRIM', 'resolved scale heat flux', 'W/s', 3, 'ZXY', '' ), &
407  vinfo( 'W_PRIM3', 'skewness of w', 'm3/s3', 3, 'ZXY', '' ), &
408  vinfo( 'TKE_RS', 'resolved scale TKE', 'm2/s2', 3, 'ZXY', '' ), &
409  vinfo( 'VELZ', 'velocity w at the half level', 'm/s', 3, 'ZHXY','' ), &
410  vinfo( 'VELX', 'velocity u at the half level', 'm/s', 3, 'ZXHY','' ), &
411  vinfo( 'VELY', 'velocity v at the half level', 'm/s', 3, 'ZXYH','' ), &
412  vinfo( 'Umet', 'eastward velocity', 'm/s', 3, 'ZXY', 'eastward_wind' ), &
413  vinfo( 'Vmet', 'northward velocity', 'm/s', 3, 'ZXY', 'northward_wind' ) /
414 
415  ! for history output and monitor
416  integer, private :: pv_hist_id (pv_nmax)
417  integer, private :: pv_monit_id(pv_nmax)
418  integer, private, allocatable :: qp_hist_id (:)
419  integer, private, allocatable :: qp_monit_id(:)
420  integer, private :: dv_hist_id (dv_nmax)
421  integer, private :: hist_id_gph
422 
423  integer, private, parameter :: im_qdry = 1
424  integer, private, parameter :: im_qtot = 2
425  integer, private, parameter :: im_evap = 3
426  integer, private, parameter :: im_prec = 4
427  integer, private, parameter :: im_engt = 5
428  integer, private, parameter :: im_engp = 6
429  integer, private, parameter :: im_engk = 7
430  integer, private, parameter :: im_engi = 8
431  integer, private, parameter :: im_engflxt = 9
432  integer, private, parameter :: im_engsfc_sh = 10
433  integer, private, parameter :: im_engsfc_lh = 11
434  integer, private, parameter :: im_engsfc_evap = 12
435  integer, private, parameter :: im_engsfc_prec = 13
436  integer, private, parameter :: im_engsfc_rd = 14
437  integer, private, parameter :: im_engtom_rd = 15
438  integer, private, parameter :: im_engsfc_lw_up = 16
439  integer, private, parameter :: im_engsfc_lw_dn = 17
440  integer, private, parameter :: im_engsfc_sw_up = 18
441  integer, private, parameter :: im_engsfc_sw_dn = 19
442  integer, private, parameter :: im_engtom_lw_up = 20
443  integer, private, parameter :: im_engtom_lw_dn = 21
444  integer, private, parameter :: im_engtom_sw_up = 22
445  integer, private, parameter :: im_engtom_sw_dn = 23
446  integer, private, parameter :: dvm_nmax = 23
447  integer, private :: dv_monit_id(dvm_nmax)
448 
449 
450  logical, private :: moist
451  real(rp), private, target, allocatable :: zero(:,:,:)
452 
453 
454  ! for restart
455  integer, private :: restart_fid = -1 ! file ID
456  logical, private :: atmos_restart_in_check_coordinates = .true.
457 
458 
459  real(rp), private, allocatable :: work3d(:,:,:)
460  real(rp), private, allocatable :: work2d(:,:)
461  real(rp), private, allocatable :: work1d(:)
462 
463  !-----------------------------------------------------------------------------
464 contains
465  !-----------------------------------------------------------------------------
467  subroutine atmos_vars_setup
468  use scale_const, only: &
469  undef => const_undef
470  use scale_prc, only: &
471  prc_abort
472  use scale_file_history, only: &
474  use scale_monitor, only: &
476  use scale_atmos_hydrometeor, only: &
478  n_hyd, &
479  i_qv, &
480  i_hc, &
481  i_hr, &
482  i_hi, &
483  i_hs, &
484  i_hg, &
485  i_hh
486  use mod_atmos_admin, only: &
488  use mod_atmos_dyn_vars, only: &
490  use mod_atmos_phy_mp_vars, only: &
492  use mod_atmos_phy_ae_vars, only: &
494  use mod_atmos_phy_ch_vars, only: &
496  use mod_atmos_phy_rd_vars, only: &
498  use mod_atmos_phy_sf_vars, only: &
500  use mod_atmos_phy_tb_vars, only: &
502  use mod_atmos_phy_bl_vars, only: &
504  use mod_atmos_phy_cp_vars, only: &
506  use mod_atmos_phy_lt_vars, only: &
508  implicit none
509 
510  namelist / param_atmos_vars / &
514  atmos_restart_in_check_coordinates, &
524  atmos_vars_checkrange, &
525  atmos_vars_checkcfl_soft, &
526  atmos_vars_checkcfl_hard
527 
528  integer :: ierr
529  integer :: iv, iq
530  !---------------------------------------------------------------------------
531 
532  log_newline
533  log_info("ATMOS_vars_setup",*) 'Setup'
534 
535  allocate( dens(ka,ia,ja) )
536  allocate( momz(ka,ia,ja) )
537  allocate( momx(ka,ia,ja) )
538  allocate( momy(ka,ia,ja) )
539  allocate( rhot(ka,ia,ja) )
540  allocate( qtrc(ka,ia,ja,max(qa,1)) )
541  !$acc enter data create(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC)
542 
543  if ( atmos_use_average ) then
544  allocate( dens_avw(ka,ia,ja) )
545  allocate( momz_avw(ka,ia,ja) )
546  allocate( momx_avw(ka,ia,ja) )
547  allocate( momy_avw(ka,ia,ja) )
548  allocate( rhot_avw(ka,ia,ja) )
549  allocate( qtrc_avw(ka,ia,ja,max(qa,1)) )
550  !$acc enter data create(DENS_avw, MOMZ_avw, MOMX_avw, MOMY_avw, RHOT_avw, QTRC_avw)
551 
552  dens_av => dens_avw
553  momz_av => momz_avw
554  momx_av => momx_avw
555  momy_av => momy_avw
556  rhot_av => rhot_avw
557  qtrc_av => qtrc_avw
558  else
559  dens_av => dens
560  momz_av => momz
561  momx_av => momx
562  momy_av => momy
563  rhot_av => rhot
564  qtrc_av => qtrc
565  endif
566 
567  allocate( dens_tp(ka,ia,ja) )
568  allocate( momz_tp(ka,ia,ja) )
569  allocate( rhou_tp(ka,ia,ja) )
570  allocate( rhov_tp(ka,ia,ja) )
571  allocate( rhot_tp(ka,ia,ja) )
572  allocate( rhoh_p(ka,ia,ja) )
573  allocate( rhoq_tp(ka,ia,ja,max(qa,1)) )
574  !$acc enter data create(DENS_tp, MOMZ_tp, RHOU_tp, RHOV_tp, RHOT_tp, RHOH_p, RHOQ_tp)
575 
576  allocate( w(ka,ia,ja) )
577  allocate( u(ka,ia,ja) )
578  allocate( v(ka,ia,ja) )
579  !$omp parallel workshare
580  w(:,:,:) = undef
581  u(:,:,:) = undef
582  v(:,:,:) = undef
583  !$omp end parallel workshare
584  !$acc enter data create(W, U, V)
585 
586  allocate( pott(ka,ia,ja) )
587  allocate( temp(ka,ia,ja) )
588  allocate( pres(ka,ia,ja) )
589  allocate( exner(ka,ia,ja) )
590  allocate( phyd(ka,ia,ja) )
591  allocate( phydh(0:ka,ia,ja) )
592  !$omp parallel workshare
593  pott(:,:,:) = undef
594  temp(:,:,:) = undef
595  pres(:,:,:) = undef
596  exner(:,:,:) = undef
597  phyd(:,:,:) = undef
598  phydh(:,:,:) = undef
599  !$omp end parallel workshare
600  !$acc enter data create(POTT, TEMP, PRES, EXNER, PHYD, PHYDH)
601 
602  allocate( qdry(ka,ia,ja) )
603  allocate( rtot(ka,ia,ja) )
604  allocate( cvtot(ka,ia,ja) )
605  allocate( cptot(ka,ia,ja) )
606  !$omp parallel workshare
607  qdry(:,:,:) = undef
608  rtot(:,:,:) = undef
609  cvtot(:,:,:) = undef
610  cptot(:,:,:) = undef
611  !$omp end parallel workshare
612  !$acc enter data create(Qdry, Rtot, CVtot, CPtot)
613 
614  allocate( prec(ia,ja) )
615  allocate( prec_engi(ia,ja) )
616  !$omp parallel workshare
617  prec(:,:) = undef
618  prec_engi(:,:) = undef
619  !$omp end parallel workshare
620  !$acc enter data create(PREC, PREC_ENGI)
621 
622  ! obsolute
623  allocate( momx_tp(ka,ia,ja) )
624  allocate( momy_tp(ka,ia,ja) )
625  !$acc enter data create(MOMX_tp, MOMY_tp)
626 
627 
628  !$omp parallel workshare
629  !$acc kernels
630  momz(1:ks-1,:,:) = 0.0_rp
631  !$acc end kernels
632  !$acc kernels
633  momz(ke:ka,:,:) = 0.0_rp
634  !$acc end kernels
635  !$omp end parallel workshare
636 
637  allocate( work3d(ka,ia,ja) )
638  allocate( work2d( ia,ja) )
639  allocate( work1d(ka ) )
640  !$acc enter data create(WORK3D, WORK2D, WORK1D)
641 
642 
643  !--- read namelist
644  rewind(io_fid_conf)
645  read(io_fid_conf,nml=param_atmos_vars,iostat=ierr)
646  if( ierr < 0 ) then !--- missing
647  log_info("ATMOS_vars_setup",*) 'Not found namelist. Default used.'
648  elseif( ierr > 0 ) then !--- fatal error
649  log_error("ATMOS_vars_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_VARS. Check!'
650  call prc_abort
651  endif
652  log_nml(param_atmos_vars)
653 
654  log_newline
655  log_info("ATMOS_vars_setup",*) 'List of prognostic variables (ATMOS) '
656  log_info_cont('(1x,A,A24,A,A48,A,A12,A)') &
657  ' |', 'VARNAME ','|', &
658  'DESCRIPTION ', '[', 'UNIT ', ']'
659  do iv = 1, pv_nmax
660  log_info_cont('(1x,A,I3,A,A24,A,A48,A,A12,A)') &
661  'NO.',iv,'|',pv_info(iv)%NAME,'|', pv_info(iv)%DESC,'[', pv_info(iv)%UNIT,']'
662  enddo
663  do iq = 1, qa
664  log_info_cont('(1x,A,I3,A,A24,A,A48,A,A12,A)') &
665  'NO.',5+iq,'|',tracer_name(iq),'|', tracer_desc(iq),'[', tracer_unit(iq),']'
666  enddo
667 
668  log_newline
669  if ( atmos_restart_in_basename /= '' ) then
670  log_info("ATMOS_vars_setup",*) 'Restart input? : YES, file = ', trim(atmos_restart_in_basename)
671  log_info("ATMOS_vars_setup",*) 'Add timelabel? : ', atmos_restart_in_postfix_timelabel
672  else
673  log_info("ATMOS_vars_setup",*) 'Restart input? : NO'
674  endif
675  if ( atmos_restart_output &
676  .AND. atmos_restart_out_basename /= '' ) then
677  log_info("ATMOS_vars_setup",*) 'Restart output? : YES, file = ', trim(atmos_restart_out_basename)
678  log_info("ATMOS_vars_setup",*) 'Add timelabel? : ', atmos_restart_out_postfix_timelabel
679  else
680  log_info("ATMOS_vars_setup",*) 'Restart output? : NO'
681  atmos_restart_output = .false.
682  endif
683 
684  if ( atmos_restart_check_basename == '' ) then
685  atmos_restart_check = .false.
686  endif
687 
688  if ( atmos_vars_checkcfl_hard > 0.0_rp ) then
689  atmos_vars_checkcfl_soft = min( atmos_vars_checkcfl_soft, atmos_vars_checkcfl_hard )
690  endif
691 
692  log_newline
693  log_info("ATMOS_vars_setup",*) 'Check restart consistency? : ', atmos_restart_check
694  log_info("ATMOS_vars_setup",*) 'Check value range of variables? : ', atmos_vars_checkrange
695  if ( atmos_vars_checkcfl_soft > 0.0_rp ) then
696  log_info("ATMOS_vars_setup",*) 'Threshold of Courant number to warn : ', atmos_vars_checkcfl_soft
697  else
698  log_info("ATMOS_vars_setup",*) 'Threshold of Courant number to warn : disabled'
699  endif
700  if ( atmos_vars_checkcfl_hard > 0.0_rp ) then
701  log_info("ATMOS_vars_setup",*) 'Threshold of Courant number to stop : ', atmos_vars_checkcfl_hard
702  else
703  log_info("ATMOS_vars_setup",*) 'Threshold of Courant number to stop : disabled'
704  endif
705 
716 
717 
718  ! water content
719  if ( atmos_hydrometeor_dry ) then
720  allocate( zero(ka,ia,ja) )
721  !$acc enter data create(ZERO)
722  !$omp parallel workshare
723  !$acc kernels
724 !OCL XFILL
725  zero(:,:,:) = 0.0_rp
726  !$acc end kernels
727  !$omp end parallel workshare
728 
729  qv => zero
730  qc => zero
731  qr => zero
732  qi => zero
733  qs => zero
734  qg => zero
735  qh => zero
736 
737  moist = .false.
738  else
739  allocate( qe(ka,ia,ja,n_hyd) )
740  !$omp parallel workshare
741 !OCL XFILL
742  qe(:,:,:,:) = undef
743  !$omp end parallel workshare
744  !$acc enter data create(Qe)
745 
746  qv => qtrc_av(:,:,:,i_qv)
747  qc => qe(:,:,:,i_hc)
748  qr => qe(:,:,:,i_hr)
749  qi => qe(:,:,:,i_hi)
750  qs => qe(:,:,:,i_hs)
751  qg => qe(:,:,:,i_hg)
752  qh => qe(:,:,:,i_hh)
753 
754  moist = .true.
755  end if
756 
757 
758  dv_calculated(dv_nmax) = .false.
759 
760  !-----< history output setup >-----
761  allocate( qp_hist_id( max(qa,1) ) )
762  allocate( qp_monit_id( max(qa,1) ) )
763  pv_hist_id(:) = -1
764  pv_monit_id(:) = -1
765  qp_hist_id(:) = -1
766  qp_monit_id(:) = -1
767  dv_hist_id(:) = -1
768  dv_monit_id(:) = -1
769 
770 
771  do iv = 1, pv_nmax
772  call file_history_reg( pv_info(iv)%NAME, pv_info(iv)%DESC, pv_info(iv)%UNIT, pv_hist_id(iv), dim_type=pv_info(iv)%dim_type, standard_name=pv_info(iv)%STDNAME )
773  end do
774 
775  do iq = 1, qa
776  call file_history_reg( tracer_name(iq), tracer_desc(iq), tracer_unit(iq), qp_hist_id(iq), dim_type='ZXY' )
777  enddo
778 
779  do iv = 1, dv_nmax
780  call file_history_reg( dv_info(iv)%NAME, dv_info(iv)%DESC, dv_info(iv)%UNIT, dv_hist_id(iv), dim_type=dv_info(iv)%dim_type, standard_name=dv_info(iv)%STDNAME )
781  end do
782 
783  call file_history_reg( "GPH", "geopotential height", "m", hist_id_gph, dim_type='ZXY', standard_name="geopotential_height" )
784 
785 
786  !-----< monitor output setup >-----
787  do iv = 1, pv_nmax
788  call monitor_reg( pv_info(iv)%NAME, pv_info(iv)%DESC, trim(pv_info(iv)%UNIT)//"*m3", & ! (in)
789  pv_monit_id(iv), & ! (out)
790  dim_type=pv_info(iv)%dim_type, is_tendency=.false. ) ! (in)
791  end do
792  do iq = 1, qa
793  call monitor_reg( tracer_name(iq), tracer_desc(iq), tracer_unit(iq)//"*kg", & ! (in)
794  qp_monit_id(iq), & ! (out)
795  dim_type='ZXY', is_tendency=.false. ) ! (in)
796  enddo
797 
798  call monitor_reg( 'QDRY', 'dry air mass', 'kg', & ! (in)
799  dv_monit_id(im_qdry), & ! (out)
800  dim_type='ZXY', is_tendency=.false. ) ! (in)
801  call monitor_reg( 'QTOT', 'water mass', 'kg', & ! (in)
802  dv_monit_id(im_qtot), & ! (out)
803  dim_type='ZXY', is_tendency=.false. ) ! (in)
804  call monitor_reg( 'EVAP', 'evaporation at the surface', 'kg', & ! (in)
805  dv_monit_id(im_evap), & ! (out)
806  dim_type='XY', is_tendency=.true. ) ! (in)
807  call monitor_reg( 'PREC', 'precipitation', 'kg', & ! (in)
808  dv_monit_id(im_prec), & ! (out)
809  dim_type='XY', is_tendency=.true. ) ! (in)
810 
811  call monitor_reg( 'ENGT', 'total energy', 'J', & ! (in)
812  dv_monit_id(im_engt), & ! (out)
813  dim_type='ZXY', is_tendency=.false. ) ! (in)
814  call monitor_reg( 'ENGP', 'potential energy', 'J', & ! (in)
815  dv_monit_id(im_engp), & ! (out)
816  dim_type='ZXY', is_tendency=.false. ) ! (in)
817  call monitor_reg( 'ENGK', 'kinetic energy', 'J', & ! (in)
818  dv_monit_id(im_engk), & ! (out)
819  dim_type='ZXY', is_tendency=.false. ) ! (in)
820  call monitor_reg( 'ENGI', 'internal energy', 'J', & ! (in)
821  dv_monit_id(im_engi), & ! (out)
822  dim_type='ZXY', is_tendency=.false. ) ! (in)
823 
824  call monitor_reg( 'ENGFLXT', 'total energy flux convergence', 'J', & ! (in)
825  dv_monit_id(im_engflxt), & ! (out)
826  dim_type='XY', is_tendency=.true. ) ! (in)
827  call monitor_reg( 'ENGSFC_SH', 'SFC sensible heat flux', 'J', & ! (in)
828  dv_monit_id(im_engsfc_sh), & ! (out)
829  dim_type='XY', is_tendency=.true. ) ! (in)
830  call monitor_reg( 'ENGSFC_LH', 'SFC latent heat flux', 'J', & ! (in)
831  dv_monit_id(im_engsfc_lh), & ! (out)
832  dim_type='XY', is_tendency=.true. ) ! (in)
833  call monitor_reg( 'ENGSFC_EVAP', 'SFC internal energy flux of the evapolation', 'J', & ! (in)
834  dv_monit_id(im_engsfc_evap), & ! (out)
835  dim_type='XY', is_tendency=.true. ) ! (in)
836  call monitor_reg( 'ENGSFC_PREC', 'SFC internal energy flux of the precipitation', 'J', & ! (in)
837  dv_monit_id(im_engsfc_prec), & ! (out)
838  dim_type='XY', is_tendency=.true. ) ! (in)
839  call monitor_reg( 'ENGSFC_RD', 'SFC net radiation flux', 'J', & ! (in)
840  dv_monit_id(im_engsfc_rd), & ! (out)
841  dim_type='XY', is_tendency=.true. ) ! (in)
842  call monitor_reg( 'ENGTOM_RD', 'TOM net radiation flux', 'J', & ! (in)
843  dv_monit_id(im_engtom_rd), & ! (out)
844  dim_type='XY', is_tendency=.true. ) ! (in)
845 
846  call monitor_reg( 'ENGSFC_LW_up', 'SFC LW upward flux', 'J', & ! (in)
847  dv_monit_id(im_engsfc_lw_up), & ! (out)
848  dim_type='XY', is_tendency=.true. ) ! (in)
849  call monitor_reg( 'ENGSFC_LW_dn', 'SFC LW downward flux', 'J', & ! (in)
850  dv_monit_id(im_engsfc_lw_dn), & ! (out)
851  dim_type='XY', is_tendency=.true. ) ! (in)
852  call monitor_reg( 'ENGSFC_SW_up', 'SFC SW upward flux', 'J', & ! (in)
853  dv_monit_id(im_engsfc_sw_up), & ! (out)
854  dim_type='XY', is_tendency=.true. ) ! (in)
855  call monitor_reg( 'ENGSFC_SW_dn', 'SFC SW downward flux', 'J', & ! (in)
856  dv_monit_id(im_engsfc_sw_dn), & ! (out)
857  dim_type='XY', is_tendency=.true. ) ! (in)
858 
859  call monitor_reg( 'ENGTOM_LW_up', 'TOM LW upward flux', 'J', & ! (in)
860  dv_monit_id(im_engtom_lw_up), & ! (out)
861  dim_type='XY', is_tendency=.true. ) ! (in)
862  call monitor_reg( 'ENGTOM_LW_dn', 'TOM LW downward flux', 'J', & ! (in)
863  dv_monit_id(im_engtom_lw_dn), & ! (out)
864  dim_type='XY', is_tendency=.true. ) ! (in)
865  call monitor_reg( 'ENGTOM_SW_up', 'TOM SW upward flux', 'J', & ! (in)
866  dv_monit_id(im_engtom_sw_up), & ! (out)
867  dim_type='XY', is_tendency=.true. ) ! (in)
868  call monitor_reg( 'ENGTOM_SW_dn', 'TOM SW downward flux', 'J', & ! (in)
869  dv_monit_id(im_engtom_sw_dn), & ! (out)
870  dim_type='XY', is_tendency=.true. ) ! (in)
871 
872  return
873  end subroutine atmos_vars_setup
874 
875  !-----------------------------------------------------------------------------
877  subroutine atmos_vars_fillhalo( &
878  FILL_BND )
879  use scale_comm_cartesc, only: &
880  comm_vars8, &
881  comm_wait
882  implicit none
883 
884  logical, intent(in), optional :: fill_bnd
885 
886  logical :: fill_bnd_
887  integer :: i, j, iq
888  !---------------------------------------------------------------------------
889 
890  fill_bnd_ = .false.
891  if ( present(fill_bnd) ) fill_bnd_ = fill_bnd
892 
893 #ifdef QUICKDEBUG
894  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
895  !$acc kernels
896  do j = jsb, jeb
897  do i = isb, ieb
898  dens( 1:ks-1,i,j) = dens(ks,i,j)
899  momz( 1:ks-2,i,j) = momz(ks-1,i,j)
900  momx( 1:ks-1,i,j) = momx(ks,i,j)
901  momy( 1:ks-1,i,j) = momy(ks,i,j)
902  rhot( 1:ks-1,i,j) = rhot(ks,i,j)
903  dens(ke+1:ka, i,j) = dens(ke,i,j)
904  momz(ke+1:ka, i,j) = momz(ke,i,j)
905  momx(ke+1:ka, i,j) = momx(ke,i,j)
906  momy(ke+1:ka, i,j) = momy(ke,i,j)
907  rhot(ke+1:ka, i,j) = rhot(ke,i,j)
908  enddo
909  enddo
910  !$acc end kernels
911 
912  !$omp parallel do private(i,j,iq) OMP_SCHEDULE_ collapse(3)
913  !$acc kernels
914  do iq = 1, qa
915  do j = jsb, jeb
916  do i = isb, ieb
917  qtrc( 1:ks-1,i,j,iq) = qtrc(ks,i,j,iq)
918  qtrc(ke+1:ka, i,j,iq) = qtrc(ke,i,j,iq)
919  enddo
920  enddo
921  enddo
922  !$acc end kernels
923 #endif
924 
925  call comm_vars8( dens(:,:,:), 1 )
926  call comm_vars8( momz(:,:,:), 2 )
927  call comm_vars8( momx(:,:,:), 3 )
928  call comm_vars8( momy(:,:,:), 4 )
929  call comm_vars8( rhot(:,:,:), 5 )
930  call comm_wait ( dens(:,:,:), 1, fill_bnd_ )
931  call comm_wait ( momz(:,:,:), 2, fill_bnd_ )
932  call comm_wait ( momx(:,:,:), 3, fill_bnd_ )
933  call comm_wait ( momy(:,:,:), 4, fill_bnd_ )
934  call comm_wait ( rhot(:,:,:), 5, fill_bnd_ )
935 
936  do iq = 1, qa
937  call comm_vars8( qtrc(:,:,:,iq), iq )
938  enddo
939  do iq = 1, qa
940  call comm_wait ( qtrc(:,:,:,iq), iq, fill_bnd_ )
941  enddo
942 
943  return
944  end subroutine atmos_vars_fillhalo
945 
946  !-----------------------------------------------------------------------------
948  subroutine atmos_vars_restart_open
949  use scale_prc, only: &
950  prc_abort
951  use scale_const, only: &
952  grav => const_grav
953  use scale_time, only: &
955  use scale_file_cartesc, only: &
957  file_cartesc_check_coordinates
958  use mod_atmos_admin, only: &
960  atmos_sw_dyn, &
961  atmos_sw_phy_mp, &
962  atmos_sw_phy_ae, &
963  atmos_sw_phy_ch, &
964  atmos_sw_phy_rd, &
965  atmos_sw_phy_sf, &
966  atmos_sw_phy_tb, &
967  atmos_sw_phy_bl, &
968  atmos_sw_phy_cp, &
970  use mod_atmos_dyn_vars, only: &
972  use mod_atmos_phy_mp_vars, only: &
974  use mod_atmos_phy_ae_vars, only: &
976  use mod_atmos_phy_ch_vars, only: &
978  use mod_atmos_phy_rd_vars, only: &
980  use mod_atmos_phy_sf_vars, only: &
982  use mod_atmos_phy_tb_vars, only: &
984  use mod_atmos_phy_bl_vars, only: &
986  use mod_atmos_phy_cp_vars, only: &
988  use mod_atmos_phy_lt_vars, only: &
990  use mod_cpl_admin, only: &
991  cpl_sw
992  implicit none
993 
994  character(len=19) :: timelabel
995  character(len=H_LONG) :: basename
996  !---------------------------------------------------------------------------
997 
998  log_newline
999  log_info("ATMOS_vars_restart_open",*) 'Open restart file (ATMOS) '
1000 
1001  call prof_rapstart('ATM_Restart', 1)
1002 
1003  if ( atmos_restart_in_basename /= '' ) then
1004 
1006  call time_gettimelabel( timelabel )
1007  basename = trim(atmos_restart_in_basename)//'_'//trim(timelabel)
1008  else
1009  basename = trim(atmos_restart_in_basename)
1010  endif
1011 
1012  log_info("ATMOS_vars_restart_open",*) 'basename: ', trim(basename)
1013 
1014  call file_cartesc_open( basename, restart_fid, aggregate=atmos_restart_in_aggregate )
1015 
1016  if ( atmos_restart_in_check_coordinates ) then
1017  call file_cartesc_check_coordinates( restart_fid, atmos=.true. )
1018  end if
1019 
1020  else
1021  log_error("ATMOS_vars_restart_open",*) 'restart file for atmosphere is not specified. STOP!'
1022  call prc_abort
1023  endif
1024 
1025  if ( atmos_use_average ) then
1026  !$omp workshare
1027  !$acc kernels
1028  dens_av(:,:,:) = dens(:,:,:)
1029  !$acc end kernels
1030  !$acc kernels
1031  momz_av(:,:,:) = momz(:,:,:)
1032  !$acc end kernels
1033  !$acc kernels
1034  momx_av(:,:,:) = momx(:,:,:)
1035  !$acc end kernels
1036  !$acc kernels
1037  momy_av(:,:,:) = momy(:,:,:)
1038  !$acc end kernels
1039  !$acc kernels
1040  rhot_av(:,:,:) = rhot(:,:,:)
1041  !$acc end kernels
1042  !$acc kernels
1043  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
1044  !$acc end kernels
1045  !$omp end workshare
1046  endif
1047 
1053  if( atmos_sw_phy_sf .and. (.not. cpl_sw) ) call atmos_phy_sf_vars_restart_open
1058 
1059  call prof_rapend('ATM_Restart', 1)
1060 
1061  return
1062  end subroutine atmos_vars_restart_open
1063 
1064  !-----------------------------------------------------------------------------
1066  subroutine atmos_vars_restart_read
1067  use scale_prc, only: &
1068  prc_abort
1069  use scale_file, only: &
1071  use scale_file_cartesc, only: &
1072  file_cartesc_read, &
1074  use mod_atmos_admin, only: &
1076  atmos_sw_dyn, &
1077  atmos_sw_phy_mp, &
1078  atmos_sw_phy_ae, &
1079  atmos_sw_phy_ch, &
1080  atmos_sw_phy_rd, &
1081  atmos_sw_phy_sf, &
1082  atmos_sw_phy_tb, &
1083  atmos_sw_phy_bl, &
1084  atmos_sw_phy_cp, &
1086  use mod_atmos_dyn_vars, only: &
1088  use mod_atmos_phy_mp_vars, only: &
1090  use mod_atmos_phy_ae_vars, only: &
1092  use mod_atmos_phy_ch_vars, only: &
1094  use mod_atmos_phy_rd_vars, only: &
1096  use mod_atmos_phy_sf_vars, only: &
1098  use mod_atmos_phy_tb_vars, only: &
1100  use mod_atmos_phy_bl_vars, only: &
1102  use mod_atmos_phy_cp_vars, only: &
1104  use mod_atmos_phy_lt_vars, only: &
1106  use mod_cpl_admin, only: &
1107  cpl_sw
1108  implicit none
1109 
1110  integer :: i, j, iq
1111  !---------------------------------------------------------------------------
1112 
1113  call prof_rapstart('ATM_Restart', 1)
1114 
1115  if ( restart_fid /= -1 ) then
1116  log_newline
1117  log_info("ATMOS_vars_restart_read",*) 'Read from restart file (ATMOS) '
1118 
1119  call file_cartesc_read( restart_fid, pv_info(i_dens)%NAME, 'ZXY', & ! [IN]
1120  dens(:,:,:) ) ! [OUT]
1121  call file_cartesc_read( restart_fid, pv_info(i_momz)%NAME, 'ZHXY', & ! [IN]
1122  momz(:,:,:) ) ! [OUT]
1123  call file_cartesc_read( restart_fid, pv_info(i_momx)%NAME, 'ZXHY', & ! [IN]
1124  momx(:,:,:) ) ! [OUT]
1125  call file_cartesc_read( restart_fid, pv_info(i_momy)%NAME, 'ZXYH', & ! [IN]
1126  momy(:,:,:) ) ! [OUT]
1127  call file_cartesc_read( restart_fid, pv_info(i_rhot)%NAME, 'ZXY', & ! [IN]
1128  rhot(:,:,:) ) ! [OUT]
1129 
1130  do iq = 1, qa
1131  call file_cartesc_read( restart_fid, tracer_name(iq), 'ZXY', & ! [IN]
1132  qtrc(:,:,:,iq) ) ! [OUT]
1133  enddo
1134 
1135  if ( file_get_aggregate(restart_fid) ) then
1136  call file_cartesc_flush( restart_fid ) ! X/Y halos have been read from file
1137  !$acc update device(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC)
1138 #ifdef QUICKDEBUG
1139  ! fill k halos
1140  !$acc kernels
1141  do j = 1, ja
1142  do i = 1, ia
1143  dens( 1:ks-1,i,j) = dens(ks,i,j)
1144  momz( 1:ks-2,i,j) = momz(ks-1,i,j)
1145  momx( 1:ks-1,i,j) = momx(ks,i,j)
1146  momy( 1:ks-1,i,j) = momy(ks,i,j)
1147  rhot( 1:ks-1,i,j) = rhot(ks,i,j)
1148  dens(ke+1:ka, i,j) = dens(ke,i,j)
1149  momz(ke+1:ka, i,j) = momz(ke,i,j)
1150  momx(ke+1:ka, i,j) = momx(ke,i,j)
1151  momy(ke+1:ka, i,j) = momy(ke,i,j)
1152  rhot(ke+1:ka, i,j) = rhot(ke,i,j)
1153  do iq = 1, qa
1154  qtrc( 1:ks-1,i,j,iq) = qtrc(ks,i,j,iq)
1155  qtrc(ke+1:ka ,i,j,iq) = qtrc(ke,i,j,iq)
1156  end do
1157  enddo
1158  enddo
1159  !$acc end kernels
1160 #endif
1161  else
1162  call atmos_vars_fillhalo
1163  end if
1164 
1166  call atmos_vars_check( force = .true. )
1167  else
1168  log_error("ATMOS_vars_restart_read",*) 'invalid restart file ID for atmosphere. STOP!'
1169  call prc_abort
1170  endif
1171 
1172  if ( atmos_use_average ) then
1173  !$omp parallel workshare
1174  !$acc kernels
1175  dens_av(:,:,:) = dens(:,:,:)
1176  !$acc end kernels
1177  !$acc kernels
1178  momz_av(:,:,:) = momz(:,:,:)
1179  !$acc end kernels
1180  !$acc kernels
1181  momx_av(:,:,:) = momx(:,:,:)
1182  !$acc end kernels
1183  !$acc kernels
1184  momy_av(:,:,:) = momy(:,:,:)
1185  !$acc end kernels
1186  !$acc kernels
1187  rhot_av(:,:,:) = rhot(:,:,:)
1188  !$acc end kernels
1189  !$acc kernels
1190  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
1191  !$acc end kernels
1192  !$omp end parallel workshare
1193  endif
1194 
1200  if ( atmos_sw_phy_sf .and. (.not. cpl_sw) ) call atmos_phy_sf_vars_restart_read
1205 
1206  call prof_rapend('ATM_Restart', 1)
1207 
1208  return
1209  end subroutine atmos_vars_restart_read
1210 
1211  !-----------------------------------------------------------------------------
1213  subroutine atmos_vars_history_setpres
1215  real_fz => atmos_grid_cartesc_real_fz
1216  use scale_atmos_bottom, only: &
1217  bottom_estimate => atmos_bottom_estimate
1218  use scale_file_history_cartesc, only: &
1220  use mod_atmos_phy_sf_vars, only: &
1221  sfc_temp => atmos_phy_sf_sfc_temp
1222  implicit none
1223 
1224  real(rp) :: sfc_dens(ia,ja)
1225  real(rp) :: sfc_pres(ia,ja)
1226  !---------------------------------------------------------------------------
1227 
1228  call prof_rapstart('ATM_History', 1)
1229 
1230  !$acc data create(SFC_TEMP, SFC_DENS, SFC_PRES)
1231  call bottom_estimate( ka, ks, ke, ia, isb, ieb, ja, jsb, jeb, &
1232  dens_av(:,:,:), pres(:,:,:), qv(:,:,:), & ! [IN]
1233  sfc_temp(:,:), & ! [IN]
1234  real_fz(:,:,:), & ! [IN]
1235  sfc_dens(:,:), sfc_pres(:,:) ) ! [OUT]
1236 
1237  call file_history_cartesc_set_pres( phyd(:,:,:), & ! [IN]
1238  phydh(:,:,:), & ! [IN]
1239  sfc_pres(:,:) ) ! [IN]
1240 
1241  !$acc end data
1242 
1243  call prof_rapend('ATM_History', 1)
1244 
1245  return
1246  end subroutine atmos_vars_history_setpres
1247 
1248  !-----------------------------------------------------------------------------
1250  subroutine atmos_vars_restart_check
1251  use scale_prc, only: &
1252  prc_myrank
1253  use scale_file, only: &
1255  use scale_file_cartesc, only: &
1257  file_cartesc_read, &
1260  implicit none
1261 
1262  real(rp) :: dens_check(ka,ia,ja) ! Density [kg/m3]
1263  real(rp) :: momz_check(ka,ia,ja) ! momentum z [kg/s/m2]
1264  real(rp) :: momx_check(ka,ia,ja) ! momentum x [kg/s/m2]
1265  real(rp) :: momy_check(ka,ia,ja) ! momentum y [kg/s/m2]
1266  real(rp) :: rhot_check(ka,ia,ja) ! DENS * POTT [K*kg/m3]
1267  real(rp) :: qtrc_check(ka,ia,ja,qa) ! tracer mixing ratio [kg/kg]
1268 
1269  character(len=H_LONG) :: basename
1270 
1271  logical :: datacheck
1272  integer :: k, i, j, iq
1273  integer :: fid
1274  !---------------------------------------------------------------------------
1275 
1276  call prof_rapstart('Debug')
1277 
1278  log_info("ATMOS_vars_restart_check",*) 'Compare last Data with ', trim(atmos_restart_check_basename), 'on PE=', prc_myrank
1279  log_info("ATMOS_vars_restart_check",*) 'criterion = ', atmos_restart_check_criterion
1280  datacheck = .true.
1281 
1282  basename = atmos_restart_check_basename
1283 
1284  call file_cartesc_open( basename, fid )
1285 
1286  call file_cartesc_read( fid, 'DENS', 'ZXY' , dens_check(:,:,:) )
1287  call file_cartesc_read( fid, 'MOMZ', 'ZHXY', momz_check(:,:,:) )
1288  call file_cartesc_read( fid, 'MOMX', 'ZXHY', momx_check(:,:,:) )
1289  call file_cartesc_read( fid, 'MOMY', 'ZXYH', momy_check(:,:,:) )
1290  call file_cartesc_read( fid, 'RHOT', 'ZXY' , rhot_check(:,:,:) )
1291  do iq = 1, qa
1292  call file_cartesc_read( fid, tracer_name(iq), 'ZXY', qtrc_check(:,:,:,iq) )
1293  end do
1294  if ( file_get_aggregate(fid) ) call file_cartesc_flush( fid )
1295 
1296  call file_cartesc_close( fid ) ! [IN]
1297 
1298  !$acc update host(DENS)
1299  do k = ks, ke
1300  do j = js, je
1301  do i = is, ie
1302  if ( abs( dens(k,i,j)-dens_check(k,i,j) ) > atmos_restart_check_criterion ) then
1303  log_error("ATMOS_vars_restart_check",*) 'there is the difference : ', dens(k,i,j)-dens_check(k,i,j)
1304  log_error_cont(*) 'at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'DENS'
1305  datacheck = .false.
1306  endif
1307  enddo
1308  enddo
1309  enddo
1310 
1311  !$acc update host(MOMZ)
1312  do k = ks-1, ke
1313  do j = js, je
1314  do i = is, ie
1315  if ( abs( momz(k,i,j)-momz_check(k,i,j) ) > atmos_restart_check_criterion ) then
1316  log_error("ATMOS_vars_restart_check",*) 'there is the difference : ', momz(k,i,j)-momz_check(k,i,j)
1317  log_error_cont(*) 'at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'MOMZ'
1318  datacheck = .false.
1319  endif
1320  enddo
1321  enddo
1322  enddo
1323 
1324  !$acc update host(MOMX)
1325  do k = ks, ke
1326  do j = js, je
1327  do i = is, ie
1328  if ( abs( momx(k,i,j)-momx_check(k,i,j) ) > atmos_restart_check_criterion ) then
1329  log_error("ATMOS_vars_restart_check",*) 'there is the difference : ', momx(k,i,j)-momx_check(k,i,j)
1330  log_error_cont(*) 'at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'MOMX'
1331  datacheck = .false.
1332  endif
1333  enddo
1334  enddo
1335  enddo
1336 
1337  !$acc update host(MOMY)
1338  do k = ks, ke
1339  do j = js, je
1340  do i = is, ie
1341  if ( abs( momy(k,i,j)-momy_check(k,i,j) ) > atmos_restart_check_criterion ) then
1342  log_error("ATMOS_vars_restart_check",*) 'there is the difference : ', momy(k,i,j)-momy_check(k,i,j)
1343  log_error_cont(*) 'at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'MOMY'
1344  datacheck = .false.
1345  endif
1346  enddo
1347  enddo
1348  enddo
1349 
1350  !$acc update host(RHOT)
1351  do k = ks, ke
1352  do j = js, je
1353  do i = is, ie
1354  if ( abs( rhot(k,i,j)-rhot_check(k,i,j) ) > atmos_restart_check_criterion ) then
1355  log_error("ATMOS_vars_restart_check",*) 'there is the difference : ', rhot(k,i,j)-rhot_check(k,i,j)
1356  log_error_cont(*) 'at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'RHOT'
1357  datacheck = .false.
1358  endif
1359  enddo
1360  enddo
1361  enddo
1362 
1363  !$acc update host(QTRC)
1364  do iq = 1, qa
1365  do k = ks, ke
1366  do j = js, je
1367  do i = is, ie
1368  if ( abs( qtrc(k,i,j,iq)-qtrc_check(k,i,j,iq) ) > atmos_restart_check_criterion ) then
1369  log_error("ATMOS_vars_restart_check",*) 'there is the difference : ', qtrc(k,i,j,iq)-qtrc_check(k,i,j,iq)
1370  log_error_cont(*) 'at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, tracer_name(iq)
1371  datacheck = .false.
1372  endif
1373  enddo
1374  enddo
1375  enddo
1376  enddo
1377 
1378  if (datacheck) then
1379  log_info("ATMOS_vars_restart_check",*) 'Data Check Clear.'
1380  else
1381  log_info("ATMOS_vars_restart_check",*) 'Data Check Failed. See std. output.'
1382  log_error("ATMOS_vars_restart_check",*) 'Data Check Failed.'
1383  endif
1384 
1385  call prof_rapend('Debug')
1386 
1387  return
1388  end subroutine atmos_vars_restart_check
1389 
1390  !-----------------------------------------------------------------------------
1392  subroutine atmos_vars_history
1393  use scale_file_history, only: &
1394  file_history_query, &
1395  file_history_put
1396  use scale_atmos_grid_cartesc_real, only: &
1397  real_cz => atmos_grid_cartesc_real_cz
1398  use mod_atmos_phy_mp_vars, only: &
1400  use mod_atmos_phy_ae_vars, only: &
1402  implicit none
1403 
1404  logical :: do_put
1405  integer :: iq, iv
1406  !---------------------------------------------------------------------------
1407 
1408  call prof_rapstart('ATM_History', 1)
1409 
1410  ! history output of prognostic variables
1411  call file_history_put ( pv_hist_id(i_dens), dens(:,:,:) )
1412  call file_history_put ( pv_hist_id(i_momz), momz(:,:,:) )
1413  call file_history_put ( pv_hist_id(i_momx), momx(:,:,:) )
1414  call file_history_put ( pv_hist_id(i_momy), momy(:,:,:) )
1415  call file_history_put ( pv_hist_id(i_rhot), rhot(:,:,:) )
1416  do iq = 1, qa
1417  call file_history_put ( qp_hist_id(iq), qtrc(:,:,:,iq) )
1418  enddo
1419 
1420 
1421  ! history output of diagnostic variables
1422  call file_history_put ( dv_hist_id(i_w ), w(:,:,:) )
1423  call file_history_put ( dv_hist_id(i_u ), u(:,:,:) )
1424  call file_history_put ( dv_hist_id(i_v ), v(:,:,:) )
1425  call file_history_put ( dv_hist_id(i_pott ), pott(:,:,:) )
1426  call file_history_put ( dv_hist_id(i_temp ), temp(:,:,:) )
1427  call file_history_put ( dv_hist_id(i_pres ), pres(:,:,:) )
1428 
1429  call file_history_put ( dv_hist_id(i_exner), exner(:,:,:) )
1430  call file_history_put ( dv_hist_id(i_phyd ), phyd(:,:,:) )
1431 
1432  call file_history_put ( dv_hist_id(i_qdry ), qdry(:,:,:) )
1433  call file_history_put ( dv_hist_id(i_rtot ), rtot(:,:,:) )
1434  call file_history_put ( dv_hist_id(i_cvtot), cvtot(:,:,:) )
1435  call file_history_put ( dv_hist_id(i_cptot), cptot(:,:,:) )
1436 
1437  do iv = i_cptot+1, dv_nmax
1438  if ( dv_hist_id(iv) > 0 ) then
1439  call file_history_query( dv_hist_id(iv), do_put )
1440 
1441  if ( do_put ) then
1442  select case( dv_info(iv)%ndims )
1443  case( 3 )
1444  call atmos_vars_get_diagnostic( dv_info(iv)%NAME, work3d(:,:,:) )
1445  call file_history_put( dv_hist_id(iv), work3d(:,:,:) )
1446  case( 2 )
1447  call atmos_vars_get_diagnostic( dv_info(iv)%NAME, work2d(:,:) )
1448  call file_history_put( dv_hist_id(iv), work2d(:,:) )
1449  case( 1 )
1450  call atmos_vars_get_diagnostic( dv_info(iv)%NAME, work1d(:) )
1451  call file_history_put( dv_hist_id(iv), work1d(:) )
1452  end select
1453  endif
1454  endif
1455  enddo
1456 
1457  call file_history_put( hist_id_gph, real_cz(:,:,:) )
1458 
1459 
1460  if ( moist ) &
1461  call atmos_phy_mp_vars_history( dens_av(:,:,:), temp(:,:,:), qtrc_av(:,:,:,:) )
1462 ! if ( .false. ) then
1463 ! call ATMOS_vars_get_diagnostic( "RH", WORK3D(:,:,:) )
1464 ! call ATMOS_PHY_AE_vars_history( QTRC_av(:,:,:,:), WORK3D(:,:,:) )
1465 ! end if
1466 
1467  call prof_rapend ('ATM_History', 1)
1468 
1469  return
1470  end subroutine atmos_vars_history
1471 
1472  !-----------------------------------------------------------------------------
1474  subroutine atmos_vars_check( force )
1475  use scale_prc, only: &
1476  prc_myrank, &
1477  prc_abort
1478  use scale_prc_cartesc, only: &
1479  prc_twod
1480  use scale_const, only: &
1481  grav => const_grav, &
1482  cvdry => const_cvdry
1483  use scale_statistics, only: &
1485  statistics_total, &
1486  statistics_detail
1487  use scale_atmos_grid_cartesc_real, only: &
1496  use mod_atmos_admin, only: &
1498  use scale_atmos_grid_cartesc, only: &
1499  rfdx => atmos_grid_cartesc_rfdx, &
1500  rfdy => atmos_grid_cartesc_rfdy
1501  use scale_atmos_grid_cartesc_real, only: &
1502  real_cz => atmos_grid_cartesc_real_cz
1503  use scale_atmos_grid_cartesc_metric, only: &
1505  use scale_time, only: &
1507  implicit none
1508  logical, intent(in), optional :: force
1509 
1510  real(rp) :: rhoq(ka,ia,ja)
1511 
1512  real(rp) :: work (ka,ia,ja,3)
1513  character(len=H_SHORT) :: wname(3)
1514  real(rp) :: cflmax
1515 
1516  integer :: k, i, j, iq
1517  logical :: check
1518  !---------------------------------------------------------------------------
1519 
1520  !$acc data create(RHOQ, WORK)
1521 
1522  if ( present(force) ) then
1523  check = force
1524  else
1525  check = atmos_vars_checkrange
1526  end if
1527 
1528  ! value check for prognostic variables
1529  if ( check ) then
1530  call valcheck( ka, ks, ke, ia, is, ie, ja, js, je, &
1531  dens(:,:,:), 0.0_rp, 2.0_rp, pv_info(i_dens)%NAME, __file__, __line__ )
1532  call valcheck( ka, ks, ke, ia, is, ie, ja, js, je, &
1533  momz(:,:,:), -200.0_rp, 200.0_rp, pv_info(i_momz)%NAME, __file__, __line__ )
1534  call valcheck( ka, ks, ke, ia, is, ie, ja, js, je, &
1535  momx(:,:,:), -200.0_rp, 200.0_rp, pv_info(i_momx)%NAME, __file__, __line__ )
1536  call valcheck( ka, ks, ke, ia, is, ie, ja, js, je, &
1537  momy(:,:,:), -200.0_rp, 200.0_rp, pv_info(i_momy)%NAME, __file__, __line__ )
1538  call valcheck( ka, ks, ke, ia, is, ie, ja, js, je, &
1539  rhot(:,:,:), 0.0_rp, 1000.0_rp, pv_info(i_rhot)%NAME, __file__, __line__ )
1540 
1541  !$omp parallel workshare
1542  !$acc kernels
1543 !OCL XFILL
1544  work(:,:,:,1) = w(:,:,:)
1545  !$acc end kernels
1546  !$acc kernels
1547 !OCL XFILL
1548  work(:,:,:,2) = u(:,:,:)
1549  !$acc end kernels
1550  !$acc kernels
1551 !OCL XFILL
1552  work(:,:,:,3) = v(:,:,:)
1553  !$acc end kernels
1554  !$omp end parallel workshare
1555 
1556  wname(1) = "W"
1557  wname(2) = "U"
1558  wname(3) = "V"
1559 
1560  call statistics_detail( ka, ks, ke, ia, is, ie, ja, js, je, 3, &
1561  wname(:), work(:,:,:,:) )
1562  endif
1563 
1564  if ( present(force) ) then
1565  check = force
1566  else
1568  end if
1569 
1570  if ( check ) then
1571 
1572  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1573  dens(:,:,:), pv_info(i_dens)%NAME, & ! (in)
1574  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1575  atmos_grid_cartesc_real_totvol ) ! (in)
1576  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1577  momz(:,:,:), pv_info(i_momz)%NAME, & ! (in)
1578  atmos_grid_cartesc_real_volwxy(:,:,:), & ! (in)
1579  atmos_grid_cartesc_real_totvolwxy ) ! (in)
1580  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1581  momx(:,:,:), pv_info(i_momx)%NAME, & ! (in)
1582  atmos_grid_cartesc_real_volzuy(:,:,:), & ! (in)
1583  atmos_grid_cartesc_real_totvolzuy ) ! (in)
1584  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1585  momy(:,:,:), pv_info(i_momy)%NAME, & ! (in)
1586  atmos_grid_cartesc_real_volzxv(:,:,:), & ! (in)
1587  atmos_grid_cartesc_real_totvolzxv ) ! (in)
1588  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1589  rhot(:,:,:), pv_info(i_rhot)%NAME, & ! (in)
1590  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1591  atmos_grid_cartesc_real_totvol ) ! (in)
1592 
1593  do iq = 1, qa
1594  !$acc kernels
1595  rhoq(ks:ke,is:ie,js:je) = dens(ks:ke,is:ie,js:je) * qtrc(ks:ke,is:ie,js:je,iq)
1596  !$acc end kernels
1597 
1598  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1599  rhoq(:,:,:), tracer_name(iq), & ! (in)
1600  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1601  atmos_grid_cartesc_real_totvol ) ! (in)
1602  enddo
1603 
1605 
1606 
1607  !$acc kernels
1608  rhoq(ks:ke,is:ie,js:je) = dens(ks:ke,is:ie,js:je) * qdry(ks:ke,is:ie,js:je)
1609  !$acc end kernels
1610  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1611  rhoq(:,:,:), 'QDRY', & ! (in)
1612  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1613  atmos_grid_cartesc_real_totvol ) ! (in)
1614 
1615  !$acc kernels
1616  rhoq(ks:ke,is:ie,js:je) = dens(ks:ke,is:ie,js:je) * ( 1.0_rp - qdry(ks:ke,is:ie,js:je) ) ! Qtotal
1617  !$acc end kernels
1618  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1619  rhoq(:,:,:), 'QTOT', & ! (in)
1620  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1621  atmos_grid_cartesc_real_totvol ) ! (in)
1622 
1623 
1624  call atmos_vars_get_diagnostic( 'ENGT', work3d(:,:,:) )
1625  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1626  work3d(:,:,:), 'ENGT', & ! (in)
1627  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1628  atmos_grid_cartesc_real_totvol ) ! (in)
1629  call atmos_vars_get_diagnostic( 'ENGP', work3d(:,:,:) )
1630  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1631  work3d(:,:,:), 'ENGP', & ! (in)
1632  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1633  atmos_grid_cartesc_real_totvol ) ! (in)
1634  call atmos_vars_get_diagnostic( 'ENGK', work3d(:,:,:) )
1635  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1636  work3d(:,:,:), 'ENGK', & ! (in)
1637  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1638  atmos_grid_cartesc_real_totvol ) ! (in)
1639  call atmos_vars_get_diagnostic( 'ENGI', work3d(:,:,:) )
1640  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1641  work3d(:,:,:), 'ENGI', & ! (in)
1642  atmos_grid_cartesc_real_vol(:,:,:), & ! (in)
1643  atmos_grid_cartesc_real_totvol ) ! (in)
1644 
1645  end if
1646 
1647 
1648  ! CFL condition check
1649  if ( ( atmos_dyn_type /= 'OFF' .AND. atmos_dyn_type /= 'NONE' ) &
1650  .AND. ( atmos_vars_checkcfl_soft > 0.0_rp .OR. atmos_vars_checkcfl_hard > 0.0_rp ) ) then
1651  !$omp parallel workshare
1652  !$acc kernels
1653 !OCL XFILL
1654  work(:,:,:,:) = 0.0_rp
1655  !$acc end kernels
1656  !$omp end parallel workshare
1657 
1658  !$omp parallel do
1659  !$acc kernels
1660  do j = js, je
1661  do i = is, ie
1662  do k = ks, ke
1663  work(k,i,j,1) = 0.5_rp * abs(momz_av(k,i,j)) / ( dens_av(k+1,i,j) + dens_av(k,i,j) ) &
1664  * time_dtsec_atmos_dyn / ( real_cz(k+1,i,j) - real_cz(k,i,j) )
1665  work(k,i,j,3) = 0.5_rp * abs(momy_av(k,i,j)) / ( dens_av(k,i,j+1) + dens_av(k,i,j) ) &
1666  * time_dtsec_atmos_dyn * rfdy(j) * mapf(i,j,2,i_xv)
1667  enddo
1668  enddo
1669  enddo
1670  !$acc end kernels
1671  if ( prc_twod ) then
1672  !$omp parallel do
1673  !$acc kernels
1674  do j = js, je
1675  do k = ks, ke
1676  work(k,is,j,2) = 0.0_rp
1677  enddo
1678  enddo
1679  !$acc end kernels
1680  else
1681  !$omp parallel do
1682  !$acc kernels
1683  do j = js, je
1684  do i = is, ie
1685  do k = ks, ke
1686  work(k,i,j,2) = 0.5_rp * abs(momx_av(k,i,j)) / ( dens_av(k,i+1,j) + dens_av(k,i,j) ) &
1687  * time_dtsec_atmos_dyn * rfdx(i) * mapf(i,j,1,i_uy)
1688  enddo
1689  enddo
1690  enddo
1691  !$acc end kernels
1692  end if
1693 
1694  cflmax = -999e10_rp
1695  !$omp parallel do collapse(3) reduction(max:CFLMAX)
1696  !$acc kernels loop reduction(max:CFLMAX)
1697  do iq = 1, 3
1698  do j = 1, ja
1699  do i = 1, ia
1700  do k = 1, ka
1701  cflmax = max( cflmax, work(k,i,j,iq) )
1702  end do
1703  end do
1704  end do
1705  end do
1706  !$acc end kernels
1707 
1708  if ( atmos_vars_checkcfl_hard > 0.0_rp .AND. cflmax > atmos_vars_checkcfl_hard ) then
1709  log_info("ATMOS_vars_check",*) "Courant number =", cflmax, " exceeded the hard limit =", atmos_vars_checkcfl_hard
1710  log_error("ATMOS_vars_check",*)"Courant number =", cflmax, " exceeded the hard limit =", atmos_vars_checkcfl_hard
1711  log_error_cont(*) "Rank =", prc_myrank
1712  log_error_cont(*) "Please set ATMOS_VARS_CHECKCFL_HARD in the namelist PARAM_ATMOS_VARS when you want to change the limit."
1713 
1714  wname(1) = "Courant num. Z"
1715  wname(2) = "Courant num. X"
1716  wname(3) = "Courant num. Y"
1717  call statistics_detail( ka, ks, ke, ia, is, ie, ja, js, je, 3, &
1718  wname(:), work(:,:,:,:), &
1719  local=.true. )
1720 
1721  call prc_abort
1722  endif
1723 
1724  if ( atmos_vars_checkcfl_soft > 0.0_rp .AND. cflmax > atmos_vars_checkcfl_soft ) then
1725  log_info("ATMOS_vars_check",*) "Courant number =", cflmax, " exceeded the soft limit =", atmos_vars_checkcfl_soft
1726  log_error("ATMOS_vars_check",*) "Courant number =", cflmax, " exceeded the soft limit =", atmos_vars_checkcfl_soft
1727  log_error_cont(*) "Rank =", prc_myrank
1728 
1729  wname(1) = "Courant num. Z"
1730  wname(2) = "Courant num. X"
1731  wname(3) = "Courant num. Y"
1732  call statistics_detail( ka, ks, ke, ia, is, ie, ja, js, je, 3, &
1733  wname(:), work(:,:,:,:), &
1734  local=.true. )
1735  endif
1736 
1737  endif
1738 
1739  !$acc end data
1740 
1741  return
1742  end subroutine atmos_vars_check
1743 
1744  !-----------------------------------------------------------------------------
1746  subroutine atmos_vars_calc_diagnostics
1748  real_cz => atmos_grid_cartesc_real_cz, &
1749  real_fz => atmos_grid_cartesc_real_fz
1750  use scale_atmos_thermodyn, only: &
1751  atmos_thermodyn_specific_heat
1752  use scale_atmos_diagnostic, only: &
1755  use scale_atmos_diagnostic_cartesc, only: &
1757  use scale_comm_cartesc, only: &
1758  comm_vars8, &
1759  comm_wait
1760  use scale_atmos_hydrometeor, only: &
1761  n_hyd
1762  use mod_atmos_phy_mp_vars, only: &
1765  use mod_atmos_phy_ae_vars, only: &
1767  implicit none
1768 
1769  integer :: iq
1770 
1771  call prof_rapstart('ATM_Diag', 1)
1772 
1773  call atmos_thermodyn_specific_heat( &
1774  ka, ks, ke, ia, 1, ia, ja, 1, ja, qa, &
1775  qtrc_av(:,:,:,:), & ! (in)
1776  tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), & ! (in)
1777  qdry(:,:,:), rtot(:,:,:), cvtot(:,:,:), cptot(:,:,:) ) ! (out)
1778 
1780  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
1781  dens_av(:,:,:), momz_av(:,:,:), momx_av(:,:,:), momy_av(:,:,:), & ! (in)
1782  w(:,:,:), u(:,:,:), v(:,:,:) ) ! (out)
1783 
1785  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
1786  dens_av(:,:,:), rhot_av(:,:,:), & ! (in)
1787  rtot(:,:,:), cvtot(:,:,:), cptot(:,:,:), & ! (in)
1788  pott(:,:,:), temp(:,:,:), pres(:,:,:), exner(:,:,:) ) ! (out)
1789 
1791  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
1792  dens_av(:,:,:), pres(:,:,:), & ! (in)
1793  real_cz(:,:,:), real_fz(:,:,:), & ! (in)
1794  phyd(:,:,:), phydh(:,:,:) ) ! (out)
1795 
1798 
1799  if ( moist ) then
1801  dens_av(:,:,:), temp(:,:,:), qtrc_av(:,:,:,:), & ! (in)
1802  qe=qe(:,:,:,:) ) ! (out)
1803  do iq = 1, n_hyd
1804  call comm_vars8(qe(:,:,:,iq), iq)
1805  end do
1806  do iq = 1, n_hyd
1807  call comm_wait (qe(:,:,:,iq), iq)
1808  end do
1809  end if
1810 
1811  ! reset diagnostic variables
1812  dv_calculated(:) = .false.
1813 
1814  call prof_rapend('ATM_Diag', 1)
1815 
1816  return
1817  end subroutine atmos_vars_calc_diagnostics
1818 
1819  !-----------------------------------------------------------------------------
1821  recursive subroutine atmos_vars_get_diagnostic_3d( &
1822  vname, &
1823  var )
1824  use scale_const, only: &
1825  grav => const_grav, &
1826  rvap => const_rvap, &
1827  cpdry => const_cpdry, &
1828  cvdry => const_cvdry
1829  use scale_prc, only: &
1830  prc_abort
1831  use scale_prc_cartesc, only: &
1832  prc_twod
1833  use scale_atmos_grid_cartesc, only: &
1834  rcdx => atmos_grid_cartesc_rcdx, &
1835  rcdy => atmos_grid_cartesc_rcdy
1836  use scale_atmos_grid_cartesc_real, only: &
1837  real_cz => atmos_grid_cartesc_real_cz, &
1838  real_fz => atmos_grid_cartesc_real_fz, &
1840  use scale_atmos_grid_cartesc_metric, only: &
1842  use scale_comm_cartesc, only: &
1843  comm_vars8, &
1844  comm_wait
1845  use scale_atmos_hydrometeor, only: &
1846  lhvc => lhv, &
1847  lhfc => lhf, &
1848  atmos_hydrometeor_lhv, &
1849  atmos_hydrometeor_lhf, &
1850  atmos_hydrometeor_lhs
1851  use scale_atmos_saturation, only: &
1852  atmos_saturation_dens2qsat_all, &
1853  atmos_saturation_psat_all, &
1854  atmos_saturation_psat_liq, &
1855  atmos_saturation_psat_ice, &
1856  atmos_saturation_tdew_liq, &
1857  atmos_saturation_pote
1858  use scale_atmos_diagnostic, only: &
1862  implicit none
1863  character(len=*), intent(in) :: vname
1864  real(rp), intent(out) :: var(:,:,:)
1865 
1866  real(rp) :: uh (ka,ia,ja)
1867  real(rp) :: vh (ka,ia,ja)
1868 
1869  real(rp) :: work(ka,ia,ja)
1870 
1871  integer :: k, i, j, iq
1872 
1873  !$acc data copyout(var) create(UH, VH, WORK)
1874 
1875  select case ( vname )
1876  case ( 'W' )
1877  !$omp parallel workshare
1878  !$acc kernels
1879  var(:,:,:) = w(:,:,:)
1880  !$acc end kernels
1881  !$omp end parallel workshare
1882 
1883  case ( 'U' )
1884  !$omp parallel workshare
1885  !$acc kernels
1886  var(:,:,:) = u(:,:,:)
1887  !$acc end kernels
1888  !$omp end parallel workshare
1889 
1890  case ( 'V' )
1891  !$omp parallel workshare
1892  !$acc kernels
1893  var(:,:,:) = v(:,:,:)
1894  !$acc end kernels
1895  !$omp end parallel workshare
1896 
1897  case ( 'PT' )
1898  !$omp parallel workshare
1899  !$acc kernels
1900  var(:,:,:) = pott(:,:,:)
1901  !$acc end kernels
1902  !$omp end parallel workshare
1903 
1904  case ( 'T' )
1905  !$omp parallel workshare
1906  !$acc kernels
1907  var(:,:,:) = temp(:,:,:)
1908  !$acc end kernels
1909  !$omp end parallel workshare
1910 
1911  case ( 'EXNER' )
1912  !$omp parallel workshare
1913  !$acc kernels
1914  var(:,:,:) = exner(:,:,:)
1915  !$acc end kernels
1916  !$omp end parallel workshare
1917 
1918  case ( 'PHYD' )
1919  !$omp parallel workshare
1920  !$acc kernels
1921  var(:,:,:) = phyd(:,:,:)
1922  !$acc end kernels
1923  !$omp end parallel workshare
1924 
1925  case ( 'QDRY' )
1926  !$omp parallel workshare
1927  !$acc kernels
1928  var(:,:,:) = qdry(:,:,:)
1929  !$acc end kernels
1930  !$omp end parallel workshare
1931 
1932  case ( 'RTOT' )
1933  !$omp parallel workshare
1934  !$acc kernels
1935  var(:,:,:) = rtot(:,:,:)
1936  !$acc end kernels
1937  !$omp end parallel workshare
1938 
1939  case ( 'CVTOT' )
1940  !$omp parallel workshare
1941  !$acc kernels
1942  var(:,:,:) = cvtot(:,:,:)
1943  !$acc end kernels
1944  !$omp end parallel workshare
1945 
1946  case ( 'CPTOT' )
1947  !$omp parallel workshare
1948  !$acc kernels
1949  var(:,:,:) = cptot(:,:,:)
1950  !$acc end kernels
1951  !$omp end parallel workshare
1952 
1953  case ( 'LHV' )
1954  if ( .not. dv_calculated(i_lhv) ) then
1955  call allocate_3d( lhv )
1956  call atmos_hydrometeor_lhv( &
1957  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
1958  temp(:,:,:), & ! (in)
1959  lhv(:,:,:) ) ! (out)
1960  dv_calculated(i_lhv) = .true.
1961  end if
1962  !$omp parallel workshare
1963  !$acc kernels
1964  var(ks:ke,:,:) = lhv(ks:ke,:,:)
1965  !$acc end kernels
1966  !$omp end parallel workshare
1967 
1968  case ( 'LHS' )
1969  if ( .not. dv_calculated(i_lhs) ) then
1970  call allocate_3d( lhs )
1971  call atmos_hydrometeor_lhs( &
1972  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
1973  temp(:,:,:), & ! (in)
1974  lhs(:,:,:) ) ! (out)
1975  dv_calculated(i_lhs) = .true.
1976  end if
1977  !$omp parallel workshare
1978  !$acc kernels
1979  var(ks:ke,:,:) = lhs(ks:ke,:,:)
1980  !$acc end kernels
1981  !$omp end parallel workshare
1982 
1983  case ( 'LHF' )
1984  if ( .not. dv_calculated(i_lhf) ) then
1985  call allocate_3d( lhf )
1986  call atmos_hydrometeor_lhf( &
1987  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
1988  temp(:,:,:), & ! (in)
1989  lhf(:,:,:) ) ! (out)
1990  dv_calculated(i_lhf) = .true.
1991  end if
1992  !$omp parallel workshare
1993  !$acc kernels
1994  var(ks:ke,:,:) = lhf(ks:ke,:,:)
1995  !$acc end kernels
1996  !$omp end parallel workshare
1997 
1998  case ( 'POTV' )
1999  if ( .not. dv_calculated(i_potv) ) then
2000  call allocate_3d( potv )
2002  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2003  pott(:,:,:), rtot(:,:,:), & ! (in)
2004  potv(:,:,:) ) ! (out)
2005  dv_calculated(i_potv) = .true.
2006  end if
2007  !$omp parallel workshare
2008  !$acc kernels
2009  var(ks:ke,:,:) = potv(ks:ke,:,:)
2010  !$acc end kernels
2011  !$omp end parallel workshare
2012 
2013  case ( 'TEML' )
2014  if ( .not. dv_calculated(i_teml) ) then
2015  call allocate_3d( teml )
2016  call atmos_vars_get_diagnostic( 'LHV', work3d(:,:,:) )
2017  call atmos_vars_get_diagnostic( 'LHS', work3d(:,:,:) )
2018  call atmos_vars_get_diagnostic( 'QLIQ', work3d(:,:,:) )
2019  call atmos_vars_get_diagnostic( 'QICE', work3d(:,:,:) )
2021  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2022  temp(:,:,:), lhv(:,:,:), lhs(:,:,:), & ! (in)
2023  qc(:,:,:), qi(:,:,:), cptot(:,:,:), & ! (in)
2024  teml(:,:,:) ) ! (out)
2025  dv_calculated(i_teml) = .true.
2026  end if
2027  !$omp parallel workshare
2028  !$acc kernels
2029  var(ks:ke,:,:) = teml(ks:ke,:,:)
2030  !$acc end kernels
2031 
2032  !$omp end parallel workshare
2033 
2034  case ( 'POTL' )
2035  if ( .not. dv_calculated(i_potl) ) then
2036  call allocate_3d( potl )
2037  call atmos_vars_get_diagnostic( 'TEML', work3d(:,:,:) )
2038 !OCL XFILL
2039  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2040  !$omp private(i,j,k) &
2041  !$omp shared(POTL,TEML,EXNER) &
2042  !$omp shared(KS,KE,IA,JA)
2043  !$acc kernels
2044  !$acc loop collapse(3) independent
2045  do j = 1, ja
2046  do i = 1, ia
2047  do k = ks, ke
2048  potl(k,i,j) = teml(k,i,j) / exner(k,i,j)
2049  enddo
2050  enddo
2051  enddo
2052  !$acc end kernels
2053  dv_calculated(i_potl) = .true.
2054  end if
2055  !$omp parallel workshare
2056  !$acc kernels
2057  var(ks:ke,:,:) = potl(ks:ke,:,:)
2058  !$acc end kernels
2059  !$omp end parallel workshare
2060 
2061  case ( 'POTE' )
2062  if ( .not. dv_calculated(i_pote) ) then
2063  call allocate_3d( pote )
2064  call atmos_saturation_pote( &
2065  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2066  dens(:,:,:), pott(:,:,:), temp(:,:,:), qv(:,:,:), & ! [IN]
2067  pote(:,:,:) ) ! [OUT]
2068  end if
2069  !$omp parallel workshare
2070  !$acc kernels
2071  var(ks:ke,:,:) = pote(ks:ke,:,:)
2072  !$acc end kernels
2073  !$omp end parallel workshare
2074 
2075  case ( 'QTOT' )
2076  if ( .not. dv_calculated(i_qtot) ) then
2077  call allocate_3d( qtot )
2078  if ( moist ) then
2079  call atmos_vars_get_diagnostic( 'QHYD', work3d(:,:,:) )
2080 !OCL XFILL
2081  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2082  !$omp private(i,j,k) &
2083  !$omp shared(QTOT,QV,QHYD) &
2084  !$omp shared(KS,KE,IA,JA)
2085  !$acc kernels
2086  !$acc loop collapse(3) independent
2087  do j = 1, ja
2088  do i = 1, ia
2089  do k = ks, ke
2090  qtot(k,i,j) = qv(k,i,j) + qhyd(k,i,j)
2091  enddo
2092  enddo
2093  enddo
2094  !$acc end kernels
2095  else
2096 !OCL XFILL
2097  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2098  !$omp private(i,j,k) &
2099  !$omp shared(QTOT) &
2100  !$omp shared(KS,KE,IA,JA)
2101  !$acc kernels
2102  do j = 1, ja
2103  do i = 1, ia
2104  do k = ks, ke
2105  qtot(k,i,j) = 0.0_rp
2106  enddo
2107  enddo
2108  enddo
2109  !$acc end kernels
2110  end if
2111  dv_calculated(i_qtot) = .true.
2112  end if
2113  !$omp parallel workshare
2114  !$acc kernels
2115  var(ks:ke,:,:) = qtot(ks:ke,:,:)
2116  !$acc end kernels
2117  !$omp end parallel workshare
2118 
2119  case ( 'QHYD' )
2120  if ( .not. dv_calculated(i_qhyd) ) then
2121  call allocate_3d( qhyd )
2122  if ( moist ) then
2123  call atmos_vars_get_diagnostic( 'QLIQ', work3d(:,:,:) )
2124  call atmos_vars_get_diagnostic( 'QICE', work3d(:,:,:) )
2125 !OCL XFILL
2126  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2127  !$omp private(i,j,k) &
2128  !$omp shared(QHYD,QLIQ,QICE) &
2129  !$omp shared(KS,KE,IA,JA)
2130  !$acc kernels
2131  !$acc loop collapse(3) independent
2132  do j = 1, ja
2133  do i = 1, ia
2134  do k = ks, ke
2135  qhyd(k,i,j) = qliq(k,i,j) + qice(k,i,j)
2136  enddo
2137  enddo
2138  enddo
2139  !$acc end kernels
2140  else
2141 !OCL XFILL
2142  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2143  !$omp private(i,j,k) &
2144  !$omp shared(QHYD) &
2145  !$omp shared(KS,KE,IA,JA)
2146  !$acc kernels
2147  do j = 1, ja
2148  do i = 1, ia
2149  do k = ks, ke
2150  qhyd(k,i,j) = 0.0_rp
2151  enddo
2152  enddo
2153  enddo
2154  !$acc end kernels
2155  end if
2156  dv_calculated(i_qhyd) = .true.
2157  end if
2158  !$omp parallel workshare
2159  !$acc kernels
2160  var(ks:ke,:,:) = qhyd(ks:ke,:,:)
2161  !$acc end kernels
2162  !$omp end parallel workshare
2163 
2164  case ( 'QLIQ' )
2165  if ( .not. dv_calculated(i_qliq) ) then
2166  call allocate_3d( qliq )
2167 !OCL XFILL
2168  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2169  !$omp private(i,j,k) &
2170  !$omp shared(QLIQ,QC,QR) &
2171  !$omp shared(KS,KE,IA,JA)
2172  !$acc kernels
2173  !$acc loop collapse(3) independent
2174  do j = 1, ja
2175  do i = 1, ia
2176  do k = ks, ke
2177  qliq(k,i,j) = qc(k,i,j) + qr(k,i,j)
2178  enddo
2179  enddo
2180  enddo
2181  !$acc end kernels
2182  dv_calculated(i_qliq) = .true.
2183  end if
2184  !$omp parallel workshare
2185  !$acc kernels
2186  var(ks:ke,:,:) = qliq(ks:ke,:,:)
2187  !$acc end kernels
2188  !$omp end parallel workshare
2189 
2190  case ( 'QICE' )
2191  if ( .not. dv_calculated(i_qice) ) then
2192  call allocate_3d( qice )
2193 !OCL XFILL
2194  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2195  !$omp private(i,j,k) &
2196  !$omp shared(QICE,QI,QS,QG,QH) &
2197  !$omp shared(KS,KE,IA,JA)
2198  !$acc kernels
2199  !$acc loop collapse(3) independent
2200  do j = 1, ja
2201  do i = 1, ia
2202  do k = ks, ke
2203  qice(k,i,j) = qi(k,i,j) + qs(k,i,j) + qg(k,i,j) + qh(k,i,j)
2204  enddo
2205  enddo
2206  enddo
2207  !$acc end kernels
2208  dv_calculated(i_qice) = .true.
2209  end if
2210  !$omp parallel workshare
2211  !$acc kernels
2212  var(ks:ke,:,:) = qice(ks:ke,:,:)
2213  !$acc end kernels
2214  !$omp end parallel workshare
2215 
2216  case ( 'QSAT' )
2217  if ( .not. dv_calculated(i_qsat) ) then
2218  call allocate_3d( qsat )
2219  call atmos_saturation_dens2qsat_all( &
2220  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2221  temp(:,:,:), dens_av(:,:,:), & ! (in)
2222  qsat(:,:,:) ) ! (out)
2223  dv_calculated(i_qsat) = .true.
2224  end if
2225  !$omp parallel workshare
2226  !$acc kernels
2227  var(ks:ke,:,:) = qsat(ks:ke,:,:)
2228  !$acc end kernels
2229  !$omp end parallel workshare
2230 
2231  case ( 'RHA' )
2232  if ( .not. dv_calculated(i_rha) ) then
2233  call allocate_3d( rha )
2234  if ( moist ) then
2235  call atmos_saturation_psat_all( &
2236  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2237  temp(:,:,:), & ! (in)
2238  work(:,:,:) ) ! (out)
2239 !OCL XFILL
2240  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2241  !$omp private(i,j,k) &
2242  !$omp shared(RHA,DENS_av,QV,WORK,TEMP) &
2243  !$omp shared(KS,KE,IA,JA)
2244  !$acc kernels
2245  !$acc loop collapse(3) independent
2246  do j = 1, ja
2247  do i = 1, ia
2248  do k = ks, ke
2249  rha(k,i,j) = dens_av(k,i,j) * qv(k,i,j) &
2250  / work(k,i,j) * rvap * temp(k,i,j) &
2251  * 100.0_rp
2252  enddo
2253  enddo
2254  enddo
2255  !$acc end kernels
2256  else
2257 !OCL XFILL
2258  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2259  !$omp private(i,j,k) &
2260  !$omp shared(RHA) &
2261  !$omp shared(KS,KE,IA,JA)
2262  !$acc kernels
2263  do j = 1, ja
2264  do i = 1, ia
2265  do k = ks, ke
2266  rha(k,i,j) = 0.0_rp
2267  enddo
2268  enddo
2269  enddo
2270  !$acc end kernels
2271  end if
2272  dv_calculated(i_rha) = .true.
2273  end if
2274  !$omp parallel workshare
2275  !$acc kernels
2276  var(ks:ke,:,:) = rha(ks:ke,:,:)
2277  !$acc end kernels
2278  !$omp end parallel workshare
2279 
2280  case ( 'RHL', 'RH' )
2281  if ( .not. dv_calculated(i_rhl) ) then
2282  call allocate_3d( rhl )
2283  if ( moist ) then
2284  call atmos_saturation_psat_liq( &
2285  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2286  temp(:,:,:), & ! (in)
2287  work(:,:,:) ) ! (out)
2288 !OCL XFILL
2289  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2290  !$omp private(i,j,k) &
2291  !$omp shared(RHL,DENS_av,QV,WORK,TEMP) &
2292  !$omp shared(KS,KE,IA,JA)
2293  !$acc kernels
2294  !$acc loop collapse(3) independent
2295  do j = 1, ja
2296  do i = 1, ia
2297  do k = ks, ke
2298  rhl(k,i,j) = dens_av(k,i,j) * qv(k,i,j) &
2299  / work(k,i,j) * rvap * temp(k,i,j) &
2300  * 100.0_rp
2301  enddo
2302  enddo
2303  enddo
2304  !$acc end kernels
2305  else
2306 !OCL XFILL
2307  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2308  !$omp private(i,j,k) &
2309  !$omp shared(RHL) &
2310  !$omp shared(KS,KE,IA,JA)
2311  !$acc kernels
2312  do j = 1, ja
2313  do i = 1, ia
2314  do k = ks, ke
2315  rhl(k,i,j) = 0.0_rp
2316  enddo
2317  enddo
2318  enddo
2319  !$acc end kernels
2320  end if
2321  dv_calculated(i_rhl) = .true.
2322  end if
2323  !$omp parallel workshare
2324  !$acc kernels
2325  var(ks:ke,:,:) = rhl(ks:ke,:,:)
2326  !$acc end kernels
2327  !$omp end parallel workshare
2328 
2329  case ( 'RHI' )
2330  if ( .not. dv_calculated(i_rhi) ) then
2331  call allocate_3d( rhi )
2332  if ( moist ) then
2333  call atmos_saturation_psat_ice( &
2334  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2335  temp(:,:,:), & ! (int)
2336  work(:,:,:) ) ! (out)
2337 !OCL XFILL
2338  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2339  !$omp private(i,j,k) &
2340  !$omp shared(RHI,DENS_av,QV,WORK,TEMP) &
2341  !$omp shared(KS,KE,IA,JA)
2342  !$acc kernels
2343  !$acc loop collapse(3) independent
2344  do j = 1, ja
2345  do i = 1, ia
2346  do k = ks, ke
2347  rhi(k,i,j) = dens_av(k,i,j) * qv(k,i,j) &
2348  / work(k,i,j) * rvap * temp(k,i,j) &
2349  * 100.0_rp
2350  enddo
2351  enddo
2352  enddo
2353  !$acc end kernels
2354  else
2355 !OCL XFILL
2356  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
2357  !$omp private(i,j,k) &
2358  !$omp shared(RHI) &
2359  !$omp shared(KS,KE,IA,JA)
2360  !$acc kernels
2361  do j = 1, ja
2362  do i = 1, ia
2363  do k = ks, ke
2364  rhi(k,i,j) = 0.0_rp
2365  enddo
2366  enddo
2367  enddo
2368  !$acc end kernels
2369  end if
2370  dv_calculated(i_rhi) = .true.
2371  end if
2372  !$omp parallel workshare
2373  !$acc kernels
2374  var(ks:ke,:,:) = rhi(ks:ke,:,:)
2375  !$acc end kernels
2376  !$omp end parallel workshare
2377 
2378  case ( 'VOR' )
2379  if ( .not. dv_calculated(i_vor) ) then
2380  call allocate_3d( vor )
2381  !!! to move to grid !!!
2382  ! at x, v, layer
2383  ! at u, y, layer
2384  if ( prc_twod ) then
2385  !$omp parallel do private(j,k) OMP_SCHEDULE_
2386  !$acc kernels
2387 !OCL XFILL
2388  do j = 1, ja-1
2389  do k = ks, ke
2390  uh(k,is,j) = ( momx_av(k,is,j) + momx_av(k,is,j+1) ) &
2391  / ( dens_av(k,is,j) + dens_av(k,is,j+1) )
2392  enddo
2393  enddo
2394  !$acc end kernels
2395  !$omp parallel do private(j,k) OMP_SCHEDULE_
2396  !$acc kernels
2397 !OCL XFILL
2398  do j = 2, ja-1
2399  do k = ks, ke
2400  vor(k,is,j) = - ( uh(k,is,j) - uh(k,is,j-1) ) * rcdy(j)
2401  enddo
2402  enddo
2403  !$acc end kernels
2404  !$omp parallel do private(k) OMP_SCHEDULE_
2405  !$acc kernels
2406  do k = ks, ke
2407  vor(k,is,1 ) = vor(k,is,2 )
2408  vor(k,is,ja) = vor(k,is,ja-1)
2409  enddo
2410  !$acc end kernels
2411  else
2412  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2413  !$acc kernels
2414 !OCL XFILL
2415  do j = 1, ja-1
2416  do i = 2, ia
2417  do k = ks, ke
2418  uh(k,i,j) = 0.5_rp * ( momx_av(k,i ,j) + momx_av(k,i ,j+1) &
2419  + momx_av(k,i-1,j) + momx_av(k,i-1,j+1) ) &
2420  / ( dens_av(k,i,j) + dens_av(k,i,j+1) )
2421  enddo
2422  enddo
2423  enddo
2424  !$acc end kernels
2425  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2426  !$acc kernels
2427 !OCL XFILL
2428  do j = 2, ja
2429  do i = 1, ia-1
2430  do k = ks, ke
2431  vh(k,i,j) = 0.5_rp * ( momy_av(k,i,j ) + momy_av(k,i+1,j ) &
2432  + momy_av(k,i,j-1) + momy_av(k,i+1,j-1) ) &
2433  / ( dens_av(k,i,j) + dens_av(k,i+1,j) )
2434  enddo
2435  enddo
2436  enddo
2437  !$acc end kernels
2438  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2439  !$acc kernels
2440 !OCL XFILL
2441  do j = 2, ja-1
2442  do i = 2, ia-1
2443  do k = ks, ke
2444  vor(k,i,j) = ( vh(k,i,j ) - vh(k,i-1,j ) ) * rcdx(i) &
2445  - ( uh(k,i ,j) - uh(k,i ,j-1) ) * rcdy(j)
2446  enddo
2447  enddo
2448  enddo
2449  !$acc end kernels
2450  !$omp parallel do private(j,k) OMP_SCHEDULE_
2451  !$acc kernels
2452  do j = 1, ja
2453  do k = ks, ke
2454  vor(k,1 ,j) = vor(k,2 ,j)
2455  vor(k,ia,j) = vor(k,ia-1,j)
2456  enddo
2457  enddo
2458  !$acc end kernels
2459  !$omp parallel do private(i,k) OMP_SCHEDULE_
2460  !$acc kernels
2461  do i = 1, ia
2462  do k = ks, ke
2463  vor(k,i,1 ) = vor(k,i,2 )
2464  vor(k,i,ja) = vor(k,i,ja-1)
2465  enddo
2466  enddo
2467  !$acc end kernels
2468  end if
2469  dv_calculated(i_vor) = .true.
2470  end if
2471  call comm_vars8( vor(:,:,:), 1 )
2472  call comm_wait ( vor(:,:,:), 1, .false. )
2473  !$omp parallel workshare
2474  !$acc kernels
2475  var(ks:ke,:,:) = vor(ks:ke,:,:)
2476  !$acc end kernels
2477  !$omp end parallel workshare
2478 
2479  case ( 'DIV' )
2480  if ( .not. dv_calculated(i_div) ) then
2481  call allocate_3d( div )
2482  call atmos_vars_get_diagnostic( 'HDIV', work3d(:,:,:) )
2483  !!!! to move to grid !!!!
2484  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2485  !$acc kernels
2486 !OCL XFILL
2487  !$acc loop collapse(3) independent
2488  do j = 1, ja
2489  do i = 1, ia
2490  do k = ks, ke
2491  div(k,i,j) = ( momz_av(k,i,j) - momz_av(k-1,i ,j ) ) / ( real_fz(k,i,j)-real_fz(k-1,i,j) ) &
2492  + hdiv(k,i,j)
2493  enddo
2494  enddo
2495  enddo
2496  !$acc end kernels
2497  dv_calculated(i_div) = .true.
2498  end if
2499  !$omp parallel workshare
2500  !$acc kernels
2501  var(ks:ke,:,:) = div(ks:ke,:,:)
2502  !$acc end kernels
2503  !$omp end parallel workshare
2504 
2505  case ( 'HDIV' )
2506  if ( .not. dv_calculated(i_hdiv) ) then
2507  call allocate_3d( hdiv )
2508  !!!! to move to grid !!!!
2509  if ( prc_twod ) then
2510  !$omp parallel do private(j,k) OMP_SCHEDULE_
2511  !$acc kernels
2512 !OCL XFILL
2513  !$acc loop collapse(2) independent
2514  do j = 2, ja
2515  do k = ks, ke
2516  hdiv(k,is,j) = ( momy_av(k,is,j) - momy_av(k ,is,j-1) ) * rcdy(j)
2517  enddo
2518  enddo
2519  !$acc end kernels
2520  !$omp parallel do private(k) OMP_SCHEDULE_
2521  !$acc kernels
2522  do k = ks, ke
2523  hdiv(k,is,1) = hdiv(k,is,2)
2524  enddo
2525  !$acc end kernels
2526  else
2527  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2528  !$acc kernels
2529 !OCL XFILL
2530  !$acc loop collapse(3) independent
2531  do j = 2, ja
2532  do i = 2, ia
2533  do k = ks, ke
2534  hdiv(k,i,j) = ( momx_av(k,i,j) - momx_av(k ,i-1,j ) ) * rcdx(i) &
2535  + ( momy_av(k,i,j) - momy_av(k ,i ,j-1) ) * rcdy(j)
2536  enddo
2537  enddo
2538  enddo
2539  !$acc end kernels
2540  !$omp parallel do private(i,k) OMP_SCHEDULE_
2541  !$acc kernels
2542  do i = 1, ia
2543  do k = ks, ke
2544  hdiv(k,i,1) = hdiv(k,i,2)
2545  enddo
2546  enddo
2547  !$acc end kernels
2548  !$omp parallel do private(j,k) OMP_SCHEDULE_
2549  !$acc kernels
2550  do j = 1, ja
2551  do k = ks, ke
2552  hdiv(k,1,j) = hdiv(k,2,j)
2553  enddo
2554  enddo
2555  !$acc end kernels
2556  end if
2557  call comm_vars8( hdiv(:,:,:), 1 )
2558  call comm_wait ( hdiv(:,:,:), 1, .false. )
2559  dv_calculated(i_hdiv) = .true.
2560  end if
2561  !$omp parallel workshare
2562  !$acc kernels
2563  var(ks:ke,:,:) = hdiv(ks:ke,:,:)
2564  !$acc end kernels
2565  !$omp end parallel workshare
2566 
2567  case ( 'Uabs' )
2568  if ( .not. dv_calculated(i_uabs) ) then
2569  call allocate_3d( uabs )
2570 !OCL XFILL
2571  !$omp parallel do private(k,i,j) OMP_SCHEDULE_ collapse(2)
2572  !$acc kernels
2573  !$acc loop collapse(3) independent
2574  do j = 1, ja
2575  do i = 1, ia
2576  do k = ks, ke
2577  uabs(k,i,j) = sqrt( w(k,i,j)**2 + u(k,i,j)**2 + v(k,i,j)**2 )
2578  enddo
2579  enddo
2580  enddo
2581  !$acc end kernels
2582  dv_calculated(i_uabs) = .true.
2583  end if
2584  !$omp parallel workshare
2585  !$acc kernels
2586  var(ks:ke,:,:) = uabs(ks:ke,:,:)
2587  !$acc end kernels
2588  !$omp end parallel workshare
2589 
2590  case ( 'N2' )
2591  if ( .not. dv_calculated(i_n2) ) then
2592  call allocate_3d( n2 )
2593  call atmos_diagnostic_get_n2( &
2594  ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2595  pott(:,:,:), rtot(:,:,:), & !(in)
2596  real_cz(:,:,:), real_fz(:,:,:), & !(in)
2597  f2h(:,:,:,:), & !(in)
2598  n2(:,:,:) ) ! (out)
2599  dv_calculated(i_n2) = .true.
2600  end if
2601  !$omp parallel workshare
2602  !$acc kernels
2603  var(ks:ke,:,:) = n2(ks:ke,:,:)
2604  !$acc end kernels
2605  !$omp end parallel workshare
2606 
2607  case ( 'MSE' )
2608  if ( .not. dv_calculated(i_mse) ) then
2609  call allocate_3d( mse )
2610  call atmos_vars_get_diagnostic( 'LHV', work3d(:,:,:) )
2611 !OCL XFILL
2612  !$omp parallel do private(k,i,j) OMP_SCHEDULE_ collapse(2)
2613  !$acc kernels
2614  !$acc loop collapse(3) independent
2615  do j = 1, ja
2616  do i = 1, ia
2617  do k = ks, ke
2618  mse(k,i,j) = cptot(k,i,j) * temp(k,i,j) &
2619  + grav * ( real_cz(k,i,j) - real_fz(ks-1,i,j) ) &
2620  + lhv(k,i,j) * qv(k,i,j)
2621  enddo
2622  enddo
2623  enddo
2624  !$acc end kernels
2625  dv_calculated(i_mse) = .true.
2626  end if
2627  !$omp parallel workshare
2628  !$acc kernels
2629  var(ks:ke,:,:) = mse(ks:ke,:,:)
2630  !$acc end kernels
2631  !$omp end parallel workshare
2632 
2633  case ( 'TDEW' )
2634  if ( .not. dv_calculated(i_tdew) ) then
2635  call allocate_3d( tdew )
2636  call atmos_saturation_tdew_liq( ka, ks, ke, ia, 1, ia, ja, 1, ja, &
2637  dens(:,:,:), temp(:,:,:), qv(:,:,:), & ! [IN]
2638  tdew(:,:,:) ) ! [OUT]
2639  dv_calculated(i_tdew) = .true.
2640  end if
2641  !$omp parallel workshare
2642  !$acc kernels
2643  var(ks:ke,:,:) = tdew(ks:ke,:,:)
2644  !$acc end kernels
2645  !$omp end parallel workshare
2646 
2647  case ( 'ENGP' )
2648  if ( .not. dv_calculated(i_engp) ) then
2649  call allocate_3d( engp )
2650  !$omp parallel do private(k,i,j) OMP_SCHEDULE_ collapse(2)
2651  !$acc kernels
2652  !$acc loop collapse(3) independent
2653  do j = 1, ja
2654  do i = 1, ia
2655  do k = ks, ke
2656  engp(k,i,j) = dens_av(k,i,j) * grav * real_cz(k,i,j)
2657  end do
2658  end do
2659  end do
2660  !$acc end kernels
2661  dv_calculated(i_engp) = .true.
2662  end if
2663  !$omp parallel workshare
2664  !$acc kernels
2665  var(ks:ke,:,:) = engp(ks:ke,:,:)
2666  !$acc end kernels
2667  !$omp end parallel workshare
2668 
2669  case ( 'ENGK' )
2670  if ( .not. dv_calculated(i_engk) ) then
2671  call allocate_3d( engk )
2672  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2673  !$acc kernels
2674  !$acc loop collapse(3) independent
2675  do j = 1, ja
2676  do i = 1, ia
2677  do k = ks, ke
2678  engk(k,i,j) = 0.5_rp * dens_av(k,i,j) &
2679  * ( w(k,i,j)**2 + u(k,i,j)**2 + v(k,i,j)**2 )
2680  end do
2681  end do
2682  end do
2683  !$acc end kernels
2684  dv_calculated(i_engk) = .true.
2685  end if
2686  !$omp parallel workshare
2687  !$acc kernels
2688  var(ks:ke,:,:) = engk(ks:ke,:,:)
2689  !$acc end kernels
2690  !$omp end parallel workshare
2691 
2692  case ( 'ENGI' )
2693  if ( .not. dv_calculated(i_engi) ) then
2694  call allocate_3d( engi )
2695  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2696  !$acc kernels
2697  !$acc loop collapse(3) independent
2698  do j = 1, ja
2699  do i = 1, ia
2700  do k = ks, ke
2701  engi(k,i,j) = dens_av(k,i,j) * qdry(k,i,j) * temp(k,i,j) * cvdry
2702  !$acc loop seq
2703  do iq = 1, qa
2704  engi(k,i,j) = engi(k,i,j) &
2705  + dens_av(k,i,j) * qtrc_av(k,i,j,iq) * ( temp(k,i,j) * tracer_cv(iq) + tracer_engi0(iq) )
2706  enddo
2707  end do
2708  end do
2709  end do
2710  !$acc end kernels
2711  dv_calculated(i_engi) = .true.
2712  end if
2713  !$omp parallel workshare
2714  !$acc kernels
2715  var(ks:ke,:,:) = engi(ks:ke,:,:)
2716  !$acc end kernels
2717  !$omp end parallel workshare
2718 
2719  case ( 'ENGT' )
2720  if ( .not. dv_calculated(i_engt) ) then
2721  call allocate_3d( engt )
2722  call atmos_vars_get_diagnostic( 'ENGP', work3d(:,:,:) )
2723  call atmos_vars_get_diagnostic( 'ENGK', work3d(:,:,:) )
2724  call atmos_vars_get_diagnostic( 'ENGI', work3d(:,:,:) )
2725  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2726  !$acc kernels
2727  !$acc loop collapse(3) independent
2728  do j = 1, ja
2729  do i = 1, ia
2730  do k = ks, ke
2731  engt(k,i,j) = engp(k,i,j) + engk(k,i,j) + engi(k,i,j)
2732  enddo
2733  enddo
2734  enddo
2735  !$acc end kernels
2736  dv_calculated(i_engt) = .true.
2737  end if
2738  !$omp parallel workshare
2739  !$acc kernels
2740  var(ks:ke,:,:) = engt(ks:ke,:,:)
2741  !$acc end kernels
2742  !$omp end parallel workshare
2743 
2744  case ( 'DENS_PRIM' )
2745  if ( .not. dv_calculated(i_dens_prim) ) then
2746  call allocate_3d( dens_prim )
2747  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
2748 !OCL XFILL
2749  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2750  !$acc kernels
2751  !$acc loop collapse(3) independent
2752  do j = 1, ja
2753  do i = 1, ia
2754  do k = ks, ke
2755  dens_prim(k,i,j) = dens_av(k,i,j) - dens_mean(k)
2756  enddo
2757  enddo
2758  enddo
2759  !$acc end kernels
2760  dv_calculated(i_dens_prim) = .true.
2761  end if
2762  !$omp parallel workshare
2763  !$acc kernels
2764  var(ks:ke,:,:) = dens_prim(ks:ke,:,:)
2765  !$acc end kernels
2766  !$omp end parallel workshare
2767 
2768  case ( 'W_PRIM' )
2769  if ( .not. dv_calculated(i_w_prim) ) then
2770  call allocate_3d( w_prim )
2771  call atmos_vars_get_diagnostic( 'W_MEAN', work1d(:) )
2772 !OCL XFILL
2773  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2774  !$acc kernels
2775  !$acc loop collapse(3) independent
2776  do j = 1, ja
2777  do i = 1, ia
2778  do k = ks, ke
2779  w_prim(k,i,j) = w(k,i,j) - w_mean(k)
2780  enddo
2781  enddo
2782  enddo
2783  !$acc end kernels
2784  dv_calculated(i_w_prim) = .true.
2785  end if
2786  !$omp parallel workshare
2787  !$acc kernels
2788  var(ks:ke,:,:) = w_prim(ks:ke,:,:)
2789  !$acc end kernels
2790  !$omp end parallel workshare
2791 
2792  case ( 'U_PRIM' )
2793  if ( .not. dv_calculated(i_u_prim) ) then
2794  call allocate_3d( u_prim )
2795  call atmos_vars_get_diagnostic( 'U_MEAN', work1d(:) )
2796 !OCL XFILL
2797  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2798  !$acc kernels
2799  !$acc loop collapse(3) independent
2800  do j = 1, ja
2801  do i = 1, ia
2802  do k = ks, ke
2803  u_prim(k,i,j) = u(k,i,j) - u_mean(k)
2804  enddo
2805  enddo
2806  enddo
2807  !$acc end kernels
2808  dv_calculated(i_u_prim) = .true.
2809  end if
2810  !$omp parallel workshare
2811  !$acc kernels
2812  var(ks:ke,:,:) = u_prim(ks:ke,:,:)
2813  !$acc end kernels
2814  !$omp end parallel workshare
2815 
2816  case ( 'V_PRIM' )
2817  if ( .not. dv_calculated(i_v_prim) ) then
2818  call allocate_3d( v_prim )
2819  call atmos_vars_get_diagnostic( 'V_MEAN', work1d(:) )
2820 !OCL XFILL
2821  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2822  !$acc kernels
2823  !$acc loop collapse(3) independent
2824  do j = 1, ja
2825  do i = 1, ia
2826  do k = ks, ke
2827  v_prim(k,i,j) = v(k,i,j) - v_mean(k)
2828  enddo
2829  enddo
2830  enddo
2831  !$acc end kernels
2832  dv_calculated(i_v_prim) = .true.
2833  end if
2834  !$omp parallel workshare
2835  !$acc kernels
2836  var(ks:ke,:,:) = v_prim(ks:ke,:,:)
2837  !$acc end kernels
2838  !$omp end parallel workshare
2839 
2840  case ( 'PT_PRIM' )
2841  if ( .not. dv_calculated(i_pt_prim) ) then
2842  call allocate_3d( pt_prim )
2843  call atmos_vars_get_diagnostic( 'PT_MEAN', work1d(:) )
2844 !OCL XFILL
2845  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2846  !$acc kernels
2847  !$acc loop collapse(3) independent
2848  do j = 1, ja
2849  do i = 1, ia
2850  do k = ks, ke
2851  pt_prim(k,i,j) = pott(k,i,j) - pt_mean(k)
2852  enddo
2853  enddo
2854  enddo
2855  !$acc end kernels
2856  dv_calculated(i_pt_prim) = .true.
2857  end if
2858  !$omp parallel workshare
2859  !$acc kernels
2860  var(ks:ke,:,:) = pt_prim(ks:ke,:,:)
2861  !$acc end kernels
2862  !$omp end parallel workshare
2863 
2864  case ( 'W_PRIM2' )
2865  if ( .not. dv_calculated(i_w_prim2) ) then
2866  call allocate_3d( w_prim2 )
2867  call atmos_vars_get_diagnostic( 'W_PRIM', work3d(:,:,:) )
2868 !OCL XFILL
2869  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2870  !$acc kernels
2871  !$acc loop collapse(3) independent
2872  do j = 1, ja
2873  do i = 1, ia
2874  do k = ks, ke
2875  w_prim2(k,i,j) = w_prim(k,i,j)**2
2876  enddo
2877  enddo
2878  enddo
2879  !$acc end kernels
2880  dv_calculated(i_w_prim2) = .true.
2881  end if
2882  !$omp parallel workshare
2883  !$acc kernels
2884  var(ks:ke,:,:) = w_prim2(ks:ke,:,:)
2885  !$acc end kernels
2886  !$omp end parallel workshare
2887 
2888  case ( 'PT_W_PRIM' )
2889  if ( .not. dv_calculated(i_pt_w_prim) ) then
2890  call allocate_3d( pt_w_prim )
2891  call atmos_vars_get_diagnostic( 'W_PRIM', work3d(:,:,:) )
2892  call atmos_vars_get_diagnostic( 'PT_PRIM', work3d(:,:,:) )
2893  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2894  !$acc kernels
2895  !$acc loop collapse(3) independent
2896  do j = 1, ja
2897  do i = 1, ia
2898  do k = ks, ke
2899  pt_w_prim(k,i,j) = w_prim(k,i,j) * pt_prim(k,i,j) * dens_av(k,i,j) * cpdry
2900  enddo
2901  enddo
2902  enddo
2903  !$acc end kernels
2904  dv_calculated(i_pt_w_prim) = .true.
2905  end if
2906  !$omp parallel workshare
2907  !$acc kernels
2908  var(ks:ke,:,:) = pt_w_prim(ks:ke,:,:)
2909  !$acc end kernels
2910  !$omp end parallel workshare
2911 
2912  case ( 'W_PRIM3' )
2913  if ( .not. dv_calculated(i_w_prim3) ) then
2914  call allocate_3d( w_prim3 )
2915  call atmos_vars_get_diagnostic( 'W_PRIM', work3d(:,:,:) )
2916 !OCL XFILL
2917  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2918  !$acc kernels
2919  !$acc loop collapse(3) independent
2920  do j = 1, ja
2921  do i = 1, ia
2922  do k = ks, ke
2923  w_prim3(k,i,j) = w_prim(k,i,j)**3
2924  enddo
2925  enddo
2926  enddo
2927  !$acc end kernels
2928  dv_calculated(i_w_prim3) = .true.
2929  end if
2930  !$omp parallel workshare
2931  !$acc kernels
2932  var(ks:ke,:,:) = w_prim3(ks:ke,:,:)
2933  !$acc end kernels
2934  !$omp end parallel workshare
2935 
2936  case ( 'TKE_RS' )
2937  if ( .not. dv_calculated(i_tke_rs) ) then
2938  call allocate_3d( tke_rs )
2939  call atmos_vars_get_diagnostic( 'W_PRIM', work3d(:,:,:) )
2940  call atmos_vars_get_diagnostic( 'U_PRIM', work3d(:,:,:) )
2941  call atmos_vars_get_diagnostic( 'V_PRIM', work3d(:,:,:) )
2942 !OCL XFILL
2943  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2944  !$acc kernels
2945  !$acc loop collapse(3) independent
2946  do j = 1, ja
2947  do i = 1, ia
2948  do k = ks, ke
2949  tke_rs(k,i,j) = 0.5_rp * ( w_prim(k,i,j)**2 + u_prim(k,i,j)**2 + v_prim(k,i,j)**2 )
2950  enddo
2951  enddo
2952  enddo
2953  !$acc end kernels
2954  dv_calculated(i_tke_rs) = .true.
2955  end if
2956  !$omp parallel workshare
2957  !$acc kernels
2958  var(ks:ke,:,:) = tke_rs(ks:ke,:,:)
2959  !$acc end kernels
2960  !$omp end parallel workshare
2961 
2962  case ( 'VELZ' )
2963  if ( .not. dv_calculated(i_velz) ) then
2964  call allocate_3d( velz )
2965 !OCL XFILL
2966  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2967  !$acc kernels
2968  !$acc loop collapse(2) independent
2969  do j = 1, ja
2970  do i = 1, ia
2971  velz(ks-1,i,j) = 0.0_rp
2972  !$acc loop independent
2973  do k = ks, ke-1
2974  velz(k,i,j) = momz(k,i,j) * 2.0_rp / ( dens(k,i,j) + dens(k+1,i,j) )
2975  end do
2976  velz(ke,i,j) = 0.0_rp
2977  enddo
2978  enddo
2979  !$acc end kernels
2980  dv_calculated(i_velz) = .true.
2981  end if
2982  !$omp parallel workshare
2983  !$acc kernels
2984  var(ks-1:ke,:,:) = velz(ks-1:ke,:,:)
2985  !$acc end kernels
2986  !$omp end parallel workshare
2987 
2988  case ( 'VELX' )
2989  if ( .not. dv_calculated(i_velx) ) then
2990  call allocate_3d( velx )
2991  if ( prc_twod ) then
2992 !OCL XFILL
2993  !$omp parallel do private(j,k) OMP_SCHEDULE_
2994  !$acc kernels
2995  !$acc loop collapse(2) independent
2996  do j = 1, ja
2997  do k = ks, ke
2998  velx(k,is,j) = momx(k,is,j) / dens(k,is,j)
2999  enddo
3000  enddo
3001  !$acc end kernels
3002  else
3003  !OCL XFILL
3004  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3005  !$acc kernels
3006  !$acc loop collapse(3) independent
3007  do j = 1, ja
3008  do i = 1, ia-1
3009  do k = ks, ke
3010  velx(k,i,j) = momx(k,i,j) * 2.0_rp / ( dens(k,i,j) + dens(k,i+1,j) )
3011  enddo
3012  enddo
3013  enddo
3014  !$acc end kernels
3015 !OCL XFILL
3016  !$omp parallel do private(j,k) OMP_SCHEDULE_
3017  !$acc kernels
3018  !$acc loop collapse(2) independent
3019  do j = 1, ja
3020  do k = ks, ke
3021  velx(k,ia,j) = momx(k,ia,j) / dens(k,ia,j)
3022  enddo
3023  enddo
3024  !$acc end kernels
3025  call comm_vars8( velx(:,:,:), 1 )
3026  call comm_wait ( velx(:,:,:), 1, .false. )
3027  end if
3028  dv_calculated(i_velx) = .true.
3029  end if
3030  !$omp parallel workshare
3031  !$acc kernels
3032  var(ks:ke,:,:) = velx(ks:ke,:,:)
3033  !$acc end kernels
3034  !$omp end parallel workshare
3035 
3036  case ( 'VELY' )
3037  if ( .not. dv_calculated(i_vely) ) then
3038  call allocate_3d( vely )
3039 !OCL XFILL
3040  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3041  !$acc kernels
3042  !$acc loop collapse(3) independent
3043  do j = 1, ja-1
3044  do i = 1, ia
3045  do k = ks, ke
3046  vely(k,i,j) = momy(k,i,j) * 2.0_rp / ( dens(k,i,j) + dens(k,i,j+1) )
3047  enddo
3048  enddo
3049  enddo
3050  !$acc end kernels
3051 !OCL XFILL
3052  !$omp parallel do private(i,k) OMP_SCHEDULE_
3053  !$acc kernels
3054  !$acc loop collapse(2) independent
3055  do i = 1, ia
3056  do k = ks, ke
3057  vely(k,i,ja) = momy(k,i,ja) / dens(k,i,ja)
3058  enddo
3059  enddo
3060  !$acc end kernels
3061  call comm_vars8( vely(:,:,:), 1 )
3062  call comm_wait ( vely(:,:,:), 1, .false. )
3063  dv_calculated(i_vely) = .true.
3064  end if
3065  !$omp parallel workshare
3066  !$acc kernels
3067  var(ks:ke,:,:) = vely(ks:ke,:,:)
3068  !$acc end kernels
3069  !$omp end parallel workshare
3070 
3071  case ( 'Umet' )
3072  if ( .not. dv_calculated(i_umet) ) then
3073  call allocate_3d( umet )
3074 !OCL XFILL
3075  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3076  !$acc kernels
3077  !$acc loop collapse(3) independent
3078  do j = 1, ja
3079  do i = 1, ia
3080  do k = ks, ke
3081  umet(k,i,j) = u(k,i,j) * rotc(i,j,1) - v(k,i,j) * rotc(i,j,2)
3082  end do
3083  end do
3084  end do
3085  !$acc end kernels
3086  dv_calculated(i_umet) = .true.
3087  end if
3088  !$omp parallel workshare
3089  !$acc kernels
3090  var(ks:ke,:,:) = umet(ks:ke,:,:)
3091  !$acc end kernels
3092  !$omp end parallel workshare
3093 
3094  case ( 'Vmet' )
3095  if ( .not. dv_calculated(i_vmet) ) then
3096  call allocate_3d( vmet )
3097 !OCL XFILL
3098  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3099  !$acc kernels
3100  !$acc loop collapse(3) independent
3101  do j = 1, ja
3102  do i = 1, ia
3103  do k = ks, ke
3104  vmet(k,i,j) = u(k,i,j) * rotc(i,j,2) + v(k,i,j) * rotc(i,j,1)
3105  end do
3106  end do
3107  end do
3108  !$acc end kernels
3109  dv_calculated(i_vmet) = .true.
3110  end if
3111  !$omp parallel workshare
3112  !$acc kernels
3113  var(ks:ke,:,:) = vmet(ks:ke,:,:)
3114  !$acc end kernels
3115  !$omp end parallel workshare
3116 
3117  case default
3118  log_error("ATMOS_vars_calc_diagnostics",*) 'name is invalid for ATMOS_vars_get_diagnostic_3D: ', trim(vname)
3119  call prc_abort
3120  end select
3121 
3122  !$acc end data
3123 
3124  return
3125  end subroutine atmos_vars_get_diagnostic_3d
3126 
3127  !-----------------------------------------------------------------------------
3129  recursive subroutine atmos_vars_get_diagnostic_2d( &
3130  vname, &
3131  var )
3132  use scale_prc, only: &
3133  prc_abort
3134  use scale_atmos_grid_cartesc_real, only: &
3135  real_cz => atmos_grid_cartesc_real_cz, &
3136  real_fz => atmos_grid_cartesc_real_fz
3137  use scale_atmos_adiabat, only: &
3138  atmos_adiabat_cape
3139  use mod_atmos_phy_mp_vars, only: &
3140  sflx_rain_mp => atmos_phy_mp_sflx_rain, &
3141  sflx_snow_mp => atmos_phy_mp_sflx_snow
3142  use mod_atmos_phy_cp_vars, only: &
3143  sflx_rain_cp => atmos_phy_cp_sflx_rain, &
3144  sflx_snow_cp => atmos_phy_cp_sflx_snow
3145  implicit none
3146 
3147  character(len=*), intent(in) :: vname
3148  real(rp), intent(out) :: var(:,:)
3149 
3150  real(rp) :: fact, sum
3151  integer :: k, i, j
3152  !---------------------------------------------------------------------------
3153 
3154  !$acc data copyout(var)
3155 
3156  select case ( vname )
3157  case ( 'LWP' )
3158  if ( .not. dv_calculated(i_lwp) ) then
3159  call allocate_2d( lwp )
3160  call atmos_vars_get_diagnostic( 'QLIQ', work3d(:,:,:) )
3161  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
3162  !$omp private(i,j,k,sum) &
3163  !$omp shared(LWP,QLIQ,DENS_av,REAL_FZ) &
3164  !$omp shared(KS,KE,IA,JA)
3165  !$acc kernels
3166  !$acc loop collapse(2) independent
3167  do j = 1, ja
3168  do i = 1, ia
3169  sum = 0.0_rp
3170  do k = ks, ke
3171  sum = sum &
3172  + qliq(k,i,j) * dens_av(k,i,j) * ( real_fz(k,i,j)-real_fz(k-1,i,j) ) * 1.e3_rp ! [kg/m2->g/m2]
3173  enddo
3174  lwp(i,j) = sum
3175  enddo
3176  enddo
3177  !$acc end kernels
3178  dv_calculated(i_lwp) = .true.
3179  end if
3180  !$omp parallel workshare
3181  !$acc kernels
3182  var(:,:) = lwp(:,:)
3183  !$acc end kernels
3184  !$omp end parallel workshare
3185 
3186  case ( 'IWP' )
3187  if ( .not. dv_calculated(i_iwp) ) then
3188  call allocate_2d( iwp )
3189  call atmos_vars_get_diagnostic( 'QICE', work3d(:,:,:) )
3190  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
3191  !$omp private(i,j,k,sum) &
3192  !$omp shared(IWP,QICE,DENS_av,REAL_FZ) &
3193  !$omp shared(KS,KE,IA,JA)
3194  !$acc kernels
3195  !$acc loop collapse(2) independent
3196  do j = 1, ja
3197  do i = 1, ia
3198  sum = 0.0_rp
3199  !$acc loop reduction(+:sum)
3200  do k = ks, ke
3201  sum = sum &
3202  + qice(k,i,j) * dens_av(k,i,j) * ( real_fz(k,i,j)-real_fz(k-1,i,j) ) * 1.e3_rp ! [kg/m2->g/m2]
3203  enddo
3204  iwp(i,j) = sum
3205  enddo
3206  enddo
3207  !$acc end kernels
3208  dv_calculated(i_iwp) = .true.
3209  end if
3210  !$omp parallel workshare
3211  !$acc kernels
3212  var(:,:) = iwp(:,:)
3213  !$acc end kernels
3214  !$omp end parallel workshare
3215 
3216  case ( 'PW' )
3217  if ( .not. dv_calculated(i_pw) ) then
3218  call allocate_2d( pw )
3219  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
3220  !$omp private(i,j,k,sum) &
3221  !$omp shared(PW,QV,DENS_av,REAL_FZ) &
3222  !$omp shared(KS,KE,IA,JA)
3223  !$acc kernels
3224  !$acc loop collapse(2) independent
3225  do j = 1, ja
3226  do i = 1, ia
3227  sum = 0.0_rp
3228  !$acc loop reduction(+:sum)
3229  do k = ks, ke
3230  sum = sum &
3231  + qv(k,i,j) * dens_av(k,i,j) * ( real_fz(k,i,j)-real_fz(k-1,i,j) ) * 1.e3_rp ! [kg/m2->g/m2]
3232  enddo
3233  pw(i,j) = sum
3234  enddo
3235  enddo
3236  !$acc end kernels
3237  dv_calculated(i_pw) = .true.
3238  end if
3239  !$omp parallel workshare
3240  !$acc kernels
3241  var(:,:) = pw(:,:)
3242  !$acc end kernels
3243  !$omp end parallel workshare
3244 
3245  case ( 'PBLH' )
3246  if ( .not. dv_calculated(i_pblh) ) then
3247  call allocate_2d( pblh )
3248  call atmos_vars_get_diagnostic( 'POTV', work3d(:,:,:) )
3249  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
3250  !$omp private(k,i,j) &
3251  !$omp private(fact) &
3252  !$omp shared(PBLH,POTV,REAL_CZ,REAL_FZ) &
3253  !$omp shared(KS,KE,IA,JA)
3254  !$acc kernels
3255  !$acc loop independent
3256  do j = 1, ja
3257  !$acc loop independent
3258  do i = 1, ia
3259  pblh(i,j) = real_cz(ks,i,j) - real_fz(ks-1,i,j)
3260  do k = ks+1, ke
3261  if ( potv(k,i,j) > potv(ks,i,j) ) then
3262  fact = ( potv(ks,i,j) - potv(k-1,i,j) ) &
3263  / ( potv(k,i,j) - potv(k-1,i,j) )
3264  pblh(i,j) = real_cz(k-1,i,j) - real_fz(ks-1,i,j) &
3265  + fact * ( real_cz(k,i,j) - real_cz(k-1,i,j) )
3266 
3267  exit
3268  endif
3269  enddo
3270  enddo
3271  enddo
3272  !$acc end kernels
3273  dv_calculated(i_pblh) = .true.
3274  end if
3275  !$omp parallel workshare
3276  !$acc kernels
3277  var(:,:) = pblh(:,:)
3278  !$acc end kernels
3279  !$omp end parallel workshare
3280 
3281  case ( 'CAPE', 'CIN', 'LCL', 'LFC', 'LNB' )
3282  if ( .not. dv_calculated(i_cape) ) then
3283  call allocate_2d( cape )
3284  call allocate_2d( cin )
3285  call allocate_2d( lcl )
3286  call allocate_2d( lfc )
3287  call allocate_2d( lnb )
3288  call atmos_adiabat_cape( &
3289  ka, ks, ke, ia, is, ie, ja, js, je, &
3290  ks, & ! (in)
3291  dens_av(:,:,:), temp(:,:,:), pres(:,:,:), & ! (in)
3292  qv(:,:,:), qc(:,:,:), qdry(:,:,:), & ! (in)
3293  rtot(:,:,:), cptot(:,:,:), & ! (in)
3294  real_cz(:,:,:), real_fz(:,:,:), & ! (in)
3295  cape(:,:), cin(:,:), lcl(:,:), lfc(:,:), lnb(:,:) ) ! (out)
3296  dv_calculated(i_cape) = .true.
3297  end if
3298  select case ( vname )
3299  case ( 'CAPE' )
3300  !$omp parallel do private(i,j) OMP_SCHEDULE_
3301  !$acc kernels
3302  do j = js, je
3303  do i = is, ie
3304  var(i,j) = cape(i,j)
3305  end do
3306  end do
3307  !$acc end kernels
3308  case ( 'CIN' )
3309  !$omp parallel do private(i,j) OMP_SCHEDULE_
3310  !$acc kernels
3311  do j = js, je
3312  do i = is, ie
3313  var(i,j) = cin(i,j)
3314  end do
3315  end do
3316  !$acc end kernels
3317  case ( 'LCL' )
3318  !$omp parallel do private(i,j) OMP_SCHEDULE_
3319  !$acc kernels
3320  do j = js, je
3321  do i = is, ie
3322  var(i,j) = lcl(i,j)
3323  end do
3324  end do
3325  !$acc end kernels
3326  case ( 'LFC' )
3327  !$omp parallel do private(i,j) OMP_SCHEDULE_
3328  !$acc kernels
3329  do j = js, je
3330  do i = is, ie
3331  var(i,j) = lfc(i,j)
3332  end do
3333  end do
3334  !$acc end kernels
3335  case ( 'LNB' )
3336  !$omp parallel do private(i,j) OMP_SCHEDULE_
3337  !$acc kernels
3338  do j = js, je
3339  do i = is, ie
3340  var(i,j) = lnb(i,j)
3341  end do
3342  end do
3343  !$acc end kernels
3344  end select
3345 
3346  case ( 'PREC' )
3347  !$omp parallel do private(i,j) OMP_SCHEDULE_
3348  !$acc kernels
3349  do j = js, je
3350  do i = is, ie
3351  var(i,j) = prec(i,j)
3352  end do
3353  end do
3354  !$acc end kernels
3355 
3356  case ( 'RAIN', 'SNOW' )
3357  if ( .not. dv_calculated(i_rain) ) then
3358  call allocate_2d( rain )
3359  call allocate_2d( snow )
3360  !$omp parallel do private(i,j) OMP_SCHEDULE_
3361  !$acc kernels
3362  !$acc loop collapse(2) independent
3363  do j = js, je
3364  do i = is, ie
3365  rain(i,j) = sflx_rain_mp(i,j) + sflx_rain_cp(i,j)
3366  snow(i,j) = sflx_snow_mp(i,j) + sflx_snow_cp(i,j)
3367  enddo
3368  enddo
3369  !$acc end kernels
3370  dv_calculated(i_rain) = .true.
3371  end if
3372  select case (vname)
3373  case ( 'RAIN' )
3374  !$omp parallel do private(i,j) OMP_SCHEDULE_
3375  !$acc kernels
3376  do j = js, je
3377  do i = is, ie
3378  var(i,j) = rain(i,j)
3379  end do
3380  end do
3381  !$acc end kernels
3382  case ( 'SNOW' )
3383  !$omp parallel do private(i,j) OMP_SCHEDULE_
3384  !$acc kernels
3385  do j = js, je
3386  do i = is, ie
3387  var(i,j) = snow(i,j)
3388  end do
3389  end do
3390  !$acc end kernels
3391  end select
3392 
3393  case default
3394  log_error("ATMOS_vars_calc_diagnostics",*) 'name is invalid for ATMOS_vars_get_diagnostic_2D: ', trim(vname)
3395  call prc_abort
3396  end select
3397 
3398  !$acc end data
3399 
3400  return
3401  end subroutine atmos_vars_get_diagnostic_2d
3402 
3403  !-----------------------------------------------------------------------------
3405  recursive subroutine atmos_vars_get_diagnostic_1d( &
3406  vname, &
3407  var )
3408  use scale_const, only: &
3409  cpdry => const_cpdry
3410  use scale_prc, only: &
3411  prc_abort
3412  use scale_statistics, only: &
3413  statistics_horizontal_mean
3414  use scale_atmos_grid_cartesc_real, only: &
3416  implicit none
3417 
3418  character(len=*), intent(in) :: vname
3419  real(rp), intent(out) :: var(:)
3420 
3421  real(rp) :: work(ka,ia,ja)
3422  integer :: k, i, j
3423  !---------------------------------------------------------------------------
3424 
3425  !$acc data copyout(var) create(WORK)
3426 
3427  select case ( vname )
3428  case ( 'DENS_MEAN' )
3429  if ( .not. dv_calculated(i_dens_mean) ) then
3430  call allocate_1d( dens_mean )
3431  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3432  dens(:,:,:), area(:,:), dens_mean(:) )
3433  dv_calculated(i_dens_mean) = .true.
3434  end if
3435  !$acc kernels
3436  var(:) = dens_mean(:)
3437  !$acc end kernels
3438 
3439  case ( 'W_MEAN' )
3440  if ( .not. dv_calculated(i_w_mean) ) then
3441  call allocate_1d( w_mean )
3442  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3443 !OCL XFILL
3444  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3445  !$acc kernels
3446  do j = jsb, jeb
3447  do i = isb, ieb
3448  do k = ks, ke
3449  work(k,i,j) = w(k,i,j) * dens_av(k,i,j)
3450  enddo
3451  enddo
3452  enddo
3453  !$acc end kernels
3454  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3455  work(:,:,:), area(:,:), w_mean(:) )
3456  !$acc kernels
3457  !$acc loop independent
3458  do k = ks, ke
3459  w_mean(k) = w_mean(k) / dens_mean(k)
3460  enddo
3461  !$acc end kernels
3462  dv_calculated(i_w_mean) = .true.
3463  end if
3464  !$acc kernels
3465  var(:) = w_mean(:)
3466  !$acc end kernels
3467 
3468  case ( 'U_MEAN' )
3469  if ( .not. dv_calculated(i_u_mean) ) then
3470  call allocate_1d( u_mean )
3471  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3472 !OCL XFILL
3473  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3474  !$acc kernels
3475  do j = jsb, jeb
3476  do i = isb, ieb
3477  do k = ks, ke
3478  work(k,i,j) = u(k,i,j) * dens_av(k,i,j)
3479  enddo
3480  enddo
3481  enddo
3482  !$acc end kernels
3483  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3484  work(:,:,:), area(:,:), u_mean(:) )
3485  !$acc kernels
3486  !$acc loop independent
3487  do k = ks, ke
3488  u_mean(k) = u_mean(k) / dens_mean(k)
3489  enddo
3490  !$acc end kernels
3491  dv_calculated(i_u_mean) = .true.
3492  end if
3493  !$acc kernels
3494  var(:) = u_mean(:)
3495  !$acc end kernels
3496 
3497  case ( 'V_MEAN' )
3498  if ( .not. dv_calculated(i_v_mean) ) then
3499  call allocate_1d( v_mean )
3500  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3501 !OCL XFILL
3502  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3503  !$acc kernels
3504  do j = jsb, jeb
3505  do i = isb, ieb
3506  do k = ks, ke
3507  work(k,i,j) = v(k,i,j) * dens_av(k,i,j)
3508  enddo
3509  enddo
3510  enddo
3511  !$acc end kernels
3512  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3513  work(:,:,:), area(:,:), v_mean(:) )
3514  !$acc kernels
3515  !$acc loop independent
3516  do k = ks, ke
3517  v_mean(k) = v_mean(k) / dens_mean(k)
3518  enddo
3519  !$acc end kernels
3520  dv_calculated(i_v_mean) = .true.
3521  end if
3522  !$acc kernels
3523  var(:) = v_mean(:)
3524  !$acc end kernels
3525 
3526  case ( 'PT_MEAN' )
3527  if ( .not. dv_calculated(i_pt_mean) ) then
3528  call allocate_1d( pt_mean )
3529  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3530  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3531  rhot(:,:,:), area(:,:), pt_mean(:) )
3532  !$acc kernels
3533  !$acc loop independent
3534  do k = ks, ke
3535  pt_mean(k) = pt_mean(k) / dens_mean(k)
3536  enddo
3537  !$acc end kernels
3538  dv_calculated(i_pt_mean) = .true.
3539  end if
3540  !$acc kernels
3541  var(:) = pt_mean(:)
3542  !$acc end kernels
3543 
3544  case ( 'T_MEAN' )
3545  if ( .not. dv_calculated(i_t_mean) ) then
3546  call allocate_1d( t_mean )
3547  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3548 !OCL XFILL
3549  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3550  !$acc kernels
3551  do j = jsb, jeb
3552  do i = isb, ieb
3553  do k = ks, ke
3554  work(k,i,j) = temp(k,i,j) * dens_av(k,i,j)
3555  enddo
3556  enddo
3557  enddo
3558  !$acc end kernels
3559  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3560  work(:,:,:), area(:,:), t_mean(:) )
3561  !$acc kernels
3562  !$acc loop independent
3563  do k = ks, ke
3564  t_mean(k) = t_mean(k) / dens_mean(k)
3565  enddo
3566  !$acc end kernels
3567  dv_calculated(i_t_mean) = .true.
3568  end if
3569  !$acc kernels
3570  var(:) = t_mean(:)
3571  !$acc end kernels
3572 
3573  case ( 'QV_MEAN' )
3574  if ( .not. dv_calculated(i_qv_mean) ) then
3575  call allocate_1d( qv_mean )
3576  if ( moist ) then
3577  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3578 !OCL XFILL
3579  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3580  !$acc kernels
3581  do j = jsb, jeb
3582  do i = isb, ieb
3583  do k = ks, ke
3584  work(k,i,j) = qv(k,i,j) * dens_av(k,i,j)
3585  enddo
3586  enddo
3587  enddo
3588  !$acc end kernels
3589  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3590  work(:,:,:), area(:,:), qv_mean(:) )
3591  !$acc kernels
3592  !$acc loop independent
3593  do k = ks, ke
3594  qv_mean(k) = qv_mean(k) / dens_mean(k)
3595  enddo
3596  !$acc end kernels
3597  else
3598  !$omp parallel do private(k) OMP_SCHEDULE_
3599  !$acc kernels
3600  do k = ks, ke
3601  qv_mean(k) = 0.0_rp
3602  enddo
3603  !$acc end kernels
3604  end if
3605  dv_calculated(i_qv_mean) = .true.
3606  end if
3607  !$acc kernels
3608  var(:) = qv_mean(:)
3609  !$acc end kernels
3610 
3611  case ( 'QHYD_MEAN' )
3612  if ( .not. dv_calculated(i_qhyd_mean) ) then
3613  call allocate_1d( qhyd_mean )
3614  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3615  call atmos_vars_get_diagnostic( 'QHYD', work3d(:,:,:) )
3616 !OCL XFILL
3617  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3618  !$acc kernels
3619  do j = jsb, jeb
3620  do i = isb, ieb
3621  do k = ks, ke
3622  work(k,i,j) = qhyd(k,i,j) * dens_av(k,i,j)
3623  enddo
3624  enddo
3625  enddo
3626  !$acc end kernels
3627  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3628  work(:,:,:), area(:,:), qhyd_mean(:) )
3629  !$acc kernels
3630  !$acc loop independent
3631  do k = ks, ke
3632  qhyd_mean(k) = qhyd_mean(k) / dens_mean(k)
3633  enddo
3634  !$acc end kernels
3635  dv_calculated(i_qhyd_mean) = .true.
3636  end if
3637  !$acc kernels
3638  var(:) = qhyd_mean(:)
3639  !$acc end kernels
3640 
3641  case ( 'QLIQ_MEAN' )
3642  if ( .not. dv_calculated(i_qliq_mean) ) then
3643  call allocate_1d( qliq_mean )
3644  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3645  call atmos_vars_get_diagnostic( 'QLIQ', work3d(:,:,:) )
3646 !OCL XFILL
3647  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3648  !$acc kernels
3649  do j = jsb, jeb
3650  do i = isb, ieb
3651  do k = ks, ke
3652  work(k,i,j) = qliq(k,i,j) * dens_av(k,i,j)
3653  enddo
3654  enddo
3655  enddo
3656  !$acc end kernels
3657  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3658  work(:,:,:), area(:,:), qliq_mean(:) )
3659  !$acc kernels
3660  !$acc loop independent
3661  do k = ks, ke
3662  qliq_mean(k) = qliq_mean(k) / dens_mean(k)
3663  enddo
3664  !$acc end kernels
3665  dv_calculated(i_qliq_mean) = .true.
3666  end if
3667  !$acc kernels
3668  var(:) = qliq_mean(:)
3669  !$acc end kernels
3670 
3671  case ( 'QICE_MEAN' )
3672  if ( .not. dv_calculated(i_qice_mean) ) then
3673  call allocate_1d( qice_mean )
3674  call atmos_vars_get_diagnostic( 'DENS_MEAN', work1d(:) )
3675  call atmos_vars_get_diagnostic( 'QICE', work3d(:,:,:) )
3676 !OCL XFILL
3677  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3678  !$acc kernels
3679  do j = jsb, jeb
3680  do i = isb, ieb
3681  do k = ks, ke
3682  work(k,i,j) = qice(k,i,j) * dens_av(k,i,j)
3683  enddo
3684  enddo
3685  enddo
3686  !$acc end kernels
3687  call statistics_horizontal_mean( ka, ks, ke, ia, is, ie, ja, js, je, &
3688  work(:,:,:), area(:,:), qice_mean(:) )
3689  !$acc kernels
3690  !$acc loop independent
3691  do k = ks, ke
3692  qice_mean(k) = qice_mean(k) / dens_mean(k)
3693  enddo
3694  !$acc end kernels
3695  dv_calculated(i_qice_mean) = .true.
3696  end if
3697  !$acc kernels
3698  var(:) = qice_mean(:)
3699  !$acc end kernels
3700 
3701  case default
3702  log_error("ATMOS_vars_calc_diagnostics",*) 'name is invalid for ATMOS_vars_get_diagnostic_1D: ', trim(vname)
3703  call prc_abort
3704  end select
3705 
3706  !$acc end data
3707 
3708  return
3709  end subroutine atmos_vars_get_diagnostic_1d
3710 
3711  !-----------------------------------------------------------------------------
3713  subroutine atmos_vars_monitor
3714  use scale_monitor, only: &
3715  monitor_put
3716  use scale_atmos_hydrometeor, only: &
3717  i_qv
3718  use mod_atmos_phy_rd_vars, only: &
3719  sflx_lw_up => atmos_phy_rd_sflx_lw_up, &
3720  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
3721  sflx_sw_up => atmos_phy_rd_sflx_sw_up, &
3722  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn, &
3723  tomflx_lw_up => atmos_phy_rd_tomflx_lw_up, &
3724  tomflx_lw_dn => atmos_phy_rd_tomflx_lw_dn, &
3725  tomflx_sw_up => atmos_phy_rd_tomflx_sw_up, &
3726  tomflx_sw_dn => atmos_phy_rd_tomflx_sw_dn
3727  use mod_atmos_phy_sf_vars, only: &
3728  sflx_sh => atmos_phy_sf_sflx_sh, &
3729  sflx_lh => atmos_phy_sf_sflx_lh, &
3730  sflx_engi => atmos_phy_sf_sflx_engi, &
3731  sflx_qtrc => atmos_phy_sf_sflx_qtrc
3732  implicit none
3733 
3734  real(rp) :: rhoq(ka,ia,ja)
3735 
3736  real(rp) :: engflxt (ia,ja) ! total flux [J/m2/s]
3737  real(rp) :: sflx_rd_net(ia,ja) ! net SFC radiation flux [J/m2/s]
3738  real(rp) :: tflx_rd_net(ia,ja) ! net TOM radiation flux [J/m2/s]
3739 
3740  integer :: k, i, j, iq
3741  !---------------------------------------------------------------------------
3742 
3743  !$acc data create(RHOQ, ENGFLXT, SFLX_RD_net, TFLX_RD_net)
3744 
3745  call monitor_put( pv_monit_id(i_dens), dens(:,:,:) )
3746  call monitor_put( pv_monit_id(i_momz), momz(:,:,:) )
3747  call monitor_put( pv_monit_id(i_momx), momx(:,:,:) )
3748  call monitor_put( pv_monit_id(i_momy), momy(:,:,:) )
3749  call monitor_put( pv_monit_id(i_rhot), rhot(:,:,:) )
3750 
3751  !##### Mass Budget #####
3752 
3753  do iq = 1, qa
3754  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3755  !$acc kernels
3756 !OCL XFILL
3757  do j = js, je
3758  do i = is, ie
3759  do k = ks, ke
3760  rhoq(k,i,j) = dens_av(k,i,j) * qtrc_av(k,i,j,iq)
3761  enddo
3762  enddo
3763  enddo
3764  !$acc end kernels
3765 
3766  call monitor_put( qp_monit_id(iq), rhoq(:,:,:) )
3767  enddo
3768 
3769  ! total dry airmass
3770  if ( dv_monit_id(im_qdry) > 0 ) then
3771  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3772  !$acc kernels
3773 !OCL XFILL
3774  do j = js, je
3775  do i = is, ie
3776  do k = ks, ke
3777  rhoq(k,i,j) = dens(k,i,j) * qdry(k,i,j)
3778  enddo
3779  enddo
3780  enddo
3781  !$acc end kernels
3782  call monitor_put( dv_monit_id(im_qdry), rhoq(:,:,:) )
3783  end if
3784 
3785  ! total vapor,liquid,solid tracers
3786  if ( dv_monit_id(im_qtot) > 0 ) then
3787  call atmos_vars_get_diagnostic( 'QTOT', work3d(:,:,:) )
3788  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
3789  !$acc kernels
3790 !OCL XFILL
3791  do j = js, je
3792  do i = is, ie
3793  do k = ks, ke
3794  rhoq(k,i,j) = dens(k,i,j) * qtot(k,i,j)
3795  enddo
3796  enddo
3797  enddo
3798  !$acc end kernels
3799  call monitor_put( dv_monit_id(im_qtot), rhoq(:,:,:) )
3800  end if
3801 
3802  ! total evapolation
3803  if ( moist ) call monitor_put( dv_monit_id(im_evap), sflx_qtrc(:,:,i_qv) )
3804 
3805  ! total precipitation
3806  if ( dv_monit_id(im_prec) > 0 ) then
3807  call monitor_put( dv_monit_id(im_prec), prec(:,:) )
3808  end if
3809 
3810 
3811  !##### Energy Budget #####
3812 
3813  if ( dv_monit_id(im_engt) > 0 ) then
3814  call atmos_vars_get_diagnostic( 'ENGT', work3d(:,:,:) )
3815  call monitor_put( dv_monit_id(im_engt), work3d(:,:,:) )
3816  end if
3817  if ( dv_monit_id(im_engp) > 0 ) then
3818  call atmos_vars_get_diagnostic( 'ENGP', work3d(:,:,:) )
3819  call monitor_put( dv_monit_id(im_engp), work3d(:,:,:) )
3820  end if
3821  if ( dv_monit_id(im_engk) > 0 ) then
3822  call atmos_vars_get_diagnostic( 'ENGK', work3d(:,:,:) )
3823  call monitor_put( dv_monit_id(im_engk), work3d(:,:,:) )
3824  end if
3825  if ( dv_monit_id(im_engi) > 0 ) then
3826  call atmos_vars_get_diagnostic( 'ENGI', work3d(:,:,:) )
3827  call monitor_put( dv_monit_id(im_engi), work3d(:,:,:) )
3828  end if
3829 
3830  ! radiation flux
3831 !OCL XFILL
3832  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3833  !$acc kernels
3834  do j = js, je
3835  do i = is, ie
3836  sflx_rd_net(i,j) = ( sflx_lw_up(i,j) - sflx_lw_dn(i,j) ) &
3837  + ( sflx_sw_up(i,j) - sflx_sw_dn(i,j) )
3838 
3839  tflx_rd_net(i,j) = ( tomflx_lw_up(i,j) - tomflx_lw_dn(i,j) ) &
3840  + ( tomflx_sw_up(i,j) - tomflx_sw_dn(i,j) )
3841 
3842  engflxt(i,j) = sflx_sh(i,j) &
3843  + sflx_engi(i,j) - prec_engi(i,j) &
3844  + sflx_rd_net(i,j) - tflx_rd_net(i,j)
3845  enddo
3846  enddo
3847  !$acc end kernels
3848 
3849  call monitor_put( dv_monit_id(im_engflxt), engflxt(:,:) )
3850 
3851  call monitor_put( dv_monit_id(im_engsfc_sh), sflx_sh(:,:) )
3852  call monitor_put( dv_monit_id(im_engsfc_lh), sflx_lh(:,:) )
3853  call monitor_put( dv_monit_id(im_engsfc_evap), sflx_engi(:,:) )
3854  call monitor_put( dv_monit_id(im_engsfc_prec), prec_engi(:,:) )
3855  call monitor_put( dv_monit_id(im_engsfc_rd), sflx_rd_net(:,:) )
3856  call monitor_put( dv_monit_id(im_engtom_rd), tflx_rd_net(:,:) )
3857 
3858  call monitor_put( dv_monit_id(im_engsfc_lw_up), sflx_lw_up(:,:) )
3859  call monitor_put( dv_monit_id(im_engsfc_lw_dn), sflx_lw_dn(:,:) )
3860  call monitor_put( dv_monit_id(im_engsfc_sw_up), sflx_sw_up(:,:) )
3861  call monitor_put( dv_monit_id(im_engsfc_sw_dn), sflx_sw_dn(:,:) )
3862 
3863  call monitor_put( dv_monit_id(im_engtom_lw_up), tomflx_lw_up(:,:) )
3864  call monitor_put( dv_monit_id(im_engtom_lw_dn), tomflx_lw_dn(:,:) )
3865  call monitor_put( dv_monit_id(im_engtom_sw_up), tomflx_sw_up(:,:) )
3866  call monitor_put( dv_monit_id(im_engtom_sw_dn), tomflx_sw_dn(:,:) )
3867 
3868  !$acc end data
3869 
3870  return
3871  end subroutine atmos_vars_monitor
3872 
3873  !-----------------------------------------------------------------------------
3875  subroutine atmos_vars_finalize
3878  use mod_atmos_admin, only: &
3880  use mod_atmos_dyn_vars, only: &
3882  use mod_atmos_phy_mp_vars, only: &
3884  use mod_atmos_phy_ae_vars, only: &
3886  use mod_atmos_phy_ch_vars, only: &
3888  use mod_atmos_phy_rd_vars, only: &
3890  use mod_atmos_phy_sf_vars, only: &
3892  use mod_atmos_phy_tb_vars, only: &
3894  use mod_atmos_phy_bl_vars, only: &
3896  use mod_atmos_phy_cp_vars, only: &
3898  use mod_atmos_phy_lt_vars, only: &
3900  implicit none
3901  !---------------------------------------------------------------------------
3902 
3903  log_newline
3904  log_info("ATMOS_vars_finalize",*) 'Finalize'
3905 
3906  if ( atmos_use_average ) then
3907  !$acc exit data delete(DENS_avw, MOMZ_avw, MOMX_avw, MOMY_avw, RHOT_avw, QTRC_avw)
3908  deallocate( dens_avw )
3909  deallocate( momz_avw )
3910  deallocate( momx_avw )
3911  deallocate( momy_avw )
3912  deallocate( rhot_avw )
3913  deallocate( qtrc_avw )
3914  endif
3915 
3916  !$acc exit data delete(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC)
3917  deallocate( dens )
3918  deallocate( momz )
3919  deallocate( momx )
3920  deallocate( momy )
3921  deallocate( rhot )
3922  deallocate( qtrc )
3923 
3924  !$acc exit data delete(DENS_tp, MOMZ_tp, RHOU_tp, RHOV_tp, RHOT_tp, RHOH_p, RHOQ_tp)
3925  deallocate( dens_tp )
3926  deallocate( momz_tp )
3927  deallocate( rhou_tp )
3928  deallocate( rhov_tp )
3929  deallocate( rhot_tp )
3930  deallocate( rhoh_p )
3931  deallocate( rhoq_tp )
3932  ! obsolute
3933  !$acc exit data delete(MOMX_tp, MOMY_tp)
3934  deallocate( momx_tp )
3935  deallocate( momy_tp )
3936 
3937 
3938  !$acc exit data delete(W, U, V)
3939  deallocate( w )
3940  deallocate( u )
3941  deallocate( v )
3942 
3943  !$acc exit data delete(POTT, TEMP, PRES, EXNER, PHYD, PHYDH)
3944  deallocate( pott )
3945  deallocate( temp )
3946  deallocate( pres )
3947  deallocate( exner )
3948  deallocate( phyd )
3949  deallocate( phydh )
3950 
3951  !$acc exit data delete(Qdry, Rtot, CVtot, CPtot)
3952  deallocate( qdry )
3953  deallocate( rtot )
3954  deallocate( cvtot)
3955  deallocate( cptot)
3956 
3957  !$acc exit data delete(PREC, PREC_ENGI)
3958  deallocate( prec )
3959  deallocate( prec_engi )
3960 
3961  !$acc exit data delete(WORK3D, WORK2D, WORK1D)
3962  deallocate( work3d )
3963  deallocate( work2d )
3964  deallocate( work1d )
3965 
3976 
3977 
3978  ! water content
3979  if ( atmos_hydrometeor_dry ) then
3980  !$acc exit data delete(ZERO)
3981  deallocate( zero )
3982 
3983  else
3984  !$acc exit data delete(Qe)
3985  deallocate( qe )
3986  end if
3987 
3988  !-----< history output finalize >-----
3989  deallocate( qp_hist_id )
3990  deallocate( qp_monit_id )
3991 
3992  return
3993  end subroutine atmos_vars_finalize
3994 
3995  !-----------------------------------------------------------------------------
3997  subroutine atmos_vars_restart_create
3998  use scale_time, only: &
4000  use scale_file_cartesc, only: &
4002  use mod_atmos_admin, only: &
4003  atmos_sw_dyn, &
4004  atmos_sw_phy_mp, &
4005  atmos_sw_phy_ae, &
4006  atmos_sw_phy_ch, &
4007  atmos_sw_phy_rd, &
4008  atmos_sw_phy_sf, &
4009  atmos_sw_phy_tb, &
4010  atmos_sw_phy_bl, &
4011  atmos_sw_phy_cp, &
4013  use mod_cpl_admin, only: &
4014  cpl_sw
4015  use mod_atmos_dyn_vars, only: &
4017  use mod_atmos_phy_mp_vars, only: &
4019  use mod_atmos_phy_ae_vars, only: &
4021  use mod_atmos_phy_ch_vars, only: &
4023  use mod_atmos_phy_rd_vars, only: &
4025  use mod_atmos_phy_sf_vars, only: &
4027  use mod_atmos_phy_tb_vars, only: &
4029  use mod_atmos_phy_bl_vars, only: &
4031  use mod_atmos_phy_cp_vars, only: &
4033  use mod_atmos_phy_lt_vars, only: &
4035 #ifdef SDM
4036  use scale_atmos_phy_mp_sdm, only: &
4037  sd_rest_flg_out, &
4038  atmos_phy_mp_sdm_restart_create
4039  use scale_time, only: &
4040  nowdaysec => time_nowdaysec
4041 #endif
4042  implicit none
4043 
4044  character(len=19) :: timelabel
4045  character(len=H_LONG) :: basename
4046  !---------------------------------------------------------------------------
4047 
4048  call prof_rapstart('ATM_Restart', 1)
4049 
4050 #ifdef SDM
4051  if( sd_rest_flg_out ) then
4052  log_info("ATMOS_vars_restart_create",*) 'Output random number for SDM '
4053  call atmos_phy_mp_sdm_restart_create(nowdaysec)
4054  endif
4055 #endif
4056 
4057  if ( atmos_restart_out_basename /= '' ) then
4058 
4059  log_newline
4060  log_info("ATMOS_vars_restart_create",*) 'Create restart file (ATMOS) '
4061 
4063  call time_gettimelabel( timelabel )
4064  basename = trim(atmos_restart_out_basename)//'_'//trim(timelabel)
4065  else
4066  basename = trim(atmos_restart_out_basename)
4067  endif
4068 
4069  log_info("ATMOS_vars_restart_create",*) 'basename: ', trim(basename)
4070 
4071  call file_cartesc_create( &
4073  restart_fid, & ! [OUT]
4074  aggregate=atmos_restart_out_aggregate ) ! [IN]
4075 
4076  allocate( pv_id(pv_nmax+qa) )
4077  endif
4078 
4089 
4090  call prof_rapend('ATM_Restart', 1)
4091 
4092  return
4093  end subroutine atmos_vars_restart_create
4094 
4095  !-----------------------------------------------------------------------------
4097  subroutine atmos_vars_restart_enddef
4098  use scale_file_cartesc, only: &
4100  use mod_atmos_admin, only: &
4101  atmos_sw_dyn, &
4102  atmos_sw_phy_mp, &
4103  atmos_sw_phy_ae, &
4104  atmos_sw_phy_ch, &
4105  atmos_sw_phy_rd, &
4106  atmos_sw_phy_sf, &
4107  atmos_sw_phy_tb, &
4108  atmos_sw_phy_bl, &
4109  atmos_sw_phy_cp, &
4111  use mod_cpl_admin, only: &
4112  cpl_sw
4113  use mod_atmos_dyn_vars, only: &
4115  use mod_atmos_phy_mp_vars, only: &
4117  use mod_atmos_phy_ae_vars, only: &
4119  use mod_atmos_phy_ch_vars, only: &
4121  use mod_atmos_phy_rd_vars, only: &
4123  use mod_atmos_phy_sf_vars, only: &
4125  use mod_atmos_phy_tb_vars, only: &
4127  use mod_atmos_phy_bl_vars, only: &
4129  use mod_atmos_phy_cp_vars, only: &
4131  use mod_atmos_phy_lt_vars, only: &
4133 #ifdef SDM
4134  use scale_atmos_phy_mp_sdm, only: &
4135  sd_rest_flg_out, &
4136  atmos_phy_mp_sdm_restart_enddef
4137 #endif
4138  implicit none
4139 
4140  !---------------------------------------------------------------------------
4141 
4142  call prof_rapstart('ATM_Restart', 1)
4143 
4144 #ifdef SDM
4145  if( sd_rest_flg_out ) then
4146  call atmos_phy_mp_sdm_restart_enddef
4147  endif
4148 #endif
4149 
4150  if ( restart_fid /= -1 ) then
4151  call file_cartesc_enddef( restart_fid ) ! [IN]
4152  endif
4153 
4164 
4165  call prof_rapend('ATM_Restart', 1)
4166 
4167  return
4168  end subroutine atmos_vars_restart_enddef
4169 
4170  !-----------------------------------------------------------------------------
4172  subroutine atmos_vars_restart_close
4173  use scale_file_cartesc, only: &
4175  use mod_atmos_admin, only: &
4176  atmos_sw_dyn, &
4177  atmos_sw_phy_mp, &
4178  atmos_sw_phy_ae, &
4179  atmos_sw_phy_ch, &
4180  atmos_sw_phy_rd, &
4181  atmos_sw_phy_sf, &
4182  atmos_sw_phy_tb, &
4183  atmos_sw_phy_bl, &
4184  atmos_sw_phy_cp, &
4186  use mod_cpl_admin, only: &
4187  cpl_sw
4188  use mod_atmos_dyn_vars, only: &
4190  use mod_atmos_phy_mp_vars, only: &
4192  use mod_atmos_phy_ae_vars, only: &
4194  use mod_atmos_phy_ch_vars, only: &
4196  use mod_atmos_phy_rd_vars, only: &
4198  use mod_atmos_phy_sf_vars, only: &
4200  use mod_atmos_phy_tb_vars, only: &
4202  use mod_atmos_phy_bl_vars, only: &
4204  use mod_atmos_phy_cp_vars, only: &
4206  use mod_atmos_phy_lt_vars, only: &
4208 #ifdef SDM
4209  use scale_atmos_phy_mp_sdm, only: &
4210  sd_rest_flg_out, &
4211  atmos_phy_mp_sdm_restart_close
4212 #endif
4213  implicit none
4214  !---------------------------------------------------------------------------
4215 
4216  call prof_rapstart('ATM_Restart', 1)
4217 
4218 #ifdef SDM
4219  if( sd_rest_flg_out ) then
4220  call atmos_phy_mp_sdm_restart_close
4221  endif
4222 #endif
4223 
4224  if ( restart_fid /= -1 ) then
4225  log_newline
4226  log_info("ATMOS_vars_restart_close",*) 'Close restart file (ATMOS) '
4227 
4228  call file_cartesc_close( restart_fid ) ! [IN]
4229 
4230  restart_fid = -1
4231 
4232  if ( allocated(pv_id) ) deallocate( pv_id )
4233  endif
4234 
4240  if( atmos_sw_phy_sf .and. (.not. cpl_sw) ) call atmos_phy_sf_vars_restart_close
4245 
4246  call prof_rapend('ATM_Restart', 1)
4247 
4248  return
4249  end subroutine atmos_vars_restart_close
4250 
4251  !-----------------------------------------------------------------------------
4253  subroutine atmos_vars_restart_def_var
4254  use scale_file_cartesc, only: &
4256  use mod_atmos_admin, only: &
4257  atmos_sw_dyn, &
4258  atmos_sw_phy_mp, &
4259  atmos_sw_phy_ae, &
4260  atmos_sw_phy_ch, &
4261  atmos_sw_phy_rd, &
4262  atmos_sw_phy_sf, &
4263  atmos_sw_phy_tb, &
4264  atmos_sw_phy_bl, &
4265  atmos_sw_phy_cp, &
4267  use mod_cpl_admin, only: &
4268  cpl_sw
4269  use mod_atmos_dyn_vars, only: &
4271  use mod_atmos_phy_mp_vars, only: &
4273  use mod_atmos_phy_ae_vars, only: &
4275  use mod_atmos_phy_ch_vars, only: &
4277  use mod_atmos_phy_rd_vars, only: &
4279  use mod_atmos_phy_sf_vars, only: &
4281  use mod_atmos_phy_tb_vars, only: &
4283  use mod_atmos_phy_bl_vars, only: &
4285  use mod_atmos_phy_cp_vars, only: &
4287  use mod_atmos_phy_lt_vars, only: &
4289 #ifdef SDM
4290  use scale_atmos_phy_mp_sdm, only: &
4291  sd_rest_flg_out, &
4292  atmos_phy_mp_sdm_restart_def_var
4293 #endif
4294  implicit none
4295 
4296  integer iq
4297  !---------------------------------------------------------------------------
4298 
4299  call prof_rapstart('ATM_Restart', 1)
4300 
4301 #ifdef SDM
4302  if( sd_rest_flg_out ) then
4303  call atmos_phy_mp_sdm_restart_def_var
4304  endif
4305 #endif
4306 
4307  if ( restart_fid /= -1 ) then
4308 
4309  call file_cartesc_def_var( restart_fid, pv_info(i_dens)%NAME, pv_info(i_dens)%DESC, pv_info(i_dens)%UNIT, 'ZXY', atmos_restart_out_dtype, &
4310  pv_id(i_dens), &
4311  standard_name=pv_info(i_dens)%STDNAME )
4312  call file_cartesc_def_var( restart_fid, pv_info(i_momz)%NAME, pv_info(i_momz)%DESC, pv_info(i_momz)%UNIT, 'ZHXY', atmos_restart_out_dtype, &
4313  pv_id(i_momz), &
4314  standard_name=pv_info(i_momz)%STDNAME )
4315  call file_cartesc_def_var( restart_fid, pv_info(i_momx)%NAME, pv_info(i_momx)%DESC, pv_info(i_momx)%UNIT, 'ZXHY', atmos_restart_out_dtype, &
4316  pv_id(i_momx), &
4317  standard_name=pv_info(i_momx)%STDNAME )
4318  call file_cartesc_def_var( restart_fid, pv_info(i_momy)%NAME, pv_info(i_momy)%DESC, pv_info(i_momy)%UNIT, 'ZXYH', atmos_restart_out_dtype, &
4319  pv_id(i_momy), &
4320  standard_name=pv_info(i_momy)%STDNAME )
4321  call file_cartesc_def_var( restart_fid, pv_info(i_rhot)%NAME, pv_info(i_rhot)%DESC, pv_info(i_rhot)%UNIT, 'ZXY', atmos_restart_out_dtype, &
4322  pv_id(i_rhot), &
4323  standard_name=pv_info(i_rhot)%STDNAME )
4324  do iq = 1, qa
4325  call file_cartesc_def_var( restart_fid, tracer_name(iq), tracer_desc(iq), tracer_unit(iq), 'ZXY', atmos_restart_out_dtype, &
4326  pv_id(pv_nmax+iq) )
4327  enddo
4328 
4329  endif
4330 
4341 
4342  call prof_rapend('ATM_Restart', 1)
4343 
4344  return
4345  end subroutine atmos_vars_restart_def_var
4346 
4347  !-----------------------------------------------------------------------------
4349  subroutine atmos_vars_restart_write
4350  use scale_file_cartesc, only: &
4351  file_cartesc_write_var
4352  use mod_atmos_admin, only: &
4353  atmos_sw_dyn, &
4354  atmos_sw_phy_mp, &
4355  atmos_sw_phy_ae, &
4356  atmos_sw_phy_ch, &
4357  atmos_sw_phy_rd, &
4358  atmos_sw_phy_sf, &
4359  atmos_sw_phy_tb, &
4360  atmos_sw_phy_bl, &
4361  atmos_sw_phy_cp, &
4363  use mod_cpl_admin, only: &
4364  cpl_sw
4365  use mod_atmos_dyn_vars, only: &
4367  use mod_atmos_phy_mp_vars, only: &
4369  use mod_atmos_phy_ae_vars, only: &
4371  use mod_atmos_phy_ch_vars, only: &
4373  use mod_atmos_phy_rd_vars, only: &
4375  use mod_atmos_phy_sf_vars, only: &
4377  use mod_atmos_phy_tb_vars, only: &
4379  use mod_atmos_phy_bl_vars, only: &
4381  use mod_atmos_phy_cp_vars, only: &
4383  use mod_atmos_phy_lt_vars, only: &
4385 #ifdef SDM
4386  use scale_atmos_phy_mp_sdm, only: &
4387  sd_rest_flg_out, &
4388  atmos_phy_mp_sdm_restart_write
4389 #endif
4390  implicit none
4391 
4392  integer iq
4393  !---------------------------------------------------------------------------
4394 
4395  call prof_rapstart('ATM_Restart', 1)
4396 
4397 #ifdef SDM
4398  if( sd_rest_flg_out ) then
4399  call atmos_phy_mp_sdm_restart_write
4400  endif
4401 #endif
4402 
4403  if ( restart_fid /= -1 ) then
4404 
4405  call atmos_vars_fillhalo
4406 
4407  call atmos_vars_check( force = .true. )
4408 
4409  call file_cartesc_write_var( restart_fid, pv_id(i_dens), dens(:,:,:), pv_info(i_dens)%NAME, 'ZXY' ) ! [IN]
4410  call file_cartesc_write_var( restart_fid, pv_id(i_momz), momz(:,:,:), pv_info(i_momz)%NAME, 'ZHXY' ) ! [IN]
4411  call file_cartesc_write_var( restart_fid, pv_id(i_momx), momx(:,:,:), pv_info(i_momx)%NAME, 'ZXHY' ) ! [IN]
4412  call file_cartesc_write_var( restart_fid, pv_id(i_momy), momy(:,:,:), pv_info(i_momy)%NAME, 'ZXYH' ) ! [IN]
4413  call file_cartesc_write_var( restart_fid, pv_id(i_rhot), rhot(:,:,:), pv_info(i_rhot)%NAME, 'ZXY' ) ! [IN]
4414 
4415  do iq = 1, qa
4416  call file_cartesc_write_var( restart_fid, pv_id(pv_nmax+iq), qtrc(:,:,:,iq), tracer_name(iq), 'ZXY' ) ! [IN]
4417  enddo
4418 
4419  endif
4420 
4426  if( atmos_sw_phy_sf .and. (.not. cpl_sw) ) call atmos_phy_sf_vars_restart_write
4431 
4432  call prof_rapend('ATM_Restart', 1)
4433 
4434  return
4435  end subroutine atmos_vars_restart_write
4436 
4437 
4438  ! private
4439  subroutine allocate_3d( ary )
4440  use scale_const, only: &
4441  undef => const_undef
4442  real(RP), intent(inout), allocatable :: ary(:,:,:)
4443 
4444  if ( .not. allocated(ary) ) then
4445  allocate( ary(ka,ia,ja) )
4446  ary(:,:,:) = undef
4447  !$acc enter data create(ary)
4448  end if
4449 
4450  return
4451  end subroutine allocate_3d
4452 
4453  subroutine allocate_2d( ary )
4454  use scale_const, only: &
4455  undef => const_undef
4456  real(RP), intent(inout), allocatable :: ary(:,:)
4457 
4458  if ( .not. allocated(ary) ) then
4459  allocate( ary(ia,ja) )
4460  ary(:,:) = undef
4461  !$acc enter data create(ary)
4462  end if
4463 
4464  return
4465  end subroutine allocate_2d
4466 
4467  subroutine allocate_1d( ary )
4468  use scale_const, only: &
4469  undef => const_undef
4470  real(RP), intent(inout), allocatable :: ary(:)
4471 
4472  if ( .not. allocated(ary) ) then
4473  allocate( ary(ka) )
4474  ary(:) = undef
4475  !$acc enter data create(ary)
4476  end if
4477 
4478  return
4479  end subroutine allocate_1d
4480 
4481 end module mod_atmos_vars
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_write
subroutine, public atmos_phy_cp_vars_restart_write
Write restart.
Definition: mod_atmos_phy_cp_vars.F90:626
mod_atmos_vars::momz_av
real(rp), dimension(:,:,:), pointer, public momz_av
Definition: mod_atmos_vars.F90:91
scale_atmos_grid_cartesc_index::isb
integer, public isb
Definition: scale_atmos_grid_cartesC_index.F90:64
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_close
subroutine, public atmos_phy_tb_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_tb_vars.F90:332
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_create
subroutine, public atmos_phy_ch_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_ch_vars.F90:284
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
mod_atmos_phy_mp_vars
module Atmosphere / Physics Cloud Microphysics
Definition: mod_atmos_phy_mp_vars.F90:12
mod_atmos_vars::atmos_vars_restart_write
subroutine, public atmos_vars_restart_write
Write restart of atmospheric variables.
Definition: mod_atmos_vars.F90:4350
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_create
subroutine, public atmos_phy_rd_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_rd_vars.F90:380
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_close
subroutine, public atmos_phy_ae_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_ae_vars.F90:421
scale_atmos_grid_cartesc_index::i_uy
integer, public i_uy
Definition: scale_atmos_grid_cartesC_index.F90:100
mod_atmos_vars::qe
real(rp), dimension(:,:,:,:), allocatable, target, public qe
Definition: mod_atmos_vars.F90:105
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_close
subroutine, public atmos_phy_ch_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_ch_vars.F90:335
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
mod_atmos_vars::rhot_avw
real(rp), dimension(:,:,:), allocatable, target, public rhot_avw
Definition: mod_atmos_vars.F90:87
scale_atmos_diagnostic::atmos_diagnostic_get_teml
subroutine, public atmos_diagnostic_get_teml(KA, KS, KE, IA, IS, IE, JA, JS, JE, TEMP, LHV, LHS, QC, QI, CPtot, TEML)
ATMOS_DIAGNOSTIC_get_teml liqued water temperature.
Definition: scale_atmos_diagnostic.F90:338
scale_time::time_nowdaysec
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:72
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_open
subroutine, public atmos_phy_ch_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_ch_vars.F90:203
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_restart_write
subroutine, public atmos_phy_lt_vars_restart_write
Write restart.
Definition: mod_atmos_phy_lt_vars.F90:408
scale_atmos_grid_cartesc_index::i_xv
integer, public i_xv
Definition: scale_atmos_grid_cartesC_index.F90:101
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
mod_atmos_dyn_vars
module Atmosphere / Dynamics
Definition: mod_atmos_dyn_vars.F90:12
scale_atmos_hydrometeor::i_hr
integer, parameter, public i_hr
liquid water rain
Definition: scale_atmos_hydrometeor.F90:98
scale_index::i_momz
integer, parameter, public i_momz
Definition: scale_index.F90:29
mod_atmos_admin::atmos_sw_phy_tb
logical, public atmos_sw_phy_tb
Definition: mod_atmos_admin.F90:57
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_open
subroutine, public atmos_phy_sf_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_sf_vars.F90:470
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_enddef
subroutine, public atmos_phy_ch_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_ch_vars.F90:321
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic
subroutine, public atmos_phy_mp_vars_get_diagnostic(DENS, TEMP, QTRC, CLDFRAC, Re, Qe, Ne)
Definition: mod_atmos_phy_mp_vars.F90:629
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
mod_atmos_vars::qr
real(rp), dimension(:,:,:), pointer, public qr
Definition: mod_atmos_vars.F90:99
mod_atmos_vars::momx_av
real(rp), dimension(:,:,:), pointer, public momx_av
Definition: mod_atmos_vars.F90:92
scale_index::i_momx
integer, parameter, public i_momx
Definition: scale_index.F90:30
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzxv
real(rp), public atmos_grid_cartesc_real_totvolzxv
total volume (zxv, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:91
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:39
mod_atmos_vars::rhoq_tp
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
Definition: mod_atmos_vars.F90:121
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_restart_enddef
subroutine, public atmos_phy_lt_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_lt_vars.F90:358
scale_tracer::tracer_desc
character(len=h_mid), dimension(qa_max), public tracer_desc
Definition: scale_tracer.F90:40
scale_tracer::tracer_unit
character(len=h_short), dimension(qa_max), public tracer_unit
Definition: scale_tracer.F90:41
scale_index
module Index
Definition: scale_index.F90:11
mod_atmos_vars::pott
real(rp), dimension(:,:,:), allocatable, target, public pott
Definition: mod_atmos_vars.F90:133
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_close
subroutine, public atmos_phy_rd_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_rd_vars.F90:431
mod_atmos_dyn_vars::atmos_dyn_vars_finalize
subroutine, public atmos_dyn_vars_finalize
Finalize.
Definition: mod_atmos_dyn_vars.F90:177
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_enddef
subroutine, public atmos_phy_tb_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_tb_vars.F90:318
scale_atmos_hydrometeor::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:100
scale_tracer::tracer_mass
real(rp), dimension(qa_max), public tracer_mass
Definition: scale_tracer.F90:47
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdx
reciprocal of face-dx
Definition: scale_atmos_grid_cartesC.F90:68
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volwxy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volwxy
control volume (wxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:85
mod_atmos_vars::atmos_vars_fillhalo
subroutine, public atmos_vars_fillhalo(FILL_BND)
HALO Communication.
Definition: mod_atmos_vars.F90:879
scale_file_cartesc::file_cartesc_enddef
subroutine, public file_cartesc_enddef(fid)
Exit netCDF file define mode.
Definition: scale_file_cartesC.F90:964
mod_atmos_admin::atmos_use_average
logical, public atmos_use_average
Definition: mod_atmos_admin.F90:49
scale_file_cartesc::file_cartesc_def_var
subroutine, public file_cartesc_def_var(fid, varname, desc, unit, dim_type, datatype, vid, standard_name, timeintv, nsteps, cell_measures)
Define a variable to file.
Definition: scale_file_cartesC.F90:3360
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_reset_diagnostics
subroutine, public atmos_phy_ae_vars_reset_diagnostics
Definition: mod_atmos_phy_ae_vars.F90:598
mod_atmos_vars::qtrc_av
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
Definition: mod_atmos_vars.F90:95
scale_atmos_diagnostic::atmos_diagnostic_get_potv
subroutine, public atmos_diagnostic_get_potv(KA, KS, KE, IA, IS, IE, JA, JS, JE, POTT, Rtot, POTV)
ATMOS_DIAGNOSTIC_get_potv virtual potential temperature.
Definition: scale_atmos_diagnostic.F90:293
mod_atmos_vars::allocate_1d
subroutine allocate_1d(ary)
Definition: mod_atmos_vars.F90:4468
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_restart_def_var
subroutine, public atmos_phy_bl_vars_restart_def_var
Write restart.
Definition: mod_atmos_phy_bl_vars.F90:377
scale_precision
module PRECISION
Definition: scale_precision.F90:14
mod_atmos_vars::rhov_tp
real(rp), dimension(:,:,:), allocatable, public rhov_tp
Definition: mod_atmos_vars.F90:118
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_def_var
subroutine, public atmos_phy_sf_vars_restart_def_var
Write restart.
Definition: mod_atmos_phy_sf_vars.F90:628
mod_atmos_phy_rd_vars::atmos_phy_rd_tomflx_lw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_tomflx_lw_dn
Definition: mod_atmos_phy_rd_vars.F90:68
mod_atmos_vars::atmos_vars_restart_def_var
subroutine, public atmos_vars_restart_def_var
Define atmospheric variables in restart file.
Definition: mod_atmos_vars.F90:4254
scale_file_history_cartesc::file_history_cartesc_set_pres
subroutine, public file_history_cartesc_set_pres(PRES, PRESH, SFC_PRES)
set hydrostatic pressure for pressure coordinate
Definition: scale_file_history_cartesC.F90:278
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_def_var
subroutine, public atmos_phy_ch_vars_restart_def_var
Write restart.
Definition: mod_atmos_phy_ch_vars.F90:355
mod_atmos_vars::qdry
real(rp), dimension(:,:,:), allocatable, target, public qdry
Definition: mod_atmos_vars.F90:140
mod_atmos_phy_lt_vars
module Atmosphere / Physics Chemistry
Definition: mod_atmos_phy_lt_vars.F90:12
mod_atmos_vars::momx_avw
real(rp), dimension(:,:,:), allocatable, target, public momx_avw
Definition: mod_atmos_vars.F90:85
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_create
subroutine, public atmos_phy_ae_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_ae_vars.F90:370
scale_const::const_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:68
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdx
reciprocal of center-dx
Definition: scale_atmos_grid_cartesC.F90:66
mod_atmos_vars::cptot
real(rp), dimension(:,:,:), allocatable, target, public cptot
Definition: mod_atmos_vars.F90:143
mod_atmos_admin::atmos_sw_dyn
logical, public atmos_sw_dyn
Definition: mod_atmos_admin.F90:51
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_open
subroutine, public atmos_phy_cp_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_cp_vars.F90:383
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_create
subroutine, public atmos_phy_cp_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_cp_vars.F90:508
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_read
subroutine, public atmos_phy_tb_vars_restart_read
Read restart.
Definition: mod_atmos_phy_tb_vars.F90:239
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolwxy
real(rp), public atmos_grid_cartesc_real_totvolwxy
total volume (wxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:89
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:174
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
mod_atmos_admin::atmos_sw_phy_mp
logical, public atmos_sw_phy_mp
Definition: mod_atmos_admin.F90:52
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_open
subroutine, public atmos_phy_rd_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_rd_vars.F90:305
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_finalize
subroutine, public atmos_phy_tb_vars_finalize
Finalize.
Definition: mod_atmos_phy_tb_vars.F90:163
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_setup
subroutine, public atmos_phy_lt_vars_setup
Setup.
Definition: mod_atmos_phy_lt_vars.F90:92
mod_atmos_vars::rhot_av
real(rp), dimension(:,:,:), pointer, public rhot_av
Definition: mod_atmos_vars.F90:94
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_setup
subroutine, public atmos_phy_ch_vars_setup
Setup.
Definition: mod_atmos_phy_ch_vars.F90:91
mod_atmos_vars::qh
real(rp), dimension(:,:,:), pointer, public qh
Definition: mod_atmos_vars.F90:103
scale_atmos_diagnostic
module atmosphere / diagnostic
Definition: scale_atmos_diagnostic.F90:12
mod_atmos_vars::phyd
real(rp), dimension(:,:,:), allocatable, target, public phyd
Definition: mod_atmos_vars.F90:137
mod_atmos_phy_rd_vars
module Atmosphere / Physics Radiation
Definition: mod_atmos_phy_rd_vars.F90:12
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
scale_atmos_diagnostic::atmos_diagnostic_get_phyd
subroutine, public atmos_diagnostic_get_phyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, CZ, FZ, PHYD, PHYDH)
ATMOS_DIAGNOSTIC_get_phyd hydrostatic pressure.
Definition: scale_atmos_diagnostic.F90:165
mod_atmos_vars::atmos_vars_restart_enddef
subroutine, public atmos_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_vars.F90:4098
mod_atmos_dyn_vars::atmos_dyn_vars_restart_read
subroutine, public atmos_dyn_vars_restart_read
Read restart.
Definition: mod_atmos_dyn_vars.F90:262
mod_atmos_phy_sf_vars
module ATMOSPHERIC Surface Variables
Definition: mod_atmos_phy_sf_vars.F90:12
mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_rain
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_rain
Definition: mod_atmos_phy_mp_vars.F90:74
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_enddef
subroutine, public atmos_phy_cp_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_cp_vars.F90:545
scale_atmos_hydrometeor::i_hh
integer, parameter, public i_hh
hail
Definition: scale_atmos_hydrometeor.F90:102
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_restart_create
subroutine, public atmos_phy_lt_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_lt_vars.F90:321
mod_atmos_vars::prec_engi
real(rp), dimension(:,:), allocatable, public prec_engi
Definition: mod_atmos_vars.F90:146
scale_index::i_dens
integer, parameter, public i_dens
Definition: scale_index.F90:28
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:114
mod_atmos_phy_rd_vars::atmos_phy_rd_tomflx_sw_up
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_tomflx_sw_up
Definition: mod_atmos_phy_rd_vars.F90:69
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_def_var
subroutine, public atmos_phy_mp_vars_restart_def_var
Define variables in restart file.
Definition: mod_atmos_phy_mp_vars.F90:474
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:80
mod_atmos_vars::atmos_restart_check_criterion
real(rp), public atmos_restart_check_criterion
Definition: mod_atmos_vars.F90:73
mod_atmos_vars::atmos_vars_check
subroutine, public atmos_vars_check(force)
Check variables for atmosphere.
Definition: mod_atmos_vars.F90:1475
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_create
subroutine, public atmos_phy_tb_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_tb_vars.F90:282
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_finalize
subroutine, public atmos_phy_mp_vars_finalize
Finalize.
Definition: mod_atmos_phy_mp_vars.F90:274
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_temp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_temp
Definition: mod_atmos_phy_sf_vars.F90:65
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_finalize
subroutine, public atmos_phy_rd_vars_finalize
Finalize.
Definition: mod_atmos_phy_rd_vars.F90:223
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:81
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_tracer::tracer_engi0
real(rp), dimension(qa_max), public tracer_engi0
Definition: scale_tracer.F90:45
scale_file
module file
Definition: scale_file.F90:15
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_write
subroutine, public atmos_phy_tb_vars_restart_write
Write restart.
Definition: mod_atmos_phy_tb_vars.F90:372
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzuy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzuy
control volume (zuy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:86
mod_atmos_vars::qv_ref
real(rp), dimension(:,:,:), allocatable, public qv_ref
Definition: mod_atmos_vars.F90:112
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_create
subroutine, public atmos_phy_mp_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_mp_vars.F90:403
scale_tracer::unit
character(len=h_short), public unit
Definition: scale_tracer.F90:41
mod_atmos_vars::atmos_restart_check
logical, public atmos_restart_check
Check value consistency?
Definition: mod_atmos_vars.F90:71
scale_atmos_grid_cartesc_index::jeb
integer, public jeb
Definition: scale_atmos_grid_cartesC_index.F90:67
scale_index::i_momy
integer, parameter, public i_momy
Definition: scale_index.F90:31
mod_atmos_vars::j
real(rp), allocatable, target, public j
Definition: mod_atmos_vars.F90:141
mod_atmos_vars::atmos_vars_get_diagnostic_1d
recursive subroutine atmos_vars_get_diagnostic_1d(vname, var)
get diagnostic variable 1D
Definition: mod_atmos_vars.F90:3408
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_atmos_vars::atmos_restart_out_basename
character(len=h_long), public atmos_restart_out_basename
Basename of the output file.
Definition: mod_atmos_vars.F90:65
mod_atmos_phy_cp_vars::atmos_phy_cp_sflx_rain
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_sflx_rain
Definition: mod_atmos_phy_cp_vars.F90:63
mod_atmos_admin::atmos_sw_phy_ae
logical, public atmos_sw_phy_ae
Definition: mod_atmos_admin.F90:53
mod_atmos_dyn_vars::atmos_dyn_vars_restart_open
subroutine, public atmos_dyn_vars_restart_open
Open restart file for read.
Definition: mod_atmos_dyn_vars.F90:225
mod_atmos_vars::atmos_vars_monitor
subroutine, public atmos_vars_monitor
monitor output
Definition: mod_atmos_vars.F90:3714
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_history
subroutine, public atmos_phy_ae_vars_history(QTRC, RH)
Definition: mod_atmos_phy_ae_vars.F90:477
mod_atmos_vars::rhou_tp
real(rp), dimension(:,:,:), allocatable, public rhou_tp
Definition: mod_atmos_vars.F90:117
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_mapf
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_mapf
map factor
Definition: scale_atmos_grid_cartesC_metric.F90:35
mod_atmos_vars::atmos_vars_get_diagnostic_2d
recursive subroutine atmos_vars_get_diagnostic_2d(vname, var)
get diagnostic variable 2D
Definition: mod_atmos_vars.F90:3132
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_diagnostic_cartesc::atmos_diagnostic_cartesc_get_vel
subroutine, public atmos_diagnostic_cartesc_get_vel(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, MOMZ, MOMX, MOMY, W, U, V)
ATMOS_DIAGNOSTIC_CARTESC_get_vel W, U, V.
Definition: scale_atmos_diagnostic_cartesC.F90:58
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:99
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
mod_atmos_phy_cp_vars::atmos_phy_cp_sflx_snow
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_sflx_snow
Definition: mod_atmos_phy_cp_vars.F90:64
mod_atmos_vars::prec
real(rp), dimension(:,:), allocatable, target, public prec
Definition: mod_atmos_vars.F90:145
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_setup
subroutine, public atmos_phy_tb_vars_setup
Setup.
Definition: mod_atmos_phy_tb_vars.F90:84
scale_io
module STDIO
Definition: scale_io.F90:10
mod_atmos_vars::atmos_restart_out_title
character(len=h_mid), public atmos_restart_out_title
Title of the output file.
Definition: mod_atmos_vars.F90:68
mod_atmos_vars::atmos_vars_restart_read
subroutine, public atmos_vars_restart_read
Read restart of atmospheric variables.
Definition: mod_atmos_vars.F90:1067
mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_lw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_dn
Definition: mod_atmos_phy_rd_vars.F90:62
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdy
reciprocal of center-dy
Definition: scale_atmos_grid_cartesC.F90:67
mod_atmos_vars::atmos_restart_out_aggregate
logical, public atmos_restart_out_aggregate
Switch to use aggregate file.
Definition: mod_atmos_vars.F90:66
mod_atmos_phy_bl_vars
module atmosphere / physics / PBL
Definition: mod_atmos_phy_bl_vars.F90:12
mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_sw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
Definition: mod_atmos_phy_rd_vars.F90:64
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:76
mod_atmos_dyn_vars::atmos_dyn_vars_restart_create
subroutine, public atmos_dyn_vars_restart_create
Create restart file.
Definition: mod_atmos_dyn_vars.F90:317
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_reset_diagnostics
subroutine, public atmos_phy_mp_vars_reset_diagnostics
Definition: mod_atmos_phy_mp_vars.F90:872
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_write
subroutine, public atmos_phy_sf_vars_restart_write
Write variables to restart file.
Definition: mod_atmos_phy_sf_vars.F90:653
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_close
subroutine, public atmos_phy_sf_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_sf_vars.F90:608
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzxv
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzxv
control volume (zxv) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:87
mod_atmos_vars::allocate_2d
subroutine allocate_2d(ary)
Definition: mod_atmos_vars.F90:4454
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:84
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_enddef
subroutine, public atmos_phy_sf_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_sf_vars.F90:594
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_finalize
subroutine, public atmos_phy_sf_vars_finalize
Finalize.
Definition: mod_atmos_phy_sf_vars.F90:364
mod_atmos_vars::atmos_vars_setup
subroutine, public atmos_vars_setup
Setup.
Definition: mod_atmos_vars.F90:468
mod_atmos_vars::atmos_vars_history_setpres
subroutine, public atmos_vars_history_setpres
Set pressure for history output.
Definition: mod_atmos_vars.F90:1214
mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_sw_up
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_up
Definition: mod_atmos_phy_rd_vars.F90:63
mod_atmos_vars::atmos_vars_restart_check
subroutine, public atmos_vars_restart_check
Check and compare between last data and sample data.
Definition: mod_atmos_vars.F90:1251
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:59
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_restart_read
subroutine, public atmos_phy_lt_vars_restart_read
Read restart.
Definition: mod_atmos_phy_lt_vars.F90:259
mod_atmos_vars::temp_ref
real(rp), dimension(:,:,:), allocatable, public temp_ref
Definition: mod_atmos_vars.F90:110
mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_lw_up
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_up
Definition: mod_atmos_phy_rd_vars.F90:61
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:77
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_setup
subroutine, public atmos_phy_sf_vars_setup
Setup.
Definition: mod_atmos_phy_sf_vars.F90:191
scale_atmos_hydrometeor::lhv
real(rp), public lhv
latent heat of vaporization for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:144
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_def_var
subroutine, public atmos_phy_rd_vars_restart_def_var
Define variables in restart file.
Definition: mod_atmos_phy_rd_vars.F90:451
mod_atmos_admin::atmos_sw_phy_lt
logical, public atmos_sw_phy_lt
Definition: mod_atmos_admin.F90:60
scale_atmos_adiabat
module atmosphere / adiabat
Definition: scale_atmos_adiabat.F90:11
mod_atmos_phy_tb_vars
module Atmosphere / Physics Turbulence
Definition: mod_atmos_phy_tb_vars.F90:12
mod_atmos_vars::pres_ref
real(rp), dimension(:,:,:), allocatable, public pres_ref
Definition: mod_atmos_vars.F90:111
scale_const::const_cvdry
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:61
mod_atmos_vars::atmos_restart_out_dtype
character(len=h_short), public atmos_restart_out_dtype
REAL4 or REAL8.
Definition: mod_atmos_vars.F90:69
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_close
subroutine, public atmos_phy_cp_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_cp_vars.F90:559
mod_atmos_vars::momy_av
real(rp), dimension(:,:,:), pointer, public momy_av
Definition: mod_atmos_vars.F90:93
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_finalize
subroutine, public atmos_phy_ae_vars_finalize
Finalize.
Definition: mod_atmos_phy_ae_vars.F90:225
scale_tracer::tracer_cv
real(rp), dimension(qa_max), public tracer_cv
Definition: scale_tracer.F90:42
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_write
subroutine, public atmos_phy_rd_vars_restart_write
Write variables to restart file.
Definition: mod_atmos_phy_rd_vars.F90:475
scale_file_cartesc::file_cartesc_close
subroutine, public file_cartesc_close(fid)
Close a netCDF file.
Definition: scale_file_cartesC.F90:1044
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_read
subroutine, public atmos_phy_ae_vars_restart_read
Read restart.
Definition: mod_atmos_phy_ae_vars.F90:322
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_restart_def_var
subroutine, public atmos_phy_lt_vars_restart_def_var
Write restart.
Definition: mod_atmos_phy_lt_vars.F90:392
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_finalize
subroutine, public atmos_phy_bl_vars_finalize
Finalize.
Definition: mod_atmos_phy_bl_vars.F90:193
mod_atmos_vars::v
real(rp), dimension(:,:,:), allocatable, target, public v
Definition: mod_atmos_vars.F90:131
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_restart_create
subroutine, public atmos_phy_bl_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_bl_vars.F90:307
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_setup
subroutine, public atmos_phy_cp_vars_setup
Setup.
Definition: mod_atmos_phy_cp_vars.F90:130
scale_tracer::tracer_name
character(len=h_short), dimension(qa_max), public tracer_name
Definition: scale_tracer.F90:39
mod_atmos_vars::w
real(rp), dimension(:,:,:), allocatable, target, public w
Definition: mod_atmos_vars.F90:129
mod_atmos_vars::qtrc_avw
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc_avw
Definition: mod_atmos_vars.F90:88
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_area
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_area
horizontal area ( xy, normal z) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:66
mod_atmos_phy_ch_vars
module Atmosphere / Physics Chemistry
Definition: mod_atmos_phy_ch_vars.F90:12
mod_atmos_vars::momy_avw
real(rp), dimension(:,:,:), allocatable, target, public momy_avw
Definition: mod_atmos_vars.F90:86
mod_atmos_vars::momz_tp
real(rp), dimension(:,:,:), allocatable, public momz_tp
Definition: mod_atmos_vars.F90:116
scale_index::i_rhot
integer, parameter, public i_rhot
Definition: scale_index.F90:32
mod_atmos_vars::momx
real(rp), dimension(:,:,:), allocatable, target, public momx
Definition: mod_atmos_vars.F90:78
mod_atmos_vars::exner
real(rp), dimension(:,:,:), allocatable, target, public exner
Definition: mod_atmos_vars.F90:136
mod_atmos_phy_rd_vars::atmos_phy_rd_tomflx_sw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_tomflx_sw_dn
Definition: mod_atmos_phy_rd_vars.F90:70
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_create
subroutine, public atmos_phy_sf_vars_restart_create
Create restart file.
Definition: mod_atmos_phy_sf_vars.F90:558
mod_atmos_admin::atmos_sw_phy_cp
logical, public atmos_sw_phy_cp
Definition: mod_atmos_admin.F90:59
mod_atmos_admin::atmos_sw_phy_ch
logical, public atmos_sw_phy_ch
Definition: mod_atmos_admin.F90:54
scale_monitor::monitor_reg
subroutine, public monitor_reg(name, desc, unit, itemid, ndims, dim_type, is_tendency)
Search existing item, or matching check between requested and registered item.
Definition: scale_monitor.F90:243
mod_atmos_vars::temp
real(rp), dimension(:,:,:), allocatable, target, public temp
Definition: mod_atmos_vars.F90:134
scale_time::time_dtsec_atmos_dyn
real(dp), public time_dtsec_atmos_dyn
time interval of dynamics [sec]
Definition: scale_time.F90:35
mod_atmos_vars::dens_tp
real(rp), dimension(:,:,:), allocatable, public dens_tp
Definition: mod_atmos_vars.F90:115
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_finalize
subroutine, public atmos_phy_ch_vars_finalize
Finalize.
Definition: mod_atmos_phy_ch_vars.F90:165
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzuy
real(rp), public atmos_grid_cartesc_real_totvolzuy
total volume (zuy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:90
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_write
subroutine, public atmos_phy_mp_vars_restart_write
Write restart.
Definition: mod_atmos_phy_mp_vars.F90:496
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_atmos_vars::momy
real(rp), dimension(:,:,:), allocatable, target, public momy
Definition: mod_atmos_vars.F90:79
scale_time
module TIME
Definition: scale_time.F90:11
mod_atmos_vars::qv
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
Definition: mod_atmos_vars.F90:97
mod_atmos_dyn_vars::atmos_dyn_vars_restart_write
subroutine, public atmos_dyn_vars_restart_write
Write variables to restart file.
Definition: mod_atmos_dyn_vars.F90:412
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, target, public atmos_phy_sf_sflx_qtrc
Definition: mod_atmos_phy_sf_vars.F90:86
scale_atmos_hydrometeor::i_hc
integer, parameter, public i_hc
liquid water cloud
Definition: scale_atmos_hydrometeor.F90:97
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_restart_open
subroutine, public atmos_phy_bl_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_bl_vars.F90:236
mod_atmos_vars::phydh
real(rp), dimension(:,:,:), allocatable, target, public phydh
Definition: mod_atmos_vars.F90:138
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_def_var
subroutine, public atmos_phy_tb_vars_restart_def_var
Write restart.
Definition: mod_atmos_phy_tb_vars.F90:352
mod_atmos_vars::atmos_vars_restart_open
subroutine, public atmos_vars_restart_open
Open restart file for reading atmospheric variables.
Definition: mod_atmos_vars.F90:949
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_setup
subroutine, public atmos_phy_rd_vars_setup
Setup.
Definition: mod_atmos_phy_rd_vars.F90:122
scale_atmos_diagnostic_cartesc
module atmosphere / diagnostic / CartesianC
Definition: scale_atmos_diagnostic_cartesC.F90:12
scale_atmos_diagnostic::atmos_diagnostic_get_therm_rhot
subroutine, public atmos_diagnostic_get_therm_rhot(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, RHOT, Rtot, CVtot, CPtot, POTT, TEMP, PRES, EXNER)
ATMOS_DIAGNOSTIC_get_therm_rhot potential temperature, temperature, pressure.
Definition: scale_atmos_diagnostic.F90:59
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_write
subroutine, public atmos_phy_ae_vars_restart_write
Write restart.
Definition: mod_atmos_phy_ae_vars.F90:457
mod_atmos_vars::pres
real(rp), dimension(:,:,:), allocatable, target, public pres
Definition: mod_atmos_vars.F90:135
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_def_var
subroutine, public atmos_phy_ae_vars_restart_def_var
Write restart.
Definition: mod_atmos_phy_ae_vars.F90:441
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_engi
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_engi
Definition: mod_atmos_phy_sf_vars.F90:87
mod_atmos_admin::atmos_sw_phy_rd
logical, public atmos_sw_phy_rd
Definition: mod_atmos_admin.F90:55
mod_atmos_vars::dens_av
real(rp), dimension(:,:,:), pointer, public dens_av
Definition: mod_atmos_vars.F90:90
mod_atmos_dyn_vars::atmos_dyn_vars_restart_enddef
subroutine, public atmos_dyn_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_dyn_vars.F90:356
mod_atmos_vars::rtot
real(rp), dimension(:,:,:), allocatable, target, public rtot
Definition: mod_atmos_vars.F90:141
mod_atmos_phy_rd_vars::atmos_phy_rd_tomflx_lw_up
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_tomflx_lw_up
Definition: mod_atmos_phy_rd_vars.F90:67
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_read
subroutine, public atmos_phy_rd_vars_restart_read
Read restart.
Definition: mod_atmos_phy_rd_vars.F90:340
mod_atmos_vars::u
real(rp), dimension(:,:,:), allocatable, target, public u
Definition: mod_atmos_vars.F90:130
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_f2h
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_real_f2h
coefficient for interpolation from full to half levels
Definition: scale_atmos_grid_cartesC_real.F90:47
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_write
subroutine, public atmos_phy_ch_vars_restart_write
Write restart.
Definition: mod_atmos_phy_ch_vars.F90:371
mod_atmos_vars::atmos_restart_out_postfix_timelabel
logical, public atmos_restart_out_postfix_timelabel
Add timelabel to the basename of output file?
Definition: mod_atmos_vars.F90:67
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_setup
subroutine, public atmos_phy_mp_vars_setup
Setup.
Definition: mod_atmos_phy_mp_vars.F90:127
scale_debug
module DEBUG
Definition: scale_debug.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_lh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
Definition: mod_atmos_phy_sf_vars.F90:81
scale_tracer::tracer_cp
real(rp), dimension(qa_max), public tracer_cp
Definition: scale_tracer.F90:43
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_read
subroutine, public atmos_phy_mp_vars_restart_read
Read restart.
Definition: mod_atmos_phy_mp_vars.F90:368
mod_atmos_dyn_vars::atmos_dyn_vars_restart_close
subroutine, public atmos_dyn_vars_restart_close
Close restart file.
Definition: mod_atmos_dyn_vars.F90:370
scale_file_history_cartesc
module file / history_cartesC
Definition: scale_file_history_cartesC.F90:12
mod_atmos_vars::atmos_vars_restart_close
subroutine, public atmos_vars_restart_close
Close restart file.
Definition: mod_atmos_vars.F90:4173
mod_atmos_vars::atmos_restart_in_postfix_timelabel
logical, public atmos_restart_in_postfix_timelabel
Add timelabel to the basename of input file?
Definition: mod_atmos_vars.F90:64
scale_file_cartesc::file_cartesc_create
subroutine, public file_cartesc_create(basename, title, datatype, fid, date, subsec, haszcoord, append, aggregate, single)
Create/open a netCDF file.
Definition: scale_file_cartesC.F90:796
mod_atmos_vars::rhoh_p
real(rp), dimension(:,:,:), allocatable, public rhoh_p
Definition: mod_atmos_vars.F90:120
mod_atmos_dyn_vars::atmos_dyn_vars_restart_def_var
subroutine, public atmos_dyn_vars_restart_def_var
Define variables in restart file.
Definition: mod_atmos_dyn_vars.F90:390
mod_atmos_admin::atmos_dyn_type
character(len=h_short), public atmos_dyn_type
Definition: mod_atmos_admin.F90:35
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:109
scale_file_cartesc::file_cartesc_flush
subroutine, public file_cartesc_flush(fid)
Flush all pending requests to a netCDF file (PnetCDF only)
Definition: scale_file_cartesC.F90:1018
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_setup
subroutine, public atmos_phy_ae_vars_setup
Setup.
Definition: mod_atmos_phy_ae_vars.F90:110
scale_time::time_gettimelabel
subroutine, public time_gettimelabel(timelabel)
generate time label
Definition: scale_time.F90:93
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
mod_atmos_vars::atmos_vars_finalize
subroutine, public atmos_vars_finalize
finalize
Definition: mod_atmos_vars.F90:3876
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_restart_read
subroutine, public atmos_phy_bl_vars_restart_read
Read restart.
Definition: mod_atmos_phy_bl_vars.F90:271
scale_atmos_bottom
module atmosphere / bottom boundary extrapolation
Definition: scale_atmos_bottom.F90:12
mod_atmos_vars::qi
real(rp), dimension(:,:,:), pointer, public qi
Definition: mod_atmos_vars.F90:100
mod_atmos_vars::atmos_vars_history
subroutine, public atmos_vars_history
History output set for atmospheric variables.
Definition: mod_atmos_vars.F90:1393
scale_tracer::tracer_r
real(rp), dimension(qa_max), public tracer_r
Definition: scale_tracer.F90:44
mod_atmos_vars::rhot_tp
real(rp), dimension(:,:,:), allocatable, public rhot_tp
Definition: mod_atmos_vars.F90:119
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
mod_atmos_vars::atmos_vars_restart_create
subroutine, public atmos_vars_restart_create
Create atmospheric restart file.
Definition: mod_atmos_vars.F90:3998
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_enddef
subroutine, public atmos_phy_ae_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_ae_vars.F90:407
mod_atmos_vars::atmos_restart_in_aggregate
logical, public atmos_restart_in_aggregate
Switch to use aggregate file.
Definition: mod_atmos_vars.F90:63
mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_snow
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_snow
Definition: mod_atmos_phy_mp_vars.F90:75
mod_atmos_vars::atmos_restart_output
logical, public atmos_restart_output
Output restart file?
Definition: mod_atmos_vars.F90:60
mod_atmos_vars::momy_tp
real(rp), dimension(:,:,:), allocatable, public momy_tp
Definition: mod_atmos_vars.F90:125
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
Definition: mod_atmos_phy_sf_vars.F90:80
mod_atmos_admin::atmos_sw_phy_bl
logical, public atmos_sw_phy_bl
Definition: mod_atmos_admin.F90:58
mod_atmos_phy_tb_vars::atmos_phy_tb_vars_restart_open
subroutine, public atmos_phy_tb_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_tb_vars.F90:204
scale_atmos_grid_cartesc_index::ieb
integer, public ieb
Definition: scale_atmos_grid_cartesC_index.F90:65
scale_atmos_diagnostic::atmos_diagnostic_get_n2
subroutine, public atmos_diagnostic_get_n2(KA, KS, KE, IA, IS, IE, JA, JS, JE, POTT, Rtot, CZ, FZ, F2H, N2)
ATMOS_DIAGNOSTIC_get_n2 N^2.
Definition: scale_atmos_diagnostic.F90:231
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:146
scale_file::file_get_aggregate
logical function, public file_get_aggregate(fid)
Definition: scale_file.F90:6316
mod_atmos_phy_ae_vars
module ATMOSPHERE / Physics Aerosol Microphysics
Definition: mod_atmos_phy_ae_vars.F90:12
scale_file_cartesc::file_cartesc_open
subroutine, public file_cartesc_open(basename, fid, single, aggregate)
open a netCDF file for read
Definition: scale_file_cartesC.F90:760
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_finalize
subroutine, public atmos_phy_cp_vars_finalize
Setup.
Definition: mod_atmos_phy_cp_vars.F90:277
mod_atmos_vars::momx_tp
real(rp), dimension(:,:,:), allocatable, public momx_tp
Definition: mod_atmos_vars.F90:124
scale_atmos_thermodyn
module atmosphere / thermodyn
Definition: scale_atmos_thermodyn.F90:11
mod_atmos_vars::atmos_vars_calc_diagnostics
subroutine, public atmos_vars_calc_diagnostics
Calc diagnostic variables.
Definition: mod_atmos_vars.F90:1747
scale_file_history::file_history_reg
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
Definition: scale_file_history.F90:685
mod_atmos_vars::atmos_vars_get_diagnostic_3d
recursive subroutine atmos_vars_get_diagnostic_3d(vname, var)
get diagnostic variable 3D
Definition: mod_atmos_vars.F90:1824
mod_atmos_phy_ae_vars::atmos_phy_ae_vars_restart_open
subroutine, public atmos_phy_ae_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_ae_vars.F90:287
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
mod_atmos_vars::allocate_3d
subroutine allocate_3d(ary)
Definition: mod_atmos_vars.F90:4440
mod_atmos_vars::qc
real(rp), dimension(:,:,:), pointer, public qc
Definition: mod_atmos_vars.F90:98
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_restart_write
subroutine, public atmos_phy_bl_vars_restart_write
Write restart.
Definition: mod_atmos_phy_bl_vars.F90:400
mod_atmos_vars::atmos_restart_check_basename
character(len=h_long), public atmos_restart_check_basename
Definition: mod_atmos_vars.F90:72
mod_atmos_phy_rd_vars::atmos_phy_rd_vars_restart_enddef
subroutine, public atmos_phy_rd_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_rd_vars.F90:417
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_restart_close
subroutine, public atmos_phy_bl_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_bl_vars.F90:357
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_open
subroutine, public atmos_phy_mp_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_mp_vars.F90:333
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_finalize
subroutine, public atmos_phy_lt_vars_finalize
Finalize.
Definition: mod_atmos_phy_lt_vars.F90:172
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:43
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_setup
subroutine, public atmos_phy_bl_vars_setup
Setup.
Definition: mod_atmos_phy_bl_vars.F90:101
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol
real(rp), public atmos_grid_cartesc_real_totvol
total volume (zxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:88
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
mod_atmos_vars::pott_ref
real(rp), dimension(:,:,:), allocatable, public pott_ref
Definition: mod_atmos_vars.F90:109
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_restart_open
subroutine, public atmos_phy_lt_vars_restart_open
Open restart file for read.
Definition: mod_atmos_phy_lt_vars.F90:224
mod_atmos_vars::dens_avw
real(rp), dimension(:,:,:), allocatable, target, public dens_avw
Definition: mod_atmos_vars.F90:83
mod_atmos_vars::momz_avw
real(rp), dimension(:,:,:), allocatable, target, public momz_avw
Definition: mod_atmos_vars.F90:84
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdy
reciprocal of face-dy
Definition: scale_atmos_grid_cartesC.F90:69
scale_atmos_grid_cartesc_index::jsb
integer, public jsb
Definition: scale_atmos_grid_cartesC_index.F90:66
mod_atmos_phy_sf_vars::atmos_phy_sf_vars_restart_read
subroutine, public atmos_phy_sf_vars_restart_read
Read restart.
Definition: mod_atmos_phy_sf_vars.F90:506
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_close
subroutine, public atmos_phy_mp_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_mp_vars.F90:454
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
mod_atmos_phy_lt_vars::atmos_phy_lt_vars_restart_close
subroutine, public atmos_phy_lt_vars_restart_close
Close restart file.
Definition: mod_atmos_phy_lt_vars.F90:372
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_history
subroutine, public atmos_phy_mp_vars_history(DENS, TEMP, QTRC)
Definition: mod_atmos_phy_mp_vars.F90:518
mod_atmos_phy_cp_vars
module Atmosphere / Physics Cumulus
Definition: mod_atmos_phy_cp_vars.F90:12
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:56
scale_file_cartesc
module file / cartesianC
Definition: scale_file_cartesC.F90:11
mod_cpl_admin::cpl_sw
logical, public cpl_sw
Definition: mod_cpl_admin.F90:33
scale_tracer::name
character(len=h_short), public name
Definition: scale_tracer.F90:39
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95
mod_atmos_vars::atmos_restart_in_basename
character(len=h_long), public atmos_restart_in_basename
Basename of the input file.
Definition: mod_atmos_vars.F90:62
mod_atmos_vars::dens_ref
real(rp), dimension(:,:,:), allocatable, public dens_ref
Definition: mod_atmos_vars.F90:108
mod_atmos_vars::cvtot
real(rp), dimension(:,:,:), allocatable, target, public cvtot
Definition: mod_atmos_vars.F90:142
scale_atmos_bottom::atmos_bottom_estimate
subroutine, public atmos_bottom_estimate(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, QV, SFC_TEMP, FZ, SFC_DENS, SFC_PRES)
Calc bottom boundary of atmosphere (just above surface)
Definition: scale_atmos_bottom.F90:51
mod_atmos_vars::qg
real(rp), dimension(:,:,:), pointer, public qg
Definition: mod_atmos_vars.F90:102
mod_atmos_phy_ch_vars::atmos_phy_ch_vars_restart_read
subroutine, public atmos_phy_ch_vars_restart_read
Read restart.
Definition: mod_atmos_phy_ch_vars.F90:238
mod_atmos_phy_bl_vars::atmos_phy_bl_vars_restart_enddef
subroutine, public atmos_phy_bl_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_bl_vars.F90:343
mod_cpl_admin
module Coupler admin
Definition: mod_cpl_admin.F90:11
mod_atmos_phy_mp_vars::atmos_phy_mp_vars_restart_enddef
subroutine, public atmos_phy_mp_vars_restart_enddef
Exit netCDF define mode.
Definition: mod_atmos_phy_mp_vars.F90:440
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_def_var
subroutine, public atmos_phy_cp_vars_restart_def_var
Write restart.
Definition: mod_atmos_phy_cp_vars.F90:579
mod_atmos_dyn_vars::atmos_dyn_vars_setup
subroutine, public atmos_dyn_vars_setup
Setup.
Definition: mod_atmos_dyn_vars.F90:81
mod_atmos_vars::qs
real(rp), dimension(:,:,:), pointer, public qs
Definition: mod_atmos_vars.F90:101
scale_atmos_hydrometeor::i_hg
integer, parameter, public i_hg
graupel
Definition: scale_atmos_hydrometeor.F90:101
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
Definition: scale_atmos_grid_cartesC_metric.F90:36
mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_read
subroutine, public atmos_phy_cp_vars_restart_read
Read restart.
Definition: mod_atmos_phy_cp_vars.F90:418
mod_atmos_admin::atmos_sw_phy_sf
logical, public atmos_sw_phy_sf
Definition: mod_atmos_admin.F90:56
scale_monitor
module MONITOR
Definition: scale_monitor.F90:12