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  !
43 
44  !-----------------------------------------------------------------------------
45  !
46  !++ Public parameters & variables
47  !
48  !-----------------------------------------------------------------------------
49  !
50  !++ Private procedure
51  !
52  !-----------------------------------------------------------------------------
53  !
54  !++ Private parameters & variables
55  !
56  real(RP), private, allocatable :: DENS_RK1(:,:,:) ! prognostic variables (registers for stage2)
57  real(RP), private, allocatable :: MOMZ_RK1(:,:,:)
58  real(RP), private, allocatable :: MOMX_RK1(:,:,:)
59  real(RP), private, allocatable :: MOMY_RK1(:,:,:)
60  real(RP), private, allocatable :: RHOT_RK1(:,:,:)
61  real(RP), private, allocatable :: PROG_RK1(:,:,:,:)
62  real(RP), private, allocatable :: DENS_RK2(:,:,:) ! prognostic variables (registers for stage3)
63  real(RP), private, allocatable :: MOMZ_RK2(:,:,:)
64  real(RP), private, allocatable :: MOMX_RK2(:,:,:)
65  real(RP), private, allocatable :: MOMY_RK2(:,:,:)
66  real(RP), private, allocatable :: RHOT_RK2(:,:,:)
67  real(RP), private, allocatable :: PROG_RK2(:,:,:,:)
68 
69  ! for communication
70  integer :: I_COMM_DENS_RK1
71  integer :: I_COMM_MOMZ_RK1
72  integer :: I_COMM_MOMX_RK1
73  integer :: I_COMM_MOMY_RK1
74  integer :: I_COMM_RHOT_RK1
75  integer, allocatable :: I_COMM_PROG_RK1(:)
76 
77  integer :: I_COMM_DENS_RK2
78  integer :: I_COMM_MOMZ_RK2
79  integer :: I_COMM_MOMX_RK2
80  integer :: I_COMM_MOMY_RK2
81  integer :: I_COMM_RHOT_RK2
82  integer, allocatable :: I_COMM_PROG_RK2(:)
83 
84  logical :: FLAG_WS2002
85  real(RP) :: fact_dt1
86  real(RP) :: fact_dt2
87 
88  !-----------------------------------------------------------------------------
89 contains
90  !-----------------------------------------------------------------------------
93  tinteg_type )
94  use scale_prc, only: &
95  prc_abort
96  use scale_const, only: &
97  undef => const_undef
98  use scale_comm_cartesc, 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  log_info("ATMOS_DYN_Tinteg_short_rk3_setup",*) "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  log_info("ATMOS_DYN_Tinteg_short_rk3_setup",*) "RK3: Wicker and Skamarock (2002) is used"
122  ! Wicker 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  ! Unlike to Heun's RK3 method, the memory arrays are not needed in this case.
128  flag_ws2002 = .true.
129  fact_dt1 = 1.0_rp / 3.0_rp
130  fact_dt2 = 1.0_rp / 2.0_rp
131  case default
132  log_error("ATMOS_DYN_Tinteg_short_rk3_setup",*) 'TINTEG_TYPE is not RK3. Check!'
133  call prc_abort
134  end select
135 
136  allocate( dens_rk1(ka,ia,ja) )
137  allocate( momz_rk1(ka,ia,ja) )
138  allocate( momx_rk1(ka,ia,ja) )
139  allocate( momy_rk1(ka,ia,ja) )
140  allocate( rhot_rk1(ka,ia,ja) )
141 
142  allocate( dens_rk2(ka,ia,ja) )
143  allocate( momz_rk2(ka,ia,ja) )
144  allocate( momx_rk2(ka,ia,ja) )
145  allocate( momy_rk2(ka,ia,ja) )
146  allocate( rhot_rk2(ka,ia,ja) )
147 
148  allocate( prog_rk1(ka,ia,ja,max(va,1)) )
149  allocate( prog_rk2(ka,ia,ja,max(va,1)) )
150  allocate( i_comm_prog_rk1(max(va,1)) )
151  allocate( i_comm_prog_rk2(max(va,1)) )
152 
153  i_comm_dens_rk1 = 1
154  i_comm_momz_rk1 = 2
155  i_comm_momx_rk1 = 3
156  i_comm_momy_rk1 = 4
157  i_comm_rhot_rk1 = 5
158  call comm_vars8_init( 'DENS_RK1', dens_rk1, i_comm_dens_rk1 )
159  call comm_vars8_init( 'MOMZ_RK1', momz_rk1, i_comm_momz_rk1 )
160  call comm_vars8_init( 'MOMX_RK1', momx_rk1, i_comm_momx_rk1 )
161  call comm_vars8_init( 'MOMY_RK1', momy_rk1, i_comm_momy_rk1 )
162  call comm_vars8_init( 'RHOT_RK1', rhot_rk1, i_comm_rhot_rk1 )
163  do iv = 1, va
164  i_comm_prog_rk1(iv) = 5 + iv
165  call comm_vars8_init( 'PROG_RK1', prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
166  enddo
167 
168 
169  i_comm_dens_rk2 = 1
170  i_comm_momz_rk2 = 2
171  i_comm_momx_rk2 = 3
172  i_comm_momy_rk2 = 4
173  i_comm_rhot_rk2 = 5
174  call comm_vars8_init( 'DENS_RK2', dens_rk2, i_comm_dens_rk2 )
175  call comm_vars8_init( 'MOMZ_RK2', momz_rk2, i_comm_momz_rk2 )
176  call comm_vars8_init( 'MOMX_RK2', momx_rk2, i_comm_momx_rk2 )
177  call comm_vars8_init( 'MOMY_RK2', momy_rk2, i_comm_momy_rk2 )
178  call comm_vars8_init( 'RHOT_RK2', rhot_rk2, i_comm_rhot_rk2 )
179  do iv = 1, va
180  i_comm_prog_rk2(iv) = 5 + iv
181  call comm_vars8_init( 'PROG_RK2', prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
182  enddo
183 
184  dens_rk1(:,:,:) = undef
185  momz_rk1(:,:,:) = undef
186  momx_rk1(:,:,:) = undef
187  momy_rk1(:,:,:) = undef
188  rhot_rk1(:,:,:) = undef
189  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
190 
191  dens_rk2(:,:,:) = undef
192  momz_rk2(:,:,:) = undef
193  momx_rk2(:,:,:) = undef
194  momy_rk2(:,:,:) = undef
195  rhot_rk2(:,:,:) = undef
196  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
197 
198  return
199  end subroutine atmos_dyn_tinteg_short_rk3_setup
200 
201  !-----------------------------------------------------------------------------
204 
205  deallocate( dens_rk1 )
206  deallocate( momz_rk1 )
207  deallocate( momx_rk1 )
208  deallocate( momy_rk1 )
209  deallocate( rhot_rk1 )
210 
211  deallocate( dens_rk2 )
212  deallocate( momz_rk2 )
213  deallocate( momx_rk2 )
214  deallocate( momy_rk2 )
215  deallocate( rhot_rk2 )
216 
217  deallocate( prog_rk1 )
218  deallocate( prog_rk2 )
219  deallocate( i_comm_prog_rk1 )
220  deallocate( i_comm_prog_rk2 )
221 
222  return
224  !-----------------------------------------------------------------------------
226  subroutine atmos_dyn_tinteg_short_rk3( &
227  DENS, MOMZ, MOMX, MOMY, RHOT, PROG, &
228  mflx_hi, tflx_hi, &
229  DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, &
230  DPRES0, CVtot, CORIOLI, &
231  num_diff, wdamp_coef, divdmp_coef, DDIV, &
232  FLAG_FCT_MOMENTUM, FLAG_FCT_T, &
233  FLAG_FCT_ALONG_STREAM, &
234  CDZ, FDZ, FDX, FDY, &
235  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
236  PHI, GSQRT, J13G, J23G, J33G, MAPF, &
237  REF_pres, REF_dens, &
238  BND_W, BND_E, BND_S, BND_N, TwoD, &
239  dt )
240  use scale_comm_cartesc, only: &
241  comm_vars8, &
242  comm_wait
243  use scale_atmos_dyn_tstep_short, only: &
244  atmos_dyn_tstep => atmos_dyn_tstep_short
245  use scale_atmos_dyn_common, only: &
247  implicit none
248 
249  real(rp), intent(inout) :: dens(ka,ia,ja)
250  real(rp), intent(inout) :: momz(ka,ia,ja)
251  real(rp), intent(inout) :: momx(ka,ia,ja)
252  real(rp), intent(inout) :: momy(ka,ia,ja)
253  real(rp), intent(inout) :: rhot(ka,ia,ja)
254  real(rp), intent(inout) :: prog(ka,ia,ja,va)
255 
256  real(rp), intent(inout) :: mflx_hi(ka,ia,ja,3)
257  real(rp), intent(out ) :: tflx_hi(ka,ia,ja,3)
258 
259  real(rp), intent(in) :: dens_t(ka,ia,ja)
260  real(rp), intent(in) :: momz_t(ka,ia,ja)
261  real(rp), intent(in) :: momx_t(ka,ia,ja)
262  real(rp), intent(in) :: momy_t(ka,ia,ja)
263  real(rp), intent(in) :: rhot_t(ka,ia,ja)
264 
265  real(rp), intent(in) :: dpres0(ka,ia,ja)
266  real(rp), intent(in) :: cvtot(ka,ia,ja)
267  real(rp), intent(in) :: corioli(ia,ja)
268  real(rp), intent(in) :: num_diff(ka,ia,ja,5,3)
269  real(rp), intent(in) :: wdamp_coef(ka)
270  real(rp), intent(in) :: divdmp_coef
271  real(rp), intent(in) :: ddiv(ka,ia,ja)
272 
273  logical, intent(in) :: flag_fct_momentum
274  logical, intent(in) :: flag_fct_t
275  logical, intent(in) :: flag_fct_along_stream
276 
277  real(rp), intent(in) :: cdz (ka)
278  real(rp), intent(in) :: fdz (ka-1)
279  real(rp), intent(in) :: fdx (ia-1)
280  real(rp), intent(in) :: fdy (ja-1)
281  real(rp), intent(in) :: rcdz(ka)
282  real(rp), intent(in) :: rcdx(ia)
283  real(rp), intent(in) :: rcdy(ja)
284  real(rp), intent(in) :: rfdz(ka-1)
285  real(rp), intent(in) :: rfdx(ia-1)
286  real(rp), intent(in) :: rfdy(ja-1)
287 
288  real(rp), intent(in) :: phi (ka,ia,ja)
289  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
290  real(rp), intent(in) :: j13g (ka,ia,ja,7)
291  real(rp), intent(in) :: j23g (ka,ia,ja,7)
292  real(rp), intent(in) :: j33g
293  real(rp), intent(in) :: mapf (ia,ja,2,4)
294 
295  real(rp), intent(in) :: ref_pres(ka,ia,ja)
296  real(rp), intent(in) :: ref_dens(ka,ia,ja)
297 
298  logical, intent(in) :: bnd_w
299  logical, intent(in) :: bnd_e
300  logical, intent(in) :: bnd_s
301  logical, intent(in) :: bnd_n
302  logical, intent(in) :: twod
303 
304  real(rp), intent(in) :: dt
305 
306  real(rp) :: dens0(ka,ia,ja)
307  real(rp) :: momz0(ka,ia,ja)
308  real(rp) :: momx0(ka,ia,ja)
309  real(rp) :: momy0(ka,ia,ja)
310  real(rp) :: rhot0(ka,ia,ja)
311  real(rp) :: prog0(ka,ia,ja,va)
312 
313  real(rp) :: mflx_hi_rk(ka,ia,ja,3,2)
314  real(rp) :: tflx_hi_rk(ka,ia,ja,3,2)
315 
316  real(rp) :: dtrk
317 
318  integer :: i, j, k, iv, n
319  !---------------------------------------------------------------------------
320 
321  call prof_rapstart("DYN_RK3_Prep",3)
322 
323 #ifdef DEBUG
324  dens_rk1(:,:,:) = undef
325  momz_rk1(:,:,:) = undef
326  momx_rk1(:,:,:) = undef
327  momy_rk1(:,:,:) = undef
328  rhot_rk1(:,:,:) = undef
329  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
330 
331  dens_rk2(:,:,:) = undef
332  momz_rk2(:,:,:) = undef
333  momx_rk2(:,:,:) = undef
334  momy_rk2(:,:,:) = undef
335  rhot_rk2(:,:,:) = undef
336  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
337 
338  mflx_hi_rk(:,:,:,:,:) = undef
339  tflx_hi_rk(:,:,:,:,:) = undef
340 #endif
341 
342 #ifdef QUICKDEBUG
343  mflx_hi( 1:ks-1,:,:,:) = undef
344  mflx_hi(ke+1:ka ,:,:,:) = undef
345 #endif
346 
347 !OCL XFILL
348  dens0 = dens
349 !OCL XFILL
350  momz0 = momz
351 !OCL XFILL
352  momx0 = momx
353 !OCL XFILL
354  momy0 = momy
355 !OCL XFILL
356  rhot0 = rhot
357 !OCL XFILL
358  if ( va > 0 ) prog0 = prog
359 
360  if ( bnd_w .and. (.not. twod) ) then
361  do j = js, je
362  do k = ks, ke
363  mflx_hi_rk(k,is-1,j,2,:) = mflx_hi(k,is-1,j,2)
364  end do
365  end do
366  end if
367  if ( bnd_e .and. (.not. twod) ) then
368  do j = js, je
369  do k = ks, ke
370  mflx_hi_rk(k,ie,j,2,:) = mflx_hi(k,ie,j,2)
371  end do
372  end do
373  end if
374  if ( bnd_s ) then
375  do i = is, ie
376  do k = ks, ke
377  mflx_hi_rk(k,i,js-1,3,:) = mflx_hi(k,i,js-1,3)
378  end do
379  end do
380  end if
381  if ( bnd_n ) then
382  do i = is, ie
383  do k = ks, ke
384  mflx_hi_rk(k,i,je,3,:) = mflx_hi(k,i,je,3)
385  end do
386  end do
387  end if
388 
389  call prof_rapend ("DYN_RK3_Prep",3)
390 
391  !------------------------------------------------------------------------
392  ! Start RK
393  !------------------------------------------------------------------------
394 
395  !##### RK1 : PROG0,PROG->PROG_RK1 #####
396 
397  call prof_rapstart("DYN_RK3",3)
398 
399  dtrk = dt * fact_dt1
400 
401  call atmos_dyn_tstep( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [OUT]
402  prog_rk1, & ! [OUT]
403  mflx_hi_rk(:,:,:,:,1), tflx_hi_rk(:,:,:,:,1), & ! [INOUT,OUT]
404  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
405  dens, momz, momx, momy, rhot, & ! [IN]
406  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
407  prog0, prog, & ! [IN]
408  dpres0, cvtot, corioli, & ! [IN]
409  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
410  flag_fct_momentum, flag_fct_t, & ! [IN]
411  flag_fct_along_stream, & ! [IN]
412  cdz, fdz, fdx, fdy, & ! [IN]
413  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
414  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
415  ref_pres, ref_dens, & ! [IN]
416  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
417  dtrk, .false. ) ! [IN]
418 
419  call prof_rapend ("DYN_RK3",3)
420  call prof_rapstart("DYN_RK3_BND",3)
421 
422  call atmos_dyn_copy_boundary( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [INOUT]
423  prog_rk1, & ! [INOUT]
424  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
425  prog0, & ! [IN]
426  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
427 
428  call prof_rapend ("DYN_RK3_BND",3)
429 
430  call comm_vars8( dens_rk1(:,:,:), i_comm_dens_rk1 )
431  call comm_vars8( momz_rk1(:,:,:), i_comm_momz_rk1 )
432  call comm_vars8( momx_rk1(:,:,:), i_comm_momx_rk1 )
433  call comm_vars8( momy_rk1(:,:,:), i_comm_momy_rk1 )
434  call comm_vars8( rhot_rk1(:,:,:), i_comm_rhot_rk1 )
435  do iv = 1, va
436  call comm_vars8( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
437  enddo
438 
439  call comm_wait ( dens_rk1(:,:,:), i_comm_dens_rk1, .false. )
440  call comm_wait ( momz_rk1(:,:,:), i_comm_momz_rk1, .false. )
441  call comm_wait ( momx_rk1(:,:,:), i_comm_momx_rk1, .false. )
442  call comm_wait ( momy_rk1(:,:,:), i_comm_momy_rk1, .false. )
443  call comm_wait ( rhot_rk1(:,:,:), i_comm_rhot_rk1, .false. )
444  do iv = 1, va
445  call comm_wait ( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv), .false. )
446  enddo
447 
448  !##### RK2 : PROG0,PROG_RK1->PROG_RK2 #####
449 
450  call prof_rapstart("DYN_RK3",3)
451 
452  dtrk = dt * fact_dt2
453 
454  call atmos_dyn_tstep( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [OUT]
455  prog_rk2, & ! [OUT]
456  mflx_hi_rk(:,:,:,:,2), tflx_hi_rk(:,:,:,:,2), & ! [INOUT,OUT]
457  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
458  dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [IN]
459  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
460  prog0, prog_rk1, & ! [IN]
461  dpres0, cvtot, corioli, & ! [IN]
462  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
463  flag_fct_momentum, flag_fct_t, & ! [IN]
464  flag_fct_along_stream, & ! [IN]
465  cdz, fdz, fdx, fdy, & ! [IN]
466  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
467  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
468  ref_pres, ref_dens, & ! [IN]
469  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
470  dtrk, .false. ) ! [IN]
471 
472  call prof_rapend ("DYN_RK3",3)
473  call prof_rapstart("DYN_RK3_BND",3)
474 
475  call atmos_dyn_copy_boundary( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [INOUT]
476  prog_rk2, & ! [INOUT]
477  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
478  prog0, & ! [IN]
479  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
480 
481  call prof_rapend ("DYN_RK3_BND",3)
482 
483  call comm_vars8( dens_rk2(:,:,:), i_comm_dens_rk2 )
484  call comm_vars8( momz_rk2(:,:,:), i_comm_momz_rk2 )
485  call comm_vars8( momx_rk2(:,:,:), i_comm_momx_rk2 )
486  call comm_vars8( momy_rk2(:,:,:), i_comm_momy_rk2 )
487  call comm_vars8( rhot_rk2(:,:,:), i_comm_rhot_rk2 )
488  do iv = 1, va
489  call comm_vars8( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
490  enddo
491 
492  call comm_wait ( dens_rk2(:,:,:), i_comm_dens_rk2, .false. )
493  call comm_wait ( momz_rk2(:,:,:), i_comm_momz_rk2, .false. )
494  call comm_wait ( momx_rk2(:,:,:), i_comm_momx_rk2, .false. )
495  call comm_wait ( momy_rk2(:,:,:), i_comm_momy_rk2, .false. )
496  call comm_wait ( rhot_rk2(:,:,:), i_comm_rhot_rk2, .false. )
497  do iv = 1, va
498  call comm_wait ( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv), .false. )
499  enddo
500 
501  !##### RK3 : PROG0,PROG_RK2->PROG #####
502 
503  call prof_rapstart("DYN_RK3",3)
504 
505  dtrk = dt
506 
507  call atmos_dyn_tstep( dens, momz, momx, momy, rhot, & ! [OUT]
508  prog, & ! [OUT]
509  mflx_hi, tflx_hi, & ! [INOUT,OUT]
510  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
511  dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [IN]
512  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
513  prog0, prog_rk2, & ! [IN]
514  dpres0, cvtot, corioli, & ! [IN]
515  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
516  flag_fct_momentum, flag_fct_t, & ! [IN]
517  flag_fct_along_stream, & ! [IN]
518  cdz, fdz, fdx, fdy, & ! [IN]
519  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
520  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
521  ref_pres, ref_dens, & ! [IN]
522  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
523  dtrk, .true. ) ! [IN]
524 
525  if ( .NOT. flag_ws2002 ) then
526  do j = js, je
527  do i = is, ie
528  do k = ks, ke
529  dens(k,i,j) = ( dens_rk1(k,i,j) * 3.0_rp &
530  + dens(k,i,j) * 3.0_rp &
531  - dens0(k,i,j) * 2.0_rp ) / 4.0_rp
532  enddo
533  enddo
534  enddo
535 
536  do j = js, je
537  do i = is, ie
538  do k = ks, ke-1
539  momz(k,i,j) = ( momz_rk1(k,i,j) * 3.0_rp &
540  + momz(k,i,j) * 3.0_rp &
541  - momz0(k,i,j) * 2.0_rp ) / 4.0_rp
542  enddo
543  enddo
544  enddo
545 
546  do j = js, je
547  do i = is, ie
548  do k = ks, ke
549  momx(k,i,j) = ( momx_rk1(k,i,j) * 3.0_rp &
550  + momx(k,i,j) * 3.0_rp &
551  - momx0(k,i,j) * 2.0_rp ) / 4.0_rp
552  enddo
553  enddo
554  enddo
555 
556  do j = js, je
557  do i = is, ie
558  do k = ks, ke
559  momy(k,i,j) = ( momy_rk1(k,i,j) * 3.0_rp &
560  + momy(k,i,j) * 3.0_rp &
561  - momy0(k,i,j) * 2.0_rp ) / 4.0_rp
562  enddo
563  enddo
564  enddo
565 
566  do j = js, je
567  do i = is, ie
568  do k = ks, ke
569  rhot(k,i,j) = ( rhot_rk1(k,i,j) * 3.0_rp &
570  + rhot(k,i,j) * 3.0_rp &
571  - rhot0(k,i,j) * 2.0_rp ) / 4.0_rp
572  enddo
573  enddo
574  enddo
575 
576  do iv = 1, va
577  do j = js, je
578  do i = is, ie
579  do k = ks, ke
580  prog(k,i,j,iv) = ( prog_rk1(k,i,j,iv) * 3.0_rp &
581  + prog(k,i,j,iv) * 3.0_rp &
582  - prog0(k,i,j,iv) * 2.0_rp ) / 4.0_rp
583  enddo
584  enddo
585  enddo
586  enddo
587 
588  do n = 1, 3
589  do j = js, je
590  do i = is, ie
591  do k = ks, ke
592  mflx_hi(k,i,j,n) = ( mflx_hi_rk(k,i,j,n,1) &
593  + mflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
594  enddo
595  enddo
596  enddo
597  enddo
598 
599  do n = 1, 3
600  do j = js, je
601  do i = is, ie
602  do k = ks, ke
603  tflx_hi(k,i,j,n) = ( tflx_hi_rk(k,i,j,n,1) &
604  + tflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
605  enddo
606  enddo
607  enddo
608  enddo
609 
610  endif
611 
612  call prof_rapend ("DYN_RK3",3)
613 
614  return
615  end subroutine atmos_dyn_tinteg_short_rk3
616 
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:350
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:40
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:174
scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3_finalize
subroutine, public atmos_dyn_tinteg_short_rk3_finalize
finalize
Definition: scale_atmos_dyn_tinteg_short_rk3.F90:204
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:45
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_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:240
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:59
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:94
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:246
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_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:215
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56