10 #if !defined(KRATOS_MANAGE_TIME_STEP_PROCESS_H_INCLUDED)
11 #define KRATOS_MANAGE_TIME_STEP_PROCESS_H_INCLUDED
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)
56 : mrModelPart(rModelPart)
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
78 mAdaptiveTimeStep = false;
79 if( CustomParameters.
Has(
"adaptive_time_step") )
80 mAdaptiveTimeStep =
true;
87 bool restarted =
false;
88 if( rCurrentProcessInfo.
Has(IS_RESTARTED) ){
89 if( rCurrentProcessInfo[IS_RESTARTED] ==
true ){
95 rCurrentProcessInfo.
SetValue(DELTA_TIME, CustomParameters[
"time_step"].GetDouble());
96 rCurrentProcessInfo.
SetValue(TIME, CustomParameters[
"start_time"].GetDouble());
99 mTime = rCurrentProcessInfo[TIME];
100 mStep = rCurrentProcessInfo[STEP];
101 mEndTime = CustomParameters[
"end_time"].
GetDouble();
103 if( mAdaptiveTimeStep )
104 SetAdaptiveTimeParameters(CustomParameters[
"adaptive_time_step"]);
177 if( mAdaptiveTimeStep )
180 double mTime = rCurrentProcessInfo[TIME] + rCurrentProcessInfo[DELTA_TIME];
182 mrModelPart.ProcessInfo[STEP] = (++mStep);
227 std::string
Info()
const override
229 return "ManageTimeStepProcess";
235 rOStream <<
"ManageTimeStepProcess";
287 bool mAdaptiveTimStep;
293 double mMinDeltaTime;
294 double mMaxDeltaTime;
296 double mReductionFactor;
297 double mIncreaseFactor;
299 double mErrorTolerance;
301 unsigned int mMinIterations;
302 unsigned int mMaxIterations;
303 unsigned int mNumberOfConstantSteps;
312 void SetAdaptiveTimeParameters(
Parameters CustomParameters)
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
332 mMinDeltaTime = CustomParameters[
"minimum_time_step"].
GetDouble();
333 mMaxDeltaTime = CustomParameters[
"maximum_time_step"].
GetDouble();
335 mReductionFactor = CustomParameters[
"reduction_factor"].
GetDouble();
336 mIncreaseFactor = CustomParameters[
"increase_factor"].
GetDouble();
338 mErrorTolerance = CustomParameters[
"error_tolerance"].
GetDouble();
340 mMinIterations = CustomParameters[
"minimum_iterations"].
GetInt();
341 mMaxIterations = CustomParameters[
"maximum_iterations"].
GetInt();
342 mNumberOfConstantSteps = CustomParameters[
"number_constant_steps"].
GetInt();
386 rOStream << std::endl;
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
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