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_vector_field_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_VECTOR_FIELD_TO_PFEM_ENTITIES_PROCESS_H_INCLUDED)
11 #define KRATOS_ASSIGN_VECTOR_FIELD_TO_PFEM_ENTITIES_PROCESS_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
22 
25 
27 
30 {
31 public:
34 
37 
39 
44  pybind11::object &rPyObject,
45  const std::string &rPyMethodName,
46  const bool SpatialFieldFunction,
47  Parameters rParameters) : BaseType(rModelPart, rPyObject, rPyMethodName, SpatialFieldFunction)
48  {
50 
51  Parameters default_parameters(R"(
52  {
53  "model_part_name":"MODEL_PART_NAME",
54  "variable_name": "VARIABLE_NAME",
55  "entity_type": "NODES",
56  "value" : [0.0, 0.0, 0.0],
57  "local_axes" : {},
58  "compound_assignment": "direct"
59  } )");
60 
61  // Validate against defaults -- this ensures no type mismatch
62  rParameters.ValidateAndAssignDefaults(default_parameters);
63 
64  this->mvariable_name = rParameters["variable_name"].GetString();
65 
66  // Admissible values for local axes, are "empty" or
67  //"local_axes" :{
68  // "origin" : [0.0, 0.0, 0.0]
69  // "axes" : [ [1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0] ]
70  // }
71 
72  this->mHasLocalOrigin = false;
73  if (rParameters["local_axes"].Has("origin"))
74  {
75  this->mHasLocalOrigin = true;
76  this->mLocalOrigin.resize(3, false);
77  for (unsigned int i = 0; i < 3; i++)
78  this->mLocalOrigin[i] = rParameters["local_axes"]["origin"][i].GetDouble();
79  }
80 
81  this->mHasLocalAxes = false;
82  if (rParameters["local_axes"].Has("axes"))
83  {
84  this->mHasLocalAxes = true;
85  this->mTransformationMatrix.resize(3, 3, false);
86  for (unsigned int i = 0; i < 3; i++)
87  for (unsigned int j = 0; j < 3; j++)
88  this->mTransformationMatrix(i, j) = rParameters["local_axes"]["axes"][i][j].GetDouble();
89  }
90 
91  if (rParameters["entity_type"].GetString() == "NODES")
92  {
93  this->mEntity = EntityType::NODES;
94  }
95  else if (rParameters["entity_type"].GetString() == "CONDITIONS")
96  {
98  }
99  else
100  {
101  KRATOS_ERROR << " Entity type " << rParameters["entity_type"].GetString() << " is not supported " << std::endl;
102  }
103 
104  if (this->mEntity == EntityType::CONDITIONS)
105  {
106 
107  if (KratosComponents<Variable<Vector>>::Has(this->mvariable_name) == false) //case of vector variable
108  {
109  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
110  }
111  else if (KratosComponents<Variable<array_1d<double, 3>>>::Has(this->mvariable_name)) //case of array_1d variable
112  {
113  KRATOS_ERROR << "trying to set a variable that is not in the model_part - variable name is " << mvariable_name << std::endl;
114  }
115  }
116  else
117  {
118  KRATOS_ERROR << " Assignment to " << rParameters["entity_type"].GetString() << " not implemented " << std::endl;
119  }
120 
121  mvector_value[0] = rParameters["value"][0].GetDouble();
122  mvector_value[1] = rParameters["value"][1].GetDouble();
123  mvector_value[2] = rParameters["value"][2].GetDouble();
124 
125  this->SetAssignmentType(rParameters["compound_assignment"].GetString(), this->mAssignment);
126 
127  KRATOS_CATCH("");
128  }
129 
132 
136 
138  void operator()()
139  {
140  Execute();
141  }
142 
146 
148  void Execute() override
149  {
150 
151  KRATOS_TRY
152 
153  ProcessInfo &rCurrentProcessInfo = mrModelPart.GetProcessInfo();
154 
155  const double &rCurrentTime = rCurrentProcessInfo[TIME];
156 
157  if (KratosComponents<Variable<Vector>>::Has(this->mvariable_name)) //case of vector variable
158  {
159 
160  Vector Value;
161  //KratosComponents< Variable<Vector> >, Vector
162  AssignValueToConditions<>(KratosComponents<Variable<Vector>>::Get(this->mvariable_name), Value, rCurrentTime);
163  }
164  else if (KratosComponents<Variable<array_1d<double, 3>>>::Has(this->mvariable_name)) //case of array_1d variable
165  {
166 
167  array_1d<double, 3> Value;
168  //KratosComponents< Variable<array_1d<double,3> > >,array_1d<double,3>
169  AssignValueToConditions<>(KratosComponents<Variable<array_1d<double, 3>>>::Get(this->mvariable_name), Value, rCurrentTime);
170  }
171  else
172  {
173  KRATOS_ERROR << "Not able to set the variable. Attempting to set variable:" << this->mvariable_name << std::endl;
174  }
175 
176  KRATOS_CATCH("");
177  }
178 
181  void ExecuteInitialize() override
182  {
183  }
184 
188  {
189  }
190 
193  {
194  }
195 
198  {
199  }
200 
202  void ExecuteBeforeOutputStep() override
203  {
204  }
205 
207  void ExecuteAfterOutputStep() override
208  {
209  }
210 
213  void ExecuteFinalize() override
214  {
215 
216  KRATOS_TRY
217 
218  if (this->mEntity == EntityType::CONDITIONS)
219  {
220 
222  if (KratosComponents<Variable<Vector>>::Has(this->mvariable_name)) //case of vector variable
223  {
224 
225  Vector Value;
226  noalias(Value) = ZeroVector(3);
227  BaseType::AssignValueToConditions<>(KratosComponents<Variable<Vector>>::Get(this->mvariable_name), Value);
228  }
229  else if (KratosComponents<Variable<array_1d<double, 3>>>::Has(this->mvariable_name)) //case of array_1d variable
230  {
231 
232  array_1d<double, 3> Value;
233  Value.clear();
234  BaseType::AssignValueToConditions<>(KratosComponents<Variable<array_1d<double, 3>>>::Get(this->mvariable_name), Value);
235  }
236  else
237  {
238  KRATOS_ERROR << "Not able to set the variable. Attempting to set variable:" << this->mvariable_name << std::endl;
239  }
240  }
241 
242  KRATOS_CATCH("");
243  }
244 
248 
252 
256 
258  std::string Info() const override
259  {
260  return "AssignVectorFieldToPfemEntitiesProcess";
261  }
262 
264  void PrintInfo(std::ostream &rOStream) const override
265  {
266  rOStream << "AssignVectorFieldToPfemEntitiesProcess";
267  }
268 
270  void PrintData(std::ostream &rOStream) const override
271  {
272  }
273 
278 
279 protected:
288 
291 
302 
303 private:
309 
310  array_1d<double, 3> mvector_value;
311 
315 
316  template <class TDataType>
317  void CallFunction(const Condition::Pointer &pCondition, const double &time, TDataType &rValue)
318  {
319 
320  Condition::GeometryType &rConditionGeometry = pCondition->GetGeometry();
321  unsigned int size = rConditionGeometry.size();
322  double value = 0;
323  unsigned int counter = 0;
324 
325  rValue.resize(size * 3, false);
326 
327  if (mIsSpatialField)
328  {
329 
330  double x = 0.0, y = 0.0, z = 0.0;
331 
332  for (unsigned int i = 0; i < size; i++)
333  {
334  this->LocalAxesTransform(rConditionGeometry[i].X(), rConditionGeometry[i].Y(), rConditionGeometry[i].Z(), x, y, z);
335 
336  value = mPyObject.attr(this->mPyMethodName.c_str())(x, y, z, time).cast<double>();
337 
338  for (unsigned int j = 0; j < 3; j++)
339  {
340  rValue[counter] = value * mvector_value[j];
341  counter++;
342  }
343  }
344  }
345  else
346  {
347 
348  value = mPyObject.attr(this->mPyMethodName.c_str())(0.0, 0.0, 0.0, time).cast<double>();
349  for (unsigned int i = 0; i < size; i++)
350  {
351  for (unsigned int j = 0; j < 3; j++)
352  {
353  rValue[counter] = value * mvector_value[j];
354  counter++;
355  }
356  }
357  }
358  }
359 
360  template <class TVarType, class TDataType>
361  void AssignValueToConditions(TVarType &rVariable, TDataType &Value, const double &rTime)
362  {
363 
364  typedef void (BaseType::*AssignmentMethodPointer)(ModelPart::ConditionType &, const TVarType &, const TDataType &);
365 
366  AssignmentMethodPointer AssignmentMethod = this->GetAssignmentMethod<AssignmentMethodPointer>();
367 
368  const int nconditions = mrModelPart.GetMesh().Conditions().size();
369 
370  if (nconditions != 0)
371  {
372  ModelPart::ConditionsContainerType::iterator it_begin = mrModelPart.GetMesh().ConditionsBegin();
373 
374  //#pragma omp parallel for //it does not work in parallel
375  for (int i = 0; i < nconditions; i++)
376  {
377  ModelPart::ConditionsContainerType::iterator it = it_begin + i;
378 
379  this->CallFunction<TDataType>(*(it.base()), rTime, Value);
380 
381  (this->*AssignmentMethod)(*it, rVariable, Value);
382  }
383  }
384  }
385 
389 
393 
396 
406 
407 }; // Class AssignVectorFieldToPfemEntitiesProcess
408 
410 
413 
417 
419 inline std::istream &operator>>(std::istream &rIStream,
421 
423 inline std::ostream &operator<<(std::ostream &rOStream,
425 {
426  rThis.PrintInfo(rOStream);
427  rOStream << std::endl;
428  rThis.PrintData(rOStream);
429 
430  return rOStream;
431 }
433 
434 } // namespace Kratos.
435 
436 #endif // KRATOS_ASSIGN_VECTOR_FIELD_TO_PFEM_ENTITIES_PROCESS_H_INCLUDED defined
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_scalar_field_to_pfem_entities_process.hpp:32
bool mIsSpatialField
Definition: assign_scalar_field_to_pfem_entities_process.hpp:328
void LocalAxesTransform(const double &rX_global, const double &rY_global, const double &rZ_global, double &rx_local, double &ry_local, double &rz_local)
Definition: assign_scalar_field_to_pfem_entities_process.hpp:344
std::string mPyMethodName
Definition: assign_scalar_field_to_pfem_entities_process.hpp:323
pybind11::object mPyObject
Definition: assign_scalar_field_to_pfem_entities_process.hpp:322
Matrix mTransformationMatrix
Definition: assign_scalar_field_to_pfem_entities_process.hpp:326
Vector mLocalOrigin
Definition: assign_scalar_field_to_pfem_entities_process.hpp:325
bool mHasLocalOrigin
Definition: assign_scalar_field_to_pfem_entities_process.hpp:330
bool mHasLocalAxes
Definition: assign_scalar_field_to_pfem_entities_process.hpp:331
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
ModelPart & mrModelPart
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:250
EntityType mEntity
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:253
void SetAssignmentType(std::string method, AssignmentType &rAssignment)
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:272
std::string mvariable_name
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:251
AssignmentType mAssignment
Definition: assign_scalar_variable_to_pfem_entities_process.hpp:254
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: assign_vector_field_to_pfem_entities_process.hpp:30
AssignScalarFieldToPfemEntitiesProcess BaseType
Definition: assign_vector_field_to_pfem_entities_process.hpp:38
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: assign_vector_field_to_pfem_entities_process.hpp:202
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: assign_vector_field_to_pfem_entities_process.hpp:207
std::string Info() const override
Turn back information as a string.
Definition: assign_vector_field_to_pfem_entities_process.hpp:258
~AssignVectorFieldToPfemEntitiesProcess() override
Destructor.
Definition: assign_vector_field_to_pfem_entities_process.hpp:131
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_vector_field_to_pfem_entities_process.hpp:264
void Execute() override
Execute method is used to execute the AssignVectorFieldToPfemEntitiesProcess algorithms.
Definition: assign_vector_field_to_pfem_entities_process.hpp:148
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: assign_vector_field_to_pfem_entities_process.hpp:192
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: assign_vector_field_to_pfem_entities_process.hpp:197
AssignVectorFieldToPfemEntitiesProcess(ModelPart &rModelPart, pybind11::object &rPyObject, const std::string &rPyMethodName, const bool SpatialFieldFunction, Parameters rParameters)
Definition: assign_vector_field_to_pfem_entities_process.hpp:43
void ExecuteBeforeSolutionLoop() override
Definition: assign_vector_field_to_pfem_entities_process.hpp:187
void ExecuteFinalize() override
Definition: assign_vector_field_to_pfem_entities_process.hpp:213
AssignVectorFieldToPfemEntitiesProcess(AssignVectorFieldToPfemEntitiesProcess const &rOther)
Copy constructor.
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_vector_field_to_pfem_entities_process.hpp:138
void ExecuteInitialize() override
Definition: assign_vector_field_to_pfem_entities_process.hpp:181
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_vector_field_to_pfem_entities_process.hpp:270
KRATOS_CLASS_POINTER_DEFINITION(AssignVectorFieldToPfemEntitiesProcess)
Pointer definition of AssignVectorFieldToPfemEntitiesProcess.
Geometry base class.
Definition: geometry.h:71
SizeType size() const
Definition: geometry.h:518
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
ConditionsContainerType & Conditions()
Definition: mesh.h:691
ConditionIterator ConditionsBegin()
Definition: mesh.h:671
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
MeshType & GetMesh(IndexType ThisIndex=0)
Definition: model_part.h:1791
Condition ConditionType
Definition: model_part.h:121
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
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
size_type size() const
Returns the number of elements in the container.
Definition: pointer_vector_set.h:502
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
BOOST_UBLAS_INLINE void clear()
Definition: array_1d.h:325
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_ERROR
Definition: exception.h:161
z
Definition: GenerateWind.py:163
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
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
time
Definition: face_heat.py:85
y
Other simbols definition.
Definition: generate_axisymmetric_navier_stokes_element.py:54
int j
Definition: quadrature.py:648
int counter
Definition: script_THERMAL_CORRECT.py:218
x
Definition: sensitivityMatrix.py:49
integer i
Definition: TensorModule.f:17