SCALE-RM
Functions/Subroutines | Variables
scale_atmos_dyn_tinteg_short_rk11s8o Module Reference

module Atmosphere / Dyn Tinteg More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_short_rk11s8o_setup (tinteg_type)
 Setup. More...
 
subroutine, public atmos_dyn_tinteg_short_rk11s8o (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)
 RK7s6o. More...
 

Variables

real(rp), dimension(11, 11), parameter, public rkcoef_a_11s8o_cooperverner1972 = reshape( (/ ZERO, 0.5_RP, 0.25_RP, RC41, RC51, RC61, RC71, RC81, RC91, RC101, ZERO, ZERO, ZERO, 0.25_RP, RC42, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC43, RC53, RC63, RC73, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC54, RC64, RC74, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC65, RC75, RC85, RC95, RC105, RC115, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC76, RC86, RC96, RC106, RC116, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC87, RC97, RC107, RC117, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC98, RC108, RC118, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC109, RC119, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC1110, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO /), shape(RKCoef_a_11s8o_CooperVerner1972) )
 
real(rp), dimension(11), parameter, public rkcoef_b_11s8o_cooperverner1972 = (/ 0.05_RP, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, 49.0_RP/180.0_RP, 16.0_RP/45.0_RP, 49.0_RP/180.0_RP, 0.05_RP /)
 

Detailed Description

module Atmosphere / Dyn Tinteg

Description
Temporal integration in Dynamical core for Atmospheric process 11 stage Runge-Kutta scheme with 8th order accuracy
Author
Team SCALE

This module provides two type of 11 stage Runge-Kutta scheme with 8th order accuracy proposed by Cooper and Verner (1972). See this source file for the detail of RK coffecients.

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_short_rk11s8o_setup()

subroutine, public scale_atmos_dyn_tinteg_short_rk11s8o::atmos_dyn_tinteg_short_rk11s8o_setup ( character(len=*)  tinteg_type)

Setup.

Definition at line 149 of file scale_atmos_dyn_tinteg_short_rk11s8o.F90.

149  use scale_prc, only: &
150  prc_abort
151  use scale_const, only: &
152  undef => const_undef
153  use scale_comm_cartesc, only: &
155  implicit none
156 
157  character(len=*) :: tinteg_type
158 
159  integer :: iv, d
160  character(H_SHORT) :: dynvar_name_list(5)
161  character(H_SHORT) :: prgvar_name_list(VA)
162  character(H_SHORT) :: flux_name_list(3)
163 
164  real(RP) :: RKCoef_a(RK_nstage,RK_nstage)
165  real(RP) :: RKCoef_b(RK_nstage)
166 
167  !---------------------------------------------------------------------------
168 
169  select case( trim(tinteg_type) )
170  case ('RK11s8o', 'RK11s8oCooperVerner1972')
171  rkcoef_a(:,:) = rkcoef_a_11s8o_cooperverner1972
172  rkcoef_b(:) = rkcoef_b_11s8o_cooperverner1972
173  case default
174  log_error("ATMOS_DYN_Tinteg_short_rk11s8o_setup",*) 'The specified TINTEG_TYPE is invalid. Check!', tinteg_type
175  call prc_abort
176  end select
177 
178  dynvar_name_list(1) = 'DENS'
179  dynvar_name_list(2) = 'MOMZ'
180  dynvar_name_list(3) = 'MOMX'
181  dynvar_name_list(4) = 'MOMY'
182  dynvar_name_list(5) = 'RHOT'
183  call rkcommon_setup( rk_dynvar, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, dynvar_name_list )
184 
185  do iv = 1, va
186  flux_name_list(iv) = 'PROG'
187  end do
188  call rkcommon_setup( rk_prgvar, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b,prgvar_name_list, comm_id_offset=5 )
189 
190  do d = 1, 3
191  flux_name_list(iv) = 'mflx_hi'
192  end do
193  call rkcommon_setup( rk_mflx_hi, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, flux_name_list, is_type_flux=.true. )
194 
195 
196  do d = 1, 3
197  flux_name_list(iv) = 'tflx_hi'
198  end do
199  call rkcommon_setup( rk_tflx_hi, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, flux_name_list, is_type_flux=.true. )
200 
201  !----------------------------------------------------
202 
203  return

