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_finalize
 finalize 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 124 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

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

subroutine, public scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve_finalize

finalize

Definition at line 398 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

398 
399  !$acc exit data delete(DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, RHOQ_t)
400  deallocate( dens_t )
401  deallocate( momz_t )
402  deallocate( momx_t )
403  deallocate( momy_t )
404  deallocate( rhot_t )
405  deallocate( rhoq_t )
406 
407  !$acc exit data delete(DENS_damp, mflx)
408  deallocate( dens_damp )
409 
410  deallocate( mflx )
411 
412  deallocate( i_comm_prog )
413  deallocate( i_comm_qtrc )
414  deallocate( i_comm_rhoq_t )
415 
416  !$acc exit data delete(ZERO)
417  deallocate( zero )
418 
419  deallocate( hist_qflx )
420  deallocate( hist_phys_qtrc )
421  deallocate( hist_damp_qtrc )
422 
423  !$acc exit data delete(zero_x, zero_y)
424  deallocate( zero_x )
425  deallocate( zero_y )
426 
427  return

Referenced by scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup().

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

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

1923  use scale_monitor, only: &
1924  monitor_put
1925  implicit none
1926 
1927  integer, intent(in) :: MONIT_ID
1928  logical, intent(in) :: BND_flag
1929  real(RP), target, intent(in) :: flx(:,:)
1930  real(RP), target, intent(in) :: flx_zero(:,:)
1931 
1932  real(RP), pointer :: flx_ptr(:,:)
1933  !------------------------------------------
1934  if ( bnd_flag ) then
1935  flx_ptr => flx
1936  else
1937  flx_ptr => flx_zero
1938  end if
1939  call monitor_put( monit_id, flx_ptr(:,:) )
1940 
1941  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 1945 of file scale_atmos_dyn_tstep_large_fvm_heve.F90.

1945  implicit none
1946  integer, intent(in) :: I_DIR
1947  real(RP), intent(inout) :: flx(KA,IA,JA,3)
1948  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
1949  real(RP), intent(in) :: MAPF (IA,JA,2,4)
1950 
1951  integer :: i, j, k
1952  !---------------------------------------------------
1953 
1954  !$acc data copy(flx) copyin(GSQRT, MAPF)
1955 
1956  if (i_dir == zdir) then
1957  !$omp parallel do
1958  !$acc kernels
1959  do j = js, je
1960  do i = is, ie
1961  do k = ks-1, ke
1962  flx(k,i,j,zdir) = flx(k,i,j,zdir) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy)
1963  end do
1964  end do
1965  end do
1966  !$acc end kernels
1967  end if
1968 
1969  if (i_dir == xdir) then
1970  !$omp parallel do
1971  !$acc kernels
1972  do j = js, je
1973  do i = isb, ieb
1974  do k = ks, ke
1975  flx(k,i,j,xdir) = flx(k,i,j,xdir) * mapf(i,j,2,i_uy) / gsqrt(k,i,j,i_uyz)
1976  end do
1977  end do
1978  end do
1979  !$acc end kernels
1980  end if
1981 
1982  if (i_dir == ydir) then
1983  !$omp parallel do
1984  !$acc kernels
1985  do j = jsb, jeb
1986  do i = is, ie
1987  do k = ks, ke
1988  flx(k,i,j,ydir) = flx(k,i,j,ydir) * mapf(i,j,1,i_xv) / gsqrt(k,i,j,i_xvz)
1989  end do
1990  end do
1991  end do
1992  !$acc end kernels
1993  end if
1994 
1995  !$acc end data
1996 
1997  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::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:64
scale_atmos_grid_cartesc_index::i_uy
integer, public i_uy
Definition: scale_atmos_grid_cartesC_index.F90:100
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:101
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer
procedure(tinteg), pointer, public atmos_dyn_tinteg_tracer
Definition: scale_atmos_dyn_tinteg_tracer.F90:76
scale_spnudge::spnudge_uv_lm
integer, public spnudge_uv_lm
Definition: scale_spnudge.F90:27
scale_spnudge::spnudge_pt_lm
integer, public spnudge_pt_lm
Definition: scale_spnudge.F90:31
scale_spnudge::spnudge_pt
logical, public spnudge_pt
Definition: scale_spnudge.F90:30
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:428
scale_prc_cartesc::prc_has_s
logical, public prc_has_s
Definition: scale_prc_cartesC.F90:51
scale_atmos_dyn_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
scale_spnudge::spnudge_pt_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_pt_alpha
Definition: scale_spnudge.F90:40
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:270
scale_spnudge::spnudge_uv_divfree
logical, public spnudge_uv_divfree
Definition: scale_spnudge.F90:26
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:1980
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
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:35
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:49
scale_spnudge::spnudge_qv
logical, public spnudge_qv
Definition: scale_spnudge.F90:34
scale_atmos_grid_cartesc_index::jeb
integer, public jeb
Definition: scale_atmos_grid_cartesC_index.F90:67
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:50
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:32
scale_atmos_grid_cartesc_index::i_xy
integer, public i_xy
Definition: scale_atmos_grid_cartesC_index.F90:99
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:95
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_spnudge::spnudge_u_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_u_alpha
Definition: scale_spnudge.F90:38
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:96
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:461
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:144
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:36
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:28
scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short
procedure(short), pointer, public atmos_dyn_tinteg_short
Definition: scale_atmos_dyn_tinteg_short.F90:112
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:243
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:1302
scale_dft::dft_g2g
subroutine, public dft_g2g(KA, KS, KE, IA, IS, IE, JA, JS, JE, LM, MM, f)
Definition: scale_dft.F90:257
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:1674
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:442
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:65
scale_spnudge::spnudge_qv_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_qv_alpha
Definition: scale_spnudge.F90:41
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:685
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:583
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid, gid)
Register variables.
Definition: scale_comm_cartesC.F90:811
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_grid_cartesc_index::jsb
integer, public jsb
Definition: scale_atmos_grid_cartesC_index.F90:66
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:48
scale_const::const_ohm
real(rp), public const_ohm
angular velocity of the planet [1/s]
Definition: scale_const.F90:48
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:56
scale_spnudge::spnudge_uv
logical, public spnudge_uv
Definition: scale_spnudge.F90:25
scale_spnudge::spnudge_v_alpha
real(rp), dimension(:,:,:), allocatable, public spnudge_v_alpha
Definition: scale_spnudge.F90:39
scale_monitor
module MONITOR
Definition: scale_monitor.F90:12