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_FVM_FLUX_TYPE, DYN_FVM_FLUX_TYPE_TRACER, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, CDZ, CDX, CDY, FDZ, FDX, FDY, 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_CV, 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...
 

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

Setup.

Definition at line 117 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_comm::comm_vars8_init(), scale_const::const_ohm, scale_const::const_undef, scale_grid_index::ia, scale_grid_index::ja, scale_grid_index::ka, scale_rm_process::prc_has_e, scale_rm_process::prc_has_n, scale_rm_process::prc_has_s, scale_rm_process::prc_has_w, scale_process::prc_mpistop(), scale_tracer::qa, and scale_index::va.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver_setup().

117  use scale_process, only: &
119  use scale_rm_process, only: &
120  prc_has_e, &
121  prc_has_w, &
122  prc_has_n, &
123  prc_has_s
124  use scale_const, only: &
125  ohm => const_ohm, &
126  undef => const_undef
127  use scale_comm, only: &
129  use scale_atmos_dyn_common, only: &
131  use scale_atmos_dyn_tinteg_short, only: &
133  use scale_atmos_dyn_tinteg_tracer, only: &
135  use scale_atmos_dyn_tinteg_large, only: &
137  use scale_atmos_dyn_tstep_short, only: &
139  use scale_atmos_dyn_tstep_tracer, only: &
141  use scale_atmos_dyn_tstep_large, only: &
143  use scale_atmos_dyn_fvm_flux, only: &
145  implicit none
146 
147  character(len=*), intent(in) :: dyn_tinteg_short_type
148  character(len=*), intent(in) :: dyn_tinteg_tracer_type
149  character(len=*), intent(in) :: dyn_tinteg_large_type
150  character(len=*), intent(in) :: dyn_tstep_tracer_type
151  character(len=*), intent(in) :: dyn_tstep_large_type
152  character(len=*), intent(in) :: dyn_fvm_flux_type
153  character(len=*), intent(in) :: dyn_fvm_flux_type_tracer
154  real(RP), intent(inout) :: dens(ka,ia,ja) ! MPI_RECV_INIT requires intent(inout)
155  real(RP), intent(inout) :: momz(ka,ia,ja)
156  real(RP), intent(inout) :: momx(ka,ia,ja)
157  real(RP), intent(inout) :: momy(ka,ia,ja)
158  real(RP), intent(inout) :: rhot(ka,ia,ja)
159  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
160  real(RP), intent(inout) :: prog(ka,ia,ja,va)
161  real(RP), intent(in) :: cdz(ka)
162  real(RP), intent(in) :: cdx(ia)
163  real(RP), intent(in) :: cdy(ja)
164  real(RP), intent(in) :: fdz(ka-1)
165  real(RP), intent(in) :: fdx(ia-1)
166  real(RP), intent(in) :: fdy(ja-1)
167  logical, intent(in) :: enable_coriolis
168  real(RP), intent(in) :: lat(ia,ja)
169  logical, optional, intent(in) :: none
170 
171  integer :: iv, iq
172  !---------------------------------------------------------------------------
173 
174  if ( present(none) ) then
175  dyn_none = none
176  endif
177 
178  allocate( i_comm_prog(max(va,1)) )
179  allocate( i_comm_qtrc(qa) )
180 
181  if ( .NOT. dyn_none ) then
182  bnd_w = .NOT. prc_has_w
183  bnd_e = .NOT. prc_has_e
184  bnd_s = .NOT. prc_has_s
185  bnd_n = .NOT. prc_has_n
186 
187  allocate( corioli(ia,ja) )
188  allocate( mflx_hi(ka,ia,ja,3) )
189  allocate( num_diff(ka,ia,ja,5,3) )
190  allocate( num_diff_q(ka,ia,ja,3) )
191  mflx_hi(:,:,:,:) = 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  ! coriolis parameter
216  if ( enable_coriolis ) then
217  corioli(:,:) = 2.0_rp * ohm * sin( lat(:,:) )
218  else
219  corioli(:,:) = 0.0_rp
220  endif
221 
222  else
223 
224  call comm_vars8_init( dens, i_comm_dens )
225  call comm_vars8_init( momz, i_comm_momz )
226  call comm_vars8_init( momx, i_comm_momx )
227  call comm_vars8_init( momy, i_comm_momy )
228  call comm_vars8_init( rhot, i_comm_rhot )
229 
230  do iv = 1, va
231  i_comm_prog(iv) = 5 + iv
232  call comm_vars8_init( prog(:,:,:,iv), i_comm_prog(iv) )
233  enddo
234 
235  do iq = 1, qa
236  i_comm_qtrc(iq) = 5 + va + iq
237  call comm_vars8_init( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
238  enddo
239 
240  endif
241 
242  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.
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
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 comm_vars8_init(var, vid)
Register variables.
Definition: scale_comm.F90:328
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(qqa), intent(in)  AQ_CV,
real(rp), dimension(ka,ia,ja), intent(in)  REF_dens,
real(rp), dimension(ka,ia,ja), intent(in)  REF_pott,
real(rp), dimension (ka,ia,ja), intent(in)  REF_qv,
real(rp), dimension(ka,ia,ja), intent(in)  REF_pres,
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 268 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, 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().

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