VM2D 1.14
Vortex methods for 2D flows simulation
Loading...
Searching...
No Matches
Velocity2DBarnesHut.cpp
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: Velocity2DBarnesHut.cpp |
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#include "Velocity2DBarnesHut.h"
41
42#include "Airfoil2D.h"
43#include "Boundary2D.h"
44#include "MeasureVP2D.h"
45#include "Mechanics2D.h"
46#include "StreamParser.h"
47#include "Wake2D.h"
48#include "World2D.h"
49
50#include "BarnesHut.h"
51
52using namespace VM2D;
53
59
64
65
66void VelocityBarnesHut::CalcConvVeloToSetOfPointsFromWake(const WakeDataBase& pointsDb, std::vector<Point2D>& velo, std::vector<double>& domainRadius, bool calcVelo, bool calcRadius)
67{
68 BH::params prm;
71
72 int nCores = omp_get_max_threads();
75
76 //Выбор неочевиден
77 //&&&
78 prm.NumOfLevelsVortex = (int)log2(pointsDb.vtx.size());
79
80 BH::BarnesHut BH(prm, pointsDb.vtx);
81 double timeBuildTree = 0.0, timeRect = 0.0, timeSummarization = 0.0, timeComputation = 0.0;
82 BH.BuildOneTree(BH.treeVrt, prm.NumOfLevelsVortex, BH.pointsCopyVrt, timeBuildTree);
83 BH.BuildEnclosingRectangle(timeRect);
84 BH.InfluenceComputation(velo, domainRadius, timeSummarization, timeComputation, calcRadius);
85 //std::cout << timeBuildTree << " " << timeRect << " " << timeSummarization << " " << timeComputation << std::endl;
86 //std::cout << "Tot_time = " << timeBuildTree + timeRect + timeSummarization + timeComputation << std::endl;
87
88#ifdef USE_CUDA
89 if (calcVelo)
90 W.getCuda().CopyMemToDev<double, 2>(velo.size(), (double*)velo.data(), pointsDb.devVelPtr);
91 if (calcRadius)
92 W.getCuda().CopyMemToDev<double, 1>(domainRadius.size(), domainRadius.data(), pointsDb.devRadPtr);
93#endif
94
95
96}
97
98
99
100//Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа
101void VelocityBarnesHut::CalcConvVPVeloToSetOfPointsFromWake(const WakeDataBase& pointsDb, std::vector<Point2D>& velo, std::vector<double>& domainRadius, bool calcVelo, bool calcRadius)
102{
103 std::vector<Point2D> selfVelo(pointsDb.vtx.size());
104 domainRadius.resize(pointsDb.vtx.size());
105
106 double cft = IDPI;
107
108#pragma warning (push)
109#pragma warning (disable: 4101)
110 //Локальные переменные для цикла
111 Point2D velI;
112 Point2D tempVel;
113 double dst2eps, dst2;
114#pragma warning (pop)
115
117
118 if (calcVelo)
119 {
120#pragma omp parallel for default(none) shared(selfVelo, cft, calcVelo, calcRadius, eps2, pointsDb, domainRadius) private(tempVel, velI, dst2, dst2eps) schedule(dynamic, DYN_SCHEDULE)
121 for (int i = 0; i < pointsDb.vtx.size(); ++i)
122 {
123 double ee2[3] = { 10000.0, 10000.0, 10000.0 };
124
125 velI.toZero();
126
127 const Point2D& posI = pointsDb.vtx[i].r();
128
129 for (size_t j = 0; j < W.getWake().vtx.size(); ++j)
130 {
131 const Point2D& posJ = W.getWake().vtx[j].r();
132
133 dst2 = (posI - posJ).length2();
134
135// //Модифицируем массив квадратов расстояний до ближайших вихрей из wake
136//#ifndef TESTONLYVELO
137// if (calcRadius)
138// VMlib::ModifyE2(ee2, dst2);
139//#endif //!TESTONLYVELO
140
141 const double& gamJ = W.getWake().vtx[j].g();
142
143 tempVel.toZero();
144 dst2eps = VMlib::boundDenom(dst2, eps2); //Сглаживать надо!!!
145
146 tempVel = { -posI[1] + posJ[1], posI[0] - posJ[0] };
147 tempVel *= (gamJ / dst2eps);
148 velI += tempVel;
149 }
150
151 for (size_t j = 0; j < W.getSource().vtx.size(); ++j)
152 {
153 const Point2D& posJ = W.getSource().vtx[j].r();
154 const double& gamJ = W.getPassport().physicalProperties.accelCft(W.getCurrentTime()) * W.getSource().vtx[j].g();
155
156 tempVel.toZero();
157
158 dst2 = dist2(posI, posJ);
159 dst2eps = VMlib::boundDenom(dst2, W.getPassport().wakeDiscretizationProperties.eps2); //Сглаживать надо!!!
160
161 tempVel = { posI[0] - posJ[0], posI[1] - posJ[1] };
162 tempVel *= (gamJ / dst2eps);
163 velI += tempVel;
164 }
165//#ifndef TESTONLYVELO
166// if (calcRadius)
167// {
168// for (size_t s = 0; s < W.getNumberOfBoundary(); ++s)
169// {
170// const auto& bou = W.getBoundary(s);
171// //Модифицируем массив квадратов расстояний до ближайших вихрей из virtualWake
172// for (size_t j = 0; j < bou.virtualWake.vtx.size(); ++j)
173// ModifyE2(ee2, dist2(posI, bou.virtualWake.vtx[j].r()));
174// }
175// }
176//#endif
177
178 velI *= cft;
179 selfVelo[i] = velI;
180
181//#ifndef TESTONLYVELO
182// if (calcRadius)
183// domainRadius[i] = 1.0 * sqrt((ee2[0] + ee2[1] + ee2[2]) / 3.0);
184//#endif
185 }
186 }
187// else if (calcRadius)
188// {
189//#pragma omp parallel for default(none) shared(selfVelo, cft, calcVelo, calcRadius, pointsDb, domainRadius) private(tempVel, velI, dst2) schedule(dynamic, DYN_SCHEDULE)
190// for (int i = 0; i < pointsDb.vtx.size(); ++i)
191// {
192// double ee2[3] = { 10000.0, 10000.0, 10000.0 };
193//
194// velI.toZero();
195//
196// const Point2D& posI = pointsDb.vtx[i].r();
197//
198// for (size_t j = 0; j < W.getWake().vtx.size(); ++j)
199// {
200// const Point2D& posJ = W.getWake().vtx[j].r();
201//
202// dst2 = dist2(posI, posJ);
203//
204// //Модифицируем массив квадратов расстояний до ближайших вихрей из wake
205// VMlib::ModifyE2(ee2, dst2);
206// }
207//
208// for (size_t s = 0; s < W.getNumberOfBoundary(); ++s)
209// {
210// for (size_t j = 0; j < W.getBoundary(s).virtualWake.vtx.size(); ++j)
211// {
212// const Point2D& posJ = W.getBoundary(s).virtualWake.vtx[j].r();
213// dst2 = dist2(posI, posJ);
214//
215// //Модифицируем массив квадратов расстояний до ближайших вихрей из virtualWake
216// ModifyE2(ee2, dst2);
217// }
218// }
219//
220// domainRadius[i] = 1.0 * sqrt((ee2[0] + ee2[1] + ee2[2]) / 3.0);
221// }
222// } //else
223
224
225 if (calcVelo)
226 for (size_t i = 0; i < velo.size(); ++i)
227 velo[i] += selfVelo[i];
228
229#ifdef USE_CUDA
230 if (calcVelo)
231 W.getCuda().CopyMemToDev<double, 2>(velo.size(), (double*)velo.data(), pointsDb.devVelPtr);
232// if (calcRadius)
233// W.getCuda().CopyMemToDev<double, 1>(domainRadius.size(), domainRadius.data(), pointsDb.devRadPtr);
234#endif
235
236}//CalcConvVPVeloToSetOfPointsFromWake(...)
237
238
239
240
241#if defined(USE_CUDA)
242void VelocityBarnesHut::GPUCalcConvVeloToSetOfPointsFromWake(std::unique_ptr<BHcu::CudaTreeInfo>& cntrTree, const WakeDataBase& pointsDb, std::vector<Point2D>& velo, std::vector<double>& domainRadius, bool calcVelo, bool calcRadius)
243{
244 const double& eps2 = W.getPassport().wakeDiscretizationProperties.eps2;
245 const double& theta = W.getPassport().numericalSchemes.nbodyTheta;
247
248 int npt = cntrTree->nObject;
249
250 float tBLD = 0.0f, tUPW = 0.0f, tDNV;
251
252 auto& inflTree = *W.getCuda().inflTreeWake;
253
254 tDNV = inflTree.DownwardTraversalVorticesToPoints(*cntrTree, (Point2D*)pointsDb.devVelPtr, pointsDb.devRadPtr, eps2, theta, order, calcRadius);
255
256 //std::cout << "nvt = " << npt << std::endl;
257 //std::cout << "tBLD = " << tBLD << std::endl;
258 //std::cout << "tUPW = " << tUPW << std::endl;
259 //std::cout << "tDNV = " << tDNV << std::endl;
260
261 std::vector<Point2D> Vel;
262 if (calcVelo)
263 {
264 Vel.resize(npt);
265 W.getCuda().CopyMemFromDev<double, 2>(npt, (double*)pointsDb.devVelPtr, (double*)Vel.data(), 20);
266
267 for (size_t q = 0; q < npt; ++q)
268 velo[q] += Vel[q];
269 }
270
271 if (calcRadius)
272 W.getCuda().CopyMemFromDev<double, 1>(npt, pointsDb.devRadPtr, domainRadius.data(), 212);
273}
274
275
276void VelocityBarnesHut::GPUCalcConvVelocityToSetOfPointsFromSheets(std::unique_ptr<BHcu::CudaTreeInfo>& cntrTree, const WakeDataBase& pointsDb, std::vector<Point2D>& velo) const
277{
278 double*& velD = pointsDb.devVelPtr;
279
281 const double& theta = W.getPassport().numericalSchemes.nbodyTheta;
283
284 std::vector<Point2D> newV(velo.size());
285
286 int npt = cntrTree->nObject;
287
288 if (npt > 0)
289 {
290 int nTotPan = (int)W.getAirfoil(0).getNumberOfPanels();
291 for (size_t q = 1; q < W.getNumberOfAirfoil(); ++q)
292 nTotPan += (int)W.getAirfoil(q).getNumberOfPanels();
293 auto& afl = W.getAirfoil(0);
294
295 //Влияние вихревых слоев
296 W.getCuda().inflTreePnlVortex->DownwardTraversalPanelsToPoints(*cntrTree, (Point2D*)velD, eps2, theta, order);
297 W.getCuda().CopyMemFromDev<double, 2>(npt, velD, (double*)newV.data());
298
299 for (size_t q = 0; q < velo.size(); ++q)
300 velo[q] += newV[q];
301
303 {
304 //Влияние слоев источников
305 W.getCuda().inflTreePnlSource->DownwardTraversalPanelsToPoints(*cntrTree, (Point2D*)velD, eps2, theta, order);
306 W.getCuda().CopyMemFromDev<double, 2>(npt, velD, (double*)newV.data());
307
308 for (size_t q = 0; q < velo.size(); ++q)
309 velo[q] += newV[q];
310 }
311 }
312}//GPUCalcConvVelocityToSetOfPointsFromSheets(...)
313#endif
314
315
316
317
318
319
320
321
322
Заголовочный файл с описанием класса Airfoil.
Заголовок основного класса BarnesHut.
Заголовочный файл с описанием класса Boundary.
Заголовочный файл с описанием класса MeasureVP.
Заголовочный файл с описанием класса Mechanics.
Заголовочный файл с описанием класса StreamParser.
Заголовочный файл с описанием класса VelocityBarnesHut.
Заголовочный файл с описанием класса Wake.
Заголовочный файл с описанием класса World2D.
Класс, определяющий основной алгоритм модификации метода Барнса - Хата
Definition BarnesHut.h:59
Класс, содержащий параметры метода Баонса - Хата для CPU.
Definition ParamsBH.h:65
double eps
Радиус вихревого элемента
Definition ParamsBH.h:68
int NumOfLevelsVortex
Максимальное количество уровней дерева
Definition ParamsBH.h:81
double theta
Параметр точности
Definition ParamsBH.h:84
int order
Точность расчета скоростей:
Definition ParamsBH.h:78
double eps2
Квадрат радиуса вихревого элемента
Definition ParamsBH.h:71
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition Airfoil2D.h:163
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition Passport2D.h:289
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition Passport2D.h:292
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition Passport2D.h:295
virtual void CalcConvVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius) override
Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа
VelocityBarnesHut(const World2D &W_)
Конструктор
virtual ~VelocityBarnesHut()
Деструктор
virtual void CalcConvVPVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius) override
Абстрактный класс, определяющий способ вычисления скоростей
Definition Velocity2D.h:105
const World2D & W
Константная ссылка на решаемую задачу
Definition Velocity2D.h:108
Класс, опеделяющий набор вихрей
std::vector< Vortex2D > vtx
Список вихревых элементов
Класс, опеделяющий текущую решаемую задачу
Definition World2D.h:74
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition World2D.h:174
bool isAnyMovableOrDeformable() const
Возврат признака того, что хотя бы один из профилей подвижный или деформируемый
Definition World2D.cpp:1765
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition World2D.h:226
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition World2D.h:157
const WakeDataBase & getSource() const
Возврат константной ссылки на источники в области течения
Definition World2D.h:236
const Gpu & getCuda() const
Возврат константной ссылки на объект, связанный с видеокартой (GPU)
Definition World2D.h:261
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition World2D.h:251
double getCurrentTime() const
Definition WorldGen.h:100
numvector< T, n > & toZero(P val=0)
Установка всех компонент вектора в константу (по умолчанию — нуль)
Definition numvector.h:527
const double IDPI
Число .
Definition defs.h:85
double boundDenom(double r2, double eps2)
Способ сглаживания скорости вихря (вихрь Рэнкина или вихрь Ламба)
Definition defs.h:563
double accelCft(double currentTime) const
Функция-множитель, позволяющая моделировать разгон
double eps
Радиус вихря
Definition Passport2D.h:126
double eps2
Квадрат радиуса вихря
Definition Passport2D.h:129