References scale_comm_cartesc::comm_vars8_init(), scale_const::const_undef, scale_prc::prc_abort(), rkcoef_a_11s8o_cooperverner1972, rkcoef_b_11s8o_cooperverner1972, and scale_index::va.

Referenced by scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_dyn_tinteg_short_rk11s8o()

subroutine, public scale_atmos_dyn_tinteg_short_rk11s8o::atmos_dyn_tinteg_short_rk11s8o ( 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(out)  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)  DPRES0,
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,
logical, intent(in)  TwoD,
real(rp), intent(in)  dt 
)

RK7s6o.

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 222 of file scale_atmos_dyn_tinteg_short_rk11s8o.F90.

222  use scale_comm_cartesc, only: &
223  comm_vars8, &
224  comm_wait
225  use scale_atmos_dyn_tstep_short, only: &
226  atmos_dyn_tstep => atmos_dyn_tstep_short
227  use scale_atmos_dyn_common, only: &
229  implicit none
230 
231  real(RP), intent(inout) :: DENS(KA,IA,JA)
232  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
233  real(RP), intent(inout) :: MOMX(KA,IA,JA)
234  real(RP), intent(inout) :: MOMY(KA,IA,JA)
235  real(RP), intent(inout) :: RHOT(KA,IA,JA)
236  real(RP), intent(inout) :: PROG(KA,IA,JA,VA)
237 
238  real(RP), intent(inout) :: mflx_hi(KA,IA,JA,3)
239  real(RP), intent(out) :: tflx_hi(KA,IA,JA,3)
240 
241  real(RP), intent(in) :: DENS_t(KA,IA,JA)
242  real(RP), intent(in) :: MOMZ_t(KA,IA,JA)
243  real(RP), intent(in) :: MOMX_t(KA,IA,JA)
244  real(RP), intent(in) :: MOMY_t(KA,IA,JA)
245  real(RP), intent(in) :: RHOT_t(KA,IA,JA)
246 
247  real(RP), intent(in) :: DPRES0(KA,IA,JA)
248  real(RP), intent(in) :: CVtot(KA,IA,JA)
249  real(RP), intent(in) :: CORIOLI(IA,JA)
250  real(RP), intent(in) :: num_diff(KA,IA,JA,5,3)
251  real(RP), intent(in) :: wdamp_coef(KA)
252  real(RP), intent(in) :: divdmp_coef
253  real(RP), intent(in) :: DDIV(KA,IA,JA)
254 
255  logical, intent(in) :: FLAG_FCT_MOMENTUM
256  logical, intent(in) :: FLAG_FCT_T
257  logical, intent(in) :: FLAG_FCT_ALONG_STREAM
258 
259  real(RP), intent(in) :: CDZ (KA)
260  real(RP), intent(in) :: FDZ (KA-1)
261  real(RP), intent(in) :: FDX (IA-1)
262  real(RP), intent(in) :: FDY (JA-1)
263  real(RP), intent(in) :: RCDZ(KA)
264  real(RP), intent(in) :: RCDX(IA)
265  real(RP), intent(in) :: RCDY(JA)
266  real(RP), intent(in) :: RFDZ(KA-1)
267  real(RP), intent(in) :: RFDX(IA-1)
268  real(RP), intent(in) :: RFDY(JA-1)
269 
270  real(RP), intent(in) :: PHI (KA,IA,JA)
271  real(RP), intent(in) :: GSQRT(KA,IA,JA,7)
272  real(RP), intent(in) :: J13G (KA,IA,JA,7)
273  real(RP), intent(in) :: J23G (KA,IA,JA,7)
274  real(RP), intent(in) :: J33G
275  real(RP), intent(in) :: MAPF (IA,JA,2,4)
276 
277  real(RP), intent(in) :: REF_pres(KA,IA,JA)
278  real(RP), intent(in) :: REF_dens(KA,IA,JA)
279 
280  logical, intent(in) :: BND_W
281  logical, intent(in) :: BND_E
282  logical, intent(in) :: BND_S
283  logical, intent(in) :: BND_N
284  logical, intent(in) :: TwoD
285 
286  real(RP), intent(in) :: dt
287 
288  integer :: i, j, k, iv, n, s
289  integer :: stage
290 
291  integer, parameter :: ko_dynvar(5) = (/ 0, 1, 0, 0, 0 /)
292  integer, parameter :: io_dynvar(5) = (/ 0, 0, 1, 0, 0 /)
293  integer, parameter :: jo_dynvar(5) = (/ 0, 0, 0, 1, 0 /)
294  integer :: ko_prgvar(VA)
295  integer :: io_prgvar(VA)
296  integer :: jo_prgvar(VA)
297 
298  !---------------------------------------------------------------------------
299 
300  call prof_rapstart("DYN_RK11s8o_Prep",3)
301 
302  call rkcommon_rkwork_alloc( rk_dynvar )
303  call rkcommon_rkwork_alloc( rk_prgvar )
304 
305 #ifdef QUICKDEBUG
306  rk_mflx_hi%buf( 1:ks-1,:,:,:) = undef
307  rk_mflx_hi%buf(ke+1:ka ,:,:,:) = undef
308 #endif
309 
310  !$omp parallel
311 
312  !$omp workshare
313 !OCL XFILL
314  rk_dynvar%work0(:,:,:,i_rk_dens) = dens(:,:,:)
315 !OCL XFILL
316  rk_dynvar%work0(:,:,:,i_rk_momz) = momz(:,:,:)
317 !OCL XFILL
318  rk_dynvar%work0(:,:,:,i_rk_momx) = momx(:,:,:)
319 !OCL XFILL
320  rk_dynvar%work0(:,:,:,i_rk_momy) = momy(:,:,:)
321 !OCL XFILL
322  rk_dynvar%work0(:,:,:,i_rk_rhot) = rhot(:,:,:)
323 !OCL XFILL
324  rk_dynvar%buf(:,:,:,:) = rk_dynvar%work0(:,:,:,:)
325  !$omp end workshare
326 
327  if (va > 0) then
328  !$omp workshare
329 !OCL XFILL
330  rk_prgvar%work0(:,:,:,:) = prog
331 !OCL XFILL
332  rk_prgvar%buf(:,:,:,:) = rk_prgvar%work0(:,:,:,:)
333  !$omp end workshare
334  end if
335  !$omp end parallel
336 
337  if ( bnd_w ) then
338  do j = js, je
339  do k = ks, ke
340  rk_mflx_hi%buf(k,is-1,j,2) = mflx_hi(k,is-1,j,2)
341  end do
342  end do
343  end if
344  if ( bnd_e ) then
345  do j = js, je
346  do k = ks, ke
347  rk_mflx_hi%buf(k,ie,j,2) = mflx_hi(k,ie,j,2)
348  end do
349  end do
350  end if
351  if ( bnd_s ) then
352  do i = is, ie
353  do k = ks, ke
354  rk_mflx_hi%buf(k,i,js-1,3) = mflx_hi(k,i,js-1,3)
355  end do
356  end do
357  end if
358  if ( bnd_n ) then
359  do i = is, ie
360  do k = ks, ke
361  rk_mflx_hi%buf(k,i,je,3) = mflx_hi(k,i,je,3)
362  end do
363  end do
364  end if
365 
366  io_prgvar(:) = 0; jo_prgvar(:) = 0; ko_prgvar(:) = 0
367 
368  call prof_rapend ("DYN_RK11s8o_Prep",3)
369 
370  !------------------------------------------------------------------------
371  ! Start RK
372  !------------------------------------------------------------------------
373 
374  do stage = 1, rk_nstage
375 
376  if ( stage > 1) then
377  call prof_rapstart("DYN_RK11s8o_BND",3)
378  call atmos_dyn_copy_boundary( rk_dynvar%buf(:,:,:,i_rk_dens), & ! [INOUT]
379  rk_dynvar%buf(:,:,:,i_rk_momz), & ! [INOUT]
380  rk_dynvar%buf(:,:,:,i_rk_momx), & ! [INOUT]
381  rk_dynvar%buf(:,:,:,i_rk_momy), & ! [INOUT]
382  rk_dynvar%buf(:,:,:,i_rk_rhot), & ! [INOUT]
383  rk_prgvar%buf(:,:,:,: ), & ! [INOUT]
384  rk_dynvar%work0(:,:,:,i_rk_dens), & ! [IN]
385  rk_dynvar%work0(:,:,:,i_rk_momz), & ! [IN]
386  rk_dynvar%work0(:,:,:,i_rk_momx), & ! [IN]
387  rk_dynvar%work0(:,:,:,i_rk_momy), & ! [IN]
388  rk_dynvar%work0(:,:,:,i_rk_rhot), & ! [IN]
389  rk_prgvar%work0(:,:,:,:), & ! [IN]
390  bnd_w, bnd_e, bnd_s, bnd_n, twod ) ! [IN]
391  call prof_rapend ("DYN_RK11s8o_BND",3)
392 
393  call rkcommon_comm( rk_dynvar )
394  call rkcommon_comm( rk_prgvar )
395  call rkcommon_comm_wait( rk_dynvar )
396  call rkcommon_comm_wait( rk_prgvar )
397  end if
398  !--
399 
400  call prof_rapstart("DYN_RK11s8o",3)
401 
402  call atmos_dyn_tstep( &
403  rk_dynvar%work(:,:,:,i_rk_dens,stage), & ! [OUT]
404  rk_dynvar%work(:,:,:,i_rk_momz,stage), & ! [OUT]
405  rk_dynvar%work(:,:,:,i_rk_momx,stage), & ! [OUT]
406  rk_dynvar%work(:,:,:,i_rk_momy,stage), & ! [OUT]
407  rk_dynvar%work(:,:,:,i_rk_rhot,stage), & ! [OUT]
408  rk_prgvar%work(:,:,:,: ,stage), & ! [OUT]
409  rk_mflx_hi%buf(:,:,:,:), rk_tflx_hi%buf(:,:,:,:), & ! [INOUT,OUT]
410  rk_dynvar%work0(:,:,:,i_rk_dens), & ! [IN]
411  rk_dynvar%work0(:,:,:,i_rk_momz), & ! [IN]
412  rk_dynvar%work0(:,:,:,i_rk_momx), & ! [IN]
413  rk_dynvar%work0(:,:,:,i_rk_momy), & ! [IN]
414  rk_dynvar%work0(:,:,:,i_rk_rhot), & ! [IN]
415  rk_dynvar%buf(:,:,:,i_rk_dens), & ! [IN]
416  rk_dynvar%buf(:,:,:,i_rk_momz), & ! [IN]
417  rk_dynvar%buf(:,:,:,i_rk_momx), & ! [IN]
418  rk_dynvar%buf(:,:,:,i_rk_momy), & ! [IN]
419  rk_dynvar%buf(:,:,:,i_rk_rhot), & ! [IN]
420  dens_t, momz_t, momx_t, momy_t, rhot_t, & ! [IN]
421  rk_prgvar%work0, rk_prgvar%buf(:,:,:,:), & ! [IN]
422  dpres0, cvtot, corioli, & ! [IN]
423  num_diff, wdamp_coef, divdmp_coef, ddiv, & ! [IN]
424  flag_fct_momentum, flag_fct_t, & ! [IN]
425  flag_fct_along_stream, & ! [IN]
426  cdz, fdz, fdx, fdy, & ! [IN]
427  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
428  phi, gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
429  ref_pres, ref_dens, & ! [IN]
430  bnd_w, bnd_e, bnd_s, bnd_n, twod, & ! [IN]
431  1.0_rp, .false. ) ! [IN]
432 
433  if ( stage < rk_nstage) then
434  call rkcommon_nextstage( rk_dynvar, stage, io_dynvar, jo_dynvar, ko_dynvar, dt )
435  call rkcommon_nextstage( rk_prgvar, stage, io_prgvar, jo_prgvar, ko_prgvar, dt )
436  else
437  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_dens, i_rk_dens, dt, dens )
438  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_momz, i_rk_momz, dt, momz )
439  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_momx, i_rk_momx, dt, momx )
440  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_momy, i_rk_momy, dt, momy )
441  call rkcommon_updatevar( rk_dynvar, io_dynvar, jo_dynvar, ko_dynvar, i_rk_rhot, i_rk_rhot, dt, rhot )
442  call rkcommon_updatevar( rk_prgvar, io_prgvar, jo_prgvar, ko_prgvar, 1, va, dt, prog )
443  end if
444  call rkcommon_updateflux( rk_mflx_hi, stage, 0, 0, 0, 3, mflx_hi )
445  call rkcommon_updateflux( rk_tflx_hi, stage, 0, 0, 0, 3, tflx_hi )
446 
447  call prof_rapend ("DYN_RK11s8o",3)
448  end do
449 
450  call rkcommon_rkwork_dealloc( rk_dynvar )
451  call rkcommon_rkwork_dealloc( rk_prgvar )
452 
453  return

