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

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

#include <Airfoil2DRect.h>

Inheritance diagram for VM2D::AirfoilRect:
Collaboration diagram for VM2D::AirfoilRect:

Public Member Functions

 AirfoilRect (const World2D &W_, const size_t numberInPassport_)
 Конструктор More...
 
 AirfoilRect (const Airfoil &afl)
 
virtual ~AirfoilRect ()
 Деструктор More...
 
void CalcNrmTauLen ()
 Вычисление нормалей, касательных и длин панелей по текущему положению вершин More...
 
virtual void GetGabarits (double gap=0.02) override
 Вычисляет габаритный прямоугольник профиля More...
 
virtual void ReadFromFile (const std::string &dir) override
 Считывание профиля из файла More...
 
virtual void Rotate (double alpha) override
 Поворот профиля More...
 
virtual void Scale (const Point2D &) override
 Масштабирование профиля More...
 
virtual void Move (const Point2D &dr) override
 Перемещение профиля More...
 
virtual std::vector< double > getA (size_t p, size_t i, const Airfoil &airfoil, size_t j) const override
 Вычисление коэффициентов матрицы A для расчета влияния панели на панель More...
 
virtual void calcIQ (size_t p, const Airfoil &otherAirfoil, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &matrPair) const override
 Вычисление коэффициентов матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2. More...
 
virtual bool IsPointInAirfoil (const Point2D &point) const override
 Определяет, находится ли точка с радиус-вектором \( \vec r \) внутри профиля More...
 
virtual void GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses (const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3) override
 Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек, обусловленных геометрией профиля, и вычисление вязкого трения More...
 
virtual void GetDiffVelocityI0I3ToWakeAndViscousStresses (const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3) override
 
