VM2D  1.12
Vortex methods for 2D flows simulation
Boundary2DVortexCollocN.cpp
Go to the documentation of this file.
1 /*--------------------------------*- VM2D -*-----------------*---------------*\
2 | ## ## ## ## #### ##### | | Version 1.12 |
3 | ## ## ### ### ## ## ## ## | VM2D: Vortex Method | 2024/01/14 |
4 | ## ## ## # ## ## ## ## | for 2D Flow Simulation *----------------*
5 | #### ## ## ## ## ## | Open Source Code |
6 | ## ## ## ###### ##### | https://www.github.com/vortexmethods/VM2D |
7 | |
8 | Copyright (C) 2017-2024 I. Marchevsky, K. Sokol, E. Ryatina, A. Kolganova |
9 *-----------------------------------------------------------------------------*
10 | File name: Boundary2DVortexCollocN.cpp |
11 | Info: Source code of VM2D |
12 | |
13 | This file is part of VM2D. |
14 | VM2D is free software: you can redistribute it and/or modify it |
15 | under the terms of the GNU General Public License as published by |
16 | the Free Software Foundation, either version 3 of the License, or |
17 | (at your option) any later version. |
18 | |
19 | VM2D is distributed in the hope that it will be useful, but WITHOUT |
20 | ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
21 | FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
22 | for more details. |
23 | |
24 | You should have received a copy of the GNU General Public License |
25 | along with VM2D. If not, see <http://www.gnu.org/licenses/>. |
26 \*---------------------------------------------------------------------------*/
27 
28 
43 
44 #include "Airfoil2D.h"
45 #include "MeasureVP2D.h"
46 #include "Mechanics2D.h"
47 #include "Passport2D.h"
48 #include "StreamParser.h"
49 #include "Velocity2D.h"
50 #include "Wake2D.h"
51 #include "World2D.h"
52 
53 using namespace VM2D;
54 
55 //Конструктор
56 BoundaryVortexCollocN::BoundaryVortexCollocN(const World2D& W_, size_t numberInPassport_) :
57  Boundary(W_, numberInPassport_, 1)
58  {
59  c.reserve(afl.getNumberOfPanels());
60  for (size_t i=0; i < afl.getNumberOfPanels(); ++i)
61  c.push_back(0.5*(afl.getR(i)+afl.getR(i+1)));
62  };
63 
64 
65 //Пересчет решения на интенсивность вихревого слоя
67 {
68  Vortex2D virtVort;
69  Point2D midNorm;
70 
71  size_t np = afl.getNumberOfPanels();
72 
74 
76 
77  //Очистка и резервирование памяти
78  virtualWake.vecHalfGamma.clear();
79  virtualWake.vecHalfGamma.reserve(np * nVortPerPan);
80 
81  //Очистка и резервирование памяти
82  virtualWake.aflPan.clear();
83  virtualWake.aflPan.reserve(np * nVortPerPan);
84 
85  //Резервирование памяти
86  virtualWake.vtx.clear();
87  virtualWake.vtx.reserve(np * nVortPerPan);
88 
89  //Очистка и резервирование памяти
90  vortexBeginEnd.clear();
91  vortexBeginEnd.reserve(np);
92 
94 
95  std::pair<int, int> pair;
96 
97  for (size_t i = 0; i < np; ++i)
98  {
99  midNorm = afl.nrm[i] * delta;
100  size_t iNext = (i < np-1) ? i + 1 : 0;
101  size_t iPrev = (i > 0) ? i - 1 : np - 1;
102 
103  double Gamma = (sol(i) / (afl.len[iPrev] + afl.len[i]) + sol(iNext) / (afl.len[iNext] + afl.len[i])) * afl.len[i];
104 
105  size_t NEWnVortPerPan = std::max(static_cast<int>(std::ceil(fabs(Gamma) / W.getPassport().wakeDiscretizationProperties.maxGamma)), nVortPerPan);
106 
107  pair.first = (int)virtualWake.vtx.size();
108 
109 
110  Point2D dr = 1.0 / NEWnVortPerPan * (afl.getR(i + 1) - afl.getR(i));
111  virtVort.g() = Gamma / NEWnVortPerPan;
112 
113  for (size_t j = 0; j < NEWnVortPerPan; ++j)
114  {
115  virtVort.r() = afl.getR(i) + (0.5 + j) * dr + delta * afl.nrm[i];
116 
117  virtualWake.vtx.push_back(virtVort);
118 
119  virtualWake.vecHalfGamma.push_back(0.5 * Gamma / afl.len[i] * afl.tau[i]);
120  virtualWake.aflPan.push_back({ numberInPassport, i });
121  }
122 
123  pair.second = (int)virtualWake.vtx.size();
124  vortexBeginEnd.push_back(pair);
125 
126  sheets.freeVortexSheet(i, 0) = Gamma / afl.len[i];
127  }
128 
129 
130 }//SolutionToFreeVortexSheetAndVirtualVortex(...)
131 
132 
133 //Генерация блока матрицы
134 void BoundaryVortexCollocN::FillMatrixSelf(Eigen::MatrixXd& matr, Eigen::VectorXd& lastLine, Eigen::VectorXd& lactCol)
135 {
136  size_t np = afl.getNumberOfPanels();
137 
138  for (size_t i = 0; i < np; ++i)
139  {
140  lactCol(i) = 1.0;
141  lastLine(i) = 1.0;
142  }
143 
144 
145 #pragma omp parallel for shared(matr)
146  for (int i = 0; i < np; ++i)
147  for (size_t j = 0; j < np; ++j)
148  {
149  matr(i, j) = -(afl.tau[i] & (c[i]-afl.getR(j)) ) * IDPI / (c[i]-afl.getR(j)).length2();
150  }
151 }//FillMatrixSelf(...)
152 
153 void BoundaryVortexCollocN::FillIQSelf(std::pair<Eigen::MatrixXd, Eigen::MatrixXd>& IQ)
154 {
155  //afl.calcIQ(1, afl, IQ);
156 }//FillIQSelf(...)
157 
158 //Генерация блока матрицы влияния от другого профиля того же типа
159 void BoundaryVortexCollocN::FillMatrixFromOther(const Boundary& otherBoundary, Eigen::MatrixXd& matr)
160 {
161 #pragma omp parallel for shared(matr)
162  for (int i = 0; i < afl.getNumberOfPanels(); ++i)
163  for (size_t j = 0; j < otherBoundary.afl.getNumberOfPanels(); ++j)
164  {
165 
166  matr(i, j) = -(afl.tau[i] & (c[i]-otherBoundary.afl.getR(j)) ) * IDPI / (c[i]-otherBoundary.afl.getR(j)).length2();
167  }
168 }//FillMatrixFromOther(...)
169 
170 
171 void BoundaryVortexCollocN::FillIQFromOther(const Boundary& otherBoundary, std::pair<Eigen::MatrixXd, Eigen::MatrixXd>& IQ)
172 {
173  //afl.calcIQ(1, otherBoundary.afl, IQ);
174 }//FillIQFromOther(...)
175 
176 
177 //Вычисление скоростей в наборе точек, вызываемых наличием слоев вихрей и источников на профиле
178 void BoundaryVortexCollocN::CalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase& pointsDb, std::vector<Point2D>& velo) const
179 {
180  std::vector<Point2D> selfVelo(pointsDb.vtx.size());
181 
182  double cft = IDPI;
183 
184 #pragma warning (push)
185 #pragma warning (disable: 4101)
186  //Локальные переменные для цикла
187  Point2D velI;
188  Point2D tempVel;
189  //double dst2eps, dst2;
190 #pragma warning (pop)
191 
192 #pragma omp parallel for default(none) shared(selfVelo, pointsDb, cft, std::cout) private(velI, tempVel)
193  for (int i = 0; i < pointsDb.vtx.size(); ++i)
194  {
195  velI.toZero();
196 
197  const Point2D& posI = pointsDb.vtx[i].r();
198 
201  for (size_t j = 0; j < sheets.getSheetSize(); ++j)
202  {
203  Point2D dj = afl.getR(j + 1) - afl.getR(j);
204  Point2D tauj = dj.unit();
205 
206  Point2D s = posI - afl.getR(j);
207  Point2D p = posI - afl.getR(j + 1);
208 
209  double a = VMlib::Alpha(p, s);
210 
211  double lambda;
212  if ( (s.length2() > 1e-16) && (p.length2() > 1e-16) )
213  lambda = VMlib::Lambda(p, s);
214  else
215  lambda = 0.0;
216 
217  Point2D skos = -a * tauj.kcross() + lambda * tauj;
218 
219  velI += sheets.freeVortexSheet(j, 0) * skos.kcross();
220  velI += sheets.attachedVortexSheet(j, 0) * skos.kcross();
221  velI += sheets.attachedSourceSheet(j, 0) * skos;
222  }//for j
223 
224  velI *= cft;
225  selfVelo[i] = velI;
226 
227  }//for i
228 
229  for (size_t i = 0; i < velo.size(); ++i)
230  velo[i] += selfVelo[i];
231 }//CalcConvVelocityToSetOfPointsFromSheets(...)
232 
233 
234 #if defined(USE_CUDA)
235 void BoundaryVortexCollocN::GPUCalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase& pointsDb, std::vector<Point2D>& velo) const
236 {
237  if (afl.numberInPassport == 0)
238  {
239  const size_t npt = pointsDb.vtx.size();
240  double*& dev_ptr_pt = pointsDb.devVtxPtr;
241 
242  size_t npnl = afl.getNumberOfPanels();
243  for (size_t q = 1; q < W.getNumberOfAirfoil(); ++q)
244  npnl += W.getAirfoil(q).getNumberOfPanels();
245 
246  double*& dev_ptr_r = afl.devRPtr;
247  double*& dev_ptr_freeVortexSheet = afl.devFreeVortexSheetPtr;
248  double*& dev_ptr_attachedVortexSheet = afl.devAttachedVortexSheetPtr;
249  double*& dev_ptr_attachedSourceSheet = afl.devAttachedSourceSheetPtr;
250 
251  std::vector<Point2D>& Vel = velo;
252  std::vector<Point2D> newV(npt);
253  double*& dev_ptr_vel = pointsDb.devVelPtr;
255 
256  //Явная синхронизация слоев не нужна, т.к. она выполняется в Gpu::RefreshAfls()
257  if (npt > 0)
258  {
259  //double tt1 = omp_get_wtime();
260  cuCalculateConvVeloWakeFromVirtual(npt, dev_ptr_pt, npnl, dev_ptr_r, dev_ptr_freeVortexSheet, dev_ptr_attachedVortexSheet, dev_ptr_attachedSourceSheet, dev_ptr_vel, eps2);
261  //double tt2 = omp_get_wtime();
262  //std::cout << "SHEET: " << tt2 - tt1 << std::endl;
263 
264  W.getCuda().CopyMemFromDev<double, 2>(npt, dev_ptr_vel, (double*)newV.data());
265 
266 
267  for (size_t q = 0; q < Vel.size(); ++q)
268  Vel[q] += newV[q];
269  }
270  }
271 }
272 //GPUCalcConvVelocityToSetOfPointsFromSheets(...)
273 #endif
274 
275 
276 void BoundaryVortexCollocN::CalcConvVelocityAtVirtualVortexes(std::vector<Point2D>& velo) const
277 {
278  std::vector<Point2D>& Vel = velo;
279 
280  //Скорости виртуальных вихрей
281 
282 #pragma omp parallel for default(none) shared(Vel)
283  for (int i = 0; i < (int)Vel.size(); ++i)
284  Vel[i] = afl.getV(virtualWake.aflPan[i].second) \
286  - W.getPassport().physicalProperties.V0(); // V0 потом прибавляется ко всем скоростям в функции MoveVortexes
287 }//CalcConvVelocityAtVirtualVortexes(...)
288 
289 
290 //Вычисление интенсивностей присоединенного вихревого слоя и присоединенного слоя источников
292 {
293 
294  for (size_t i = 0; i < sheets.getSheetSize(); ++i)
295  {
298  }
299 
300  for (size_t i = 0; i < sheets.getSheetSize(); ++i)
301  {
302  sheets.attachedVortexSheet(i, 0) = 0.5 * (afl.getV(i) + afl.getV(i + 1)) & afl.tau[i];
303  sheets.attachedSourceSheet(i, 0) = 0.5 * (afl.getV(i) + afl.getV(i + 1)) & afl.nrm[i];
304  }
305 }//ComputeAttachedSheetsIntensity()
306 
307 
308 //Вычисляет влияния части подряд идущих вихрей из вихревого следа на прямолинейную панель для правой части
309 void BoundaryVortexCollocN::GetInfluenceFromVorticesToRectPanel(size_t panel, const Vortex2D* ptr, ptrdiff_t count, std::vector<double>& wakeRhs) const
310 {
311  double& velI = wakeRhs[0];
312 
313  const Point2D& posI = c[panel];
314 
315  for (size_t it = 0; it != count; ++it)
316  {
317  const Vortex2D& vt = ptr[it];
318 
319  const Point2D& posJ = vt.r();
320  const double& gamJ = vt.g();
321 
322  Point2D d = posI - posJ;
323  //std::cout << "d = " << d << std::endl;
324 
325  velI -= (gamJ * afl.len[panel] / std::max(d.length2(), W.getPassport().wakeDiscretizationProperties.eps2)) * (afl.tau[panel] & d);
326  }
327 }// GetInfluenceFromVorticesToRectPanel(...)
328 
329 
330 
331 //Вычисляет влияния части подряд идущих источников в области течения на прямолинейную панель для правой части
332 void BoundaryVortexCollocN::GetInfluenceFromSourcesToRectPanel(size_t panel, const Vortex2D* ptr, ptrdiff_t count, std::vector<double>& wakeRhs) const
333 {
334  double& velI = wakeRhs[0];
335 
336  const Point2D& posI = 0.5 * (afl.getR(panel) + afl.getR(panel + 1));
337 
338  for (size_t it = 0; it != count; ++it)
339  {
340  const Vortex2D& vt = ptr[it];
341 
342  const Point2D& posJ = vt.r();
343  const double& gamJ = vt.g();
344 
345  Point2D d = posI - posJ;
346 
347  velI += (gamJ * afl.len[panel] / std::max(d.length2(), W.getPassport().wakeDiscretizationProperties.eps2)) * (afl.nrm[panel] & d);
348  }
349 }// GetInfluenceFromSourcesToRectPanel(...)
350 
351 
352 //Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения
354 {
355  vel.toZero();
356 
357  const Point2D& posI = ptr.r();
358 
359  Point2D dj = afl.getR(panel + 1) - afl.getR(panel);
360  Point2D tauj = dj.unit();
361 
362  Point2D s = posI - afl.getR(panel);
363  Point2D p = posI - afl.getR(panel + 1);
364 
365  double a = VMlib::Alpha(p, s);
366 
367  double lambda;
368  if ((s.length2() > 1e-16) && (p.length2() > 1e-16))
369  lambda = VMlib::Lambda(p, s);
370  else
371  lambda = 0.0;
372 
373  vel += sheets.attachedSourceSheet(panel, 0) * (-a * tauj.kcross() + lambda * tauj);
374 }// GetInfluenceFromSourceSheetAtRectPanelToVortex(...)
375 
376 //Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вихрь в области течения
378 {
379  vel.toZero();
380 
381  const Point2D& posI = ptr.r();
382 
383  Point2D dj = afl.getR(panel + 1) - afl.getR(panel);
384  Point2D tauj = dj.unit();
385 
386  Point2D s = posI - afl.getR(panel);
387  Point2D p = posI - afl.getR(panel + 1);
388  double a = VMlib::Alpha(p, s);
389 
390  double lambda;
391  if ((s.length2() > 1e-16) && (p.length2() > 1e-16))
392  lambda = VMlib::Lambda(p, s);
393  else
394  lambda = 0.0;
395 
396  Point2D skos = -a * tauj.kcross() + lambda * tauj;
397 
398  vel += sheets.freeVortexSheet(panel, 0) * skos.kcross();
399  vel += sheets.attachedVortexSheet(panel, 0) * skos.kcross();
400 
401 }// GetInfluenceFromVortexSheetAtRectPanelToVortex(...)
402 
403 
404 //Вычисляет влияния набегающего потока на прямолинейную панель для правой части
405 void BoundaryVortexCollocN::GetInfluenceFromVInfToRectPanel(std::vector<double>& vInfRhs) const
406 {
407  size_t np = afl.getNumberOfPanels();
408 
409 #pragma omp parallel for default(none) shared(vInfRhs, np)
410  for (int i = 0; i < np; ++i)
411  vInfRhs[i] = afl.nrm[i] & W.getPassport().physicalProperties.V0();
412 
413 }// GetInfluenceFromVInfToRectPanel(...)
414 
const double & freeVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:99
virtual void CalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase &pointsDb, std::vector< Point2D > &velo) const override
Вычисление конвективных скоростей в наборе точек, вызываемых наличием слоев вихрей и источников на пр...
Заголовочный файл с описанием класса Passport (двумерный) и cоответствующими структурами ...
virtual void FillIQFromOther(const Boundary &otherBoundary, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &IQ) override
Генерация блока матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2, влияние одного профиля на другой...
const double & attachedVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:104
virtual void FillMatrixSelf(Eigen::MatrixXd &matr, Eigen::VectorXd &lastLine, Eigen::VectorXd &lactCol) override
Генерация блока матрицы
Заголовочный файл с описанием класса Wake.
virtual void GetInfluenceFromVorticesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const override
Вычисление влияния части подряд идущих вихрей из вихревого следа на прямолинейную панель для правой ч...
double delta
Расстояние, на которое рождаемый вихрь отодвигается от профиля
Definition: Passport2D.h:151
Заголовочный файл с описанием класса World2D.
virtual void GetInfluenceFromVInfToRectPanel(std::vector< double > &vInfRhs) const override
Вычисление влияния набегающего потока на прямолинейную панель для правой части
Заголовочный файл с описанием класса Airfoil.
double maxGamma
Максимально допустимая циркуляция вихря
Definition: Passport2D.h:157
const double IDPI
Число .
Definition: defs.h:76
Абстрактный класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле ...
Definition: Boundary2D.h:63
std::vector< Point2D > vecHalfGamma
Скорость вихрей виртуального следа конкретного профиля (равна Gamma/2) используется для расчета давле...
Definition: VirtualWake2D.h:70
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
Заголовочный файл с описанием класса Mechanics.
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
const Point2D & getV(size_t q) const
Возврат константной ссылки на скорость вершины профиля
Definition: Airfoil2D.h:113
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
HD double & g()
Функция для доступа к циркуляции вихря
Definition: Vortex2D.h:93
double Alpha(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления угла между векторами
Definition: defs.cpp:259
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition: numvector.h:383
Заголовочный файл с описанием класса BoundaryVortexCollocN.
BoundaryVortexCollocN(const World2D &W_, size_t numberInPassport_)
Конструктор
virtual void FillMatrixFromOther(const Boundary &otherBoundary, Eigen::MatrixXd &matr) override
Генерация блока матрицы влияния от другого профиля того же типа
std::vector< std::pair< size_t, size_t > > aflPan
Пара чисел: номер профиля и номер панели, на которой рожден виртуальный вихрь
Definition: VirtualWake2D.h:73
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition: defs.cpp:265
const size_t numberInPassport
Номер профиля в паспорте
Definition: Boundary2D.h:70
virtual void GetInfluenceFromVortexSheetAtRectPanelToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const override
Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вих...
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
virtual void CalcConvVelocityAtVirtualVortexes(std::vector< Point2D > &velo) const override
Вычисление конвективной скорости в точках расположения виртуальных вихрей
virtual void SolutionToFreeVortexSheetAndVirtualVortex(const Eigen::VectorXd &sol) override
Пересчет решения на интенсивность вихревого слоя и на рождаемые вихри на конкретном профиле ...
Definition: Airfoil2D.h:45
Заголовочный файл с описанием класса StreamParser.
std::vector< Vortex2D > vtx
Список вихревых элементов
size_t getSheetSize() const
Definition: Sheet2D.h:94
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
auto unit(P newlen=1) const -> numvector< typename std::remove_const< decltype(this->data[0]*newlen)>::type, n >
Вычисление орта вектора или вектора заданной длины, коллинеарного данному
Definition: numvector.h:399
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
Заголовочный файл с описанием класса MeasureVP.
const double & attachedSourceSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:109
Sheet oldSheets
Слои на профиле с предыдущего шага
Definition: Boundary2D.h:98
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
virtual void ComputeAttachedSheetsIntensity() override
Вычисление интенсивностей присоединенного вихревого слоя и присоединенного слоя источников ...
const Gpu & getCuda() const
Возврат константной ссылки на объект, связанный с видеокартой (GPU)
Definition: World2D.h:227
std::vector< Point2D > c
Контрольные точки - центры панелей
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
int minVortexPerPanel
Минимальное число вихрей, рождаемых на каждой панели профииля
Definition: Passport2D.h:154
virtual void GetInfluenceFromSourcesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const override
Вычисление влияния части подряд источников из области течения на прямолинейную панель для правой част...
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
virtual void FillIQSelf(std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &IQ) override
Генерация блока матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2, влияние профиля самого на себя ...
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
Point2D V0() const
Функция скорости набегающего потока с учетом разгона
Definition: Passport2D.h:93
Класс, опеделяющий текущую решаемую задачу
Definition: World2D.h:68
Класс, опеделяющий набор вихрей
numvector< T, n > & toZero(P val=0)
Установка всех компонент вектора в константу (по умолчанию — нуль)
Definition: numvector.h:524
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
Заголовочный файл с описанием класса Velocity.
std::vector< std::pair< int, int > > vortexBeginEnd
Номера первого и последнего вихрей, рождаемых на каждой панели профиля (формируется после решения СЛА...
Definition: Boundary2D.h:82
const World2D & W
Константная ссылка на решаемую задачу
Definition: Boundary2D.h:67
virtual void GetInfluenceFromSourceSheetAtRectPanelToVortex(size_t panel, const Vortex2D &ptr, Point2D &vel) const override
Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения ...
double eps2
Квадрат радиуса вихря
Definition: Passport2D.h:142
const Airfoil & afl
Definition: Boundary2D.h:76
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition: numvector.h:507