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_properties_to_nodes_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: July 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_ASSIGN_PROPERTIES_TO_NODES_PROCESS_H_INCLUDED)
11 #define KRATOS_ASSIGN_PROPERTIES_TO_NODES_PROCESS_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 #include "includes/model_part.h"
20 #include "processes/process.h"
21 
23 
24 namespace Kratos
25 {
26 
29 
31 
34 {
35 public:
38 
40  typedef Node NodeType;
41 
43  typedef typename PropertiesContainerType::Pointer PropertiesContainerPointerType;
44 
48 
53  Parameters rParameters
54  ) : Process() , mrModelPart(model_part)
55  {
57 
58  Parameters default_parameters( R"(
59  {
60  "fluid_mixture": false,
61  "solid_mixture": false
62  } )" );
63 
64 
65  // Validate against defaults -- this ensures no type mismatch
66  rParameters.ValidateAndAssignDefaults(default_parameters);
67 
68  mFluidMixture = rParameters["fluid_mixture"].GetBool();
69  mSolidMixture = rParameters["solid_mixture"].GetBool();
70 
71  mpProperties = mrModelPart.pProperties();
72 
73  KRATOS_CATCH("");
74  }
75 
76 
79 
80 
84 
86  void operator()()
87  {
88  Execute();
89  }
90 
91 
95 
96 
98  void Execute() override
99  {
100  }
101 
104  void ExecuteInitialize() override
105  {
106  KRATOS_TRY
107 
108  this->AssignPropertiesToNodes();
109 
110  KRATOS_CATCH("")
111  }
112 
116  {
117  }
118 
119 
122  {
123  KRATOS_TRY
124 
125  this->AssignMaterialPercentageToNodes();
126 
127  KRATOS_CATCH("")
128  }
129 
132  {
133  }
134 
135 
137  void ExecuteBeforeOutputStep() override
138  {
139  }
140 
141 
143  void ExecuteAfterOutputStep() override
144  {
145  }
146 
147 
150  void ExecuteFinalize() override
151  {
152  }
153 
154 
158 
159 
163 
164 
168 
170  std::string Info() const override
171  {
172  return "AssignPropertiesToNodesProcess";
173  }
174 
176  void PrintInfo(std::ostream& rOStream) const override
177  {
178  rOStream << "AssignPropertiesToNodesProcess";
179  }
180 
182  void PrintData(std::ostream& rOStream) const override
183  {
184  }
185 
186 
191 
192 protected:
193 
202 
205 
219 
220 private:
221 
227 
228  ModelPart& mrModelPart;
229 
230  bool mFluidMixture;
231 
232  bool mSolidMixture;
233 
234  PropertiesContainerPointerType mpProperties;
235 
239 
240 
241  void AssignPropertiesToNodes()
242  {
243  const int nnodes = mrModelPart.GetMesh().Nodes().size();
244 
245  if(nnodes != 0)
246  {
247  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh().NodesBegin();
248 
249  //#pragma omp parallel for
250  for(int i = 0; i<nnodes; i++)
251  {
252  ModelPart::NodesContainerType::iterator it = it_begin + i;
253  it->SetValue(PROPERTIES_VECTOR,mpProperties);
254  }
255  }
256  }
257 
258  void AssignMaterialPercentageToNodes()
259  {
260  const int nnodes = mrModelPart.GetMesh().Nodes().size();
261 
262  if(nnodes != 0)
263  {
264  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh().NodesBegin();
265 
266  //#pragma omp parallel for
267  for(int i = 0; i<nnodes; i++)
268  {
269  ModelPart::NodesContainerType::iterator it = it_begin + i;
270  Vector MaterialPercentage;
271  this->CalculateMaterialPercentage(*it, MaterialPercentage);
272  it->SetValue(MATERIAL_PERCENTAGE, MaterialPercentage);
273  }
274  }
275  }
276 
277 
278  void CalculateMaterialPercentage(NodeType& rNode, Vector& MaterialPercentage)
279  {
280  KRATOS_TRY
281 
282  unsigned int size = mpProperties->size();
283  MaterialPercentage.resize(size,false);
284  noalias(MaterialPercentage) = ZeroVector(size);
285 
286  double counter = 0;
287  if( rNode.Is(FLUID) && mFluidMixture ){
288 
289  ElementWeakPtrVectorType& nElements = rNode.GetValue(NEIGHBOUR_ELEMENTS);
290 
291  for(auto& i_nelem : nElements)
292  {
293  if(i_nelem.Is(FLUID)){
294  unsigned int id = i_nelem.GetProperties().Id();
295  if( id < size ){
296  MaterialPercentage[id] += 1;
297  ++counter;
298  }
299  }
300  }
301  }
302  else if( rNode.Is(SOLID) && mSolidMixture ){
303 
304  ElementWeakPtrVectorType& nElements = rNode.GetValue(NEIGHBOUR_ELEMENTS);
305  for(auto& i_nelem : nElements)
306  {
307  if(i_nelem.Is(SOLID)){
308  unsigned int id = i_nelem.GetProperties().Id();
309  if( id < size ){
310  MaterialPercentage[id] += 1;
311  ++counter;
312  }
313  }
314  }
315 
316  }
317 
318  double divider = 1.0;
319  if( counter != 0 )
320  divider = 1.0/counter;
321 
322  for(unsigned int i=0; i<size; ++i)
323  MaterialPercentage[i] *= divider;
324 
325  KRATOS_CATCH("")
326  }
327 
334 
337 
338 
348 
349 }; // Class AssignPropertiesToNodesProcess
350 
352 
355 
356 
360 
361 
363 inline std::istream& operator >> (std::istream& rIStream,
365 
367 inline std::ostream& operator << (std::ostream& rOStream,
368  const AssignPropertiesToNodesProcess& rThis)
369 {
370  rThis.PrintInfo(rOStream);
371  rOStream << std::endl;
372  rThis.PrintData(rOStream);
373 
374  return rOStream;
375 }
377 
378 
379 } // namespace Kratos.
380 
381 #endif // KRATOS_ASSIGN_PROPERTIES_TO_NODES_PROCESS_H_INCLUDED defined
The base class for fixing scalar variable Dof or array_1d component Dof processes in Kratos.
Definition: assign_properties_to_nodes_process.hpp:34
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: assign_properties_to_nodes_process.hpp:176
void ExecuteInitializeSolutionStep() override
this function will be executed at every time step BEFORE performing the solve phase
Definition: assign_properties_to_nodes_process.hpp:121
void ExecuteBeforeOutputStep() override
this function will be executed at every time step BEFORE writing the output
Definition: assign_properties_to_nodes_process.hpp:137
void ExecuteFinalizeSolutionStep() override
this function will be executed at every time step AFTER performing the solve phase
Definition: assign_properties_to_nodes_process.hpp:131
void ExecuteBeforeSolutionLoop() override
Definition: assign_properties_to_nodes_process.hpp:115
void ExecuteInitialize() override
Definition: assign_properties_to_nodes_process.hpp:104
void ExecuteFinalize() override
Definition: assign_properties_to_nodes_process.hpp:150
PointerVectorSet< Properties, IndexedObject > PropertiesContainerType
Definition: assign_properties_to_nodes_process.hpp:42
void ExecuteAfterOutputStep() override
this function will be executed at every time step AFTER writing the output
Definition: assign_properties_to_nodes_process.hpp:143
AssignPropertiesToNodesProcess(AssignPropertiesToNodesProcess const &rOther)
Copy constructor.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: assign_properties_to_nodes_process.hpp:182
GlobalPointersVector< Element > ElementWeakPtrVectorType
Definition: assign_properties_to_nodes_process.hpp:45
void Execute() override
Execute method is used to execute the AssignPropertiesToNodesProcess algorithms.
Definition: assign_properties_to_nodes_process.hpp:98
virtual ~AssignPropertiesToNodesProcess()
Destructor.
Definition: assign_properties_to_nodes_process.hpp:78
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: assign_properties_to_nodes_process.hpp:86
Node NodeType
NodeType.
Definition: assign_properties_to_nodes_process.hpp:40
PropertiesContainerType::Pointer PropertiesContainerPointerType
Definition: assign_properties_to_nodes_process.hpp:43
AssignPropertiesToNodesProcess(ModelPart &model_part, Parameters rParameters)
Definition: assign_properties_to_nodes_process.hpp:52
std::string Info() const override
Turn back information as a string.
Definition: assign_properties_to_nodes_process.hpp:170
KRATOS_CLASS_POINTER_DEFINITION(AssignPropertiesToNodesProcess)
Pointer definition of AssignPropertiesToNodesProcess.
bool Is(Flags const &rOther) const
Definition: flags.h:274
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
PropertiesContainerType::Pointer pProperties(IndexType ThisIndex=0)
Definition: model_part.h:1008
This class defines the node.
Definition: node.h:65
TVariableType::Type & GetValue(const TVariableType &rThisVariable)
Definition: node.h:466
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
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
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
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
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
model_part
Definition: face_heat.py:14
int counter
Definition: script_THERMAL_CORRECT.py:218
int nnodes
Definition: sensitivityMatrix.py:24
integer i
Definition: TensorModule.f:17