VM2D  1.12
Vortex methods for 2D flows simulation
VM2D::BoundaryLinLayerAver Class Reference

Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле More...

#include <Boundary2DLinLayerAver.h>

Inheritance diagram for VM2D::BoundaryLinLayerAver:
Collaboration diagram for VM2D::BoundaryLinLayerAver:

Public Member Functions

 BoundaryLinLayerAver (const World2D &W_, size_t numberInPassport_)
 Конструктор More...
 
virtual ~BoundaryLinLayerAver ()
 Деструктор More...
 
virtual void FillMatrixSelf (Eigen::MatrixXd &matr, Eigen::VectorXd &lastLine, Eigen::VectorXd &lactCol) override
 Генерация блока матрицы More...
 
virtual void FillIQSelf (std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &IQ) override
 Генерация блока матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2, влияние профиля самого на себя More...
 
virtual void FillMatrixFromOther (const Boundary &otherBoundary, Eigen::MatrixXd &matr) override
 Генерация блока матрицы влияния от другого профиля того же типа More...
 
virtual void FillIQFromOther (const Boundary &otherBoundary, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &IQ) override
 Генерация блока матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2, влияние одного профиля на другой More...
 
virtual void SolutionToFreeVortexSheetAndVirtualVortex (const Eigen::VectorXd &sol) override
 Пересчет решения на интенсивность вихревого слоя и на рождаемые вихри на конкретном профиле More...
 
virtual void CalcConvVelocityToSetOfPointsFromSheets (const WakeDataBase &pointsDb, std::vector< Point2D > &velo) const override
 Вычисление конвективных скоростей в наборе точек, вызываемых наличием слоев вихрей и источников на профиле More...
 
virtual void CalcConvVelocityAtVirtualVortexes (std::vector< Point2D > &velo) const override
 Вычисление конвективной скорости в точках расположения виртуальных вихрей More...
 
virtual void ComputeAttachedSheetsIntensity () override
 Вычисление интенсивностей присоединенного вихревого слоя и присоединенного слоя источников More...
 
