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.
integration_point_to_node_transformation_utility.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: Michael Andre, https://github.com/msandre
11 //
12 
13 #if !defined(KRATOS_INTEGRATION_POINT_TO_NODE_TRANSFORMATION_UTILITY_H_INCLUDED)
14 #define KRATOS_INTEGRATION_POINT_TO_NODE_TRANSFORMATION_UTILITY_H_INCLUDED
15 
16 // System includes
17 #include <string>
18 #include <iostream>
19 
20 // External includes
21 
22 // Project includes
23 #include "includes/define.h"
24 #include "includes/element.h"
25 #include "includes/model_part.h"
26 #include "utilities/openmp_utils.h"
27 
28 // Application includes
30 
31 namespace Kratos
32 {
33 
36 
39 
50  template<unsigned int TDim, unsigned int TNumNodes = TDim + 1>
52 
53  public:
54 
57 
60 
61  template<class TVariableType>
63  ModelPart& rModelPart) const
64  {
65 #pragma omp parallel
66  {
67  ModelPart::NodeIterator NodesBegin;
68  ModelPart::NodeIterator NodesEnd;
69  OpenMPUtils::PartitionedIterators(rModelPart.Nodes(),NodesBegin,NodesEnd);
70 
71  for (ModelPart::NodeIterator itNode = NodesBegin; itNode != NodesEnd; ++itNode)
72  {
73  itNode->FastGetSolutionStepValue(rVariable) = rVariable.Zero();
74  itNode->FastGetSolutionStepValue(NODAL_AREA) = 0.0;
75  }
76  }
77 
78 #pragma omp parallel
79  {
82  OpenMPUtils::PartitionedIterators(rModelPart.Elements(),ElemBegin,ElemEnd);
83  std::vector<TVariableType> ValuesOnIntPoint;
84 
85  for (ModelPart::ElementIterator itElem = ElemBegin; itElem != ElemEnd; ++itElem)
86  {
87  const auto& r_process_info = rModelPart.GetProcessInfo();
88  itElem->CalculateOnIntegrationPoints(rVariable,ValuesOnIntPoint,r_process_info);
89  Element::GeometryType& rGeom = itElem->GetGeometry();
90  const double Weight = rGeom.Volume() / (double) TNumNodes;
91  for (unsigned int iNode = 0; iNode < rGeom.size(); iNode++)
92  {
93  rGeom[iNode].SetLock();
94  rGeom[iNode].FastGetSolutionStepValue(rVariable) += Weight * ValuesOnIntPoint[0];
95  rGeom[iNode].FastGetSolutionStepValue(NODAL_AREA) += Weight;
96  rGeom[iNode].UnSetLock();
97  }
98  }
99  }
100 
101  rModelPart.GetCommunicator().AssembleCurrentData(rVariable);
102  rModelPart.GetCommunicator().AssembleCurrentData(NODAL_AREA);
103 
104 #pragma omp parallel
105  {
106  ModelPart::NodeIterator NodesBegin;
107  ModelPart::NodeIterator NodesEnd;
108  OpenMPUtils::PartitionedIterators(rModelPart.Nodes(),NodesBegin,NodesEnd);
109 
110  for (ModelPart::NodeIterator itNode = NodesBegin; itNode != NodesEnd; ++itNode)
111  {
112  const double NodalArea = itNode->FastGetSolutionStepValue(NODAL_AREA);
113  itNode->FastGetSolutionStepValue(rVariable) /= NodalArea;
114  }
115  }
116  }
117 
118  }; // class IntegrationPointToNodalDataTransformationUtility
119 
121 
123 
124 } // namespace Kratos
125 
126 #endif // KRATOS_INTEGRATION_POINT_TO_NODAL_DATA_TRANSFORMATION_UTILITY_H_INCLUDED defined
virtual bool AssembleCurrentData(Variable< int > const &ThisVariable)
Definition: communicator.cpp:502
Geometry base class.
Definition: geometry.h:71
SizeType size() const
Definition: geometry.h:518
virtual double Volume() const
This method calculate and return volume of this geometry.
Definition: geometry.h:1358
A utility for transforming values on integration points to nodes.
Definition: integration_point_to_node_transformation_utility.h:51
KRATOS_CLASS_POINTER_DEFINITION(IntegrationPointToNodeTransformationUtility)
Pointer definition of IntegrationPointToNodeTransformationUtility.
void TransformFromIntegrationPointsToNodes(const Variable< TVariableType > &rVariable, ModelPart &rModelPart) const
Definition: integration_point_to_node_transformation_utility.h:62
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Communicator & GetCommunicator()
Definition: model_part.h:1821
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
MeshType::ElementIterator ElementIterator
Definition: model_part.h:174
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
static void PartitionedIterators(TVector &rVector, typename TVector::iterator &rBegin, typename TVector::iterator &rEnd)
Generate a partition for an std::vector-like array, providing iterators to the begin and end position...
Definition: openmp_utils.h:179
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
const TDataType & Zero() const
This method returns the zero value of the variable type.
Definition: variable.h:346
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
TABLE_NUMBER_ANGULAR_VELOCITY TABLE_NUMBER_MOMENT I33 BEAM_INERTIA_ROT_UNIT_LENGHT_Y KRATOS_DEFINE_APPLICATION_VARIABLE(DEM_APPLICATION, double, BEAM_INERTIA_ROT_UNIT_LENGHT_Z) typedef std double
Definition: DEM_application_variables.h:182