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.
poromechanics_explicit_nonlocal_strategy.hpp
Go to the documentation of this file.
1 
2 // | / |
3 // ' / __| _` | __| _ \ __|
4 // . \ | ( | | ( |\__ `
5 // _|\_\_| \__,_|\__|\___/ ____/
6 // Multi-Physics
7 //
8 // License: BSD License
9 // Kratos default license: kratos/license.txt
10 //
11 // Main authors: Ignasi de Pouplana
12 //
13 
14 
15 #if !defined(KRATOS_POROMECHANICS_EXPLICIT_NONLOCAL_STRATEGY)
16 #define KRATOS_POROMECHANICS_EXPLICIT_NONLOCAL_STRATEGY
17 
18 // Project includes
23 
24 // Application includes
26 
27 namespace Kratos {
28 
29 template <class TSparseSpace,
30  class TDenseSpace,
31  class TLinearSolver
32  >
34  : public PoromechanicsExplicitStrategy<TSparseSpace, TDenseSpace, TLinearSolver> {
35 public:
36 
38 
52 
54  typedef typename Node::DofType DofType;
55  typedef typename DofType::Pointer DofPointerType;
56 
61 
62 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
63 
67  typename TSchemeType::Pointer pScheme,
68  Parameters& rParameters,
69  bool CalculateReactions = false,
70  bool ReformDofSetAtEachStep = false,
71  bool MoveMeshFlag = false
72  ) : PoromechanicsExplicitStrategy<TSparseSpace, TDenseSpace, TLinearSolver>(model_part, pScheme, rParameters,
73  CalculateReactions, ReformDofSetAtEachStep, MoveMeshFlag)
74  {
76  mSearchNeighboursAtEachStep = rParameters["search_neighbours_step"].GetBool();
77  }
78 
79  //------------------------------------------------------------------------------------
80 
83 
84 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
85 
89  void Initialize() override
90  {
92 
93  if (mInitializeWasPerformed == false) {
95 
96  if(mNonlocalDamageIsInitialized == false) {
97  if(BaseType::GetModelPart().GetProcessInfo()[DOMAIN_SIZE]==2) {
99  } else {
101  }
103 
105  }
106  }
107 
108  KRATOS_CATCH( "" )
109  }
110 
111 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
112 
117  void InitializeSolutionStep() override {
118 
119  KRATOS_TRY
120 
122 
123  if(mNonlocalDamageIsInitialized == false) {
124  if(BaseType::GetModelPart().GetProcessInfo()[DOMAIN_SIZE]==2){
126  } else {
128  }
130 
132  }
133 
134  KRATOS_CATCH("")
135  }
136 
137 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
138 
142  bool SolveSolutionStep() override
143  {
144  ModelPart& r_model_part = BaseType::GetModelPart();
145 
146  // Some dummy sets and matrices
147  DofsArrayType dof_set_dummy;
151 
152  // Initialize the non linear iteration
153  mpScheme->InitializeNonLinIteration(BaseType::GetModelPart(), rA, rDx, rb);
155 
156  mpScheme->Predict(r_model_part, dof_set_dummy, rA, rDx, rb);
157 
158  // Move the mesh if needed
161 
162  mpScheme->InitializeNonLinIteration(BaseType::GetModelPart(), rA, rDx, rb);
164 
165  // Explicitly integrates the equation of motion.
166  mpScheme->Update(r_model_part, dof_set_dummy, rA, rDx, rb);
167 
168  // CONVERGENCE CHECK
169  this->CheckConvergence(r_model_part);
170 
171  // Move the mesh if needed
174 
175  // Finalize the non linear iteration
176  mpScheme->FinalizeNonLinIteration(BaseType::GetModelPart(), rA, rDx, rb);
178 
179  // Calculate reactions if required
181  this->CalculateReactions(mpScheme, r_model_part);
182  }
183 
184  return true;
185  }
186 
187 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
188 
189  void FinalizeSolutionStep() override
190  {
191  KRATOS_TRY
192 
194 
195  if(mSearchNeighboursAtEachStep == true) {
198  }
199 
200  KRATOS_CATCH("")
201  }
202 
203 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
204 
205  void Clear() override
206  {
207  KRATOS_TRY
208 
210 
211  if(mSearchNeighboursAtEachStep == false) {
214  }
215 
216  KRATOS_CATCH( "" )
217  }
218 
219 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
220 
221 protected:
222 
227 
228 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
229 
230 
231 }; // Class PoromechanicsExplicitNonlocalStrategy
232 
233 } // namespace Kratos
234 
235 #endif // KRATOS_POROMECHANICS_EXPLICIT_NONLOCAL_STRATEGY defined
Dof represents a degree of freedom (DoF).
Definition: dof.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
BaseType::TSystemVectorType TSystemVectorType
Definition: implicit_solving_strategy.h:72
Scheme< TSparseSpace, TDenseSpace > TSchemeType
Definition: implicit_solving_strategy.h:82
BaseType::TSystemMatrixType TSystemMatrixType
Definition: implicit_solving_strategy.h:70
This strategy is used for the explicit time integration.
Definition: mechanical_explicit_strategy.hpp:54
TSchemeType::Pointer mpScheme
Definition: mechanical_explicit_strategy.hpp:587
bool mInitializeWasPerformed
Flag telling if the initialize was performed.
Definition: mechanical_explicit_strategy.hpp:608
void Clear() override
Clears the internal storage.
Definition: mechanical_explicit_strategy.hpp:405
bool mCalculateReactionsFlag
Flag telling if it is needed or not to compute the reactions.
Definition: mechanical_explicit_strategy.hpp:602
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Definition: nonlocal_damage_2D_utilities.hpp:26
Definition: nonlocal_damage_3D_utilities.hpp:26
Definition: nonlocal_damage_utilities.hpp:36
virtual void SearchGaussPointsNeighbours(Parameters *pParameters, ModelPart &rModelPart)
Definition: nonlocal_damage_utilities.hpp:82
void CalculateNonlocalEquivalentStrain(Parameters *pParameters, const ProcessInfo &CurrentProcessInfo)
Definition: nonlocal_damage_utilities.hpp:89
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
bool GetBool() const
This method returns the boolean contained in the current Parameter.
Definition: kratos_parameters.cpp:675
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
Definition: poromechanics_explicit_nonlocal_strategy.hpp:34
TSchemeType::Pointer mpScheme
Definition: mechanical_explicit_strategy.hpp:587
BaseType::TSchemeType TSchemeType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:42
BaseType::TSystemVectorPointerType TSystemVectorPointerType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:47
void Clear() override
Clears the internal storage.
Definition: poromechanics_explicit_nonlocal_strategy.hpp:205
ImplicitSolvingStrategy< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:39
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:51
MechanicalExplicitStrategy< TSparseSpace, TDenseSpace, TLinearSolver > GrandMotherType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:40
~PoromechanicsExplicitNonlocalStrategy() override
Destructor.
Definition: poromechanics_explicit_nonlocal_strategy.hpp:82
void FinalizeSolutionStep() override
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: poromechanics_explicit_nonlocal_strategy.hpp:189
PoromechanicsExplicitNonlocalStrategy(ModelPart &model_part, typename TSchemeType::Pointer pScheme, Parameters &rParameters, bool CalculateReactions=false, bool ReformDofSetAtEachStep=false, bool MoveMeshFlag=false)
Constructor.
Definition: poromechanics_explicit_nonlocal_strategy.hpp:65
PoromechanicsExplicitStrategy< TSparseSpace, TDenseSpace, TLinearSolver > MotherType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:41
bool mInitializeWasPerformed
Flag telling if the initialize was performed.
Definition: mechanical_explicit_strategy.hpp:608
bool mCalculateReactionsFlag
Flag telling if it is needed or not to compute the reactions.
Definition: mechanical_explicit_strategy.hpp:602
Node::DofType DofType
DoF types definition.
Definition: poromechanics_explicit_nonlocal_strategy.hpp:54
bool mNonlocalDamageIsInitialized
Definition: poromechanics_explicit_nonlocal_strategy.hpp:225
BaseType::TSystemMatrixPointerType TSystemMatrixPointerType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:46
NonlocalDamageUtilities * mpNonlocalDamageUtility
Member Variables.
Definition: poromechanics_explicit_nonlocal_strategy.hpp:224
bool SolveSolutionStep() override
Solves the current step. This function returns true if a solution has been found, false otherwise.
Definition: poromechanics_explicit_nonlocal_strategy.hpp:142
KRATOS_CLASS_POINTER_DEFINITION(PoromechanicsExplicitNonlocalStrategy)
Parameters * mpParameters
Member Variables.
Definition: poromechanics_explicit_strategy.hpp:281
BaseType::NodesArrayType NodesArrayType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:48
BaseType::TSystemVectorType TSystemVectorType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:45
DofType::Pointer DofPointerType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:55
void InitializeSolutionStep() override
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: poromechanics_explicit_nonlocal_strategy.hpp:117
BaseType::ElementsArrayType ElementsArrayType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:49
BaseType::TSystemMatrixType TSystemMatrixType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:44
void Initialize() override
Initialization of member variables and prior operations.
Definition: poromechanics_explicit_nonlocal_strategy.hpp:89
BaseType::DofsArrayType DofsArrayType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:43
BaseType::ConditionsArrayType ConditionsArrayType
Definition: poromechanics_explicit_nonlocal_strategy.hpp:50
bool mSearchNeighboursAtEachStep
Definition: poromechanics_explicit_nonlocal_strategy.hpp:226
Definition: poromechanics_explicit_strategy.hpp:34
void InitializeSolutionStep() override
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: poromechanics_explicit_strategy.hpp:174
virtual void CheckConvergence(ModelPart &rModelPart)
Definition: poromechanics_explicit_strategy.hpp:342
void Initialize() override
Initialization of member variables and prior operations.
Definition: poromechanics_explicit_strategy.hpp:120
virtual void CalculateReactions(typename TSchemeType::Pointer pScheme, ModelPart &rModelPart)
Name of the nodal variable associated to every SubModelPart.
Definition: poromechanics_explicit_strategy.hpp:295
Parameters * mpParameters
Member Variables.
Definition: poromechanics_explicit_strategy.hpp:281
void FinalizeSolutionStep() override
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: poromechanics_explicit_strategy.hpp:260
TSparseSpace::VectorPointerType TSystemVectorPointerType
Definition: solving_strategy.h:77
ModelPart::NodesContainerType NodesArrayType
Definition: solving_strategy.h:89
ModelPart & GetModelPart()
Operations to get the pointer to the model.
Definition: solving_strategy.h:350
ModelPart::ElementsContainerType ElementsArrayType
Definition: solving_strategy.h:91
TSparseSpace::MatrixType TSystemMatrixType
Definition: solving_strategy.h:71
TSparseSpace::MatrixPointerType TSystemMatrixPointerType
Definition: solving_strategy.h:75
ModelPart::ConditionsContainerType ConditionsArrayType
Definition: solving_strategy.h:93
bool MoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:290
TSparseSpace::VectorType TSystemVectorType
Definition: solving_strategy.h:73
virtual void MoveMesh()
This function is designed to move the mesh.
Definition: solving_strategy.h:330
TDenseSpace::VectorType LocalSystemVectorType
Definition: solving_strategy.h:81
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
ProcessInfo & GetProcessInfo(ModelPart &rModelPart)
Definition: add_model_part_to_python.cpp:54
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
model_part
Definition: face_heat.py:14