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.
integration_values_extrapolation_to_nodes_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 #if !defined(KRATOS_INTEGRATION_VALUES_EXTRAPOLATION_TO_NODES_PROCESS )
14 #define KRATOS_INTEGRATION_VALUES_EXTRAPOLATION_TO_NODES_PROCESS
15 
16 // System includes
17 #include <unordered_map>
18 
19 // External includes
20 
21 // Project includes
22 #include "processes/process.h"
23 #include "containers/model.h"
24 #include "includes/key_hash.h"
26 
27 namespace Kratos
28 {
31 
35 
39 
43 
47 
57 class KRATOS_API(KRATOS_CORE) IntegrationValuesExtrapolationToNodesProcess
58  : public Process
59 {
60 public:
63 
64  // Node type definition
65  typedef Node NodeType;
66 
69 
71  typedef std::size_t SizeType;
72 
74  typedef std::size_t IndexType;
75 
78 
82 
86 
93  Model& rModel,
94  Parameters ThisParameters = Parameters(R"({})")
95  );
96 
103  ModelPart& rMainModelPart,
104  Parameters ThisParameters = Parameters(R"({})")
105  );
106 
109 
113 
114  void operator()()
115  {
116  Execute();
117  }
118 
122 
126  void Execute() override;
127 
131  void ExecuteBeforeSolutionLoop() override;
132 
136  void ExecuteFinalizeSolutionStep() override;
137 
141  void ExecuteFinalize() override;
142 
146  const Parameters GetDefaultParameters() const override;
147 
151 
155 
159 
160  /************************************ GET INFO *************************************/
161  /***********************************************************************************/
162 
163  std::string Info() const override
164  {
165  return "IntegrationValuesExtrapolationToNodesProcess";
166  }
167 
168  /************************************ PRINT INFO ***********************************/
169  /***********************************************************************************/
170 
171  void PrintInfo(std::ostream& rOStream) const override
172  {
173  rOStream << Info();
174  }
175 
179 
181 
182 protected:
183 
186 
190 
194 
198 
202 
206 
210 
212 
213 private:
216 
220 
221  ModelPart& mrModelPart;
222 
223  bool mExtrapolateNonHistorical;
224  bool mAreaAverage;
225 
226  std::vector<const Variable<double>*> mDoubleVariable;
227  std::vector<const Variable<array_1d<double, 3>>*> mArrayVariable;
228  std::vector<const Variable<Vector>*> mVectorVariable;
229  std::vector<const Variable<Matrix>*> mMatrixVariable;
230 
231  std::unordered_map<const Variable<Vector>*, SizeType, pVariableHasher, pVariableComparator> mSizeVectors;
232  std::unordered_map<const Variable<Matrix>*, std::pair<SizeType, SizeType>, pVariableHasher, pVariableComparator> mSizeMatrixes;
233 
234  const Variable<double>* mpAverageVariable;
235 
236  SizeType mEchoLevel;
237 
241 
245 
249  void InitializeMaps();
250 
254  void InitializeVariables();
255 
259 
263 
267 
269 
270 }; // Class IntegrationValuesExtrapolationToNodesProcess
271 
273 
276 
277 
281 
282 /****************************** INPUT STREAM FUNCTION ******************************/
283 /***********************************************************************************/
284 
285 inline std::istream& operator >> (std::istream& rIStream,
287 
288 /***************************** OUTPUT STREAM FUNCTION ******************************/
289 /***********************************************************************************/
290 
291 inline std::ostream& operator << (std::ostream& rOStream,
293 {
294  return rOStream;
295 }
296 
298 
299 } // namespace Kratos.
300 
301 #endif // KRATOS_INTEGRATION_VALUES_EXTRAPOLATION_TO_NODES_PROCESS defined
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: periodic_interface_process.hpp:55
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Geometry base class.
Definition: geometry.h:71
This process extrapolates vales from the integration points to the nodes.
Definition: integration_values_extrapolation_to_nodes_process.h:59
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: integration_values_extrapolation_to_nodes_process.h:171
std::size_t IndexType
Defining the index type.
Definition: integration_values_extrapolation_to_nodes_process.h:74
std::string Info() const override
Turn back information as a string.
Definition: integration_values_extrapolation_to_nodes_process.h:163
Node NodeType
Definition: integration_values_extrapolation_to_nodes_process.h:65
KRATOS_CLASS_POINTER_DEFINITION(IntegrationValuesExtrapolationToNodesProcess)
Pointer definition of IntegrationValuesExtrapolationToNodesProcess.
std::size_t SizeType
Defining the size type.
Definition: integration_values_extrapolation_to_nodes_process.h:71
Geometry< NodeType > GeometryType
Geometry type definition.
Definition: integration_values_extrapolation_to_nodes_process.h:68
void operator()()
Definition: integration_values_extrapolation_to_nodes_process.h:114
~IntegrationValuesExtrapolationToNodesProcess() override=default
Destructor.
This class aims to manage different model parts across multi-physics simulations.
Definition: model.h:60
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
The base class for all processes in Kratos.
Definition: process.h:49
void InitializeVariables(TContainerType &rContainer, const Variable< TDataType > &rOutputVariable, const Variable< TDataType > &rReferenceVariable)
Definition: temporal_method_utilities.h:36
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
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
This is a key comparer between two variables pointers.
Definition: key_hash.h:196
This is a hasher for variables pointers.
Definition: key_hash.h:183