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.
manage_time_step_process.h
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: August 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_MANAGE_TIME_STEP_PROCESS_H_INCLUDED)
11 #define KRATOS_MANAGE_TIME_STEP_PROCESS_H_INCLUDED
12 
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
19 #include "includes/model_part.h"
21 #include "processes/process.h"
22 
23 namespace Kratos
24 {
25 
28 
30 
33 {
34 public:
37 
40 
44 
45 
47  double& rMinDeltaTime, double& rMaxDeltaTime,
48  double& rReductionFactor, double& rIncreaseFactor,
49  double& rErrorTolerance, unsigned int& rMinIterations,
50  unsigned int& rMaxIterations, unsigned int& rNumberOfConstantSteps) : Process(Flags()), mrModelPart(rModelPart), mMinDeltaTime(rMinDeltaTime), mMaxDeltaTime(rMaxDeltaTime), mReductionFactor(rReductionFactor), mIncreaseFactor(rIncreaseFactor), mErrorTolerance(rErrorTolerance), mMinIterations(rMinIterations), mMaxIterations(rMaxIterations), mNumberOfConstantSteps(rNumberOfConstantSteps)
51  {
52  }
53 
55  Parameters CustomParameters )
56  : mrModelPart(rModelPart)
57  {
59 
60  Parameters DefaultParameters( R"(
61  {
62  "time_step": 1.0,
63  "start_time": 0.0,
64  "end_time": 1.0,
65  "adaptive_time_step":{
66  "minimum_time_step": 0.1,
67  "maximum_time_step": 1.0,
68  "reduction_factor": 2.0,
69  "increase_factor": 1.0,
70  "error_tolerance": 1e-4,
71  "minimum_iterations": 2,
72  "maximum_iterations": 10,
73  "number_constant_steps": 4
74  }
75  } )" );
76 
77 
78  mAdaptiveTimeStep = false;
79  if( CustomParameters.Has("adaptive_time_step") )
80  mAdaptiveTimeStep = true;
81 
82  //validate against defaults -- this also ensures no type mismatch
83  CustomParameters.ValidateAndAssignDefaults(DefaultParameters);
84 
85  ProcessInfo& rCurrentProcessInfo = mrModelPart.GetProcessInfo();
86 
87  bool restarted = false;
88  if( rCurrentProcessInfo.Has(IS_RESTARTED) ){
89  if( rCurrentProcessInfo[IS_RESTARTED] == true ){
90  restarted = true;
91  }
92  }
93 
94  if( !restarted ){
95  rCurrentProcessInfo.SetValue(DELTA_TIME, CustomParameters["time_step"].GetDouble());
96  rCurrentProcessInfo.SetValue(TIME, CustomParameters["start_time"].GetDouble());
97  }
98 
99  mTime = rCurrentProcessInfo[TIME];
100  mStep = rCurrentProcessInfo[STEP];
101  mEndTime = CustomParameters["end_time"].GetDouble();
102 
103  if( mAdaptiveTimeStep )
104  SetAdaptiveTimeParameters(CustomParameters["adaptive_time_step"]);
105 
106  KRATOS_CATCH(" ")
107 
108  }
109 
112 
113 
117 
119  void operator()()
120  {
121  Execute();
122  }
123 
124 
128 
129 
131  void Execute() override
132  {
133  KRATOS_TRY;
134 
135 // const int nelements = mrModelPart.Elements().size();
136 
137 // if (nelements != 0)
138 // {
139 // ModelPart::ElementsContainerType::iterator it_begin =
140 // mrModelPart.ElementsBegin();
141 
142 // #pragma omp parallel for
143 // for (int i = 0; i < nelements; i++)
144 // {
145 // ModelPart::ElementsContainerType::iterator it = it_begin + i;
146 
147 // if (this->MatchTransferFlags(*(it.base())))
148 // {
149 // this->AssignFlags(*(it.base()));
150 // }
151 // }
152 // }
153 
154 
155 
156  KRATOS_CATCH("");
157  }
158 
161  void ExecuteInitialize() override
162  {
163  }
164 
168  {
169  }
170 
171 
174  {
175  ProcessInfo& rCurrentProcessInfo = mrModelPart.GetProcessInfo();
176 
177  if( mAdaptiveTimeStep )
178  PredictTimeStep();
179 
180  double mTime = rCurrentProcessInfo[TIME] + rCurrentProcessInfo[DELTA_TIME];
181 
182  mrModelPart.ProcessInfo[STEP] = (++mStep);
183 
184  mrModelPart.CloneTimeStep(mTime);
185  }
186 
189  {
190  }
191 
192 
194  void ExecuteBeforeOutputStep() override
195  {
196  }
197 
198 
200  void ExecuteAfterOutputStep() override
201  {
202  }
203 
204 
207  void ExecuteFinalize() override
208  {
209  }
210 
211 
215 
216 
220 
221 
225 
227  std::string Info() const override
228  {
229  return "ManageTimeStepProcess";
230  }
231 
233  void PrintInfo(std::ostream& rOStream) const override
234  {
235  rOStream << "ManageTimeStepProcess";
236  }
237 
239  void PrintData(std::ostream& rOStream) const override
240  {
241  }
242 
243 
248 
249 protected:
250 
259 
262 
276 
277 private:
278 
284 
285  ModelPart& mrModelPart;
286 
287  bool mAdaptiveTimStep;
288 
289  double mTime;
290  double mStep;
291  double mEndTime;
292 
293  double mMinDeltaTime;
294  double mMaxDeltaTime;
295 
296  double mReductionFactor;
297  double mIncreaseFactor;
298 
299  double mErrorTolerance;
300 
301  unsigned int mMinIterations;
302  unsigned int mMaxIterations;
303  unsigned int mNumberOfConstantSteps;
304 
311 
312  void SetAdaptiveTimeParameters(Parameters CustomParameters)
313  {
314  KRATOS_TRY
315 
316  Parameters DefaultParameters( R"(
317  {
318  "minimum_time_step": 0.1,
319  "maximum_time_step": 1.0,
320  "reduction_factor": 2.0,
321  "increase_factor": 1.0,
322  "error_tolerance": 1e-4,
323  "minimum_iterations": 2,
324  "maximum_iterations": 10,
325  "number_constant_steps": 4
326  } )" );
327 
328  //validate against defaults -- this also ensures no type mismatch
329  CustomParameters.ValidateAndAssignDefaults(DefaultParameters);
330 
331 
332  mMinDeltaTime = CustomParameters["minimum_time_step"].GetDouble();
333  mMaxDeltaTime = CustomParameters["maximum_time_step"].GetDouble();
334 
335  mReductionFactor = CustomParameters["reduction_factor"].GetDouble();
336  mIncreaseFactor = CustomParameters["increase_factor"].GetDouble();
337 
338  mErrorTolerance = CustomParameters["error_tolerance"].GetDouble();
339 
340  mMinIterations = CustomParameters["minimum_iterations"].GetInt();
341  mMaxIterations = CustomParameters["maximum_iterations"].GetInt();
342  mNumberOfConstantSteps = CustomParameters["number_constant_steps"].GetInt();
343 
344  KRATOS_CATCH(" ")
345  }
349 
351  ManageTimeStepProcess& operator=(ManageTimeStepProcess const& rOther);
352 
362 
363 }; // Class ManageTimeStepProcess
364 
365 
367 
370 
371 
375 
376 
378 inline std::istream& operator >> (std::istream& rIStream,
379  ManageTimeStepProcess& rThis);
380 
382 inline std::ostream& operator << (std::ostream& rOStream,
383  const ManageTimeStepProcess& rThis)
384 {
385  rThis.PrintInfo(rOStream);
386  rOStream << std::endl;
387  rThis.PrintData(rOStream);
388 
389  return rOStream;
390 }
392 
393 
394 } // namespace Kratos.
395 
396 #endif // KRATOS_MANAGE_TIME_STEP_PROCESS_H_INCLUDED defined
bool Has(const Variable< TDataType > &rThisVariable) const
Checks if the data container has a value associated with a given variable.
Definition: data_value_container.h:382
void SetValue(const Variable< TDataType > &rThisVariable, TDataType const &rValue)
Sets the value for a given variable.
Definition: data_value_container.h:320
Definition: flags.h:58
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: manage_time_step_process.h:33
void ExecuteBeforeSolutionLoop() override
Definition: manage_time_step_process.h:167
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: manage_time_step_process.h:233
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: manage_time_step_process.h:194
virtual ~ManageTimeStepProcess()
Destructor.
Definition: manage_time_step_process.h:111
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: manage_time_step_process.h:173
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: manage_time_step_process.h:119
ManageTimeStepProcess(ManageTimeStepProcess const &rOther)
Copy constructor.
void ExecuteFinalize() override
Definition: manage_time_step_process.h:207
void Execute() override
Execute method is used to execute the ManageTimeStepProcess algorithms.
Definition: manage_time_step_process.h:131
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: manage_time_step_process.h:188
std::string Info() const override
Turn back information as a string.
Definition: manage_time_step_process.h:227
void ExecuteInitialize() override
Definition: manage_time_step_process.h:161
ManageTimeStepProcess(ModelPart &rModelPart, Parameters CustomParameters)
Definition: manage_time_step_process.h:54
ManageTimeStepProcess(ModelPart &rModelPart, double &rMinDeltaTime, double &rMaxDeltaTime, double &rReductionFactor, double &rIncreaseFactor, double &rErrorTolerance, unsigned int &rMinIterations, unsigned int &rMaxIterations, unsigned int &rNumberOfConstantSteps)
Definition: manage_time_step_process.h:46
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: manage_time_step_process.h:200
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: manage_time_step_process.h:239
KRATOS_CLASS_POINTER_DEFINITION(ManageTimeStepProcess)
Pointer definition of ManageTimeStepProcess.
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
IndexType CloneTimeStep()
Definition: model_part.cpp:143
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
double GetDouble() const
This method returns the double contained in the current Parameter.
Definition: kratos_parameters.cpp:657
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 Has(const std::string &rEntry) const
This method checks if the Parameter contains a certain entry.
Definition: kratos_parameters.cpp:520
The base class for all processes in Kratos.
Definition: process.h:49
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.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
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