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