FE-Project
Loading...
Searching...
No Matches
scale_timeint_rk::timeint_rk Type Reference

Derived type to provide RK scheme. More...

Public Member Functions

procedure, public init (this, rk_scheme_name, dt, var_num, ndim, size_each_var)
 Initialize a object to provide RK scheme.
procedure, public final (this)
 Finalize a object to provide RK scheme.
procedure, public get_implicit_diagfac (this, nowstage)
 Get a diagonal component of coefficient matrix at a stage for implicit RK scheme.
procedure, public get_deltime (this)
 Get timestep.
procedure, public advance1d (this, nowstage, q, varid, is, ie)
 Advance a variable at the current stage with explicit RK part.
procedure, public advance1d_varlist (this, nowstage, var_list, varids, is, ie)
 Advance variables at the current stage with explicit RK part.
procedure, public advance_trcvar_1d (this, nowstage, q, varid, is, ie, ddens, ddens0, dens_hyd)
 Advance tracer variable data at the current stage with explicit RK part.
procedure, public storevar0_1d (this, q, varid, is, ie)
 Store variable data at the time level n for the case of IMEX RK scheme.
procedure, public storeimplicit1d (this, nowstage, q, varid, is, ie)
 Store variable data after the implicit part of current stage in IMEX RK scheme.
procedure, public advance2d (this, nowstage, q, varid, is, ie, js, je)
 Advance a variable at the current stage with explicit RK part.
procedure, public advance2d_varlist (this, nowstage, var_list, varids, is, ie, js, je)
 Advance variables at the current stage with explicit RK part.
procedure, public advance_trcvar_2d (this, nowstage, q, varid, is, ie, js, je, ddens, ddens0, dens_hyd)
 Advance tracer variable data at the current stage with explicit RK part.
procedure, public storevar0_2d (this, q, varid, is, ie, js, je)
 Store variable data at the time level n for the case of IMEX RK scheme.
procedure, public storeimplicit2d (this, nowstage, q, varid, is, ie, js, je)
 Store variable data after the implicit part of current stage in IMEX RK scheme.
procedure, public advance3d (this, nowstage, q, varid, is, ie, js, je, ks, ke)
 Advance a variable at the current stage with explicit RK part.
procedure, public advance3d_varlist (this, nowstage, var_list, varids, is, ie, js, je, ks, ke)
 Advance variables at the current stage with explicit RK part.
procedure, public advance_trcvar_3d (this, nowstage, q, varid, is, ie, js, je, ks, ke, ddens, ddens0, dens_hyd)
 Advance tracer variable data at the current stage with explicit RK part.
procedure, public storevar0_3d (this, q, varid, is, ie, js, je, ks, ke)
 Store variable data at the time level n for the case of IMEX RK scheme.
procedure, public storeimplicit3d (this, nowstage, q, varid, is, ie, js, je, ks, ke)
 Store variable data after the implicit part of current stage in IMEX RK scheme.
generic, public advance advance1d
generic, public advance advance2d
generic, public advance advance3d
generic, public advance_varlist advance1d_varlist
generic, public advance_varlist advance2d_varlist
generic, public advance_varlist advance3d_varlist
generic, public advance_trcvar advance_trcvar_1d
generic, public advance_trcvar advance_trcvar_2d
generic, public advance_trcvar advance_trcvar_3d
generic, public storevar0 storevar0_1d
generic, public storevar0 storevar0_2d
generic, public storevar0 storevar0_3d
generic, public storeimplicit storeimplicit1d
generic, public storeimplicit storeimplicit2d
generic, public storeimplicit storeimplicit3d

Public Attributes

integer, public nstage
 Number of stage in RK scheme.
integer, public tend_buf_size
 Buffer size to store tendencies.
real(rp), dimension(:,:), allocatable, public coef_a_ex
 Coefficients A in the Butcher table (explicit part)
real(rp), dimension(:), allocatable, public coef_b_ex
 Coefficients b in the Butcher table (explicit part)
real(rp), dimension(:), allocatable, public coef_c_ex
 Coefficients c in the Butcher table (explicit part)
real(rp), dimension(:,:), allocatable, public coef_sig_ex
real(rp), dimension(:,:), allocatable, public coef_gam_ex
real(rp), dimension(:,:), allocatable, public coef_a_im
 Coefficients A in the Butcher table (implicit part)
real(rp), dimension(:), allocatable, public coef_b_im
 Coefficients b in the Butcher table (implicit part)
real(rp), dimension(:), allocatable, public coef_c_im
 Coefficients c in the Butcher table (implicit part)
integer, dimension(:), allocatable tend_buf_indmap
integer var_num
real(rp), dimension(:,:,:), allocatable, public var_buf1d_ex
real(rp), dimension(:,:,:), allocatable, public tend_buf1d_ex
real(rp), dimension(:,:,:), allocatable, public tend_buf1d_im
real(rp), dimension(:,:), allocatable, public var0_1d
real(rp), dimension(:,:,:,:), allocatable, public var_buf2d_ex
real(rp), dimension(:,:,:,:), allocatable, public tend_buf2d_ex
real(rp), dimension(:,:,:,:), allocatable, public tend_buf2d_im
real(rp), dimension(:,:,:), allocatable, public var0_2d
real(rp), dimension(:,:,:,:,:), allocatable, public var_buf3d_ex
real(rp), dimension(:,:,:,:,:), allocatable, public tend_buf3d_ex
real(rp), dimension(:,:,:,:,:), allocatable, public tend_buf3d_im
real(rp), dimension(:,:,:,:), allocatable, public var0_3d
logical, public imex_flag
integer, public ndim

Detailed Description

Derived type to provide RK scheme.

Definition at line 35 of file scale_timeint_rk.F90.

Member Function/Subroutine Documentation

◆ init()

procedure, public scale_timeint_rk::timeint_rk::init ( class(timeint_rk), intent(inout) this,
character(*), intent(in) rk_scheme_name,
real(rp), intent(in) dt,
integer, intent(in) var_num,
integer, intent(in) ndim,
integer, dimension(ndim), intent(in) size_each_var )

Initialize a object to provide RK scheme.

Parameters
rk_scheme_nameName of RK scheme
dtTimestep
ndimNumber of spatial Dimension
var_numNumber of variables
size_each_varArray to store size of each dimension

Definition at line 81 of file scale_timeint_rk.F90.

◆ final()

procedure, public scale_timeint_rk::timeint_rk::final ( class(timeint_rk), intent(inout) this)

Finalize a object to provide RK scheme.

Definition at line 82 of file scale_timeint_rk.F90.

◆ get_implicit_diagfac()

procedure, public scale_timeint_rk::timeint_rk::get_implicit_diagfac ( class(timeint_rk), intent(in) this,
integer, intent(in) nowstage )

Get a diagonal component of coefficient matrix at a stage for implicit RK scheme.

Parameters
nowstageCurrent stage

Definition at line 83 of file scale_timeint_rk.F90.

◆ get_deltime()

procedure, public scale_timeint_rk::timeint_rk::get_deltime ( class(timeint_rk), intent(in) this)

Get timestep.

Definition at line 84 of file scale_timeint_rk.F90.

◆ advance1d()

procedure, public scale_timeint_rk::timeint_rk::advance1d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie )

Advance a variable at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
qArray to store variable data
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction

Definition at line 85 of file scale_timeint_rk.F90.

◆ advance1d_varlist()

procedure, public scale_timeint_rk::timeint_rk::advance1d_varlist ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
type(timeint_rk_var), dimension(size(varids)), intent(inout) var_list,
integer, dimension(:), intent(in) varids,
integer, intent(in) is,
integer, intent(in) ie )

Advance variables at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
var_listArray of timeint_rk_var object which saves pointer to variable data
varIDsIndices of the targeting variables
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction

Definition at line 86 of file scale_timeint_rk.F90.

◆ advance_trcvar_1d()

procedure, public scale_timeint_rk::timeint_rk::advance_trcvar_1d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
real(rp), dimension (:), intent(in) ddens,
real(rp), dimension(:), intent(in) ddens0,
real(rp), dimension(:), intent(in) dens_hyd )

