Class providing helper methods and assemble and solve routines for general 3D HHO schemes.
More...
|
| typedef std::function< Eigen::MatrixXd(Cell *, ElementQuad &)> | MatrixFType |
| | type for the global operator as a function of a Cell and an ElementQuad
|
| |
| typedef std::function< Eigen::VectorXd(Cell *, ElementQuad &)> | VectorFType |
| | type for the load vector as a function of a Cell and an ElementQuad
|
| |
|
| | HHO3D::HHO3D (HybridCore &, const size_t, const size_t, const bool use_threads=true, size_t doeT=0, size_t doeF=0) |
| | Class constructor: initialises the model by providing a HybridCore object, and the exact solution and boundary conditions of the model.
|
| |
| void | HHO3D::assemble () |
| | A general assemble routine that calculates the statically condensed matrices required by solve.
|
| |
| UVector | HHO3D::solve () |
| | Solves the statically condensed system.
|
| |
| UVector | HHO3D::neumann_solve () |
| | Solves the system when the model is ill posed (not yet running)
|
| |
| double | HHO3D::energy_norm (const UVector) |
| | Returns the energy norm of a given UVector.
|
| |
| void | HHO3D::set_global_operator (const MatrixFType &) |
| | Set the global operator.
|
| |
| void | HHO3D::set_load_vector (const VectorFType &) |
| | Set the load vector.
|
| |
| void | HHO3D::plot (const std::string, const UVector &, const FType< double > &) |
| | Plot the numerical and exact solutions.
|
| |
| VectorFType | HHO3D::standard_load_vector (const CellFType< double > &) |
| | Returns the standard load vector (f, v_T)_T with no Neumann boundary conditions.
|
| |
| VectorFType | HHO3D::standard_load_vector (const CellFType< double > &, const CellFType< VectorRd > &, const BoundaryConditions &) |
| | Returns the standard load vector (f, v_T)_T.
|
| |
| VectorFType | HHO3D::standard_load_vector (const CellFType< double > &, const FType< double > &, const BoundaryConditions &) |
| | Returns the standard load vector (f, v_T)_T.
|
| |
| void | HHO3D::set_dirichlet (const FType< double > &, const size_t) |
| | Set the Dirichlet boundary conditions.
|
| |
| void | HHO3D::set_dirichlet (const size_t) |
| | Set the Dirichlet boundary condition to zero.
|
| |
| Eigen::SparseMatrix< double > | HHO3D::get_SysMat () |
| | Return the (statically condensed) matrix system.
|
| |
| double | HHO3D::get_assembly_time () const |
| | CPU time to assemble the scheme.
|
| |
| double | HHO3D::get_solving_time () const |
| | CPU time to solve the scheme.
|
| |
| double | HHO3D::get_solving_error () const |
| | Residual after solving the scheme.
|
| |
Class providing helper methods and assemble and solve routines for general 3D HHO schemes.
◆ MatrixFType
type for the global operator as a function of a Cell and an ElementQuad
The HHO3D class contains methods applicable to a general 3D HHO model. It contains assemble and solve routines usable by any scheme. The global operator and load vector are passed to the model through setters and are then accessible by assemble and solve. The class also contains helper methods for writing HHO schemes.
◆ VectorFType
type for the load vector as a function of a Cell and an ElementQuad
◆ assemble()
A general assemble routine that calculates the statically condensed matrices required by solve.
◆ energy_norm()
Returns the energy norm of a given UVector.
◆ get_assembly_time()
| double HHO3D::get_assembly_time |
( |
| ) |
const |
|
inline |
CPU time to assemble the scheme.
◆ get_solving_error()
| double HHO3D::get_solving_error |
( |
| ) |
const |
|
inline |
Residual after solving the scheme.
◆ get_solving_time()
| double HHO3D::get_solving_time |
( |
| ) |
const |
|
inline |
CPU time to solve the scheme.
◆ get_SysMat()
| Eigen::SparseMatrix< double > HHO3D::get_SysMat |
( |
| ) |
|
|
inline |
Return the (statically condensed) matrix system.
◆ HHO3D()
Class constructor: initialises the model by providing a HybridCore object, and the exact solution and boundary conditions of the model.
- Parameters
-
| hho | A reference to the HybridCore object containing mesh data |
| L | Cell polynomial degree |
| K | Face polynomial degree |
| use_threads | Optional argument to indicate if threads should be used |
| doeT | Optional argument to set cell quadrature degree. Default is L + K + 1 |
| doeF | Optional argument to set face quadrature degree. Default is 2 * K + 1 |
◆ neumann_solve()
Solves the system when the model is ill posed (not yet running)
◆ plot()
Plot the numerical and exact solutions.
- Parameters
-
| plot_file | Plot file |
| sol | Numerical solution |
| exact_sol | Exact solution |
◆ set_dirichlet() [1/2]
Set the Dirichlet boundary conditions.
- Parameters
-
| f | Function on the Dirichlet BC faces |
| n_dir_faces | Number of Dirichlet faces |
◆ set_dirichlet() [2/2]
Set the Dirichlet boundary condition to zero.
- Parameters
-
| n_dir_faces | Number of Dirichlet faces |
◆ set_global_operator()
◆ set_load_vector()
◆ solve()
Solves the statically condensed system.
◆ standard_load_vector() [1/3]
Returns the standard load vector (f, v_T)_T with no Neumann boundary conditions.
- Parameters
-
◆ standard_load_vector() [2/3]
Returns the standard load vector (f, v_T)_T.
- Parameters
-
| source | Source term |
| f | Gradient to be dotted with the normal vector on each Neumann face |
| BC | A reference to the BC to determine Neumann faces |
◆ standard_load_vector() [3/3]
Returns the standard load vector (f, v_T)_T.
- Parameters
-
| source | Source term |
| f | Function on the Neumann BC faces |
| BC | A reference to the BC to determine Neumann faces |