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, coriolis_type, coriolis_f0, coriolis_beta, coriolis_y0, CY, 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, BND_QA, 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, 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...
 

Variables

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

Detailed Description

module Atmosphere / Dynamics FENT + FCT

Description
Dynamical core for Atmospheric process Full explicit, no terrain + tracer FCT limiter
Author
Team SCALE
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 (w/ CHECK_MASS) kg/m2/s mflx_hi
MFLXY momentum flux of y-direction (w/ CHECK_MASS) kg/m2/s mflx_hi
MFLXZ momentum flux of z-direction (w/ CHECK_MASS) kg/m2/s mflx_hi
TFLXX potential temperature flux of x-direction (w/ CHECK_MASS) K*kg/m2/s tflx_hi
TFLXY potential temperature flux of y-direction (w/ CHECK_MASS) K*kg/m2/s tflx_hi
TFLXZ potential temperature flux of z-direction (w/ CHECK_MASS) 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,
character(len=*), intent(in)  coriolis_type,
real(rp), intent(in)  coriolis_f0,
real(rp), intent(in)  coriolis_beta,
real(rp), intent(in)  coriolis_y0,
real(rp), dimension(ja), intent(in)  CY,
real(rp), dimension(ia,ja), intent(in)  lat,
logical, intent(in), optional  none 
)

Setup.

Definition at line 103 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_cartesc::comm_vars8_init(), scale_const::const_ohm, scale_const::const_undef, coriolis, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::ka, 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_tracer::qa, and scale_index::va.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver_setup().

103  use scale_prc, only: &
104  prc_abort
105  use scale_prc_cartesc, only: &
106  prc_has_e, &
107  prc_has_w, &
108  prc_has_n, &
109  prc_has_s
110  use scale_const, only: &
111  ohm => const_ohm, &
112  undef => const_undef
113  use scale_comm_cartesc, only: &
115  use scale_atmos_dyn_common, only: &
118  use scale_atmos_dyn_tinteg_short, only: &
120  use scale_atmos_dyn_tinteg_tracer, only: &
122  use scale_atmos_dyn_tinteg_large, only: &
124  use scale_atmos_dyn_tstep_short, only: &
126  use scale_atmos_dyn_tstep_tracer, only: &
128  use scale_atmos_dyn_tstep_large, only: &
130  use scale_atmos_dyn_fvm_flux, only: &
132  implicit none
133 
134  character(len=*), intent(in) :: dyn_tinteg_short_type
135  character(len=*), intent(in) :: dyn_tinteg_tracer_type
136  character(len=*), intent(in) :: dyn_tinteg_large_type
137  character(len=*), intent(in) :: dyn_tstep_tracer_type
138  character(len=*), intent(in) :: dyn_tstep_large_type
139  character(len=*), intent(in) :: dyn_fvm_flux_type
140  character(len=*), intent(in) :: dyn_fvm_flux_type_tracer
141  real(RP), intent(inout) :: dens(ka,ia,ja) ! MPI_RECV_INIT requires intent(inout)
142  real(RP), intent(inout) :: momz(ka,ia,ja)
143  real(RP), intent(inout) :: momx(ka,ia,ja)
144  real(RP), intent(inout) :: momy(ka,ia,ja)
145  real(RP), intent(inout) :: rhot(ka,ia,ja)
146  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
147  real(RP), intent(inout) :: prog(ka,ia,ja,va)
148  real(RP), intent(in) :: cdz(ka)
149  real(RP), intent(in) :: cdx(ia)
150  real(RP), intent(in) :: cdy(ja)
151  real(RP), intent(in) :: fdz(ka-1)
152  real(RP), intent(in) :: fdx(ia-1)
153  real(RP), intent(in) :: fdy(ja-1)
154  real(RP), intent(in) :: wdamp_tau
155  real(RP), intent(in) :: wdamp_height
156  real(RP), intent(in) :: fz(0:ka)
157  character(len=*), intent(in) :: coriolis_type
158  real(RP), intent(in) :: coriolis_f0
159  real(RP), intent(in) :: coriolis_beta
160  real(RP), intent(in) :: coriolis_y0
161  real(RP), intent(in) :: cy(ja)
162  real(RP), intent(in) :: lat(ia,ja)
163  logical, optional, intent(in) :: none
164 
165  integer :: j
166  integer :: iv, iq
167  !---------------------------------------------------------------------------
168 
169  dyn_none = .false.
170 
171  if ( present(none) ) then
172  dyn_none = none
173  endif
174 
175  allocate( i_comm_prog(max(va,1)) )
176  allocate( i_comm_qtrc(qa) )
177 
178  if ( .NOT. dyn_none ) then
179  bnd_w = .NOT. prc_has_w
180  bnd_e = .NOT. prc_has_e
181  bnd_s = .NOT. prc_has_s
182  bnd_n = .NOT. prc_has_n
183 
184  allocate( coriolis(ia,ja) )
185  allocate( mflx_hi(ka,ia,ja,3) )
186  allocate( num_diff(ka,ia,ja,5,3) )
187  allocate( num_diff_q(ka,ia,ja,3) )
188  allocate( wdamp_coef(ka) )
189  mflx_hi(:,:,:,:) = undef
190  num_diff(:,:,:,:,:) = undef
191  num_diff_q(:,:,:,:) = undef
192 
193  call atmos_dyn_fvm_flux_setup ( dyn_fvm_flux_type, & ! [IN]
194  dyn_fvm_flux_type_tracer ) ! [IN]
195 
197 
198  call atmos_dyn_tstep_tracer_setup ( dyn_tstep_tracer_type ) ! [IN]
199 
200  call atmos_dyn_tstep_large_setup ( dyn_tstep_large_type, & ! [IN]
201  dens, momz, momx, momy, rhot, & ! [INOUT]
202  qtrc, prog, & ! [INOUT]
203  mflx_hi ) ! [INOUT]
204 
205  call atmos_dyn_tinteg_short_setup ( dyn_tinteg_short_type ) ! [IN]
206 
207  call atmos_dyn_tinteg_tracer_setup( dyn_tinteg_tracer_type ) ! [IN]
208 
209  call atmos_dyn_tinteg_large_setup ( dyn_tinteg_large_type ) ! [IN]
210 
211  ! numerical diffusion
212  call atmos_dyn_filter_setup( num_diff, num_diff_q, & ! [INOUT]
213  cdz, cdx, cdy, fdz, fdx, fdy ) ! [IN]
214 
215  ! numerical diffusion
216  call atmos_dyn_wdamp_setup( wdamp_coef(:), & ! [INOUT]
217  wdamp_tau, wdamp_height, & ! [IN]
218  fz(:) ) ! [IN]
219 
220  ! coriolis parameter
221  select case ( coriolis_type )
222  case ( 'PLANE' )
223  do j = 1, ja
224  coriolis(:,j) = coriolis_f0 + coriolis_beta * ( cy(j) - coriolis_y0 )
225  end do
226  case ( 'SPHERE' )
227  coriolis(:,:) = 2.0_rp * ohm * sin( lat(:,:) )
228  case default
229  log_error("ATMOS_DYN_setup",*) 'Coriolis type is invalid: ', trim(coriolis_type)
230  log_error_cont(*) 'The type must be PLANE or SPHERE'
231  call prc_abort
232  end select
233 
234  else
235 
236  call comm_vars8_init( 'DENS', dens, i_comm_dens )
237  call comm_vars8_init( 'MOMZ', momz, i_comm_momz )
238  call comm_vars8_init( 'MOMX', momx, i_comm_momx )
239  call comm_vars8_init( 'MOMY', momy, i_comm_momy )
240  call comm_vars8_init( 'RHOT', rhot, i_comm_rhot )
241 
242  do iv = 1, va
243  i_comm_prog(iv) = 5 + iv
244  call comm_vars8_init( 'PROG', prog(:,:,:,iv), i_comm_prog(iv) )
245  enddo
246 
247  do iq = 1, qa
248  i_comm_qtrc(iq) = 5 + va + iq
249  call comm_vars8_init( 'QTRC', qtrc(:,:,:,iq), i_comm_qtrc(iq) )
250  enddo
251 
252  endif
253 
254  return
subroutine, public atmos_dyn_tinteg_large_setup(ATMOS_DYN_Tinteg_large_TYPE)
Register.
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.
module process / cartesC
module Atmosphere / Dynamical scheme
module Atmosphere / Dynamics Temporal integration
logical, public prc_has_s
logical, public prc_has_n
logical, public prc_has_e
subroutine, public atmos_dyn_wdamp_setup(wdamp_coef, wdamp_tau, wdamp_height, FZ)
Setup.
real(rp), public const_undef
Definition: scale_const.F90:41
module COMMUNICATION
real(rp), public const_ohm
angular velocity of the planet [1/s]
Definition: scale_const.F90:45
module PROCESS
Definition: scale_prc.F90:11
subroutine, public atmos_dyn_filter_setup(num_diff, num_diff_q, CDZ, CDX, CDY, FDZ, FDX, FDY)
Setup.
module Atmosphere / Dynamics common
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
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.
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
logical, public prc_has_w
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,
integer, intent(in)  BND_QA,
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), 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 283 of file scale_atmos_dyn.F90.

References scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol, scale_comm_cartesc::comm_datatype, scale_comm_cartesc::comm_world, coriolis, scale_atmos_grid_cartesc_index::i_xy, scale_atmos_grid_cartesc_index::i_xyz, 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_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::kmax, scale_atmos_grid_cartesc_index::ks, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, scale_atmos_grid_cartesc_index::xdir, scale_atmos_grid_cartesc_index::ydir, and scale_atmos_grid_cartesc_index::zdir.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver().

283  use scale_comm_cartesc, only: &
284  comm_vars8, &
285  comm_wait
286  use scale_atmos_dyn_tinteg_large, only: &
288  implicit none
289 
290  real(RP), intent(inout) :: dens(ka,ia,ja)
291  real(RP), intent(inout) :: momz(ka,ia,ja)
292  real(RP), intent(inout) :: momx(ka,ia,ja)
293  real(RP), intent(inout) :: momy(ka,ia,ja)
294  real(RP), intent(inout) :: rhot(ka,ia,ja)
295  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
296  real(RP), intent(inout) :: prog(ka,ia,ja,va)
297 
298  real(RP), intent(inout) :: dens_av(ka,ia,ja)
299  real(RP), intent(inout) :: momz_av(ka,ia,ja)
300  real(RP), intent(inout) :: momx_av(ka,ia,ja)
301  real(RP), intent(inout) :: momy_av(ka,ia,ja)
302  real(RP), intent(inout) :: rhot_av(ka,ia,ja)
303  real(RP), intent(inout) :: qtrc_av(ka,ia,ja,qa)
304 
305  real(RP), intent(in) :: dens_tp(ka,ia,ja)
306  real(RP), intent(in) :: momz_tp(ka,ia,ja)
307  real(RP), intent(in) :: momx_tp(ka,ia,ja)
308  real(RP), intent(in) :: momy_tp(ka,ia,ja)
309  real(RP), intent(in) :: rhot_tp(ka,ia,ja)
310  real(RP), intent(in) :: rhoq_tp(ka,ia,ja,qa)
311 
312  real(RP), intent(in) :: cdz (ka)
313  real(RP), intent(in) :: cdx (ia)
314  real(RP), intent(in) :: cdy (ja)
315  real(RP), intent(in) :: fdz (ka-1)
316  real(RP), intent(in) :: fdx (ia-1)
317  real(RP), intent(in) :: fdy (ja-1)
318  real(RP), intent(in) :: rcdz(ka)
319  real(RP), intent(in) :: rcdx(ia)
320  real(RP), intent(in) :: rcdy(ja)
321  real(RP), intent(in) :: rfdz(ka-1)
322  real(RP), intent(in) :: rfdx(ia-1)
323  real(RP), intent(in) :: rfdy(ja-1)
324 
325  real(RP), intent(in) :: phi (ka,ia,ja)
326  real(RP), intent(in) :: gsqrt(ka,ia,ja,7)
327  real(RP), intent(in) :: j13g (ka,ia,ja,7)
328  real(RP), intent(in) :: j23g (ka,ia,ja,7)
329  real(RP), intent(in) :: j33g
330  real(RP), intent(in) :: mapf (ia,ja,2,4)
331 
332  real(RP), intent(in) :: aq_r (qa)
333  real(RP), intent(in) :: aq_cv (qa)
334  real(RP), intent(in) :: aq_cp (qa)
335  real(RP), intent(in) :: aq_mass(qa)
336 
337  real(RP), intent(in) :: ref_dens(ka,ia,ja)
338  real(RP), intent(in) :: ref_pott(ka,ia,ja)
339  real(RP), intent(in) :: ref_qv (ka,ia,ja)
340  real(RP), intent(in) :: ref_pres(ka,ia,ja)
341  real(RP), intent(in) :: nd_coef
342  real(RP), intent(in) :: nd_coef_q
343  integer, intent(in) :: nd_order
344  real(RP), intent(in) :: nd_sfc_fact
345  logical, intent(in) :: nd_use_rs
346 
347  integer, intent(in) :: bnd_qa
348  real(RP), intent(in) :: bnd_smoother_fact
349 
350  real(RP), intent(in) :: damp_dens(ka,ia,ja)
351  real(RP), intent(in) :: damp_velz(ka,ia,ja)
352  real(RP), intent(in) :: damp_velx(ka,ia,ja)
353  real(RP), intent(in) :: damp_vely(ka,ia,ja)
354  real(RP), intent(in) :: damp_pott(ka,ia,ja)
355  real(RP), intent(in) :: damp_qtrc(ka,ia,ja,bnd_qa)
356 
357  real(RP), intent(in) :: damp_alpha_dens(ka,ia,ja)
358  real(RP), intent(in) :: damp_alpha_velz(ka,ia,ja)
359  real(RP), intent(in) :: damp_alpha_velx(ka,ia,ja)
360  real(RP), intent(in) :: damp_alpha_vely(ka,ia,ja)
361  real(RP), intent(in) :: damp_alpha_pott(ka,ia,ja)
362  real(RP), intent(in) :: damp_alpha_qtrc(ka,ia,ja,bnd_qa)
363 
364  real(RP), intent(in) :: divdmp_coef
365 
366  logical, intent(in) :: flag_tracer_split_tend
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  integer, intent(in) :: i_qv
375 
376  real(DP), intent(in) :: dtsec
377  real(DP), intent(in) :: dtsec_dyn
378 
379  ! for time integartion
380  real(RP) :: dens00 (ka,ia,ja) ! saved density before update
381 
382  ! For tracer advection
383  real(RP) :: tflx_hi(ka,ia,ja,3) ! rho * theta * vel(x,y,z) @ (u,v,w)-face high order
384  real(RP) :: dt
385 
386  integer :: i, j, k, iq
387  !---------------------------------------------------------------------------
388 
389  log_progress(*) 'atmosphere / dynamics'
390 
391  dt = real(dtsec, kind=rp)
392 
393  if ( dyn_none ) then
394 
395  call prof_rapstart("DYN_Tinteg", 2)
396 
397  do j = js, je
398  do i = is, ie
399  do k = ks, ke
400  dens00(k,i,j) = dens(k,i,j) ! save previous value
401 
402  dens(k,i,j) = dens(k,i,j) + dens_tp(k,i,j) * dt
403  enddo
404  enddo
405  enddo
406 
407  do j = js, je
408  do i = is, ie
409  do k = ks, ke
410  momz(k,i,j) = momz(k,i,j) + momz_tp(k,i,j) * dt
411  enddo
412  enddo
413  enddo
414 
415  do j = js, je
416  do i = is, ie
417  do k = ks, ke
418  momx(k,i,j) = momx(k,i,j) + momx_tp(k,i,j) * dt
419  enddo
420  enddo
421  enddo
422 
423  do j = js, je
424  do i = is, ie
425  do k = ks, ke
426  momy(k,i,j) = momy(k,i,j) + momy_tp(k,i,j) * dt
427  enddo
428  enddo
429  enddo
430 
431  do j = js, je
432  do i = is, ie
433  do k = ks, ke
434  rhot(k,i,j) = rhot(k,i,j) + rhot_tp(k,i,j) * dt
435  enddo
436  enddo
437  enddo
438 
439  do iq = 1, qa
440  do j = js, je
441  do i = is, ie
442  do k = ks, ke
443  qtrc(k,i,j,iq) = max( 0.0_rp, &
444  qtrc(k,i,j,iq) * dens00(k,i,j) + rhoq_tp(k,i,j,iq) * dt &
445  ) / dens(k,i,j)
446  enddo
447  enddo
448  enddo
449  enddo
450 
451  call comm_vars8( dens(:,:,:), i_comm_dens )
452  call comm_vars8( momz(:,:,:), i_comm_momz )
453  call comm_vars8( momx(:,:,:), i_comm_momx )
454  call comm_vars8( momy(:,:,:), i_comm_momy )
455  call comm_vars8( rhot(:,:,:), i_comm_rhot )
456  do iq = 1, qa
457  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
458  enddo
459 
460  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
461  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
462  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
463  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
464  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
465  do iq = 1, qa
466  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
467  enddo
468 
469  if ( use_average ) then
470  dens_av(:,:,:) = dens(:,:,:)
471  momz_av(:,:,:) = momz(:,:,:)
472  momx_av(:,:,:) = momx(:,:,:)
473  momy_av(:,:,:) = momy(:,:,:)
474  rhot_av(:,:,:) = rhot(:,:,:)
475  qtrc_av(:,:,:,:) = qtrc(:,:,:,:)
476  endif
477 
478  call prof_rapend("DYN_Tinteg", 2)
479 
480  return
481  endif ! if DYN_NONE == .true.
482 
483  call prof_rapstart("DYN_Tinteg", 2)
484 
485  call atmos_dyn_tinteg_large( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
486  prog, & ! [INOUT]
487  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
488  mflx_hi, tflx_hi, & ! [OUT]
489  num_diff, num_diff_q, & ! [OUT;WORK]
490  dens_tp, momz_tp, momx_tp, momy_tp, rhot_tp, rhoq_tp, & ! [IN]
491  coriolis, & ! [IN]
492  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
493  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
494  phi, gsqrt, & ! [IN]
495  j13g, j23g, j33g, mapf, & ! [IN]
496  aq_r, aq_cv, aq_cp, aq_mass, & ! [IN]
497  ref_dens, ref_pott, ref_qv, ref_pres, & ! [IN]
498  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
499  nd_coef, nd_coef_q, nd_order, nd_sfc_fact, nd_use_rs, & ! [IN]
500  bnd_qa, bnd_smoother_fact, & ! [IN]
501  damp_dens, damp_velz, damp_velx, & ! [IN]
502  damp_vely, damp_pott, damp_qtrc, & ! [IN]
503  damp_alpha_dens, damp_alpha_velz, damp_alpha_velx, & ! [IN]
504  damp_alpha_vely, damp_alpha_pott, damp_alpha_qtrc, & ! [IN]
505  wdamp_coef, & ! [IN]
506  divdmp_coef, & ! [IN]
507  flag_tracer_split_tend, & ! [IN]
508  flag_fct_momentum, flag_fct_t, flag_fct_tracer, & ! [IN]
509  flag_fct_along_stream, & ! [IN]
510  use_average, & ! [IN]
511  i_qv, & ! [IN]
512  dtsec, dtsec_dyn ) ! [IN]
513 
514  call prof_rapend ("DYN_Tinteg", 2)
515 
516 #ifdef CHECK_MASS
517  call check_mass( dens, damp_dens, & ! [IN]
518  mflx_hi, tflx_hi, & ! [IN]
519  gsqrt, mapf, & ! [IN]
520  rcdx, rcdy, & ! [IN]
521  dt, & ! [IN]
522  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
523 #endif
524 
525  return
module Atmosphere / Dynamics Temporal integration
module COMMUNICATION
procedure(large), pointer, public atmos_dyn_tinteg_large
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ coriolis

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

Definition at line 45 of file scale_atmos_dyn.F90.

Referenced by atmos_dyn(), and atmos_dyn_setup().

45  real(RP), public, allocatable :: coriolis(:,:) ! coriolis term