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_finalize
 finalize 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 95 of file scale_atmos_dyn.F90.

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

subroutine, public scale_atmos_dyn::atmos_dyn_finalize

finalize

Definition at line 236 of file scale_atmos_dyn.F90.

236  use scale_atmos_dyn_tstep_large, only: &
238  use scale_atmos_dyn_tinteg_short, only: &
240  use scale_atmos_dyn_tinteg_tracer, only: &
242  use scale_atmos_dyn_fvm_numfilter, only: &
244  use scale_spnudge, only: &
246 
247  deallocate( i_comm_prog )
248  deallocate( i_comm_qtrc )
249 
250  if ( .NOT. dyn_none ) then
251 
253 
255 
257 
259 
260  call spnudge_finalize
261 
262  !$acc exit data delete(num_diff, num_diff_q, wdamp_coef)
263  deallocate( num_diff )
264  deallocate( num_diff_q )
265  deallocate( wdamp_coef )
266  end if
267 
268  return

References scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_finalize(), scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_finalize, scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_finalize, scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_finalize, and scale_spnudge::spnudge_finalize().

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver_finalize().

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,ia,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 299 of file scale_atmos_dyn.F90.

299  use scale_prc_cartesc, only: &
300  prc_twod
301  use scale_comm_cartesc, only: &
302  comm_vars8, &
303  comm_wait
304  use scale_atmos_dyn_tinteg_large, only: &
306  implicit none
307 
308  real(RP), intent(inout) :: DENS(KA,IA,JA)
309  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
310  real(RP), intent(inout) :: MOMX(KA,IA,JA)
311  real(RP), intent(inout) :: MOMY(KA,IA,JA)
312  real(RP), intent(inout) :: RHOT(KA,IA,JA)
313  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
314  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
315 
316  real(RP), intent(inout) :: DENS_av(KA,IA,JA)
317  real(RP), intent(inout) :: MOMZ_av(KA,IA,JA)
318  real(RP), intent(inout) :: MOMX_av(KA,IA,JA)
319  real(RP), intent(inout) :: MOMY_av(KA,IA,JA)
320  real(RP), intent(inout) :: RHOT_av(KA,IA,JA)
321  real(RP), intent(inout) :: QTRC_av(KA,IA,JA,QA)
322 
323  real(RP), intent(in) :: DENS_tp(KA,IA,JA)
324  real(RP), intent(in) :: MOMZ_tp(KA,IA,JA)
325  real(RP), intent(in) :: MOMX_tp(KA,IA,JA)
326  real(RP), intent(in) :: MOMY_tp(KA,IA,JA)
327  real(RP), intent(in) :: RHOT_tp(KA,IA,JA)
328  real(RP), intent(in) :: RHOQ_tp(KA,IA,JA,QA)
329 
330  real(RP), intent(in) :: CORIOLIS(IA,JA)
331 
332  real(RP), intent(in) :: CDZ (KA)
333  real(RP), intent(in) :: CDX (IA)
334  real(RP), intent(in) :: CDY (JA)
335  real(RP), intent(in) :: FDZ (KA-1)
336  real(RP), intent(in) :: FDX (IA-1)
337  real(RP), intent(in) :: FDY (JA-1)
338  real(RP), intent(in) :: RCDZ(KA)
339  real(RP), intent(in) :: RCDX(IA)
340  real(RP), intent(in) :: RCDY(JA)
341  real(RP), intent(in) :: RFDZ(KA-1)
342  real(RP), intent(in) :: RFDX(IA-1)
343  real(RP), intent(in) :: RFDY(JA-1)
344 
345  real(RP), intent(in) :: PHI (KA,IA,JA)
346  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
347  real(RP), intent(in) :: J13G (KA,IA,JA,7)
348  real(RP), intent(in) :: J23G (KA,IA,JA,7)
349  real(RP), intent(in) :: J33G
350  real(RP), intent(in) :: MAPF (IA,JA,2,4)
351 
352  real(RP), intent(in) :: AQ_R (QA)
353  real(RP), intent(in) :: AQ_CV (QA)
354  real(RP), intent(in) :: AQ_CP (QA)
355  real(RP), intent(in) :: AQ_MASS(QA)
356 
357  real(RP), intent(in) :: REF_dens(KA,IA,JA)
358  real(RP), intent(in) :: REF_pott(KA,IA,JA)
359  real(RP), intent(in) :: REF_qv (KA,IA,JA)
360  real(RP), intent(in) :: REF_pres(KA,IA,JA)
361  real(RP), intent(in) :: ND_COEF
362  real(RP), intent(in) :: ND_COEF_Q
363  integer, intent(in) :: ND_LAPLACIAN_NUM
364  real(RP), intent(in) :: ND_SFC_FACT
365  logical, intent(in) :: ND_USE_RS
366 
367  integer, intent(in) :: BND_QA
368  integer, intent(in) :: BND_IQ(QA)
369  real(RP), intent(in) :: BND_SMOOTHER_FACT
370 
371  real(RP), intent(in) :: DAMP_DENS(KA,IA,JA)
372  real(RP), intent(in) :: DAMP_VELZ(KA,IA,JA)
373  real(RP), intent(in) :: DAMP_VELX(KA,IA,JA)
374  real(RP), intent(in) :: DAMP_VELY(KA,IA,JA)
375  real(RP), intent(in) :: DAMP_POTT(KA,IA,JA)
376  real(RP), intent(in) :: DAMP_QTRC(KA,IA,JA,BND_QA)
377 
378  real(RP), intent(in) :: DAMP_alpha_DENS(KA,IA,JA)
379  real(RP), intent(in) :: DAMP_alpha_VELZ(KA,IA,JA)
380  real(RP), intent(in) :: DAMP_alpha_VELX(KA,IA,JA)
381  real(RP), intent(in) :: DAMP_alpha_VELY(KA,IA,JA)
382  real(RP), intent(in) :: DAMP_alpha_POTT(KA,IA,JA)
383  real(RP), intent(in) :: DAMP_alpha_QTRC(KA,IA,JA,BND_QA)
384  real(RP), intent(in) :: MFLUX_OFFSET_X(KA,JA,2)
385  real(RP), intent(in) :: MFLUX_OFFSET_Y(KA,IA,2)
386 
387  real(RP), intent(in) :: divdmp_coef
388 
389  logical, intent(in) :: FLAG_TRACER_SPLIT_TEND
390  logical, intent(in) :: FLAG_FCT_MOMENTUM
391  logical, intent(in) :: FLAG_FCT_T
392  logical, intent(in) :: FLAG_FCT_TRACER
393  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
394 
395  logical, intent(in) :: USE_AVERAGE
396 
397  integer, intent(in) :: I_QV
398 
399  real(DP), intent(in) :: DTSEC
400  real(DP), intent(in) :: DTSEC_DYN
401 
402  ! for time integartion
403  real(RP) :: DENS00 (KA,IA,JA) ! saved density before update
404 
405  ! For tracer advection
406  real(RP) :: dt
407 
408  integer :: i, j, k, iq
409  !---------------------------------------------------------------------------
410 
411  log_progress(*) 'atmosphere / dynamics'
412 
413  dt = real(dtsec, kind=rp)
414 
415  if ( dyn_none ) then
416 
417  call prof_rapstart("DYN_Tinteg", 2)
418 
419  !$acc data create(DENS00)
420 
421  !$acc kernels
422  do j = js, je
423  do i = is, ie
424  do k = ks, ke
425  dens00(k,i,j) = dens(k,i,j) ! save previous value
426 
427  dens(k,i,j) = dens(k,i,j) + dens_tp(k,i,j) * dt
428  enddo
429  enddo
430  enddo
431  !$acc end kernels
432 
433  !$acc kernels
434  do j = js, je
435  do i = is, ie
436  do k = ks, ke
437  momz(k,i,j) = momz(k,i,j) + momz_tp(k,i,j) * dt
438  enddo
439  enddo
440  enddo
441  !$acc end kernels
442 
443  !$acc kernels
444  do j = js, je
445  do i = is, ie
446  do k = ks, ke
447  momx(k,i,j) = momx(k,i,j) + momx_tp(k,i,j) * dt
448  enddo
449  enddo
450  enddo
451  !$acc end kernels
452 
453  !$acc kernels
454  do j = js, je
455  do i = is, ie
456  do k = ks, ke
457  momy(k,i,j) = momy(k,i,j) + momy_tp(k,i,j) * dt
458  enddo
459  enddo
460  enddo
461  !$acc end kernels
462 
463  !$acc kernels
464  do j = js, je
465  do i = is, ie
466  do k = ks, ke
467  rhot(k,i,j) = rhot(k,i,j) + rhot_tp(k,i,j) * dt
468  enddo
469  enddo
470  enddo
471  !$acc end kernels
472 
473  !$acc kernels
474  do iq = 1, qa
475  do j = js, je
476  do i = is, ie
477  do k = ks, ke
478  qtrc(k,i,j,iq) = max( 0.0_rp, &
479  qtrc(k,i,j,iq) * dens00(k,i,j) + rhoq_tp(k,i,j,iq) * dt &
480  ) / dens(k,i,j)
481  enddo
482  enddo
483  enddo
484  enddo
485  !$acc end kernels
486 
487  !$acc end data
488 
489  call comm_vars8( dens(:,:,:), i_comm_dens )
490  call comm_vars8( momz(:,:,:), i_comm_momz )
491  call comm_vars8( momx(:,:,:), i_comm_momx )
492  call comm_vars8( momy(:,:,:), i_comm_momy )
493  call comm_vars8( rhot(:,:,:), i_comm_rhot )
494  do iq = 1, qa
495  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
496  enddo
497 
498  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
499  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
500  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
501  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
502  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
503  do iq = 1, qa
504  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
505  enddo
506 
507  if ( use_average ) then
508  !$acc kernels
509  dens_av(:,:,:) = dens(:,:,:)
510  !$acc end kernels
511  !$acc kernels
512  momz_av(:,:,:) = momz(:,:,:)
513  !$acc end kernels
514  !$acc kernels
515  momx_av(:,:,:) = momx(:,:,:)
516  !$acc end kernels
517  !$acc kernels
518  momy_av(:,:,:) = momy(:,:,:)
519  !$acc end kernels
520  !$acc kernels
521  rhot_av(:,:,:) = rhot(:,:,:)
522  !$acc end kernels
523  !$acc kernels
524  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
525  !$acc end kernels
526  endif
527 
528  call prof_rapend("DYN_Tinteg", 2)
529 
530  return
531  endif ! if DYN_NONE == .true.
532 
533  call prof_rapstart("DYN_Tinteg", 2)
534 
535  !$acc data copy(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG) &
536  !$acc copyin(DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, &
537  !$acc CORIOLIS, &
538  !$acc CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
539  !$acc PHI, GSQRT, J13G, J23G, MAPF, &
540  !$acc AQ_R, AQ_CV, AQ_CP, AQ_MASS, &
541  !$acc REF_dens, REF_pott, REF_qv, REF_pres, &
542  !$acc BND_IQ, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, &
543  !$acc DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, &
544  !$acc MFLUX_OFFSET_X, MFLUX_OFFSET_Y)
545 
546  !$acc data copy(DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av) if(USE_AVERAGE)
547 
548 
549  call atmos_dyn_tinteg_large( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
550  prog, & ! [INOUT]
551  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
552  num_diff, num_diff_q, & ! [OUT;WORK]
553  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! [IN]
554  coriolis, & ! [IN]
555  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
556  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
557  phi, gsqrt, & ! [IN]
558  j13g, j23g, j33g, mapf, & ! [IN]
559  aq_r, aq_cv, aq_cp, aq_mass, & ! [IN]
560  ref_dens, ref_pott, ref_qv, ref_pres, & ! [IN]
561  bnd_w, bnd_e, bnd_s, bnd_n, prc_twod, & ! [IN]
562  nd_coef, nd_coef_q, nd_laplacian_num, & ! [IN]
563  nd_sfc_fact, nd_use_rs, & ! [IN]
564  bnd_qa, bnd_iq, bnd_smoother_fact, & ! [IN]
565  damp_dens, damp_velz, damp_velx, & ! [IN]
566  damp_vely, damp_pott, damp_qtrc, & ! [IN]
567  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! [IN]
568  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! [IN]
569  mflux_offset_x, mflux_offset_y, & ! [IN]
570  wdamp_coef, divdmp_coef, & ! [IN]
571  flag_tracer_split_tend, & ! [IN]
572  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! [IN]
573  flag_fct_along_stream, & ! [IN]
574  use_average, & ! [IN]
575  i_qv, & ! [IN]
576  dtsec, dtsec_dyn ) ! [IN]
577 
578  call prof_rapend ("DYN_Tinteg", 2)
579 
580  !$acc end data
581  !$acc end data
582 
583  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:51
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
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:51
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:135
scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_finalize
procedure(finalize), pointer, public atmos_dyn_tinteg_tracer_finalize
Definition: scale_atmos_dyn_tinteg_tracer.F90:78
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:49
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:205
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:50
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_finalize
procedure(finalize), pointer, public atmos_dyn_tstep_large_finalize
Definition: scale_atmos_dyn_tstep_large.F90:183
scale_spnudge::spnudge_finalize
subroutine, public spnudge_finalize
Definition: scale_spnudge.F90:268
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:104
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_tinteg_short::atmos_dyn_tinteg_short_finalize
procedure(finalize), pointer, public atmos_dyn_tinteg_short_finalize
Definition: scale_atmos_dyn_tinteg_short.F90:114
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid, gid)
Register variables.
Definition: scale_comm_cartesC.F90:811
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:306
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
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:99
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:48
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:56
scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_finalize
subroutine, public atmos_dyn_fvm_numfilter_finalize
finalize
Definition: scale_atmos_dyn_fvm_numfilter.F90:396