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.
fractional_step_settings.h
Go to the documentation of this file.
1 #ifndef KRATOS_FRACTIONAL_STEP_SETTINGS_H
2 #define KRATOS_FRACTIONAL_STEP_SETTINGS_H
3 
4 // System includes
5 #include <string>
6 #include <iostream>
7 
8 
9 // External includes
10 
11 
12 // Project includes
13 #include "includes/define.h"
23 #include "processes/process.h"
24 
25 // Application includes
27 #include "custom_utilities/solver_settings.h"
28 
29 
30 namespace Kratos
31 {
34 
37 
41 
45 
49 
53 
55 template< class TSparseSpace,
56  class TDenseSpace,
57  class TLinearSolver
58  >
59 class FractionalStepSettings: public SolverSettings<TSparseSpace,TDenseSpace,TLinearSolver>
60 {
61 public:
64 
67 
72 
75 
79 
82  const unsigned int ThisDomainSize,
83  const unsigned int ThisTimeOrder,
84  const bool UseSlip,
85  const bool MoveMeshFlag,
86  const bool ReformDofSet):
87  BaseType(rModelPart,ThisDomainSize,ThisTimeOrder,UseSlip,MoveMeshFlag,ReformDofSet)
88  {}
89 
92 
96 
97 
101 
105 
106  void SetStrategy(StrategyLabel const& rStrategyLabel,
107  typename TLinearSolver::Pointer pLinearSolver,
108  const double Tolerance,
109  const unsigned int MaxIter) override
110  {
111  KRATOS_TRY;
112 
113  // pointer types for solution strategy construcion
114  typedef typename Scheme< TSparseSpace, TDenseSpace >::Pointer SchemePointerType;
115 // typedef typename ConvergenceCriteria< TSparseSpace, TDenseSpace >::Pointer ConvergenceCriteriaPointerType;
116  typedef typename BuilderAndSolver<TSparseSpace, TDenseSpace, TLinearSolver>::Pointer BuilderSolverTypePointer;
117 
118  // Default, fixed flags
119  bool CalculateReactions = false;
120  bool CalculateNormDxFlag = true;
121 
122  ModelPart& rModelPart = BaseType::GetModelPart();
123  bool UseSlip = BaseType::UseSlipConditions();
124  // Modification of the DofSet is managed by the fractional step strategy, not the auxiliary velocity and pressure strategies.
125  bool ReformDofSet = false; //BaseType::GetReformDofSet();
126  unsigned int EchoLevel = BaseType::GetEchoLevel();
127  unsigned int StrategyEchoLevel = (EchoLevel > 0) ? (EchoLevel-1) : 0;
128 
129  if ( rStrategyLabel == BaseType::Velocity )
130  {
131  // Velocity Builder and Solver
132  BuilderSolverTypePointer pBuildAndSolver = BuilderSolverTypePointer(new ResidualBasedBlockBuilderAndSolver<TSparseSpace, TDenseSpace, TLinearSolver >
133  (pLinearSolver));
134 
135  SchemePointerType pScheme;
136  //initializing fractional velocity solution step
137  if (UseSlip)
138  {
139  unsigned int DomainSize = BaseType::GetDomainSize();
140  SchemePointerType Temp = SchemePointerType(new ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace > (DomainSize,DomainSize));
141  pScheme.swap(Temp);
142  }
143  else
144  {
145  SchemePointerType Temp = SchemePointerType(new ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace > ());
146  pScheme.swap(Temp);
147  }
148 
149  // Strategy
151  rModelPart,
152  pScheme,
153  pBuildAndSolver,
154  CalculateReactions,
155  ReformDofSet,
156  CalculateNormDxFlag));
157 
158  }
159  else if ( rStrategyLabel == BaseType::Pressure )
160  {
161  // Pressure Builder and Solver
162  BuilderSolverTypePointer pBuildAndSolver = BuilderSolverTypePointer(new ResidualBasedBlockBuilderAndSolver<TSparseSpace, TDenseSpace, TLinearSolver > (pLinearSolver));
163  SchemePointerType pScheme = SchemePointerType(new ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace > ());
164 
165  // Strategy
167  rModelPart,
168  pScheme,
169  pBuildAndSolver,
170  CalculateReactions,
171  ReformDofSet,
172  CalculateNormDxFlag));
173  }
174  else
175  {
176  KRATOS_THROW_ERROR(std::runtime_error,"Error in FractionalStepSettings: Unknown strategy label.","");
177  }
178 
179  BaseType::mTolerances[rStrategyLabel] = Tolerance;
180 
181  BaseType::mMaxIter[rStrategyLabel] = MaxIter;
182 
183  BaseType::mStrategies[rStrategyLabel]->SetEchoLevel(StrategyEchoLevel);
184 
185  KRATOS_CATCH("");
186  }
187 
188  void SetTurbulenceModel(TurbulenceModelLabel const& rTurbulenceModel,
189  typename TLinearSolver::Pointer pLinearSolver,
190  const double Tolerance,
191  const unsigned int MaxIter) override
192  {
193  KRATOS_TRY;
194 
196 
197  ModelPart& rModelPart = BaseType::GetModelPart();
198  unsigned int DomainSize = BaseType::GetDomainSize();
199  unsigned int TimeOrder = BaseType::GetTimeOrder();
200  bool ReformDofSet = BaseType::GetReformDofSet();
201 
202  if (rTurbulenceModel == BaseType::SpalartAllmaras)
203  {
205  (rModelPart,pLinearSolver,DomainSize,Tolerance,MaxIter,ReformDofSet,TimeOrder));
206  }
207  else
208  {
209  KRATOS_THROW_ERROR(std::runtime_error,"Error in FractionalStepSettings: Unknown turbulence model label.","");
210  }
211 
212  KRATOS_CATCH("");
213  }
214 
215  void SetTurbulenceModel(ProcessPointerType pTurbulenceModel) override
216  {
217  BaseType::SetTurbulenceModel(pTurbulenceModel);
218  }
219 
223 
227 
229  std::string Info() const override
230  {
231  std::stringstream buffer;
232  buffer << "FractionalStepSettings" ;
233  return buffer.str();
234  }
235 
237  void PrintInfo(std::ostream& rOStream) const override {rOStream << "FractionalStepSettings";}
238 
240  void PrintData(std::ostream& rOStream) const override {}
241 
242 
246 
247 
249 
250 protected:
253 
254 
258 
259 
263 
264 
268 
269 
273 
274 
278 
279 
283 
284 
286 
287 private:
290 
291 
295 
296 
300 
301 
305 
306 
310 
311 
315 
316 
320 
323 
325  FractionalStepSettings& operator=(FractionalStepSettings const& rOther){}
326 
328  FractionalStepSettings(FractionalStepSettings const& rOther){}
329 
330 
332 
333 }; // Class FractionalStepSettings
334 
336 
339 
340 
344 
345 
347 template< class TDenseSpace, class TSparseSpace, class TLinearSolver >
348 inline std::istream& operator >> (std::istream& rIStream,
350 {
351  return rIStream;
352 }
353 
355 template< class TDenseSpace, class TSparseSpace, class TLinearSolver >
356 inline std::ostream& operator << (std::ostream& rOStream,
358 {
359  rThis.PrintInfo(rOStream);
360  rOStream << std::endl;
361  rThis.PrintData(rOStream);
362 
363  return rOStream;
364 }
366 
368 
369 } // namespace Kratos.
370 
371 #endif // KRATOS_FRACTIONAL_STEP_SETTINGS_H
Current class provides an implementation for the base builder and solving operations.
Definition: builder_and_solver.h:64
Helper class to define solution strategies for FS_Strategy.
Definition: fractional_step_settings.h:60
BaseType::StrategyLabel StrategyLabel
Definition: fractional_step_settings.h:73
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: fractional_step_settings.h:237
BaseType::TurbulenceModelLabel TurbulenceModelLabel
Definition: fractional_step_settings.h:74
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: fractional_step_settings.h:240
void SetStrategy(StrategyLabel const &rStrategyLabel, typename TLinearSolver::Pointer pLinearSolver, const double Tolerance, const unsigned int MaxIter) override
Definition: fractional_step_settings.h:106
KRATOS_CLASS_POINTER_DEFINITION(FractionalStepSettings)
Pointer definition of FractionalStepSettings.
BaseType::StrategyType StrategyType
Definition: fractional_step_settings.h:69
void SetTurbulenceModel(TurbulenceModelLabel const &rTurbulenceModel, typename TLinearSolver::Pointer pLinearSolver, const double Tolerance, const unsigned int MaxIter) override
Definition: fractional_step_settings.h:188
BaseType::StrategyPointerType StrategyPointerType
Definition: fractional_step_settings.h:70
std::string Info() const override
Turn back information as a string.
Definition: fractional_step_settings.h:229
void SetTurbulenceModel(ProcessPointerType pTurbulenceModel) override
Definition: fractional_step_settings.h:215
BaseType::ProcessPointerType ProcessPointerType
Definition: fractional_step_settings.h:71
FractionalStepSettings(ModelPart &rModelPart, const unsigned int ThisDomainSize, const unsigned int ThisTimeOrder, const bool UseSlip, const bool MoveMeshFlag, const bool ReformDofSet)
Constructor.
Definition: fractional_step_settings.h:81
SolverSettings< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: fractional_step_settings.h:68
~FractionalStepSettings() override
Destructor.
Definition: fractional_step_settings.h:91
Implicit solving strategy base class This is the base class from which we will derive all the implici...
Definition: implicit_solving_strategy.h:61
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Current class provides an implementation for standard builder and solving operations.
Definition: residualbased_block_builder_and_solver.h:82
This class provides the implementation of a static scheme.
Definition: residualbased_incrementalupdate_static_scheme.h:57
Scheme for the solution of problems involving a slip condition.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:70
This is a very simple strategy to solve linearly the problem.
Definition: residualbased_linear_strategy.h:64
This class provides the implementation of the basic tasks that are needed by the solution strategy.
Definition: scheme.h:56
Helper class to define solution strategies for FS_Strategy.
Definition: solver_settings.h:53
StrategyLabel
Definition: solver_settings.h:66
@ Velocity
Definition: solver_settings.h:66
@ Pressure
Definition: solver_settings.h:66
virtual unsigned int GetDomainSize() const
Definition: solver_settings.h:146
ModelPart & GetModelPart()
Definition: solver_settings.h:285
virtual void SetTurbulenceModel(TurbulenceModelLabel const &rTurbulenceModel, typename TLinearSolver::Pointer pLinearSolver, const double Tolerance, const unsigned int MaxIter)
Definition: solver_settings.h:124
bool GetReformDofSet()
Definition: solver_settings.h:226
StrategyType::Pointer StrategyPointerType
Definition: solver_settings.h:62
virtual bool UseSlipConditions() const
Definition: solver_settings.h:156
std::map< StrategyLabel, double > mTolerances
Definition: solver_settings.h:292
std::map< StrategyLabel, StrategyPointerType > mStrategies
Definition: solver_settings.h:290
std::map< StrategyLabel, unsigned int > mMaxIter
Definition: solver_settings.h:294
Process::Pointer ProcessPointerType
Definition: solver_settings.h:63
virtual unsigned int GetEchoLevel()
Definition: solver_settings.h:221
TurbulenceModelLabel
Definition: solver_settings.h:68
@ SpalartAllmaras
Definition: solver_settings.h:68
bool mHaveTurbulenceModel
Definition: solver_settings.h:298
ProcessPointerType mpTurbulenceModel
Definition: solver_settings.h:296
virtual unsigned int GetTimeOrder() const
Definition: solver_settings.h:151
An impelementation of the Spalart-Allmaras turbulence model for incompressible flows.
Definition: spalart_allmaras_turbulence_model.h:77
#define KRATOS_THROW_ERROR(ExceptionType, ErrorMessage, MoreInfo)
Definition: define.h:77
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
static int EchoLevel
Definition: co_sim_EMPIRE_API.h:42
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
Temp
Definition: PecletTest.py:105