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_for_chimera.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 //
8 // License: BSD License
9 // Kratos default license: kratos/license.txt
10 //
11 // Authors: Aditya Ghantasala, https://github.com/adityaghantasala
12 // Navaneeth K Narayanan
13 // Rishith Ellath Meethal
14 //
15 
16 
17 #ifndef KRATOS_FRACTIONAL_STEP_SETTINGS_FOR_CHIMERA_H
18 #define KRATOS_FRACTIONAL_STEP_SETTINGS_FOR_CHIMERA_H
19 
20 // System includes
21 #include <string>
22 #include <iostream>
23 
24 
25 // External includes
26 
27 
28 // Project includes
29 #include "includes/define.h"
30 
38 #include "processes/process.h"
40 #include "containers/model.h"
41 
42 // Application includes
43 #include "custom_utilities/solver_settings.h"
46 
47 
48 namespace Kratos
49 {
52 
55 
59 
63 
67 
71 
73 template< class TSparseSpace,
74  class TDenseSpace,
75  class TLinearSolver
76  >
77 class FractionalStepSettingsForChimera: public SolverSettings<TSparseSpace,TDenseSpace,TLinearSolver>
78 {
79 public:
82 
85 
92  typedef typename ResidualBasedBlockBuilderAndSolverType::Pointer ResidualBasedBlockBuilderAndSolverPointerType;
93 
94 
98 
101  const std::size_t ThisDomainSize,
102  const std::size_t ThisTimeOrder,
103  const bool use_slip,
104  const bool MoveMeshFlag,
105  const bool reform_dof_set):
106  BaseType(r_model_part,ThisDomainSize,ThisTimeOrder,use_slip,MoveMeshFlag,reform_dof_set)
107  {}
108 
111 
112 
115 
118 
121 
122 
126 
127 
131 
135 
136  void SetStrategy(StrategyLabel const& rStrategyLabel,
137  typename TLinearSolver::Pointer pLinearSolver,
138  const double Tolerance,
139  const unsigned int MaxIter) override
140  {
141  KRATOS_TRY;
142 
143  // pointer types for solution strategy construction
144  typedef typename Scheme< TSparseSpace, TDenseSpace >::Pointer SchemePointerType;
145  typedef ResidualBasedBlockBuilderAndSolverWithConstraintsForChimera<TSparseSpace, TDenseSpace, TLinearSolver > ResidualBasedBlockBuilderAndSolverWithConstraintsForChimeraType;
146 
147  // Default, fixed flags
148  bool calculate_reactions = false;
149  bool calculate_norm_dx_flag = true;
150 
151  ModelPart& r_model_part = BaseType::GetModelPart();
152 
153  bool use_slip = BaseType::UseSlipConditions();
154  // Modification of the DofSet is managed by the fractional step strategy, not the auxiliary velocity and pressure strategies.
155  bool reform_dof_set = BaseType::GetReformDofSet();
156  std::size_t echo_level = BaseType::GetEchoLevel();
157  std::size_t strategy_echo_level = (echo_level > 0) ? (echo_level-1) : 0;
158 
159  if ( rStrategyLabel == BaseType::Velocity )
160  {
161  std::string fs_vel_mp_name = r_model_part.Name()+"fs_velocity_model_part";
162  ModelPart &r_fs_velocity_model_part = r_model_part.CreateSubModelPart(fs_vel_mp_name);
163  FastTransferBetweenModelPartsProcess(r_fs_velocity_model_part, r_model_part).Execute();
164  // Velocity Builder and Solver
165  ResidualBasedBlockBuilderAndSolverPointerType p_build_and_solver = Kratos::make_shared<ResidualBasedBlockBuilderAndSolverWithConstraintsForChimeraType>(pLinearSolver);
166  p_build_and_solver->SetEchoLevel(strategy_echo_level);
167 
168  SchemePointerType p_scheme;
169  //initializing fractional velocity solution step
170  if (use_slip)
171  {
172  std::size_t domain_size = BaseType::GetDomainSize();
173  SchemePointerType p_temp = Kratos::make_shared<ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace>> (domain_size,domain_size);
174  p_scheme.swap(p_temp);
175  }
176  else
177  {
178  SchemePointerType p_temp = Kratos::make_shared<ResidualBasedIncrementalUpdateStaticScheme<TSparseSpace, TDenseSpace >>();
179  p_scheme.swap(p_temp);
180  }
181 
182  // Strategy
183  BaseType::mStrategies[rStrategyLabel] = Kratos::make_shared< ResidualBasedLinearStrategy<TSparseSpace, TDenseSpace, TLinearSolver >>
184  (r_fs_velocity_model_part,
185  p_scheme,
186  p_build_and_solver,
187  calculate_reactions,
188  reform_dof_set,
189  calculate_norm_dx_flag);
190  }
191  else if ( rStrategyLabel == BaseType::Pressure )
192  {
193  std::string fs_pressure_mp_name = r_model_part.Name()+"fs_pressure_model_part";
194  ModelPart &r_fs_pressure_model_part = r_model_part.CreateSubModelPart(fs_pressure_mp_name);
195  FastTransferBetweenModelPartsProcess(r_fs_pressure_model_part, r_model_part).Execute();
196  // Pressure Builder and Solver
197  ResidualBasedBlockBuilderAndSolverPointerType p_build_and_solver = Kratos::make_shared<ResidualBasedBlockBuilderAndSolverWithConstraintsForChimeraType>(pLinearSolver);
198  SchemePointerType p_scheme = Kratos::make_shared<ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >> ();
199  p_build_and_solver->SetEchoLevel(strategy_echo_level);
200 
201  // Strategy
202  BaseType::mStrategies[rStrategyLabel] = Kratos::make_shared<ResidualBasedLinearStrategy<TSparseSpace, TDenseSpace, TLinearSolver >>(
203  r_fs_pressure_model_part,
204  p_scheme,
205  p_build_and_solver,
206  calculate_reactions,
207  reform_dof_set,
208  calculate_norm_dx_flag);
209  }
210  else
211  {
212  KRATOS_ERROR << "Error in FractionalStepSettingsForChimera: Unknown strategy label."<<std::endl;
213  }
214 
215  BaseType::mTolerances[rStrategyLabel] = Tolerance;
216 
217  BaseType::mMaxIter[rStrategyLabel] = MaxIter;
218 
219  BaseType::mStrategies[rStrategyLabel]->SetEchoLevel(strategy_echo_level);
220 
221  KRATOS_CATCH("");
222  }
223 
224 
225  bool FindStrategy(StrategyLabel const& rStrategyLabel,
226  StrategyPointerType& pThisStrategy) override
227  {
228  pThisStrategy = BaseType::mStrategies[rStrategyLabel];
229  if(pThisStrategy != nullptr)
230  return true;
231  return false;
232  }
233 
237 
241 
243  std::string Info() const override
244  {
245  std::stringstream buffer;
246  buffer << "FractionalStepSettingsForChimera" ;
247  return buffer.str();
248  }
249 
251  void PrintInfo(std::ostream& rOStream) const override {rOStream << "FractionalStepSettingsForChimera";}
252 
254  void PrintData(std::ostream& rOStream) const override {}
255 
256 
260 
261 
263 
264 
265 }; // Class FractionalStepSettingsForChimera
266 
268 
271 
272 
276 
277 
279 template< class TDenseSpace, class TSparseSpace, class TLinearSolver >
280 inline std::istream& operator >> (std::istream& rIStream,
282 {
283  return rIStream;
284 }
285 
287 template< class TDenseSpace, class TSparseSpace, class TLinearSolver >
288 inline std::ostream& operator << (std::ostream& rOStream,
290 {
291  rThis.PrintInfo(rOStream);
292  rOStream << std::endl;
293  rThis.PrintData(rOStream);
294 
295  return rOStream;
296 }
298 
300 
301 } // namespace Kratos.
302 
303 #endif // KRATOS_FRACTIONAL_STEP_SETTINGS_FOR_CHIMERA_H
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
void Execute() override
Execute method is used to execute the FastTransferBetweenModelPartsProcess algorithms.
Definition: fast_transfer_between_model_parts_process.cpp:60
Helper class to define solution strategies for FS_Strategy.
Definition: fractional_step_settings_for_chimera.h:78
void SetStrategy(StrategyLabel const &rStrategyLabel, typename TLinearSolver::Pointer pLinearSolver, const double Tolerance, const unsigned int MaxIter) override
Definition: fractional_step_settings_for_chimera.h:136
std::string Info() const override
Turn back information as a string.
Definition: fractional_step_settings_for_chimera.h:243
FractionalStepSettingsForChimera(ModelPart &r_model_part, const std::size_t ThisDomainSize, const std::size_t ThisTimeOrder, const bool use_slip, const bool MoveMeshFlag, const bool reform_dof_set)
Constructor.
Definition: fractional_step_settings_for_chimera.h:100
BaseType::StrategyType StrategyType
Definition: fractional_step_settings_for_chimera.h:87
BaseType::StrategyPointerType StrategyPointerType
Definition: fractional_step_settings_for_chimera.h:88
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: fractional_step_settings_for_chimera.h:254
bool FindStrategy(StrategyLabel const &rStrategyLabel, StrategyPointerType &pThisStrategy) override
Definition: fractional_step_settings_for_chimera.h:225
BaseType::ProcessPointerType ProcessPointerType
Definition: fractional_step_settings_for_chimera.h:89
ResidualBasedBlockBuilderAndSolver< TSparseSpace, TDenseSpace, TLinearSolver > ResidualBasedBlockBuilderAndSolverType
Definition: fractional_step_settings_for_chimera.h:91
FractionalStepSettingsForChimera()=delete
Default constructor.
ResidualBasedBlockBuilderAndSolverType::Pointer ResidualBasedBlockBuilderAndSolverPointerType
Definition: fractional_step_settings_for_chimera.h:92
FractionalStepSettingsForChimera & operator=(FractionalStepSettingsForChimera const &rOther)=delete
Assignment operator.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: fractional_step_settings_for_chimera.h:251
FractionalStepSettingsForChimera(FractionalStepSettingsForChimera const &rOther)=delete
Copy constructor.
BaseType::StrategyLabel StrategyLabel
Definition: fractional_step_settings_for_chimera.h:90
~FractionalStepSettingsForChimera()=default
Destructor.
KRATOS_CLASS_POINTER_DEFINITION(FractionalStepSettingsForChimera)
Pointer definition of FractionalStepSettingsForChimera.
SolverSettings< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: fractional_step_settings_for_chimera.h:86
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
std::string & Name()
Definition: model_part.h:1811
ModelPart & CreateSubModelPart(std::string const &NewSubModelPartName)
Definition: model_part.cpp:2000
Current class provides an implementation for standard builder and solving operations.
Definition: residualbased_block_builder_and_solver.h:82
Current class provides an implementation for applying the chimera constraints that is enforcing conti...
Definition: residualbased_block_builder_and_solver_with_constraints_for_chimera.h:73
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
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
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
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
int domain_size
Definition: face_heat.py:4
echo_level
Definition: script.py:68