|
HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
Implementations of the Locally Enriched Polytopal Non-Conforming method. More...
Classes | |
| class | HArDCore2D::LEPNC_diffusion |
| The LEPNC_diffusion class provides an implementation of the LEPNC method for the stationnary diffusion problem \(-div(K\nabla u)=f\). More... | |
| class | HArDCore2D::LEPNC_StefanPME |
| LEPNC scheme for diffusion equation \(u - \div(K \nabla(\zeta(u))) = f\). More... | |
| class | HArDCore2D::LEPNC_StefanPME_Transient |
| LEPNC scheme for diffusion equation \(\partial_t u - \div(K \nabla(\zeta(u))) = f\). More... | |
| class | HArDCore2D::LEPNCCore |
Typedefs | |
| using | HArDCore2D::LEPNC_diffusion::solution_function_type = std::function< double(double, double)> |
| type for solution | |
| using | HArDCore2D::LEPNC_diffusion::source_function_type = std::function< double(double, double, Cell *)> |
| type for source | |
| using | HArDCore2D::LEPNC_diffusion::grad_function_type = std::function< VectorRd(double, double, Cell *)> |
| type for gradient | |
| using | HArDCore2D::LEPNC_diffusion::tensor_function_type = std::function< Eigen::Matrix2d(double, double, Cell *)> |
| type for diffusion tensor | |
| using | HArDCore2D::LEPNC_StefanPME::solution_function_type = std::function< double(double, double)> |
| type for solution | |
| using | HArDCore2D::LEPNC_StefanPME::source_function_type = std::function< double(double, double, Cell *)> |
| type for source | |
| using | HArDCore2D::LEPNC_StefanPME::grad_function_type = std::function< VectorRd(double, double, Cell *)> |
| type for gradient | |
| using | HArDCore2D::LEPNC_StefanPME::tensor_function_type = std::function< Eigen::Matrix2d(double, double, Cell *)> |
| type for diffusion tensor | |
| using | HArDCore2D::LEPNC_StefanPME_Transient::solution_function_type = std::function< double(const double &, const VectorRd &)> |
| type for solution | |
| using | HArDCore2D::LEPNC_StefanPME_Transient::source_function_type = std::function< double(const double &, const VectorRd &, const Cell *)> |
| type for source | |
| using | HArDCore2D::LEPNC_StefanPME_Transient::grad_function_type = std::function< VectorRd(const double &, const VectorRd &, const Cell *)> |
| type for gradient | |
| using | HArDCore2D::LEPNC_StefanPME_Transient::tensor_function_type = std::function< Eigen::Matrix2d(const double, const double, const Cell *)> |
| type for diffusion tensor (does not depend on time) | |
| using | HArDCore2D::LEPNCCore::cell_basis_type = std::function< double(double, double)> |
| type for cell basis | |
| using | HArDCore2D::LEPNCCore::cell_gradient_type = std::function< VectorRd(double, double)> |
| type for gradients of cell basis | |
Functions | |
| HArDCore2D::LEPNC_diffusion::LEPNC_diffusion (LEPNCCore &nc, tensor_function_type kappa, size_t deg_kappa, source_function_type source, BoundaryConditions BC, solution_function_type exact_solution, grad_function_type grad_exact_solution, std::string solver_type, std::ostream &output=std::cout) | |
| Constructor of the class. | |
| Eigen::VectorXd | HArDCore2D::LEPNC_diffusion::solve () |
| Assemble and solve the scheme. | |
| double | HArDCore2D::LEPNC_diffusion::EnergyNorm (const Eigen::VectorXd Xh) const |
| Discrete energy norm (associated to the diffusion operator) | |
| double | HArDCore2D::LEPNC_diffusion::get_assembly_time () const |
| cpu time to assemble the scheme | |
| double | HArDCore2D::LEPNC_diffusion::get_solving_time () const |
| cpu time to solve the scheme | |
| double | HArDCore2D::LEPNC_diffusion::get_solving_error () const |
| residual after solving the scheme | |
| double | HArDCore2D::LEPNC_diffusion::get_itime (size_t idx) const |
| various intermediate assembly times | |
| HArDCore2D::LEPNC_StefanPME::LEPNC_StefanPME (LEPNCCore &nc, tensor_function_type kappa, size_t deg_kappa, source_function_type source, BoundaryConditions BC, solution_function_type exact_solution, grad_function_type grad_exact_solution, TestCaseNonLinearity::nonlinearity_function_type zeta, double weight, std::string solver_type, std::ostream &output=std::cout) | |
| Constructor of the class. | |
| Eigen::VectorXd | HArDCore2D::LEPNC_StefanPME::solve () |
| Assemble and solve the scheme. | |
| Eigen::VectorXd | HArDCore2D::LEPNC_StefanPME::apply_nonlinearity (const Eigen::VectorXd &Y, const std::string type) const |
| Compute non-linearity on vector (depends if weight=0, weight=1 or weight\in (0,1) ) | |
| double | HArDCore2D::LEPNC_StefanPME::L2_MassLumped (const Eigen::VectorXd Xh) const |
| Mass-lumped L2 norm of a function given by a vector. | |
| double | HArDCore2D::LEPNC_StefanPME::EnergyNorm (const Eigen::VectorXd Xh) const |
| Discrete energy norm (associated to the diffusion operator) | |
| double | HArDCore2D::LEPNC_StefanPME::get_assembly_time () const |
| cpu time to assemble the scheme | |
| double | HArDCore2D::LEPNC_StefanPME::get_solving_time () const |
| cpu time to solve the scheme | |
| double | HArDCore2D::LEPNC_StefanPME::get_solving_error () const |
| residual after solving the scheme | |
| double | HArDCore2D::LEPNC_StefanPME::get_itime (size_t idx) const |
| various intermediate assembly times | |
| HArDCore2D::LEPNC_StefanPME_Transient::LEPNC_StefanPME_Transient (LEPNCCore &nc, tensor_function_type kappa, size_t deg_kappa, source_function_type source, BoundaryConditions BC, solution_function_type exact_solution, grad_function_type grad_exact_solution, TestCaseNonLinearity::nonlinearity_function_type zeta, double weight, std::string solver_type, std::ostream &output=std::cout) | |
| Constructor of the class. | |
| Eigen::VectorXd | HArDCore2D::LEPNC_StefanPME_Transient::iterate (const double tps, const double dt, const Eigen::VectorXd Xn) |
| Execute one time iteration. | |
| Eigen::VectorXd | HArDCore2D::LEPNC_StefanPME_Transient::apply_nonlinearity (const Eigen::VectorXd &Y, const std::string type) const |
| Compute non-linearity on vector (depends if weight=0, weight=1 or weight\in (0,1) ) | |
| double | HArDCore2D::LEPNC_StefanPME_Transient::L2_MassLumped (const Eigen::VectorXd Xh) const |
| Mass-lumped L2 norm of a function given by a vector. | |
| double | HArDCore2D::LEPNC_StefanPME_Transient::Lp_MassLumped (const Eigen::VectorXd Xh, double p) const |
| Mass-lumped Lp norm of a function given by a vector. | |
| double | HArDCore2D::LEPNC_StefanPME_Transient::EnergyNorm (const Eigen::VectorXd Xh) const |
| Discrete energy norm (associated to the diffusion operator) | |
| double | HArDCore2D::LEPNC_StefanPME_Transient::get_assembly_time () const |
| cpu time to assemble the scheme | |
| double | HArDCore2D::LEPNC_StefanPME_Transient::get_solving_time () const |
| cpu time to solve the scheme | |
| double | HArDCore2D::LEPNC_StefanPME_Transient::get_itime (size_t idx) const |
| various intermediate assembly times | |
| double | HArDCore2D::LEPNC_StefanPME_Transient::get_solving_error () const |
| residual after solving the scheme | |
| size_t | HArDCore2D::LEPNC_StefanPME_Transient::get_nb_newton () const |
| number of Newton iterations | |
| Eigen::MatrixXd | HArDCore2D::LEPNC_StefanPME_Transient::get_MassT (size_t iT) const |
| Mass matrix in cell iT. | |
| HArDCore2D::LEPNCCore::LEPNCCore (const Mesh *mesh_ptr, const size_t K, const size_t L, std::ostream &output=std::cout) | |
| Class constructor: initialise the LEPNCCore class, and creates non-conforming basis functions (and gradients) | |
| const cell_basis_type & | HArDCore2D::LEPNCCore::nc_basis (size_t iT, size_t i) const |
| Return a reference to the i'th non-conforming basis function of the cell iT. | |
| const cell_gradient_type & | HArDCore2D::LEPNCCore::nc_gradient (size_t iT, size_t i) const |
| Return a reference to the gradient of the i'th non-conforming basis function of the cell iT. | |
| const VectorRd | HArDCore2D::LEPNCCore::ml_node (size_t iT, size_t i) const |
| Return the i'th nodal point in cell iT (for mass lumping) | |
| const std::vector< Eigen::ArrayXd > | HArDCore2D::LEPNCCore::nc_basis_quad (const size_t iT, const QuadratureRule quad) const |
| Computes non-conforming basis functions at the given quadrature nodes. | |
| const std::vector< Eigen::ArrayXXd > | HArDCore2D::LEPNCCore::grad_nc_basis_quad (const size_t iT, const QuadratureRule quad) const |
| Compute \((\nabla \phi_i^{nc})_{i\in I}\) at the quadrature nodes, where \((\phi_i^{nc})_{i\in I}\) are the cell basis functions. | |
| Eigen::MatrixXd | HArDCore2D::LEPNCCore::gram_matrix (const std::vector< Eigen::ArrayXd > &f_quad, const std::vector< Eigen::ArrayXd > &g_quad, const size_t &nrows, const size_t &ncols, const QuadratureRule &quad, const bool &sym, std::vector< double > L2weight={}) const |
| Eigen::MatrixXd | HArDCore2D::LEPNCCore::gram_matrix (const std::vector< Eigen::ArrayXXd > &F_quad, const std::vector< Eigen::ArrayXXd > &G_quad, const size_t &nrows, const size_t &ncols, const QuadratureRule &quad, const bool &sym, std::vector< Eigen::Matrix2d > L2Weight={}) const |
| Overloaded version of the previous one for vector-valued functions: the functions (F_i) and (G_j) are vector-valued functions. | |
| template<typename Function > | |
| Eigen::VectorXd | HArDCore2D::LEPNCCore::nc_interpolate_moments (const Function &f, size_t doe) const |
| Interpolates a continuous function on the degrees of freedom, using the moments on the basis functions. The first ones are the cell DOFs (DimPoly<Cell>(1) for each cell), the last ones are the edge DOFs (one for each edge) | |
| template<typename Function > | |
| Eigen::VectorXd | HArDCore2D::LEPNCCore::nc_interpolate_ml (const Function &f, size_t doe) const |
| Interpolates a continuous function on the degrees of freedom, using the moments on the basis functions associated to the edges and the nodal values (corresponding to mass-lumping) on the cell basis functions. The first ones are the cell DOFs (DimPoly<Cell>(1) for each cell), the last ones are the edge DOFs (one for each edge) | |
| Eigen::VectorXd | HArDCore2D::LEPNCCore::nc_restr (const Eigen::VectorXd &Xh, size_t iT) const |
| Extract from a global vector Xh of unknowns the non-conforming unknowns corresponding to cell iT. | |
| double | HArDCore2D::LEPNCCore::nc_L2norm (const Eigen::VectorXd &Xh) const |
| Compute L2 norm of a discrete function (given by coefficients on the basis functions) | |
| double | HArDCore2D::LEPNCCore::nc_L2norm_ml (const Eigen::VectorXd &Xh) const |
| Compute L2 norm of the mass-lumped discrete function (given by coefficients on the basis functions) | |
| double | HArDCore2D::LEPNCCore::nc_H1norm (const Eigen::VectorXd &Xh) const |
| Compute broken H1 norm of a discrete function (given by coefficients on the basis functions) | |
| double | HArDCore2D::LEPNCCore::nc_evaluate_in_cell (const Eigen::VectorXd XTF, size_t iT, double x, double y) const |
| Evaluates a non-conforming discrete function in the cell iT at point (x,y) | |
| Eigen::VectorXd | HArDCore2D::LEPNCCore::nc_VertexValues (const Eigen::VectorXd Xh, const double weight=0) |
| From a non-conforming discrete function, computes a vector of values at the vertices of the mesh. | |
Implementations of the Locally Enriched Polytopal Non-Conforming method.
Locally Enriched Polytopal Non-Conforming method.
| using HArDCore2D::LEPNCCore::cell_basis_type = std::function<double(double, double)> |
type for cell basis
| using HArDCore2D::LEPNCCore::cell_gradient_type = std::function<VectorRd(double, double)> |
type for gradients of cell basis
| using HArDCore2D::LEPNC_diffusion::grad_function_type = std::function<VectorRd(double,double,Cell*)> |
type for gradient
| using HArDCore2D::LEPNC_StefanPME::grad_function_type = std::function<VectorRd(double,double,Cell*)> |
type for gradient
| using HArDCore2D::LEPNC_StefanPME_Transient::grad_function_type = std::function<VectorRd(const double&, const VectorRd&, const Cell*)> |
type for gradient
| using HArDCore2D::LEPNC_diffusion::solution_function_type = std::function<double(double,double)> |
type for solution
| using HArDCore2D::LEPNC_StefanPME::solution_function_type = std::function<double(double,double)> |
type for solution
| using HArDCore2D::LEPNC_StefanPME_Transient::solution_function_type = std::function<double(const double&, const VectorRd&)> |
type for solution
| using HArDCore2D::LEPNC_diffusion::source_function_type = std::function<double(double,double,Cell*)> |
type for source
| using HArDCore2D::LEPNC_StefanPME::source_function_type = std::function<double(double,double,Cell*)> |
type for source
| using HArDCore2D::LEPNC_StefanPME_Transient::source_function_type = std::function<double(const double&, const VectorRd&, const Cell*)> |
type for source
| using HArDCore2D::LEPNC_diffusion::tensor_function_type = std::function<Eigen::Matrix2d(double,double,Cell*)> |
type for diffusion tensor
| using HArDCore2D::LEPNC_StefanPME::tensor_function_type = std::function<Eigen::Matrix2d(double,double,Cell*)> |
type for diffusion tensor
| using HArDCore2D::LEPNC_StefanPME_Transient::tensor_function_type = std::function<Eigen::Matrix2d(const double, const double, const Cell*)> |
type for diffusion tensor (does not depend on time)
| Eigen::VectorXd HArDCore2D::LEPNC_StefanPME::apply_nonlinearity | ( | const Eigen::VectorXd & | Y, |
| const std::string | type | ||
| ) | const |
Compute non-linearity on vector (depends if weight=0, weight=1 or weight\in (0,1) )
| Eigen::VectorXd HArDCore2D::LEPNC_StefanPME_Transient::apply_nonlinearity | ( | const Eigen::VectorXd & | Y, |
| const std::string | type | ||
| ) | const |
Compute non-linearity on vector (depends if weight=0, weight=1 or weight\in (0,1) )
| double HArDCore2D::LEPNC_diffusion::EnergyNorm | ( | const Eigen::VectorXd | Xh | ) | const |
Discrete energy norm (associated to the diffusion operator)
| double HArDCore2D::LEPNC_StefanPME::EnergyNorm | ( | const Eigen::VectorXd | Xh | ) | const |
Discrete energy norm (associated to the diffusion operator)
| double HArDCore2D::LEPNC_StefanPME_Transient::EnergyNorm | ( | const Eigen::VectorXd | Xh | ) | const |
Discrete energy norm (associated to the diffusion operator)
| double HArDCore2D::LEPNC_diffusion::get_assembly_time | ( | ) | const |
cpu time to assemble the scheme
| double HArDCore2D::LEPNC_StefanPME::get_assembly_time | ( | ) | const |
cpu time to assemble the scheme
|
inline |
cpu time to assemble the scheme
| double HArDCore2D::LEPNC_diffusion::get_itime | ( | size_t | idx | ) | const |
various intermediate assembly times
| double HArDCore2D::LEPNC_StefanPME::get_itime | ( | size_t | idx | ) | const |
various intermediate assembly times
|
inline |
various intermediate assembly times
|
inline |
Mass matrix in cell iT.
|
inline |
number of Newton iterations
| double HArDCore2D::LEPNC_diffusion::get_solving_error | ( | ) | const |
residual after solving the scheme
| double HArDCore2D::LEPNC_StefanPME::get_solving_error | ( | ) | const |
residual after solving the scheme
|
inline |
residual after solving the scheme
| double HArDCore2D::LEPNC_diffusion::get_solving_time | ( | ) | const |
cpu time to solve the scheme
| double HArDCore2D::LEPNC_StefanPME::get_solving_time | ( | ) | const |
cpu time to solve the scheme
|
inline |
cpu time to solve the scheme
| const std::vector< Eigen::ArrayXXd > HArDCore2D::LEPNCCore::grad_nc_basis_quad | ( | const size_t | iT, |
| const QuadratureRule | quad | ||
| ) | const |
Compute \((\nabla \phi_i^{nc})_{i\in I}\) at the quadrature nodes, where \((\phi_i^{nc})_{i\in I}\) are the cell basis functions.
| iT | global index of the cell |
| quad | quadrature rules in the cell |
| Eigen::MatrixXd HArDCore2D::LEPNCCore::gram_matrix | ( | const std::vector< Eigen::ArrayXd > & | f_quad, |
| const std::vector< Eigen::ArrayXd > & | g_quad, | ||
| const size_t & | nrows, | ||
| const size_t & | ncols, | ||
| const QuadratureRule & | quad, | ||
| const bool & | sym, | ||
| std::vector< double > | L2weight = {} |
||
| ) | const |
Create the matrix of L2 products of two families (f_i) and (g_j) of functions (this is not really a Gram matrix, unless the two families are the same)
Create the matrix of L2 products of two families (f_i) and (g_j) of functions (this is not really a Gram matrix, unless the two families are the same)
| f_quad | Values of functions (f1,f2,...) at the quadrature nodes |
| g_quad | Values of functions (g1,g2,...) at the quadrature nodes |
| nrows | Number of rows of the matrix - typically number of functions f_i (but could be less) |
| ncols | Number of columns of the matrix - typically number of functions g_j (but could be less) |
| quad | Quadrature nodes for integration |
| sym | True if the matrix is pseudo-symmetric (that is, #f<=#g and f_i=g_i if i<=#f) |
| L2weight | Optional weight for the L2 product. If provided, should be a std::vector<double> of the weight at the quadrature nodes |
| Eigen::MatrixXd HArDCore2D::LEPNCCore::gram_matrix | ( | const std::vector< Eigen::ArrayXXd > & | F_quad, |
| const std::vector< Eigen::ArrayXXd > & | G_quad, | ||
| const size_t & | nrows, | ||
| const size_t & | ncols, | ||
| const QuadratureRule & | quad, | ||
| const bool & | sym, | ||
| std::vector< Eigen::Matrix2d > | L2Weight = {} |
||
| ) | const |
Overloaded version of the previous one for vector-valued functions: the functions (F_i) and (G_j) are vector-valued functions.
| F_quad | Values of functions (F1,F2,...) at the quadrature nodes |
| G_quad | Values of functions (G1,G2,...) at the quadrature nodes |
| nrows | Number of rows of the matrix - typically number of functions F_i (but could be less) |
| ncols | Number of rows of the matrix - typically number of functions G_j (but could be less) |
| quad | Quadrature nodes for integration |
| sym | True if the matrix is pseudo-symmetric (that is, #F<=#G and F_i=G_i if i<=#F) |
| L2Weight | Optional weight for the L2 product. If provided, should be a std::vector<Eigen::Matrix2d> of the weight at the quadrature nodes |
| Eigen::VectorXd HArDCore2D::LEPNC_StefanPME_Transient::iterate | ( | const double | tps, |
| const double | dt, | ||
| const Eigen::VectorXd | Xn | ||
| ) |
Execute one time iteration.
| tps | current time (for computing source) |
| dt | time step |
| Xn | Unkown at previous time step |
| double HArDCore2D::LEPNC_StefanPME::L2_MassLumped | ( | const Eigen::VectorXd | Xh | ) | const |
Mass-lumped L2 norm of a function given by a vector.
| double HArDCore2D::LEPNC_StefanPME_Transient::L2_MassLumped | ( | const Eigen::VectorXd | Xh | ) | const |
Mass-lumped L2 norm of a function given by a vector.
| HArDCore2D::LEPNC_diffusion::LEPNC_diffusion | ( | LEPNCCore & | nc, |
| tensor_function_type | kappa, | ||
| size_t | deg_kappa, | ||
| source_function_type | source, | ||
| BoundaryConditions | BC, | ||
| solution_function_type | exact_solution, | ||
| grad_function_type | grad_exact_solution, | ||
| std::string | solver_type, | ||
| std::ostream & | output = std::cout |
||
| ) |
Constructor of the class.
| nc | instance of the core nonconforming class (with basis functions, etc.) |
| kappa | diffusion tensor |
| deg_kappa | polynomial degree of the diffusion tensor |
| source | source term |
| BC | type of boundary conditions |
| exact_solution | exact solution |
| grad_exact_solution | gradient of the exact solution |
| solver_type | type of solver to use for the global system (bicgstab at the moment) |
| output | optional argument for output of messages |
| HArDCore2D::LEPNC_StefanPME::LEPNC_StefanPME | ( | LEPNCCore & | nc, |
| tensor_function_type | kappa, | ||
| size_t | deg_kappa, | ||
| source_function_type | source, | ||
| BoundaryConditions | BC, | ||
| solution_function_type | exact_solution, | ||
| grad_function_type | grad_exact_solution, | ||
| TestCaseNonLinearity::nonlinearity_function_type | zeta, | ||
| double | weight, | ||
| std::string | solver_type, | ||
| std::ostream & | output = std::cout |
||
| ) |
Constructor of the class.
| nc | instance of the core nonconforming class (with basis functions, etc.) |
| kappa | diffusion tensor |
| deg_kappa | polynomial degree of the diffusion tensor |
| source | source term |
| BC | type of boundary conditions |
| exact_solution | exact solution |
| grad_exact_solution | gradient of the exact solution |
| zeta | function describing the nonlinearity |
| weight | proportion of weight (between 0 and 1) of mass-lumping on the edges |
| solver_type | type of solver to use for the global system (bicgstab at the moment) |
| output | optional argument for output of messages |
| HArDCore2D::LEPNC_StefanPME_Transient::LEPNC_StefanPME_Transient | ( | LEPNCCore & | nc, |
| tensor_function_type | kappa, | ||
| size_t | deg_kappa, | ||
| source_function_type | source, | ||
| BoundaryConditions | BC, | ||
| solution_function_type | exact_solution, | ||
| grad_function_type | grad_exact_solution, | ||
| TestCaseNonLinearity::nonlinearity_function_type | zeta, | ||
| double | weight, | ||
| std::string | solver_type, | ||
| std::ostream & | output = std::cout |
||
| ) |
Constructor of the class.
| nc | instance of the core nonconforming class (with basis functions, etc.) |
| kappa | diffusion tensor |
| deg_kappa | polynomial degree of the diffusion tensor |
| source | source term |
| BC | type of boundary conditions |
| exact_solution | exact solution |
| grad_exact_solution | gradient of the exact solution |
| zeta | function describing the nonlinearity |
| weight | proportion of weight (between 0 and 1) of mass-lumping on the edges |
| solver_type | type of solver to use for the global system (bicgstab at the moment) |
| output | optional argument for output of messages |
| HArDCore2D::LEPNCCore::LEPNCCore | ( | const Mesh * | mesh_ptr, |
| const size_t | K, | ||
| const size_t | L, | ||
| std::ostream & | output = std::cout |
||
| ) |
Class constructor: initialise the LEPNCCore class, and creates non-conforming basis functions (and gradients)
| mesh_ptr | A pointer to the loaded mesh |
| K | The degree of the edge polynomials |
| L | The degree of the cell polynomials |
| output | Optional argument for specifying outputs of messages. |
| double HArDCore2D::LEPNC_StefanPME_Transient::Lp_MassLumped | ( | const Eigen::VectorXd | Xh, |
| double | p | ||
| ) | const |
Mass-lumped Lp norm of a function given by a vector.
|
inline |
Return the i'th nodal point in cell iT (for mass lumping)
| iT | The global region number of the cell |
| i | The index of the desired nodal point |
|
inline |
Return a reference to the i'th non-conforming basis function of the cell iT.
For i=DimPoly<Cell>(1) to nedge-1, nc_basis(iT,i) is the basis function associated with edge i (that is, product of all distances to all other edges, scaled to have an integral over edge i equal to 1). For i=0, ..., DimPoly<Cell>(1)-1 the basis function corresponds to 1, x, y, to which we subtract a linear combination of the basis functions corresponding to the edge in order to obtain basis functions with a zero integral on each edge
| iT | The global region number of the cell |
| i | The index of the desired basis function |
| const std::vector< Eigen::ArrayXd > HArDCore2D::LEPNCCore::nc_basis_quad | ( | const size_t | iT, |
| const QuadratureRule | quad | ||
| ) | const |
Computes non-conforming basis functions at the given quadrature nodes.
| iT | global index of the cell/edge |
| quad | quadrature nodes and weights on the cell/edge |
| double HArDCore2D::LEPNCCore::nc_evaluate_in_cell | ( | const Eigen::VectorXd | XTF, |
| size_t | iT, | ||
| double | x, | ||
| double | y | ||
| ) | const |
Evaluates a non-conforming discrete function in the cell iT at point (x,y)
|
inline |
Return a reference to the gradient of the i'th non-conforming basis function of the cell iT.
The gradient functions are indexed the same as the basis functions.
| iT | The global region number of the cell |
| i | The index of the desired basis function |
| double HArDCore2D::LEPNCCore::nc_H1norm | ( | const Eigen::VectorXd & | Xh | ) | const |
Compute broken H1 norm of a discrete function (given by coefficients on the basis functions)
| Eigen::VectorXd HArDCore2D::LEPNCCore::nc_interpolate_ml | ( | const Function & | f, |
| size_t | doe | ||
| ) | const |
Interpolates a continuous function on the degrees of freedom, using the moments on the basis functions associated to the edges and the nodal values (corresponding to mass-lumping) on the cell basis functions. The first ones are the cell DOFs (DimPoly<Cell>(1) for each cell), the last ones are the edge DOFs (one for each edge)
| Eigen::VectorXd HArDCore2D::LEPNCCore::nc_interpolate_moments | ( | const Function & | f, |
| size_t | doe | ||
| ) | const |
Interpolates a continuous function on the degrees of freedom, using the moments on the basis functions. The first ones are the cell DOFs (DimPoly<Cell>(1) for each cell), the last ones are the edge DOFs (one for each edge)
| double HArDCore2D::LEPNCCore::nc_L2norm | ( | const Eigen::VectorXd & | Xh | ) | const |
Compute L2 norm of a discrete function (given by coefficients on the basis functions)
| double HArDCore2D::LEPNCCore::nc_L2norm_ml | ( | const Eigen::VectorXd & | Xh | ) | const |
Compute L2 norm of the mass-lumped discrete function (given by coefficients on the basis functions)
| Eigen::VectorXd HArDCore2D::LEPNCCore::nc_restr | ( | const Eigen::VectorXd & | Xh, |
| size_t | iT | ||
| ) | const |
Extract from a global vector Xh of unknowns the non-conforming unknowns corresponding to cell iT.
| Eigen::VectorXd HArDCore2D::LEPNCCore::nc_VertexValues | ( | const Eigen::VectorXd | Xh, |
| const double | weight = 0 |
||
| ) |
From a non-conforming discrete function, computes a vector of values at the vertices of the mesh.
| Xh | non-conforming function (coefficients on basis) |
| weight | weight put on the edge values when evaluating the functions at the vertices |
| Eigen::VectorXd HArDCore2D::LEPNC_diffusion::solve | ( | ) |
Assemble and solve the scheme.
| Eigen::VectorXd HArDCore2D::LEPNC_StefanPME::solve | ( | ) |
Assemble and solve the scheme.