SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tinteg_short_rk3 Module Reference

module Atmosphere / Dyn Tinteg More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_short_rk3_setup (tinteg_type)
 Setup. More...
 
subroutine, public atmos_dyn_tinteg_short_rk3 (DENS, MOMZ, MOMX, MOMY, RHOT, PROG, mflx_hi, tflx_hi, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, Rtot, CVtot, CORIOLI, num_diff, divdmp_coef, DDIV, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_ALONG_STREAM, CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, REF_pres, REF_dens, BND_W, BND_E, BND_S, BND_N, dt)
 RK3. More...
 

Detailed Description

module Atmosphere / Dyn Tinteg

Description
Temporal integration in Dynamical core for Atmospheric process three step Runge-Kutta scheme
Author
Team SCALE
History

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_short_rk3_setup()

subroutine, public scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3_setup ( character(len=*)  tinteg_type)

Setup.

Definition at line 94 of file scale_atmos_dyn_tinteg_short_rk3.F90.

References scale_comm::comm_vars8_init(), scale_const::const_undef, scale_grid_index::ia, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::ja, scale_grid_index::ka, scale_process::prc_mpistop(), and scale_index::va.

Referenced by scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup().

94  use scale_process, only: &
96  use scale_const, only: &
97  undef => const_undef
98  use scale_comm, only: &
100  implicit none
101 
102  character(len=*) :: tinteg_type
103 
104  integer :: iv
105  !---------------------------------------------------------------------------
106 
107  select case ( tinteg_type )
108  case ( 'RK3' )
109  if( io_l ) write(io_fid_log,*) "*** RK3: Heun's method is used"
110  ! Heun's method
111  ! k1 = f(\phi_n); r1 = \phi_n + k1 * dt / 3
112  ! k2 = f(r1); r2 = \phi_n + k2 * dt * 2 / 3
113  ! k3 = f(r2); r3 = \phi_n + k3 * dt
114  ! \phi_{n+1} = \phi_n + ( k1 + 3 * k3 ) dt / 4
115  ! = \phi_n + ( (r1-\phi_n) * 3 + (r3-\phi_n) * 3 ) / 4
116  ! = ( r1 * 3 + r3 * 3 - \phi_n * 2 ) / 4
117  flag_ws2002 = .false.
118  fact_dt1 = 1.0_rp / 3.0_rp
119  fact_dt2 = 2.0_rp / 3.0_rp
120  case ( 'RK3WS2002' )
121  if( io_l ) write(io_fid_log,*) "*** RK3: Wichere and Skamarock (2002) is used"
122  ! Wicher and Skamarock (2002) RK3 scheme
123  ! k1 = f(\phi_n); r1 = \phi_n + k1 * dt / 3
124  ! k2 = f(r1); r2 = \phi_n + k2 * dt / 2
125  ! k3 = f(r2); r3 = \phi_n + k3 * dt
126  ! \phi_{n+1} = r3
127  flag_ws2002 = .true.
128  fact_dt1 = 1.0_rp / 3.0_rp
129  fact_dt2 = 1.0_rp / 2.0_rp
130  case default
131  write(*,*) 'xxx TINTEG_TYPE is not RK3. Check!'
132  call prc_mpistop
133  end select
134 
135  allocate( dens_rk1(ka,ia,ja) )
136  allocate( momz_rk1(ka,ia,ja) )
137  allocate( momx_rk1(ka,ia,ja) )
138  allocate( momy_rk1(ka,ia,ja) )
139  allocate( rhot_rk1(ka,ia,ja) )
140 
141  allocate( dens_rk2(ka,ia,ja) )
142  allocate( momz_rk2(ka,ia,ja) )
143  allocate( momx_rk2(ka,ia,ja) )
144  allocate( momy_rk2(ka,ia,ja) )
145  allocate( rhot_rk2(ka,ia,ja) )
146 
147  allocate( prog_rk1(ka,ia,ja,max(va,1)) )
148  allocate( prog_rk2(ka,ia,ja,max(va,1)) )
149  allocate( i_comm_prog_rk1(max(va,1)) )
150  allocate( i_comm_prog_rk2(max(va,1)) )
151 
152  call comm_vars8_init( dens_rk1, i_comm_dens_rk1 )
153  call comm_vars8_init( momz_rk1, i_comm_momz_rk1 )
154  call comm_vars8_init( momx_rk1, i_comm_momx_rk1 )
155  call comm_vars8_init( momy_rk1, i_comm_momy_rk1 )
156  call comm_vars8_init( rhot_rk1, i_comm_rhot_rk1 )
157  do iv = 1, va
158  i_comm_prog_rk1(iv) = 5 + iv
159  call comm_vars8_init( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
160  enddo
161 
162  call comm_vars8_init( dens_rk2, i_comm_dens_rk2 )
163  call comm_vars8_init( momz_rk2, i_comm_momz_rk2 )
164  call comm_vars8_init( momx_rk2, i_comm_momx_rk2 )
165  call comm_vars8_init( momy_rk2, i_comm_momy_rk2 )
166  call comm_vars8_init( rhot_rk2, i_comm_rhot_rk2 )
167  do iv = 1, va
168  i_comm_prog_rk2(iv) = 5 + iv
169  call comm_vars8_init( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
170  enddo
171 
172  dens_rk1(:,:,:) = undef
173  momz_rk1(:,:,:) = undef
174  momx_rk1(:,:,:) = undef
175  momy_rk1(:,:,:) = undef
176  rhot_rk1(:,:,:) = undef
177  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
178 
179  dens_rk2(:,:,:) = undef
180  momz_rk2(:,:,:) = undef
181  momx_rk2(:,:,:) = undef
182  momy_rk2(:,:,:) = undef
183  rhot_rk2(:,:,:) = undef
184  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
185 
186  return
subroutine, public prc_mpistop
Abort MPI.
integer, public va
Definition: scale_index.F90:38
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), public const_undef
Definition: scale_const.F90:43
integer, public ia
of x whole cells (local, with HALO)
subroutine, public comm_vars8_init(var, vid)
Register variables.
Definition: scale_comm.F90:328
integer, public ka
of z whole cells (local, with HALO)
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tinteg_short_rk3()

subroutine, public scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3 ( 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,va), intent(inout)  PROG,
real(rp), dimension(ka,ia,ja,3), intent(inout)  mflx_hi,
real(rp), dimension(ka,ia,ja,3), intent(inout)  tflx_hi,
real(rp), dimension(ka,ia,ja), intent(in)  DENS_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY_t,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT_t,
real(rp), dimension(ka,ia,ja), intent(in)  Rtot,
real(rp), dimension(ka,ia,ja), intent(in)  CVtot,
real(rp), dimension(ia,ja), intent(in)  CORIOLI,
real(rp), dimension(ka,ia,ja,5,3), intent(in)  num_diff,
real(rp), intent(in)  divdmp_coef,
real(rp), dimension(ka,ia,ja), intent(in)  DDIV,
logical, intent(in)  FLAG_FCT_MOMENTUM,
logical, intent(in)  FLAG_FCT_T,
logical, intent(in)  FLAG_FCT_ALONG_STREAM,
real(rp), dimension (ka), intent(in)  CDZ,
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(ka,ia,ja), intent(in)  REF_pres,
real(rp), dimension(ka,ia,ja), intent(in)  REF_dens,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
real(rp), intent(in)  dt 
)

