SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tstep_large_fvm_heve Module Reference

module Atmosphere / Dynamics More...

Functions/Subroutines

subroutine, public atmos_dyn_tstep_large_fvm_heve_setup (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG)
 Setup. More...
 
subroutine, public atmos_dyn_tstep_large_fvm_heve (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, num_diff, num_diff_q, QTRC0, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLI, 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, BND_W, BND_E, BND_S, BND_N, TwoD, 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, wdamp_coef, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTLS, DTSS, Llast)
 Dynamical Process. More...
 
subroutine monitor_put_lateral_flux (MONIT_ID, BND_flag, flx, flx_zero)
 
subroutine multiply_flux_by_metric_xyz (I_DIR, flx, GSQRT, MAPF)
 

Detailed Description

module Atmosphere / Dynamics

Description
HEVE FVM scheme for large time step in Atmospheric dynamical process
Author
Team SCALE
NAMELIST
  • ATMOS_DYN_TSTEP_LARGE_FVM_HEVE
    nametypedefault valuecomment
    EVAL_TYPE_NUMFILTER character(len=H_SHORT) 'TENDENCY'

History Output
namedescriptionunitvariable
XFLX_{TRACER_NAME} {TRACER_NAME} flux of x-direction;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s {'XFLX_'//trim}
YFLX_{TRACER_NAME} {TRACER_NAME} flux of y-direction;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s {'YFLX_'//trim}
ZFLX_{TRACER_NAME} {TRACER_NAME} flux of z-direction;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s {'ZFLX_'//trim}
DENS_t_damp tendency of dencity due to damping kg/m3/s DENS_t_damp
DENS_t_phys tendency of dencity due to physics kg/m3/s DENS_t_phys
MOMX_t_damp tendency of momentum x due to damping kg/m2/s2 MOMX_t_damp
MOMX_t_phys tendency of momentum x due to physics kg/m2/s2 MOMX_t_phys
MOMY_t_damp tendency of momentum y due to damping kg/m2/s2 MOMY_t_damp
MOMY_t_phys tendency of momentum y due to physics kg/m2/s2 MOMY_t_phys
MOMZ_t_damp tendency of momentum z due to damping kg/m2/s2 MOMZ_t_damp
MOMZ_t_phys tendency of momentum z due to physics kg/m2/s2 MOMZ_t_phys
RHOT_t_damp tendency of rho*theta temperature due to damping K kg/m3/s RHOT_t_damp
RHOT_t_phys tendency of rho*theta temperature due to physics K*kg/m3/s RHOT_t_phys
XFLX_MOM momentum flux of x-direction kg/m2/s XFLX_MOM
XFLX_RHOT potential temperature flux of x-direction K*kg/m2/s XFLX_RHOT
YFLX_MOM momentum flux of y-direction kg/m2/s YFLX_MOM
YFLX_RHOT potential temperature flux of y-direction K*kg/m2/s YFLX_RHOT
ZFLX_MOM momentum flux of z-direction kg/m2/s ZFLX_MOM
ZFLX_RHOT potential temperature flux of z-direction K*kg/m2/s ZFLX_RHOT
{TRACER_NAME}_t_damp tendency of {TRACER_NAME} due to damping;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
kg/m3/s trim(TRACER_NAME(iq))//'_t_damp
{TRACER_NAME}_t_phys tendency of {TRACER_NAME} due to physics;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
kg/m3/s trim(TRACER_NAME(iq))//'_t_phys

Function/Subroutine Documentation

◆ atmos_dyn_tstep_large_fvm_heve_setup()

subroutine, public scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve_setup ( 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 
)

Setup.

Definition at line 123 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

123  use scale_prc, only: &
124  prc_abort
125  use scale_prc_cartesc, only: &
126  prc_twod, &
127  prc_has_e, &
128  prc_has_w, &
129  prc_has_n, &
130  prc_has_s
131  use scale_const, only: &
132  ohm => const_ohm, &
133  undef => const_undef
134  use scale_comm_cartesc, only: &
136  use scale_file_history, only: &
138  file_history_put
139  use scale_monitor, only: &
140  monitor_reg, &
141  monitor_put
142  implicit none
143 
144  ! MPI_RECV_INIT requires intent(inout)
145  real(RP), intent(inout) :: DENS(KA,IA,JA)
146  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
147  real(RP), intent(inout) :: MOMX(KA,IA,JA)
148  real(RP), intent(inout) :: MOMY(KA,IA,JA)
149  real(RP), intent(inout) :: RHOT(KA,IA,JA)
150  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
151  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
152 
153  integer :: iv, iq
154 
155  namelist /atmos_dyn_tstep_large_fvm_heve/ &
156  eval_type_numfilter
157 
158  integer :: ierr
159  !---------------------------------------------------------------------------
160 
161  !--- read namelist
162  rewind(io_fid_conf)
163  read(io_fid_conf,nml=atmos_dyn_tstep_large_fvm_heve,iostat=ierr)
164  if( ierr < 0 ) then !--- missing
165  log_info("ATMOS_DYN_Tstep_large_fvm_heve_setup",*) 'Not found namelist. Default used.'
166  elseif( ierr > 0 ) then !--- fatal error
167  log_error("ATMOS_DYN_Tstep_large_fvm_heve_setup",*) 'Not appropriate names in namelist ATMOS_DYN_TSTEP_LARGE_FVM_HEVE. Check!'
168  call prc_abort
169  endif
171 
172  select case( eval_type_numfilter )
173  case( 'TENDENCY', 'FILTER' )
174  case default
175  log_error("ATMOS_DYN_Tstep_large_fvm_heve_setup",*) 'The specfied value of EVAL_TYPE_NUMFILTER is not appropriate. Check!'
176  call prc_abort
177  end select
178  !--
179 
180  allocate( dens_t(ka,ia,ja) )
181  allocate( momz_t(ka,ia,ja) )
182  allocate( momx_t(ka,ia,ja) )
183  allocate( momy_t(ka,ia,ja) )
184  allocate( rhot_t(ka,ia,ja) )
185  allocate( rhoq_t(ka,ia,ja,qa) )
186 
187  allocate( dens_damp(ka,ia,ja) )
188 
189  allocate( mflx(ka,ia,ja,3) )
190 
191  allocate( i_comm_prog(max(va,1)) )
192  allocate( i_comm_qtrc(qa) )
193  allocate( i_comm_rhoq_t(qa) )
194 
195  call comm_vars8_init( 'DENS', dens, i_comm_dens )
196  call comm_vars8_init( 'MOMZ', momz, i_comm_momz )
197  call comm_vars8_init( 'MOMX', momx, i_comm_momx )
198  call comm_vars8_init( 'MOMY', momy, i_comm_momy )
199  call comm_vars8_init( 'RHOT', rhot, i_comm_rhot )
200  do iv = 1, va
201  i_comm_prog(iv) = 5 + iv
202  call comm_vars8_init( 'PROG', prog(:,:,:,iv), i_comm_prog(iv) )
203  end do
204 
205  call comm_vars8_init( 'DENS_t', dens_t, i_comm_dens_t )
206  call comm_vars8_init( 'MOMZ_t', momz_t, i_comm_momz_t )
207  call comm_vars8_init( 'MOMX_t', momx_t, i_comm_momx_t )
208  call comm_vars8_init( 'MOMY_t', momy_t, i_comm_momy_t )
209  call comm_vars8_init( 'RHOT_t', rhot_t, i_comm_rhot_t )
210 
211  do iq = 1, qa
212  i_comm_rhoq_t(iq) = 5 + va + iq
213  i_comm_qtrc(iq) = 5 + va + iq
214 
215  call comm_vars8_init( 'RHOQ_t', rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq) )
216  call comm_vars8_init( 'QTRC', qtrc(:,:,:,iq), i_comm_qtrc(iq) )
217  end do
218 
219  call comm_vars8_init( 'DENS_t', dens_t, i_comm_dens_damp )
220 
221  call comm_vars8_init( 'mflx_Z', mflx(:,:,:,zdir), i_comm_mflx_z )
222  call comm_vars8_init( 'mflx_X', mflx(:,:,:,xdir), i_comm_mflx_x )
223  call comm_vars8_init( 'mflx_Y', mflx(:,:,:,ydir), i_comm_mflx_y )
224 
225  allocate( zero(ka,ia,ja) )
226  zero(:,:,:) = 0.0_rp
227 
228  mflx(:,:,:,:) = undef
229  if ( prc_twod ) then
230  mflx(:,:,:,xdir) = 0.0_rp
231  end if
232 
233 
234  ! history
235  call file_history_reg( 'ZFLX_MOM', 'momentum flux of z-direction', 'kg/m2/s', & ! [IN]
236  hist_mflx(1), & ! [OUT]
237  dim_type='ZHXY' ) ! [IN]
238  call file_history_reg( 'XFLX_MOM', 'momentum flux of x-direction', 'kg/m2/s', & ! [IN]
239  hist_mflx(2), & ! [OUT]
240  dim_type='ZXHY' ) ! [IN]
241  call file_history_reg( 'YFLX_MOM', 'momentum flux of y-direction', 'kg/m2/s', & ! [IN]
242  hist_mflx(3), & ! [OUT]
243  dim_type='ZXYH' ) ! [IN]
244 
245  call file_history_reg( 'ZFLX_RHOT', 'potential temperature flux of z-direction', 'K*kg/m2/s', & ! [IN]
246  hist_tflx(1), & ! [OUT]
247  dim_type='ZHXY' ) ! [IN]
248  call file_history_reg( 'XFLX_RHOT', 'potential temperature flux of x-direction', 'K*kg/m2/s', & ! [IN]
249  hist_tflx(2), & ! [OUT]
250  dim_type='ZXHY' ) ! [IN]
251  call file_history_reg( 'YFLX_RHOT', 'potential temperature flux of y-direction', 'K*kg/m2/s', & ! [IN]
252  hist_tflx(3), & ! [OUT]
253  dim_type='ZXYH' ) ! [IN]
254 
255  call file_history_reg( 'DENS_t_phys', 'tendency of dencity due to physics', 'kg/m3/s', & ! [IN]
256  hist_phys(1) ) ! [OUT]
257  call file_history_reg( 'MOMZ_t_phys', 'tendency of momentum z due to physics', 'kg/m2/s2', & ! [IN]
258  hist_phys(2), & ! [OUT]
259  dim_type='ZHXY' ) ! [IN]
260  call file_history_reg( 'MOMX_t_phys', 'tendency of momentum x due to physics', 'kg/m2/s2', & ! [IN]
261  hist_phys(3), & ! [OUT]
262  dim_type='ZXHY' ) ! [IN]
263  call file_history_reg( 'MOMY_t_phys', 'tendency of momentum y due to physics', 'kg/m2/s2', & ! [IN]
264  hist_phys(4), & ! [OUT]
265  dim_type='ZXYH' ) ! [IN]
266  call file_history_reg( 'RHOT_t_phys', 'tendency of rho*theta temperature due to physics', 'K*kg/m3/s', & ! [IN]
267  hist_phys(5) ) ! [OUT]
268 
269  call file_history_reg( 'DENS_t_damp', 'tendency of dencity due to damping', 'kg/m3/s', & ! [IN]
270  hist_damp(1) ) ! [OUT]
271  call file_history_reg( 'MOMZ_t_damp', 'tendency of momentum z due to damping', 'kg/m2/s2', & ! [IN]
272  hist_damp(2), & ! [OUT]
273  dim_type='ZHXY' ) ! [IN]
274  call file_history_reg( 'MOMX_t_damp', 'tendency of momentum x due to damping', 'kg/m2/s2', & ! [IN]
275  hist_damp(3), & ! [OUT]
276  dim_type='ZXHY' ) ! [IN]
277  call file_history_reg( 'MOMY_t_damp', 'tendency of momentum y due to damping', 'kg/m2/s2', & ! [IN]
278  hist_damp(4), & ! [OUT]
279  dim_type='ZXYH' ) ! [IN]
280  call file_history_reg( 'RHOT_t_damp', 'tendency of rho*theta temperature due to damping', 'K kg/m3/s', & ! [IN]
281  hist_damp(5) ) ! [OUT]
282 
283  allocate( hist_qflx(3,qa) )
284  allocate( hist_phys_qtrc(qa) )
285  allocate( hist_damp_qtrc(qa) )
286  do iq = 1, qa
287  call file_history_reg( 'ZFLX_'//trim(tracer_name(iq)), trim(tracer_name(iq))//' flux of z-direction', 'kg/m2/s', & ! [IN]
288  hist_qflx(1,iq), & ! [OUT]
289  dim_type='ZHXY' ) ! [IN]
290  call file_history_reg( 'XFLX_'//trim(tracer_name(iq)), trim(tracer_name(iq))//' flux of x-direction', 'kg/m2/s', & ! [IN]
291  hist_qflx(2,iq), & ! [OUT]
292  dim_type='ZXHY' ) ! [IN]
293  call file_history_reg( 'YFLX_'//trim(tracer_name(iq)), trim(tracer_name(iq))//' flux of y-direction', 'kg/m2/s', & ! [IN]
294  hist_qflx(3,iq), & ! [OUT]
295  dim_type='ZXYH' ) ! [IN]
296  call file_history_reg( trim(tracer_name(iq))//'_t_phys', 'tendency of '//trim(tracer_name(iq))//' due to physics', 'kg/m3/s', &
297  hist_phys_qtrc(iq) )
298  call file_history_reg( trim(tracer_name(iq))//'_t_damp', 'tendency of '//trim(tracer_name(iq))//' due to damping', 'kg/m3/s', &
299  hist_damp_qtrc(iq) )
300  end do
301 
302  ! for history at t=0
303  do iv = 1, 3
304  call file_history_put( hist_mflx(iv), zero(:,:,:) )
305  call file_history_put( hist_tflx(iv), zero(:,:,:) )
306  end do
307  do iv = 1, 5
308  call file_history_put( hist_phys(iv), zero(:,:,:) )
309  call file_history_put( hist_damp(iv), zero(:,:,:) )
310  end do
311  do iq = 1, qa
312  do iv = 1, 3
313  call file_history_put( hist_qflx(iv,iq), zero(:,:,:) )
314  end do
315  call file_history_put( hist_phys_qtrc(iq), zero(:,:,:) )
316  call file_history_put( hist_damp_qtrc(iq), zero(:,:,:) )
317  end do
318 
319 
320  ! for monitor
321  allocate( zero_x(ka,ja) )
322  allocate( zero_y(ka,ia) )
323  zero_x(:,:) = 0.0_rp
324  zero_y(:,:) = 0.0_rp
325 
326  call monitor_reg( "MASSTND_DAMP", "mass tendency by the damping", "kg", & ! [IN]
327  monit_damp_mass, & ! [OUT]
328  is_tendency=.true. ) ! [IN]
329 
330  call monitor_reg( "MASSFLX_WEST", "mass flux at the western boundary", "kg", & ! [IN]
331  monit_mflx_west, & ! [OUT]
332  dim_type="ZY-W", is_tendency=.true. ) ! [IN]
333  call monitor_reg( "MASSFLX_EAST", "mass flux at the eastern boundary", "kg", & ! [IN]
334  monit_mflx_east, & ! [OUT]
335  dim_type="ZY-E", is_tendency=.true. ) ! [IN]
336  call monitor_reg( "MASSFLX_SOUTH", "mass flux at the southern boundary", "kg", & ! [IN]
337  monit_mflx_south, & ! [OUT]
338  dim_type="ZX-S", is_tendency=.true. ) ! [IN]
339  call monitor_reg( "MASSFLX_NORTH", "mass flux at the northern boundary", "kg", & ! [IN]
340  monit_mflx_north, & ! [OUT]
341  dim_type="ZX-N", is_tendency=.true. ) ! [IN]
342 
343  call monitor_reg( "QTOTTND_DAMP", "water mass tendency by the damping", "kg", & ! [IN]
344  monit_damp_qtot, & ! [OUT]
345  is_tendency=.true. ) ! [IN]
346 
347  call monitor_reg( "QTOTFLX_WEST", "water mass flux at the western boundary", "kg", & ! [IN]
348  monit_qflx_west, & ! [OUT]
349  dim_type="ZY-W", is_tendency=.true. ) ! [IN]
350  call monitor_reg( "QTOTFLX_EAST", "water mass flux at the eastern boundary", "kg", & ! [IN]
351  monit_qflx_east, & ! [OUT]
352  dim_type="ZY-E", is_tendency=.true. ) ! [IN]
353  call monitor_reg( "QTOTFLX_SOUTH", "water mass flux at the southern boundary", "kg", & ! [IN]
354  monit_qflx_south, & ! [OUT]
355  dim_type="ZX-S", is_tendency=.true. ) ! [IN]
356  call monitor_reg( "QTOTFLX_NORTH", "water mass flux at the northern boundary", "kg", & ! [IN]
357  monit_qflx_north, & ! [OUT]
358  dim_type="ZX-N", is_tendency=.true. ) ! [IN]
359 
360  ! at t=0
361  call monitor_put( monit_damp_mass, zero(:,:,:) )
362  call monitor_put( monit_mflx_west, zero_x(:,:) )
363  call monitor_put( monit_mflx_east, zero_x(:,:) )
364  call monitor_put( monit_mflx_south, zero_y(:,:) )
365  call monitor_put( monit_mflx_north, zero_y(:,:) )
366 
367  call monitor_put( monit_damp_qtot, zero(:,:,:) )
368  call monitor_put( monit_qflx_west, zero_x(:,:) )
369  call monitor_put( monit_qflx_east, zero_x(:,:) )
370  call monitor_put( monit_qflx_south, zero_y(:,:) )
371  call monitor_put( monit_qflx_north, zero_y(:,:) )
372 
373  return

References atmos_dyn_tstep_large_fvm_heve(), scale_comm_cartesc::comm_vars8_init(), scale_const::const_ohm, scale_const::const_undef, scale_file_history::file_history_reg(), scale_atmos_grid_cartesc_index::ia, scale_io::io_fid_conf, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::ka, scale_monitor::monitor_reg(), 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_tracer::tracer_name, scale_index::va, scale_atmos_grid_cartesc_index::xdir, scale_atmos_grid_cartesc_index::ydir, and scale_atmos_grid_cartesc_index::zdir.

Referenced by scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup().

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

◆ atmos_dyn_tstep_large_fvm_heve()

subroutine, public scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve ( 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,5,3), intent(out)  num_diff,
real(rp), dimension(ka,ia,ja,3), intent(out)  num_diff_q,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC0,
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)  CORIOLI,
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,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
logical, intent(in)  TwoD,
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), dimension(ka), intent(in)  wdamp_coef,
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)  DTLS,
real(dp), intent(in)  DTSS,
logical, intent(in)  Llast 
)

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 407 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

407  use scale_const, only: &
408  eps => const_eps
409  use scale_comm_cartesc, only: &
410  comm_vars8, &
411  comm_wait
412  use scale_atmos_dyn_common, only: &
416  use scale_atmos_dyn_fvm_numfilter, only: &
420  use scale_file_history, only: &
421  file_history_query, &
422  file_history_put, &
424  use scale_monitor, only: &
425  monitor_put
426  use scale_atmos_dyn_tinteg_short, only: &
428  use scale_atmos_dyn_tinteg_tracer, only: &
430  use scale_spnudge, only: &
431  spnudge_uv, &
433  spnudge_pt, &
434  spnudge_qv, &
435  spnudge_u_alpha, &
436  spnudge_v_alpha, &
439  spnudge_uv_lm, &
440  spnudge_uv_mm, &
441  spnudge_pt_lm, &
442  spnudge_pt_mm, &
443  spnudge_qv_lm, &
445  use scale_dft, only: &
446  dft_g2g_divfree, &
447  dft_g2g
448  implicit none
449 
450  real(RP), intent(inout) :: DENS(KA,IA,JA)
451  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
452  real(RP), intent(inout) :: MOMX(KA,IA,JA)
453  real(RP), intent(inout) :: MOMY(KA,IA,JA)
454  real(RP), intent(inout) :: RHOT(KA,IA,JA)
455  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
456  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
457 
458  real(RP), intent(inout) :: DENS_av(KA,IA,JA)
459  real(RP), intent(inout) :: MOMZ_av(KA,IA,JA)
460  real(RP), intent(inout) :: MOMX_av(KA,IA,JA)
461  real(RP), intent(inout) :: MOMY_av(KA,IA,JA)
462  real(RP), intent(inout) :: RHOT_av(KA,IA,JA)
463  real(RP), intent(inout) :: QTRC_av(KA,IA,JA,QA)
464 
465  real(RP), intent(out) :: num_diff(KA,IA,JA,5,3)
466  real(RP), intent(out) :: num_diff_q(KA,IA,JA,3)
467 
468  real(RP), intent(in) :: QTRC0(KA,IA,JA,QA)
469 
470  real(RP), intent(in) :: DENS_tp(KA,IA,JA)
471  real(RP), intent(in) :: MOMZ_tp(KA,IA,JA)
472  real(RP), intent(in) :: MOMX_tp(KA,IA,JA)
473  real(RP), intent(in) :: MOMY_tp(KA,IA,JA)
474  real(RP), intent(in) :: RHOT_tp(KA,IA,JA)
475  real(RP), intent(in) :: RHOQ_tp(KA,IA,JA,QA)
476 
477  real(RP), intent(in) :: CORIOLI(IA,JA)
478 
479  real(RP), intent(in) :: CDZ (KA)
480  real(RP), intent(in) :: CDX (IA)
481  real(RP), intent(in) :: CDY (JA)
482  real(RP), intent(in) :: FDZ (KA-1)
483  real(RP), intent(in) :: FDX (IA-1)
484  real(RP), intent(in) :: FDY (JA-1)
485  real(RP), intent(in) :: RCDZ(KA)
486  real(RP), intent(in) :: RCDX(IA)
487  real(RP), intent(in) :: RCDY(JA)
488  real(RP), intent(in) :: RFDZ(KA-1)
489  real(RP), intent(in) :: RFDX(IA-1)
490  real(RP), intent(in) :: RFDY(JA-1)
491 
492  real(RP), intent(in) :: PHI (KA,IA,JA)
493  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
494  real(RP), intent(in) :: J13G (KA,IA,JA,7)
495  real(RP), intent(in) :: J23G (KA,IA,JA,7)
496  real(RP), intent(in) :: J33G
497  real(RP), intent(in) :: MAPF (IA,JA,2,4)
498 
499  real(RP), intent(in) :: AQ_R (QA)
500  real(RP), intent(in) :: AQ_CV (QA)
501  real(RP), intent(in) :: AQ_CP (QA)
502  real(RP), intent(in) :: AQ_MASS(QA)
503 
504  real(RP), intent(in) :: REF_dens(KA,IA,JA)
505  real(RP), intent(in) :: REF_pott(KA,IA,JA)
506  real(RP), intent(in) :: REF_qv (KA,IA,JA)
507  real(RP), intent(in) :: REF_pres(KA,IA,JA)
508 
509  logical, intent(in) :: BND_W
510  logical, intent(in) :: BND_E
511  logical, intent(in) :: BND_S
512  logical, intent(in) :: BND_N
513  logical, intent(in) :: TwoD
514 
515  real(RP), intent(in) :: ND_COEF
516  real(RP), intent(in) :: ND_COEF_Q
517  integer, intent(in) :: ND_LAPLACIAN_NUM
518  real(RP), intent(in) :: ND_SFC_FACT
519  logical, intent(in) :: ND_USE_RS
520 
521  integer, intent(in) :: BND_QA
522  integer, intent(in) :: BND_IQ(QA)
523  real(RP), intent(in) :: BND_SMOOTHER_FACT
524 
525  real(RP), intent(in) :: DAMP_DENS(KA,IA,JA)
526  real(RP), intent(in) :: DAMP_VELZ(KA,IA,JA)
527  real(RP), intent(in) :: DAMP_VELX(KA,IA,JA)
528  real(RP), intent(in) :: DAMP_VELY(KA,IA,JA)
529  real(RP), intent(in) :: DAMP_POTT(KA,IA,JA)
530  real(RP), intent(in) :: DAMP_QTRC(KA,IA,JA,BND_QA)
531 
532  real(RP), intent(in) :: DAMP_alpha_DENS(KA,IA,JA)
533  real(RP), intent(in) :: DAMP_alpha_VELZ(KA,IA,JA)
534  real(RP), intent(in) :: DAMP_alpha_VELX(KA,IA,JA)
535  real(RP), intent(in) :: DAMP_alpha_VELY(KA,IA,JA)
536  real(RP), intent(in) :: DAMP_alpha_POTT(KA,IA,JA)
537  real(RP), intent(in) :: DAMP_alpha_QTRC(KA,IA,JA,BND_QA)
538  real(RP), intent(in) :: MFLUX_OFFSET_X(KA,JA,2)
539  real(RP), intent(in) :: MFLUX_OFFSET_Y(KA,IA,2)
540 
541  real(RP), intent(in) :: wdamp_coef(KA)
542  real(RP), intent(in) :: divdmp_coef
543 
544  logical, intent(in) :: FLAG_TRACER_SPLIT_TEND
545  logical, intent(in) :: FLAG_FCT_MOMENTUM
546  logical, intent(in) :: FLAG_FCT_T
547  logical, intent(in) :: FLAG_FCT_TRACER
548  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
549 
550  logical, intent(in) :: USE_AVERAGE
551 
552  integer, intent(in) :: I_QV
553 
554  real(DP), intent(in) :: DTLS
555  real(DP), intent(in) :: DTSS
556  logical , intent(in) :: Llast
557 
558  ! for time integartion
559  real(RP) :: DENS00 (KA,IA,JA) ! saved density before small step loop
560  real(RP) :: qflx (KA,IA,JA,3)
561 
562  ! diagnostic variables
563  real(RP) :: DDIV (KA,IA,JA) ! 3 dimensional divergence
564  real(RP) :: DPRES0 (KA,IA,JA) ! pressure deviation
565  real(RP) :: RT2P (KA,IA,JA) ! factor of RHOT to PRES
566  real(RP) :: REF_rhot(KA,IA,JA) ! reference of RHOT
567 
568  real(RP) :: DENS_tq(KA,IA,JA)
569  real(RP) :: diff(KA,IA,JA), diff2(KA,IA,JA), diff3(KA,IA,JA)
570  real(RP) :: damp
571  real(RP) :: damp_t_DENS(KA,IA,JA)
572  real(RP) :: damp_t_MOMZ(KA,IA,JA)
573  real(RP) :: damp_t_MOMX(KA,IA,JA)
574  real(RP) :: damp_t_MOMY(KA,IA,JA)
575  real(RP) :: damp_t_RHOT(KA,IA,JA)
576  real(RP) :: damp_t_QTRC(KA,IA,JA)
577 
578  real(RP) :: tflx(KA,IA,JA,3)
579 
580  ! For tracer advection
581  real(RP) :: mflx_av(KA,IA,JA,3) ! rho * vel(x,y,z) @ (u,v,w)-face average
582 
583  real(RP) :: dtl
584  real(RP) :: dts
585  integer :: nstep
586 
587  ! for history
588  logical :: do_put
589 
590  ! for monitor
591  real(RP), target :: qflx_west (KA,JA)
592  real(RP), target :: qflx_east (KA,JA)
593  real(RP), target :: qflx_south(KA,IA)
594  real(RP), target :: qflx_north(KA,IA)
595  logical :: MONIT_lateral_flag(3)
596 
597  integer :: i, j, k, iq, iqb, step
598  integer :: iv
599  integer :: n
600  !---------------------------------------------------------------------------
601 
602  call prof_rapstart("DYN_Large_Preparation", 2)
603 
604  dts = real(dtss, kind=rp) ! short time step
605  dtl = real(dtls, kind=rp) ! large time step
606  nstep = ceiling( ( dtl - eps ) / dts )
607  dts = dtl / nstep ! dts is divisor of dtl and smaller or equal to dtss
608 
609  monit_lateral_flag(zdir) = .false.
610  monit_lateral_flag(xdir) = monit_mflx_west > 0 .or. monit_mflx_east > 0
611  monit_lateral_flag(ydir) = monit_mflx_south > 0 .or. monit_mflx_north > 0
612 
613 #ifdef DEBUG
614  log_info("ATMOS_DYN_Tstep_large_fvm_heve",*) 'Dynamics large time step'
615  log_info_cont('(1x,A,F0.2,A,F0.2,A,I0)') &
616  '-> DT_large, DT_small, DT_large/DT_small : ', dtl, ', ', dts, ', ', nstep
617 
618  dens00(:,:,:) = undef
619  num_diff(:,:,:,:,:) = undef
620 #endif
621 
622 !OCL XFILL
623  dens00(:,:,:) = dens(:,:,:)
624 
625  if ( use_average ) then
626 !OCL XFILL
627  !$omp parallel do
628  do j = jsb, jeb
629  do i = isb, ieb
630  do k = ks, ke
631  dens_av(k,i,j) = 0.0_rp
632  momz_av(k,i,j) = 0.0_rp
633  momx_av(k,i,j) = 0.0_rp
634  momy_av(k,i,j) = 0.0_rp
635  rhot_av(k,i,j) = 0.0_rp
636  end do
637  end do
638  end do
639  endif
640 
641 !OCL XFILL
642  mflx_av(:,:,:,:) = 0.0_rp
643 
644 
645 !OCL XFILL
646  !$omp parallel do
647  do j = js, je
648  do k = ks, ke
649  qflx_west(k,j) = 0.0_rp
650  qflx_east(k,j) = 0.0_rp
651  end do
652  end do
653 !OCL XFILL
654  !$omp parallel do
655  do i = is, ie
656  do k = ks, ke
657  qflx_south(k,i) = 0.0_rp
658  qflx_north(k,i) = 0.0_rp
659  end do
660  end do
661 
662  !- prepare some variables for pressure linearization
664  dpres0, rt2p, ref_rhot, & ! (out)
665  rhot, qtrc, ref_pres, aq_r, aq_cv, aq_cp, aq_mass ) ! (in)
666 
667  call prof_rapend ("DYN_Large_Preparation", 2)
668 
669  !###########################################################################
670  ! Update DENS,MONZ,MOMX,MOMY,MOMZ,RHOT
671  !###########################################################################
672 
673  call prof_rapstart("DYN_Large_Tendency", 2)
674 
675 !OCL XFILL
676  dens_tq(:,:,:) = 0.0_rp
677 
678  do iq = 1, qa
679 
680  iqb = bnd_iq(iq)
681 
682  if ( iqb > 0 ) then
683 
684  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
685 !OCL XFILL
686  do j = js-1, je+1
687  do i = max(is-1,1), min(ie+1,ia)
688  do k = ks, ke
689  diff(k,i,j) = qtrc(k,i,j,iq) - damp_qtrc(k,i,j,iqb)
690  enddo
691  enddo
692  enddo
693 
694  call file_history_query( hist_damp_qtrc(iq), do_put )
695  if ( twod ) then
696  !$omp parallel do default(none) OMP_SCHEDULE_ &
697  !$omp private(j,k,damp) &
698  !$omp shared(IS,JS,JE,KS,KE,iq,iqb) &
699  !$omp shared(RHOQ_t,RHOQ_tp,DENS_tq,DAMP_alpha_QTRC,diff,BND_SMOOTHER_FACT,DENS00,TRACER_MASS,I_QV) &
700  !$omp shared(damp_t_QTRC,do_put)
701 !OCL XFILL
702  do j = js, je
703  do k = ks, ke
704  damp = - damp_alpha_qtrc(k,is,j,iqb) &
705  * ( diff(k,is,j) & ! rayleigh damping
706  - ( diff(k,is,j-1) + diff(k,is,j+1) - diff(k,is,j)*2.0_rp ) &
707  * 0.25_rp * bnd_smoother_fact ) ! horizontal smoother
708  damp = damp * dens00(k,is,j)
709  if ( do_put ) damp_t_qtrc(k,is,j) = damp
710  rhoq_t(k,is,j,iq) = rhoq_tp(k,is,j,iq) + damp
711  dens_tq(k,is,j) = dens_tq(k,is,j) + damp * tracer_mass(iq) ! only for mass tracer
712  enddo
713  enddo
714  else
715  if( iq == i_qv .and. spnudge_qv ) then
716 
717  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
718 !OCL XFILL
719  do j = js, je
720  do i = is, ie
721  do k = ks, ke
722  diff2(k,i,j) = diff(k,i,j)
723  enddo
724  enddo
725  enddo
726 
728 
729  else
730 
731  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
732 !OCL XFILL
733  do j = js, je
734  do i = is, ie
735  do k = ks, ke
736  diff2(k,i,j) = 0
737  enddo
738  enddo
739  enddo
740 
741  endif
742 
743  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
744  !$omp private(i,j,k,damp) &
745  !$omp shared(JS,JE,IS,IE,KS,KE,iq,iqb) &
746  !$omp shared(RHOQ_t,RHOQ_tp,DENS_tq,DAMP_alpha_QTRC,diff,diff2,BND_SMOOTHER_FACT,SPNUDGE_qv_alpha,DENS00,TRACER_MASS) &
747  !$omp shared(damp_t_QTRC,do_put)
748 !OCL XFILL
749  do j = js, je
750  do i = is, ie
751  do k = ks, ke
752  damp = - damp_alpha_qtrc(k,i,j,iqb) &
753  * ( diff(k,i,j) & ! rayleigh damping
754  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
755  * 0.125_rp * bnd_smoother_fact ) & ! horizontal smoother
756  - spnudge_qv_alpha(k,i,j) * diff2(k,i,j)
757  damp = damp * dens00(k,i,j)
758  if ( do_put ) damp_t_qtrc(k,i,j) = damp
759  rhoq_t(k,i,j,iq) = rhoq_tp(k,i,j,iq) + damp
760  dens_tq(k,i,j) = dens_tq(k,i,j) + damp * tracer_mass(iq) ! only for mass tracer
761  enddo
762  enddo
763  enddo
764  end if
765 
766  if ( llast ) then
767  if ( do_put ) call file_history_put( hist_damp_qtrc(iq), damp_t_qtrc(:,:,:) )
768  call file_history_put( hist_phys_qtrc(iq), rhoq_tp(:,:,:,iq) )
769  end if
770 
771  call atmos_dyn_fill_halo( rhoq_t(:,:,:,iq), 0.0_rp, .false., .true. )
772  call comm_vars8( rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq) )
773  call comm_wait ( rhoq_t(:,:,:,iq), i_comm_rhoq_t(iq), .false. )
774 
775  else
776 
777  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
778 !OCL XFILL
779  do j = 1, ja
780  do i = 1, ia
781  do k = 1, ka
782  rhoq_t(k,i,j,iq) = rhoq_tp(k,i,j,iq)
783  enddo
784  enddo
785  enddo
786 
787  end if
788 
789  end do
790 
791  call prof_rapend ("DYN_Large_Tendency", 2)
792 
793  call prof_rapstart("DYN_Large_Boundary", 2)
794 
795  if ( bnd_w .and. (.not. twod) ) then
796  !$omp parallel do private(j,k) OMP_SCHEDULE_
797  do j = js, je
798  do k = ks, ke
799  mflx(k,is-1,j,xdir) = ( momx(k,is-1,j) + mflux_offset_x(k,j,1) ) &
800  * gsqrt(k,is-1,j,i_uyz) / mapf(is-1,j,2,i_uy)
801  enddo
802  enddo
803  end if
804  if ( bnd_e .and. (.not. twod) ) then
805  !$omp parallel do private(j,k) OMP_SCHEDULE_
806  do j = js, je
807  do k = ks, ke
808  mflx(k,ie,j,xdir) = ( momx(k,ie,j) + mflux_offset_x(k,j,2) ) &
809  * gsqrt(k,ie,j,i_uyz) / mapf(ie,j,2,i_uy)
810  enddo
811  enddo
812  end if
813  if ( bnd_s ) then
814  !$omp parallel do private(i,k) OMP_SCHEDULE_
815  do i = is, ie
816  do k = ks, ke
817  mflx(k,i,js-1,ydir) = ( momy(k,i,js-1) + mflux_offset_y(k,i,1) ) &
818  * gsqrt(k,i,js-1,i_xvz) / mapf(i,js-1,1,i_xv)
819  enddo
820  enddo
821  end if
822  if ( bnd_n ) then
823  !$omp parallel do private(i,k) OMP_SCHEDULE_
824  do i = is, ie
825  do k = ks, ke
826  mflx(k,i,je,ydir) = ( momy(k,i,je) + mflux_offset_y(k,i,2) ) &
827  * gsqrt(k,i,je,i_xvz) / mapf(i,je,1,i_xv)
828  enddo
829  enddo
830  end if
831 
832  call prof_rapend ("DYN_Large_Boundary", 2)
833 
834 !OCL XFILL
835  !$omp parallel do collapse(2)
836  do j = 1, ja
837  do i = 1, ia
838  do k = 1, ka
839  damp_t_dens(k,i,j) = 0.0_rp
840  damp_t_momz(k,i,j) = 0.0_rp
841  damp_t_momx(k,i,j) = 0.0_rp
842  damp_t_momy(k,i,j) = 0.0_rp
843  damp_t_rhot(k,i,j) = 0.0_rp
844  end do
845  end do
846  end do
847 
848  call atmos_dyn_fill_halo( dens_damp, 0.0_rp, .true., .true. )
849  call atmos_dyn_fill_halo( dens_t, 0.0_rp, .false., .true. )
850 
851  do step = 1, nstep
852 
853  !-----< prepare tendency >-----
854 
855  call prof_rapstart("DYN_Large_Tendency", 2)
856 
857  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
858 !OCL XFILL
859  do j = js-1, je+1
860  do i = max(is-1,1), min(ie+1,ia)
861  do k = ks, ke
862  diff(k,i,j) = dens(k,i,j) - damp_dens(k,i,j)
863  enddo
864  enddo
865  enddo
866 
867  call file_history_query( hist_damp(1), do_put )
868  if ( twod ) then
869  !$omp parallel do default(none) OMP_SCHEDULE_ &
870  !$omp private(j,k) &
871  !$omp shared(IS,JS,JE,KS,KE) &
872  !$omp shared(DAMP_alpha_DENS,diff,DENS_tq,DENS_t,DENS_tp,BND_SMOOTHER_FACT,EPS) &
873  !$omp shared(damp_t_DENS,DENS_damp,do_put,nstep)
874 !OCL XFILL
875  do j = js, je
876  do k = ks, ke
877  dens_damp(k,is,j) = - damp_alpha_dens(k,is,j) &
878  * ( diff(k,is,j) & ! rayleigh damping
879  - ( diff(k,is,j-1) + diff(k,is,j+1) - diff(k,is,j)*2.0_rp ) &
880  * 0.25_rp * bnd_smoother_fact ) & ! horizontal smoother
881  + dens_tq(k,is,j) * ( 0.5_rp - sign( 0.5_rp, damp_alpha_dens(k,is,j)-eps ) ) ! dencity change due to rayleigh damping for tracers
882  dens_t(k,is,j) = dens_tp(k,is,j) & ! tendency from physical step
883  + dens_damp(k,is,j)
884  if ( do_put ) damp_t_dens(k,is,j) = damp_t_dens(k,is,j) + dens_damp(k,is,j) / nstep
885  enddo
886  enddo
887  else
888  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
889  !$omp private(i,j,k) &
890  !$omp shared(JS,JE,IS,IE,KS,KE) &
891  !$omp shared(DAMP_alpha_DENS,diff,DENS_tq,DENS_t,DENS_tp,BND_SMOOTHER_FACT,EPS) &
892  !$omp shared(damp_t_DENS,DENS_damp,do_put,nstep)
893 !OCL XFILL
894  do j = js, je
895  do i = is, ie
896  do k = ks, ke
897  dens_damp(k,i,j) = - damp_alpha_dens(k,i,j) &
898  * ( diff(k,i,j) & ! rayleigh damping
899  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
900  * 0.125_rp * bnd_smoother_fact ) & ! horizontal smoother
901  + dens_tq(k,i,j) * ( 0.5_rp - sign( 0.5_rp, damp_alpha_dens(k,i,j)-eps ) ) ! dencity change due to rayleigh damping for tracers
902  dens_t(k,i,j) = dens_tp(k,i,j) & ! tendency from physical step
903  + dens_damp(k,i,j)
904  if ( do_put ) damp_t_dens(k,i,j) = damp_t_dens(k,i,j) + dens_damp(k,i,j) / nstep
905  enddo
906  enddo
907  enddo
908  end if
909 
910  call comm_vars8( dens_damp(:,:,:), i_comm_dens_damp )
911  call comm_vars8( dens_t(:,:,:), i_comm_dens_t )
912 
913  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
914 !OCL XFILL
915  do j = js-1, je+1
916  do i = max(is-1,1), min(ie+1,ia)
917  do k = ks, ke-1
918  diff(k,i,j) = momz(k,i,j) - damp_velz(k,i,j) * ( dens(k,i,j)+dens(k+1,i,j) ) * 0.5_rp
919  enddo
920  enddo
921  enddo
922 
923  call file_history_query( hist_damp(2), do_put )
924  if ( twod ) then
925  !$omp parallel do default(none) OMP_SCHEDULE_ &
926  !$omp private(j,k,damp) &
927  !$omp shared(IS,JS,JE,KS,KE) &
928  !$omp shared(DENS_damp,DENS,MOMZ) &
929  !$omp shared(DAMP_alpha_VELZ,diff,BND_SMOOTHER_FACT,MOMZ_t,MOMZ_tp) &
930  !$omp shared(damp_t_MOMZ,do_put,nstep)
931 !OCL XFILL
932  do j = js, je
933  do k = ks, ke-1
934  damp = - damp_alpha_velz(k,is,j) &
935  * ( diff(k,is,j) & ! rayleigh damping
936  - ( diff(k,is,j-1) + diff(k,is,j+1) - diff(k,is,j)*2.0_rp ) &
937  * 0.25_rp * bnd_smoother_fact ) ! horizontal smoother
938 
939  momz_t(k,is,j) = momz_tp(k,is,j) & ! tendency from physical step
940  + damp &
941  + ( dens_damp(k,is,j) + dens_damp(k+1,is,j) ) * momz(k,is,j) / ( dens(k,is,j) + dens(k,is,j) )
942  if ( do_put ) damp_t_momz(k,is,j) = damp_t_momz(k,is,j) + damp / nstep
943  enddo
944  enddo
945  else
946  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
947  !$omp private(i,j,k,damp) &
948  !$omp shared(JS,JE,IS,IE,KS,KE) &
949  !$omp shared(DENS_damp,DENS,MOMZ) &
950  !$omp shared(DAMP_alpha_VELZ,diff,BND_SMOOTHER_FACT,MOMZ_t,MOMZ_tp) &
951  !$omp shared(damp_t_MOMZ,do_put,nstep)
952 !OCL XFILL
953  do j = js, je
954  do i = is, ie
955  do k = ks, ke-1
956  damp = - damp_alpha_velz(k,i,j) &
957  * ( diff(k,i,j) & ! rayleigh damping
958  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
959  * 0.125_rp * bnd_smoother_fact ) ! horizontal smoother
960 
961  momz_t(k,i,j) = momz_tp(k,i,j) & ! tendency from physical step
962  + damp &
963  + ( dens_damp(k,i,j) + dens_damp(k+1,i,j) ) * momz(k,i,j) / ( dens(k,i,j) + dens(k,i,j) )
964  if ( do_put ) damp_t_momz(k,i,j) = damp_t_momz(k,i,j) + damp / nstep
965  enddo
966  enddo
967  enddo
968  end if
969 !OCL XFILL
970  do j = js, je
971  do i = is, ie
972  momz_t( 1:ks-1,i,j) = 0.0_rp
973  momz_t(ke:ka ,i,j) = 0.0_rp
974  enddo
975  enddo
976  call comm_vars8( momz_t(:,:,:), i_comm_momz_t )
977 
978  call comm_wait( dens_damp(:,:,:), i_comm_dens_damp )
979 
980  if ( twod ) then
981  !$omp parallel do private(j,k) OMP_SCHEDULE_
982 !OCL XFILL
983  do j = js-1, je+1
984  do k = ks, ke
985  diff(k,is,j) = momx(k,is,j) - damp_velx(k,is,j) * dens(k,is,j)
986  enddo
987  enddo
988  else
989  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
990 !OCL XFILL
991  do j = js-1, je+1
992  do i = is-1, ie+1
993  do k = ks, ke
994  diff(k,i,j) = momx(k,i,j) - damp_velx(k,i,j) * ( dens(k,i,j)+dens(k,i+1,j) ) * 0.5_rp
995  enddo
996  enddo
997  enddo
998  end if
999 
1000 
1001  call file_history_query( hist_damp(3), do_put )
1002  if ( twod ) then
1003 !OCL XFILL
1004  !$omp parallel do default(none) OMP_SCHEDULE_ &
1005  !$omp private(j,k,damp) &
1006  !$omp shared(IS,JS,JE,KS,KE) &
1007  !$omp shared(DENS_damp,DENS,MOMX) &
1008  !$omp shared(DAMP_alpha_VELX,diff,BND_SMOOTHER_FACT,MOMX_tp,MOMX_t) &
1009  !$omp shared(damp_t_MOMX,do_put,nstep)
1010  do j = js, je
1011  do k = ks, ke
1012  damp = - damp_alpha_velx(k,is,j) &
1013  * ( diff(k,is,j) & ! rayleigh damping
1014  - ( diff(k,is,j-1) + diff(k,is,j+1) - diff(k,is,j)*2.0_rp ) &
1015  * 0.25_rp * bnd_smoother_fact ) ! horizontal smoother
1016  momx_t(k,is,j) = momx_tp(k,is,j) & ! tendency from physical step
1017  + damp &
1018  + dens_damp(k,is,j) * momx(k,is,j) / dens(k,is,j)
1019  if ( do_put ) damp_t_momx(k,is,j) = damp_t_momx(k,is,j) + damp / nstep
1020  enddo
1021  enddo
1022  else
1023  if( spnudge_uv ) then
1024  if( spnudge_uv_divfree ) then
1025 
1026  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1027 !OCL XFILL
1028  do j = js-1, je+1
1029  do i = is-1, ie+1
1030  do k = ks, ke
1031  diff3(k,i,j) = momy(k,i,j) - damp_vely(k,i,j) * ( dens(k,i,j)+dens(k,i,j+1) ) * 0.5_rp
1032  enddo
1033  enddo
1034  enddo
1035 
1036  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1037 !OCL XFILL
1038  do j = js, je
1039  do i = is, ie
1040  do k = ks, ke
1041  diff2(k,i,j) = diff(k,i,j) / ( ( dens(k,i,j)+dens(k,i+1,j) ) * 0.5_rp )
1042  enddo
1043  enddo
1044  enddo
1045 
1046  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1047 !OCL XFILL
1048  do j = js, je
1049  do i = is, ie
1050  do k = ks, ke
1051  diff3(k,i,j) = diff3(k,i,j) / ( ( dens(k,i,j)+dens(k,i,j+1) ) * 0.5_rp )
1052  enddo
1053  enddo
1054  enddo
1055 
1057 
1058  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1059 !OCL XFILL
1060  do j = js, je
1061  do i = is, ie
1062  do k = ks, ke
1063  diff2(k,i,j) = diff2(k,i,j) * ( dens(k,i,j)+dens(k,i+1,j) ) * 0.5_rp
1064  enddo
1065  enddo
1066  enddo
1067 
1068  else
1069 
1070  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1071 !OCL XFILL
1072  do j = js, je
1073  do i = is, ie
1074  do k = ks, ke
1075  diff2(k,i,j) = diff(k,i,j) / ( ( dens(k,i,j)+dens(k,i+1,j) ) * 0.5_rp )
1076  enddo
1077  enddo
1078  enddo
1079 
1081 
1082  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1083  do j = js, je
1084  do i = is, ie
1085  do k = ks, ke
1086  diff2(k,i,j) = diff2(k,i,j) * ( dens(k,i,j)+dens(k,i+1,j) ) * 0.5_rp
1087  enddo
1088  enddo
1089  enddo
1090 
1091  endif
1092  else
1093 
1094  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1095 !OCL XFILL
1096  do j = js, je
1097  do i = is, ie
1098  do k = ks, ke
1099  diff2(k,i,j) = 0
1100  enddo
1101  enddo
1102  enddo
1103 
1104  endif
1105 
1106 !OCL XFILL
1107  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1108  !$omp private(i,j,k,damp) &
1109  !$omp shared(JS,JE,IS,IE,KS,KE) &
1110  !$omp shared(DENS_damp,DENS,MOMX) &
1111  !$omp shared(DAMP_alpha_VELX,diff,diff2,BND_SMOOTHER_FACT,SPNUDGE_u_alpha,MOMX_tp,MOMX_t) &
1112  !$omp shared(damp_t_MOMX,do_put,nstep)
1113  do j = js, je
1114  do i = is, ie
1115  do k = ks, ke
1116  damp = - damp_alpha_velx(k,i,j) &
1117  * ( diff(k,i,j) & ! rayleigh damping
1118  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
1119  * 0.125_rp * bnd_smoother_fact ) & ! horizontal smoother
1120  - spnudge_u_alpha(k,i,j) * diff2(k,i,j)
1121  momx_t(k,i,j) = momx_tp(k,i,j) & ! tendency from physical step
1122  + damp &
1123  + ( dens_damp(k,i,j) + dens_damp(k,i+1,j) ) * momx(k,i,j) / ( dens(k,i,j) + dens(k,i+1,j) )
1124  if ( do_put ) damp_t_momx(k,i,j) = damp_t_momx(k,i,j) + damp / nstep
1125  enddo
1126  enddo
1127  enddo
1128  end if
1129 
1130  call atmos_dyn_fill_halo( momx_t, 0.0_rp, .false., .true. )
1131  call comm_vars8( momx_t(:,:,:), i_comm_momx_t )
1132 
1133  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1134 !OCL XFILL
1135  do j = js-1, je+1
1136  do i = max(is-1,1), min(ie+1,ia)
1137  do k = ks, ke
1138  diff(k,i,j) = momy(k,i,j) - damp_vely(k,i,j) * ( dens(k,i,j)+dens(k,i,j+1) ) * 0.5_rp
1139  enddo
1140  enddo
1141  enddo
1142 
1143  call file_history_query( hist_damp(4), do_put )
1144  if ( twod ) then
1145 !OCL XFILL
1146  !$omp parallel do default(none) OMP_SCHEDULE_ &
1147  !$omp private(j,k,damp) &
1148  !$omp shared(IS,JS,JE,KS,KE) &
1149  !$omp shared(DENS_damp,DENS,MOMY) &
1150  !$omp shared(DAMP_alpha_VELY,diff,BND_SMOOTHER_FACT,MOMY_tp,MOMY_t) &
1151  !$omp shared(damp_t_MOMY,do_put,nstep)
1152  do j = js, je
1153  do k = ks, ke
1154  damp = - damp_alpha_vely(k,is,j) &
1155  * ( diff(k,is,j) & ! rayleigh damping
1156  - ( diff(k,is,j-1) + diff(k,is,j+1) - diff(k,is,j)*2.0_rp ) &
1157  * 0.25_rp * bnd_smoother_fact ) ! horizontal smoother
1158  momy_t(k,is,j) = momy_tp(k,is,j) & ! tendency from physical step
1159  + damp &
1160  + ( dens_damp(k,is,j) + dens_damp(k,is,j+1) ) * momy(k,is,j) / ( dens(k,is,j) + dens(k,is,j+1) )
1161  if ( do_put ) damp_t_momy(k,is,j) = damp_t_momy(k,is,j) + damp / nstep
1162  enddo
1163  enddo
1164  else
1165  if( spnudge_uv ) then
1166  if( spnudge_uv_divfree ) then
1167 
1168  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1169 !OCL XFILL
1170  do j = js, je
1171  do i = is, ie
1172  do k = ks, ke
1173  diff2(k,i,j) = diff3(k,i,j) * ( dens(k,i,j)+dens(k,i,j+1) ) * 0.5_rp
1174  enddo
1175  enddo
1176  enddo
1177 
1178  else
1179  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1180 !OCL XFILL
1181  do j = js, je
1182  do i = is, ie
1183  do k = ks, ke
1184  diff2(k,i,j) = diff(k,i,j) / ( ( dens(k,i,j)+dens(k,i,j+1) ) * 0.5_rp )
1185  enddo
1186  enddo
1187  enddo
1188 
1190 
1191  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1192 !OCL XFILL
1193  do j = js, je
1194  do i = is, ie
1195  do k = ks, ke
1196  diff2(k,i,j) = diff2(k,i,j) * ( dens(k,i,j)+dens(k,i,j+1) ) * 0.5_rp
1197  enddo
1198  enddo
1199  enddo
1200  endif
1201  endif
1202 
1203 !OCL XFILL
1204  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1205  !$omp private(i,j,k,damp) &
1206  !$omp shared(JS,JE,IS,IE,KS,KE) &
1207  !$omp shared(DENS_damp,DENS,MOMY) &
1208  !$omp shared(DAMP_alpha_VELY,diff,diff2,BND_SMOOTHER_FACT,SPNUDGE_v_alpha,MOMY_tp,MOMY_t) &
1209  !$omp shared(damp_t_MOMY,do_put,nstep)
1210  do j = js, je
1211  do i = is, ie
1212  do k = ks, ke
1213  damp = - damp_alpha_vely(k,i,j) &
1214  * ( diff(k,i,j) & ! rayleigh damping
1215  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
1216  * 0.125_rp * bnd_smoother_fact ) & ! horizontal smoother
1217  - spnudge_v_alpha(k,i,j) * diff2(k,i,j)
1218  momy_t(k,i,j) = momy_tp(k,i,j) & ! tendency from physical step
1219  + damp &
1220  + ( dens_damp(k,i,j) + dens_damp(k,i,j+1) ) * momy(k,i,j) / ( dens(k,i,j) + dens(k,i,j+1) )
1221  if ( do_put ) damp_t_momy(k,i,j) = damp_t_momy(k,i,j) + damp / nstep
1222  enddo
1223  enddo
1224  enddo
1225  end if
1226 !OCL XFILL
1227 
1228  call atmos_dyn_fill_halo( momy_t, 0.0_rp, .false., .true. )
1229  call comm_vars8( momy_t(:,:,:), i_comm_momy_t )
1230 
1231  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1232 !OCL XFILL
1233  do j = js-1, je+2
1234  do i = max(is-1,1), min(ie+2,ia)
1235  do k = ks, ke
1236  diff(k,i,j) = rhot(k,i,j) - damp_pott(k,i,j) * dens(k,i,j)
1237  enddo
1238  enddo
1239  enddo
1240 
1241  call file_history_query( hist_damp(5), do_put )
1242  if ( twod ) then
1243 !OCL XFILL
1244  !$omp parallel do default(none) OMP_SCHEDULE_ &
1245  !$omp private(j,k,damp) &
1246  !$omp shared(IS,JS,JE,KS,KE) &
1247  !$omp shared(DENS_damp,DENS,RHOT) &
1248  !$omp shared(DAMP_alpha_POTT,diff,BND_SMOOTHER_FACT,RHOT_t,RHOT_tp) &
1249  !$omp shared(damp_t_RHOT,do_put,nstep)
1250  do j = js, je
1251  do k = ks, ke
1252  damp = - damp_alpha_pott(k,is,j) &
1253  * ( diff(k,is,j) & ! rayleigh damping
1254  - ( diff(k,is,j-1) + diff(k,is,j+1) - diff(k,is,j)*2.0_rp ) &
1255  * 0.25_rp * bnd_smoother_fact ) ! horizontal smoother
1256  rhot_t(k,is,j) = rhot_tp(k,is,j) & ! tendency from physical step
1257  + damp &
1258  + dens_damp(k,is,j) * rhot(k,is,j) / dens(k,is,j)
1259  if ( do_put ) damp_t_rhot(k,is,j) = damp_t_rhot(k,is,j) + damp / nstep
1260  enddo
1261  enddo
1262  else
1263  if( spnudge_pt ) then
1264 
1265  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1266 !OCL XFILL
1267  do j = js, je
1268  do i = is, ie
1269  do k = ks, ke
1270  diff2(k,i,j) = diff(k,i,j) / dens(k,i,j)
1271  enddo
1272  enddo
1273  enddo
1274 
1276 
1277  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1278 !OCL XFILL
1279  do j = js, je
1280  do i = is, ie
1281  do k = ks, ke
1282  diff2(k,i,j) = diff2(k,i,j) * dens(k,i,j)
1283  enddo
1284  enddo
1285  enddo
1286 
1287  else
1288 
1289  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1290 !OCL XFILL
1291  do j = js, je
1292  do i = is, ie
1293  do k = ks, ke
1294  diff2(k,i,j) = 0
1295  enddo
1296  enddo
1297  enddo
1298 
1299  endif
1300 
1301 !OCL XFILL
1302  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1303  !$omp private(i,j,k,damp) &
1304  !$omp shared(JS,JE,IS,IE,KS,KE) &
1305  !$omp shared(DENS_damp,DENS,RHOT) &
1306  !$omp shared(DAMP_alpha_POTT,diff,diff2,BND_SMOOTHER_FACT,SPNUDGE_pt_alpha,RHOT_t,RHOT_tp) &
1307  !$omp shared(damp_t_RHOT,do_put,nstep)
1308  do j = js, je
1309  do i = is, ie
1310  do k = ks, ke
1311  damp = - damp_alpha_pott(k,i,j) &
1312  * ( diff(k,i,j) & ! rayleigh damping
1313  - ( diff(k,i-1,j) + diff(k,i+1,j) + diff(k,i,j-1) + diff(k,i,j+1) - diff(k,i,j)*4.0_rp ) &
1314  * 0.125_rp * bnd_smoother_fact ) & ! horizontal smoother
1315  - spnudge_pt_alpha(k,i,j) * diff2(k,i,j)
1316  rhot_t(k,i,j) = rhot_tp(k,i,j) & ! tendency from physical step
1317  + damp &
1318  + dens_damp(k,i,j) * rhot(k,i,j) / dens(k,i,j)
1319  if ( do_put ) damp_t_rhot(k,i,j) = damp_t_rhot(k,i,j) + damp / nstep
1320  enddo
1321  enddo
1322  enddo
1323  end if
1324 !OCL XFILL
1325 
1326  call atmos_dyn_fill_halo( rhot_t, 0.0_rp, .false., .true. )
1327  call comm_vars8( rhot_t(:,:,:), i_comm_rhot_t )
1328 
1329  call comm_wait ( dens_t(:,:,:), i_comm_dens_t, .false. )
1330  call comm_wait ( momz_t(:,:,:), i_comm_momz_t, .false. )
1331  call comm_wait ( momx_t(:,:,:), i_comm_momx_t, .false. )
1332  call comm_wait ( momy_t(:,:,:), i_comm_momy_t, .false. )
1333  call comm_wait ( rhot_t(:,:,:), i_comm_rhot_t, .false. )
1334 
1335  call prof_rapend ("DYN_Large_Tendency", 2)
1336 
1337  call prof_rapstart("DYN_Large_Numfilter", 2)
1338 
1339  !-----< prepare the fluxes of explicit numerical diffusion >-----
1340  if ( nd_coef == 0.0_rp .or. eval_type_numfilter == 'FILTER' ) then
1341 !OCL XFILL
1342  num_diff(:,:,:,:,:) = 0.0_rp
1343  else
1344  call atmos_dyn_fvm_numfilter_flux( num_diff(:,:,:,:,:), & ! [OUT]
1345  dens, momz, momx, momy, rhot, & ! [IN]
1346  cdz, cdx, cdy, fdz, fdx, fdy, twod, dts, & ! [IN]
1347  ref_dens, ref_pott, & ! [IN]
1348  nd_coef, nd_laplacian_num, nd_sfc_fact, nd_use_rs ) ! [IN]
1349  endif
1350 
1351  call prof_rapend ("DYN_Large_Numfilter", 2)
1352 
1353  !-----< calculate the divegence term >-----
1354 
1355  if ( divdmp_coef > 0.0_rp ) then
1356 
1357  call atmos_dyn_divergence( ddiv, & ! (out)
1358  momz, momx, momy, & ! (in)
1359  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
1360  twod, & ! (in)
1361  rcdz, rcdx, rcdy, rfdz, fdz ) ! (in)
1362 
1363  else
1364 
1365 !XFILL
1366  ddiv = 0.0_rp
1367 
1368  end if
1369 
1370  !------------------------------------------------------------------------
1371  ! Start short time integration
1372  !------------------------------------------------------------------------
1373 
1374  call file_history_set_disable( .not. ( llast .and. ( step==nstep ) ) )
1375 
1376  call prof_rapstart("DYN_Short_Tinteg", 2)
1377 
1378  call atmos_dyn_tinteg_short( dens, momz, momx, momy, rhot, prog, & ! (inout)
1379  mflx, tflx, & ! (inout, out)
1380  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! (in)
1381  dpres0, rt2p, corioli, & ! (in)
1382  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! (in)
1383  flag_fct_momentum, flag_fct_t, & ! (in)
1384  flag_fct_along_stream, & ! (in)
1385  cdz, fdz, fdx, fdy, & ! (in)
1386  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
1387  phi, gsqrt, j13g, j23g, j33g, mapf, & ! (in)
1388  ref_dens, ref_rhot, & ! (in)
1389  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! (in)
1390  dts ) ! (in)
1391 
1392  call prof_rapend ("DYN_Short_Tinteg", 2)
1393 
1394 
1395  call file_history_set_disable( .false. )
1396 
1397  if ( nd_coef > 0.0_rp .and. eval_type_numfilter == 'FILTER' ) then
1398  call comm_vars8( dens(:,:,:), i_comm_dens )
1399  call comm_vars8( momz(:,:,:), i_comm_momz )
1400  call comm_vars8( momx(:,:,:), i_comm_momx )
1401  call comm_vars8( momy(:,:,:), i_comm_momy )
1402  call comm_vars8( rhot(:,:,:), i_comm_rhot )
1403  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
1404  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
1405  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
1406  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
1407  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
1408 
1410  num_diff, & ! (out)
1411  dens, momz, momx, momy, rhot, & ! (inout)
1412  cdz, cdx, cdy, fdz, fdx, fdy, & ! (in)
1413  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! (in)
1414  twod, dts, & ! (in)
1415  gsqrt, mapf, ref_dens, ref_pott, & ! (in)
1416  nd_coef, nd_laplacian_num, nd_sfc_fact, nd_use_rs ) ! (in)
1417  endif
1418 
1419  !$omp parallel do default(none) private(i,j,iv) OMP_SCHEDULE_ collapse(2) &
1420  !$omp shared(JSB,JEB,ISB,IEB,KS,KA,DENS,MOMZ,MOMX,MOMY,RHOT,VA,PROG,KE)
1421  do j = jsb, jeb
1422  do i = isb, ieb
1423  dens( 1:ks-1,i,j) = dens(ks,i,j)
1424  momz( 1:ks-1,i,j) = 0.0_rp
1425  momx( 1:ks-1,i,j) = momx(ks,i,j)
1426  momy( 1:ks-1,i,j) = momy(ks,i,j)
1427  rhot( 1:ks-1,i,j) = rhot(ks,i,j)
1428  do iv = 1, va
1429  prog( 1:ks-1,i,j,iv) = prog(ks,i,j,iv)
1430  end do
1431  dens(ke+1:ka, i,j) = dens(ke,i,j)
1432  momz(ke+1:ka, i,j) = 0.0_rp
1433  momx(ke+1:ka, i,j) = momx(ke,i,j)
1434  momy(ke+1:ka, i,j) = momy(ke,i,j)
1435  rhot(ke+1:ka, i,j) = rhot(ke,i,j)
1436  do iv = 1, va
1437  prog(ke+1:ka, i,j,iv) = prog(ke,i,j,iv)
1438  end do
1439  enddo
1440  enddo
1441 
1442  call comm_vars8( dens(:,:,:), i_comm_dens )
1443  call comm_vars8( momz(:,:,:), i_comm_momz )
1444  call comm_vars8( momx(:,:,:), i_comm_momx )
1445  call comm_vars8( momy(:,:,:), i_comm_momy )
1446  call comm_vars8( rhot(:,:,:), i_comm_rhot )
1447  do iv = 1, va
1448  call comm_vars8( prog(:,:,:,iv), i_comm_prog(iv) )
1449  end do
1450  call comm_wait ( dens(:,:,:), i_comm_dens, .false. )
1451  call comm_wait ( momz(:,:,:), i_comm_momz, .false. )
1452  call comm_wait ( momx(:,:,:), i_comm_momx, .false. )
1453  call comm_wait ( momy(:,:,:), i_comm_momy, .false. )
1454  call comm_wait ( rhot(:,:,:), i_comm_rhot, .false. )
1455  do iv = 1, va
1456  call comm_wait ( prog(:,:,:,iv), i_comm_prog(iv), .false. )
1457  end do
1458 
1459  if ( use_average ) then
1460  !$omp parallel do
1461  do j = jsb, jeb
1462  do i = isb, ieb
1463  do k = ks, ke
1464  dens_av(k,i,j) = dens_av(k,i,j) + dens(k,i,j) / nstep
1465  momz_av(k,i,j) = momz_av(k,i,j) + momz(k,i,j) / nstep
1466  momx_av(k,i,j) = momx_av(k,i,j) + momx(k,i,j) / nstep
1467  momy_av(k,i,j) = momy_av(k,i,j) + momy(k,i,j) / nstep
1468  rhot_av(k,i,j) = rhot_av(k,i,j) + rhot(k,i,j) / nstep
1469  end do
1470  end do
1471  end do
1472  endif
1473 
1474  !$omp parallel
1475  do n = 1, 3
1476  !$omp do
1477  do j = jsb, jeb
1478  do i = isb, ieb
1479  do k = ks, ke
1480  mflx_av(k,i,j,n) = mflx_av(k,i,j,n) + mflx(k,i,j,n)
1481  end do
1482  end do
1483  end do
1484  !$omp end do nowait
1485  end do
1486  !$omp end parallel
1487 
1488  enddo ! dynamical steps
1489 
1490 
1491  !###########################################################################
1492  ! Update Tracers
1493  !###########################################################################
1494 
1495  !$omp parallel
1496  do n = 1, 3
1497 !OCL XFILL
1498  !$omp do
1499  do j = jsb, jeb
1500  do i = isb, ieb
1501  do k = ks, ke
1502  mflx(k,i,j,n) = mflx_av(k,i,j,n) / nstep
1503  end do
1504  end do
1505  end do
1506  !$omp end do nowait
1507  end do
1508  !$omp end parallel
1509 
1510  call comm_vars8( mflx(:,:,:,zdir), i_comm_mflx_z )
1511  if ( .not. twod ) call comm_vars8( mflx(:,:,:,xdir), i_comm_mflx_x )
1512  call comm_vars8( mflx(:,:,:,ydir), i_comm_mflx_y )
1513  call comm_wait ( mflx(:,:,:,zdir), i_comm_mflx_z, .false. )
1514  if ( .not. twod ) call comm_wait ( mflx(:,:,:,xdir), i_comm_mflx_x, .false. )
1515  call comm_wait ( mflx(:,:,:,ydir), i_comm_mflx_y, .false. )
1516 
1517 #ifndef DRY
1518 
1519  !------------------------------------------------------------------------
1520  ! Update each tracer
1521  !------------------------------------------------------------------------
1522 
1523  do iq = 1, qa
1524 
1525  if ( tracer_advc(iq) ) then
1526 
1527  call prof_rapstart("DYN_Large_Numfilter", 2)
1528 
1529  if ( nd_coef_q == 0.0_rp ) then
1530 !OCL XFILL
1531  num_diff_q(:,:,:,:) = 0.0_rp
1532  else
1533  call atmos_dyn_fvm_numfilter_flux_q( num_diff_q(:,:,:,:), & ! [OUT]
1534  dens00, qtrc(:,:,:,iq), iq==i_qv, & ! [IN]
1535  cdz, cdx, cdy, twod, dtl, & ! [IN]
1536  ref_qv, iq, & ! [IN]
1537  nd_coef_q, nd_laplacian_num, nd_sfc_fact, nd_use_rs ) ! [IN]
1538  endif
1539 
1540  call prof_rapend ("DYN_Large_Numfilter", 2)
1541 
1542  call prof_rapstart("DYN_Tracer_Tinteg", 2)
1543 
1544  if ( flag_tracer_split_tend ) then
1545  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1546  !$omp shared(KA,IA,JA,iq,QTRC,RHOQ_t,DENS,dtl)
1547  do j = 1, ja
1548  do i = 1, ia
1549  do k = 1, ka
1550  qtrc(k,i,j,iq) = qtrc(k,i,j,iq) &
1551  + rhoq_t(k,i,j,iq) * dtl / dens(k,i,j)
1552  end do
1553  end do
1554  end do
1555  rhoq_tn => zero
1556  else
1557  rhoq_tn => rhoq_t(:,:,:,iq)
1558  end if
1559 
1560  call atmos_dyn_tinteg_tracer( &
1561  qtrc(:,:,:,iq), & ! (inout)
1562  qflx(:,:,:,:), & ! (out)
1563  qtrc0(:,:,:,iq), rhoq_tn, & ! (in)
1564  dens00, dens, & ! (in)
1565  mflx, num_diff_q, & ! (in)
1566  gsqrt, mapf(:,:,:,i_xy), & ! (in)
1567  cdz, rcdz, rcdx, rcdy, & ! (in)
1568  bnd_w, bnd_e, bnd_s, bnd_n, & ! (in)
1569  twod, & ! (in)
1570  dtl, & ! (in)
1571  llast .AND. flag_fct_tracer, & ! (in)
1572  flag_fct_along_stream ) ! (in)
1573 
1574  call prof_rapend ("DYN_Tracer_Tinteg", 2)
1575 
1576  if ( llast ) then
1577  do iv = 1, 3
1578  call file_history_query( hist_qflx(iv,iq), do_put )
1579  if ( do_put .or. monit_lateral_flag(iv) ) then
1580  call multiply_flux_by_metric_xyz( iv, qflx, gsqrt, mapf )
1581  call file_history_put( hist_qflx(iv,iq), qflx(:,:,:,iv) )
1582  end if
1583  end do
1584 
1585  if ( tracer_mass(iq) == 1.0_rp ) then
1586  if ( bnd_w .and. monit_qflx_west > 0 ) then
1587  !$omp parallel do
1588  do j = js, je
1589  do k = ks, ke
1590  qflx_west(k,j) = qflx_west(k,j) + qflx(k,is-1,j,xdir)
1591  end do
1592  end do
1593  end if
1594  if ( bnd_e .and. monit_qflx_east > 0 ) then
1595  !$omp parallel do
1596  do j = js, je
1597  do k = ks, ke
1598  qflx_east(k,j) = qflx_east(k,j) + qflx(k,ie,j,xdir)
1599  end do
1600  end do
1601  end if
1602  if ( bnd_s .and. monit_qflx_south > 0 ) then
1603  !$omp parallel do
1604  do i = is, ie
1605  do k = ks, ke
1606  qflx_south(k,i) = qflx_south(k,i) + qflx(k,i,js-1,ydir)
1607  end do
1608  end do
1609  end if
1610  if ( bnd_n .and. monit_qflx_north > 0 ) then
1611  !$omp parallel do
1612  do i = is, ie
1613  do k = ks, ke
1614  qflx_north(k,i) = qflx_north(k,i) + qflx(k,i,je,ydir)
1615  end do
1616  end do
1617  end if
1618 
1619  end if
1620 
1621  end if
1622 
1623  else
1624 
1625  !$omp parallel do
1626  do j = js, je
1627  do i = is, ie
1628  do k = ks, ke
1629  qtrc(k,i,j,iq) = ( qtrc0(k,i,j,iq) * dens00(k,i,j) &
1630  + rhoq_t(k,i,j,iq) * dtl ) / dens(k,i,j)
1631  end do
1632  end do
1633  end do
1634 
1635  end if
1636 
1637  call comm_vars8( qtrc(:,:,:,iq), i_comm_qtrc(iq) )
1638 
1639  if ( use_average ) then
1640  do j = jsb, jeb
1641  do i = isb, ieb
1642  do k = ks, ke
1643  qtrc_av(k,i,j,iq) = qtrc(k,i,j,iq)
1644  end do
1645  end do
1646  end do
1647  endif
1648 
1649  enddo ! scalar quantities loop
1650 
1651  do iq = 1, qa
1652  call comm_wait ( qtrc(:,:,:,iq), i_comm_qtrc(iq), .false. )
1653  enddo
1654 #endif
1655 
1656  if ( llast ) then
1657  !- history ---------------------------------------
1658  call file_history_put( hist_phys(1), dens_tp )
1659  call file_history_put( hist_phys(2), momz_tp )
1660  call file_history_put( hist_phys(3), momx_tp )
1661  call file_history_put( hist_phys(4), momy_tp )
1662  call file_history_put( hist_phys(5), rhot_tp )
1663 
1664  call file_history_put( hist_damp(1), damp_t_dens )
1665  call file_history_put( hist_damp(2), damp_t_momz )
1666  call file_history_put( hist_damp(3), damp_t_momx )
1667  call file_history_put( hist_damp(4), damp_t_momy )
1668  call file_history_put( hist_damp(5), damp_t_rhot )
1669 
1670  do iv = 1, 3
1671  call file_history_query( hist_mflx(iv), do_put )
1672  if ( do_put .or. monit_lateral_flag(iv) ) then
1673  call multiply_flux_by_metric_xyz( iv, mflx, gsqrt, mapf )
1674  call file_history_put( hist_mflx(iv), mflx(:,:,:,iv) )
1675  end if
1676  end do
1677 
1678  do iv = 1, 3
1679  call file_history_query( hist_tflx(iv), do_put )
1680  if ( do_put ) then
1681  call multiply_flux_by_metric_xyz( iv, tflx, gsqrt, mapf )
1682  call file_history_put( hist_tflx(iv), tflx(:,:,:,iv) )
1683  end if
1684  end do
1685 
1686  !- monitor mass budget ------------------------------------
1687  call monitor_put( monit_damp_mass, damp_t_dens(:,:,:) )
1688  call monitor_put_lateral_flux( monit_mflx_west, bnd_w .and. monit_mflx_west > 0, mflx(:,is-1,:,xdir), zero_x )
1689  call monitor_put_lateral_flux( monit_mflx_east, bnd_e .and. monit_mflx_east > 0, mflx(:,ie,:,xdir), zero_x )
1690  call monitor_put_lateral_flux( monit_mflx_south, bnd_s .and. monit_mflx_south > 0, mflx(:,:,js-1,ydir), zero_y )
1691  call monitor_put_lateral_flux( monit_mflx_north, bnd_n .and. monit_mflx_north > 0, mflx(:,:,je,ydir), zero_y )
1692 
1693  call monitor_put( monit_damp_qtot, dens_tq(:,:,:) )
1694  call monitor_put_lateral_flux( monit_qflx_west, bnd_w, qflx_west, zero_x )
1695  call monitor_put_lateral_flux( monit_qflx_east, bnd_e, qflx_east, zero_x )
1696  call monitor_put_lateral_flux( monit_qflx_south, bnd_s, qflx_south, zero_y )
1697  call monitor_put_lateral_flux( monit_qflx_north, bnd_n, qflx_north, zero_y )
1698  end if
1699 
1700  return

