13 #if !defined(KRATOS_GENERATE_EMBEDDED_SKIN_UTILITY_H_INCLUDED )
14 #define KRATOS_GENERATE_EMBEDDED_SKIN_UTILITY_H_INCLUDED
64 template<std::
size_t TDim>
74 typedef std::unordered_map<
76 std::tuple< const Element::Pointer, const unsigned int >,
98 const std::string LevelSetType =
"continuous",
99 const std::vector<std::string>& InterpolatedSkinVariables = {}) :
100 mrModelPart(rModelPart),
101 mrSkinModelPart(rSkinModelPart),
102 mLevelSetType(LevelSetType ==
"continuous" ? Continuous : Discontinuous),
104 mInterpolatedSkinVariables(InterpolatedSkinVariables) {};
132 void InterpolateMeshVariableToSkin(
142 void InterpolateMeshVariableToSkin(
154 void InterpolateDiscontinuousMeshVariableToSkin(
157 const std::string &rInterfaceSide);
167 void InterpolateDiscontinuousMeshVariableToSkin(
170 const std::string &rInterfaceSide);
207 std::vector<std::string> mInterpolatedSkinVariables;
238 void ComputeElementSkin(
239 const Element::Pointer pElement,
240 const Vector &rNodalDistances,
241 unsigned int &rTempNodeId,
242 unsigned int &rTempCondId,
243 Properties::Pointer pCondProp,
255 bool inline ElementIsSplit(
257 const Vector &rNodalDistances);
270 template<
class TVarType>
271 void InterpolateMeshVariableToSkinSpecialization(
274 const std::string &rInterfaceSide =
"positive")
277 KRATOS_ERROR_IF((mrModelPart.NodesBegin())->SolutionStepsDataHas(rMeshVariable) ==
false)
278 <<
"Mesh model part solution step data missing variable: " << rMeshVariable << std::endl;
280 <<
"Generated skin model part solution step data missing variable: " << rSkinVariable << std::endl;
283 KRATOS_ERROR_IF(mrModelPart.NumberOfElements() == 0) <<
"Mesh model part has no elements.";
287 Element::Pointer p_elem;
288 #pragma omp parallel for private (i_edge, p_elem)
289 for (
int i_node = 0; i_node < static_cast<int>(mrSkinModelPart.
NumberOfNodes()); ++i_node) {
291 auto it_node = mrSkinModelPart.
NodesBegin() + i_node;
292 Node::Pointer p_node = &(*it_node);
295 const auto i_node_info = mEdgeNodesMap.find(p_node);
296 if (i_node_info != mEdgeNodesMap.end()){
298 std::tie(p_elem, i_edge) = std::get<1>(*i_node_info);
301 const auto p_elem_geom = p_elem->pGetGeometry();
302 const auto elem_dist = this->SetDistancesVector(*p_elem);
303 const auto p_mod_sh_func = pCreateModifiedShapeFunctions(p_elem_geom, elem_dist);
306 const auto edge_sh_func = this->GetModifiedShapeFunctionsValuesOnEdge(
311 const auto edge_N =
row(edge_sh_func, i_edge);
312 const auto &r_elem_geom = p_elem->GetGeometry();
313 auto &r_value = it_node->FastGetSolutionStepValue(rSkinVariable);
314 r_value = rSkinVariable.
Zero();
315 for (
unsigned int i_elem_node = 0; i_elem_node < r_elem_geom.PointsNumber(); ++i_elem_node) {
316 r_value += edge_N[i_elem_node] * r_elem_geom[i_elem_node].FastGetSolutionStepValue(rMeshVariable);
319 KRATOS_ERROR <<
"Intersected edge node " << it_node->Id() <<
" not found in intersected edges nodes map" << std::endl;
331 void RenumberAndAddSkinEntities(
342 const Vector SetDistancesVector(
const Element &rElement);
350 DivideGeometry<Node>::Pointer SetDivideGeometryUtility(
351 const Geometry<Node> &rGeometry,
352 const Vector &rNodalDistances);
360 Geometry< Node >::Pointer pCreateNewConditionGeometry(
375 Condition::Pointer pCreateNewCondition(
378 const unsigned int &rConditionId,
379 const Properties::Pointer pConditionProperties);
387 static const std::string GetConditionType();
395 Properties::Pointer SetSkinEntitiesProperties();
405 ModifiedShapeFunctions::UniquePointer pCreateModifiedShapeFunctions(
406 const Geometry<Node>::Pointer pGeometry,
407 const Vector& rNodalDistances);
417 Matrix GetModifiedShapeFunctionsValues(
418 const ModifiedShapeFunctions::UniquePointer &rpModifiedShapeFunctions,
419 const std::string &rInterfaceSide)
const;
429 Matrix GetModifiedShapeFunctionsValuesOnEdge(
430 const ModifiedShapeFunctions::UniquePointer &rpModifiedShapeFunctions,
431 const std::string &rInterfaceSide)
const;
448 EmbeddedSkinUtility&
operator=(EmbeddedSkinUtility
const& rOther) =
delete;
451 EmbeddedSkinUtility(EmbeddedSkinUtility
const& rOther) =
delete;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Conditions.
Definition: condition.h:59
Geometry< NodeType >::PointsArrayType NodesArrayType
definition of nodes container type, redefined from GeometryType
Definition: condition.h:86
Utility to compute the skin representation from a distance function.
Definition: embedded_skin_utility.h:66
EmbeddedSkinUtility(ModelPart &rModelPart, ModelPart &rSkinModelPart, const std::string LevelSetType="continuous", const std::vector< std::string > &InterpolatedSkinVariables={})
Default constructor.
Definition: embedded_skin_utility.h:95
KRATOS_CLASS_POINTER_DEFINITION(EmbeddedSkinUtility)
Pointer definition of EmbeddedSkinUtility.
std::unordered_map< Node::Pointer, std::tuple< const Element::Pointer, const unsigned int >, SharedPointerHasher< Node::Pointer >, SharedPointerComparator< Node::Pointer > > EdgeNodesMapType
Definition: embedded_skin_utility.h:78
LevelSetTypeEnum
Definition: embedded_skin_utility.h:85
virtual ~EmbeddedSkinUtility()=default
Destructor.
KratosGeometryType
Definition: geometry_data.h:110
Geometry base class.
Definition: geometry.h:71
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
SizeType NumberOfNodes(IndexType ThisIndex=0) const
Definition: model_part.h:341
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
const TDataType & Zero() const
This method returns the zero value of the variable type.
Definition: variable.h:346
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
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
AMatrix::MatrixRow< const TExpressionType > row(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression, std::size_t RowIndex)
Definition: amatrix_interface.h:649
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