SCALE-RM
Functions/Subroutines | Variables
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_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, enable_coriolis, lat, 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, 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_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, DTSEC, DTSEC_DYN)
 Dynamical Process. More...
 

Variables

real(rp), dimension(:,:), allocatable, public corioli
 

Detailed Description

module Atmosphere / Dynamics FENT + FCT

Description
Dynamical core for Atmospheric process Full explicit, no terrain + tracer FCT limiter
Author
Team SCALE
History
  • 2011-11-11 (H.Yashiro) [new] Imported from SCALE-LES ver.2
  • 2011-11-11 (H.Yashiro) [mod] Merged with Y.Miyamoto's
  • 2011-12-11 (H.Yashiro) [mod] Use reference state
  • 2011-12-26 (Y.Miyamoto) [mod] Add numerical diffusion into mass flux calc
  • 2012-01-04 (H.Yashiro) [mod] Nonblocking communication (Y.Ohno)
  • 2012-01-25 (H.Yashiro) [fix] Bugfix (Y.Miyamoto)
  • 2011-01-25 (H.Yashiro) [mod] sprit as "full" FCT (Y.Miyamoto)
  • 2012-02-14 (H.Yashiro) [mod] Cache tiling
  • 2012-03-14 (H.Yashiro) [mod] Bugfix (Y.Miyamoto)
  • 2012-03-23 (H.Yashiro) [mod] Explicit index parameter inclusion
  • 2012-04-09 (H.Yashiro) [mod] Integrate RDMA communication
  • 2012-06-10 (Y.Miyamoto) [mod] large-scale divergence (from H.Yashiro's)
  • 2012-07-13 (H.Yashiro) [mod] prevent conditional branching in FCT
  • 2012-07-27 (Y.Miyamoto) [mod] divegence damping option
  • 2012-08-16 (S.Nishizawa) [mod] use FCT for momentum and temperature
  • 2012-09-21 (Y.Sato) [mod] merge DYCOMS-II experimental set
  • 2013-03-26 (Y.Sato) [mod] modify Large scale forcing and corioli forcing
  • 2013-04-04 (Y.Sato) [mod] modify Large scale forcing
  • 2013-06-14 (S.Nishizawa) [mod] enable to change order of numerical diffusion
  • 2013-06-18 (S.Nishizawa) [mod] split part of RK to other files
  • 2013-06-20 (S.Nishizawa) [mod] split large scale sining to other file
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
ALLMOM_lb_hz horizontally total momentum flux from lateral boundary kg/m2/s allmflx_lb_horizontal
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
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

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_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)  enable_coriolis,
real(rp), dimension(ia,ja), intent(in)  lat,
logical, intent(in), optional  none 
)

Setup.

Definition at line 122 of file scale_atmos_dyn.F90.

