|
HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
Implementation of the HHO scheme with full gradient for variable diffusion equations. More...
Classes | |
| struct | HArDCore2D::FullGradientDiffusion |
| Assemble a diffusion problem. More... | |
Typedefs | |
| typedef Eigen::SparseMatrix< double > | HArDCore2D::FullGradientDiffusion::SystemMatrixType |
| typedef std::function< double(const VectorRd &)> | HArDCore2D::FullGradientDiffusion::ForcingTermType |
| typedef std::function< double(const VectorRd &)> | HArDCore2D::FullGradientDiffusion::SolutionType |
| typedef std::function< VectorRd(const VectorRd &)> | HArDCore2D::FullGradientDiffusion::SolutionGradientType |
| typedef IntegralWeight | HArDCore2D::FullGradientDiffusion::PermeabilityType |
Functions | |
| HArDCore2D::FullGradientDiffusion::FullGradientDiffusion (const HHOSpace &hho_space, const BoundaryConditions &BC, bool use_threads, std::ostream &output=std::cout) | |
| Constructor. | |
| void | HArDCore2D::FullGradientDiffusion::assembleLinearSystem (const ForcingTermType &f, const PermeabilityType &kappa, const SolutionType &u, Eigen::VectorXd &UDir) |
| Assemble the global system | |
| size_t | HArDCore2D::FullGradientDiffusion::numSCDOFs () const |
| Returns the number of statically condensed DOFs (the cell DOFs) | |
| size_t | HArDCore2D::FullGradientDiffusion::numDirDOFs () const |
| Returns the number of Dirichlet DOFs. | |
| size_t | HArDCore2D::FullGradientDiffusion::numSkeletalDOFs () const |
| Returns the number of DOFs after SC but before eliminating Dirichlet DOFs. | |
| size_t | HArDCore2D::FullGradientDiffusion::sizeSystem () const |
| Returns the size of the final system, after application of SC and removal of Dirichlet BCs. | |
| const HHOSpace & | HArDCore2D::FullGradientDiffusion::hhospace () const |
| Returns the HHO space. | |
| const Mesh & | HArDCore2D::FullGradientDiffusion::mesh () const |
| Returns the mesh. | |
| const SystemMatrixType & | HArDCore2D::FullGradientDiffusion::systemMatrix () const |
| Returns the linear system matrix. | |
| SystemMatrixType & | HArDCore2D::FullGradientDiffusion::systemMatrix () |
| Returns the linear system matrix. | |
| const Eigen::VectorXd & | HArDCore2D::FullGradientDiffusion::systemVector () const |
| Returns the linear system right-hand side vector. | |
| Eigen::VectorXd & | HArDCore2D::FullGradientDiffusion::systemVector () |
| Returns the linear system right-hand side vector. | |
| const double & | HArDCore2D::FullGradientDiffusion::stabilizationParameter () const |
| Returns the stabilization parameter (scaling) | |
| double & | HArDCore2D::FullGradientDiffusion::stabilizationParameter () |
| Returns the stabilization parameter. | |
| const SystemMatrixType & | HArDCore2D::FullGradientDiffusion::scMatrix () const |
| Returns the static condensation recovery operator. | |
| Eigen::VectorXd & | HArDCore2D::FullGradientDiffusion::scVector () |
| Returns the static condensation rhs. | |
| std::vector< double > | HArDCore2D::FullGradientDiffusion::computeEnergyNorms (const std::vector< Eigen::VectorXd > &list_dofs) const |
| Compute the discrete energy norm of a family of vectors representing the dofs. | |
Variables | |
| static const double | HArDCore2D::PI = boost::math::constants::pi<double>() |
| static const VectorRd | HArDCore2D::vec_a = VectorRd(1.,2.) |
| static FullGradientDiffusion::SolutionType | HArDCore2D::linear_u |
| static FullGradientDiffusion::SolutionGradientType | HArDCore2D::linear_gradu |
| static FullGradientDiffusion::ForcingTermType | HArDCore2D::linear_f |
| static FullGradientDiffusion::PermeabilityType | HArDCore2D::linear_kappa = FullGradientDiffusion::PermeabilityType(1.) |
| static FullGradientDiffusion::SolutionType | HArDCore2D::trigonometric_u |
| static FullGradientDiffusion::SolutionGradientType | HArDCore2D::trigonometric_gradu |
| static FullGradientDiffusion::ForcingTermType | HArDCore2D::trigonometric_f |
| static FullGradientDiffusion::PermeabilityType | HArDCore2D::trigonometric_kappa = FullGradientDiffusion::PermeabilityType(1.) |
Implementation of the HHO scheme with full gradient for variable diffusion equations.
| typedef std::function<double(const VectorRd &)> HArDCore2D::FullGradientDiffusion::ForcingTermType |
| typedef std::function<VectorRd(const VectorRd &)> HArDCore2D::FullGradientDiffusion::SolutionGradientType |
| typedef std::function<double(const VectorRd &)> HArDCore2D::FullGradientDiffusion::SolutionType |
| typedef Eigen::SparseMatrix<double> HArDCore2D::FullGradientDiffusion::SystemMatrixType |
| void FullGradientDiffusion::assembleLinearSystem | ( | const ForcingTermType & | f, |
| const PermeabilityType & | kappa, | ||
| const SolutionType & | u, | ||
| Eigen::VectorXd & | UDir | ||
| ) |
Assemble the global system
| f | Forcing term |
| kappa | Permeability |
| u | Exact solution for boundary condition |
| UDir | Vector filled in by Dirichlets BCs |
| std::vector< double > FullGradientDiffusion::computeEnergyNorms | ( | const std::vector< Eigen::VectorXd > & | list_dofs | ) | const |
Compute the discrete energy norm of a family of vectors representing the dofs.
| list_dofs | The list of vectors representing the dofs |
| FullGradientDiffusion::FullGradientDiffusion | ( | const HHOSpace & | hho_space, |
| const BoundaryConditions & | BC, | ||
| bool | use_threads, | ||
| std::ostream & | output = std::cout |
||
| ) |
Constructor.
| hho_space | HHO space and operators |
| BC | Boundary conditions |
| use_threads | True for parallel execution, false for sequential execution |
| output | Output stream to print status messages |
|
inline |
Returns the HHO space.
|
inline |
Returns the mesh.
|
inline |
Returns the number of Dirichlet DOFs.
|
inline |
Returns the number of statically condensed DOFs (the cell DOFs)
|
inline |
Returns the number of DOFs after SC but before eliminating Dirichlet DOFs.
|
inline |
Returns the static condensation recovery operator.
|
inline |
Returns the static condensation rhs.
|
inline |
Returns the size of the final system, after application of SC and removal of Dirichlet BCs.
|
inline |
Returns the stabilization parameter.
|
inline |
Returns the stabilization parameter (scaling)
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system right-hand side vector.
|
inline |
Returns the linear system right-hand side vector.
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |