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_solving_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: Ruben Zorrilla
11 //
12 //
13 
14 #if !defined(KRATOS_EXPLICIT_SOLVING_STRATEGY)
15 #define KRATOS_EXPLICIT_SOLVING_STRATEGY
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "solving_strategy.h"
23 #include "includes/define.h"
24 #include "includes/model_part.h"
30 
31 namespace Kratos
32 {
33 
36 
40 
44 
48 
52 
56 template <class TSparseSpace, class TDenseSpace>
57 class ExplicitSolvingStrategy : public SolvingStrategy<TSparseSpace, TDenseSpace>
58 {
59 public:
62 
63  // The explicit builder and solver definition
65 
66  // The explicit builder and solver pointer definition
67  typedef typename ExplicitBuilderType::Pointer ExplicitBuilderPointerType;
68 
69  // The DOF type from the explicit builder and solver class
71 
74 
77 
80 
84 
89  {
90  }
91 
98  ModelPart &rModelPart,
99  Parameters ThisParameters)
100  : BaseType(rModelPart, ThisParameters)
101  {
102  mpExplicitBuilder = Kratos::make_unique<ExplicitBuilder<TSparseSpace, TDenseSpace>>();
103  }
104 
113  ModelPart &rModelPart,
114  typename ExplicitBuilderType::Pointer pExplicitBuilder,
115  bool MoveMeshFlag = false,
116  int RebuildLevel = 0)
117  : BaseType(rModelPart, MoveMeshFlag)
118  , mpExplicitBuilder(pExplicitBuilder)
119  {
120  SetRebuildLevel(RebuildLevel);
121  }
122 
130  ModelPart &rModelPart,
131  bool MoveMeshFlag = false,
132  int RebuildLevel = 0)
133  : BaseType(rModelPart, MoveMeshFlag)
134  {
135  SetRebuildLevel(RebuildLevel);
136  mpExplicitBuilder = Kratos::make_unique<ExplicitBuilder<TSparseSpace, TDenseSpace>>();
137  }
138 
142 
146 
152  typename BaseType::Pointer Create(
153  ModelPart& rModelPart,
154  Parameters ThisParameters
155  ) const override
156  {
157  return Kratos::make_shared<ClassType>(rModelPart, ThisParameters);
158  }
159 
163 
167 
171  void Initialize() override
172  {
173  // Initialize elements, conditions and constraints
175 
176  // Set the explicit DOFs rebuild level
177  if (mRebuildLevel != 0) {
178  mpExplicitBuilder->SetResetDofSetFlag(true);
179  }
180 
181  // If the mesh is updated at each step, we require to accordingly update the lumped mass at each step
183  mpExplicitBuilder->SetResetLumpedMassVectorFlag(true);
184  }
185 
186  // Call the explicit builder and solver initialize (Set up DOF set and lumped mass vector)
187  mpExplicitBuilder->Initialize(BaseType::GetModelPart());
188 
189  }
190 
194  void Clear() override
195  {
196  // This clears the DOF set and lumped mass vector
197  mpExplicitBuilder->Clear();
198  }
199 
204  void InitializeSolutionStep() override
205  {
206  // InitializeSolutionStep elements, conditions and constraints
208 
209  // Call the builder and solver initialize solution step
210  mpExplicitBuilder->InitializeSolutionStep(BaseType::GetModelPart());
211  }
212 
217  void FinalizeSolutionStep() override
218  {
219  // FinalizeSolutionStep elements, conditions and constraints
221 
222  // Call the builder and solver finalize solution step (the reactions are computed in here)
223  mpExplicitBuilder->FinalizeSolutionStep(BaseType::GetModelPart());
224  }
225 
230  bool SolveSolutionStep() override
231  {
232  // Call the initialize non-linear iteration
234 
235  // Apply constraints
236  if(BaseType::GetModelPart().MasterSlaveConstraints().size() != 0) {
237  mpExplicitBuilder->ApplyConstraints(BaseType::GetModelPart());
238  }
239 
240  // Solve the problem assuming that a lumped mass matrix is used
242 
243  // If required, update the mesh with the obtained solution
246  }
247 
248  // Call the finalize non-linear iteration
250 
251  return true;
252  }
253 
263  void SetRebuildLevel(int Level) override
264  {
265  mRebuildLevel = Level;
266  }
267 
277  int GetRebuildLevel() const override
278  {
279  return mRebuildLevel;
280  }
281 
287  {
288  return mpExplicitBuilder;
289  };
290 
296  {
297  KRATOS_ERROR_IF(mpExplicitBuilder == nullptr) << "Asking for builder and solver when it is empty" << std::endl;
298  return *mpExplicitBuilder;
299  };
300 
306  {
307  KRATOS_ERROR_IF(mpExplicitBuilder == nullptr) << "Asking for builder and solver when it is empty" << std::endl;
308  return *mpExplicitBuilder;
309  };
310 
315  double GetResidualNorm() override
316  {
317  // Get the required data from the explicit builder and solver
318  auto& r_explicit_bs = GetExplicitBuilder();
319  auto& r_dof_set = r_explicit_bs.GetDofSet();
320 
321  // Calculate the explicit residual
322  r_explicit_bs.BuildRHS(BaseType::GetModelPart());
323 
324  // Calculate the residual norm
325  double res_norm = 0.0;
326  res_norm = block_for_each<SumReduction<double>>(
327  r_dof_set,
328  [](DofType &rDof){return rDof.GetSolutionStepReactionValue();}
329  );
330 
331  return res_norm;
332  }
333 
339  {
340  Parameters default_parameters = Parameters(R"(
341  {
342  "explicit_solving_strategy" : "explicit_solving_strategy",
343  "rebuild_level" : 0,
344  "explicit_builder_settings" : {
345  "name": "explicit_builder"
346  }
347  })");
348 
349  // Getting base class default parameters
350  const Parameters base_default_parameters = BaseType::GetDefaultParameters();
351  default_parameters.RecursivelyAddMissingParameters(base_default_parameters);
352 
353  return default_parameters;
354  }
355 
360  static std::string Name()
361  {
362  return "explicit_solving_strategy";
363  }
364 
368 
370  std::string Info() const override
371  {
372  return "ExplicitSolvingStrategy";
373  }
374 
376 protected:
379 
380  // Settings for the rebuilding of the DOF set
382 
386 
390 
394 
401  {
402  KRATOS_ERROR << "Calling the base ExplicitSolvingStrategy SolveWithLumpedMassMatrix(). Implement the specific explicit scheme solution update in a derived class" << std::endl;
403  }
404 
410  virtual inline double GetDeltaTime()
411  {
412  return BaseType::GetModelPart().GetProcessInfo().GetValue(DELTA_TIME);
413  }
414 
419  void AssignSettings(const Parameters ThisParameters) override
420  {
421  // Add base strategy settings
422  BaseType::AssignSettings(ThisParameters);
423 
424  const bool rebuild_level = ThisParameters["rebuild_level"].GetInt();
425  SetRebuildLevel(rebuild_level);
426  }
427 
431 
435 
439 
440 private:
444 
445 
449 
450  ExplicitBuilderPointerType mpExplicitBuilder = nullptr;
451 
455 
459 
463 
467 
471 
473 }; /* Class NewExplicitSolvingStrategy */
474 
476 
479 
481 
482 } /* namespace Kratos.*/
483 
484 #endif /* KRATOS_EXPLICIT_SOLVING_STRATEGY defined */
TDataType & GetValue(const Variable< TDataType > &rThisVariable)
Gets the value associated with a given variable.
Definition: data_value_container.h:268
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
Current class provides an implementation for the base explicit builder and solving operations.
Definition: explicit_builder.h:68
Explicit solving strategy base class.
Definition: explicit_solving_strategy.h:58
SolvingStrategy< TSparseSpace, TDenseSpace > BaseType
The definition of the base class.
Definition: explicit_solving_strategy.h:73
int GetRebuildLevel() const override
This returns the build level.
Definition: explicit_solving_strategy.h:277
Parameters GetDefaultParameters() const override
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: explicit_solving_strategy.h:338
ExplicitBuilderType::Pointer ExplicitBuilderPointerType
Definition: explicit_solving_strategy.h:67
virtual double GetDeltaTime()
Get the Delta Time object This method returns the DELTA_TIME from the ProcessInfo container.
Definition: explicit_solving_strategy.h:410
ExplicitSolvingStrategy()
Default constructor. (empty)
Definition: explicit_solving_strategy.h:88
void InitializeSolutionStep() override
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: explicit_solving_strategy.h:204
KRATOS_CLASS_POINTER_DEFINITION(ExplicitSolvingStrategy)
ExplicitBuilder< TSparseSpace, TDenseSpace > ExplicitBuilderType
Definition: explicit_solving_strategy.h:64
std::string Info() const override
Turn back information as a string.
Definition: explicit_solving_strategy.h:370
BaseType::Pointer Create(ModelPart &rModelPart, Parameters ThisParameters) const override
Create method.
Definition: explicit_solving_strategy.h:152
void Clear() override
Clears the internal storage.
Definition: explicit_solving_strategy.h:194
double GetResidualNorm() override
Operations to get the residual norm.
Definition: explicit_solving_strategy.h:315
int mRebuildLevel
Definition: explicit_solving_strategy.h:381
ExplicitSolvingStrategy< TSparseSpace, TDenseSpace > ClassType
The definition of the current class.
Definition: explicit_solving_strategy.h:76
virtual void SolveWithLumpedMassMatrix()
Calculate the explicit update This method is intended to implement the explicit update calculation No...
Definition: explicit_solving_strategy.h:400
void FinalizeSolutionStep() override
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: explicit_solving_strategy.h:217
ExplicitBuilderType::DofType DofType
Definition: explicit_solving_strategy.h:70
ExplicitBuilderPointerType pGetExplicitBuilder()
Operations to get the pointer to the explicit builder and solver.
Definition: explicit_solving_strategy.h:286
void AssignSettings(const Parameters ThisParameters) override
This method assigns settings to member variables.
Definition: explicit_solving_strategy.h:419
ExplicitSolvingStrategy(ModelPart &rModelPart, Parameters ThisParameters)
Default constructor. (with parameters)
Definition: explicit_solving_strategy.h:97
void Initialize() override
Initialization of member variables and prior operations.
Definition: explicit_solving_strategy.h:171
void SetRebuildLevel(int Level) override
Definition: explicit_solving_strategy.h:263
const ExplicitBuilderType & GetExplicitBuilder() const
Operations to get the explicit builder and solver.
Definition: explicit_solving_strategy.h:305
ExplicitBuilderType & GetExplicitBuilder()
Operations to get the explicit builder and solver.
Definition: explicit_solving_strategy.h:295
virtual ~ExplicitSolvingStrategy()
Definition: explicit_solving_strategy.h:145
ExplicitSolvingStrategy(const ExplicitSolvingStrategy &Other)=delete
ExplicitSolvingStrategy(ModelPart &rModelPart, bool MoveMeshFlag=false, int RebuildLevel=0)
Default constructor.
Definition: explicit_solving_strategy.h:129
ExplicitSolvingStrategy(ModelPart &rModelPart, typename ExplicitBuilderType::Pointer pExplicitBuilder, bool MoveMeshFlag=false, int RebuildLevel=0)
Default constructor.
Definition: explicit_solving_strategy.h:112
static std::string Name()
Returns the name of the class as used in the settings (snake_case format)
Definition: explicit_solving_strategy.h:360
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 class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
int GetInt() const
This method returns the integer contained in the current Parameter.
Definition: kratos_parameters.cpp:666
void RecursivelyAddMissingParameters(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing contain at least all parameters...
Definition: kratos_parameters.cpp:1457
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 GetMoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:299
virtual void AssignSettings(const Parameters ThisParameters)
This method assigns settings to member variables.
Definition: solving_strategy.h:519
bool MoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:290
virtual void MoveMesh()
This function is designed to move the mesh.
Definition: solving_strategy.h:330
virtual Parameters GetDefaultParameters() const
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: solving_strategy.h:397
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
void InitializeAllEntities(ModelPart &rModelPart)
This method initializes all the active entities (conditions, elements, constraints)
Definition: entities_utilities.cpp:199
void FinalizeNonLinearIterationAllEntities(ModelPart &rModelPart)
This method calls FinalizeNonLinearIteration for all the entities (conditions, elements,...
Definition: entities_utilities.cpp:243
void InitializeSolutionStepAllEntities(ModelPart &rModelPart)
This method calls InitializeSolution for all the entities (conditions, elements, constraints)
Definition: entities_utilities.cpp:210
void FinalizeSolutionStepAllEntities(ModelPart &rModelPart)
This method calls FinalizeSolutionStep for all the entities (conditions, elements,...
Definition: entities_utilities.cpp:221
void InitializeNonLinearIterationAllEntities(ModelPart &rModelPart)
This method calls InitializeNonLinearIteration for all the entities (conditions, elements,...
Definition: entities_utilities.cpp:232
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21