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.
iterative_solver_strategy.h
Go to the documentation of this file.
1 //
2 // Author: Miquel Santasusana msantasusana@cimne.upc.edu
3 //
4 
5 #if !defined(KRATOS_ITERATIVE_SOLVER_STRATEGY)
6 #define KRATOS_ITERATIVE_SOLVER_STRATEGY
8 
9 namespace Kratos
10 {
12  {
13  public:
14 
16 
21 
25 
28 
31 
33  ExplicitSolverSettings& settings,
34  const double max_delta_time,
35  const double n_step_search,
36  const double safety_factor,
37  const int delta_option,
38  ParticleCreatorDestructor::Pointer p_creator_destructor,
39  DEM_FEM_Search::Pointer p_dem_fem_search,
40  SpatialSearch::Pointer pSpSearch,
41  Parameters strategy_parameters)
42  :ExplicitSolverStrategy(settings, max_delta_time, n_step_search, safety_factor, delta_option, p_creator_destructor, p_dem_fem_search, pSpSearch, strategy_parameters)
43  {
44  BaseType::GetParticleCreatorDestructor() = p_creator_destructor;
45  }
46 
49  {
50  //Timer::SetOuputFile("TimesPartialRelease");
51  //Timer::PrintTimingInformation();
52  }
53 
54 
55  virtual void Initialize() override
56  {
57 
59  ModelPart& r_model_part = BaseType::GetModelPart();
62  BaseType::ForceOperations(r_model_part);
64  //SchemeInitialize();
65 
66  KRATOS_CATCH("")
67  }// Initialize()
68 
69 
71  {
73  //BaseType::GetScheme()->Calculate(BaseType::GetModelPart(),0);
74  //BaseType::GetScheme()->Calculate(BaseType::GetClusterModelPart(),0);
75  }
76 
77 
79  {
81  //BaseType::GetScheme()->Calculate(BaseType::GetModelPart(),1); //TODO: better call Predict function (would be empty in general)
82  //BaseType::GetScheme()->Calculate(BaseType::GetClusterModelPart(),1);
83  }
84 
86  {
88  //BaseType::GetScheme()->Calculate(BaseType::GetModelPart(),2); //TODO: better call Correct function (normal operations would be in that method)
89  //BaseType::GetScheme()->Calculate(BaseType::GetClusterModelPart(),2);
90  }
91 
92  virtual double SolveSolutionStep() override
93  {
94 
96 
97  ModelPart& r_model_part = BaseType::GetModelPart();
98 
99  SchemePredict();
100  BaseType::SearchDEMOperations(r_model_part);
101  BaseType::SearchFEMOperations(r_model_part);
102  BaseType::ForceOperations(r_model_part);
103  SchemeCorrect();
104 
105  return 0.00;
106 
107  KRATOS_CATCH("")
108 
109  }//SolveSolutionStep()
110 
111  };//ClassIterativeSolverStrategy
112 
113 
114 
115 
116 } // namespace Kratos.
117 
118 #endif // KRATOS_FILENAME_H_INCLUDED defined
Definition: explicit_solver_strategy.h:54
Definition: explicit_solver_strategy.h:70
ModelPart & GetModelPart()
Definition: explicit_solver_strategy.h:256
void SearchDEMOperations(ModelPart &r_model_part, bool has_mpi=true)
Definition: explicit_solver_strategy.cpp:423
virtual void FinalizeSolutionStep()
Definition: explicit_solver_strategy.cpp:673
ModelPart::ElementsContainerType ElementsArrayType
Definition: explicit_solver_strategy.h:74
void InitializeSolutionStep()
Definition: explicit_solver_strategy.cpp:621
virtual void Initialize()
Definition: explicit_solver_strategy.cpp:125
virtual void ForceOperations(ModelPart &r_model_part)
Definition: explicit_solver_strategy.cpp:495
ModelPart::NodesContainerType NodesArrayType
Definition: explicit_solver_strategy.h:73
virtual void PerformTimeIntegrationOfMotion(int StepFlag=0)
Definition: explicit_solver_strategy.cpp:557
ElementsArrayType::iterator ElementsIterator
Definition: explicit_solver_strategy.h:75
ParticleCreatorDestructor::Pointer & GetParticleCreatorDestructor()
Definition: explicit_solver_strategy.h:272
void SearchFEMOperations(ModelPart &r_model_part, bool has_mpi=true)
Definition: explicit_solver_strategy.cpp:469
ModelPart::ConditionsContainerType ConditionsArrayType
Definition: explicit_solver_strategy.h:76
typename TContainerType::iterator ptr_iterator
Definition: global_pointers_vector.h:85
Definition: iterative_solver_strategy.h:12
BaseType::NodesArrayType NodesArrayType
Definition: iterative_solver_strategy.h:17
IterativeSolverStrategy()
Default constructor.
Definition: iterative_solver_strategy.h:30
BaseType::ConditionsArrayType ConditionsArrayType
Definition: iterative_solver_strategy.h:20
virtual void Initialize() override
Definition: iterative_solver_strategy.h:55
virtual double SolveSolutionStep() override
Definition: iterative_solver_strategy.h:92
BaseType::ElementsIterator ElementsIterator
Definition: iterative_solver_strategy.h:19
ParticleWeakVectorType::ptr_iterator ParticleWeakIteratorType_ptr
Definition: iterative_solver_strategy.h:24
BaseType::ElementsArrayType ElementsArrayType
Definition: iterative_solver_strategy.h:18
IterativeSolverStrategy(ExplicitSolverSettings &settings, const double max_delta_time, const double n_step_search, const double safety_factor, const int delta_option, ParticleCreatorDestructor::Pointer p_creator_destructor, DEM_FEM_Search::Pointer p_dem_fem_search, SpatialSearch::Pointer pSpSearch, Parameters strategy_parameters)
Definition: iterative_solver_strategy.h:32
KRATOS_CLASS_POINTER_DEFINITION(IterativeSolverStrategy)
Pointer definition of ExplicitSolverStrategy.
virtual ~IterativeSolverStrategy()
Destructor.
Definition: iterative_solver_strategy.h:48
void SchemeInitialize()
Definition: iterative_solver_strategy.h:70
ExplicitSolverStrategy BaseType
Definition: iterative_solver_strategy.h:15
void SchemePredict()
Definition: iterative_solver_strategy.h:78
void SchemeCorrect()
Definition: iterative_solver_strategy.h:85
GlobalPointersVector< Element > ParticleWeakVectorType
Definition: iterative_solver_strategy.h:22
GlobalPointersVector< Element >::iterator ParticleWeakIteratorType
Definition: iterative_solver_strategy.h:23
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
#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
safety_factor
Definition: edgebased_PureConvection.py:110