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.
apply_constant_scalarvalue_process.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Riccardo Rossi
11 //
12 //
13 
14 #if !defined(KRATOS_APPLY_CONSTANT_VALUE_PROCESS_H_INCLUDED )
15 #define KRATOS_APPLY_CONSTANT_VALUE_PROCESS_H_INCLUDED
16 
17 
18 
19 // System includes
20 #include <string>
21 #include <iostream>
22 
23 
24 // External includes
25 
26 
27 // Project includes
28 #include "includes/define.h"
29 #include "includes/kratos_flags.h"
31 #include "processes/process.h"
33 
34 namespace Kratos
35 {
36 
39 
41 
43 class KRATOS_API(KRATOS_CORE) ApplyConstantScalarValueProcess : public Process
44 {
45 public:
48  KRATOS_DEFINE_LOCAL_FLAG(VARIABLE_IS_FIXED);
49 
50 
53 
58  Parameters rParameters
59  ) : Process(Flags()) , mr_model_part(model_part)
60  {
62 
63 //only include validation with c++11 since raw_literals do not exist in c++03
64 
65 
66  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
67  // So that an error is thrown if they don't exist
68  rParameters["value"];
69  rParameters["variable_name"];
70  rParameters["model_part_name"];
71 
72  // Now validate agains defaults -- this also ensures no type mismatch
73 
74  rParameters.ValidateAndAssignDefaults(GetDefaultParameters());
75 
76  mmesh_id = rParameters["mesh_id"].GetInt();
77  mvariable_name = rParameters["variable_name"].GetString();
78  this->Set( VARIABLE_IS_FIXED, rParameters["is_fixed"].GetBool());
79 
80  if( KratosComponents< Variable<double> >::Has( mvariable_name ) ) //case of double variable
81  {
82  mdouble_value = rParameters["value"].GetDouble();
83 
84  if( model_part.GetNodalSolutionStepVariablesList().Has( KratosComponents< Variable<double> >::Get( mvariable_name ) ) == false )
85  {
86  KRATOS_THROW_ERROR(std::runtime_error,"trying to fix a variable that is not in the model_part - variable name is ",mvariable_name);
87  }
88  }
89  else if( KratosComponents< Variable<int> >::Has( mvariable_name ) ) //case of int variable
90  {
91  mint_value = rParameters["value"].GetInt();
92 
93  if( model_part.GetNodalSolutionStepVariablesList().Has( KratosComponents< Variable<int> >::Get( mvariable_name ) ) == false )
94  {
95  KRATOS_THROW_ERROR(std::runtime_error,"trying to fix a variable that is not in the model_part - variable name is ",mvariable_name);
96  }
97 
98  if(this->Is(VARIABLE_IS_FIXED))
99  {
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","");
101  }
102  }
103  else if( KratosComponents< Variable<bool> >::Has( mvariable_name ) ) //case of bool variable
104  {
105  mbool_value = rParameters["value"].GetBool();
106 
107  if( model_part.GetNodalSolutionStepVariablesList().Has( KratosComponents< Variable<bool> >::Get( mvariable_name ) ) == false )
108  {
109  KRATOS_THROW_ERROR(std::runtime_error,"trying to fix a variable that is not in the model_part - variable name is ",mvariable_name);
110  }
111 
112  if(this->Is(VARIABLE_IS_FIXED))
113  {
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","");
115  }
116  }
117 
118  KRATOS_CATCH("");
119  }
120 
122  const Variable<double>& rVariable,
123  const double double_value,
124  std::size_t mesh_id,
125  const Flags options
126  ) : Process(options) , mr_model_part(model_part),mdouble_value(double_value), mint_value(0), mbool_value(false),mmesh_id(mesh_id)
127  {
128  KRATOS_TRY;
129 
130  if(this->IsDefined(VARIABLE_IS_FIXED) == false )
131  {
132  KRATOS_THROW_ERROR(std::runtime_error,"please specify if the variable is to be fixed or not (flag VARIABLE_IS_FIXED)","");
133  }
134 
135  if( model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) == false )
136  {
137  KRATOS_THROW_ERROR(std::runtime_error,"trying to fix a variable that is not in the model_part - variable name is ",rVariable);
138  }
139 
140  mvariable_name = rVariable.Name();
141 
142  KRATOS_CATCH("");
143  }
144 
146  const Variable< int >& rVariable,
147  const int int_value,
148  std::size_t mesh_id,
149  const Flags options
150  ) : Process(options) , mr_model_part(model_part),mdouble_value(0.0), mint_value(int_value), mbool_value(false),mmesh_id(mesh_id)
151  {
152  KRATOS_TRY;
153 
154  if(this->IsDefined(VARIABLE_IS_FIXED) == false )
155  {
156  KRATOS_THROW_ERROR(std::runtime_error,"Please specify if the variable is to be fixed or not (flag VARIABLE_IS_FIXED)","");
157  }
158  if(this->Is(VARIABLE_IS_FIXED))
159  {
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","");
161  }
162 
163  if( model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) == false )
164  {
165  KRATOS_THROW_ERROR(std::runtime_error,"Trying to fix a variable that is not in the model_part - variable name is ",rVariable);
166  }
167 
168  mvariable_name = rVariable.Name();
169 
170  KRATOS_CATCH("");
171  }
172 
174  const Variable< bool >& rVariable,
175  const bool bool_value,
176  std::size_t mesh_id,
177  const Flags options
178  ) : Process(options) , mr_model_part(model_part),mdouble_value(0.0), mint_value(0), mbool_value(bool_value),mmesh_id(mesh_id)
179  {
180  KRATOS_TRY;
181 
182  if(this->IsDefined(VARIABLE_IS_FIXED) == false )
183  {
184  KRATOS_THROW_ERROR(std::runtime_error,"Please specify if the variable is to be fixed or not (flag VARIABLE_IS_FIXED)","");
185  }
186  if(this->Is(VARIABLE_IS_FIXED))
187  {
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","");
189  }
190 
191  if( model_part.GetNodalSolutionStepVariablesList().Has( rVariable ) == false )
192  {
193  KRATOS_THROW_ERROR(std::runtime_error,"Trying to fix a variable that is not in the model_part - variable name is ",rVariable);
194  }
195 
196  mvariable_name = rVariable.Name();
197 
198  KRATOS_CATCH("");
199  }
200 
201 
204 
205 
209 
211  void operator()()
212  {
213  Execute();
214  }
215 
216  const Parameters GetDefaultParameters() const override
217  {
218  const Parameters default_parameters( R"(
219  {
220  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
221  "mesh_id": 0,
222  "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
223  "is_fixed": false,
224  "value" : 1.0
225  } )" );
226  return default_parameters;
227  }
228 
232 
233 
235  void Execute() override {}
236 
239  void ExecuteInitialize() override
240  {
241  KRATOS_TRY;
242  const bool is_fixed = this->Is(VARIABLE_IS_FIXED);
243 
244  if( KratosComponents< Variable<double> >::Has( mvariable_name ) ) //case of double variable
245  {
246  InternalApplyValue<>(KratosComponents< Variable<double> >::Get(mvariable_name) , is_fixed, mdouble_value);
247  }
248  else if( KratosComponents< Variable<int> >::Has( mvariable_name ) ) //case of int variable
249  {
250  InternalApplyValueWithoutFixing<>(KratosComponents< Variable<int> >::Get(mvariable_name) , mint_value);
251  }
252  else if( KratosComponents< Variable<bool> >::Has( mvariable_name ) ) //case of bool variable
253  {
254  InternalApplyValueWithoutFixing<>(KratosComponents< Variable<bool> >::Get(mvariable_name), mbool_value);
255  }
256  else
257  {
258  KRATOS_THROW_ERROR(std::logic_error, "Not able to fix the variable. Attempting to fix variable:",mvariable_name);
259  }
260 
261  KRATOS_CATCH("");
262  }
263 
267  {
268  }
269 
270 
273  {
274  }
275 
278  {
279  }
280 
281 
283  void ExecuteBeforeOutputStep() override
284  {
285  }
286 
287 
289  void ExecuteAfterOutputStep() override
290  {
291  }
292 
293 
296  void ExecuteFinalize() override
297  {
298  }
299 
300 
304 
305 
309 
310 
314 
316  std::string Info() const override
317  {
318  return "ApplyConstantScalarValueProcess";
319  }
320 
322  void PrintInfo(std::ostream& rOStream) const override
323  {
324  rOStream << "ApplyConstantScalarValueProcess";
325  }
326 
328  void PrintData(std::ostream& rOStream) const override
329  {
330  }
331 
332 
336 
337 
339 protected:
340 
342  std::string mvariable_name;
346  std::size_t mmesh_id;
347 
348 private:
351  template< class TVarType, class TDataType >
352  void InternalApplyValue(const TVarType& rVar, const bool to_be_fixed, const TDataType value)
353  {
354  const int nnodes = mr_model_part.GetMesh(mmesh_id).Nodes().size();
355 
356  if(nnodes != 0)
357  {
358  block_for_each(mr_model_part.GetMesh(mmesh_id).Nodes(), [&](Node& rNode){
359  if(to_be_fixed)
360  {
361  rNode.Fix(rVar);
362  }
363  rNode.FastGetSolutionStepValue(rVar) = value;
364  });
365  }
366  }
367 
368  template< class TVarType, class TDataType >
369  void InternalApplyValueWithoutFixing(const TVarType& rVar, const TDataType value)
370  {
371  const int nnodes = mr_model_part.GetMesh(mmesh_id).Nodes().size();
372 
373  if(nnodes != 0)
374  {
375  VariableUtils().SetVariable(rVar, value, mr_model_part.GetMesh(mmesh_id).Nodes());
376  }
377  }
378 
382 
384  ApplyConstantScalarValueProcess& operator=(ApplyConstantScalarValueProcess const& rOther);
385 
387  //ApplyConstantScalarValueProcess(ApplyConstantScalarValueProcess const& rOther);
388 
389 
391 
392 }; // Class ApplyConstantScalarValueProcess
393 
395 
398 
399 
403 
404 
406 inline std::istream& operator >> (std::istream& rIStream,
408 
410 inline std::ostream& operator << (std::ostream& rOStream,
411  const ApplyConstantScalarValueProcess& rThis)
412 {
413  rThis.PrintInfo(rOStream);
414  rOStream << std::endl;
415  rThis.PrintData(rOStream);
416 
417  return rOStream;
418 }
420 
421 
422 } // namespace Kratos.
423 
424 #endif // KRATOS_APPLY_CONSTANT_VALUE_PROCESS_H_INCLUDED defined
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
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
Definition: flags.h:58
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