SCALE-RM
scale_atmos_dyn_tinteg_short_rk3.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
15 #include "inc_openmp.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
25  use scale_index
26  use scale_tracer
27 #ifdef DEBUG
28  use scale_debug, only: &
29  check
30  use scale_const, only: &
31  undef => const_undef, &
32  iundef => const_undef2
33 #endif
34  !-----------------------------------------------------------------------------
35  implicit none
36  private
37  !-----------------------------------------------------------------------------
38  !
39  !++ Public procedure
40  !
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 (+1/3 step)
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 (+2/3 step)
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 = 1
71  integer :: i_comm_momz_rk1 = 2
72  integer :: i_comm_momx_rk1 = 3
73  integer :: i_comm_momy_rk1 = 4
74  integer :: i_comm_rhot_rk1 = 5
75  integer, allocatable :: i_comm_prog_rk1(:)
76 
77  integer :: i_comm_dens_rk2 = 1
78  integer :: i_comm_momz_rk2 = 2
79  integer :: i_comm_momx_rk2 = 3
80  integer :: i_comm_momy_rk2 = 4
81  integer :: i_comm_rhot_rk2 = 5
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_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
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  Rtot, CVtot, CORIOLI, &
196  num_diff, 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, &
204  dt )
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
576  end subroutine atmos_dyn_tinteg_short_rk3
577 
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
integer, public va
Definition: scale_index.F90:38
module Atmosphere / Dynamical scheme
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
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.
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:58
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
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)
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
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module Atmosphere / Dynamics common
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
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)
subroutine, public atmos_dyn_tinteg_short_rk3_setup(tinteg_type)
Setup.
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
integer, public ja
of y whole cells (local, with HALO)