virtual void GetInfluenceFromVorticesToPanel (size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const override
 Вычисление влияния части подряд идущих вихрей из вихревого следа на панель для правой части More...
 
virtual void GetInfluenceFromSourcesToPanel (size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const override
 Вычисление влияния части подряд идущих источников из области течения на панель для правой части More...
 
virtual void GetInfluenceFromSourceSheetToVortex (size_t panel, const Vortex2D &vtx, Point2D &vel) const override
 Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения More...
 
virtual void GetInfluenceFromVortexSheetToVortex (size_t panel, const Vortex2D &vtx, Point2D &vel) const override
 Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вихрь в области течения More...
 
virtual void GetInfluenceFromVInfToPanel (std::vector< double > &vInfRhs) const override
 Вычисление влияния набегающего потока на панель для правой части More...
 
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...
 
bool isAfter (size_t i, size_t j) const
 Проверка, идет ли вершина i следом за вершиной j. More...
 
bool isInsideGabarits (const Point2D &r) const
 Определяет, находится ли точка с радиус-вектором \( \vec r \) внутри габаритного прямоугольника профиля More...
 
bool isOutsideGabarits (const Point2D &r) const
 Определяет, находится ли точка с радиус-вектором \( \vec r \) вне габаритного прямоугольника профиля More...
 
void calcMeanEpsOverPanel ()
 Вычисление средних значений eps на панелях 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
Марчевский Илья Константинович
Сокол Ксения Сергеевна
Рятина Евгения Павловна
Колганова Александра Олеговна
\Version 1.12
\date 14 января 2024 г.

Definition at line 64 of file Airfoil2DRect.h.

Constructor & Destructor Documentation

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

Конструктор

Definition at line 70 of file Airfoil2DRect.h.

71  :Airfoil(W_, numberInPassport_)
72  { };
Airfoil(const World2D &W_, const size_t numberInPassport_)
Definition: Airfoil2D.cpp:55
VM2D::AirfoilRect::AirfoilRect ( const Airfoil afl)
inline

Definition at line 74 of file Airfoil2DRect.h.

74 : Airfoil(afl) {};
Airfoil(const World2D &W_, const size_t numberInPassport_)
Definition: Airfoil2D.cpp:55
virtual VM2D::AirfoilRect::~AirfoilRect ( )
inlinevirtual

Деструктор

Definition at line 77 of file Airfoil2DRect.h.

77 { };

Here is the call graph for this function:

Member Function Documentation

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

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

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

Implements VM2D::Airfoil.

Definition at line 708 of file Airfoil2DRect.cpp.

709 {
710  bool self = (&otherAirfoil == this);
711 
712  size_t aflNSelf = numberInPassport;
713  size_t aflNOther = otherAirfoil.numberInPassport;
714 
715  size_t npI;
716  size_t npJ;
717 
718  numvector<double, 3> alpha, lambda;
719 
720  //auxillary vectors
721  Point2D p1, s1, p2, s2, di, dj, i00, i01, i10, i11;
722  numvector<Point2D, 3> v00, v11;
723  numvector<Point2D, 2> v01, v10;
724 
725 #pragma omp parallel for \
726  default(none) \
727  shared(otherAirfoil, self, aflNSelf, aflNOther, matrPair, p, IDPI, IQPI) \
728  private(npI, npJ, alpha, lambda, p1, s1, p2, s2, di, dj, i00, i01, i10, i11, v00, v11, v01, v10) schedule(dynamic, DYN_SCHEDULE)
729  for(int i = 0; i < (int)getNumberOfPanels(); ++i)
730  for (int j = 0; j < (int)otherAirfoil.getNumberOfPanels(); ++j)
731  {
732  npI = getNumberOfPanels();
733  npJ = otherAirfoil.getNumberOfPanels();
734 
735  if ((i == j) && self)
736  {
737 
738  matrPair.first(i, j) = 0.0;
739  matrPair.second(i, j) = 0.0;
740 
741  if (p == 2)
742  {
743  matrPair.first(i, npJ + j) = 0.0;
744  matrPair.first(npI + i, j) = 0.0;
745 
746  matrPair.second(i, npJ + j) = -IQPI;
747  matrPair.second(npI + i, j) = IQPI;
748 
749  matrPair.first(npI + i, npJ + j) = 0.0;
750  matrPair.second(npI + i, npJ + j) = 0.0;
751 
752  }
753  if (p > 2)
754  throw (-42);
755  }//if i==j
756  else
757  {
758  const Point2D& taui = tau[i];
759  const Point2D& tauj = otherAirfoil.tau[j];
760 
761  p1 = getR(i + 1) - otherAirfoil.getR(j + 1);
762  s1 = getR(i + 1) - otherAirfoil.getR(j);
763  p2 = getR(i) - otherAirfoil.getR(j + 1);
764  s2 = getR(i) - otherAirfoil.getR(j);
765  di = getR(i + 1) - getR(i);
766  dj = otherAirfoil.getR(j + 1) - otherAirfoil.getR(j);
767 
768  alpha = { \
769  (self && isAfter(j, i)) ? 0.0 : VMlib::Alpha(s2, s1), \
770  VMlib::Alpha(s2, p1), \
771  (self && isAfter(i, j)) ? 0.0 : VMlib::Alpha(p1, p2) \
772  };
773 
774  lambda = { \
775  (self && isAfter(j, i)) ? 0.0 : VMlib::Lambda(s2, s1), \
776  VMlib::Lambda(s2, p1), \
777  (self && isAfter(i, j)) ? 0.0 : VMlib::Lambda(p1, p2) \
778  };
779 
780  v00 = {
781  VMlib::Omega(s1, taui, tauj),
782  -VMlib::Omega(di, taui, tauj),
783  VMlib::Omega(p2, taui, tauj)
784  };
785 
786  i00 = IDPI / len[i] * (-(alpha[0] * v00[0] + alpha[1] * v00[1] + alpha[2] * v00[2]).kcross() \
787  + (lambda[0] * v00[0] + lambda[1] * v00[1] + lambda[2] * v00[2]));
788 
789  matrPair.first(i, j) = i00 & nrm[i];
790  matrPair.second(i, j) = i00 & taui;
791 
792 
793  if (p == 2)
794  {
795  v01 = {
796  0.5 / (dj.length()) * (((p1 + s1) & tauj) * VMlib::Omega(s1, taui, tauj) - s1.length2() * taui),
797  -0.5 * di.length() / dj.length() * VMlib::Omega(s1 + p2, tauj, tauj)
798  };
799 
800  i01 = IDPI / len[i] * (-((alpha[0] + alpha[2]) * v01[0] + (alpha[1] + alpha[2]) * v01[1]).kcross() \
801  + ((lambda[0] + lambda[2]) * v01[0] + (lambda[1] + lambda[2]) * v01[1]) - 0.5 * di.length() * tauj);
802 
803  matrPair.first(i, npJ + j) = i01 & nrm[i];
804  matrPair.second(i, npJ + j) = i01 & taui;
805 
806 
807  v10 = {
808  -0.5 / di.length() * (((s1 + s2) & taui) * VMlib::Omega(s1, taui, tauj) - s1.length2() * tauj),
809  0.5 * dj.length() / di.length() * VMlib::Omega(s1 + p2, taui, taui)
810  };
811 
812  i10 = IDPI / len[i] * (-((alpha[0] + alpha[2]) * v10[0] + alpha[2] * v10[1]).kcross() \
813  + ((lambda[0] + lambda[2]) * v10[0] + lambda[2] * v10[1]) + 0.5 * dj.length() * taui);
814 
815  matrPair.first(npI + i, j) = i10 & nrm[i];
816  matrPair.second(npI + i, j) = i10 & taui;
817 
818 
819  v11 = {
820  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),
821  -di.length() / (12.0 * dj.length()) * VMlib::Omega(di, tauj, tauj),
822  -dj.length() / (12.0 * di.length()) * VMlib::Omega(dj, taui, taui)
823  };
824 
825  i11 = IDPI / len[i] * (-((alpha[0] + alpha[2]) * v11[0] + (alpha[1] + alpha[2]) * v11[1] + alpha[2] * v11[2]).kcross()\
826  + (lambda[0] + lambda[2]) * v11[0] + (lambda[1] + lambda[2]) * v11[1] + lambda[2] * v11[2] \
827  + 1.0 / 12.0 * (dj.length() * taui + di.length() * tauj - 2.0 * VMlib::Omega(s1, taui, tauj)));
828 
829  matrPair.first(npI + i, npJ + j) = i11 & nrm[i];
830  matrPair.second(npI + i, npJ + j) = i11 & taui;
831  }
832 
833 
834  if (p > 2)
835  throw (-42);
836  }//else(i == j)
837 
838  }//for(...)
839 }//getIQ(...)
Шаблонный класс, определяющий вектор фиксированной длины Фактически представляет собой массив...
Definition: numvector.h:95
const double IDPI
Число .
Definition: defs.h:76
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
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
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition: numvector.h:383
const double IQPI
Число .
Definition: defs.h:82
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
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
bool isAfter(size_t i, size_t j) const
Проверка, идет ли вершина i следом за вершиной j.
Definition: Airfoil2D.cpp:62

Here is the call graph for this function:

Here is the caller graph for this function:

void Airfoil::calcMeanEpsOverPanel ( )
inherited

Вычисление средних значений 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:

void AirfoilRect::CalcNrmTauLen ( )

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

Definition at line 616 of file Airfoil2DRect.cpp.

617 {
618  if (nrm.size() != r_.size())
619  {
620  nrm.resize(r_.size());
621  tau.resize(r_.size());
622  len.resize(r_.size());
623  }
624 
625  Point2D rpan;
626 
627  for (size_t i = 0; i < r_.size(); ++i)
628  {
629  rpan = (getR(i + 1) - getR(i));
630  len[i] = rpan.length();
631  tau[i] = rpan.unit();
632 
633  nrm[i] = { tau[i][1], -tau[i][0] };
634  }
635 
636  //закольцовываем
637  nrm.push_back(nrm[0]);
638  tau.push_back(tau[0]);
639  len.push_back(len[0]);
640 }//CalcNrmTauLen()
P length() const
Вычисление 2-нормы (длины) вектора
Definition: numvector.h:371
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Класс, опеделяющий двумерный вектор
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
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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

Implements VM2D::Airfoil.

Definition at line 665 of file Airfoil2DRect.cpp.

666 {
667  std::vector<double> res(p*p, 0.0);
668 
669  if ( (i == j) && (&airfoil == this))
670  {
671  res[0] = 0.5 * (tau[i] ^ nrm[i]);
672  if (p == 1)
673  return res;
674 
675  res[3] = (1.0 / 12.0) * res[0];
676  if (p == 2)
677  return res;
678 
679  if (p > 2)
680  throw (-42);
681 
682  }//if i==j
683 
684  const auto& miq = W.getIQ(numberInPassport, airfoil.numberInPassport);
685 
686  res[0] = miq.first(i, j);
687 
688  if (p == 1)
689  return res;
690 
691  res[1] = miq.first(i, airfoil.getNumberOfPanels() + j);
692 
693  res[2] = miq.first(getNumberOfPanels() + i, j);
694 
695  res[3] = miq.first(getNumberOfPanels() + i, airfoil.getNumberOfPanels() + j);
696 
697  if (p == 2)
698  return res;
699 
700  if (p > 2)
701  throw (-42);
702 
703  //Хотя сюда никогда и не попадем
704  return res;
705 }//getA(...)
const std::pair< Eigen::MatrixXd, Eigen::MatrixXd > & getIQ(size_t i, size_t j) const
Возврат константной ссылки на объект, связанный с матрицей интегралов от (r-xi)/|r-xi|^2.
Definition: World2D.h:237
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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

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

Implements VM2D::Airfoil.

Definition at line 303 of file Airfoil2DRect.cpp.

304 {
305  std::vector<double> selfI0;
306  std::vector<Point2D> selfI3;
307 
308  std::vector<double> currViscousStress;
309  currViscousStress.resize(r_.size(), 0.0);
310 
311 
312  selfI0.resize(pointsDb.vtx.size(), 0.0);
313  selfI3.resize(pointsDb.vtx.size(), { 0.0, 0.0 });
314 
315 #pragma warning (push)
316 #pragma warning (disable: 4101)
317  //Локальные переменные для цикла
318  Point2D q, xi, xi_m, v0;
319  double lxi, lxi_m, lenj_m;
320 
321  Point2D mn;
322  int new_n;
323  Point2D h;
324 
325  Point2D vec;
326  double s, d;
327 
328  double vs;
329  double iDDomRad, domRad, expon;
330 #pragma warning (pop)
331 
332 #pragma omp parallel for \
333  default(none) \
334  shared(domainRadius, pointsDb, selfI0, selfI3, currViscousStress, PI) \
335  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)
336  for (int i = 0; i < pointsDb.vtx.size(); ++i)
337  {
338  const Vortex2D& vtxI = pointsDb.vtx[i];
339 
340  domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
341  iDDomRad = 1.0 / domRad;
342 
343  for (size_t j = 0; j < r_.size(); ++j)
344  {
345  vs = 0.0;
346  q = vtxI.r() - 0.5 * (getR(j) + getR(j + 1));
347  vec = tau[j];
348 
349  s = q & vec;
350  d = fabs(q & nrm[j]);
351 
352  if ((d < 50.0 * len[j]) && (fabs(s) < 50.0 * len[j]))
353  {
354  v0 = vec * len[j];
355 
356  if ((d > 5.0 * len[j]) || (fabs(s) > 5.0 * len[j]))
357  {
358  xi = q * iDDomRad;
359  lxi = xi.length();
360 
361  expon = exp(-lxi) * len[j];
362  mn = nrm[j] * expon;
363 
364  if (selfI0[i] != -PI * domRad)
365  {
366  selfI0[i] += (xi & mn) * (lxi + 1.0) / (lxi * lxi);
367  selfI3[i] += mn;
368  }
369 
370  vs = vtxI.g() * expon / (PI * sqr(meanEpsOverPanel[j]));
371  }
372  else if ((d >= 0.1 * len[j]) || (fabs(s) > 0.5 * len[j]))
373  {
374  vs = 0.0;
375  double den = (fabs(s) < 0.5 * len[j]) ? d : (fabs(s) + d - 0.5 * len[j]);
376 
377  new_n = std::max(1, static_cast<int>(ceil(5.0 * len[j] / den)));
378  //new_n = 100;
379  h = v0 * (1.0 / new_n);
380 
381  for (int m = 0; m < new_n; ++m)
382  {
383  xi_m = (vtxI.r() - (getR(j) + h * (m + 0.5))) * iDDomRad;
384  lxi_m = xi_m.length();
385 
386  lenj_m = len[j] / new_n;
387  expon = exp(-lxi_m) * lenj_m;
388 
389  mn = nrm[j] * expon;
390  if (selfI0[i] != -PI * domRad)
391  {
392  selfI0[i] += (xi_m & mn) * (lxi_m + 1.0) / (lxi_m * lxi_m);
393  selfI3[i] += mn;
394  }
395  vs += expon;
396  }//for m
397  vs *= vtxI.g() / (PI * sqr(meanEpsOverPanel[j]));
398  }
399  else
400  {
401 
402  selfI0[i] = -PI * domRad;
403  double mnog = 2.0 * domRad * (1.0 - exp(-len[j] * iDDomRad / 2.0) * cosh(fabs(s) * iDDomRad));
404  selfI3[i] = nrm[j] * mnog;
405  vs = mnog * vtxI.g() / (PI * sqr(meanEpsOverPanel[j]));
406  }
407  }//if d<50 len
408 
409 #pragma omp atomic
410  currViscousStress[j] += vs;
411  }//for j
412  }//for i
413 
414 
415  if (&pointsDb == &(W.getWake()))
416  for (size_t i = 0; i < viscousStress.size(); ++i)
417  {
418  viscousStress[i] += currViscousStress[i];
419  }
420 
421  for (size_t i = 0; i < I0.size(); ++i)
422  {
423  domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
424 
425  if (I0[i] != -PI * domRad)
426  {
427  if (selfI0[i] == -PI * domRad)
428  {
429  I0[i] = selfI0[i];
430  I3[i] = selfI3[i];
431  }
432  else
433  {
434  I0[i] += selfI0[i];
435  I3[i] += selfI3[i];
436  }
437  }
438  }
439 }; //GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(...)
std::vector< double > viscousStress
Касательные напряжения на панелях профиля
Definition: Airfoil2D.h:188
const double PI
Число .
Definition: defs.h:73
double m
Масса профиля
Definition: Airfoil2D.h:86
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
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
T sqr(T x)
Возведение числа в квадрат
Definition: defs.h:430
std::vector< Vortex2D > vtx
Список вихревых элементов
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177
std::vector< double > meanEpsOverPanel
Средние значения Eps на панелях
Definition: Airfoil2D.h:92
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
double getMinEpsAst() const
Функция минимально возможного значения для epsAst.
Definition: Passport2D.h:166

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements VM2D::Airfoil.

