SCALE-RM
mod_atmos_vars.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
15 !-------------------------------------------------------------------------------
16 #include "inc_openmp.h"
18  !-----------------------------------------------------------------------------
19  !
20  !++ used modules
21  !
22  use scale_precision
23  use scale_stdio
24  use scale_prof
25  use scale_debug
27  use scale_index
28  use scale_tracer
29  !-----------------------------------------------------------------------------
30  implicit none
31  private
32  !-----------------------------------------------------------------------------
33  !
34  !++ Public procedure
35  !
36  public :: atmos_vars_setup
37  public :: atmos_vars_fillhalo
38  public :: atmos_vars_restart_read
39  public :: atmos_vars_restart_write
40  public :: atmos_vars_restart_check
41  public :: atmos_vars_history
42  public :: atmos_vars_total
43  public :: atmos_vars_diagnostics
44  public :: atmos_vars_monitor
45 
46  !-----------------------------------------------------------------------------
47  !
48  !++ Public parameters & variables
49  !
50  logical, public :: atmos_restart_output = .false.
51  logical, public :: atmos_restart_check = .false.
52 
53  character(len=H_LONG), public :: atmos_restart_in_basename = ''
54  character(len=H_LONG), public :: atmos_restart_out_basename = ''
55  character(len=H_MID), public :: atmos_restart_out_title = 'ATMOS restart'
56  character(len=H_MID), public :: atmos_restart_out_dtype = 'DEFAULT'
57  logical, public :: atmos_restart_in_allowmissingq = .false.
58 
59  character(len=H_LONG), public :: atmos_restart_check_basename = 'restart_check'
60  real(RP), public :: atmos_restart_check_criterion = 1.e-6_rp
61 
62  ! prognostic variables
63  real(RP), public, target, allocatable :: dens(:,:,:) ! Density [kg/m3]
64  real(RP), public, target, allocatable :: momz(:,:,:) ! momentum z [kg/m2/s]
65  real(RP), public, target, allocatable :: momx(:,:,:) ! momentum x [kg/m2/s]
66  real(RP), public, target, allocatable :: momy(:,:,:) ! momentum y [kg/m2/s]
67  real(RP), public, target, allocatable :: rhot(:,:,:) ! DENS * POTT [K*kg/m3]
68  real(RP), public, target, allocatable :: qtrc(:,:,:,:) ! ratio of mass of tracer to total mass[kg/kg]
69 
70  real(RP), public, target, allocatable :: dens_avw(:,:,:)
71  real(RP), public, target, allocatable :: momz_avw(:,:,:)
72  real(RP), public, target, allocatable :: momx_avw(:,:,:)
73  real(RP), public, target, allocatable :: momy_avw(:,:,:)
74  real(RP), public, target, allocatable :: rhot_avw(:,:,:)
75  real(RP), public, target, allocatable :: qtrc_avw(:,:,:,:)
76 
77  real(RP), public, pointer :: dens_av(:,:,:)
78  real(RP), public, pointer :: momz_av(:,:,:)
79  real(RP), public, pointer :: momx_av(:,:,:)
80  real(RP), public, pointer :: momy_av(:,:,:)
81  real(RP), public, pointer :: rhot_av(:,:,:)
82  real(RP), public, pointer :: qtrc_av(:,:,:,:)
83 
84  ! tendency by physical processes
85  real(RP), public, allocatable :: dens_tp(:,:,:)
86  real(RP), public, allocatable :: momz_tp(:,:,:)
87  real(RP), public, allocatable :: momx_tp(:,:,:)
88  real(RP), public, allocatable :: momy_tp(:,:,:)
89  real(RP), public, allocatable :: rhot_tp(:,:,:)
90  real(RP), public, allocatable :: rhoq_tp(:,:,:,:)
91 
92  ! diagnostic variables
93  real(RP), public, allocatable :: temp(:,:,:) ! temperature [K]
94  real(RP), public, allocatable :: pres(:,:,:) ! pressure [Pa=J/m3]
95  real(RP), public, allocatable :: w (:,:,:) ! velocity w [m/s]
96  real(RP), public, allocatable :: u (:,:,:) ! velocity u [m/s]
97  real(RP), public, allocatable :: v (:,:,:) ! velocity v [m/s]
98  real(RP), public, allocatable :: pott(:,:,:) ! potential temperature [K]
99 
100  !-----------------------------------------------------------------------------
101  !
102  !++ Private procedure
103  !
104  !-----------------------------------------------------------------------------
105  !
106  !++ Private parameters & variables
107  !
108  logical, private :: atmos_vars_checkrange = .false.
109  real(RP), private :: atmos_vars_checkcfl = 0.0_rp
110 
111  integer, private, parameter :: vmax = 5
112 
113  character(len=H_SHORT), private :: var_name(vmax)
114  character(len=H_MID), private :: var_desc(vmax)
115  character(len=H_SHORT), private :: var_unit(vmax)
116 
117  data var_name / 'DENS', &
118  'MOMZ', &
119  'MOMX', &
120  'MOMY', &
121  'RHOT' /
122  data var_desc / 'density', &
123  'momentum z', &
124  'momentum x', &
125  'momentum y', &
126  'rho * theta' /
127  data var_unit / 'kg/m3', &
128  'kg/m2/s', &
129  'kg/m2/s', &
130  'kg/m2/s', &
131  'kg/m3*K' /
132 
133  ! history output of prognostic variables
134  integer, private :: var_hist_id(vmax)
135 
136  integer, private, allocatable :: aq_hist_id(:)
137 
138  ! history & monitor output of diagnostic variables
139  integer, private, parameter :: ad_nmax = 70 ! number of diagnostic variables for history output
140 
141  integer, private, parameter :: i_w = 1 ! velocity w at cell center
142  integer, private, parameter :: i_u = 2 ! velocity u at cell center
143  integer, private, parameter :: i_v = 3 ! velocity v at cell center
144  integer, private, parameter :: i_pott = 4 ! potential temperature
145 
146  integer, private, parameter :: i_qdry = 5 ! ratio of dry air to total mass
147  integer, private, parameter :: i_qtot = 6 ! ratio of total tracer to total mass
148  integer, private, parameter :: i_qhyd = 7 ! ratio of total hydrometeor to total mass
149  integer, private, parameter :: i_qliq = 8 ! ratio of total liquid water to total mass
150  integer, private, parameter :: i_qice = 9 ! ratio of total ice water to total mass
151 
152  integer, private, parameter :: i_lwp = 10 ! liquid water path
153  integer, private, parameter :: i_iwp = 11 ! ice water path
154  integer, private, parameter :: i_pw = 12 ! ice water path
155 
156  integer, private, parameter :: i_rtot = 13 ! total gas constant
157  integer, private, parameter :: i_cptot = 14 ! total heat capacity (constant pressure)
158  integer, private, parameter :: i_pres = 15 ! pressure
159  integer, private, parameter :: i_temp = 16 ! temperature
160 
161  integer, private, parameter :: i_potl = 17 ! liquid water potential temperature
162  integer, private, parameter :: i_rha = 18 ! relative humidity (liquid+ice)
163  integer, private, parameter :: i_rhl = 19 ! relative humidity against to liquid
164  integer, private, parameter :: i_rhi = 20 ! relative humidity against to ice
165 
166  integer, private, parameter :: i_vor = 21 ! vertical vorticity
167  integer, private, parameter :: i_div = 22 ! divergence
168  integer, private, parameter :: i_hdiv = 23 ! horizontal divergence
169 
170  integer, private, parameter :: i_dens_prim = 24 ! prime term of density
171  integer, private, parameter :: i_w_prim = 25 ! prime term of w
172  integer, private, parameter :: i_u_prim = 26 ! prime term of u
173  integer, private, parameter :: i_v_prim = 27 ! prime term of v
174  integer, private, parameter :: i_pott_prim = 28 ! prime term of potential temperature
175  integer, private, parameter :: i_w_prim2 = 29 ! variance of w
176  integer, private, parameter :: i_pt_w_prim = 30 ! resolved scale heat flux
177  integer, private, parameter :: i_w_prim3 = 31 ! skewness of w
178  integer, private, parameter :: i_tke_rs = 32 ! resolved scale TKE
179 
180  integer, private, parameter :: i_engp = 33 ! potential energy
181  integer, private, parameter :: i_engk = 34 ! kinetic energy
182  integer, private, parameter :: i_engi = 35 ! internal energy
183  integer, private, parameter :: i_engt = 36 ! total energy
184 
185  integer, private, parameter :: i_engsfc_sh = 37
186  integer, private, parameter :: i_engsfc_lh = 38
187  integer, private, parameter :: i_engsfc_rd = 39
188  integer, private, parameter :: i_engtoa_rd = 40
189 
190  integer, private, parameter :: i_engsfc_lw_up = 41
191  integer, private, parameter :: i_engsfc_lw_dn = 42
192  integer, private, parameter :: i_engsfc_sw_up = 43
193  integer, private, parameter :: i_engsfc_sw_dn = 44
194 
195  integer, private, parameter :: i_engtoa_lw_up = 45
196  integer, private, parameter :: i_engtoa_lw_dn = 46
197  integer, private, parameter :: i_engtoa_sw_up = 47
198  integer, private, parameter :: i_engtoa_sw_dn = 48
199 
200  integer, private, parameter :: i_engflxt = 49
201 
202  integer, private, parameter :: i_evap = 50
203  integer, private, parameter :: i_prcp = 51
204 
205  integer, private, parameter :: i_dens_mean = 52
206  integer, private, parameter :: i_w_mean = 53
207  integer, private, parameter :: i_u_mean = 54
208  integer, private, parameter :: i_v_mean = 55
209  integer, private, parameter :: i_pott_mean = 56
210  integer, private, parameter :: i_t_mean = 57
211 
212  integer, private, parameter :: i_qv_mean = 58
213  integer, private, parameter :: i_qhyd_mean = 59
214  integer, private, parameter :: i_qliq_mean = 60
215  integer, private, parameter :: i_qice_mean = 61
216 
217  integer, private, parameter :: i_qsat = 62
218 
219  integer, private, parameter :: i_uabs = 63
220 
221  integer, private, parameter :: i_cape = 64
222  integer, private, parameter :: i_cin = 65
223  integer, private, parameter :: i_lcl = 66
224  integer, private, parameter :: i_lfc = 67
225  integer, private, parameter :: i_lnb = 68
226 
227  integer, private, parameter :: i_pblh = 69
228  integer, private, parameter :: i_mse = 70
229 
230  integer, private :: ad_hist_id (ad_nmax)
231  integer, private :: ad_prep_sw (ad_nmax)
232  integer, private :: ad_monit_id(ad_nmax)
233 
234  !-----------------------------------------------------------------------------
235 contains
236  !-----------------------------------------------------------------------------
238  subroutine atmos_vars_setup
239  use scale_process, only: &
241  use scale_history, only: &
242  hist_reg
243  use scale_monitor, only: &
244  monit_reg
245  use mod_atmos_admin, only: &
247  use mod_atmos_dyn_vars, only: &
249  use mod_atmos_phy_mp_vars, only: &
251  use mod_atmos_phy_ae_vars, only: &
253  use mod_atmos_phy_ch_vars, only: &
255  use mod_atmos_phy_rd_vars, only: &
257  use mod_atmos_phy_sf_vars, only: &
259  use mod_atmos_phy_tb_vars, only: &
261  use mod_atmos_phy_cp_vars, only: &
263  implicit none
264 
265  namelist / param_atmos_vars / &
275  atmos_vars_checkrange, &
276  atmos_vars_checkcfl
277 
278  logical :: zinterp ! dummy
279  integer :: ierr
280  integer :: iv, iq
281  !---------------------------------------------------------------------------
282 
283  if( io_l ) write(io_fid_log,*)
284  if( io_l ) write(io_fid_log,*) '++++++ Module[VARS] / Categ[ATMOS] / Origin[SCALE-RM]'
285 
286  allocate( dens(ka,ia,ja) )
287  allocate( momz(ka,ia,ja) )
288  allocate( momx(ka,ia,ja) )
289  allocate( momy(ka,ia,ja) )
290  allocate( rhot(ka,ia,ja) )
291  allocate( qtrc(ka,ia,ja,qa) )
292 
293  if ( atmos_use_average ) then
294  allocate( dens_avw(ka,ia,ja) )
295  allocate( momz_avw(ka,ia,ja) )
296  allocate( momx_avw(ka,ia,ja) )
297  allocate( momy_avw(ka,ia,ja) )
298  allocate( rhot_avw(ka,ia,ja) )
299  allocate( qtrc_avw(ka,ia,ja,qa) )
300 
301  dens_av => dens_avw
302  momz_av => momz_avw
303  momx_av => momx_avw
304  momy_av => momy_avw
305  rhot_av => rhot_avw
306  qtrc_av => qtrc_avw
307  else
308  dens_av => dens
309  momz_av => momz
310  momx_av => momx
311  momy_av => momy
312  rhot_av => rhot
313  qtrc_av => qtrc
314  endif
315 
316  allocate( dens_tp(ka,ia,ja) )
317  allocate( momz_tp(ka,ia,ja) )
318  allocate( momx_tp(ka,ia,ja) )
319  allocate( momy_tp(ka,ia,ja) )
320  allocate( rhot_tp(ka,ia,ja) )
321  allocate( rhoq_tp(ka,ia,ja,qa) )
322 
323  allocate( temp(ka,ia,ja) )
324  allocate( pres(ka,ia,ja) )
325  allocate( w(ka,ia,ja) )
326  allocate( u(ka,ia,ja) )
327  allocate( v(ka,ia,ja) )
328  allocate( pott(ka,ia,ja) )
329 
330  momz(1:ks-1,:,:) = 0.0_rp
331  momz(ke:ka,:,:) = 0.0_rp
332 
333  !--- read namelist
334  rewind(io_fid_conf)
335  read(io_fid_conf,nml=param_atmos_vars,iostat=ierr)
336  if( ierr < 0 ) then !--- missing
337  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
338  elseif( ierr > 0 ) then !--- fatal error
339  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_VARS. Check!'
340  call prc_mpistop
341  endif
342  if( io_lnml ) write(io_fid_log,nml=param_atmos_vars)
343 
344  if( io_l ) write(io_fid_log,*)
345  if( io_l ) write(io_fid_log,*) '*** List of prognostic variables (ATMOS) ***'
346  if( io_l ) write(io_fid_log,'(1x,A,A15,A,A32,3(A))') &
347  '*** |','VARNAME ','|', 'DESCRIPTION ','[', 'UNIT ',']'
348  do iv = 1, vmax
349  if( io_l ) write(io_fid_log,'(1x,A,i3,A,A15,A,A32,3(A))') &
350  '*** NO.',iv,'|',var_name(iv),'|', var_desc(iv),'[', var_unit(iv),']'
351  enddo
352  do iq = 1, qa
353  if( io_l ) write(io_fid_log,'(1x,A,i3,A,A15,A,A32,3(A))') &
354  '*** NO.',5+iq,'|',aq_name(iq),'|', aq_desc(iq),'[', aq_unit(iq),']'
355  enddo
356 
357  if( io_l ) write(io_fid_log,*)
358  if ( atmos_restart_in_basename /= '' ) then
359  if( io_l ) write(io_fid_log,*) '*** Restart input? : ', trim(atmos_restart_in_basename)
360  else
361  if( io_l ) write(io_fid_log,*) '*** Restart input? : NO'
362  endif
363  if ( atmos_restart_output &
364  .AND. atmos_restart_out_basename /= '' ) then
365  if( io_l ) write(io_fid_log,*) '*** Restart output? : ', trim(atmos_restart_out_basename)
366  else
367  if( io_l ) write(io_fid_log,*) '*** Restart output? : NO'
368  atmos_restart_output = .false.
369  endif
370 
371  if ( atmos_restart_check_basename == '' ) then
372  atmos_restart_check = .false.
373  endif
374 
375  if( io_l ) write(io_fid_log,*)
376  if( io_l ) write(io_fid_log,*) '*** Check restart consistency? : ', atmos_restart_check
377  if( io_l ) write(io_fid_log,*) '*** Check value range of variables? : ', atmos_vars_checkrange
378  if ( atmos_vars_checkcfl > 0.0_rp ) then
379  if( io_l ) write(io_fid_log,*) '*** Check CFL condition? : YES'
380  if( io_l ) write(io_fid_log,*) '*** Limit of Courant number : ', atmos_vars_checkcfl
381  else
382  if( io_l ) write(io_fid_log,*) '*** Check CFL condition? : NO'
383  endif
384 
393 
394 
395 
396  !##### todo: the part below should be moved to the mod_atmos_diag #####
397 
398  allocate( aq_hist_id(qa))
399 
400  var_hist_id(:) = -1
401  aq_hist_id(:) = -1
402  ad_hist_id(:) = -1
403  ad_monit_id(:) = -1
404  ad_prep_sw(:) = -1
405 
406  call hist_reg( var_hist_id(i_dens), zinterp, var_name(i_dens), var_desc(i_dens), var_unit(i_dens), ndim=3 )
407  call hist_reg( var_hist_id(i_momz), zinterp, var_name(i_momz), var_desc(i_momz), var_unit(i_momz), ndim=3, zdim='half' )
408  call hist_reg( var_hist_id(i_momx), zinterp, var_name(i_momx), var_desc(i_momx), var_unit(i_momx), ndim=3, xdim='half' )
409  call hist_reg( var_hist_id(i_momy), zinterp, var_name(i_momy), var_desc(i_momy), var_unit(i_momy), ndim=3, ydim='half' )
410  call hist_reg( var_hist_id(i_rhot), zinterp, var_name(i_rhot), var_desc(i_rhot), var_unit(i_rhot), ndim=3 )
411  do iq = 1, qa
412  call hist_reg( aq_hist_id(iq), zinterp, aq_name(iq), aq_desc(iq), aq_unit(iq), ndim=3 )
413  enddo
414 
415  call hist_reg( ad_hist_id(i_w ) , zinterp, 'W', 'velocity w', 'm/s', ndim=3 )
416  call hist_reg( ad_hist_id(i_u ) , zinterp, 'U', 'velocity u', 'm/s', ndim=3 )
417  call hist_reg( ad_hist_id(i_v ) , zinterp, 'V', 'velocity v', 'm/s', ndim=3 )
418  call hist_reg( ad_hist_id(i_pott) , zinterp, 'PT', 'potential temp.', 'K', ndim=3 )
419 
420  call hist_reg( ad_hist_id(i_qdry) , zinterp, 'QDRY', 'dry air', 'kg/kg', ndim=3 )
421  call hist_reg( ad_hist_id(i_qtot) , zinterp, 'QTOT', 'total water', 'kg/kg', ndim=3 )
422  call hist_reg( ad_hist_id(i_qhyd) , zinterp, 'QHYD', 'total hydrometeors', 'kg/kg', ndim=3 )
423  call hist_reg( ad_hist_id(i_qliq) , zinterp, 'QLIQ', 'total liquid water', 'kg/kg', ndim=3 )
424  call hist_reg( ad_hist_id(i_qice) , zinterp, 'QICE', 'total ice water', 'kg/kg', ndim=3 )
425 
426  call hist_reg( ad_hist_id(i_lwp) , zinterp, 'LWP', 'liquid water path', 'g/m2', ndim=2 )
427  call hist_reg( ad_hist_id(i_iwp) , zinterp, 'IWP', 'ice water path', 'g/m2', ndim=2 )
428  call hist_reg( ad_hist_id(i_pw ) , zinterp, 'PW', 'precipitable water', 'g/m2', ndim=2 )
429 
430  call hist_reg( ad_hist_id(i_rtot) , zinterp, 'RTOT', 'Total gas constant', 'J/kg/K', ndim=3 )
431  call hist_reg( ad_hist_id(i_cptot), zinterp, 'CPTOT', 'Total heat capacity', 'J/kg/K', ndim=3 )
432  call hist_reg( ad_hist_id(i_pres) , zinterp, 'PRES', 'pressure', 'Pa', ndim=3 )
433  call hist_reg( ad_hist_id(i_temp) , zinterp, 'T', 'temperature', 'K', ndim=3 )
434 
435  call hist_reg( ad_hist_id(i_potl) , zinterp, 'LWPT', 'liq. potential temp.', 'K', ndim=3 )
436  call hist_reg( ad_hist_id(i_rha) , zinterp, 'RHA', 'relative humidity(liq+ice)', '%', ndim=3 )
437  call hist_reg( ad_hist_id(i_rhl) , zinterp, 'RH', 'relative humidity(liq)', '%', ndim=3 )
438  call hist_reg( ad_hist_id(i_rhi) , zinterp, 'RHI', 'relative humidity(ice)', '%', ndim=3 )
439 
440  call hist_reg( ad_hist_id(i_vor) , zinterp, 'VOR', 'vertical vorticity', '1/s', ndim=3 )
441  call hist_reg( ad_hist_id(i_div) , zinterp, 'DIV', 'divergence', '1/s', ndim=3 )
442  call hist_reg( ad_hist_id(i_hdiv) , zinterp, 'HDIV', 'horizontal divergence', '1/s', ndim=3 )
443  call hist_reg( ad_hist_id(i_uabs) , zinterp, 'Uabs', 'absolute velocity', 'm/s', ndim=3 )
444 
445  call hist_reg( ad_hist_id(i_cape) , zinterp, 'CAPE', 'convection avail. pot. energy', 'm2/s2', ndim=2 )
446  call hist_reg( ad_hist_id(i_cin) , zinterp, 'CIN', 'convection inhibition', 'm2/s2', ndim=2 )
447  call hist_reg( ad_hist_id(i_lcl) , zinterp, 'LCL', 'lifted condensation level', 'm', ndim=2 )
448  call hist_reg( ad_hist_id(i_lfc) , zinterp, 'LFC', 'level of free convection', 'm', ndim=2 )
449  call hist_reg( ad_hist_id(i_lnb) , zinterp, 'LNB', 'level of neutral buoyancy', 'm', ndim=2 )
450 
451  call hist_reg( ad_hist_id(i_pblh) , zinterp, 'PBLH', 'PBL height', 'm', ndim=2 )
452  call hist_reg( ad_hist_id(i_mse) , zinterp, 'MSE', 'moist static energy', 'm2/s2', ndim=3 )
453 
454  call hist_reg( ad_hist_id(i_dens_mean), zinterp, 'DENS_MEAN', 'horiz. mean of density', 'kg/m3', ndim=1 )
455  call hist_reg( ad_hist_id(i_w_mean), zinterp, 'W_MEAN', 'horiz. mean of w', 'm/s', ndim=1 )
456  call hist_reg( ad_hist_id(i_u_mean), zinterp, 'U_MEAN', 'horiz. mean of u', 'm/s', ndim=1 )
457  call hist_reg( ad_hist_id(i_v_mean), zinterp, 'V_MEAN', 'horiz. mean of v', 'm/s', ndim=1 )
458  call hist_reg( ad_hist_id(i_pott_mean), zinterp, 'PT_MEAN', 'horiz. mean of pot.', 'K', ndim=1 )
459  call hist_reg( ad_hist_id(i_t_mean), zinterp, 'T_MEAN', 'horiz. mean of t', 'K', ndim=1 )
460  call hist_reg( ad_hist_id(i_qv_mean), zinterp, 'QV_MEAN', 'horiz. mean of QV', '1', ndim=1 )
461  call hist_reg( ad_hist_id(i_qhyd_mean), zinterp, 'QHYD_MEAN', 'horiz. mean of QHYD', '1', ndim=1 )
462  call hist_reg( ad_hist_id(i_qliq_mean), zinterp, 'QLIQ_MEAN', 'horiz. mean of QLIQ', '1', ndim=1 )
463  call hist_reg( ad_hist_id(i_qice_mean), zinterp, 'QICE_MEAN', 'horiz. mean of QICE', '1', ndim=1 )
464 
465  call hist_reg( ad_hist_id(i_dens_prim), zinterp, 'DENS_PRIM', 'horiz. deviation of density', 'kg/m3', ndim=3 )
466  call hist_reg( ad_hist_id(i_w_prim ), zinterp, 'W_PRIM', 'horiz. deviation of w', 'm/s', ndim=3 )
467  call hist_reg( ad_hist_id(i_u_prim ), zinterp, 'U_PRIM', 'horiz. deviation of u', 'm/s', ndim=3 )
468  call hist_reg( ad_hist_id(i_v_prim ), zinterp, 'V_PRIM', 'horiz. deviation of v', 'm/s', ndim=3 )
469  call hist_reg( ad_hist_id(i_pott_prim), zinterp, 'PT_PRIM', 'horiz. deviation of pot. temp.', 'K', ndim=3 )
470  call hist_reg( ad_hist_id(i_w_prim2 ), zinterp, 'W_PRIM2', 'variance of w', 'm2/s2', ndim=3 )
471  call hist_reg( ad_hist_id(i_pt_w_prim), zinterp, 'PT_W_PRIM', 'resolved scale heat flux', 'W/s', ndim=3 )
472  call hist_reg( ad_hist_id(i_w_prim3 ), zinterp, 'W_PRIM3', 'skewness of w', 'm3/s3', ndim=3 )
473  call hist_reg( ad_hist_id(i_tke_rs ), zinterp, 'TKE_RS', 'resolved scale TKE', 'm2/s2', ndim=3 )
474 
475  call hist_reg( ad_hist_id(i_engt) , zinterp, 'ENGT', 'total energy', 'J/m3', ndim=3 )
476  call hist_reg( ad_hist_id(i_engp) , zinterp, 'ENGP', 'potential energy', 'J/m3', ndim=3 )
477  call hist_reg( ad_hist_id(i_engk) , zinterp, 'ENGK', 'kinetic energy', 'J/m3', ndim=3 )
478  call hist_reg( ad_hist_id(i_engi) , zinterp, 'ENGI', 'internal energy', 'J/m3', ndim=3 )
479 
480  !-----< monitor output setup >-----
481 
482  call monit_reg( ad_monit_id(i_qdry), 'QDRY', 'dry air mass', 'kg', ndim=3, isflux=.false. )
483  call monit_reg( ad_monit_id(i_qtot), 'QTOT', 'water mass', 'kg', ndim=3, isflux=.false. )
484  call monit_reg( ad_monit_id(i_evap), 'EVAP', 'evaporation', 'kg', ndim=2, isflux=.true. )
485  call monit_reg( ad_monit_id(i_prcp), 'PRCP', 'precipitation', 'kg', ndim=2, isflux=.true. )
486 
487  call monit_reg( ad_monit_id(i_engt), 'ENGT', 'total energy', 'J', ndim=3, isflux=.false. )
488  call monit_reg( ad_monit_id(i_engp), 'ENGP', 'potential energy', 'J', ndim=3, isflux=.false. )
489  call monit_reg( ad_monit_id(i_engk), 'ENGK', 'kinetic energy', 'J', ndim=3, isflux=.false. )
490  call monit_reg( ad_monit_id(i_engi), 'ENGI', 'internal energy', 'J', ndim=3, isflux=.false. )
491 
492  call monit_reg( ad_monit_id(i_engflxt), 'ENGFLXT', 'total energy flux', 'J', ndim=2, isflux=.true. )
493 
494  call monit_reg( ad_monit_id(i_engsfc_sh), 'ENGSFC_SH', 'total energy', 'J', ndim=2, isflux=.true. )
495  call monit_reg( ad_monit_id(i_engsfc_lh), 'ENGSFC_LH', 'potential energy', 'J', ndim=2, isflux=.true. )
496  call monit_reg( ad_monit_id(i_engsfc_rd), 'ENGSFC_RD', 'kinetic energy', 'J', ndim=2, isflux=.true. )
497  call monit_reg( ad_monit_id(i_engtoa_rd), 'ENGTOA_RD', 'internal energy', 'J', ndim=2, isflux=.true. )
498 
499  call monit_reg( ad_monit_id(i_engsfc_lw_up), 'ENGSFC_LW_up', 'total energy', 'J', ndim=2, isflux=.true. )
500  call monit_reg( ad_monit_id(i_engsfc_lw_dn), 'ENGSFC_LW_dn', 'potential energy', 'J', ndim=2, isflux=.true. )
501  call monit_reg( ad_monit_id(i_engsfc_sw_up), 'ENGSFC_SW_up', 'kinetic energy', 'J', ndim=2, isflux=.true. )
502  call monit_reg( ad_monit_id(i_engsfc_sw_dn), 'ENGSFC_SW_dn', 'internal energy', 'J', ndim=2, isflux=.true. )
503 
504  call monit_reg( ad_monit_id(i_engtoa_lw_up), 'ENGTOA_LW_up', 'total energy', 'J', ndim=2, isflux=.true. )
505  call monit_reg( ad_monit_id(i_engtoa_lw_dn), 'ENGTOA_LW_dn', 'potential energy', 'J', ndim=2, isflux=.true. )
506  call monit_reg( ad_monit_id(i_engtoa_sw_up), 'ENGTOA_SW_up', 'kinetic energy', 'J', ndim=2, isflux=.true. )
507  call monit_reg( ad_monit_id(i_engtoa_sw_dn), 'ENGTOA_SW_dn', 'internal energy', 'J', ndim=2, isflux=.true. )
508 
509  if ( ad_hist_id(i_w) > 0 ) then
510  ad_prep_sw(i_w) = 1
511  endif
512  if ( ad_hist_id(i_u) > 0 ) then
513  ad_prep_sw(i_u) = 1
514  endif
515  if ( ad_hist_id(i_v) > 0 ) then
516  ad_prep_sw(i_v) = 1
517  endif
518  if ( ad_hist_id(i_pott) > 0 ) then
519  ad_prep_sw(i_pott) = 1
520  endif
521 
522  if ( ad_hist_id(i_qdry) > 0 &
523  .OR. ad_monit_id(i_qdry) > 0 ) then
524  ad_prep_sw(i_qdry) = 1
525  endif
526  if ( ad_hist_id(i_qtot) > 0 &
527  .OR. ad_monit_id(i_qtot) > 0 ) then
528  ad_prep_sw(i_qdry) = 1
529  ad_prep_sw(i_qtot) = 1
530  endif
531  if ( ad_hist_id(i_qhyd) > 0 ) then
532  ad_prep_sw(i_qhyd) = 1
533  endif
534  if ( ad_hist_id(i_qliq) > 0 ) then
535  ad_prep_sw(i_qliq) = 1
536  endif
537  if ( ad_hist_id(i_qice) > 0 ) then
538  ad_prep_sw(i_qice) = 1
539  endif
540 
541  if ( ad_hist_id(i_lwp) > 0 ) then
542  ad_prep_sw(i_qliq) = 1
543  ad_prep_sw(i_lwp) = 1
544  endif
545  if ( ad_hist_id(i_iwp) > 0 ) then
546  ad_prep_sw(i_qice) = 1
547  ad_prep_sw(i_iwp) = 1
548  endif
549  if ( ad_hist_id(i_pw) > 0 ) then
550  ad_prep_sw(i_pw) = 1
551  endif
552 
553  if ( ad_hist_id(i_rtot) > 0 ) then
554  ad_prep_sw(i_qdry) = 1
555  ad_prep_sw(i_rtot) = 1
556  endif
557  if ( ad_hist_id(i_cptot) > 0 ) then
558  ad_prep_sw(i_qdry) = 1
559  ad_prep_sw(i_cptot) = 1
560  endif
561  if ( ad_hist_id(i_pres) > 0 ) then
562  ad_prep_sw(i_qdry) = 1
563  ad_prep_sw(i_rtot) = 1
564  ad_prep_sw(i_cptot) = 1
565  ad_prep_sw(i_pres) = 1
566  endif
567  if ( ad_hist_id(i_temp) > 0 ) then
568  ad_prep_sw(i_qdry) = 1
569  ad_prep_sw(i_rtot) = 1
570  ad_prep_sw(i_cptot) = 1
571  ad_prep_sw(i_pres) = 1
572  ad_prep_sw(i_temp) = 1
573  endif
574 
575  if ( ad_hist_id(i_potl) > 0 ) then
576  ad_prep_sw(i_pott) = 1
577  ad_prep_sw(i_qdry) = 1
578  ad_prep_sw(i_rtot) = 1
579  ad_prep_sw(i_cptot) = 1
580  ad_prep_sw(i_pres) = 1
581  ad_prep_sw(i_temp) = 1
582  ad_prep_sw(i_potl) = 1
583  endif
584  if ( ad_hist_id(i_rha) > 0 &
585  .OR. ad_hist_id(i_rhl) > 0 &
586  .OR. ad_hist_id(i_rhi) > 0 ) then
587  ad_prep_sw(i_qdry) = 1
588  ad_prep_sw(i_rtot) = 1
589  ad_prep_sw(i_cptot) = 1
590  ad_prep_sw(i_pres) = 1
591  ad_prep_sw(i_temp) = 1
592  ad_prep_sw(i_qsat) = 1
593  endif
594 
595 
596  if ( ad_hist_id(i_vor) > 0 ) then
597  ad_prep_sw(i_vor) = 1
598  endif
599 
600  if ( ad_prep_sw(i_div) > 0 ) then
601  ad_prep_sw(i_hdiv) = 1
602  endif
603 
604  if ( ad_hist_id(i_uabs) > 0 ) then
605  ad_prep_sw(i_u) = 1
606  ad_prep_sw(i_v) = 1
607  ad_prep_sw(i_uabs) = 1
608  endif
609 
610  if ( ad_hist_id(i_cape) > 0 &
611  .OR. ad_hist_id(i_cin) > 0 &
612  .OR. ad_hist_id(i_lcl) > 0 &
613  .OR. ad_hist_id(i_lfc) > 0 &
614  .OR. ad_hist_id(i_lnb) > 0 ) then
615  ad_prep_sw(i_cape) = 1
616  ad_prep_sw(i_cin) = 1
617  ad_prep_sw(i_lcl) = 1
618  ad_prep_sw(i_lfc) = 1
619  ad_prep_sw(i_lnb) = 1
620  endif
621 
622  if ( ad_hist_id(i_pblh) > 0 ) then
623  ad_prep_sw(i_pott) = 1
624  ad_prep_sw(i_pblh) = 1
625  endif
626 
627  if ( ad_hist_id(i_mse) > 0 ) then
628  ad_prep_sw(i_cptot) = 1
629  ad_prep_sw(i_temp) = 1
630  ad_prep_sw(i_mse) = 1
631  endif
632 
633  if ( ad_hist_id(i_dens_prim) > 0 ) then
634  ad_prep_sw(i_dens_prim) = 1
635  ad_prep_sw(i_dens_mean) = 1
636  endif
637 
638  if ( ad_hist_id(i_w_prim) > 0 ) then
639  ad_prep_sw(i_w) = 1
640  ad_prep_sw(i_w_prim) = 1
641  ad_prep_sw(i_dens_mean) = 1
642  ad_prep_sw(i_w_mean) = 1
643  endif
644 
645  if ( ad_hist_id(i_u_prim) > 0 ) then
646  ad_prep_sw(i_u) = 1
647  ad_prep_sw(i_u_prim) = 1
648  ad_prep_sw(i_dens_mean) = 1
649  ad_prep_sw(i_u_mean) = 1
650  endif
651 
652  if ( ad_hist_id(i_v_prim) > 0 ) then
653  ad_prep_sw(i_v) = 1
654  ad_prep_sw(i_v_prim) = 1
655  ad_prep_sw(i_dens_mean) = 1
656  ad_prep_sw(i_v_mean) = 1
657  endif
658 
659  if ( ad_hist_id(i_pott_prim) > 0 ) then
660  ad_prep_sw(i_pott) = 1
661  ad_prep_sw(i_pott_prim) = 1
662  ad_prep_sw(i_dens_mean) = 1
663  ad_prep_sw(i_pott_mean) = 1
664  endif
665 
666  if ( ad_hist_id(i_dens_mean) > 0 ) then
667  ad_prep_sw(i_dens_mean) = 1
668  endif
669 
670  if ( ad_hist_id(i_w_mean) > 0 ) then
671  ad_prep_sw(i_w_mean) = 1
672  endif
673 
674  if ( ad_hist_id(i_u_mean) > 0 ) then
675  ad_prep_sw(i_u_mean) = 1
676  endif
677 
678  if ( ad_hist_id(i_v_mean) > 0 ) then
679  ad_prep_sw(i_v_mean) = 1
680  endif
681 
682  if ( ad_hist_id(i_pott_mean) > 0 ) then
683  ad_prep_sw(i_pott_mean) = 1
684  end if
685 
686  if ( ad_hist_id(i_t_mean) > 0 ) then
687  ad_prep_sw(i_t_mean) = 1
688  end if
689 
690  if ( ad_hist_id(i_qv_mean) > 0 ) then
691  ad_prep_sw(i_qv_mean) = 1
692  end if
693 
694  if ( ad_hist_id(i_qhyd_mean) > 0 ) then
695  ad_prep_sw(i_qhyd) = 1
696  ad_prep_sw(i_qhyd_mean) = 1
697  end if
698 
699  if ( ad_hist_id(i_qliq_mean) > 0 ) then
700  ad_prep_sw(i_qliq) = 1
701  ad_prep_sw(i_qliq_mean) = 1
702  end if
703 
704  if ( ad_hist_id(i_qice_mean) > 0 ) then
705  ad_prep_sw(i_qice) = 1
706  ad_prep_sw(i_qice_mean) = 1
707  end if
708 
709  if ( ad_hist_id(i_w_prim2) > 0 ) then
710  ad_prep_sw(i_w) = 1
711  ad_prep_sw(i_w_prim) = 1
712  ad_prep_sw(i_w_prim2) = 1
713  ad_prep_sw(i_dens_mean) = 1
714  ad_prep_sw(i_w_mean) = 1
715  endif
716 
717  if ( ad_hist_id(i_pt_w_prim) > 0 ) then
718  ad_prep_sw(i_w) = 1
719  ad_prep_sw(i_w_prim) = 1
720  ad_prep_sw(i_pott) = 1
721  ad_prep_sw(i_pott_prim) = 1
722  ad_prep_sw(i_pt_w_prim) = 1
723  ad_prep_sw(i_dens_mean) = 1
724  ad_prep_sw(i_w_mean) = 1
725  ad_prep_sw(i_pott_mean) = 1
726  endif
727 
728  if ( ad_hist_id(i_w_prim3) > 0 ) then
729  ad_prep_sw(i_w) = 1
730  ad_prep_sw(i_w_prim) = 1
731  ad_prep_sw(i_w_prim3) = 1
732  ad_prep_sw(i_dens_mean) = 1
733  ad_prep_sw(i_w_mean) = 1
734  endif
735 
736  if ( ad_hist_id(i_tke_rs) > 0 ) then
737  ad_prep_sw(i_w) = 1
738  ad_prep_sw(i_u) = 1
739  ad_prep_sw(i_v) = 1
740  ad_prep_sw(i_w_prim) = 1
741  ad_prep_sw(i_u_prim) = 1
742  ad_prep_sw(i_v_prim) = 1
743  ad_prep_sw(i_tke_rs) = 1
744  ad_prep_sw(i_dens_mean) = 1
745  ad_prep_sw(i_w_mean) = 1
746  ad_prep_sw(i_u_mean) = 1
747  ad_prep_sw(i_v_mean) = 1
748  endif
749 
750  if ( ad_hist_id(i_engp) > 0 &
751  .OR. ad_monit_id(i_engp) > 0 ) then
752  ad_prep_sw(i_engp) = 1
753  endif
754  if ( ad_hist_id(i_engk) > 0 &
755  .OR. ad_monit_id(i_engk) > 0 ) then
756  ad_prep_sw(i_w) = 1
757  ad_prep_sw(i_u) = 1
758  ad_prep_sw(i_v) = 1
759  ad_prep_sw(i_engk) = 1
760  endif
761  if ( ad_hist_id(i_engi) > 0 &
762  .OR. ad_monit_id(i_engi) > 0 ) then
763  ad_prep_sw(i_qdry) = 1
764  ad_prep_sw(i_rtot) = 1
765  ad_prep_sw(i_cptot) = 1
766  ad_prep_sw(i_pres) = 1
767  ad_prep_sw(i_temp) = 1
768  ad_prep_sw(i_engi) = 1
769  endif
770  if ( ad_hist_id(i_engt) > 0 &
771  .OR. ad_monit_id(i_engt) > 0 ) then
772  ad_prep_sw(i_engp) = 1
773  ad_prep_sw(i_w) = 1
774  ad_prep_sw(i_u) = 1
775  ad_prep_sw(i_v) = 1
776  ad_prep_sw(i_engk) = 1
777  ad_prep_sw(i_qdry) = 1
778  ad_prep_sw(i_rtot) = 1
779  ad_prep_sw(i_cptot) = 1
780  ad_prep_sw(i_pres) = 1
781  ad_prep_sw(i_temp) = 1
782  ad_prep_sw(i_engi) = 1
783  ad_prep_sw(i_engt) = 1
784  endif
785 
786  return
787  end subroutine atmos_vars_setup
788 
789  !-----------------------------------------------------------------------------
791  subroutine atmos_vars_fillhalo( &
792  FILL_BND )
793  use scale_comm, only: &
794  comm_vars8, &
795  comm_wait
796  implicit none
797 
798  logical, intent(in), optional :: FILL_BND
799 
800  logical :: FILL_BND_
801  integer :: i, j, iq
802  !---------------------------------------------------------------------------
803 
804  fill_bnd_ = .false.
805  if ( present(fill_bnd) ) fill_bnd_ = fill_bnd
806 
807  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
808  do j = jsb, jeb
809  do i = isb, ieb
810  dens( 1:ks-1,i,j) = dens(ks,i,j)
811  momz( 1:ks-1,i,j) = momz(ks,i,j)
812  momx( 1:ks-1,i,j) = momx(ks,i,j)
813  momy( 1:ks-1,i,j) = momy(ks,i,j)
814  rhot( 1:ks-1,i,j) = rhot(ks,i,j)
815  dens(ke+1:ka, i,j) = dens(ke,i,j)
816  momz(ke+1:ka, i,j) = momz(ke,i,j)
817  momx(ke+1:ka, i,j) = momx(ke,i,j)
818  momy(ke+1:ka, i,j) = momy(ke,i,j)
819  rhot(ke+1:ka, i,j) = rhot(ke,i,j)
820  enddo
821  enddo
822 
823  !$omp parallel do private(i,j,iq) OMP_SCHEDULE_ collapse(3)
824  do iq = 1, qa
825  do j = jsb, jeb
826  do i = isb, ieb
827  qtrc( 1:ks-1,i,j,iq) = qtrc(ks,i,j,iq)
828  qtrc(ke+1:ka, i,j,iq) = qtrc(ke,i,j,iq)
829  enddo
830  enddo
831  enddo
832 
833  call comm_vars8( dens(:,:,:), 1 )
834  call comm_vars8( momz(:,:,:), 2 )
835  call comm_vars8( momx(:,:,:), 3 )
836  call comm_vars8( momy(:,:,:), 4 )
837  call comm_vars8( rhot(:,:,:), 5 )
838  call comm_wait ( dens(:,:,:), 1, fill_bnd_ )
839  call comm_wait ( momz(:,:,:), 2, fill_bnd_ )
840  call comm_wait ( momx(:,:,:), 3, fill_bnd_ )
841  call comm_wait ( momy(:,:,:), 4, fill_bnd_ )
842  call comm_wait ( rhot(:,:,:), 5, fill_bnd_ )
843 
844  do iq = 1, qa
845  call comm_vars8( qtrc(:,:,:,iq), iq )
846  enddo
847  do iq = 1, qa
848  call comm_wait ( qtrc(:,:,:,iq), iq, fill_bnd_ )
849  enddo
850 
851  return
852  end subroutine atmos_vars_fillhalo
853 
854  !-----------------------------------------------------------------------------
856  subroutine atmos_vars_restart_read
857  use scale_process, only: &
859  use scale_const, only: &
860  grav => const_grav
861  use scale_fileio, only: &
862  fileio_read
863  use scale_atmos_thermodyn, only: &
864  thermodyn_qd => atmos_thermodyn_qd, &
865  thermodyn_temp_pres => atmos_thermodyn_temp_pres
866  use mod_atmos_admin, only: &
868  atmos_sw_dyn, &
869  atmos_sw_phy_mp, &
870  atmos_sw_phy_ae, &
871  atmos_sw_phy_ch, &
872  atmos_sw_phy_rd, &
873  atmos_sw_phy_sf, &
874  atmos_sw_phy_tb, &
876  use mod_atmos_dyn_vars, only: &
878  use mod_atmos_phy_mp_vars, only: &
880  use mod_atmos_phy_ae_vars, only: &
882  use mod_atmos_phy_ch_vars, only: &
884  use mod_atmos_phy_rd_vars, only: &
886  use mod_atmos_phy_sf_vars, only: &
888  use mod_atmos_phy_tb_vars, only: &
890  use mod_atmos_phy_cp_vars, only: &
892  implicit none
893 
894  integer :: iq
895  !---------------------------------------------------------------------------
896 
897  if( io_l ) write(io_fid_log,*)
898  if( io_l ) write(io_fid_log,*) '*** Input restart file (ATMOS) ***'
899 
900  if ( atmos_restart_in_basename /= '' ) then
901  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(atmos_restart_in_basename)
902 
903  call fileio_read( dens(:,:,:), & ! [OUT]
904  atmos_restart_in_basename, var_name(1), 'ZXY', step=1 ) ! [IN]
905  call fileio_read( momz(:,:,:), & ! [OUT]
906  atmos_restart_in_basename, var_name(2), 'ZXY', step=1 ) ! [IN]
907  call fileio_read( momx(:,:,:), & ! [OUT]
908  atmos_restart_in_basename, var_name(3), 'ZXY', step=1 ) ! [IN]
909  call fileio_read( momy(:,:,:), & ! [OUT]
910  atmos_restart_in_basename, var_name(4), 'ZXY', step=1 ) ! [IN]
911  call fileio_read( rhot(:,:,:), & ! [OUT]
912  atmos_restart_in_basename, var_name(5), 'ZXY', step=1 ) ! [IN]
913 
914  do iq = 1, qa
915  call fileio_read( qtrc(:,:,:,iq), & ! [OUT]
916  atmos_restart_in_basename, aq_name(iq), 'ZXY', step=1 ) ! [IN]
917  enddo
918 
920 
921  call atmos_vars_total
922  else
923  write(*,*) '*** restart file for atmosphere is not specified. STOP!'
924  call prc_mpistop
925  endif
926 
927  if ( atmos_use_average ) then
928  dens_av(:,:,:) = dens(:,:,:)
929  momz_av(:,:,:) = momz(:,:,:)
930  momx_av(:,:,:) = momx(:,:,:)
931  momy_av(:,:,:) = momy(:,:,:)
932  rhot_av(:,:,:) = rhot(:,:,:)
933  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
934  endif
935 
944 
945  return
946  end subroutine atmos_vars_restart_read
947 
948  !-----------------------------------------------------------------------------
950  subroutine atmos_vars_restart_write
951  use scale_time, only: &
953  use scale_fileio, only: &
954  fileio_write
955  use mod_atmos_admin, only: &
956  atmos_sw_dyn, &
957  atmos_sw_phy_mp, &
958  atmos_sw_phy_ae, &
959  atmos_sw_phy_ch, &
960  atmos_sw_phy_rd, &
961  atmos_sw_phy_sf, &
962  atmos_sw_phy_tb, &
964  use mod_atmos_dyn_vars, only: &
966  use mod_atmos_phy_mp_vars, only: &
968  use mod_atmos_phy_ae_vars, only: &
970  use mod_atmos_phy_ch_vars, only: &
972  use mod_atmos_phy_rd_vars, only: &
974  use mod_atmos_phy_sf_vars, only: &
976  use mod_atmos_phy_tb_vars, only: &
978  use mod_atmos_phy_cp_vars, only: &
980 #ifdef _SDM
981  use scale_atmos_phy_mp_sdm, only: &
982  sd_rest_flg_out, &
983  atmos_phy_mp_sdm_restart_out
984  use scale_time, only: &
985  nowsec => time_nowsec
986 #endif
987  implicit none
988 
989  character(len=20) :: timelabel
990  character(len=H_LONG) :: basename
991 
992  integer :: iq
993  !---------------------------------------------------------------------------
994 
995 #ifdef _SDM
996  if( sd_rest_flg_out ) then
997  if( io_l ) write(io_fid_log,*) '*** Output random number for SDM ***'
998  call atmos_phy_mp_sdm_restart_out(nowsec)
999  endif
1000 #endif
1001 
1002  if ( atmos_restart_out_basename /= '' ) then
1003 
1004  call time_gettimelabel( timelabel )
1005  write(basename,'(A,A,A)') trim(atmos_restart_out_basename), '_', trim(timelabel)
1006 
1007  if( io_l ) write(io_fid_log,*)
1008  if( io_l ) write(io_fid_log,*) '*** Output restart file (ATMOS) ***'
1009  if( io_l ) write(io_fid_log,*) '*** basename: ', trim(basename)
1010 
1011  call atmos_vars_fillhalo
1012 
1013  call atmos_vars_total
1014 
1015  call fileio_write( dens(:,:,:), basename, atmos_restart_out_title, & ! [IN]
1016  var_name(i_dens), var_desc(i_dens), var_unit(i_dens), 'ZXY', atmos_restart_out_dtype ) ! [IN]
1017  call fileio_write( momz(:,:,:), basename, atmos_restart_out_title, & ! [IN]
1018  var_name(i_momz), var_desc(i_momz), var_unit(i_momz), 'ZHXY', atmos_restart_out_dtype ) ! [IN]
1019  call fileio_write( momx(:,:,:), basename, atmos_restart_out_title, & ! [IN]
1020  var_name(i_momx), var_desc(i_momx), var_unit(i_momx), 'ZXHY', atmos_restart_out_dtype ) ! [IN]
1021  call fileio_write( momy(:,:,:), basename, atmos_restart_out_title, & ! [IN]
1022  var_name(i_momy), var_desc(i_momy), var_unit(i_momy), 'ZXYH', atmos_restart_out_dtype ) ! [IN]
1023  call fileio_write( rhot(:,:,:), basename, atmos_restart_out_title, & ! [IN]
1024  var_name(i_rhot), var_desc(i_rhot), var_unit(i_rhot), 'ZXY', atmos_restart_out_dtype ) ! [IN]
1025 
1026  do iq = 1, qa
1027  call fileio_write( qtrc(:,:,:,iq), basename, atmos_restart_out_title, & ! [IN]
1028  aq_name(iq), aq_desc(iq), aq_unit(iq), 'ZXY', atmos_restart_out_dtype ) ! [IN]
1029  enddo
1030 
1031  endif
1032 
1041 
1042  return
1043  end subroutine atmos_vars_restart_write
1044 
1045  !-----------------------------------------------------------------------------
1047  subroutine atmos_vars_restart_check
1048  use scale_process, only: &
1049  prc_myrank
1050  use gtool_file, only: &
1051  fileread
1052  implicit none
1053 
1054  real(RP) :: DENS_check(ka,ia,ja) ! Density [kg/m3]
1055  real(RP) :: MOMZ_check(ka,ia,ja) ! momentum z [kg/s/m2]
1056  real(RP) :: MOMX_check(ka,ia,ja) ! momentum x [kg/s/m2]
1057  real(RP) :: MOMY_check(ka,ia,ja) ! momentum y [kg/s/m2]
1058  real(RP) :: RHOT_check(ka,ia,ja) ! DENS * POTT [K*kg/m3]
1059  real(RP) :: QTRC_check(ka,ia,ja,qa) ! tracer mixing ratio [kg/kg]
1060 
1061  real(RP) :: restart_atmos(kmax,imax,jmax)
1062 
1063  character(len=H_LONG) :: basename
1064 
1065  logical :: datacheck
1066  integer :: k, i, j, iq
1067  !---------------------------------------------------------------------------
1068 
1069  call prof_rapstart('Debug')
1070 
1071  write(*,*) 'Compare last Data with ', trim(atmos_restart_check_basename), 'on PE=', prc_myrank
1072  write(*,*) '*** criterion = ', atmos_restart_check_criterion
1073  datacheck = .true.
1074 
1075  basename = atmos_restart_check_basename
1076 
1077  call fileread( restart_atmos(:,:,:), basename, 'DENS', 1, prc_myrank )
1078  dens_check(ks:ke,is:ie,js:je) = restart_atmos(1:kmax,1:imax,1:jmax)
1079  do k = ks, ke
1080  do j = js, je
1081  do i = is, ie
1082  if ( abs( dens(k,i,j)-dens_check(k,i,j) ) > atmos_restart_check_criterion ) then
1083  write(*,*) 'xxx there is the difference : ', dens(k,i,j)-dens_check(k,i,j)
1084  write(*,*) 'xxx at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'DENS'
1085  datacheck = .false.
1086  endif
1087  enddo
1088  enddo
1089  enddo
1090 
1091  call fileread( restart_atmos(:,:,:), basename, 'MOMZ', 1, prc_myrank )
1092  momz_check(ks:ke,is:ie,js:je) = restart_atmos(1:kmax,1:imax,1:jmax)
1093  do k = ks, ke
1094  do j = js, je
1095  do i = is, ie
1096  if ( abs( momz(k,i,j)-momz_check(k,i,j) ) > atmos_restart_check_criterion ) then
1097  write(*,*) 'xxx there is the difference : ', momz(k,i,j)-momz_check(k,i,j)
1098  write(*,*) 'xxx at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'MOMZ'
1099  datacheck = .false.
1100  endif
1101  enddo
1102  enddo
1103  enddo
1104 
1105  call fileread( restart_atmos(:,:,:), basename, 'MOMX', 1, prc_myrank )
1106  momx_check(ks:ke,is:ie,js:je) = restart_atmos(1:kmax,1:imax,1:jmax)
1107  do k = ks, ke
1108  do j = js, je
1109  do i = is, ie
1110  if ( abs( momx(k,i,j)-momx_check(k,i,j) ) > atmos_restart_check_criterion ) then
1111  write(*,*) 'xxx there is the difference : ', momx(k,i,j)-momx_check(k,i,j)
1112  write(*,*) 'xxx at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'MOMX'
1113  datacheck = .false.
1114  endif
1115  enddo
1116  enddo
1117  enddo
1118 
1119  call fileread( restart_atmos(:,:,:), basename, 'MOMY', 1, prc_myrank )
1120  momy_check(ks:ke,is:ie,js:je) = restart_atmos(1:kmax,1:imax,1:jmax)
1121  do k = ks, ke
1122  do j = js, je
1123  do i = is, ie
1124  if ( abs( momy(k,i,j)-momy_check(k,i,j) ) > atmos_restart_check_criterion ) then
1125  write(*,*) 'xxx there is the difference : ', momy(k,i,j)-momy_check(k,i,j)
1126  write(*,*) 'xxx at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'MOMY'
1127  datacheck = .false.
1128  endif
1129  enddo
1130  enddo
1131  enddo
1132 
1133  call fileread( restart_atmos(:,:,:), basename, 'RHOT', 1, prc_myrank )
1134  rhot_check(ks:ke,is:ie,js:je) = restart_atmos(1:kmax,1:imax,1:jmax)
1135  do k = ks, ke
1136  do j = js, je
1137  do i = is, ie
1138  if ( abs( rhot(k,i,j)-rhot_check(k,i,j) ) > atmos_restart_check_criterion ) then
1139  write(*,*) 'xxx there is the difference : ', rhot(k,i,j)-rhot_check(k,i,j)
1140  write(*,*) 'xxx at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, 'RHOT'
1141  datacheck = .false.
1142  endif
1143  enddo
1144  enddo
1145  enddo
1146 
1147  do iq = 1, qa
1148  call fileread( restart_atmos(:,:,:), basename, aq_name(iq), 1, prc_myrank )
1149  qtrc_check(ks:ke,is:ie,js:je,iq) = restart_atmos(1:kmax,1:imax,1:jmax)
1150  do k = ks, ke
1151  do j = js, je
1152  do i = is, ie
1153  if ( abs( qtrc(k,i,j,iq)-qtrc_check(k,i,j,iq) ) > atmos_restart_check_criterion ) then
1154  write(*,*) 'xxx there is the difference : ', qtrc(k,i,j,iq)-qtrc_check(k,i,j,iq)
1155  write(*,*) 'xxx at (PE-id,k,i,j,varname) : ', prc_myrank, k, i, j, aq_name(iq)
1156  datacheck = .false.
1157  endif
1158  enddo
1159  enddo
1160  enddo
1161  enddo
1162 
1163  if (datacheck) then
1164  if( io_l ) write(io_fid_log,*) 'Data Check Clear.'
1165  write(*,*) 'Data Check Clear.'
1166  else
1167  if( io_l ) write(io_fid_log,*) 'Data Check Failed. See std. output.'
1168  write(*,*) 'Data Check Failed.'
1169  endif
1170 
1171  call prof_rapend('Debug')
1172 
1173  return
1174  end subroutine atmos_vars_restart_check
1175 
1176  !-----------------------------------------------------------------------------
1178  subroutine atmos_vars_history
1179  use scale_const, only: &
1180  grav => const_grav, &
1181  rdry => const_rdry, &
1182  rvap => const_rvap, &
1183  cpdry => const_cpdry, &
1184  cvdry => const_cvdry, &
1185  lhv => const_lhv, &
1186  lhf => const_lhf, &
1187  p00 => const_pre00
1188  use scale_grid, only: &
1189  rcdx => grid_rcdx, &
1190  rcdy => grid_rcdy
1191  use scale_grid_real, only: &
1192  real_cz, &
1193  real_fz
1194  use scale_comm, only: &
1196  use scale_history, only: &
1197  hist_in
1198  use scale_atmos_thermodyn, only: &
1199  thermodyn_qd => atmos_thermodyn_qd, &
1200  thermodyn_templhv => atmos_thermodyn_templhv, &
1201  cpw => aq_cp, &
1202  cvw => aq_cv
1203  use scale_atmos_saturation, only: &
1204  saturation_psat_all => atmos_saturation_psat_all, &
1205  saturation_psat_liq => atmos_saturation_psat_liq, &
1206  saturation_psat_ice => atmos_saturation_psat_ice
1207  use scale_atmos_adiabat, only: &
1208  adiabat_cape => atmos_adiabat_cape
1209  implicit none
1210 
1211  real(RP) :: QDRY (ka,ia,ja) ! dry air [kg/kg]
1212  real(RP) :: QTOT (ka,ia,ja) ! total water [kg/kg]
1213  real(RP) :: QHYD (ka,ia,ja) ! total hydrometeor [kg/kg]
1214  real(RP) :: QLIQ (ka,ia,ja) ! total liquid water [kg/kg]
1215  real(RP) :: QICE (ka,ia,ja) ! total ice water [kg/kg]
1216  real(RP) :: RHOQ (ka,ia,ja)
1217 
1218  real(RP) :: LWP (ia,ja) ! liquid water path [g/m2]
1219  real(RP) :: IWP (ia,ja) ! ice water path [g/m2]
1220  real(RP) :: PW (ia,ja) ! precipitable water [g/m2]
1221 
1222  real(RP) :: RTOT (ka,ia,ja) ! Total gas constant [J/kg/K]
1223  real(RP) :: CPTOT (ka,ia,ja) ! Total heat capacity [J/kg/K]
1224  real(RP) :: CPovCV(ka,ia,ja) ! Cp/Cv
1225 
1226  real(RP) :: POTL (ka,ia,ja) ! liquid water potential temperature [K]
1227  real(RP) :: RHA (ka,ia,ja) ! relative humidity (liquid+ice) [%]
1228  real(RP) :: RHL (ka,ia,ja) ! relative humidity against to liquid [%]
1229  real(RP) :: RHI (ka,ia,ja) ! relative humidity against to ice [%]
1230 
1231  real(RP) :: VOR (ka,ia,ja) ! vertical vorticity [1/s]
1232  real(RP) :: DIV (ka,ia,ja) ! divergence [1/s]
1233  real(RP) :: HDIV (ka,ia,ja) ! horizontal divergence [1/s]
1234  real(RP) :: Uabs (ka,ia,ja) ! absolute velocity [m/s]
1235 
1236  real(RP) :: CAPE (ia,ja) ! CAPE [m2/s2]
1237  real(RP) :: CIN (ia,ja) ! CIN [m2/s2]
1238  real(RP) :: LCL (ia,ja) ! LCL height [m]
1239  real(RP) :: LFC (ia,ja) ! LFC height [m]
1240  real(RP) :: LNB (ia,ja) ! LNB height [m]
1241 
1242  real(RP) :: PBLH (ia,ja) ! PBL height [m]
1243  real(RP) :: POTTv (ka,ia,ja) ! vertual potential temperature [K]
1244  real(RP) :: fact
1245 
1246  real(RP) :: MSE (ka,ia,ja) ! MSE [m2/s2]
1247  real(RP) :: LHvap (ka,ia,ja) ! latent heat for vaporization [m2/s2]
1248 
1249  real(RP) :: DENS_PRIM(ka,ia,ja) ! horiz. deviation of density [kg/m3]
1250  real(RP) :: W_PRIM (ka,ia,ja) ! horiz. deviation of w [m/s]
1251  real(RP) :: U_PRIM (ka,ia,ja) ! horiz. deviation of u [m/s]
1252  real(RP) :: V_PRIM (ka,ia,ja) ! horiz. deviation of v [m/s]
1253  real(RP) :: POTT_PRIM(ka,ia,ja) ! horiz. deviation of pot. temp. [K]
1254  real(RP) :: W_PRIM2 (ka,ia,ja) ! variance of w [m2/s2]
1255  real(RP) :: PT_W_PRIM(ka,ia,ja) ! resolved scale heat flux [W/s]
1256  real(RP) :: W_PRIM3 (ka,ia,ja) ! skewness of w [m3/s3]
1257  real(RP) :: TKE_RS (ka,ia,ja) ! resolved scale TKE [m2/s2]
1258  real(RP) :: DENS_MEAN(ka) ! horiz. mean of density [kg/m3]
1259  real(RP) :: W_MEAN (ka) ! horiz. mean of w [m/s]
1260  real(RP) :: U_MEAN (ka) ! horiz. mean of u [m/s]
1261  real(RP) :: V_MEAN (ka) ! horiz. mean of v [m/s]
1262  real(RP) :: PT_MEAN (ka) ! horiz. mean of pot. [K]
1263  real(RP) :: T_MEAN (ka) ! horiz. mean of t [K]
1264  real(RP) :: QV_MEAN (ka) ! horiz. mean of QV
1265  real(RP) :: QHYD_MEAN(ka) ! horiz. mean of QHYD
1266  real(RP) :: QLIQ_MEAN(ka) ! horiz. mean of QLIQ
1267  real(RP) :: QICE_MEAN(ka) ! horiz. mean of QICE
1268 
1269  real(RP) :: ENGT (ka,ia,ja) ! total energy [J/m3]
1270  real(RP) :: ENGP (ka,ia,ja) ! potential energy [J/m3]
1271  real(RP) :: ENGK (ka,ia,ja) ! kinetic energy [J/m3]
1272  real(RP) :: ENGI (ka,ia,ja) ! internal energy [J/m3]
1273 
1274  real(RP) :: PSAT (ka,ia,ja)
1275  real(RP) :: UH (ka,ia,ja)
1276  real(RP) :: VH (ka,ia,ja)
1277 
1278  integer :: k, i, j, iq
1279  !---------------------------------------------------------------------------
1280 
1281  ! value check for prognostic variables
1282  if ( atmos_vars_checkrange ) then
1283  call valcheck( dens(:,:,:), 0.0_rp, 2.0_rp, var_name(i_dens), __file__, __line__ )
1284  call valcheck( momz(:,:,:), -200.0_rp, 200.0_rp, var_name(i_momz), __file__, __line__ )
1285  call valcheck( momx(:,:,:), -200.0_rp, 200.0_rp, var_name(i_momx), __file__, __line__ )
1286  call valcheck( momy(:,:,:), -200.0_rp, 200.0_rp, var_name(i_momy), __file__, __line__ )
1287  call valcheck( rhot(:,:,:), 0.0_rp, 1000.0_rp, var_name(i_rhot), __file__, __line__ )
1288  endif
1289 
1290  ! history output of prognostic variables
1291  call hist_in( dens(:,:,:), var_name(i_dens), var_desc(i_dens), var_unit(i_dens) )
1292  call hist_in( momz(:,:,:), var_name(i_momz), var_desc(i_momz), var_unit(i_momz) )
1293  call hist_in( momx(:,:,:), var_name(i_momx), var_desc(i_momx), var_unit(i_momx) )
1294  call hist_in( momy(:,:,:), var_name(i_momy), var_desc(i_momy), var_unit(i_momy) )
1295  call hist_in( rhot(:,:,:), var_name(i_rhot), var_desc(i_rhot), var_unit(i_rhot) )
1296  do iq = 1, qa
1297  call hist_in( qtrc(:,:,:,iq), aq_name(iq), aq_desc(iq), aq_unit(iq) )
1298  enddo
1299 
1300  ! prepare and history output of diagnostic variables
1301 
1302 ! if ( AD_PREP_sw(I_W) > 0 ) then
1303 ! !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1304 ! do j = JS, JE
1305 ! do i = IS, IE
1306 ! do k = KS, KE
1307 ! W(k,i,j) = 0.5_RP * ( MOMZ(k-1,i,j)+MOMZ(k,i,j) ) / DENS(k,i,j)
1308 ! enddo
1309 ! enddo
1310 ! enddo
1311 ! endif
1312 !
1313 ! if ( AD_PREP_sw(I_U) > 0 ) then
1314 ! !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1315 ! do j = JS, JE
1316 ! do i = IS, IE
1317 ! do k = KS, KE
1318 ! U(k,i,j) = 0.5_RP * ( MOMX(k,i-1,j)+MOMX(k,i,j) ) / DENS(k,i,j)
1319 ! enddo
1320 ! enddo
1321 ! enddo
1322 ! endif
1323 !
1324 ! if ( AD_PREP_sw(I_V) > 0 ) then
1325 ! !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1326 ! do j = JS, JE
1327 ! do i = IS, IE
1328 ! do k = KS, KE
1329 ! V(k,i,j) = 0.5_RP * ( MOMY(k,i,j-1)+MOMY(k,i,j) ) / DENS(k,i,j)
1330 ! enddo
1331 ! enddo
1332 ! enddo
1333 ! endif
1334 !
1335 ! if ( AD_PREP_sw(I_POTT) > 0 ) then
1336 ! !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1337 ! do j = JS, JE
1338 ! do i = IS, IE
1339 ! do k = KS, KE
1340 ! POTT(k,i,j) = RHOT(k,i,j) / DENS(k,i,j)
1341 ! enddo
1342 ! enddo
1343 ! enddo
1344 ! endif
1345 
1346  if ( ad_prep_sw(i_qdry) > 0 ) then
1347  call thermodyn_qd( qdry(:,:,:), & ! [OUT]
1348  qtrc(:,:,:,:) ) ! [IN]
1349  endif
1350 
1351  if ( ad_prep_sw(i_qtot) > 0 ) then
1352 !OCL XFILL
1353  do j = jsb, jeb
1354  do i = isb, ieb
1355  do k = ks, ke
1356  qtot(k,i,j) = 1.0_rp - qdry(k,i,j)
1357  enddo
1358  enddo
1359  enddo
1360  endif
1361 
1362  if ( ad_prep_sw(i_qhyd) > 0 ) then
1363 !OCL XFILL
1364  qhyd(:,:,:) = 0.0_rp
1365 
1366  do iq = qws, qwe
1367  qhyd(:,:,:) = qhyd(:,:,:) + qtrc(:,:,:,iq)
1368  enddo
1369  do iq = qis, qie
1370  qhyd(:,:,:) = qhyd(:,:,:) + qtrc(:,:,:,iq)
1371  enddo
1372  endif
1373 
1374  if ( ad_prep_sw(i_qliq) > 0 ) then
1375 !OCL XFILL
1376  qliq(:,:,:) = 0.0_rp
1377 
1378  do iq = qws, qwe
1379  qliq(:,:,:) = qliq(:,:,:) + qtrc(:,:,:,iq)
1380  enddo
1381  endif
1382 
1383  if ( ad_prep_sw(i_qice) > 0 ) then
1384 !OCL XFILL
1385  qice(:,:,:) = 0.0_rp
1386 
1387  do iq = qis, qie
1388  qice(:,:,:) = qice(:,:,:) + qtrc(:,:,:,iq)
1389  enddo
1390  endif
1391 
1392  if ( ad_prep_sw(i_lwp) > 0 ) then
1393  do j = jsb, jeb
1394  do i = isb, ieb
1395  lwp(i,j) = 0.0_rp
1396  do k = ks, ke
1397  lwp(i,j) = lwp(i,j) &
1398  + qliq(k,i,j) * dens(k,i,j) * ( real_fz(k,i,j)-real_fz(k-1,i,j) ) * 1.e3_rp ! [kg/m2->g/m2]
1399  enddo
1400  enddo
1401  enddo
1402  endif
1403 
1404  if ( ad_prep_sw(i_iwp) > 0 ) then
1405  do j = jsb, jeb
1406  do i = isb, ieb
1407  iwp(i,j) = 0.0_rp
1408  do k = ks, ke
1409  iwp(i,j) = iwp(i,j) &
1410  + qice(k,i,j) * dens(k,i,j) * ( real_fz(k,i,j)-real_fz(k-1,i,j) ) * 1.e3_rp ! [kg/m2->g/m2]
1411  enddo
1412  enddo
1413  enddo
1414  endif
1415 
1416  if ( ad_prep_sw(i_pw) > 0 ) then
1417  do j = jsb, jeb
1418  do i = isb, ieb
1419  pw(i,j) = 0.0_rp
1420  do k = ks, ke
1421  pw(i,j) = pw(i,j) &
1422  + qtrc(k,i,j,i_qv) * dens(k,i,j) * ( real_fz(k,i,j)-real_fz(k-1,i,j) ) * 1.e3_rp ! [kg/m2->g/m2]
1423  enddo
1424  enddo
1425  enddo
1426  endif
1427 
1428  if ( ad_prep_sw(i_rtot) > 0 ) then
1429  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1430 !OCL XFILL
1431  do j = jsb, jeb
1432  do i = isb, ieb
1433  do k = ks, ke
1434  rtot(k,i,j) = rdry * qdry(k,i,j) + rvap * qtrc(k,i,j,i_qv)
1435  enddo
1436  enddo
1437  enddo
1438  endif
1439 
1440  if ( ad_prep_sw(i_cptot) > 0 ) then
1441  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1442 !OCL XFILL
1443  do j = jsb, jeb
1444  do i = isb, ieb
1445  do k = ks, ke
1446  cptot(k,i,j) = cpdry * qdry(k,i,j)
1447  enddo
1448  enddo
1449  enddo
1450 
1451  do iq = qqs, qqe
1452  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1453  do j = jsb, jeb
1454  do i = isb, ieb
1455  do k = ks, ke
1456  cptot(k,i,j) = cptot(k,i,j) + qtrc(k,i,j,iq) * cpw(iq)
1457  enddo
1458  enddo
1459  enddo
1460  enddo
1461 
1462  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1463 !OCL XFILL
1464  do j = jsb, jeb
1465  do i = isb, ieb
1466  do k = ks, ke
1467  cpovcv(k,i,j) = cptot(k,i,j) / ( cptot(k,i,j) - rtot(k,i,j) )
1468  enddo
1469  enddo
1470  enddo
1471  endif
1472 
1473 ! if ( AD_PREP_sw(I_PRES) > 0 ) then
1474 ! !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1475 ! do j = JS, JE
1476 ! do i = IS, IE
1477 ! do k = KS, KE
1478 ! PRES(k,i,j) = P00 * ( RHOT(k,i,j) * RTOT(k,i,j) / P00 )**CPovCV(k,i,j)
1479 ! enddo
1480 ! enddo
1481 ! enddo
1482 ! endif
1483 !
1484 ! if ( AD_PREP_sw(I_TEMP) > 0 ) then
1485 ! !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1486 ! do j = JS, JE
1487 ! do i = IS, IE
1488 ! do k = KS, KE
1489 ! TEMP(k,i,j) = PRES(k,i,j) / ( DENS(k,i,j) * RTOT(k,i,j) )
1490 ! enddo
1491 ! enddo
1492 ! enddo
1493 ! endif
1494 
1495  if ( ad_prep_sw(i_potl) > 0 ) then
1496  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1497 !OCL XFILL
1498  do j = jsb, jeb
1499  do i = isb, ieb
1500  do k = ks, ke
1501  potl(k,i,j) = pott(k,i,j) &
1502  - lhv / cpdry * qliq(k,i,j) * pott(k,i,j) / temp(k,i,j)
1503  enddo
1504  enddo
1505  enddo
1506  endif
1507 
1508  if ( ad_prep_sw(i_qsat) > 0 ) then
1509 ! call SATURATION_dens2qsat_all( QSAT(:,:,:), & ! [OUT]
1510 ! TEMP(:,:,:), & ! [IN]
1511 ! DENS(:,:,:) ) ! [IN]
1512  end if
1513 
1514  if ( ad_hist_id(i_rha) > 0 ) then
1515  call saturation_psat_all( psat(:,:,:), & ! [OUT]
1516  temp(:,:,:) ) ! [IN]
1517  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1518 !OCL XFILL
1519  do j = jsb, jeb
1520  do i = isb, ieb
1521  do k = ks, ke
1522  rha(k,i,j) = dens(k,i,j) * qtrc(k,i,j,i_qv) &
1523  / psat(k,i,j) * rvap * temp(k,i,j) &
1524  * 100.0_rp
1525  enddo
1526  enddo
1527  enddo
1528  endif
1529 
1530  if ( ad_hist_id(i_rhl) > 0 ) then
1531  call saturation_psat_liq( psat(:,:,:), & ! [OUT]
1532  temp(:,:,:) ) ! [IN]
1533  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1534 !OCL XFILL
1535  do j = jsb, jeb
1536  do i = isb, ieb
1537  do k = ks, ke
1538  rhl(k,i,j) = dens(k,i,j) * qtrc(k,i,j,i_qv) &
1539  / psat(k,i,j) * rvap * temp(k,i,j) &
1540  * 100.0_rp
1541  enddo
1542  enddo
1543  enddo
1544  endif
1545 
1546  if ( ad_hist_id(i_rhi) > 0 ) then
1547  call saturation_psat_ice( psat(:,:,:), & ! [OUT]
1548  temp(:,:,:) ) ! [IN]
1549  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1550 !OCL XFILL
1551  do j = jsb, jeb
1552  do i = isb, ieb
1553  do k = ks, ke
1554  rhi(k,i,j) = dens(k,i,j) * qtrc(k,i,j,i_qv) &
1555  / psat(k,i,j) * rvap * temp(k,i,j) &
1556  * 100.0_rp
1557  enddo
1558  enddo
1559  enddo
1560  endif
1561 
1562  if ( ad_prep_sw(i_vor) > 0 ) then
1563  ! at x, v, layer
1564  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1565 !OCL XFILL
1566  do j = 1, ja-1
1567  do i = 2, ia
1568  do k = ks, ke
1569  uh(k,i,j) = 0.5_rp * ( momx(k,i,j)+momx(k,i,j+1)+momx(k,i-1,j)+momx(k,i-1,j+1) ) &
1570  / ( dens(k,i,j)+dens(k,i,j+1) )
1571  enddo
1572  enddo
1573  enddo
1574 
1575  ! at u, y, layer
1576  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1577 !OCL XFILL
1578  do j = 2, ja
1579  do i = 1, ia-1
1580  do k = ks, ke
1581  vh(k,i,j) = 0.5_rp * ( momy(k,i,j)+momy(k,i+1,j)+momy(k,i,j-1)+momy(k,i+1,j-1) ) &
1582  / ( dens(k,i,j)+dens(k,i+1,j) )
1583  enddo
1584  enddo
1585  enddo
1586 
1587  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1588 !OCL XFILL
1589  do j = 2, ja-1
1590  do i = 2, ia-1
1591  do k = ks, ke
1592  vor(k,i,j) = ( vh(k,i,j ) - vh(k,i-1,j ) ) * rcdx(i) &
1593  - ( uh(k,i ,j) - uh(k,i ,j-1) ) * rcdy(j)
1594  enddo
1595  enddo
1596  enddo
1597 
1598  !$omp parallel do private(j,k) OMP_SCHEDULE_
1599  do j = 1, ja
1600  do k = ks, ke
1601  vor(k,1 ,j) = vor(k,2 ,j)
1602  vor(k,ia,j) = vor(k,ia-1,j)
1603  enddo
1604  enddo
1605 
1606  !$omp parallel do private(i,k) OMP_SCHEDULE_
1607  do i = 1, ia
1608  do k = ks, ke
1609  vor(k,i,1 ) = vor(k,i,2 )
1610  vor(k,i,ja) = vor(k,i,ja-1)
1611  enddo
1612  enddo
1613  endif
1614 
1615  if ( ad_prep_sw(i_hdiv) > 0 ) then
1616  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1617 !OCL XFILL
1618  do j = 2, ja
1619  do i = 2, ia
1620  do k = ks, ke
1621  hdiv(k,i,j) = ( momx(k,i,j) - momx(k ,i-1,j ) ) * rcdx(i) &
1622  + ( momy(k,i,j) - momy(k ,i ,j-1) ) * rcdy(j)
1623  enddo
1624  enddo
1625  enddo
1626  !$omp parallel do private(i,k) OMP_SCHEDULE_
1627  do i = 1, ia
1628  do k = ks, ke
1629  hdiv(k,i,1) = hdiv(k,i,2)
1630  enddo
1631  enddo
1632  !$omp parallel do private(j,k) OMP_SCHEDULE_
1633  do j = 1, ja
1634  do k = ks, ke
1635  hdiv(k,1,j) = hdiv(k,2,j)
1636  enddo
1637  enddo
1638  endif
1639 
1640  if ( ad_prep_sw(i_div) > 0 ) then
1641  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1642 !OCL XFILL
1643  do j = 1, ja
1644  do i = 1, ia
1645  do k = ks, ke
1646  div(k,i,j) = ( momz(k,i,j) - momz(k-1,i ,j ) ) * ( real_fz(k,i,j)-real_fz(k-1,i,j) ) &
1647  + hdiv(k,i,j)
1648  enddo
1649  enddo
1650  enddo
1651  endif
1652 
1653  if ( ad_prep_sw(i_uabs) > 0 ) then
1654 !OCL XFILL
1655  do j = 1, ja
1656  do i = 1, ia
1657  do k = ks, ke
1658  uabs(k,i,j) = sqrt( u(k,i,j) * u(k,i,j) &
1659  + v(k,i,j) * v(k,i,j) )
1660  enddo
1661  enddo
1662  enddo
1663  endif
1664 
1665  if ( ad_prep_sw(i_dens_mean) > 0 ) then
1666  call comm_horizontal_mean( dens_mean(:), dens(:,:,:) )
1667  end if
1668 
1669  if ( ad_prep_sw(i_dens_prim) > 0 ) then
1670 !OCL XFILL
1671  do j = 1, ja
1672  do i = 1, ia
1673  do k = ks, ke
1674  dens_prim(k,i,j) = dens(k,i,j) - dens_mean(k)
1675  enddo
1676  enddo
1677  enddo
1678  endif
1679 
1680  if ( ad_prep_sw(i_w_mean) > 0 ) then
1681 !OCL XFILL
1682  do j = 1, ja
1683  do i = 1, ia
1684  do k = ks, ke
1685  w_prim(k,i,j) = w(k,i,j) * dens(k,i,j)
1686  enddo
1687  enddo
1688  enddo
1689  call comm_horizontal_mean( w_mean(:), w_prim(:,:,:) )
1690  do k = ks, ke
1691  w_mean(k) = w_mean(k) / dens_mean(k)
1692  enddo
1693  end if
1694  if ( ad_prep_sw(i_w_prim) > 0 ) then
1695 !OCL XFILL
1696  do j = 1, ja
1697  do i = 1, ia
1698  do k = ks, ke
1699  w_prim(k,i,j) = w(k,i,j) - w_mean(k)
1700  enddo
1701  enddo
1702  enddo
1703  endif
1704 
1705  if ( ad_prep_sw(i_u_mean) > 0 ) then
1706 !OCL XFILL
1707  do j = 1, ja
1708  do i = 1, ia
1709  do k = ks, ke
1710  u_prim(k,i,j) = u(k,i,j) * dens(k,i,j)
1711  enddo
1712  enddo
1713  enddo
1714  call comm_horizontal_mean( u_mean(:), u_prim(:,:,:) )
1715  do k = ks, ke
1716  u_mean(k) = u_mean(k) / dens_mean(k)
1717  enddo
1718  end if
1719  if ( ad_prep_sw(i_u_prim) > 0 ) then
1720 !OCL XFILL
1721  do j = 1, ja
1722  do i = 1, ia
1723  do k = ks, ke
1724  u_prim(k,i,j) = u(k,i,j) - u_mean(k)
1725  enddo
1726  enddo
1727  enddo
1728  endif
1729 
1730  if ( ad_prep_sw(i_v_mean) > 0 ) then
1731 !OCL XFILL
1732  do j = 1, ja
1733  do i = 1, ia
1734  do k = ks, ke
1735  v_prim(k,i,j) = v(k,i,j) * dens(k,i,j)
1736  enddo
1737  enddo
1738  enddo
1739  call comm_horizontal_mean( v_mean(:), v_prim(:,:,:) )
1740  do k = ks, ke
1741  v_mean(k) = v_mean(k) / dens_mean(k)
1742  enddo
1743  end if
1744  if ( ad_prep_sw(i_v_prim) > 0 ) then
1745 !OCL XFILL
1746  do j = 1, ja
1747  do i = 1, ia
1748  do k = ks, ke
1749  v_prim(k,i,j) = v(k,i,j) - v_mean(k)
1750  enddo
1751  enddo
1752  enddo
1753  endif
1754 
1755  if ( ad_prep_sw(i_t_mean) > 0 ) then
1756 !OCL XFILL
1757  do j = 1, ja
1758  do i = 1, ia
1759  do k = ks, ke
1760  pott_prim(k,i,j) = temp(k,i,j) * dens(k,i,j)
1761  enddo
1762  enddo
1763  enddo
1764  call comm_horizontal_mean( t_mean(:), pott_prim(:,:,:) )
1765  do k = ks, ke
1766  t_mean(k) = t_mean(k) / dens_mean(k)
1767  enddo
1768  end if
1769 
1770  if ( ad_prep_sw(i_pott_mean) > 0 ) then
1771  call comm_horizontal_mean( pt_mean(:), rhot(:,:,:) )
1772  do k = ks, ke
1773  pt_mean(k) = pt_mean(k) / dens_mean(k)
1774  enddo
1775  end if
1776  if ( ad_prep_sw(i_pott_prim) > 0 ) then
1777 !OCL XFILL
1778  do j = 1, ja
1779  do i = 1, ia
1780  do k = ks, ke
1781  pott_prim(k,i,j) = pott(k,i,j) - pt_mean(k)
1782  enddo
1783  enddo
1784  enddo
1785  endif
1786 
1787  if ( ad_prep_sw(i_qv_mean) > 0 ) then
1788 !OCL XFILL
1789  do j = 1, ja
1790  do i = 1, ia
1791  do k = ks, ke
1792  rhoq(k,i,j) = qtrc(k,i,j,i_qv) * dens(k,i,j)
1793  enddo
1794  enddo
1795  enddo
1796  call comm_horizontal_mean( qv_mean(:), rhoq(:,:,:) )
1797  do k = ks, ke
1798  qv_mean(k) = qv_mean(k) / dens_mean(k)
1799  enddo
1800  end if
1801 
1802  if ( ad_prep_sw(i_qhyd_mean) > 0 ) then
1803 !OCL XFILL
1804  do j = 1, ja
1805  do i = 1, ia
1806  do k = ks, ke
1807  rhoq(k,i,j) = qhyd(k,i,j) * dens(k,i,j)
1808  enddo
1809  enddo
1810  enddo
1811  call comm_horizontal_mean( qhyd_mean(:), rhoq(:,:,:) )
1812  do k = ks, ke
1813  qhyd_mean(k) = qhyd_mean(k) / dens_mean(k)
1814  enddo
1815  end if
1816 
1817  if ( ad_prep_sw(i_qliq_mean) > 0 ) then
1818 !OCL XFILL
1819  do j = 1, ja
1820  do i = 1, ia
1821  do k = ks, ke
1822  rhoq(k,i,j) = qliq(k,i,j) * dens(k,i,j)
1823  enddo
1824  enddo
1825  enddo
1826  call comm_horizontal_mean( qliq_mean(:), rhoq(:,:,:) )
1827  do k = ks, ke
1828  qliq_mean(k) = qliq_mean(k) / dens_mean(k)
1829  enddo
1830  end if
1831 
1832  if ( ad_prep_sw(i_qice_mean) > 0 ) then
1833 !OCL XFILL
1834  do j = 1, ja
1835  do i = 1, ia
1836  do k = ks, ke
1837  rhoq(k,i,j) = qice(k,i,j) * dens(k,i,j)
1838  enddo
1839  enddo
1840  enddo
1841  call comm_horizontal_mean( qice_mean(:), rhoq(:,:,:) )
1842  do k = ks, ke
1843  qice_mean(k) = qice_mean(k) / dens_mean(k)
1844  enddo
1845  end if
1846 
1847  if ( ad_prep_sw(i_w_prim2) > 0 ) then
1848 !OCL XFILL
1849  do j = 1, ja
1850  do i = 1, ia
1851  do k = ks, ke
1852  w_prim2(k,i,j) = w_prim(k,i,j) * w_prim(k,i,j)
1853  enddo
1854  enddo
1855  enddo
1856  endif
1857 
1858  if ( ad_prep_sw(i_pt_w_prim) > 0 ) then
1859 !OCL XFILL
1860  do j = 1, ja
1861  do i = 1, ia
1862  do k = ks, ke
1863  pt_w_prim(k,i,j) = w_prim(k,i,j) * pott_prim(k,i,j) * dens(k,i,j) * cpdry
1864  enddo
1865  enddo
1866  enddo
1867  endif
1868 
1869  if ( ad_prep_sw(i_w_prim3) > 0 ) then
1870 !OCL XFILL
1871  do j = 1, ja
1872  do i = 1, ia
1873  do k = ks, ke
1874  w_prim3(k,i,j) = w_prim(k,i,j) * w_prim(k,i,j) * w_prim(k,i,j)
1875  enddo
1876  enddo
1877  enddo
1878  endif
1879 
1880  if ( ad_prep_sw(i_tke_rs) > 0 ) then
1881 !OCL XFILL
1882  do j = 1, ja
1883  do i = 1, ia
1884  do k = ks, ke
1885  tke_rs(k,i,j) = 0.5_rp * ( w_prim(k,i,j) * w_prim(k,i,j) &
1886  + u_prim(k,i,j) * u_prim(k,i,j) &
1887  + v_prim(k,i,j) * v_prim(k,i,j) )
1888  enddo
1889  enddo
1890  enddo
1891  endif
1892 
1893  if ( ad_prep_sw(i_engp) > 0 ) then
1894  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1895 !OCL XFILL
1896  do j = 1, ja
1897  do i = 1, ia
1898  do k = ks, ke
1899  engp(k,i,j) = dens(k,i,j) * grav * real_cz(k,i,j)
1900  enddo
1901  enddo
1902  enddo
1903  endif
1904 
1905  if ( ad_prep_sw(i_engk) > 0 ) then
1906  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1907 !OCL XFILL
1908  do j = 1, ja
1909  do i = 1, ia
1910  do k = ks, ke
1911  engk(k,i,j) = 0.5_rp * dens(k,i,j) * ( w(k,i,j)**2 &
1912  + u(k,i,j)**2 &
1913  + v(k,i,j)**2 )
1914  enddo
1915  enddo
1916  enddo
1917  endif
1918 
1919  if ( ad_prep_sw(i_engi) > 0 ) then
1920  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1921 !OCL XFILL
1922  do j = jsb, jeb
1923  do i = isb, ieb
1924  do k = ks, ke
1925  engi(k,i,j) = dens(k,i,j) * qdry(k,i,j) * temp(k,i,j) * cvdry
1926  enddo
1927  enddo
1928  enddo
1929 
1930  do iq = qqs, qqe
1931  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
1932  do j = jsb, jeb
1933  do i = isb, ieb
1934  do k = ks, ke
1935  engi(k,i,j) = engi(k,i,j) &
1936  + dens(k,i,j) * qtrc(k,i,j,iq) * temp(k,i,j) * cvw(iq)
1937  enddo
1938  enddo
1939  enddo
1940  enddo
1941 
1942  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
1943  do j = jsb, jeb
1944  do i = isb, ieb
1945  do k = ks, ke
1946  engi(k,i,j) = engi(k,i,j) &
1947  + dens(k,i,j) * qtrc(k,i,j,i_qv) * lhv ! Latent Heat [vapor->liquid]
1948  enddo
1949  enddo
1950  enddo
1951 
1952  do iq = qis, qie
1953  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
1954  do j = jsb, jeb
1955  do i = isb, ieb
1956  do k = ks, ke
1957  engi(k,i,j) = engi(k,i,j) &
1958  - dens(k,i,j) * qtrc(k,i,j,iq) * lhf ! Latent Heat [ice->liquid]
1959  enddo
1960  enddo
1961  enddo
1962  enddo
1963  endif
1964 
1965  if ( ad_prep_sw(i_engt) > 0 ) then
1966  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1967 !OCL XFILL
1968  do j = jsb, jeb
1969  do i = isb, ieb
1970  do k = ks, ke
1971  engt(k,i,j) = engp(k,i,j) + engk(k,i,j) + engi(k,i,j)
1972  enddo
1973  enddo
1974  enddo
1975  endif
1976 
1977  call hist_in( w(:,:,:), 'W', 'velocity w', 'm/s' )
1978  call hist_in( u(:,:,:), 'U', 'velocity u', 'm/s' )
1979  call hist_in( v(:,:,:), 'V', 'velocity v', 'm/s' )
1980  call hist_in( pott(:,:,:), 'PT', 'potential temp.', 'K' )
1981 
1982  call hist_in( qdry(:,:,:), 'QDRY', 'dry air', 'kg/kg' )
1983  call hist_in( qtot(:,:,:), 'QTOT', 'total water', 'kg/kg' )
1984  call hist_in( qhyd(:,:,:), 'QHYD', 'total hydrometeors', 'kg/kg' )
1985  call hist_in( qliq(:,:,:), 'QLIQ', 'total liquid water', 'kg/kg' )
1986  call hist_in( qice(:,:,:), 'QICE', 'total ice water', 'kg/kg' )
1987 
1988  call hist_in( lwp(:,:), 'LWP', 'liquid water path', 'g/m2' )
1989  call hist_in( iwp(:,:), 'IWP', 'ice water path', 'g/m2' )
1990  call hist_in( pw(:,:), 'PW', 'precipitable water', 'g/m2' )
1991 
1992  call hist_in( rtot(:,:,:), 'RTOT', 'Total gas constant', 'J/kg/K' )
1993  call hist_in( cptot(:,:,:), 'CPTOT', 'Total heat capacity', 'J/kg/K' )
1994  call hist_in( pres(:,:,:), 'PRES', 'pressure', 'Pa' )
1995  call hist_in( temp(:,:,:), 'T', 'temperature', 'K' )
1996 
1997  call hist_in( potl(:,:,:), 'LWPT', 'liq. potential temp.', 'K' )
1998  call hist_in( rha(:,:,:), 'RHA', 'relative humidity(liq+ice)','%' )
1999  call hist_in( rhl(:,:,:), 'RH' , 'relative humidity(liq)', '%' )
2000  call hist_in( rhi(:,:,:), 'RHI', 'relative humidity(ice)', '%' )
2001 
2002  call hist_in( vor(:,:,:), 'VOR', 'vertical vorticity', '1/s' )
2003  call hist_in( div(:,:,:), 'DIV', 'divergence', '1/s' )
2004  call hist_in( hdiv(:,:,:), 'HDIV', 'horizontal divergence', '1/s' )
2005  call hist_in( uabs(:,:,:), 'Uabs', 'absolute velocity', 'm/s' )
2006 
2007  call hist_in( dens_mean(:), 'DENS_MEAN', 'horiz. mean of density', 'kg/m3' )
2008  call hist_in( w_mean(:), 'W_MEAN', 'horiz. mean of w', 'm/s' )
2009  call hist_in( u_mean(:), 'U_MEAN', 'horiz. mean of u', 'm/s' )
2010  call hist_in( v_mean(:), 'V_MEAN', 'horiz. mean of v', 'm/s' )
2011  call hist_in( pt_mean(:), 'PT_MEAN', 'horiz. mean of pot.', 'K' )
2012  call hist_in( t_mean(:), 'T_MEAN', 'horiz. mean of t', 'K' )
2013  call hist_in( qv_mean(:), 'QV_MEAN', 'horiz. mean of QV', '1' )
2014  call hist_in( qhyd_mean(:), 'QHYD_MEAN', 'horiz. mean of QHYD', '1' )
2015  call hist_in( qliq_mean(:), 'QLIQ_MEAN', 'horiz. mean of QLIQ', '1' )
2016  call hist_in( qice_mean(:), 'QICE_MEAN', 'horiz. mean of QICE', '1' )
2017 
2018  call hist_in( dens_prim(:,:,:), 'DENS_PRIM', 'horiz. deviation of density', 'kg/m3' )
2019  call hist_in( w_prim(:,:,:), 'W_PRIM', 'horiz. deviation of w', 'm/s' )
2020  call hist_in( u_prim(:,:,:), 'U_PRIM', 'horiz. deviation of u', 'm/s' )
2021  call hist_in( v_prim(:,:,:), 'V_PRIM', 'horiz. deviation of v', 'm/s' )
2022  call hist_in( pott_prim(:,:,:), 'PT_PRIM', 'horiz. deviation of pot. temp.', 'K' )
2023  call hist_in( w_prim2(:,:,:), 'W_PRIM2', 'variance of w', 'm2/s2' )
2024  call hist_in( pt_w_prim(:,:,:), 'PT_W_PRIM', 'resolved scale heat flux', 'W/s' )
2025  call hist_in( w_prim3(:,:,:), 'W_PRIM3', 'skewness of w', 'm3/s3' )
2026  call hist_in( tke_rs(:,:,:), 'TKE_RS', 'resolved scale TKE', 'm2/s2' )
2027 
2028  call hist_in( engt(:,:,:), 'ENGT', 'total energy', 'J/m3' )
2029  call hist_in( engp(:,:,:), 'ENGP', 'potential energy', 'J/m3' )
2030  call hist_in( engk(:,:,:), 'ENGK', 'kinetic energy', 'J/m3' )
2031  call hist_in( engi(:,:,:), 'ENGI', 'internal energy', 'J/m3' )
2032 
2033  if ( ad_prep_sw(i_cape) > 0 &
2034  .OR. ad_prep_sw(i_cin) > 0 &
2035  .OR. ad_prep_sw(i_lcl) > 0 &
2036  .OR. ad_prep_sw(i_lfc) > 0 &
2037  .OR. ad_prep_sw(i_lnb) > 0 ) then
2038 
2039  call adiabat_cape( ks, & ! [IN]
2040  dens(:,:,:), & ! [IN]
2041  temp(:,:,:), & ! [IN]
2042  pres(:,:,:), & ! [IN]
2043  qtrc(:,:,:,:), & ! [IN]
2044  real_cz(:,:,:), & ! [IN]
2045  real_fz(:,:,:), & ! [IN]
2046  cape(:,:), & ! [OUT]
2047  cin(:,:), & ! [OUT]
2048  lcl(:,:), & ! [OUT]
2049  lfc(:,:), & ! [OUT]
2050  lnb(:,:) ) ! [OUT]
2051 
2052  endif
2053 
2054  call hist_in( cape(:,:), 'CAPE', 'convection avail. pot. energy', 'm2/s2' )
2055  call hist_in( cin(:,:), 'CIN', 'convection inhibition', 'm2/s2' )
2056  call hist_in( lcl(:,:), 'LCL', 'lifted condensation level', 'm' )
2057  call hist_in( lfc(:,:), 'LFC', 'level of free convection', 'm' )
2058  call hist_in( lnb(:,:), 'LNB', 'level of neutral buoyancy', 'm' )
2059 
2060  if ( ad_prep_sw(i_pblh) > 0 ) then
2061  do j = js, je
2062  do i = is, ie
2063  do k = ks, ke
2064  pottv(k,i,j) = pott(k,i,j) * ( 1.0_rp + 0.61_rp * qtrc(k,i,j,i_qv) - 1.61 * qtrc(k,i,j,i_qc) )
2065  enddo
2066  enddo
2067  enddo
2068 
2069  do j = js, je
2070  do i = is, ie
2071  pblh(i,j) = real_cz(ks,i,j) - real_fz(ks-1,i,j)
2072 
2073  do k = ks+1, ke
2074  if ( pottv(k,i,j) > pottv(ks,i,j) ) then
2075  fact = ( pottv(ks,i,j) - pottv(k-1,i,j) ) &
2076  / ( pottv(k,i,j) - pottv(k-1,i,j) )
2077 
2078  pblh(i,j) = real_cz(k-1,i,j) - real_fz(ks-1,i,j) &
2079  + fact * ( real_cz(k,i,j) - real_cz(k-1,i,j) )
2080 
2081  exit
2082  endif
2083  enddo
2084  enddo
2085  enddo
2086  endif
2087  call hist_in( pblh(:,:), 'PBLH', 'PBL height', 'm' )
2088 
2089  if ( ad_prep_sw(i_mse) > 0 ) then
2090  call thermodyn_templhv( lhvap(:,:,:), & ! [OUT]
2091  temp(:,:,:) ) ! [IN]
2092 
2093  do j = js, je
2094  do i = is, ie
2095  do k = ks, ke
2096  mse(k,i,j) = cptot(k,i,j) * temp(k,i,j) &
2097  + grav * ( real_cz(k,i,j) - real_fz(ks-1,i,j) ) &
2098  + lhvap(k,i,j) * qtrc(k,i,j,i_qv)
2099  enddo
2100  enddo
2101  enddo
2102  endif
2103  call hist_in( mse(:,:,:), 'MSE', 'moist static energy', 'm2/s2' )
2104 
2105  return
2106  end subroutine atmos_vars_history
2107 
2108  !-----------------------------------------------------------------------------
2110  subroutine atmos_vars_total
2111  use scale_const, only: &
2112  grav => const_grav, &
2113  cvdry => const_cvdry, &
2114  lhv => const_lhv, &
2115  lhf => const_lhf
2116  use scale_grid_real, only: &
2117  real_cz
2118  use scale_rm_statistics, only: &
2120  stat_total
2121  use scale_atmos_thermodyn, only: &
2122  thermodyn_qd => atmos_thermodyn_qd, &
2123  thermodyn_temp_pres => atmos_thermodyn_temp_pres, &
2124  cvw => aq_cv
2125  implicit none
2126 
2127  real(RP) :: W (ka,ia,ja) ! velocity w at cell center [m/s]
2128  real(RP) :: U (ka,ia,ja) ! velocity u at cell center [m/s]
2129  real(RP) :: V (ka,ia,ja) ! velocity v at cell center [m/s]
2130 
2131  real(RP) :: QDRY(ka,ia,ja) ! dry air [kg/kg]
2132  real(RP) :: PRES(ka,ia,ja) ! pressure [Pa]
2133  real(RP) :: TEMP(ka,ia,ja) ! temperature [K]
2134 
2135  real(RP) :: ENGT(ka,ia,ja) ! total energy [J/m3]
2136  real(RP) :: ENGP(ka,ia,ja) ! potential energy [J/m3]
2137  real(RP) :: ENGK(ka,ia,ja) ! kinetic energy [J/m3]
2138  real(RP) :: ENGI(ka,ia,ja) ! internal energy [J/m3]
2139 
2140  real(RP) :: RHOQ(ka,ia,ja)
2141 
2142  real(RP) :: total ! dummy
2143  integer :: i, j, k, iq
2144  !---------------------------------------------------------------------------
2145 
2146  if ( statistics_checktotal ) then
2147 
2148  call stat_total( total, dens(:,:,:), var_name(i_dens) )
2149  call stat_total( total, momz(:,:,:), var_name(i_momz) )
2150  call stat_total( total, momx(:,:,:), var_name(i_momx) )
2151  call stat_total( total, momy(:,:,:), var_name(i_momy) )
2152  call stat_total( total, rhot(:,:,:), var_name(i_rhot) )
2153 
2154  do iq = 1, qa
2155  rhoq(:,:,:) = dens(:,:,:) * qtrc(:,:,:,iq)
2156 
2157  call stat_total( total, rhoq(:,:,:), aq_name(iq) )
2158  enddo
2159 
2160  call thermodyn_qd( qdry(:,:,:), & ! [OUT]
2161  qtrc(:,:,:,:) ) ! [IN]
2162 
2163  call thermodyn_temp_pres( temp(:,:,:), & ! [OUT]
2164  pres(:,:,:), & ! [OUT]
2165  dens(:,:,:), & ! [IN]
2166  rhot(:,:,:), & ! [IN]
2167  qtrc(:,:,:,:) ) ! [IN]
2168 
2169  rhoq(ks:ke,is:ie,js:je) = dens(ks:ke,is:ie,js:je) * qdry(ks:ke,is:ie,js:je)
2170 
2171  call stat_total( total, rhoq(:,:,:), 'QDRY' )
2172 
2173  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
2174 
2175  call stat_total( total, rhoq(:,:,:), 'QTOT' )
2176 
2177  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2178  do j = js, je
2179  do i = is, ie
2180  do k = ks, ke
2181  w(k,i,j) = 0.5_rp * ( momz(k-1,i,j)+momz(k,i,j) ) / dens(k,i,j)
2182  u(k,i,j) = 0.5_rp * ( momx(k,i-1,j)+momx(k,i,j) ) / dens(k,i,j)
2183  v(k,i,j) = 0.5_rp * ( momy(k,i,j-1)+momy(k,i,j) ) / dens(k,i,j)
2184 
2185  engp(k,i,j) = dens(k,i,j) * grav * real_cz(k,i,j)
2186 
2187  engk(k,i,j) = 0.5_rp * dens(k,i,j) * ( w(k,i,j)**2 &
2188  + u(k,i,j)**2 &
2189  + v(k,i,j)**2 )
2190 
2191  engi(k,i,j) = dens(k,i,j) * qdry(k,i,j) * temp(k,i,j) * cvdry
2192  do iq = qqs, qqe
2193  engi(k,i,j) = engi(k,i,j) &
2194  + dens(k,i,j) * qtrc(k,i,j,iq) * temp(k,i,j) * cvw(iq)
2195  enddo
2196 
2197  engi(k,i,j) = engi(k,i,j) + dens(k,i,j) * qtrc(k,i,j,i_qv) * lhv ! Latent Heat [vapor->liquid]
2198 
2199  do iq = qis, qie
2200  engi(k,i,j) = engi(k,i,j) - dens(k,i,j) * qtrc(k,i,j,iq) * lhf ! Latent Heat [ice->liquid]
2201  enddo
2202 
2203  engt(k,i,j) = engp(k,i,j) + engk(k,i,j) + engi(k,i,j)
2204  enddo
2205  enddo
2206  enddo
2207 
2208  call stat_total( total, engp(:,:,:), 'ENGP' )
2209  call stat_total( total, engk(:,:,:), 'ENGK' )
2210  call stat_total( total, engi(:,:,:), 'ENGI' )
2211  call stat_total( total, engt(:,:,:), 'ENGT' )
2212 
2213  endif
2214 
2215  return
2216  end subroutine atmos_vars_total
2217 
2218  !-----------------------------------------------------------------------------
2220  subroutine atmos_vars_diagnostics
2221  use scale_comm, only: &
2222  comm_vars8, &
2223  comm_wait
2224  use scale_atmos_thermodyn, only: &
2225  thermodyn_temp_pres => atmos_thermodyn_temp_pres
2226  implicit none
2227 
2228  integer :: k, i, j
2229  !---------------------------------------------------------------------------
2230 
2231  if( io_l ) write(io_fid_log,*) '*** Calc diagnostics'
2232 
2233  call thermodyn_temp_pres( temp(:,:,:), & ! [OUT]
2234  pres(:,:,:), & ! [OUT]
2235  dens(:,:,:), & ! [IN]
2236  rhot(:,:,:), & ! [IN]
2237  qtrc(:,:,:,:) ) ! [IN]
2238 
2239 !OCL XFILL
2240  do j = 1, ja
2241  do i = 1, ia
2242  do k = ks+1, ke-1
2243  w(k,i,j) = 0.5_rp * ( momz(k-1,i,j)+momz(k,i,j) ) / dens(k,i,j)
2244  enddo
2245  enddo
2246  enddo
2247 !OCL XFILL
2248  do j = 1, ja
2249  do i = 1, ia
2250  w(ks,i,j) = 0.5_rp * ( momz(ks,i,j) ) / dens(ks,i,j)
2251  enddo
2252  enddo
2253 !OCL XFILL
2254  do j = 1, ja
2255  do i = 1, ia
2256  w(ke,i,j) = 0.5_rp * ( momz(ke-1,i,j) ) / dens(ke,i,j)
2257  enddo
2258  enddo
2259 
2260 !OCL XFILL
2261  do j = 1, ja
2262  do i = 2, ia
2263  do k = ks, ke
2264  u(k,i,j) = 0.5_rp * ( momx(k,i-1,j)+momx(k,i,j) ) / dens(k,i,j)
2265  enddo
2266  enddo
2267  enddo
2268 !OCL XFILL
2269  do j = 1, ja
2270  do k = ks, ke
2271  u(k,1,j) = momx(k,1,j) / dens(k,1,j)
2272  enddo
2273  enddo
2274 
2275 !OCL XFILL
2276  do j = 2, ja
2277  do i = 1, ia
2278  do k = ks, ke
2279  v(k,i,j) = 0.5_rp * ( momy(k,i,j-1)+momy(k,i,j) ) / dens(k,i,j)
2280  enddo
2281  enddo
2282  enddo
2283 !OCL XFILL
2284  do i = 1, ia
2285  do k = ks, ke
2286  v(k,i,1) = momy(k,i,1) / dens(k,i,1)
2287  enddo
2288  enddo
2289 
2290  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
2291  do j = 1, ja
2292  do i = 1, ia
2293  w( 1:ks-1,i,j) = w(ks,i,j)
2294  u( 1:ks-1,i,j) = u(ks,i,j)
2295  v( 1:ks-1,i,j) = v(ks,i,j)
2296  w(ke+1:ka, i,j) = w(ke,i,j)
2297  u(ke+1:ka, i,j) = u(ke,i,j)
2298  v(ke+1:ka, i,j) = v(ke,i,j)
2299  enddo
2300  enddo
2301 
2302  call comm_vars8( u(:,:,:), 1 )
2303  call comm_vars8( v(:,:,:), 2 )
2304  call comm_wait ( u(:,:,:), 1, .false. )
2305  call comm_wait ( v(:,:,:), 2, .false. )
2306 
2307 !OCL XFILL
2308  do j = 1, ja
2309  do i = 1, ia
2310  do k = ks, ke
2311  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
2312  enddo
2313  enddo
2314  enddo
2315 
2316  return
2317  end subroutine atmos_vars_diagnostics
2318 
2319  !-----------------------------------------------------------------------------
2321  subroutine atmos_vars_monitor
2322  use scale_process, only: &
2323  prc_myrank
2324  use scale_const, only: &
2325  grav => const_grav, &
2326  cvdry => const_cvdry, &
2327  lhv => const_lhv, &
2328  lhf => const_lhf
2329  use scale_grid, only: &
2330  rfdx => grid_rfdx, &
2331  rfdy => grid_rfdy
2332  use scale_grid_real, only: &
2333  real_cz, &
2334  real_fz
2335  use scale_gridtrans, only: &
2336  mapf => gtrans_mapf, &
2337  i_uy, &
2338  i_xv
2339  use scale_comm, only: &
2340  comm_vars8, &
2341  comm_wait
2342  use scale_rm_statistics, only: &
2344  stat_total, &
2345  stat_detail
2346  use scale_monitor, only: &
2347  monit_put, &
2348  monit_in
2349  use scale_time, only: &
2351  use scale_atmos_thermodyn, only: &
2352  thermodyn_qd => atmos_thermodyn_qd, &
2353  thermodyn_temp_pres => atmos_thermodyn_temp_pres, &
2354  cvw => aq_cv
2355  use mod_atmos_phy_mp_vars, only: &
2356  sflx_rain => atmos_phy_mp_sflx_rain, &
2357  sflx_snow => atmos_phy_mp_sflx_snow
2358  use mod_atmos_phy_rd_vars, only: &
2359  sflx_lw_up => atmos_phy_rd_sflx_lw_up, &
2360  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
2361  sflx_sw_up => atmos_phy_rd_sflx_sw_up, &
2362  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn, &
2363  toaflx_lw_up => atmos_phy_rd_toaflx_lw_up, &
2364  toaflx_lw_dn => atmos_phy_rd_toaflx_lw_dn, &
2365  toaflx_sw_up => atmos_phy_rd_toaflx_sw_up, &
2366  toaflx_sw_dn => atmos_phy_rd_toaflx_sw_dn
2367  use mod_atmos_phy_sf_vars, only: &
2368  sflx_sh => atmos_phy_sf_sflx_sh, &
2369  sflx_lh => atmos_phy_sf_sflx_lh, &
2370  sflx_qtrc => atmos_phy_sf_sflx_qtrc
2371  implicit none
2372 
2373  real(RP) :: QDRY(ka,ia,ja) ! dry air [kg/kg]
2374  real(RP) :: RHOQ(ka,ia,ja) ! DENS * tracer [kg/m3]
2375  real(RP) :: PRCP(ia,ja) ! rain + snow [kg/m2/s]
2376 
2377  real(RP) :: ENGT(ka,ia,ja) ! total energy [J/m3]
2378  real(RP) :: ENGP(ka,ia,ja) ! potential energy [J/m3]
2379  real(RP) :: ENGK(ka,ia,ja) ! kinetic energy [J/m3]
2380  real(RP) :: ENGI(ka,ia,ja) ! internal energy [J/m3]
2381 
2382  real(RP) :: ENGFLXT (ia,ja) ! total flux [J/m2/s]
2383  real(RP) :: SFLX_RD_net (ia,ja) ! net SFC radiation flux [J/m2/s]
2384  real(RP) :: TOAFLX_RD_net(ia,ja) ! net TOA radiation flux [J/m2/s]
2385 
2386  real(RP) :: WORK (ka,ia,ja,3)
2387  character(len=H_SHORT) :: WNAME(3)
2388  real(RP) :: CFLMAX
2389 
2390  integer :: k, i, j, iq
2391  !---------------------------------------------------------------------------
2392 
2393  if( io_l ) write(io_fid_log,*) '*** Monitor'
2394 
2395  call monit_in( dens(:,:,:), var_name(i_dens), var_desc(i_dens), var_unit(i_dens), ndim=3, isflux=.false. )
2396  call monit_in( momz(:,:,:), var_name(i_momz), var_desc(i_momz), var_unit(i_momz), ndim=3, isflux=.false. )
2397  call monit_in( momx(:,:,:), var_name(i_momx), var_desc(i_momx), var_unit(i_momx), ndim=3, isflux=.false. )
2398  call monit_in( momy(:,:,:), var_name(i_momy), var_desc(i_momy), var_unit(i_momy), ndim=3, isflux=.false. )
2399  call monit_in( rhot(:,:,:), var_name(i_rhot), var_desc(i_rhot), var_unit(i_rhot), ndim=3, isflux=.false. )
2400 
2401  !##### Mass Budget #####
2402  do iq = 1, qa
2403  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2404 !OCL XFILL
2405  do j = js, je
2406  do i = is, ie
2407  do k = ks, ke
2408  rhoq(k,i,j) = dens(k,i,j) * qtrc(k,i,j,iq)
2409  enddo
2410  enddo
2411  enddo
2412 
2413  call monit_in( rhoq(:,:,:), aq_name(iq), aq_desc(iq), aq_unit(iq), ndim=3, isflux=.false. )
2414  enddo
2415 
2416  ! total dry airmass
2417 
2418  call thermodyn_qd( qdry(:,:,:), & ! [OUT]
2419  qtrc(:,:,:,:) ) ! [IN]
2420 
2421  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2422 !OCL XFILL
2423  do j = js, je
2424  do i = is, ie
2425  do k = ks, ke
2426  rhoq(k,i,j) = dens(k,i,j) * qdry(k,i,j)
2427  enddo
2428  enddo
2429  enddo
2430  call monit_put( ad_monit_id(i_qdry), rhoq(:,:,:) )
2431 
2432  ! total vapor,liquid,solid tracers
2433  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2434 !OCL XFILL
2435  do j = js, je
2436  do i = is, ie
2437  do k = ks, ke
2438  rhoq(k,i,j) = dens(k,i,j) * ( 1.0_rp - qdry(k,i,j) )
2439  enddo
2440  enddo
2441  enddo
2442  call monit_put( ad_monit_id(i_qtot), rhoq(:,:,:) )
2443 
2444  ! total evapolation
2445  call monit_put( ad_monit_id(i_evap), sflx_qtrc(:,:,i_qv) )
2446 
2447  ! total precipitation
2448 !OCL XFILL
2449  do j = js, je
2450  do i = is, ie
2451  prcp(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
2452  enddo
2453  enddo
2454  call monit_put( ad_monit_id(i_prcp), prcp(:,:) )
2455 
2456  !##### Energy Budget #####
2457 
2458  call thermodyn_temp_pres( temp(:,:,:), & ! [OUT]
2459  pres(:,:,:), & ! [OUT]
2460  dens(:,:,:), & ! [IN]
2461  rhot(:,:,:), & ! [IN]
2462  qtrc(:,:,:,:) ) ! [IN]
2463 
2464  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2465  do j = js, je
2466  do i = is, ie
2467  do k = ks, ke
2468  engp(k,i,j) = dens(k,i,j) * grav * real_cz(k,i,j)
2469 
2470  engk(k,i,j) = 0.5_rp * dens(k,i,j) * ( w(k,i,j)**2 &
2471  + u(k,i,j)**2 &
2472  + v(k,i,j)**2 )
2473 
2474  engi(k,i,j) = dens(k,i,j) * qdry(k,i,j) * temp(k,i,j) * cvdry
2475  do iq = qqs, qqe
2476  engi(k,i,j) = engi(k,i,j) &
2477  + dens(k,i,j) * qtrc(k,i,j,iq) * temp(k,i,j) * cvw(iq)
2478  enddo
2479 
2480  engi(k,i,j) = engi(k,i,j) + dens(k,i,j) * qtrc(k,i,j,i_qv) * lhv ! Latent Heat [vapor->liquid]
2481 
2482  do iq = qis, qie
2483  engi(k,i,j) = engi(k,i,j) - dens(k,i,j) * qtrc(k,i,j,iq) * lhf ! Latent Heat [ice->liquid]
2484  enddo
2485  enddo
2486  enddo
2487  enddo
2488 
2489  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2490 !OCL XFILL
2491  do j = js, je
2492  do i = is, ie
2493  do k = ks, ke
2494  engt(k,i,j) = engp(k,i,j) + engk(k,i,j) + engi(k,i,j)
2495  enddo
2496  enddo
2497  enddo
2498 
2499 !OCL XFILL
2500  do j = js, je
2501  do i = is, ie
2502  sflx_rd_net(i,j) = ( sflx_lw_up(i,j) - sflx_lw_dn(i,j) ) &
2503  + ( sflx_sw_up(i,j) - sflx_sw_dn(i,j) )
2504 
2505  toaflx_rd_net(i,j) = ( toaflx_lw_up(i,j) - toaflx_lw_dn(i,j) ) &
2506  + ( toaflx_sw_up(i,j) - toaflx_sw_dn(i,j) )
2507  enddo
2508  enddo
2509 
2510  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
2511 !OCL XFILL
2512  do j = js, je
2513  do i = is, ie
2514  engflxt(i,j) = sflx_sh(i,j) + sflx_lh(i,j) &
2515  + sflx_rd_net(i,j) - toaflx_rd_net(i,j)
2516  enddo
2517  enddo
2518 
2519  call monit_put( ad_monit_id(i_engt), engt(:,:,:) )
2520  call monit_put( ad_monit_id(i_engp), engp(:,:,:) )
2521  call monit_put( ad_monit_id(i_engk), engk(:,:,:) )
2522  call monit_put( ad_monit_id(i_engi), engi(:,:,:) )
2523 
2524  call monit_put( ad_monit_id(i_engflxt), engflxt(:,:) )
2525 
2526 
2527  call monit_put( ad_monit_id(i_engsfc_sh), sflx_sh(:,:) )
2528  call monit_put( ad_monit_id(i_engsfc_lh), sflx_lh(:,:) )
2529  call monit_put( ad_monit_id(i_engsfc_rd), sflx_rd_net(:,:) )
2530  call monit_put( ad_monit_id(i_engtoa_rd), toaflx_rd_net(:,:) )
2531 
2532  call monit_put( ad_monit_id(i_engsfc_lw_up), sflx_lw_up(:,:) )
2533  call monit_put( ad_monit_id(i_engsfc_lw_dn), sflx_lw_dn(:,:) )
2534  call monit_put( ad_monit_id(i_engsfc_sw_up), sflx_sw_up(:,:) )
2535  call monit_put( ad_monit_id(i_engsfc_sw_dn), sflx_sw_dn(:,:) )
2536 
2537  call monit_put( ad_monit_id(i_engtoa_lw_up), toaflx_lw_up(:,:) )
2538  call monit_put( ad_monit_id(i_engtoa_lw_dn), toaflx_lw_dn(:,:) )
2539  call monit_put( ad_monit_id(i_engtoa_sw_up), toaflx_sw_up(:,:) )
2540  call monit_put( ad_monit_id(i_engtoa_sw_dn), toaflx_sw_dn(:,:) )
2541 
2542  if ( atmos_vars_checkrange ) then
2543 !OCL XFILL
2544  work(:,:,:,1) = w(:,:,:)
2545 !OCL XFILL
2546  work(:,:,:,2) = u(:,:,:)
2547 !OCL XFILL
2548  work(:,:,:,3) = v(:,:,:)
2549 
2550  wname(1) = "W"
2551  wname(2) = "U"
2552  wname(3) = "V"
2553 
2554  call stat_detail( work(:,:,:,:), wname(:) )
2555  endif
2556 
2557  if ( atmos_vars_checkcfl > 0.0_rp ) then
2558 !OCL XFILL
2559  work(:,:,:,:) = 0.0_rp
2560 
2561  do j = js, je
2562  do i = is, ie
2563  do k = ks, ke
2564  work(k,i,j,1) = 0.5_rp * abs(momz(k,i,j)) / ( dens(k+1,i,j) + dens(k,i,j) ) &
2565  * time_dtsec_atmos_dyn / ( real_cz(k+1,i,j) - real_cz(k,i,j) )
2566  work(k,i,j,2) = 0.5_rp * abs(momx(k,i,j)) / ( dens(k,i+1,j) + dens(k,i,j) ) &
2567  * time_dtsec_atmos_dyn * rfdx(i) * mapf(i,j,1,i_uy)
2568  work(k,i,j,3) = 0.5_rp * abs(momy(k,i,j)) / ( dens(k,i,j+1) + dens(k,i,j) ) &
2569  * time_dtsec_atmos_dyn * rfdy(j) * mapf(i,j,2,i_xv)
2570  enddo
2571  enddo
2572  enddo
2573 
2574  cflmax = maxval( work(:,:,:,:) )
2575  if ( cflmax > atmos_vars_checkcfl ) then
2576  if( io_l ) write(io_fid_log,*) "*** [ATMOS_vars_monitor] Courant number exceeded the upper limit. : ", cflmax
2577  write(*,*) "*** [ATMOS_vars_monitor] Courant number exceeded the upper limit. : ", cflmax, &
2578  ", rank = ", prc_myrank
2579 
2580  wname(1) = "Courant num. Z"
2581  wname(2) = "Courant num. X"
2582  wname(3) = "Courant num. Y"
2583 
2584  call stat_detail( work(:,:,:,:), wname(:), supress_globalcomm=.true. )
2585  endif
2586  endif
2587 
2588  return
2589  end subroutine atmos_vars_monitor
2590 
2591 end module mod_atmos_vars
integer, public qie
module ATMOS admin
integer, public imax
of computational cells: x
real(rp), dimension(:,:,:), allocatable, public dens_tp
integer, parameter, public i_rhot
Definition: scale_index.F90:35
module ATMOSPHERE / Adiabatic process
character(len=h_mid), dimension(:), allocatable, public aq_desc
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
logical, public atmos_sw_phy_cp
logical, public statistics_checktotal
calc&report variable totals to logfile?
module GTOOL_FILE
Definition: gtool_file.f90:17
real(rp), dimension(:,:,:), allocatable, target, public momz
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:59
integer, public je
end point of inner domain: y, local
module Atmosphere / Physics Cumulus
real(rp), dimension(:), allocatable, public grid_rcdy
reciprocal of center-dy
subroutine, public atmos_phy_sf_vars_restart_read
Read restart.
subroutine, public atmos_phy_tb_vars_restart_read
Read restart.
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
subroutine, public hist_reg(itemid, zinterp, item, desc, unit, ndim, xdim, ydim, zdim)
Register/Append variable to history file.
module ATMOSPHERE / Saturation adjustment
subroutine, public atmos_phy_ch_vars_setup
Setup.
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sflx_qtrc
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
logical, public atmos_sw_phy_rd
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:,:), allocatable, public momy_tp
real(dp), public time_nowsec
subday part of current time [sec]
Definition: scale_time.F90:68
real(rp), dimension(:), allocatable, public aq_cp
CP for each hydrometeors [J/kg/K].
subroutine, public atmos_vars_fillhalo(FILL_BND)
HALO Communication.
integer, public qqe
module Atmosphere / Physics Cloud Microphysics
integer, parameter, public i_momx
Definition: scale_index.F90:33
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public atmos_phy_cp_vars_restart_write
Write restart.
real(rp), dimension(:), allocatable, public grid_rcdx
reciprocal of center-dx
real(rp), dimension(:,:,:), allocatable, public pres
module Atmosphere / Dynamics
subroutine, public atmos_phy_tb_vars_restart_write
Write restart.
integer, parameter, public i_momz
Definition: scale_index.F90:32
integer, public qwe
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
real(rp), dimension(:,:,:), allocatable, target, public momx
subroutine, public atmos_vars_restart_write
Write restart of atmospheric variables.
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public rhot_tp
integer, public qa
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), public atmos_restart_check_criterion
integer, public qws
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
logical, public atmos_sw_phy_ae
real(rp), dimension(:), allocatable, public grid_rfdy
reciprocal of face-dy
real(rp), dimension(:,:,:), allocatable, target, public dens
subroutine, public atmos_dyn_vars_restart_write
Write restart.
integer, public qis
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_up
integer, parameter, public i_dens
Definition: scale_index.F90:31
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_rain
character(len=h_long), public atmos_restart_check_basename
integer, parameter, public i_momy
Definition: scale_index.F90:34
module FILE I/O (netcdf)
subroutine, public atmos_phy_sf_vars_setup
Setup.
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_sw_up
module Statistics
integer, public ieb
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_dn
logical, public atmos_restart_output
output restart file?
module Atmosphere / Physics Radiation
subroutine, public atmos_adiabat_cape(Kstr, DENS, TEMP, PRES, QTRC, CZ, FZ, CAPE, CIN, LCL, LFC, LNB)
Calc CAPE and CIN Type of parcel method: Pseudo-adiabatic ascend from lowermost layer of the model Re...
module grid index
real(rp), dimension(:,:,:), allocatable, target, public momz_avw
module ATMOSPHERIC Surface Variables
subroutine, public atmos_phy_sf_vars_restart_write
Write restart.
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc_avw
subroutine, public atmos_phy_cp_vars_setup
Setup.
real(rp), dimension(:,:,:), pointer, public momx_av
logical, public atmos_sw_phy_tb
real(rp), dimension(:,:,:,:), allocatable, public gtrans_mapf
map factor
module TRACER
module Index
Definition: scale_index.F90:14
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_up
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:,:,:), allocatable, public temp
character(len=h_mid), public atmos_restart_out_dtype
REAL4 or REAL8.
subroutine, public comm_horizontal_mean(varmean, var)
calculate horizontal mean (global total with communication)
Definition: scale_comm.F90:516
module GRIDTRANS
module GRID (real space)
subroutine, public time_gettimelabel(timelabel)
generate time label
Definition: scale_time.F90:90
subroutine, public monit_reg(itemid, item, desc, unit, ndim, isflux)
Search existing item, or matching check between requested and registered item.
integer, public ka
of z whole cells (local, with HALO)
integer, public i_uy
integer, public i_qv
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:93
module MONITOR
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
integer, public kmax
of computational cells: z
subroutine, public atmos_vars_setup
Setup.
logical, public atmos_sw_dyn
subroutine, public atmos_vars_total
Budget monitor for atmosphere.
character(len=h_short), dimension(:), allocatable, public aq_name
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_lh
module COMMUNICATION
Definition: scale_comm.F90:23
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
integer, public js
start point of inner domain: y, local
subroutine, public atmos_phy_tb_vars_setup
Setup.
module TIME
Definition: scale_time.F90:15
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_lw_dn
subroutine, public atmos_dyn_vars_restart_read
Read restart.
real(dp), public time_dtsec_atmos_dyn
time interval of dynamics [sec]
Definition: scale_time.F90:38
module PROCESS
logical, public atmos_sw_phy_sf
subroutine, public atmos_phy_rd_vars_setup
Setup.
logical, public atmos_sw_phy_ch
module Atmosphere / Physics Turbulence
subroutine, public stat_detail(var, varname, supress_globalcomm)
Search global maximum & minimum value.
real(rp), dimension(:,:,:), allocatable, target, public momx_avw
real(rp), dimension(:,:,:), pointer, public dens_av
real(rp), public const_lhv
latent heat of vaporizaion for use
Definition: scale_const.F90:77
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
logical, public atmos_sw_phy_mp
real(rp), dimension(:), allocatable, public aq_cv
CV for each hydrometeors [J/kg/K].
character(len=h_short), dimension(:), allocatable, public aq_unit
subroutine, public atmos_phy_rd_vars_restart_write
Write restart.
integer, public ks
start point of inner domain: z, local
subroutine, public atmos_dyn_vars_setup
Setup.
subroutine, public atmos_phy_ae_vars_setup
Setup.
subroutine, public atmos_phy_ae_vars_restart_read
Read restart.
integer, public prc_myrank
process num in local communicator
real(rp), dimension(:,:,:), allocatable, public momx_tp
real(rp), dimension(:,:,:), allocatable, target, public momy
module GRID (cartesian)
module Atmosphere / Physics Chemistry
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
real(rp), public const_lhf
latent heat of fusion for use
Definition: scale_const.F90:79
subroutine, public atmos_vars_monitor
monitor output
real(rp), dimension(:,:,:), allocatable, public pott
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_snow
real(rp), dimension(:,:,:), allocatable, public momz_tp
integer, public i_xv
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:,:), allocatable, target, public dens_avw
real(rp), dimension(:,:,:), allocatable, public v
subroutine, public atmos_phy_mp_vars_setup
Setup.
module profiler
Definition: scale_prof.F90:10
character(len=h_mid), public atmos_restart_out_title
title of the output file
integer, public ie
end point of inner domain: x, local
logical, public atmos_restart_check
check value consistency?
subroutine, public atmos_phy_mp_vars_restart_write
Write restart.
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public atmos_phy_mp_vars_restart_read
Read restart.
real(rp), dimension(:,:,:), allocatable, public w
module ATMOSPHERE / Thermodynamics
integer, public qqs
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_sw_dn
character(len=h_long), public atmos_restart_in_basename
basename of the restart file
module PRECISION
character(len=h_long), public atmos_restart_out_basename
basename of the output file
subroutine, public atmos_phy_ae_vars_restart_write
Write restart.
subroutine, public atmos_phy_rd_vars_restart_read
Read restart.
real(rp), dimension(:,:,:), allocatable, public u
module HISTORY
real(rp), dimension(:,:,:), pointer, public momz_av
real(rp), dimension(:), allocatable, public grid_rfdx
reciprocal of face-dx
integer, public isb
subroutine, public atmos_phy_ch_vars_restart_write
Write restart.
subroutine, public atmos_phy_cp_vars_restart_read
Read restart.
real(rp), dimension(:,:,:), allocatable, target, public rhot_avw
real(rp), dimension(:,:,:), pointer, public rhot_av
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
subroutine, public atmos_vars_history
History output set for atmospheric variables.
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), dimension(:,:,:), pointer, public momy_av
integer, public jsb
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
subroutine, public atmos_vars_restart_check
Check and compare between last data and sample data.
real(rp), dimension(:,:,:), allocatable, target, public momy_avw
subroutine, public atmos_phy_ch_vars_restart_read
Read restart.
module ATMOSPHERE / Physics Aerosol Microphysics
logical, public atmos_restart_in_allowmissingq
integer, public jmax
of computational cells: y
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_lw_up
subroutine, public atmos_vars_restart_read
Read restart of atmospheric variables.
logical, public atmos_use_average
subroutine, public atmos_vars_diagnostics
Calc diagnostic variables.
integer, public i_qc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
integer, public ja
of y whole cells (local, with HALO)