v0.15.0
Loading...
Searching...
No Matches
LegendrePolynomial.hpp
Go to the documentation of this file.
1/** \file LegendrePolynomial.hpp
2\brief Implementation of Legendre polynomials
3
4This file provides functions and classes for calculating Legendre polynomial
5basis functions and their derivatives for use in finite element approximations.
6
7*/
8
9
10
11#ifndef __LEGENDREPOLYNOMIALS_HPP__
12#define __LEGENDREPOLYNOMIALS_HPP__
13
14namespace MoFEM {
15
16/**
17\brief Calculate Legendre approximation basis
18
19\ingroup mofem_base_functions
20
21This function computes Legendre polynomials and their derivatives for use in
22finite element approximations. Legendre polynomials form a complete orthogonal
23system on the interval [-1,1] and are widely used in high-order finite element
24methods due to their excellent numerical properties.
25
26The Legendre polynomials are defined by the recurrence relation:
27\f[
28P_0(s)=1;\quad P_1(s) = s
29\f]
30and the following terms are generated inductively using Bonnet's recursion formula:
31\f[
32P_{n+1}(s)=\frac{2n+1}{n+1}sP_n(s)-\frac{n}{n+1}P_{n-1}(s)
33\f]
34
35The parameter \f$s\f$ is the normalized coordinate:
36\f[
37s\in[-1,1] \quad \textrm{and}\; s=s(\xi_0,\xi_1,\xi_2)
38\f]
39where \f$\xi_i\f$ are barycentric coordinates of the element.
40
41For more information on Legendre polynomials, see:
42https://en.wikipedia.org/wiki/Legendre_polynomials
43
44\param p approximation order (degree of highest polynomial)
45\param s position parameter \f$s\in[-1,1]\f$ at which to evaluate polynomials
46\param diff_s derivatives of coordinate transformation, i.e. \f$\frac{\partial s}{\partial \xi_i}\f$
47\param L [out] array of Legendre polynomial values \f$P_0(s), P_1(s), \ldots, P_p(s)\f$
48\param diffL [out] array of polynomial derivatives, i.e. \f$\frac{\partial P_k}{\partial \xi_i}\f$
49\param dim spatial dimension for derivative calculations
50\return PetscErrorCode indicating success or failure
51
52\note The output arrays L and diffL must be pre-allocated with sufficient size
53\note For complex-valued evaluations, use the overloaded version with complex arguments
54
55*/
56PetscErrorCode Legendre_polynomials(int p, double s,
57 double *diff_s, double *L,
58 double *diffL, int dim);
59
60
61/**
62 * \brief Complex-valued Legendre polynomial evaluation
63 *
64 * This overloaded version computes Legendre polynomials for complex arguments,
65 * which is useful for analytical continuation and complex finite element methods.
66 * The recurrence relations and mathematical properties are identical to the
67 * real-valued version.
68 *
69 * \param p approximation order (degree of highest polynomial)
70 * \param s complex position parameter at which to evaluate polynomials
71 * \param diff_s derivatives of coordinate transformation (real-valued)
72 * \param L [out] array of complex Legendre polynomial values
73 * \param diffL [out] array of complex polynomial derivatives
74 * \param dim spatial dimension for derivative calculations
75 * \return PetscErrorCode indicating success or failure
76 *
77 * \see Legendre_polynomials(int, double, double*, double*, double*, int) for real version
78 */
79PetscErrorCode Legendre_polynomials(
80 int p, std::complex<double> s, std::complex<double>* diff_s,
81 std::complex<double>* L, std::complex<double>* diffL, int dim);
82
83/**
84 * \brief Context class for Legendre basis function arguments
85 * \ingroup mofem_base_functions
86 *
87 * This class stores parameters and data required for evaluating Legendre
88 * polynomial basis functions and their derivatives.
89 */
91
92 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
93 UnknownInterface **iface) const;
94
95 int P; ///< Approximation order
96 double *diffS; ///< Derivatives of s parameter
97 int dIm; ///< Dimension
98 boost::shared_ptr<MatrixDouble> baseFunPtr; ///< Pointer to base function matrix
99 boost::shared_ptr<MatrixDouble> baseDiffFunPtr; ///< Pointer to base function derivatives matrix
100
101 /// Function pointer for polynomial evaluation
102 PetscErrorCode (*basePolynomialsType0)(int p, double s, double *diff_s,
103 double *L, double *diffL,
104 const int dim);
105
106 /**
107 * \brief Constructor
108 * \param p approximation order
109 * \param diff_s derivatives of s parameter
110 * \param dim dimension
111 * \param base_fun_ptr pointer to base function matrix
112 * \param base_diff_fun_ptr pointer to base function derivatives matrix
113 */
114 LegendrePolynomialCtx(int p, double *diff_s, int dim,
115 boost::shared_ptr<MatrixDouble> base_fun_ptr,
116 boost::shared_ptr<MatrixDouble> base_diff_fun_ptr)
117 : P(p), diffS(diff_s), dIm(dim), baseFunPtr(base_fun_ptr),
118 baseDiffFunPtr(base_diff_fun_ptr),
121};
122
123/**
124 * \brief Class for calculating Legendre basis functions
125 * \ingroup mofem_base_functions
126 *
127 * This class implements the BaseFunction interface to provide Legendre
128 * polynomial evaluation capabilities for finite element approximations.
129 */
131
132 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
133 UnknownInterface **iface) const;
134
137
138 /**
139 * \brief Evaluate Legendre basis functions at given points
140 * \param pts matrix of evaluation points
141 * \param ctx_ptr context containing evaluation parameters
142 * \return error code
143 */
145 boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
146};
147
148
149
150} // namespace MoFEM
151
152#endif //__LEGENDREPOLYNOMIALS_HPP__
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, int dim)
Calculate Legendre approximation basis.
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
Base class used to exchange data between element data structures and class calculating base functions...
Base class if inherited used to calculate base functions.
Context class for Legendre basis function arguments.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
LegendrePolynomialCtx(int p, double *diff_s, int dim, boost::shared_ptr< MatrixDouble > base_fun_ptr, boost::shared_ptr< MatrixDouble > base_diff_fun_ptr)
Constructor.
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Function pointer for polynomial evaluation.
double * diffS
Derivatives of s parameter.
boost::shared_ptr< MatrixDouble > baseDiffFunPtr
Pointer to base function derivatives matrix.
boost::shared_ptr< MatrixDouble > baseFunPtr
Pointer to base function matrix.
Class for calculating Legendre 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 Legendre basis functions at given points.
base class for all interface classes