Definition at line 442 of file Airfoil2DRect.cpp.

443 {
444  if (&pointsDb == &(W.getWake()))
445  {
446  viscousStress.clear();
447  viscousStress.resize(r_.size(), 0.0);
448  }
449 
450  GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(pointsDb, domainRadius, I0, I3);
451 }//GetDiffVelocityI0I3ToWakeAndViscousStresses(...)
std::vector< double > viscousStress
Касательные напряжения на панелях профиля
Definition: Airfoil2D.h:188
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
virtual void GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3) override
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек, обусловленных геометрией профиля, и вычисление вязкого трения
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64

Here is the call graph for this function:

Here is the caller graph for this function:

void AirfoilRect::GetGabarits ( double  gap = 0.02)
overridevirtual

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

Implements VM2D::Airfoil.

Definition at line 643 of file Airfoil2DRect.cpp.

644 {
645  lowLeft = { 1E+10, 1E+10 };
646  upRight = { -1E+10, -1E+10 };
647 
648  for (size_t i = 0; i < r_.size(); ++i)
649  {
650  lowLeft[0] = std::min(lowLeft[0], r_[i][0]);
651  lowLeft[1] = std::min(lowLeft[1], r_[i][1]);
652 
653  upRight[0] = std::max(upRight[0], r_[i][0]);
654  upRight[1] = std::max(upRight[1], r_[i][1]);
655  }
656 
657  Point2D size = upRight - lowLeft;
658  lowLeft -= size*gap;
659  upRight += size*gap;
660 }//GetGabarits(...)
Point2D lowLeft
Левый нижний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:190
Point2D upRight
Правый верхний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:191
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75