virtual void GetInfluenceFromVorticesToRectPanel (size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const override
 Вычисление влияния части подряд идущих вихрей из вихревого следа на прямолинейную панель для правой части More...
 
virtual void GetInfluenceFromSourcesToRectPanel (size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const override
 Вычисление влияния части подряд источников из области течения на прямолинейную панель для правой части More...
 
virtual void GetInfluenceFromVortexSheetAtRectPanelToVortex (size_t panel, const Vortex2D &vtx, Point2D &vel) const override
 Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вихрь в области течения ///. More...
 
virtual void GetInfluenceFromSourceSheetAtRectPanelToVortex (size_t panel, const Vortex2D &vtx, Point2D &vel) const override
 Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения More...
 
virtual void GetInfluenceFromVInfToRectPanel (std::vector< double > &vInfRhs) const override
 Вычисление влияния набегающего потока на прямолинейную панель для правой части More...
 
size_t GetUnknownsSize () const
 Возврат размерности вектора решения More...
 

Public Attributes

const Airfoilafl
 
int minVortexPerPanel
 Минимальное число вихрей, рождаемых на панели профиля и формирующих виртуальный вихревой след More...
 
std::vector< std::pair< int, int > > vortexBeginEnd
 Номера первого и последнего вихрей, рождаемых на каждой панели профиля (формируется после решения СЛАУ) More...
 
VirtualWake virtualWake
 Виртуальный вихревой след конкретного профиля More...
 
size_t sheetDim
 Размерность параметров каждого из слоев на каждой из панелей More...
 
Sheet sheets
 Слои на профиле More...
 
Sheet oldSheets
 Слои на профиле с предыдущего шага More...
 

Protected Attributes

const World2DW
 Константная ссылка на решаемую задачу More...
 
const size_t numberInPassport
 Номер профиля в паспорте More...
 

Detailed Description

Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле

Способ удовлетворения граничного условия:

  • вихревой слой представляется кусочно-линейном виде с разрывами в каждом узле профиля;
  • условие ортогональности невязки граничного условия константным и линейным проекционным функциям.
Author
Марчевский Илья Константинович
Сокол Ксения Сергеевна
Рятина Евгения Павловна
Колганова Александра Олеговна
\Version 1.12
\date 14 января 2024 г.

Definition at line 65 of file Boundary2DLinLayerAver.h.

Constructor & Destructor Documentation

VM2D::BoundaryLinLayerAver::BoundaryLinLayerAver ( const World2D W_,
size_t  numberInPassport_ 
)
inline

Конструктор

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

Definition at line 74 of file Boundary2DLinLayerAver.h.

74  :
75  Boundary(W_, numberInPassport_, 2)
76  { };
Boundary(const World2D &W_, size_t numberInPassport_, int sheetDim_)
Конструктор
Definition: Boundary2D.cpp:53
virtual VM2D::BoundaryLinLayerAver::~BoundaryLinLayerAver ( )
inlinevirtual

Деструктор

Definition at line 79 of file Boundary2DLinLayerAver.h.

79 {};

Here is the call graph for this function:

Member Function Documentation

void BoundaryLinLayerAver::CalcConvVelocityAtVirtualVortexes ( std::vector< Point2D > &  velo) const
overridevirtual

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

Вычисление производится в точках расположения виртуальных вихрей. Скорости находятся как скорости соответствующих точек профиля плюс половина интенсивности присоединенного вихревого слоя, умноженная на касательную к профилю на данной панели.

Parameters
[out]veloссылка на заполняемый список скоростей
Warning
Массив velo заполняется путем присвоения, а не прибавления значений

Implements VM2D::Boundary.

Definition at line 355 of file Boundary2DLinLayerAver.cpp.

356 {
357  std::vector<Point2D>& Vel = velo;
358 
359  //Скорости виртуальных вихрей
360 
361 #pragma omp parallel for default(none) shared(Vel)
362  for (int i = 0; i < (int)Vel.size(); ++i)
363  Vel[i] = afl.getV(virtualWake.aflPan[i].second) \
366 }//CalcConvVelocityAtVirtualVortexes(...)
std::vector< Point2D > vecHalfGamma
Скорость вихрей виртуального следа конкретного профиля (равна Gamma/2) используется для расчета давле...
Definition: VirtualWake2D.h:70
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
const Point2D & getV(size_t q) const
Возврат константной ссылки на скорость вершины профиля
Definition: Airfoil2D.h:113
std::vector< std::pair< size_t, size_t > > aflPan
Пара чисел: номер профиля и номер панели, на которой рожден виртуальный вихрь
Definition: VirtualWake2D.h:73
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
Point2D V0() const
Функция скорости набегающего потока с учетом разгона
Definition: Passport2D.h:93
const World2D & W
Константная ссылка на решаемую задачу
Definition: Boundary2D.h:67
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::CalcConvVelocityToSetOfPointsFromSheets ( const WakeDataBase pointsDb,
std::vector< Point2D > &  velo 
) const
overridevirtual

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

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

Parameters
[in]pointsDbконстантная ссылка на базу данных вихрей, в точках которых вычисляются скорости
[out]veloссылка на вектор скоростей, которые приобретают точки из-за влияния слоев вихрей и источников на профиле
Warning
velo — накапливается!
Todo:
Тут надо разобраться, как должно быть...
Todo:
сделать if(move || deform)
Todo:
почему не sheets.freeVortexSheet(j, 0)?

Implements VM2D::Boundary.

Definition at line 268 of file Boundary2DLinLayerAver.cpp.

269 {
270  std::vector<Point2D> selfVelo(pointsDb.vtx.size());
271 
272 #pragma warning (push)
273 #pragma warning (disable: 4101)
274  //Локальные переменные для цикла
275  Point2D velI;
276  Point2D tempVel;
277 #pragma warning (pop)
278 
279 #pragma omp parallel for default(none) shared(selfVelo, pointsDb, std::cout) private(velI, tempVel)
280  for (int i = 0; i < pointsDb.vtx.size(); ++i)
281  {
282  velI.toZero();
283 
284  const Point2D& posI = pointsDb.vtx[i].r();
285 
288  for (size_t j = 0; j < afl.getNumberOfPanels(); ++j)
289  {
290  Point2D dj = afl.getR(j + 1) - afl.getR(j);
291  Point2D tauj = dj.unit();
292 
293  Point2D s = posI - afl.getR(j);
294  Point2D p = posI - afl.getR(j + 1);
295 
296  double a = VMlib::Alpha(p, s);
297 
298  double lambda = VMlib::Lambda(p, s);
299 
300  Point2D u1 = 0.5 / dj.length() * VMlib::Omega(p + s, tauj, tauj) ;
301 
302  Point2D skos0 = -a * tauj.kcross() + lambda * tauj;
303  Point2D skos1 = -a * u1.kcross() + lambda * u1 - tauj;
304 
306  velI += sheets.freeVortexSheet(j, 0) * skos0.kcross() + sheets.freeVortexSheet(j, 1) * skos1.kcross();
307  velI += sheets.attachedVortexSheet(j, 0) * skos0.kcross() + sheets.attachedVortexSheet(j, 1) * skos1.kcross();
308  velI += sheets.attachedSourceSheet(j, 0) * skos0 + sheets.attachedSourceSheet(j, 1) * skos1;
309  }//for j
310 
311  velI *= IDPI;
312  selfVelo[i] = velI;
313  }//for i
314 
315  for (size_t i = 0; i < velo.size(); ++i)
316  velo[i] += selfVelo[i];
317 }//CalcConvVelocityToSetOfPointsFromSheets(...)
const double & freeVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:99
const double & attachedVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:104
const double IDPI
Число .
Definition: defs.h:76
P length() const
Вычисление 2-нормы (длины) вектора
Definition: numvector.h:371
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
double Alpha(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления угла между векторами
Definition: defs.cpp:259
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
Point2D Omega(const Point2D &a, const Point2D &b, const Point2D &c)
Вспомогательная функция вычисления величины .
Definition: defs.cpp:271
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition: defs.cpp:265
std::vector< Vortex2D > vtx
Список вихревых элементов
Класс, опеделяющий двумерный вектор
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
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
const double & attachedSourceSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:109
const Airfoil & afl
Definition: Boundary2D.h:76
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition: numvector.h:507

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::ComputeAttachedSheetsIntensity ( )
overridevirtual

Вычисление интенсивностей присоединенного вихревого слоя и присоединенного слоя источников

Implements VM2D::Boundary.

Definition at line 370 of file Boundary2DLinLayerAver.cpp.

371 {
372  for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
373  {
378  }
379 
380  const Airfoil* oldAfl = (W.getCurrentStep() == 0) ? &afl : &W.getOldAirfoil(numberInPassport);
381 
382  for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
383  {
384  sheets.attachedVortexSheet(i, 0) = 0.5 * (afl.getV(i + 1) + afl.getV(i)) & afl.tau[i];
386 
387  sheets.attachedSourceSheet(i, 0) = 0.5 * (afl.getV(i + 1) + afl.getV(i)) & afl.nrm[i];
389 
390 
391  }
392 }//ComputeAttachedSheetsIntensity()
const double & attachedVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:104
const Point2D & getV(size_t q) const
Возврат константной ссылки на скорость вершины профиля
Definition: Airfoil2D.h:113
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
double dt
Шаг по времени
Definition: PassportGen.h:65
const Airfoil & getOldAirfoil(size_t i) const
Возврат константной ссылки на объект старого профиля
Definition: World2D.h:136
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
const size_t numberInPassport
Номер профиля в паспорте
Definition: Boundary2D.h:70
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
const double & attachedSourceSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:109
Sheet oldSheets
Слои на профиле с предыдущего шага
Definition: Boundary2D.h:98
double phiAfl
Поворот профиля
Definition: Airfoil2D.h:80
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177
Абстрактный класс, определяющий обтекаемый профиль
Definition: Airfoil2D.h:60
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition: WorldGen.h:91
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
const World2D & W
Константная ссылка на решаемую задачу
Definition: Boundary2D.h:67
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::FillIQFromOther ( const Boundary otherBoundary,
std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &  IQ 
)
overridevirtual

Генерация блока матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2, влияние одного профиля на другой

Генерирует блок матрицы влияния от другого профиля того же типа

Parameters
[in]otherBoundaryконстантная ссылка на граничное условие на втором профиле
[out]IQссылка на пару матриц, выражающих взаимные влияния (касательные и нормальные) панелей профиля
Todo:
Пока считается, что граничные условия одинаковые

Implements VM2D::Boundary.

Definition at line 261 of file Boundary2DLinLayerAver.cpp.

262 {
263  afl.calcIQ(2, otherBoundary.afl, IQ);
264 }//FillIQFromOther(...)
virtual void calcIQ(size_t p, const Airfoil &otherAirfoil, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &matrPair) const =0
Вычисление коэффициентов матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2.
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::FillIQSelf ( std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &  IQ)
overridevirtual

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

Parameters
[out]IQссылка на генерируемую матрицу

Implements VM2D::Boundary.

Definition at line 235 of file Boundary2DLinLayerAver.cpp.

236 {
237  afl.calcIQ(2, afl, IQ);
238 }//FillIQSelf(...)
virtual void calcIQ(size_t p, const Airfoil &otherAirfoil, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &matrPair) const =0
Вычисление коэффициентов матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2.
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::FillMatrixFromOther ( const Boundary otherBoundary,
Eigen::MatrixXd &  matr 
)
overridevirtual

Генерация блока матрицы влияния от другого профиля того же типа

Генерирует блок матрицы влияния от другого профиля того же типа

Parameters
[in]otherBoundaryконстантная ссылка на граничное условие на втором профиле
[out]matrссылка на генерируемый блок матрицы
Todo:
Пока считается, что граничные условия одинаковые

Implements VM2D::Boundary.

Definition at line 242 of file Boundary2DLinLayerAver.cpp.

243 {
244  size_t np = afl.getNumberOfPanels();
245  size_t npOther = otherBoundary.afl.getNumberOfPanels();
246 
247  std::vector<double> res(4, 0.0);
248 
249  for (size_t i = 0; i < np; ++i)
250  for (size_t j = 0; j < npOther; ++j)
251  {
252  res = afl.getA(2, i, otherBoundary.afl, j);
253 
254  matr(i, j) = res[0];
255  matr(i, npOther + j) = res[1];
256  matr(np + i, j) = res[2];
257  matr(np + i, npOther + j) = res[3];
258  }
259 }//FillMatrixFromOther(...)
virtual std::vector< double > getA(size_t p, size_t i, const Airfoil &airfoil, size_t j) const =0
Вычисление коэффициентов матрицы A для расчета влияния панели на панель
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::FillMatrixSelf ( Eigen::MatrixXd &  matr,
Eigen::VectorXd &  lastLine,
Eigen::VectorXd &  lactCol 
)
overridevirtual

Генерация блока матрицы

Генерирует следующие компоненты матрицы:

  • диагональный блок матрицы — влияние данного профиля на самого себя;
  • нижнюю строку для матрицы для данного профиля;
  • правый столбец матрицы для данного профиля.
Parameters
[out]matrссылка на генерируемую матрицу
[out]lastLineссылка на нижнюю строку
[out]lactColссылка на правый столбец

Implements VM2D::Boundary.

Definition at line 209 of file Boundary2DLinLayerAver.cpp.

210 {
211  size_t np = afl.getNumberOfPanels();
212 
213  std::vector<double> res(4, 0.0);
214 
215  for (size_t i = 0; i < np; ++i)
216  {
217  lastCol(i) = 1.0;
218  lastLine(i) = afl.len[i];
219  lastCol(np + i) = 0.0;
220  lastLine(np + i) = 0.0;
221  }
222 
223  for (size_t i = 0; i < np; ++i)
224  for (size_t j = 0; j < np; ++j)
225  {
226  res = afl.getA(2, i, afl, j);
227  matr(i, j) = res[0];
228  matr(i, np + j) = res[1];
229  matr(np + i, j) = res[2];
230  matr(np + i, np + j) = res[3];
231  }
232 }//FillMatrixSelf(...)
virtual std::vector< double > getA(size_t p, size_t i, const Airfoil &airfoil, size_t j) const =0
Вычисление коэффициентов матрицы A для расчета влияния панели на панель
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::GetInfluenceFromSourceSheetAtRectPanelToVortex ( size_t  panel,
const Vortex2D vtx,
Point2D vel 
) const
overridevirtual

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

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

Implements VM2D::Boundary.

Definition at line 459 of file Boundary2DLinLayerAver.cpp.

460 {
461  vel.toZero();
462 
463  const Point2D& posI = ptr.r();
464 
465  Point2D dj = afl.getR(panel + 1) - afl.getR(panel);
466  Point2D tauj = dj.unit();
467 
468  Point2D s = posI - afl.getR(panel);
469  Point2D p = posI - afl.getR(panel + 1);
470 
471  Point2D u1 = 0.5 / dj.length()* VMlib::Omega(p + s, tauj, tauj);
472 
473  double a = VMlib::Alpha(p, s);
474 
475  double lambda;
476  if ((s.length2() > 1e-16) && (p.length2() > 1e-16))
477  lambda = VMlib::Lambda(p, s);
478  else
479  lambda = 0.0;
480 
481  vel += sheets.attachedSourceSheet(panel, 0) * (-a * u1.kcross() + lambda * u1 - tauj);
482  vel *= IDPI;
483 }// GetInfluenceFromSourceSheetAtRectPanelToVortex(...)
const double IDPI
Число .
Definition: defs.h:76
P length() const
Вычисление 2-нормы (длины) вектора
Definition: numvector.h:371
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
double Alpha(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления угла между векторами
Definition: defs.cpp:259
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition: numvector.h:383
Point2D Omega(const Point2D &a, const Point2D &b, const Point2D &c)
Вспомогательная функция вычисления величины .
Definition: defs.cpp:271
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition: defs.cpp:265
Класс, опеделяющий двумерный вектор
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
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
const double & attachedSourceSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:109
numvector< T, n > & toZero(P val=0)
Установка всех компонент вектора в константу (по умолчанию — нуль)
Definition: numvector.h:524
const Airfoil & afl
Definition: Boundary2D.h:76
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition: numvector.h:507

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::GetInfluenceFromSourcesToRectPanel ( size_t  panel,
const Vortex2D ptr,
ptrdiff_t  count,
std::vector< double > &  wakeRhs 
) const
overridevirtual

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

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

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

Implements VM2D::Boundary.

Definition at line 427 of file Boundary2DLinLayerAver.cpp.

428 {
429  double& velI = wakeRhs[0];
430  double& velILin = wakeRhs[1];
431 
432  const Point2D& posI0 = afl.getR(panel);
433  const Point2D& posI1 = afl.getR(panel + 1);
434  Point2D di = posI1 - posI0;
435  const Point2D& taui = afl.tau[panel];
436 
437  for (size_t it = 0; it != count; ++it)
438  {
439  const Vortex2D& vt = ptr[it];
440  const Point2D& posJ = vt.r();
441  const double& gamJ = vt.g();
442 
443  Point2D s = posJ - posI0;
444  Point2D p = posJ - posI1;
445 
446  Point2D u1 = 0.5 / di.length() * VMlib::Omega(p + s, taui, taui);
447 
448  double alpha = VMlib::Alpha(p, s);
449  double lambda = VMlib::Lambda(p, s);
450 
451  velI -= gamJ * lambda;
452 
453  velILin -= gamJ * (alpha * (u1 & taui.kcross()) + lambda * (u1 & taui) - 1.0);
454  }
455 
456 }//GetInfluenceFromSourcesToRectPanel(...)
P length() const
Вычисление 2-нормы (длины) вектора
Definition: numvector.h:371
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
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
Point2D Omega(const Point2D &a, const Point2D &b, const Point2D &c)
Вспомогательная функция вычисления величины .
Definition: defs.cpp:271
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition: defs.cpp:265
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
const Airfoil & afl
Definition: Boundary2D.h:76
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition: numvector.h:507

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::GetInfluenceFromVInfToRectPanel ( std::vector< double > &  vInfRhs) const
overridevirtual

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

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

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

Implements VM2D::Boundary.

Definition at line 512 of file Boundary2DLinLayerAver.cpp.

513 {
514  size_t np = afl.getNumberOfPanels();
515 
516  vInfRhs.resize(2 * np);
517 
518 #pragma omp parallel for default(none) shared(vInfRhs, np)
519  for (int i = 0; i < np; ++i)
520  {
521  vInfRhs[i] = afl.tau[i] & W.getPassport().physicalProperties.V0();
522  vInfRhs[np + i] = 0;
523  }
524 }
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Point2D V0() const
Функция скорости набегающего потока с учетом разгона
Definition: Passport2D.h:93
const World2D & W
Константная ссылка на решаемую задачу
Definition: Boundary2D.h:67
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::GetInfluenceFromVortexSheetAtRectPanelToVortex ( size_t  panel,
const Vortex2D vtx,
Point2D vel 
) const
overridevirtual

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

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

Implements VM2D::Boundary.

Definition at line 486 of file Boundary2DLinLayerAver.cpp.

487 {
488  vel.toZero();
489 
490  const Point2D& posI = ptr.r();
491 
492  Point2D dj = afl.getR(panel + 1) - afl.getR(panel);
493  Point2D tauj = dj.unit();
494 
495  Point2D s = posI - afl.getR(panel);
496  Point2D p = posI - afl.getR(panel + 1);
497 
498  Point2D u1 = 0.5 / dj.length()* VMlib::Omega(p + s, tauj, tauj);
499 
500  double lambda;
501  if ((s.length2() > 1e-16) && (p.length2() > 1e-16))
502  lambda = VMlib::Lambda(p, s);
503  else
504  lambda = 0.0;
505 
506  vel += sheets.freeVortexSheet(panel, 0) * (lambda * u1.kcross() - tauj.kcross());
507  vel += sheets.attachedVortexSheet(panel, 0) * (lambda * u1.kcross() - tauj.kcross());
508  vel *= IDPI;
509 }// GetInfluenceFromVortexSheetAtRectPanelToVortex(...)
const double & freeVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:99
const double & attachedVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:104
const double IDPI
Число .
Definition: defs.h:76
P length() const
Вычисление 2-нормы (длины) вектора
Definition: numvector.h:371
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition: numvector.h:383
Point2D Omega(const Point2D &a, const Point2D &b, const Point2D &c)
Вспомогательная функция вычисления величины .
Definition: defs.cpp:271
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition: defs.cpp:265
Класс, опеделяющий двумерный вектор
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
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
numvector< T, n > & toZero(P val=0)
Установка всех компонент вектора в константу (по умолчанию — нуль)
Definition: numvector.h:524
const Airfoil & afl
Definition: Boundary2D.h:76
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition: numvector.h:507

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::GetInfluenceFromVorticesToRectPanel ( size_t  panel,
const Vortex2D ptr,
ptrdiff_t  count,
std::vector< double > &  wakeRhs 
) const
overridevirtual

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

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

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

Implements VM2D::Boundary.

Definition at line 395 of file Boundary2DLinLayerAver.cpp.

396 {
397  double& velI = wakeRhs[0];
398  double& velILin = wakeRhs[1];
399 
400  const Point2D& posI0 = afl.getR(panel);
401  const Point2D& posI1 = afl.getR(panel + 1);
402  Point2D di = posI1 - posI0;
403  const Point2D& taui = afl.tau[panel];
404 
405  for (size_t it = 0; it != count; ++it)
406  {
407  const Vortex2D& vt = ptr[it];
408  const Point2D& posJ = vt.r();
409  const double& gamJ = vt.g();
410 
411  Point2D s = posJ - posI0;
412  Point2D p = posJ - posI1;
413 
414  Point2D u1 = 0.5 / di.length() * VMlib::Omega(p + s, taui, taui);
415 
416  double alpha = VMlib::Alpha(p, s);
417  double lambda = VMlib::Lambda(p, s);
418 
419  velI -= gamJ * alpha;
420 
421  velILin -= gamJ * (alpha * (u1 & taui) + lambda * (u1 & (-taui.kcross())));
422  }
423 }//GetInfluenceFromVorticesToRectPanel(...)
P length() const
Вычисление 2-нормы (длины) вектора
Definition: numvector.h:371
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
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
Point2D Omega(const Point2D &a, const Point2D &b, const Point2D &c)
Вспомогательная функция вычисления величины .
Definition: defs.cpp:271
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition: defs.cpp:265
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
const Airfoil & afl
Definition: Boundary2D.h:76
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition: numvector.h:507

Here is the call graph for this function:

Here is the caller graph for this function:

size_t Boundary::GetUnknownsSize ( ) const
inherited

Возврат размерности вектора решения

(без учета регуляризирующей переменной)

Returns
размерность вектора решения

Definition at line 71 of file Boundary2D.cpp.

72 {
73  return sheetDim * afl.getNumberOfPanels();
74 }
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
size_t sheetDim
Размерность параметров каждого из слоев на каждой из панелей
Definition: Boundary2D.h:92
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void BoundaryLinLayerAver::SolutionToFreeVortexSheetAndVirtualVortex ( const Eigen::VectorXd &  sol)
overridevirtual

Пересчет решения на интенсивность вихревого слоя и на рождаемые вихри на конкретном профиле

1) Приводит решение к интенсивности вихревого слоя и записывает его в sheets.freeVortexSheet:

  • если неизвестное — интенсивность вихря, то он "размазывается" по панели;
  • если неизвестное — интенсивность слоя, то она передается непосредственно.

2) Приводит интенсивность вихревого слоя к рождаемым вихрям, а также вычисляет их положения

Parameters
[in]solвектор решения СЛАУ

Implements VM2D::Boundary.

Definition at line 57 of file Boundary2DLinLayerAver.cpp.

58 {
59  Vortex2D virtVort;
60  Point2D midNorm;
61 
63 
65 
66  //Очистка и резервирование памяти
67  virtualWake.vecHalfGamma.clear();
68  virtualWake.vecHalfGamma.reserve(afl.getNumberOfPanels() * nVortPerPan);
69 
70  //Очистка и резервирование памяти
71  virtualWake.aflPan.clear();
72  virtualWake.aflPan.reserve(afl.getNumberOfPanels() * nVortPerPan);
73 
74  //Очистка и резервирование памяти
75  virtualWake.vtx.clear();
76  virtualWake.vtx.reserve(afl.getNumberOfPanels() * nVortPerPan);
77 
78  //Очистка и резервирование памяти
79  vortexBeginEnd.clear();
81 
83 
84  std::pair<int, int> pair;
85 
86  for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
87  {
88  midNorm = afl.nrm[i] * delta;
89 
90  //число участков, на которые разбивается панель для сброса вихрей
91  //определяется через наибольшее значение решения на профиле, т.е. в крайней точке
92 
93 
94  pair.first = (int)virtualWake.vtx.size();
95 
96  double a = sol(i) - 0.5 * sol(afl.getNumberOfPanels() + i);
97  double b = sol(i) + 0.5 * sol(afl.getNumberOfPanels() + i);
98 
99  if (fabs(a - b) < 1e-10)
100  {
101  size_t NEWnVortPerPan = (size_t)std::max((int)std::ceil(fabs(sol(i)*afl.len[i]) / maxG), nVortPerPan);
102  Point2D dr = 1.0 / NEWnVortPerPan * (afl.getR(i + 1) - afl.getR(i));
103 
104  for (size_t j = 0; j < NEWnVortPerPan; ++j)
105  {
106  virtVort.r() = afl.getR(i) + dr * (j * 1.0 + 0.5) + midNorm;
107  virtVort.g() = sol(i) * afl.len[i] / NEWnVortPerPan;
108  virtualWake.vtx.push_back(virtVort);
109 
110  virtualWake.vecHalfGamma.push_back(0.5 * sol(i) * afl.tau[i]);
111  virtualWake.aflPan.push_back({ numberInPassport, i });
112  }
113  }
114  else if (a * b >= 0.0)
115  {
116  size_t NEWnVortPerPan = (size_t)std::max((int)std::ceil(fabs(0.5*(a+b)*afl.len[i]) / maxG), nVortPerPan);
117  std::vector<double> ds(NEWnVortPerPan+1, 0.0);
118 
119  for (size_t k = 1; k <= NEWnVortPerPan; ++k)
120  {
121  if ((a > 0) || (b > 0))
122  ds[k] = afl.len[i] / (a - b) * (a - sqrt((b*b * k + a * a*(NEWnVortPerPan - k)) / NEWnVortPerPan));
123  else
124  ds[k] = afl.len[i] / (a - b) * (a + sqrt((b*b * k + a * a*(NEWnVortPerPan - k)) / NEWnVortPerPan));
125  }
126 
127  for (size_t j = 0; j < NEWnVortPerPan; ++j)
128  {
129  virtVort.r() = afl.getR(i) + 0.5*(ds[j]+ds[j+1])*afl.tau[i] + midNorm;
130  virtVort.g() = 0.5 * (a+b) * afl.len[i] / NEWnVortPerPan;
131  virtualWake.vtx.push_back(virtVort);
132 
133  virtualWake.vecHalfGamma.push_back(0.5 * (a + 0.5*(ds[j] + ds[j + 1]) * (b-a)/ afl.len[i]) * afl.tau[i]);
134  virtualWake.aflPan.push_back({ numberInPassport, i });
135  }
136  }
137  else
138  {
139  double sast = -a * afl.len[i] / (b - a);
140 
141  //from 0 to sast
142  {
143  size_t NEWnVortPerPan = (size_t)std::max((int)std::ceil(fabs(0.5*(a + 0.0)*sast) / maxG), (int)std::ceil(nVortPerPan * sast / afl.len[i]));
144  std::vector<double> ds(NEWnVortPerPan + 1, 0.0);
145 
146  for (size_t k = 1; k <= NEWnVortPerPan; ++k)
147  {
148  if (a > 0)
149  ds[k] = sast / (a - 0.0) * (a - sqrt((a * a*(NEWnVortPerPan - k)) / NEWnVortPerPan));
150  else
151  ds[k] = sast / (a - 0.0) * (a + sqrt((a * a*(NEWnVortPerPan - k)) / NEWnVortPerPan));
152  }
153 
154  for (size_t j = 0; j < NEWnVortPerPan; ++j)
155  {
156  virtVort.r() = afl.getR(i) + 0.5*(ds[j] + ds[j + 1])*afl.tau[i] + midNorm;
157  virtVort.g() = 0.5 * (a + 0.0) * sast / NEWnVortPerPan;
158  virtualWake.vtx.push_back(virtVort);
159 
160  virtualWake.vecHalfGamma.push_back(0.5 * (a + 0.5*(ds[j] + ds[j + 1]) * (0.0 - a) / sast) * afl.tau[i]);
161  virtualWake.aflPan.push_back({ numberInPassport, i });
162  }
163  }
164 
165  double sastast = afl.len[i] - sast;
166 
167  //from sats to len
168  {
169  size_t NEWnVortPerPan = (size_t)std::max((int)std::ceil(fabs(0.5*(0.0 + b)*sastast) / maxG), (int)std::ceil(nVortPerPan * sastast / afl.len[i]));
170  std::vector<double> ds(NEWnVortPerPan + 1, 0.0);
171 
172  for (size_t k = 1; k <= NEWnVortPerPan; ++k)
173  {
174  if (b > 0)
175  ds[k] = sastast / (0.0 - b) * (0.0 - sqrt((b*b * k) / NEWnVortPerPan));
176  else
177  ds[k] = sastast / (0.0 - b) * (0.0 + sqrt((b*b * k) / NEWnVortPerPan));
178  }
179 
180  for (size_t j = 0; j < NEWnVortPerPan; ++j)
181  {
182  virtVort.r() = afl.getR(i) + sast * afl.tau[i] + 0.5*(ds[j] + ds[j + 1])*afl.tau[i] + midNorm;
183  virtVort.g() = 0.5 * (0.0 + b) * sastast / NEWnVortPerPan;
184  virtualWake.vtx.push_back(virtVort);
185 
186  virtualWake.vecHalfGamma.push_back(0.5 * (0.0 + 0.5*(ds[j] + ds[j + 1]) * (b - 0.0) / sastast) * afl.tau[i]);
187  virtualWake.aflPan.push_back({ numberInPassport, i });
188  }
189 
190  }
191  }
192 
193  pair.second = (int)virtualWake.vtx.size();
194  vortexBeginEnd.push_back(pair);
195  }
196 
197 
198  for (size_t j = 0; j < afl.getNumberOfPanels(); ++j)
199  {
200  sheets.freeVortexSheet(j, 0) = sol(j);
201  sheets.freeVortexSheet(j, 1) = sol(afl.getNumberOfPanels() + j);
202  }
203 
204 }//SolutionToFreeVortexSheetAndVirtualVortex(...)
const double & freeVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:99
double delta
Расстояние, на которое рождаемый вихрь отодвигается от профиля
Definition: Passport2D.h:151
double maxGamma
Максимально допустимая циркуляция вихря
Definition: Passport2D.h:157
std::vector< Point2D > vecHalfGamma
Скорость вихрей виртуального следа конкретного профиля (равна Gamma/2) используется для расчета давле...
Definition: VirtualWake2D.h:70
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
HD double & g()
Функция для доступа к циркуляции вихря
Definition: Vortex2D.h:93
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::vector< std::pair< size_t, size_t > > aflPan
Пара чисел: номер профиля и номер панели, на которой рожден виртуальный вихрь
Definition: VirtualWake2D.h:73
const size_t numberInPassport
Номер профиля в паспорте
Definition: Boundary2D.h:70
std::vector< Vortex2D > vtx
Список вихревых элементов
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
int minVortexPerPanel
Минимальное число вихрей, рождаемых на каждой панели профииля
Definition: Passport2D.h:154
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
std::vector< std::pair< int, int > > vortexBeginEnd
Номера первого и последнего вихрей, рождаемых на каждой панели профиля (формируется после решения СЛА...
Definition: Boundary2D.h:82
const World2D & W
Константная ссылка на решаемую задачу
Definition: Boundary2D.h:67
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

const Airfoil& VM2D::Boundary::afl
inherited

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

Definition at line 76 of file Boundary2D.h.

int VM2D::Boundary::minVortexPerPanel
inherited

Минимальное число вихрей, рождаемых на панели профиля и формирующих виртуальный вихревой след

Definition at line 79 of file Boundary2D.h.

const size_t VM2D::Boundary::numberInPassport
protectedinherited

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

Definition at line 70 of file Boundary2D.h.

Sheet VM2D::Boundary::oldSheets
inherited

Слои на профиле с предыдущего шага

Definition at line 98 of file Boundary2D.h.

size_t VM2D::Boundary::sheetDim
inherited

Размерность параметров каждого из слоев на каждой из панелей

Указывает, сколькими числами задается интенсивность каждого из слоев на каждой панели:

  • 1 — одно число — задается только среднее значение;
  • 2 — два числа — задается среднее значение и "наклон".

Definition at line 92 of file Boundary2D.h.

Sheet VM2D::Boundary::sheets
inherited

Слои на профиле

Definition at line 95 of file Boundary2D.h.

VirtualWake VM2D::Boundary::virtualWake
inherited

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

Definition at line 85 of file Boundary2D.h.

std::vector<std::pair<int, int> > VM2D::Boundary::vortexBeginEnd
inherited

Номера первого и последнего вихрей, рождаемых на каждой панели профиля (формируется после решения СЛАУ)

Definition at line 82 of file Boundary2D.h.

const World2D& VM2D::Boundary::W
protectedinherited

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

Definition at line 67 of file Boundary2D.h.


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