10 #if !defined(KRATOS_ASSIGN_FLAGS_TO_MODEL_PART_ENTITIES_PROCESS_H_INCLUDED)
11 #define KRATOS_ASSIGN_FLAGS_TO_MODEL_PART_ENTITIES_PROCESS_H_INCLUDED
47 const std::vector<Flags>& rAssignFlags) :
Process(
Flags()) , mrModelPart(rModelPart), mEntityType(EntityType), mrTransferFlags(
std::vector<
Flags>()), mrAssignFlags(rAssignFlags)
56 const std::vector<Flags>& rTransferFlags) :
Process(
Flags()) , mrModelPart(rModelPart), mEntityType(EntityType), mrTransferFlags(rTransferFlags), mrAssignFlags(rAssignFlags)
89 if (mEntityType ==
"Nodes")
95 ModelPart::NodesContainerType::iterator it_begin =
98 #pragma omp parallel for
101 ModelPart::NodesContainerType::iterator it = it_begin +
i;
103 if (this->MatchTransferFlags(*(it.base())))
105 this->AssignFlags(*(it.base()));
110 else if (mEntityType ==
"Elements")
112 const int nelements = mrModelPart.
Elements().size();
116 ModelPart::ElementsContainerType::iterator it_begin =
119 #pragma omp parallel for
120 for (
int i = 0;
i < nelements;
i++)
122 ModelPart::ElementsContainerType::iterator it = it_begin +
i;
124 if (this->MatchTransferFlags(*(it.base())))
126 this->AssignFlags(*(it.base()));
131 else if (mEntityType ==
"Conditions")
133 const int nconditions = mrModelPart.
Conditions().size();
135 if (nconditions != 0)
137 ModelPart::ConditionsContainerType::iterator it_begin =
141 for (
int i = 0;
i < nconditions;
i++)
143 ModelPart::ConditionsContainerType::iterator it = it_begin +
i;
145 if (this->MatchTransferFlags(*(it.base())))
147 this->AssignFlags(*(it.base()));
214 std::string
Info()
const override
216 return "AssignFlagsToModelPartEntitiesProcess";
222 rOStream <<
"AssignFlagsToModelPartEntitiesProcess";
274 const std::string mEntityType;
276 const std::vector<Flags> mrTransferFlags;
277 const std::vector<Flags> mrAssignFlags;
283 bool MatchTransferFlags(
const Node::Pointer& pNode)
286 for(
unsigned int i = 0;
i<mrTransferFlags.size();
i++)
288 if( pNode->IsNot(mrTransferFlags[
i]) )
296 void AssignFlags(
const Node::Pointer& pNode)
299 for(
unsigned int i = 0;
i<mrAssignFlags.size();
i++)
300 pNode->
Set(mrAssignFlags[
i]);
305 bool MatchTransferFlags(
const Element::Pointer& pElement)
308 for(
unsigned int i = 0;
i<mrTransferFlags.size();
i++)
310 if( pElement->IsNot(mrTransferFlags[
i]) )
318 void AssignFlags(
const Element::Pointer& pElement)
321 for(
unsigned int i = 0;
i<mrAssignFlags.size();
i++)
322 pElement->Set(mrAssignFlags[
i]);
327 bool MatchTransferFlags(
const Condition::Pointer& pCondition)
330 for(
unsigned int i = 0;
i<mrTransferFlags.size();
i++)
332 if( pCondition->IsNot(mrTransferFlags[
i]) )
340 void AssignFlags(
const Condition::Pointer& pCondition)
343 for(
unsigned int i = 0;
i<mrAssignFlags.size();
i++)
344 pCondition->Set(mrAssignFlags[
i]);
392 rOStream << std::endl;
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_flags_to_model_part_entities_process.hpp:33
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: assign_flags_to_model_part_entities_process.hpp:187
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: assign_flags_to_model_part_entities_process.hpp:170
AssignFlagsToModelPartEntitiesProcess(ModelPart &rModelPart, const std::string EntityType, const std::vector< Flags > &rAssignFlags, const std::vector< Flags > &rTransferFlags)
Definition: assign_flags_to_model_part_entities_process.hpp:55
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: assign_flags_to_model_part_entities_process.hpp:175
void ExecuteBeforeSolutionLoop() override
Definition: assign_flags_to_model_part_entities_process.hpp:164
void ExecuteFinalize() override
Definition: assign_flags_to_model_part_entities_process.hpp:194
AssignFlagsToModelPartEntitiesProcess(AssignFlagsToModelPartEntitiesProcess const &rOther)
Copy constructor.
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: assign_flags_to_model_part_entities_process.hpp:181
void Execute() override
Execute method is used to execute the AssignFlagsToModelPartEntitiesProcess algorithms.
Definition: assign_flags_to_model_part_entities_process.hpp:85
~AssignFlagsToModelPartEntitiesProcess() override
Destructor.
Definition: assign_flags_to_model_part_entities_process.hpp:65
std::string Info() const override
Turn back information as a string.
Definition: assign_flags_to_model_part_entities_process.hpp:214
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_flags_to_model_part_entities_process.hpp:226
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_flags_to_model_part_entities_process.hpp:73
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_flags_to_model_part_entities_process.hpp:220
AssignFlagsToModelPartEntitiesProcess(ModelPart &rModelPart, const std::string EntityType, const std::vector< Flags > &rAssignFlags)
Definition: assign_flags_to_model_part_entities_process.hpp:46
void ExecuteInitialize() override
Definition: assign_flags_to_model_part_entities_process.hpp:158
KRATOS_CLASS_POINTER_DEFINITION(AssignFlagsToModelPartEntitiesProcess)
Pointer definition of AssignFlagsToModelPartEntitiesProcess.
void Set(const Flags ThisFlag)
Definition: flags.cpp:33
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ElementIterator ElementsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1169
ConditionIterator ConditionsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1361
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ConditionsContainerType & Conditions(IndexType ThisIndex=0)
Definition: model_part.h:1381
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
The base class for all processes in Kratos.
Definition: process.h:49
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
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 nnodes
Definition: sensitivityMatrix.py:24
namespace
Definition: array_1d.h:793
integer i
Definition: TensorModule.f:17