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_vectorvalue_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 
15 #if !defined(KRATOS_APPLY_CONSTANT_VECTORVALUE_PROCESS_H_INCLUDED )
16 #define KRATOS_APPLY_CONSTANT_VECTORVALUE_PROCESS_H_INCLUDED
17 
18 
19 
20 // System includes
21 #include <string>
22 #include <iostream>
23 
24 
25 // External includes
26 
27 
28 // Project includes
29 #include "includes/define.h"
30 #include "includes/kratos_flags.h"
32 #include "processes/process.h"
33 
34 namespace Kratos
35 {
36 
39 
41 
45 {
46 public:
49  KRATOS_DEFINE_LOCAL_FLAG(X_COMPONENT_FIXED);
50  KRATOS_DEFINE_LOCAL_FLAG(Y_COMPONENT_FIXED);
51  KRATOS_DEFINE_LOCAL_FLAG(Z_COMPONENT_FIXED);
52 
55 
59 
64  {
66 
67 
68  Parameters default_parameters( R"(
69  {
70  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
71  "mesh_id": 0,
72  "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
73  "is_fixed_x": false,
74  "is_fixed_y": false,
75  "is_fixed_z": false,
76  "modulus" : 1.0,
77  "direction": [1.0, 0.0, 0.0]
78  } )" );
79 
80  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
81  // So that an error is thrown if they don't exist
82  if(parameters["direction"].IsArray() == true && parameters["direction"].size() != 3)
83  {
84  KRATOS_THROW_ERROR(std::runtime_error,"direction vector is not a vector or it does not have size 3. Direction vector currently passed",parameters.PrettyPrintJsonString());
85  }
86  if(parameters["modulus"].IsNumber() == false)
87  {
88  KRATOS_THROW_ERROR(std::runtime_error,"modulus shall be a number. Parameter list in which is included is :", parameters.PrettyPrintJsonString());
89  }
90  if(parameters["variable_name"].IsString() == false)
91  {
92  KRATOS_THROW_ERROR(std::runtime_error,"vairbale_name shall be a String. Parameter list in which is included is :", parameters.PrettyPrintJsonString());
93  }
94  if(parameters["model_part_name"].IsString() == false)
95  {
96  KRATOS_THROW_ERROR(std::runtime_error,"model_part_name shall be a String. Parameter list in which is included is :", parameters.PrettyPrintJsonString());
97  }
98 
99  //now validate agains defaults -- this also ensures no type mismatch
100  parameters.ValidateAndAssignDefaults(default_parameters);
101 
102  // Read from the parameters and assign to the values
103  mmesh_id = parameters["mesh_id"].GetInt();
104 
105  this->Set(X_COMPONENT_FIXED, parameters["is_fixed_x"].GetBool());
106  this->Set(Y_COMPONENT_FIXED, parameters["is_fixed_y"].GetBool());
107  this->Set(Z_COMPONENT_FIXED, parameters["is_fixed_z"].GetBool());
108 
109  // Get the modulus and variable name
110  mvariable_name = parameters["variable_name"].GetString();
111  mmodulus = parameters["modulus"].GetDouble();
112 // mvalue = parameters["value"].GetDouble();
113 
114  mdirection.resize(3,false);
115  mdirection[0] = parameters["direction"][0].GetDouble();
116  mdirection[1] = parameters["direction"][1].GetDouble();
117  mdirection[2] = parameters["direction"][2].GetDouble();
118 
119  const double dim_norm = norm_2(mdirection);
120  if(dim_norm < 1e-20)
121  {
122  KRATOS_THROW_ERROR(std::runtime_error," Norm of direction given is approximately zero. Please give a direction vector with a non zero norm : current value of direction vector = ",mdirection);
123  }
124 
125  // Normalize the direction
126  mdirection /= dim_norm;
127 
129  {
130  KRATOS_THROW_ERROR(std::runtime_error,"Not defined the variable ",mvariable_name);
131  }
133 
134 
135  if(mmesh_id >= model_part.NumberOfMeshes())
136  {
137  KRATOS_THROW_ERROR(std::runtime_error,"Mesh does not exist in model_part: mesh id is --> ",mmesh_id);
138  }
139 
140  if( model_part.GetNodalSolutionStepVariablesList().Has(rVariable) == false )
141  {
142  std::string err_msg = std::string("Trying to fix a variable that is not in the model_part - variable: ")+mvariable_name;
143  KRATOS_THROW_ERROR(std::runtime_error,err_msg,mvariable_name);
144  }
145 
146  if(mdirection.size() != 3)
147  {
148  KRATOS_THROW_ERROR(std::runtime_error,"Direction vector is expected to have size 3. Direction vector currently passed",mdirection);
149  }
150 
151  typedef Variable<double> component_type;
152  if(KratosComponents< component_type >::Has(mvariable_name+std::string("_X")) == false)
153  {
154  KRATOS_THROW_ERROR(std::runtime_error,"Not defined the variable ",mvariable_name+std::string("_X"));
155  }
156  if(KratosComponents< component_type >::Has(mvariable_name+std::string("_Y")) == false)
157  {
158  KRATOS_THROW_ERROR(std::runtime_error,"Not defined the variable ",mvariable_name+std::string("_Y"));
159  }
160  if(KratosComponents< component_type >::Has(mvariable_name+std::string("_Z")) == false)
161  {
162  KRATOS_THROW_ERROR(std::runtime_error,"Not defined the variable ",mvariable_name+std::string("_Z"));
163  }
164 
165  KRATOS_CATCH("");
166  }
167 
168 
170  const Variable< array_1d<double, 3 > >& rVariable,
171  const double modulus,
172  const Vector& direction,
173  std::size_t mesh_id,
174  const Flags options
175  ) : Process(options) , mr_model_part(model_part), mmodulus(modulus),mdirection(direction),mmesh_id(mesh_id)
176  {
177  KRATOS_TRY;
178 
179  if(mesh_id >= model_part.NumberOfMeshes())
180  {
181  KRATOS_THROW_ERROR(std::runtime_error,"Mesh does not exist in model_part: mesh id is --> ",mesh_id);
182  }
183 
184  if(this->IsDefined(X_COMPONENT_FIXED) == false )
185  {
186  KRATOS_THROW_ERROR(std::runtime_error,"Please specify if component x is to be fixed or not (flag X_COMPONENT_FIXED)","");
187  }
188  if(this->IsDefined(Y_COMPONENT_FIXED) == false )
189  {
190  KRATOS_THROW_ERROR(std::runtime_error,"Please specify if component y is to be fixed or not (flag Y_COMPONENT_FIXED)","");
191  }
192  if(this->IsDefined(Z_COMPONENT_FIXED) == false )
193  {
194  KRATOS_THROW_ERROR(std::runtime_error,"Please specify if the variable is to be fixed or not (flag Z_COMPONENT_FIXED)","");
195  }
196 
197  mvariable_name = rVariable.Name();
198 
199  if( model_part.GetNodalSolutionStepVariablesList().Has(rVariable) == false )
200  {
201  std::string err_msg = std::string("Trying to fix a variable that is not in the model_part - variable: ")+mvariable_name;
202  KRATOS_THROW_ERROR(std::runtime_error,err_msg,mvariable_name);
203  }
204 
205  if(direction.size() != 3)
206  {
207  KRATOS_THROW_ERROR(std::runtime_error,"Direction vector is expected to have size 3. Direction vector currently passed",mdirection);
208  }
209 
210  typedef Variable<double> component_type;
211  if(KratosComponents< component_type >::Has(mvariable_name+std::string("_X")) == false)
212  {
213  KRATOS_THROW_ERROR(std::runtime_error,"Not defined the variable ",mvariable_name+std::string("_X"));
214  }
215  if(KratosComponents< component_type >::Has(mvariable_name+std::string("_Y")) == false)
216  {
217  KRATOS_THROW_ERROR(std::runtime_error,"Not defined the variable ",mvariable_name+std::string("_Y"));
218  }
219  if(KratosComponents< component_type >::Has(mvariable_name+std::string("_Z")) == false)
220  {
221  KRATOS_THROW_ERROR(std::runtime_error,"Not defined the variable ",mvariable_name+std::string("_Z"));
222  }
223 
224  KRATOS_CATCH("");
225  }
226 
227 
228 
231 
232 
236 
238  void operator()()
239  {
240  Execute();
241  }
242 
243 
247 
248 
250  void Execute() override {}
251 
254  void ExecuteInitialize() override
255  {
256  //compute the value to be applied
258  typedef Variable<double> component_type;
259 
260  const component_type& varx = KratosComponents< component_type >::Get(mvariable_name+std::string("_X"));
261  const component_type& vary = KratosComponents< component_type >::Get(mvariable_name+std::string("_Y"));
262  const component_type& varz = KratosComponents< component_type >::Get(mvariable_name+std::string("_Z"));
263 
264  InternalApplyValue<component_type >(varx, this->Is(X_COMPONENT_FIXED), value[0]);
265  InternalApplyValue<component_type >(vary, this->Is(Y_COMPONENT_FIXED), value[1]);
266  InternalApplyValue<component_type >(varz, this->Is(Z_COMPONENT_FIXED), value[2]);
267  }
268 
272  {
273  }
274 
275 
278  {
279  }
280 
283  {
284  }
285 
286 
288  void ExecuteBeforeOutputStep() override
289  {
290  }
291 
292 
294  void ExecuteAfterOutputStep() override
295  {
296  }
297 
298 
301  void ExecuteFinalize() override
302  {
303  }
304 
305 
309 
310 
314 
315 
319 
321  std::string Info() const override
322  {
323  return "ApplyConstantVectorValueProcess";
324  }
325 
327  void PrintInfo(std::ostream& rOStream) const override
328  {
329  rOStream << "ApplyConstantVectorValueProcess";
330  }
331 
333  void PrintData(std::ostream& rOStream) const override
334  {
335  }
336 
337 
341 
342 
344 protected:
345 
347  std::string mvariable_name;
348  double mmodulus;
350  std::size_t mmesh_id;
351 
352 private:
355 
356  template< class TVarType >
357  void InternalApplyValue(const TVarType& rVar, const bool to_be_fixed, const double value)
358  {
359  const int nnodes = mr_model_part.GetMesh(mmesh_id).Nodes().size();
360 
361  if(nnodes != 0)
362  {
363  ModelPart::NodesContainerType::iterator it_begin = mr_model_part.GetMesh(mmesh_id).NodesBegin();
364 // ModelPart::NodesContainerType::iterator it_end = mr_model_part.GetMesh(mmesh_id).NodesEnd();
365 
366  //check if the dofs are there (on the first node)
367  if(to_be_fixed && (it_begin->HasDofFor(rVar) == false) )
368  {
369  KRATOS_THROW_ERROR(std::runtime_error, " Trying to fix a dofs which was not allocated. Variable is --> ",rVar.Name() );
370  }
371 
373  if(to_be_fixed)
374  {
375  rNode.Fix(rVar);
376  }
377  rNode.FastGetSolutionStepValue(rVar) = value;
378  });
379  }
380  }
381 
385 
387  ApplyConstantVectorValueProcess& operator=(ApplyConstantVectorValueProcess const& rOther);
388 
390  //ApplyConstantVectorValueProcess(ApplyConstantVectorValueProcess const& rOther);
391 
392 
394 
395 }; // Class ApplyConstantVectorValueProcess
396 
400 
401 
403 
406 
407 
411 
412 
414 inline std::istream& operator >> (std::istream& rIStream,
416 
418 inline std::ostream& operator << (std::ostream& rOStream,
419  const ApplyConstantVectorValueProcess& rThis)
420 {
421  rThis.PrintInfo(rOStream);
422  rOStream << std::endl;
423  rThis.PrintData(rOStream);
424 
425  return rOStream;
426 }
428 
429 
430 } // namespace Kratos.
431 
432 #endif // KRATOS_APPLY_CONSTANT_VECTORVALUE_PROCESS_H_INCLUDED defined
433 
434 
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
The base class for all processes in Kratos.
Definition: apply_constant_vectorvalue_process.h:45
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: apply_constant_vectorvalue_process.h:333
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: apply_constant_vectorvalue_process.h:294
std::string Info() const override
Turn back information as a string.
Definition: apply_constant_vectorvalue_process.h:321
ApplyConstantVectorValueProcess(ModelPart &model_part, const Variable< array_1d< double, 3 > > &rVariable, const double modulus, const Vector &direction, std::size_t mesh_id, const Flags options)
Definition: apply_constant_vectorvalue_process.h:169
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: apply_constant_vectorvalue_process.h:282
KRATOS_CLASS_POINTER_DEFINITION(ApplyConstantVectorValueProcess)
Pointer definition of ApplyConstantVectorValueProcess.
void ExecuteInitialize() override
Definition: apply_constant_vectorvalue_process.h:254
double mmodulus
Definition: apply_constant_vectorvalue_process.h:348
std::size_t mmesh_id
Definition: apply_constant_vectorvalue_process.h:350
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: apply_constant_vectorvalue_process.h:288
ApplyConstantVectorValueProcess(ModelPart &model_part, Parameters parameters)
Default constructor.
Definition: apply_constant_vectorvalue_process.h:61
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: apply_constant_vectorvalue_process.h:238
~ApplyConstantVectorValueProcess() override
Destructor.
Definition: apply_constant_vectorvalue_process.h:230
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: apply_constant_vectorvalue_process.h:277
void ExecuteFinalize() override
Definition: apply_constant_vectorvalue_process.h:301
void Execute() override
Execute method is used to execute the ApplyConstantVectorValueProcess algorithms.
Definition: apply_constant_vectorvalue_process.h:250
ModelPart & mr_model_part
Definition: apply_constant_vectorvalue_process.h:346
Vector mdirection
Definition: apply_constant_vectorvalue_process.h:349
std::string mvariable_name
Definition: apply_constant_vectorvalue_process.h:347
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: apply_constant_vectorvalue_process.h:327
void ExecuteBeforeSolutionLoop() override
Definition: apply_constant_vectorvalue_process.h:271
Definition: flags.h:58
void Set(const Flags ThisFlag)
Definition: flags.cpp:33
bool IsDefined(Flags const &rOther) const
Definition: flags.h:279
bool Is(Flags const &rOther) const
Definition: flags.h:274
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
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
NodesContainerType & Nodes()
Definition: mesh.h:346
NodeIterator NodesBegin()
Definition: mesh.h:326
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
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#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
#define KRATOS_CREATE_LOCAL_FLAG(class_name, name, position)
Definition: define.h:685
std::istream & operator>>(std::istream &rIStream, ApplyConstantVectorValueProcess &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const ApplyConstantVectorValueProcess &rThis)
output stream function
Definition: apply_constant_vectorvalue_process.h:418
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
TExpressionType::data_type norm_2(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:625
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
model_part
Definition: face_heat.py:14
parameters
Definition: fluid_chimera_analysis.py:35
string err_msg
Definition: fluid_chimera_analysis.py:21
int nnodes
Definition: sensitivityMatrix.py:24
e
Definition: run_cpp_mpi_tests.py:31