References scale_atmos_dyn_common::atmos_dyn_copy_boundary(), scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_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().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ rkcoef_a_11s8o_cooperverner1972

real(rp), dimension(11,11), parameter, public scale_atmos_dyn_tinteg_short_rk11s8o::rkcoef_a_11s8o_cooperverner1972 = reshape( (/ ZERO, 0.5_RP, 0.25_RP, RC41, RC51, RC61, RC71, RC81, RC91, RC101, ZERO, ZERO, ZERO, 0.25_RP, RC42, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC43, RC53, RC63, RC73, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC54, RC64, RC74, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC65, RC75, RC85, RC95, RC105, RC115, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC76, RC86, RC96, RC106, RC116, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC87, RC97, RC107, RC117, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC98, RC108, RC118, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC109, RC119, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, RC1110, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO /), shape(RKCoef_a_11s8o_CooperVerner1972) )

Definition at line 102 of file scale_atmos_dyn_tinteg_short_rk11s8o.F90.

102  real(RP), parameter, public :: RKCoef_a_11s8o_CooperVerner1972(11,11) = reshape( &
103  (/ zero, 0.5_rp, 0.25_rp, rc41, rc51, rc61, rc71, rc81, rc91, rc101, zero, &
104  zero, zero, 0.25_rp, rc42, zero, zero, zero, zero, zero, zero, zero, &
105  zero, zero, zero, rc43, rc53, rc63, rc73, zero, zero, zero, zero, &
106  zero, zero, zero, zero, rc54, rc64, rc74, zero, zero, zero, zero, &
107  zero, zero, zero, zero, zero, rc65, rc75, rc85, rc95, rc105, rc115, &
108  zero, zero, zero, zero, zero, zero, rc76, rc86, rc96, rc106, rc116, &
109  zero, zero, zero, zero, zero, zero, zero, rc87, rc97, rc107, rc117, &
110  zero, zero, zero, zero, zero, zero, zero, zero, rc98, rc108, rc118, &
111  zero, zero, zero, zero, zero, zero, zero, zero, zero, rc109, rc119, &
112  zero, zero, zero, zero, zero, zero, zero, zero, zero, zero, rc1110, &
113  zero, zero, zero, zero, zero, zero, zero, zero, zero, zero, zero /), &
114  shape(rkcoef_a_11s8o_cooperverner1972) )

