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_field_to_entities_process.h
Go to the documentation of this file.
1 //
2 // | / |
3 // ' / __| _` | __| _ \ __|
4 // . \ | ( | | ( |\__ `
5 // _|\_\_| \__,_|\__|\___/ ____/
6 // Multi-Physics
7 //
8 // License: BSD License
9 // Kratos default license: kratos/license.txt
10 //
11 // Main authors: Riccardo Rossi
12 // Josep Maria Carbonell
13 // Vicente Mataix Ferrandiz
14 //
15 
16 #if !defined(KRATOS_ASSIGN_SCALAR_FIELD_TO_ENTITIES_PROCESS_H_INCLUDED )
17 #define KRATOS_ASSIGN_SCALAR_FIELD_TO_ENTITIES_PROCESS_H_INCLUDED
18 
19 // System includes
20 
21 // External includes
22 
23 // Project includes
25 #include "includes/model_part.h"
27 #include "processes/process.h"
28 
29 namespace Kratos
30 {
31 
34 
48 template<class TEntity>
49 class KRATOS_API(KRATOS_CORE) AssignScalarFieldToEntitiesProcess
50  : public Process
51 {
52 public:
55 
57  typedef Node NodeType;
58 
61 
63  typedef std::size_t IndexType;
64 
66  typedef std::size_t SizeType;
67 
70 
73 
77 
84  ModelPart& rModelPart,
85  Parameters rParameters
86  );
87 
90 
94 
96  void operator()()
97  {
98  Execute();
99  }
100 
104 
108  void Execute() override;
109 
114  {
115  Execute();
116  }
117 
121  const Parameters GetDefaultParameters() const override;
122 
126 
127 
131 
132 
136 
138  std::string Info() const override
139  {
140  return "AssignScalarFieldToEntitiesProcess";
141  }
142 
144  void PrintInfo(std::ostream& rOStream) const override
145  {
146  rOStream << "AssignScalarFieldToEntitiesProcess";
147  }
148 
150  void PrintData(std::ostream& rOStream) const override
151  {
152  }
153 
158 
159 protected:
160 
169 
184 
185 private:
186 
192 
193  ModelPart& mrModelPart;
194 
196 
197  std::string mVariableName;
198 
199  std::size_t mMeshId = 0;
200 
204 
211  void CallFunction(
212  const typename TEntity::Pointer& pEntity,
213  const double Time,
214  Vector& rValue
215  );
216 
223  void CallFunctionComponents(
224  const typename TEntity::Pointer& pEntity,
225  const double Time,
226  double& rValue
227  );
228 
235  void CallFunctionLocalSystem(
236  const typename TEntity::Pointer& pEntity,
237  const double Time,
238  Vector& rValue
239  );
240 
247  void CallFunctionLocalSystemComponents(
248  const typename TEntity::Pointer& pEntity,
249  const double Time,
250  double& rValue
251  );
252 
260  void AssignTimeDependentValue(
261  const typename TEntity::Pointer& pEntity,
262  const double Time,
263  Vector& rValue,
264  const double Value
265  );
266 
272  template< class TVarType >
273  void InternalAssignValueVector(
274  TVarType& rVar,
275  const double Time
276  )
277  {
278  auto& r_entities_array = GetEntitiesContainer();
279  const SizeType number_of_entities = r_entities_array.size();
280 
281  Vector Value;
282 
283  if(number_of_entities != 0) {
284  auto it_begin = r_entities_array.begin();
285 
286  if(mpFunction->DependsOnSpace()) {
287  if(mpFunction->UseLocalSystem()) {
288  // WARNING: do not parallelize with openmp. python GIL prevents it
289  for(IndexType i = 0; i<number_of_entities; ++i) {
290  auto it_entity = it_begin + i;
291  this->CallFunctionLocalSystem(*(it_entity.base()), Time, Value);
292  it_entity->SetValue(rVar, Value);
293  }
294  } else {
295  // WARNING: do not parallelize with openmp. python GIL prevents it
296  for(IndexType i = 0; i<number_of_entities; ++i) {
297  auto it_entity = it_begin + i;
298  this->CallFunction(*(it_entity.base()), Time, Value);
299  it_entity->SetValue(rVar, Value);
300  }
301  }
302  } else { // only varies in time
303  const double TimeValue = mpFunction->CallFunction(0.0, 0.0, 0.0, Time);
304  // WARNING: do not parallelize with openmp. python GIL prevents it
305  for(IndexType i = 0; i<number_of_entities; ++i) {
306  auto it_entity = it_begin + i;
307  this->AssignTimeDependentValue(*(it_entity.base()), Time, Value, TimeValue);
308  it_entity->SetValue(rVar, Value);
309  }
310  }
311  }
312  }
313 
319  template< class TVarType >
320  void InternalAssignValueScalar(
321  TVarType& rVar,
322  const double Time
323  )
324  {
325  auto& r_entities_array = GetEntitiesContainer();
326  const SizeType number_of_entities = r_entities_array.size();
327 
328  if(number_of_entities != 0) {
329  auto it_begin = r_entities_array.begin();
330 
331  if(mpFunction->DependsOnSpace()) {
332  double Value;
333 
334  if(mpFunction->UseLocalSystem()) {
335  // WARNING: do not parallelize with openmp. python GIL prevents it
336  for(IndexType i = 0; i<number_of_entities; ++i) {
337  auto it_entity = it_begin + i;
338  this->CallFunctionLocalSystemComponents(*(it_entity.base()), Time, Value);
339  it_entity->SetValue(rVar, Value);
340  }
341  } else {
342  // WARNING: do not parallelize with openmp. python GIL prevents it
343  for(IndexType i = 0; i<number_of_entities; ++i) {
344  auto it_entity = it_begin + i;
345  this->CallFunctionComponents(*(it_entity.base()), Time, Value);
346  it_entity->SetValue(rVar, Value);
347  }
348  }
349  } else { // only varies in time
350  const double TimeValue = mpFunction->CallFunction(0.0, 0.0, 0.0, Time);
351  // WARNING: do not parallelize with openmp. python GIL prevents it
352  for(IndexType i = 0; i<number_of_entities; ++i) {
353  auto it_entity = it_begin + i;
354  it_entity->SetValue(rVar, TimeValue);
355  }
356 
357  }
358  }
359  }
360 
365  EntityContainerType& GetEntitiesContainer();
366 
373 
375  AssignScalarFieldToEntitiesProcess& operator=(AssignScalarFieldToEntitiesProcess const& rOther);
376 
377 
387 
388 }; // Class AssignScalarFieldToEntitiesProcess
389 
390 
392 
395 
396 
400 
401 
403 template<class TEntity>
404 inline std::istream& operator >> (std::istream& rIStream,
406 
408 template<class TEntity>
409 inline std::ostream& operator << (std::ostream& rOStream,
411 {
412  rThis.PrintInfo(rOStream);
413  rOStream << std::endl;
414  rThis.PrintData(rOStream);
415 
416  return rOStream;
417 }
419 
420 
421 } // namespace Kratos.
422 
423 #endif // KRATOS_ASSIGN_SCALAR_FIELD_TO_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_field_to_entities_process.h:51
std::size_t SizeType
The sizeType definition.
Definition: assign_scalar_field_to_entities_process.h:66
KRATOS_CLASS_POINTER_DEFINITION(AssignScalarFieldToEntitiesProcess)
Pointer definition of AssignScalarFieldToEntitiesProcess.
Geometry< NodeType > GeometryType
The definitionof the geometry.
Definition: assign_scalar_field_to_entities_process.h:60
Node NodeType
The definitionof the node.
Definition: assign_scalar_field_to_entities_process.h:57
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_scalar_field_to_entities_process.h:96
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_scalar_field_to_entities_process.h:144
PointerVectorSet< TEntity, IndexedObject > EntityContainerType
The container of the entities.
Definition: assign_scalar_field_to_entities_process.h:69
void ExecuteInitializeSolutionStep() override
This function will be executed at every time step BEFORE performing the solve phase.
Definition: assign_scalar_field_to_entities_process.h:113
std::string Info() const override
Turn back information as a string.
Definition: assign_scalar_field_to_entities_process.h:138
std::size_t IndexType
The IndexType definition.
Definition: assign_scalar_field_to_entities_process.h:63
~AssignScalarFieldToEntitiesProcess() override
Destructor.
Definition: assign_scalar_field_to_entities_process.h:89
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_scalar_field_to_entities_process.h:150
Geometry base class.
Definition: geometry.h:71
iterator begin()
Definition: amatrix_interface.h:241
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
The base class for all processes in Kratos.
Definition: process.h:49
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::unique_ptr< T > unique_ptr
Definition: smart_pointers.h:33
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
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