SCALE-RM
scale_atmos_dyn_tinteg_short_rk3.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
14 #include "scalelib.h"
16  !-----------------------------------------------------------------------------
17  !
18  !++ used modules
19  !
20  use scale_precision
21  use scale_io
22  use scale_prof
24  use scale_index
25  use scale_tracer
26 #if defined DEBUG || defined QUICKDEBUG
27  use scale_debug, only: &
28  check
29  use scale_const, only: &
30  undef => const_undef, &
31  iundef => const_undef2
32 #endif
33  !-----------------------------------------------------------------------------
34  implicit none
35  private
36  !-----------------------------------------------------------------------------
37  !
38  !++ Public procedure
39  !
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Public parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private procedure
50  !
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private parameters & variables
54  !
55  real(RP), private, allocatable :: DENS_RK1(:,:,:) ! prognostic variables (registers for stage2)
56  real(RP), private, allocatable :: MOMZ_RK1(:,:,:)
57  real(RP), private, allocatable :: MOMX_RK1(:,:,:)
58  real(RP), private, allocatable :: MOMY_RK1(:,:,:)
59  real(RP), private, allocatable :: RHOT_RK1(:,:,:)
60  real(RP), private, allocatable :: PROG_RK1(:,:,:,:)
61  real(RP), private, allocatable :: DENS_RK2(:,:,:) ! prognostic variables (registers for stage3)
62  real(RP), private, allocatable :: MOMZ_RK2(:,:,:)
63  real(RP), private, allocatable :: MOMX_RK2(:,:,:)
64  real(RP), private, allocatable :: MOMY_RK2(:,:,:)
65  real(RP), private, allocatable :: RHOT_RK2(:,:,:)
66  real(RP), private, allocatable :: PROG_RK2(:,:,:,:)
67 
68  ! for communication
69  integer :: I_COMM_DENS_RK1 = 1
70  integer :: I_COMM_MOMZ_RK1 = 2
71  integer :: I_COMM_MOMX_RK1 = 3
72  integer :: I_COMM_MOMY_RK1 = 4
73  integer :: I_COMM_RHOT_RK1 = 5
74  integer, allocatable :: I_COMM_PROG_RK1(:)
75 
76  integer :: I_COMM_DENS_RK2 = 1
77  integer :: I_COMM_MOMZ_RK2 = 2
78  integer :: I_COMM_MOMX_RK2 = 3
79  integer :: I_COMM_MOMY_RK2 = 4
80  integer :: I_COMM_RHOT_RK2 = 5
81  integer, allocatable :: I_COMM_PROG_RK2(:)
82 
83  logical :: FLAG_WS2002
84  real(RP) :: fact_dt1
85  real(RP) :: fact_dt2
86 
87  !-----------------------------------------------------------------------------
88 contains
89  !-----------------------------------------------------------------------------
92  tinteg_type )
93  use scale_prc, only: &
94  prc_abort
95  use scale_const, only: &
96  undef => const_undef
97  use scale_comm_cartesc, only: &
99  implicit none
100 
101  character(len=*) :: tinteg_type
102 
103  integer :: iv
104  !---------------------------------------------------------------------------
105 
106  select case( tinteg_type )
107  case( 'RK3' )
108  log_info("ATMOS_DYN_Tinteg_short_rk3_setup",*) "RK3: Heun's method is used"
109  ! Heun's method
110  ! k1 = f(\phi_n); r1 = \phi_n + k1 * dt / 3
111  ! k2 = f(r1); r2 = \phi_n + k2 * dt * 2 / 3
112  ! k3 = f(r2); r3 = \phi_n + k3 * dt
113  ! \phi_{n+1} = \phi_n + ( k1 + 3 * k3 ) dt / 4
114  ! = \phi_n + ( (r1-\phi_n) * 3 + (r3-\phi_n) * 3 ) / 4
115  ! = ( r1 * 3 + r3 * 3 - \phi_n * 2 ) / 4
116  flag_ws2002 = .false.
117  fact_dt1 = 1.0_rp / 3.0_rp
118  fact_dt2 = 2.0_rp / 3.0_rp
119  case( 'RK3WS2002' )
120  log_info("ATMOS_DYN_Tinteg_short_rk3_setup",*) "RK3: Wicker and Skamarock (2002) is used"
121  ! Wicker and Skamarock (2002) RK3 scheme
122  ! k1 = f(\phi_n); r1 = \phi_n + k1 * dt / 3
123  ! k2 = f(r1); r2 = \phi_n + k2 * dt / 2
124  ! k3 = f(r2); r3 = \phi_n + k3 * dt
125  ! \phi_{n+1} = r3
126  ! Unlike to Heun's RK3 method, the memory arrays are not needed in this case.
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  log_error("ATMOS_DYN_Tinteg_short_rk3_setup",*) 'TINTEG_TYPE is not RK3. Check!'
132  call prc_abort
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', dens_rk1, i_comm_dens_rk1 )
153  call comm_vars8_init( 'MOMZ_RK1', momz_rk1, i_comm_momz_rk1 )
154  call comm_vars8_init( 'MOMX_RK1', momx_rk1, i_comm_momx_rk1 )
155  call comm_vars8_init( 'MOMY_RK1', momy_rk1, i_comm_momy_rk1 )
156  call comm_vars8_init( 'RHOT_RK1', 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', prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
160  enddo
161 
162  call comm_vars8_init( 'DENS_RK2', dens_rk2, i_comm_dens_rk2 )
163  call comm_vars8_init( 'MOMZ_RK2', momz_rk2, i_comm_momz_rk2 )
164  call comm_vars8_init( 'MOMX_RK2', momx_rk2, i_comm_momx_rk2 )
165  call comm_vars8_init( 'MOMY_RK2', momy_rk2, i_comm_momy_rk2 )
166  call comm_vars8_init( 'RHOT_RK2', 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', 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
187  end subroutine atmos_dyn_tinteg_short_rk3_setup
188 
189  !-----------------------------------------------------------------------------
191  subroutine atmos_dyn_tinteg_short_rk3( &
192  DENS, MOMZ, MOMX, MOMY, RHOT, PROG, &
193  mflx_hi, tflx_hi, &
194  DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, &
195  DPRES0, CVtot, CORIOLI, &
196  num_diff, wdamp_coef, divdmp_coef, DDIV, &
197  FLAG_FCT_MOMENTUM, FLAG_FCT_T, &
198  FLAG_FCT_ALONG_STREAM, &
199  CDZ, FDZ, FDX, FDY, &
200  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
201  PHI, GSQRT, J13G, J23G, J33G, MAPF, &
202  REF_pres, REF_dens, &
203  BND_W, BND_E, BND_S, BND_N, TwoD, &
204  dt )
205  use scale_comm_cartesc, 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(out ) :: 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) :: dpres0(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) :: wdamp_coef(ka)
235  real(rp), intent(in) :: divdmp_coef
236  real(rp), intent(in) :: ddiv(ka,ia,ja)
237 
238  logical, intent(in) :: flag_fct_momentum
239  logical, intent(in) :: flag_fct_t
240  logical, intent(in) :: flag_fct_along_stream
241 
242  real(rp), intent(in) :: cdz (ka)
243  real(rp), intent(in) :: fdz (ka-1)
244  real(rp), intent(in) :: fdx (ia-1)
245  real(rp), intent(in) :: fdy (ja-1)
246  real(rp), intent(in) :: rcdz(ka)
247  real(rp), intent(in) :: rcdx(ia)
248  real(rp), intent(in) :: rcdy(ja)
249  real(rp), intent(in) :: rfdz(ka-1)
250  real(rp), intent(in) :: rfdx(ia-1)
251  real(rp), intent(in) :: rfdy(ja-1)
252 
253  real(rp), intent(in) :: phi (ka,ia,ja)
254  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
255  real(rp), intent(in) :: j13g (ka,ia,ja,7)
256  real(rp), intent(in) :: j23g (ka,ia,ja,7)
257  real(rp), intent(in) :: j33g
258  real(rp), intent(in) :: mapf (ia,ja,2,4)
259 
260  real(rp), intent(in) :: ref_pres(ka,ia,ja)
261  real(rp), intent(in) :: ref_dens(ka,ia,ja)
262 
263  logical, intent(in) :: bnd_w
264  logical, intent(in) :: bnd_e
265  logical, intent(in) :: bnd_s
266  logical, intent(in) :: bnd_n
267  logical, intent(in) :: twod
268 
269  real(rp), intent(in) :: dt
270 
271  real(rp) :: dens0(ka,ia,ja)
272  real(rp) :: momz0(ka,ia,ja)
273  real(rp) :: momx0(ka,ia,ja)
274  real(rp) :: momy0(ka,ia,ja)
275  real(rp) :: rhot0(ka,ia,ja)
276  real(rp) :: prog0(ka,ia,ja,va)
277 
278  real(rp) :: mflx_hi_rk(ka,ia,ja,3,2)
279  real(rp) :: tflx_hi_rk(ka,ia,ja,3,2)
280 
281  real(rp) :: dtrk
282 
283  integer :: i, j, k, iv, n
284  !---------------------------------------------------------------------------
285 
286  call prof_rapstart("DYN_RK3_Prep",3)
287 
288 #ifdef DEBUG
289  dens_rk1(:,:,:) = undef
290  momz_rk1(:,:,:) = undef
291  momx_rk1(:,:,:) = undef
292  momy_rk1(:,:,:) = undef
293  rhot_rk1(:,:,:) = undef
294  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
295 
296  dens_rk2(:,:,:) = undef
297  momz_rk2(:,:,:) = undef
298  momx_rk2(:,:,:) = undef
299  momy_rk2(:,:,:) = undef
300  rhot_rk2(:,:,:) = undef
301  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
302 
303  mflx_hi_rk(:,:,:,:,:) = undef
304  tflx_hi_rk(:,:,:,:,:) = undef
305 #endif
306 
307 #ifdef QUICKDEBUG
308  mflx_hi( 1:ks-1,:,:,:) = undef
309  mflx_hi(ke+1:ka ,:,:,:) = undef
310 #endif
311 
312 !OCL XFILL
313  dens0 = dens
314 !OCL XFILL
315  momz0 = momz
316 !OCL XFILL
317  momx0 = momx
318 !OCL XFILL
319  momy0 = momy
320 !OCL XFILL
321  rhot0 = rhot
322 !OCL XFILL
323  if ( va > 0 ) prog0 = prog
324 
325  if ( bnd_w .and. (.not. twod) ) then
326  do j = js, je
327  do k = ks, ke
328  mflx_hi_rk(k,is-1,j,2,:) = mflx_hi(k,is-1,j,2)
329  end do
330  end do
331  end if
332  if ( bnd_e .and. (.not. twod) ) then
333  do j = js, je
334  do k = ks, ke
335  mflx_hi_rk(k,ie,j,2,:) = mflx_hi(k,ie,j,2)
336  end do
337  end do
338  end if
339  if ( bnd_s ) then
340  do i = is, ie
341  do k = ks, ke
342  mflx_hi_rk(k,i,js-1,3,:) = mflx_hi(k,i,js-1,3)
343  end do
344  end do
345  end if
346  if ( bnd_n ) then
347  do i = is, ie
348  do k = ks, ke
349  mflx_hi_rk(k,i,je,3,:) = mflx_hi(k,i,je,3)
350  end do
351  end do
352  end if
353 
354  call prof_rapend ("DYN_RK3_Prep",3)
355 
356  !------------------------------------------------------------------------
357  ! Start RK
358  !------------------------------------------------------------------------
359 
360  !##### RK1 : PROG0,PROG->PROG_RK1 #####
361 
362  call prof_rapstart("DYN_RK3",3)
363 
364  dtrk = dt * fact_dt1
365 
366  call atmos_dyn_tstep( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [OUT]
367  prog_rk1, & ! [OUT]
368  mflx_hi_rk(:,:,:,:,1), tflx_hi_rk(:,:,:,:,1), & ! [INOUT,OUT]
369  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
370  dens, momz, momx, momy, rhot, & ! [IN]
371  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
372  prog0, prog, & ! [IN]
373  dpres0, cvtot, corioli, & ! [IN]
374  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
375  flag_fct_momentum, flag_fct_t, & ! [IN]
376  flag_fct_along_stream, & ! [IN]
377  cdz, fdz, fdx, fdy, & ! [IN]
378  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
379  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
380  ref_pres, ref_dens, & ! [IN]
381  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
382  dtrk, .false. ) ! [IN]
383 
384  call prof_rapend ("DYN_RK3",3)
385  call prof_rapstart("DYN_RK3_BND",3)
386 
387  call atmos_dyn_copy_boundary( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [INOUT]
388  prog_rk1, & ! [INOUT]
389  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
390  prog0, & ! [IN]
391  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
392 
393  call prof_rapend ("DYN_RK3_BND",3)
394 
395  call comm_vars8( dens_rk1(:,:,:), i_comm_dens_rk1 )
396  call comm_vars8( momz_rk1(:,:,:), i_comm_momz_rk1 )
397  call comm_vars8( momx_rk1(:,:,:), i_comm_momx_rk1 )
398  call comm_vars8( momy_rk1(:,:,:), i_comm_momy_rk1 )
399  call comm_vars8( rhot_rk1(:,:,:), i_comm_rhot_rk1 )
400  do iv = 1, va
401  call comm_vars8( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
402  enddo
403 
404  call comm_wait ( dens_rk1(:,:,:), i_comm_dens_rk1, .false. )
405  call comm_wait ( momz_rk1(:,:,:), i_comm_momz_rk1, .false. )
406  call comm_wait ( momx_rk1(:,:,:), i_comm_momx_rk1, .false. )
407  call comm_wait ( momy_rk1(:,:,:), i_comm_momy_rk1, .false. )
408  call comm_wait ( rhot_rk1(:,:,:), i_comm_rhot_rk1, .false. )
409  do iv = 1, va
410  call comm_wait ( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv), .false. )
411  enddo
412 
413  !##### RK2 : PROG0,PROG_RK1->PROG_RK2 #####
414 
415  call prof_rapstart("DYN_RK3",3)
416 
417  dtrk = dt * fact_dt2
418 
419  call atmos_dyn_tstep( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [OUT]
420  prog_rk2, & ! [OUT]
421  mflx_hi_rk(:,:,:,:,2), tflx_hi_rk(:,:,:,:,2), & ! [INOUT,OUT]
422  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
423  dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [IN]
424  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
425  prog0, prog_rk1, & ! [IN]
426  dpres0, cvtot, corioli, & ! [IN]
427  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
428  flag_fct_momentum, flag_fct_t, & ! [IN]
429  flag_fct_along_stream, & ! [IN]
430  cdz, fdz, fdx, fdy, & ! [IN]
431  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
432  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
433  ref_pres, ref_dens, & ! [IN]
434  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
435  dtrk, .false. ) ! [IN]
436 
437  call prof_rapend ("DYN_RK3",3)
438  call prof_rapstart("DYN_RK3_BND",3)
439 
440  call atmos_dyn_copy_boundary( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [INOUT]
441  prog_rk2, & ! [INOUT]
442  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
443  prog0, & ! [IN]
444  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
445 
446  call prof_rapend ("DYN_RK3_BND",3)
447 
448  call comm_vars8( dens_rk2(:,:,:), i_comm_dens_rk2 )
449  call comm_vars8( momz_rk2(:,:,:), i_comm_momz_rk2 )
450  call comm_vars8( momx_rk2(:,:,:), i_comm_momx_rk2 )
451  call comm_vars8( momy_rk2(:,:,:), i_comm_momy_rk2 )
452  call comm_vars8( rhot_rk2(:,:,:), i_comm_rhot_rk2 )
453  do iv = 1, va
454  call comm_vars8( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
455  enddo
456 
457  call comm_wait ( dens_rk2(:,:,:), i_comm_dens_rk2, .false. )
458  call comm_wait ( momz_rk2(:,:,:), i_comm_momz_rk2, .false. )
459  call comm_wait ( momx_rk2(:,:,:), i_comm_momx_rk2, .false. )
460  call comm_wait ( momy_rk2(:,:,:), i_comm_momy_rk2, .false. )
461  call comm_wait ( rhot_rk2(:,:,:), i_comm_rhot_rk2, .false. )
462  do iv = 1, va
463  call comm_wait ( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv), .false. )
464  enddo
465 
466  !##### RK3 : PROG0,PROG_RK2->PROG #####
467 
468  call prof_rapstart("DYN_RK3",3)
469 
470  dtrk = dt
471 
472  call atmos_dyn_tstep( dens, momz, momx, momy, rhot, & ! [OUT]
473  prog, & ! [OUT]
474  mflx_hi, tflx_hi, & ! [INOUT,OUT]
475  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
476  dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [IN]
477  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
478  prog0, prog_rk2, & ! [IN]
479  dpres0, cvtot, corioli, & ! [IN]
480  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
481  flag_fct_momentum, flag_fct_t, & ! [IN]
482  flag_fct_along_stream, & ! [IN]
483  cdz, fdz, fdx, fdy, & ! [IN]
484  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
485  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
486  ref_pres, ref_dens, & ! [IN]
487  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
488  dtrk, .true. ) ! [IN]
489 
490  if ( .NOT. flag_ws2002 ) then
491  do j = js, je
492  do i = is, ie
493  do k = ks, ke
494  dens(k,i,j) = ( dens_rk1(k,i,j) * 3.0_rp &
495  + dens(k,i,j) * 3.0_rp &
496  - dens0(k,i,j) * 2.0_rp ) / 4.0_rp
497  enddo
498  enddo
499  enddo
500 
501  do j = js, je
502  do i = is, ie
503  do k = ks, ke-1
504  momz(k,i,j) = ( momz_rk1(k,i,j) * 3.0_rp &
505  + momz(k,i,j) * 3.0_rp &
506  - momz0(k,i,j) * 2.0_rp ) / 4.0_rp
507  enddo
508  enddo
509  enddo
510 
511  do j = js, je
512  do i = is, ie
513  do k = ks, ke
514  momx(k,i,j) = ( momx_rk1(k,i,j) * 3.0_rp &
515  + momx(k,i,j) * 3.0_rp &
516  - momx0(k,i,j) * 2.0_rp ) / 4.0_rp
517  enddo
518  enddo
519  enddo
520 
521  do j = js, je
522  do i = is, ie
523  do k = ks, ke
524  momy(k,i,j) = ( momy_rk1(k,i,j) * 3.0_rp &
525  + momy(k,i,j) * 3.0_rp &
526  - momy0(k,i,j) * 2.0_rp ) / 4.0_rp
527  enddo
528  enddo
529  enddo
530 
531  do j = js, je
532  do i = is, ie
533  do k = ks, ke
534  rhot(k,i,j) = ( rhot_rk1(k,i,j) * 3.0_rp &
535  + rhot(k,i,j) * 3.0_rp &
536  - rhot0(k,i,j) * 2.0_rp ) / 4.0_rp
537  enddo
538  enddo
539  enddo
540 
541  do iv = 1, va
542  do j = js, je
543  do i = is, ie
544  do k = ks, ke
545  prog(k,i,j,iv) = ( prog_rk1(k,i,j,iv) * 3.0_rp &
546  + prog(k,i,j,iv) * 3.0_rp &
547  - prog0(k,i,j,iv) * 2.0_rp ) / 4.0_rp
548  enddo
549  enddo
550  enddo
551  enddo
552 
553  do n = 1, 3
554  do j = js, je
555  do i = is, ie
556  do k = ks, ke
557  mflx_hi(k,i,j,n) = ( mflx_hi_rk(k,i,j,n,1) &
558  + mflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
559  enddo
560  enddo
561  enddo
562  enddo
563 
564  do n = 1, 3
565  do j = js, je
566  do i = is, ie
567  do k = ks, ke
568  tflx_hi(k,i,j,n) = ( tflx_hi_rk(k,i,j,n,1) &
569  + tflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
570  enddo
571  enddo
572  enddo
573  enddo
574 
575  endif
576 
577  call prof_rapend ("DYN_RK3",3)
578 
579  return
580  end subroutine atmos_dyn_tinteg_short_rk3
581 
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_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_index
module Index
Definition: scale_index.F90:11
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
scale_atmos_dyn_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:159
scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short
procedure(short), pointer, public atmos_dyn_tstep_short
Definition: scale_atmos_dyn_tstep_short.F90:135
scale_index::va
integer, public va
Definition: scale_index.F90:35
scale_prc
module PROCESS
Definition: scale_prc.F90:11
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_io
module STDIO
Definition: scale_io.F90:10
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:44
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3
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, DPRES0, CVtot, CORIOLI, num_diff, wdamp_coef, 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, TwoD, dt)
RK3.
Definition: scale_atmos_dyn_tinteg_short_rk3.F90:205
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
scale_prof
module profiler
Definition: scale_prof.F90:11
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_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_tstep_short
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_short.F90:12
scale_atmos_dyn_tinteg_short_rk3
module Atmosphere / Dyn Tinteg
Definition: scale_atmos_dyn_tinteg_short_rk3.F90:15
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
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_short_rk3::atmos_dyn_tinteg_short_rk3_setup
subroutine, public atmos_dyn_tinteg_short_rk3_setup(tinteg_type)
Setup.
Definition: scale_atmos_dyn_tinteg_short_rk3.F90:93
scale_debug
module DEBUG
Definition: scale_debug.F90:11
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_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:217
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_atmos_dyn_common::atmos_dyn_copy_boundary
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, TwoD)
Definition: scale_atmos_dyn_common.F90:188
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56