VM2D 1.14
Vortex methods for 2D flows simulation
Loading...
Searching...
No Matches
Velocity2D.h
Go to the documentation of this file.
1/*--------------------------------*- VM2D -*-----------------*---------------*\
2| ## ## ## ## #### ##### | | Version 1.14 |
3| ## ## ### ### ## ## ## ## | VM2D: Vortex Method | 2026/03/06 |
4| ## ## ## # ## ## ## ## | for 2D Flow Simulation *----------------*
5| #### ## ## ## ## ## | Open Source Code |
6| ## ## ## ###### ##### | https://www.github.com/vortexmethods/VM2D |
7| |
8| Copyright (C) 2017-2026 I. Marchevsky, K. Sokol, E. Ryatina, A. Kolganova |
9*-----------------------------------------------------------------------------*
10| File name: Velocity2D.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 VELOCITY2D_H
41#define VELOCITY2D_H
42
43#include <memory>
44
45#include "defs.h"
46#include "cudaTreeInfo.h"
47#include "OptimizedVelocity2D.h"
48
49namespace VM2D
50{
51
52 class Airfoil;
53 class Boundary;
54 class WakeDataBase;
55 class World2D;
56
57
70 {
72 std::vector<Point2D> convVelo;
73
75 std::vector<Point2D> diffVelo;
76
78 std::vector<double> I0;
79
81 std::vector<double> I1;
82
84 std::vector<Point2D> I2;
85
87 std::vector<Point2D> I3;
88
90 std::vector<double> epsastWake;
91 };
92
105 {
106 protected:
108 const World2D& W;
109
110 public:
113
115 std::vector<VortexesParams> virtualVortexesParams;
116
118
122 Velocity(const World2D& W_) :
123 W(W_)
124 {
125 virtualVortexesParams.resize(0);
126 };
127
131 void CalcConvVelo();
132
140 virtual void CalcConvVeloToSetOfPointsFromWake(const WakeDataBase& pointsDb, std::vector<Point2D>& velo, std::vector<double>& domainRadius, bool calcVelo, bool calcRadius) = 0;
141 virtual void CalcConvVPVeloToSetOfPointsFromWake(const WakeDataBase& pointsDb, std::vector<Point2D>& velo, std::vector<double>& domainRadius, bool calcVelo, bool calcRadius) {};
142
143#if defined (USE_CUDA)
144 virtual void GPUCalcConvVeloToSetOfPointsFromWake(std::unique_ptr<BHcu::CudaTreeInfo>& cntrTree, const WakeDataBase& pointsDb, std::vector<Point2D>& velo, std::vector<double>& domainRadius, bool calcVelo, bool calcRadius) = 0;
145 virtual void GPUCalcConvVelocityToSetOfPointsFromSheets(std::unique_ptr<BHcu::CudaTreeInfo>& cntrTree, const WakeDataBase& pointsDb, std::vector<Point2D>& velo) const = 0;
146#endif
147
149 virtual void CalcVeloToWakeVP();
150
158 void CalcDiffVeloI1I2ToSetOfPointsFromWake(const WakeDataBase& pointsDb, const std::vector<double>& domainRadius, const WakeDataBase& vorticesDb, std::vector<double>& I1, std::vector<Point2D>& I2);
159 void CalcDiffVeloI1I2ToSetOfPointsFromSheets(const WakeDataBase& pointsDb, const std::vector<double>& domainRadius, const Boundary& bnd, std::vector<double>& I1, std::vector<Point2D>& I2);
160
161 void CalcDiffVeloI1I2ToWakeFromSheets(const WakeDataBase& pointsDb, const std::vector<double>& domainRadius, const Boundary& bnd, std::vector<double>& I1, std::vector<Point2D>& I2);
162 void CalcDiffVeloI1I2ToWakeFromWake(const WakeDataBase& pointsDb, const std::vector<double>& domainRadius, const WakeDataBase& vorticesDb, std::vector<double>& I1, std::vector<Point2D>& I2);
163#if defined(USE_CUDA)
164 void GPUCalcDiffVeloI1I2ToSetOfPointsFromWake(const WakeDataBase& pointsDb, const std::vector<double>& domainRadius, const WakeDataBase& vorticesDb, std::vector<double>& I1, std::vector<Point2D>& I2, bool useMesh = false);
165 void GPUCalcDiffVeloI1I2ToSetOfPointsFromSheets(const WakeDataBase& pointsDb, const std::vector<double>& domainRadius, const Boundary& bnd, std::vector<double>& I1, std::vector<Point2D>& I2, bool useMesh = false);
166 void GPUDiffVeloFAST(const std::vector<double>& domainRadius, std::vector<double>& I1, std::vector<Point2D>& I2);
167#endif
168
174 void CalcDiffVeloI1I2();
175 void CalcDiffVeloI0I3();
176
177
181 void LimitDiffVelo(std::vector<Point2D>& diffVel);
182
186 void CalcDiffVelo();
187
191 void ResizeAndZero();
192
196 void SaveVisStress();
197
198
205 void GetWakeInfluenceToRhs(const Airfoil& afl, std::vector<double>& wakeRhs) const;
206#if defined(USE_CUDA)
207 void GPUGetWakeInfluenceToRhs(const Airfoil& afl, std::vector<double>& wakeVelo) const;
208 void GPUFASTGetWakeInfluenceToRhs(const Airfoil& afl, std::vector<double>& wakeVelo) const;
209#endif
210 void FillRhs(/*Eigen::VectorXd& rhs, Eigen::VectorXd& rhsSkos,*/ Eigen::VectorXd& rhsReord) const;
211
212
214 virtual ~Velocity() { };
215 };
216
217}//namespace VM2D
218
219#endif
Абстрактный класс, определяющий обтекаемый профиль
Definition Airfoil2D.h:182
Абстрактный класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле
Definition Boundary2D.h:65
Абстрактный класс, определяющий способ вычисления скоростей
Definition Velocity2D.h:105
void CalcConvVelo()
Вычисление конвективных скоростей вихрей и виртуальных вихрей в вихревом следе, а также в точках wake...
void CalcDiffVelo()
Вычисление диффузионных скоростей
void LimitDiffVelo(std::vector< Point2D > &diffVel)
Контроль больших значений диффузионных скоростей
void GetWakeInfluenceToRhs(const Airfoil &afl, std::vector< double > &wakeRhs) const
Генерация вектора влияния вихревого следа на профиль
virtual void CalcConvVPVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius)
Definition Velocity2D.h:141
void CalcDiffVeloI1I2ToWakeFromSheets(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)
void CalcDiffVeloI1I2ToWakeFromWake(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition Velocity2D.h:115
virtual void CalcConvVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius)=0
Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition Velocity2D.h:112
virtual void CalcVeloToWakeVP()
Вычисление скоростей в точках wakeVP.
void FillRhs(Eigen::VectorXd &rhsReord) const
void CalcDiffVeloI0I3()
void CalcDiffVeloI1I2ToSetOfPointsFromWake(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек
OptimizedVelocity optimizedVelocity
Definition Velocity2D.h:117
void CalcDiffVeloI1I2ToSetOfPointsFromSheets(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)
void ResizeAndZero()
Очистка старых массивов под хранение скоростей, выделение новой памяти и обнуление
void SaveVisStress()
Сохранение вязких напряжений
Velocity(const World2D &W_)
Конструктор
Definition Velocity2D.h:122
virtual ~Velocity()
Деструктор
Definition Velocity2D.h:214
void CalcDiffVeloI1I2()
Вычисление диффузионных скоростей вихрей и виртуальных вихрей в вихревом следе
const World2D & W
Константная ссылка на решаемую задачу
Definition Velocity2D.h:108
Класс, опеделяющий набор вихрей
Класс, опеделяющий текущую решаемую задачу
Definition World2D.h:74
Описание базовых вспомогательных функций
Структура, определяющая параметры виртуальных вихрей для отдельного профиля
Definition Velocity2D.h:70
std::vector< Point2D > convVelo
Вектор конвективных скоростей вихрей
Definition Velocity2D.h:72
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition Velocity2D.h:90
std::vector< Point2D > I2
Вектор числителей (I2) диффузионных скоростей вихрей (обусловленных завихренностью)
Definition Velocity2D.h:84
std::vector< double > I0
Вектор знаменателей (I0) диффузионных скоростей вихрей (обусловленных профилем)
Definition Velocity2D.h:78
std::vector< Point2D > diffVelo
Вектор диффузионных скоростей вихрей
Definition Velocity2D.h:75
std::vector< Point2D > I3
Вектор числителей (I3) диффузионных скоростей вихрей (обусловленных профилем)
Definition Velocity2D.h:87
std::vector< double > I1
Вектор знаменателей (I1) диффузионных скоростей вихрей (обусловленных завихренностью)
Definition Velocity2D.h:81