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.
potential_flow_utilities.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: Inigo Lopez
11 
12 #if !defined(KRATOS_POTENTIAL_FLOW_UTILITIES_H_INCLUDED )
13 #define KRATOS_POTENTIAL_FLOW_UTILITIES_H_INCLUDED
14 
15 
16 // Project includes
17 #include "containers/array_1d.h"
19 #include "utilities/geometry_utilities.h"
20 
21 namespace Kratos
22 {
23  // forward-declaring
24  class Element;
25  class ProcessInfo;
26  class ModelPart;
27  template< class TDataType >
28  class GlobalPointersVector;
29  template< class TDataType >
30  class Dof;
31  class Node;
32  template< class TPointType >
33  class Geometry;
34 
35 namespace PotentialFlowUtilities
36 {
37 template <unsigned int TNumNodes, unsigned int TDim>
39 {
40  template<typename TGeometryType>
41  ElementalData(const TGeometryType& rGeometry)
42  {
44  }
45 
47  double vol;
48 
51 };
52 
53 typedef Node NodeType;
55 
56 template <int Dim, int NumNodes>
58 
59 template <int Dim, int NumNodes>
61 
62 template <int Dim, int NumNodes>
64  const Element& rElement, const array_1d<double, NumNodes>& rDistances);
65 
66 template <int Dim, int NumNodes>
68  const Element& rElement, const array_1d<double, NumNodes>& rDistances);
69 
70 template <int Dim, int NumNodes>
72  const Element& rElement, const array_1d<double, NumNodes>& rDistances);
73 
74 template <int Dim, int NumNodes>
76 
77 template <int Dim, int NumNodes>
79 
80 template <int Dim, int NumNodes>
82 
83 template <int Dim, int NumNodes>
85 
86 template <int Dim, int NumNodes>
87 array_1d<double, Dim> ComputePerturbedVelocity(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
88 
89 template <int Dim, int NumNodes>
90 array_1d<double, Dim> ComputePerturbedVelocityLowerElement(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
91 
92 template <int Dim, int NumNodes>
93 double ComputeMaximumVelocitySquared(const ProcessInfo& rCurrentProcessInfo);
94 
95 double ComputeVacuumVelocitySquared(const ProcessInfo& rCurrentProcessInfo);
96 
97 template <int Dim, int NumNodes>
98 double ComputeClampedVelocitySquared(const array_1d<double, Dim>& rVelocity, const ProcessInfo& rCurrentProcessInfo);
99 
100 template <int Dim, int NumNodes>
101 double ComputeVelocityMagnitude(const double localMachNumberSquared, const ProcessInfo& rCurrentProcessInfo);
102 
103 template <int Dim, int NumNodes>
104 double ComputeIncompressiblePressureCoefficient(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
105 
106 template <int Dim, int NumNodes>
107 double ComputePerturbationIncompressiblePressureCoefficient(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
108 
109 template <int Dim, int NumNodes>
110 double ComputeCompressiblePressureCoefficient(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
111 
112 template <int Dim, int NumNodes>
113 double ComputePerturbationCompressiblePressureCoefficient(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
114 
115 template <int Dim, int NumNodes>
116 double ComputeLocalSpeedOfSound(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
117 
118 template <int Dim, int NumNodes>
119 double ComputeLocalSpeedofSoundSquared(const array_1d<double, Dim>& rVelocity,const ProcessInfo& rCurrentProcessInfo);
120 
121 template <int Dim, int NumNodes>
122 double ComputeSquaredSpeedofSoundFactor(const double localVelocitySquared, const ProcessInfo& rCurrentProcessInfo);
123 
124 template <int Dim, int NumNodes>
125 double ComputePerturbationLocalSpeedOfSound(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
126 
127 template <int Dim, int NumNodes>
128 double ComputeLocalMachNumber(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
129 
130 template <int Dim, int NumNodes>
131 double ComputeLocalMachNumberSquared(const array_1d<double, Dim>& rVelocity, const ProcessInfo& rCurrentProcessInfo);
132 
133 template <int Dim, int NumNodes>
134 double ComputeDerivativeLocalMachSquaredWRTVelocitySquared(const array_1d<double, Dim>& rVelocity, const double localMachNumberSquared,const ProcessInfo& rCurrentProcessInfo);
135 
136 template <int Dim, int NumNodes>
137 double ComputePerturbationLocalMachNumber(const Element& rElement, const ProcessInfo& rCurrentProcessInfo);
138 
139 template <int Dim, int NumNodes>
140 double ComputeUpwindFactor(double localMachNumberSquared,const ProcessInfo& rCurrentProcessInfo);
141 
142 template <int Dim, int NumNodes>
143 double SelectMaxUpwindFactor(const array_1d<double, Dim>& rCurrentVelocity, const array_1d<double, Dim>& rUpwindVelocity, const ProcessInfo& rCurrentProcessInfo);
144 
145 template <int Dim, int NumNodes>
146 size_t ComputeUpwindFactorCase(array_1d<double, 3>& rUpwindFactorOptions);
147 
148 template <int Dim, int NumNodes>
149 double ComputeUpwindFactorDerivativeWRTMachSquared(const double localMachNumberSquared,const ProcessInfo& rCurrentProcessInfo);
150 
151 template <int Dim, int NumNodes>
152 double ComputeUpwindFactorDerivativeWRTVelocitySquared(const array_1d<double, Dim>& rVelocity,const ProcessInfo& rCurrentProcessInfo);
153 
154 template <int Dim, int NumNodes>
155 double ComputeDensity(const double localMachNumberSquared, const ProcessInfo& rCurrentProcessInfo);
156 
157 template <int Dim, int NumNodes>
158 double ComputeUpwindedDensity(const array_1d<double, Dim>& rCurrentVelocity, const array_1d<double, Dim>& rUpwindVelocity, const ProcessInfo& rCurrentProcessInfo);
159 
160 template <int Dim, int NumNodes>
161 double ComputeDensityDerivativeWRTVelocitySquared(const double localMachNumberSquared, const ProcessInfo& rCurrentProcessInfo);
162 
163 template <int Dim, int NumNodes>
164 double ComputeUpwindedDensityDerivativeWRTVelocitySquaredSupersonicAccelerating(const array_1d<double, Dim>& rCurrentVelocity, const double currentMachNumberSquared, const double upwindMachNumberSquared, const ProcessInfo& rCurrentProcessInfo);
165 
166 template <int Dim, int NumNodes>
167 double ComputeUpwindedDensityDerivativeWRTVelocitySquaredSupersonicDeaccelerating(const double currentMachNumberSquared, const double upwindMachNumberSquared, const ProcessInfo& rCurrentProcessInfo);
168 
169 template <int Dim, int NumNodes>
170 double ComputeUpwindedDensityDerivativeWRTUpwindVelocitySquaredSupersonicAccelerating(const double currentMachNumberSquared, const double upwindMachNumberSquared, const ProcessInfo& rCurrentProcessInfo);
171 
172 template <int Dim, int NumNodes>
173 double ComputeUpwindedDensityDerivativeWRTUpwindVelocitySquaredSupersonicDeaccelerating(const array_1d<double, Dim>& rUpwindVelocity, const double currentMachNumberSquared, const double upwindMachNumberSquared, const ProcessInfo& rCurrentProcessInfo);
174 
175 template <int Dim, int NumNodes>
177 
178 bool CheckIfElementIsTrailingEdge(const Element& rElement);
179 
180 template <int Dim>
181 void CheckIfWakeConditionsAreFulfilled(const ModelPart& rWakeModelPart, const double& rTolerance, const int& rEchoLevel);
182 
183 template <int Dim, int NumNodes>
184 bool CheckWakeCondition(const Element& rElement, const double& rTolerance, const int& rEchoLevel);
185 
186 template <int Dim, int NumNodes>
187 void GetSortedIds(std::vector<size_t>& Ids, const GeometryType& rGeom);
188 
189 template <int Dim, int NumNodes>
191 
192 template<int Dim>
193 Vector ComputeKuttaNormal(const double angle);
194 
195 template <class TContainerType>
196 double CalculateArea(TContainerType& rContainer);
197 
198 template <int Dim, int NumNodes>
199 void ComputePotentialJump(ModelPart& rWakeModelPart);
200 
201 template <int Dim, int NumNodes>
202 void AddKuttaConditionPenaltyTerm(const Element& rElement,
203  Matrix& rLeftHandSideMatrix,
204  Vector& rRightHandSideVector,
205  const ProcessInfo& rCurrentProcessInfo);
206 
207 template <int Dim, int NumNodes>
209  Vector& rRightHandSideVector,
210  const ProcessInfo& rCurrentProcessInfo);
211 
212 template <int Dim, int NumNodes>
214  Matrix& rLeftHandSideMatrix,
215  const ProcessInfo& rCurrentProcessInfo);
216 
217 template <int Dim, int NumNodes>
219  Matrix& rLeftHandSideMatrix,
220  Vector& rRightHandSideVector,
221  const ProcessInfo& rCurrentProcessInfo);
222 
223 } // namespace PotentialFlow
224 } // namespace Kratos
225 
226 #endif // KRATOS_POTENTIAL_FLOW_UTILITIES_H_INCLUDED defined
Base class for all Elements.
Definition: element.h:60
Geometry base class.
Definition: geometry.h:71
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
This class defines the node.
Definition: node.h:65
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Kratos::ModelPart ModelPart
Definition: kratos_wrapper.h:31
void AddPotentialGradientStabilizationTerm(Element &rElement, Matrix &rLeftHandSideMatrix, Vector &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:1158
bool CheckWakeCondition(const Element &rElement, const double &rTolerance, const int &rEchoLevel)
Definition: potential_flow_utilities.cpp:897
double CalculateArea(TContainerType &rContainer)
Definition: potential_flow_utilities.cpp:970
array_1d< double, Dim > ComputePerturbedVelocityLowerElement(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:138
double ComputeUpwindFactorDerivativeWRTVelocitySquared(const array_1d< double, Dim > &rVelocity, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:671
double ComputePerturbationLocalMachNumber(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:568
double ComputeIncompressiblePressureCoefficient(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:304
double ComputeUpwindFactorDerivativeWRTMachSquared(const double localMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:651
bool CheckIfElementIsTrailingEdge(const Element &rElement)
Definition: potential_flow_utilities.cpp:864
double ComputeUpwindedDensityDerivativeWRTUpwindVelocitySquaredSupersonicAccelerating(const double currentMachNumberSquared, const double upwindMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:805
double ComputePerturbationLocalSpeedOfSound(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:468
Geometry< NodeType > GeometryType
Definition: potential_flow_utilities.h:54
void AddKuttaConditionPenaltyTerm(const Element &rElement, Matrix &rLeftHandSideMatrix, Vector &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:1011
void AddKuttaConditionPenaltyPerturbationLHS(const Element &rElement, Matrix &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:1060
double ComputeLocalSpeedofSoundSquared(const array_1d< double, Dim > &rVelocity, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:424
array_1d< double, Dim > ComputeVelocityUpperWakeElement(const Element &rElement)
Definition: potential_flow_utilities.cpp:280
size_t ComputeUpwindFactorCase(array_1d< double, 3 > &rUpwindFactorOptions)
Definition: potential_flow_utilities.cpp:629
double ComputeLocalMachNumberSquared(const array_1d< double, Dim > &rVelocity, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:510
double ComputeDerivativeLocalMachSquaredWRTVelocitySquared(const array_1d< double, Dim > &rVelocity, const double localMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:529
double ComputeDensity(const double localMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:686
double ComputePerturbationIncompressiblePressureCoefficient(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:321
double ComputeUpwindedDensity(const array_1d< double, Dim > &rCurrentVelocity, const array_1d< double, Dim > &rUpwindVelocity, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:712
bool CheckIfElementIsCutByDistance(const BoundedVector< double, NumNodes > &rNodalDistances)
Definition: potential_flow_utilities.cpp:843
double ComputePerturbationCompressiblePressureCoefficient(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:368
void CheckIfWakeConditionsAreFulfilled(const ModelPart &rWakeModelPart, const double &rTolerance, const int &rEchoLevel)
Definition: potential_flow_utilities.cpp:880
double ComputeLocalSpeedOfSound(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:398
void GetNodeNeighborElementCandidates(GlobalPointersVector< Element > &ElementCandidates, const GeometryType &rGeom)
Definition: potential_flow_utilities.cpp:933
double ComputeClampedVelocitySquared(const array_1d< double, Dim > &rVelocity, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:215
double ComputeCompressiblePressureCoefficient(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:342
double ComputeSquaredSpeedofSoundFactor(const double localVelocitySquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:447
Node NodeType
Definition: potential_flow_utilities.h:53
BoundedVector< double, NumNodes > GetPotentialOnUpperWakeElement(const Element &rElement, const array_1d< double, NumNodes > &rDistances)
Definition: potential_flow_utilities.cpp:74
double ComputeUpwindedDensityDerivativeWRTUpwindVelocitySquaredSupersonicDeaccelerating(const array_1d< double, Dim > &rUpwindVelocity, const double currentMachNumberSquared, const double upwindMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:821
array_1d< double, Dim > ComputeVelocity(const Element &rElement)
Definition: potential_flow_utilities.cpp:112
double ComputeUpwindFactor(double localMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:585
double ComputeMaximumVelocitySquared(const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:153
BoundedVector< double, NumNodes > GetPotentialOnNormalElement(const Element &rElement)
Definition: potential_flow_utilities.cpp:28
void ComputePotentialJump(ModelPart &rWakeModelPart)
Definition: potential_flow_utilities.cpp:980
array_1d< double, Dim > ComputeVelocityLowerWakeElement(const Element &rElement)
Definition: potential_flow_utilities.cpp:292
BoundedVector< double, 2 *NumNodes > GetPotentialOnWakeElement(const Element &rElement, const array_1d< double, NumNodes > &rDistances)
Definition: potential_flow_utilities.cpp:55
double ComputeUpwindedDensityDerivativeWRTVelocitySquaredSupersonicDeaccelerating(const double currentMachNumberSquared, const double upwindMachSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:787
double ComputeDensityDerivativeWRTVelocitySquared(const double localMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:733
array_1d< double, NumNodes > GetWakeDistances(const Element &rElement)
Definition: potential_flow_utilities.cpp:22
array_1d< double, Dim > ComputePerturbedVelocity(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:123
double ComputeVacuumVelocitySquared(const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:186
double ComputeUpwindedDensityDerivativeWRTVelocitySquaredSupersonicAccelerating(const array_1d< double, Dim > &rCurrentVelocity, const double currentMachNumberSquared, const double upwindMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:765
BoundedVector< double, NumNodes > GetPotentialOnLowerWakeElement(const Element &rElement, const array_1d< double, NumNodes > &rDistances)
Definition: potential_flow_utilities.cpp:93
void GetSortedIds(std::vector< size_t > &Ids, const GeometryType &rGeom)
Definition: potential_flow_utilities.cpp:921
double ComputeLocalMachNumber(const Element &rElement, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:497
void AddKuttaConditionPenaltyPerturbationRHS(const Element &rElement, Vector &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:1103
array_1d< double, Dim > ComputeVelocityNormalElement(const Element &rElement)
Definition: potential_flow_utilities.cpp:270
double SelectMaxUpwindFactor(const array_1d< double, Dim > &rCurrentVelocity, const array_1d< double, Dim > &rUpwindVelocity, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:608
Vector ComputeKuttaNormal(const double angle)
double ComputeVelocityMagnitude(const double localMachNumberSquared, const ProcessInfo &rCurrentProcessInfo)
Definition: potential_flow_utilities.cpp:238
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ProcessInfo
Definition: edgebased_PureConvection.py:116
Definition: potential_flow_utilities.h:39
double vol
Definition: potential_flow_utilities.h:47
array_1d< double, TNumNodes > distances
Definition: potential_flow_utilities.h:46
array_1d< double, TNumNodes > potentials
Definition: potential_flow_utilities.h:46
BoundedMatrix< double, TNumNodes, TDim > DN_DX
Definition: potential_flow_utilities.h:49
ElementalData(const TGeometryType &rGeometry)
Definition: potential_flow_utilities.h:41
array_1d< double, TNumNodes > N
Definition: potential_flow_utilities.h:50