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.
assign_scalar_variable_to_pfem_entities_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemFluidDynamicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: August 2016 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_ASSIGN_SCALAR_VARIABLES_TO_PFEM_ENTITIES_PROCESS_H_INCLUDED)
11 #define KRATOS_ASSIGN_SCALAR_VARIABLES_TO_PFEM_ENTITIES_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 KRATOS_API(PFEM_FLUID_DYNAMICS_APPLICATION) AssignScalarVariableToPfemEntitiesProcess : public Process
34 {
35 public:
38 
41 
42  enum class EntityType { NODES, CONDITIONS, ELEMENTS };
43 
44  enum class AssignmentType { DIRECT, ADDITION, SUBTRACTION, MULTIPLICATION, DIVISION };
45 
46 
50  AssignScalarVariableToPfemEntitiesProcess(ModelPart& rModelPart) : Process(Flags()), mrModelPart(rModelPart) {}
51 
53  Parameters rParameters) : Process(Flags()) , mrModelPart(rModelPart)
54  {
56 
57  Parameters default_parameters( R"(
58  {
59  "model_part_name":"MODEL_PART_NAME",
60  "variable_name": "VARIABLE_NAME",
61  "entity_type": "NODES",
62  "value" : 1.0,
63  "compound_assignment": "direct"
64  } )" );
65 
66 
67  // Validate against defaults -- this ensures no type mismatch
68  rParameters.ValidateAndAssignDefaults(default_parameters);
69 
70  mvariable_name = rParameters["variable_name"].GetString();
71 
72  if( rParameters["entity_type"].GetString() == "NODES" ){
73  mEntity = EntityType::NODES;
74  }
75  else if( rParameters["entity_type"].GetString() == "CONDITIONS" ){
76  mEntity = EntityType::CONDITIONS;
77  }
78  else if( rParameters["entity_type"].GetString() == "ELEMENTS" ){
79  mEntity = EntityType::ELEMENTS;
80  }
81  else{
82  KRATOS_ERROR <<" Entity type "<< rParameters["entity_type"].GetString() << " is not supported " << std::endl;
83  }
84 
85  if( mEntity == EntityType::NODES ){
86 
87  if( KratosComponents< Variable<double> >::Has( mvariable_name ) ) //case of double variable
88  {
89  if( rModelPart.GetNodalSolutionStepVariablesList().Has( KratosComponents< Variable<double> >::Get( mvariable_name ) ) == false )
90  {
91  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
92  }
93 
94  mdouble_value = rParameters["value"].GetDouble();
95  }
96  else if( KratosComponents< Variable<int> >::Has( mvariable_name ) ) //case of int variable
97  {
98  if( rModelPart.GetNodalSolutionStepVariablesList().Has( KratosComponents< Variable<int> >::Get( mvariable_name ) ) == false )
99  {
100  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
101  }
102 
103  mint_value = rParameters["value"].GetInt();
104 
105  }
106  else if( KratosComponents< Variable<bool> >::Has( mvariable_name ) ) //case of bool variable
107  {
108  if( rModelPart.GetNodalSolutionStepVariablesList().Has( KratosComponents< Variable<bool> >::Get( mvariable_name ) ) == false )
109  {
110  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
111  }
112 
113  mbool_value = rParameters["value"].GetBool();
114  }
115 
116  }
117  else if( mEntity == EntityType::CONDITIONS || mEntity == EntityType::ELEMENTS ){
118 
119  if( KratosComponents< Variable<double> >::Has( mvariable_name ) ) //case of double variable
120  {
121  mdouble_value = rParameters["value"].GetDouble();
122  }
123  else if( KratosComponents< Variable<int> >::Has( mvariable_name ) ) //case of int variable
124  {
125  mint_value = rParameters["value"].GetInt();
126  }
127  else if( KratosComponents< Variable<bool> >::Has( mvariable_name ) ) //case of bool variable
128  {
129  mbool_value = rParameters["value"].GetBool();
130  }
131  else{
132  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
133  }
134 
135  }
136 
137  this->SetAssignmentType(rParameters["compound_assignment"].GetString(), mAssignment);
138 
139  KRATOS_CATCH("")
140  }
141 
142 
145 
146 
150 
152  void operator()()
153  {
154  Execute();
155  }
156 
157 
161 
162 
164  void Execute() override;
165 
166 
169  void ExecuteInitialize() override
170  {
171  }
172 
176  {
177  }
178 
179 
182  {
183  }
184 
187  {
188  }
189 
190 
192  void ExecuteBeforeOutputStep() override
193  {
194  }
195 
196 
198  void ExecuteAfterOutputStep() override
199  {
200  }
201 
202 
205  void ExecuteFinalize() override;
206 
207 
211 
212 
216 
217 
221 
223  std::string Info() const override
224  {
225  return "AssignScalarVariableToPfemEntitiesProcess";
226  }
227 
229  void PrintInfo(std::ostream& rOStream) const override
230  {
231  rOStream << "AssignScalarVariableToPfemEntitiesProcess";
232  }
233 
235  void PrintData(std::ostream& rOStream) const override
236  {
237  }
238 
239 
244 
245 protected:
246 
249 
251  std::string mvariable_name;
252 
255 
262 
265 
269 
270  // set assignment method
271 
272  void SetAssignmentType(std::string method, AssignmentType& rAssignment)
273  {
274  //compound_assignment:
275 
276  //implemented:
277  // = direct
278  // += addition
279  // -= subtraction
280  // *= multiplication
281  // /= division
282 
283  if( method == "direct" ){
284  rAssignment = AssignmentType::DIRECT;
285  }
286  else if( method == "addition" ){
287  rAssignment = AssignmentType::ADDITION;
288  }
289  else if( method == "subtraction" ){
290  rAssignment = AssignmentType::SUBTRACTION;
291  }
292  else if( method == "multiplication" ){
293  rAssignment = AssignmentType::MULTIPLICATION;
294  }
295  else if( method == "division" ){
296  rAssignment = AssignmentType::DIVISION;
297  }
298  else{
299  KRATOS_ERROR <<" Assignment type "<< method << " is not supported " << std::endl;
300  }
301 
302  }
303 
304  // nodes
305 
306  template< class TVarType, class TDataType >
307  void DirectAssignValue(ModelPart::NodeType& rNode, const TVarType& rVariable, const TDataType& value)
308  {
309  rNode.FastGetSolutionStepValue(rVariable) = value;
310  }
311 
312  template< class TVarType, class TDataType >
313  void AddAssignValue(ModelPart::NodeType& rNode, const TVarType& rVariable, const TDataType& value)
314  {
315  rNode.FastGetSolutionStepValue(rVariable) += value;
316  }
317 
318  template< class TVarType, class TDataType >
319  void SubtractAssignValue(ModelPart::NodeType& rNode, const TVarType& rVariable, const TDataType& value)
320  {
321  rNode.FastGetSolutionStepValue(rVariable) -= value;
322  }
323 
324  template< class TVarType, class TDataType >
325  void MultiplyAssignValue(ModelPart::NodeType& rNode, const TVarType& rVariable, const TDataType& value)
326  {
327  rNode.FastGetSolutionStepValue(rVariable) *= value;
328  }
329 
330  template< class TVarType, class TDataType >
331  void DivideAssignValue(ModelPart::NodeType& rNode, const TVarType& rVariable, const TDataType& value)
332  {
333  if(value!=0)
334  rNode.FastGetSolutionStepValue(rVariable) /= value;
335  }
336 
337 
338  // override for the bool type (only direct assign)
339  // void AddAssignValue(ModelPart::NodeType& rNode, const Variable<bool>& rVariable, const bool& value)
340  // {
341  // rNode.FastGetSolutionStepValue(rVariable) = value;
342  // }
343 
344  // void SubtractAssignValue(ModelPart::NodeType& rNode, const Variable<bool>& rVariable, const bool& value)
345  // {
346  // rNode.FastGetSolutionStepValue(rVariable) = value;
347  // }
348 
349  // void MultiplyAssignValue(ModelPart::NodeType& rNode, const Variable<bool>& rVariable, const bool& value)
350  // {
351  // rNode.FastGetSolutionStepValue(rVariable) = value;
352  // }
353 
354  // void DivideAssignValue(ModelPart::NodeType& rNode, const Variable<bool>& rVariable, const bool& value)
355  // {
356  // rNode.FastGetSolutionStepValue(rVariable) = value;
357  // }
358 
359  // elements and conditions
360 
361  template< class TEntityType, class TVarType, class TDataType >
362  void DirectAssignValue(TEntityType& rEntity, const TVarType& rVariable, const TDataType& value)
363  {
364  rEntity.SetValue(rVariable,value);
365  }
366 
367  template< class TEntityType, class TVarType, class TDataType >
368  void AddAssignValue(TEntityType& rEntity, const TVarType& rVariable, const TDataType& value)
369  {
370  TDataType AddedValue = rEntity.GetValue(rVariable)+value;
371  rEntity.SetValue(rVariable,AddedValue);
372  }
373 
374  template< class TEntityType, class TVarType, class TDataType >
375  void SubtractAssignValue(TEntityType& rEntity, const TVarType& rVariable, const TDataType& value)
376  {
377  TDataType SubtractedValue = rEntity.GetValue(rVariable)-value;
378  rEntity.SetValue(rVariable,SubtractedValue);
379  }
380 
381  template< class TEntityType, class TVarType, class TDataType >
382  void MultiplyAssignValue(TEntityType& rEntity, const TVarType& rVariable, const TDataType value)
383  {
384  TDataType MultipliedValue = rEntity.GetValue(rVariable)*value;
385  rEntity.SetValue(rVariable,MultipliedValue);
386  }
387 
388  template< class TEntityType, class TVarType, class TDataType >
389  void DivideAssignValue(TEntityType& rEntity, const TVarType& rVariable, const TDataType& value)
390  {
391  TDataType DividedValue = rEntity.GetValue(rVariable);
392  if(value!=0)
393  DividedValue/=value;
394  rEntity.SetValue(rVariable,DividedValue);
395  }
396 
397  template< class TEntityType >
398  void MultiplyAssignValue(TEntityType& rEntity, const Variable<Vector>& rVariable, const Vector& value)
399  {
400  Vector MultipliedValue = rEntity.GetValue(rVariable);
401  for(unsigned int i=0; i<MultipliedValue.size(); ++i)
402  {
403  MultipliedValue[i]*=value[i];
404  }
405  rEntity.SetValue(rVariable,MultipliedValue);
406  }
407 
408  template< class TEntityType >
409  void DivideAssignValue(TEntityType& rEntity, const Variable<Vector>& rVariable, const Vector& value)
410  {
411  Vector DividedValue = rEntity.GetValue(rVariable);
412  for(unsigned int i=0; i<DividedValue.size(); ++i)
413  {
414  if(value[i]!=0)
415  DividedValue[i]/=value[i];
416  }
417  rEntity.SetValue(rVariable,DividedValue);
418  }
419 
420  template< class TEntityType >
421  void MultiplyAssignValue(TEntityType& rEntity, const Variable<array_1d<double,3>>& rVariable, const array_1d<double,3>& value)
422  {
423  Vector MultipliedValue = rEntity.GetValue(rVariable);
424  for(unsigned int i=0; i<3; ++i)
425  {
426  MultipliedValue[i]*=value[i];
427  }
428  rEntity.SetValue(rVariable,MultipliedValue);
429  }
430 
431  template< class TEntityType >
432  void DivideAssignValue(TEntityType& rEntity, const Variable<array_1d<double,3>>& rVariable, const array_1d<double,3>& value)
433  {
434  Vector DividedValue = rEntity.GetValue(rVariable);
435  for(unsigned int i=0; i<3; ++i)
436  {
437  if(value[i]!=0)
438  DividedValue[i]/=value[i];
439  }
440  rEntity.SetValue(rVariable,DividedValue);
441  }
442 
443  // override for the bool type (only direct assign)
444  // template<class TEntityType>
445  // void AddAssignValue(TEntityType& rEntity, const Variable<bool>& rVariable, const bool& value)
446  // {
447  // this->DirectAssignValue(rEntity,rVariable,value);
448  // }
449 
450  // template<class TEntityType>
451  // void SubtractAssignValue(TEntityType& rEntity, const Variable<bool>& rVariable, const bool& value)
452  // {
453  // this->DirectAssignValue(rEntity,rVariable,value);
454  // }
455 
456  // template<class TEntityType>
457  // void MultiplyAssignValue(TEntityType& rEntity, const Variable<bool>& rVariable, const bool value)
458  // {
459  // this->DirectAssignValue(rEntity,rVariable,value);
460  // }
461 
462  // template<class TEntityType>
463  // void DivideAssignValue(TEntityType& rEntity, const Variable<bool>& rVariable, const bool& value)
464  // {
465  // this->DirectAssignValue(rEntity,rVariable,value);
466  // }
467 
468  template< class TMethodPointerType >
469  TMethodPointerType GetAssignmentMethod()
470  {
471  TMethodPointerType AssignmentMethod = nullptr;
472  switch( mAssignment )
473  {
474  case AssignmentType::DIRECT:
476  break;
477  case AssignmentType::ADDITION:
479  break;
480  case AssignmentType::SUBTRACTION:
482  break;
483  case AssignmentType::MULTIPLICATION:
485  break;
486  case AssignmentType::DIVISION:
488  break;
489  default:
490  KRATOS_ERROR << "Unexpected value for Assignment method " << std::endl;
491  }
492  return AssignmentMethod;
493  }
494 
505 
506 private:
507 
513 
514  bool mbool_value;
515  int mint_value;
516  double mdouble_value;
517 
521 
522  template< class TVarType, class TDataType >
523  void AssignValueToNodes(const TVarType& rVariable, const TDataType value);
524 
525  template< class TVarType, class TDataType >
526  void AssignValueToConditions(const TVarType& rVariable, const TDataType value);
527 
528  template< class TVarType, class TDataType >
529  void AssignValueToElements(const TVarType& rVariable, const TDataType value);
530 
537 
540 
541 
551 
552 }; // Class AssignScalarVariableToPfemEntitiesProcess
553 
554 
556 
559 
563 
564 
566 inline std::istream& operator >> (std::istream& rIStream,
568 
570 inline std::ostream& operator << (std::ostream& rOStream,
572 {
573  rThis.PrintInfo(rOStream);
574  rOStream << std::endl;
575  rThis.PrintData(rOStream);
576 
577  return rOStream;
578 }
580 
581 
582 } // namespace Kratos.
583 
584 #endif // KRATOS_ASSIGN_SCALAR_VARIABLES_TO_PFEM_ENTITIES_PROCESS_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:34
AssignScalarVariableToPfemEntitiesProcess(ModelPart &rModelPart, Parameters rParameters)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:52
TMethodPointerType GetAssignmentMethod()
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:469
void DivideAssignValue(TEntityType &rEntity, const Variable< array_1d< double, 3 >> &rVariable, const array_1d< double, 3 > &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:432
void SubtractAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:319
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:186
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:229
std::string Info() const override
Turn back information as a string.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:223
void ExecuteInitialize() override
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:169
AssignScalarVariableToPfemEntitiesProcess(AssignScalarVariableToPfemEntitiesProcess const &rOther)
Copy constructor.
ModelPart & mrModelPart
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:250
void MultiplyAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:382
EntityType mEntity
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:253
void ExecuteBeforeSolutionLoop() override
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:175
void DirectAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:362
void SetAssignmentType(std::string method, AssignmentType &rAssignment)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:272
void MultiplyAssignValue(TEntityType &rEntity, const Variable< array_1d< double, 3 >> &rVariable, const array_1d< double, 3 > &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:421
void DivideAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:389
void MultiplyAssignValue(TEntityType &rEntity, const Variable< Vector > &rVariable, const Vector &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:398
~AssignScalarVariableToPfemEntitiesProcess() override
Destructor.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:144
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:152
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:181
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:192
void AddAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:368
KRATOS_CLASS_POINTER_DEFINITION(AssignScalarVariableToPfemEntitiesProcess)
Pointer definition of AssignScalarVariableToPfemEntitiesProcess.
EntityType
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:42
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:235
AssignmentType
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:44
std::string mvariable_name
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:251
AssignScalarVariableToPfemEntitiesProcess(ModelPart &rModelPart)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:50
void AddAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:313
void SubtractAssignValue(TEntityType &rEntity, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:375
void DivideAssignValue(TEntityType &rEntity, const Variable< Vector > &rVariable, const Vector &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:409
void DivideAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:331
void MultiplyAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:325
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:198
void DirectAssignValue(ModelPart::NodeType &rNode, const TVarType &rVariable, const TDataType &value)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:307
AssignmentType mAssignment
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:254
Definition: flags.h:58
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
VariablesList & GetNodalSolutionStepVariablesList()
Definition: model_part.h:549
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
The base class for all processes in Kratos.
Definition: process.h:49
bool Has(const VariableData &rThisVariable) const
Definition: variables_list.h:372
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
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
integer i
Definition: TensorModule.f:17