KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
transfer_between_model_parts_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosDelaunayMeshingApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: August 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED)
11 #define KRATOS_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED
12 
13 
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
20 #include "includes/model_part.h"
22 #include "processes/process.h"
23 
24 namespace Kratos
25 {
26 
29 
31 
34 {
35 public:
38 
41 
45 
47  ModelPart& rOriginModelPart,
48  const std::string EntityType
49  ) : Process(Flags()) , mrDestinationModelPart(rDestinationModelPart), mrOriginModelPart(rOriginModelPart), mEntityType(EntityType), mrTransferFlags(std::vector<Flags>()), mrAssignFlags(std::vector<Flags>())
50  {
52 
53  KRATOS_CATCH("");
54  }
55 
56 
58  ModelPart& rOriginModelPart,
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>() )
62  {
64 
65 
66  KRATOS_CATCH("");
67  }
68 
69 
71  ModelPart& rOriginModelPart,
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)
76  {
78 
79  KRATOS_CATCH("");
80  }
81 
82 
85 
86 
90 
92  void operator()()
93  {
94  Execute();
95  }
96 
97 
101 
102 
104  void Execute() override
105  {
106  KRATOS_TRY;
107 
108  if (mEntityType == "Nodes")
109  {
110  const int nnodes = mrOriginModelPart.Nodes().size();
111 
112  if (nnodes != 0)
113  {
114  ModelPart::NodesContainerType::iterator it_begin =
115  mrOriginModelPart.NodesBegin();
116  //#pragma omp parallel for //some nodes are not added in
117  //parallel
118  for (int i = 0; i < nnodes; i++)
119  {
120  ModelPart::NodesContainerType::iterator it = it_begin + i;
121 
122  if (this->MatchTransferFlags(*(it.base())))
123  {
124  // mrDestinationModelPart.AddNode(*(it.base()));
125  this->AssignFlags(*(it.base()));
126  mrDestinationModelPart.Nodes().push_back(*(it.base()));
127  // std::cout<<" Node Inserted "<<it->Id()<<std::endl;
128  }
129  }
130  mrDestinationModelPart.Nodes().Unique();
131  }
132  }
133  else if (mEntityType == "Elements")
134  {
135  const int nelements = mrOriginModelPart.Elements().size();
136 
137  if (nelements != 0)
138  {
139  ModelPart::ElementsContainerType::iterator it_begin =
140  mrOriginModelPart.ElementsBegin();
141  //#pragma omp parallel for //some elements are not added in
142  // parallel
143  for (int i = 0; i < nelements; i++)
144  {
145  ModelPart::ElementsContainerType::iterator it = it_begin + i;
146 
147  if (this->MatchTransferFlags(*(it.base())))
148  {
149  // mrDestinationModelPart.AddNode(*(it.base()));
150  this->AssignFlags(*(it.base()));
151  mrDestinationModelPart.Elements().push_back(*(it.base()));
152  // std::cout<<" Element Inserted "<<it->Id()<<std::endl;
153  }
154  }
155  mrDestinationModelPart.Elements().Unique();
156  }
157  }
158  else if (mEntityType == "Conditions")
159  {
160  const int nconditions = mrOriginModelPart.Conditions().size();
161 
162  if (nconditions != 0)
163  {
164  ModelPart::ConditionsContainerType::iterator it_begin =
165  mrOriginModelPart.ConditionsBegin();
166  //#pragma omp parallel for //some elements are not added in
167  // parallel
168  for (int i = 0; i < nconditions; i++)
169  {
170  ModelPart::ConditionsContainerType::iterator it = it_begin + i;
171 
172  if (this->MatchTransferFlags(*(it.base())))
173  {
174  // mrDestinationModelPart.AddNode(*(it.base()));
175  this->AssignFlags(*(it.base()));
176  mrDestinationModelPart.Conditions().push_back(*(it.base()));
177  // std::cout<<" Condition Inserted
178  // "<<it->Id()<<std::endl;
179  }
180  }
181  mrDestinationModelPart.Conditions().Unique();
182  }
183  }
184 
185  KRATOS_CATCH("");
186  }
187 
190  void ExecuteInitialize() override
191  {
192  }
193 
197  {
198  }
199 
200 
203  {
204  }
205 
208  {
209  }
210 
211 
213  void ExecuteBeforeOutputStep() override
214  {
215  }
216 
217 
219  void ExecuteAfterOutputStep() override
220  {
221  }
222 
223 
226  void ExecuteFinalize() override
227  {
228  }
229 
230 
234 
235 
239 
240 
244 
246  std::string Info() const override
247  {
248  return "TransferBetweenModelPartsProcess";
249  }
250 
252  void PrintInfo(std::ostream& rOStream) const override
253  {
254  rOStream << "TransferBetweenModelPartsProcess";
255  }
256 
258  void PrintData(std::ostream& rOStream) const override
259  {
260  }
261 
262 
267 
268 protected:
269 
278 
281 
295 
296 private:
297 
303 
304  ModelPart& mrDestinationModelPart;
305  ModelPart& mrOriginModelPart;
306 
307  const std::string mEntityType;
308 
309  const std::vector<Flags> mrTransferFlags;
310  const std::vector<Flags> mrAssignFlags;
311 
315 
316  bool MatchTransferFlags(const Node::Pointer& pNode)
317  {
318 
319  for(unsigned int i = 0; i<mrTransferFlags.size(); i++)
320  {
321  if( pNode->IsNot(mrTransferFlags[i]) )
322  return false;
323  }
324 
325  return true;
326 
327  }
328 
329  void AssignFlags(const Node::Pointer& pNode)
330  {
331 
332  for(unsigned int i = 0; i<mrAssignFlags.size(); i++)
333  pNode->Set(mrAssignFlags[i]);
334 
335  }
336 
337 
338  bool MatchTransferFlags(const Element::Pointer& pElement)
339  {
340 
341  for(unsigned int i = 0; i<mrTransferFlags.size(); i++)
342  {
343  if( pElement->IsNot(mrTransferFlags[i]) )
344  return false;
345  }
346 
347  return true;
348 
349  }
350 
351  void AssignFlags(const Element::Pointer& pElement)
352  {
353 
354  for(unsigned int i = 0; i<mrAssignFlags.size(); i++)
355  pElement->Set(mrAssignFlags[i]);
356 
357  }
358 
359 
360  bool MatchTransferFlags(const Condition::Pointer& pCondition)
361  {
362 
363  for(unsigned int i = 0; i<mrTransferFlags.size(); i++)
364  {
365  if( pCondition->IsNot(mrTransferFlags[i]) )
366  return false;
367  }
368 
369  return true;
370 
371  }
372 
373  void AssignFlags(const Condition::Pointer& pCondition)
374  {
375 
376  for(unsigned int i = 0; i<mrAssignFlags.size(); i++)
377  pCondition->Set(mrAssignFlags[i]);
378 
379  }
380 
387 
390 
391 
401 
402 }; // Class TransferBetweenModelPartsProcess
403 
404 
406 
409 
410 
414 
415 
417 inline std::istream& operator >> (std::istream& rIStream,
419 
421 inline std::ostream& operator << (std::ostream& rOStream,
423 {
424  rThis.PrintInfo(rOStream);
425  rOStream << std::endl;
426  rThis.PrintData(rOStream);
427 
428  return rOStream;
429 }
431 
432 
433 } // namespace Kratos.
434 
435 #endif // KRATOS_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED defined
Definition: flags.h:58
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