SCALE-RM
scale_atmos_dyn_tinteg_short_rk3.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
12 #include "scalelib.h"
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_io
20  use scale_prof
22  use scale_index
23  use scale_tracer
24 #if defined DEBUG || defined QUICKDEBUG
25  use scale_debug, only: &
26  check
27  use scale_const, only: &
28  undef => const_undef, &
29  iundef => const_undef2
30 #endif
31  !-----------------------------------------------------------------------------
32  implicit none
33  private
34  !-----------------------------------------------------------------------------
35  !
36  !++ Public procedure
37  !
40 
41  !-----------------------------------------------------------------------------
42  !
43  !++ Public parameters & variables
44  !
45  !-----------------------------------------------------------------------------
46  !
47  !++ Private procedure
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private parameters & variables
52  !
53  real(RP), private, allocatable :: dens_rk1(:,:,:) ! prognostic variables (+1/3 step)
54  real(RP), private, allocatable :: momz_rk1(:,:,:)
55  real(RP), private, allocatable :: momx_rk1(:,:,:)
56  real(RP), private, allocatable :: momy_rk1(:,:,:)
57  real(RP), private, allocatable :: rhot_rk1(:,:,:)
58  real(RP), private, allocatable :: prog_rk1(:,:,:,:)
59  real(RP), private, allocatable :: dens_rk2(:,:,:) ! prognostic variables (+2/3 step)
60  real(RP), private, allocatable :: momz_rk2(:,:,:)
61  real(RP), private, allocatable :: momx_rk2(:,:,:)
62  real(RP), private, allocatable :: momy_rk2(:,:,:)
63  real(RP), private, allocatable :: rhot_rk2(:,:,:)
64  real(RP), private, allocatable :: prog_rk2(:,:,:,:)
65 
66  ! for communication
67  integer :: i_comm_dens_rk1 = 1
68  integer :: i_comm_momz_rk1 = 2
69  integer :: i_comm_momx_rk1 = 3
70  integer :: i_comm_momy_rk1 = 4
71  integer :: i_comm_rhot_rk1 = 5
72  integer, allocatable :: i_comm_prog_rk1(:)
73 
74  integer :: i_comm_dens_rk2 = 1
75  integer :: i_comm_momz_rk2 = 2
76  integer :: i_comm_momx_rk2 = 3
77  integer :: i_comm_momy_rk2 = 4
78  integer :: i_comm_rhot_rk2 = 5
79  integer, allocatable :: i_comm_prog_rk2(:)
80 
81  logical :: flag_ws2002
82  real(RP) :: fact_dt1
83  real(RP) :: fact_dt2
84 
85  !-----------------------------------------------------------------------------
86 contains
87  !-----------------------------------------------------------------------------
90  tinteg_type )
91  use scale_prc, only: &
92  prc_abort
93  use scale_const, only: &
94  undef => const_undef
95  use scale_comm_cartesc, only: &
97  implicit none
98 
99  character(len=*) :: tinteg_type
100 
101  integer :: iv
102  !---------------------------------------------------------------------------
103 
104  select case( tinteg_type )
105  case( 'RK3' )
106  log_info("ATMOS_DYN_Tinteg_short_rk3_setup",*) "RK3: Heun's method is used"
107  ! Heun's method
108  ! k1 = f(\phi_n); r1 = \phi_n + k1 * dt / 3
109  ! k2 = f(r1); r2 = \phi_n + k2 * dt * 2 / 3
110  ! k3 = f(r2); r3 = \phi_n + k3 * dt
111  ! \phi_{n+1} = \phi_n + ( k1 + 3 * k3 ) dt / 4
112  ! = \phi_n + ( (r1-\phi_n) * 3 + (r3-\phi_n) * 3 ) / 4
113  ! = ( r1 * 3 + r3 * 3 - \phi_n * 2 ) / 4
114  flag_ws2002 = .false.
115  fact_dt1 = 1.0_rp / 3.0_rp
116  fact_dt2 = 2.0_rp / 3.0_rp
117  case( 'RK3WS2002' )
118  log_info("ATMOS_DYN_Tinteg_short_rk3_setup",*) "RK3: Wichere and Skamarock (2002) is used"
119  ! Wicher and Skamarock (2002) RK3 scheme
120  ! k1 = f(\phi_n); r1 = \phi_n + k1 * dt / 3
121  ! k2 = f(r1); r2 = \phi_n + k2 * dt / 2
122  ! k3 = f(r2); r3 = \phi_n + k3 * dt
123  ! \phi_{n+1} = r3
124  flag_ws2002 = .true.
125  fact_dt1 = 1.0_rp / 3.0_rp
126  fact_dt2 = 1.0_rp / 2.0_rp
127  case default
128  log_error("ATMOS_DYN_Tinteg_short_rk3_setup",*) 'TINTEG_TYPE is not RK3. Check!'
129  call prc_abort
130  end select
131 
132  allocate( dens_rk1(ka,ia,ja) )
133  allocate( momz_rk1(ka,ia,ja) )
134  allocate( momx_rk1(ka,ia,ja) )
135  allocate( momy_rk1(ka,ia,ja) )
136  allocate( rhot_rk1(ka,ia,ja) )
137 
138  allocate( dens_rk2(ka,ia,ja) )
139  allocate( momz_rk2(ka,ia,ja) )
140  allocate( momx_rk2(ka,ia,ja) )
141  allocate( momy_rk2(ka,ia,ja) )
142  allocate( rhot_rk2(ka,ia,ja) )
143 
144  allocate( prog_rk1(ka,ia,ja,max(va,1)) )
145  allocate( prog_rk2(ka,ia,ja,max(va,1)) )
146  allocate( i_comm_prog_rk1(max(va,1)) )
147  allocate( i_comm_prog_rk2(max(va,1)) )
148 
149  call comm_vars8_init( 'DENS_RK1', dens_rk1, i_comm_dens_rk1 )
150  call comm_vars8_init( 'MOMZ_RK1', momz_rk1, i_comm_momz_rk1 )
151  call comm_vars8_init( 'MOMX_RK1', momx_rk1, i_comm_momx_rk1 )
152  call comm_vars8_init( 'MOMY_RK1', momy_rk1, i_comm_momy_rk1 )
153  call comm_vars8_init( 'RHOT_RK1', rhot_rk1, i_comm_rhot_rk1 )
154  do iv = 1, va
155  i_comm_prog_rk1(iv) = 5 + iv
156  call comm_vars8_init( 'PROG_RK1', prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
157  enddo
158 
159  call comm_vars8_init( 'DENS_RK2', dens_rk2, i_comm_dens_rk2 )
160  call comm_vars8_init( 'MOMZ_RK2', momz_rk2, i_comm_momz_rk2 )
161  call comm_vars8_init( 'MOMX_RK2', momx_rk2, i_comm_momx_rk2 )
162  call comm_vars8_init( 'MOMY_RK2', momy_rk2, i_comm_momy_rk2 )
163  call comm_vars8_init( 'RHOT_RK2', rhot_rk2, i_comm_rhot_rk2 )
164  do iv = 1, va
165  i_comm_prog_rk2(iv) = 5 + iv
166  call comm_vars8_init( 'PROG_RK2', prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
167  enddo
168 
169  dens_rk1(:,:,:) = undef
170  momz_rk1(:,:,:) = undef
171  momx_rk1(:,:,:) = undef
172  momy_rk1(:,:,:) = undef
173  rhot_rk1(:,:,:) = undef
174  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
175 
176  dens_rk2(:,:,:) = undef
177  momz_rk2(:,:,:) = undef
178  momx_rk2(:,:,:) = undef
179  momy_rk2(:,:,:) = undef
180  rhot_rk2(:,:,:) = undef
181  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
182 
183  return
184  end subroutine atmos_dyn_tinteg_short_rk3_setup
185 
186  !-----------------------------------------------------------------------------
188  subroutine atmos_dyn_tinteg_short_rk3( &
189  DENS, MOMZ, MOMX, MOMY, RHOT, PROG, &
190  mflx_hi, tflx_hi, &
191  DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, &
192  DPRES0, CVtot, CORIOLI, &
193  num_diff, wdamp_coef, divdmp_coef, DDIV, &
194  FLAG_FCT_MOMENTUM, FLAG_FCT_T, &
195  FLAG_FCT_ALONG_STREAM, &
196  CDZ, FDZ, FDX, FDY, &
197  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
198  PHI, GSQRT, J13G, J23G, J33G, MAPF, &
199  REF_pres, REF_dens, &
200  BND_W, BND_E, BND_S, BND_N, &
201  dt )
202  use scale_comm_cartesc, only: &
203  comm_vars8, &
204  comm_wait
205  use scale_atmos_dyn_tstep_short, only: &
206  atmos_dyn_tstep => atmos_dyn_tstep_short
207  use scale_atmos_dyn_common, only: &
209  implicit none
210 
211  real(RP), intent(inout) :: DENS(ka,ia,ja)
212  real(RP), intent(inout) :: MOMZ(ka,ia,ja)
213  real(RP), intent(inout) :: MOMX(ka,ia,ja)
214  real(RP), intent(inout) :: MOMY(ka,ia,ja)
215  real(RP), intent(inout) :: RHOT(ka,ia,ja)
216  real(RP), intent(inout) :: PROG(ka,ia,ja,va)
217 
218  real(RP), intent(inout) :: mflx_hi(ka,ia,ja,3)
219  real(RP), intent(inout) :: tflx_hi(ka,ia,ja,3)
220 
221  real(RP), intent(in) :: DENS_t(ka,ia,ja)
222  real(RP), intent(in) :: MOMZ_t(ka,ia,ja)
223  real(RP), intent(in) :: MOMX_t(ka,ia,ja)
224  real(RP), intent(in) :: MOMY_t(ka,ia,ja)
225  real(RP), intent(in) :: RHOT_t(ka,ia,ja)
226 
227  real(RP), intent(in) :: DPRES0(ka,ia,ja)
228  real(RP), intent(in) :: CVtot(ka,ia,ja)
229  real(RP), intent(in) :: CORIOLI(ia,ja)
230  real(RP), intent(in) :: num_diff(ka,ia,ja,5,3)
231  real(RP), intent(in) :: wdamp_coef(ka)
232  real(RP), intent(in) :: divdmp_coef
233  real(RP), intent(in) :: DDIV(ka,ia,ja)
234 
235  logical, intent(in) :: FLAG_FCT_MOMENTUM
236  logical, intent(in) :: FLAG_FCT_T
237  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
238 
239  real(RP), intent(in) :: CDZ (ka)
240  real(RP), intent(in) :: FDZ (ka-1)
241  real(RP), intent(in) :: FDX (ia-1)
242  real(RP), intent(in) :: FDY (ja-1)
243  real(RP), intent(in) :: RCDZ(ka)
244  real(RP), intent(in) :: RCDX(ia)
245  real(RP), intent(in) :: RCDY(ja)
246  real(RP), intent(in) :: RFDZ(ka-1)
247  real(RP), intent(in) :: RFDX(ia-1)
248  real(RP), intent(in) :: RFDY(ja-1)
249 
250  real(RP), intent(in) :: PHI (ka,ia,ja)
251  real(RP), intent(in) :: GSQRT(ka,ia,ja,7)
252  real(RP), intent(in) :: J13G (ka,ia,ja,7)
253  real(RP), intent(in) :: J23G (ka,ia,ja,7)
254  real(RP), intent(in) :: J33G
255  real(RP), intent(in) :: MAPF (ia,ja,2,4)
256 
257  real(RP), intent(in) :: REF_pres(ka,ia,ja)
258  real(RP), intent(in) :: REF_dens(ka,ia,ja)
259 
260  logical, intent(in) :: BND_W
261  logical, intent(in) :: BND_E
262  logical, intent(in) :: BND_S
263  logical, intent(in) :: BND_N
264 
265  real(RP), intent(in) :: dt
266 
267  real(RP) :: DENS0(ka,ia,ja)
268  real(RP) :: MOMZ0(ka,ia,ja)
269  real(RP) :: MOMX0(ka,ia,ja)
270  real(RP) :: MOMY0(ka,ia,ja)
271  real(RP) :: RHOT0(ka,ia,ja)
272  real(RP) :: PROG0(ka,ia,ja,va)
273 
274  real(RP) :: mflx_hi_RK(ka,ia,ja,3,2)
275  real(RP) :: tflx_hi_RK(ka,ia,ja,3,2)
276 
277  real(RP) :: dtrk
278 
279  integer :: i, j, k, iv, n
280  !---------------------------------------------------------------------------
281 
282  call prof_rapstart("DYN_RK3_Prep",3)
283 
284 #ifdef DEBUG
285  dens_rk1(:,:,:) = undef
286  momz_rk1(:,:,:) = undef
287  momx_rk1(:,:,:) = undef
288  momy_rk1(:,:,:) = undef
289  rhot_rk1(:,:,:) = undef
290  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
291 
292  dens_rk2(:,:,:) = undef
293  momz_rk2(:,:,:) = undef
294  momx_rk2(:,:,:) = undef
295  momy_rk2(:,:,:) = undef
296  rhot_rk2(:,:,:) = undef
297  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
298 
299  mflx_hi_rk(:,:,:,:,:) = undef
300  tflx_hi_rk(:,:,:,:,:) = undef
301 #endif
302 
303 #ifdef QUICKDEBUG
304  mflx_hi( 1:ks-1,:,:,:) = undef
305  mflx_hi(ke+1:ka ,:,:,:) = 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  dpres0, cvtot, corioli, & ! [IN]
370  num_diff, wdamp_coef, 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, .false. ) ! [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  dpres0, cvtot, corioli, & ! [IN]
423  num_diff, wdamp_coef, 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, .false. ) ! [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  dpres0, cvtot, corioli, & ! [IN]
476  num_diff, wdamp_coef, 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, .true. ) ! [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 
module DEBUG
Definition: scale_debug.F90:11
integer, public va
Definition: scale_index.F90:35
module Atmosphere / Dynamical scheme
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
real(rp), public const_undef
Definition: scale_const.F90:41
module COMMUNICATION
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
procedure(short), pointer, public atmos_dyn_tstep_short
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, 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, 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, dt)
RK3.
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
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:11
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
module STDIO
Definition: scale_io.F90:10
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210