VM2D 1.14
Vortex methods for 2D flows simulation
Loading...
Searching...
No Matches
VM2D::Airfoil Class Reference

Абстрактный класс, определяющий обтекаемый профиль More...

#include <Airfoil2D.h>

Inheritance diagram for VM2D::Airfoil:
Collaboration diagram for VM2D::Airfoil:

Public Member Functions

void lightningTest ()
 Тест на "отвещенность".
 
 IFCUDA (mutable double *devRPtr)
 Указатель на девайсе, где хранятся вершины профиля
 
 IFCUDA (mutable double *devPsnPtr)
 Указатель на девайсе, где хранятся псевдонормали профиля
 
 IFCUDA (mutable double *devRhsPtr)
 Указатель на девайсе, где хранится правая часть (константная) матрицы
 
 IFCUDA (mutable double *devRhsLinPtr)
 Указатель на девайсе, где хранится правая часть (линейная) матрицы
 
 IFCUDA (mutable std::vector< double > tmpRhs)
 Указатель на хосте, где хранится временная часть матрицы, полученная с девайса
 
 IFCUDA (mutable double *devFreeVortexSheetPtr)
 Указатель на девайсе, где хранятся интенсивности (константные) свободного вихревого слоя на панелях
 
 IFCUDA (mutable double *devFreeVortexSheetLinPtr)
 Указатель на девайсе, где хранятся интенсивности (линейные) свободного вихревого слоя на панелях
 
 IFCUDA (mutable double *devAttachedVortexSheetPtr)
 Указатель на девайсе, где хранятся интенсивности (константные) присоединенного вихревого слоя на панелях
 
 IFCUDA (mutable double *devAttachedVortexSheetLinPtr)
 Указатель на девайсе, где хранятся интенсивности (линейные) присоединенного вихревого слоя на панелях
 
 IFCUDA (mutable double *devAttachedSourceSheetPtr)
 Указатель на девайсе, где хранятся интенсивности (константные) присоединенного слоя источников на панелях
 
 IFCUDA (mutable double *devAttachedSourceSheetLinPtr)
 Указатель на девайсе, где хранятся интенсивности (линейные) присоединенного слоя источников на панелях
 
 IFCUDA (mutable double *devMeanEpsOverPanelPtr)
 Указатель на девайсе, где хранятся средние eps на панелях
 
 IFCUDA (mutable double *devViscousStressesPtr)
 Указатель на девайсе, где хранится вектор (по панелям) для силы вязкого трения
 
 IFCUDA (mutable std::vector< double > tmpViscousStresses)
 Указатель на хосте, где хранится временная часть вектора (по панелям) для силы вязкого трения
 
 Airfoil (const World2D &W_, const size_t numberInPassport_)
 
virtual ~Airfoil ()
 Деструктор
 
bool isAfter (size_t i, size_t j) const
 Проверка, идет ли вершина i следом за вершиной j.
 
void CalcNrmTauLen ()
 Вычисление нормалей, касательных и длин панелей по текущему положению вершин
 
virtual void Rotate (double alpha)
 Поворот профиля
 
virtual void Scale (const Point2D &)
 Масштабирование профиля
 
virtual void Move (const Point2D &dr)
 Перемещение профиля
 
virtual void GetGabarits (double gap=0.02)
 Вычисляет габаритный прямоугольник профиля
 
bool isInsideGabarits (const Point2D &r) const
 Определяет, находится ли точка с радиус-вектором \( \vec r \) внутри габаритного прямоугольника профиля
 
bool isOutsideGabarits (const Point2D &r) const
 Определяет, находится ли точка с радиус-вектором \( \vec r \) вне габаритного прямоугольника профиля
 
virtual bool IsPointInAirfoil (const Point2D &point) const
 Определяет, находится ли точка с радиус-вектором \( \vec r \) внутри профиля
 
virtual void ReadFromFile (const std::string &dir)
 Считывание профиля из файла
 
virtual std::vector< double > getA (size_t p, size_t i, const Airfoil &airfoil, size_t j) const
 Вычисление коэффициентов матрицы A для расчета влияния панели на панель
 
virtual void calcIQ (size_t p, const Airfoil &otherAirfoil, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &matrPair) const
 Вычисление коэффициентов матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2.
 
void calcMeanEpsOverPanel ()
 Вычисление средних значений eps на панелях
 
virtual void GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses (const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)
 Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек, обусловленных геометрией профиля, и вычисление вязкого трения
 
virtual void GetDiffVelocityI0I3ToWakeAndViscousStresses (const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)
 