Advance tracer variable data at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
qArray to store tracer variable data
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
DDENSArray to store density perturbation data at the time level n+1
DDENS0Array to store density perturbation data at the time level n
DENS_hydArray to hydrostatic part of density

Definition at line 87 of file scale_timeint_rk.F90.

◆ storevar0_1d()

procedure, public scale_timeint_rk::timeint_rk::storevar0_1d ( class(timeint_rk), intent(inout) this,
real(rp), dimension(:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie )

Store variable data at the time level n for the case of IMEX RK scheme.

Parameters
qArray to store variable data at the time level n
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction

Definition at line 88 of file scale_timeint_rk.F90.

◆ storeimplicit1d()

procedure, public scale_timeint_rk::timeint_rk::storeimplicit1d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie )

Store variable data after the implicit part of current stage in IMEX RK scheme.

Parameters
nowstageCurrent stage
qArray to store variable data after the implicit part of current stage
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction

Definition at line 89 of file scale_timeint_rk.F90.

◆ advance2d()

procedure, public scale_timeint_rk::timeint_rk::advance2d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je )

Advance a variable at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
qArray to store variable data
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction

Definition at line 90 of file scale_timeint_rk.F90.

◆ advance2d_varlist()

procedure, public scale_timeint_rk::timeint_rk::advance2d_varlist ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
type(timeint_rk_var), dimension(size(varids)), intent(inout) var_list,
integer, dimension(:), intent(in) varids,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je )

Advance variables at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
var_listArray of timeint_rk_var object which saves pointer to variable data
varIDsIndices of the targeting variables
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction

Definition at line 91 of file scale_timeint_rk.F90.

◆ advance_trcvar_2d()

procedure, public scale_timeint_rk::timeint_rk::advance_trcvar_2d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je,
real(rp), dimension (:,:), intent(in) ddens,
real(rp), dimension(:,:), intent(in) ddens0,
real(rp), dimension(:,:), intent(in) dens_hyd )

Advance tracer variable data at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
qArray to store tracer variable data
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction
DDENSArray to store density perturbation data at the time level n+1
DDENS0Array to store density perturbation data at the time level n
DENS_hydArray to hydrostatic part of density

Definition at line 92 of file scale_timeint_rk.F90.

◆ storevar0_2d()

procedure, public scale_timeint_rk::timeint_rk::storevar0_2d ( class(timeint_rk), intent(inout) this,
real(rp), dimension(:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je )

Store variable data at the time level n for the case of IMEX RK scheme.

Parameters
qArray to store variable data at the time level n
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction

Definition at line 93 of file scale_timeint_rk.F90.

◆ storeimplicit2d()

procedure, public scale_timeint_rk::timeint_rk::storeimplicit2d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je )

Store variable data after the implicit part of current stage in IMEX RK scheme.

Parameters
nowstageCurrent stage
qArray to store variable data after the implicit part of current stage
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction

Definition at line 94 of file scale_timeint_rk.F90.

◆ advance3d()

procedure, public scale_timeint_rk::timeint_rk::advance3d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:,:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je,
integer, intent(in) ks,
integer, intent(in) ke )

Advance a variable at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
qArray to store variable data
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction
ksIndex at which begins the loop for the corresponding direction
keIndex at which finishes the loop for the corresponding direction

Definition at line 95 of file scale_timeint_rk.F90.

◆ advance3d_varlist()

procedure, public scale_timeint_rk::timeint_rk::advance3d_varlist ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
type(timeint_rk_var), dimension(size(varids)), intent(inout) var_list,
integer, dimension(:), intent(in) varids,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je,
integer, intent(in) ks,
integer, intent(in) ke )

Advance variables at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
var_listArray of timeint_rk_var object which saves pointer to variable data
varIDsIndices of the targeting variables
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction
ksIndex at which begins the loop for the corresponding direction
keIndex at which finishes the loop for the corresponding direction

Definition at line 96 of file scale_timeint_rk.F90.

◆ advance_trcvar_3d()

