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.
body_normal_calculation_utils.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 #include "includes/model_part.h"
16 
17 #if !defined(KRATOS_BODY_NORMAL_CALCULATION_UTILS )
18 #define KRATOS_BODY_NORMAL_CALCULATION_UTILS
19 
20 
21 /* System includes */
22 
23 
24 /* External includes */
25 
26 
27 /* Project includes */
28 //#include "utilities/math_utils.h"
29 #include "utilities/geometry_utilities.h"
30 
31 
32 namespace Kratos
33 {
34 
84 {
85 public:
112  //***********************************************************************
113  //***********************************************************************
116 
125  ModelPart& r_model_part,
126  int dimension)
127  {
128  KRATOS_TRY
129 
130  ModelPart::ElementsContainerType& rElements = r_model_part.Elements();
131 
132  //resetting the normals - only for the nodes on which we will do the calculate
134 
135  for(ModelPart::NodesContainerType::iterator it = r_model_part.NodesBegin();
136  it !=r_model_part.NodesEnd(); it++)
137  {
138  noalias(it->FastGetSolutionStepValue(NORMAL)) = zero;
139  }
140 
141 // for(ModelPart::ElementsContainerType::iterator it = rElements.begin();
142 // it !=rElements.end(); it++)
143 // {
144 // Element::GeometryType& rNodes = it->GetGeometry();
145 // for(unsigned int in = 0; in<rNodes.size(); in++)
146 // noalias((rNodes[in]).FastGetSolutionStepValue(NORMAL)) = zero;
147 // }
148 
149 
150  //calculating the normals and storing on the conditions
152  if(dimension == 2)
153  {
156  double Volume;
157  for(ModelPart::ElementsContainerType::iterator it = rElements.begin(); it !=rElements.end(); it++)
158  {
159  Element::GeometryType& geom = it->GetGeometry();
160  GeometryUtils::CalculateGeometryData(geom, DN_DX, N, Volume);
161 
162  for(unsigned int i = 0; i<geom.size(); i++)
163  {
164  array_1d<double,3>& normal = geom[i].FastGetSolutionStepValue(NORMAL);
165  for(unsigned int j=0; j<2; j++)
166  {
167  normal[j] += Volume*DN_DX(i,j);
168  }
169  }
170  }
171  }
172  else if(dimension == 3)
173  {
176  double Volume;
177  for(ModelPart::ElementsContainerType::iterator it = rElements.begin(); it !=rElements.end(); it++)
178  {
179  Element::GeometryType& geom = it->GetGeometry();
180  GeometryUtils::CalculateGeometryData(geom, DN_DX, N, Volume);
181 
182  for(unsigned int i = 0; i<geom.size(); i++)
183  {
184  array_1d<double,3>& normal = geom[i].FastGetSolutionStepValue(NORMAL);
185  for(unsigned int j=0; j<3; j++)
186  {
187  normal[j] += Volume*DN_DX(i,j);
188  }
189  }
190  }
191  }
192 
193  r_model_part.GetCommunicator().AssembleCurrentData(NORMAL);
194 // r_model_part.GetCommunicator().AssembleCurrentData(PRESSURE);
195 
196 
197  KRATOS_CATCH("")
198 
199  }
200 
201 
202 
220 private:
228  //this function adds the Contribution of one of the geometries
229  //to the corresponding nodes
230 
255  //BodyNormalCalculationUtils(void);
256 
257  //BodyNormalCalculationUtils(BodyNormalCalculationUtils& rSource);
258 
259 
262 }; /* Class ClassName */
263 
272 } /* namespace Kratos.*/
273 
274 #endif /* KRATOS_BODY_NORMAL_CALCULATION_UTILS defined */
275 
Definition: body_normal_calculation_utils.h:84
void CalculateBodyNormals(ModelPart &r_model_part, int dimension)
Definition: body_normal_calculation_utils.h:124
ModelPart::ElementsContainerType ElementsArrayType
Definition: body_normal_calculation_utils.h:89
ModelPart::NodesContainerType NodesArrayType
Definition: body_normal_calculation_utils.h:88
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
static void CalculateGeometryData(const GeometryType &rGeometry, BoundedMatrix< double, 4, 3 > &rDN_DX, array_1d< double, 4 > &rN, double &rVolume)
This function is designed to compute the shape function derivatives, shape functions and volume in 3D...
Definition: geometry_utilities.h:176
Definition: amatrix_interface.h:41
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Communicator & GetCommunicator()
Definition: model_part.h:1821
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
NodeIterator NodesEnd(IndexType ThisIndex=0)
Definition: model_part.h:497
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
boost::indirect_iterator< typename TContainerType::iterator > iterator
Definition: pointer_vector_set.h:95
#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
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
int dimension
Definition: isotropic_damage_automatic_differentiation.py:123
int j
Definition: quadrature.py:648
N
Definition: sensitivityMatrix.py:29
integer i
Definition: TensorModule.f:17
zero
Definition: test_pureconvectionsolver_benchmarking.py:94