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 #if defined DEBUG || defined QUICKDEBUG
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', 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  Rtot, 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, &
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) :: 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 
268  real(RP), intent(in) :: dt
269 
270  real(RP) :: dens0(ka,ia,ja)
271  real(RP) :: momz0(ka,ia,ja)
272  real(RP) :: momx0(ka,ia,ja)
273  real(RP) :: momy0(ka,ia,ja)
274  real(RP) :: rhot0(ka,ia,ja)
275  real(RP) :: prog0(ka,ia,ja,va)
276 
277  real(RP) :: mflx_hi_rk(ka,ia,ja,3,2)
278  real(RP) :: tflx_hi_rk(ka,ia,ja,3,2)
279 
280  real(RP) :: dtrk
281 
282  integer :: i, j, k, iv, n
283  !---------------------------------------------------------------------------
284 
285  call prof_rapstart("DYN_RK3_Prep",3)
286 
287 #ifdef DEBUG
288  dens_rk1(:,:,:) = undef
289  momz_rk1(:,:,:) = undef
290  momx_rk1(:,:,:) = undef
291  momy_rk1(:,:,:) = undef
292  rhot_rk1(:,:,:) = undef
293  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
294 
295  dens_rk2(:,:,:) = undef
296  momz_rk2(:,:,:) = undef
297  momx_rk2(:,:,:) = undef
298  momy_rk2(:,:,:) = undef
299  rhot_rk2(:,:,:) = undef
300  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
301 
302  mflx_hi_rk(:,:,:,:,:) = undef
303  tflx_hi_rk(:,:,:,:,:) = undef
304 #endif
305 
306 #ifdef QUICKDEBUG
307  mflx_hi( 1:ks-1,:,:,:) = undef
308  mflx_hi(ke+1:ka ,:,:,:) = undef
309 #endif
310 
311 !OCL XFILL
312  dens0 = dens
313 !OCL XFILL
314  momz0 = momz
315 !OCL XFILL
316  momx0 = momx
317 !OCL XFILL
318  momy0 = momy
319 !OCL XFILL
320  rhot0 = rhot
321 !OCL XFILL
322  if ( va > 0 ) prog0 = prog
323 
324  if ( bnd_w ) then
325  do j = js, je
326  do k = ks, ke
327  mflx_hi_rk(k,is-1,j,2,:) = mflx_hi(k,is-1,j,2)
328  end do
329  end do
330  end if
331  if ( bnd_e ) then
332  do j = js, je
333  do k = ks, ke
334  mflx_hi_rk(k,ie,j,2,:) = mflx_hi(k,ie,j,2)
335  end do
336  end do
337  end if
338  if ( bnd_s ) then
339  do i = is, ie
340  do k = ks, ke
341  mflx_hi_rk(k,i,js-1,3,:) = mflx_hi(k,i,js-1,3)
342  end do
343  end do
344  end if
345  if ( bnd_n ) then
346  do i = is, ie
347  do k = ks, ke
348  mflx_hi_rk(k,i,je,3,:) = mflx_hi(k,i,je,3)
349  end do
350  end do
351  end if
352 
353  call prof_rapend ("DYN_RK3_Prep",3)
354 
355  !------------------------------------------------------------------------
356  ! Start RK
357  !------------------------------------------------------------------------
358 
359  !##### RK1 : PROG0,PROG->PROG_RK1 #####
360 
361  call prof_rapstart("DYN_RK3",3)
362 
363  dtrk = dt * fact_dt1
364 
365  call atmos_dyn_tstep( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [OUT]
366  prog_rk1, & ! [OUT]
367  mflx_hi_rk(:,:,:,:,1), tflx_hi_rk(:,:,:,:,1), & ! [INOUT]
368  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
369  dens, momz, momx, momy, rhot, & ! [IN]
370  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
371  prog0, prog, & ! [IN]
372  rtot, cvtot, corioli, & ! [IN]
373  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
374  flag_fct_momentum, flag_fct_t, & ! [IN]
375  flag_fct_along_stream, & ! [IN]
376  cdz, fdz, fdx, fdy, & ! [IN]
377  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
378  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
379  ref_pres, ref_dens, & ! [IN]
380  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
381  dtrk, .false. ) ! [IN]
382 
383  call prof_rapend ("DYN_RK3",3)
384  call prof_rapstart("DYN_RK3_BND",3)
385 
386  call atmos_dyn_copy_boundary( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [INOUT]
387  prog_rk1, & ! [INOUT]
388  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
389  prog0, & ! [IN]
390  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
391 
392  call prof_rapend ("DYN_RK3_BND",3)
393 
394  call comm_vars8( dens_rk1(:,:,:), i_comm_dens_rk1 )
395  call comm_vars8( momz_rk1(:,:,:), i_comm_momz_rk1 )
396  call comm_vars8( momx_rk1(:,:,:), i_comm_momx_rk1 )
397  call comm_vars8( momy_rk1(:,:,:), i_comm_momy_rk1 )
398  call comm_vars8( rhot_rk1(:,:,:), i_comm_rhot_rk1 )
399  do iv = 1, va
400  call comm_vars8( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
401  enddo
402 
403  call comm_wait ( dens_rk1(:,:,:), i_comm_dens_rk1, .false. )
404  call comm_wait ( momz_rk1(:,:,:), i_comm_momz_rk1, .false. )
405  call comm_wait ( momx_rk1(:,:,:), i_comm_momx_rk1, .false. )
406  call comm_wait ( momy_rk1(:,:,:), i_comm_momy_rk1, .false. )
407  call comm_wait ( rhot_rk1(:,:,:), i_comm_rhot_rk1, .false. )
408  do iv = 1, va
409  call comm_wait ( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv), .false. )
410  enddo
411 
412  !##### RK2 : PROG0,PROG_RK1->PROG_RK2 #####
413 
414  call prof_rapstart("DYN_RK3",3)
415 
416  dtrk = dt * fact_dt2
417 
418  call atmos_dyn_tstep( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [OUT]
419  prog_rk2, & ! [OUT]
420  mflx_hi_rk(:,:,:,:,2), tflx_hi_rk(:,:,:,:,2), & ! [INOUT]
421  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
422  dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [IN]
423  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
424  prog0, prog_rk1, & ! [IN]
425  rtot, cvtot, corioli, & ! [IN]
426  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
427  flag_fct_momentum, flag_fct_t, & ! [IN]
428  flag_fct_along_stream, & ! [IN]
429  cdz, fdz, fdx, fdy, & ! [IN]
430  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
431  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
432  ref_pres, ref_dens, & ! [IN]
433  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
434  dtrk, .false. ) ! [IN]
435 
436  call prof_rapend ("DYN_RK3",3)
437  call prof_rapstart("DYN_RK3_BND",3)
438 
439  call atmos_dyn_copy_boundary( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [INOUT]
440  prog_rk2, & ! [INOUT]
441  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
442  prog0, & ! [IN]
443  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
444 
445  call prof_rapend ("DYN_RK3_BND",3)
446 
447  call comm_vars8( dens_rk2(:,:,:), i_comm_dens_rk2 )
448  call comm_vars8( momz_rk2(:,:,:), i_comm_momz_rk2 )
449  call comm_vars8( momx_rk2(:,:,:), i_comm_momx_rk2 )
450  call comm_vars8( momy_rk2(:,:,:), i_comm_momy_rk2 )
451  call comm_vars8( rhot_rk2(:,:,:), i_comm_rhot_rk2 )
452  do iv = 1, va
453  call comm_vars8( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
454  enddo
455 
456  call comm_wait ( dens_rk2(:,:,:), i_comm_dens_rk2, .false. )
457  call comm_wait ( momz_rk2(:,:,:), i_comm_momz_rk2, .false. )
458  call comm_wait ( momx_rk2(:,:,:), i_comm_momx_rk2, .false. )
459  call comm_wait ( momy_rk2(:,:,:), i_comm_momy_rk2, .false. )
460  call comm_wait ( rhot_rk2(:,:,:), i_comm_rhot_rk2, .false. )
461  do iv = 1, va
462  call comm_wait ( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv), .false. )
463  enddo
464 
465  !##### RK3 : PROG0,PROG_RK2->PROG #####
466 
467  call prof_rapstart("DYN_RK3",3)
468 
469  dtrk = dt
470 
471  call atmos_dyn_tstep( dens, momz, momx, momy, rhot, & ! [OUT]
472  prog, & ! [OUT]
473  mflx_hi, tflx_hi, & ! [INOUT]
474  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
475  dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [IN]
476  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
477  prog0, prog_rk2, & ! [IN]
478  rtot, cvtot, corioli, & ! [IN]
479  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
480  flag_fct_momentum, flag_fct_t, & ! [IN]
481  flag_fct_along_stream, & ! [IN]
482  cdz, fdz, fdx, fdy, & ! [IN]
483  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
484  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
485  ref_pres, ref_dens, & ! [IN]
486  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
487  dtrk, .true. ) ! [IN]
488 
489  if ( .NOT. flag_ws2002 ) then
490  do j = js, je
491  do i = is, ie
492  do k = ks, ke
493  dens(k,i,j) = ( dens_rk1(k,i,j) * 3.0_rp &
494  + dens(k,i,j) * 3.0_rp &
495  - dens0(k,i,j) * 2.0_rp ) / 4.0_rp
496  enddo
497  enddo
498  enddo
499 
500  do j = js, je
501  do i = is, ie
502  do k = ks, ke-1
503  momz(k,i,j) = ( momz_rk1(k,i,j) * 3.0_rp &
504  + momz(k,i,j) * 3.0_rp &
505  - momz0(k,i,j) * 2.0_rp ) / 4.0_rp
506  enddo
507  enddo
508  enddo
509 
510  do j = js, je
511  do i = is, ie
512  do k = ks, ke
513  momx(k,i,j) = ( momx_rk1(k,i,j) * 3.0_rp &
514  + momx(k,i,j) * 3.0_rp &
515  - momx0(k,i,j) * 2.0_rp ) / 4.0_rp
516  enddo
517  enddo
518  enddo
519 
520  do j = js, je
521  do i = is, ie
522  do k = ks, ke
523  momy(k,i,j) = ( momy_rk1(k,i,j) * 3.0_rp &
524  + momy(k,i,j) * 3.0_rp &
525  - momy0(k,i,j) * 2.0_rp ) / 4.0_rp
526  enddo
527  enddo
528  enddo
529 
530  do j = js, je
531  do i = is, ie
532  do k = ks, ke
533  rhot(k,i,j) = ( rhot_rk1(k,i,j) * 3.0_rp &
534  + rhot(k,i,j) * 3.0_rp &
535  - rhot0(k,i,j) * 2.0_rp ) / 4.0_rp
536  enddo
537  enddo
538  enddo
539 
540  do iv = 1, va
541  do j = js, je
542  do i = is, ie
543  do k = ks, ke
544  prog(k,i,j,iv) = ( prog_rk1(k,i,j,iv) * 3.0_rp &
545  + prog(k,i,j,iv) * 3.0_rp &
546  - prog0(k,i,j,iv) * 2.0_rp ) / 4.0_rp
547  enddo
548  enddo
549  enddo
550  enddo
551 
552  do n = 1, 3
553  do j = js, je
554  do i = is, ie
555  do k = ks, ke
556  mflx_hi(k,i,j,n) = ( mflx_hi_rk(k,i,j,n,1) &
557  + mflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
558  enddo
559  enddo
560  enddo
561  enddo
562 
563  do n = 1, 3
564  do j = js, je
565  do i = is, ie
566  do k = ks, ke
567  tflx_hi(k,i,j,n) = ( tflx_hi_rk(k,i,j,n,1) &
568  + tflx_hi(k,i,j,n ) * 3.0_rp ) / 4.0_rp
569  enddo
570  enddo
571  enddo
572  enddo
573 
574  endif
575 
576  call prof_rapend ("DYN_RK3",3)
577 
578  return
579  end subroutine atmos_dyn_tinteg_short_rk3
580 
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
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm.F90:313
module Atmosphere / Dynamical scheme
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
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 whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
procedure(short), pointer, public atmos_dyn_tstep_short
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, 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, dt)
RK3.
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:156
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:204
integer, public ja
of whole cells: y, local, with HALO