Here is the caller graph for this function:

void AirfoilRect::GetInfluenceFromSourceSheetToVortex ( size_t  panel,
const Vortex2D vtx,
Point2D vel 
) const
overridevirtual

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

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

Implements VM2D::Airfoil.

Definition at line 855 of file Airfoil2DRect.cpp.

856 {
858 }//GetInfluenceFromSourceSheetToVortex(...)
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
virtual void GetInfluenceFromSourceSheetAtRectPanelToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения ...

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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

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

Implements VM2D::Airfoil.

Definition at line 849 of file Airfoil2DRect.cpp.

850 {
852 }//GetInfluenceFromSourcesToPanel(...)
virtual void GetInfluenceFromSourcesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const =0
Вычисление влияния части подряд источников из области течения на прямолинейную панель для правой част...
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153

Here is the call graph for this function:

Here is the caller graph for this function:

void AirfoilRect::GetInfluenceFromVInfToPanel ( std::vector< double > &  vInfRhs) const
overridevirtual

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

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

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

Implements VM2D::Airfoil.

Definition at line 866 of file Airfoil2DRect.cpp.

867 {
869 }//GetInfluenceFromVInfToPanel(...)
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
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:

void AirfoilRect::GetInfluenceFromVortexSheetToVortex ( size_t  panel,
const Vortex2D vtx,
Point2D vel 
) const
overridevirtual

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

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

