SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tstep_large_fvm_heve Module Reference

module Atmosphere / Dynamics More...

Functions/Subroutines

subroutine, public atmos_dyn_tstep_large_fvm_heve_setup (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, mflx_hi)
 Setup. More...
 
subroutine, public atmos_dyn_tstep_large_fvm_heve (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, mflx_hi, tflx_hi, num_diff, num_diff_q, QTRC0, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLI, CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, AQ_CV, REF_dens, REF_pott, REF_qv, REF_pres, BND_W, BND_E, BND_S, BND_N, ND_COEF, ND_COEF_Q, ND_ORDER, ND_SFC_FACT, ND_USE_RS, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, divdmp_coef, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, DTLS, DTSS, Llast)
 Dynamical Process. More...
 
subroutine check_mass (DENS, DAMP_DENS, mflx_hi, tflx_hi, GSQRT, MAPF, RCDX, RCDY, dt, step, BND_W, BND_E, BND_S, BND_N)
 

Detailed Description

module Atmosphere / Dynamics

Description
HEVE FVM scheme for large time step in Atmospheric dynamical process
Author
Team SCALE
History
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
ALLMOM_lb_hz horizontally total momentum flux from lateral boundary kg/m2/s allmflx_lb_horizontal
DENS_t_damp tendency of dencity due to rayleigh damping kg/m3/s damp_t
DENS_t_phys tendency of dencity due to physics kg/m3/s DENS_tp
MFLXX momentum flux of x-direction kg/m2/s mflx_hi
MFLXY momentum flux of y-direction kg/m2/s mflx_hi
MFLXZ momentum flux of z-direction kg/m2/s mflx_hi
MOMX_t_damp tendency of momentum x due to rayleigh damping kg/m2/s2 damp_t
MOMX_t_phys tendency of momentum x due to physics kg/m2/s2 MOMX_tp
MOMY_t_damp tendency of momentum y due to rayleigh damping kg/m2/s2 damp_t
MOMY_t_phys tendency of momentum y due to physics kg/m2/s2 MOMY_tp
MOMZ_t_damp tendency of momentum z due to rayleigh damping kg/m2/s2 damp_t
MOMZ_t_phys tendency of momentum z due to physics kg/m2/s2 MOMZ_tp
RHOT_t_damp tendency of rho*theta temperature due to rayleigh damping K kg/m3/s damp_t
RHOT_t_phys tendency of rho*theta temperature due to physics K kg/m3/s RHOT_tp
TFLXX potential temperature flux of x-direction K*kg/m2/s tflx_hi
TFLXY potential temperature flux of y-direction K*kg/m2/s tflx_hi
TFLXZ potential temperature flux of z-direction K*kg/m2/s tflx_hi
trim(AQ_NAME(iq))//'_t_damp' tendency of '//trim(AQ_NAME(iq))//' due to rayleigh damping kg/kg/s damp_t
trim(AQ_NAME(iq))//'_t_phys' tendency of '//trim(AQ_NAME(iq))//' due to physics kg/kg/s RHOQ_tp

Function/Subroutine Documentation

◆ atmos_dyn_tstep_large_fvm_heve_setup()

subroutine, public scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve_setup ( real(rp), dimension(ka,ia,ja), intent(inout)  DENS,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMX,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMY,
real(rp), dimension(ka,ia,ja), intent(inout)  RHOT,
real(rp), dimension(ka,ia,ja,qa), intent(inout)  QTRC,
real(rp), dimension(ka,ia,ja,va), intent(inout)  PROG,
real(rp), dimension(ka,ia,ja,3), intent(inout)  mflx_hi 
)

Setup.

Definition at line 95 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

References scale_comm::comm_vars8_init(), scale_const::const_ohm, scale_const::const_undef, scale_grid_index::ia, scale_grid_index::ja, scale_grid_index::ka, scale_rm_process::prc_has_e, scale_rm_process::prc_has_n, scale_rm_process::prc_has_s, scale_rm_process::prc_has_w, scale_process::prc_mpistop(), scale_tracer::qa, scale_index::va, scale_grid_index::xdir, scale_grid_index::ydir, and scale_grid_index::zdir.

Referenced by scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup().

95  use scale_process, only: &
97  use scale_rm_process, only: &
98  prc_has_e, &
99  prc_has_w, &
100  prc_has_n, &
101  prc_has_s
102  use scale_const, only: &
103  ohm => const_ohm, &
104  undef => const_undef
105  use scale_comm, only: &
107  implicit none
108 
109  ! MPI_RECV_INIT requires intent(inout)
110  real(RP), intent(inout) :: dens(ka,ia,ja)
111  real(RP), intent(inout) :: momz(ka,ia,ja)
112  real(RP), intent(inout) :: momx(ka,ia,ja)
113  real(RP), intent(inout) :: momy(ka,ia,ja)
114  real(RP), intent(inout) :: rhot(ka,ia,ja)
115  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
116  real(RP), intent(inout) :: prog(ka,ia,ja,va)
117  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3)
118 
119  integer :: iv, iq
120  !---------------------------------------------------------------------------
121 
122  allocate( dens_t(ka,ia,ja) )
123  allocate( momz_t(ka,ia,ja) )
124  allocate( momx_t(ka,ia,ja) )
125  allocate( momy_t(ka,ia,ja) )
126  allocate( rhot_t(ka,ia,ja) )
127  allocate( rhoq_t(ka,ia,ja,qa) )
128 
129  allocate( i_comm_prog(max(va,1)) )
130  allocate( i_comm_qtrc(qa) )
131  allocate( i_comm_rhoq_t(qa) )
132 
133  call comm_vars8_init( dens, i_comm_dens )
134  call comm_vars8_init( momz, i_comm_momz )
135  call comm_vars8_init( momx, i_comm_momx )
136  call comm_vars8_init( momy, i_comm_momy )
137  call comm_vars8_init( rhot, i_comm_rhot )
138  do iv = 1, va
139  i_comm_prog(iv) = 5 + iv
140  call comm_vars8_init( prog(:,:,:,iv), i_comm_prog(iv) )
141  end do
142 
143  call comm_vars8_init( dens_t, i_comm_dens_t )
144  call comm_vars8_init( momz_t, i_comm_momz_t )
145  call comm_vars8_init( momx_t, i_comm_momx_t )
146  call comm_vars8_init( momy_t, i_comm_momy_t )
147  call comm_vars8_init( rhot_t, i_comm_rhot_t )
148 
149  do iq = 1, qa
150  i_comm_rhoq_t(iq) = 5 + va + iq
151  i_comm_qtrc(iq) = 5 + va + iq
152 
153  call comm_vars8_init( rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq) )
154  call comm_vars8_init( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
155  end do
156 
157  call comm_vars8_init( mflx_hi(:,:,:,zdir), i_comm_mflx_z )
158  call comm_vars8_init( mflx_hi(:,:,:,xdir), i_comm_mflx_x )
159  call comm_vars8_init( mflx_hi(:,:,:,ydir), i_comm_mflx_y )
160 
161  mflx_hi(:,:,:,:) = undef
162 
163  return
logical, public prc_has_n
subroutine, public prc_mpistop
Abort MPI.
integer, public va
Definition: scale_index.F90:38
integer, parameter, public zdir
logical, public prc_has_e
integer, parameter, public ydir
integer, parameter, public xdir
logical, public prc_has_s
real(rp), public const_undef
Definition: scale_const.F90:43
real(rp), public const_ohm
angular velocity of the planet [1/s]
Definition: scale_const.F90:47
integer, public ia
of x whole cells (local, with HALO)
subroutine, public comm_vars8_init(var, vid)
Register variables.
Definition: scale_comm.F90:328
integer, public ka
of z whole cells (local, with HALO)
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module RM PROCESS
logical, public prc_has_w
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tstep_large_fvm_heve()

subroutine, public scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve ( real(rp), dimension(ka,ia,ja), intent(inout)  DENS,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMX,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMY,
real(rp), dimension(ka,ia,ja), intent(inout)  RHOT,
real(rp), dimension(ka,ia,ja,qa), intent(inout)  QTRC,
real(rp), dimension(ka,ia,ja,va), intent(inout)  PROG,
real(rp), dimension(ka,ia,ja), intent(inout)  DENS_av,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMZ_av,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMX_av,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMY_av,
real(rp), dimension(ka,ia,ja), intent(inout)  RHOT_av,
real(rp), dimension(ka,ia,ja,qa), intent(inout)  QTRC_av,
real(rp), dimension(ka,ia,ja,3), intent(out)  mflx_hi,
real(rp), dimension(ka,ia,ja,3), intent(out)  tflx_hi,
real(rp), dimension(ka,ia,ja,5,3), intent(out)  num_diff,
real(rp), dimension(ka,ia,ja,3), intent(out)  num_diff_q,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC0,
real(rp), dimension(ka,ia,ja), intent(in)  DENS_tp,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ_tp,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX_tp,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY_tp,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT_tp,
real(rp), dimension(ka,ia,ja,qa), intent(in)  RHOQ_tp,
real(rp), dimension(ia,ja), intent(in)  CORIOLI,
real(rp), dimension (ka), intent(in)  CDZ,
real(rp), dimension (ia), intent(in)  CDX,
real(rp), dimension (ja), intent(in)  CDY,
real(rp), dimension (ka-1), intent(in)  FDZ,
real(rp), dimension (ia-1), intent(in)  FDX,
real(rp), dimension (ja-1), intent(in)  FDY,
real(rp), dimension(ka), intent(in)  RCDZ,
real(rp), dimension(ia), intent(in)  RCDX,
real(rp), dimension(ja), intent(in)  RCDY,
real(rp), dimension(ka-1), intent(in)  RFDZ,
real(rp), dimension(ia-1), intent(in)  RFDX,
real(rp), dimension(ja-1), intent(in)  RFDY,
real(rp), dimension (ka,ia,ja), intent(in)  PHI,
real(rp), dimension(ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension (ka,ia,ja,7), intent(in)  J13G,
real(rp), dimension (ka,ia,ja,7), intent(in)  J23G,
real(rp), intent(in)  J33G,
real(rp), dimension (ia,ja,2,4), intent(in)  MAPF,
real(rp), dimension(qqa), intent(in)  AQ_CV,
real(rp), dimension(ka,ia,ja), intent(in)  REF_dens,
real(rp), dimension(ka,ia,ja), intent(in)  REF_pott,
real(rp), dimension (ka,ia,ja), intent(in)  REF_qv,
real(rp), dimension(ka,ia,ja), intent(in)  REF_pres,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
real(rp), intent(in)  ND_COEF,
real(rp), intent(in)  ND_COEF_Q,
integer, intent(in)  ND_ORDER,
real(rp), intent(in)  ND_SFC_FACT,
logical, intent(in)  ND_USE_RS,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_DENS,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_VELZ,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_VELX,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_VELY,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_POTT,
real(rp), dimension(ka,ia,ja,bnd_qa), intent(in)  DAMP_QTRC,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_alpha_DENS,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_alpha_VELZ,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_alpha_VELX,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_alpha_VELY,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_alpha_POTT,
real(rp), dimension(ka,ia,ja,bnd_qa), intent(in)  DAMP_alpha_QTRC,
real(rp), intent(in)  divdmp_coef,
logical, intent(in)  FLAG_FCT_MOMENTUM,
logical, intent(in)  FLAG_FCT_T,
logical, intent(in)  FLAG_FCT_TRACER,
logical, intent(in)  FLAG_FCT_ALONG_STREAM,
logical, intent(in)  USE_AVERAGE,
real(dp), intent(in)  DTLS,
real(dp), intent(in)  DTSS,
logical, intent(in)  Llast 
)

Dynamical Process.

Parameters
[in]phigeopotential
[in]gsqrtvertical metrics {G}^1/2
[in]j13g(1,3) element of Jacobian matrix
[in]j23g(2,3) element of Jacobian matrix
[in]j33g(3,3) element of Jacobian matrix
[in]mapfmap factor
[in]ref_presreference pressure

Definition at line 193 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

References scale_tracer::aq_name, scale_atmos_boundary::atmos_boundary_smoother_fact, scale_atmos_dyn_common::atmos_dyn_divergence(), scale_atmos_dyn_common::atmos_dyn_fct(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xyz, scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxx_xyz_ud1(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xyz, scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxy_xyz_ud1(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xyz, scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_xyz_ud1(), scale_atmos_dyn_common::atmos_dyn_numfilter_coef(), scale_atmos_dyn_common::atmos_dyn_numfilter_coef_q(), scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short, scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer, scale_atmos_boundary::bnd_qa, check_mass(), scale_const::const_cvdry, scale_const::const_eps, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_rvap, scale_tracer::i_qv, scale_gridtrans::i_uy, scale_gridtrans::i_uyz, scale_gridtrans::i_xv, scale_gridtrans::i_xvz, scale_gridtrans::i_xy, scale_gridtrans::i_xyw, scale_gridtrans::i_xyz, scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, scale_tracer::qqe, scale_tracer::qqs, scale_index::va, scale_grid_index::xdir, scale_grid_index::ydir, and scale_grid_index::zdir.

Referenced by scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup().

193  use scale_const, only: &
194  eps => const_eps, &
195  p0 => const_pre00, &
196  rdry => const_rdry, &
197  rvap => const_rvap, &
198  cvdry => const_cvdry
199  use scale_comm, only: &
200  comm_vars8, &
201  comm_wait
202  use scale_gridtrans, only: &
203  i_xyz, &
204  i_xyw, &
205  i_uyz, &
206  i_xvz, &
207  i_xy, &
208  i_uy, &
209  i_xv
210  use scale_atmos_dyn_common, only: &
215  use scale_atmos_dyn_fvm_flux_ud1, only: &
219  use scale_atmos_dyn_fvm_flux, only: &
223  use scale_atmos_boundary, only: &
224  bnd_qa, &
225  bnd_smoother_fact => atmos_boundary_smoother_fact
226  use scale_history, only: &
227 #ifdef hist_tend
228  hist_in, &
229 #endif
230  hist_switch
231  use scale_atmos_dyn_tinteg_short, only: &
233  use scale_atmos_dyn_tinteg_tracer, only: &
235  implicit none
236 
237  real(RP), intent(inout) :: dens(ka,ia,ja)
238  real(RP), intent(inout) :: momz(ka,ia,ja)
239  real(RP), intent(inout) :: momx(ka,ia,ja)
240  real(RP), intent(inout) :: momy(ka,ia,ja)
241  real(RP), intent(inout) :: rhot(ka,ia,ja)
242  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
243  real(RP), intent(inout) :: prog(ka,ia,ja,va)
244 
245  real(RP), intent(inout) :: dens_av(ka,ia,ja)
246  real(RP), intent(inout) :: momz_av(ka,ia,ja)
247  real(RP), intent(inout) :: momx_av(ka,ia,ja)
248  real(RP), intent(inout) :: momy_av(ka,ia,ja)
249  real(RP), intent(inout) :: rhot_av(ka,ia,ja)
250  real(RP), intent(inout) :: qtrc_av(ka,ia,ja,qa)
251 
252  real(RP), intent(out) :: mflx_hi(ka,ia,ja,3)
253  real(RP), intent(out) :: tflx_hi(ka,ia,ja,3)
254 
255  real(RP), intent(out) :: num_diff(ka,ia,ja,5,3)
256  real(RP), intent(out) :: num_diff_q(ka,ia,ja,3)
257 
258  real(RP), intent(in) :: qtrc0(ka,ia,ja,qa)
259 
260  real(RP), intent(in) :: dens_tp(ka,ia,ja)
261  real(RP), intent(in) :: momz_tp(ka,ia,ja)
262  real(RP), intent(in) :: momx_tp(ka,ia,ja)
263  real(RP), intent(in) :: momy_tp(ka,ia,ja)
264  real(RP), intent(in) :: rhot_tp(ka,ia,ja)
265  real(RP), intent(in) :: rhoq_tp(ka,ia,ja,qa)
266 
267  real(RP), intent(in) :: corioli(ia,ja)
268 
269  real(RP), intent(in) :: cdz (ka)
270  real(RP), intent(in) :: cdx (ia)
271  real(RP), intent(in) :: cdy (ja)
272  real(RP), intent(in) :: fdz (ka-1)
273  real(RP), intent(in) :: fdx (ia-1)
274  real(RP), intent(in) :: fdy (ja-1)
275  real(RP), intent(in) :: rcdz(ka)
276  real(RP), intent(in) :: rcdx(ia)
277  real(RP), intent(in) :: rcdy(ja)
278  real(RP), intent(in) :: rfdz(ka-1)
279  real(RP), intent(in) :: rfdx(ia-1)
280  real(RP), intent(in) :: rfdy(ja-1)
281 
282  real(RP), intent(in) :: phi (ka,ia,ja)
283  real(RP), intent(in) :: gsqrt(ka,ia,ja,7)
284  real(RP), intent(in) :: j13g (ka,ia,ja,7)
285  real(RP), intent(in) :: j23g (ka,ia,ja,7)
286  real(RP), intent(in) :: j33g
287  real(RP), intent(in) :: mapf (ia,ja,2,4)
288 
289  real(RP), intent(in) :: aq_cv(qqa)
290 
291  real(RP), intent(in) :: ref_dens(ka,ia,ja)
292  real(RP), intent(in) :: ref_pott(ka,ia,ja)
293  real(RP), intent(in) :: ref_qv (ka,ia,ja)
294  real(RP), intent(in) :: ref_pres(ka,ia,ja)
295 
296  logical, intent(in) :: bnd_w
297  logical, intent(in) :: bnd_e
298  logical, intent(in) :: bnd_s
299  logical, intent(in) :: bnd_n
300 
301  real(RP), intent(in) :: nd_coef
302  real(RP), intent(in) :: nd_coef_q
303  integer, intent(in) :: nd_order
304  real(RP), intent(in) :: nd_sfc_fact
305  logical, intent(in) :: nd_use_rs
306 
307  real(RP), intent(in) :: damp_dens(ka,ia,ja)
308  real(RP), intent(in) :: damp_velz(ka,ia,ja)
309  real(RP), intent(in) :: damp_velx(ka,ia,ja)
310  real(RP), intent(in) :: damp_vely(ka,ia,ja)
311  real(RP), intent(in) :: damp_pott(ka,ia,ja)
312  real(RP), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
313 
314  real(RP), intent(in) :: damp_alpha_dens(ka,ia,ja)
315  real(RP), intent(in) :: damp_alpha_velz(ka,ia,ja)
316  real(RP), intent(in) :: damp_alpha_velx(ka,ia,ja)
317  real(RP), intent(in) :: damp_alpha_vely(ka,ia,ja)
318  real(RP), intent(in) :: damp_alpha_pott(ka,ia,ja)
319  real(RP), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
320 
321  real(RP), intent(in) :: divdmp_coef
322 
323  logical, intent(in) :: flag_fct_momentum
324  logical, intent(in) :: flag_fct_t
325  logical, intent(in) :: flag_fct_tracer
326  logical, intent(in) :: flag_fct_along_stream
327 
328  logical, intent(in) :: use_average
329 
330  real(DP), intent(in) :: dtls
331  real(DP), intent(in) :: dtss
332  logical , intent(in) :: llast
333 
334  ! for time integartion
335  real(RP) :: dens00 (ka,ia,ja) ! saved density before small step loop
336 
337  ! diagnostic variables
338  real(RP) :: ddiv (ka,ia,ja) ! 3 dimensional divergence
339  real(RP) :: dpres0 (ka,ia,ja) ! pressure deviation
340  real(RP) :: rt2p (ka,ia,ja) ! factor of RHOT to PRES
341  real(RP) :: ref_rhot(ka,ia,ja) ! reference of RHOT
342 
343  real(RP) :: qdry ! dry air
344  real(RP) :: rtot ! total R
345  real(RP) :: cvtot ! total CV
346  real(RP) :: cptot ! total CP
347  real(RP) :: pres ! pressure
348 
349  real(RP) :: dens_tq(ka,ia,ja)
350  real(RP) :: diff(ka,ia,ja)
351  real(RP) :: damp
352 #ifdef HIST_TEND
353  real(RP) :: damp_t(ka,ia,ja)
354 #endif
355 
356  ! For tracer advection
357  real(RP) :: mflx_av (ka,ia,ja,3) ! rho * vel(x,y,z) @ (u,v,w)-face average
358  real(RP) :: qflx_hi (ka,ia,ja,3) ! rho * vel(x,y,z) * phi @ (u,v,w)-face high order
359  real(RP) :: qflx_lo (ka,ia,ja,3) ! rho * vel(x,y,z) * phi, monotone flux
360  real(RP) :: qflx_anti(ka,ia,ja,3) ! anti-diffusive flux
361 
362  real(RP) :: dtl
363  real(RP) :: dts
364  integer :: nstep
365 
366  integer :: iis, iie
367  integer :: jjs, jje
368  integer :: i, j, k, iq, step
369  integer :: iv
370 
371  real(RP) :: diff_coef
372  !---------------------------------------------------------------------------
373 
374  call prof_rapstart("DYN_Large_Preparation", 2)
375 
376  dts = real(DTSS, kind=RP) ! short time step
377  dtl = real(DTLS, kind=RP) ! large time step
378  nstep = ceiling( ( dtl - eps ) / dts )
379  dts = dtl / nstep ! dts is divisor of dtl and smaller or equal to dtss
380 
381 #ifdef DEBUG
382  if( io_l ) write(io_fid_log,*) '*** Dynamics large time step'
383  if( io_l ) write(io_fid_log,'(1x,A,F0.2,A,F0.2,A,I0)') &
384  '*** -> DT_large, DT_small, DT_large/DT_small : ', dtl, ', ', dts, ', ', nstep
385 
386  dens00(:,:,:) = undef
387 
388  num_diff(:,:,:,:,:) = undef
389 
390  mflx_hi(:,:,:,:) = undef
391  tflx_hi(:,:,:,:) = undef
392 
393  qflx_hi(:,:,:,:) = undef
394  qflx_lo(:,:,:,:) = undef
395 #endif
396 
397 !OCL XFILL
398  dens00(:,:,:) = dens(:,:,:)
399 
400  if ( use_average ) then
401 !OCL XFILL
402  dens_av(:,:,:) = 0.0_rp
403 !OCL XFILL
404  momz_av(:,:,:) = 0.0_rp
405 !OCL XFILL
406  momx_av(:,:,:) = 0.0_rp
407 !OCL XFILL
408  momy_av(:,:,:) = 0.0_rp
409 !OCL XFILL
410  rhot_av(:,:,:) = 0.0_rp
411  endif
412 
413 #ifndef DRY
414 !OCL XFILL
415  mflx_av(:,:,:,:) = 0.0_rp
416 
417 #endif
418 
419  !------------------------------------------------------------------------
420  ! prepare thermodynamical data
421  ! specific heat
422  ! pressure data ( linearization )
423  !
424  ! pres = P0 * ( R * rhot / P0 )**(CP/CV)
425  ! d pres / d rhot ~ CP*R/CV * ( R * rhot / P0 )**(R/CV)
426  ! = CP*R/CV * ( pres / P0 )**(R/CP)
427  ! = CP*R/CV * temp / pott
428  ! = CP/CV * pres / rhot
429  ! pres ~ P0 * ( R * rhot0 / P0 ) ** (CP/CV) + CV*R/CP * ( pres / P0 )**(R/CP) * rhot'
430  !------------------------------------------------------------------------
431 !OCL XFILL
432  do j = 1, ja
433  do i = 1, ia
434  do k = ks, ke
435 #ifdef DRY
436  pres = p0 * ( rdry * rhot(k,i,j) / p0 )**cpovcv
437  rt2p(k,i,j) = cpovcv * pres / rhot(k,i,j)
438 #else
439  cvtot = 0.0_rp
440  qdry = 1.0_rp
441  do iq = qqs, qqe
442  cvtot = cvtot + aq_cv(iq) * qtrc(k,i,j,iq)
443  qdry = qdry - qtrc(k,i,j,iq)
444  enddo
445  cvtot = cvdry * qdry + cvtot
446  rtot = rdry * qdry + rvap * qtrc(k,i,j,i_qv)
447  cptot = cvtot + rtot
448  pres = p0 * ( rtot * rhot(k,i,j) / p0 )**( cptot / cvtot )
449  rt2p(k,i,j) = cptot / cvtot * pres / rhot(k,i,j)
450 #endif
451  dpres0(k,i,j) = pres - ref_pres(k,i,j)
452  ref_rhot(k,i,j) = rhot(k,i,j)
453  end do
454  dpres0(ks-1,i,j) = dpres0(ks+1,i,j) + ( ref_pres(ks+1,i,j) - ref_pres(ks-1,i,j) )
455  dpres0(ke+1,i,j) = dpres0(ke-1,i,j) + ( ref_pres(ke-1,i,j) - ref_pres(ke+1,i,j) )
456  end do
457  end do
458 
459  call prof_rapend ("DYN_Large_Preparation", 2)
460 
461  !###########################################################################
462  ! Update DENS,MONZ,MOMX,MOMY,MOMZ,RHOT
463  !###########################################################################
464 
465  call prof_rapstart("DYN_Large_Tendency", 2)
466 
467 #ifdef HIST_TEND
468 !OCL XFILL
469  damp_t(:,:,:) = 0.0_rp
470 #endif
471 
472 !OCL XFILL
473  dens_tq(:,:,:) = 0.0_rp
474 
475  do iq = 1, bnd_qa
476 
477  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
478 !OCL XFILL
479  do j = js-1, je+2
480  do i = is-1, ie+2
481  do k = ks, ke
482  diff(k,i,j) = qtrc(k,i,j,iq) - damp_qtrc(k,i,j,iq)
483  enddo
484  enddo
485  enddo
486  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
487 !OCL XFILL
488  do j = js, je
489  do i = is, ie
490  do k = ks, ke
491  damp = - damp_alpha_qtrc(k,i,j,iq) &
492  * ( diff(k,i,j) & ! rayleigh damping
493  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
494  * 0.125_rp * bnd_smoother_fact ) ! horizontal smoother
495  rhoq_t(k,i,j,iq) = rhoq_tp(k,i,j,iq) + damp * dens00(k,i,j)
496 #ifdef HIST_TEND
497  damp_t(k,i,j) = damp
498 #endif
499  enddo
500  enddo
501  enddo
502 #ifdef HIST_TEND
503  call hist_in(rhoq_tp(:,:,:,iq), trim(aq_name(iq))//'_t_phys', &
504  'tendency of '//trim(aq_name(iq))//' due to physics', 'kg/kg/s' )
505  call hist_in(damp_t, trim(aq_name(iq))//'_t_damp', &
506  'tendency of '//trim(aq_name(iq))//' due to rayleigh damping', 'kg/kg/s' )
507 #endif
508 !OCL XFILL
509  do j = js, je
510  do i = is, ie
511  rhoq_t( 1:ks-1,i,j,iq) = 0.0_rp
512  rhoq_t(ke+1:ka ,i,j,iq) = 0.0_rp
513  enddo
514  enddo
515 
516  call comm_vars8( rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq) )
517 
518  do j = js, je
519  do i = is, ie
520  do k = ks, ke
521  dens_tq(k,i,j) = dens_tq(k,i,j) + rhoq_t(k,i,j,iq)
522  end do
523  end do
524  end do
525 
526  call comm_wait ( rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq), .false. )
527 
528  end do
529 
530  !$omp parallel do private(i,j,k,iq) OMP_SCHEDULE_ collapse(3)
531 !OCL XFILL
532  do iq = bnd_qa+1, qa
533  do j = 1, ja
534  do i = 1, ia
535  do k = 1, ka
536  rhoq_t(k,i,j,iq) = rhoq_tp(k,i,j,iq)
537  enddo
538  enddo
539  enddo
540  end do
541 
542  call prof_rapend ("DYN_Large_Tendency", 2)
543 
544  call prof_rapstart("DYN_Large_Boundary", 2)
545 
546  if ( bnd_w ) then
547  !$omp parallel do private(j,k) OMP_SCHEDULE_ collapse(2)
548  do j = js, je
549  do k = ks, ke
550  mflx_hi(k,is-1,j,xdir) = gsqrt(k,is-1,j,i_uyz) * momx(k,is-1,j) / mapf(is-1,j,2,i_uy)
551  enddo
552  enddo
553  end if
554  if ( bnd_e ) then
555  !$omp parallel do private(j,k) OMP_SCHEDULE_ collapse(2)
556  do j = js, je
557  do k = ks, ke
558  mflx_hi(k,ie,j,xdir) = gsqrt(k,ie,j,i_uyz) * momx(k,ie,j) / mapf(ie,j,2,i_uy)
559  enddo
560  enddo
561  end if
562  if ( bnd_s ) then
563  !$omp parallel do private(i,k) OMP_SCHEDULE_ collapse(2)
564  do i = is, ie
565  do k = ks, ke
566  mflx_hi(k,i,js-1,ydir) = gsqrt(k,i,js-1,i_xvz) * momy(k,i,js-1) / mapf(i,js-1,1,i_xv)
567  enddo
568  enddo
569  end if
570  if ( bnd_n ) then
571  !$omp parallel do private(i,k) OMP_SCHEDULE_ collapse(2)
572  do i = is, ie
573  do k = ks, ke
574  mflx_hi(k,i,je,ydir) = gsqrt(k,i,je,i_xvz) * momy(k,i,je) / mapf(i,je,1,i_xv)
575  enddo
576  enddo
577  end if
578 
579  call prof_rapend ("DYN_Large_Boundary", 2)
580 
581 
582  do step = 1, nstep
583 
584  call hist_switch( llast .AND. step == nstep )
585 
586  !-----< prepare tendency >-----
587 
588  call prof_rapstart("DYN_Large_Tendency", 2)
589 
590  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
591 !OCL XFILL
592  do j = js-1, je+1
593  do i = is-1, ie+1
594  do k = ks, ke
595  diff(k,i,j) = dens(k,i,j) - damp_dens(k,i,j)
596  enddo
597  enddo
598  enddo
599  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
600 !OCL XFILL
601  do j = js, je
602  do i = is, ie
603  do k = ks, ke
604  damp = - damp_alpha_dens(k,i,j) &
605  * ( diff(k,i,j) & ! rayleigh damping
606  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
607  * 0.125_rp * bnd_smoother_fact ) & ! horizontal smoother
608  + dens_tq(k,i,j) ! dencity change due to rayleigh damping for tracers
609  dens_t(k,i,j) = dens_tp(k,i,j) & ! tendency from physical step
610  + damp
611 #ifdef HIST_TEND
612  damp_t(k,i,j) = damp
613 #endif
614  enddo
615  enddo
616  enddo
617 !OCL XFILL
618  do j = js, je
619  do i = is, ie
620  dens_t( 1:ks-1,i,j) = 0.0_rp
621  dens_t(ke+1:ka ,i,j) = 0.0_rp
622  enddo
623  enddo
624  call comm_vars8( dens_t(:,:,:), i_comm_dens_t )
625 #ifdef HIST_TEND
626  call hist_in(dens_tp, 'DENS_t_phys', 'tendency of dencity due to physics', 'kg/m3/s' )
627  call hist_in(damp_t, 'DENS_t_damp', 'tendency of dencity due to rayleigh damping', 'kg/m3/s' )
628 #endif
629 
630  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
631 !OCL XFILL
632  do j = js-1, je+1
633  do i = is-1, ie+1
634  do k = ks, ke-1
635  diff(k,i,j) = momz(k,i,j) - damp_velz(k,i,j) * ( dens(k,i,j)+dens(k+1,i,j) ) * 0.5_rp
636  enddo
637  enddo
638  enddo
639  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
640 !OCL XFILL
641  do j = js, je
642  do i = is, ie
643  do k = ks, ke-1
644  damp = - damp_alpha_velz(k,i,j) &
645  * ( diff(k,i,j) & ! rayleigh damping
646  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
647  * 0.125_rp * bnd_smoother_fact ) ! horizontal smoother
648  momz_t(k,i,j) = momz_tp(k,i,j) & ! tendency from physical step
649  + damp
650 #ifdef HIST_TEND
651  damp_t(k,i,j) = damp
652 #endif
653  enddo
654  enddo
655  enddo
656 !OCL XFILL
657  do j = js, je
658  do i = is, ie
659  momz_t( 1:ks-1,i,j) = 0.0_rp
660  momz_t(ke:ka ,i,j) = 0.0_rp
661  enddo
662  enddo
663  call comm_vars8( momz_t(:,:,:), i_comm_momz_t )
664 #ifdef HIST_TEND
665  call hist_in(momz_tp, 'MOMZ_t_phys', 'tendency of momentum z due to physics', 'kg/m2/s2', zdim='half' )
666  call hist_in(damp_t, 'MOMZ_t_damp', 'tendency of momentum z due to rayleigh damping', 'kg/m2/s2', zdim='half' )
667 #endif
668 
669  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
670 !OCL XFILL
671  do j = js-1, je+1
672  do i = is-1, ie+1
673  do k = ks, ke
674  diff(k,i,j) = momx(k,i,j) - damp_velx(k,i,j) * ( dens(k,i,j)+dens(k,i+1,j) ) * 0.5_rp
675  enddo
676  enddo
677  enddo
678 !OCL XFILL
679  do j = js, je
680  do i = is, ie
681  do k = ks, ke
682  damp = - damp_alpha_velx(k,i,j) &
683  * ( diff(k,i,j) & ! rayleigh damping
684  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
685  * 0.125_rp * bnd_smoother_fact ) ! horizontal smoother
686  momx_t(k,i,j) = momx_tp(k,i,j) & ! tendency from physical step
687  + damp
688 #ifdef HIST_TEND
689  damp_t(k,i,j) = damp
690 #endif
691  enddo
692  enddo
693  enddo
694 !OCL XFILL
695  do j = js, je
696  do i = is, ie
697  momx_t( 1:ks-1,i,j) = 0.0_rp
698  momx_t(ke+1:ka ,i,j) = 0.0_rp
699  enddo
700  enddo
701  call comm_vars8( momx_t(:,:,:), i_comm_momx_t )
702 #ifdef HIST_TEND
703  call hist_in(momx_tp, 'MOMX_t_phys', 'tendency of momentum x due to physics', 'kg/m2/s2', xdim='half' )
704  call hist_in(damp_t, 'MOMX_t_damp', 'tendency of momentum x due to rayleigh damping', 'kg/m2/s2', xdim='half' )
705 #endif
706 
707  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
708 !OCL XFILL
709  do j = js-1, je+1
710  do i = is-1, ie+1
711  do k = ks, ke
712  diff(k,i,j) = momy(k,i,j) - damp_vely(k,i,j) * ( dens(k,i,j)+dens(k,i,j+1) ) * 0.5_rp
713  enddo
714  enddo
715  enddo
716 !OCL XFILL
717  do j = js, je
718  do i = is, ie
719  do k = ks, ke
720  damp = - damp_alpha_vely(k,i,j) &
721  * ( diff(k,i,j) & ! rayleigh damping
722  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
723  * 0.125_rp * bnd_smoother_fact ) ! horizontal smoother
724  momy_t(k,i,j) = momy_tp(k,i,j) & ! tendency from physical step
725  + damp
726 #ifdef HIST_TEND
727  damp_t(k,i,j) = damp
728 #endif
729  enddo
730  enddo
731  enddo
732 !OCL XFILL
733  do j = js, je
734  do i = is, ie
735  momy_t( 1:ks-1,i,j) = 0.0_rp
736  momy_t(ke+1:ka ,i,j) = 0.0_rp
737  enddo
738  enddo
739  call comm_vars8( momy_t(:,:,:), i_comm_momy_t )
740 #ifdef HIST_TEND
741  call hist_in(momy_tp, 'MOMY_t_phys', 'tendency of momentum y due to physics', 'kg/m2/s2', ydim='half' )
742  call hist_in(damp_t, 'MOMY_t_damp', 'tendency of momentum y due to rayleigh damping', 'kg/m2/s2', ydim='half' )
743 #endif
744 
745  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
746 !OCL XFILL
747  do j = js-1, je+2
748  do i = is-1, ie+2
749  do k = ks, ke
750  diff(k,i,j) = rhot(k,i,j) - damp_pott(k,i,j) * dens(k,i,j)
751  enddo
752  enddo
753  enddo
754 !OCL XFILL
755  do j = js, je
756  do i = is, ie
757  do k = ks, ke
758  damp = - damp_alpha_pott(k,i,j) &
759  * ( diff(k,i,j) & ! rayleigh damping
760  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
761  * 0.125_rp * bnd_smoother_fact ) ! horizontal smoother
762  rhot_t(k,i,j) = rhot_tp(k,i,j) & ! tendency from physical step
763  + damp
764 #ifdef HIST_TEND
765  damp_t(k,i,j) = damp
766 #endif
767  enddo
768  enddo
769  enddo
770 !OCL XFILL
771  do j = js, je
772  do i = is, ie
773  rhot_t( 1:ks-1,i,j) = 0.0_rp
774  rhot_t(ke+1:ka ,i,j) = 0.0_rp
775  enddo
776  enddo
777  call comm_vars8( rhot_t(:,:,:), i_comm_rhot_t )
778 #ifdef HIST_TEND
779  call hist_in(rhot_tp, 'RHOT_t_phys', 'tendency of rho*theta temperature due to physics', 'K kg/m3/s' )
780  call hist_in(damp_t, 'RHOT_t_damp', 'tendency of rho*theta temperature due to rayleigh damping', 'K kg/m3/s' )
781 #endif
782 
783  call comm_wait ( dens_t(:,:,:), i_comm_dens_t, .false. )
784  call comm_wait ( momz_t(:,:,:), i_comm_momz_t, .false. )
785  call comm_wait ( momx_t(:,:,:), i_comm_momx_t, .false. )
786  call comm_wait ( momy_t(:,:,:), i_comm_momy_t, .false. )
787  call comm_wait ( rhot_t(:,:,:), i_comm_rhot_t, .false. )
788 
789  call prof_rapend ("DYN_Large_Tendency", 2)
790 
791  call prof_rapstart("DYN_Large_Numfilter", 2)
792 
793  !-----< prepare numerical diffusion coefficient >-----
794 
795  if ( nd_coef == 0.0_rp ) then
796 !OCL XFILL
797  num_diff(:,:,:,:,:) = 0.0_rp
798  else
799  call atmos_dyn_numfilter_coef( num_diff(:,:,:,:,:), & ! [OUT]
800  dens, momz, momx, momy, rhot, & ! [IN]
801  cdz, cdx, cdy, fdz, fdx, fdy, dts, & ! [IN]
802  ref_dens, ref_pott, & ! [IN]
803  nd_coef, nd_order, nd_sfc_fact, nd_use_rs ) ! [IN]
804  endif
805 
806  call prof_rapend ("DYN_Large_Numfilter", 2)
807 
808  if ( divdmp_coef > 0.0_rp ) then
809 
810  call atmos_dyn_divergence( ddiv, & ! (out)
811  momz, momx, momy, & ! (in)
812  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
813  rcdz, rcdx, rcdy, rfdz, fdz ) ! (in)
814 
815  else
816 
817 !XFILL
818  ddiv = 0.0_rp
819 
820  end if
821 
822  !------------------------------------------------------------------------
823  ! Start short time integration
824  !------------------------------------------------------------------------
825 
826  call prof_rapstart("DYN_Short_Tinteg", 2)
827 
828  call atmos_dyn_tinteg_short( dens, momz, momx, momy, rhot, prog, & ! (inout)
829  mflx_hi, tflx_hi, & ! (inout)
830  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! (in)
831  dpres0, rt2p, corioli, & ! (in)
832  num_diff, divdmp_coef, ddiv, & ! (in)
833  flag_fct_momentum, flag_fct_t, & ! (in)
834  flag_fct_along_stream, & ! (in)
835  cdz, fdz, fdx, fdy, & ! (in)
836  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
837  phi, gsqrt, j13g, j23g, j33g, mapf, & ! (in)
838  ref_dens, ref_rhot, & ! (in)
839  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
840  dts ) ! (in)
841 
842  call prof_rapend ("DYN_Short_Tinteg", 2)
843 
844 #ifdef CHECK_MASS
845  call check_mass( &
846  dens, damp_dens, &
847  mflx_hi, tflx_hi, &
848  gsqrt, mapf, &
849  rcdx, rcdy, &
850  dts, step, &
851  bnd_w, bnd_e, bnd_s, bnd_n )
852 #endif
853 
854  do j = js, je
855  do i = is, ie
856  dens( 1:ks-1,i,j) = dens(ks,i,j)
857  momz( 1:ks-1,i,j) = momz(ks,i,j)
858  momx( 1:ks-1,i,j) = momx(ks,i,j)
859  momy( 1:ks-1,i,j) = momy(ks,i,j)
860  rhot( 1:ks-1,i,j) = rhot(ks,i,j)
861  do iv = 1, va
862  prog( 1:ks-1,i,j,iv) = prog(ks,i,j,iv)
863  end do
864  dens(ke+1:ka, i,j) = dens(ke,i,j)
865  momz(ke+1:ka, i,j) = momz(ke,i,j)
866  momx(ke+1:ka, i,j) = momx(ke,i,j)
867  momy(ke+1:ka, i,j) = momy(ke,i,j)
868  rhot(ke+1:ka, i,j) = rhot(ke,i,j)
869  do iv = 1, va
870  prog(ke+1:ka, i,j,iv) = prog(ke,i,j,iv)
871  end do
872  enddo
873  enddo
874 
875  call comm_vars8( dens(:,:,:), i_comm_dens )
876  call comm_vars8( momz(:,:,:), i_comm_momz )
877  call comm_vars8( momx(:,:,:), i_comm_momx )
878  call comm_vars8( momy(:,:,:), i_comm_momy )
879  call comm_vars8( rhot(:,:,:), i_comm_rhot )
880  do iv = 1, va
881  call comm_vars8( prog(:,:,:,iv), i_comm_prog(iv) )
882  end do
883  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
884  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
885  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
886  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
887  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
888  do iv = 1, va
889  call comm_wait ( prog(:,:,:,iv), i_comm_prog(iv), .false. )
890  end do
891 
892  if ( use_average ) then
893  dens_av(:,:,:) = dens_av(:,:,:) + dens(:,:,:)
894  momz_av(:,:,:) = momz_av(:,:,:) + momz(:,:,:)
895  momx_av(:,:,:) = momx_av(:,:,:) + momx(:,:,:)
896  momy_av(:,:,:) = momy_av(:,:,:) + momy(:,:,:)
897  rhot_av(:,:,:) = rhot_av(:,:,:) + rhot(:,:,:)
898  endif
899 
900 #ifndef DRY
901  mflx_av(:,:,:,:) = mflx_av(:,:,:,:) + mflx_hi(:,:,:,:)
902 #endif
903 
904  enddo ! dynamical steps
905 
906  if ( use_average ) then
907  dens_av(:,:,:) = dens_av(:,:,:) / nstep
908  momz_av(:,:,:) = momz_av(:,:,:) / nstep
909  momx_av(:,:,:) = momx_av(:,:,:) / nstep
910  momy_av(:,:,:) = momy_av(:,:,:) / nstep
911  rhot_av(:,:,:) = rhot_av(:,:,:) / nstep
912  endif
913 
914 #ifndef DRY
915  !###########################################################################
916  ! Update Tracers
917  !###########################################################################
918 
919 !OCL XFILL
920  mflx_hi(:,:,:,:) = mflx_av(:,:,:,:) / nstep
921 
922  call comm_vars8( mflx_hi(:,:,:,zdir), i_comm_mflx_z )
923  call comm_vars8( mflx_hi(:,:,:,xdir), i_comm_mflx_x )
924  call comm_vars8( mflx_hi(:,:,:,ydir), i_comm_mflx_y )
925  call comm_wait ( mflx_hi(:,:,:,zdir), i_comm_mflx_z, .false. )
926  call comm_wait ( mflx_hi(:,:,:,xdir), i_comm_mflx_x, .false. )
927  call comm_wait ( mflx_hi(:,:,:,ydir), i_comm_mflx_y, .false. )
928 
929  if ( use_average ) then
930 !OCL XFILL
931  qtrc_av(:,:,:,:) = 0.0_rp
932  endif
933 
934  !------------------------------------------------------------------------
935  ! Update each tracer
936  !------------------------------------------------------------------------
937 
938 #ifdef _SDM
939  do iq = 1, i_qv
940 #else
941  do iq = 1, qa
942 #endif
943 
944  call prof_rapstart("DYN_Large_Numfilter", 2)
945 
946  if ( nd_coef_q == 0.0_rp ) then
947 !OCL XFILL
948  num_diff_q(:,:,:,:) = 0.0_rp
949  else
950  call atmos_dyn_numfilter_coef_q( num_diff_q(:,:,:,:), & ! [OUT]
951  dens00, qtrc(:,:,:,iq), & ! [IN]
952  cdz, cdx, cdy, dtl, & ! [IN]
953  ref_qv, iq, & ! [IN]
954  nd_coef_q, nd_order, nd_sfc_fact, nd_use_rs ) ! [IN]
955  endif
956 
957  call prof_rapend ("DYN_Large_Numfilter", 2)
958 
959  call prof_rapstart("DYN_Tracer_Tinteg", 2)
960 
962  qtrc(:,:,:,iq), & ! (inout)
963  qtrc0(:,:,:,iq), rhoq_t(:,:,:,iq), &! (in)
964  dens00, dens, & ! (in)
965  mflx_hi, num_diff_q, & ! (in)
966  gsqrt, mapf(:,:,:,i_xy), & ! (in)
967  cdz, rcdz, rcdx, rcdy, & ! (in)
968  dtl, & ! (in)
969  llast .AND. flag_fct_tracer, & ! (in)
970  flag_fct_along_stream ) ! (in)
971 
972 
973  call prof_rapend ("DYN_Tracer_Tinteg", 2)
974 
975  if ( use_average ) then
976  qtrc_av(:,:,:,iq) = qtrc(:,:,:,iq)
977  endif
978 
979  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
980 
981  enddo ! scalar quantities loop
982 
983  do iq = 1, qa
984  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
985  enddo
986 #endif
987 
988  return
integer, public is
start point of inner domain: x, local
integer, public i_xvz
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
subroutine, public atmos_dyn_numfilter_coef_q(num_diff_q, DENS, QTRC, CDZ, CDX, CDY, dt, REF_qv, iq, ND_COEF, ND_ORDER, ND_SFC_FACT, ND_USE_RS)
Calc coefficient of numerical filter.
integer, public va
Definition: scale_index.F90:38
module Atmosphere / Dynamics Temporal integration
subroutine, public atmos_dyn_numfilter_coef(num_diff, DENS, MOMZ, MOMX, MOMY, RHOT, CDZ, CDX, CDY, FDZ, FDX, FDY, DT, REF_dens, REF_pott, ND_COEF, ND_ORDER, ND_SFC_FACT, ND_USE_RS)
Calc coefficient of numerical filter.
subroutine, public atmos_dyn_fvm_fluxy_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, parameter, public zdir
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxx_xyz
integer, parameter, public ydir
integer, public ke
end point of inner domain: z, local
integer, parameter, public xdir
integer, public i_xy
module Atmosphere / Dynamics Temporal integration
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
integer, public ia
of x whole cells (local, with HALO)
module GRIDTRANS
integer, public ka
of z whole cells (local, with HALO)
integer, public i_uy
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:93
integer, public i_xyw
module COMMUNICATION
Definition: scale_comm.F90:23
subroutine, public atmos_dyn_fvm_fluxz_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
integer, public js
start point of inner domain: y, local
procedure(tinteg), pointer, public atmos_dyn_tinteg_tracer
module Atmosphere / Dynamics common
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
integer, public ks
start point of inner domain: z, local
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
integer, public i_xv
integer, public ie
end point of inner domain: x, local
real(rp), public const_eps
small number
Definition: scale_const.F90:36
module scale_atmos_dyn_fvm_flux
integer, public i_uyz
subroutine, public atmos_dyn_fct(qflx_anti, phi_in, DENS0, DENS, qflx_hi, qflx_lo, mflx_hi, rdz, rdx, rdy, GSQRT, MAPF, dt, flag_vect)
Flux Correction Transport Limiter.
subroutine, public atmos_dyn_divergence(DDIV, MOMZ, MOMX, MOMY, GSQRT, J13G, J23G, J33G, MAPF, RCDZ, RCDX, RCDY, RFDZ, FDZ)
module ATMOSPHERE / Boundary treatment
subroutine, public atmos_dyn_fvm_fluxx_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
module HISTORY
procedure(short), pointer, public atmos_dyn_tinteg_short
integer, public i_xyz
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxz_xyz
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxy_xyz
module scale_atmos_dyn_fvm_flux_ud1
real(rp), public atmos_boundary_smoother_fact
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ check_mass()

subroutine scale_atmos_dyn_tstep_large_fvm_heve::check_mass ( real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension(ka,ia,ja), intent(in)  DAMP_DENS,
real(rp), dimension (ka,ia,ja,3), intent(in)  mflx_hi,
real(rp), dimension (ka,ia,ja,3), intent(in)  tflx_hi,
real(rp), dimension (ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension ( ia,ja,2,7), intent(in)  MAPF,
real(rp), dimension(ia), intent(in)  RCDX,
real(rp), dimension(ja), intent(in)  RCDY,
real(rp), intent(in)  dt,
integer, intent(in)  step,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N 
)

Definition at line 1000 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

References scale_comm::comm_datatype, scale_comm::comm_world, scale_gridtrans::i_xy, scale_gridtrans::i_xyz, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, scale_grid_real::real_vol, scale_grid_index::xdir, scale_grid_index::ydir, and scale_grid_index::zdir.

Referenced by atmos_dyn_tstep_large_fvm_heve().

1000  use mpi
1001  use scale_grid_real, only: &
1002  vol => real_vol
1003  use scale_comm, only: &
1004  comm_datatype, &
1005  comm_world
1006  use scale_history, only: &
1007  hist_in
1008  use scale_gridtrans, only: &
1009  i_xyz, &
1010  i_xy
1011  implicit none
1012  real(RP), intent(in) :: dens (ka,ia,ja)
1013  real(RP), intent(in) :: damp_dens(ka,ia,ja)
1014  real(RP), intent(in) :: mflx_hi (ka,ia,ja,3)
1015  real(RP), intent(in) :: tflx_hi (ka,ia,ja,3)
1016  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
1017  real(RP), intent(in) :: mapf ( ia,ja,2,7)
1018  real(RP), intent(in) :: rcdx(ia)
1019  real(RP), intent(in) :: rcdy(ja)
1020  real(RP), intent(in) :: dt
1021  integer, intent(in) :: step
1022  logical, intent(in) :: bnd_w
1023  logical, intent(in) :: bnd_e
1024  logical, intent(in) :: bnd_s
1025  logical, intent(in) :: bnd_n
1026 
1027  ! lateral boundary flux
1028  real(RP) :: mflx_lb_horizontal(ka)
1029  real(RP) :: allmflx_lb_horizontal(ka)
1030  real(RP) :: mflx_lb_total
1031  real(RP) :: mass_total
1032  real(RP) :: mass_total2
1033  real(RP) :: allmflx_lb_total
1034  real(RP) :: allmass_total
1035  real(RP) :: allmass_total2
1036 
1037  integer :: k, i, j
1038  integer :: ierr
1039 
1040 
1041  call hist_in(mflx_hi(:,:,:,zdir), 'MFLXZ', 'momentum flux of z-direction', 'kg/m2/s', zdim='half' )
1042  call hist_in(mflx_hi(:,:,:,xdir), 'MFLXX', 'momentum flux of x-direction', 'kg/m2/s', xdim='half' )
1043  call hist_in(mflx_hi(:,:,:,ydir), 'MFLXY', 'momentum flux of y-direction', 'kg/m2/s', ydim='half' )
1044 
1045  call hist_in(tflx_hi(:,:,:,zdir), 'TFLXZ', 'potential temperature flux of z-direction', 'K*kg/m2/s', zdim='half' )
1046  call hist_in(tflx_hi(:,:,:,xdir), 'TFLXX', 'potential temperature flux of x-direction', 'K*kg/m2/s', xdim='half' )
1047  call hist_in(tflx_hi(:,:,:,ydir), 'TFLXY', 'potential temperature flux of y-direction', 'K*kg/m2/s', ydim='half' )
1048 
1049  mflx_lb_total = 0.0_rp
1050  mflx_lb_horizontal(:) = 0.0_rp
1051  allmflx_lb_horizontal(:) = 0.0_rp
1052 
1053  if ( bnd_w ) then ! for western boundary
1054  i = is
1055  do j = js, je
1056  do k = ks, ke
1057  mflx_lb_total = mflx_lb_total + mflx_hi(k,i-1,j,xdir) * rcdx(i) * vol(k,i,j) &
1058  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1059  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) + mflx_hi(k,i-1,j,xdir) * rcdx(i) * vol(k,i,j) &
1060  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1061 
1062  end do
1063  end do
1064  end if
1065  if ( bnd_e ) then ! for eastern boundary
1066  i = ie
1067  do j = js, je
1068  do k = ks, ke
1069  mflx_lb_total = mflx_lb_total - mflx_hi(k,i,j,xdir) * rcdx(i) * vol(k,i,j) &
1070  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1071  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) - mflx_hi(k,i,j,xdir) * rcdx(i) * vol(k,i,j) &
1072  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1073  end do
1074  end do
1075  end if
1076  if ( bnd_s ) then ! for sourthern boundary
1077  j = js
1078  do i = is, ie
1079  do k = ks, ke
1080  mflx_lb_total = mflx_lb_total + mflx_hi(k,i,j-1,ydir) * rcdy(j) * vol(k,i,j) &
1081  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1082  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) + mflx_hi(k,i,j-1,ydir) * rcdy(j) * vol(k,i,j) &
1083  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1084  end do
1085  end do
1086  end if
1087  if ( bnd_n ) then ! for northern boundary
1088  j = je
1089  do i = is, ie
1090  do k = ks, ke
1091  mflx_lb_total = mflx_lb_total - mflx_hi(k,i,j,ydir) * rcdy(j) * vol(k,i,j) &
1092  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1093  mflx_lb_horizontal(k) = mflx_lb_horizontal(k) - mflx_hi(k,i,j,ydir) * rcdy(j) * vol(k,i,j) &
1094  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) * dt
1095  end do
1096  end do
1097  end if
1098 
1099  mass_total = 0.0_rp
1100  mass_total2 = 0.0_rp
1101 
1102  ! check total mass in the inner region
1103  do j = js, je
1104  do i = is, ie
1105  do k = ks, ke
1106  mass_total = mass_total + dens(k,i,j) * vol(k,i,j)
1107  mass_total2 = mass_total2 + damp_dens(k,i,j) * vol(k,i,j)
1108  end do
1109  end do
1110  end do
1111 
1112  call mpi_allreduce( mflx_lb_total, &
1113  allmflx_lb_total, &
1114  1, &
1115  comm_datatype, &
1116  mpi_sum, &
1117  comm_world, &
1118  ierr )
1119 
1120  if( io_l ) write(io_fid_log,'(A,1x,I1,1x,ES24.17)') 'total mflx_lb:', step, allmflx_lb_total
1121 
1122  call mpi_allreduce( mass_total, &
1123  allmass_total, &
1124  1, &
1125  comm_datatype, &
1126  mpi_sum, &
1127  comm_world, &
1128  ierr )
1129 
1130  if( io_l ) write(io_fid_log,'(A,1x,I1,1x,ES24.17)') 'total mass :', step, allmass_total
1131 
1132  call mpi_allreduce( mass_total2, &
1133  allmass_total2, &
1134  1, &
1135  comm_datatype, &
1136  mpi_sum, &
1137  comm_world, &
1138  ierr )
1139 
1140  if( io_l ) write(io_fid_log,'(A,1x,I1,1x,ES24.17)') 'total mass2 :', step, allmass_total2
1141 
1142  call mpi_allreduce( mflx_lb_horizontal(ks:ke), &
1143  allmflx_lb_horizontal(ks:ke), &
1144  kmax, &
1145  comm_datatype, &
1146  mpi_sum, &
1147  comm_world, &
1148  ierr )
1149 
1150  call hist_in(allmflx_lb_horizontal(:), 'ALLMOM_lb_hz', &
1151  'horizontally total momentum flux from lateral boundary', 'kg/m2/s' )
1152 
1153  return
integer, public is
start point of inner domain: x, local
integer, public comm_datatype
datatype of variable
Definition: scale_comm.F90:117
integer, public je
end point of inner domain: y, local
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, parameter, public zdir
integer, parameter, public ydir
integer, public ke
end point of inner domain: z, local
integer, parameter, public xdir
integer, public i_xy
integer, public ia
of x whole cells (local, with HALO)
module GRIDTRANS
module GRID (real space)
integer, public ka
of z whole cells (local, with HALO)
integer, public comm_world
communication world ID
Definition: scale_comm.F90:118
integer, public kmax
of computational cells: z
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
real(rp), dimension(:,:,:), allocatable, public real_vol
control volume [m3]
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
module HISTORY
integer, public i_xyz
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function: