SCALE-RM
Functions/Subroutines
scale_atmos_dyn Module Reference

module Atmosphere / Dynamics FENT + FCT More...

Functions/Subroutines

subroutine, public atmos_dyn_setup (DYN_Tinteg_Short_TYPE, DYN_Tinteg_Tracer_TYPE, DYN_Tinteg_Large_TYPE, DYN_Tstep_Tracer_TYPE, DYN_Tstep_Large_TYPE, DYN_Tstep_Short_TYPE, DYN_FVM_FLUX_TYPE, DYN_FVM_FLUX_TYPE_TRACER, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, CDZ, CDX, CDY, FDZ, FDX, FDY, wdamp_tau, wdamp_height, FZ, none)
 Setup. More...
 
subroutine, public atmos_dyn (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLIS, CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, AQ_R, AQ_CV, AQ_CP, AQ_MASS, REF_dens, REF_pott, REF_qv, REF_pres, ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, ND_SFC_FACT, ND_USE_RS, BND_QA, BND_IQ, BND_SMOOTHER_FACT, 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, MFLUX_OFFSET_X, MFLUX_OFFSET_Y, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTSEC, DTSEC_DYN)
 Dynamical Process. More...
 

Detailed Description

module Atmosphere / Dynamics FENT + FCT

Description
Dynamical core for Atmospheric process Full explicit, no terrain + tracer FCT limiter
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_dyn_setup()

subroutine, public scale_atmos_dyn::atmos_dyn_setup ( character(len=*), intent(in)  DYN_Tinteg_Short_TYPE,
character(len=*), intent(in)  DYN_Tinteg_Tracer_TYPE,
character(len=*), intent(in)  DYN_Tinteg_Large_TYPE,
character(len=*), intent(in)  DYN_Tstep_Tracer_TYPE,
character(len=*), intent(in)  DYN_Tstep_Large_TYPE,
character(len=*), intent(in)  DYN_Tstep_Short_TYPE,
character(len=*), intent(in)  DYN_FVM_FLUX_TYPE,
character(len=*), intent(in)  DYN_FVM_FLUX_TYPE_TRACER,
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), 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), intent(in)  wdamp_tau,
real(rp), intent(in)  wdamp_height,
real(rp), dimension(0:ka), intent(in)  FZ,
logical, intent(in), optional  none 
)

Setup.

Definition at line 94 of file scale_atmos_dyn.F90.

94  use scale_prc, only: &
95  prc_abort
96  use scale_prc_cartesc, only: &
97  prc_twod, &
98  prc_has_e, &
99  prc_has_w, &
100  prc_has_n, &
101  prc_has_s
102  use scale_const, only: &
103  undef => const_undef
104  use scale_comm_cartesc, only: &
106  use scale_atmos_dyn_common, only: &
108  use scale_atmos_dyn_fvm_numfilter, only: &
110  use scale_atmos_dyn_tinteg_short, only: &
112  use scale_atmos_dyn_tinteg_tracer, only: &
114  use scale_atmos_dyn_tinteg_large, only: &
116  use scale_atmos_dyn_tstep_short, only: &
118  use scale_atmos_dyn_tstep_tracer, only: &
120  use scale_atmos_dyn_tstep_large, only: &
122  use scale_atmos_dyn_fvm_flux, only: &
124  use scale_spnudge, only: &
126  implicit none
127 
128  character(len=*), intent(in) :: DYN_Tinteg_Short_TYPE
129  character(len=*), intent(in) :: DYN_Tinteg_Tracer_TYPE
130  character(len=*), intent(in) :: DYN_Tinteg_Large_TYPE
131  character(len=*), intent(in) :: DYN_Tstep_Tracer_TYPE
132  character(len=*), intent(in) :: DYN_Tstep_Large_TYPE
133  character(len=*), intent(in) :: DYN_Tstep_Short_TYPE
134  character(len=*), intent(in) :: DYN_FVM_FLUX_TYPE
135  character(len=*), intent(in) :: DYN_FVM_FLUX_TYPE_TRACER
136  real(RP), intent(inout) :: DENS(KA,IA,JA) ! MPI_RECV_INIT requires intent(inout)
137  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
138  real(RP), intent(inout) :: MOMX(KA,IA,JA)
139  real(RP), intent(inout) :: MOMY(KA,IA,JA)
140  real(RP), intent(inout) :: RHOT(KA,IA,JA)
141  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
142  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
143  real(RP), intent(in) :: CDZ(KA)
144  real(RP), intent(in) :: CDX(IA)
145  real(RP), intent(in) :: CDY(JA)
146  real(RP), intent(in) :: FDZ(KA-1)
147  real(RP), intent(in) :: FDX(IA-1)
148  real(RP), intent(in) :: FDY(JA-1)
149  real(RP), intent(in) :: wdamp_tau
150  real(RP), intent(in) :: wdamp_height
151  real(RP), intent(in) :: FZ(0:KA)
152  logical, optional, intent(in) :: none
153 
154  integer :: iv, iq
155  !---------------------------------------------------------------------------
156 
157  dyn_none = .false.
158 
159  if ( present(none) ) then
160  dyn_none = none
161  endif
162 
163  allocate( i_comm_prog(max(va,1)) )
164  allocate( i_comm_qtrc(qa) )
165 
166  if ( .NOT. dyn_none ) then
167  bnd_w = ( .NOT. prc_has_w ) .and. ( .NOT. prc_twod )
168  bnd_e = ( .NOT. prc_has_e ) .and. ( .NOT. prc_twod )
169  bnd_s = .NOT. prc_has_s
170  bnd_n = .NOT. prc_has_n
171 
172  allocate( num_diff(ka,ia,ja,5,3) )
173  allocate( num_diff_q(ka,ia,ja,3) )
174  allocate( wdamp_coef(ka) )
175  num_diff(:,:,:,:,:) = undef
176  num_diff_q(:,:,:,:) = undef
177 
178  call atmos_dyn_fvm_flux_setup ( dyn_fvm_flux_type, & ! [IN]
179  dyn_fvm_flux_type_tracer ) ! [IN]
180 
182 
183  call atmos_dyn_tstep_tracer_setup ( dyn_tstep_tracer_type ) ! [IN]
184 
185  call atmos_dyn_tstep_large_setup ( dyn_tstep_large_type, & ! [IN]
186  dens, momz, momx, momy, rhot, & ! [INOUT]
187  qtrc, prog ) ! [INOUT]
188 
189  call atmos_dyn_tinteg_short_setup ( dyn_tinteg_short_type, dyn_tstep_short_type ) ! [IN]
190 
191  call atmos_dyn_tinteg_tracer_setup( dyn_tinteg_tracer_type ) ! [IN]
192 
193  call atmos_dyn_tinteg_large_setup ( dyn_tinteg_large_type ) ! [IN]
194 
195  ! numerical diffusion
196  call atmos_dyn_fvm_numfilter_setup( num_diff, num_diff_q, & ! [INOUT]
197  cdz, cdx, cdy, fdz, fdx, fdy ) ! [IN]
198 
199  ! sponge layer
200  call atmos_dyn_wdamp_setup( wdamp_coef(:), & ! [INOUT]
201  wdamp_tau, wdamp_height, & ! [IN]
202  fz(:) ) ! [IN]
203 
204  ! setup for spectral nudging
205  call spnudge_setup( ka, ks, ke, ia, is, ie, ja, js, je )
206 
207 
208  else
209 
210  call comm_vars8_init( 'DENS', dens, i_comm_dens )
211  call comm_vars8_init( 'MOMZ', momz, i_comm_momz )
212  call comm_vars8_init( 'MOMX', momx, i_comm_momx )
213  call comm_vars8_init( 'MOMY', momy, i_comm_momy )
214  call comm_vars8_init( 'RHOT', rhot, i_comm_rhot )
215 
216  do iv = 1, va
217  i_comm_prog(iv) = 5 + iv
218  call comm_vars8_init( 'PROG', prog(:,:,:,iv), i_comm_prog(iv) )
219  enddo
220 
221  do iq = 1, qa
222  i_comm_qtrc(iq) = 5 + va + iq
223  call comm_vars8_init( 'QTRC', qtrc(:,:,:,iq), i_comm_qtrc(iq) )
224  enddo
225 
226  endif
227 
228  return

References scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_setup(), scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_setup(), scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large_setup(), scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup(), scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup(), scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup(), scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short_setup, scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup(), scale_atmos_dyn_common::atmos_dyn_wdamp_setup(), scale_comm_cartesc::comm_vars8_init(), scale_const::const_undef, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_prc::prc_abort(), scale_prc_cartesc::prc_has_e, scale_prc_cartesc::prc_has_n, scale_prc_cartesc::prc_has_s, scale_prc_cartesc::prc_has_w, scale_prc_cartesc::prc_twod, scale_tracer::qa, scale_spnudge::spnudge_setup(), and scale_index::va.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn()

subroutine, public scale_atmos_dyn::atmos_dyn ( 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), 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)  CORIOLIS,
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 (qa), intent(in)  AQ_R,
real(rp), dimension (qa), intent(in)  AQ_CV,
real(rp), dimension (qa), intent(in)  AQ_CP,
real(rp), dimension(qa), intent(in)  AQ_MASS,
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,
real(rp), intent(in)  ND_COEF,
real(rp), intent(in)  ND_COEF_Q,
integer, intent(in)  ND_LAPLACIAN_NUM,
real(rp), intent(in)  ND_SFC_FACT,
logical, intent(in)  ND_USE_RS,
integer, intent(in)  BND_QA,
integer, dimension(qa), intent(in)  BND_IQ,
real(rp), intent(in)  BND_SMOOTHER_FACT,
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), dimension(ka,ja,2), intent(in)  MFLUX_OFFSET_X,
real(rp), dimension(ka,ja,2), intent(in)  MFLUX_OFFSET_Y,
real(rp), intent(in)  divdmp_coef,
logical, intent(in)  FLAG_TRACER_SPLIT_TEND,
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,
integer, intent(in)  I_QV,
real(dp), intent(in)  DTSEC,
real(dp), intent(in)  DTSEC_DYN 
)

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 260 of file scale_atmos_dyn.F90.

