10 #if !defined(KRATOS_GENERATE_NEW_CONTACT_CONDITIONS_MESHER_PROCESS_H_INCLUDED )
11 #define KRATOS_GENERATE_NEW_CONTACT_CONDITIONS_MESHER_PROCESS_H_INCLUDED
19 #include <boost/timer.hpp>
82 : mrModelPart(rModelPart),
83 mrRemesh(rRemeshingParameters)
112 if( mEchoLevel > 0 ){
113 std::cout<<
" [ GENERATE NEW CONTACT ELEMENTS: "<<std::endl;
114 std::cout<<
" Total Conditions BEFORE: ["<<mrModelPart.
Conditions().size()<<
"] ];"<<std::endl;
118 std::cout<<
" ModelPart Supplied do not corresponds to the Meshing Domain: ("<<mrModelPart.
Name()<<
" != "<<mrRemesh.
SubModelPartName<<
")"<<std::endl;
124 std::cout<<
" ERROR : no selection of elements performed before building the elements "<<std::endl;
148 std::cout<<
" [START contact Element Generation "<<std::endl;
156 ModelPart::NodesContainerType::iterator nodes_begin = mrModelPart.
NodesBegin();
158 for(
int el = 0;
el<OutNumberOfElements; ++
el)
163 for(
unsigned int i=0;
i<nds; ++
i)
180 Vertices.
push_back(*(nodes_begin + OutElementList[
el*nds+
i]-1).base());
182 Vertices.
back().Set(CONTACT);
187 Condition::Pointer pContactCondition = rReferenceCondition.
Create(
id, Vertices, pProperties);
193 bool condition_found=
false;
197 if(!condition_found){
199 std::cout<<
" MASTER CONDITION NOT FOUND: Contact Element Release "<<std::endl;
210 pContactCondition->SetValue(MASTER_CONDITION, pMasterCondition);
211 pContactCondition->SetValue(MASTER_ELEMENTS, pMasterCondition->GetValue(MASTER_ELEMENTS) );
212 pContactCondition->SetValue(MASTER_NODES, pMasterCondition->GetValue(MASTER_NODES) );
214 if( pContactCondition->Is(SELECTED) ){
221 std::vector<bool> edge_nodes(4);
222 std::fill(edge_nodes.begin(), edge_nodes.end(),
false);
228 if(rGeometry[
j].Id()==rMasterGeometry[
i].Id()){
229 edge_nodes[
i] =
true;
235 for(
unsigned int i=0;
i<4; ++
i)
237 if(!edge_nodes[
i] && rMasterGeometry[
i].Id() != rMasterNode.
Id())
238 pContactCondition->GetValue(MASTER_NODES).push_back( rMasterGeometry(
i) );
243 pContactCondition->SetValue(NORMAL, pMasterCondition->GetValue(NORMAL) );
244 pContactCondition->Set(CONTACT);
249 pContactCondition->Set(ACTIVE);
260 for(ModelPart::NodesContainerType::iterator in = mrModelPart.
NodesBegin() ; in != mrModelPart.
NodesEnd(); ++in)
265 if( mEchoLevel > 0 ){
266 std::cout<<
" [END contact Elements Generation ["<<
id-previous_id<<
"] ]"<<std::endl;
268 std::cout<<
" Total Conditions AFTER: ["<<mrModelPart.
Conditions().size()<<
"] ];"<<std::endl;
272 std::cout<<
" [Contact Candidates:"<<
id-previous_id<<
"]"<<std::endl;
294 std::string
Info()
const override
296 return "GenerateNewContactConditionsMesherProcess";
302 rOStream <<
"GenerateNewContactConditionsMesherProcess";
424 rOStream << std::endl;
Base class for all Conditions.
Definition: condition.h:59
virtual Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new condition pointer.
Definition: condition.h:205
Base class for all Elements.
Definition: element.h:60
TVariableType::Type & GetValue(const TVariableType &rThisVariable)
Definition: geometrical_object.h:248
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
Geometry base class.
Definition: geometry.h:71
SizeType PointsNumber() const
Definition: geometry.h:528
void push_back(PointPointerType x)
Definition: geometry.h:548
SizeType size() const
Definition: geometry.h:518
PointReferenceType back()
Definition: geometry.h:507
The base class for processes passed to the solution scheme.
Definition: mesher_process.hpp:37
Short class definition.
Definition: mesher_utilities.hpp:49
Condition::Pointer FindMasterCondition(Condition::Pointer &pCondition, ModelPart::ConditionsContainerType &rModelConditions, bool &condition_found)
Definition: mesher_utilities.cpp:1642
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
std::string & Name()
Definition: model_part.h:1811
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ConditionsContainerType & Conditions(IndexType ThisIndex=0)
Definition: model_part.h:1381
void AddCondition(ConditionType::Pointer pNewCondition, IndexType ThisIndex=0)
Definition: model_part.cpp:1436
NodeIterator NodesEnd(IndexType ThisIndex=0)
Definition: model_part.h:497
ModelPart & GetParentModelPart()
Definition: model_part.cpp:2124
ConditionIterator ConditionsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1371
This class defines the node.
Definition: node.h:65
IndexType Id() const
Definition: node.h:262
TVariableType::Type & GetValue(const TVariableType &rThisVariable)
Definition: node.h:466
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
Refine Mesh Elements Process 2D and 3D.
Definition: select_elements_mesher_process.hpp:47
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: select_elements_mesher_process.hpp:95
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
static int EchoLevel
Definition: co_sim_EMPIRE_API.h:42
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
int j
Definition: quadrature.py:648
el
Definition: read_stl.py:25
integer i
Definition: TensorModule.f:17
int * GetElementList()
Definition: mesher_utilities.hpp:178
int & GetNumberOfElements()
Definition: mesher_utilities.hpp:183
Definition: mesher_utilities.hpp:631
Condition const & GetReferenceCondition()
Definition: mesher_utilities.hpp:847
std::vector< int > PreservedElements
Definition: mesher_utilities.hpp:669
bool MeshElementsSelectedFlag
Definition: mesher_utilities.hpp:662
MeshContainer OutMesh
Definition: mesher_utilities.hpp:675
PropertiesType::Pointer GetProperties()
Definition: mesher_utilities.hpp:838
std::string SubModelPartName
Definition: mesher_utilities.hpp:642
std::vector< int > NodalPreIds
Definition: mesher_utilities.hpp:668