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.
embedded_skin_visualization_process.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: Ruben Zorrilla
11 //
12 //
13 
14 #ifndef KRATOS_EMBEDDED_SKIN_VISUALIZATION_PROCESS_H
15 #define KRATOS_EMBEDDED_SKIN_VISUALIZATION_PROCESS_H
16 
17 // System includes
18 #include <string>
19 #include <iostream>
20 #include <unordered_map>
21 
22 // External includes
23 #include <boost/functional/hash.hpp> //TODO: remove this dependence when Kratos has en internal one
24 #include <boost/unordered_map.hpp> //TODO: remove this dependence when Kratos has en internal one
25 
26 // Project includes
27 #include "includes/define.h"
28 #include "includes/key_hash.h"
30 #include "processes/process.h"
33 
34 // Application includes
35 
36 
37 namespace Kratos
38 {
41 
44 
48 
52 
56 
57 
61 
75 class KRATOS_API(FLUID_DYNAMICS_APPLICATION) EmbeddedSkinVisualizationProcess : public Process
76 {
77 public:
80 
87  enum class LevelSetType
88  {
89  Continuous,
90  Discontinuous
91  };
92 
99  enum class ShapeFunctionsType
100  {
101  Ausas,
102  Standard
103  };
104 
105  struct Hash{
106  std::size_t operator()(const std::pair<unsigned int,bool>& k) const{
107  std::size_t h1 = std::hash<unsigned int>()(std::get<0>(k));
108  std::size_t h2 = std::hash<bool>()(std::get<1>(k));
109  return h1 ^ (h2 << 1);
110  }
111  };
112 
113  struct KeyEqual{
114  bool operator()(const std::pair<unsigned int,bool>& lhs, const std::pair<unsigned int,bool>& rhs) const{
115  return ((std::get<0>(lhs) == std::get<0>(rhs)) && (std::get<1>(lhs) == std::get<1>(rhs)));
116  }
117  };
118 
121 
122  typedef std::unordered_map<
123  Node::Pointer,
124  std::tuple< const Node::Pointer, const Node::Pointer, const double, const double >,
127 
131 
139  static ModelPart& CreateAndPrepareVisualizationModelPart(
140  Model& rModel,
141  const Parameters rParameters
142  );
143 
150  static void CheckAndSetLevelSetType(
151  const Parameters rParameters,
152  LevelSetType& rLevelSetType);
153 
160  static void CheckAndSetShapeFunctionsType(
161  const Parameters rParameters,
162  ShapeFunctionsType& rShapeFunctionsType);
163 
171  static const std::string CheckAndReturnDistanceVariableName(
172  const Parameters rParameters,
173  const LevelSetType& rLevelSetType);
174 
182  template<class TDataType>
183  static void FillVariablesList(
184  const Parameters rParameters,
185  std::vector<const Variable<TDataType>*>& rVariablesList);
186 
188 
200  ModelPart& rModelPart,
201  ModelPart& rVisualizationModelPart,
202  const std::vector<const Variable< double>* >& rVisualizationScalarVariables,
203  const std::vector<const Variable< array_1d<double, 3> >* >& rVisualizationVectorVariables,
204  const std::vector<const Variable< double>* >& rVisualizationNonHistoricalScalarVariables,
205  const std::vector<const Variable< array_1d<double, 3> >* >& rVisualizationNonHistoricalVectorVariables,
206  const LevelSetType& rLevelSetType,
207  const ShapeFunctionsType& rShapeFunctionsType,
208  const bool ReformModelPartAtEachTimeStep = false);
209 
217  ModelPart& rModelPart,
218  ModelPart& rVisualizationModelPart,
219  Parameters rParameters);
220 
227  Model& rModel,
228  Parameters rParameters);
229 
232 
236 
240 
241  void ExecuteBeforeSolutionLoop() override;
242 
243  void ExecuteBeforeOutputStep() override;
244 
245  void ExecuteAfterOutputStep() override;
246 
247  int Check() override;
248 
254  static Parameters StaticGetDefaultParameters();
255 
259  const Parameters GetDefaultParameters() const override
260  {
261  return StaticGetDefaultParameters();
262  }
263 
267 
268 
272 
276 
278  std::string Info() const override
279  {
280  std::stringstream buffer;
281  buffer << "EmbeddedSkinVisualizationProcess" ;
282  return buffer.str();
283  }
284 
286  void PrintInfo(std::ostream& rOStream) const override {rOStream << "EmbeddedSkinVisualizationProcess";}
287 
289  void PrintData(std::ostream& rOStream) const override {}
290 
294 
296 
297 private:
300 
301 
305 
306  // Unordered map to relate the newly created nodes and the origin mesh ones
307  CutNodesMapType mCutNodesMap;
308 
309  // Container with the splitting newly created elements
310  ModelPart::ElementsContainerType mNewElementsPointers;
311 
312  // Reference to the origin model part
313  ModelPart& mrModelPart;
314 
315  // Reference to the visualization model part
316  ModelPart& mrVisualizationModelPart;
317 
318  // Level set type. Current available options continuous and discontinuous
319  const LevelSetType mLevelSetType;
320 
321  // Shape functions type. Current available options ausas and standard
322  const ShapeFunctionsType mShapeFunctionsType;
323 
324  // If true, the visualization model part is created each time step (required in case the level set function is not constant)
325  const bool mReformModelPartAtEachTimeStep;
326 
327  // Pointer to the variable that stores the nodal level set function
328  const Variable<double>* mpNodalDistanceVariable;
329 
330  // Pointer to the variable that stores the elemental level set function
331  const Variable<Vector>* mpElementalDistanceVariable;
332 
333  // Vector containing the scalar variables to be interpolated in the visualization mesh
334  const std::vector<const Variable<double>*> mVisualizationScalarVariables;
335 
336  // Vector containing the vector variables to be interpolated in the visualization mesh
337  const std::vector<const Variable<array_1d<double, 3>>*> mVisualizationVectorVariables;
338 
339  // Vector containing the non-historical scalar variables to be interpolated in the visualization mesh
340  const std::vector<const Variable<double>*> mVisualizationNonHistoricalScalarVariables;
341 
342  // Vector containing the non-historical vector variables to be interpolated in the visualization mesh
343  const std::vector<const Variable<array_1d<double, 3>>*> mVisualizationNonHistoricalVectorVariables;
344 
348 
349 
353 
363  template<bool IsHistorical>
364  double& AuxiliaryGetValue(
365  Node& rNode,
366  const Variable<double>& rVariable);
367 
377  template<bool IsHistorical>
378  array_1d<double,3>& AuxiliaryGetValue(
379  Node& rNode,
380  const Variable<array_1d<double,3>>& rVariable);
381 
386  void CreateVisualizationMesh();
387 
391  void ComputeNewNodesInterpolation();
392 
405  template<class TDataType, bool IsHistorical>
406  void InterpolateVariablesListValues(
407  const Node::Pointer& rpNode,
408  const Node::Pointer& rpNodeI,
409  const Node::Pointer& rpNodeJ,
410  const double WeightI,
411  const double WeightJ,
412  const std::vector<const Variable<TDataType>*>& rVariablesList);
413 
420  template<const bool IsDistributed>
421  void CopyOriginNodes();
422 
426  void CopyOriginNodalValues();
427 
438  template<class TDataType, bool IsHistorical>
439  void CopyVariablesListValues(
440  const ModelPart::NodeIterator& rItOriginNode,
441  ModelPart::NodeIterator& rItVisualizationNode,
442  const std::vector<const Variable<TDataType>*>& rVariablesList);
443 
447  void CreateVisualizationGeometries();
448 
455  template<class TDataType>
456  void InitializeNonHistoricalVariables(const std::vector<const Variable<TDataType>*>& rNonHistoricalVariablesVector);
457 
464  bool ElementIsSplit(
465  const Geometry<Node>::Pointer pGeometry,
466  const Vector &rNodalDistances);
467 
474  bool ElementIsIncised(const Vector &rEdgeDistancesExtrapolated);
475 
482  bool ElementIsPositive(
483  Geometry<Node>::Pointer pGeometry,
484  const Vector &rNodalDistances);
485 
492  const Vector SetDistancesVector(ModelPart::ElementIterator ItElem);
493 
499  const inline Vector SetEdgeDistancesExtrapolatedVector(const Element& rElem);
500 
508  ModifiedShapeFunctions::Pointer SetModifiedShapeFunctionsUtility(
509  const Geometry<Node>::Pointer pGeometry,
510  const Vector &rNodalDistances);
511 
519  ModifiedShapeFunctions::Pointer SetAusasIncisedModifiedShapeFunctionsUtility(
520  const Geometry<Node>::Pointer pGeometry,
521  const Vector &rNodalDistancesWithExtra,
522  const Vector &rEdgeDistancesExtrapolated);
523 
530  Geometry< Node >::Pointer SetNewConditionGeometry(
531  const GeometryData::KratosGeometryType &rOriginGeometryType,
532  const Condition::NodesArrayType &rNewNodesArray);
533 
540  std::tuple< Properties::Pointer , Properties::Pointer > SetVisualizationProperties();
541 
546  void RemoveVisualizationProperties();
547 
548  template<const bool IsDistributed>
549  static void SetPartitionIndexFromOriginNode(
550  const Node& rOriginNode,
551  Node& rVisualizationNode);
552 
553  template<const bool IsDistributed>
554  static void SetPartitionIndex(
555  const int PartitionIndex,
556  Node& rVisualizationNode);
557 
561 
565 
569 
572 
575 
578 
580 
581 }; // Class EmbeddedSkinVisualizationProcess
582 
586 
590 
592 
594 
595 }; // namespace Kratos.
596 
597 #endif // KRATOS_EMBEDDED_SKIN_VISUALIZATION_PROCESS_H
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Elements.
Definition: element.h:60
This process saves the intersected elements in a different model part for its visualization.
Definition: embedded_skin_visualization_process.h:76
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: embedded_skin_visualization_process.h:286
KRATOS_CLASS_POINTER_DEFINITION(EmbeddedSkinVisualizationProcess)
Pointer definition of EmbeddedSkinVisualizationProcess.
const Parameters GetDefaultParameters() const override
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: embedded_skin_visualization_process.h:259
std::unordered_map< Node::Pointer, std::tuple< const Node::Pointer, const Node::Pointer, const double, const double >, SharedPointerHasher< Node::Pointer >, SharedPointerComparator< Node::Pointer > > CutNodesMapType
Definition: embedded_skin_visualization_process.h:126
LevelSetType
Enum class with the available level set types Auxiliary enum class to store the available level set t...
Definition: embedded_skin_visualization_process.h:88
std::string Info() const override
Turn back information as a string.
Definition: embedded_skin_visualization_process.h:278
ShapeFunctionsType
Enum class with the available shape functions type Auxiliary enum class to store the available shape ...
Definition: embedded_skin_visualization_process.h:100
~EmbeddedSkinVisualizationProcess() override
Destructor.
Definition: embedded_skin_visualization_process.h:231
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: embedded_skin_visualization_process.h:289
KratosGeometryType
Definition: geometry_data.h:110
Geometry base class.
Definition: geometry.h:71
This class aims to manage different model parts across multi-physics simulations.
Definition: model.h:60
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
MeshType::ElementIterator ElementIterator
Definition: model_part.h:174
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
The base class for all processes in Kratos.
Definition: process.h:49
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
rhs
Definition: generate_frictional_mortar_condition.py:297
lhs
Definition: generate_frictional_mortar_condition.py:297
int k
Definition: quadrature.py:595
Definition: embedded_skin_visualization_process.h:105
std::size_t operator()(const std::pair< unsigned int, bool > &k) const
Definition: embedded_skin_visualization_process.h:106
Definition: embedded_skin_visualization_process.h:113
bool operator()(const std::pair< unsigned int, bool > &lhs, const std::pair< unsigned int, bool > &rhs) const
Definition: embedded_skin_visualization_process.h:114
This is a key comparer between two shared pointers.
Definition: key_hash.h:312
This is a hasher for shared pointers.
Definition: key_hash.h:294