procedure, public scale_timeint_rk::timeint_rk::advance_trcvar_3d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:,:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je,
integer, intent(in) ks,
integer, intent(in) ke,
real(rp), dimension (:,:,:), intent(in) ddens,
real(rp), dimension(:,:,:), intent(in) ddens0,
real(rp), dimension(:,:,:), intent(in) dens_hyd )

Advance tracer variable data at the current stage with explicit RK part.

Parameters
nowstageCurrent stage
qArray to store tracer variable data
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction
ksIndex at which begins the loop for the corresponding direction
keIndex at which finishes the loop for the corresponding direction
DDENSArray to store density perturbation data at the time level n+1
DDENS0Array to store density perturbation data at the time level n
DENS_hydArray to hydrostatic part of density

Definition at line 97 of file scale_timeint_rk.F90.

◆ storevar0_3d()

procedure, public scale_timeint_rk::timeint_rk::storevar0_3d ( class(timeint_rk), intent(inout) this,
real(rp), dimension(:,:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je,
integer, intent(in) ks,
integer, intent(in) ke )

Store variable data at the time level n for the case of IMEX RK scheme.

Parameters
qArray to store variable data at the time level n
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction
ksIndex at which begins the loop for the corresponding direction
keIndex at which finishes the loop for the corresponding direction

Definition at line 98 of file scale_timeint_rk.F90.

◆ storeimplicit3d()

procedure, public scale_timeint_rk::timeint_rk::storeimplicit3d ( class(timeint_rk), intent(inout) this,
integer, intent(in) nowstage,
real(rp), dimension(:,:,:), intent(inout) q,
integer, intent(in) varid,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) js,
integer, intent(in) je,
integer, intent(in) ks,
integer, intent(in) ke )

Store variable data after the implicit part of current stage in IMEX RK scheme.

Parameters
nowstageCurrent stage
qArray to store variable data after the implicit part of current stage
varIDIndex of the targeting variable
isIndex at which begins the loop for the corresponding direction
ieIndex at which finishes the loop for the corresponding direction
jsIndex at which begins the loop for the corresponding direction
jeIndex at which finishes the loop for the corresponding direction
ksIndex at which begins the loop for the corresponding direction
keIndex at which finishes the loop for the corresponding direction

Definition at line 99 of file scale_timeint_rk.F90.

◆ advance() [1/3]

generic, public scale_timeint_rk::timeint_rk::advance

Definition at line 100 of file scale_timeint_rk.F90.

◆ advance() [2/3]

generic, public scale_timeint_rk::timeint_rk::advance

Definition at line 100 of file scale_timeint_rk.F90.

◆ advance() [3/3]

generic, public scale_timeint_rk::timeint_rk::advance

Definition at line 100 of file scale_timeint_rk.F90.

◆ advance_varlist() [1/3]

generic, public scale_timeint_rk::timeint_rk::advance_varlist

Definition at line 101 of file scale_timeint_rk.F90.

◆ advance_varlist() [2/3]

generic, public scale_timeint_rk::timeint_rk::advance_varlist

Definition at line 101 of file scale_timeint_rk.F90.

◆ advance_varlist() [3/3]

generic, public scale_timeint_rk::timeint_rk::advance_varlist

Definition at line 101 of file scale_timeint_rk.F90.

◆ advance_trcvar() [1/3]

generic, public scale_timeint_rk::timeint_rk::advance_trcvar

Definition at line 102 of file scale_timeint_rk.F90.

◆ advance_trcvar() [2/3]

generic, public scale_timeint_rk::timeint_rk::advance_trcvar

Definition at line 102 of file scale_timeint_rk.F90.

◆ advance_trcvar() [3/3]

generic, public scale_timeint_rk::timeint_rk::advance_trcvar

Definition at line 102 of file scale_timeint_rk.F90.

◆ storevar0() [1/3]

generic, public scale_timeint_rk::timeint_rk::storevar0

Definition at line 103 of file scale_timeint_rk.F90.

◆ storevar0() [2/3]

generic, public scale_timeint_rk::timeint_rk::storevar0

Definition at line 103 of file scale_timeint_rk.F90.

◆ storevar0() [3/3]

generic, public scale_timeint_rk::timeint_rk::storevar0

