v0.15.0
Loading...
Searching...
No Matches
LobattoPolynomial.hpp
Go to the documentation of this file.
1/** \file LobattoPolynomial.hpp
2\brief Implementation of Lobatto polynomials
3
4This file provides functions and classes for calculating Lobatto polynomial
5basis functions and their derivatives for use in finite element approximations.
6Lobatto polynomials are particularly useful for spectral and high-order methods.
7
8*/
9
10
11
12#ifndef __LOBATTOPOLYNOMIALS_HPP__
13#define __LOBATTOPOLYNOMIALS_HPP__
14
15namespace MoFEM {
16
17
18/**
19 \brief Calculate Lobatto basis functions \cite FUENTES2015353.
20
21 \ingroup mofem_base_functions
22
23 Lobatto polynomials are hierarchical basis functions that include vertex modes
24 and edge/internal modes. They are constructed from Legendre polynomials and
25 provide excellent conditioning properties for high-order finite element methods.
26
27 The first function corresponds to order 2 and the sequence continues up to order p.
28 These polynomials are particularly useful for spectral element methods due to
29 their orthogonality properties on the reference interval [-1,1].
30
31 For more information on Lobatto polynomials, see:
32 https://en.wikipedia.org/wiki/Gauss%E2%80%93Lobatto_quadrature
33
34 \param p approximation order (maximum polynomial degree)
35 \param s coordinate mapping from edge to \f$[-1, 1]\f$, i.e., \f$s(\xi_1,\cdot,\xi_{dim})\in[-1,1]\f$
36 \param diff_s [in] Jacobian of the coordinate transformation, i.e. \f$\frac{\partial s}{\partial \xi_i}\f$
37 \param L [out] values of basis functions evaluated at coordinate s
38 \param diffL [out] derivatives of basis functions at s, i.e. \f$\frac{\partial L}{\partial \xi_i}\f$
39 \param dim spatial dimension
40 \return PetscErrorCode indicating success or failure
41
42 \note Output arrays L and diffL must be pre-allocated with sufficient size
43*/
44PetscErrorCode Lobatto_polynomials(int p, double s, double *diff_s, double *L,
45 double *diffL, const int dim);
46
47/**
48 * \brief Context class for Lobatto basis function arguments
49 * \ingroup mofem_base_functions
50 *
51 * This class extends LegendrePolynomialCtx to provide context for evaluating
52 * Lobatto polynomial basis functions. It inherits all functionality from the
53 * Legendre context and specializes the polynomial evaluation function.
54 */
56
57 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
58 UnknownInterface **iface) const;
59
60 /**
61 * \brief Constructor
62 * \param p approximation order
63 * \param diff_s derivatives of coordinate parameter
64 * \param dim spatial dimension
65 * \param base_fun_ptr pointer to base function matrix
66 * \param base_diff_fun_ptr pointer to base function derivatives matrix
67 */
68 LobattoPolynomialCtx(int p, double *diff_s, int dim,
69 boost::shared_ptr<MatrixDouble> base_fun_ptr,
70 boost::shared_ptr<MatrixDouble> base_diff_fun_ptr)
71 : LegendrePolynomialCtx(p, diff_s, dim, base_fun_ptr, base_diff_fun_ptr) {
73 }
75};
76
77/**
78 * \brief Class for calculating Lobatto basis functions
79 * \ingroup mofem_base_functions
80 *
81 * This class extends LegendrePolynomial to provide Lobatto polynomial
82 * evaluation capabilities. Lobatto polynomials are hierarchical basis
83 * functions that include both vertex and higher-order modes, making them
84 * ideal for spectral and high-order finite element methods.
85 */
87
88 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
89 UnknownInterface **iface) const;
90
93
94 /**
95 * \brief Evaluate Lobatto basis functions at given points
96 * \param pts matrix of evaluation points
97 * \param ctx_ptr context containing evaluation parameters
98 * \return error code
99 */
101 boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
102};
103
104/**
105 * \brief Context class for Kernel Lobatto basis function arguments
106 * \ingroup mofem_base_functions
107 *
108 * This class provides context for evaluating Kernel Lobatto polynomial basis
109 * functions. Kernel Lobatto polynomials are a variant that may have different
110 * orthogonality or boundary conditions compared to standard Lobatto polynomials.
111 */
113
114 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
115 UnknownInterface **iface) const;
116
117 /**
118 * \brief Constructor
119 * \param p approximation order
120 * \param diff_s derivatives of coordinate parameter
121 * \param dim spatial dimension
122 * \param base_fun_ptr pointer to base function matrix
123 * \param base_diff_fun_ptr pointer to base function derivatives matrix
124 */
125 KernelLobattoPolynomialCtx(int p, double *diff_s, int dim,
126 boost::shared_ptr<MatrixDouble> base_fun_ptr,
127 boost::shared_ptr<MatrixDouble> base_diff_fun_ptr)
128 : LegendrePolynomialCtx(p, diff_s, dim, base_fun_ptr, base_diff_fun_ptr) {
130 }
132};
133
134/**
135 * \brief Class for calculating Kernel Lobatto basis functions
136 * \ingroup mofem_base_functions
137 *
138 * This class extends LegendrePolynomial to provide Kernel Lobatto polynomial
139 * evaluation capabilities. Kernel Lobatto polynomials are a specialized variant
140 * that may be used for specific boundary conditions or integration schemes.
141 */
143
144 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
145 UnknownInterface **iface) const;
146
149
150 /**
151 * \brief Evaluate Kernel Lobatto basis functions at given points
152 * \param pts matrix of evaluation points
153 * \param ctx_ptr context containing evaluation parameters
154 * \return error code
155 */
157 boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
158};
159
160} // namespace MoFEM
161
162#endif
PetscErrorCode Lobatto_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Lobatto basis functions .
PetscErrorCode LobattoKernel_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Kernel Lobatto base functions.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
Context class for Kernel Lobatto basis function arguments.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
KernelLobattoPolynomialCtx(int p, double *diff_s, int dim, boost::shared_ptr< MatrixDouble > base_fun_ptr, boost::shared_ptr< MatrixDouble > base_diff_fun_ptr)
Constructor.
Class for calculating Kernel Lobatto basis functions.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Evaluate Kernel Lobatto basis functions at given points.
Context class for Legendre basis function arguments.
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Function pointer for polynomial evaluation.
Class for calculating Legendre basis functions.
Context class for Lobatto basis function arguments.
LobattoPolynomialCtx(int p, double *diff_s, int dim, boost::shared_ptr< MatrixDouble > base_fun_ptr, boost::shared_ptr< MatrixDouble > base_diff_fun_ptr)
Constructor.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Class for calculating Lobatto basis functions.
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Evaluate Lobatto basis functions at given points.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
base class for all interface classes