VM2D  1.12
Vortex methods for 2D flows simulation
Boundary2D.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: Boundary2D.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 BOUNDARY_H
41 #define BOUNDARY_H
42 
43 #include <memory>
44 #include "Sheet2D.h"
45 #include "VirtualWake2D.h"
46 
47 namespace VM2D
48 {
49 
50  class World2D;
51  class Airfoil;
52 
63  class Boundary
64  {
65  protected:
67  const World2D& W;
68 
70  const size_t numberInPassport;
71 
72  public:
73 
76  const Airfoil& afl;
77 
80 
82  std::vector<std::pair<int, int>> vortexBeginEnd;
83 
86 
92  size_t sheetDim;
93 
96 
99 
105  Boundary(const World2D& W_, size_t numberInPassport_, int sheetDim_);
106 
108  virtual ~Boundary() { };
109 
120  virtual void FillMatrixSelf(Eigen::MatrixXd& matr, Eigen::VectorXd& lastLine, Eigen::VectorXd& lactCol) = 0;
121 
125  virtual void FillIQSelf(std::pair<Eigen::MatrixXd, Eigen::MatrixXd>& IQ) = 0;
126 
134  virtual void FillMatrixFromOther(const Boundary& otherBoundary, Eigen::MatrixXd& matr) = 0;
135 
143  virtual void FillIQFromOther(const Boundary& otherBoundary, std::pair<Eigen::MatrixXd, Eigen::MatrixXd>& IQ) = 0;
144 
145 
154  virtual void CalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase& pointsDb, std::vector<Point2D>& velo) const = 0;
155 #if defined(USE_CUDA)
156  virtual void GPUCalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase& pointsDb, std::vector<Point2D>& velo) const = 0;
157 #endif
158 
167  virtual void CalcConvVelocityAtVirtualVortexes(std::vector<Point2D>& velo) const = 0;
168 
169 
175  size_t GetUnknownsSize() const;
176 
177 
187  virtual void SolutionToFreeVortexSheetAndVirtualVortex(const Eigen::VectorXd& sol) = 0;
188 
190  virtual void ComputeAttachedSheetsIntensity() = 0;
191 
192 
201 
202  virtual void GetInfluenceFromVorticesToRectPanel(size_t panel, const Vortex2D* ptr, ptrdiff_t count, std::vector<double>& wakeRhs) const = 0;
203 
212  virtual void GetInfluenceFromSourcesToRectPanel(size_t panel, const Vortex2D* ptr, ptrdiff_t count, std::vector<double>& wakeRhs) const = 0;
213 
219  virtual void GetInfluenceFromVortexSheetAtRectPanelToVortex(size_t panel, const Vortex2D& vtx, Point2D& vel) const = 0;
220 
226  virtual void GetInfluenceFromSourceSheetAtRectPanelToVortex(size_t panel, const Vortex2D& vtx, Point2D& vel) const = 0;
227 
233  virtual void GetInfluenceFromVInfToRectPanel(std::vector<double>& vInfRhs) const = 0;
234 
235 };
236 
237 }//namespace VM2D
238 
239 #endif
virtual void GetInfluenceFromSourcesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const =0
Вычисление влияния части подряд источников из области течения на прямолинейную панель для правой част...
virtual void ComputeAttachedSheetsIntensity()=0
Вычисление интенсивностей присоединенного вихревого слоя и присоединенного слоя источников ...
virtual void CalcConvVelocityAtVirtualVortexes(std::vector< Point2D > &velo) const =0
Вычисление конвективной скорости в точках расположения виртуальных вихрей
virtual void GetInfluenceFromVorticesToRectPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &wakeRhs) const =0
Вычисление влияния части подряд идущих вихрей из вихревого следа на прямолинейную панель для правой ч...
Заголовочный файл с описанием класса VirtualWake.
virtual void FillMatrixFromOther(const Boundary &otherBoundary, Eigen::MatrixXd &matr)=0
Генерация блока матрицы влияния от другого профиля того же типа
Абстрактный класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле ...
Definition: Boundary2D.h:63
Boundary(const World2D &W_, size_t numberInPassport_, int sheetDim_)
Конструктор
Definition: Boundary2D.cpp:53
virtual void GetInfluenceFromVortexSheetAtRectPanelToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния вихревых слоев (свободный + присоединенный) конкретной прямолинейной панели на вих...
int minVortexPerPanel
Минимальное число вихрей, рождаемых на панели профиля и формирующих виртуальный вихревой след ...
Definition: Boundary2D.h:79
Класс, опеделяющий слои на поверхности обтекаемого профиля
Definition: Sheet2D.h:61
Класс, опеделяющий вихревой след (пелену)
Definition: VirtualWake2D.h:60
virtual ~Boundary()
Деструктор
Definition: Boundary2D.h:108
virtual void FillIQFromOther(const Boundary &otherBoundary, std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &IQ)=0
Генерация блока матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2, влияние одного профиля на другой...
Заголовочный файл с описанием класса SheetV.
const size_t numberInPassport
Номер профиля в паспорте
Definition: Boundary2D.h:70
virtual void FillMatrixSelf(Eigen::MatrixXd &matr, Eigen::VectorXd &lastLine, Eigen::VectorXd &lactCol)=0
Генерация блока матрицы
virtual void CalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase &pointsDb, std::vector< Point2D > &velo) const =0
Вычисление конвективных скоростей в наборе точек, вызываемых наличием слоев вихрей и источников на пр...
Definition: Airfoil2D.h:45
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
virtual void GetInfluenceFromVInfToRectPanel(std::vector< double > &vInfRhs) const =0
Вычисление влияния набегающего потока на прямолинейную панель для правой части
Sheet oldSheets
Слои на профиле с предыдущего шага
Definition: Boundary2D.h:98
virtual void SolutionToFreeVortexSheetAndVirtualVortex(const Eigen::VectorXd &sol)=0
Пересчет решения на интенсивность вихревого слоя и на рождаемые вихри на конкретном профиле ...
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
virtual void FillIQSelf(std::pair< Eigen::MatrixXd, Eigen::MatrixXd > &IQ)=0
Генерация блока матрицы, состоящей из интегралов от (r-xi)/|r-xi|^2, влияние профиля самого на себя ...
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
Абстрактный класс, определяющий обтекаемый профиль
Definition: Airfoil2D.h:60
virtual void GetInfluenceFromSourceSheetAtRectPanelToVortex(size_t panel, const Vortex2D &vtx, Point2D &vel) const =0
Вычисление влияния слоя источников конкретной прямолинейной панели на вихрь в области течения ...
size_t sheetDim
Размерность параметров каждого из слоев на каждой из панелей
Definition: Boundary2D.h:92
Класс, опеделяющий текущую решаемую задачу
Definition: World2D.h:68
Класс, опеделяющий набор вихрей
std::vector< std::pair< int, int > > vortexBeginEnd
Номера первого и последнего вихрей, рождаемых на каждой панели профиля (формируется после решения СЛА...
Definition: Boundary2D.h:82
const World2D & W
Константная ссылка на решаемую задачу
Definition: Boundary2D.h:67
size_t GetUnknownsSize() const
Возврат размерности вектора решения
Definition: Boundary2D.cpp:71
const Airfoil & afl
Definition: Boundary2D.h:76