Inherits ProblemStatBase, and StandardProblemIterationAdaptor< ProblemStat< Traits > >.
Public Types | |
using | GlobalBasis = typename Traits::GlobalBasis |
using | GridView = typename GlobalBasis::GridView |
using | Grid = AdaptiveGrid_t< typename GridView::Grid > |
using | Element = typename GridView::template Codim< 0 >::Entity |
using | WorldVector = typename Element::Geometry::GlobalCoordinate |
using | WorldMatrix = FieldMatrix< typename WorldVector::field_type, WorldVector::dimension, WorldVector::dimension > |
using | value_type = typename Traits::CoefficientType |
using | Mat = typename Traits::Backend::template Matrix< GlobalBasis, GlobalBasis >::template Impl< value_type > |
using | Vec = typename Traits::Backend::template Vector< GlobalBasis >::template Impl< value_type > |
using | LinearSolverInterface = AMDiS::LinearSolverInterface< Mat, Vec, Vec > |
using | LinearSolver = AMDiS::LinearSolver< Mat, Vec, Vec > |
using | SystemMatrix = BiLinearForm< GlobalBasis, GlobalBasis, value_type, typename Traits::Backend > |
using | SystemVector = LinearForm< GlobalBasis, value_type, typename Traits::Backend > |
using | SolutionVector = DOFVector< GlobalBasis, value_type, typename Traits::Backend > |
Public Member Functions | |
ProblemStat (std::string const &name) | |
Constructor. Takes the name of the problem that is used to access values corresponding to this problem in the parameter file. | |
template<class Grid_ > | |
ProblemStat (std::string const &name, Grid_ &&grid) | |
template<class Grid_ , class Basis_ , class B_ = Underlying_t<Basis_>, REQUIRES(Concepts::GlobalBasis< B_ >) > | |
ProblemStat (std::string const &name, Grid_ &&grid, Basis_ &&globalBasis) | |
Constructor taking a grid and basis. Wraps both in shared pointers. | |
template<class Grid_ , class PBF_ , class GV_ = typename Underlying_t<Grid_>::LeafGridView, REQUIRES(Concepts::PreBasisFactory< PBF_, GV_ >) > | |
ProblemStat (std::string const &name, Grid_ &&grid, PBF_ const &preBasisFactory) | |
Constructor taking a grid and pre-basis factory to create a global basis on the fly. | |
void | initialize (Flag initFlag, Self *adoptProblem=nullptr, Flag adoptFlag=INIT_NOTHING) |
Initialisation of the problem. More... | |
void | restore (Flag initFlag) |
Read the grid and solution from backup files and initialize the problem. More... | |
void | addConstraint (BoundaryCondition< SystemMatrix, SolutionVector, SystemVector > constraint) |
Flag | oneIteration (AdaptInfo &adaptInfo, Flag toDo=FULL_ITERATION) override |
Implementation of StandardProblemIteration::oneIteration. More... | |
void | buildAfterAdapt (AdaptInfo &adaptInfo, Flag flag, bool asmMatrix=true, bool asmVector=true) override |
Implementation of ProblemStatBase::buildAfterCoarse. More... | |
void | assemble (AdaptInfo &adaptInfo) |
Assemble the linear system by calling buildAfterAdapt with asmMatrix and asmVector set to true. | |
void | solve (AdaptInfo &adaptInfo, bool createMatrixData=true, bool storeMatrixData=false) override |
Implementation of ProblemStatBase::solve. More... | |
void | estimate (AdaptInfo &) override |
Implementation of ProblemStatBase::estimate. More... | |
Flag | adaptGrid (AdaptInfo &adaptInfo) override |
Implementation of ProblemStatBase::refineMesh. More... | |
Flag | markElements (AdaptInfo &adaptInfo) override |
Implementation of ProblemStatBase::markElements. More... | |
Flag | globalCoarsen (int n) override |
Uniform global grid coarsening by up to n level. More... | |
Flag | globalRefine (int n) override |
Uniform global refinement by n level. More... | |
void | writeFiles (AdaptInfo &adaptInfo, bool force=false) |
Writes output files. If force=true write even if timestep out of write rhythm. | |
std::string const & | name () const override |
Implementation of ProblemStatBase::name. More... | |
std::shared_ptr< Grid > | grid () |
Return the grid_. | |
std::shared_ptr< Grid const > | grid () const |
GridView | gridView () const |
Return the gridView of the basis. | |
std::shared_ptr< BoundaryManager< Grid > > | boundaryManager () |
Return the boundary manager to identify boundary segments. | |
std::shared_ptr< BoundaryManager< Grid > const > | boundaryManager () const |
std::shared_ptr< GlobalBasis > | globalBasis () |
Return the globalBasis_. | |
std::shared_ptr< GlobalBasis const > | globalBasis () const |
std::shared_ptr< LinearSolverInterface > | solver () |
Return a reference to the linear solver, linearSolver. | |
std::shared_ptr< LinearSolverInterface const > | solver () const |
std::shared_ptr< SystemMatrix > | systemMatrix () |
Returns a reference to system-matrix, systemMatrix_. | |
std::shared_ptr< SystemMatrix const > | systemMatrix () const |
std::shared_ptr< SolutionVector > | solutionVector () |
Returns a reference to the solution vector, solution_. | |
std::shared_ptr< SolutionVector const > | solutionVector () const |
std::shared_ptr< SystemVector > | rhsVector () |
Return a reference to the rhs system-vector, rhs. | |
std::shared_ptr< SystemVector const > | rhsVector () const |
template<class Range = void, class... Indices> | |
auto | solution (Indices... ii) |
Return a mutable view to a solution component. More... | |
template<class Range = void, class... Indices> | |
auto | solution (Indices... ii) const |
Return a const view to a solution component. More... | |
template<class Solver_ > | |
void | setSolver (Solver_ &&solver) |
Set a new linear solver for the problem. | |
template<class Grid_ > | |
void | setGrid (Grid_ &&grid) |
template<class Marker_ > | |
void | addMarker (Marker_ &&m) |
Store the shared_ptr and the name of the marker in the problem. More... | |
void | removeMarker (std::string name) |
Remove a marker with the given name from the problem. | |
void | removeMarker (Marker< Grid > const &marker) |
Remove a marker from the problem. | |
template<class FileWriter_ > | |
void | addFileWriter (FileWriter_ &&f) |
Add another filewriter to the problem. | |
void | clearFileWriter () |
Deletes all filewriters. | |
template<class Operator , class RowTreePath = RootTreePath, class ColTreePath = RootTreePath> | |
void | addMatrixOperator (Operator const &op, RowTreePath row={}, ColTreePath col={}) |
Add an operator to A. More... | |
template<class Operator , class RowTreePath = RootTreePath, class ColTreePath = RootTreePath> | |
void | addMatrixOperator (BoundaryType b, Operator const &op, RowTreePath row={}, ColTreePath col={}) |
Operator evaluated on the boundary of the domain with boundary index b More... | |
template<class Operator , class TreePath = RootTreePath> | |
void | addVectorOperator (Operator const &op, TreePath path={}) |
Add an operator to rhs. More... | |
template<class Operator , class TreePath = RootTreePath> | |
void | addVectorOperator (BoundaryType b, Operator const &op, TreePath path={}) |
Operator evaluated on the boundary of the domain with boundary index b More... | |
template<class Predicate , class RowTreePath , class ColTreePath , class Values > | |
void | addDirichletBC (Predicate const &predicate, RowTreePath row, ColTreePath col, Values const &values) |
Add boundary conditions to the system. More... | |
template<class RowTreePath , class ColTreePath , class Values > | |
void | addDirichletBC (BoundaryType id, RowTreePath row, ColTreePath col, Values const &values) |
template<class Identifier , class Values > | |
void | addDirichletBC (Identifier &&id, Values &&values) |
void | addPeriodicBC (BoundaryType id, WorldMatrix const &A, WorldVector const &b) |
virtual Flag | markElements (AdaptInfo &adaptInfo)=0 |
Marks mesh elements for refinement and coarsening. More... | |
virtual void | buildAfterAdapt (AdaptInfo &adaptInfo, Flag flag, bool assembleMatrix, bool assembleVector)=0 |
Assembling of system matrices and vectors after coarsening. By the last two parameters, assembling can be restricted to either matrices or vectors only. More... | |
virtual Flag | adaptGrid (AdaptInfo &adaptInfo)=0 |
Refinement/coarsening of the grid. More... | |
virtual Flag | globalCoarsen (int n)=0 |
virtual Flag | globalRefine (int n)=0 |
virtual void | solve (AdaptInfo &adaptInfo, bool createMatrixData=true, bool storeMatrixData=false)=0 |
Solves the assembled system. The result is an approximative solution. The last two boolean arguments can be used to controll successive solutions of systems with the same matrix. More... | |
virtual void | estimate (AdaptInfo &adaptInfo)=0 |
A posteriori error estimation of the calculated solution. Should store a local error estimation at each elements leaf data and return the total error sum. More... | |
virtual std::string const & | name () const =0 |
Returns the name of the problem. More... | |
Public Member Functions inherited from StandardProblemIterationAdaptor< ProblemStat< Traits > > | |
StandardProblemIterationAdaptor (StandardProblemIterationAdaptor const &) | |
StandardProblemIterationAdaptor (StandardProblemIterationAdaptor &&) | |
Public Member Functions inherited from StandardProblemIteration | |
StandardProblemIteration (ProblemStatBase &prob) | |
constructor | |
void | beginIteration (AdaptInfo &adaptInfo) override |
Implementation of ProblemIterationIterface::beginIteration() More... | |
Flag | oneIteration (AdaptInfo &adaptInfo, Flag toDo) override |
Implementation of ProblemIterationInterface::oneIteration() More... | |
void | endIteration (AdaptInfo &adaptInfo) override |
Implementation of ProblemIterationInterface::endIteration() More... | |
std::string const & | name () const override |
Returns the name of the problem. More... | |
int | numProblems () const override |
Returns number of managed problems. More... | |
ProblemStatBase & | problem (int number=0) override |
Return the managed ProblemStat problem, by number. More... | |
ProblemStatBase & | problem (std::string const &name) override |
Return the managed ProblemStat problem, by name. More... | |
virtual void | beginIteration (AdaptInfo &) |
Called before each adaption loop iteration. More... | |
virtual Flag | oneIteration (AdaptInfo &adaptInfo, Flag toDo=FULL_ITERATION)=0 |
Determines the execution order of the single adaption steps. If adapt is true, mesh adaption will be performed. This allows to avoid mesh adaption, e.g. in timestep adaption loops of timestep adaptive strategies. More... | |
virtual void | endIteration (AdaptInfo &) |
Called after each adaption loop iteration. More... | |
virtual int | numProblems () const =0 |
Returns number of managed problems. More... | |
virtual ProblemStatBase & | problem (int number=0)=0 |
Returns the problem with the given number. If only one problem is managed by this master problem, the number hasn't to be given. More... | |
virtual ProblemStatBase & | problem (std::string const &name)=0 |
Returns the problem with the given name. More... | |
virtual std::string const & | name () const =0 |
Returns the name of the problem. More... | |
Static Public Attributes | |
static constexpr int | dim = Grid::dimension |
Dimension of the grid. | |
static constexpr int | dow = Grid::dimensionworld |
Dimension of the world. | |
Protected Member Functions | |
void | createGlobalBasis () |
void | createGrid () |
void | createMatricesAndVectors () |
void | createSolver () |
void | createMarker () |
void | createFileWriter () |
void | adoptGlobalBasis (std::shared_ptr< GlobalBasis > globalBasis) |
void | adoptGrid (std::shared_ptr< Grid > const &grid, std::shared_ptr< BoundaryManager< Grid > > const &boundaryManager) |
void | adoptGrid (std::shared_ptr< Grid > const &grid) |
void | adoptGrid (std::shared_ptr< typename Grid::HostGrid > const &hostGrid) |
Protected Member Functions inherited from StandardProblemIteration | |
Flag | buildAndAdapt (AdaptInfo &adaptInfo, Flag toDo) |
Nested assemblage and mesh adaption. | |
Protected Attributes | |
std::string | name_ |
Name of this problem. | |
std::shared_ptr< Grid > | grid_ |
Grid of this problem. | |
std::string | gridName_ = "mesh" |
Name of the grid. | |
std::shared_ptr< BoundaryManager< Grid > > | boundaryManager_ |
Management of boundary conditions. | |
std::shared_ptr< GlobalBasis > | globalBasis_ |
FE space of this problem. | |
std::list< std::shared_ptr< FileWriterInterface > > | filewriter_ |
A FileWriter object. | |
std::map< std::string, std::shared_ptr< Marker< Grid > > > | marker_ |
Pointer to the adaptation markers. | |
std::shared_ptr< LinearSolverInterface > | linearSolver_ |
Pointer to the estimators for this problem. More... | |
std::shared_ptr< SystemMatrix > | systemMatrix_ |
Matrix that is filled during assembling. | |
std::shared_ptr< SolutionVector > | solution_ |
Vector with the solution components. | |
std::shared_ptr< SystemVector > | rhs_ |
std::map< std::string, std::vector< double > > | estimates_ |
std::list< BoundaryCondition< SystemMatrix, SolutionVector, SystemVector > > | constraints_ |
List of constraints to apply to matrix, solution and rhs. | |
Protected Attributes inherited from StandardProblemIteration | |
ProblemStatBase & | problem_ |
The problem to solve. | |
Friends | |
class | ProblemInstat< Traits > |
|
inline |
Constructor taking additionally a grid that is used instead of the default created grid, ProblemStat
References ProblemStat< Traits >::grid().
Implementation of ProblemStatBase::refineMesh.
Implements ProblemStatBase.
void addDirichletBC | ( | Predicate const & | predicate, |
RowTreePath | row, | ||
ColTreePath | col, | ||
Values const & | values | ||
) |
Add boundary conditions to the system.
Dirichlet boundary condition Enforce Dirichlet boundary values for the solution vector on boundary regions identified by the predicate.
predicate | Functor bool(WorldVector) returning true for all DOFs on the boundary that should be assigned a value. |
row | TreePath identifying the sub-basis in the global basis tree corresponding to the row basis. |
col | TreePath identifying the sub-basis in the global basis tree corresponding to the column basis. |
values | Functor Range(WorldVector) or any GridFunction that is evaluated in the DOFs identified by the predicate. |
Example:
References AMDiS::Concepts::Functor, AMDiS::makeGridFunction(), and AMDiS::Concepts::Predicate.
|
inline |
Store the shared_ptr and the name of the marker in the problem.
Note: multiple markers can be added but must have different names
|
inline |
Operator evaluated on the boundary of the domain with boundary index b
Adds an operator to the list of boundary operators to be assembled in quadrature points on the boundary intersections.
b | Boundary identifier where to assemble this operator. Can be constructed from an integer. |
op | A (pre-) local operator, |
row | TreePath identifying the sub-basis in the global basis tree corresponding to the row basis. |
col | TreePath identifying the sub-basis in the global basis tree corresponding to the column basis. |
Example:
|
inline |
Add an operator to A.
Operator evaluated on the whole element Adds an operator to the list of element operators to be assembled in quadrature points inside the element.
op | A (pre-) local operator, |
row | TreePath identifying the sub-basis in the global basis tree corresponding to the row basis. |
col | TreePath identifying the sub-basis in the global basis tree corresponding to the column basis. |
Example:
void addPeriodicBC | ( | BoundaryType | id, |
WorldMatrix const & | A, | ||
WorldVector const & | b | ||
) |
Add a periodic boundary conditions to the system, by specifying a face transformation y = A*x + b of coordinates. We assume, that A is orthonormal.
|
inline |
Operator evaluated on the boundary of the domain with boundary index b
Adds an operator to the list of boundary operators to be assembled in quadrature points on the boundary intersections.
b | Boundary identifier where to assemble this operator. Can be constructed from an integer. |
op | A (pre-) local operator, |
path | TreePath identifying the sub-basis in the global basis tree corresponding to the row basis. |
Example:
|
inline |
Add an operator to rhs.
Operator evaluated on the whole element Adds an operator to the list of element operators to be assembled in quadrature points inside the element.
op | A (pre-) local operator, |
path | TreePath identifying the sub-basis in the global basis tree corresponding to the row basis. |
Example:
|
overridevirtual |
Implementation of ProblemStatBase::buildAfterCoarse.
Implements ProblemStatBase.
|
inlineoverridevirtual |
Implementation of ProblemStatBase::estimate.
Implements ProblemStatBase.
Uniform global grid coarsening by up to n level.
Implements ProblemStatBase.
Uniform global refinement by n level.
Implements ProblemStatBase.
References ProblemStat< Traits >::globalRefine().
Referenced by ProblemStat< Traits >::globalRefine().
Initialisation of the problem.
Parameters read in initialize()
[GRID_NAME]->global refinements
: nr of initial global refinements References ProblemStat< Traits >::boundaryManager_, ProblemStat< Traits >::estimates_, ProblemStat< Traits >::globalBasis_, ProblemStat< Traits >::grid_, Flag::isSet(), ProblemStat< Traits >::linearSolver_, ProblemStat< Traits >::marker_, ProblemStat< Traits >::rhs_, ProblemStat< Traits >::solution_, and ProblemStat< Traits >::systemMatrix_.
Implementation of ProblemStatBase::markElements.
Implements ProblemStatBase.
|
inlineoverridevirtual |
Implementation of ProblemStatBase::name.
Implements ProblemStatBase.
Implementation of StandardProblemIteration::oneIteration.
Implements ProblemIterationInterface.
void restore | ( | Flag | initFlag | ) |
Read the grid and solution from backup files and initialize the problem.
Parameters read in restore() for problem with name 'PROB'
[PROB]->restore->grid
: name of the grid backup file[PROB]->restore->solution
: name of the solution backup file References Flag::isSet().
|
inline |
Set the grid. Stores pointer and initializes feSpaces matrices and vectors, as well as markers and file-writers. If grid is given as reference, wrap it into a non-destroying shared_ptr
|
inline |
Return a mutable view to a solution component.
Range | The range type return by evaluating the view in coordinates. If not specified, it is automatically selected using RangeType_t template. |
|
inline |
Return a const view to a solution component.
Range | The range type return by evaluating the view in coordinates. If not specified, it is automatically selected using RangeType_t template. |
|
overridevirtual |
Implementation of ProblemStatBase::solve.
Implements ProblemStatBase.
|
protected |
A vector with the local element error estimates for each node in the basis tree, indexed by [to_string(treePath)][element index]
Referenced by ProblemStat< Traits >::initialize().
|
protected |
Pointer to the estimators for this problem.
An object of the linearSolver Interface
Referenced by ProblemStat< Traits >::initialize().
|
protected |
Vector (load-vector) corresponding to the right-hand side of the equation, filled during assembling
Referenced by ProblemStat< Traits >::initialize().