14 #ifndef KRATOS_EMBEDDED_SKIN_VISUALIZATION_PROCESS_H
15 #define KRATOS_EMBEDDED_SKIN_VISUALIZATION_PROCESS_H
20 #include <unordered_map>
23 #include <boost/functional/hash.hpp>
24 #include <boost/unordered_map.hpp>
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);
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)));
122 typedef std::unordered_map<
124 std::tuple< const Node::Pointer, const Node::Pointer, const double, const double >,
139 static ModelPart& CreateAndPrepareVisualizationModelPart(
150 static void CheckAndSetLevelSetType(
160 static void CheckAndSetShapeFunctionsType(
171 static const std::string CheckAndReturnDistanceVariableName(
182 template<
class TDataType>
183 static void FillVariablesList(
204 const std::vector<
const Variable< double>* >& rVisualizationNonHistoricalScalarVariables,
208 const bool ReformModelPartAtEachTimeStep =
false);
241 void ExecuteBeforeSolutionLoop()
override;
243 void ExecuteBeforeOutputStep()
override;
245 void ExecuteAfterOutputStep()
override;
247 int Check()
override;
254 static Parameters StaticGetDefaultParameters();
261 return StaticGetDefaultParameters();
278 std::string
Info()
const override
280 std::stringstream buffer;
281 buffer <<
"EmbeddedSkinVisualizationProcess" ;
286 void PrintInfo(std::ostream& rOStream)
const override {rOStream <<
"EmbeddedSkinVisualizationProcess";}
289 void PrintData(std::ostream& rOStream)
const override {}
307 CutNodesMapType mCutNodesMap;
319 const LevelSetType mLevelSetType;
322 const ShapeFunctionsType mShapeFunctionsType;
325 const bool mReformModelPartAtEachTimeStep;
334 const std::vector<const Variable<double>*> mVisualizationScalarVariables;
337 const std::vector<const Variable<array_1d<double, 3>>*> mVisualizationVectorVariables;
340 const std::vector<const Variable<double>*> mVisualizationNonHistoricalScalarVariables;
343 const std::vector<const Variable<array_1d<double, 3>>*> mVisualizationNonHistoricalVectorVariables;
363 template<
bool IsHistorical>
364 double& AuxiliaryGetValue(
377 template<
bool IsHistorical>
386 void CreateVisualizationMesh();
391 void ComputeNewNodesInterpolation();
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,
420 template<const
bool IsDistributed>
421 void CopyOriginNodes();
426 void CopyOriginNodalValues();
438 template<
class TDataType,
bool IsHistorical>
439 void CopyVariablesListValues(
447 void CreateVisualizationGeometries();
455 template<
class TDataType>
456 void InitializeNonHistoricalVariables(
const std::vector<
const Variable<TDataType>*>& rNonHistoricalVariablesVector);
466 const Vector &rNodalDistances);
474 bool ElementIsIncised(
const Vector &rEdgeDistancesExtrapolated);
482 bool ElementIsPositive(
484 const Vector &rNodalDistances);
499 const inline Vector SetEdgeDistancesExtrapolatedVector(
const Element& rElem);
508 ModifiedShapeFunctions::Pointer SetModifiedShapeFunctionsUtility(
510 const Vector &rNodalDistances);
519 ModifiedShapeFunctions::Pointer SetAusasIncisedModifiedShapeFunctionsUtility(
521 const Vector &rNodalDistancesWithExtra,
522 const Vector &rEdgeDistancesExtrapolated);
540 std::tuple< Properties::Pointer , Properties::Pointer > SetVisualizationProperties();
546 void RemoveVisualizationProperties();
548 template<const
bool IsDistributed>
549 static void SetPartitionIndexFromOriginNode(
550 const Node& rOriginNode,
551 Node& rVisualizationNode);
553 template<const
bool IsDistributed>
554 static void SetPartitionIndex(
555 const int PartitionIndex,
556 Node& rVisualizationNode);
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