VM2D  1.12
Vortex methods for 2D flows simulation
VM2D::Airfoil Class Referenceabstract

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

#include <Airfoil2D.h>

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

Public Member Functions

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

Public Attributes

const World2DW
 Константная ссылка на решаемую задачу More...
 
const size_t numberInPassport
 Номер профиля в паспорте More...
 
Point2D rcm
 Положение центра масс профиля More...
 
double phiAfl
 Поворот профиля More...
 
double area
 Площадь профиля More...
 
double m
 Масса профиля More...
 
double J
 Полярный момент инерции профиля относительно центра масс More...
 
std::vector< double > meanEpsOverPanel
 Средние значения Eps на панелях More...
 
bool inverse
 Признак разворота нормалей (для расчета внутренних течений) More...
 
std::vector< Point2Dnrm
 Нормали к панелям профиля More...
 
std::vector< Point2Dtau
 Касательные к панелям профиля More...
 
std::vector< double > len
 Длины панелей профиля More...
 
std::vector< double > viscousStress
 Касательные напряжения на панелях профиля More...
 
Point2D lowLeft
 Левый нижний угол габаритного прямоугольника профиля More...
 
Point2D upRight
 Правый верхний угол габаритного прямоугольника профиля More...
 
std::vector< double > gammaThrough
 Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге More...
 

Protected Attributes

std::vector< Point2Dr_
 Координаты начал панелей More...
 
std::vector< Point2Dv_
 Скорости начал панелей More...
 

Detailed Description

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

Author
Марчевский Илья Константинович
Сокол Ксения Сергеевна
Рятина Евгения Павловна
Колганова Александра Олеговна 1.12
Date
14 января 2024 г.

Definition at line 60 of file Airfoil2D.h.

Constructor & Destructor Documentation

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

Конструктор

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

Definition at line 55 of file Airfoil2D.cpp.

56  : W(W_), numberInPassport(numberInPassport_)
57 { }
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
virtual VM2D::Airfoil::~Airfoil ( )
inlinevirtual

Деструктор

Definition at line 205 of file Airfoil2D.h.

205 { };

Here is the call graph for this function:

Member Function Documentation

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

void Airfoil::calcMeanEpsOverPanel ( )

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

Definition at line 82 of file Airfoil2D.cpp.

