VM2D  1.12
Vortex methods for 2D flows simulation
Airfoil2D.h
Go to the documentation of this file.
1 /*--------------------------------*- VM2D -*-----------------*---------------*\
2 | ## ## ## ## #### ##### | | Version 1.12 |
3 | ## ## ### ### ## ## ## ## | VM2D: Vortex Method | 2024/01/14 |
4 | ## ## ## # ## ## ## ## | for 2D Flow Simulation *----------------*
5 | #### ## ## ## ## ## | Open Source Code |
6 | ## ## ## ###### ##### | https://www.github.com/vortexmethods/VM2D |
7 | |
8 | Copyright (C) 2017-2024 I. Marchevsky, K. Sokol, E. Ryatina, A. Kolganova |
9 *-----------------------------------------------------------------------------*
10 | File name: Airfoil2D.h |
11 | Info: Source code of VM2D |
12 | |
13 | This file is part of VM2D. |
14 | VM2D is free software: you can redistribute it and/or modify it |
15 | under the terms of the GNU General Public License as published by |
16 | the Free Software Foundation, either version 3 of the License, or |
17 | (at your option) any later version. |
18 | |
19 | VM2D is distributed in the hope that it will be useful, but WITHOUT |
20 | ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
21 | FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
22 | for more details. |
23 | |
24 | You should have received a copy of the GNU General Public License |
25 | along with VM2D. If not, see <http://www.gnu.org/licenses/>. |
26 \*---------------------------------------------------------------------------*/
27 
28 
40 #ifndef AIRFOIL_H
41 #define AIRFOIL_H
42 
43 #include "defs.h"
44 
45 namespace VM2D
46 {
47 
48  class World2D;
49  class WakeDataBase;
50 
60  class Airfoil
61  {
62  protected:
64  std::vector<Point2D> r_;
65 
67  std::vector<Point2D> v_;
68 
69  public:
71  const World2D& W;
72 
74  const size_t numberInPassport;
75 
78 
80  double phiAfl;
81 
83  double area;
84 
86  double m;
87 
89  double J;
90 
92  std::vector<double> meanEpsOverPanel;
93 
101  const Point2D& getR(size_t q) const
102  {
103  return (q < r_.size()) ? r_[q] : r_[0];
104  };
105 
113  const Point2D& getV(size_t q) const
114  {
115  return (q < v_.size()) ? v_[q] : v_[0];
116  };
117 
121  void setV(const Point2D& vel)
122  {
123  v_.clear();
124  v_.resize(r_.size(), vel);
125  }
126 
130  void setV(const std::vector<Point2D>& vel)
131  {
132  v_.clear();
133  v_.insert(v_.end(), vel.begin(), vel.end());
134  }
135 
139  size_t getNumberOfPanels() const { return r_.size(); };
140 
142  bool inverse;
143 
145  IFCUDA(mutable double* devRPtr);
146 
148  IFCUDA(mutable double* devRhsPtr);
149 
151  IFCUDA(mutable double* devRhsLinPtr);
152 
154  IFCUDA(mutable std::vector<double> tmpRhs);
155 
157  IFCUDA(mutable double* devFreeVortexSheetPtr);
158 
160  IFCUDA(mutable double* devAttachedVortexSheetPtr);
161 
163  IFCUDA(mutable double* devAttachedSourceSheetPtr);
164 
166  IFCUDA(mutable double* devMeanEpsOverPanelPtr);
167 
169  IFCUDA(mutable double* devViscousStressesPtr);
170 
172  IFCUDA(mutable std::vector<double> tmpViscousStresses);
173 
177  std::vector<Point2D> nrm;
178 
182  std::vector<Point2D> tau;
183 
185  std::vector<double> len;
186 
188  std::vector<double> viscousStress;
189 
192 
196  std::vector<double> gammaThrough;
197 
201  Airfoil(const World2D& W_, const size_t numberInPassport_);
202 
203 
205  virtual ~Airfoil() { };
206 
207 
213  bool isAfter(size_t i, size_t j) const;
214 
219  virtual void Rotate(double alpha) = 0;
220 
225  virtual void Scale(const Point2D&) = 0;
226 
232  virtual void Move(const Point2D& dr) = 0;
233 
239  virtual void GetGabarits(double gap = 0.02) = 0;
240 
245  bool isInsideGabarits(const Point2D& r) const;
246 
251  bool isOutsideGabarits(const Point2D& r) const;
252 
257  virtual bool IsPointInAirfoil(const Point2D& point) const = 0;
258 
266  virtual void ReadFromFile(const std::string& dir) = 0;
267 
275  virtual std::vector<double> getA(size_t p, size_t i, const Airfoil& airfoil, size_t j) const = 0;
276 
283  virtual void calcIQ(size_t p, const Airfoil& otherAirfoil, std::pair<Eigen::MatrixXd, Eigen::MatrixXd>& matrPair) const = 0;
284 
285 
287  void calcMeanEpsOverPanel();
288 
289 
300  virtual void GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(const WakeDataBase& pointsDb, std::vector<double>& domainRadius, std::vector<double>& I0, std::vector<Point2D>& I3) = 0;
301 #if defined(USE_CUDA)
302  virtual void GPUGetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(const WakeDataBase& pointsDb, std::vector<double>& domainRadius, std::vector<double>& I0, std::vector<Point2D>& I3) = 0;
303 #endif
304  virtual void GetDiffVelocityI0I3ToWakeAndViscousStresses(const WakeDataBase& pointsDb, std::vector<double>& domainRadius, std::vector<double>& I0, std::vector<Point2D>& I3) = 0;
305 
314  virtual void GetInfluenceFromVorticesToPanel(size_t panel, const Vortex2D* ptr, ptrdiff_t count, std::vector<double>& panelRhs) const = 0;
315 
316 
325  virtual void GetInfluenceFromSourcesToPanel(size_t panel, const Vortex2D* ptr, ptrdiff_t count, std::vector<double>& panelRhs) const = 0;
326 
327 
333  virtual void GetInfluenceFromSourceSheetToVortex(size_t panel, const Vortex2D& vtx, Point2D& vel) const = 0;
334 
340  virtual void GetInfluenceFromVortexSheetToVortex(size_t panel, const Vortex2D& vtx, Point2D& vel) const = 0;
341 
347  virtual void GetInfluenceFromVInfToPanel(std::vector<double>& vInfRhs) const = 0;
348 
349 
350 };
351 
352 }//namespace VM2D
353 
354 #endif
virtual void GetDiffVelocityI0I3ToWakeAndViscousStresses(const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)=0
std::vector< double > viscousStress
Касательные напряжения на панелях профиля
Definition: Airfoil2D.h:188
void setV(const std::vector< Point2D > &vel)
Установка скоростей всех вершин профиля
Definition: Airfoil2D.h:130
Airfoil(const World2D &W_, const size_t numberInPassport_)
Definition: Airfoil2D.cpp:55
double m
Масса профиля
Definition: Airfoil2D.h:86
virtual void GetInfluenceFromSourceSheetToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения ...
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
virtual void Scale(const Point2D &)=0
Масштабирование профиля
double J
Полярный момент инерции профиля относительно центра масс
Definition: Airfoil2D.h:89
Point2D lowLeft
Левый нижний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:190
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
Описание базовых вспомогательных функций
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
const Point2D & getV(size_t q) const
Возврат константной ссылки на скорость вершины профиля
Definition: Airfoil2D.h:113
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
Point2D upRight
Правый верхний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:191
virtual ~Airfoil()
Деструктор
Definition: Airfoil2D.h:205
virtual void GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)=0
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек, обусловленных геометрией профиля, и вычисление вязкого трения
IFCUDA(mutable double *devRPtr)
Указатель на девайсе, где хранятся вершины профиля
bool isInsideGabarits(const Point2D &r) const
Определяет, находится ли точка с радиус-вектором внутри габаритного прямоугольника профиля ...
Definition: Airfoil2D.cpp:69
Definition: Airfoil2D.h:45
virtual void GetInfluenceFromVorticesToPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const =0
Вычисление влияния части подряд идущих вихрей из вихревого следа на панель для правой части ...
std::vector< Point2D > r_
Координаты начал панелей
Definition: Airfoil2D.h:64
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
virtual void ReadFromFile(const std::string &dir)=0
Считывание профиля из файла
double phiAfl
Поворот профиля
Definition: Airfoil2D.h:80
virtual void Rotate(double alpha)=0
Поворот профиля
virtual void GetInfluenceFromVInfToPanel(std::vector< double > &vInfRhs) const =0
Вычисление влияния набегающего потока на панель для правой части
virtual void calcIQ(size_t p, const Airfoil &otherAirfoil, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &matrPair) const =0
Вычисление коэффициентов матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2.
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
Point2D rcm
Положение центра масс профиля
Definition: Airfoil2D.h:77
virtual void GetGabarits(double gap=0.02)=0
Вычисляет габаритный прямоугольник профиля
void calcMeanEpsOverPanel()
Вычисление средних значений eps на панелях
Definition: Airfoil2D.cpp:82
void setV(const Point2D &vel)
Установка постоянной скорости всех вершин профиля
Definition: Airfoil2D.h:121
virtual bool IsPointInAirfoil(const Point2D &point) const =0
Определяет, находится ли точка с радиус-вектором внутри профиля
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition: Airfoil2D.h:177
Абстрактный класс, определяющий обтекаемый профиль
Definition: Airfoil2D.h:60
std::vector< double > meanEpsOverPanel
Средние значения Eps на панелях
Definition: Airfoil2D.h:92
Класс, опеделяющий текущую решаемую задачу
Definition: World2D.h:68
Класс, опеделяющий набор вихрей
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
bool isAfter(size_t i, size_t j) const
Проверка, идет ли вершина i следом за вершиной j.
Definition: Airfoil2D.cpp:62
virtual void GetInfluenceFromVortexSheetToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вих...
virtual void GetInfluenceFromSourcesToPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const =0
Вычисление влияния части подряд идущих источников из области течения на панель для правой части ...
virtual void Move(const Point2D &dr)=0
Перемещение профиля
bool isOutsideGabarits(const Point2D &r) const
Определяет, находится ли точка с радиус-вектором вне габаритного прямоугольника профиля ...
Definition: Airfoil2D.cpp:76