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.
stokes_initialization_process.h
Go to the documentation of this file.
1 #ifndef KRATOS_STOKES_INITIALIZATION_PROCESS_H
2 #define KRATOS_STOKES_INITIALIZATION_PROCESS_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"
14 #include "containers/model.h"
15 #include "includes/model_part.h"
16 #include "processes/process.h"
22 
23 namespace Kratos
24 {
27 
30 
34 
38 
42 
46 
48 
53 template<class TSparseSpace,
54  class TDenseSpace,
55  class TLinearSolver
56  >
58 {
59 public:
62 
65 
69 
71  typename TLinearSolver::Pointer pLinearSolver,
72  unsigned int DomainSize,
73  const Variable<int>& PeriodicPairIndicesVar):
74  Process(),
75  mrReferenceModelPart(rModelPart),
76  mpLinearSolver(pLinearSolver),
77  mDomainSize(DomainSize)
78  {
79  KRATOS_TRY;
80 
81  if(!mrReferenceModelPart.GetModel().HasModelPart("StokesModelPart"))
82  mrReferenceModelPart.GetModel().DeleteModelPart("StokesModelPart");
83 
84  ModelPart& r_stokes_part = mrReferenceModelPart.GetModel().CreateModelPart("StokesModelPart");
85 
87  r_stokes_part.SetBufferSize(1);
88  r_stokes_part.SetNodes( mrReferenceModelPart.pNodes() );
91 
92  // Retrieve Stokes element model
93  std::string ElementName;
94  if (mDomainSize == 2)
95  ElementName = std::string("StationaryStokes2D");
96  else
97  ElementName = std::string("StationaryStokes3D");
98 
99  const Element& rReferenceElement = KratosComponents<Element>::Get(ElementName);
100 
101  // Generate Stokes elements
102  for (ModelPart::ElementsContainerType::const_iterator itElem = mrReferenceModelPart.ElementsBegin(); itElem != mrReferenceModelPart.ElementsEnd(); itElem++)
103  {
104  Element::Pointer pElem = rReferenceElement.Create(itElem->Id(), itElem->GetGeometry(), itElem->pGetProperties() );
105  r_stokes_part.Elements().push_back(pElem);
106  }
107 
108  // pointer types for the solution strategy construcion
109  typedef typename Scheme< TSparseSpace, TDenseSpace >::Pointer SchemePointerType;
110  typedef typename BuilderAndSolver<TSparseSpace, TDenseSpace, TLinearSolver>::Pointer BuilderSolverTypePointer;
112 
113  // Solution scheme: Linear static scheme
114  SchemePointerType pScheme = SchemePointerType( new ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace > () );
115 
116  // Builder and solver
117  BuilderSolverTypePointer pBuildAndSolver;
119  PeriodicPairIndicesVar));
120 
121 
122  // Strategy
123  bool ReactionFlag = false;
124  bool ReformDofSetFlag = false;
125  bool CalculateNormDxFlag = false;
126  bool MoveMeshFlag = false;
128  r_stokes_part,
129  pScheme,
130  pBuildAndSolver,
131  ReactionFlag,
132  ReformDofSetFlag,
133  CalculateNormDxFlag,
134  MoveMeshFlag) );
135  mpSolutionStrategy->SetEchoLevel(0);
136  mpSolutionStrategy->Check();
137 
138  mIsCleared = false;
139 
140  KRATOS_CATCH("");
141  }
142 
144 
146  {
147  mrReferenceModelPart.GetModel().DeleteModelPart("StokesModelPart");
148  // mpSolutionStrategy->Clear();
149  }
150 
151 
155 
156 
160 
162  void Execute() override
163  {
164  KRATOS_TRY;
165 
166  if (!mIsCleared)
167  {
168  // Solve Stokes problem.
169  this->Solve();
170 
171  // Destroy the auxiliary ModelPart and solution structures.
172  this->Clear();
173  }
174  else
175  {
176  KRATOS_THROW_ERROR(std::runtime_error,"Stokes initialization process was called twice","");
177  }
178 
179  KRATOS_CATCH("");
180  }
181 
185 
186  void SetConditions(ModelPart& rStokesPart, ModelPart::ConditionsContainerType::Pointer pConditions)
187  {
188 
189  ModelPart& r_stokes_part = mrReferenceModelPart.GetModel().GetModelPart("StokesModelPart");
190  rStokesPart.SetConditions(pConditions);
191  r_stokes_part.GetCommunicator().LocalMesh().SetConditions(pConditions);
192  }
193 
197 
198 
202 
204 
205  std::string Info() const override
206  {
207  std::stringstream buffer;
208  buffer << " StokesInitializationProcess";
209  return buffer.str();
210  }
211 
213 
214  void PrintInfo(std::ostream& rOStream) const override
215  {
216  rOStream << " StokesInitializationProcess";
217  }
218 
220 
221  void PrintData(std::ostream& rOStream) const override
222  {
223  }
224 
225 
229 
230 
232 
233 protected:
236 
237 
241 
243 
244  typename TLinearSolver::Pointer mpLinearSolver;
245 
246  unsigned int mDomainSize;
247 
249 
251 
255 
258  typename TLinearSolver::Pointer pLinearSolver,
259  unsigned int DomainSize,
260  const StokesInitializationProcess* pThis):
261  Process(),
262  mrReferenceModelPart(rModelPart),
263  mpLinearSolver(pLinearSolver),
264  mDomainSize(DomainSize)
265  {}
266 
270 
271 
275 
276 
277 
278  // Solve Stokes problem.
279  virtual void Solve()
280  {
281  KRATOS_TRY;
282 
283  mpSolutionStrategy->Solve();
284 
285  KRATOS_CATCH("");
286  }
287 
288 
289  // Liberate memory.
290  void Clear() override
291  {
292  mpSolutionStrategy->Clear();
293  mIsCleared = true;
294  }
295 
299 
300 
304 
305 
309 
311 
312 private:
315 
316 
320 
321 
325 
326 
330 
331 
335 
336 
340 
341 
345 
347 
348  StokesInitializationProcess & operator=( StokesInitializationProcess const& rOther)
349  {
350  return *this;
351  }
352 
354 
356  { }
357 
358 
360 
361 }; // Class StokesInitializationProcess
362 
364 
367 
368 
372 
373 
375 
376 template<class TSparseSpace,
377  class TDenseSpace,
378  class TLinearSolver
379  >
380 inline std::istream & operator >>(std::istream& rIStream,
382 {
383  return rIStream;
384 }
385 
387 
388 template<class TSparseSpace,
389  class TDenseSpace,
390  class TLinearSolver
391  >
392 inline std::ostream & operator <<(std::ostream& rOStream,
394 {
395  rThis.PrintInfo(rOStream);
396  rOStream << std::endl;
397  rThis.PrintData(rOStream);
398 
399  return rOStream;
400 }
402 
404 
405 } // namespace Kratos.
406 
407 
408 #endif // KRATOS_STOKES_INITIALIZATION_PROCESS_H
Current class provides an implementation for the base builder and solving operations.
Definition: builder_and_solver.h:64
MeshType & LocalMesh()
Returns the reference to the mesh storing all local entities.
Definition: communicator.cpp:245
Base class for all Elements.
Definition: element.h:60
virtual Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const
It creates a new element pointer.
Definition: element.h:202
Implicit solving strategy base class This is the base class from which we will derive all the implici...
Definition: implicit_solving_strategy.h:61
static const TComponentType & Get(const std::string &rName)
Retrieves a component with the specified name.
Definition: kratos_components.h:114
void SetConditions(typename ConditionsContainerType::Pointer pOtherConditions)
Definition: mesh.h:706
ModelPart & GetModelPart(const std::string &rFullModelPartName)
This method returns a model part given a certain name.
Definition: model.cpp:107
ModelPart & CreateModelPart(const std::string &ModelPartName, IndexType NewBufferSize=1)
This method creates a new model part contained in the current Model with a given name and buffer size...
Definition: model.cpp:37
void DeleteModelPart(const std::string &ModelPartName)
This method deletes a modelpart with a given name.
Definition: model.cpp:64
bool HasModelPart(const std::string &rFullModelPartName) const
This method checks if a certain a model part exists given a certain name.
Definition: model.cpp:178
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ElementIterator ElementsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1169
void SetNodes(NodesContainerType::Pointer pOtherNodes, IndexType ThisIndex=0)
Definition: model_part.h:522
void SetProcessInfo(ProcessInfo::Pointer pNewProcessInfo)
Definition: model_part.h:1766
void SetProperties(PropertiesContainerType::Pointer pOtherProperties, IndexType ThisIndex=0)
Definition: model_part.h:1013
ProcessInfo::Pointer pGetProcessInfo()
Definition: model_part.h:1756
Communicator & GetCommunicator()
Definition: model_part.h:1821
void SetBufferSize(IndexType NewBufferSize)
This method sets the suffer size of the model part database.
Definition: model_part.cpp:2171
NodesContainerType::Pointer pNodes(IndexType ThisIndex=0)
Definition: model_part.h:517
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
PropertiesContainerType::Pointer pProperties(IndexType ThisIndex=0)
Definition: model_part.h:1008
ElementIterator ElementsEnd(IndexType ThisIndex=0)
Definition: model_part.h:1179
void SetConditions(ConditionsContainerType::Pointer pOtherConditions, IndexType ThisIndex=0)
Definition: model_part.h:1396
Model & GetModel()
Definition: model_part.h:323
VariablesList & GetNodalSolutionStepVariablesList()
Definition: model_part.h:549
The base class for all processes in Kratos.
Definition: process.h:49
Definition: residualbased_block_builder_and_solver_periodic.h:69
This class provides the implementation of a static scheme.
Definition: residualbased_incrementalupdate_static_scheme.h:57
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
A process to provide initial values for Navier-Stokes problems.
Definition: stokes_initialization_process.h:58
~StokesInitializationProcess() override
Destructor.
Definition: stokes_initialization_process.h:145
void Clear() override
This method clears the assignation of the conditions.
Definition: stokes_initialization_process.h:290
StokesInitializationProcess(ModelPart &rModelPart, typename TLinearSolver::Pointer pLinearSolver, unsigned int DomainSize, const StokesInitializationProcess *pThis)
Protected constructor to be used by derived classes.
Definition: stokes_initialization_process.h:257
void Execute() override
Solve a stationary Stokes problem.
Definition: stokes_initialization_process.h:162
ImplicitSolvingStrategy< TSparseSpace, TDenseSpace, TLinearSolver >::Pointer mpSolutionStrategy
Definition: stokes_initialization_process.h:248
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: stokes_initialization_process.h:221
StokesInitializationProcess(ModelPart &rModelPart, typename TLinearSolver::Pointer pLinearSolver, unsigned int DomainSize, const Variable< int > &PeriodicPairIndicesVar)
Definition: stokes_initialization_process.h:70
virtual void Solve()
Definition: stokes_initialization_process.h:279
std::string Info() const override
Turn back information as a string.
Definition: stokes_initialization_process.h:205
TLinearSolver::Pointer mpLinearSolver
Definition: stokes_initialization_process.h:244
ModelPart & mrReferenceModelPart
Definition: stokes_initialization_process.h:242
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: stokes_initialization_process.h:214
bool mIsCleared
Definition: stokes_initialization_process.h:250
KRATOS_CLASS_POINTER_DEFINITION(StokesInitializationProcess)
Pointer definition of StokesInitializationProcess.
void SetConditions(ModelPart &rStokesPart, ModelPart::ConditionsContainerType::Pointer pConditions)
Definition: stokes_initialization_process.h:186
unsigned int mDomainSize
Definition: stokes_initialization_process.h:246
#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
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