SCALE-RM
scale_atmos_dyn_tinteg_short_rk4.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
18 !-------------------------------------------------------------------------------
19 #include "scalelib.h"
21  !-----------------------------------------------------------------------------
22  !
23  !++ used modules
24  !
25  use scale_precision
26  use scale_io
27  use scale_prof
29  use scale_index
30  use scale_tracer
31 #if defined DEBUG || defined QUICKDEBUG
32  use scale_debug, only: &
33  check
34  use scale_const, only: &
35  undef => const_undef, &
36  iundef => const_undef2
37 #endif
38  !-----------------------------------------------------------------------------
39  implicit none
40  private
41  !-----------------------------------------------------------------------------
42  !
43  !++ Public procedure
44  !
48 
49  !-----------------------------------------------------------------------------
50  !
51  !++ Public parameters & variables
52  !
53  !-----------------------------------------------------------------------------
54  !
55  !++ Private procedure
56  !
57  !-----------------------------------------------------------------------------
58  !
59  !++ Private parameters & variables
60  !
61  real(RP), private, allocatable :: DENS_RK1(:,:,:) ! prognostic variables (registers for stage2)
62  real(RP), private, allocatable :: MOMZ_RK1(:,:,:)
63  real(RP), private, allocatable :: MOMX_RK1(:,:,:)
64  real(RP), private, allocatable :: MOMY_RK1(:,:,:)
65  real(RP), private, allocatable :: RHOT_RK1(:,:,:)
66  real(RP), private, allocatable :: PROG_RK1(:,:,:,:)
67  real(RP), private, allocatable :: DENS_RK2(:,:,:) ! prognostic variables (registers for stage3)
68  real(RP), private, allocatable :: MOMZ_RK2(:,:,:)
69  real(RP), private, allocatable :: MOMX_RK2(:,:,:)
70  real(RP), private, allocatable :: MOMY_RK2(:,:,:)
71  real(RP), private, allocatable :: RHOT_RK2(:,:,:)
72  real(RP), private, allocatable :: PROG_RK2(:,:,:,:)
73  real(RP), private, allocatable :: DENS_RK3(:,:,:) ! prognostic variables (registers for stage4)
74  real(RP), private, allocatable :: MOMZ_RK3(:,:,:)
75  real(RP), private, allocatable :: MOMX_RK3(:,:,:)
76  real(RP), private, allocatable :: MOMY_RK3(:,:,:)
77  real(RP), private, allocatable :: RHOT_RK3(:,:,:)
78  real(RP), private, allocatable :: PROG_RK3(:,:,:,:)
79 
80  ! for communication
81  integer :: I_COMM_DENS_RK1 = 1
82  integer :: I_COMM_MOMZ_RK1 = 2
83  integer :: I_COMM_MOMX_RK1 = 3
84  integer :: I_COMM_MOMY_RK1 = 4
85  integer :: I_COMM_RHOT_RK1 = 5
86  integer, allocatable :: I_COMM_PROG_RK1(:)
87 
88  integer :: I_COMM_DENS_RK2 = 1
89  integer :: I_COMM_MOMZ_RK2 = 2
90  integer :: I_COMM_MOMX_RK2 = 3
91  integer :: I_COMM_MOMY_RK2 = 4
92  integer :: I_COMM_RHOT_RK2 = 5
93  integer, allocatable :: I_COMM_PROG_RK2(:)
94 
95  integer :: I_COMM_DENS_RK3 = 1
96  integer :: I_COMM_MOMZ_RK3 = 2
97  integer :: I_COMM_MOMX_RK3 = 3
98  integer :: I_COMM_MOMY_RK3 = 4
99  integer :: I_COMM_RHOT_RK3 = 5
100  integer, allocatable :: I_COMM_PROG_RK3(:)
101 
102  !-----------------------------------------------------------------------------
103 contains
104  !-----------------------------------------------------------------------------
107  tinteg_type )
108  use scale_prc, only: &
109  prc_abort
110  use scale_const, only: &
111  undef => const_undef
112  use scale_comm_cartesc, only: &
114  implicit none
115 
116  character(len=*) :: tinteg_type
117 
118  integer :: iv
119  !---------------------------------------------------------------------------
120 
121  if ( tinteg_type /= 'RK4' ) then
122  log_error("ATMOS_DYN_Tinteg_short_rk4_setup",*) 'TINTEG_TYPE is not RK4. Check!'
123  call prc_abort
124  end if
125 
126  allocate( dens_rk1(ka,ia,ja) )
127  allocate( momz_rk1(ka,ia,ja) )
128  allocate( momx_rk1(ka,ia,ja) )
129  allocate( momy_rk1(ka,ia,ja) )
130  allocate( rhot_rk1(ka,ia,ja) )
131  !$acc enter data create(DENS_RK1, MOMZ_RK1, MOMX_RK1, MOMY_RK1, RHOT_RK1)
132 
133  allocate( dens_rk2(ka,ia,ja) )
134  allocate( momz_rk2(ka,ia,ja) )
135  allocate( momx_rk2(ka,ia,ja) )
136  allocate( momy_rk2(ka,ia,ja) )
137  allocate( rhot_rk2(ka,ia,ja) )
138  !$acc enter data create(DENS_RK2, MOMZ_RK2, MOMX_RK2, MOMY_RK2, RHOT_RK2)
139 
140  allocate( dens_rk3(ka,ia,ja) )
141  allocate( momz_rk3(ka,ia,ja) )
142  allocate( momx_rk3(ka,ia,ja) )
143  allocate( momy_rk3(ka,ia,ja) )
144  allocate( rhot_rk3(ka,ia,ja) )
145  !$acc enter data create(DENS_RK3, MOMZ_RK3, MOMX_RK3, MOMY_RK3, RHOT_RK3)
146 
147  allocate( prog_rk1(ka,ia,ja,max(va,1)) )
148  allocate( prog_rk2(ka,ia,ja,max(va,1)) )
149  allocate( prog_rk3(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  allocate( i_comm_prog_rk3(max(va,1)) )
153  !$acc enter data create(PROG_RK1, PROG_RK2, PROG_RK3)
154 
155  i_comm_dens_rk1 = 1
156  i_comm_momz_rk1 = 2
157  i_comm_momx_rk1 = 3
158  i_comm_momy_rk1 = 4
159  i_comm_rhot_rk1 = 5
160  call comm_vars8_init( 'DENS_RK1', dens_rk1, i_comm_dens_rk1 )
161  call comm_vars8_init( 'MOMZ_RK1', momz_rk1, i_comm_momz_rk1 )
162  call comm_vars8_init( 'MOMX_RK1', momx_rk1, i_comm_momx_rk1 )
163  call comm_vars8_init( 'MOMY_RK1', momy_rk1, i_comm_momy_rk1 )
164  call comm_vars8_init( 'RHOT_RK1', rhot_rk1, i_comm_rhot_rk1 )
165  do iv = 1, va
166  i_comm_prog_rk1(iv) = 5 + iv
167  call comm_vars8_init( 'PROG_RK1', prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
168  enddo
169 
170  i_comm_dens_rk2 = 1
171  i_comm_momz_rk2 = 2
172  i_comm_momx_rk2 = 3
173  i_comm_momy_rk2 = 4
174  i_comm_rhot_rk2 = 5
175  call comm_vars8_init( 'DENS_RK2', dens_rk2, i_comm_dens_rk2 )
176  call comm_vars8_init( 'MOMZ_RK2', momz_rk2, i_comm_momz_rk2 )
177  call comm_vars8_init( 'MOMX_RK2', momx_rk2, i_comm_momx_rk2 )
178  call comm_vars8_init( 'MOMY_RK2', momy_rk2, i_comm_momy_rk2 )
179  call comm_vars8_init( 'RHOT_RK2', rhot_rk2, i_comm_rhot_rk2 )
180  do iv = 1, va
181  i_comm_prog_rk2(iv) = 5 + iv
182  call comm_vars8_init( 'PROG_RK2', prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
183  enddo
184 
185  i_comm_dens_rk3 = 1
186  i_comm_momz_rk3 = 2
187  i_comm_momx_rk3 = 3
188  i_comm_momy_rk3 = 4
189  i_comm_rhot_rk3 = 5
190  call comm_vars8_init( 'DENS_RK3', dens_rk3, i_comm_dens_rk3 )
191  call comm_vars8_init( 'MOMZ_RK3', momz_rk3, i_comm_momz_rk3 )
192  call comm_vars8_init( 'MOMX_RK3', momx_rk3, i_comm_momx_rk3 )
193  call comm_vars8_init( 'MOMY_RK3', momy_rk3, i_comm_momy_rk3 )
194  call comm_vars8_init( 'RHOT_RK3', rhot_rk3, i_comm_rhot_rk3 )
195  do iv = 1, va
196  i_comm_prog_rk3(iv) = 5 + iv
197  call comm_vars8_init( 'PROG_RK3', prog_rk3(:,:,:,iv), i_comm_prog_rk3(iv) )
198  enddo
199 
200  dens_rk1(:,:,:) = undef
201  momz_rk1(:,:,:) = undef
202  momx_rk1(:,:,:) = undef
203  momy_rk1(:,:,:) = undef
204  rhot_rk1(:,:,:) = undef
205  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
206 
207  dens_rk2(:,:,:) = undef
208  momz_rk2(:,:,:) = undef
209  momx_rk2(:,:,:) = undef
210  momy_rk2(:,:,:) = undef
211  rhot_rk2(:,:,:) = undef
212  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
213 
214  dens_rk3(:,:,:) = undef
215  momz_rk3(:,:,:) = undef
216  momx_rk3(:,:,:) = undef
217  momy_rk3(:,:,:) = undef
218  rhot_rk3(:,:,:) = undef
219  if ( va > 0 ) prog_rk3(:,:,:,:) = undef
220 
221  return
222  end subroutine atmos_dyn_tinteg_short_rk4_setup
223 
224  !-----------------------------------------------------------------------------
227 
228  !$acc exit data delete(DENS_RK1, MOMZ_RK1, MOMX_RK1, MOMY_RK1, RHOT_RK1)
229  deallocate( dens_rk1 )
230  deallocate( momz_rk1 )
231  deallocate( momx_rk1 )
232  deallocate( momy_rk1 )
233  deallocate( rhot_rk1 )
234 
235  !$acc exit data delete(DENS_RK2, MOMZ_RK2, MOMX_RK2, MOMY_RK2, RHOT_RK2)
236  deallocate( dens_rk2 )
237  deallocate( momz_rk2 )
238  deallocate( momx_rk2 )
239  deallocate( momy_rk2 )
240  deallocate( rhot_rk2 )
241 
242  !$acc exit data delete(DENS_RK3, MOMZ_RK3, MOMX_RK3, MOMY_RK3, RHOT_RK3)
243  deallocate( dens_rk3 )
244  deallocate( momz_rk3 )
245  deallocate( momx_rk3 )
246  deallocate( momy_rk3 )
247  deallocate( rhot_rk3 )
248 
249  !$acc exit data delete(PROG_RK1, PROG_RK2, PROG_RK3)
250  deallocate( prog_rk1 )
251  deallocate( prog_rk2 )
252  deallocate( prog_rk3 )
253  deallocate( i_comm_prog_rk1 )
254  deallocate( i_comm_prog_rk2 )
255  deallocate( i_comm_prog_rk3 )
256 
257  return
259  !-----------------------------------------------------------------------------
261  subroutine atmos_dyn_tinteg_short_rk4( &
262  DENS, MOMZ, MOMX, MOMY, RHOT, PROG, &
263  mflx_hi, tflx_hi, &
264  DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, &
265  DPRES0, CVtot, CORIOLI, &
266  num_diff, wdamp_coef, divdmp_coef, DDIV, &
267  FLAG_FCT_MOMENTUM, FLAG_FCT_T, &
268  FLAG_FCT_ALONG_STREAM, &
269  CDZ, FDZ, FDX, FDY, &
270  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
271  PHI, GSQRT, J13G, J23G, J33G, MAPF, &
272  REF_pres, REF_dens, &
273  BND_W, BND_E, BND_S, BND_N, TwoD, &
274  dt )
275  use scale_comm_cartesc, only: &
276  comm_vars8, &
277  comm_wait
278  use scale_atmos_dyn_tstep_short, only: &
279  atmos_dyn_tstep => atmos_dyn_tstep_short
280  use scale_atmos_dyn_common, only: &
282  implicit none
283 
284  real(rp), intent(inout) :: dens(ka,ia,ja)
285  real(rp), intent(inout) :: momz(ka,ia,ja)
286  real(rp), intent(inout) :: momx(ka,ia,ja)
287  real(rp), intent(inout) :: momy(ka,ia,ja)
288  real(rp), intent(inout) :: rhot(ka,ia,ja)
289  real(rp), intent(inout) :: prog(ka,ia,ja,va)
290 
291  real(rp), intent(inout) :: mflx_hi(ka,ia,ja,3)
292  real(rp), intent(out) :: tflx_hi(ka,ia,ja,3)
293 
294  real(rp), intent(in) :: dens_t(ka,ia,ja)
295  real(rp), intent(in) :: momz_t(ka,ia,ja)
296  real(rp), intent(in) :: momx_t(ka,ia,ja)
297  real(rp), intent(in) :: momy_t(ka,ia,ja)
298  real(rp), intent(in) :: rhot_t(ka,ia,ja)
299 
300  real(rp), intent(in) :: dpres0(ka,ia,ja)
301  real(rp), intent(in) :: cvtot(ka,ia,ja)
302  real(rp), intent(in) :: corioli(ia,ja)
303  real(rp), intent(in) :: num_diff(ka,ia,ja,5,3)
304  real(rp), intent(in) :: wdamp_coef(ka)
305  real(rp), intent(in) :: divdmp_coef
306  real(rp), intent(in) :: ddiv(ka,ia,ja)
307 
308  logical, intent(in) :: flag_fct_momentum
309  logical, intent(in) :: flag_fct_t
310  logical, intent(in) :: flag_fct_along_stream
311 
312  real(rp), intent(in) :: cdz (ka)
313  real(rp), intent(in) :: fdz (ka-1)
314  real(rp), intent(in) :: fdx (ia-1)
315  real(rp), intent(in) :: fdy (ja-1)
316  real(rp), intent(in) :: rcdz(ka)
317  real(rp), intent(in) :: rcdx(ia)
318  real(rp), intent(in) :: rcdy(ja)
319  real(rp), intent(in) :: rfdz(ka-1)
320  real(rp), intent(in) :: rfdx(ia-1)
321  real(rp), intent(in) :: rfdy(ja-1)
322 
323  real(rp), intent(in) :: phi (ka,ia,ja)
324  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
325  real(rp), intent(in) :: j13g (ka,ia,ja,7)
326  real(rp), intent(in) :: j23g (ka,ia,ja,7)
327  real(rp), intent(in) :: j33g
328  real(rp), intent(in) :: mapf (ia,ja,2,4)
329 
330  real(rp), intent(in) :: ref_pres(ka,ia,ja)
331  real(rp), intent(in) :: ref_dens(ka,ia,ja)
332 
333  logical, intent(in) :: bnd_w
334  logical, intent(in) :: bnd_e
335  logical, intent(in) :: bnd_s
336  logical, intent(in) :: bnd_n
337  logical, intent(in) :: twod
338 
339  real(rp), intent(in) :: dt
340 
341  real(rp) :: dens0(ka,ia,ja)
342  real(rp) :: momz0(ka,ia,ja)
343  real(rp) :: momx0(ka,ia,ja)
344  real(rp) :: momy0(ka,ia,ja)
345  real(rp) :: rhot0(ka,ia,ja)
346  real(rp) :: prog0(ka,ia,ja,va)
347 
348  real(rp) :: mflx_hi_rk(ka,ia,ja,3,3)
349  real(rp) :: tflx_hi_rk(ka,ia,ja,3,3)
350 
351  real(rp) :: dtrk
352 
353  integer :: i, j, k, iv, n
354  !---------------------------------------------------------------------------
355 
356  call prof_rapstart("DYN_RK4_Prep",3)
357 
358  !$acc data copy(DENS, MOMZ, MOMX, MOMY, RHOT, PROG, mflx_hi) &
359  !$acc copyout(tflx_hi) &
360  !$acc copyin(DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, &
361  !$acc DPRES0, CVtot, CORIOLI, num_diff, wdamp_coef, DDIV, &
362  !$acc CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
363  !$acc PHI, GSQRT, J13G, J23G, MAPF, &
364  !$acc REF_pres, REF_dens) &
365  !$acc create(DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, PROG0, mflx_hi_RK, tflx_hi_RK)
366 
367 #ifdef DEBUG
368  !$acc kernels
369  dens_rk1(:,:,:) = undef
370  momz_rk1(:,:,:) = undef
371  momx_rk1(:,:,:) = undef
372  momy_rk1(:,:,:) = undef
373  rhot_rk1(:,:,:) = undef
374  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
375 
376  dens_rk2(:,:,:) = undef
377  momz_rk2(:,:,:) = undef
378  momx_rk2(:,:,:) = undef
379  momy_rk2(:,:,:) = undef
380  rhot_rk2(:,:,:) = undef
381  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
382 
383  dens_rk3(:,:,:) = undef
384  momz_rk3(:,:,:) = undef
385  momx_rk3(:,:,:) = undef
386  momy_rk3(:,:,:) = undef
387  rhot_rk3(:,:,:) = undef
388  if ( va > 0 ) prog_rk3(:,:,:,:) = undef
389 
390  mflx_hi_rk(:,:,:,:,:) = undef
391  tflx_hi_rk(:,:,:,:,:) = undef
392  !$acc end kernels
393 #endif
394 
395 #ifdef QUICKDEBUG
396  !$acc kernels
397  mflx_hi( 1:ks-1,:,:,:) = undef
398  mflx_hi(ke+1:ka ,:,:,:) = undef
399  !$acc end kernels
400 #endif
401 
402  !$acc kernels
403 !OCL XFILL
404  dens0(:,:,:) = dens(:,:,:)
405  !$acc end kernels
406  !$acc kernels
407 !OCL XFILL
408  momz0(:,:,:) = momz(:,:,:)
409  !$acc end kernels
410  !$acc kernels
411 !OCL XFILL
412  momx0(:,:,:) = momx(:,:,:)
413  !$acc end kernels
414  !$acc kernels
415 !OCL XFILL
416  momy0(:,:,:) = momy(:,:,:)
417  !$acc end kernels
418  !$acc kernels
419 !OCL XFILL
420  rhot0(:,:,:) = rhot(:,:,:)
421  !$acc end kernels
422  if ( va > 0 ) then
423  !$acc kernels
424 !OCL XFILL
425  prog0(:,:,:,:) = prog(:,:,:,:)
426  !$acc end kernels
427  end if
428 
429  if ( bnd_w .and. (.not. twod) ) then
430  !$acc kernels
431  do j = js, je
432  do k = ks, ke
433  mflx_hi_rk(k,is-1,j,2,:) = mflx_hi(k,is-1,j,2)
434  end do
435  end do
436  !$acc end kernels
437  end if
438  if ( bnd_e .and. (.not. twod) ) then
439  !$acc kernels
440  do j = js, je
441  do k = ks, ke
442  mflx_hi_rk(k,ie,j,2,:) = mflx_hi(k,ie,j,2)
443  end do
444  end do
445  !$acc end kernels
446  end if
447  if ( bnd_s ) then
448  !$acc kernels
449  do i = is, ie
450  do k = ks, ke
451  mflx_hi_rk(k,i,js-1,3,:) = mflx_hi(k,i,js-1,3)
452  end do
453  end do
454  !$acc end kernels
455  end if
456  if ( bnd_n ) then
457  !$acc kernels
458  do i = is, ie
459  do k = ks, ke
460  mflx_hi_rk(k,i,je,3,:) = mflx_hi(k,i,je,3)
461  end do
462  end do
463  !$acc end kernels
464  end if
465 
466  call prof_rapend ("DYN_RK4_Prep",3)
467 
468  !------------------------------------------------------------------------
469  ! Start RK
470  !------------------------------------------------------------------------
471 
472  !##### RK1 : PROG0,PROG->PROG_RK1 #####
473 
474  call prof_rapstart("DYN_RK4",3)
475 
476  dtrk = dt / 2.0_rp
477 
478  call atmos_dyn_tstep( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [OUT]
479  prog_rk1, & ! [OUT]
480  mflx_hi_rk(:,:,:,:,1), tflx_hi_rk(:,:,:,:,1), & ! [INOUT,OUT]
481  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
482  dens, momz, momx, momy, rhot, & ! [IN]
483  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
484  prog0, prog, & ! [IN]
485  dpres0, cvtot, corioli, & ! [IN]
486  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
487  flag_fct_momentum, flag_fct_t, & ! [IN]
488  flag_fct_along_stream, & ! [IN]
489  cdz, fdz, fdx, fdy, & ! [IN]
490  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
491  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
492  ref_pres, ref_dens, & ! [IN]
493  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
494  dtrk, .false. ) ! [IN]
495 
496  call prof_rapend ("DYN_RK4",3)
497  call prof_rapstart("DYN_RK4_BND",3)
498 
499  call atmos_dyn_copy_boundary( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [INOUT]
500  prog_rk1, & ! [INOUT]
501  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
502  prog0, & ! [IN]
503  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
504 
505  call prof_rapend ("DYN_RK4_BND",3)
506 
507  call comm_vars8( dens_rk1(:,:,:), i_comm_dens_rk1 )
508  call comm_vars8( momz_rk1(:,:,:), i_comm_momz_rk1 )
509  call comm_vars8( momx_rk1(:,:,:), i_comm_momx_rk1 )
510  call comm_vars8( momy_rk1(:,:,:), i_comm_momy_rk1 )
511  call comm_vars8( rhot_rk1(:,:,:), i_comm_rhot_rk1 )
512  do iv = 1, va
513  call comm_vars8( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
514  enddo
515 
516  call comm_wait ( dens_rk1(:,:,:), i_comm_dens_rk1, .false. )
517  call comm_wait ( momz_rk1(:,:,:), i_comm_momz_rk1, .false. )
518  call comm_wait ( momx_rk1(:,:,:), i_comm_momx_rk1, .false. )
519  call comm_wait ( momy_rk1(:,:,:), i_comm_momy_rk1, .false. )
520  call comm_wait ( rhot_rk1(:,:,:), i_comm_rhot_rk1, .false. )
521  do iv = 1, va
522  call comm_wait ( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv), .false. )
523  enddo
524 
525  !##### RK2 : PROG0,PROG_RK1->PROG_RK2 #####
526 
527  call prof_rapstart("DYN_RK4",3)
528 
529  dtrk = dt / 2.0_rp
530 
531  call atmos_dyn_tstep( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [OUT]
532  prog_rk2, & ! [OUT]
533  mflx_hi_rk(:,:,:,:,2), tflx_hi_rk(:,:,:,:,2), & ! [INOUT,OUT]
534  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
535  dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [IN]
536  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
537  prog0, prog_rk1, & ! [IN]
538  dpres0, cvtot, corioli, & ! [IN]
539  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
540  flag_fct_momentum, flag_fct_t, & ! [IN]
541  flag_fct_along_stream, & ! [IN]
542  cdz, fdz, fdx, fdy, & ! [IN]
543  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
544  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
545  ref_pres, ref_dens, & ! [IN]
546  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
547  dtrk, .false. ) ! [IN]
548 
549  call prof_rapend ("DYN_RK4",3)
550  call prof_rapstart("DYN_RK4_BND",3)
551 
552  call atmos_dyn_copy_boundary( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [INOUT]
553  prog_rk2, & ! [INOUT]
554  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
555  prog0, & ! [IN]
556  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
557 
558  call prof_rapend ("DYN_RK4_BND",3)
559 
560  call comm_vars8( dens_rk2(:,:,:), i_comm_dens_rk2 )
561  call comm_vars8( momz_rk2(:,:,:), i_comm_momz_rk2 )
562  call comm_vars8( momx_rk2(:,:,:), i_comm_momx_rk2 )
563  call comm_vars8( momy_rk2(:,:,:), i_comm_momy_rk2 )
564  call comm_vars8( rhot_rk2(:,:,:), i_comm_rhot_rk2 )
565  do iv = 1, va
566  call comm_vars8( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
567  enddo
568 
569  call comm_wait ( dens_rk2(:,:,:), i_comm_dens_rk2, .false. )
570  call comm_wait ( momz_rk2(:,:,:), i_comm_momz_rk2, .false. )
571  call comm_wait ( momx_rk2(:,:,:), i_comm_momx_rk2, .false. )
572  call comm_wait ( momy_rk2(:,:,:), i_comm_momy_rk2, .false. )
573  call comm_wait ( rhot_rk2(:,:,:), i_comm_rhot_rk2, .false. )
574  do iv = 1, va
575  call comm_wait ( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv), .false. )
576  enddo
577 
578  !##### RK3 : PROG0,PROG_RK2->PROG_RK3 #####
579 
580  call prof_rapstart("DYN_RK4",3)
581 
582  dtrk = dt
583 
584  call atmos_dyn_tstep( dens_rk3, momz_rk3, momx_rk3, momy_rk3, rhot_rk3, & ! [OUT]
585  prog_rk3, & ! [OUT]
586  mflx_hi_rk(:,:,:,:,3), tflx_hi_rk(:,:,:,:,3), & ! [INOUT,OUT]
587  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
588  dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [IN]
589  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
590  prog0, prog_rk2, & ! [IN]
591  dpres0, cvtot, corioli, & ! [IN]
592  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
593  flag_fct_momentum, flag_fct_t, & ! [IN]
594  flag_fct_along_stream, & ! [IN]
595  cdz, fdz, fdx, fdy, & ! [IN]
596  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
597  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
598  ref_pres, ref_dens, & ! [IN]
599  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
600  dtrk, .false. ) ! [IN]
601 
602  call prof_rapend ("DYN_RK4",3)
603  call prof_rapstart("DYN_RK4_BND",3)
604 
605  call atmos_dyn_copy_boundary( dens_rk3, momz_rk3, momx_rk3, momy_rk3, rhot_rk3, & ! [INOUT]
606  prog_rk3, & ! [INOUT]
607  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
608  prog0, & ! [IN]
609  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
610 
611  call prof_rapend ("DYN_RK4_BND",3)
612 
613  call comm_vars8( dens_rk3(:,:,:), i_comm_dens_rk3 )
614  call comm_vars8( momz_rk3(:,:,:), i_comm_momz_rk3 )
615  call comm_vars8( momx_rk3(:,:,:), i_comm_momx_rk3 )
616  call comm_vars8( momy_rk3(:,:,:), i_comm_momy_rk3 )
617  call comm_vars8( rhot_rk3(:,:,:), i_comm_rhot_rk3 )
618  do iv = 1, va
619  call comm_vars8( prog_rk3(:,:,:,iv), i_comm_prog_rk3(iv) )
620  enddo
621 
622  call comm_wait ( dens_rk3(:,:,:), i_comm_dens_rk3, .false. )
623  call comm_wait ( momz_rk3(:,:,:), i_comm_momz_rk3, .false. )
624  call comm_wait ( momx_rk3(:,:,:), i_comm_momx_rk3, .false. )
625  call comm_wait ( momy_rk3(:,:,:), i_comm_momy_rk3, .false. )
626  call comm_wait ( rhot_rk3(:,:,:), i_comm_rhot_rk3, .false. )
627  do iv = 1, va
628  call comm_wait ( prog_rk3(:,:,:,iv), i_comm_prog_rk3(iv), .false. )
629  enddo
630 
631  !##### RK4 : PROG0,PROG_RK3->PROG #####
632 
633  call prof_rapstart("DYN_RK4",3)
634 
635  dtrk = dt
636 
637  call atmos_dyn_tstep( dens, momz, momx, momy, rhot, & ! [OUT]
638  prog, & ! [OUT]
639  mflx_hi, tflx_hi, & ! [INOUT,OUT]
640  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
641  dens_rk3, momz_rk3, momx_rk3, momy_rk3, rhot_rk3, & ! [IN]
642  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
643  prog0, prog_rk3, & ! [IN]
644  dpres0, cvtot, corioli, & ! [IN]
645  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
646  flag_fct_momentum, flag_fct_t, & ! [IN]
647  flag_fct_along_stream, & ! [IN]
648  cdz, fdz, fdx, fdy, & ! [IN]
649  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
650  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
651  ref_pres, ref_dens, & ! [IN]
652  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
653  dtrk, .true. ) ! [IN]
654 
655  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
656  !$omp shared(JS,JE,IS,IE,KS,KE,DENS,DENS_RK1,DENS_RK2,DENS_RK3,DENS0)
657  !$acc kernels
658  do j = js, je
659  do i = is, ie
660  do k = ks, ke
661  dens(k,i,j) = ( dens_rk1(k,i,j) * 2.0_rp &
662  + dens_rk2(k,i,j) * 4.0_rp &
663  + dens_rk3(k,i,j) * 2.0_rp &
664  + dens(k,i,j) &
665  - dens0(k,i,j) * 3.0_rp ) / 6.0_rp
666  enddo
667  enddo
668  enddo
669  !$acc end kernels
670 
671  !$acc kernels
672  do j = js, je
673  do i = is, ie
674  do k = ks, ke-1
675  momz(k,i,j) = ( momz_rk1(k,i,j) * 2.0_rp &
676  + momz_rk2(k,i,j) * 4.0_rp &
677  + momz_rk3(k,i,j) * 2.0_rp &
678  + momz(k,i,j) &
679  - momz0(k,i,j) * 3.0_rp ) / 6.0_rp
680  enddo
681  enddo
682  enddo
683  !$acc end kernels
684 
685  !$acc kernels
686  do j = js, je
687  do i = is, ie
688  do k = ks, ke
689  momx(k,i,j) = ( momx_rk1(k,i,j) * 2.0_rp &
690  + momx_rk2(k,i,j) * 4.0_rp &
691  + momx_rk3(k,i,j) * 2.0_rp &
692  + momx(k,i,j) &
693  - momx0(k,i,j) * 3.0_rp ) / 6.0_rp
694  enddo
695  enddo
696  enddo
697  !$acc end kernels
698 
699  !$acc kernels
700  do j = js, je
701  do i = is, ie
702  do k = ks, ke
703  momy(k,i,j) = ( momy_rk1(k,i,j) * 2.0_rp &
704  + momy_rk2(k,i,j) * 4.0_rp &
705  + momy_rk3(k,i,j) * 2.0_rp &
706  + momy(k,i,j) &
707  - momy0(k,i,j) * 3.0_rp ) / 6.0_rp
708  enddo
709  enddo
710  enddo
711  !$acc end kernels
712 
713  !$omp parallel do default(none) &
714  !$omp shared(JS,JE,IS,IE,KS,KE,RHOT,RHOT_RK1,RHOT_RK2,RHOT_RK3,RHOT0) &
715  !$omp private(i,j,k) OMP_SCHEDULE_ collapse(2)
716  !$acc kernels
717  do j = js, je
718  do i = is, ie
719  do k = ks, ke
720  rhot(k,i,j) = ( rhot_rk1(k,i,j) * 2.0_rp &
721  + rhot_rk2(k,i,j) * 4.0_rp &
722  + rhot_rk3(k,i,j) * 2.0_rp &
723  + rhot(k,i,j) &
724  - rhot0(k,i,j) * 3.0_rp ) / 6.0_rp
725  enddo
726  enddo
727  enddo
728  !$acc end kernels
729 
730  !$acc kernels
731  do iv = 1, va
732  do j = js, je
733  do i = is, ie
734  do k = ks, ke
735  prog(k,i,j,iv) = ( prog_rk1(k,i,j,iv) * 2.0_rp &
736  + prog_rk2(k,i,j,iv) * 4.0_rp &
737  + prog_rk3(k,i,j,iv) * 2.0_rp &
738  + prog(k,i,j,iv) &
739  - prog0(k,i,j,iv) * 3.0_rp ) / 6.0_rp
740  enddo
741  enddo
742  enddo
743  enddo
744  !$acc end kernels
745 
746 
747  call atmos_dyn_copy_boundary( dens, momz, momx, momy, rhot, & ! [INOUT]
748  prog, & ! [INOUT]
749  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
750  prog0, & ! [IN]
751  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
752 
753  !$acc kernels
754  do n = 1, 3
755  do j = js, je
756  do i = is, ie
757  do k = ks, ke
758  mflx_hi(k,i,j,n) = ( mflx_hi_rk(k,i,j,n,1) &
759  + mflx_hi_rk(k,i,j,n,2) * 2.0_rp &
760  + mflx_hi_rk(k,i,j,n,3) * 2.0_rp &
761  + mflx_hi(k,i,j,n ) ) / 6.0_rp
762  enddo
763  enddo
764  enddo
765  enddo
766  !$acc end kernels
767 
768  !$acc kernels
769  do n = 1, 3
770  do j = js, je
771  do i = is, ie
772  do k = ks, ke
773  tflx_hi(k,i,j,n) = ( tflx_hi_rk(k,i,j,n,1) &
774  + tflx_hi_rk(k,i,j,n,2) * 2.0_rp &
775  + tflx_hi_rk(k,i,j,n,3) * 2.0_rp &
776  + tflx_hi(k,i,j,n ) ) / 6.0_rp
777  enddo
778  enddo
779  enddo
780  enddo
781  !$acc end kernels
782 
783  !$acc end data
784 
785  call prof_rapend ("DYN_RK4",3)
786 
787  return
788  end subroutine atmos_dyn_tinteg_short_rk4
789 
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_rk4::atmos_dyn_tinteg_short_rk4
subroutine, public atmos_dyn_tinteg_short_rk4(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_rk4.F90:275
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_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4_setup
subroutine, public atmos_dyn_tinteg_short_rk4_setup(tinteg_type)
Setup.
Definition: scale_atmos_dyn_tinteg_short_rk4.F90:108
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_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_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_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_atmos_dyn_tinteg_short_rk4
module Atmosphere / Dyn Tinteg
Definition: scale_atmos_dyn_tinteg_short_rk4.F90:20
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_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4_finalize
subroutine, public atmos_dyn_tinteg_short_rk4_finalize
finalize
Definition: scale_atmos_dyn_tinteg_short_rk4.F90:227
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