KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
add_dofs_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: June 2017 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_ADD_DOFS_PROCESS_H_INCLUDED)
11 #define KRATOS_ADD_DOFS_PROCESS_H_INCLUDED
12 
13 
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
20 #include "includes/model_part.h"
22 #include "processes/process.h"
23 
24 namespace Kratos
25 {
26 
29 
31 
33 class AddDofsProcess : public Process
34 {
35 public:
38 
42 
45 
50  Parameters rParameters
51  ) : Process() , mrModelPart(model_part)
52  {
54 
55  Parameters default_parameters( R"(
56  {
57  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
58  "variables_list": [],
59  "reactions_list": []
60 
61  } )" );
62 
63 
64  // Validate against defaults -- this ensures no type mismatch
65  rParameters.ValidateAndAssignDefaults(default_parameters);
66 
67  // Check variables vs reactions consistency
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;
70 
71 
72  for(unsigned int i=0; i<rParameters["variables_list"].size(); i++)
73  {
74  if( !rParameters["variables_list"][i].IsString() )
75  KRATOS_ERROR << "variables_list contains a non-string variable name "<<std::endl;
76 
77  std::string variable_name = rParameters["variables_list"][i].GetString();
78 
79  bool supplied_reaction = true;
80  if(rParameters["reactions_list"][i].IsNull())
81  supplied_reaction = false;
82 
83  if( KratosComponents< VectorVariableType >::Has( variable_name ) ){ //case of array_1d (vector with components) variable
84 
85  const VectorVariableType& VectorVariable = KratosComponents< VectorVariableType >::Get(variable_name);
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;
88  }
89  else{
90  for(unsigned int j=0; j<3; j++)
91  {
92  std::string component_name = variable_name;
93  component_name += ms_components[j];
94  const ComponentType& ComponentVariable = KratosComponents< ComponentType >::Get(component_name);
95 
96  if(supplied_reaction){
97  std::string reaction_component_name = rParameters["reactions_list"][i].GetString();
98  reaction_component_name += ms_components[j];
99  const ComponentType& ReactionComponentVariable = KratosComponents< ComponentType >::Get(reaction_component_name);
100  m_component_variables_list.push_back(&ComponentVariable);
101  m_component_reactions_list.push_back(&ReactionComponentVariable);
102  }
103  else{
104  m_component_variables_no_reaction_list.push_back(&ComponentVariable);
105  }
106 
107  }
108  }
109  }
110  else if( KratosComponents< ComponentType >::Has(variable_name) ){ //case of component variable
111 
112  const ComponentType& ComponentVariable = KratosComponents< ComponentType >::Get(variable_name);
113 
114  if( model_part.GetNodalSolutionStepVariablesList().Has( ComponentVariable.GetSourceVariable() ) == false ){
115 
116  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
117  }
118  else{
119 
120  if(supplied_reaction){
121  std::string reaction_name = rParameters["reactions_list"][i].GetString();
122  const ComponentType& ReactionComponentVariable = KratosComponents< ComponentType >::Get(reaction_name);
123  m_component_variables_list.push_back(&ComponentVariable);
124  m_component_reactions_list.push_back(&ReactionComponentVariable);
125  }
126  else{
127  m_component_variables_no_reaction_list.push_back(&ComponentVariable);
128  }
129 
130  }
131 
132 
133  }
134  else if( KratosComponents< ScalarVariableType >::Has( variable_name ) ){ //case of double variable
135 
136  const ScalarVariableType& ScalarVariable = KratosComponents< ScalarVariableType >::Get( variable_name );
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;
139  }
140  else{
141 
142  if(supplied_reaction){
143  std::string reaction_name = rParameters["reactions_list"][i].GetString();
144  const ScalarVariableType& ReactionVariable = KratosComponents< ScalarVariableType >::Get(reaction_name);
145  m_scalar_variables_list.push_back(&ScalarVariable);
146  m_scalar_reactions_list.push_back(&ReactionVariable);
147  }
148  else{
149  m_scalar_variables_no_reaction_list.push_back(&ScalarVariable);
150  }
151 
152  }
153 
154  }
155  else{
156  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
157  }
158  }
159 
160 
161  KRATOS_CATCH("")
162  }
163 
164 
166  const pybind11::list& rVariablesList,
167  const pybind11::list& rReactionsList
168  ) : Process(), mrModelPart(model_part)
169  {
170  KRATOS_TRY
171 
172  unsigned int number_variables = len(rVariablesList);
173  unsigned int number_reactions = len(rReactionsList);
174 
175  // Check variables vs reactions consistency
176  if( number_variables != number_reactions )
177  KRATOS_ERROR << "variables_list and reactions_list has not the same number of components "<<std::endl;
178 
179  for(unsigned int i=0; i<number_variables; i++)
180  {
181 
182  //std::string variable_name = boost::python::extract<std::string>(rVariablesList[i]);
183  //std::string reaction_name = boost::python::extract<std::string>(rReactionsList[i]);
184  std::string variable_name = pybind11::cast<std::string>(rVariablesList[i]);
185  std::string reaction_name = pybind11::cast<std::string>(rReactionsList[i]);
186 
187  bool supplied_reaction = true;
188  if(reaction_name == "NOT_DEFINED")
189  supplied_reaction = false;
190 
191  if( KratosComponents< VectorVariableType >::Has( variable_name ) ){ //case of array_1d (vector with components) variable
192 
193  const VectorVariableType& VectorVariable = KratosComponents< VectorVariableType >::Get(variable_name);
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;
196  }
197  else{
198  for(unsigned int j=0; j<3; j++)
199  {
200  std::string component_name = variable_name;
201  component_name += ms_components[j];
202  const ComponentType& ComponentVariable = KratosComponents< ComponentType >::Get(component_name);
203 
204  if(supplied_reaction){
205  std::string reaction_component_name = reaction_name;
206  reaction_component_name += ms_components[j];
207  const ComponentType& ReactionComponentVariable = KratosComponents< ComponentType >::Get(reaction_component_name);
208  m_component_variables_list.push_back(&ComponentVariable);
209  m_component_reactions_list.push_back(&ReactionComponentVariable);
210  }
211  else{
212  m_component_variables_no_reaction_list.push_back(&ComponentVariable);
213  }
214 
215  }
216  }
217  }
218  else if( KratosComponents< ComponentType >::Has(variable_name) ){ //case of component variable
219 
220  const ComponentType& ComponentVariable = KratosComponents< ComponentType >::Get(variable_name);
221 
222  if( model_part.GetNodalSolutionStepVariablesList().Has( ComponentVariable.GetSourceVariable() ) == false ){
223 
224  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
225  }
226  else{
227 
228  if(supplied_reaction){
229  const ComponentType& ReactionComponentVariable = KratosComponents< ComponentType >::Get(reaction_name);
230  m_component_variables_list.push_back(&ComponentVariable);
231  m_component_reactions_list.push_back(&ReactionComponentVariable);
232  }
233  else{
234  m_component_variables_list.push_back(&ComponentVariable);
235  }
236 
237  }
238 
239  }
240  else if( KratosComponents< ScalarVariableType >::Has( variable_name ) ){ //case of double variable
241 
242  const ScalarVariableType& ScalarVariable = KratosComponents< ScalarVariableType >::Get( variable_name );
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;
245  }
246  else{
247 
248  if(supplied_reaction){
249  const ScalarVariableType& ReactionVariable = KratosComponents< ScalarVariableType >::Get(reaction_name);
250  m_scalar_variables_list.push_back(&ScalarVariable);
251  m_scalar_reactions_list.push_back(&ReactionVariable);
252  }
253  else{
254  m_scalar_variables_no_reaction_list.push_back(&ScalarVariable);
255  }
256 
257  }
258 
259  }
260  else{
261  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is "<<variable_name<<std::endl;
262  }
263  }
264 
265  KRATOS_CATCH("")
266  }
267 
268 
270  ~AddDofsProcess() override {}
271 
272 
276 
278  void operator()()
279  {
280  Execute();
281  }
282 
283 
287 
288 
290  void Execute() override
291  {
292 
293  KRATOS_TRY;
294 
295  int number_of_nodes = mrModelPart.NumberOfNodes();
296  ModelPart::NodeConstantIterator nodes_begin = mrModelPart.NodesBegin();
297 
298  /*
299  //1nd way: (fastest) generating the dofs for the initial node and add to others (still fails if a variable or a dof is set when mdpa is read)
300  AddNodalDofs(nodes_begin);
301  ModelPart::NodeType::DofsContainerType& reference_dofs = nodes_begin->GetDofs();
302  #pragma omp parallel for
303  for (int k=0; k<number_of_nodes; k++)
304  {
305  ModelPart::NodeConstantIterator it = nodes_begin + k;
306 
307  for(ModelPart::NodeType::DofsContainerType::iterator iii = reference_dofs.begin(); iii != reference_dofs.end(); iii++)
308  {
309  it->pAddDof( *iii );
310  }
311  }
312  */
313 
314  //2nd way: (faster)
315  // #pragma omp parallel for
316  for (int k=0; k<number_of_nodes; k++)
317  {
318  ModelPart::NodeConstantIterator it = nodes_begin + k;
319  AddNodalDofs(it);
320  }
321 
322 
323  /*
324  //3rt way: add dofs in the standard way one by one to all nodes (slower)
325  AddNodalDofs();
326  */
327 
328  //CheckNodalData(nodes_begin);
329 
330  KRATOS_CATCH("");
331 
332  }
333 
336  void ExecuteInitialize() override
337  {
338  }
339 
343  {
344  }
345 
346 
349  {
350  }
351 
354  {
355  }
356 
357 
359  void ExecuteBeforeOutputStep() override
360  {
361  }
362 
363 
365  void ExecuteAfterOutputStep() override
366  {
367  }
368 
369 
372  void ExecuteFinalize() override
373  {
374  }
375 
376 
380 
381 
385 
386 
390 
392  std::string Info() const override
393  {
394  return "AddDofsProcess";
395  }
396 
398  void PrintInfo(std::ostream& rOStream) const override
399  {
400  rOStream << "AddDofsProcess";
401  }
402 
404  void PrintData(std::ostream& rOStream) const override
405  {
406  }
407 
408 
413 
414 protected:
415 
424 
427 
441 
442 private:
443 
449 
450  ModelPart& mrModelPart;
451 
452  const std::vector<std::string> ms_components {"_X", "_Y", "_Z"};
453 
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;
457 
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;
461 
462 
466 
467 
468  void AddNodalDofs()
469  {
470  KRATOS_TRY
471 
472  int number_of_nodes = mrModelPart.NumberOfNodes();
473  ModelPart::NodeConstantIterator nodes_begin = mrModelPart.NodesBegin();
474 
475  for( unsigned int i=0; i < m_component_variables_list.size(); i++ )
476  {
477  #pragma omp parallel for
478  for (int k=0; k<number_of_nodes; k++)
479  {
480  ModelPart::NodeConstantIterator it = nodes_begin + k;
481  it->pAddDof(*m_component_variables_list[i],*m_component_reactions_list[i]);
482  }
483  }
484 
485  for( unsigned int j=0; j < m_component_variables_no_reaction_list.size(); j++ )
486  {
487  #pragma omp parallel for
488  for (int k=0; k<number_of_nodes; k++)
489  {
490  ModelPart::NodeConstantIterator it = nodes_begin + k;
491  it->pAddDof(*m_component_variables_no_reaction_list[j]);
492  }
493  }
494 
495  for( unsigned int l=0; l < m_scalar_variables_list.size(); l++ )
496  {
497  #pragma omp parallel for
498  for (int k=0; k<number_of_nodes; k++)
499  {
500  ModelPart::NodeConstantIterator it = nodes_begin + k;
501  it->pAddDof(*m_scalar_variables_list[l],*m_scalar_reactions_list[l]);
502  }
503  }
504 
505  for( unsigned int m=0; m < m_scalar_variables_no_reaction_list.size(); m++ )
506  {
507  #pragma omp parallel for
508  for (int k=0; k<number_of_nodes; k++)
509  {
510  ModelPart::NodeConstantIterator it = nodes_begin + k;
511  it->pAddDof(*m_scalar_variables_no_reaction_list[m]);
512  }
513  }
514 
515  KRATOS_CATCH(" ")
516  }
517 
518 
519  void AddNodalDofs( ModelPart::NodeConstantIterator& node_it )
520  {
521  KRATOS_TRY
522 
523  for( unsigned int i=0; i < m_component_variables_list.size(); i++ )
524  {
525  node_it->pAddDof(*m_component_variables_list[i],*m_component_reactions_list[i]);
526  }
527 
528  for( unsigned int j=0; j < m_component_variables_no_reaction_list.size(); j++ )
529  {
530  node_it->pAddDof(*m_component_variables_no_reaction_list[j]);
531  }
532 
533  for( unsigned int l=0; l < m_scalar_variables_list.size(); l++ )
534  {
535  node_it->pAddDof(*m_scalar_variables_list[l],*m_scalar_reactions_list[l]);
536  }
537 
538  for( unsigned int m=0; m < m_scalar_variables_no_reaction_list.size(); m++ )
539  {
540  node_it->pAddDof(*m_scalar_variables_no_reaction_list[m]);
541  }
542 
543  KRATOS_CATCH(" ")
544  }
545 
546 
547  void CheckNodalData( ModelPart::NodeConstantIterator& node_it )
548  {
549  KRATOS_TRY
550 
551  std::cout<<" CHECK VARIABLES LIST KEYS "<<std::endl;
552 
553  VariablesListDataValueContainer VariablesList = (node_it)->SolutionStepData();
554 
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;
558 
559  KRATOS_CATCH(" ")
560  }
561 
568 
570  AddDofsProcess& operator=(AddDofsProcess const& rOther);
571 
572 
582 
583 }; // Class AddDofsProcess
584 
586 
589 
590 
594 
595 
597 inline std::istream& operator >> (std::istream& rIStream,
598  AddDofsProcess& rThis);
599 
601 inline std::ostream& operator << (std::ostream& rOStream,
602  const AddDofsProcess& rThis)
603 {
604  rThis.PrintInfo(rOStream);
605  rOStream << std::endl;
606  rThis.PrintData(rOStream);
607 
608  return rOStream;
609 }
611 
612 
613 } // namespace Kratos.
614 
615 #endif // KRATOS_ADD_DOFS_PROCESS_H_INCLUDED defined
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