Referenced by atmos_dyn_tinteg_short_rk11s8o_setup().

◆ rkcoef_b_11s8o_cooperverner1972

real(rp), dimension(11), parameter, public scale_atmos_dyn_tinteg_short_rk11s8o::rkcoef_b_11s8o_cooperverner1972 = (/ 0.05_RP, ZERO, ZERO, ZERO, ZERO, ZERO, ZERO, 49.0_RP/180.0_RP, 16.0_RP/45.0_RP, 49.0_RP/180.0_RP, 0.05_RP /)

Definition at line 116 of file scale_atmos_dyn_tinteg_short_rk11s8o.F90.

116  real(RP), parameter, public :: RKCoef_b_11s8o_CooperVerner1972(11) = &
117  (/ 0.05_rp, zero, zero, zero, zero, zero, zero, 49.0_rp/180.0_rp, 16.0_rp/45.0_rp, &
118  49.0_rp/180.0_rp, 0.05_rp /)

Referenced by atmos_dyn_tinteg_short_rk11s8o_setup().

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_atmos_dyn_common
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_common.F90:12
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
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_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_const
module CONSTANT
Definition: scale_const.F90:11
scale_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
Definition: scale_comm_cartesC.F90:294
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_dyn_tstep_short
module Atmosphere / Dynamical scheme
Definition: scale_atmos_dyn_tstep_short.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_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_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