260  use scale_prc_cartesc, only: &
261  prc_twod
262  use scale_comm_cartesc, only: &
263  comm_vars8, &
264  comm_wait
265  use scale_atmos_dyn_tinteg_large, only: &
267  implicit none
268 
269  real(RP), intent(inout) :: DENS(KA,IA,JA)
270  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
271  real(RP), intent(inout) :: MOMX(KA,IA,JA)
272  real(RP), intent(inout) :: MOMY(KA,IA,JA)
273  real(RP), intent(inout) :: RHOT(KA,IA,JA)
274  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
275  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
276 
277  real(RP), intent(inout) :: DENS_av(KA,IA,JA)
278  real(RP), intent(inout) :: MOMZ_av(KA,IA,JA)
279  real(RP), intent(inout) :: MOMX_av(KA,IA,JA)
280  real(RP), intent(inout) :: MOMY_av(KA,IA,JA)
281  real(RP), intent(inout) :: RHOT_av(KA,IA,JA)
282  real(RP), intent(inout) :: QTRC_av(KA,IA,JA,QA)
283 
284  real(RP), intent(in) :: DENS_tp(KA,IA,JA)
285  real(RP), intent(in) :: MOMZ_tp(KA,IA,JA)
286  real(RP), intent(in) :: MOMX_tp(KA,IA,JA)
287  real(RP), intent(in) :: MOMY_tp(KA,IA,JA)
288  real(RP), intent(in) :: RHOT_tp(KA,IA,JA)
289  real(RP), intent(in) :: RHOQ_tp(KA,IA,JA,QA)
290 
291  real(RP), intent(in) :: CORIOLIS(IA,JA)
292 
293  real(RP), intent(in) :: CDZ (KA)
294  real(RP), intent(in) :: CDX (IA)
295  real(RP), intent(in) :: CDY (JA)
296  real(RP), intent(in) :: FDZ (KA-1)
297  real(RP), intent(in) :: FDX (IA-1)
298  real(RP), intent(in) :: FDY (JA-1)
299  real(RP), intent(in) :: RCDZ(KA)
300  real(RP), intent(in) :: RCDX(IA)
301  real(RP), intent(in) :: RCDY(JA)
302  real(RP), intent(in) :: RFDZ(KA-1)
303  real(RP), intent(in) :: RFDX(IA-1)
304  real(RP), intent(in) :: RFDY(JA-1)
305 
306  real(RP), intent(in) :: PHI (KA,IA,JA)
307  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
308  real(RP), intent(in) :: J13G (KA,IA,JA,7)
309  real(RP), intent(in) :: J23G (KA,IA,JA,7)
310  real(RP), intent(in) :: J33G
311  real(RP), intent(in) :: MAPF (IA,JA,2,4)
312 
313  real(RP), intent(in) :: AQ_R (QA)
314  real(RP), intent(in) :: AQ_CV (QA)
315  real(RP), intent(in) :: AQ_CP (QA)
316  real(RP), intent(in) :: AQ_MASS(QA)
317 
318  real(RP), intent(in) :: REF_dens(KA,IA,JA)
319  real(RP), intent(in) :: REF_pott(KA,IA,JA)
320  real(RP), intent(in) :: REF_qv (KA,IA,JA)
321  real(RP), intent(in) :: REF_pres(KA,IA,JA)
322  real(RP), intent(in) :: ND_COEF
323  real(RP), intent(in) :: ND_COEF_Q
324  integer, intent(in) :: ND_LAPLACIAN_NUM
325  real(RP), intent(in) :: ND_SFC_FACT
326  logical, intent(in) :: ND_USE_RS
327 
328  integer, intent(in) :: BND_QA
329  integer, intent(in) :: BND_IQ(QA)
330  real(RP), intent(in) :: BND_SMOOTHER_FACT
331 
332  real(RP), intent(in) :: DAMP_DENS(KA,IA,JA)
333  real(RP), intent(in) :: DAMP_VELZ(KA,IA,JA)
334  real(RP), intent(in) :: DAMP_VELX(KA,IA,JA)
335  real(RP), intent(in) :: DAMP_VELY(KA,IA,JA)
336  real(RP), intent(in) :: DAMP_POTT(KA,IA,JA)
337  real(RP), intent(in) :: DAMP_QTRC(KA,IA,JA,BND_QA)
338 
339  real(RP), intent(in) :: DAMP_alpha_DENS(KA,IA,JA)
340  real(RP), intent(in) :: DAMP_alpha_VELZ(KA,IA,JA)
341  real(RP), intent(in) :: DAMP_alpha_VELX(KA,IA,JA)
342  real(RP), intent(in) :: DAMP_alpha_VELY(KA,IA,JA)
343  real(RP), intent(in) :: DAMP_alpha_POTT(KA,IA,JA)
344  real(RP), intent(in) :: DAMP_alpha_QTRC(KA,IA,JA,BND_QA)
345  real(RP), intent(in) :: MFLUX_OFFSET_X(KA,JA,2)
346  real(RP), intent(in) :: MFLUX_OFFSET_Y(KA,JA,2)
347 
348  real(RP), intent(in) :: divdmp_coef
349 
350  logical, intent(in) :: FLAG_TRACER_SPLIT_TEND
351  logical, intent(in) :: FLAG_FCT_MOMENTUM
352  logical, intent(in) :: FLAG_FCT_T
353  logical, intent(in) :: FLAG_FCT_TRACER
354  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
355 
356  logical, intent(in) :: USE_AVERAGE
357 
358  integer, intent(in) :: I_QV
359 
360  real(DP), intent(in) :: DTSEC
361  real(DP), intent(in) :: DTSEC_DYN
362 
363  ! for time integartion
364  real(RP) :: DENS00 (KA,IA,JA) ! saved density before update
365 
366  ! For tracer advection
367  real(RP) :: dt
368 
369  integer :: i, j, k, iq
370  !---------------------------------------------------------------------------
371 
372  log_progress(*) 'atmosphere / dynamics'
373 
374  dt = real(dtsec, kind=rp)
375 
376  if ( dyn_none ) then
377 
378  call prof_rapstart("DYN_Tinteg", 2)
379 
380  do j = js, je
381  do i = is, ie
382  do k = ks, ke
383  dens00(k,i,j) = dens(k,i,j) ! save previous value
384 
385  dens(k,i,j) = dens(k,i,j) + dens_tp(k,i,j) * dt
386  enddo
387  enddo
388  enddo
389 
390  do j = js, je
391  do i = is, ie
392  do k = ks, ke
393  momz(k,i,j) = momz(k,i,j) + momz_tp(k,i,j) * dt
394  enddo
395  enddo
396  enddo
397 
398  do j = js, je
399  do i = is, ie
400  do k = ks, ke
401  momx(k,i,j) = momx(k,i,j) + momx_tp(k,i,j) * dt
402  enddo
403  enddo
404  enddo
405 
406  do j = js, je
407  do i = is, ie
408  do k = ks, ke
409  momy(k,i,j) = momy(k,i,j) + momy_tp(k,i,j) * dt
410  enddo
411  enddo
412  enddo
413 
414  do j = js, je
415  do i = is, ie
416  do k = ks, ke
417  rhot(k,i,j) = rhot(k,i,j) + rhot_tp(k,i,j) * dt
418  enddo
419  enddo
420  enddo
421 
422  do iq = 1, qa
423  do j = js, je
424  do i = is, ie
425  do k = ks, ke
426  qtrc(k,i,j,iq) = max( 0.0_rp, &
427  qtrc(k,i,j,iq) * dens00(k,i,j) + rhoq_tp(k,i,j,iq) * dt &
428  ) / dens(k,i,j)
429  enddo
430  enddo
431  enddo
432  enddo
433 
434  call comm_vars8( dens(:,:,:), i_comm_dens )
435  call comm_vars8( momz(:,:,:), i_comm_momz )
436  call comm_vars8( momx(:,:,:), i_comm_momx )
437  call comm_vars8( momy(:,:,:), i_comm_momy )
438  call comm_vars8( rhot(:,:,:), i_comm_rhot )
439  do iq = 1, qa
440  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
441  enddo
442 
443  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
444  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
445  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
446  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
447  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
448  do iq = 1, qa
449  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
450  enddo
451 
452  if ( use_average ) then
453  dens_av(:,:,:) = dens(:,:,:)
454  momz_av(:,:,:) = momz(:,:,:)
455  momx_av(:,:,:) = momx(:,:,:)
456  momy_av(:,:,:) = momy(:,:,:)
457  rhot_av(:,:,:) = rhot(:,:,:)
458  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
459  endif
460 
461  call prof_rapend("DYN_Tinteg", 2)
462 
463  return
464  endif ! if DYN_NONE == .true.
465 
466  call prof_rapstart("DYN_Tinteg", 2)
467 
468  call atmos_dyn_tinteg_large( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
469  prog, & ! [INOUT]
470  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
471  num_diff, num_diff_q, & ! [OUT;WORK]
472  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! [IN]
473  coriolis, & ! [IN]
474  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
475  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
476  phi, gsqrt, & ! [IN]
477  j13g, j23g, j33g, mapf, & ! [IN]
478  aq_r, aq_cv, aq_cp, aq_mass, & ! [IN]
479  ref_dens, ref_pott, ref_qv, ref_pres, & ! [IN]
480  bnd_w, bnd_e, bnd_s, bnd_n, prc_twod, & ! [IN]
481  nd_coef, nd_coef_q, nd_laplacian_num, & ! [IN]
482  nd_sfc_fact, nd_use_rs, & ! [IN]
483  bnd_qa, bnd_iq, bnd_smoother_fact, & ! [IN]
484  damp_dens, damp_velz, damp_velx, & ! [IN]
485  damp_vely, damp_pott, damp_qtrc, & ! [IN]
486  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! [IN]
487  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! [IN]
488  mflux_offset_x, mflux_offset_y, & ! [IN]
489  wdamp_coef, divdmp_coef, & ! [IN]
490  flag_tracer_split_tend, & ! [IN]
491  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! [IN]
492  flag_fct_along_stream, & ! [IN]
493  use_average, & ! [IN]
494  i_qv, & ! [IN]
495  dtsec, dtsec_dyn ) ! [IN]
496 
497  call prof_rapend ("DYN_Tinteg", 2)
498 
499  return

References scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_prc_cartesc::prc_twod, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_tracer::qa.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_spnudge::spnudge_setup
subroutine, public spnudge_setup(KA, KS, KE, IA, IS, IE, JA, JS, JE)
Definition: scale_spnudge.F90:50
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large_setup
subroutine, public atmos_dyn_tinteg_large_setup(ATMOS_DYN_Tinteg_large_TYPE)
Register.
Definition: scale_atmos_dyn_tinteg_large.F90:190
scale_prc_cartesc::prc_has_s
logical, public prc_has_s
Definition: scale_prc_cartesC.F90:50
scale_atmos_dyn_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
scale_atmos_dyn_fvm_flux
module scale_atmos_dyn_fvm_flux
Definition: scale_atmos_dyn_fvm_flux.F90:13
scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large
procedure(large), pointer, public atmos_dyn_tinteg_large
Definition: scale_atmos_dyn_tinteg_large.F90:169
scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup
subroutine, public atmos_dyn_tinteg_short_setup(ATMOS_DYN_Tinteg_short_TYPE, ATMOS_DYN_Tstep_short_TYPE)
Register.
Definition: scale_atmos_dyn_tinteg_short.F90:130
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:48
scale_atmos_dyn_tinteg_short
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_short.F90:11
scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup
subroutine, public atmos_dyn_tstep_large_setup(Tstep_large_type, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG)
Register.
Definition: scale_atmos_dyn_tstep_large.F90:201
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:49
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
scale_atmos_dyn_common::atmos_dyn_wdamp_setup
subroutine, public atmos_dyn_wdamp_setup(wdamp_coef, wdamp_tau, wdamp_height, FZ)
Setup.
Definition: scale_atmos_dyn_common.F90:69
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_atmos_dyn_tstep_tracer
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_tracer.F90:12
scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_setup
subroutine, public atmos_dyn_fvm_numfilter_setup(num_diff, num_diff_q, CDZ, CDX, CDY, FDZ, FDX, FDY)
Setup.
Definition: scale_atmos_dyn_fvm_numfilter.F90:103
scale_atmos_dyn_fvm_numfilter
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_fvm_numfilter.F90:12
scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup
subroutine, public atmos_dyn_tstep_tracer_setup(ATMOS_DYN_TSTEP_TRACER_TYPE)
Register.
Definition: scale_atmos_dyn_tstep_tracer.F90:93
scale_atmos_dyn_tstep_short
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_short.F90:12
scale_spnudge
Definition: scale_spnudge.F90:2
scale_atmos_dyn_tinteg_tracer
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_tracer.F90:12
scale_atmos_dyn_tstep_large
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_large.F90:12
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_dyn_tinteg_large
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_large.F90:12
scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short_setup
procedure(short_setup), pointer, public atmos_dyn_tstep_short_setup
Definition: scale_atmos_dyn_tstep_short.F90:133
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_setup
subroutine, public atmos_dyn_fvm_flux_setup(scheme, scheme_tracer)
setup
Definition: scale_atmos_dyn_fvm_flux.F90:281
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup
subroutine, public atmos_dyn_tinteg_tracer_setup(ATMOS_DYN_Tinteg_tracer_TYPE)
Register.
Definition: scale_atmos_dyn_tinteg_tracer.F90:95
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:47
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:55