14 #if !defined(KRATOS_AUTOMATIC_DT_PROCESS )
15 #define KRATOS_AUTOMATIC_DT_PROCESS
47 mrModelPart(rModelPart)
53 "model_part_name":"MODEL_PART_NAME",
54 "correction_factor" : 5.0e-3
60 mCorrectionFactor = rParameters[
"correction_factor"].
GetDouble();
79 void ExecuteBeforeSolutionLoop()
override;
82 std::string
Info()
const override
84 return "AutomaticDTProcess";
88 void PrintInfo(std::ostream& rOStream)
const override
90 rOStream <<
"AutomaticDTProcess";
94 void PrintData(std::ostream& rOStream)
const override
128 rOStream << std::endl;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: automatic_dt_process.hpp:34
KRATOS_CLASS_POINTER_DEFINITION(AutomaticDTProcess)
std::string Info() const override
Turn back information as a string.
Definition: automatic_dt_process.hpp:82
AutomaticDTProcess(ModelPart &rModelPart, Parameters rParameters)
Constructor.
Definition: automatic_dt_process.hpp:43
void Execute() override
Execute method is used to execute the AutomaticDTProcess algorithms.
Definition: automatic_dt_process.hpp:73
ModelPart & mrModelPart
Member Variables.
Definition: automatic_dt_process.hpp:104
double mCorrectionFactor
Definition: automatic_dt_process.hpp:105
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: automatic_dt_process.hpp:88
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: automatic_dt_process.hpp:94
~AutomaticDTProcess() override
Destructor.
Definition: automatic_dt_process.hpp:68
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
double GetDouble() const
This method returns the double contained in the current Parameter.
Definition: kratos_parameters.cpp:657
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
The base class for all processes in Kratos.
Definition: process.h:49
#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