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