KASKADE 7 development version
|
#include <boost/timer/timer.hpp>
#include "dune/grid/config.h"
#include "dune/istl/operators.hh"
#include "dune/istl/solvers.hh"
#include "dune/common/dynvector.hh"
#include "linalg/dynamicMatrix.hh"
#include "utilities/timing.hh"
#include <cassert>
Go to the source code of this file.
Classes | |
class | Kaskade::SDCTimeGrid |
Abstract base class of time grids for (block) spectral defect correction methods. More... | |
class | Kaskade::SDCGridImplementationBase |
A convenience class simplifying the implementation of different SDCTimeGrid derived classes. More... | |
class | Kaskade::LobattoTimeGrid |
spectral time grid for defect correction methods with Lobatto points on \( [a,b]\) More... | |
class | Kaskade::RadauTimeGrid |
spectral time grid for defect correction methods with Radau points on \( [a,b]\). More... | |
class | Kaskade::GaussTimeGrid |
spectral time grid for defect correction methods with Gauss points More... | |
class | Kaskade::Sdc< Vector > |
Encapsulates the state of an SDC fixed point iteration and provides stepping. More... | |
Namespaces | |
namespace | Kaskade |
Functions | |
template<class Matrix , class Vectors , class ReactionDerivatives , class Solver > | |
Matrix::field_type | Kaskade::sdcIterationStep (SDCTimeGrid const &grid, SDCTimeGrid::RealMatrix const &Shat, Solver const &solve, Matrix const &M, Matrix const &A, Vectors const &r, ReactionDerivatives const &f_u, Vectors const &u, Vectors &du) |
A single spectral defect correction iteration sweep. More... | |
template<class Matrix , class Vectors , class Solver > | |
Matrix::field_type | Kaskade::sdcIMEXStep (SDCTimeGrid const &grid, SDCTimeGrid::RealMatrix const &Shat, Solver const &solve, Matrix const &M, Matrix const &A, Vectors const &r, Vectors const &u, Vectors &du) |
A single IMEX SDIRK spectral deferred correction sweep. More... | |
template<class Matrix , class Vectors , class ReactionDerivatives , class Solver > | |
Matrix::field_type | Kaskade::sdcIterationStep2 (SDCTimeGrid const &grid, SDCTimeGrid::RealMatrix const &Shat, Solver const &solve, Matrix const &M, Matrix const &A, Vectors const &rUi, ReactionDerivatives const &f_u, Vectors const &Mdu, Vectors &du) |
A single spectral defect correction iteration sweep. More... | |
template<class Matrix , class Vectors , class ReactionDerivatives > | |
double | Kaskade::waveformRelaxationStep2 (SDCTimeGrid const &grid, Matrix const &M, Matrix const &A, Vectors const &rUi, ReactionDerivatives const &f_u, Vectors &du) |
A single waveform relaxation iteration. More... | |
Approximate integration matrices | |
The convergence properties of SDC are essentially defined by the selection of an approximate integration matrix \( \hat S \) acting as a lower triangular preconditioner for the spectral integration matrix \( S \). We provide several choices. | |
void | Kaskade::eulerIntegrationMatrix (SDCTimeGrid const &grid, SDCTimeGrid::RealMatrix &Shat) |
Triangular approximate integration matrix \( \hat S \) corresponding to the implicit Euler integrator. More... | |
SDCTimeGrid::RealMatrix | Kaskade::eulerIntegrationMatrix (SDCTimeGrid const &grid) |
Triangular approximate integration matrix \( \hat S \) corresponding to the implicit Euler integrator. More... | |
void | Kaskade::luIntegrationMatrix (SDCTimeGrid const &grid, SDCTimeGrid::RealMatrix &Shat) |
Triangular approximate integration matrix \( \hat S \) corresponding to a specialized DIRK sweep. More... | |
SDCTimeGrid::RealMatrix | Kaskade::luIntegrationMatrix (SDCTimeGrid const &grid) |
Triangular approximate integration matrix \( \hat S \) corresponding to a specialized DIRK sweep. More... | |
SDCTimeGrid::RealMatrix | Kaskade::singleEulerIntegrationMatrix (SDCTimeGrid const &grid) |
Triangular approximate integration matrix \( \hat S \) corresponding to an Euler sweep with constant "step size". More... | |
SDCTimeGrid::RealMatrix | Kaskade::SDIRKIntegrationMatrix (SDCTimeGrid const &grid) |
Triangular approximate integration matrix \( \hat S \) corresponding to an SDIRK sweep (i.e. with constant "step size"). More... | |