![]() |
SG++-Doxygen-Documentation
|
This class provides a simple-to-use solver of the multi dimensional Poisson Equation on Sparse Grids. More...
#include <PoissonEquationSolver.hpp>
Public Member Functions | |
| void | constructGrid (sgpp::base::BoundingBox &myBoundingBox, size_t level) |
| Use this routine the construct a regular grid to solve a PDE. More... | |
| void | initGridWithExpHeat (sgpp::base::DataVector &alpha, double factor=1.0) |
| Inits the grid with a heat distribution based on the e-function. More... | |
| void | initGridWithExpHeatFullDomain (sgpp::base::DataVector &alpha, double factor=1.0) |
| Inits the grid with a heat distribution based on the e-function. More... | |
| void | initGridWithSmoothHeat (sgpp::base::DataVector &alpha, double mu, double sigma, double factor) |
| Inits the grid with a smooth heat distribution (based on a std-normal distribution) on its boundaries. More... | |
| void | initGridWithSmoothHeatFullDomain (sgpp::base::DataVector &alpha, double mu, double sigma, double factor) |
| Inits the grid with a smooth heat distribution (based on a std-normal distribution) on its boundaries. More... | |
| void | initScreen () |
| Inits the screen object. More... | |
| PoissonEquationSolver () | |
| Std-Constructor of the solver. More... | |
| void | solvePDE (sgpp::base::DataVector &alpha, sgpp::base::DataVector &rhs, size_t maxCGIterations, double epsilonCG, bool verbose=false) |
| abstract method to solve an elliptic PDE. More... | |
| void | storeInnerRHS (sgpp::base::DataVector &alpha, std::string tFilename) |
| Routine to export the RHS of the inner system which has to be solved in order to solve the Poisson equation. More... | |
| void | storeInnerSolution (sgpp::base::DataVector &alpha, size_t maxCGIterations, double epsilonCG, std::string tFilename) |
| Routine to export the solution of the inner system which has been calculated by Up/Down scheme. More... | |
| virtual | ~PoissonEquationSolver () |
| Std-Destructor of the solver. More... | |
Public Member Functions inherited from sgpp::pde::EllipticPDESolver | |
| EllipticPDESolver () | |
| Std-Constructor of the solver. More... | |
| virtual | ~EllipticPDESolver () |
| Std-Destructor of the solver. More... | |
Public Member Functions inherited from sgpp::pde::PDESolver | |
| void | coarsenInitialGridSurplus (sgpp::base::DataVector &alpha, double dThreshold) |
| Coarsens a grid by taking the grid's coefficients into account. More... | |
| void | deleteGrid () |
| deletes the grid created within that solver More... | |
| void | evaluateCuboid (sgpp::base::DataVector &alpha, sgpp::base::DataVector &FunctionValues, sgpp::base::DataMatrix &EvaluationPoints) |
| Evaluates the sparse grid's function given by the stored grid and the alpha coefficients. More... | |
| double | evaluatePoint (sgpp::base::DataVector &evalPoint, sgpp::base::DataVector &alpha) |
| Determines the value of the function in the d-dimensional space. More... | |
| std::string | getGrid () const |
| gets the a string the describes the grid which is currently used to solve More... | |
| size_t | getNumberDimensions () const |
| use this the determine the number of dimensions that are currently used in the solver. More... | |
| size_t | getNumberGridPoints () const |
| use this to determine the number of grid points, used to solve the current problem More... | |
| size_t | getNumberInnerGridPoints () const |
| use this to determine the number of inner grid points, used to solve the current problem More... | |
| PDESolver () | |
| Std-Constructor of the solver. More... | |
| virtual void | printGrid (sgpp::base::DataVector &alpha, size_t PointesPerDimension, std::string tfilename) const |
| This is some kind of debug functionality. More... | |
| virtual void | printGridDomain (sgpp::base::DataVector &alpha, size_t PointesPerDimension, sgpp::base::BoundingBox &GridArea, std::string tfilename) const |
| This is some kind of debug functionality. More... | |
| virtual void | printLevelIndexGrid (std::string tfilename) const |
| Prints the level,index pairs of the grid for each Gridpoint to a file. More... | |
| virtual void | printSparseGrid (sgpp::base::DataVector &alpha, std::string tfilename, bool bSurplus) const |
| Prints the sgpp::base::Grid Points of the Sparse sgpp::base::Grid either with their node basis value or their hierarchical surplus. More... | |
| virtual void | printSparseGridExpTransform (sgpp::base::DataVector &alpha, std::string tfilename, bool bSurplus) const |
| Prints the sgpp::base::Grid Points of the Sparse sgpp::base::Grid either with their node basis value or their hierarchical surplus. More... | |
| void | refineInitialGridSurplus (sgpp::base::DataVector &alpha, int numRefinePoints, double dThreshold) |
| Refines a grid by taking the grid's coefficients into account. More... | |
| void | refineInitialGridSurplusSubDomain (sgpp::base::DataVector &alpha, int numRefinePoints, double dThreshold, std::vector< double > &norm_mu, std::vector< double > &norm_sigma) |
| Refines a grid by taking the grid's coefficients into account. More... | |
| void | refineInitialGridSurplusToMaxLevel (sgpp::base::DataVector &alpha, double dThreshold, sgpp::base::level_t maxLevel) |
| Refines a grid by taking the grid's coefficients into account. More... | |
| void | refineInitialGridSurplusToMaxLevelSubDomain (sgpp::base::DataVector &alpha, double dThreshold, sgpp::base::level_t maxLevel, std::vector< double > &norm_mu, std::vector< double > &norm_sigma) |
| Refines a grid by taking the grid's coefficients into account. More... | |
| void | setGrid (const std::string &serializedGrid) |
| Sets the grid used in this BlackScholes Solver by an given serialized string of the grid. More... | |
| virtual | ~PDESolver () |
| Std-Destructor of the solver. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from sgpp::pde::PDESolver | |
| virtual void | getGridNormalDistribution (sgpp::base::DataVector &alpha, std::vector< double > &norm_mu, std::vector< double > &norm_sigma) |
| This function calculates for every grid point the value of a normal distribution given by norm_mu and norm_sigma. More... | |
Protected Attributes inherited from sgpp::pde::PDESolver | |
| bool | bGridConstructed |
| stores if the grid was created inside the solver More... | |
| size_t | dim |
| the dimension of the grid More... | |
| int | levels |
| the number of levels used for an regular grid More... | |
| sgpp::base::BoundingBox * | myBoundingBox |
| Stores Pointer to the sgpp::base::Grid's Bounding Box. More... | |
| sgpp::base::Grid * | myGrid |
| The Sparse sgpp::base::Grid needed in this classificator. More... | |
| sgpp::base::GridStorage * | myGridStorage |
| Stores Pointer to the Girs's Storage. More... | |
This class provides a simple-to-use solver of the multi dimensional Poisson Equation on Sparse Grids.
The class's aim is, to hide all complex details of solving the Poisson Equation on Sparse Grids!
| sgpp::pde::PoissonEquationSolver::PoissonEquationSolver | ( | ) |
Std-Constructor of the solver.
References sgpp::pde::PDESolver::bGridConstructed.
|
virtual |
Std-Destructor of the solver.
|
virtual |
Use this routine the construct a regular grid to solve a PDE.
| myBoundingBox | reference to a bounding box that describes the grid |
| level | number of the regular's grid levels |
Implements sgpp::pde::PDESolver.
References sgpp::pde::PDESolver::bGridConstructed, sgpp::pde::PDESolver::dim, sgpp::base::Grid::getBoundingBox(), sgpp::base::BoundingBox::getDimension(), sgpp::base::Grid::getGenerator(), sgpp::base::Grid::getStorage(), level, sgpp::pde::PDESolver::levels, sgpp::pde::PDESolver::myBoundingBox, sgpp::pde::PDESolver::myGrid, sgpp::pde::PDESolver::myGridStorage, and sgpp::base::GridGenerator::regular().
| void sgpp::pde::PoissonEquationSolver::initGridWithExpHeat | ( | sgpp::base::DataVector & | alpha, |
| double | factor = 1.0 |
||
| ) |
Inits the grid with a heat distribution based on the e-function.
The e-function is shifted in that way the right boundary values becomes 1 (in case of factor = 1)
| alpha | reference to the coefficient's vector |
| factor | a constant factor used to enlarge the exp-functions input parameter |
References sgpp::pde::PDESolver::bGridConstructed, sgpp::op_factory::createOperationHierarchisation(), sgpp::pde::PDESolver::dim, sgpp::base::OperationHierarchisation::doHierarchisation(), sgpp::base::BoundingBox::getBoundary(), sgpp::base::Grid::getBoundingBox(), sgpp::base::HashGridStorage::getCoordinates(), sgpp::base::HashGridStorage::getPoint(), sgpp::base::Grid::getSize(), python.statsfileInfo::i, python.utils.statsfile2gnuplot::j, sgpp::base::BoundingBox1D::leftBoundary, sgpp::pde::PDESolver::myBoundingBox, sgpp::pde::PDESolver::myGrid, sgpp::pde::PDESolver::myGridStorage, sgpp::base::BoundingBox1D::rightBoundary, and analyse_erg::tmp.
| void sgpp::pde::PoissonEquationSolver::initGridWithExpHeatFullDomain | ( | sgpp::base::DataVector & | alpha, |
| double | factor = 1.0 |
||
| ) |
Inits the grid with a heat distribution based on the e-function.
The e-function is shifted in that way the right boundary values becomes 1 (in case of factor = 1)
| alpha | reference to the coefficient's vector |
| factor | a constant factor used to enlarge the exp-functions input parameter |
References sgpp::pde::PDESolver::bGridConstructed, sgpp::op_factory::createOperationHierarchisation(), sgpp::pde::PDESolver::dim, sgpp::base::OperationHierarchisation::doHierarchisation(), sgpp::base::BoundingBox::getBoundary(), sgpp::base::Grid::getBoundingBox(), sgpp::base::HashGridStorage::getCoordinates(), sgpp::base::HashGridStorage::getPoint(), sgpp::base::Grid::getSize(), python.statsfileInfo::i, python.utils.statsfile2gnuplot::j, sgpp::pde::PDESolver::myGrid, sgpp::pde::PDESolver::myGridStorage, and analyse_erg::tmp.
| void sgpp::pde::PoissonEquationSolver::initGridWithSmoothHeat | ( | sgpp::base::DataVector & | alpha, |
| double | mu, | ||
| double | sigma, | ||
| double | factor | ||
| ) |
Inits the grid with a smooth heat distribution (based on a std-normal distribution) on its boundaries.
Coefficients of inner grid points are set to zero since an elliptical PDE is solved
| alpha | reference to the coefficients vector |
| mu | the exspected value of the normal distribution |
| sigma | the sigma of the normal distribution |
| factor | a factor that is used to stretch the function values |
References sgpp::pde::PDESolver::bGridConstructed, sgpp::op_factory::createOperationHierarchisation(), sgpp::pde::PDESolver::dim, sgpp::base::OperationHierarchisation::doHierarchisation(), sgpp::base::BoundingBox::getBoundary(), sgpp::base::HashGridStorage::getCoordinates(), sgpp::base::HashGridStorage::getPoint(), sgpp::base::Grid::getSize(), python.statsfileInfo::i, python.utils.statsfile2gnuplot::j, sgpp::base::BoundingBox1D::leftBoundary, mu, sgpp::pde::PDESolver::myBoundingBox, sgpp::pde::PDESolver::myGrid, sgpp::pde::PDESolver::myGridStorage, sgpp::base::BoundingBox1D::rightBoundary, and analyse_erg::tmp.
| void sgpp::pde::PoissonEquationSolver::initGridWithSmoothHeatFullDomain | ( | sgpp::base::DataVector & | alpha, |
| double | mu, | ||
| double | sigma, | ||
| double | factor | ||
| ) |
Inits the grid with a smooth heat distribution (based on a std-normal distribution) on its boundaries.
Coefficients of inner grid points aren't set to zero since they are used to hint an adaptive refinement of the grid BEFORE solving the PDE.
| alpha | reference to the coefficients vector |
| mu | the exspected value of the normal distribution |
| sigma | the sigma of the normal distribution |
| factor | a factor that is used to stretch the function values |
References sgpp::pde::PDESolver::bGridConstructed, sgpp::op_factory::createOperationHierarchisation(), sgpp::pde::PDESolver::dim, sgpp::base::OperationHierarchisation::doHierarchisation(), sgpp::base::HashGridStorage::getCoordinates(), sgpp::base::HashGridStorage::getPoint(), sgpp::base::Grid::getSize(), python.statsfileInfo::i, python.utils.statsfile2gnuplot::j, mu, sgpp::pde::PDESolver::myGrid, sgpp::pde::PDESolver::myGridStorage, and analyse_erg::tmp.
| void sgpp::pde::PoissonEquationSolver::initScreen | ( | ) |
Inits the screen object.
References sgpp::base::ScreenOutput::writeTitle().
|
virtual |
abstract method to solve an elliptic PDE.
All solver of elliptic PDEs have to implement this method.
| alpha | the coefficients of the Sparse Gird's basis functions will be in this vector after solving |
| rhs | the right hand side of the SLE |
| maxCGIterations | the maximum of interation in the CG solver |
| epsilonCG | the epsilon used in the CG |
| verbose | enables verbose output during solving |
Implements sgpp::pde::EllipticPDESolver.
References sgpp::pde::OperationEllipticPDESolverSystemDirichlet::generateRHS(), sgpp::pde::OperationEllipticPDESolverSystemDirichlet::getGridCoefficientsForCG(), sgpp::pde::OperationEllipticPDESolverSystem::getNumGridPointsComplete(), sgpp::pde::OperationEllipticPDESolverSystem::getNumGridPointsInner(), sgpp::pde::OperationEllipticPDESolverSystemDirichlet::getSolutionBoundGrid(), sgpp::pde::PDESolver::myGrid, sgpp::solver::ConjugateGradients::solve(), sgpp::base::SGppStopwatch::start(), and sgpp::base::SGppStopwatch::stop().
| void sgpp::pde::PoissonEquationSolver::storeInnerRHS | ( | sgpp::base::DataVector & | alpha, |
| std::string | tFilename | ||
| ) |
Routine to export the RHS of the inner system which has to be solved in order to solve the Poisson equation.
| alpha | the start solution |
| tFilename | file into which the rhs is written |
References sgpp::pde::OperationEllipticPDESolverSystemDirichlet::generateRHS(), sgpp::base::DataVector::get(), sgpp::base::DataVector::getSize(), python.statsfileInfo::i, sgpp::pde::PDESolver::myGrid, sgpp::base::SGppStopwatch::start(), and sgpp::base::SGppStopwatch::stop().
| void sgpp::pde::PoissonEquationSolver::storeInnerSolution | ( | sgpp::base::DataVector & | alpha, |
| size_t | maxCGIterations, | ||
| double | epsilonCG, | ||
| std::string | tFilename | ||
| ) |
Routine to export the solution of the inner system which has been calculated by Up/Down scheme.
| alpha | the start solution |
| maxCGIterations | the maximum of interation in the CG solver |
| epsilonCG | the epsilon used in the C |
| tFilename | file into which the rhs is written |
References sgpp::pde::OperationEllipticPDESolverSystemDirichlet::generateRHS(), sgpp::base::DataVector::get(), sgpp::pde::OperationEllipticPDESolverSystemDirichlet::getGridCoefficientsForCG(), sgpp::base::DataVector::getSize(), python.statsfileInfo::i, sgpp::pde::PDESolver::myGrid, and sgpp::solver::ConjugateGradients::solve().