v0.9.0
LegendrePolynomial.hpp
Go to the documentation of this file.
1 /** \file LegendrePolynomial.hpp
2 \brief Implementation of Legendre polynomial
3 
4 */
5 
6 /* This file is part of MoFEM.
7  * MoFEM is free software: you can redistribute it and/or modify it under
8  * the terms of the GNU Lesser General Public License as published by the
9  * Free Software Foundation, either version 3 of the License, or (at your
10  * option) any later version.
11  *
12  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
13  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  * License for more details.
16  *
17  * You should have received a copy of the GNU Lesser General Public
18  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
19 
20 #ifndef __LEGENDREPOLYNOMIALS_HPP__
21 #define __LEGENDREPOLYNOMIALS_HPP__
22 
23 namespace MoFEM {
24 
27 
28 /**
29  * \brief Class used to give arguments to Legendre base functions
30  * \ingroup mofem_base_functions
31  */
33 
35  BaseFunctionUnknownInterface **iface) const;
36 
37  int P;
38  double *diffS;
39  int dIm;
40  boost::shared_ptr<MatrixDouble> baseFunPtr;
41  boost::shared_ptr<MatrixDouble> baseDiffFunPtr;
42 
43  PetscErrorCode (*basePolynomialsType0)(int p, double s, double *diff_s,
44  double *L, double *diffL,
45  const int dim);
46 
47  LegendrePolynomialCtx(int p, double *diff_s, int dim,
48  boost::shared_ptr<MatrixDouble> base_fun_ptr,
49  boost::shared_ptr<MatrixDouble> base_diff_fun_ptr)
50  : P(p), diffS(diff_s), dIm(dim), baseFunPtr(base_fun_ptr),
51  baseDiffFunPtr(base_diff_fun_ptr),
54 };
55 
56 /**
57  * \brief Calculating Legendre base functions
58  * \ingroup mofem_base_functions
59  */
61 
63  BaseFunctionUnknownInterface **iface) const;
64 
67 
69  boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
70 };
71 
72 } // namespace MoFEM
73 
74 #endif //__LEGENDREPOLYNOMIALS_HPP__
Calculating Legendre base functions.
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
MoFEM interface unique ID.
boost::shared_ptr< MatrixDouble > baseFunPtr
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:74
static const MOFEMuuid IDD_LEGENDRE_BASE_FUNCTION
std::bitset< BITINTERFACEUID_SIZE > BitIntefaceId
Definition: Types.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
MoFEMErrorCode query_interface(const MOFEMuuid &uuid, BaseFunctionUnknownInterface **iface) const
Base class if inherited used to calculate base functions.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:66
boost::shared_ptr< MatrixDouble > baseDiffFunPtr
MoFEMErrorCode query_interface(const MOFEMuuid &uuid, BaseFunctionUnknownInterface **iface) const
Base class used to exchange data between element data structures and class calculating base functions...
LegendrePolynomialCtx(int p, double *diff_s, int dim, boost::shared_ptr< MatrixDouble > base_fun_ptr, boost::shared_ptr< MatrixDouble > base_diff_fun_ptr)
Class used to give arguments to Legendre base functions.