References scale_atmos_dyn_common::atmos_dyn_divergence(), scale_atmos_dyn_common::atmos_dyn_fill_halo(), scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_apply_numfilter(), scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_flux(), scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_flux_q(), scale_atmos_dyn_common::atmos_dyn_prep_pres_linearization(), scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short, scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer, scale_const::const_eps, scale_dft::dft_g2g(), scale_dft::dft_g2g_divfree(), scale_file_history::file_history_set_disable(), scale_atmos_grid_cartesc_index::i_uy, scale_atmos_grid_cartesc_index::i_uyz, scale_atmos_grid_cartesc_index::i_xv, scale_atmos_grid_cartesc_index::i_xvz, scale_atmos_grid_cartesc_index::i_xy, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::ieb, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::jsb, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, monitor_put_lateral_flux(), multiply_flux_by_metric_xyz(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, scale_spnudge::spnudge_pt, scale_spnudge::spnudge_pt_alpha, scale_spnudge::spnudge_pt_lm, scale_spnudge::spnudge_pt_mm, scale_spnudge::spnudge_qv, scale_spnudge::spnudge_qv_alpha, scale_spnudge::spnudge_qv_lm, scale_spnudge::spnudge_qv_mm, scale_spnudge::spnudge_u_alpha, scale_spnudge::spnudge_uv, scale_spnudge::spnudge_uv_divfree, scale_spnudge::spnudge_uv_lm, scale_spnudge::spnudge_uv_mm, scale_spnudge::spnudge_v_alpha, scale_tracer::tracer_advc, scale_tracer::tracer_mass, scale_index::va, scale_atmos_grid_cartesc_index::xdir, scale_atmos_grid_cartesc_index::ydir, and scale_atmos_grid_cartesc_index::zdir.

Referenced by atmos_dyn_tstep_large_fvm_heve_setup(), and scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup().

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

◆ monitor_put_lateral_flux()

subroutine scale_atmos_dyn_tstep_large_fvm_heve::monitor_put_lateral_flux ( integer, intent(in)  MONIT_ID,
logical, intent(in)  BND_flag,
real(rp), dimension(:,:), intent(in), target  flx,
real(rp), dimension(:,:), intent(in), target  flx_zero 
)

Definition at line 1706 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

1706  use scale_monitor, only: &
1707  monitor_put
1708  implicit none
1709 
1710  integer, intent(in) :: MONIT_ID
1711  logical, intent(in) :: BND_flag
1712  real(RP), target, intent(in) :: flx(:,:)
1713  real(RP), target, intent(in) :: flx_zero(:,:)
1714 
1715  real(RP), pointer :: flx_ptr(:,:)
1716  !------------------------------------------
1717  if ( bnd_flag ) then
1718  flx_ptr => flx
1719  else
1720  flx_ptr => flx_zero
1721  end if
1722  call monitor_put( monit_id, flx_ptr(:,:) )
1723 
1724  return

Referenced by atmos_dyn_tstep_large_fvm_heve().

Here is the caller graph for this function:

◆ multiply_flux_by_metric_xyz()

subroutine scale_atmos_dyn_tstep_large_fvm_heve::multiply_flux_by_metric_xyz ( integer, intent(in)  I_DIR,
real(rp), dimension(ka,ia,ja,3), intent(inout)  flx,
real(rp), dimension(ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension (ia,ja,2,4), intent(in)  MAPF 
)

Definition at line 1728 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

1728  implicit none
1729  integer, intent(in) :: I_DIR
1730  real(RP), intent(inout) :: flx(KA,IA,JA,3)
1731  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
1732  real(RP), intent(in) :: MAPF (IA,JA,2,4)
1733 
1734  integer :: i, j, k
1735  !---------------------------------------------------
1736 
1737  if (i_dir == zdir) then
1738  !$omp parallel do
1739  do j = js, je
1740  do i = is, ie
1741  do k = ks-1, ke
1742  flx(k,i,j,zdir) = flx(k,i,j,zdir) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyw)
1743  end do
1744  end do
1745  end do
1746  end if
1747 
1748  if (i_dir == xdir) then
1749  !$omp parallel do
1750  do j = js, je
1751  do i = isb, ieb
1752  do k = ks, ke
1753  flx(k,i,j,xdir) = flx(k,i,j,xdir) * mapf(i,j,2,i_uy) / gsqrt(k,i,j,i_uyz)
1754  end do
1755  end do
1756  end do
1757  end if
1758 
1759  if (i_dir == ydir) then
1760  !$omp parallel do
1761  do j = jsb, jeb
1762  do i = is, ie
1763  do k = ks, ke
1764  flx(k,i,j,ydir) = flx(k,i,j,ydir) * mapf(i,j,1,i_xv) / gsqrt(k,i,j,i_xvz)
1765  end do
1766  end do
1767  end do
1768  end if
1769 
1770  return

References scale_atmos_grid_cartesc_index::i_uy, scale_atmos_grid_cartesc_index::i_uyz, scale_atmos_grid_cartesc_index::i_xv, scale_atmos_grid_cartesc_index::i_xvz, scale_atmos_grid_cartesc_index::i_xy, scale_atmos_grid_cartesc_index::i_xyw, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::ieb, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::jsb, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_atmos_grid_cartesc_index::xdir, scale_atmos_grid_cartesc_index::ydir, and scale_atmos_grid_cartesc_index::zdir.

Referenced by atmos_dyn_tstep_large_fvm_heve().

Here is the caller graph for this function:
scale_atmos_grid_cartesc_index::isb
integer, public isb
Definition: scale_atmos_grid_cartesC_index.F90:63
scale_atmos_grid_cartesc_index::i_uy
integer, public i_uy
Definition: scale_atmos_grid_cartesC_index.F90:99
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_atmos_grid_cartesc_index::i_xv
integer, public i_xv
Definition: scale_atmos_grid_cartesC_index.F90:100
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer
procedure(tinteg), pointer, public atmos_dyn_tinteg_tracer
Definition: scale_atmos_dyn_tinteg_tracer.F90:74
scale_spnudge::spnudge_uv_lm
integer, public spnudge_uv_lm
Definition: scale_spnudge.F90:26
scale_spnudge::spnudge_pt_lm
integer, public spnudge_pt_lm
Definition: scale_spnudge.F90:30
scale_spnudge::spnudge_pt
logical, public spnudge_pt
Definition: scale_spnudge.F90:29
scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_flux
subroutine, public atmos_dyn_fvm_numfilter_flux(num_diff, DENS, MOMZ, MOMX, MOMY, RHOT, CDZ, CDX, CDY, FDZ, FDX, FDY, TwoD, DT, REF_dens, REF_pott, ND_COEF, ND_LAPLACIAN_NUM, ND_SFC_FACT, ND_USE_RS)
Calculate fluxes with numerical filter for prognostic variables of dynamical core.
Definition: scale_atmos_dyn_fvm_numfilter.F90:390
scale_prc_cartesc::prc_has_s
logical, public prc_has_s
Definition: scale_prc_cartesC.F90:50
scale_atmos_dyn_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
scale_spnudge::spnudge_pt_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_pt_alpha
Definition: scale_spnudge.F90:39
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_dft::dft_g2g_divfree
subroutine, public dft_g2g_divfree(KA, KS, KE, IA, IS, IE, JA, JS, JE, LM, MM, u, v)
Definition: scale_dft.F90:227
scale_spnudge::spnudge_uv_divfree
logical, public spnudge_uv_divfree
Definition: scale_spnudge.F90:25
scale_file_history::file_history_set_disable
subroutine, public file_history_set_disable(switch)
set switch to turn on/off history
Definition: scale_file_history.F90:1712
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_spnudge::spnudge_qv_lm
integer, public spnudge_qv_lm
Definition: scale_spnudge.F90:34
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:48
scale_spnudge::spnudge_qv
logical, public spnudge_qv
Definition: scale_spnudge.F90:33
scale_atmos_grid_cartesc_index::jeb
integer, public jeb
Definition: scale_atmos_grid_cartesC_index.F90:66
scale_atmos_dyn_tinteg_short
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_short.F90:11
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:49
scale_index::va
integer, public va
Definition: scale_index.F90:35
scale_dft
Definition: scale_dft.F90:3
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_spnudge::spnudge_pt_mm
integer, public spnudge_pt_mm
Definition: scale_spnudge.F90:31
scale_atmos_grid_cartesc_index::i_xy
integer, public i_xy
Definition: scale_atmos_grid_cartesC_index.F90:98
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_atmos_grid_cartesc_index::i_uyz
integer, public i_uyz
Definition: scale_atmos_grid_cartesC_index.F90:94
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
scale_spnudge::spnudge_u_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_u_alpha
Definition: scale_spnudge.F90:37
scale_atmos_grid_cartesc_index::zdir
integer, parameter, public zdir
Definition: scale_atmos_grid_cartesC_index.F90:32
scale_atmos_grid_cartesc_index::i_xvz
integer, public i_xvz
Definition: scale_atmos_grid_cartesC_index.F90:95
scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve
subroutine, public atmos_dyn_tstep_large_fvm_heve(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, num_diff, num_diff_q, QTRC0, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLI, 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, BND_W, BND_E, BND_S, BND_N, TwoD, 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, wdamp_coef, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTLS, DTSS, Llast)
Dynamical Process.
Definition: scale_atmos_dyn_tstep_large_fvm_heve.F90:407
scale_atmos_dyn_common::atmos_dyn_fill_halo
subroutine, public atmos_dyn_fill_halo(var, fill_constval, lateral_halo, top_bottom_halo)
Definition: scale_atmos_dyn_common.F90:129
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_spnudge::spnudge_qv_mm
integer, public spnudge_qv_mm
Definition: scale_spnudge.F90:35
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_dyn_fvm_numfilter
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_fvm_numfilter.F90:12
scale_spnudge::spnudge_uv_mm
integer, public spnudge_uv_mm
Definition: scale_spnudge.F90:27
scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short
procedure(short), pointer, public atmos_dyn_tinteg_short
Definition: scale_atmos_dyn_tinteg_short.F90:109
scale_monitor::monitor_reg
subroutine, public monitor_reg(name, desc, unit, itemid, ndims, dim_type, is_tendency)
Search existing item, or matching check between requested and registered item.
Definition: scale_monitor.F90:241
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_flux_q
subroutine, public atmos_dyn_fvm_numfilter_flux_q(num_diff_q, DENS, QTRC, is_qv, CDZ, CDX, CDY, TwoD, dt, REF_qv, iq, ND_COEF, ND_LAPLACIAN_NUM, ND_SFC_FACT, ND_USE_RS)
Calculate fluxes with numerical filter for tracer variables.
Definition: scale_atmos_dyn_fvm_numfilter.F90:1160
scale_dft::dft_g2g
subroutine, public dft_g2g(KA, KS, KE, IA, IS, IE, JA, JS, JE, LM, MM, f)
Definition: scale_dft.F90:214
scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_apply_numfilter
subroutine, public atmos_dyn_fvm_apply_numfilter(num_diff, DENS, MOMZ, MOMX, MOMY, RHOT, CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, TwoD, DT, GSQRT, MAPF, REF_dens, REF_pott, ND_COEF, ND_LAPLACIAN_NUM, ND_SFC_FACT, ND_USE_RS)
Definition: scale_atmos_dyn_fvm_numfilter.F90:1490
scale_spnudge
Definition: scale_spnudge.F90:2
scale_atmos_grid_cartesc_index::xdir
integer, parameter, public xdir
Definition: scale_atmos_grid_cartesC_index.F90:33
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_atmos_dyn_tinteg_tracer
module Atmosphere / Dynamics Temporal integration
Definition: scale_atmos_dyn_tinteg_tracer.F90:12
scale_atmos_dyn_common::atmos_dyn_divergence
subroutine, public atmos_dyn_divergence(DDIV, MOMZ, MOMX, MOMY, GSQRT, J13G, J23G, J33G, MAPF, TwoD, RCDZ, RCDX, RCDY, RFDZ, FDZ)
Definition: scale_atmos_dyn_common.F90:382
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc_index::ydir
integer, parameter, public ydir
Definition: scale_atmos_grid_cartesC_index.F90:34
scale_atmos_grid_cartesc_index::ieb
integer, public ieb
Definition: scale_atmos_grid_cartesC_index.F90:64
scale_spnudge::spnudge_qv_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_qv_alpha
Definition: scale_spnudge.F90:40
scale_file_history::file_history_reg
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
Definition: scale_file_history.F90:650
scale_atmos_dyn_common::atmos_dyn_prep_pres_linearization
subroutine, public atmos_dyn_prep_pres_linearization(DPRES, RT2P, REF_rhot, RHOT, QTRC, REF_pres, AQ_R, AQ_CV, AQ_CP, AQ_MASS)
Definition: scale_atmos_dyn_common.F90:510
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_atmos_grid_cartesc_index::jsb
integer, public jsb
Definition: scale_atmos_grid_cartesC_index.F90:65
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:47
scale_const::const_ohm
real(rp), public const_ohm
angular velocity of the planet [1/s]
Definition: scale_const.F90:45
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:55
scale_spnudge::spnudge_uv
logical, public spnudge_uv
Definition: scale_spnudge.F90:24
scale_spnudge::spnudge_v_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_v_alpha
Definition: scale_spnudge.F90:38
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:91
scale_monitor
module MONITOR
Definition: scale_monitor.F90:12