Definition at line 103 of file scale_timeint_rk.F90.

◆ storeimplicit() [1/3]

generic, public scale_timeint_rk::timeint_rk::storeimplicit

Definition at line 104 of file scale_timeint_rk.F90.

◆ storeimplicit() [2/3]

generic, public scale_timeint_rk::timeint_rk::storeimplicit

Definition at line 104 of file scale_timeint_rk.F90.

◆ storeimplicit() [3/3]

generic, public scale_timeint_rk::timeint_rk::storeimplicit

Definition at line 104 of file scale_timeint_rk.F90.

Member Data Documentation

◆ nstage

integer, public scale_timeint_rk::timeint_rk::nstage

Number of stage in RK scheme.

Definition at line 42 of file scale_timeint_rk.F90.

42 integer, public :: nstage !< Number of stage in RK scheme

◆ tend_buf_size

integer, public scale_timeint_rk::timeint_rk::tend_buf_size

Buffer size to store tendencies.

Definition at line 43 of file scale_timeint_rk.F90.

43 integer, public :: tend_buf_size !< Buffer size to store tendencies

◆ coef_a_ex

real(rp), dimension(:,:), allocatable, public scale_timeint_rk::timeint_rk::coef_a_ex

Coefficients A in the Butcher table (explicit part)

Definition at line 46 of file scale_timeint_rk.F90.

46 real(RP), public, allocatable :: coef_a_ex(:,:) !< Coefficients A in the Butcher table (explicit part)

◆ coef_b_ex

real(rp), dimension(:), allocatable, public scale_timeint_rk::timeint_rk::coef_b_ex

Coefficients b in the Butcher table (explicit part)

Definition at line 47 of file scale_timeint_rk.F90.

47 real(RP), public, allocatable :: coef_b_ex(:) !< Coefficients b in the Butcher table (explicit part)

◆ coef_c_ex

real(rp), dimension(:), allocatable, public scale_timeint_rk::timeint_rk::coef_c_ex

Coefficients c in the Butcher table (explicit part)

Definition at line 48 of file scale_timeint_rk.F90.

48 real(RP), public, allocatable :: coef_c_ex(:) !< Coefficients c in the Butcher table (explicit part)

◆ coef_sig_ex

real(rp), dimension(:,:), allocatable, public scale_timeint_rk::timeint_rk::coef_sig_ex

Definition at line 51 of file scale_timeint_rk.F90.

51 real(RP), public, allocatable :: coef_sig_ex(:,:)

◆ coef_gam_ex

real(rp), dimension(:,:), allocatable, public scale_timeint_rk::timeint_rk::coef_gam_ex

Definition at line 52 of file scale_timeint_rk.F90.

52 real(RP), public, allocatable :: coef_gam_ex(:,:)

◆ coef_a_im

real(rp), dimension(:,:), allocatable, public scale_timeint_rk::timeint_rk::coef_a_im

Coefficients A in the Butcher table (implicit part)

Definition at line 55 of file scale_timeint_rk.F90.

55 real(RP), public, allocatable :: coef_a_im(:,:) !< Coefficients A in the Butcher table (implicit part)

◆ coef_b_im

real(rp), dimension(:), allocatable, public scale_timeint_rk::timeint_rk::coef_b_im

Coefficients b in the Butcher table (implicit part)

Definition at line 56 of file scale_timeint_rk.F90.

56 real(RP), public, allocatable :: coef_b_im(:) !< Coefficients b in the Butcher table (implicit part)

◆ coef_c_im

real(rp), dimension(:), allocatable, public scale_timeint_rk::timeint_rk::coef_c_im

Coefficients c in the Butcher table (implicit part)

Definition at line 57 of file scale_timeint_rk.F90.

57 real(RP), public, allocatable :: coef_c_im(:) !< Coefficients c in the Butcher table (implicit part)

◆ tend_buf_indmap

integer, dimension(:), allocatable scale_timeint_rk::timeint_rk::tend_buf_indmap

Definition at line 59 of file scale_timeint_rk.F90.

59 integer, allocatable :: tend_buf_indmap(:)

