v0.14.0
Loading...
Searching...
No Matches
LobattoPolynomial.hpp
Go to the documentation of this file.
1/** \file LobattoPolynomial.hpp
2\brief Implementation of Lobatto polynomial
3
4*/
5
6
7
8#ifndef __LOBATTOPOLYNOMIALS_HPP__
9#define __LOBATTOPOLYNOMIALS_HPP__
10
11namespace MoFEM {
12
13/**
14 * \brief Class used to give arguments to Lobatto base functions
15 * \ingroup mofem_base_functions
16 */
18
19 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
20 UnknownInterface **iface) const;
21
22 LobattoPolynomialCtx(int p, double *diff_s, int dim,
23 boost::shared_ptr<MatrixDouble> base_fun_ptr,
24 boost::shared_ptr<MatrixDouble> base_diff_fun_ptr)
25 : LegendrePolynomialCtx(p, diff_s, dim, base_fun_ptr, base_diff_fun_ptr) {
27 }
29};
30
31/**
32 * \brief Calculating Lobatto base functions
33 * \ingroup mofem_base_functions
34 */
36
37 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
38 UnknownInterface **iface) const;
39
42
44 boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
45};
46
47/**
48 * \brief Class used to give arguments to Kernel Lobatto base functions
49 * \ingroup mofem_base_functions
50 */
52
53 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
54 UnknownInterface **iface) const;
55
56 KernelLobattoPolynomialCtx(int p, double *diff_s, int dim,
57 boost::shared_ptr<MatrixDouble> base_fun_ptr,
58 boost::shared_ptr<MatrixDouble> base_diff_fun_ptr)
59 : LegendrePolynomialCtx(p, diff_s, dim, base_fun_ptr, base_diff_fun_ptr) {
61 }
63};
64
65/**
66 * \brief Calculating Lobatto base functions
67 * \ingroup mofem_base_functions
68 */
70
71 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
72 UnknownInterface **iface) const;
73
76
78 boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
79};
80
81} // namespace MoFEM
82
83#endif
static Index< 'p', 3 > p
const int dim
PetscErrorCode Lobatto_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Lobatto base 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
Class used to give arguments to Kernel Lobatto base functions.
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)
Calculating Lobatto base functions.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Class used to give arguments to Legendre base functions.
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculating Legendre base functions.
Class used to give arguments to Lobatto base functions.
LobattoPolynomialCtx(int p, double *diff_s, int dim, boost::shared_ptr< MatrixDouble > base_fun_ptr, boost::shared_ptr< MatrixDouble > base_diff_fun_ptr)
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Calculating Lobatto base functions.
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
base class for all interface classes