8 #if !defined(KRATOS_APPLY_COMPONENT_TABLE_PROCESS )
9 #define KRATOS_APPLY_COMPONENT_TABLE_PROCESS
43 "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
44 "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
57 rParameters[
"variable_name"];
58 rParameters[
"model_part_name"];
67 unsigned int TableId = rParameters[
"table"].
GetInt();
103 #pragma omp parallel for
106 ModelPart::NodesContainerType::iterator it = it_begin +
i;
110 it->Fix(var_component);
129 double value =
mpTable->GetValue(Time);
137 #pragma omp parallel for
140 ModelPart::NodesContainerType::iterator it = it_begin +
i;
142 it->FastGetSolutionStepValue(var_component) = value;
150 std::string
Info()
const override
152 return "ApplyComponentTableProcessDam";
158 rOStream <<
"ApplyComponentTableProcessDam";
200 rOStream << std::endl;
Definition: apply_component_table_process.hpp:22
ApplyComponentTableProcessDam(ModelPart &model_part, Parameters rParameters)
Constructor.
Definition: apply_component_table_process.hpp:34
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: apply_component_table_process.hpp:121
KRATOS_CLASS_POINTER_DEFINITION(ApplyComponentTableProcessDam)
~ApplyComponentTableProcessDam() override
Destructor.
Definition: apply_component_table_process.hpp:77
double minitial_value
Definition: apply_component_table_process.hpp:175
TableType::Pointer mpTable
Definition: apply_component_table_process.hpp:176
std::string mvariable_name
Definition: apply_component_table_process.hpp:173
bool mis_fixed
Definition: apply_component_table_process.hpp:174
Table< double, double > TableType
Defining a table with double argument and result type as table type.
Definition: apply_component_table_process.hpp:29
ModelPart & mr_model_part
Member Variables.
Definition: apply_component_table_process.hpp:172
void Execute() override
Execute method is used to execute the ApplyComponentTableProcessDam algorithms.
Definition: apply_component_table_process.hpp:82
void ExecuteInitialize() override
Definition: apply_component_table_process.hpp:90
std::string Info() const override
Turn back information as a string.
Definition: apply_component_table_process.hpp:150
double mTimeUnitConverter
Definition: apply_component_table_process.hpp:177
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: apply_component_table_process.hpp:162
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: apply_component_table_process.hpp:156
static const TComponentType & Get(const std::string &rName)
Retrieves a component with the specified name.
Definition: kratos_components.h:114
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
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
std::string GetString() const
This method returns the string contained in the current Parameter.
Definition: kratos_parameters.cpp:684
bool GetBool() const
This method returns the boolean contained in the current Parameter.
Definition: kratos_parameters.cpp:675
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
model_part
Definition: face_heat.py:14
int nnodes
Definition: sensitivityMatrix.py:24
integer i
Definition: TensorModule.f:17