10 #if !defined(KRATOS_ADD_DOFS_PROCESS_H_INCLUDED)
11 #define KRATOS_ADD_DOFS_PROCESS_H_INCLUDED
57 "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
68 if( rParameters[
"variables_list"].size() != rParameters[
"reactions_list"].size() )
69 KRATOS_ERROR <<
"variables_list and reactions_list has not the same number of components "<<std::endl;
72 for(
unsigned int i=0;
i<rParameters[
"variables_list"].
size();
i++)
74 if( !rParameters[
"variables_list"][
i].IsString() )
75 KRATOS_ERROR <<
"variables_list contains a non-string variable name "<<std::endl;
77 std::string variable_name = rParameters[
"variables_list"][
i].
GetString();
79 bool supplied_reaction =
true;
80 if(rParameters[
"reactions_list"][
i].IsNull())
81 supplied_reaction =
false;
86 if(
model_part.GetNodalSolutionStepVariablesList().Has( VectorVariable ) ==
false ){
87 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
90 for(
unsigned int j=0;
j<3;
j++)
92 std::string component_name = variable_name;
93 component_name += ms_components[
j];
96 if(supplied_reaction){
97 std::string reaction_component_name = rParameters[
"reactions_list"][
i].
GetString();
98 reaction_component_name += ms_components[
j];
100 m_component_variables_list.push_back(&ComponentVariable);
101 m_component_reactions_list.push_back(&ReactionComponentVariable);
104 m_component_variables_no_reaction_list.push_back(&ComponentVariable);
116 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
120 if(supplied_reaction){
121 std::string reaction_name = rParameters[
"reactions_list"][
i].
GetString();
123 m_component_variables_list.push_back(&ComponentVariable);
124 m_component_reactions_list.push_back(&ReactionComponentVariable);
127 m_component_variables_no_reaction_list.push_back(&ComponentVariable);
137 if(
model_part.GetNodalSolutionStepVariablesList().Has( ScalarVariable ) ==
false ){
138 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
142 if(supplied_reaction){
143 std::string reaction_name = rParameters[
"reactions_list"][
i].
GetString();
145 m_scalar_variables_list.push_back(&ScalarVariable);
146 m_scalar_reactions_list.push_back(&ReactionVariable);
149 m_scalar_variables_no_reaction_list.push_back(&ScalarVariable);
156 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
166 const pybind11::list& rVariablesList,
167 const pybind11::list& rReactionsList
172 unsigned int number_variables = len(rVariablesList);
173 unsigned int number_reactions = len(rReactionsList);
176 if( number_variables != number_reactions )
177 KRATOS_ERROR <<
"variables_list and reactions_list has not the same number of components "<<std::endl;
179 for(
unsigned int i=0;
i<number_variables;
i++)
184 std::string variable_name = pybind11::cast<std::string>(rVariablesList[
i]);
185 std::string reaction_name = pybind11::cast<std::string>(rReactionsList[
i]);
187 bool supplied_reaction =
true;
188 if(reaction_name ==
"NOT_DEFINED")
189 supplied_reaction =
false;
194 if(
model_part.GetNodalSolutionStepVariablesList().Has( VectorVariable ) ==
false ){
195 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
198 for(
unsigned int j=0;
j<3;
j++)
200 std::string component_name = variable_name;
201 component_name += ms_components[
j];
204 if(supplied_reaction){
205 std::string reaction_component_name = reaction_name;
206 reaction_component_name += ms_components[
j];
208 m_component_variables_list.push_back(&ComponentVariable);
209 m_component_reactions_list.push_back(&ReactionComponentVariable);
212 m_component_variables_no_reaction_list.push_back(&ComponentVariable);
224 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
228 if(supplied_reaction){
230 m_component_variables_list.push_back(&ComponentVariable);
231 m_component_reactions_list.push_back(&ReactionComponentVariable);
234 m_component_variables_list.push_back(&ComponentVariable);
243 if(
model_part.GetNodalSolutionStepVariablesList().Has( ScalarVariable ) ==
false ){
244 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
248 if(supplied_reaction){
250 m_scalar_variables_list.push_back(&ScalarVariable);
251 m_scalar_reactions_list.push_back(&ReactionVariable);
254 m_scalar_variables_no_reaction_list.push_back(&ScalarVariable);
261 KRATOS_ERROR <<
"trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
316 for (
int k=0;
k<number_of_nodes;
k++)
392 std::string
Info()
const override
394 return "AddDofsProcess";
400 rOStream <<
"AddDofsProcess";
452 const std::vector<std::string> ms_components {
"_X",
"_Y",
"_Z"};
454 std::vector<ComponentType const *> m_component_variables_list;
455 std::vector<ComponentType const *> m_component_reactions_list;
456 std::vector<ComponentType const *> m_component_variables_no_reaction_list;
458 std::vector<ScalarVariableType const *> m_scalar_variables_list;
459 std::vector<ScalarVariableType const *> m_scalar_reactions_list;
460 std::vector<ScalarVariableType const *> m_scalar_variables_no_reaction_list;
472 int number_of_nodes = mrModelPart.NumberOfNodes();
475 for(
unsigned int i=0;
i < m_component_variables_list.size();
i++ )
477 #pragma omp parallel for
478 for (
int k=0;
k<number_of_nodes;
k++)
481 it->pAddDof(*m_component_variables_list[
i],*m_component_reactions_list[
i]);
485 for(
unsigned int j=0;
j < m_component_variables_no_reaction_list.size();
j++ )
487 #pragma omp parallel for
488 for (
int k=0;
k<number_of_nodes;
k++)
491 it->pAddDof(*m_component_variables_no_reaction_list[
j]);
495 for(
unsigned int l=0;
l < m_scalar_variables_list.size();
l++ )
497 #pragma omp parallel for
498 for (
int k=0;
k<number_of_nodes;
k++)
501 it->pAddDof(*m_scalar_variables_list[
l],*m_scalar_reactions_list[
l]);
505 for(
unsigned int m=0;
m < m_scalar_variables_no_reaction_list.size();
m++ )
507 #pragma omp parallel for
508 for (
int k=0;
k<number_of_nodes;
k++)
511 it->pAddDof(*m_scalar_variables_no_reaction_list[
m]);
523 for(
unsigned int i=0;
i < m_component_variables_list.size();
i++ )
525 node_it->pAddDof(*m_component_variables_list[
i],*m_component_reactions_list[
i]);
528 for(
unsigned int j=0;
j < m_component_variables_no_reaction_list.size();
j++ )
530 node_it->pAddDof(*m_component_variables_no_reaction_list[
j]);
533 for(
unsigned int l=0;
l < m_scalar_variables_list.size();
l++ )
535 node_it->pAddDof(*m_scalar_variables_list[
l],*m_scalar_reactions_list[
l]);
538 for(
unsigned int m=0;
m < m_scalar_variables_no_reaction_list.size();
m++ )
540 node_it->pAddDof(*m_scalar_variables_no_reaction_list[
m]);
551 std::cout<<
" CHECK VARIABLES LIST KEYS "<<std::endl;
553 VariablesListDataValueContainer VariablesList = (node_it)->SolutionStepData();
555 std::cout<<
" list size "<<VariablesList.pGetVariablesList()->size()<<std::endl;
556 std::cout<<
" Variable: "<<(*VariablesList.pGetVariablesList())[0]<<std::endl;
557 std::cout<<
" end "<<std::endl;
605 rOStream << std::endl;
The base class for fixing scalar variable Dof or array_1d component Dof processes in Kratos.
Definition: add_dofs_process.hpp:34
void ExecuteFinalize() override
Definition: add_dofs_process.hpp:372
Variable< double > ScalarVariableType
Definition: add_dofs_process.hpp:40
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: add_dofs_process.hpp:353
Variable< array_1d< double, 3 > > VectorVariableType
Definition: add_dofs_process.hpp:39
void ExecuteInitialize() override
Definition: add_dofs_process.hpp:336
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: add_dofs_process.hpp:365
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: add_dofs_process.hpp:359
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: add_dofs_process.hpp:398
KRATOS_CLASS_POINTER_DEFINITION(AddDofsProcess)
Pointer definition of AddDofsProcess.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: add_dofs_process.hpp:404
AddDofsProcess(ModelPart &model_part, Parameters rParameters)
Definition: add_dofs_process.hpp:49
std::string Info() const override
Turn back information as a string.
Definition: add_dofs_process.hpp:392
Variable< double > ComponentType
Definition: add_dofs_process.hpp:41
AddDofsProcess(AddDofsProcess const &rOther)
Copy constructor.
void Execute() override
Execute method is used to execute the AddDofsProcess algorithms.
Definition: add_dofs_process.hpp:290
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: add_dofs_process.hpp:348
~AddDofsProcess() override
Destructor.
Definition: add_dofs_process.hpp:270
AddDofsProcess(ModelPart &model_part, const pybind11::list &rVariablesList, const pybind11::list &rReactionsList)
Definition: add_dofs_process.hpp:165
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: add_dofs_process.hpp:278
void ExecuteBeforeSolutionLoop() override
Definition: add_dofs_process.hpp:342
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
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
SizeType NumberOfNodes(IndexType ThisIndex=0) const
Definition: model_part.h:341
MeshType::NodeConstantIterator NodeConstantIterator
Definition: model_part.h:140
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
SizeType size() const
This method returns the total size of the current array parameter.
Definition: kratos_parameters.cpp:993
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 VariableData & GetSourceVariable() const
Definition: variable_data.h:232
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
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 k
Definition: quadrature.py:595
int j
Definition: quadrature.py:648
int m
Definition: run_marine_rain_substepping.py:8
integer i
Definition: TensorModule.f:17
integer l
Definition: TensorModule.f:17