10 #if !defined(KRATOS_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED)
11 #define KRATOS_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED
48 const std::string EntityType
49 ) :
Process(
Flags()) , mrDestinationModelPart(rDestinationModelPart), mrOriginModelPart(rOriginModelPart), mEntityType(EntityType), mrTransferFlags(
std::vector<
Flags>()), mrAssignFlags(
std::vector<
Flags>())
59 const std::string EntityType,
60 const std::vector<Flags>& rTransferFlags
61 ) :
Process(
Flags()) , mrDestinationModelPart(rDestinationModelPart), mrOriginModelPart(rOriginModelPart), mEntityType(EntityType), mrTransferFlags(rTransferFlags), mrAssignFlags(
std::vector<
Flags>() )
72 const std::string EntityType,
73 const std::vector<Flags>& rTransferFlags,
74 const std::vector<Flags>& rAssignFlags
75 ) :
Process(
Flags()) , mrDestinationModelPart(rDestinationModelPart), mrOriginModelPart(rOriginModelPart), mEntityType(EntityType), mrTransferFlags(rTransferFlags), mrAssignFlags(rAssignFlags)
108 if (mEntityType ==
"Nodes")
110 const int nnodes = mrOriginModelPart.
Nodes().size();
114 ModelPart::NodesContainerType::iterator it_begin =
120 ModelPart::NodesContainerType::iterator it = it_begin +
i;
122 if (this->MatchTransferFlags(*(it.base())))
125 this->AssignFlags(*(it.base()));
126 mrDestinationModelPart.
Nodes().push_back(*(it.base()));
130 mrDestinationModelPart.
Nodes().Unique();
133 else if (mEntityType ==
"Elements")
135 const int nelements = mrOriginModelPart.
Elements().size();
139 ModelPart::ElementsContainerType::iterator it_begin =
143 for (
int i = 0;
i < nelements;
i++)
145 ModelPart::ElementsContainerType::iterator it = it_begin +
i;
147 if (this->MatchTransferFlags(*(it.base())))
150 this->AssignFlags(*(it.base()));
151 mrDestinationModelPart.
Elements().push_back(*(it.base()));
155 mrDestinationModelPart.
Elements().Unique();
158 else if (mEntityType ==
"Conditions")
160 const int nconditions = mrOriginModelPart.
Conditions().size();
162 if (nconditions != 0)
164 ModelPart::ConditionsContainerType::iterator it_begin =
168 for (
int i = 0;
i < nconditions;
i++)
170 ModelPart::ConditionsContainerType::iterator it = it_begin +
i;
172 if (this->MatchTransferFlags(*(it.base())))
175 this->AssignFlags(*(it.base()));
176 mrDestinationModelPart.
Conditions().push_back(*(it.base()));
246 std::string
Info()
const override
248 return "TransferBetweenModelPartsProcess";
254 rOStream <<
"TransferBetweenModelPartsProcess";
307 const std::string mEntityType;
309 const std::vector<Flags> mrTransferFlags;
310 const std::vector<Flags> mrAssignFlags;
316 bool MatchTransferFlags(
const Node::Pointer& pNode)
319 for(
unsigned int i = 0;
i<mrTransferFlags.size();
i++)
321 if( pNode->IsNot(mrTransferFlags[
i]) )
329 void AssignFlags(
const Node::Pointer& pNode)
332 for(
unsigned int i = 0;
i<mrAssignFlags.size();
i++)
333 pNode->
Set(mrAssignFlags[
i]);
338 bool MatchTransferFlags(
const Element::Pointer& pElement)
341 for(
unsigned int i = 0;
i<mrTransferFlags.size();
i++)
343 if( pElement->IsNot(mrTransferFlags[
i]) )
351 void AssignFlags(
const Element::Pointer& pElement)
354 for(
unsigned int i = 0;
i<mrAssignFlags.size();
i++)
355 pElement->Set(mrAssignFlags[
i]);
360 bool MatchTransferFlags(
const Condition::Pointer& pCondition)
363 for(
unsigned int i = 0;
i<mrTransferFlags.size();
i++)
365 if( pCondition->IsNot(mrTransferFlags[
i]) )
373 void AssignFlags(
const Condition::Pointer& pCondition)
376 for(
unsigned int i = 0;
i<mrAssignFlags.size();
i++)
377 pCondition->Set(mrAssignFlags[
i]);
425 rOStream << std::endl;
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
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: transfer_between_model_parts_process.hpp:34
KRATOS_CLASS_POINTER_DEFINITION(TransferBetweenModelPartsProcess)
Pointer definition of TransferBetweenModelPartsProcess.
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: transfer_between_model_parts_process.hpp:219
~TransferBetweenModelPartsProcess() override
Destructor.
Definition: transfer_between_model_parts_process.hpp:84
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: transfer_between_model_parts_process.hpp:213
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: transfer_between_model_parts_process.hpp:202
std::string Info() const override
Turn back information as a string.
Definition: transfer_between_model_parts_process.hpp:246
TransferBetweenModelPartsProcess(ModelPart &rDestinationModelPart, ModelPart &rOriginModelPart, const std::string EntityType, const std::vector< Flags > &rTransferFlags, const std::vector< Flags > &rAssignFlags)
Definition: transfer_between_model_parts_process.hpp:70
TransferBetweenModelPartsProcess(ModelPart &rDestinationModelPart, ModelPart &rOriginModelPart, const std::string EntityType, const std::vector< Flags > &rTransferFlags)
Definition: transfer_between_model_parts_process.hpp:57
TransferBetweenModelPartsProcess(ModelPart &rDestinationModelPart, ModelPart &rOriginModelPart, const std::string EntityType)
Definition: transfer_between_model_parts_process.hpp:46
TransferBetweenModelPartsProcess(TransferBetweenModelPartsProcess const &rOther)
Copy constructor.
void ExecuteFinalize() override
Definition: transfer_between_model_parts_process.hpp:226
void Execute() override
Execute method is used to execute the TransferBetweenModelPartsProcess algorithms.
Definition: transfer_between_model_parts_process.hpp:104
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: transfer_between_model_parts_process.hpp:92
void ExecuteInitialize() override
Definition: transfer_between_model_parts_process.hpp:190
void ExecuteBeforeSolutionLoop() override
Definition: transfer_between_model_parts_process.hpp:196
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: transfer_between_model_parts_process.hpp:207
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: transfer_between_model_parts_process.hpp:252
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: transfer_between_model_parts_process.hpp:258
#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