DUBeat
1.0.1
High-order discontinuous Galerkin methods and applications to cardiac electrophysiology
|
Class used to discretize analytical solutions as linear combinations of Lagrangian or Dubiner basis. More...
#include <DUB_FEM_handler.hpp>
Public Member Functions | |
DUBFEMHandler (const unsigned int degree, const DoFHandlerDG< basis > &dof_hand) | |
Constructor. More... | |
DUBFEMHandler (DUBFEMHandler< basis > &DUBFEMHandler)=default | |
Default copy constructor. More... | |
DUBFEMHandler (const DUBFEMHandler< basis > &DUBFEMHandler)=default | |
Default const copy constructor. More... | |
DUBFEMHandler (DUBFEMHandler< basis > &&DUBFEMHandler)=default | |
Default move constructor. More... | |
lifex::LinAlg::MPI::Vector | dubiner_to_fem (const lifex::LinAlg::MPI::Vector &dub_solution) const |
Conversion of a discretized solution vector from Dubiner coefficients to FEM coefficients. More... | |
lifex::LinAlg::MPI::Vector | dubiner_to_fem (const lifex::LinAlg::MPI::Vector &dub_solution, const std::string &FEM_mesh_path, const std::string &subsection, const MPI_Comm &mpi_comm_, const unsigned int degree_fem=1, const double scaling_factor=1) const |
Same as dubiner_to_fem but allows to choose the grid refinement and polynomial order for the FE evaluations. More... | |
lifex::LinAlg::MPI::Vector | fem_to_dubiner (const lifex::LinAlg::MPI::Vector &fem_solution) const |
Conversion of a discretized solution vector from FEM coefficients to Dubiner coefficients. More... | |
lifex::LinAlg::MPI::Vector | fem_to_dubiner (const lifex::LinAlg::MPI::Vector &fem_solution, const std::string &FEM_mesh_path, const std::string &subsection, const MPI_Comm &mpi_comm_, const unsigned int degree_fem=1, const double scaling_factor=1) const |
As fem_to_dubiner but it works with a FEM different grid refinement. More... | |
lifex::LinAlg::MPI::Vector | analytical_to_dubiner (lifex::LinAlg::MPI::Vector dub_solution, const std::shared_ptr< dealii::Function< lifex::dim >> &u_analytical) const |
Conversion of an analytical solution to a vector of Dubiner coefficients. More... | |
bool | is_in_unit_cell (const dealii::Point< lifex::dim > &unit_q, const double tol=0) const |
To check if a point is inside the reference cell, needed in other methods. More... | |
Public Member Functions inherited from DUBValues< lifex::dim > | |
DUBValues (const unsigned int degree) | |
Constructor. More... | |
DUBValues (DUBValues< dim > &DUBValues)=default | |
Default copy constructor. More... | |
DUBValues (const DUBValues< dim > &DUBValues)=default | |
Default const copy constructor. More... | |
DUBValues (DUBValues< dim > &&DUBValues)=default | |
Default move constructor. More... | |
unsigned int | get_dofs_per_cell () const |
Return the number of degrees of freedom per element. More... | |
unsigned int | get_dofs_per_cell (const unsigned int degree) const |
Same as get_dofs_per_cell() but with variable space degree. More... | |
double | shape_value (const unsigned int function_no, const dealii::Point< dim > &quadrature_point) const |
Evaluation of the Dubiner basis functions. More... | |
double | shape_value (const unsigned int function_no, const dealii::Point< 2 > &quadrature_point) const |
Evaluation of the function basis in two dimensions. More... | |
double | shape_value (const unsigned int function_no, const dealii::Point< 3 > &quadrature_point) const |
Evaluation of the function basis in three dimensions. More... | |
dealii::Tensor< 1, dim > | shape_grad (const unsigned int function_no, const dealii::Point< dim > &quadrature_point) const |
Evaluation of the Dubiner basis functions gradients. More... | |
dealii::Tensor< 1, 2 > | shape_grad (const unsigned int function_no, const dealii::Point< 2 > &quadrature_point) const |
Evaluation of the gradient of the function basis in two dimensions. More... | |
dealii::Tensor< 1, 3 > | shape_grad (const unsigned int function_no, const dealii::Point< 3 > &quadrature_point) const |
Evaluation of the gradient of the function basis in three dimensions. More... | |
virtual | ~DUBValues ()=default |
Default destructor. More... | |
Private Attributes | |
const DoFHandlerDG< basis > & | dof_handler |
Dof handler object of the problem. More... | |
const unsigned int | n_quad_points |
Number of quadrature points in the volume element. More... | |
Additional Inherited Members | |
Public Attributes inherited from DUBValues< lifex::dim > | |
unsigned int | dofs_per_cell |
Number of degrees of freedom. More... | |
Protected Member Functions inherited from DUBValues< lifex::dim > | |
double | eval_jacobi_polynomial (const unsigned int n, const unsigned int alpha, const unsigned int beta, const double eval_point) const |
Evaluation of the \((n,\alpha,\beta)\)-Jacobi polynomial (used to evaluate Dubiner basis). More... | |
std::array< unsigned int, dim > | fun_coeff_conversion (const unsigned int n) const |
Conversion from the FEM basis taxonomy ( \(1^{st}\) basis function, \(2^{nd}\) basis function \(\ldots\)) to the Dubiner basis taxonomy (where every basis function is indexed by a tuple of indexes \((i,j)\) in 2D and three indexes \((i,j,k)\) in 3D). More... | |
std::array< unsigned int, 2 > | fun_coeff_conversion (const unsigned int n) const |
Specialization in two dimensions, i.e. More... | |
std::array< unsigned int, 3 > | fun_coeff_conversion (const unsigned int n) const |
Specialization in three dimensions, i.e. More... | |
Protected Attributes inherited from DUBValues< lifex::dim > | |
const unsigned int | poly_degree |
Polynomial degree used. More... | |
const double | tol |
Default tolerance. More... | |
Class used to discretize analytical solutions as linear combinations of Lagrangian or Dubiner basis.
It also provides conversions for solution vectors with respect to the Lagrangian basis to the Dubiner basis and viceversa. This class is necessary, for instance, to discretize initial analytical solutions for time-dependent problems or to obtain solution vectors for contour plots at the end of the system solving.
Definition at line 60 of file DUB_FEM_handler.hpp.
|
inline |
Constructor.
Definition at line 68 of file DUB_FEM_handler.hpp.
|
default |
Default copy constructor.
|
default |
Default const copy constructor.
|
default |
Default move constructor.
lifex::LinAlg::MPI::Vector DUBFEMHandler< basis >::analytical_to_dubiner | ( | lifex::LinAlg::MPI::Vector | dub_solution, |
const std::shared_ptr< dealii::Function< lifex::dim >> & | u_analytical | ||
) | const |
Conversion of an analytical solution to a vector of Dubiner coefficients.
Definition at line 426 of file DUB_FEM_handler.hpp.
lifex::LinAlg::MPI::Vector DUBFEMHandler< basis >::dubiner_to_fem | ( | const lifex::LinAlg::MPI::Vector & | dub_solution | ) | const |
Conversion of a discretized solution vector from Dubiner coefficients to FEM coefficients.
The output FEM vector belongs to a space of order at most 2 due to the deal.II current availabilities. See dof_handler_DG.hpp description for more information.
Definition at line 135 of file DUB_FEM_handler.hpp.
lifex::LinAlg::MPI::Vector DUBFEMHandler< basis >::dubiner_to_fem | ( | const lifex::LinAlg::MPI::Vector & | dub_solution, |
const std::string & | FEM_mesh_path, | ||
const std::string & | subsection, | ||
const MPI_Comm & | mpi_comm_, | ||
const unsigned int | degree_fem = 1 , |
||
const double | scaling_factor = 1 |
||
) | const |
Same as dubiner_to_fem but allows to choose the grid refinement and polynomial order for the FE evaluations.
This generalization overcomes the issue of a contour visualization that might be less refined than the numerical solution itself.
Definition at line 189 of file DUB_FEM_handler.hpp.
lifex::LinAlg::MPI::Vector DUBFEMHandler< basis >::fem_to_dubiner | ( | const lifex::LinAlg::MPI::Vector & | fem_solution | ) | const |
Conversion of a discretized solution vector from FEM coefficients to Dubiner coefficients.
Definition at line 272 of file DUB_FEM_handler.hpp.
lifex::LinAlg::MPI::Vector DUBFEMHandler< basis >::fem_to_dubiner | ( | const lifex::LinAlg::MPI::Vector & | fem_solution, |
const std::string & | FEM_mesh_path, | ||
const std::string & | subsection, | ||
const MPI_Comm & | mpi_comm_, | ||
const unsigned int | degree_fem = 1 , |
||
const double | scaling_factor = 1 |
||
) | const |
As fem_to_dubiner but it works with a FEM different grid refinement.
Definition at line 322 of file DUB_FEM_handler.hpp.
bool DUBFEMHandler< basis >::is_in_unit_cell | ( | const dealii::Point< lifex::dim > & | unit_q, |
const double | tol = 0 |
||
) | const |
To check if a point is inside the reference cell, needed in other methods.
Definition at line 468 of file DUB_FEM_handler.hpp.
|
private |
Dof handler object of the problem.
Definition at line 64 of file DUB_FEM_handler.hpp.
|
private |
Number of quadrature points in the volume element.
By default: \((degree+2)^{dim}\).
Definition at line 68 of file DUB_FEM_handler.hpp.