SCALE-RM
Functions/Subroutines
scale_atmos_dyn_tinteg_short_rk4 Module Reference

module Atmosphere / Dyn Tinteg More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_short_rk4_setup (tinteg_type)
 Setup. More...
 
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, Rtot, CVtot, CORIOLI, num_diff, wdamp_coef, divdmp_coef, DDIV, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_ALONG_STREAM, CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, REF_pres, REF_dens, BND_W, BND_E, BND_S, BND_N, dt)
 RK3. More...
 

Detailed Description

module Atmosphere / Dyn Tinteg

Description
Temporal integration in Dynamical core for Atmospheric process four step Runge-Kutta scheme
Author
Team SCALE
History
  • 2016-04-30 (S.Nishizawa) [new]

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_short_rk4_setup()

subroutine, public scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4_setup ( character(len=*)  tinteg_type)

Setup.

Definition at line 103 of file scale_atmos_dyn_tinteg_short_rk4.F90.

References scale_comm::comm_vars8_init(), scale_const::const_undef, scale_grid_index::ia, scale_grid_index::ja, scale_grid_index::ka, scale_process::prc_mpistop(), and scale_index::va.

Referenced by scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup().

103  use scale_process, only: &
105  use scale_const, only: &
106  undef => const_undef
107  use scale_comm, only: &
109  implicit none
110 
111  character(len=*) :: tinteg_type
112 
113  integer :: iv
114  !---------------------------------------------------------------------------
115 
116  if ( tinteg_type /= 'RK4' ) then
117  write(*,*) 'xxx TINTEG_TYPE is not RK4. Check!'
118  call prc_mpistop
119  end if
120 
121  allocate( dens_rk1(ka,ia,ja) )
122  allocate( momz_rk1(ka,ia,ja) )
123  allocate( momx_rk1(ka,ia,ja) )
124  allocate( momy_rk1(ka,ia,ja) )
125  allocate( rhot_rk1(ka,ia,ja) )
126 
127  allocate( dens_rk2(ka,ia,ja) )
128  allocate( momz_rk2(ka,ia,ja) )
129  allocate( momx_rk2(ka,ia,ja) )
130  allocate( momy_rk2(ka,ia,ja) )
131  allocate( rhot_rk2(ka,ia,ja) )
132 
133  allocate( dens_rk3(ka,ia,ja) )
134  allocate( momz_rk3(ka,ia,ja) )
135  allocate( momx_rk3(ka,ia,ja) )
136  allocate( momy_rk3(ka,ia,ja) )
137  allocate( rhot_rk3(ka,ia,ja) )
138 
139  allocate( prog_rk1(ka,ia,ja,max(va,1)) )
140  allocate( prog_rk2(ka,ia,ja,max(va,1)) )
141  allocate( prog_rk3(ka,ia,ja,max(va,1)) )
142  allocate( i_comm_prog_rk1(max(va,1)) )
143  allocate( i_comm_prog_rk2(max(va,1)) )
144  allocate( i_comm_prog_rk3(max(va,1)) )
145 
146  call comm_vars8_init( 'DENS_RK1', dens_rk1, i_comm_dens_rk1 )
147  call comm_vars8_init( 'MOMZ_RK1', momz_rk1, i_comm_momz_rk1 )
148  call comm_vars8_init( 'MOMX_RK1', momx_rk1, i_comm_momx_rk1 )
149  call comm_vars8_init( 'MOMY_RK1', momy_rk1, i_comm_momy_rk1 )
150  call comm_vars8_init( 'RHOT_RK1', rhot_rk1, i_comm_rhot_rk1 )
151  do iv = 1, va
152  i_comm_prog_rk1(iv) = 5 + iv
153  call comm_vars8_init( 'PROG_RK1', prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
154  enddo
155 
156  call comm_vars8_init( 'DENS_RK2', dens_rk2, i_comm_dens_rk2 )
157  call comm_vars8_init( 'MOMZ_RK2', momz_rk2, i_comm_momz_rk2 )
158  call comm_vars8_init( 'MOMX_RK2', momx_rk2, i_comm_momx_rk2 )
159  call comm_vars8_init( 'MOMY_RK2', momy_rk2, i_comm_momy_rk2 )
160  call comm_vars8_init( 'RHOT_RK2', rhot_rk2, i_comm_rhot_rk2 )
161  do iv = 1, va
162  i_comm_prog_rk2(iv) = 5 + iv
163  call comm_vars8_init( 'PROG_RK2', prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
164  enddo
165 
166  call comm_vars8_init( 'DENS_RK3', dens_rk3, i_comm_dens_rk3 )
167  call comm_vars8_init( 'MOMZ_RK3', momz_rk3, i_comm_momz_rk3 )
168  call comm_vars8_init( 'MOMX_RK3', momx_rk3, i_comm_momx_rk3 )
169  call comm_vars8_init( 'MOMY_RK3', momy_rk3, i_comm_momy_rk3 )
170  call comm_vars8_init( 'RHOT_RK3', rhot_rk3, i_comm_rhot_rk3 )
171  do iv = 1, va
172  i_comm_prog_rk3(iv) = 5 + iv
173  call comm_vars8_init( 'PROG_RK3', prog_rk3(:,:,:,iv), i_comm_prog_rk3(iv) )
174  enddo
175 
176  dens_rk1(:,:,:) = undef
177  momz_rk1(:,:,:) = undef
178  momx_rk1(:,:,:) = undef
179  momy_rk1(:,:,:) = undef
180  rhot_rk1(:,:,:) = undef
181  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
182 
183  dens_rk2(:,:,:) = undef
184  momz_rk2(:,:,:) = undef
185  momx_rk2(:,:,:) = undef
186  momy_rk2(:,:,:) = undef
187  rhot_rk2(:,:,:) = undef
188  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
189 
190  dens_rk3(:,:,:) = undef
191  momz_rk3(:,:,:) = undef
192  momx_rk3(:,:,:) = undef
193  momy_rk3(:,:,:) = undef
194  rhot_rk3(:,:,:) = undef
195  if ( va > 0 ) prog_rk3(:,:,:,:) = undef
196 
197  return
subroutine, public prc_mpistop
Abort MPI.
integer, public va
Definition: scale_index.F90:38
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm.F90:313
real(rp), public const_undef
Definition: scale_const.F90:43
integer, public ia
of whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public ja
of whole cells: y, local, with HALO
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tinteg_short_rk4()

subroutine, public scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4 ( real(rp), dimension(ka,ia,ja), intent(inout)  DENS,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMX,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMY,
real(rp), dimension(ka,ia,ja), intent(inout)  RHOT,
real(rp), dimension(ka,ia,ja,va), intent(inout)  PROG,
real(rp), dimension(ka,ia,ja,3), intent(inout)  mflx_hi,
real(rp), dimension(ka,ia,ja,3), intent(inout)  tflx_hi,
real(rp), dimension(ka,ia,ja), intent(in)  DENS_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX_t,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY_t,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT_t,
real(rp), dimension(ka,ia,ja), intent(in)  Rtot,
real(rp), dimension(ka,ia,ja), intent(in)  CVtot,
real(rp), dimension(ia,ja), intent(in)  CORIOLI,
real(rp), dimension(ka,ia,ja,5,3), intent(in)  num_diff,
real(rp), dimension(ka), intent(in)  wdamp_coef,
real(rp), intent(in)  divdmp_coef,
real(rp), dimension(ka,ia,ja), intent(in)  DDIV,
logical, intent(in)  FLAG_FCT_MOMENTUM,
logical, intent(in)  FLAG_FCT_T,
logical, intent(in)  FLAG_FCT_ALONG_STREAM,
real(rp), dimension (ka), intent(in)  CDZ,
real(rp), dimension (ka-1), intent(in)  FDZ,
real(rp), dimension (ia-1), intent(in)  FDX,
real(rp), dimension (ja-1), intent(in)  FDY,
real(rp), dimension(ka), intent(in)  RCDZ,
real(rp), dimension(ia), intent(in)  RCDX,
real(rp), dimension(ja), intent(in)  RCDY,
real(rp), dimension(ka-1), intent(in)  RFDZ,
real(rp), dimension(ia-1), intent(in)  RFDX,
real(rp), dimension(ja-1), intent(in)  RFDY,
real(rp), dimension (ka,ia,ja), intent(in)  PHI,
real(rp), dimension(ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension (ka,ia,ja,7), intent(in)  J13G,
real(rp), dimension (ka,ia,ja,7), intent(in)  J23G,
real(rp), intent(in)  J33G,
real(rp), dimension (ia,ja,2,4), intent(in)  MAPF,
real(rp), dimension(ka,ia,ja), intent(in)  REF_pres,
real(rp), dimension(ka,ia,ja), intent(in)  REF_dens,
logical, intent(in)  BND_W,
logical, intent(in)  BND_E,
logical, intent(in)  BND_S,
logical, intent(in)  BND_N,
real(rp), intent(in)  dt 
)

RK3.

Parameters
[in]phigeopotential
[in]gsqrtvertical metrics {G}^1/2
[in]j13g(1,3) element of Jacobian matrix
[in]j23g(2,3) element of Jacobian matrix
[in]j33g(3,3) element of Jacobian matrix
[in]mapfmap factor
[in]ref_presreference pressure

Definition at line 216 of file scale_atmos_dyn_tinteg_short_rk4.F90.

References scale_atmos_dyn_common::atmos_dyn_copy_boundary(), scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_index::va.

Referenced by scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup().

216  use scale_comm, only: &
217  comm_vars8, &
218  comm_wait
219  use scale_atmos_dyn_tstep_short, only: &
220  atmos_dyn_tstep => atmos_dyn_tstep_short
221  use scale_atmos_dyn_common, only: &
223  implicit none
224 
225  real(RP), intent(inout) :: DENS(KA,IA,JA)
226  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
227  real(RP), intent(inout) :: MOMX(KA,IA,JA)
228  real(RP), intent(inout) :: MOMY(KA,IA,JA)
229  real(RP), intent(inout) :: RHOT(KA,IA,JA)
230  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
231 
232  real(RP), intent(inout) :: mflx_hi(KA,IA,JA,3)
233  real(RP), intent(inout) :: tflx_hi(KA,IA,JA,3)
234 
235  real(RP), intent(in) :: DENS_t(KA,IA,JA)
236  real(RP), intent(in) :: MOMZ_t(KA,IA,JA)
237  real(RP), intent(in) :: MOMX_t(KA,IA,JA)
238  real(RP), intent(in) :: MOMY_t(KA,IA,JA)
239  real(RP), intent(in) :: RHOT_t(KA,IA,JA)
240 
241  real(RP), intent(in) :: Rtot(KA,IA,JA)
242  real(RP), intent(in) :: CVtot(KA,IA,JA)
243  real(RP), intent(in) :: CORIOLI(IA,JA)
244  real(RP), intent(in) :: num_diff(KA,IA,JA,5,3)
245  real(RP), intent(in) :: wdamp_coef(KA)
246  real(RP), intent(in) :: divdmp_coef
247  real(RP), intent(in) :: DDIV(KA,IA,JA)
248 
249  logical, intent(in) :: FLAG_FCT_MOMENTUM
250  logical, intent(in) :: FLAG_FCT_T
251  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
252 
253  real(RP), intent(in) :: CDZ (KA)
254  real(RP), intent(in) :: FDZ (KA-1)
255  real(RP), intent(in) :: FDX (IA-1)
256  real(RP), intent(in) :: FDY (JA-1)
257  real(RP), intent(in) :: RCDZ(KA)
258  real(RP), intent(in) :: RCDX(IA)
259  real(RP), intent(in) :: RCDY(JA)
260  real(RP), intent(in) :: RFDZ(KA-1)
261  real(RP), intent(in) :: RFDX(IA-1)
262  real(RP), intent(in) :: RFDY(JA-1)
263 
264  real(RP), intent(in) :: PHI (KA,IA,JA)
265  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
266  real(RP), intent(in) :: J13G (KA,IA,JA,7)
267  real(RP), intent(in) :: J23G (KA,IA,JA,7)
268  real(RP), intent(in) :: J33G
269  real(RP), intent(in) :: MAPF (IA,JA,2,4)
270 
271  real(RP), intent(in) :: REF_pres(KA,IA,JA)
272  real(RP), intent(in) :: REF_dens(KA,IA,JA)
273 
274  logical, intent(in) :: BND_W
275  logical, intent(in) :: BND_E
276  logical, intent(in) :: BND_S
277  logical, intent(in) :: BND_N
278 
279  real(RP), intent(in) :: dt
280 
281  real(RP) :: DENS0(KA,IA,JA)
282  real(RP) :: MOMZ0(KA,IA,JA)
283  real(RP) :: MOMX0(KA,IA,JA)
284  real(RP) :: MOMY0(KA,IA,JA)
285  real(RP) :: RHOT0(KA,IA,JA)
286  real(RP) :: PROG0(KA,IA,JA,VA)
287 
288  real(RP) :: mflx_hi_RK(KA,IA,JA,3,3)
289  real(RP) :: tflx_hi_RK(KA,IA,JA,3,3)
290 
291  real(RP) :: dtrk
292 
293  integer :: i, j, k, iv, n
294  !---------------------------------------------------------------------------
295 
296  call prof_rapstart("DYN_RK4_Prep",3)
297 
298 #ifdef DEBUG
299  dens_rk1(:,:,:) = undef
300  momz_rk1(:,:,:) = undef
301  momx_rk1(:,:,:) = undef
302  momy_rk1(:,:,:) = undef
303  rhot_rk1(:,:,:) = undef
304  if ( va > 0 ) prog_rk1(:,:,:,:) = undef
305 
306  dens_rk2(:,:,:) = undef
307  momz_rk2(:,:,:) = undef
308  momx_rk2(:,:,:) = undef
309  momy_rk2(:,:,:) = undef
310  rhot_rk2(:,:,:) = undef
311  if ( va > 0 ) prog_rk2(:,:,:,:) = undef
312 
313  dens_rk3(:,:,:) = undef
314  momz_rk3(:,:,:) = undef
315  momx_rk3(:,:,:) = undef
316  momy_rk3(:,:,:) = undef
317  rhot_rk3(:,:,:) = undef
318  if ( va > 0 ) prog_rk3(:,:,:,:) = undef
319 
320  mflx_hi_rk(:,:,:,:,:) = undef
321  tflx_hi_rk(:,:,:,:,:) = undef
322 #endif
323 
324 #ifdef QUICKDEBUG
325  mflx_hi( 1:ks-1,:,:,:) = undef
326  mflx_hi(ke+1:ka ,:,:,:) = undef
327 #endif
328 
329 !OCL XFILL
330  dens0 = dens
331 !OCL XFILL
332  momz0 = momz
333 !OCL XFILL
334  momx0 = momx
335 !OCL XFILL
336  momy0 = momy
337 !OCL XFILL
338  rhot0 = rhot
339 !OCL XFILL
340  if ( va > 0 ) prog0 = prog
341 
342  if ( bnd_w ) then
343  do j = js, je
344  do k = ks, ke
345  mflx_hi_rk(k,is-1,j,2,:) = mflx_hi(k,is-1,j,2)
346  end do
347  end do
348  end if
349  if ( bnd_e ) then
350  do j = js, je
351  do k = ks, ke
352  mflx_hi_rk(k,ie,j,2,:) = mflx_hi(k,ie,j,2)
353  end do
354  end do
355  end if
356  if ( bnd_s ) then
357  do i = is, ie
358  do k = ks, ke
359  mflx_hi_rk(k,i,js-1,3,:) = mflx_hi(k,i,js-1,3)
360  end do
361  end do
362  end if
363  if ( bnd_n ) then
364  do i = is, ie
365  do k = ks, ke
366  mflx_hi_rk(k,i,je,3,:) = mflx_hi(k,i,je,3)
367  end do
368  end do
369  end if
370 
371  call prof_rapend ("DYN_RK4_Prep",3)
372 
373  !------------------------------------------------------------------------
374  ! Start RK
375  !------------------------------------------------------------------------
376 
377  !##### RK1 : PROG0,PROG->PROG_RK1 #####
378 
379  call prof_rapstart("DYN_RK4",3)
380 
381  dtrk = dt / 2.0_rp
382 
383  call atmos_dyn_tstep( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [OUT]
384  prog_rk1, & ! [OUT]
385  mflx_hi_rk(:,:,:,:,1), tflx_hi_rk(:,:,:,:,1), & ! [INOUT]
386  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
387  dens, momz, momx, momy, rhot, & ! [IN]
388  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
389  prog0, prog, & ! [IN]
390  rtot, cvtot, corioli, & ! [IN]
391  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
392  flag_fct_momentum, flag_fct_t, & ! [IN]
393  flag_fct_along_stream, & ! [IN]
394  cdz, fdz, fdx, fdy, & ! [IN]
395  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
396  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
397  ref_pres, ref_dens, & ! [IN]
398  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
399  dtrk, .false. ) ! [IN]
400 
401  call prof_rapend ("DYN_RK4",3)
402  call prof_rapstart("DYN_RK4_BND",3)
403 
404  call atmos_dyn_copy_boundary( dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [INOUT]
405  prog_rk1, & ! [INOUT]
406  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
407  prog0, & ! [IN]
408  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
409 
410  call prof_rapend ("DYN_RK4_BND",3)
411 
412  call comm_vars8( dens_rk1(:,:,:), i_comm_dens_rk1 )
413  call comm_vars8( momz_rk1(:,:,:), i_comm_momz_rk1 )
414  call comm_vars8( momx_rk1(:,:,:), i_comm_momx_rk1 )
415  call comm_vars8( momy_rk1(:,:,:), i_comm_momy_rk1 )
416  call comm_vars8( rhot_rk1(:,:,:), i_comm_rhot_rk1 )
417  do iv = 1, va
418  call comm_vars8( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv) )
419  enddo
420 
421  call comm_wait ( dens_rk1(:,:,:), i_comm_dens_rk1, .false. )
422  call comm_wait ( momz_rk1(:,:,:), i_comm_momz_rk1, .false. )
423  call comm_wait ( momx_rk1(:,:,:), i_comm_momx_rk1, .false. )
424  call comm_wait ( momy_rk1(:,:,:), i_comm_momy_rk1, .false. )
425  call comm_wait ( rhot_rk1(:,:,:), i_comm_rhot_rk1, .false. )
426  do iv = 1, va
427  call comm_wait ( prog_rk1(:,:,:,iv), i_comm_prog_rk1(iv), .false. )
428  enddo
429 
430  !##### RK2 : PROG0,PROG_RK1->PROG_RK2 #####
431 
432  call prof_rapstart("DYN_RK4",3)
433 
434  dtrk = dt / 2.0_rp
435 
436  call atmos_dyn_tstep( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [OUT]
437  prog_rk2, & ! [OUT]
438  mflx_hi_rk(:,:,:,:,2), tflx_hi_rk(:,:,:,:,2), & ! [INOUT]
439  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
440  dens_rk1, momz_rk1, momx_rk1, momy_rk1, rhot_rk1, & ! [IN]
441  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
442  prog0, prog_rk1, & ! [IN]
443  rtot, cvtot, corioli, & ! [IN]
444  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
445  flag_fct_momentum, flag_fct_t, & ! [IN]
446  flag_fct_along_stream, & ! [IN]
447  cdz, fdz, fdx, fdy, & ! [IN]
448  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
449  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
450  ref_pres, ref_dens, & ! [IN]
451  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
452  dtrk, .false. ) ! [IN]
453 
454  call prof_rapend ("DYN_RK4",3)
455  call prof_rapstart("DYN_RK4_BND",3)
456 
457  call atmos_dyn_copy_boundary( dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [INOUT]
458  prog_rk2, & ! [INOUT]
459  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
460  prog0, & ! [IN]
461  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
462 
463  call prof_rapend ("DYN_RK4_BND",3)
464 
465  call comm_vars8( dens_rk2(:,:,:), i_comm_dens_rk2 )
466  call comm_vars8( momz_rk2(:,:,:), i_comm_momz_rk2 )
467  call comm_vars8( momx_rk2(:,:,:), i_comm_momx_rk2 )
468  call comm_vars8( momy_rk2(:,:,:), i_comm_momy_rk2 )
469  call comm_vars8( rhot_rk2(:,:,:), i_comm_rhot_rk2 )
470  do iv = 1, va
471  call comm_vars8( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv) )
472  enddo
473 
474  call comm_wait ( dens_rk2(:,:,:), i_comm_dens_rk2, .false. )
475  call comm_wait ( momz_rk2(:,:,:), i_comm_momz_rk2, .false. )
476  call comm_wait ( momx_rk2(:,:,:), i_comm_momx_rk2, .false. )
477  call comm_wait ( momy_rk2(:,:,:), i_comm_momy_rk2, .false. )
478  call comm_wait ( rhot_rk2(:,:,:), i_comm_rhot_rk2, .false. )
479  do iv = 1, va
480  call comm_wait ( prog_rk2(:,:,:,iv), i_comm_prog_rk2(iv), .false. )
481  enddo
482 
483  !##### RK3 : PROG0,PROG_RK2->PROG_RK3 #####
484 
485  call prof_rapstart("DYN_RK4",3)
486 
487  dtrk = dt
488 
489  call atmos_dyn_tstep( dens_rk3, momz_rk3, momx_rk3, momy_rk3, rhot_rk3, & ! [OUT]
490  prog_rk3, & ! [OUT]
491  mflx_hi_rk(:,:,:,:,3), tflx_hi_rk(:,:,:,:,3), & ! [INOUT]
492  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
493  dens_rk2, momz_rk2, momx_rk2, momy_rk2, rhot_rk2, & ! [IN]
494  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
495  prog0, prog_rk2, & ! [IN]
496  rtot, cvtot, corioli, & ! [IN]
497  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
498  flag_fct_momentum, flag_fct_t, & ! [IN]
499  flag_fct_along_stream, & ! [IN]
500  cdz, fdz, fdx, fdy, & ! [IN]
501  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
502  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
503  ref_pres, ref_dens, & ! [IN]
504  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
505  dtrk, .false. ) ! [IN]
506 
507  call prof_rapend ("DYN_RK4",3)
508  call prof_rapstart("DYN_RK4_BND",3)
509 
510  call atmos_dyn_copy_boundary( dens_rk3, momz_rk3, momx_rk3, momy_rk3, rhot_rk3, & ! [INOUT]
511  prog_rk3, & ! [INOUT]
512  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
513  prog0, & ! [IN]
514  bnd_w, bnd_e, bnd_s, bnd_n ) ! [IN]
515 
516  call prof_rapend ("DYN_RK4_BND",3)
517 
518  call comm_vars8( dens_rk3(:,:,:), i_comm_dens_rk3 )
519  call comm_vars8( momz_rk3(:,:,:), i_comm_momz_rk3 )
520  call comm_vars8( momx_rk3(:,:,:), i_comm_momx_rk3 )
521  call comm_vars8( momy_rk3(:,:,:), i_comm_momy_rk3 )
522  call comm_vars8( rhot_rk3(:,:,:), i_comm_rhot_rk3 )
523  do iv = 1, va
524  call comm_vars8( prog_rk3(:,:,:,iv), i_comm_prog_rk3(iv) )
525  enddo
526 
527  call comm_wait ( dens_rk3(:,:,:), i_comm_dens_rk3, .false. )
528  call comm_wait ( momz_rk3(:,:,:), i_comm_momz_rk3, .false. )
529  call comm_wait ( momx_rk3(:,:,:), i_comm_momx_rk3, .false. )
530  call comm_wait ( momy_rk3(:,:,:), i_comm_momy_rk3, .false. )
531  call comm_wait ( rhot_rk3(:,:,:), i_comm_rhot_rk3, .false. )
532  do iv = 1, va
533  call comm_wait ( prog_rk3(:,:,:,iv), i_comm_prog_rk3(iv), .false. )
534  enddo
535 
536  !##### RK4 : PROG0,PROG_RK3->PROG #####
537 
538  call prof_rapstart("DYN_RK4",3)
539 
540  dtrk = dt
541 
542  call atmos_dyn_tstep( dens, momz, momx, momy, rhot, & ! [OUT]
543  prog, & ! [OUT]
544  mflx_hi, tflx_hi, & ! [INOUT]
545  dens0, momz0, momx0, momy0, rhot0, & ! [IN]
546  dens_rk3, momz_rk3, momx_rk3, momy_rk3, rhot_rk3, & ! [IN]
547  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
548  prog0, prog_rk3, & ! [IN]
549  rtot, cvtot, corioli, & ! [IN]
550  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
551  flag_fct_momentum, flag_fct_t, & ! [IN]
552  flag_fct_along_stream, & ! [IN]
553  cdz, fdz, fdx, fdy, & ! [IN]
554  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
555  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
556  ref_pres, ref_dens, & ! [IN]
557  bnd_w, bnd_e, bnd_s, bnd_n, & ! [IN]
558  dtrk, .true. ) ! [IN]
559 
560  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
561  !$omp shared(JS,JE,IS,IE,KS,KE,DENS,DENS_RK1,DENS_RK2,DENS_RK3,DENS0)
562  do j = js, je
563  do i = is, ie
564  do k = ks, ke
565  dens(k,i,j) = ( dens_rk1(k,i,j) * 2.0_rp &
566  + dens_rk2(k,i,j) * 4.0_rp &
567  + dens_rk3(k,i,j) * 2.0_rp &
568  + dens(k,i,j) &
569  - dens0(k,i,j) * 3.0_rp ) / 6.0_rp
570  enddo
571  enddo
572  enddo
573 
574  do j = js, je
575  do i = is, ie
576  do k = ks, ke-1
577  momz(k,i,j) = ( momz_rk1(k,i,j) * 2.0_rp &
578  + momz_rk2(k,i,j) * 4.0_rp &
579  + momz_rk3(k,i,j) * 2.0_rp &
580  + momz(k,i,j) &
581  - momz0(k,i,j) * 3.0_rp ) / 6.0_rp
582  enddo
583  enddo
584  enddo
585 
586  do j = js, je
587  do i = is, ie
588  do k = ks, ke
589  momx(k,i,j) = ( momx_rk1(k,i,j) * 2.0_rp &
590  + momx_rk2(k,i,j) * 4.0_rp &
591  + momx_rk3(k,i,j) * 2.0_rp &
592  + momx(k,i,j) &
593  - momx0(k,i,j) * 3.0_rp ) / 6.0_rp
594  enddo
595  enddo
596  enddo
597 
598  do j = js, je
599  do i = is, ie
600  do k = ks, ke
601  momy(k,i,j) = ( momy_rk1(k,i,j) * 2.0_rp &
602  + momy_rk2(k,i,j) * 4.0_rp &
603  + momy_rk3(k,i,j) * 2.0_rp &
604  + momy(k,i,j) &
605  - momy0(k,i,j) * 3.0_rp ) / 6.0_rp
606  enddo
607  enddo
608  enddo
609 
610  !$omp parallel do default(none) &
611  !$omp shared(JS,JE,IS,IE,KS,KE,RHOT,RHOT_RK1,RHOT_RK2,RHOT_RK3,RHOT0) &
612  !$omp private(i,j,k) OMP_SCHEDULE_ collapse(2)
613  do j = js, je
614  do i = is, ie
615  do k = ks, ke
616  rhot(k,i,j) = ( rhot_rk1(k,i,j) * 2.0_rp &
617  + rhot_rk2(k,i,j) * 4.0_rp &
618  + rhot_rk3(k,i,j) * 2.0_rp &
619  + rhot(k,i,j) &
620  - rhot0(k,i,j) * 3.0_rp ) / 6.0_rp
621  enddo
622  enddo
623  enddo
624 
625  do iv = 1, va
626  do j = js, je
627  do i = is, ie
628  do k = ks, ke
629  prog(k,i,j,iv) = ( prog_rk1(k,i,j,iv) * 2.0_rp &
630  + prog_rk2(k,i,j,iv) * 4.0_rp &
631  + prog_rk3(k,i,j,iv) * 2.0_rp &
632  + prog(k,i,j,iv) &
633  - prog0(k,i,j,iv) * 3.0_rp ) / 6.0_rp
634  enddo
635  enddo
636  enddo
637  enddo
638 
639 
640  do n = 1, 3
641  do j = js, je
642  do i = is, ie
643  do k = ks, ke
644  mflx_hi(k,i,j,n) = ( mflx_hi_rk(k,i,j,n,1) &
645  + mflx_hi_rk(k,i,j,n,2) * 2.0_rp &
646  + mflx_hi_rk(k,i,j,n,3) * 2.0_rp &
647  + mflx_hi(k,i,j,n ) ) / 6.0_rp
648  enddo
649  enddo
650  enddo
651  enddo
652 
653  do n = 1, 3
654  do j = js, je
655  do i = is, ie
656  do k = ks, ke
657  tflx_hi(k,i,j,n) = ( tflx_hi_rk(k,i,j,n,1) &
658  + tflx_hi_rk(k,i,j,n,2) * 2.0_rp &
659  + tflx_hi_rk(k,i,j,n,3) * 2.0_rp &
660  + tflx_hi(k,i,j,n ) ) / 6.0_rp
661  enddo
662  enddo
663  enddo
664  enddo
665 
666  call prof_rapend ("DYN_RK4",3)
667 
668  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public va
Definition: scale_index.F90:38
module Atmosphere / Dynamical scheme
integer, public ke
end point of inner domain: z, local
integer, public ka
of whole cells: z, local, with HALO
procedure(short), pointer, public atmos_dyn_tstep_short
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
module Atmosphere / Dynamics common
integer, public ks
start point of inner domain: z, local
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)
integer, public ie
end point of inner domain: x, local
Here is the call graph for this function:
Here is the caller graph for this function: