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_laplacian_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_LAPLACIAN_SIMPLEX_ELEMENT_H_INCLUDED )
8 #define KRATOS_COMPUTE_LAPLACIAN_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) ComputeLaplacianSimplex : 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 
109  using BaseType::AddIntegrationPointRHSContribution;
110 
114 
115  //Constructors.
116 
118 
122  BaseType(NewId)
123  {}
124 
126 
131  BaseType(NewId, ThisNodes)
132  {}
133 
135 
139  ComputeLaplacianSimplex(IndexType NewId, GeometryType::Pointer pGeometry) :
140  BaseType(NewId, pGeometry)
141  {}
142 
144 
149  ComputeLaplacianSimplex(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) :
150  BaseType(NewId, pGeometry, pProperties)
151  {}
152 
155  {}
156 
157 
161 
162 
166 
168 
175  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes,
176  PropertiesType::Pointer pProperties) const override
177  {
178  return Element::Pointer(new ComputeLaplacianSimplex(NewId, this->GetGeometry().Create(ThisNodes), pProperties));
179  }
180 
182 
188  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
189 
191 
195  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
196 
198 
203 // virtual void CalculateOnIntegrationPoints(const Variable<array_1d<double, 3 > >& rVariable,
204 // std::vector<array_1d<double, 3 > >& rValues,
205 // const ProcessInfo& rCurrentProcessInfo)
206 // {
207 //
208 // }
209 
210 
214 
218 
220 
228  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
229 
230 
234 
235 
239 
241  virtual std::string Info() const override
242  {
243  std::stringstream buffer;
244  buffer << "ComputeLaplacianSimplex #" << this->Id();
245  return buffer.str();
246  }
247 
249  virtual void PrintInfo(std::ostream& rOStream) const override
250  {
251  rOStream << "ComputeLaplacianSimplex" << TDim << "D";
252  }
253 
254 // /// Print object's data.
255 // virtual void PrintData(std::ostream& rOStream) const;
256 
260 
261 
263 
264 protected:
267 
268 
272 
273 
277 
278 
282 
283 
287 
288 
292 
293 
297 
298 
300 
301 private:
304 
308 
312 
313  friend class Serializer;
314 
318 
319  void AddIntegrationPointRHSContribution(VectorType& F, const BoundedMatrix<double, TNumNodes, TDim>& rShapeDeriv, const double Weight);
320 
324 
325 
329 
330 
334 
335 
339 
342 
345 
347 
348 }; // Class ComputeLaplacianSimplex
349 
351 
354 
355 
359 
360 
362 template< unsigned int TDim >
363 inline std::istream& operator >>(std::istream& rIStream,
365 {
366  return rIStream;
367 }
368 
370 template< unsigned int TDim >
371 inline std::ostream& operator <<(std::ostream& rOStream,
372  const ComputeLaplacianSimplex<TDim>& rThis)
373 {
374  rThis.PrintInfo(rOStream);
375  rOStream << std::endl;
376  rThis.PrintData(rOStream);
377 
378  return rOStream;
379 }
381 
383 
384 } // namespace Kratos.
385 
386 #endif // KRATOS_COMPUTE_LAPLACIAN_SIMPLEX_ELEMENT_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
A post-processing element to recover the Laplacian from the velocity solution.
Definition: calculate_laplacian_simplex_element.h:63
std::size_t IndexType
Definition: calculate_laplacian_simplex_element.h:90
ComputeLaplacianSimplex(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: calculate_laplacian_simplex_element.h:130
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(ComputeLaplacianSimplex)
Pointer definition of ComputeLaplacianSimplex.
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new element of this type.
Definition: calculate_laplacian_simplex_element.h:175
virtual std::string Info() const override
Turn back information as a string.
Definition: calculate_laplacian_simplex_element.h:241
PointerVectorSet< Dof< double >, IndexedObject > DofsArrayType
Definition: calculate_laplacian_simplex_element.h:98
Matrix MatrixType
Matrix type for local contributions to the linear system.
Definition: calculate_laplacian_simplex_element.h:86
ComputeLaplacianSimplex(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constuctor using geometry and properties.
Definition: calculate_laplacian_simplex_element.h:149
Kratos::Vector ShapeFunctionsType
Type for shape function values container.
Definition: calculate_laplacian_simplex_element.h:101
Vector VectorType
Vector type for local contributions to the linear system.
Definition: calculate_laplacian_simplex_element.h:83
std::size_t SizeType
Definition: calculate_laplacian_simplex_element.h:92
virtual void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: calculate_laplacian_simplex_element.h:249
ComputeMaterialDerivativeSimplex< TDim, TNumNodes > BaseType
Definition: calculate_laplacian_simplex_element.h:71
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition of nodes container type, redefined from GeometryType.
Definition: calculate_laplacian_simplex_element.h:80
ComputeLaplacianSimplex(IndexType NewId=0)
Default constuctor.
Definition: calculate_laplacian_simplex_element.h:121
Kratos::Matrix ShapeFunctionDerivativesType
Type for a matrix containing the shape function gradients.
Definition: calculate_laplacian_simplex_element.h:104
Geometry< NodeType > GeometryType
Geometry type (using with given NodeType)
Definition: calculate_laplacian_simplex_element.h:77
std::vector< std::size_t > EquationIdVectorType
Definition: calculate_laplacian_simplex_element.h:94
ComputeLaplacianSimplex(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using a geometry object.
Definition: calculate_laplacian_simplex_element.h:139
Node NodeType
Node type (default is: Node)
Definition: calculate_laplacian_simplex_element.h:74
virtual ~ComputeLaplacianSimplex()
Destructor.
Definition: calculate_laplacian_simplex_element.h:154
std::vector< Dof< double >::Pointer > DofsVectorType
Definition: calculate_laplacian_simplex_element.h:96
GeometryType::ShapeFunctionsGradientsType ShapeFunctionDerivativesArrayType
Type for an array of shape function gradient matrices.
Definition: calculate_laplacian_simplex_element.h:107
Properties PropertiesType
Definition: calculate_laplacian_simplex_element.h:88
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
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