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.
calculate_component_gradient_simplex_element.h
Go to the documentation of this file.
1 //
2 // Project Name: Kratos
3 // Last Modified by: $Author: gcasas $
4 // Date: $Date: 2016-03-12
5 //
6 
7 #if !defined(KRATOS_COMPUTE_COMPONENT_GRADIENT_SIMPLEX_ELEMENT_H_INCLUDED )
8 #define KRATOS_COMPUTE_COMPONENT_GRADIENT_SIMPLEX_ELEMENT_H_INCLUDED
9 
10 // System includes
11 #include <string>
12 #include <iostream>
13 
14 
15 // External includes
16 
17 
18 // Project includes
19 #include "containers/array_1d.h"
20 #include "includes/define.h"
21 #include "includes/element.h"
22 #include "includes/serializer.h"
23 #include "geometries/geometry.h"
24 #include "utilities/math_utils.h"
25 #include "utilities/geometry_utilities.h"
26 
27 // Application includes
28 #include "includes/variables.h"
31 
32 namespace Kratos
33 {
34 
37 
40 
44 
48 
52 
56 
58 
60 template< unsigned int TDim,
61  unsigned int TNumNodes = TDim + 1 >
62 class KRATOS_API(SWIMMING_DEM_APPLICATION) ComputeComponentGradientSimplex : public ComputeMaterialDerivativeSimplex<TDim, TNumNodes>
63 {
64 public:
67 
70 
72 
74  typedef Node NodeType;
75 
78 
81 
83  typedef Vector VectorType;
84 
86  typedef Matrix MatrixType;
87 
89 
90  typedef std::size_t IndexType;
91 
92  typedef std::size_t SizeType;
93 
94  typedef std::vector<std::size_t> EquationIdVectorType;
95 
96  typedef std::vector< Dof<double>::Pointer > DofsVectorType;
97 
99 
102 
105 
108 
112 
113  //Constructors.
114 
116 
120  BaseType(NewId), mCurrentComponent('X')
121  {}
122 
124 
129  BaseType(NewId, ThisNodes), mCurrentComponent('X')
130  {}
131 
133 
137  ComputeComponentGradientSimplex(IndexType NewId, GeometryType::Pointer pGeometry) :
138  BaseType(NewId, pGeometry), mCurrentComponent('X')
139  {}
140 
142 
147  ComputeComponentGradientSimplex(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) :
148  BaseType(NewId, pGeometry, pProperties), mCurrentComponent('X')
149  {}
150 
153  {}
154 
155 
159 
160 
164 
166 
173  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes,
174  PropertiesType::Pointer pProperties) const override
175  {
176  return Element::Pointer(new ComputeComponentGradientSimplex(NewId, this->GetGeometry().Create(ThisNodes), pProperties));
177  }
178 
180  virtual void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
181  VectorType& rRightHandSideVector,
182  const ProcessInfo& rCurrentProcessInfo) override;
183 
184 
186 
192  virtual void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
193 
195 
199  virtual void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
200 
204 
208 
210 
218  virtual int Check(const ProcessInfo& rCurrentProcessInfo) const override;
219 
223 
224 
228 
230  virtual std::string Info() const override
231  {
232  std::stringstream buffer;
233  buffer << "ComputeComponentGradientSimplex #" << this->Id();
234  return buffer.str();
235  }
236 
238  virtual void PrintInfo(std::ostream& rOStream) const override
239  {
240  rOStream << "ComputeComponentGradientSimplex" << TDim << "D";
241  }
242 
243 // /// Print object's data.
244 // virtual void PrintData(std::ostream& rOStream) const;
245 
249 
250 
252 
253 protected:
256 
257 
262 
266 
267 
271 
272 
276 
277 
281 
282 
286 
287 
289 
290 private:
293 
297 
301 
302  friend class Serializer;
303 
307 
308  void AddIntegrationPointRHSContribution(VectorType& F,
309  const array_1d<double,TNumNodes>& rShapeFunc,
310  const BoundedMatrix<double, TNumNodes, TDim>& rShapeDeriv,
311  const double Weight) override;
315 
316 
320 
321 
325 
326 
330 
333 
336 
338 
339 }; // Class ComputeComponentGradientSimplex
340 
342 
345 
346 
350 
351 
353 template< unsigned int TDim >
354 inline std::istream& operator >>(std::istream& rIStream,
356 {
357  return rIStream;
358 }
359 
361 template< unsigned int TDim >
362 inline std::ostream& operator <<(std::ostream& rOStream,
364 {
365  rThis.PrintInfo(rOStream);
366  rOStream << std::endl;
367  rThis.PrintData(rOStream);
368 
369  return rOStream;
370 }
372 
374 
375 } // namespace Kratos.
376 
377 #endif // KRATOS_COMPUTE_COMPONENT_GRADIENT_SIMPLEX_ELEMENT_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
A post-processing element to recover the Laplacian from the velocity solution.
Definition: calculate_component_gradient_simplex_element.h:63
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition of nodes container type, redefined from GeometryType.
Definition: calculate_component_gradient_simplex_element.h:80
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: calculate_component_gradient_simplex_element.h:101
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new element of this type.
Definition: calculate_component_gradient_simplex_element.h:173
ComputeComponentGradientSimplex(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: calculate_component_gradient_simplex_element.h:128
std::size_t SizeType
Definition: calculate_component_gradient_simplex_element.h:92
ComputeComponentGradientSimplex(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constuctor using geometry and properties.
Definition: calculate_component_gradient_simplex_element.h:147
ComputeMaterialDerivativeSimplex< TDim, TNumNodes > BaseType
Definition: calculate_component_gradient_simplex_element.h:71
std::size_t IndexType
Definition: calculate_component_gradient_simplex_element.h:90
virtual ~ComputeComponentGradientSimplex()
Destructor.
Definition: calculate_component_gradient_simplex_element.h:152
ComputeComponentGradientSimplex(IndexType NewId=0)
Default constuctor.
Definition: calculate_component_gradient_simplex_element.h:119
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: calculate_component_gradient_simplex_element.h:86
ComputeComponentGradientSimplex(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using a geometry object.
Definition: calculate_component_gradient_simplex_element.h:137
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: calculate_component_gradient_simplex_element.h:107
Geometry< NodeType > GeometryType
Geometry type (using with given NodeType)
Definition: calculate_component_gradient_simplex_element.h:77
Node NodeType
Node type (default is: Node)
Definition: calculate_component_gradient_simplex_element.h:74
std::vector< std::size_t > EquationIdVectorType
Definition: calculate_component_gradient_simplex_element.h:94
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: calculate_component_gradient_simplex_element.h:96
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: calculate_component_gradient_simplex_element.h:98
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ComputeComponentGradientSimplex)
Pointer definition of ComputeComponentGradientSimplex.
Vector VectorType
Vector type for local contributions to the linear system.
Definition: calculate_component_gradient_simplex_element.h:83
virtual std::string Info() const override
Turn back information as a string.
Definition: calculate_component_gradient_simplex_element.h:230
Kratos::Matrix ShapeFunctionDerivativesType
Type for a matrix containing the shape function gradients.
Definition: calculate_component_gradient_simplex_element.h:104
char mCurrentComponent
Definition: calculate_component_gradient_simplex_element.h:261
Properties PropertiesType
Definition: calculate_component_gradient_simplex_element.h:88
virtual void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: calculate_component_gradient_simplex_element.h:238
A post-processing element to recover the Laplacian from the velocity solution.
Definition: calculate_mat_deriv_simplex_element.h:62
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: element.h:1135
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
Geometry base class.
Definition: geometry.h:71
This object defines an indexed object.
Definition: indexed_object.h:54
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
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.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
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
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
F
Definition: hinsberg_optimization.py:144