NuTo
Numerics Tool
ElementFem.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <vector>
8 
9 namespace NuTo
10 {
11 
13 {
14 public:
15  ElementFem(std::vector<NodeSimple*> nodes, const InterpolationSimple& interpolation)
16  : mInterpolation(interpolation)
17  , mShape(interpolation.GetShape())
18  {
19  for (NodeSimple* node : nodes)
20  mNodes.push_back(*node);
21  assert(static_cast<int>(mNodes.size()) == interpolation.GetNumNodes());
22  }
23 
24  ElementFem(std::initializer_list<std::reference_wrapper<NuTo::NodeSimple>> nodes,
25  const InterpolationSimple& interpolation)
26  : mNodes(nodes)
27  , mInterpolation(interpolation)
28  , mShape(interpolation.GetShape())
29  {
30  assert(static_cast<int>(mNodes.size()) == interpolation.GetNumNodes());
31  }
32 
33  virtual Eigen::VectorXd ExtractNodeValues(int instance = 0) const override
34  {
35  const int dim = GetDofDimension();
36  Eigen::VectorXd nodeValues(GetNumNodes() * dim);
37  for (int i = 0; i < GetNumNodes(); ++i)
38  nodeValues.segment(dim * i, dim) = GetNode(i).GetValues(instance);
39  return nodeValues;
40  }
41 
42  Eigen::MatrixXd GetNMatrix(NaturalCoords ipCoords) const override
43  {
45  }
46 
47  Eigen::VectorXd GetShapeFunctions(NaturalCoords ipCoords) const override
48  {
49  return Interpolation().GetShapeFunctions(ipCoords);
50  }
51 
52  Eigen::MatrixXd GetDerivativeShapeFunctions(NaturalCoords ipCoords) const override
53  {
54  return Interpolation().GetDerivativeShapeFunctions(ipCoords);
55  }
56 
57  int GetDofDimension() const override
58  {
59  return GetNode(0).GetNumValues();
60  }
61 
62  Eigen::VectorXi GetDofNumbering() const override
63  {
64  Eigen::VectorXi dofNumbering(GetNumNodes() * GetDofDimension());
65  int i = 0;
66  for (int iNode = 0; iNode < GetNumNodes(); ++iNode)
67  {
68  const auto& node = GetNode(iNode);
69  for (int iDof = 0; iDof < GetDofDimension(); ++iDof)
70  {
71  dofNumbering[i++] = node.GetDofNumber(iDof);
72  }
73  }
74  return dofNumbering;
75  }
76 
77  int GetNumNodes() const override
78  {
79  return mNodes.size();
80  }
81 
83  {
84  return mInterpolation;
85  }
86 
88  {
89  assert(i < static_cast<int>(mNodes.size()));
90  return mNodes[i];
91  }
92 
93 
94  const NodeSimple& GetNode(int i) const
95  {
96  assert(i < static_cast<int>(mNodes.size()));
97  return mNodes[i];
98  }
99 
100  const Shape& GetShape() const
101  {
102  return mShape;
103  }
104 
105 private:
106  std::vector<std::reference_wrapper<NodeSimple>> mNodes;
107  std::reference_wrapper<const InterpolationSimple> mInterpolation;
108  const Shape& mShape;
109 };
110 } /* NuTo */
int GetNumNodes() const override
Definition: ElementFem.h:77
Eigen::VectorXd GetShapeFunctions(NaturalCoords ipCoords) const override
Definition: ElementFem.h:47
Definition: ElementFem.h:12
virtual int GetNumNodes() const =0
returns the number of nodes
Store node values and its dof.
Definition: NodeSimple.h:11
const InterpolationSimple & Interpolation() const
Definition: ElementFem.h:82
NodeSimple & GetNode(int i)
Definition: ElementFem.h:87
Definition: ElementInterface.h:8
virtual Eigen::VectorXd GetShapeFunctions(const NaturalCoords &naturalIpCoords) const =0
calculates the shape functions
int GetNumValues() const
Definition: NodeSimple.h:82
Eigen::MatrixXd GetNMatrix(NaturalCoords ipCoords) const override
Definition: ElementFem.h:42
Base class for the interpolation. The derived classes provide information about the actual interpolat...
Definition: InterpolationSimple.h:12
const Shape & GetShape() const
Definition: ElementFem.h:100
virtual Eigen::MatrixXd GetDerivativeShapeFunctions(const NaturalCoords &naturalIpCoords) const =0
calculates the derivative shape functions
node
Definition: Brick8NCoupling.py:104
const NodeSimple & GetNode(int i) const
Definition: ElementFem.h:94
ElementFem(std::initializer_list< std::reference_wrapper< NuTo::NodeSimple >> nodes, const InterpolationSimple &interpolation)
Definition: ElementFem.h:24
Definition: Shape.h:20
virtual Eigen::VectorXd ExtractNodeValues(int instance=0) const override
extracts all node values of this element
Definition: ElementFem.h:33
Definition: Exception.h:6
const Eigen::VectorXd & GetValues(int instance=0) const
Definition: NodeSimple.h:50
Eigen::MatrixXd GetDerivativeShapeFunctions(NaturalCoords ipCoords) const override
Definition: ElementFem.h:52
constexpr int dim
Definition: ConstraintNodeToElement2D.cpp:24
int GetDofDimension() const override
Definition: ElementFem.h:57
Eigen::MatrixXd N(const Eigen::VectorXd &shapeFunctions, int numNodes, int dim)
Definition: Matrix.h:8
Eigen::VectorXi GetDofNumbering() const override
extract the dof numbers from its nodes.
Definition: ElementFem.h:62
ElementFem(std::vector< NodeSimple * > nodes, const InterpolationSimple &interpolation)
Definition: ElementFem.h:15