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.
explicit_runge_kutta_4_eulerian_convdiff_strategy.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: Riccardo Tosi
11 //
12 //
13 
14 #if !defined(KRATOS_EXPLICIT_RUNGE_KUTTA_4_EULERIAN_CONVDIFF_STRATEGY)
15 #define KRATOS_EXPLICIT_RUNGE_KUTTA_4_EULERIAN_CONVDIFF_STRATEGY
16 
17 /* System includes */
18 
19 /* External includes */
20 
21 /* Project includes */
24 
25 namespace Kratos
26 {
27 
30 
34 
38 
42 
46 
55 template <class TSparseSpace, class TDenseSpace>
57 {
58 public:
61 
62  // The base class definition
64 
65  // The explicit builder and solver definition
67 
70 
73 
77 
84  ModelPart &rModelPart,
85  Parameters ThisParameters)
86  : BaseType(rModelPart, ThisParameters)
87  {
88  }
89 
97  ModelPart &rModelPart,
98  typename ExplicitBuilderType::Pointer pExplicitBuilder,
99  bool MoveMeshFlag = false,
100  int RebuildLevel = 0)
101  : BaseType(rModelPart, pExplicitBuilder, MoveMeshFlag, RebuildLevel)
102  {
103  }
104 
111  ModelPart &rModelPart,
112  bool MoveMeshFlag = false,
113  int RebuildLevel = 0)
114  : BaseType(rModelPart, MoveMeshFlag, RebuildLevel)
115  {
116  }
117 
121 
125 
129 
130 
134 
140  void Initialize() override
141  {
142  KRATOS_TRY;
143 
144 
145  auto& r_model_part = BaseType::GetModelPart();
146  auto& r_process_info = r_model_part.GetProcessInfo();
147  ConvectionDiffusionSettings::Pointer p_settings = r_process_info[CONVECTION_DIFFUSION_SETTINGS];
148  auto& r_settings = *p_settings;
149 
150  // Call the base method
152  // If required, initialize the OSS projection variables
153  if (r_process_info[OSS_SWITCH]) {
154  for (auto& r_node : r_model_part.GetCommunicator().LocalMesh().Nodes()) {
155  r_node.SetValue(r_settings.GetProjectionVariable(), 0.0);
156  }
157  }
158 
159  KRATOS_CATCH("");
160  }
161 
165 
167  std::string Info() const override
168  {
169  return "ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion";
170  }
171 
173  void PrintInfo(std::ostream &rOStream) const override
174  {
175  rOStream << Info();
176  }
177 
179  void PrintData(std::ostream &rOStream) const override
180  {
181  rOStream << Info();
182  }
183 
185 
186 protected:
189 
190 
194 
195 
199 
203 
209  {
210  KRATOS_TRY;
211 
213 
214  auto& r_model_part = BaseType::GetModelPart();
215  auto& r_process_info = r_model_part.GetProcessInfo();
216 
217  // execute OSS step, if needed
218  if (r_process_info[OSS_SWITCH] == 1) {
219  CalculateOSSNodalProjections();
220  }
221 
222  KRATOS_CATCH("");
223  };
224 
229  virtual void InitializeRungeKuttaLastSubStep() override
230  {
231  KRATOS_TRY;
232 
234  auto& r_model_part = BaseType::GetModelPart();
235  auto& r_process_info = r_model_part.GetProcessInfo();
236 
237  // execute OSS step, if needed
238  if (r_process_info[OSS_SWITCH] == 1) {
239  CalculateOSSNodalProjections();
240  }
241 
242  KRATOS_CATCH("");
243  };
244 
249  virtual bool SolveSolutionStep() override
250  {
251  KRATOS_TRY;
252 
254  auto& r_model_part = BaseType::GetModelPart();
255  auto& r_process_info = r_model_part.GetProcessInfo();
256  if (r_process_info.GetValue(RUNGE_KUTTA_STEP) == 4) {
257  // execute OSS step, if needed
258  if (r_process_info[OSS_SWITCH] == 1) {
259  CalculateOSSNodalProjections();
260  }
261  }
262 
263  return true;
264 
265  KRATOS_CATCH("");
266  };
267 
271 
272 
276 
277 
281 
282 
284 private:
287 
288 
292 
293 
297 
298 
302 
306  virtual void CalculateOSSNodalProjections()
307  {
308  KRATOS_TRY;
309 
310  // Get the required data from the explicit builder and solver
311  const auto p_explicit_bs = BaseType::pGetExplicitBuilder();
312  const auto& r_lumped_mass_vector = p_explicit_bs->GetLumpedMassMatrixVector();
313 
314  // Get model part data
315  auto& r_model_part = BaseType::GetModelPart();
316  auto& r_process_info = r_model_part.GetProcessInfo();
317 
318  ConvectionDiffusionSettings::Pointer p_settings = r_process_info[CONVECTION_DIFFUSION_SETTINGS];
319  auto& r_settings = *p_settings;
320 
321  // Initialize the projection value
322  block_for_each(r_model_part.Nodes(), [&](Node& rNode){
323  rNode.GetValue(r_settings.GetProjectionVariable()) = 0.0;
324  });
325 
326  // Calculate the unknown projection
327  double unknown_proj;
328  block_for_each(r_model_part.Elements(), [&](ModelPart::ElementType& rElement){
329  rElement.Calculate(r_settings.GetProjectionVariable(), unknown_proj, r_process_info);
330  });
331  IndexPartition<int>(r_model_part.NumberOfNodes()).for_each(
332  [&](int i_node){
333  auto it_node = r_model_part.NodesBegin() + i_node;
334  const double mass = r_lumped_mass_vector(i_node);
335  it_node->FastGetSolutionStepValue(r_settings.GetProjectionVariable()) = it_node->GetValue(r_settings.GetProjectionVariable()) / mass;
336  }
337  );
338 
339  KRATOS_CATCH("");
340  };
341 
345 
346 
350 
351 
355 
356 
358 }; /* Class ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion */
359 
361 
364 
366 
367 } /* namespace Kratos.*/
368 
369 #endif /* KRATOS_EXPLICIT_RUNGE_KUTTA_4_EULERIAN_CONVDIFF_STRATEGY defined */
Base class for all Elements.
Definition: element.h:60
ExplicitBuilderPointerType pGetExplicitBuilder()
Operations to get the pointer to the explicit builder and solver.
Definition: explicit_solving_strategy.h:286
bool SolveSolutionStep() override
Solves the current step. The function always return true as convergence is not checked in the explici...
Definition: explicit_solving_strategy.h:230
This strategy adds the orthogonal subgrid projections computation to the base explicit runge kutta 4 ...
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:57
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:173
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(ModelPart &rModelPart, Parameters ThisParameters)
Default constructor. (with parameters)
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:83
TDenseSpace::VectorType LocalSystemVectorType
The local vector definition.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:69
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(ModelPart &rModelPart, bool MoveMeshFlag=false, int RebuildLevel=0)
Default constructor.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:110
void Initialize() override
Initialization of variables.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:140
virtual bool SolveSolutionStep() override
Finalize the Runge Kutta explicit solver step and calculate the orthogonal subscale projections if re...
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:249
ExplicitSolvingStrategyRungeKutta4< TSparseSpace, TDenseSpace > BaseType
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:63
std::string Info() const override
Turn back information as a string.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:167
BaseType::ExplicitBuilderType ExplicitBuilderType
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:66
virtual void InitializeRungeKuttaIntermediateSubStep() override
Initialize the Runge-Kutta substep.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:208
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(ModelPart &rModelPart, typename ExplicitBuilderType::Pointer pExplicitBuilder, bool MoveMeshFlag=false, int RebuildLevel=0)
Default constructor.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:96
virtual void InitializeRungeKuttaLastSubStep() override
Initialize the Runge-Kutta substep.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:229
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: explicit_runge_kutta_4_eulerian_convdiff_strategy.h:179
ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion(const ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion &Other)=delete
KRATOS_CLASS_POINTER_DEFINITION(ExplicitSolvingStrategyRungeKutta4ConvectionDiffusion)
Family of explicit Runge-Kutta schemes.
Definition: explicit_solving_strategy_runge_kutta.h:66
virtual void InitializeRungeKuttaLastSubStep()
Initialize the Runge-Kutta last substep This method is intended to implement all the operations requi...
Definition: explicit_solving_strategy_runge_kutta.h:333
typename BaseType::ExplicitBuilderType ExplicitBuilderType
The explicit builder and solver definition.
Definition: explicit_solving_strategy_runge_kutta.h:81
void Initialize() override
Initialization of member variables and prior operations.
Definition: explicit_solving_strategy_runge_kutta.h:195
virtual void InitializeRungeKuttaIntermediateSubStep()
Initialize the Runge-Kutta intermediate substep This method is intended to implement all the operatio...
Definition: explicit_solving_strategy_runge_kutta.h:321
This class is useful for index iteration over containers.
Definition: parallel_utilities.h:451
void for_each(TUnaryFunction &&f)
Definition: parallel_utilities.h:514
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
Solving strategy base class This is the base class from which we will derive all the strategies (impl...
Definition: solving_strategy.h:64
ModelPart & GetModelPart()
Operations to get the pointer to the model.
Definition: solving_strategy.h:350
bool MoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:290
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
void block_for_each(TIterator itBegin, TIterator itEnd, TFunction &&rFunction)
Execute a functor on all items of a range in parallel.
Definition: parallel_utilities.h:299