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.
fast_transfer_between_model_parts_process.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Vicente Mataix Ferrandiz
11 //
12 //
13 
14 #if !defined(KRATOS_FAST_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED )
15 #define KRATOS_FAST_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "includes/model_part.h"
23 #include "processes/process.h"
24 #include "utilities/openmp_utils.h"
25 
26 namespace Kratos
27 {
30 
34 
38 
42 
45 
54 class KRATOS_API(KRATOS_CORE) FastTransferBetweenModelPartsProcess
55  : public Process
56 {
57 public:
60 
62  typedef Node NodeType;
63 
64  // General containers type definitions
70 
72  typedef std::size_t SizeType;
73 
76 
80 
85  enum class EntityTransfered {
86  NODES = 0,
87  ELEMENTS = 1,
88  NODESANDELEMENTS = 2,
89  CONDITIONS = 3,
90  NODESANDCONDITIONS = 4,
91  CONSTRAINTS = 5,
92  NODESANDCONSTRAINTS = 6,
93  GEOMETRIES = 7,
94  NODESANDGEOMETRIES = 8,
95  ALL = 9
96  };
97 
101 
111  ModelPart& rDestinationModelPart,
112  ModelPart& rOriginModelPart,
113  const EntityTransfered Entity = EntityTransfered::ALL,
114  const Flags Flag = Flags(),
115  const bool ReplicateEntities = false
116  );
117 
120 
124 
126  void operator()();
127 
128 
132 
134  void Execute() override;
135 
139 
140 
144 
145 
149 
151  std::string Info() const override
152  {
153  return "FastTransferBetweenModelPartsProcess";
154  }
155 
157  void PrintInfo(std::ostream& rOStream) const override
158  {
159  rOStream << "FastTransferBetweenModelPartsProcess";
160  }
161 
163  void PrintData(std::ostream& rOStream) const override
164  {
165  }
166 
171 
172 protected:
173 
182 
183 // /// Copy constructor.
184 // FastTransferBetweenModelPartsProcess(FastTransferBetweenModelPartsProcess const& rOther);
185 
199 
200 private:
201 
207 
208  ModelPart& mrDestinationModelPart;
209 
210  ModelPart& mrOriginModelPart;
211 
212  const EntityTransfered mEntity;
213 
214  const Flags mFlag;
215 
219 
223 
227  void TransferWithoutFlags();
228 
232  void TransferWithFlags();
233 
238  void ReorderAllIds(ModelPart& rThhisModelPart);
239 
243  void ReplicateWithoutFlags();
244 
248  void ReplicateWithFlags();
249 
253 
256 
266 
267 }; // Class FastTransferBetweenModelPartsProcess
268 
269 
271 
274 
275 
279 
280 
282 inline std::istream& operator >> (std::istream& rIStream,
284 
286 inline std::ostream& operator << (std::ostream& rOStream,
288 {
289  rThis.PrintInfo(rOStream);
290  rOStream << std::endl;
291  rThis.PrintData(rOStream);
292 
293  return rOStream;
294 }
296 
297 
298 } // namespace Kratos.
299 
300 #endif // KRATOS_FAST_TRANSFER_BETWEEN_MODEL_PARTS_PROCESS_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: fast_transfer_between_model_parts_process.h:56
EntityTransfered
This enum helps us to identify the elements to transfer between the modelparts.
Definition: fast_transfer_between_model_parts_process.h:85
ModelPart::MasterSlaveConstraintContainerType MasterSlaveConstraintArrayType
Definition: fast_transfer_between_model_parts_process.h:68
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: fast_transfer_between_model_parts_process.h:163
ModelPart::ElementsContainerType ElementsArrayType
Definition: fast_transfer_between_model_parts_process.h:67
ModelPart::NodesContainerType NodesArrayType
Definition: fast_transfer_between_model_parts_process.h:65
ModelPart::ConditionsContainerType ConditionsArrayType
Definition: fast_transfer_between_model_parts_process.h:66
Node NodeType
The type used for the node.
Definition: fast_transfer_between_model_parts_process.h:62
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: fast_transfer_between_model_parts_process.h:157
~FastTransferBetweenModelPartsProcess() override
Destructor.
Definition: fast_transfer_between_model_parts_process.h:119
std::string Info() const override
Turn back information as a string.
Definition: fast_transfer_between_model_parts_process.h:151
std::size_t SizeType
The type used to identify the size.
Definition: fast_transfer_between_model_parts_process.h:72
KRATOS_CLASS_POINTER_DEFINITION(FastTransferBetweenModelPartsProcess)
Pointer definition of FastTransferBetweenModelPartsProcess.
ModelPart::GeometriesMapType GeometriesMapType
Definition: fast_transfer_between_model_parts_process.h:69
Definition: flags.h:58
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class defines the node.
Definition: node.h:65
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
The base class for all processes in Kratos.
Definition: process.h:49
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