10 #if !defined(KRATOS_BUILD_STRING_SKIN_PROCESS_H_INCLUDED)
11 #define KRATOS_BUILD_STRING_SKIN_PROCESS_H_INCLUDED
72 ) :
Process() , mrModelPart(rModelPart), mSides(sides), mRadius(
radius)
114 CreateSkinElements();
116 TransferSkinToOutput();
121 (
i)->
Set(RIGID,
true);
122 (
i)->
Set(ACTIVE,
true);
128 (
i)->
Set(RIGID,
true);
129 (
i)->
Set(ACTIVE,
true);
148 SetInActiveFlag(mrModelPart);
173 SetActiveFlag(mrModelPart);
211 std::string
Info()
const override
213 return "BuildStringSkinProcess";
219 rOStream <<
"BuildStringSkinProcess";
290 void CreateGeneratrix()
301 Node::Pointer Starter;
302 for(
auto i_node(mrModelPart.NodesBegin()); i_node != mrModelPart.NodesEnd(); ++i_node)
305 if( nNodes.
size() <= 1 )
306 Starter = *i_node.base();
310 unsigned int previous_id = 0;
312 Element::Pointer CurrentElement;
313 for(
unsigned int i=0;
i<mrModelPart.NumberOfNodes();
i++)
320 mGeneratrixNodes.push_back( Starter );
324 for(
auto& i_nelem : nElements)
328 bool selected =
true;
329 for(
unsigned int j = 0;
j < rGeometry.
size();
j++)
331 if( rGeometry[
j].Id() == previous_id )
339 for(
unsigned int j = 0;
j < rGeometry.
size();
j++)
343 previous_id = Starter->
Id();
344 Starter = rGeometry(
j);
377 SetActiveFlag(mrModelPart,TO_ERASE);
379 mrModelPart.RemoveNodes();
380 mrModelPart.RemoveElements();
382 std::cout<<
" [String_builder] [Defined by "<<mGeneratrixNodes.size()<<
" control points]"<<std::endl;
390 void CreateSkinNodes()
406 ModelPart::NodesContainerType::iterator nodes_begin = mGeneratrixNodes.begin();
409 for(
unsigned int i=0;
i<mGeneratrixNodes.size();
i++)
412 Point0[0] = (nodes_begin+
i)->X0();
413 Point0[1] = (nodes_begin+
i)->Y0();
414 Point0[2] = (nodes_begin+
i)->Z0();
416 Point[0] = (nodes_begin+
i)->
X();
417 Point[1] = (nodes_begin+
i)->
Y();
418 Point[2] = (nodes_begin+
i)->
Z();
425 array_1d<double,3>& NodeRotation = (nodes_begin+
i)->FastGetSolutionStepValue( ROTATION );
429 for(
unsigned int j=0;
j<3;
j++ )
431 SectionRotation[
j] = NodeRotation[
j];
445 double RadiusCorrection = 1;
447 if(
i == mGeneratrixNodes.size()-1 ){
448 BasePoint[0] = (nodes_begin+(
i-1))->X0();
449 BasePoint[1] = (nodes_begin+(
i-1))->Y0();
450 BasePoint[2] = (nodes_begin+(
i-1))->Z0();
451 DirectionZ = Point0 - BasePoint;
454 BasePoint[0] = (nodes_begin+(
i+1))->X0();
455 BasePoint[1] = (nodes_begin+(
i+1))->Y0();
456 BasePoint[2] = (nodes_begin+(
i+1))->Z0();
457 DirectionZ = BasePoint - Point0;
460 BasePoint[0] = (nodes_begin+(
i-1))->X0();
461 BasePoint[1] = (nodes_begin+(
i-1))->Y0();
462 BasePoint[2] = (nodes_begin+(
i-1))->Z0();
464 DirectionZ1 = (Point0 - BasePoint);
465 DirectionZ1 /=
norm_2(DirectionZ1);
467 BasePoint[0] = (nodes_begin+(
i+1))->X0();
468 BasePoint[1] = (nodes_begin+(
i+1))->Y0();
469 BasePoint[2] = (nodes_begin+(
i+1))->Z0();
471 DirectionZ2 = (BasePoint - Point0);
472 DirectionZ2 /=
norm_2(DirectionZ2);
474 DirectionZ = DirectionZ1 + DirectionZ2;
477 DirectionZ /=
norm_2(DirectionZ);
479 DirectionZ = DirectionZ1;
481 DirectionEllipse = DirectionZ1 - DirectionZ2;
483 if(
norm_2(DirectionEllipse))
484 DirectionEllipse/=
norm_2(DirectionEllipse);
486 RadiusCorrection =
inner_prod(DirectionZ,DirectionZ1);
487 RadiusCorrection +=
inner_prod(DirectionZ,DirectionZ2);
488 RadiusCorrection *= 0.5;
490 RadiusCorrection = 1.0/RadiusCorrection;
492 RadiusCorrection = 1.0;
498 if( DirectionZ[0] == 0 && DirectionZ[1] == 0 && DirectionZ[2] == 1){
500 DirectionX = DirectionY;
501 DirectionY = (-1) *
Temp;
504 double EllipsoidalCorrection = 1;
506 for(
unsigned int k=0;
k<mSides;
k++)
511 RotationAxis = DirectionZ *
alpha;
515 RotatedDirectionX = DirectionX;
517 RotationQuaternion.RotateVector3(RotatedDirectionX);
520 EllipsoidalCorrection =
inner_prod(DirectionEllipse,RotatedDirectionX);
521 EllipsoidalCorrection = 1 + ( RadiusCorrection - 1 ) * (EllipsoidalCorrection * EllipsoidalCorrection);
529 BasePoint = Radius * EllipsoidalCorrection * RotatedDirectionX;
531 SectionQuaternion.RotateVector3(BasePoint);
542 mrModelPart.AddNode(this->CreateNode(mrModelPart.GetParentModelPart(), BasePoint,
node_id));
562 int number_of_angles = mSides;
574 ModelPart::NodesContainerType::iterator nodes_begin = mGeneratrixNodes.begin();
576 ModelPart::NodesContainerType::iterator skin_nodes_begin = mrModelPart.NodesBegin();
582 for(
unsigned int i=0;
i<mGeneratrixNodes.size();
i++)
585 Point0[0] = (nodes_begin+
i)->X0();
586 Point0[1] = (nodes_begin+
i)->Y0();
587 Point0[2] = (nodes_begin+
i)->Z0();
589 Point[0] = (nodes_begin+
i)->
X();
590 Point[1] = (nodes_begin+
i)->
Y();
591 Point[2] = (nodes_begin+
i)->
Z();
599 PointType& NodeRotation = (nodes_begin+
i)->FastGetSolutionStepValue( ROTATION );
603 for(
unsigned int j=0;
j<3;
j++ )
605 SectionRotation[
j] = NodeRotation[
j];
618 double RadiusCorrection = 1;
620 if(
i == mGeneratrixNodes.size()-1 ){
621 BasePoint[0] = (nodes_begin+(
i-1))->X0();
622 BasePoint[1] = (nodes_begin+(
i-1))->Y0();
623 BasePoint[2] = (nodes_begin+(
i-1))->Z0();
624 DirectionZ = Point0 - BasePoint;
627 BasePoint[0] = (nodes_begin+(
i+1))->X0();
628 BasePoint[1] = (nodes_begin+(
i+1))->Y0();
629 BasePoint[2] = (nodes_begin+(
i+1))->Z0();
630 DirectionZ = BasePoint - Point0;
633 BasePoint[0] = (nodes_begin+(
i-1))->X0();
634 BasePoint[1] = (nodes_begin+(
i-1))->Y0();
635 BasePoint[2] = (nodes_begin+(
i-1))->Z0();
637 DirectionZ1 = (Point0 - BasePoint);
638 DirectionZ1 /=
norm_2(DirectionZ1);
640 BasePoint[0] = (nodes_begin+(
i+1))->X0();
641 BasePoint[1] = (nodes_begin+(
i+1))->Y0();
642 BasePoint[2] = (nodes_begin+(
i+1))->Z0();
644 DirectionZ2 = (BasePoint - Point0);
645 DirectionZ2 /=
norm_2(DirectionZ2);
647 DirectionZ = DirectionZ1 + DirectionZ2;
648 DirectionZ /=
norm_2(DirectionZ);
650 DirectionEllipse = DirectionZ1 - DirectionZ2;
652 if(
norm_2(DirectionEllipse))
653 DirectionEllipse/=
norm_2(DirectionEllipse);
655 RadiusCorrection =
inner_prod(DirectionZ,DirectionZ1);
656 RadiusCorrection +=
inner_prod(DirectionZ,DirectionZ2);
657 RadiusCorrection *= 0.5;
658 RadiusCorrection = 1.0/RadiusCorrection;
666 if( DirectionZ[0] == 0 && DirectionZ[1] == 0 && DirectionZ[2] == 1){
668 DirectionX = DirectionY;
669 DirectionY = (-1) *
Temp;
674 double EllipsoidalCorrection = 1;
676 for(
int k=0;
k<number_of_angles;
k++)
681 RotationAxis = DirectionZ *
alpha;
685 RotatedDirectionX = DirectionX;
687 RotationQuaternion.RotateVector3(RotatedDirectionX);
692 EllipsoidalCorrection =
inner_prod(DirectionEllipse,RotatedDirectionX);
693 EllipsoidalCorrection = 1 + ( RadiusCorrection - 1 ) * (EllipsoidalCorrection * EllipsoidalCorrection);
699 BasePoint = Radius * EllipsoidalCorrection * RotatedDirectionX;
701 SectionQuaternion.RotateVector3(BasePoint);
705 RadiusVector[0] = BasePoint[0];
706 RadiusVector[1] = BasePoint[1];
707 RadiusVector[2] = BasePoint[2];
711 PreviousPosition[0] = (skin_nodes_begin+
counter)->X0();
712 PreviousPosition[1] = (skin_nodes_begin+
counter)->Y0();
713 PreviousPosition[2] = (skin_nodes_begin+
counter)->Z0();
717 (skin_nodes_begin+
counter)->
X() = BasePoint[0];
718 (skin_nodes_begin+
counter)->
Y() = BasePoint[1];
719 (skin_nodes_begin+
counter)->
Z() = BasePoint[2];
722 PointType& Displacement = (skin_nodes_begin+
counter)->FastGetSolutionStepValue(DISPLACEMENT);
724 for(
int j=0;
j<3;
j++)
725 Displacement[
j] = BasePoint[
j]-PreviousPosition[
j];
727 bool DynamicVariables =
false;
729 if( DynamicVariables ){
731 array_1d<double, 3 >& Velocity = (nodes_begin+
i)->FastGetSolutionStepValue(VELOCITY);
732 array_1d<double, 3 >& Acceleration = (nodes_begin+
i)->FastGetSolutionStepValue(ACCELERATION);
733 array_1d<double, 3 >& AngularVelocity = (nodes_begin+
i)->FastGetSolutionStepValue(ANGULAR_VELOCITY);
734 array_1d<double, 3 >& AngularAcceleration = (nodes_begin+
i)->FastGetSolutionStepValue(ANGULAR_ACCELERATION);
736 Matrix SkewSymVariable(3,3);
746 Variable =
prod(SkewSymVariable,RadiusVector);
748 (skin_nodes_begin+
counter)->FastGetSolutionStepValue(VELOCITY) = Velocity + Variable;
757 AngularVariable =
prod(SkewSymVariable,Variable);
763 Variable =
prod(SkewSymVariable,RadiusVector);
765 (skin_nodes_begin+
counter)->FastGetSolutionStepValue(ACCELERATION) = Acceleration + Variable + AngularVariable;
784 void CreateSkinElements()
789 return this->CreateSkinQuadrilaterals();
799 void CreateSkinTriangles()
803 unsigned int number_of_angles = mSides;
805 unsigned int wall_nodes_number_id = mMaxId;
810 unsigned int number_of_elements = (mrModelPart.Nodes().back().Id() - wall_nodes_number_id) - (number_of_angles - 1);
813 GeometryType::Pointer pFace;
814 ConditionType::Pointer pSkinCondition;
819 Properties::Pointer pProperties = mrModelPart.GetParentModelPart().pGetProperties(0);
827 unsigned int condition_id = GetMaxConditionId(mrModelPart.GetParentModelPart());
831 std::vector<int> FaceNodesIds(3);
833 for(
unsigned int i=1;
i<number_of_elements;
i++)
837 if(
counter < number_of_angles ){
840 FaceNodesIds[0] = wall_nodes_number_id +
i ;
841 FaceNodesIds[2] = wall_nodes_number_id +
i + number_of_angles;
842 FaceNodesIds[1] = wall_nodes_number_id +
i + number_of_angles + 1;
851 for(
unsigned int j=0;
j<3;
j++)
852 FaceNodes1.push_back(mrModelPart.pGetNode(FaceNodesIds[
j]));
854 pFace = Kratos::make_shared<Triangle3DType>(FaceNodes1);
856 pSkinCondition = Kratos::make_intrusive<Condition>(condition_id, pFace, pProperties);
858 pSkinCondition->Set(ACTIVE,
false);
861 mrModelPart.AddCondition(pSkinCondition);
866 FaceNodesIds[0] = wall_nodes_number_id +
i ;
867 FaceNodesIds[2] = wall_nodes_number_id +
i + number_of_angles + 1;
868 FaceNodesIds[1] = wall_nodes_number_id +
i + 1;
878 for(
unsigned int j=0;
j<3;
j++)
879 FaceNodes2.push_back(mrModelPart.pGetNode(FaceNodesIds[
j]));
881 pFace = Kratos::make_shared<Triangle3DType>(FaceNodes2);
883 pSkinCondition = Kratos::make_intrusive<Condition>(condition_id, pFace, pProperties);
885 pSkinCondition->Set(ACTIVE,
false);
888 mrModelPart.AddCondition(pSkinCondition);
893 else if(
counter == number_of_angles ){
896 FaceNodesIds[0] = wall_nodes_number_id +
i ;
897 FaceNodesIds[1] = wall_nodes_number_id +
i + 1 - number_of_angles;
898 FaceNodesIds[2] = wall_nodes_number_id +
i + 1;
903 for(
unsigned int j=0;
j<3;
j++)
904 FaceNodes1.push_back(mrModelPart.pGetNode(FaceNodesIds[
j]));
906 pFace = Kratos::make_shared<Triangle3DType>(FaceNodes1);
908 pSkinCondition = Kratos::make_intrusive<Condition>(condition_id, pFace, pProperties);
910 pSkinCondition->Set(ACTIVE,
false);
913 mrModelPart.AddCondition(pSkinCondition);
918 FaceNodesIds[0] = wall_nodes_number_id +
i ;
919 FaceNodesIds[1] = wall_nodes_number_id +
i + 1;
920 FaceNodesIds[2] = wall_nodes_number_id +
i + number_of_angles;
925 for(
unsigned int j=0;
j<3;
j++)
926 FaceNodes2.push_back(mrModelPart.pGetNode(FaceNodesIds[
j]));
928 pFace = Kratos::make_shared<Triangle3DType>(FaceNodes2);
930 pSkinCondition = Kratos::make_intrusive<Condition>(condition_id, pFace, pProperties);
932 pSkinCondition->Set(ACTIVE,
false);
935 mrModelPart.AddCondition(pSkinCondition);
950 void CreateSkinQuadrilaterals()
954 unsigned int number_of_angles = mSides;
956 unsigned int wall_nodes_number_id = mMaxId;
961 unsigned int number_of_elements = (mrModelPart.Nodes().back().Id() - wall_nodes_number_id) - (number_of_angles - 1);
964 GeometryType::Pointer pFace;
965 ConditionType::Pointer pSkinCondition;
970 Properties::Pointer pProperties = mrModelPart.GetParentModelPart().pGetProperties(0);
978 unsigned int condition_id = GetMaxConditionId(mrModelPart.GetParentModelPart());
982 std::vector<int> FaceNodesIds(4);
985 for(
unsigned int i=1;
i<number_of_elements;
i++)
989 if(
counter < number_of_angles ){
991 FaceNodesIds[0] = wall_nodes_number_id +
i ;
992 FaceNodesIds[1] = wall_nodes_number_id +
i + number_of_angles;
993 FaceNodesIds[2] = wall_nodes_number_id +
i + number_of_angles + 1;
994 FaceNodesIds[3] = wall_nodes_number_id +
i + 1;
1003 for(
unsigned int j=0;
j<4;
j++)
1004 FaceNodes.push_back(mrModelPart.pGetNode(FaceNodesIds[
j]));
1006 pFace = Kratos::make_shared<Quadrilateral3DType>(FaceNodes);
1008 pSkinCondition = Kratos::make_intrusive<Condition>(condition_id, pFace, pProperties);
1010 pSkinCondition->Set(ACTIVE,
false);
1013 mrModelPart.AddCondition(pSkinCondition);
1018 else if(
counter == number_of_angles ){
1020 FaceNodesIds[0] = wall_nodes_number_id +
i ;
1021 FaceNodesIds[1] = wall_nodes_number_id +
i + 1 - number_of_angles;
1022 FaceNodesIds[2] = wall_nodes_number_id +
i + 1;
1023 FaceNodesIds[3] = wall_nodes_number_id +
i + number_of_angles;
1028 for(
unsigned int j=0;
j<4;
j++)
1029 FaceNodes.push_back(mrModelPart.pGetNode(FaceNodesIds[
j]));
1031 pFace = Kratos::make_shared<Quadrilateral3DType>(FaceNodes);
1033 pSkinCondition = Kratos::make_intrusive<Condition>(condition_id, pFace, pProperties);
1035 pSkinCondition->Set(ACTIVE,
false);
1038 mrModelPart.AddCondition(pSkinCondition);
1052 static inline unsigned int GetMaxNodeId(ModelPart& rModelPart)
1056 unsigned int max_id = rModelPart.Nodes().back().Id();
1058 for(ModelPart::NodesContainerType::iterator i_node = rModelPart.NodesBegin(); i_node!= rModelPart.NodesEnd(); i_node++)
1060 if(i_node->Id() > max_id)
1061 max_id = i_node->Id();
1072 static inline unsigned int GetMaxElementId(ModelPart& rModelPart)
1076 if( rModelPart.NumberOfElements() == 0 )
1079 unsigned int max_id = rModelPart.Elements().back().Id();
1081 for(ModelPart::ElementsContainerType::iterator i_elem = rModelPart.ElementsBegin(); i_elem!= rModelPart.ElementsEnd(); i_elem++)
1083 if(i_elem->Id() > max_id)
1084 max_id = i_elem->Id();
1096 static inline unsigned int GetMaxConditionId(ModelPart& rModelPart)
1100 if( rModelPart.NumberOfConditions() == 0 )
1103 unsigned int max_id = rModelPart.Conditions().back().Id();
1105 for(ModelPart::ConditionsContainerType::iterator i_cond = rModelPart.ConditionsBegin(); i_cond!= rModelPart.ConditionsEnd(); i_cond++)
1107 if(i_cond->Id() > max_id)
1108 max_id = i_cond->Id();
1119 NodeType::Pointer CreateNode (ModelPart& rModelPart,
PointType& rPoint,
const unsigned int& rNodeId)
1124 NodeType::Pointer
Node = rModelPart.CreateNewNode( rNodeId, rPoint[0], rPoint[1], rPoint[2]);
1130 for(NodeType::DofsContainerType::iterator iii = reference_dofs.begin(); iii != reference_dofs.end(); iii++)
1133 Node->pAddDof( rDof );
1139 for(NodeType::DofsContainerType::iterator iii = new_dofs.begin(); iii != new_dofs.end(); iii++)
1170 void CleanElementNeighbours()
1181 unsigned int AverageNodes = 2;
1182 unsigned int AverageElements = 2;
1185 for(
auto& i_node : rNodes)
1189 nNodes.resize(AverageNodes);
1193 nElements.resize(AverageElements);
1197 for(
auto& i_elem : rElements)
1204 nElements.resize(size);
1213 template<
class TDataType>
void AddUniquePointer
1214 (GlobalPointersVector<TDataType>&
v,
const typename TDataType::WeakPointer candidate)
1218 while (
i != endit && (
i)->Id() != (candidate)->Id())
1224 v.push_back(candidate);
1235 for(
auto i_nelem(nElements.begin()); i_nelem != nElements.end(); ++i_nelem)
1238 Geometry<Node >& nGeometry = i_nelem->GetGeometry();
1239 if(nGeometry.LocalSpaceDimension() == 1){
1240 for(
unsigned int node_i = 0; node_i < nGeometry.size(); ++node_i)
1242 if(nGeometry[node_i].Id() == Id_1)
1244 if(i_nelem->Id() != i_elem->Id())
1246 return *i_nelem.base();
1252 return *i_elem.base();
1258 void SearchNeighbours()
1269 CleanElementNeighbours();
1274 for(
auto i_elem(rElements.begin()); i_elem != rElements.end(); ++i_elem)
1277 for(
unsigned int i = 0;
i < rGeometry.size();
i++)
1279 rGeometry[
i].GetValue(NEIGHBOUR_ELEMENTS).push_back( *i_elem.base() );
1284 for(
auto& i_node : rNodes)
1287 for(
auto& i_nelem : nElements)
1290 for(
unsigned int i = 0;
i < rGeometry.size();
i++)
1292 if( rGeometry[
i].Id() != i_node.Id() )
1295 AddUniquePointer< Node >(nNodes, rGeometry(
i));
1307 for(
auto i_elem(rElements.begin()); i_elem != rElements.end(); ++i_elem)
1310 Geometry<Node >& rGeometry = i_elem->GetGeometry();
1311 if( rGeometry.FacesNumber() == 2 ){
1316 if( nElements.size() != 2 )
1321 unsigned int size = rGeometry.size();
1324 nElements(0) = CheckForNeighbourElems1D(rGeometry[0].Id(), rGeometry[0].
GetValue(NEIGHBOUR_ELEMENTS), i_elem);
1326 nElements(1) = CheckForNeighbourElems1D(rGeometry[size-1].Id(), rGeometry[size-1].
GetValue(NEIGHBOUR_ELEMENTS), i_elem);
1336 void TransferSkinToOutput()
1340 ModelPart& rOutputModelPart = GetOutputModelPart();
1342 std::vector<std::size_t> NodeIds;
1346 NodeIds.push_back(
i->Id());
1348 std::vector<std::size_t> ConditionIds;
1352 ConditionIds.push_back(
i->Id());
1355 rOutputModelPart.AddNodes(NodeIds);
1356 rOutputModelPart.AddConditions(ConditionIds);
1358 SetInActiveFlag(rOutputModelPart,TO_ERASE);
1367 std::string OutputModelPartName;
1368 ModelPart& rMainModelPart = mrModelPart.GetParentModelPart();
1371 if( i_mp->Is(ACTIVE) )
1372 OutputModelPartName = i_mp->Name();
1375 return (rMainModelPart.GetSubModelPart(OutputModelPartName));
1381 virtual void SetActiveFlag(ModelPart& rModelPart,
Flags id_flag = ACTIVE)
1388 (
i)->
Set(id_flag,
true);
1394 (
i)->
Set(id_flag,
true);
1401 virtual void SetInActiveFlag(ModelPart& rModelPart,
Flags id_flag = ACTIVE)
1409 (
i)->
Set(id_flag,
false);
1415 (
i)->
Set(id_flag,
false);
1463 rOStream << std::endl;
Definition: beam_math_utilities.hpp:31
static void CalculateLocalAxesVectors(TVector3 &rLocalX, TVector3 &rLocalY, TVector3 &rLocalZ)
Definition: beam_math_utilities.hpp:723
static void VectorToSkewSymmetricTensor(const TVector3 &rVector, TMatrix3 &rSkewSymmetricTensor)
Definition: beam_math_utilities.hpp:231
Definition: build_string_skin_process.hpp:40
array_1d< double, 3 > PointType
Definition: build_string_skin_process.hpp:45
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: build_string_skin_process.hpp:90
Quadrilateral3D4< NodeType > Quadrilateral3DType
Definition: build_string_skin_process.hpp:55
Quaternion< double > QuaternionType
Definition: build_string_skin_process.hpp:58
BuildStringSkinProcess(BuildStringSkinProcess const &rOther)
Copy constructor.
Point3D< NodeType > Point3DType
Definition: build_string_skin_process.hpp:54
BuildStringSkinProcess(ModelPart &rModelPart, unsigned int sides, double radius)
Definition: build_string_skin_process.hpp:69
~BuildStringSkinProcess() override
Destructor.
Definition: build_string_skin_process.hpp:82
ModelPart::PropertiesType PropertiesType
Definition: build_string_skin_process.hpp:51
KRATOS_CLASS_POINTER_DEFINITION(BuildStringSkinProcess)
Pointer definition of BuildStringSkinProcess.
GlobalPointersVector< Element > ElementWeakPtrVectorType
Definition: build_string_skin_process.hpp:61
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: build_string_skin_process.hpp:155
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: build_string_skin_process.hpp:143
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: build_string_skin_process.hpp:180
ModelPart::ConditionType ConditionType
Definition: build_string_skin_process.hpp:50
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: build_string_skin_process.hpp:217
GlobalPointersVector< Condition > ConditionWeakPtrVectorType
Definition: build_string_skin_process.hpp:62
void ExecuteInitialize() override
Definition: build_string_skin_process.hpp:108
BeamMathUtils< double > BeamMathUtilsType
Definition: build_string_skin_process.hpp:57
void Execute() override
Execute method is used to execute the BuildStringSkinProcess algorithms.
Definition: build_string_skin_process.hpp:102
ModelPart::NodeType NodeType
Definition: build_string_skin_process.hpp:46
ModelPart::NodesContainerType NodesContainerType
Definition: build_string_skin_process.hpp:47
void ExecuteBeforeSolutionLoop() override
Definition: build_string_skin_process.hpp:137
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: build_string_skin_process.hpp:168
ModelPart::ElementsContainerType ElementsContainerType
Definition: build_string_skin_process.hpp:48
NodesContainerType::Pointer NodesContainerPointerType
Definition: build_string_skin_process.hpp:49
Triangle3D3< NodeType > Triangle3DType
Definition: build_string_skin_process.hpp:53
std::string Info() const override
Turn back information as a string.
Definition: build_string_skin_process.hpp:211
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: build_string_skin_process.hpp:223
void ExecuteFinalize() override
Definition: build_string_skin_process.hpp:191
ConditionType::GeometryType GeometryType
Definition: build_string_skin_process.hpp:52
GlobalPointersVector< Node > NodeWeakPtrVectorType
Definition: build_string_skin_process.hpp:60
Base class for all Conditions.
Definition: condition.h:59
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
void FixDof()
Definition: dof.h:338
Geometry< NodeType > GeometryType
definition of the geometry type with given NodeType
Definition: element.h:83
void Set(const Flags ThisFlag)
Definition: flags.cpp:33
Flags()
Default constructor.
Definition: flags.h:119
Geometry base class.
Definition: geometry.h:71
SizeType size() const
Definition: geometry.h:518
IndexType const & Id() const
Id of this Geometry.
Definition: geometry.h:964
PointerVector< TPointType > PointsArrayType
Definition: geometry.h:118
virtual SizeType FacesNumber() const
Returns the number of faces of the current geometry.
Definition: geometry.h:2152
This class is a vector which stores global pointers.
Definition: global_pointers_vector.h:61
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: global_pointers_vector.h:79
void clear()
Definition: global_pointers_vector.h:361
size_type size() const
Definition: global_pointers_vector.h:307
void resize(size_type new_dim) const
Definition: global_pointers_vector.h:366
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ConditionIterator ConditionsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1361
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
SubModelPartsContainerType::iterator SubModelPartIterator
Iterator over the sub model parts of this model part.
Definition: model_part.h:258
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
NodeIterator NodesEnd(IndexType ThisIndex=0)
Definition: model_part.h:497
ModelPart & GetParentModelPart()
Definition: model_part.cpp:2124
MeshType::ConditionIterator ConditionIterator
Definition: model_part.h:189
ConditionIterator ConditionsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1371
This class defines the node.
Definition: node.h:65
std::vector< std::unique_ptr< Dof< double > >> DofsContainerType
The DoF container type definition.
Definition: node.h:92
Definition: point_3d.h:53
void reserve(size_type dim)
Definition: pointer_vector.h:319
The base class for all processes in Kratos.
Definition: process.h:49
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
A four node 3D quadrilateral geometry with bi-linear shape functions.
Definition: quadrilateral_3d_4.h:76
static Quaternion FromRotationVector(double rx, double ry, double rz)
Definition: quaternion.h:427
A three node 3D triangle geometry with linear shape functions.
Definition: triangle_3d_3.h:77
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
Kratos::ModelPart ModelPart
Definition: kratos_wrapper.h:31
constexpr double Pi
Definition: global_variables.h:25
Parameters GetValue(Parameters &rParameters, const std::string &rEntry)
Definition: add_kratos_parameters_to_python.cpp:53
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
TExpressionType::data_type norm_2(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:625
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
AMatrix::MatrixProductExpression< TExpression1Type, TExpression2Type > prod(AMatrix::MatrixExpression< TExpression1Type, TCategory1 > const &First, AMatrix::MatrixExpression< TExpression2Type, TCategory2 > const &Second)
Definition: amatrix_interface.h:568
TExpression1Type::data_type inner_prod(AMatrix::MatrixExpression< TExpression1Type, TCategory1 > const &First, AMatrix::MatrixExpression< TExpression2Type, TCategory2 > const &Second)
Definition: amatrix_interface.h:592
Element::WeakPointer ElementWeakPtrType
Definition: generate_new_conditions_mesher_process.hpp:44
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
Temp
Definition: PecletTest.py:105
v
Definition: generate_convection_diffusion_explicit_element.py:114
alpha
Definition: generate_convection_diffusion_explicit_element.py:113
tuple Q
Definition: isotropic_damage_automatic_differentiation.py:235
float radius
Definition: mesh_to_mdpa_converter.py:18
int k
Definition: quadrature.py:595
int j
Definition: quadrature.py:648
int node_id
Definition: read_stl.py:12
int counter
Definition: script_THERMAL_CORRECT.py:218
int current_id
Output settings end ####.
Definition: script.py:194
integer i
Definition: TensorModule.f:17
Configure::PointType PointType
Definition: transfer_utility.h:245