SCALE-RM
Functions/Subroutines | Variables
scale_atmos_dyn_tinteg_short_rk7s6o Module Reference

module Atmosphere / Dyn Tinteg More...

Functions/Subroutines

subroutine, public atmos_dyn_tinteg_short_rk7s6o_setup (tinteg_type)
 Setup. More...
 
subroutine, public atmos_dyn_tinteg_short_rk7s6o (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(7, 7), parameter, public rkcoef_a_7s6o_lawson1967 = reshape( (/ 0.0_RP, 3.0_RP/19.0_RP, 9.0_RP/152.0_RP, 94474764.0_RP/318611987.0_RP, -76607525678.0_RP/925997907411.0_RP, -113193410749715476.0_RP/1376008387821185625.0_RP, 510341547912673.0_RP/1709758911034368.0_RP, 0.0_RP, 0.0_RP, 27.0_RP/152.0_RP, -310753854.0_RP/318611987.0_RP, 309768324.0_RP/200562683.0_RP, 68309142.0_RP/42280325.0_RP, -3074637.0_RP/21410624.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 375818328.0_RP/318611987.0_RP, -57882086555344.0_RP/37088653028409.0_RP, -9901869473098663108168.0_RP/5940196722617929711875.0_RP, 205532548800199165.0_RP/6225256605226855824.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 643400862141470.0_RP/704684407539771.0_RP, 8947230518934447694268.0_RP/9333225588784524496875.0_RP, 32370527990426718666299.0_RP/90521226376106372167680.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -8377112295767292.0_RP/1089624335851065625.0_RP, 2610287999955961017.0_RP/236243323046620160.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -2690946369187951875.0_RP/253991013039290368.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP /), shape(RKCoef_a_7s6o_Lawson1967) )
 
real(rp), dimension(7), parameter, public rkcoef_b_7s6o_lawson1967 = (/ 119490041.0_RP/1597112640.0_RP, 0.0_RP, 55710603179056.0_RP/168638187800205.0_RP, 5739605598843081731.0_RP/28834038834414422400.0_RP, 1477688286853979.0_RP/291957783566400.0_RP, -298030839900625.0_RP/62778200252544.0_RP, 5352656.0_RP/65415735.0_RP /)
 
real(rp), dimension(7, 7), parameter, public rkcoef_a_7s6o_butcher1964 = reshape( (/ 0.0_RP, 1.0_RP/3.0_RP, 0.0_RP, 1.0_RP/12.0_RP, -1.0_RP/16.0_RP, 0.0_RP, 9.0_RP/44.0_RP, 0.0_RP, 0.0_RP, 2.0_RP/3.0_RP, 1.0_RP/3.0_RP, 9.0_RP/8.0_RP, 9.0_RP/8.0_RP, -9.0_RP/11.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -1.0_RP/12.0_RP, -3.0_RP/16.0_RP, -3_RP/8.0_RP, 63.0_RP/44.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -3.0_RP/8.0_RP, -3.0_RP/4.0_RP, 18.0_RP/11.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 1.0_RP/2.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -16.0_RP/11.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP /), shape(RKCoef_a_7s6o_Butcher1964) )
 
real(rp), dimension(7), parameter, public rkcoef_b_7s6o_butcher1964 = 1.0_RP/120.0_RP * (/ 11.0_RP, 0.0_RP, 81.0_RP, 81.0_RP, -32.0_RP, -32.0_RP, 11.0_RP /)
 

Detailed Description

module Atmosphere / Dyn Tinteg

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

This module provides two type of 7 stage Runge-Kutta scheme with 6th order accuracy:

Function/Subroutine Documentation

◆ atmos_dyn_tinteg_short_rk7s6o_setup()

subroutine, public scale_atmos_dyn_tinteg_short_rk7s6o::atmos_dyn_tinteg_short_rk7s6o_setup ( character(len=*)  tinteg_type)

Setup.

Definition at line 126 of file scale_atmos_dyn_tinteg_short_rk7s6o.F90.

126  use scale_prc, only: &
127  prc_abort
128  use scale_const, only: &
129  undef => const_undef
130  use scale_comm_cartesc, only: &
132  implicit none
133 
134  character(len=*) :: tinteg_type
135 
136  integer :: iv, d
137  character(H_SHORT) :: dynvar_name_list(5)
138  character(H_SHORT) :: prgvar_name_list(VA)
139  character(H_SHORT) :: flux_name_list(3)
140 
141  real(RP) :: RKCoef_a(RK_nstage,RK_nstage)
142  real(RP) :: RKCoef_b(RK_nstage)
143 
144  !---------------------------------------------------------------------------
145 
146  select case( trim(tinteg_type) )
147  case ('RK7s6o', 'RK7s6oLawson1967')
148  rkcoef_a(:,:) = rkcoef_a_7s6o_lawson1967
149  rkcoef_b(:) = rkcoef_b_7s6o_lawson1967
150  case ('RK7s6oButcher1964')
151  rkcoef_a(:,:) = rkcoef_a_7s6o_butcher1964
152  rkcoef_b(:) = rkcoef_b_7s6o_butcher1964
153  case default
154  log_error("ATMOS_DYN_Tinteg_short_rk7s6o_setup",*) 'The specified TINTEG_TYPE is invalid. Check!', tinteg_type
155  call prc_abort
156  end select
157 
158  dynvar_name_list(1) = 'DENS'
159  dynvar_name_list(2) = 'MOMZ'
160  dynvar_name_list(3) = 'MOMX'
161  dynvar_name_list(4) = 'MOMY'
162  dynvar_name_list(5) = 'RHOT'
163  call rkcommon_setup( rk_dynvar, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, dynvar_name_list )
164 
165  do iv = 1, va
166  flux_name_list(iv) = 'PROG'
167  end do
168  call rkcommon_setup( rk_prgvar, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b,prgvar_name_list, comm_id_offset=5 )
169 
170  do d = 1, 3
171  flux_name_list(iv) = 'mflx_hi'
172  end do
173  call rkcommon_setup( rk_mflx_hi, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, flux_name_list, is_type_flux=.true. )
174 
175 
176  do d = 1, 3
177  flux_name_list(iv) = 'tflx_hi'
178  end do
179  call rkcommon_setup( rk_tflx_hi, rk_nstage, rk_nregister, rkcoef_a, rkcoef_b, flux_name_list, is_type_flux=.true. )
180 
181  !----------------------------------------------------
182 
183  return

References scale_comm_cartesc::comm_vars8_init(), scale_const::const_undef, scale_prc::prc_abort(), rkcoef_a_7s6o_butcher1964, rkcoef_a_7s6o_lawson1967, rkcoef_b_7s6o_butcher1964, rkcoef_b_7s6o_lawson1967, 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_rk7s6o()

subroutine, public scale_atmos_dyn_tinteg_short_rk7s6o::atmos_dyn_tinteg_short_rk7s6o ( 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 202 of file scale_atmos_dyn_tinteg_short_rk7s6o.F90.

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

real(rp), dimension(7,7), parameter, public scale_atmos_dyn_tinteg_short_rk7s6o::rkcoef_a_7s6o_lawson1967 = reshape( (/ 0.0_RP, 3.0_RP/19.0_RP, 9.0_RP/152.0_RP, 94474764.0_RP/318611987.0_RP, -76607525678.0_RP/925997907411.0_RP, -113193410749715476.0_RP/1376008387821185625.0_RP, 510341547912673.0_RP/1709758911034368.0_RP, 0.0_RP, 0.0_RP, 27.0_RP/152.0_RP, -310753854.0_RP/318611987.0_RP, 309768324.0_RP/200562683.0_RP, 68309142.0_RP/42280325.0_RP, -3074637.0_RP/21410624.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 375818328.0_RP/318611987.0_RP, -57882086555344.0_RP/37088653028409.0_RP, -9901869473098663108168.0_RP/5940196722617929711875.0_RP, 205532548800199165.0_RP/6225256605226855824.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 643400862141470.0_RP/704684407539771.0_RP, 8947230518934447694268.0_RP/9333225588784524496875.0_RP, 32370527990426718666299.0_RP/90521226376106372167680.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -8377112295767292.0_RP/1089624335851065625.0_RP, 2610287999955961017.0_RP/236243323046620160.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -2690946369187951875.0_RP/253991013039290368.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP /), shape(RKCoef_a_7s6o_Lawson1967) )

Definition at line 66 of file scale_atmos_dyn_tinteg_short_rk7s6o.F90.

66  real(RP), parameter, public :: RKCoef_a_7s6o_Lawson1967(7,7) = reshape( &
67  (/ 0.0_rp, 3.0_rp/19.0_rp, 9.0_rp/152.0_rp, 94474764.0_rp/318611987.0_rp, -76607525678.0_rp/925997907411.0_rp, -113193410749715476.0_rp/1376008387821185625.0_rp, 510341547912673.0_rp/1709758911034368.0_rp, &
68  0.0_rp, 0.0_rp, 27.0_rp/152.0_rp, -310753854.0_rp/318611987.0_rp, 309768324.0_rp/200562683.0_rp, 68309142.0_rp/42280325.0_rp, -3074637.0_rp/21410624.0_rp, &
69  0.0_rp, 0.0_rp, 0.0_rp, 375818328.0_rp/318611987.0_rp, -57882086555344.0_rp/37088653028409.0_rp, -9901869473098663108168.0_rp/5940196722617929711875.0_rp, 205532548800199165.0_rp/6225256605226855824.0_rp, &
70  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 643400862141470.0_rp/704684407539771.0_rp, 8947230518934447694268.0_rp/9333225588784524496875.0_rp, 32370527990426718666299.0_rp/90521226376106372167680.0_rp, &
71  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, -8377112295767292.0_rp/1089624335851065625.0_rp, 2610287999955961017.0_rp/236243323046620160.0_rp, &
72  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, -2690946369187951875.0_rp/253991013039290368.0_rp, &
73  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp /), &
74  shape(rkcoef_a_7s6o_lawson1967) )

Referenced by atmos_dyn_tinteg_short_rk7s6o_setup().

◆ rkcoef_b_7s6o_lawson1967

real(rp), dimension(7), parameter, public scale_atmos_dyn_tinteg_short_rk7s6o::rkcoef_b_7s6o_lawson1967 = (/ 119490041.0_RP/1597112640.0_RP, 0.0_RP, 55710603179056.0_RP/168638187800205.0_RP, 5739605598843081731.0_RP/28834038834414422400.0_RP, 1477688286853979.0_RP/291957783566400.0_RP, -298030839900625.0_RP/62778200252544.0_RP, 5352656.0_RP/65415735.0_RP /)

Definition at line 76 of file scale_atmos_dyn_tinteg_short_rk7s6o.F90.

76  real(RP), parameter, public :: RKCoef_b_7s6o_Lawson1967(7) = &
77  (/ 119490041.0_rp/1597112640.0_rp, 0.0_rp, 55710603179056.0_rp/168638187800205.0_rp, &
78  5739605598843081731.0_rp/28834038834414422400.0_rp, 1477688286853979.0_rp/291957783566400.0_rp, -298030839900625.0_rp/62778200252544.0_rp, &
79  5352656.0_rp/65415735.0_rp /)

Referenced by atmos_dyn_tinteg_short_rk7s6o_setup().

◆ rkcoef_a_7s6o_butcher1964

real(rp), dimension(7,7), parameter, public scale_atmos_dyn_tinteg_short_rk7s6o::rkcoef_a_7s6o_butcher1964 = reshape( (/ 0.0_RP, 1.0_RP/3.0_RP, 0.0_RP, 1.0_RP/12.0_RP, -1.0_RP/16.0_RP, 0.0_RP, 9.0_RP/44.0_RP, 0.0_RP, 0.0_RP, 2.0_RP/3.0_RP, 1.0_RP/3.0_RP, 9.0_RP/8.0_RP, 9.0_RP/8.0_RP, -9.0_RP/11.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -1.0_RP/12.0_RP, -3.0_RP/16.0_RP, -3_RP/8.0_RP, 63.0_RP/44.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -3.0_RP/8.0_RP, -3.0_RP/4.0_RP, 18.0_RP/11.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 1.0_RP/2.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, -16.0_RP/11.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP, 0.0_RP /), shape(RKCoef_a_7s6o_Butcher1964) )

Definition at line 84 of file scale_atmos_dyn_tinteg_short_rk7s6o.F90.

84  real(RP), parameter, public :: RKCoef_a_7s6o_Butcher1964(7,7) = reshape( &
85  (/ 0.0_rp, 1.0_rp/3.0_rp, 0.0_rp, 1.0_rp/12.0_rp, -1.0_rp/16.0_rp, 0.0_rp, 9.0_rp/44.0_rp, &
86  0.0_rp, 0.0_rp, 2.0_rp/3.0_rp, 1.0_rp/3.0_rp, 9.0_rp/8.0_rp, 9.0_rp/8.0_rp, -9.0_rp/11.0_rp, &
87  0.0_rp, 0.0_rp, 0.0_rp, -1.0_rp/12.0_rp, -3.0_rp/16.0_rp, -3_rp/8.0_rp, 63.0_rp/44.0_rp, &
88  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, -3.0_rp/8.0_rp, -3.0_rp/4.0_rp, 18.0_rp/11.0_rp, &
89  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 1.0_rp/2.0_rp, 0.0_rp, &
90  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, -16.0_rp/11.0_rp, &
91  0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp, 0.0_rp /), &
92  shape(rkcoef_a_7s6o_butcher1964) )

Referenced by atmos_dyn_tinteg_short_rk7s6o_setup().

◆ rkcoef_b_7s6o_butcher1964

real(rp), dimension(7), parameter, public scale_atmos_dyn_tinteg_short_rk7s6o::rkcoef_b_7s6o_butcher1964 = 1.0_RP/120.0_RP * (/ 11.0_RP, 0.0_RP, 81.0_RP, 81.0_RP, -32.0_RP, -32.0_RP, 11.0_RP /)

Definition at line 94 of file scale_atmos_dyn_tinteg_short_rk7s6o.F90.

94  real(RP), parameter, public :: RKCoef_b_7s6o_Butcher1964(7) = &
95  1.0_rp/120.0_rp * (/ 11.0_rp, 0.0_rp, 81.0_rp, 81.0_rp, -32.0_rp, -32.0_rp, 11.0_rp /)

Referenced by atmos_dyn_tinteg_short_rk7s6o_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:350
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_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_comm_cartesc::comm_vars8_init
subroutine, public comm_vars8_init(varname, var, vid, gid)
Register variables.
Definition: scale_comm_cartesC.F90:811
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_dyn_common::atmos_dyn_copy_boundary
subroutine, public atmos_dyn_copy_boundary(DENS, MOMZ, MOMX, MOMY, RHOT, PROG, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, PROG0, BND_W, BND_E, BND_S, BND_N, TwoD)
Definition: scale_atmos_dyn_common.F90:215
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56