Implements VM2D::Airfoil.

Definition at line 861 of file Airfoil2DRect.cpp.

862 {
864 }//GetInfluenceFromVortexSheetToVortex(...)
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
virtual void GetInfluenceFromVortexSheetAtRectPanelToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вих...
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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

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

Implements VM2D::Airfoil.

Definition at line 842 of file Airfoil2DRect.cpp.

843 {
845 }//GetInfluenceFromVorticesToPanel(...)
virtual void GetInfluenceFromVorticesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const =0
Вычисление влияния части подряд идущих вихрей из вихревого следа на прямолинейную панель для правой ч...
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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
inlineinherited

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

Организовано "зацикливание" в сторону увеличения индекса, т.е. 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
inlineinherited

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

Организовано "зацикливание" в сторону увеличения индекса, т.е. 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)
inherited

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

Проверка, идет ли вершина 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
inherited

Определяет, находится ли точка с радиус-вектором \( \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
inherited

Определяет, находится ли точка с радиус-вектором \( \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:

bool AirfoilRect::IsPointInAirfoil ( const Point2D point) const
overridevirtual

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

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

Implements VM2D::Airfoil.

Definition at line 560 of file Airfoil2DRect.cpp.

561 {
562  double sumAngle = 0.0;
563 
564  Point2D v1, v2;
565 
566  for (size_t i = 0; i < r_.size(); ++i)
567  {
568  v1 = getR(i) - point;
569  v2 = getR(i+1) - point;
570  sumAngle += atan2(v1^v2, v1 & v2);
571  }
572 
573  if (fabs(sumAngle) < 0.1)
574  return false;
575  else
576  return true;
577 }//IsPointInAirfoil(...)
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75

Here is the call graph for this function:

Here is the caller graph for this function:

void AirfoilRect::Move ( const Point2D dr)
overridevirtual

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

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

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

Implements VM2D::Airfoil.

Definition at line 580 of file Airfoil2DRect.cpp.

581 {
582  for (size_t i = 0; i < r_.size(); ++i)
583  r_[i] += dr;
584  rcm += dr;
585  CalcNrmTauLen();
586  GetGabarits();
587 }//Move(...)
virtual void GetGabarits(double gap=0.02) override
Вычисляет габаритный прямоугольник профиля
void CalcNrmTauLen()
Вычисление нормалей, касательных и длин панелей по текущему положению вершин
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Point2D rcm
Положение центра масс профиля
Definition: Airfoil2D.h:77

Here is the call graph for this function:

Here is the caller graph for this function:

void AirfoilRect::ReadFromFile ( const std::string &  dir)
overridevirtual

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

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

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

Implements VM2D::Airfoil.

Definition at line 59 of file Airfoil2DRect.cpp.

60 {
62  std::string filename = dir + param.fileAirfoil;
63  std::ifstream airfoilFile;
64 
65  if (fileExistTest(filename, W.getInfo(), { "txt", "TXT" }))
66  {
67  std::stringstream airfoilFile(VMlib::Preprocessor(filename).resultString);
68 
69  VMlib::StreamParser airfoilParser(W.getInfo(), "airfoil file parser", airfoilFile);
70 
71  m = 0.0; //TODO
72  J = 0.0; //TODO
73  rcm = { 0.0, 0.0 }; //TODO
74  phiAfl = 0.0;
75 
76  inverse = param.inverse;
77 
78 
79 
80  //*
81 
83  //Проверяем отсутствие "задвоенных" подряд идущих точек
84  //airfoilParser.get("r", r_);
85  std::vector<Point2D> rFromFile;
86  airfoilParser.get("r", rFromFile);
87 
88  if (rFromFile.size() > 0)
89  {
90  if (param.requiredNPanels <= rFromFile.size()) {
91 
92  r_.reserve(rFromFile.size());
93  r_.push_back(rFromFile[0]);
94  for (size_t i = 1; i < rFromFile.size(); ++i)
95  if ((rFromFile[i] - rFromFile[i - 1]).length2() > 1e-12)
96  r_.push_back(rFromFile[i]);
97 
98  //Если первая совпадает с последней, то убираем последнюю
99  if ((r_.back() - r_.front()).length2() < 1e-12)
100  r_.resize(r_.size() - 1);
101  }
102  else
103  {
104  W.getInfo('e') << "Airfoil shape is given explicitely, it can not be automatically split!" << std::endl;
105  exit(200);
106  }
107  }
108  else
109  {
110  //std::vector<GeomPoint> geomFromFile(X.size()); //new
111  std::vector<GeomPoint> geomFromFile; //VM2D
112 
113  //for (int q = 0; q < X.size(); ++q) // new
114  //{ //
115  // geomFromFile[q].r = { X[q], Y[q] }; //
116  // geomFromFile[q].type = key[q]; //
117  //} //
118  airfoilParser.get("geometry", geomFromFile); //VM2D
119 
120  if ((geomFromFile.front().r - geomFromFile.back().r).length2() < 1e-12)
121  geomFromFile.resize(geomFromFile.size() - 1);
122 
123 
124  //size_t reqN = NumberOfPoints; //new
125  size_t reqN = param.requiredNPanels; //VM2D
126 
127  std::vector<double> L(geomFromFile.size());
128  double totalLength = 0.0;
129  std::vector<size_t> nc = {};
130  std::vector<size_t> ni = {};
131 
132  for (size_t i = 0; i < geomFromFile.size(); ++i)
133  {
134  const Point2D& p1 = geomFromFile[i].r;
135  const Point2D& p2 = geomFromFile[(i + 1) % geomFromFile.size()].r;
136 
137  L[i] = (p2 - p1).length();
138  totalLength += L[i];
139  }
140 
141 
142  //Ищем первый сплайн
143  size_t i = 0;
144  std::vector<size_t> splineStart, splineFinish;
145 
146  while (i < geomFromFile.size())
147  {
148  while (i < geomFromFile.size() && !(geomFromFile[i].type == "c"))
149  ++i;
150 
151  if (i < geomFromFile.size())
152  {
153  splineStart.push_back(i);
154  ++i;
155  //ищем второй конец
156  while (i < geomFromFile.size() && !(geomFromFile[i].type == "c"))
157  ++i;
158 
159  splineFinish.push_back(i);
160  }
161  }
162 
163  if ((splineStart.size() > 0) && (splineStart[0] != 0))
164  splineFinish.back() = splineStart[0];
165 
166  //std::vector<Point2D> r_; //new
167 
168  if (reqN == 0)
169  {
170  r_.reserve(geomFromFile.size());
171  for (size_t i = 0; i < geomFromFile.size(); ++i)
172  r_.push_back(geomFromFile[i].r);
173  }
174  else
175  {
176  double hRef = totalLength / reqN;
177 
178  std::vector<int> nPanels;
179 
180  bool cyclic = (splineStart.size() == 0);
181 
182  if (cyclic)
183  {
184  splineStart.push_back(0);
185  splineFinish.push_back(geomFromFile.size());
186  }
187 
188  for (size_t s = 0; s < splineStart.size(); ++s)
189  {
190  //Сплайн
191  std::vector<double> T, X, Y;
192 
193  size_t splineLegs = splineFinish[s] - splineStart[s];
194  if (splineFinish[s] <= splineStart[s])
195  splineLegs += geomFromFile.size();
196 
197  T.reserve(splineLegs + 1);
198  X.reserve(T.size());
199  Y.reserve(T.size());
200 
201  double tbuf = 0.0;
202 
203  for (size_t i = splineStart[s]; i <= splineStart[s] + splineLegs; ++i)
204  {
205  const Point2D& pt = geomFromFile[i % geomFromFile.size()].r;
206 
207  T.push_back(tbuf);
208  X.push_back(pt[0]);
209  Y.push_back(pt[1]);
210 
211  if ((i == splineStart[s]) && (splineFinish[s] - splineStart[s] == 1))
212  {
213  double halfL = (i < geomFromFile.size()) ? 0.5 * L[i] : 0.5 * (geomFromFile.front().r - geomFromFile.back().r).length();
214  tbuf += halfL;
215  T.push_back(tbuf);
216  Point2D nextPt = geomFromFile[(i + 1) % geomFromFile.size()].r;
217 
218  X.push_back(0.5 * (pt[0] + nextPt[0]));
219  Y.push_back(0.5 * (pt[1] + nextPt[1]));
220 
221  tbuf += halfL;
222  }
223  else
224  tbuf += (i < geomFromFile.size()) ? L[i] : (geomFromFile.front().r - geomFromFile.back().r).length();
225  }
226 
227 
228 
229  tk::spline s1, s2;
230  if (cyclic)
231  {
232  s1.set_boundary(tk::spline::cyclic);
233  s2.set_boundary(tk::spline::cyclic);
234  }
235  else
236  {
237  s1.set_boundary(tk::spline::second_deriv, 0.0, tk::spline::second_deriv, 0.0);
238  s2.set_boundary(tk::spline::second_deriv, 0.0, tk::spline::second_deriv, 0.0);
239  }
240 
241  s1.set_points(T, X);
242  s2.set_points(T, Y);
243 
244  int NumberOfPanelsSpline = (int)ceil(T.back() / hRef);
245  double hSpline = T.back() / NumberOfPanelsSpline;
246 
247  for (int i = 0; i < NumberOfPanelsSpline; ++i)
248  r_.push_back({ s1(hSpline * i), s2(hSpline * i) });
249 
250  nPanels.push_back(NumberOfPanelsSpline);
251 
252  }
253  }
254  }//else geometry
255 
256 
257 
258  //std::ofstream of("airfoil.txt");
259  //for (size_t i = 0; i < r_.size(); ++i)
260  // of << r_[i] << "," << std::endl;
261  //of.close();
262 
263  if (r_.size() == 0)
264  {
265  W.getInfo('e') << "No points on airfoil contour!" << std::endl;
266  exit(200);
267  }
268 
269  if (inverse)
270  std::reverse(r_.begin(), r_.end());
271 
272 
273  v_.resize(0);
274  for (size_t q = 0; q < r_.size(); ++q)
275  v_.push_back({ 0.0, 0.0 });
276 
277  Move(param.basePoint);
278  Scale(param.scale);
279 
280  double rotationAngle = param.angle;
282  rotationAngle = -param.angle - 0.5 * PI;
283 
284  Rotate(-rotationAngle);
285  //в конце Rotate нормали, касательные и длины вычисляются сами
286 
287  gammaThrough.clear();
288  gammaThrough.resize(r_.size(), 0.0);
289 
290  //Вычисляем площадь
291  area = 0.0;
292  for (size_t q = 0; q < r_.size(); ++q)
293  {
294  const Point2D& cntq = 0.5 * (getR(q) + getR(q + 1));
295  const Point2D& drq = getR(q + 1) - getR(q);
296  area += fabs(cntq[0] * drq[1] - cntq[1] * drq[0]);
297  }
298  area *= 0.5;
299  }
300 }//ReadFromFile(...)
virtual void Rotate(double alpha) override
Поворот профиля
bool geographicalAngles
Признак работы в "географической" системе координат
Definition: Passport2D.h:283
const double PI
Число .
Definition: defs.h:73
double m
Масса профиля
Definition: Airfoil2D.h:86
bool inverse
Признак разворота нормалей (для расчета внутреннего течения)
Definition: Passport2D.h:237
double area
Площадь профиля
Definition: Airfoil2D.h:83
bool inverse
Признак разворота нормалей (для расчета внутренних течений)
Definition: Airfoil2D.h:139
const World2D & W
Константная ссылка на решаемую задачу
Definition: Airfoil2D.h:71
std::vector< Point2D > v_
Скорости начал панелей
Definition: Airfoil2D.h:67
std::vector< double > gammaThrough
Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
Definition: Airfoil2D.h:196
double J
Полярный момент инерции профиля относительно центра масс
Definition: Airfoil2D.h:89
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
Класс, позволяющий выполнять предварительную обработку файлов
Definition: Preprocessor.h:59
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
virtual void Scale(const Point2D &) override
Масштабирование профиля
size_t requiredNPanels
Желаемое число панелей для разбиения геометрии
Definition: Passport2D.h:219
double angle
Угол поворота (угол атаки)
Definition: Passport2D.h:231
virtual void Move(const Point2D &dr) override
Перемещение профиля
Point2D basePoint
Смещение центра масс (перенос профиля)
Definition: Passport2D.h:222
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
bool fileExistTest(std::string &fileName, LogStream &info, const std::list< std::string > &extList={})
Проверка существования файла
Definition: defs.h:324
Структура, задающая параметры профиля
Definition: Passport2D.h:213
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода ...
Definition: WorldGen.h:74
double phiAfl
Поворот профиля
Definition: Airfoil2D.h:80
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Point2D rcm
Положение центра масс профиля
Definition: Airfoil2D.h:77
std::string fileAirfoil
Имя файла с начальным состоянием профилей (без полного пути)
Definition: Passport2D.h:216
std::vector< AirfoilParams > airfoilParams
Список структур с параметрами профилей
Definition: Passport2D.h:280
Point2D scale
Коэффициент масштабирования
Definition: Passport2D.h:225
Класс, позволяющий выполнять разбор файлов и строк с настройками и параметрами
Definition: StreamParser.h:151

Here is the call graph for this function:

Here is the caller graph for this function:

void AirfoilRect::Rotate ( double  alpha)
overridevirtual

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

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

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

Implements VM2D::Airfoil.

Definition at line 591 of file Airfoil2DRect.cpp.

592 {
593  phiAfl += alpha;
594  nummatrix<double, 2, 2> rotMatrix = { { cos(alpha), -sin(alpha) }, { sin(alpha), cos(alpha) } };
595 
596  for (size_t i = 0; i < r_.size(); ++i)
597  r_[i] = rcm + (rotMatrix & (r_[i] - rcm));
598 
599 
600  CalcNrmTauLen();
601  GetGabarits();
602 }//Rotate(...)
virtual void GetGabarits(double gap=0.02) override
Вычисляет габаритный прямоугольник профиля
void CalcNrmTauLen()
Вычисление нормалей, касательных и длин панелей по текущему положению вершин
Шаблонный класс, определяющий матрицу фиксированного размера Фактически представляет собой массив...
Definition: nummatrix.h:62
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
double phiAfl
Поворот профиля
Definition: Airfoil2D.h:80
Point2D rcm
Положение центра масс профиля
Definition: Airfoil2D.h:77

Here is the call graph for this function:

Here is the caller graph for this function:

void AirfoilRect::Scale ( const Point2D )
overridevirtual

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

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

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

Implements VM2D::Airfoil.

Definition at line 606 of file Airfoil2DRect.cpp.

607 {
608  for (size_t i = 0; i < r_.size(); ++i)
609  r_[i] = rcm + Point2D{ factor[0] * (r_[i] - rcm)[0], factor[1] * (r_[i] - rcm)[1] };
610 
611  CalcNrmTauLen();
612  GetGabarits();
613 }//Scale(...)
virtual void GetGabarits(double gap=0.02) override
Вычисляет габаритный прямоугольник профиля
void CalcNrmTauLen()
Вычисление нормалей, касательных и длин панелей по текущему положению вершин
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
Point2D rcm
Положение центра масс профиля
Definition: Airfoil2D.h:77

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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)
inlineinherited

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

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
inherited

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

Definition at line 83 of file Airfoil2D.h.

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

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

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

Definition at line 196 of file Airfoil2D.h.

bool VM2D::Airfoil::inverse
inherited

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

Definition at line 139 of file Airfoil2D.h.

double VM2D::Airfoil::J
inherited

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

Definition at line 89 of file Airfoil2D.h.

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

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

Definition at line 185 of file Airfoil2D.h.

Point2D VM2D::Airfoil::lowLeft
inherited

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

Definition at line 190 of file Airfoil2D.h.

double VM2D::Airfoil::m
inherited

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

Definition at line 86 of file Airfoil2D.h.

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

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

Definition at line 92 of file Airfoil2D.h.

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

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

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

Definition at line 177 of file Airfoil2D.h.

const size_t VM2D::Airfoil::numberInPassport
inherited

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

Definition at line 74 of file Airfoil2D.h.

double VM2D::Airfoil::phiAfl
inherited

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

Definition at line 80 of file Airfoil2D.h.

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

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

Definition at line 64 of file Airfoil2D.h.

Point2D VM2D::Airfoil::rcm
inherited

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

Definition at line 77 of file Airfoil2D.h.

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

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

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

Definition at line 182 of file Airfoil2D.h.

Point2D VM2D::Airfoil::upRight
inherited

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

Definition at line 191 of file Airfoil2D.h.

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

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

Definition at line 67 of file Airfoil2D.h.

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

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

Definition at line 188 of file Airfoil2D.h.

const World2D& VM2D::Airfoil::W
inherited

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

Definition at line 71 of file Airfoil2D.h.


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