NuTo
Numerics Tool
InterpolationTrussLobatto.h
Go to the documentation of this file.
1 #pragma once
2 #include <eigen3/Eigen/Core>
4 #include "nuto/math/Legendre.h"
6 #include "nuto/base/Exception.h"
7 
8 namespace NuTo
9 {
11 {
12 public:
13 
14  static Eigen::VectorXd LocalCoords(int order);
15 
16  static Eigen::VectorXd BarycentricWeights(const Eigen::VectorXd& nodes);
17 
18  static Eigen::VectorXd ShapeFunctions(const double x, const Eigen::VectorXd& nodes);
19 
20  static Eigen::VectorXd DerivativeShapeFunctions(const double x, const Eigen::VectorXd& nodes);
21 
22  InterpolationTrussLobatto(int order);
23 
24  std::unique_ptr<InterpolationSimple> Clone() const override;
25 
26  Eigen::VectorXd GetShapeFunctions(const NaturalCoords& naturalIpCoords) const override;
27 
28  Eigen::MatrixXd GetDerivativeShapeFunctions(const NaturalCoords& naturalIpCoords) const override;
29 
30  NaturalCoords GetLocalCoords(int nodeId) const override;
31 
32  int GetNumNodes() const override;
33 
34  const Shape& GetShape() const override;
35 
36 private:
37  Eigen::VectorXd mNodes;
38  Line mShape;
39 };
40 } /* NuTo */
static Eigen::VectorXd ShapeFunctions(const double x, const Eigen::VectorXd &nodes)
Definition: InterpolationTrussLobatto.cpp:43
std::unique_ptr< InterpolationSimple > Clone() const override
Definition: InterpolationTrussLobatto.cpp:110
Definition: Line.h:8
static Eigen::VectorXd BarycentricWeights(const Eigen::VectorXd &nodes)
Definition: InterpolationTrussLobatto.cpp:28
NaturalCoords GetLocalCoords(int nodeId) const override
returns the local node coordinates
Definition: InterpolationTrussLobatto.cpp:138
static Eigen::VectorXd LocalCoords(int order)
Definition: InterpolationTrussLobatto.cpp:10
Eigen::MatrixXd GetDerivativeShapeFunctions(const NaturalCoords &naturalIpCoords) const override
calculates the derivative shape functions
Definition: InterpolationTrussLobatto.cpp:126
int GetNumNodes() const override
returns the number of nodes
Definition: InterpolationTrussLobatto.cpp:145
const Shape & GetShape() const override
Definition: InterpolationTrussLobatto.cpp:150
Base class for the interpolation. The derived classes provide information about the actual interpolat...
Definition: InterpolationSimple.h:12
InterpolationTrussLobatto(int order)
Definition: InterpolationTrussLobatto.cpp:105
static Eigen::VectorXd DerivativeShapeFunctions(const double x, const Eigen::VectorXd &nodes)
Definition: InterpolationTrussLobatto.cpp:79
Definition: Shape.h:20
Definition: InterpolationTrussLobatto.h:10
Eigen::VectorXd GetShapeFunctions(const NaturalCoords &naturalIpCoords) const override
calculates the shape functions
Definition: InterpolationTrussLobatto.cpp:115
Definition: Exception.h:6