NuTo
Numerics Tool
InterpolationBrickLinear.h
Go to the documentation of this file.
1 #pragma once
4 
5 namespace NuTo
6 {
8 {
9 public:
10  static Eigen::Matrix<double, 3, 1> LocalCoords(int rNodeIndex);
11 
12  static Eigen::Matrix<double, 8, 1> ShapeFunctions(const Eigen::VectorXd& rCoordinates);
13 
14  static Eigen::Matrix<double, 8, 3> DerivativeShapeFunctions(const Eigen::VectorXd& rCoordinates);
15 
16  std::unique_ptr<InterpolationSimple> Clone() const override;
17 
18  Eigen::VectorXd GetShapeFunctions(const NaturalCoords& naturalIpCoords) const override;
19 
20  Eigen::MatrixXd GetDerivativeShapeFunctions(const NaturalCoords& naturalIpCoords) const override;
21 
22  NaturalCoords GetLocalCoords(int nodeId) const override;
23 
24  int GetNumNodes() const override;
25 
26  const Shape& GetShape() const override;
27 
28 private:
29  Hexahedron mShape;
30 };
31 } /* NuTo */
Eigen::MatrixXd GetDerivativeShapeFunctions(const NaturalCoords &naturalIpCoords) const override
calculates the derivative shape functions
Definition: InterpolationBrickLinear.cpp:17
static Eigen::Matrix< double, 8, 1 > ShapeFunctions(const Eigen::VectorXd &rCoordinates)
Definition: InterpolationBrickLinear.cpp:63
NaturalCoords GetLocalCoords(int nodeId) const override
returns the local node coordinates
Definition: InterpolationBrickLinear.cpp:22
Base class for the interpolation. The derived classes provide information about the actual interpolat...
Definition: InterpolationSimple.h:12
std::unique_ptr< InterpolationSimple > Clone() const override
Definition: InterpolationBrickLinear.cpp:6
static Eigen::Matrix< double, 8, 3 > DerivativeShapeFunctions(const Eigen::VectorXd &rCoordinates)
Definition: InterpolationBrickLinear.cpp:87
static Eigen::Matrix< double, 3, 1 > LocalCoords(int rNodeIndex)
Definition: InterpolationBrickLinear.cpp:37
int GetNumNodes() const override
returns the number of nodes
Definition: InterpolationBrickLinear.cpp:27
const Shape & GetShape() const override
Definition: InterpolationBrickLinear.cpp:32
Definition: Shape.h:20
Eigen::VectorXd GetShapeFunctions(const NaturalCoords &naturalIpCoords) const override
calculates the shape functions
Definition: InterpolationBrickLinear.cpp:11
Definition: Exception.h:6
Definition: SerializeStreamOut.h:9
Definition: Hexahedron.h:8
Definition: InterpolationBrickLinear.h:7