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.
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: Riccardo Rossi
11 // Ruben Zorrilla
12 //
13 
14 #if !defined(KRATOS_SOLVING_STRATEGY)
15 #define KRATOS_SOLVING_STRATEGY
16 
17 // System includes
18 
19 
20 // External includes
21 
22 
23 // Project includes
24 #include "includes/define.h"
25 #include "includes/model_part.h"
31 
32 namespace Kratos
33 {
34 
37 
38 
42 
43 
47 
48 
52 
53 
57 
62 template<class TSparseSpace, class TDenseSpace>
64 {
65 public:
68 
69  typedef typename TSparseSpace::DataType TDataType;
70 
72 
74 
75  typedef typename TSparseSpace::MatrixPointerType TSystemMatrixPointerType;
76 
77  typedef typename TSparseSpace::VectorPointerType TSystemVectorPointerType;
78 
80 
82 
84 
85  typedef typename ModelPart::DofType TDofType;
86 
88 
90 
92 
94 
97 
101 
105  explicit SolvingStrategy() { }
106 
112  explicit SolvingStrategy(
113  ModelPart& rModelPart,
114  Parameters ThisParameters)
115  : mpModelPart(&rModelPart)
116  {
117  // Validate and assign defaults
118  ThisParameters = this->ValidateAndAssignParameters(ThisParameters, this->GetDefaultParameters());
119  this->AssignSettings(ThisParameters);
120  }
121 
127  explicit SolvingStrategy(
128  ModelPart& rModelPart,
129  bool MoveMeshFlag = false)
130  : mpModelPart(&rModelPart)
131  {
133  }
134 
137  virtual ~SolvingStrategy(){}
138 
142 
143 
147 
153  virtual typename ClassType::Pointer Create(
154  ModelPart& rModelPart,
155  Parameters ThisParameters) const
156  {
157  return Kratos::make_shared<ClassType>(rModelPart, ThisParameters);
158  }
159 
164  virtual void Predict()
165  {
166  }
167 
171  virtual void Initialize()
172  {
173  }
174 
183  virtual double Solve()
184  {
185  Initialize();
187  Predict();
190 
191  return 0.0;
192 
193  }
194 
198  virtual void Clear()
199  {
200  }
201 
206  virtual bool IsConverged()
207  {
208  return true;
209  }
210 
216  virtual void CalculateOutputData()
217  {
218  }
219 
224  virtual void InitializeSolutionStep()
225  {
226  }
227 
232  virtual void FinalizeSolutionStep()
233  {
234  }
235 
239  virtual bool SolveSolutionStep()
240  {
241  return true;
242  }
243 
255  virtual void SetEchoLevel(const int Level)
256  {
257  mEchoLevel = Level;
258  }
259 
272  {
273  return mEchoLevel;
274  }
275 
280  void SetMoveMeshFlag(bool Flag)
281  {
282  mMoveMeshFlag = Flag;
283  }
284 
285  // TODO: DEPRECATE THIS IN FAVOR OF GetMoveMeshFlag()
291  {
292  return mMoveMeshFlag;
293  }
294 
300  {
301  return mMoveMeshFlag;
302  }
303 
310  virtual void SetRebuildLevel(int Level)
311  {
312  KRATOS_ERROR << "Accessing the strategy base class \'SetRebuildLevel\'. Please implement it in your derived class." << std::endl;
313  }
314 
321  virtual int GetRebuildLevel() const
322  {
323  KRATOS_ERROR << "Accessing the strategy base class \'GetRebuildLevel\'. Please implement it in your derived class." << std::endl;
324  }
325 
330  virtual void MoveMesh()
331  {
332  KRATOS_TRY
333 
334  KRATOS_ERROR_IF_NOT(GetModelPart().HasNodalSolutionStepVariable(DISPLACEMENT_X)) << "It is impossible to move the mesh since the DISPLACEMENT var is not in the Model Part. Either use SetMoveMeshFlag(False) or add DISPLACEMENT to the list of variables" << std::endl;
335 
336  block_for_each(GetModelPart().Nodes(), [](Node& rNode){
337  noalias(rNode.Coordinates()) = rNode.GetInitialPosition().Coordinates();
338  noalias(rNode.Coordinates()) += rNode.FastGetSolutionStepValue(DISPLACEMENT);
339  });
340 
341  KRATOS_INFO_IF("SolvingStrategy", this->GetEchoLevel() != 0) << " MESH MOVED " << std::endl;
342 
343  KRATOS_CATCH("")
344  }
345 
351  {
352  return *mpModelPart;
353  };
354 
359  const ModelPart& GetModelPart() const
360  {
361  return *mpModelPart;
362  };
363 
368  virtual double GetResidualNorm()
369  {
370  return 0.0;
371  }
372 
377  virtual int Check()
378  {
379  KRATOS_TRY
380 
381  // Check if displacement var is needed
382  if (mMoveMeshFlag) {
383  VariableUtils().CheckVariableExists<>(DISPLACEMENT, GetModelPart().Nodes());
384  }
385 
386  GetModelPart().Check();
387 
388  return 0;
389 
390  KRATOS_CATCH("")
391  }
392 
398  {
399  const Parameters default_parameters = Parameters(R"(
400  {
401  "name" : "solving_strategy",
402  "move_mesh_flag" : false,
403  "echo_level" : 1
404  })");
405  return default_parameters;
406  }
407 
413  {
414  KRATOS_TRY
415 
416  KRATOS_ERROR << "\'GetSystemMatrix\' not implemented in base \'SolvingStrategy\'" << std::endl;
417 
418  KRATOS_CATCH("");
419  }
420 
426  {
427  KRATOS_TRY
428 
429  KRATOS_ERROR << "\'GetSystemVector\' not implemented in base \'SolvingStrategy\'" << std::endl;
430 
431  KRATOS_CATCH("");
432  }
433 
439  {
440  KRATOS_TRY
441 
442  KRATOS_ERROR << "\'GetSolutionVector\' not implemented in base \'SolvingStrategy\'" << std::endl;
443 
444  KRATOS_CATCH("");
445  }
446 
451  static std::string Name()
452  {
453  return "solving_strategy";
454  }
455 
459 
461  virtual std::string Info() const
462  {
463  return "SolvingStrategy";
464  }
465 
467  virtual void PrintInfo(std::ostream& rOStream) const
468  {
469  rOStream << Info();
470  }
471 
473  virtual void PrintData(std::ostream& rOStream) const
474  {
475  rOStream << Info();
476  }
477 
479 
480 protected:
483 
484  // Level of echo for the solving strategy
486 
490 
491 
495 
496 
500 
508  Parameters ThisParameters,
509  const Parameters DefaultParameters) const
510  {
511  ThisParameters.ValidateAndAssignDefaults(DefaultParameters);
512  return ThisParameters;
513  }
514 
519  virtual void AssignSettings(const Parameters ThisParameters)
520  {
521  // By default mesh is not moved
522  mMoveMeshFlag = ThisParameters["move_mesh_flag"].GetBool();
523 
524  // Be default the minimal information is shown
525  mEchoLevel = ThisParameters["echo_level"].GetInt();
526  }
527 
531 
532 
536 
537 
541 
542 private:
543 
547 
548 
552 
553  ModelPart* mpModelPart = nullptr;
554 
555  bool mMoveMeshFlag;
556 
560 
561 
565 
566 
570 
571 
575 
576 
580 
583  SolvingStrategy(const SolvingStrategy& Other);
584 
586 
587 }; /* Class SolvingStrategy */
588 
590 
593 
594 
596 
597 } /* namespace Kratos.*/
598 
599 #endif /* KRATOS_SOLVING_STRATEGY defined */
Dof represents a degree of freedom (DoF).
Definition: dof.h:86
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
virtual int Check() const
run input validation
Definition: model_part.cpp:2204
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
This class defines the node.
Definition: node.h:65
TVariableType::Type & FastGetSolutionStepValue(const TVariableType &rThisVariable)
Definition: node.h:435
const PointType & GetInitialPosition() const
Definition: node.h:559
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 ValidateAndAssignDefaults(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing match the form prescribed by th...
Definition: kratos_parameters.cpp:1306
bool GetBool() const
This method returns the boolean contained in the current Parameter.
Definition: kratos_parameters.cpp:675
CoordinatesArrayType const & Coordinates() const
Definition: point.h:215
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
Solving strategy base class This is the base class from which we will derive all the strategies (impl...
Definition: solving_strategy.h:64
virtual bool SolveSolutionStep()
Solves the current step. This function returns true if a solution has been found, false otherwise.
Definition: solving_strategy.h:239
ModelPart::DofType TDofType
Definition: solving_strategy.h:85
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: solving_strategy.h:467
TSparseSpace::VectorPointerType TSystemVectorPointerType
Definition: solving_strategy.h:77
virtual void Predict()
Operation to predict the solution ... if it is not called a trivial predictor is used in which the va...
Definition: solving_strategy.h:164
virtual std::string Info() const
Turn back information as a string.
Definition: solving_strategy.h:461
virtual double GetResidualNorm()
Operations to get the residual norm.
Definition: solving_strategy.h:368
TDenseSpace::MatrixType LocalSystemMatrixType
Definition: solving_strategy.h:79
virtual void InitializeSolutionStep()
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: solving_strategy.h:224
ModelPart::NodesContainerType NodesArrayType
Definition: solving_strategy.h:89
virtual void SetRebuildLevel(int Level)
Set the Rebuild Level value This functions sets the rebuild level of the strategy It is only intended...
Definition: solving_strategy.h:310
virtual int GetRebuildLevel() const
Get the Rebuild Level value This function returns the rebuild level of the strategy It is only intend...
Definition: solving_strategy.h:321
virtual bool IsConverged()
This should be considered as a "post solution" convergence check which is useful for coupled analysis...
Definition: solving_strategy.h:206
TSparseSpace::DataType TDataType
Definition: solving_strategy.h:69
virtual TSystemVectorType & GetSolutionVector()
This method returns the solution vector.
Definition: solving_strategy.h:438
virtual TSystemMatrixType & GetSystemMatrix()
This method returns the LHS matrix.
Definition: solving_strategy.h:412
virtual void FinalizeSolutionStep()
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: solving_strategy.h:232
virtual void SetEchoLevel(const int Level)
This sets the level of echo for the solving strategy.
Definition: solving_strategy.h:255
ModelPart & GetModelPart()
Operations to get the pointer to the model.
Definition: solving_strategy.h:350
ModelPart::ElementsContainerType ElementsArrayType
Definition: solving_strategy.h:91
virtual void Clear()
Clears the internal storage.
Definition: solving_strategy.h:198
virtual Parameters ValidateAndAssignParameters(Parameters ThisParameters, const Parameters DefaultParameters) const
This method validate and assign default parameters.
Definition: solving_strategy.h:507
int GetEchoLevel()
This returns the level of echo for the solving strategy.
Definition: solving_strategy.h:271
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: solving_strategy.h:473
TSparseSpace::MatrixType TSystemMatrixType
Definition: solving_strategy.h:71
virtual ClassType::Pointer Create(ModelPart &rModelPart, Parameters ThisParameters) const
Create method.
Definition: solving_strategy.h:153
virtual void CalculateOutputData()
This operations should be called before printing the results when non trivial results (e....
Definition: solving_strategy.h:216
virtual void Initialize()
Initialization of member variables and prior operations.
Definition: solving_strategy.h:171
SolvingStrategy(ModelPart &rModelPart, Parameters ThisParameters)
Default constructor. (with parameters)
Definition: solving_strategy.h:112
SolvingStrategy()
Default constructor.
Definition: solving_strategy.h:105
virtual ~SolvingStrategy()
Definition: solving_strategy.h:137
bool GetMoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:299
const ModelPart & GetModelPart() const
Operations to get the pointer to the model.
Definition: solving_strategy.h:359
virtual int Check()
Function to perform expensive checks.
Definition: solving_strategy.h:377
TSparseSpace::MatrixPointerType TSystemMatrixPointerType
Definition: solving_strategy.h:75
ModelPart::ConditionsContainerType ConditionsArrayType
Definition: solving_strategy.h:93
SolvingStrategy< TSparseSpace, TDenseSpace > ClassType
Definition: solving_strategy.h:83
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 double Solve()
The problem of interest is solved.
Definition: solving_strategy.h:183
int mEchoLevel
Definition: solving_strategy.h:485
void SetMoveMeshFlag(bool Flag)
This function sets the flag that says if the mesh is moved.
Definition: solving_strategy.h:280
TSparseSpace::VectorType TSystemVectorType
Definition: solving_strategy.h:73
virtual void MoveMesh()
This function is designed to move the mesh.
Definition: solving_strategy.h:330
SolvingStrategy(ModelPart &rModelPart, bool MoveMeshFlag=false)
Default constructor.
Definition: solving_strategy.h:127
virtual TSystemVectorType & GetSystemVector()
This method returns the RHS vector.
Definition: solving_strategy.h:425
virtual Parameters GetDefaultParameters() const
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: solving_strategy.h:397
KRATOS_CLASS_POINTER_DEFINITION(SolvingStrategy)
ModelPart::DofsArrayType DofsArrayType
Definition: solving_strategy.h:87
static std::string Name()
Returns the name of the class as used in the settings (snake_case format)
Definition: solving_strategy.h:451
TDenseSpace::VectorType LocalSystemVectorType
Definition: solving_strategy.h:81
This class implements a set of auxiliar, already parallelized, methods to perform some common tasks r...
Definition: variable_utils.h:63
int CheckVariableExists(const TVarType &rVariable, const NodesContainerType &rNodes)
Checks if all the nodes of a node set has the specified variable.
Definition: variable_utils.h:1051
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF_NOT(conditional)
Definition: exception.h:163
#define KRATOS_INFO_IF(label, conditional)
Definition: logger.h:251
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
bool HasNodalSolutionStepVariable(ModelPart &rModelPart, Variable< TDataType > const &rThisVariable)
Definition: add_model_part_to_python.cpp:39
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
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484