1 #ifndef KRATOS_STOKES_INITIALIZATION_PROCESS_H
2 #define KRATOS_STOKES_INITIALIZATION_PROCESS_H
53 template<
class TSparseSpace,
71 typename TLinearSolver::Pointer pLinearSolver,
72 unsigned int DomainSize,
93 std::string ElementName;
95 ElementName = std::string(
"StationaryStokes2D");
97 ElementName = std::string(
"StationaryStokes3D");
104 Element::Pointer pElem = rReferenceElement.
Create(itElem->Id(), itElem->GetGeometry(), itElem->pGetProperties() );
105 r_stokes_part.
Elements().push_back(pElem);
117 BuilderSolverTypePointer pBuildAndSolver;
119 PeriodicPairIndicesVar));
123 bool ReactionFlag =
false;
124 bool ReformDofSetFlag =
false;
125 bool CalculateNormDxFlag =
false;
126 bool MoveMeshFlag =
false;
176 KRATOS_THROW_ERROR(std::runtime_error,
"Stokes initialization process was called twice",
"");
205 std::string
Info()
const override
207 std::stringstream buffer;
208 buffer <<
" StokesInitializationProcess";
216 rOStream <<
" StokesInitializationProcess";
258 typename TLinearSolver::Pointer pLinearSolver,
259 unsigned int DomainSize,
376 template<
class TSparseSpace,
388 template<
class TSparseSpace,
396 rOStream << std::endl;
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