The scheme describes the discretization to use in the simulation.
There are two parts that need to be configured:
spatial
discretizationtemporal
discretizationOptionally a stabilization
method may be defined for the scheme, see the
atl_stabilization_module for more details on that definition.
The kind of spatial discretization is chosen via the name
setting within
the spatial
table. The following spatial schemes are available:
The main configuration option for the spatial discretization is the polynomial
degree to use to represent the state in the DG elements.
This polynomial degree is set by the option m
in the spatial table.
It may either be a simple scalar value, defining a single polynomial degree to
be used for all elements in the domain, or a table that provides polynomial
degrees based on the local refinement level of elements.
This can be achieved either by providing the polynomial degree for each level
individually or by choosing a predefined scheme to choose the polynomial
degree for the levels.
Individual definitions take the following form:
m = {
{ level = 4, m = 4 },
{ level = 5, m = 6 }
}
A predefined scheme is offered by 'fixedfact'
, where the polynomial degree
on each level is computed by the following formula.
m(iLevel) = nint( base_order * factor**(maxLevel-iLevel)) - 1
Here the base_order
and factor
need to be defined in the configuration, where
the base_order
sets the minimal polynomial degree (+1) that is to be used on the
finest level. factor
defines the factor, by which the scheme order is to be
increased by each level. The polynomial degree definition for this case looks as
follows:
m = {
predefined = 'fixedfact',
base_order = 4,
factor = 1.5
}
The default factor for the 'fixedfact' scheme is , which allows for approximately the same time step restriction across the levels for hyperbolic equations.
Besides the polynomial degree m
it is also possible to choose the polynomial
space to use for multidimensional representation.
This modg_space
is either 'P' or 'Q'.
P indicates a multidimensional polynomial, where the sum of the mode indices
is at most equal to the configured polynomial degree m
.
Q indicates that each index in the different dimensions itself may at most be
m
. The 'Q' space is the default but requires more computational effort and memory
especially for 3D simulations.
The explicit time integration is configured by the temporal
table within scheme
.
Following schemes are available:
steps=4
steps=4
steps=2
The 'imexRungeKutta' scheme should be used when penalization terms are to be used in the flow simulations. The 'explicitRungeKuttaTaylor' is mainly intended for the solution of linear equation systems.
The time step width is controlled by a control
subtable and the time step can
either be chosend adaptively according to the CFL condition, or set as a fixed
time step. See also atl_global_time_integration_module.
A complete definition of the scheme without the optional stabilization
table
takes the following form:
scheme = {
spatial = {
name = 'modg',
modg_space = 'Q',
m = 11
},
temporal = {
name = 'explicitRungeKutta',
steps = 4,
control = {
name = 'cfl',
cfl = 0.8
}
}
}
For details on the optional stabilization
see the atl_stabilization_module.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public, | parameter | :: | atl_modg_scheme_prp | = | 6 | |
integer, | public, | parameter | :: | atl_modg_2d_scheme_prp | = | 7 | |
integer, | public, | parameter | :: | atl_modg_1d_scheme_prp | = | 8 |
Datatype to specify the timestepping method.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=rk), | private | :: | dt | The local timestep. |
type to define a one dimensional stencil for reconstructions.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | private | :: | stencil | the 1D stencil in treelm coordinates. |
|||
integer, | private | :: | nElems | the number of elements in the stencil, including the cell itself you reconstruct for. |
|||
integer, | private, | allocatable | :: | elemPos(:) | relative position of the stencil elements to the current cell. Note, that this vector has length (nElems-1) since the current cell itself is not stored here. |
||
integer, | private, | allocatable | :: | ngElemPos(:) | relative position of the stencil elements in negative direction to the current cell. |
||
integer, | private, | allocatable | :: | bnd(:,:) | for each element of the mesh we store the lowest and highest left shift that build correct stencils (i.e. correct means: not including any boundary element). The first dimension is the number of elements associated with this stencil. The second dimension is 2, the first is the lowest possible left shift index the second is the highest possible left shift index. |
type specifying all informations about the stencil for the dimension by dimension reconstruction.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
type(atl_oneDimStencil_type), | private | :: | xStencil | the stencil in x direction |
|||
type(atl_oneDimStencil_type), | private | :: | yStencil | the stencil in y direction |
|||
type(atl_oneDimStencil_type), | private | :: | zStencil | the stencil in z direction |
type containing all the informations related to the scheme, e.g.: time and space discretization, scheme order, etc.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | private | :: | scheme | integer representing the current discretization scheme. |
|||
integer, | private | :: | nDoFs | the number of degrees of freedom for the selected scheme for a single cell and a single variable of the equation. For example we have: P1PM => nDofs=4, P2PM = 10). This number includes only the degrees of freedom which will be stored. We do not include the number of reconstructed degrees of freedom here! |
|||
integer, | private | :: | nDoFsRecons | the number of reconstructed degrees of freedom for the selected scheme for a single cell and a single variable of the equation (including the reconstructed degrees of freedoms). |
|||
integer, | private | :: | nFaceDofs | The number of dofs on the faces. |
|||
type(atl_local_timestep_type), | private | :: | time | variable to specify the space integration. levelwise information of time discretization |
|||
type(atl_modg_scheme_type), | private | :: | modg | Parameters of the modal discontinuous Galerkin scheme if scheme is set to modg. |
|||
type(atl_modg_2d_scheme_type), | private | :: | modg_2d | Parameters of the modal discontinuous Galerkin scheme if scheme is set to modg 2d. |
|||
type(atl_modg_1d_scheme_type), | private | :: | modg_1d | Parameters of the modal discontinuous Galerkin scheme if scheme is set to modg 1d. |
|||
type(ply_modg_basis_type), | private | :: | modg_basis | Informations about the polynomial basis of a MODG scheme. |
|||
type(atl_stabilization_type), | private, | allocatable | :: | stabilization(:) | The stabilization(s) for the scheme. Applied one after each other. Starting with index 1, then 2, ... |
||
real(kind=rk), | private, | allocatable | :: | dl_prod(:,:) | Precomputed Scalar Products |
||
real(kind=rk), | private, | allocatable | :: | dl_prodDiff(:,:) | |||
real(kind=rk), | private, | allocatable | :: | temp_over(:,:,:) | Temp Arrays needed for evaluation of physical fluxes |
||
real(kind=rk), | private, | allocatable | :: | temp_modal(:,:,:) | |||
real(kind=rk), | private, | allocatable | :: | temp_nodal(:,:,:) |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | schemeID |
subroutine to intialize a scheme as specified by a given lua script file.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(atl_scheme_type), | intent(out) | :: | me(minlevel:maxlevel) | the scheme you want to initialize. |
||
type(flu_State), | intent(in) | :: | conf | flu binding to lua configuration file. |
||
integer, | intent(in) | :: | minLevel | The global minimum level of the mesh |
||
integer, | intent(in) | :: | maxLevel | The global maximum level of the mesh |
Subroutine do define a specific stencil for a certain scheme.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | nDims | Number of dimensions to consider in the equation. |
||
type(tem_stencilHeader_type), | intent(inout) | :: | me | the neighbor list you want to init. |
precompute the scalar products of the anstaz and test function
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=rk), | intent(out), | allocatable | :: | dl_prod(:,:) | ||
real(kind=rk), | intent(out), | allocatable | :: | dl_prodDiff(:,:) | ||
integer, | intent(in) | :: | maxPolyDegree |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(atl_local_timestep_type), | intent(inout) | :: | me | the scheme you want to initialize. |