14 #if !defined(KRATOS_APPLY_CONSTANT_VALUE_PROCESS_H_INCLUDED )
15 #define KRATOS_APPLY_CONSTANT_VALUE_PROCESS_H_INCLUDED
69 rParameters[
"variable_name"];
70 rParameters[
"model_part_name"];
76 mmesh_id = rParameters[
"mesh_id"].
GetInt();
77 mvariable_name = rParameters[
"variable_name"].
GetString();
78 this->Set( VARIABLE_IS_FIXED, rParameters[
"is_fixed"].GetBool());
82 mdouble_value = rParameters[
"value"].
GetDouble();
86 KRATOS_THROW_ERROR(std::runtime_error,
"trying to fix a variable that is not in the model_part - variable name is ",mvariable_name);
91 mint_value = rParameters[
"value"].
GetInt();
95 KRATOS_THROW_ERROR(std::runtime_error,
"trying to fix a variable that is not in the model_part - variable name is ",mvariable_name);
98 if(this->Is(VARIABLE_IS_FIXED))
100 KRATOS_THROW_ERROR(std::runtime_error,
"sorry it is not possible to fix variables of type Variable<int>. Only double variables or vector components can be fixed",
"");
105 mbool_value = rParameters[
"value"].
GetBool();
109 KRATOS_THROW_ERROR(std::runtime_error,
"trying to fix a variable that is not in the model_part - variable name is ",mvariable_name);
112 if(this->Is(VARIABLE_IS_FIXED))
114 KRATOS_THROW_ERROR(std::runtime_error,
"sorry it is not possible to fix variables of type Variable<bool>. Only double variables or vector components can be fixed",
"");
123 const double double_value,
126 ) :
Process(options) , mr_model_part(
model_part),mdouble_value(double_value), mint_value(0), mbool_value(false),mmesh_id(mesh_id)
130 if(this->IsDefined(VARIABLE_IS_FIXED) ==
false )
132 KRATOS_THROW_ERROR(std::runtime_error,
"please specify if the variable is to be fixed or not (flag VARIABLE_IS_FIXED)",
"");
135 if(
model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) ==
false )
137 KRATOS_THROW_ERROR(std::runtime_error,
"trying to fix a variable that is not in the model_part - variable name is ",rVariable);
140 mvariable_name = rVariable.
Name();
150 ) :
Process(options) , mr_model_part(
model_part),mdouble_value(0.0), mint_value(int_value), mbool_value(false),mmesh_id(mesh_id)
154 if(this->IsDefined(VARIABLE_IS_FIXED) ==
false )
156 KRATOS_THROW_ERROR(std::runtime_error,
"Please specify if the variable is to be fixed or not (flag VARIABLE_IS_FIXED)",
"");
158 if(this->Is(VARIABLE_IS_FIXED))
160 KRATOS_THROW_ERROR(std::runtime_error,
"Sorry it is not possible to fix variables of type Variable<int>. Only double variables or vector components can be fixed",
"");
163 if(
model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) ==
false )
165 KRATOS_THROW_ERROR(std::runtime_error,
"Trying to fix a variable that is not in the model_part - variable name is ",rVariable);
168 mvariable_name = rVariable.
Name();
175 const bool bool_value,
178 ) :
Process(options) , mr_model_part(
model_part),mdouble_value(0.0), mint_value(0), mbool_value(bool_value),mmesh_id(mesh_id)
182 if(this->IsDefined(VARIABLE_IS_FIXED) ==
false )
184 KRATOS_THROW_ERROR(std::runtime_error,
"Please specify if the variable is to be fixed or not (flag VARIABLE_IS_FIXED)",
"");
186 if(this->Is(VARIABLE_IS_FIXED))
188 KRATOS_THROW_ERROR(std::runtime_error,
"Sorry it is not possible to fix variables of type Variable<int>. Only double variables or vector components can be fixed",
"");
191 if(
model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) ==
false )
193 KRATOS_THROW_ERROR(std::runtime_error,
"Trying to fix a variable that is not in the model_part - variable name is ",rVariable);
196 mvariable_name = rVariable.
Name();
220 "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
222 "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
226 return default_parameters;
242 const bool is_fixed = this->Is(VARIABLE_IS_FIXED);
258 KRATOS_THROW_ERROR(std::logic_error,
"Not able to fix the variable. Attempting to fix variable:",mvariable_name);
316 std::string
Info()
const override
318 return "ApplyConstantScalarValueProcess";
324 rOStream <<
"ApplyConstantScalarValueProcess";
351 template<
class TVarType,
class TDataType >
352 void InternalApplyValue(
const TVarType& rVar,
const bool to_be_fixed,
const TDataType value)
368 template<
class TVarType,
class TDataType >
369 void InternalApplyValueWithoutFixing(
const TVarType& rVar,
const TDataType value)
371 const int nnodes = mr_model_part.GetMesh(mmesh_id).Nodes().size();
375 VariableUtils().SetVariable(rVar, value, mr_model_part.GetMesh(mmesh_id).Nodes());
384 ApplyConstantScalarValueProcess&
operator=(ApplyConstantScalarValueProcess
const& rOther);
414 rOStream << std::endl;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
The base class for all processes in Kratos.
Definition: apply_constant_scalarvalue_process.h:44
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: apply_constant_scalarvalue_process.h:211
ApplyConstantScalarValueProcess(ModelPart &model_part, const Variable< bool > &rVariable, const bool bool_value, std::size_t mesh_id, const Flags options)
Definition: apply_constant_scalarvalue_process.h:173
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: apply_constant_scalarvalue_process.h:283
void Execute() override
Execute method is used to execute the ApplyConstantScalarValueProcess algorithms.
Definition: apply_constant_scalarvalue_process.h:235
int mint_value
Definition: apply_constant_scalarvalue_process.h:344
const Parameters GetDefaultParameters() const override
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: apply_constant_scalarvalue_process.h:216
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: apply_constant_scalarvalue_process.h:272
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: apply_constant_scalarvalue_process.h:322
KRATOS_CLASS_POINTER_DEFINITION(ApplyConstantScalarValueProcess)
Pointer definition of ApplyConstantScalarValueProcess.
ApplyConstantScalarValueProcess(ModelPart &model_part, const Variable< double > &rVariable, const double double_value, std::size_t mesh_id, const Flags options)
Definition: apply_constant_scalarvalue_process.h:121
~ApplyConstantScalarValueProcess() override
Destructor.
Definition: apply_constant_scalarvalue_process.h:203
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: apply_constant_scalarvalue_process.h:277
std::string mvariable_name
Definition: apply_constant_scalarvalue_process.h:342
void ExecuteInitialize() override
Definition: apply_constant_scalarvalue_process.h:239
ApplyConstantScalarValueProcess(ModelPart &model_part, Parameters rParameters)
Definition: apply_constant_scalarvalue_process.h:57
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: apply_constant_scalarvalue_process.h:328
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: apply_constant_scalarvalue_process.h:289
KRATOS_DEFINE_LOCAL_FLAG(VARIABLE_IS_FIXED)
bool mbool_value
Definition: apply_constant_scalarvalue_process.h:345
ModelPart & mr_model_part
Definition: apply_constant_scalarvalue_process.h:341
std::string Info() const override
Turn back information as a string.
Definition: apply_constant_scalarvalue_process.h:316
double mdouble_value
Definition: apply_constant_scalarvalue_process.h:343
ApplyConstantScalarValueProcess(ModelPart &model_part, const Variable< int > &rVariable, const int int_value, std::size_t mesh_id, const Flags options)
Definition: apply_constant_scalarvalue_process.h:145
void ExecuteBeforeSolutionLoop() override
Definition: apply_constant_scalarvalue_process.h:266
std::size_t mmesh_id
Definition: apply_constant_scalarvalue_process.h:346
void ExecuteFinalize() override
Definition: apply_constant_scalarvalue_process.h:296
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
NodesContainerType & Nodes()
Definition: mesh.h:346
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType & GetMesh(IndexType ThisIndex=0)
Definition: model_part.h:1791
This class defines the node.
Definition: node.h:65
TVariableType::Type & FastGetSolutionStepValue(const TVariableType &rThisVariable)
Definition: node.h:435
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
size_type size() const
Returns the number of elements in the container.
Definition: pointer_vector_set.h:502
The base class for all processes in Kratos.
Definition: process.h:49
const std::string & Name() const
Definition: variable_data.h:201
#define KRATOS_THROW_ERROR(ExceptionType, ErrorMessage, MoreInfo)
Definition: define.h:77
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
bool Has(const std::string &ModelerName)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:24
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
void block_for_each(TIterator itBegin, TIterator itEnd, TFunction &&rFunction)
Execute a functor on all items of a range in parallel.
Definition: parallel_utilities.h:299
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