83 {
84  meanEpsOverPanel.clear();
86 
87 
88  double midEps;
89 
90  const Boundary& bnd = W.getBoundary(numberInPassport);
92 
93  for (size_t i = 0; i < getNumberOfPanels(); ++i)
94  {
95  /*
96  midEps = 0.0;
97 
98  for (int j = bnd.vortexBeginEnd[i].first; j < bnd.vortexBeginEnd[i].second; ++j)
99  midEps += virtVortParams.epsastWake[j];
100 
101  midEps /= (bnd.vortexBeginEnd[i].second - bnd.vortexBeginEnd[i].first);
102 
103  meanEpsOverPanel[i] = midEps;
104  */
105 
106  /*
107  midEps = 0.0;
108 
109  for (int j = bnd.vortexBeginEnd[i].first; j < bnd.vortexBeginEnd[i].second; ++j)
110  midEps += std::max(virtVortParams.epsastWake[j], 0.5 * len[i] / (bnd.vortexBeginEnd[i].second - bnd.vortexBeginEnd[i].first));
111 
112  midEps /= (bnd.vortexBeginEnd[i].second - bnd.vortexBeginEnd[i].first);
113  meanEpsOverPanel[i] = midEps;
114  //*/
115 
116  midEps = 0.0;
117  for (int j = bnd.vortexBeginEnd[i].first; j < bnd.vortexBeginEnd[i].second; ++j)
118  midEps += virtVortParams.epsastWake[j];
119 
120  midEps /= (bnd.vortexBeginEnd[i].second - bnd.vortexBeginEnd[i].first);
121  meanEpsOverPanel[i] = midEps;
122  }//for i
123 }
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const Velocity & getVelocity() const
Возврат константной ссылки на объект для вычисления скоростей
Definition: World2D.h:212
Абстрактный класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле ...
Definition: Boundary2D.h:63
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition: Velocity2D.h:108
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition: Velocity2D.h:85
Структура, определяющая параметры виртуальных вихрей для отдельного профиля
Definition: Velocity2D.h:64
std::vector< double > meanEpsOverPanel
Средние значения Eps на панелях
Definition: Airfoil2D.h:92
std::vector< std::pair< int, int > > vortexBeginEnd
Номера первого и последнего вихрей, рождаемых на каждой панели профиля (формируется после решения СЛА...
Definition: Boundary2D.h:82

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

size_t VM2D::Airfoil::getNumberOfPanels ( ) const
inline

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

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

Definition at line 139 of file Airfoil2D.h.

139 { return r_.size(); };
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64

Here is the caller graph for this function:

const Point2D& VM2D::Airfoil::getR ( size_t  q) const
inline

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

Организовано "зацикливание" в сторону увеличения индекса, т.е. getR[size()] = getR[0];
Это позволяет удобно обращаться к getR(i) и getR(i+1) как к началу и концу i-й панели

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

Definition at line 101 of file Airfoil2D.h.

102  {
103  return (q < r_.size()) ? r_[q] : r_[0];
104  };
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64

Here is the caller graph for this function:

const Point2D& VM2D::Airfoil::getV ( size_t  q) const
inline

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

Организовано "зацикливание" в сторону увеличения индекса, т.е. getV[size()] = getV[0];
Это позволяет удобно обращаться к getV(i) и getV(i+1) как к скоростям начала и конца i-й панели

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

Definition at line 113 of file Airfoil2D.h.

114  {
115  return (q < v_.size()) ? v_[q] : v_[0];
116  };
std::vector< Point2D > v_
Скорости начал панелей
Definition: Airfoil2D.h:67

Here is the caller graph for this function:

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

Definition at line 62 of file Airfoil2D.cpp.

63 {
64  return ((i == j + 1) || (i == 0 && j == r_.size() - 1));
65 }//isAfter(...)
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64

Here is the caller graph for this function:

bool Airfoil::isInsideGabarits ( const Point2D r) const

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

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

Definition at line 69 of file Airfoil2D.cpp.

70 {
71  return (r[0] <= upRight[0] && (r[0] >= lowLeft[0] && r[1] >= lowLeft[1] && r[1] <= upRight[1]));
72 }//isInsideGabarits(...)
Point2D lowLeft
Левый нижний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:190
Point2D upRight
Правый верхний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:191

Here is the caller graph for this function:

bool Airfoil::isOutsideGabarits ( const Point2D r) const

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

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

Definition at line 76 of file Airfoil2D.cpp.

77 {
78  return (r[0] > upRight[0] || (r[0] < lowLeft[0] || r[1] < lowLeft[1] || r[1] > upRight[1]));
79 }//isOutsideGabarits(...)
Point2D lowLeft
Левый нижний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:190
Point2D upRight
Правый верхний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:191

Here is the caller graph for this function:

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

virtual void VM2D::Airfoil::Move ( const Point2D dr)
pure virtual

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

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

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

virtual void VM2D::Airfoil::Rotate ( double  alpha)
pure virtual

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

virtual void VM2D::Airfoil::Scale ( const Point2D )
pure virtual

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

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

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

Implemented in VM2D::AirfoilRect.

Here is the caller graph for this function:

void VM2D::Airfoil::setV ( const Point2D vel)
inline

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

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

Definition at line 121 of file Airfoil2D.h.

122  {
123  v_.clear();
124  v_.resize(r_.size(), vel);
125  }
std::vector< Point2D > v_
Скорости начал панелей
Definition: Airfoil2D.h:67
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64

Here is the caller graph for this function:

void VM2D::Airfoil::setV ( const std::vector< Point2D > &  vel)
inline

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

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

Definition at line 130 of file Airfoil2D.h.

131  {
132  v_.clear();
133  v_.insert(v_.end(), vel.begin(), vel.end());
134  }
std::vector< Point2D > v_
Скорости начал панелей
Definition: Airfoil2D.h:67

Member Data Documentation

double VM2D::Airfoil::area

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

Definition at line 83 of file Airfoil2D.h.

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

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

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

Definition at line 196 of file Airfoil2D.h.

bool VM2D::Airfoil::inverse

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

Definition at line 139 of file Airfoil2D.h.

double VM2D::Airfoil::J

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

Definition at line 89 of file Airfoil2D.h.

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

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

Definition at line 185 of file Airfoil2D.h.

Point2D VM2D::Airfoil::lowLeft

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

Definition at line 190 of file Airfoil2D.h.

double VM2D::Airfoil::m

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

Definition at line 86 of file Airfoil2D.h.

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

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

Definition at line 92 of file Airfoil2D.h.

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

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

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

Definition at line 177 of file Airfoil2D.h.

const size_t VM2D::Airfoil::numberInPassport

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

Definition at line 74 of file Airfoil2D.h.

double VM2D::Airfoil::phiAfl

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

Definition at line 80 of file Airfoil2D.h.

std::vector<Point2D> VM2D::Airfoil::r_
protected

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

Definition at line 64 of file Airfoil2D.h.

Point2D VM2D::Airfoil::rcm

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

Definition at line 77 of file Airfoil2D.h.

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

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

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

Definition at line 182 of file Airfoil2D.h.

Point2D VM2D::Airfoil::upRight

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

Definition at line 191 of file Airfoil2D.h.

std::vector<Point2D> VM2D::Airfoil::v_
protected

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

Definition at line 67 of file Airfoil2D.h.

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

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

Definition at line 188 of file Airfoil2D.h.

const World2D& VM2D::Airfoil::W

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

Definition at line 71 of file Airfoil2D.h.


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