RK3.

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 205 of file scale_atmos_dyn_tinteg_short_rk3.F90.

References scale_atmos_dyn_common::atmos_dyn_copy_boundary(), scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_index::va.

Referenced by scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup().

205  use scale_comm, only: &
206  comm_vars8, &
207  comm_wait
208  use scale_atmos_dyn_tstep_short, only: &
209  atmos_dyn_tstep => atmos_dyn_tstep_short
210  use scale_atmos_dyn_common, only: &
212  implicit none
213 
214  real(RP), intent(inout) :: dens(ka,ia,ja)
215  real(RP), intent(inout) :: momz(ka,ia,ja)
216  real(RP), intent(inout) :: momx(ka,ia,ja)
217  real(RP), intent(inout) :: momy(ka,ia,ja)
218  real(RP), intent(inout) :: rhot(ka,ia,ja)
219  real(RP), intent(inout) :: prog(ka,ia,ja,va)
220 
221  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3)
222  real(RP), intent(inout) :: tflx_hi(ka,ia,ja,3)
223 
224  real(RP), intent(in) :: dens_t(ka,ia,ja)
225  real(RP), intent(in) :: momz_t(ka,ia,ja)
226  real(RP), intent(in) :: momx_t(ka,ia,ja)
227  real(RP), intent(in) :: momy_t(ka,ia,ja)
228  real(RP), intent(in) :: rhot_t(ka,ia,ja)
229 
230  real(RP), intent(in) :: rtot(ka,ia,ja)
231  real(RP), intent(in) :: cvtot(ka,ia,ja)
232  real(RP), intent(in) :: corioli(ia,ja)
233  real(RP), intent(in) :: num_diff(ka,ia,ja,5,3)
234  real(RP), intent(in) :: divdmp_coef
235  real(RP), intent(in) :: ddiv(ka,ia,ja)
236 
237  logical, intent(in) :: flag_fct_momentum
238  logical, intent(in) :: flag_fct_t
239  logical, intent(in) :: flag_fct_along_stream
240 
241  real(RP), intent(in) :: cdz (ka)
242  real(RP), intent(in) :: fdz (ka-1)
243  real(RP), intent(in) :: fdx (ia-1)
244  real(RP), intent(in) :: fdy (ja-1)
245  real(RP), intent(in) :: rcdz(ka)
246  real(RP), intent(in) :: rcdx(ia)
247  real(RP), intent(in) :: rcdy(ja)
248  real(RP), intent(in) :: rfdz(ka-1)
249  real(RP), intent(in) :: rfdx(ia-1)
250  real(RP), intent(in) :: rfdy(ja-1)
251 
252  real(RP), intent(in) :: phi (ka,ia,ja)
253  real(RP), intent(in) :: gsqrt(ka,ia,ja,7)
254  real(RP), intent(in) :: j13g (ka,ia,ja,7)
255  real(RP), intent(in) :: j23g (ka,ia,ja,7)
256  real(RP), intent(in) :: j33g
257  real(RP), intent(in) :: mapf (ia,ja,2,4)
258 
259  real(RP), intent(in) :: ref_pres(ka,ia,ja)
260  real(RP), intent(in) :: ref_dens(ka,ia,ja)
261 
262  logical, intent(in) :: bnd_w
263  logical, intent(in) :: bnd_e
264  logical, intent(in) :: bnd_s
265  logical, intent(in) :: bnd_n
266 
267  real(RP), intent(in) :: dt
268 
269  real(RP) :: dens0(ka,ia,ja)
270  real(RP) :: momz0(ka,ia,ja)
271  real(RP) :: momx0(ka,ia,ja)
272  real(RP) :: momy0(ka,ia,ja)
273  real(RP) :: rhot0(ka,ia,ja)
274  real(RP) :: prog0(ka,ia,ja,va)
275 
276  real(RP) :: mflx_hi_rk(ka,ia,ja,3,2)
277  real(RP) :: tflx_hi_rk(ka,ia,ja,3,2)
278 
279  real(RP) :: dtrk
280 
281  integer :: i, j, k, iv, n
282  !---------------------------------------------------------------------------
283 
284  call prof_rapstart("DYN_RK3_Prep",3)
285 
286 #ifdef DEBUG
287  dens_rk1(:,:,:) = undef
288  momz_rk1(:,:,:) = undef
289  momx_rk1(:,:,:) = undef
290  momy_rk1(:,:,:) = undef
291  rhot_rk1(:,:,:) = undef
292  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
293 
294  dens_rk2(:,:,:) = undef
295  momz_rk2(:,:,:) = undef
296  momx_rk2(:,:,:) = undef
297  momy_rk2(:,:,:) = undef
298  rhot_rk2(:,:,:) = undef
299  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
300 
301  mflx_hi(:,:,:,:) = undef
302  tflx_hi(:,:,:,:) = undef
303 
304  mflx_hi_rk(:,:,:,:,:) = undef
305  tflx_hi_rk(:,:,:,:,:) = undef
306 #endif
307 
308 !OCL XFILL
309  dens0 = dens
310 !OCL XFILL
311  momz0 = momz
312 !OCL XFILL
313  momx0 = momx
314 !OCL XFILL
315  momy0 = momy
316 !OCL XFILL
317  rhot0 = rhot
318 !OCL XFILL
319  if ( va > 0 ) prog0 = prog
320 
321  if ( bnd_w ) then
322  do j = js, je
323  do k = ks, ke
324  mflx_hi_rk(k,is-1,j,2,:) = mflx_hi(k,is-1,j,2)
325  end do
326  end do
327  end if
328  if ( bnd_e ) then
329  do j = js, je
330  do k = ks, ke
331  mflx_hi_rk(k,ie,j,2,:) = mflx_hi(k,ie,j,2)
332  end do
333  end do
334  end if
335  if ( bnd_s ) then
336  do i = is, ie
337  do k = ks, ke
338  mflx_hi_rk(k,i,js-1,3,:) = mflx_hi(k,i,js-1,3)
339  end do
340  end do
341  end if
342  if ( bnd_n ) then
343  do i = is, ie
344  do k = ks, ke
345  mflx_hi_rk(k,i,je,3,:) = mflx_hi(k,i,je,3)
346  end do
347  end do
348  end if
349 
350  call prof_rapend ("DYN_RK3_Prep",3)
351 
352  !------------------------------------------------------------------------
353  ! Start RK
354  !------------------------------------------------------------------------
355 
356  !##### RK1 : PROG0,PROG->PROG_RK1 #####
357 
358  call prof_rapstart("DYN_RK3",3)
359 
360  dtrk = dt * fact_dt1
361 
362  call atmos_dyn_tstep( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [OUT]
363  prog_rk1, & ! [OUT]
364  mflx_hi_rk(:,:,:,:,1), tflx_hi_rk(:,:,:,:,1), & ! [INOUT]
365  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
366  dens, momz, momx, momy, rhot, & ! [IN]
367  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
368  prog0, prog, & ! [IN]
369  rtot, cvtot, corioli, & ! [IN]
370  num_diff, divdmp_coef, ddiv, & ! [IN]
371  flag_fct_momentum, flag_fct_t, & ! [IN]
372  flag_fct_along_stream, & ! [IN]
373  cdz, fdz, fdx, fdy, & ! [IN]
374  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
375  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
376  ref_pres, ref_dens, & ! [IN]
377  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
378  dtrk, dt ) ! [IN]
379 
380  call prof_rapend ("DYN_RK3",3)
381  call prof_rapstart("DYN_RK3_BND",3)
382 
383  call atmos_dyn_copy_boundary( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [INOUT]
384  prog_rk1, & ! [INOUT]
385  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
386  prog0, & ! [IN]
387  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
388 
389  call prof_rapend ("DYN_RK3_BND",3)
390 
391  call comm_vars8( dens_rk1(:,:,:), i_comm_dens_rk1 )
392  call comm_vars8( momz_rk1(:,:,:), i_comm_momz_rk1 )
393  call comm_vars8( momx_rk1(:,:,:), i_comm_momx_rk1 )
394  call comm_vars8( momy_rk1(:,:,:), i_comm_momy_rk1 )
395  call comm_vars8( rhot_rk1(:,:,:), i_comm_rhot_rk1 )
396  do iv = 1, va
397  call comm_vars8( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
398  enddo
399 
400  call comm_wait ( dens_rk1(:,:,:), i_comm_dens_rk1, .false. )
401  call comm_wait ( momz_rk1(:,:,:), i_comm_momz_rk1, .false. )
402  call comm_wait ( momx_rk1(:,:,:), i_comm_momx_rk1, .false. )
403  call comm_wait ( momy_rk1(:,:,:), i_comm_momy_rk1, .false. )
404  call comm_wait ( rhot_rk1(:,:,:), i_comm_rhot_rk1, .false. )
405  do iv = 1, va
406  call comm_wait ( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv), .false. )
407  enddo
408 
409  !##### RK2 : PROG0,PROG_RK1->PROG_RK2 #####
410 
411  call prof_rapstart("DYN_RK3",3)
412 
413  dtrk = dt * fact_dt2
414 
415  call atmos_dyn_tstep( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [OUT]
416  prog_rk2, & ! [OUT]
417  mflx_hi_rk(:,:,:,:,2), tflx_hi_rk(:,:,:,:,2), & ! [INOUT]
418  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
419  dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [IN]
420  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
421  prog0, prog_rk1, & ! [IN]
422  rtot, cvtot, corioli, & ! [IN]
423  num_diff, divdmp_coef, ddiv, & ! [IN]
424  flag_fct_momentum, flag_fct_t, & ! [IN]
425  flag_fct_along_stream, & ! [IN]
426  cdz, fdz, fdx, fdy, & ! [IN]
427  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
428  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
429  ref_pres, ref_dens, & ! [IN]
430  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
431  dtrk, dt ) ! [IN]
432 
433  call prof_rapend ("DYN_RK3",3)
434  call prof_rapstart("DYN_RK3_BND",3)
435 
436  call atmos_dyn_copy_boundary( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [INOUT]
437  prog_rk2, & ! [INOUT]
438  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
439  prog0, & ! [IN]
440  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
441 
442  call prof_rapend ("DYN_RK3_BND",3)
443 
444  call comm_vars8( dens_rk2(:,:,:), i_comm_dens_rk2 )
445  call comm_vars8( momz_rk2(:,:,:), i_comm_momz_rk2 )
446  call comm_vars8( momx_rk2(:,:,:), i_comm_momx_rk2 )
447  call comm_vars8( momy_rk2(:,:,:), i_comm_momy_rk2 )
448  call comm_vars8( rhot_rk2(:,:,:), i_comm_rhot_rk2 )
449  do iv = 1, va
450  call comm_vars8( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
451  enddo
452 
453  call comm_wait ( dens_rk2(:,:,:), i_comm_dens_rk2, .false. )
454  call comm_wait ( momz_rk2(:,:,:), i_comm_momz_rk2, .false. )
455  call comm_wait ( momx_rk2(:,:,:), i_comm_momx_rk2, .false. )
456  call comm_wait ( momy_rk2(:,:,:), i_comm_momy_rk2, .false. )
457  call comm_wait ( rhot_rk2(:,:,:), i_comm_rhot_rk2, .false. )
458  do iv = 1, va
459  call comm_wait ( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv), .false. )
460  enddo
461 
462  !##### RK3 : PROG0,PROG_RK2->PROG #####
463 
464  call prof_rapstart("DYN_RK3",3)
465 
466  dtrk = dt
467 
468  call atmos_dyn_tstep( dens, momz, momx, momy, rhot, & ! [OUT]
469  prog, & ! [OUT]
470  mflx_hi, tflx_hi, & ! [INOUT]
471  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
472  dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [IN]
473  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
474  prog0, prog_rk2, & ! [IN]
475  rtot, cvtot, corioli, & ! [IN]
476  num_diff, divdmp_coef, ddiv, & ! [IN]
477  flag_fct_momentum, flag_fct_t, & ! [IN]
478  flag_fct_along_stream, & ! [IN]
479  cdz, fdz, fdx, fdy, & ! [IN]
480  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
481  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
482  ref_pres, ref_dens, & ! [IN]
483  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
484  dtrk, dt ) ! [IN]
485 
486  if ( .NOT. flag_ws2002 ) then
487  do j = js, je
488  do i = is, ie
489  do k = ks, ke
490  dens(k,i,j) = ( dens_rk1(k,i,j) * 3.0_rp &
491  + dens(k,i,j) * 3.0_rp &
492  - dens0(k,i,j) * 2.0_rp ) / 4.0_rp
493  enddo
494  enddo
495  enddo
496 
497  do j = js, je
498  do i = is, ie
499  do k = ks, ke-1
500  momz(k,i,j) = ( momz_rk1(k,i,j) * 3.0_rp &
501  + momz(k,i,j) * 3.0_rp &
502  - momz0(k,i,j) * 2.0_rp ) / 4.0_rp
503  enddo
504  enddo
505  enddo
506 
507  do j = js, je
508  do i = is, ie
509  do k = ks, ke
510  momx(k,i,j) = ( momx_rk1(k,i,j) * 3.0_rp &
511  + momx(k,i,j) * 3.0_rp &
512  - momx0(k,i,j) * 2.0_rp ) / 4.0_rp
513  enddo
514  enddo
515  enddo
516 
517  do j = js, je
518  do i = is, ie
519  do k = ks, ke
520  momy(k,i,j) = ( momy_rk1(k,i,j) * 3.0_rp &
521  + momy(k,i,j) * 3.0_rp &
522  - momy0(k,i,j) * 2.0_rp ) / 4.0_rp
523  enddo
524  enddo
525  enddo
526 
527  do j = js, je
528  do i = is, ie
529  do k = ks, ke
530  rhot(k,i,j) = ( rhot_rk1(k,i,j) * 3.0_rp &
531  + rhot(k,i,j) * 3.0_rp &
532  - rhot0(k,i,j) * 2.0_rp ) / 4.0_rp
533  enddo
534  enddo
535  enddo
536 
537  do iv = 1, va
538  do j = js, je
539  do i = is, ie
540  do k = ks, ke
541  prog(k,i,j,iv) = ( prog_rk1(k,i,j,iv) * 3.0_rp &
542  + prog(k,i,j,iv) * 3.0_rp &
543  - prog0(k,i,j,iv) * 2.0_rp ) / 4.0_rp
544  enddo
545  enddo
546  enddo
547  enddo
548 
549  do n = 1, 3
550  do j = js, je
551  do i = is, ie
552  do k = ks, ke
553  mflx_hi(k,i,j,n) = ( mflx_hi_rk(k,i,j,n,1) &
554  + mflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
555  enddo
556  enddo
557  enddo
558  enddo
559 
560  do n = 1, 3
561  do j = js, je
562  do i = is, ie
563  do k = ks, ke
564  tflx_hi(k,i,j,n) = ( tflx_hi_rk(k,i,j,n,1) &
565  + tflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
566  enddo
567  enddo
568  enddo
569  enddo
570 
571  endif
572 
573  call prof_rapend ("DYN_RK3",3)
574 
575  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public va
Definition: scale_index.F90:38
module Atmosphere / Dynamical scheme
integer, public ke
end point of inner domain: z, local
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
procedure(short), pointer, public atmos_dyn_tstep_short
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
subroutine, public atmos_dyn_copy_boundary(DENS, MOMZ, MOMX, MOMY, RHOT, PROG, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, PROG0, BND_W, BND_E, BND_S, BND_N)
integer, public ie
end point of inner domain: x, local
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function: