KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
convection_diffusion_reaction_element.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Suneth Warnakulasuriya
11 //
12 
13 #if !defined(KRATOS_CONVECTION_DIFFUSION_REACTION_ELEMENT_H_INCLUDED)
14 #define KRATOS_CONVECTION_DIFFUSION_REACTION_ELEMENT_H_INCLUDED
15 
16 // System includes
17 
18 // External includes
19 
20 // Project includes
21 #include "includes/define.h"
22 #include "includes/element.h"
23 
24 // Application includes
25 
26 namespace Kratos
27 {
30 
31 template <
32  unsigned int TDim,
33  unsigned int TNumNodes,
34  class TConvectionDiffusionReactionData>
36 {
37 public:
40 
41  using BaseType = Element;
42 
44  using NodeType = Node;
45 
48 
51 
53  using VectorType = Vector;
54 
56  using MatrixType = Matrix;
57 
58  using IndexType = std::size_t;
59 
60  using EquationIdVectorType = std::vector<IndexType>;
61 
62  using DofsVectorType = std::vector<Dof<double>::Pointer>;
63 
66 
68 
69  using ConvectionDiffusionReactionDataType = TConvectionDiffusionReactionData;
70 
71 
74 
79 
83 
88  IndexType NewId = 0)
89  : Element(NewId)
90  {
91  }
92 
97  IndexType NewId,
98  const NodesArrayType& ThisNodes)
99  : Element(NewId, ThisNodes)
100  {
101  }
102 
107  IndexType NewId,
108  GeometryType::Pointer pGeometry)
109  : Element(NewId, pGeometry)
110  {
111  }
112 
117  IndexType NewId,
118  GeometryType::Pointer pGeometry,
119  typename PropertiesType::Pointer pProperties)
120  : Element(NewId, pGeometry, pProperties)
121  {
122  }
123 
129  : Element(rOther)
130  {
131  }
132 
137 
141 
154  Element::Pointer Create(
155  IndexType NewId,
156  NodesArrayType const& ThisNodes,
157  typename PropertiesType::Pointer pProperties) const override
158  {
159  KRATOS_TRY
160  return Kratos::make_intrusive<CurrentElementType>(
161  NewId, this->GetGeometry().Create(ThisNodes), pProperties);
162  KRATOS_CATCH("");
163  }
164 
172  Element::Pointer Create(
173  IndexType NewId,
174  GeometryType::Pointer pGeom,
175  typename PropertiesType::Pointer pProperties) const override
176  {
177  KRATOS_TRY
178  return Kratos::make_intrusive<CurrentElementType>(NewId, pGeom, pProperties);
179  KRATOS_CATCH("");
180  }
181 
189  Element::Pointer Clone(
190  IndexType NewId,
191  NodesArrayType const& ThisNodes) const override
192  {
193  KRATOS_TRY
194  return Kratos::make_intrusive<CurrentElementType>(
195  NewId, this->GetGeometry().Create(ThisNodes), this->pGetProperties());
196  KRATOS_CATCH("");
197  }
198 
199  void EquationIdVector(
200  EquationIdVectorType& rResult,
201  const ProcessInfo& CurrentProcessInfo) const override;
202 
208  void GetDofList(
209  DofsVectorType& rElementalDofList,
210  const ProcessInfo& CurrentProcessInfo) const override;
211 
212  void GetValuesVector(
213  Vector& rValues,
214  int Step = 0) const override;
215 
217  Vector& rValues,
218  int Step = 0) const override;
219 
221  Vector& rValues,
222  int Step = 0) const override;
223 
240  MatrixType& rLeftHandSideMatrix,
241  VectorType& rRightHandSideVector,
242  const ProcessInfo& rCurrentProcessInfo) override;
243 
251  VectorType& rRightHandSideVector,
252  const ProcessInfo& rCurrentProcessInfo) override;
253 
261  MatrixType& rDampingMatrix,
262  VectorType& rRightHandSideVector,
263  const ProcessInfo& rCurrentProcessInfo) override;
264 
280  void CalculateMassMatrix(
281  MatrixType& rMassMatrix,
282  const ProcessInfo& rCurrentProcessInfo) override;
290  MatrixType& rDampingMatrix,
291  const ProcessInfo& rCurrentProcessInfo) override;
292 
302  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
303 
305 
309 
311  std::string Info() const override
312  {
313  std::stringstream buffer;
314  buffer << "ConvectionDiffusionReactionElement #" << Id();
315  return buffer.str();
316  }
318  void PrintInfo(std::ostream& rOStream) const override
319  {
320  rOStream << "CDR" << TConvectionDiffusionReactionData::GetName();
321  }
322 
324 
325 protected:
328 
335  void GetValuesArray(
337  const int Step = 0) const;
338 
346  virtual void CalculateGeometryData(
347  Vector& rGaussWeights,
348  Matrix& rNContainer,
349  ShapeFunctionDerivativesArrayType& rDN_DX) const;
350 
351  void AddLumpedMassMatrix(
352  Matrix& rMassMatrix,
353  const double Mass) const;
354 
356  Matrix& rDampingMatrix,
357  const double ReactionTerm,
358  const double EffectiveKinematicViscosity,
359  const Vector& rVelocityConvectiveTerms,
360  const double GaussWeight,
361  const Vector& rGaussShapeFunctions,
362  const Matrix& rGaussdNa_dNb) const;
363 
365 
366 private:
369 
370  friend class Serializer;
371 
372  void save(Serializer& rSerializer) const override
373  {
374  KRATOS_TRY
375 
377 
378  KRATOS_CATCH("");
379  }
380  void load(Serializer& rSerializer) override
381  {
382  KRATOS_TRY
383 
385 
386  KRATOS_CATCH("");
387  }
388 
390 
391 }; // Class ConvectionDiffusionReactionElement
392 
396 
397 template <unsigned int TDim, unsigned int TNumNodes, class TConvectionDiffusionReactionData>
398 inline std::istream& operator>>(
399  std::istream& rIStream,
401 
403 template <unsigned int TDim, unsigned int TNumNodes, class TConvectionDiffusionReactionData>
404 inline std::ostream& operator<<(
405  std::ostream& rOStream,
407 {
408  rThis.PrintInfo(rOStream);
409  rOStream << " : " << std::endl;
410  rThis.PrintData(rOStream);
411  return rOStream;
412 }
413 
415 
416 } // namespace Kratos.
417 
418 #endif // KRATOS_CONVECTION_DIFFUSION_REACTION_ELEMENT_H_INCLUDED defined
Definition: convection_diffusion_reaction_element.h:36
void AddDampingMatrixGaussPointContributions(Matrix &rDampingMatrix, const double ReactionTerm, const double EffectiveKinematicViscosity, const Vector &rVelocityConvectiveTerms, const double GaussWeight, const Vector &rGaussShapeFunctions, const Matrix &rGaussdNa_dNb) const
Definition: convection_diffusion_reaction_element.cpp:317
void AddLumpedMassMatrix(Matrix &rMassMatrix, const double Mass) const
Definition: convection_diffusion_reaction_element.cpp:307
void GetValuesVector(Vector &rValues, int Step=0) const override
Definition: convection_diffusion_reaction_element.cpp:74
Element::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
Definition: convection_diffusion_reaction_element.h:189
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, typename PropertiesType::Pointer pProperties) const override
Definition: convection_diffusion_reaction_element.h:154
TConvectionDiffusionReactionData ConvectionDiffusionReactionDataType
Definition: convection_diffusion_reaction_element.h:69
Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, typename PropertiesType::Pointer pProperties) const override
Definition: convection_diffusion_reaction_element.h:172
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: convection_diffusion_reaction_element.cpp:134
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &CurrentProcessInfo) const override
Definition: convection_diffusion_reaction_element.cpp:40
ConvectionDiffusionReactionElement(ConvectionDiffusionReactionElement const &rOther)
Definition: convection_diffusion_reaction_element.h:127
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: convection_diffusion_reaction_element.h:56
ConvectionDiffusionReactionElement(IndexType NewId=0)
Definition: convection_diffusion_reaction_element.h:87
std::string Info() const override
Turn back information as a string.
Definition: convection_diffusion_reaction_element.h:311
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ConvectionDiffusionReactionElement)
ConvectionDiffusionReactionElement(IndexType NewId, GeometryType::Pointer pGeometry, typename PropertiesType::Pointer pProperties)
Definition: convection_diffusion_reaction_element.h:116
std::vector< IndexType > EquationIdVectorType
Definition: convection_diffusion_reaction_element.h:60
ConvectionDiffusionReactionElement(IndexType NewId, const NodesArrayType &ThisNodes)
Definition: convection_diffusion_reaction_element.h:96
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: convection_diffusion_reaction_element.cpp:199
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: convection_diffusion_reaction_element.h:65
void GetSecondDerivativesVector(Vector &rValues, int Step=0) const override
Definition: convection_diffusion_reaction_element.cpp:114
void CalculateDampingMatrix(MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: convection_diffusion_reaction_element.cpp:227
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &CurrentProcessInfo) const override
Definition: convection_diffusion_reaction_element.cpp:57
virtual void CalculateGeometryData(Vector &rGaussWeights, Matrix &rNContainer, ShapeFunctionDerivativesArrayType &rDN_DX) const
Calculates shape function data for this element.
Definition: convection_diffusion_reaction_element.cpp:294
void GetValuesArray(BoundedVector< double, TNumNodes > &rValues, const int Step=0) const
Get the Values Array.
Definition: convection_diffusion_reaction_element.cpp:82
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: convection_diffusion_reaction_element.cpp:274
~ConvectionDiffusionReactionElement() override=default
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: convection_diffusion_reaction_element.h:62
void GetFirstDerivativesVector(Vector &rValues, int Step=0) const override
Definition: convection_diffusion_reaction_element.cpp:99
ConvectionDiffusionReactionElement(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: convection_diffusion_reaction_element.h:106
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: convection_diffusion_reaction_element.cpp:151
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: convection_diffusion_reaction_element.h:318
void CalculateLocalVelocityContribution(MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
CalculateLocalVelocityContribution Calculate the local contribution in terms of velocity and pressure...
Definition: convection_diffusion_reaction_element.cpp:190
GeometryData::IntegrationMethod GetIntegrationMethod() const override
Definition: convection_diffusion_reaction_element.cpp:288
Base class for all Elements.
Definition: element.h:60
Properties PropertiesType
Definition: element.h:80
Element(IndexType NewId=0)
Definition: element.h:121
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: element.h:1135
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
PropertiesType::Pointer pGetProperties()
returns the pointer to the property of the element. Does not throw an error, to allow copying of elem...
Definition: element.h:1014
std::size_t IndexType
Definition: flags.h:74
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
GeometryData::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: geometry.h:189
IndexType Id() const
Definition: indexed_object.h:107
This class defines the node.
Definition: node.h:65
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
def load(f)
Definition: ode_solve.py:307