References scale_atmos_dyn_common::atmos_dyn_filter_setup(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_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::comm_vars8_init(), scale_const::const_ohm, scale_const::const_undef, corioli, 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, and scale_index::va.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver_setup().

122  use scale_process, only: &
124  use scale_rm_process, only: &
125  prc_has_e, &
126  prc_has_w, &
127  prc_has_n, &
128  prc_has_s
129  use scale_const, only: &
130  ohm => const_ohm, &
131  undef => const_undef
132  use scale_comm, only: &
134  use scale_atmos_dyn_common, only: &
137  use scale_atmos_dyn_tinteg_short, only: &
139  use scale_atmos_dyn_tinteg_tracer, only: &
141  use scale_atmos_dyn_tinteg_large, only: &
143  use scale_atmos_dyn_tstep_short, only: &
145  use scale_atmos_dyn_tstep_tracer, only: &
147  use scale_atmos_dyn_tstep_large, only: &
149  use scale_atmos_dyn_fvm_flux, only: &
151  implicit none
152 
153  character(len=*), intent(in) :: DYN_Tinteg_Short_TYPE
154  character(len=*), intent(in) :: DYN_Tinteg_Tracer_TYPE
155  character(len=*), intent(in) :: DYN_Tinteg_Large_TYPE
156  character(len=*), intent(in) :: DYN_Tstep_Tracer_TYPE
157  character(len=*), intent(in) :: DYN_Tstep_Large_TYPE
158  character(len=*), intent(in) :: DYN_FVM_FLUX_TYPE
159  character(len=*), intent(in) :: DYN_FVM_FLUX_TYPE_TRACER
160  real(RP), intent(inout) :: DENS(KA,IA,JA) ! MPI_RECV_INIT requires intent(inout)
161  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
162  real(RP), intent(inout) :: MOMX(KA,IA,JA)
163  real(RP), intent(inout) :: MOMY(KA,IA,JA)
164  real(RP), intent(inout) :: RHOT(KA,IA,JA)
165  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
166  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
167  real(RP), intent(in) :: CDZ(KA)
168  real(RP), intent(in) :: CDX(IA)
169  real(RP), intent(in) :: CDY(JA)
170  real(RP), intent(in) :: FDZ(KA-1)
171  real(RP), intent(in) :: FDX(IA-1)
172  real(RP), intent(in) :: FDY(JA-1)
173  real(RP), intent(in) :: wdamp_tau
174  real(RP), intent(in) :: wdamp_height
175  real(RP), intent(in) :: FZ(0:KA)
176  logical, intent(in) :: enable_coriolis
177  real(RP), intent(in) :: lat(IA,JA)
178  logical, optional, intent(in) :: none
179 
180  integer :: iv, iq
181  !---------------------------------------------------------------------------
182 
183  dyn_none = .false.
184 
185  if ( present(none) ) then
186  dyn_none = none
187  endif
188 
189  allocate( i_comm_prog(max(va,1)) )
190  allocate( i_comm_qtrc(qa) )
191 
192  if ( .NOT. dyn_none ) then
193  bnd_w = .NOT. prc_has_w
194  bnd_e = .NOT. prc_has_e
195  bnd_s = .NOT. prc_has_s
196  bnd_n = .NOT. prc_has_n
197 
198  allocate( corioli(ia,ja) )
199  allocate( mflx_hi(ka,ia,ja,3) )
200  allocate( num_diff(ka,ia,ja,5,3) )
201  allocate( num_diff_q(ka,ia,ja,3) )
202  allocate( wdamp_coef(ka) )
203  mflx_hi(:,:,:,:) = undef
204 
205  call atmos_dyn_fvm_flux_setup ( dyn_fvm_flux_type, & ! [IN]
206  dyn_fvm_flux_type_tracer ) ! [IN]
207 
209 
210  call atmos_dyn_tstep_tracer_setup ( dyn_tstep_tracer_type ) ! [IN]
211 
212  call atmos_dyn_tstep_large_setup ( dyn_tstep_large_type, & ! [IN]
213  dens, momz, momx, momy, rhot, & ! [INOUT]
214  qtrc, prog, & ! [INOUT]
215  mflx_hi ) ! [INOUT]
216 
217  call atmos_dyn_tinteg_short_setup ( dyn_tinteg_short_type ) ! [IN]
218 
219  call atmos_dyn_tinteg_tracer_setup( dyn_tinteg_tracer_type ) ! [IN]
220 
221  call atmos_dyn_tinteg_large_setup ( dyn_tinteg_large_type ) ! [IN]
222 
223  ! numerical diffusion
224  call atmos_dyn_filter_setup( num_diff, num_diff_q, & ! [INOUT]
225  cdz, cdx, cdy, fdz, fdx, fdy ) ! [IN]
226 
227  ! numerical diffusion
228  call atmos_dyn_wdamp_setup( wdamp_coef(:), & ! [INOUT]
229  wdamp_tau, wdamp_height, & ! [IN]
230  fz(:) ) ! [IN]
231 
232  ! coriolis parameter
233  if ( enable_coriolis ) then
234  corioli(:,:) = 2.0_rp * ohm * sin( lat(:,:) )
235  else
236  corioli(:,:) = 0.0_rp
237  endif
238 
239  else
240 
241  call comm_vars8_init( 'DENS', dens, i_comm_dens )
242  call comm_vars8_init( 'MOMZ', momz, i_comm_momz )
243  call comm_vars8_init( 'MOMX', momx, i_comm_momx )
244  call comm_vars8_init( 'MOMY', momy, i_comm_momy )
245  call comm_vars8_init( 'RHOT', rhot, i_comm_rhot )
246 
247  do iv = 1, va
248  i_comm_prog(iv) = 5 + iv
249  call comm_vars8_init( 'PROG', prog(:,:,:,iv), i_comm_prog(iv) )
250  enddo
251 
252  do iq = 1, qa
253  i_comm_qtrc(iq) = 5 + va + iq
254  call comm_vars8_init( 'QTRC', qtrc(:,:,:,iq), i_comm_qtrc(iq) )
255  enddo
256 
257  endif
258 
259  return
logical, public prc_has_n
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_dyn_tinteg_large_setup(ATMOS_DYN_Tinteg_large_TYPE)
Register.
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm.F90:313
module Atmosphere / Dynamical scheme
module Atmosphere / Dynamics Temporal integration
module Atmosphere / Dynamics Temporal integration
subroutine, public atmos_dyn_tstep_tracer_setup(ATMOS_DYN_TSTEP_TRACER_TYPE)
Register.
logical, public prc_has_e
module Atmosphere / Dynamical scheme
module Atmosphere / Dynamics Temporal integration
logical, public prc_has_s
subroutine, public atmos_dyn_wdamp_setup(wdamp_coef, wdamp_tau, wdamp_height, FZ)
Setup.
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
subroutine, public atmos_dyn_filter_setup(num_diff, num_diff_q, CDZ, CDX, CDY, FDZ, FDX, FDY)
Setup.
module COMMUNICATION
Definition: scale_comm.F90:23
module Atmosphere / Dynamics common
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module RM PROCESS
subroutine, public atmos_dyn_tstep_large_setup(Tstep_large_type, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, mflx_hi)
Register.
subroutine, public atmos_dyn_fvm_flux_setup(scheme, scheme_tracer)
setup
module scale_atmos_dyn_fvm_flux
subroutine, public atmos_dyn_tinteg_short_setup(ATMOS_DYN_Tinteg_short_TYPE)
Register.
logical, public prc_has_w
procedure(short_setup), pointer, public atmos_dyn_tstep_short_setup
subroutine, public atmos_dyn_tinteg_tracer_setup(ATMOS_DYN_Tinteg_tracer_TYPE)
Register.
module Atmosphere / Dynamical scheme
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 (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_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)  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 285 of file scale_atmos_dyn.F90.

References scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large, scale_atmos_boundary::bnd_qa, scale_comm::comm_datatype, scale_comm::comm_world, corioli, 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_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, scale_grid_real::real_vol, scale_grid_index::xdir, scale_grid_index::ydir, and scale_grid_index::zdir.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver().

285  use scale_comm, only: &
286  comm_vars8, &
287  comm_wait
288  use scale_atmos_boundary, only: &
289  bnd_qa
290  use scale_atmos_dyn_tinteg_large, only: &
292  implicit none
293 
294  real(RP), intent(inout) :: DENS(KA,IA,JA)
295  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
296  real(RP), intent(inout) :: MOMX(KA,IA,JA)
297  real(RP), intent(inout) :: MOMY(KA,IA,JA)
298  real(RP), intent(inout) :: RHOT(KA,IA,JA)
299  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
300  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
301 
302  real(RP), intent(inout) :: DENS_av(KA,IA,JA)
303  real(RP), intent(inout) :: MOMZ_av(KA,IA,JA)
304  real(RP), intent(inout) :: MOMX_av(KA,IA,JA)
305  real(RP), intent(inout) :: MOMY_av(KA,IA,JA)
306  real(RP), intent(inout) :: RHOT_av(KA,IA,JA)
307  real(RP), intent(inout) :: QTRC_av(KA,IA,JA,QA)
308 
309  real(RP), intent(in) :: DENS_tp(KA,IA,JA)
310  real(RP), intent(in) :: MOMZ_tp(KA,IA,JA)
311  real(RP), intent(in) :: MOMX_tp(KA,IA,JA)
312  real(RP), intent(in) :: MOMY_tp(KA,IA,JA)
313  real(RP), intent(in) :: RHOT_tp(KA,IA,JA)
314  real(RP), intent(in) :: RHOQ_tp(KA,IA,JA,QA)
315 
316  real(RP), intent(in) :: CDZ (KA)
317  real(RP), intent(in) :: CDX (IA)
318  real(RP), intent(in) :: CDY (JA)
319  real(RP), intent(in) :: FDZ (KA-1)
320  real(RP), intent(in) :: FDX (IA-1)
321  real(RP), intent(in) :: FDY (JA-1)
322  real(RP), intent(in) :: RCDZ(KA)
323  real(RP), intent(in) :: RCDX(IA)
324  real(RP), intent(in) :: RCDY(JA)
325  real(RP), intent(in) :: RFDZ(KA-1)
326  real(RP), intent(in) :: RFDX(IA-1)
327  real(RP), intent(in) :: RFDY(JA-1)
328 
329  real(RP), intent(in) :: PHI (KA,IA,JA)
330  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
331  real(RP), intent(in) :: J13G (KA,IA,JA,7)
332  real(RP), intent(in) :: J23G (KA,IA,JA,7)
333  real(RP), intent(in) :: J33G
334  real(RP), intent(in) :: MAPF (IA,JA,2,4)
335 
336  real(RP), intent(in) :: AQ_R (QA)
337  real(RP), intent(in) :: AQ_CV (QA)
338  real(RP), intent(in) :: AQ_CP (QA)
339  real(RP), intent(in) :: AQ_MASS(QA)
340 
341  real(RP), intent(in) :: REF_dens(KA,IA,JA)
342  real(RP), intent(in) :: REF_pott(KA,IA,JA)
343  real(RP), intent(in) :: REF_qv (KA,IA,JA)
344  real(RP), intent(in) :: REF_pres(KA,IA,JA)
345  real(RP), intent(in) :: ND_COEF
346  real(RP), intent(in) :: ND_COEF_Q
347  integer, intent(in) :: ND_ORDER
348  real(RP), intent(in) :: ND_SFC_FACT
349  logical, intent(in) :: ND_USE_RS
350 
351  real(RP), intent(in) :: DAMP_DENS(KA,IA,JA)
352  real(RP), intent(in) :: DAMP_VELZ(KA,IA,JA)
353  real(RP), intent(in) :: DAMP_VELX(KA,IA,JA)
354  real(RP), intent(in) :: DAMP_VELY(KA,IA,JA)
355  real(RP), intent(in) :: DAMP_POTT(KA,IA,JA)
356  real(RP), intent(in) :: DAMP_QTRC(KA,IA,JA,BND_QA)
357 
358  real(RP), intent(in) :: DAMP_alpha_DENS(KA,IA,JA)
359  real(RP), intent(in) :: DAMP_alpha_VELZ(KA,IA,JA)
360  real(RP), intent(in) :: DAMP_alpha_VELX(KA,IA,JA)
361  real(RP), intent(in) :: DAMP_alpha_VELY(KA,IA,JA)
362  real(RP), intent(in) :: DAMP_alpha_POTT(KA,IA,JA)
363  real(RP), intent(in) :: DAMP_alpha_QTRC(KA,IA,JA,BND_QA)
364 
365  real(RP), intent(in) :: divdmp_coef
366 
367  logical, intent(in) :: FLAG_FCT_MOMENTUM
368  logical, intent(in) :: FLAG_FCT_T
369  logical, intent(in) :: FLAG_FCT_TRACER
370  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
371 
372  logical, intent(in) :: USE_AVERAGE
373 
374  real(DP), intent(in) :: DTSEC
375  real(DP), intent(in) :: DTSEC_DYN
376 
377  ! for time integartion
378  real(RP) :: DENS00 (KA,IA,JA) ! saved density before update
379 
380  ! For tracer advection
381  real(RP) :: tflx_hi(KA,IA,JA,3) ! rho * theta * vel(x,y,z) @ (u,v,w)-face high order
382  real(RP) :: dt
383 
384  integer :: i, j, k, iq
385  !---------------------------------------------------------------------------
386 
387  if( io_l ) write(io_fid_log,*) '*** Atmos dynamics step'
388 
389  dt = real(DTSEC, kind=rp)
390 
391  if ( dyn_none ) then
392 
393  call prof_rapstart("DYN_Tinteg", 2)
394 
395  do j = js, je
396  do i = is, ie
397  do k = ks, ke
398  dens00(k,i,j) = dens(k,i,j) ! save previous value
399 
400  dens(k,i,j) = dens(k,i,j) + dens_tp(k,i,j) * dt
401  enddo
402  enddo
403  enddo
404 
405  do j = js, je
406  do i = is, ie
407  do k = ks, ke
408  momz(k,i,j) = momz(k,i,j) + momz_tp(k,i,j) * dt
409  enddo
410  enddo
411  enddo
412 
413  do j = js, je
414  do i = is, ie
415  do k = ks, ke
416  momx(k,i,j) = momx(k,i,j) + momx_tp(k,i,j) * dt
417  enddo
418  enddo
419  enddo
420 
421  do j = js, je
422  do i = is, ie
423  do k = ks, ke
424  momy(k,i,j) = momy(k,i,j) + momy_tp(k,i,j) * dt
425  enddo
426  enddo
427  enddo
428 
429  do j = js, je
430  do i = is, ie
431  do k = ks, ke
432  rhot(k,i,j) = rhot(k,i,j) + rhot_tp(k,i,j) * dt
433  enddo
434  enddo
435  enddo
436 
437  do iq = 1, qa
438  do j = js, je
439  do i = is, ie
440  do k = ks, ke
441  qtrc(k,i,j,iq) = max( 0.0_rp, &
442  qtrc(k,i,j,iq) * dens00(k,i,j) + rhoq_tp(k,i,j,iq) * dt &
443  ) / dens(k,i,j)
444  enddo
445  enddo
446  enddo
447  enddo
448 
449  call comm_vars8( dens(:,:,:), i_comm_dens )
450  call comm_vars8( momz(:,:,:), i_comm_momz )
451  call comm_vars8( momx(:,:,:), i_comm_momx )
452  call comm_vars8( momy(:,:,:), i_comm_momy )
453  call comm_vars8( rhot(:,:,:), i_comm_rhot )
454  do iq = 1, qa
455  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
456  enddo
457 
458  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
459  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
460  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
461  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
462  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
463  do iq = 1, qa
464  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
465  enddo
466 
467  if ( use_average ) then
468  dens_av(:,:,:) = dens(:,:,:)
469  momz_av(:,:,:) = momz(:,:,:)
470  momx_av(:,:,:) = momx(:,:,:)
471  momy_av(:,:,:) = momy(:,:,:)
472  rhot_av(:,:,:) = rhot(:,:,:)
473  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
474  endif
475 
476  call prof_rapend("DYN_Tinteg", 2)
477 
478  return
479  endif ! if DYN_NONE == .true.
480 
481  call prof_rapstart("DYN_Tinteg", 2)
482 
483  call atmos_dyn_tinteg_large( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
484  prog, & ! [INOUT]
485  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
486  mflx_hi, tflx_hi, & ! [OUT]
487  num_diff, num_diff_q, & ! [OUT;WORK]
488  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! [IN]
489  corioli, & ! [IN]
490  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
491  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
492  phi, gsqrt, & ! [IN]
493  j13g, j23g, j33g, mapf, & ! [IN]
494  aq_r, aq_cv, aq_cp, aq_mass, & ! [IN]
495  ref_dens, ref_pott, ref_qv, ref_pres, & ! [IN]
496  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
497  nd_coef, nd_coef_q, nd_order, nd_sfc_fact, nd_use_rs, & ! [IN]
498  damp_dens, damp_velz, damp_velx, & ! [IN]
499  damp_vely, damp_pott, damp_qtrc, & ! [IN]
500  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! [IN]
501  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! [IN]
502  wdamp_coef, & ! [IN]
503  divdmp_coef, & ! [IN]
504  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! [IN]
505  flag_fct_along_stream, & ! [IN]
506  use_average, & ! [IN]
507  dtsec, dtsec_dyn ) ! [IN]
508 
509  call prof_rapend ("DYN_Tinteg", 2)
510 
511 #ifdef CHECK_MASS
512  call check_mass( dens, damp_dens, & ! [IN]
513  mflx_hi, tflx_hi, & ! [IN]
514  gsqrt, mapf, & ! [IN]
515  rcdx, rcdy, & ! [IN]
516  dt, & ! [IN]
517  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
518 #endif
519 
520  return
module Atmosphere / Dynamics Temporal integration
procedure(large), pointer, public atmos_dyn_tinteg_large
module COMMUNICATION
Definition: scale_comm.F90:23
module ATMOSPHERE / Boundary treatment
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ corioli

real(rp), dimension (:,:), allocatable, public scale_atmos_dyn::corioli

Definition at line 68 of file scale_atmos_dyn.F90.

Referenced by atmos_dyn(), and atmos_dyn_setup().

68  real(RP), public, allocatable :: CORIOLI (:,:) ! coriolis term