◆ var_num

integer scale_timeint_rk::timeint_rk::var_num

Definition at line 60 of file scale_timeint_rk.F90.

60 integer :: var_num

◆ var_buf1d_ex

real(rp), dimension(:,:,:), allocatable, public scale_timeint_rk::timeint_rk::var_buf1d_ex

Definition at line 61 of file scale_timeint_rk.F90.

61 real(RP), public, allocatable :: var_buf1D_ex(:,:,:)

◆ tend_buf1d_ex

real(rp), dimension(:,:,:), allocatable, public scale_timeint_rk::timeint_rk::tend_buf1d_ex

Definition at line 62 of file scale_timeint_rk.F90.

62 real(RP), public, allocatable :: tend_buf1D_ex(:,:,:)

◆ tend_buf1d_im

real(rp), dimension(:,:,:), allocatable, public scale_timeint_rk::timeint_rk::tend_buf1d_im

Definition at line 63 of file scale_timeint_rk.F90.

63 real(RP), public, allocatable :: tend_buf1D_im(:,:,:)

◆ var0_1d

real(rp), dimension(:,:), allocatable, public scale_timeint_rk::timeint_rk::var0_1d

Definition at line 64 of file scale_timeint_rk.F90.

64 real(RP), public, allocatable :: var0_1D(:,:)

◆ var_buf2d_ex

real(rp), dimension(:,:,:,:), allocatable, public scale_timeint_rk::timeint_rk::var_buf2d_ex

Definition at line 66 of file scale_timeint_rk.F90.

66 real(RP), public, allocatable :: var_buf2D_ex(:,:,:,:)

◆ tend_buf2d_ex

real(rp), dimension(:,:,:,:), allocatable, public scale_timeint_rk::timeint_rk::tend_buf2d_ex

Definition at line 67 of file scale_timeint_rk.F90.

67 real(RP), public, allocatable :: tend_buf2D_ex(:,:,:,:)

◆ tend_buf2d_im

real(rp), dimension(:,:,:,:), allocatable, public scale_timeint_rk::timeint_rk::tend_buf2d_im

Definition at line 68 of file scale_timeint_rk.F90.

68 real(RP), public, allocatable :: tend_buf2D_im(:,:,:,:)

◆ var0_2d

real(rp), dimension(:,:,:), allocatable, public scale_timeint_rk::timeint_rk::var0_2d

Definition at line 69 of file scale_timeint_rk.F90.

69 real(RP), public, allocatable :: var0_2D(:,:,:)

◆ var_buf3d_ex

real(rp), dimension(:,:,:,:,:), allocatable, public scale_timeint_rk::timeint_rk::var_buf3d_ex

Definition at line 71 of file scale_timeint_rk.F90.

71 real(RP), public, allocatable :: var_buf3D_ex(:,:,:,:,:)

◆ tend_buf3d_ex

real(rp), dimension(:,:,:,:,:), allocatable, public scale_timeint_rk::timeint_rk::tend_buf3d_ex

Definition at line 72 of file scale_timeint_rk.F90.

72 real(RP), public, allocatable :: tend_buf3D_ex(:,:,:,:,:)

◆ tend_buf3d_im

real(rp), dimension(:,:,:,:,:), allocatable, public scale_timeint_rk::timeint_rk::tend_buf3d_im

Definition at line 73 of file scale_timeint_rk.F90.

73 real(RP), public, allocatable :: tend_buf3D_im(:,:,:,:,:)

◆ var0_3d

real(rp), dimension(:,:,:,:), allocatable, public scale_timeint_rk::timeint_rk::var0_3d

Definition at line 74 of file scale_timeint_rk.F90.

74 real(RP), public, allocatable :: var0_3D(:,:,:,:)

◆ imex_flag

logical, public scale_timeint_rk::timeint_rk::imex_flag

Definition at line 78 of file scale_timeint_rk.F90.

78 logical, public :: imex_flag

◆ ndim

integer, public scale_timeint_rk::timeint_rk::ndim

Definition at line 79 of file scale_timeint_rk.F90.

79 integer, public :: ndim

The documentation for this type was generated from the following file: