Implementation of the HHO scheme with full gradient for variable diffusion equations.
More...
Implementation of the HHO scheme with full gradient for variable diffusion equations.
◆ ForcingTermType
◆ PermeabilityType
◆ SolutionGradientType
◆ SolutionType
◆ SystemMatrixType
◆ assembleLinearSystem()
Assemble the global system
- Parameters
-
| f | Forcing term |
| kappa | Permeability |
| u | Exact solution for boundary condition |
| UDir | Vector filled in by Dirichlets BCs |
◆ computeEnergyNorms()
| 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.
- Parameters
-
| list_dofs | The list of vectors representing the dofs |
◆ FullGradientDiffusion()
Constructor.
- Parameters
-
| 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 |
◆ hhospace()
| const HHOSpace & HArDCore3D::FullGradientDiffusion::hhospace |
( |
| ) |
const |
|
inline |
◆ mesh()
| const Mesh & HArDCore3D::FullGradientDiffusion::mesh |
( |
| ) |
const |
|
inline |
◆ numDirDOFs()
| size_t HArDCore3D::FullGradientDiffusion::numDirDOFs |
( |
| ) |
const |
|
inline |
Returns the number of Dirichlet DOFs.
◆ numSCDOFs()
| size_t HArDCore3D::FullGradientDiffusion::numSCDOFs |
( |
| ) |
const |
|
inline |
Returns the number of statically condensed DOFs (the cell DOFs)
◆ numSkeletalDOFs()
| size_t HArDCore3D::FullGradientDiffusion::numSkeletalDOFs |
( |
| ) |
const |
|
inline |
Returns the number of DOFs after SC but before eliminating Dirichlet DOFs.
◆ scMatrix()
Returns the static condensation recovery operator.
◆ scVector()
| Eigen::VectorXd & HArDCore3D::FullGradientDiffusion::scVector |
( |
| ) |
|
|
inline |
Returns the static condensation rhs.
◆ sizeSystem()
| size_t HArDCore3D::FullGradientDiffusion::sizeSystem |
( |
| ) |
const |
|
inline |
Returns the size of the final system, after application of SC and removal of Dirichlet BCs.
◆ stabilizationParameter() [1/2]
| double & HArDCore3D::FullGradientDiffusion::stabilizationParameter |
( |
| ) |
|
|
inline |
Returns the stabilization parameter.
◆ stabilizationParameter() [2/2]
| const double & HArDCore3D::FullGradientDiffusion::stabilizationParameter |
( |
| ) |
const |
|
inline |
Returns the stabilization parameter (scaling)
◆ systemMatrix() [1/2]
Returns the linear system matrix.
◆ systemMatrix() [2/2]
Returns the linear system matrix.
◆ systemVector() [1/2]
| Eigen::VectorXd & HArDCore3D::FullGradientDiffusion::systemVector |
( |
| ) |
|
|
inline |
Returns the linear system right-hand side vector.
◆ systemVector() [2/2]
| const Eigen::VectorXd & HArDCore3D::FullGradientDiffusion::systemVector |
( |
| ) |
const |
|
inline |
Returns the linear system right-hand side vector.
◆ linear_f
Initial value:= [](const VectorRd & x) -> double {
return 0.;
}
◆ linear_gradu
Initial value:= [](const VectorRd & x) -> VectorRd {
}
static const VectorRd vec_a
Definition hho-fullgradientdiff.hpp:205
◆ linear_kappa
◆ linear_u
Initial value:= [](const VectorRd & x) -> double {
return vec_a.dot(x) + 4.;
}
◆ PI
◆ trigonometric_f
Initial value:= [](const VectorRd & x) -> double {
}
static const double PI
Definition ddr-magnetostatics.hpp:187
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition ddr-magnetostatics.hpp:238
◆ trigonometric_gradu
Initial value:= [](const VectorRd & x) -> VectorRd {
cos(
PI*x(0)) * sin(
PI*x(1)) * sin(
PI*x(2)),
sin(
PI*x(0)) * cos(
PI*x(1)) * sin(
PI*x(2)),
sin(
PI*x(0)) * sin(
PI*x(1)) * cos(
PI*x(2))
);
}
◆ trigonometric_kappa
◆ trigonometric_u
Initial value:= [](const VectorRd & x) -> double {
return sin(
PI*x(0)) * sin(
PI*x(1)) * sin(
PI*x(2));
}
◆ vec_a
| const VectorRd HArDCore3D::vec_a = VectorRd(1.,2.,3.) |
|
static |