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.
set_mesh_velocity_for_thermal_coupling_process.hpp
Go to the documentation of this file.
1 //-------------------------------------------------------------
2 // ___ __ ___ _ _ _
3 // KRATOS| _ \/ _|___ _ __ | __| |_ _(_)__| |
4 // | _/ _/ -_) ' \| _|| | || | / _` |
5 // |_| |_| \___|_|_|_|_| |_|\_,_|_\__,_|DYNAMICS
6 //
7 // BSD License: PfemFluidDynamicsApplication/license.txt
8 //
9 // Main authors: Massimiliano Zecchetto
10 // Collaborators:
11 //
12 //-------------------------------------------------------------
13 //
14 
15 #if !defined(SET_MESH_VELOCITY_FOR_THERMAL_COUPLING_PROCESS)
16 #define SET_MESH_VELOCITY_FOR_THERMAL_COUPLING_PROCESS
17 
18 #include "includes/define.h"
19 #include "includes/model_part.h"
20 #include "processes/process.h"
21 
22 namespace Kratos
23 {
24 
32 
33 public:
34 
36 
39 
42 
43  void operator()() {
44  Execute();
45  }
46 
47  void Execute() override {
48  KRATOS_TRY;
49 
50  this->SetMeshVelocity(rModelPart);
51 
52  KRATOS_CATCH("");
53  }
54 
55  void ExecuteInitialize() override {}
56 
57  void ExecuteInitializeSolutionStep() override {}
58 
59 protected:
60 
62 
63 private:
64 
65  void SetMeshVelocity(ModelPart& rModelPart) const {
66  const auto& it_node_begin = rModelPart.NodesBegin();
67  #pragma omp parallel for
68  for (int i = 0; i < static_cast<int>(rModelPart.Nodes().size()); i++) {
69  auto it_node = it_node_begin + i;
70  noalias(it_node->FastGetSolutionStepValue(MESH_VELOCITY)) = it_node->FastGetSolutionStepValue(VELOCITY);
71  }
72  }
73 
74 }; // Class SetMeshVelocityForThermalCouplingProcess
75 
77 inline std::istream &operator>>(std::istream &rIStream,
79 
81 inline std::ostream &operator<<(std::ostream &rOStream,
83  rThis.PrintInfo(rOStream);
84  rOStream << std::endl;
85  rThis.PrintData(rOStream);
86 
87  return rOStream;
88 }
89 
90 } // namespace Kratos.
91 
92 #endif /* SET_MESH_VELOCITY_FOR_THERMAL_COUPLING_PROCESS defined */
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
The base class for all processes in Kratos.
Definition: process.h:49
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: process.h:204
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: process.h:210
This method sets the MESH_VELOCITY equal to the nodal VELOCITY.
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:31
KRATOS_CLASS_POINTER_DEFINITION(SetMeshVelocityForThermalCouplingProcess)
ModelPart & rModelPart
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:61
void ExecuteInitialize() override
This function is designed for being called at the beginning of the computations right after reading t...
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:55
~SetMeshVelocityForThermalCouplingProcess() override
Destructor.
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:41
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:47
void ExecuteInitializeSolutionStep() override
This function will be executed at every time step BEFORE performing the solve phase.
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:57
void operator()()
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:43
SetMeshVelocityForThermalCouplingProcess(ModelPart &model_part)
Constructor.
Definition: set_mesh_velocity_for_thermal_coupling_process.hpp:38
#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
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
model_part
Definition: face_heat.py:14
integer i
Definition: TensorModule.f:17