NuTo
Numerics Tool
CellIpData.h
Go to the documentation of this file.
1 #pragma once
2 
8 
9 namespace NuTo
10 {
11 
15 {
16 public:
22  CellIpData(const CellData& cellData, NuTo::Jacobian jacobian, NaturalCoords ipCoords, int ipId)
23  : mCellData(cellData)
24  , mJacobian(std::move(jacobian))
25  , mIPCoords(std::move(ipCoords))
26  , mIpId(ipId)
27  {
28  }
29 
32  Eigen::VectorXd GlobalCoordinates() const
33  {
34  return Interpolate(mCellData.Elements().CoordinateElement(), mIPCoords);
35  }
36 
39  CellIds Ids() const
40  {
41  return {mCellData.GetCellId(), mIpId};
42  }
43 
48  Eigen::VectorXd Value(DofType dofType, int instance = 0) const
49  {
50  return N(dofType) * NodeValueVector(dofType, instance);
51  }
52 
53  const NuTo::Jacobian& GetJacobian() const
54  {
55  return mJacobian;
56  }
57 
58  double Value(ScalarDofType dofType, int instance = 0) const
59  {
60  return Value(DofType(dofType), instance)[0];
61  }
62 
70  Eigen::VectorXd Apply(DofType dofType, const Nabla::Interface& b, int instance = 0) const
71  {
72  return B(dofType, b) * NodeValueVector(dofType, instance);
73  }
74 
77  const Eigen::MatrixXd& N(DofType dofType) const
78  {
79  Eigen::MatrixXd& N = mNs[dofType];
80  if (N.size() == 0) // simplest memoization using a mutable mNs to keep it const
81  N = mCellData.Elements().DofElement(dofType).GetNMatrix(mIPCoords);
82  return N;
83  }
84 
89  const Eigen::MatrixXd& B(DofType dofType, const Nabla::Interface& b) const
90  {
91  Eigen::MatrixXd& B = mBs[dofType];
92  if (B.size() == 0) // simplest memoization using a mutable mBs to keep it const
93  B = b(CalculateDerivativeShapeFunctionsGlobal(dofType));
94  return B;
95  }
96 
100  const Eigen::VectorXd& NodeValueVector(DofType dofType, int instance = 0) const
101  {
102  return mCellData.GetNodeValues(dofType, instance);
103  }
104 
105 private:
108  Eigen::MatrixXd CalculateDerivativeShapeFunctionsGlobal(DofType dofType) const
109  {
110  Eigen::MatrixXd dShapeNatural =
111  mCellData.Elements().DofElement(dofType).GetDerivativeShapeFunctions(mIPCoords);
112  return mJacobian.TransformDerivativeShapeFunctions(dShapeNatural);
113  }
114 
115  const CellData& mCellData;
116  NuTo::Jacobian mJacobian;
117  NaturalCoords mIPCoords;
118  int mIpId;
119 
120  mutable DofContainer<Eigen::MatrixXd> mBs;
121  mutable DofContainer<Eigen::MatrixXd> mNs;
122 };
123 } /* NuTo */
virtual const ElementInterface & DofElement(DofType) const =0
Definition: DifferentialOperators.h:11
int GetCellId() const
Definition: CellData.h:23
Definition: Jacobian.h:10
named pair for cellId and ipId to make the argument list shorter and avoid accidental mixup of both ...
Definition: CellIds.h:7
const NuTo::Jacobian & GetJacobian() const
Definition: CellIpData.h:53
const Eigen::MatrixXd & B(DofType dofType, const Nabla::Interface &b) const
Returns a memoized copy of the B matrix for a given dof type.
Definition: CellIpData.h:89
Eigen::VectorXd Interpolate(const ElementInterface &element, NaturalCoords ipCoords)
Definition: ElementInterface.h:31
Similar to NuTo::CellData.
Definition: CellIpData.h:14
CellIds Ids() const
Access to the cellId and ipId, compressed in CellIds.
Definition: CellIpData.h:39
Definition: DofType.h:8
Eigen::MatrixXd TransformDerivativeShapeFunctions(const Eigen::MatrixXd &global) const
Definition: Jacobian.h:100
CellIpData(const CellData &cellData, NuTo::Jacobian jacobian, NaturalCoords ipCoords, int ipId)
ctor
Definition: CellIpData.h:22
virtual Eigen::MatrixXd GetNMatrix(NaturalCoords ipCoords) const =0
virtual const ElementInterface & CoordinateElement() const =0
Eigen::VectorXd Value(DofType dofType, int instance=0) const
Calculates the value of a dof at the integration point.
Definition: CellIpData.h:48
Eigen::VectorXd Apply(DofType dofType, const Nabla::Interface &b, int instance=0) const
Calculates the gradient (derivative of the value with respect to x) for a given dof type at the integ...
Definition: CellIpData.h:70
const ElementCollection & Elements() const
Definition: CellData.h:45
virtual Eigen::MatrixXd GetDerivativeShapeFunctions(NaturalCoords ipCoords) const =0
double Value(ScalarDofType dofType, int instance=0) const
Definition: CellIpData.h:58
Definition: DofType.h:44
Definition: Exception.h:6
const Eigen::VectorXd & GetNodeValues(DofType dofType, int instance=0) const
Definition: CellData.h:33
Extracts &#39;cell data&#39; like nodal values from the cell.
Definition: CellData.h:14
Eigen::VectorXd GlobalCoordinates() const
Caluclate the global integration point coordinates.
Definition: CellIpData.h:32
const Eigen::MatrixXd & N(DofType dofType) const
Returns a memoized copy of the N matrix for a given dof type.
Definition: CellIpData.h:77
const Eigen::VectorXd & NodeValueVector(DofType dofType, int instance=0) const
Returns memoized nodal values.
Definition: CellIpData.h:100