virtual void GetInfluenceFromVorticesToPanel (size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const
 Вычисление влияния части подряд идущих вихрей из вихревого следа на панель для правой части
 
virtual void GetInfluenceFromSourcesToPanel (size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const
 Вычисление влияния части подряд идущих источников из области течения на панель для правой части
 
virtual void GetInfluenceFromSourceSheetToVortex (size_t panel, const Vortex2D &vtx, Point2D &vel) const
 Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения

 
virtual void GetInfluenceFromVortexSheetToVortex (size_t panel, const Vortex2D &vtx, Point2D &vel) const
 Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вихрь в области течения
 
virtual void GetInfluenceFromVInfToPanel (std::vector< double > &vInfRhs) const
 Вычисление влияния набегающего потока на панель для правой части
 
const Point2DgetR (size_t q) const
 Возврат константной ссылки на вершину профиля

 
Point2DsetR (size_t q)
 Возврат ссылки на вершину профиля

 
const Point2DgetV (size_t q) const
 Возврат константной ссылки на скорость вершины профиля

 
void setV (const Point2D &vel)
 Установка постоянной скорости всех вершин профиля
 
void setV (const std::vector< Point2D > &vel)
 Установка скоростей всех вершин профиля
 
size_t getNumberOfPanels () const
 Возврат количества панелей на профиле
 

Public Attributes

const World2DW
 Константная ссылка на решаемую задачу
 
const size_t numberInPassport
 Номер профиля в паспорте
 
double m
 Масса профиля
 
double J
 Полярный момент инерции профиля относительно центра масс
 
std::vector< std::vector< Point2D > > possibleWays
 Возможные пути внутри профиля от точки (0, 0) к центрам всех панелей
 
std::vector< int > wayToVertex
 Номера путей к вершинам
 
std::vector< double > meanEpsOverPanel
 Средние значения Eps на панелях
 
std::vector< double > viscousStress
 Нейросеть для коэффициентов I0 и I3 диффузионной скорости
 
Point2D lowLeft
 Левый нижний угол габаритного прямоугольника профиля
 
Point2D upRight
 Правый верхний угол габаритного прямоугольника профиля
 
std::vector< double > gammaThrough
 Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
 
bool inverse
 Признак разворота нормалей (для расчета внутренних течений)
 
std::vector< Point2Dnrm
 Нормали к панелям профиля
 
std::vector< std::pair< Point2D, Point2D > > psn
 Псевдонормали к панелям профиля
 
std::vector< Point2Dtau
 Касательные к панелям профиля
 
std::vector< double > len
 Длины панелей профиля

 
Point2D rcm
 Положение центра масс профиля
 
double phiAfl
 Поворот профиля
 
double area
 Площадь профиля
 

Protected Attributes

std::vector< Point2Dr_
 Координаты начал панелей

 
std::vector< Point2Dv_
 Скорости начал панелей

 

Detailed Description

Абстрактный класс, определяющий обтекаемый профиль

Author
Марчевский Илья Константинович
Сокол Ксения Сергеевна
Рятина Евгения Павловна
Колганова Александра Олеговна \Version 1.14
Date
6 марта 2026 г.

Definition at line 181 of file Airfoil2D.h.

Constructor & Destructor Documentation

◆ Airfoil()

Airfoil::Airfoil ( const World2D W_,
const size_t  numberInPassport_ 
)

Конструктор

Parameters
[in]W_константная ссылка на решаемую задачу
[in]numberInPassport_номер профиля в паспорте задачи

Definition at line 58 of file Airfoil2D.cpp.

59 : W(W_)
60 , numberInPassport(numberInPassport_)
61#ifdef USE_CUDA
62#ifdef USE_CUDNN
63 , NN()
64 , vtp(0)
65#endif
66#endif
67{
68
69}
const World2D & W
Константная ссылка на решаемую задачу
Definition Airfoil2D.h:185
const size_t numberInPassport
Номер профиля в паспорте
Definition Airfoil2D.h:188

◆ ~Airfoil()

virtual VM2D::Airfoil::~Airfoil ( )
inlinevirtual

Деструктор

Definition at line 285 of file Airfoil2D.h.

285{ };

Member Function Documentation

◆ calcIQ()

void Airfoil::calcIQ ( size_t  p,
const Airfoil otherAirfoil,
std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &  matrPair 
) const
virtual

Вычисление коэффициентов матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2.

Parameters
[in]pразмерность матрицы-результата
[in]otherAirfoilконстантная ссылка на профиль, от которого рассчитывается влияние
[out]matrPairссылка на пару матриц, выражающих взаимные влияния (касательные и нормальные) панелей профиля return соответствующий блок матрицы СЛАУ, вытянутый в линию

Definition at line 1180 of file Airfoil2D.cpp.

1181{
1182 bool self = (&otherAirfoil == this);
1183
1184 size_t aflNSelf = numberInPassport;
1185 size_t aflNOther = otherAirfoil.numberInPassport;
1186
1187 size_t npI;
1188 size_t npJ;
1189
1190 numvector<double, 3> alpha, lambda;
1191
1192 //auxillary vectors
1193 Point2D p1, s1, p2, s2, di, dj, i00, i01, i10, i11;
1194 numvector<Point2D, 3> v00, v11;
1195 numvector<Point2D, 2> v01, v10;
1196
1197#pragma omp parallel for \
1198 default(none) \
1199 shared(otherAirfoil, self, aflNSelf, aflNOther, matrPair, p, IDPI, IQPI) \
1200 private(npI, npJ, alpha, lambda, p1, s1, p2, s2, di, dj, i00, i01, i10, i11, v00, v11, v01, v10) schedule(dynamic, DYN_SCHEDULE)
1201 for (int i = 0; i < (int)getNumberOfPanels(); ++i)
1202 for (int j = 0; j < (int)otherAirfoil.getNumberOfPanels(); ++j)
1203 {
1204 npI = getNumberOfPanels();
1205 npJ = otherAirfoil.getNumberOfPanels();
1206
1207 if ((i == j) && self)
1208 {
1209
1210 matrPair.first(i, j) = 0.0;
1211 matrPair.second(i, j) = 0.0;
1212
1213 if (p == 2)
1214 {
1215 matrPair.first(i, npJ + j) = 0.0;
1216 matrPair.first(npI + i, j) = 0.0;
1217
1218 matrPair.second(i, npJ + j) = -IQPI;
1219 matrPair.second(npI + i, j) = IQPI;
1220
1221 matrPair.first(npI + i, npJ + j) = 0.0;
1222 matrPair.second(npI + i, npJ + j) = 0.0;
1223
1224 }
1225 if (p > 2)
1226 throw (-42);
1227 }//if i==j
1228 else
1229 {
1230 const Point2D& taui = tau[i];
1231 const Point2D& tauj = otherAirfoil.tau[j];
1232
1233 p1 = getR(i + 1) - otherAirfoil.getR(j + 1);
1234 s1 = getR(i + 1) - otherAirfoil.getR(j);
1235 p2 = getR(i) - otherAirfoil.getR(j + 1);
1236 s2 = getR(i) - otherAirfoil.getR(j);
1237 di = getR(i + 1) - getR(i);
1238 dj = otherAirfoil.getR(j + 1) - otherAirfoil.getR(j);
1239
1240 alpha = { \
1241 (self && isAfter(j, i)) ? 0.0 : VMlib::Alpha(s2, s1), \
1242 VMlib::Alpha(s2, p1), \
1243 (self && isAfter(i, j)) ? 0.0 : VMlib::Alpha(p1, p2) \
1244 };
1245
1246 lambda = { \
1247 (self && isAfter(j, i)) ? 0.0 : VMlib::Lambda(s2, s1), \
1248 VMlib::Lambda(s2, p1), \
1249 (self && isAfter(i, j)) ? 0.0 : VMlib::Lambda(p1, p2) \
1250 };
1251
1252 v00 = {
1253 VMlib::Omega(s1, taui, tauj),
1254 -VMlib::Omega(di, taui, tauj),
1255 VMlib::Omega(p2, taui, tauj)
1256 };
1257
1258 i00 = IDPI / len[i] * (-(alpha[0] * v00[0] + alpha[1] * v00[1] + alpha[2] * v00[2]).kcross() \
1259 + (lambda[0] * v00[0] + lambda[1] * v00[1] + lambda[2] * v00[2]));
1260
1261 matrPair.first(i, j) = i00 & nrm[i];
1262 matrPair.second(i, j) = i00 & taui;
1263
1264
1265 if (p == 2)
1266 {
1267 v01 = {
1268 0.5 / (dj.length()) * (((p1 + s1) & tauj) * VMlib::Omega(s1, taui, tauj) - s1.length2() * taui),
1269 -0.5 * di.length() / dj.length() * VMlib::Omega(s1 + p2, tauj, tauj)
1270 };
1271
1272 i01 = IDPI / len[i] * (-((alpha[0] + alpha[2]) * v01[0] + (alpha[1] + alpha[2]) * v01[1]).kcross() \
1273 + ((lambda[0] + lambda[2]) * v01[0] + (lambda[1] + lambda[2]) * v01[1]) - 0.5 * di.length() * tauj);
1274
1275 matrPair.first(i, npJ + j) = i01 & nrm[i];
1276 matrPair.second(i, npJ + j) = i01 & taui;
1277
1278
1279 v10 = {
1280 -0.5 / di.length() * (((s1 + s2) & taui) * VMlib::Omega(s1, taui, tauj) - s1.length2() * tauj),
1281 0.5 * dj.length() / di.length() * VMlib::Omega(s1 + p2, taui, taui)
1282 };
1283
1284 i10 = IDPI / len[i] * (-((alpha[0] + alpha[2]) * v10[0] + alpha[2] * v10[1]).kcross() \
1285 + ((lambda[0] + lambda[2]) * v10[0] + lambda[2] * v10[1]) + 0.5 * dj.length() * taui);
1286
1287 matrPair.first(npI + i, j) = i10 & nrm[i];
1288 matrPair.second(npI + i, j) = i10 & taui;
1289
1290
1291 v11 = {
1292 1.0 / (12.0 * di.length() * dj.length()) * (2.0 * (s1 & VMlib::Omega(s1 - 3.0 * p2, taui, tauj)) * VMlib::Omega(s1, taui, tauj) - s1.length2() * (s1 - 3.0 * p2)) - 0.25 * VMlib::Omega(s1, taui, tauj),
1293 -di.length() / (12.0 * dj.length()) * VMlib::Omega(di, tauj, tauj),
1294 -dj.length() / (12.0 * di.length()) * VMlib::Omega(dj, taui, taui)
1295 };
1296
1297 i11 = IDPI / len[i] * (-((alpha[0] + alpha[2]) * v11[0] + (alpha[1] + alpha[2]) * v11[1] + alpha[2] * v11[2]).kcross()\
1298 + (lambda[0] + lambda[2]) * v11[0] + (lambda[1] + lambda[2]) * v11[1] + lambda[2] * v11[2] \
1299 + 1.0 / 12.0 * (dj.length() * taui + di.length() * tauj - 2.0 * VMlib::Omega(s1, taui, tauj)));
1300
1301 matrPair.first(npI + i, npJ + j) = i11 & nrm[i];
1302 matrPair.second(npI + i, npJ + j) = i11 & taui;
1303 }
1304
1305
1306 if (p > 2)
1307 throw (-42);
1308 }//else(i == j)
1309
1310 }//for(...)
1311}//getIQ(...)
std::vector< double > len
Длины панелей профиля
Definition Airfoil2D.h:94
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition Airfoil2D.h:113
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition Airfoil2D.h:81
std::vector< Point2D > tau
Касательные к панелям профиля
Definition Airfoil2D.h:91
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition Airfoil2D.h:163
bool isAfter(size_t i, size_t j) const
Проверка, идет ли вершина i следом за вершиной j.
Definition Airfoil2D.cpp:73
Шаблонный класс, определяющий вектор фиксированной длины Фактически представляет собой массив,...
Definition numvector.h:99
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition numvector.h:389
P length() const
Вычисление 2-нормы (длины) вектора
Definition numvector.h:377
const double IQPI
Число .
Definition defs.h:91
const double IDPI
Число .
Definition defs.h:85
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition defs.cpp:268
double Alpha(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления угла между векторами
Definition defs.cpp:262
Point2D Omega(const Point2D &a, const Point2D &b, const Point2D &c)
Вспомогательная функция вычисления величины .
Definition defs.cpp:274
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calcMeanEpsOverPanel()

void Airfoil::calcMeanEpsOverPanel ( )

Вычисление средних значений eps на панелях

Definition at line 93 of file Airfoil2D.cpp.

94{
95 meanEpsOverPanel.clear();
97
98
99 double midEps;
100
103
104 for (size_t i = 0; i < getNumberOfPanels(); ++i)
105 {
106 midEps = 0.0;
107 for (int j = bnd.vortexBeginEnd[i].first; j < bnd.vortexBeginEnd[i].second; ++j)
108 midEps += virtVortParams.epsastWake[j];
109 //midEps += std::max(virtVortParams.epsastWake[j], 0.5 * len[i] / (bnd.vortexBeginEnd[i].second - bnd.vortexBeginEnd[i].first));
110
111 midEps /= (bnd.vortexBeginEnd[i].second - bnd.vortexBeginEnd[i].first);
112 meanEpsOverPanel[i] = midEps;
113 }//for i
114}//calcMeanEpsOverPanel()
std::vector< double > meanEpsOverPanel
Средние значения Eps на панелях
Definition Airfoil2D.h:209
Абстрактный класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле
Definition Boundary2D.h:65
std::vector< std::pair< int, int > > vortexBeginEnd
Номера первого и последнего вихрей, рождаемых на каждой панели профиля (формируется после решения СЛА...
Definition Boundary2D.h:83
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition Velocity2D.h:115
const Velocity & getVelocity() const
Возврат константной ссылки на объект для вычисления скоростей
Definition World2D.h:241
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition World2D.h:180
Структура, определяющая параметры виртуальных вихрей для отдельного профиля
Definition Velocity2D.h:70
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition Velocity2D.h:90
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcNrmTauLen()

void Airfoil::CalcNrmTauLen ( )

Вычисление нормалей, касательных и длин панелей по текущему положению вершин

Definition at line 1047 of file Airfoil2D.cpp.

1048{
1049 if (nrm.size() != r_.size())
1050 {
1051 nrm.resize(r_.size());
1052 psn.resize(r_.size());
1053 tau.resize(r_.size());
1054 len.resize(r_.size());
1055 }
1056
1057 Point2D rpan;
1058
1059 //считаем длины, касательные и нормали
1060#pragma omp parallel for private(rpan)
1061 for (int i = 0; i < (int)r_.size(); ++i)
1062 {
1063 rpan = (getR(i + 1) - getR(i));
1064 len[i] = rpan.length();
1065 tau[i] = rpan.unit();
1066 nrm[i] = { tau[i][1], -tau[i][0] };
1067 }
1068
1069 //считаем псевдонормали
1070#pragma omp parallel for
1071 for (int i = 0; i < (int)r_.size(); ++i)
1072 {
1073 int next = (i == (int)r_.size() - 1 ? 0 : i + 1);
1074 int prev = (i == 0 ? (int)r_.size() - 1 : i - 1);
1075
1076 psn[i] = std::make_pair( (nrm[prev] + nrm[i]).unit(), (nrm[i] + nrm[next]).unit() );
1077 }
1078
1079 //закольцовываем
1080 //nrm.push_back(nrm[0]);
1081 //tau.push_back(tau[0]);
1082 //len.push_back(len[0]);
1083}//CalcNrmTauLen()
std::vector< Point2D > r_
Координаты начал панелей
Definition Airfoil2D.h:69
std::vector< std::pair< Point2D, Point2D > > psn
Псевдонормали к панелям профиля
Definition Airfoil2D.h:86
auto unit(P newlen=1) const -> numvector< typename std::remove_const< decltype(this->data[0] *newlen)>::type, n >
Вычисление орта вектора или вектора заданной длины, коллинеарного данному
Definition numvector.h:405
Here is the call graph for this function:
Here is the caller graph for this function:

◆ getA()

std::vector< double > Airfoil::getA ( size_t  p,
size_t  i,
const Airfoil airfoil,
size_t  j 
) const
virtual

Вычисление коэффициентов матрицы A для расчета влияния панели на панель

Parameters
[in]pразмерность матрицы - результата
[in]iномер панели, на которую оказывается влияние
[in]airfoilконстантная ссылка на профиль, от которого рассчитывается влияние
[in]jномер влияющей панели return соответствующий блок матрицы СЛАУ, вытянутый в линию

Definition at line 1139 of file Airfoil2D.cpp.

1140{
1141 std::vector<double> res(p * p, 0.0);
1142
1143 if ((i == j) && (&airfoil == this))
1144 {
1145 res[0] = 0.5 * (tau[i] ^ nrm[i]);
1146 if (p == 1)
1147 return res;
1148
1149 res[3] = (1.0 / 12.0) * res[0];
1150 if (p == 2)
1151 return res;
1152
1153 if (p > 2)
1154 throw (-42);
1155
1156 }//if i==j
1157
1158 const auto& miq = W.getIQ(numberInPassport, airfoil.numberInPassport);
1159
1160 res[0] = miq.first(i, j);
1161
1162 if (p == 1)
1163 return res;
1164
1165 res[1] = miq.first(i, airfoil.getNumberOfPanels() + j);
1166 res[2] = miq.first(getNumberOfPanels() + i, j);
1167 res[3] = miq.first(getNumberOfPanels() + i, airfoil.getNumberOfPanels() + j);
1168
1169 if (p == 2)
1170 return res;
1171
1172 if (p > 2)
1173 throw (-42);
1174
1175 //Хотя сюда никогда и не попадем
1176 return res;
1177}//getA(...)
const std::pair< Eigen::MatrixXd, Eigen::MatrixXd > & getIQ(size_t i, size_t j) const
Возврат константной ссылки на объект, связанный с матрицей интегралов от (r-xi)/|r-xi|^2.
Definition World2D.h:271
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses()

void Airfoil::GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses ( const WakeDataBase pointsDb,
std::vector< double > &  domainRadius,
std::vector< double > &  I0,
std::vector< Point2D > &  I3 
)
virtual

Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек, обусловленных геометрией профиля, и вычисление вязкого трения

Вычисляет диффузионные скорости в наборе точек, которые обусловленных геометрией профиля, и вычисляет вязкое трение

Parameters
[in]pointsDbконстантная ссылка на базу данных вихрей, в которых вычисляются скорости
[in]domainRadiusссылка на радиусы вихрей
[out]I0ссылка на вектор знаменателей диффузионных скоростей, которые приобретают точки из-за влияния геометрии профиля
[out]I3ссылка на вектор числителей диффузионных скоростей, которые приобретают точки из-за влияния геометрии профиля
Warning
Векторы I0, I3 — накапливаются!

Definition at line 492 of file Airfoil2D.cpp.

493{
494 std::vector<double> selfI0;
495 std::vector<Point2D> selfI3;
496
497 std::vector<double> currViscousStress;
498 currViscousStress.resize(r_.size(), 0.0);
499
500
501 selfI0.resize(pointsDb.vtx.size(), 0.0);
502 selfI3.resize(pointsDb.vtx.size(), { 0.0, 0.0 });
503
504#pragma warning (push)
505#pragma warning (disable: 4101)
506 //Локальные переменные для цикла
507 Point2D q, xi, xi_m, v0;
508 double lxi, lxi_m, lenj_m;
509
510 Point2D mn;
511 int new_n;
512 Point2D h;
513
514 Point2D vec;
515 double s, d;
516
517 double vs;
518 double iDDomRad, domRad, expon;
519#pragma warning (pop)
520
521#pragma omp parallel for \
522 default(none) \
523 shared(domainRadius, pointsDb, selfI0, selfI3, currViscousStress, PI) \
524 private(xi, xi_m, lxi, lxi_m, lenj_m, v0, q, new_n, mn, h, d, s, vec, vs, expon, domRad, iDDomRad) schedule(dynamic, DYN_SCHEDULE)
525 for (int i = 0; i < pointsDb.vtx.size(); ++i)
526 {
527 const Vortex2D& vtxI = pointsDb.vtx[i];
528
529 domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
530 iDDomRad = 1.0 / domRad;
531
532 for (size_t j = 0; j < r_.size(); ++j)
533 {
534 vs = 0.0;
535 q = vtxI.r() - 0.5 * (getR(j) + getR(j + 1));
536 vec = tau[j];
537
538 s = q & vec;
539 d = fabs(q & nrm[j]);
540
541 if ((d < 50.0 * len[j]) && (fabs(s) < 50.0 * len[j]))
542 {
543 v0 = vec * len[j];
544
545 if ((d > 5.0 * len[j]) || (fabs(s) > 5.0 * len[j]))
546 {
547 xi = q * iDDomRad;
548 lxi = xi.length();
549
550 expon = exp(-lxi) * len[j];
551 mn = nrm[j] * expon;
552
553 //if (selfI0[i] != -PI * domRad)
554 //{
555 selfI0[i] += (xi & mn) * (lxi + 1.0) / (lxi * lxi);
556 selfI3[i] += mn;
557 //}
558
559 vs = vtxI.g() * expon / (PI * sqr(meanEpsOverPanel[j]));
560 }
561 else if ((d >= 0.1 * len[j]) || (fabs(s) > 0.5 * len[j]))
562 {
563 vs = 0.0;
564 double den = (fabs(s) < 0.5 * len[j]) ? d : (fabs(s) + d - 0.5 * len[j]);
565
566 new_n = std::max(1, static_cast<int>(ceil(5.0 * len[j] / den)));
567 //new_n = 100;
568 h = v0 * (1.0 / new_n);
569
570 for (int m = 0; m < new_n; ++m)
571 {
572 xi_m = (vtxI.r() - (getR(j) + h * (m + 0.5))) * iDDomRad;
573 lxi_m = xi_m.length();
574
575 lenj_m = len[j] / new_n;
576 expon = exp(-lxi_m) * lenj_m;
577
578 mn = nrm[j] * expon;
579 //if (selfI0[i] != -PI * domRad)
580 {
581 selfI0[i] += (xi_m & mn) * (lxi_m + 1.0) / (lxi_m * lxi_m);
582 selfI3[i] += mn;
583 }
584 vs += expon;
585 }//for m
586 vs *= vtxI.g() / (PI * sqr(meanEpsOverPanel[j]));
587 }
588 else
589 {
590
591 selfI0[i] += -PI * domRad;
592 double mnog = 2.0 * domRad * (1.0 - exp(-len[j] * iDDomRad / 2.0) * cosh(fabs(s) * iDDomRad));
593 selfI3[i] += nrm[j] * mnog;
594 vs = mnog * vtxI.g() / (PI * sqr(meanEpsOverPanel[j]));
595 }
596 }//if d<50 len
597
598#pragma omp atomic
599 currViscousStress[j] += vs;
600 }//for j
601 }//for i
602
603
604 if (&pointsDb == &(W.getWake()))
605 for (size_t i = 0; i < viscousStress.size(); ++i)
606 {
607 viscousStress[i] += currViscousStress[i];
608 }
609
610 for (size_t i = 0; i < I0.size(); ++i)
611 {
612 domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
613
614 if (I0[i] != -PI * domRad)
615 {
616 if (selfI0[i] == -PI * domRad)
617 {
618 I0[i] = selfI0[i];
619 I3[i] = selfI3[i];
620 }
621 else
622 {
623 I0[i] += selfI0[i];
624 I3[i] += selfI3[i];
625 }
626 }
627 }
628}; //GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(...)
std::vector< double > viscousStress
Нейросеть для коэффициентов I0 и I3 диффузионной скорости
Definition Airfoil2D.h:268
double m
Масса профиля
Definition Airfoil2D.h:193
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition Passport2D.h:287
std::vector< Vortex2D > vtx
Список вихревых элементов
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition World2D.h:226
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition World2D.h:251
Класс, опеделяющий двумерный вихревой элемент
Definition Vortex2D.h:59
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition Vortex2D.h:87
HD double & g()
Функция для доступа к циркуляции вихря
Definition Vortex2D.h:95
const double PI
Число .
Definition defs.h:82
T sqr(T x)
Умножение a на комплексно сопряженноe к b.
Definition defsBH.h:101
double getMinEpsAst() const
Функция минимально возможного значения для epsAst.
Definition Passport2D.h:148
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetDiffVelocityI0I3ToWakeAndViscousStresses()

void Airfoil::GetDiffVelocityI0I3ToWakeAndViscousStresses ( const WakeDataBase pointsDb,
std::vector< double > &  domainRadius,
std::vector< double > &  I0,
std::vector< Point2D > &  I3 
)
virtual

Definition at line 631 of file Airfoil2D.cpp.

632{
633 if (&pointsDb == &(W.getWake()))
634 {
635 viscousStress.clear();
636 viscousStress.resize(r_.size(), 0.0);
637 }
638
639 GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(pointsDb, domainRadius, I0, I3);
640}//GetDiffVelocityI0I3ToWakeAndViscousStresses(...)
virtual void GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек,...
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetGabarits()

void Airfoil::GetGabarits ( double  gap = 0.02)
virtual

Вычисляет габаритный прямоугольник профиля

Заполняет значения полей lowLeft и upRight по габаритному прямоугольнику профиля с учетом зазора

Parameters
[in]gapотносительная величина зазора в долях от размера габаритного прямоугольника (по умолчанию 0.02, что соответствует 2 %)

Definition at line 1027 of file Airfoil2D.cpp.

1028{
1029 lowLeft = { 1E+10, 1E+10 };
1030 upRight = { -1E+10, -1E+10 };
1031
1032 for (size_t i = 0; i < r_.size(); ++i)
1033 {
1034 lowLeft[0] = std::min(lowLeft[0], r_[i][0]);
1035 lowLeft[1] = std::min(lowLeft[1], r_[i][1]);
1036
1037 upRight[0] = std::max(upRight[0], r_[i][0]);
1038 upRight[1] = std::max(upRight[1], r_[i][1]);
1039 }
1040
1041 Point2D size = upRight - lowLeft;
1042 lowLeft -= size * gap;
1043 upRight += size * gap;
1044}//GetGabarits(...)
Point2D upRight
Правый верхний угол габаритного прямоугольника профиля
Definition Airfoil2D.h:271
Point2D lowLeft
Левый нижний угол габаритного прямоугольника профиля
Definition Airfoil2D.h:270
Here is the caller graph for this function:

◆ GetInfluenceFromSourceSheetToVortex()

void Airfoil::GetInfluenceFromSourceSheetToVortex ( size_t  panel,
const Vortex2D vtx,
Point2D vel 
) const
virtual

Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения

Parameters
[in]panelномер панели профиля, от которой считается влияние
[in]vtxссылка на вихрь
[out]velссылка на вектор полученной скорости

Definition at line 1326 of file Airfoil2D.cpp.

1327{
1329}//GetInfluenceFromSourceSheetToVortex(...)
virtual void GetInfluenceFromSourceSheetAtRectPanelToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения
Here is the call graph for this function:

◆ GetInfluenceFromSourcesToPanel()

void Airfoil::GetInfluenceFromSourcesToPanel ( size_t  panel,
const Vortex2D ptr,
ptrdiff_t  count,
std::vector< double > &  panelRhs 
) const
virtual

Вычисление влияния части подряд идущих источников из области течения на панель для правой части

Вычисляет влияния части подряд идущих источников из области течения на панель для правой части

Parameters
[in]panelномер панели профиля, на которую считается влияние
[in]ptrуказатель на начало диапазона источников
[in]countдлина диапазона источников
[out]panelRhsссылка на вектор полученного влияния для правой части СЛАУ для конкретной панели

Definition at line 1320 of file Airfoil2D.cpp.

1321{
1323}//GetInfluenceFromSourcesToPanel(...)
virtual void GetInfluenceFromSourcesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const =0
Вычисление влияния части подряд источников из области течения на прямолинейную панель для правой част...
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetInfluenceFromVInfToPanel()

void Airfoil::GetInfluenceFromVInfToPanel ( std::vector< double > &  vInfRhs) const
virtual

Вычисление влияния набегающего потока на панель для правой части

Вычисляет влияния набегающего потока на панель для правой части

Parameters
[out]vInfRhsссылка на вектор полученного влияния для правой части СЛАУ для всех панелей профиля

Definition at line 1337 of file Airfoil2D.cpp.

1338{
1340}//GetInfluenceFromVInfToPanel(...)
virtual void GetInfluenceFromVInfToRectPanel(std::vector< double > &vInfRhs) const =0
Вычисление влияния набегающего потока на прямолинейную панель для правой части
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetInfluenceFromVortexSheetToVortex()

void Airfoil::GetInfluenceFromVortexSheetToVortex ( size_t  panel,
const Vortex2D vtx,
Point2D vel 
) const
virtual

Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вихрь в области течения

Parameters
[in]panelномер панели профиля, от которой считается влияние
[in]vtxссылка на вихрь
[out]velссылка на вектор полученной скорости

Definition at line 1332 of file Airfoil2D.cpp.

1333{
1335}//GetInfluenceFromVortexSheetToVortex(...)
virtual void GetInfluenceFromVortexSheetAtRectPanelToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вих...
Here is the call graph for this function:

◆ GetInfluenceFromVorticesToPanel()

void Airfoil::GetInfluenceFromVorticesToPanel ( size_t  panel,
const Vortex2D ptr,
ptrdiff_t  count,
std::vector< double > &  panelRhs 
) const
virtual

Вычисление влияния части подряд идущих вихрей из вихревого следа на панель для правой части

Вычисляет влияния части подряд идущих вихрей из вихревого следа на панель для правой части

Parameters
[in]panelномер панели профиля, на которую считается влияние
[in]ptrуказатель на начало диапазона вихрей
[in]countдлина диапазона вихрей
[out]panelRhsссылка на вектор полученного влияния для правой части СЛАУ для конкретной панели

Definition at line 1313 of file Airfoil2D.cpp.

1314{
1316}//GetInfluenceFromVorticesToPanel(...)
virtual void GetInfluenceFromVorticesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const =0
Вычисление влияния части подряд идущих вихрей из вихревого следа на прямолинейную панель для правой ч...
Here is the call graph for this function:
Here is the caller graph for this function:

◆ getNumberOfPanels()

size_t VM2D::AirfoilGeometry::getNumberOfPanels ( ) const
inlineinherited

Возврат количества панелей на профиле

Returns
количество панелей на профиле

Definition at line 163 of file Airfoil2D.h.

163{ return r_.size(); };
Here is the caller graph for this function:

◆ getR()

const Point2D & VM2D::AirfoilGeometry::getR ( size_t  q) const
inlineinherited

Возврат константной ссылки на вершину профиля

Организовано "зацикливание" в сторону увеличения индекса, т.е. getR[size()] = getR[0];

Это позволяет удобно обращаться к getR(i) и getR(i+1) как к началу и концу i-й панели

Parameters
[in]qномер вершины профиля return константную ссылку на вершину профиля

Definition at line 113 of file Airfoil2D.h.

114 {
115 return (q < r_.size()) ? r_[q] : r_[0];
116 };
Here is the caller graph for this function:

◆ getV()

const Point2D & VM2D::AirfoilGeometry::getV ( size_t  q) const
inlineinherited

Возврат константной ссылки на скорость вершины профиля

Организовано "зацикливание" в сторону увеличения индекса, т.е. getV[size()] = getV[0];

Это позволяет удобно обращаться к getV(i) и getV(i+1) как к скоростям начала и конца i-й панели

Parameters
[in]qномер вершины профиля return константную ссылку на скорость вершины профиля

Definition at line 137 of file Airfoil2D.h.

138 {
139 return (q < v_.size()) ? v_[q] : v_[0];
140 };
std::vector< Point2D > v_
Скорости начал панелей
Definition Airfoil2D.h:72
Here is the caller graph for this function:

◆ IFCUDA() [1/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devAttachedSourceSheetLinPtr)

Указатель на девайсе, где хранятся интенсивности (линейные) присоединенного слоя источников на панелях

◆ IFCUDA() [2/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devAttachedSourceSheetPtr)

Указатель на девайсе, где хранятся интенсивности (константные) присоединенного слоя источников на панелях

◆ IFCUDA() [3/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devAttachedVortexSheetLinPtr)

Указатель на девайсе, где хранятся интенсивности (линейные) присоединенного вихревого слоя на панелях

◆ IFCUDA() [4/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devAttachedVortexSheetPtr)

Указатель на девайсе, где хранятся интенсивности (константные) присоединенного вихревого слоя на панелях

◆ IFCUDA() [5/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devFreeVortexSheetLinPtr)

Указатель на девайсе, где хранятся интенсивности (линейные) свободного вихревого слоя на панелях

◆ IFCUDA() [6/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devFreeVortexSheetPtr)

Указатель на девайсе, где хранятся интенсивности (константные) свободного вихревого слоя на панелях

◆ IFCUDA() [7/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devMeanEpsOverPanelPtr)

Указатель на девайсе, где хранятся средние eps на панелях

◆ IFCUDA() [8/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devPsnPtr)

Указатель на девайсе, где хранятся псевдонормали профиля

◆ IFCUDA() [9/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devRhsLinPtr)

Указатель на девайсе, где хранится правая часть (линейная) матрицы

◆ IFCUDA() [10/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devRhsPtr)

Указатель на девайсе, где хранится правая часть (константная) матрицы

◆ IFCUDA() [11/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devRPtr)

Указатель на девайсе, где хранятся вершины профиля

◆ IFCUDA() [12/14]

VM2D::Airfoil::IFCUDA ( mutable double *  devViscousStressesPtr)

Указатель на девайсе, где хранится вектор (по панелям) для силы вязкого трения

◆ IFCUDA() [13/14]

VM2D::Airfoil::IFCUDA ( mutable std::vector< double >  tmpRhs)

Указатель на хосте, где хранится временная часть матрицы, полученная с девайса

◆ IFCUDA() [14/14]

VM2D::Airfoil::IFCUDA ( mutable std::vector< double >  tmpViscousStresses)

Указатель на хосте, где хранится временная часть вектора (по панелям) для силы вязкого трения

◆ isAfter()

bool Airfoil::isAfter ( size_t  i,
size_t  j 
) const

Проверка, идет ли вершина i следом за вершиной j.

Parameters
[in]iпроверяемая вершина
[in]jконтрольная вершина
Returns
true, если i-я вершина следует зп j-й в порядке обхода профиля

Definition at line 73 of file Airfoil2D.cpp.

74{
75 return ((i == j + 1) || (i == 0 && j == r_.size() - 1));
76}//isAfter(...)
Here is the caller graph for this function:

◆ isInsideGabarits()

bool Airfoil::isInsideGabarits ( const Point2D r) const

Определяет, находится ли точка с радиус-вектором \( \vec r \) внутри габаритного прямоугольника профиля

Parameters
[in]rконстантная ссылка на текущее положение точки
Returns
true, если точка внутри габаритного прямоугольника

Definition at line 80 of file Airfoil2D.cpp.

81{
82 return (r[0] <= upRight[0] && (r[0] >= lowLeft[0] && r[1] >= lowLeft[1] && r[1] <= upRight[1]));
83}//isInsideGabarits(...)

◆ isOutsideGabarits()

bool Airfoil::isOutsideGabarits ( const Point2D r) const

Определяет, находится ли точка с радиус-вектором \( \vec r \) вне габаритного прямоугольника профиля

Parameters
[in]rконстантная ссылка на текущее положение точки
Returns
true, если точка вне габаритного прямоугольника

Definition at line 87 of file Airfoil2D.cpp.

88{
89 return (r[0] > upRight[0] || (r[0] < lowLeft[0] || r[1] < lowLeft[1] || r[1] > upRight[1]));
90}//isOutsideGabarits(...)
Here is the caller graph for this function:

◆ IsPointInAirfoil()

bool Airfoil::IsPointInAirfoil ( const Point2D point) const
virtual

Определяет, находится ли точка с радиус-вектором \( \vec r \) внутри профиля

Parameters
[in]pointконстантная ссылка на текущее положение точки
Returns
true, если точка внутри профиля

Definition at line 1007 of file Airfoil2D.cpp.

1008{
1009 double sumAngle = 0.0;
1010
1011 Point2D v1, v2;
1012
1013 for (size_t i = 0; i < r_.size(); ++i)
1014 {
1015 v1 = getR(i) - point;
1016 v2 = getR(i + 1) - point;
1017 sumAngle += atan2(v1 ^ v2, v1 & v2);
1018 }
1019
1020 if (fabs(sumAngle) < 0.1)
1021 return false;
1022 else
1023 return true;
1024}//IsPointInAirfoil(...)
Here is the call graph for this function:

◆ lightningTest()

void Airfoil::lightningTest ( )

Тест на "отвещенность".

Definition at line 117 of file Airfoil2D.cpp.

118{
120 {
121
122 //Тестируем на "освещенность"
123 wayToVertex.clear();
124 wayToVertex.resize(getNumberOfPanels(), -1);
125
126 for (size_t i = 0; i < getNumberOfPanels(); ++i)
127 {
128 const Point2D& dest = 0.5 * (getR(i) + getR(i + 1));
129
130 bool flag = false;
131 int wayNumber = -1;
132 size_t j;
133
134 do
135 {
136 flag = false;
137 ++wayNumber;
138 j = 0;
139 for (; j < getNumberOfPanels(); ++j)
140 {
141 const Point2D& aflRj = getR(j);
142 const Point2D& aflRj1 = getR(j + 1);
143
144 auto check = [aflRj, aflRj1](const Point2D& start, const Point2D& finish)
145 {
146 return ((((aflRj - start) ^ (finish - start)) * ((aflRj1 - start) ^ (finish - start)) <= 0) && \
147 (((start - aflRj) ^ (aflRj1 - aflRj)) * ((finish - aflRj) ^ (aflRj1 - aflRj)) <= 0));
148 };
149
150 if (i == j)
151 continue;
152
153 for (size_t q = 0; (!flag) && (q < ((wayNumber == 0) ? 1 : possibleWays[wayNumber - 1].size() + 1)); ++q)
154 {
155 Point2D start = ((q == 0) ? rcm : possibleWays[wayNumber - 1][q - 1]);
156 Point2D finish = (wayNumber == 0 || q == possibleWays[wayNumber - 1].size()) ? dest : possibleWays[wayNumber - 1][q];
157
158 flag = (flag || check(start, finish));
159 }
160
161 if (flag)
162 break;
163 }//for j
164
165 if (j == getNumberOfPanels())
166 {
167 wayToVertex[i] = wayNumber;
168 break;
169 }
170 } while ((wayToVertex[i] == -1) && (wayNumber < possibleWays.size()));
171
172 if (wayToVertex[i] == -1)
173 {
174 std::cout << "Possible way to vertex inside airfoil not found" << std::endl;
175 std::cout << "!!!" << std::endl;
176 std::cout << "dest = " << 0.5 * (getR(i) + getR(i + 1)) << std::endl;
177 std::cout << "j = " << j << std::endl;
178 //std::cout << "q = " << q << ", path[q] = " << path[q] << std::endl;
179 //std::cout << "i = " << i << std::endl;
180 std::cout << "rcm = " << rcm << std::endl;
181 std::cout << "possibleWays.size() = " << possibleWays.size() << std::endl;
182
183 for (size_t w = 0; w < possibleWays.size(); ++w)
184 {
185 std::cout << "Possible way # " << w << std::endl;
186 for (size_t p = 0; p < possibleWays[w].size(); ++p)
187 {
188 std::cout << possibleWays[w][p] << " ";
189 }
190 std::cout << std::endl;
191 }
192
193
194 std::ofstream airfoilFileStep(W.getPassport().dir + "afl" + std::to_string(W.getCurrentStep()));
195 for (size_t s = 0; s < getNumberOfPanels(); ++s)
196 airfoilFileStep << getR(s)[0] << " " << getR(s)[1] << "\n";
197 airfoilFileStep.close();
198
199 exit(767);
200 }
201 }//for i
202 }
203}//lightningTest()
Point2D rcm
Положение центра масс профиля
Definition Airfoil2D.h:97
std::vector< int > wayToVertex
Номера путей к вершинам
Definition Airfoil2D.h:202
std::vector< std::vector< Point2D > > possibleWays
Возможные пути внутри профиля от точки (0, 0) к центрам всех панелей
Definition Airfoil2D.h:199
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
std::string dir
Рабочий каталог задачи
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition WorldGen.h:99
size_t size() const
Definition numvector.h:114
int saveVPstep
Шаг вычисления и сохранения скорости и давления
Definition PassportGen.h:80
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Move()

void Airfoil::Move ( const Point2D dr)
virtual

Перемещение профиля

Плоскопараллельно перемещает профиль на вектор \( \overrightarrow{dr} \)

Parameters
[in]drконстантная ссылка на вектор перемещения

Definition at line 1086 of file Airfoil2D.cpp.

1087{
1088 for (size_t i = 0; i < r_.size(); ++i)
1089 r_[i] += dr;
1090 rcm += dr;
1091
1092 for (size_t q = 0; q < possibleWays.size(); ++q)
1093 for (Point2D& pts : possibleWays[q])
1094 pts += dr;
1095
1096 CalcNrmTauLen();
1097 GetGabarits();
1098}//Move(...)
virtual void GetGabarits(double gap=0.02)
Вычисляет габаритный прямоугольник профиля
void CalcNrmTauLen()
Вычисление нормалей, касательных и длин панелей по текущему положению вершин
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReadFromFile()

void Airfoil::ReadFromFile ( const std::string &  dir)
virtual

Считывание профиля из файла

Считывание геометрии профиля из файла, вычисление всех прочих параметров профиля
После загрузки из файла профиль поворачивается на нужный угол и масштабируется на нужный коэффициент

Warning
Сейчас масса, момент инерции и скорости вершин зануляются.
Parameters
[in]dirконстантная ссылка на строку — имя каталога, где лежит cчитываемый файл

Definition at line 209 of file Airfoil2D.cpp.

210{
212 std::string filename = dir + param.fileAirfoil;
213 std::ifstream airfoilFile;
214
215 if (fileExistTest(filename, W.getInfo(), true, { "txt", "TXT" }))
216 {
217 std::stringstream airfoilFile(VMlib::Preprocessor(filename).resultString);
218
219 VMlib::StreamParser airfoilParser(W.getInfo(), "airfoil file parser", airfoilFile);
220
221 m = 0.0; //TODO
222 J = 0.0; //TODO
223 rcm = { 0.0, 0.0 }; //TODO
224 phiAfl = 0.0;
225
226 inverse = param.inverse;
227
228
229
230 //*
231
233 //Проверяем отсутствие "задвоенных" подряд идущих точек
234 //airfoilParser.get("r", r_);
235 std::vector<Point2D> rFromFile;
236 airfoilParser.get("r", rFromFile);
237
238 if (rFromFile.size() > 0)
239 {
240 if (param.requiredNPanels <= rFromFile.size()) {
241
242 r_.reserve(rFromFile.size());
243 r_.push_back(rFromFile[0]);
244 for (size_t i = 1; i < rFromFile.size(); ++i)
245 if ((rFromFile[i] - rFromFile[i - 1]).length2() > 1e-12)
246 r_.push_back(rFromFile[i]);
247
248 //Если первая совпадает с последней, то убираем последнюю
249 if ((r_.back() - r_.front()).length2() < 1e-12)
250 r_.resize(r_.size() - 1);
251 }
252 else
253 {
254 W.getInfo('e') << "Airfoil shape is given explicitely, it can not be automatically split!" << std::endl;
255 exit(200);
256 }
257 }
258 else
259 {
260 std::vector<GeomPoint> geomFromFile;
261
262 airfoilParser.get("geometry", geomFromFile);
263
264 if ((geomFromFile.front() - geomFromFile.back()).length2() < 1e-12)
265 geomFromFile.resize(geomFromFile.size() - 1);
266
267 size_t reqN = param.requiredNPanels;
268
269 std::vector<double> L(geomFromFile.size());
270 double totalLength = 0.0;
271 std::vector<size_t> nc = {};
272 std::vector<size_t> ni = {};
273
274 for (size_t i = 0; i < geomFromFile.size(); ++i)
275 {
276 const Point2D& p1 = geomFromFile[i];
277 const Point2D& p2 = geomFromFile[(i + 1) % geomFromFile.size()];
278
279 L[i] = (p2 - p1).length();
280 totalLength += L[i];
281 }
282
283 /*
284 Point2D p1, p2, p3, cc, dd;
285 double ie = 0;
286 int numOfDivPnt = 12;
287 double alpha;
288 double rscale = 0.1;
289 std::vector<Point2D> pnt(numOfDivPnt);
290 std::vector<Point2D> newGeom(numOfDivPnt * geomFromFile.size());
291
292 for (size_t i = 0; i < geomFromFile.size(); ++i)
293 {
294 p1 = (i == 0 ? geomFromFile[(int)geomFromFile.size() - 1].r : geomFromFile[i - 1].r);
295 p2 = geomFromFile[i].r;
296 p3 = geomFromFile[(i + 1) % geomFromFile.size()].r;
297
298 Point2D dc = (p1 - p2).unit() + (p3 - p2).unit();
299 ie = dc.dist2To(dc.proj(p1 - p2));
300 if (ie < 1e-3) ie = 0;
301
302 alpha = (PI - acos(((p3 - p2) & (p1 - p2)) / (sqrt(((p3 - p2) & (p3 - p2)) * ((p1 - p2) & (p1 - p2)))))) / (numOfDivPnt - 1);
303
304 cc = p2 + (rscale * ie) * dc;
305
306 dd = p2 + (cc - p2).proj(p1 - p2) - cc;
307 for (int j = 0; j < numOfDivPnt; ++j)
308 newGeom[i * numOfDivPnt + j] = cc + dd.rotated(j * alpha);
309 }
310 */
311
312 //Ищем первый сплайн
313 size_t i = 0;
314 std::vector<size_t> splineStart, splineFinish;
315
316 while (i < geomFromFile.size())
317 {
318 while (i < geomFromFile.size() && !(geomFromFile[i].type == "c"))
319 ++i;
320
321 if (i < geomFromFile.size())
322 {
323 splineStart.push_back(i);
324 ++i;
325 //ищем второй конец
326 while (i < geomFromFile.size() && !(geomFromFile[i].type == "c"))
327 ++i;
328
329 splineFinish.push_back(i);
330 }
331 }
332
333 if ((splineStart.size() > 0) && (splineStart[0] != 0))
334 splineFinish.back() = splineStart[0];
335
336 if (reqN == 0)
337 {
338 r_.reserve(geomFromFile.size());
339 for (size_t i = 0; i < geomFromFile.size(); ++i)
340 r_.push_back(geomFromFile[i]);
341 }
342 else
343 {
344 double hRef = totalLength / reqN;
345
346 std::vector<int> nPanels;
347
348 bool cyclic = (splineStart.size() == 0);
349
350 if (cyclic)
351 {
352 splineStart.push_back(0);
353 splineFinish.push_back(geomFromFile.size());
354 }
355
356 for (size_t s = 0; s < splineStart.size(); ++s)
357 {
358 //Сплайн
359 std::vector<double> T, X, Y;
360
361 size_t splineLegs = splineFinish[s] - splineStart[s];
362 if (splineFinish[s] <= splineStart[s])
363 splineLegs += geomFromFile.size();
364
365 T.reserve(splineLegs + 1);
366 X.reserve(T.size());
367 Y.reserve(T.size());
368
369 double tbuf = 0.0;
370
371 for (size_t i = splineStart[s]; i <= splineStart[s] + splineLegs; ++i)
372 {
373 const Point2D& pt = geomFromFile[i % geomFromFile.size()];
374
375 T.push_back(tbuf);
376 X.push_back(pt[0]);
377 Y.push_back(pt[1]);
378
379 if ((i == splineStart[s]) && (splineFinish[s] - splineStart[s] == 1))
380 {
381 double halfL = (i < geomFromFile.size()) ? 0.5 * L[i] : 0.5 * (geomFromFile.front() - geomFromFile.back()).length();
382 tbuf += halfL;
383 T.push_back(tbuf);
384 Point2D nextPt = geomFromFile[(i + 1) % geomFromFile.size()];
385
386 X.push_back(0.5 * (pt[0] + nextPt[0]));
387 Y.push_back(0.5 * (pt[1] + nextPt[1]));
388
389 tbuf += halfL;
390 }
391 else
392 tbuf += (i < geomFromFile.size()) ? L[i] : (geomFromFile.front() - geomFromFile.back()).length();
393 }
394
395 tk::spline s1, s2;
396 if (cyclic)
397 {
398 s1.set_boundary(tk::spline::cyclic);
399 s2.set_boundary(tk::spline::cyclic);
400 }
401 else
402 {
403 s1.set_boundary(tk::spline::second_deriv, 0.0, tk::spline::second_deriv, 0.0);
404 s2.set_boundary(tk::spline::second_deriv, 0.0, tk::spline::second_deriv, 0.0);
405 }
406
407 s1.set_points(T, X);
408 s2.set_points(T, Y);
409
410 int NumberOfPanelsSpline = (int)ceil(T.back() / hRef);
411 double hSpline = T.back() / NumberOfPanelsSpline;
412
413 for (int i = 0; i < NumberOfPanelsSpline; ++i)
414 r_.push_back({ s1(hSpline * i), s2(hSpline * i) });
415
416 nPanels.push_back(NumberOfPanelsSpline);
417
418 }
419 }
420 }//else geometry
421
422 if (r_.size() == 0)
423 {
424 W.getInfo('e') << "No points on airfoil contour!" << std::endl;
425 exit(200);
426 }
427
428 if (inverse)
429 std::reverse(r_.begin(), r_.end());
430
431 v_.resize(0);
432 for (size_t q = 0; q < r_.size(); ++q)
433 v_.push_back({ 0.0, 0.0 });
434
435
436 int nPossibleWays;
437 int defaultNPossibleWays = 0;
438 airfoilParser.get("nPossibleWays", nPossibleWays, &defaultNPossibleWays, false);
439
440 possibleWays.resize(nPossibleWays);
441 for (int q = 0; q < nPossibleWays; ++q)
442 airfoilParser.get("possibleWay" + std::to_string(q + 1), possibleWays[q]);
443
444 //определяем начальные габаритные размеры
445 auto xMinMax = std::minmax_element(r_.begin(), r_.end(), Point2D::cmp<'x'>);
446 auto yMinMax = std::minmax_element(r_.begin(), r_.end(), Point2D::cmp<'y'>);
447
449 { Point2D({ (*xMinMax.first)[0], (*yMinMax.first)[1] }), Point2D({ (*xMinMax.second)[0], (*yMinMax.second)[1] }) };
450
451 Move(param.basePoint);
452 Scale(param.scale);
453
454 double rotationAngle = param.angle;
455
456 //Для обдува ветром, когда углы считаются по компасу
457 //if (W.getPassport().geographicalAngles)
458 // rotationAngle = -param.angle - 0.5 * PI;
459
460 Rotate(-rotationAngle);
461 //в конце Rotate нормали, касательные и длины вычисляются сами
462
463 std::string fname = W.getPassport().dir + "/airfoil-" + std::to_string(numberInPassport) + ".pfl";
464 //if (!fileExistTest(fname, W.getInfo()))
465 {
466 std::ofstream of(fname);
467 for (size_t i = 0; i < r_.size(); ++i)
468 of << r_[i][0] << " " << r_[i][1] << std::endl;
469 of.close();
470 }
471
472
473 gammaThrough.clear();
474 gammaThrough.resize(r_.size(), 0.0);
475
476 //Вычисляем площадь
477 area = 0.0;
478 for (size_t q = 0; q < r_.size(); ++q)
479 {
480 const Point2D& cntq = 0.5 * (getR(q) + getR(q + 1));
481 const Point2D& drq = getR(q + 1) - getR(q);
482 area += (cntq[0] * drq[1] - cntq[1] * drq[0]);
483 }
484 area *= 0.5;
485 area = fabs(area);
486
488 }
489}//ReadFromFile(...)
double phiAfl
Поворот профиля
Definition Airfoil2D.h:100
double area
Площадь профиля
Definition Airfoil2D.h:103
bool inverse
Признак разворота нормалей (для расчета внутренних течений)
Definition Airfoil2D.h:76
virtual void Scale(const Point2D &)
Масштабирование профиля
std::vector< double > gammaThrough
Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
Definition Airfoil2D.h:276
double J
Полярный момент инерции профиля относительно центра масс
Definition Airfoil2D.h:196
virtual void Move(const Point2D &dr)
Перемещение профиля
void lightningTest()
Тест на "отвещенность".
virtual void Rotate(double alpha)
Поворот профиля
std::vector< AirfoilParams > airfoilParams
Список структур с параметрами профилей
Definition Passport2D.h:268
Passport & getNonConstPassport() const
Возврат неконстантной ссылки на паспорт
Definition World2D.h:256
Класс, позволяющий выполнять предварительную обработку файлов
Класс, позволяющий выполнять разбор файлов и строк с настройками и параметрами
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода
Definition WorldGen.h:82
bool fileExistTest(std::string &fileName, LogStream &info, bool exitKey=false, const std::list< std::string > &extList={})
Проверка существования файла
Definition defs.h:342
Point2D basePoint
Смещение центра масс (перенос профиля)
Definition Passport2D.h:210
std::string fileAirfoil
Имя файла с начальным состоянием профилей (без полного пути)
Definition Passport2D.h:201
double angle
Угол поворота (угол атаки)
Definition Passport2D.h:219
Point2D scale
Коэффициент масштабирования
Definition Passport2D.h:213
size_t requiredNPanels
Желаемое число панелей для разбиения геометрии
Definition Passport2D.h:207
bool inverse
Признак разворота нормалей (для расчета внутреннего течения)
Definition Passport2D.h:225
Here is the call graph for this function:

◆ Rotate()

void Airfoil::Rotate ( double  alpha)
virtual

Поворот профиля

Поворачивает профиль на угол \( \alpha \) вокруг центра масс

Parameters
[in]alphaугол поворота против часовой стрелки в радианах

Definition at line 1101 of file Airfoil2D.cpp.

1102{
1103 phiAfl += alpha;
1104 nummatrix<double, 2, 2> rotMatrix = { { cos(alpha), -sin(alpha) }, { sin(alpha), cos(alpha) } };
1105
1106 for (size_t i = 0; i < r_.size(); ++i)
1107 r_[i] = rcm + (rotMatrix & (r_[i] - rcm));
1108
1109 for (size_t q = 0; q < possibleWays.size(); ++q)
1110 for (Point2D& pts : possibleWays[q])
1111 {
1112 Point2D oldPts = pts;
1113 pts = rcm + (rotMatrix & (oldPts - rcm));
1114 }
1115
1116 CalcNrmTauLen();
1117 GetGabarits();
1118}//Rotate(...)
Шаблонный класс, определяющий матрицу фиксированного размера Фактически представляет собой массив,...
Definition nummatrix.h:66
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Scale()

void Airfoil::Scale ( const Point2D factor)
virtual

Масштабирование профиля

Масштабирует профиль в factor раз относительно центра масс

Parameters
[in]factorмасштабный коэффициент

Definition at line 1122 of file Airfoil2D.cpp.

1123{
1124 for (size_t i = 0; i < r_.size(); ++i)
1125 r_[i] = rcm + Point2D{ factor[0] * (r_[i] - rcm)[0], factor[1] * (r_[i] - rcm)[1] };
1126
1127 for (size_t q = 0; q < possibleWays.size(); ++q)
1128 for (Point2D& pts : possibleWays[q])
1129 {
1130 Point2D oldPts = pts;
1131 pts = rcm + Point2D{ factor[0] * (oldPts - rcm)[0], factor[1] * (oldPts - rcm)[1] };
1132 }
1133
1134 CalcNrmTauLen();
1135 GetGabarits();
1136}//Scale(...)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setR()

Point2D & VM2D::AirfoilGeometry::setR ( size_t  q)
inlineinherited

Возврат ссылки на вершину профиля

Организовано "зацикливание" в сторону увеличения индекса, т.е. getR[size()] = getR[0];

Это позволяет удобно обращаться к getR(i) и getR(i+1) как к началу и концу i-й панели

Parameters
[in]qномер вершины профиля return ссылку на вершину профиля

Definition at line 125 of file Airfoil2D.h.

126 {
127 return (q < r_.size()) ? r_[q] : r_[0];
128 };
Here is the caller graph for this function:

◆ setV() [1/2]

void VM2D::AirfoilGeometry::setV ( const Point2D vel)
inlineinherited

Установка постоянной скорости всех вершин профиля

Parameters
[in]velконстантная ссылка на величину устанавливаемой скорости

Definition at line 145 of file Airfoil2D.h.

146 {
147 v_.clear();
148 v_.resize(r_.size(), vel);
149 }
Here is the caller graph for this function:

◆ setV() [2/2]

void VM2D::AirfoilGeometry::setV ( const std::vector< Point2D > &  vel)
inlineinherited

Установка скоростей всех вершин профиля

Parameters
[in]velконстантная ссылка на вектор величин скоростей вершин профиля

Definition at line 154 of file Airfoil2D.h.

155 {
156 v_.clear();
157 v_.insert(v_.end(), vel.begin(), vel.end());
158 }

Member Data Documentation

◆ area

double VM2D::AirfoilGeometry::area
inherited

Площадь профиля

Definition at line 103 of file Airfoil2D.h.

◆ gammaThrough

std::vector<double> VM2D::Airfoil::gammaThrough

Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге

Используются в правой части системы, чтобы компенсировать вихри, проникшие в профиль

Definition at line 276 of file Airfoil2D.h.

◆ inverse

bool VM2D::AirfoilGeometry::inverse
inherited

Признак разворота нормалей (для расчета внутренних течений)

Definition at line 76 of file Airfoil2D.h.

◆ J

double VM2D::Airfoil::J

Полярный момент инерции профиля относительно центра масс

Definition at line 196 of file Airfoil2D.h.

◆ len

std::vector<double> VM2D::AirfoilGeometry::len
inherited

Длины панелей профиля

Definition at line 94 of file Airfoil2D.h.

◆ lowLeft

Point2D VM2D::Airfoil::lowLeft

Левый нижний угол габаритного прямоугольника профиля

Definition at line 270 of file Airfoil2D.h.

◆ m

double VM2D::Airfoil::m

Масса профиля

Definition at line 193 of file Airfoil2D.h.

◆ meanEpsOverPanel

std::vector<double> VM2D::Airfoil::meanEpsOverPanel

Средние значения Eps на панелях

Definition at line 209 of file Airfoil2D.h.

◆ nrm

std::vector<Point2D> VM2D::AirfoilGeometry::nrm
inherited

Нормали к панелям профиля

Нормали задаются внешними, нормированными на единицу

Definition at line 81 of file Airfoil2D.h.

◆ numberInPassport

const size_t VM2D::Airfoil::numberInPassport

Номер профиля в паспорте

Definition at line 188 of file Airfoil2D.h.

◆ phiAfl

double VM2D::AirfoilGeometry::phiAfl
inherited

Поворот профиля

Definition at line 100 of file Airfoil2D.h.

◆ possibleWays

std::vector<std::vector<Point2D> > VM2D::Airfoil::possibleWays

Возможные пути внутри профиля от точки (0, 0) к центрам всех панелей

Definition at line 199 of file Airfoil2D.h.

◆ psn

std::vector<std::pair<Point2D, Point2D> > VM2D::AirfoilGeometry::psn
inherited

Псевдонормали к панелям профиля

Псевдонормали к началу и концу панели, нормированные на единицу

Definition at line 86 of file Airfoil2D.h.

◆ r_

std::vector<Point2D> VM2D::AirfoilGeometry::r_
protectedinherited

Координаты начал панелей

Definition at line 69 of file Airfoil2D.h.

◆ rcm

Point2D VM2D::AirfoilGeometry::rcm
inherited

Положение центра масс профиля

Definition at line 97 of file Airfoil2D.h.

◆ tau

std::vector<Point2D> VM2D::AirfoilGeometry::tau
inherited

Касательные к панелям профиля

Касательные соответствуют обходу профиля против часовой стрелки, задаются нормированными на единицу

Definition at line 91 of file Airfoil2D.h.

◆ upRight

Point2D VM2D::Airfoil::upRight

Правый верхний угол габаритного прямоугольника профиля

Definition at line 271 of file Airfoil2D.h.

◆ v_

std::vector<Point2D> VM2D::AirfoilGeometry::v_
protectedinherited

Скорости начал панелей

Definition at line 72 of file Airfoil2D.h.

◆ viscousStress

std::vector<double> VM2D::Airfoil::viscousStress

Нейросеть для коэффициентов I0 и I3 диффузионной скорости

Касательные напряжения на панелях профиля

Definition at line 268 of file Airfoil2D.h.

◆ W

const World2D& VM2D::Airfoil::W

Константная ссылка на решаемую задачу

Definition at line 185 of file Airfoil2D.h.

◆ wayToVertex

std::vector<int> VM2D::Airfoil::wayToVertex

Номера путей к вершинам

Definition at line 202 of file Airfoil2D.h.


The documentation for this class was generated from the following files: