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_model_part_elements_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemFluidDynamicsApplication $
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_MODEL_PART_ELEMENTS_PROCESS_H_INCLUDED)
11 #define KRATOS_TRANSFER_MODEL_PART_ELEMENTS_PROCESS_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 #include "includes/model_part.h"
20 #include "processes/process.h"
21 
22 namespace Kratos
23 {
24 
27 
29 
32 {
33 public:
36 
39 
44  ModelPart &rGuestModelPart) : mrHostModelPart(rHostModelPart), mrGuestModelPart(rGuestModelPart)
45  {
47 
48  KRATOS_CATCH("");
49  }
50 
53 
57 
59  void operator()()
60  {
61  Execute();
62  }
63 
67 
69  void Execute() override
70  {
71  KRATOS_TRY;
72 
73  const int nel = mrGuestModelPart.Elements().size();
74 
75  if (nel != 0)
76  {
77  ModelPart::ElementsContainerType::iterator el_begin = mrGuestModelPart.ElementsBegin();
78 
79  //#pragma omp parallel for //some nodes are not added in parallel
80  for (int i = 0; i < nel; i++)
81  {
82  ModelPart::ElementsContainerType::iterator el = el_begin + i;
83 
84  mrHostModelPart.Elements().push_back(*(el.base()));
85  }
86  }
87 
88  KRATOS_CATCH("");
89  }
90 
94 
98 
102 
104  std::string Info() const override
105  {
106  return "TransferModelPartElementsProcess";
107  }
108 
110  void PrintInfo(std::ostream &rOStream) const override
111  {
112  rOStream << "TransferModelPartElementsProcess";
113  }
114 
119 
120 protected:
129 
132 
146 
147 private:
153 
154  ModelPart &mrHostModelPart;
155  ModelPart &mrGuestModelPart;
156 
163 
166 
176 
177 }; // Class TransferModelPartElementsProcess
178 
180 
183 
187 
189 inline std::istream &operator>>(std::istream &rIStream,
191 
193 inline std::ostream &operator<<(std::ostream &rOStream,
195 {
196  rThis.PrintInfo(rOStream);
197  rOStream << std::endl;
198  rThis.PrintData(rOStream);
199 
200  return rOStream;
201 }
203 
204 } // namespace Kratos.
205 
206 #endif // KRATOS_TRANSFER_MODEL_PART_ELEMENTS_PROCESS_H_INCLUDED defined
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
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
The base class for all processes in Kratos.
Definition: process.h:49
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: process.h:210
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: transfer_model_part_elements_process.hpp:32
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: transfer_model_part_elements_process.hpp:110
std::string Info() const override
Turn back information as a string.
Definition: transfer_model_part_elements_process.hpp:104
void Execute() override
Execute method is used to execute the TransferModelPartElementsProcess algorithms.
Definition: transfer_model_part_elements_process.hpp:69
TransferModelPartElementsProcess(ModelPart &rHostModelPart, ModelPart &rGuestModelPart)
Definition: transfer_model_part_elements_process.hpp:43
virtual ~TransferModelPartElementsProcess()
Destructor.
Definition: transfer_model_part_elements_process.hpp:52
TransferModelPartElementsProcess(TransferModelPartElementsProcess const &rOther)
Copy constructor.
KRATOS_CLASS_POINTER_DEFINITION(TransferModelPartElementsProcess)
Pointer definition of TransferModelPartElementsProcess.
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: transfer_model_part_elements_process.hpp:59
#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
el
Definition: read_stl.py:25
integer i
Definition: TensorModule.f:17