10 #if !defined(KRATOS_FIX_SCALAR_DOF_PROCESS_H_INCLUDED)
11 #define KRATOS_FIX_SCALAR_DOF_PROCESS_H_INCLUDED
33 class FixScalarDofProcess :
public Process
53 "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
54 "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME"
61 mvariable_name = rParameters[
"variable_name"].
GetString();
67 KRATOS_THROW_ERROR(std::runtime_error,
"trying to set a variable that is not in the model_part - variable name is ",mvariable_name);
75 KRATOS_THROW_ERROR(std::runtime_error,
"trying to set a variable that is not in the model_part - variable name is ",mvariable_name);
83 KRATOS_THROW_ERROR(std::runtime_error,
"trying to set a variable that is not in the model_part - variable name is ",mvariable_name);
96 if(
model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) ==
false )
98 KRATOS_THROW_ERROR(std::runtime_error,
"trying to set a variable that is not in the model_part - variable name is ",rVariable);
101 mvariable_name = rVariable.
Name();
112 if(
model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) ==
false )
114 KRATOS_THROW_ERROR(std::runtime_error,
"Trying to set a variable that is not in the model_part - variable name is ",rVariable);
117 mvariable_name = rVariable.
Name();
128 if(
model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) ==
false )
130 KRATOS_THROW_ERROR(std::runtime_error,
"Trying to set a variable that is not in the model_part - variable name is ",rVariable);
133 mvariable_name = rVariable.
Name();
171 KRATOS_THROW_ERROR(std::logic_error,
"Not able to set the variable. Attempting to set variable:",mvariable_name);
236 std::string
Info()
const override
238 return "FixScalarDofProcess";
244 rOStream <<
"FixScalarDofProcess";
295 std::string mvariable_name;
301 template<
class TVarType >
302 void InternalFixDof(TVarType& rVar)
304 const int nnodes = mrModelPart.GetMesh().Nodes().size();
308 ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh().NodesBegin();
310 #pragma omp parallel for
313 ModelPart::NodesContainerType::iterator it = it_begin +
i;
356 inline std::istream&
operator >> (std::istream& rIStream,
357 FixScalarDofProcess& rThis);
360 inline std::ostream&
operator << (std::ostream& rOStream,
361 const FixScalarDofProcess& rThis)
363 rThis.PrintInfo(rOStream);
364 rOStream << std::endl;
365 rThis.PrintData(rOStream);
The base class for fixing scalar variable Dof or array_1d component Dof processes in Kratos.
Definition: fix_scalar_dof_process.hpp:38
FixScalarDofProcess(FixScalarDofProcess const &rOther)
Copy constructor.
void ExecuteBeforeSolutionLoop() override
Definition: fix_scalar_dof_process.hpp:186
std::string Info() const override
Turn back information as a string.
Definition: fix_scalar_dof_process.hpp:236
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: fix_scalar_dof_process.hpp:242
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: fix_scalar_dof_process.hpp:209
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: fix_scalar_dof_process.hpp:192
FixScalarDofProcess(ModelPart &model_part, Parameters rParameters)
Definition: fix_scalar_dof_process.hpp:45
void Execute() override
Execute method is used to execute the FixScalarDofProcess algorithms.
Definition: fix_scalar_dof_process.hpp:165
void ExecuteInitialize() override
Definition: fix_scalar_dof_process.hpp:180
FixScalarDofProcess(ModelPart &model_part, const Variable< bool > &rVariable)
Definition: fix_scalar_dof_process.hpp:122
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: fix_scalar_dof_process.hpp:248
KRATOS_CLASS_POINTER_DEFINITION(FixScalarDofProcess)
Pointer definition of FixScalarDofProcess.
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: fix_scalar_dof_process.hpp:197
FixScalarDofProcess(ModelPart &model_part, const Variable< double > &rVariable)
Definition: fix_scalar_dof_process.hpp:90
~FixScalarDofProcess() override
Destructor.
Definition: fix_scalar_dof_process.hpp:140
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: fix_scalar_dof_process.hpp:203
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: fix_scalar_dof_process.hpp:148
void ExecuteFinalize() override
Definition: fix_scalar_dof_process.hpp:216
FixScalarDofProcess(ModelPart &model_part, const Variable< int > &rVariable)
Definition: fix_scalar_dof_process.hpp:107
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
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
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
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
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