VM2D 1.14
Vortex methods for 2D flows simulation
Loading...
Searching...
No Matches
TreeBH.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: TreeBH.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
39#ifndef TREEBH_H
40#define TREEBH_H
41
42
43#include <memory>
44
45#include "defsBH.h"
46
47
48namespace BH
49{
51 extern long long op;
52
65 {
67 unsigned int key;
68
71
74
76 double g;
77 };
78
79
91 struct treeCellT
92 {
94 int parent;
95
98
101
104 bool particle;
105
107 bool leaf;
108
111
114
116 int level;
117
119 std::vector<Point2D> mom;
120
122 std::vector<Point2D> E;
123
127 std::vector<int> closeCells;
128
129 //габариты ячейки
131
132 //посчитан или нет габаритный прямоугольник
134 };
135
136
147 {
148 private:
150 const params& prm;
151
153 std::vector<PointsCopy>& pointsCopy;
154
155 public:
157 std::vector<TParticleCode> mortonCodes;
158
159 private:
166
169
175 int Delta(int i, int j) const;
176
180 int PrefixLength(int cell) const;
181
185 std::pair<unsigned int, int> Prefix(int cell) const;
186
188 void SetCellGeometry(int cell);
189
191 void BuildInternalTreeCell(int i);
192
194 static unsigned int ExpandBits(unsigned int v);
195
198 static unsigned int Morton2D(const Point2D& r);
199
201 static std::pair<Point2D, Point2D> FindEnclosingRectangleOld(const std::vector<PointsCopy>& points);
202
203 //Поиск габаритного прямоугольника системы точек
204 //работает под OpenMP 4.0
205 //static std::pair<Point2D, Point2D> FindEnclosingRectangle(const std::vector<PointsCopy>& points);
206
208 void RSort_Node3(TParticleCode* m, TParticleCode* m_temp, unsigned int n);
209
211 inline void RSort_step3(TParticleCode* source, TParticleCode* dest, unsigned int n, unsigned int* offset, unsigned char sortable_bit);
212
214 inline void RSort_Parallel(TParticleCode* m, TParticleCode* m_temp, unsigned int n, unsigned int* s);
215
216 public:
217
219 std::vector<double> iFact;
220
223 std::vector<double> binomCft;
224
226 const int maxTreeLevel;
227
234 void getStatistics(int& numParticles, int& treshold, std::pair<int, int>& partLevel, std::pair<int, int>& leafsLevel, int& lowLevelCells) const;
235
238 double& getBinom(int n, int k)
239 {
240 return binomCft[n * (prm.order + 1) + k];
241 }
242
243 //Конструктор
244 MortonTree(const params& prm_, int maxTreeLevel_, std::vector<PointsCopy>& points);
245
247 std::vector<int> mortonLowCells;
248 //std::vector<std::vector<int>> mortonLowCells;
249
257 std::vector<treeCellT> mortonTree;
258
260 void MakeRootMortonTree();
261
264
267
269 // Два варианта заполнения списка нижних вершин:
270 // первый рекурсивный, быстро работает в последовательном варианте
271 // второй нерекурсивный, хорошо распараллеливается, но неэффективен в последовательном варианте
272 void FillMortonLowCells(int cell = 0);
273 void FillMortonLowCellsA();
274
278 void CalculateMortonTreeParams(int cell, int omplevel);
279
285 void CalcLocalCoeffToLowLevel(int lowCell, std::unique_ptr<MortonTree>& treeInf, int fromWho, bool calcCloseTrees);
286
290 void CalcVeloBiotSavart(int lowCell, std::unique_ptr<MortonTree>& treeInf);//, bool calcRadius);
291
292
295 void CalcVeloTaylorExpansion(int lowCell);
296
299 void AddVelo(int lowCell);
300
301
302 //спускается до листьев рекурсивно и запускает расчет габаритного прямоугольника
303 void MinMaxCellsSearch(int numCell);
304
305 //расчет габаритного прямоугольника
306 void MinMaxCellsCalc(BH::treeCellT& cell);
307
308 //расчет габаритных прямоугольников при подъеме от листов к корню
309 void GabsPyr(int numCell);
310
311 //
312 void ShareGabs(int numCell/*, Point2D& minx_miny, Point2D& maxx_maxy*/);
313
314 };
315
316}//namespace BH
317
318
319#endif
Структура, соответствующая мортоновскому дереву
Definition TreeBH.h:147
void BuildMortonParticlesTree()
Построение верхушек дерева — отдельных частиц
Definition TreeBH.cpp:330
void MinMaxCellsSearch(int numCell)
Definition TreeBH.cpp:809
std::vector< treeCellT > mortonTree
Мортоновское дерево Размер дерева равен сумме числа внутренних вершин (число частиц без единицы) и чи...
Definition TreeBH.h:257
void BuildInternalTreeCell(int i)
Функция построения i-й внутренней ячейки дерева
Definition TreeBH.cpp:264
int Delta(int i, int j) const
Definition TreeBH.cpp:187
std::vector< PointsCopy > & pointsCopy
Ссылка на список обобщенных вихрей, по которым строится дерево
Definition TreeBH.h:153
void CalcLocalCoeffToLowLevel(int lowCell, std::unique_ptr< MortonTree > &treeInf, int fromWho, bool calcCloseTrees)
Definition TreeBH.cpp:436
void MinMaxCellsCalc(BH::treeCellT &cell)
Definition TreeBH.cpp:829
void RSort_step3(TParticleCode *source, TParticleCode *dest, unsigned int n, unsigned int *offset, unsigned char sortable_bit)
Вспомогательная функция к RSort_Node3.
Definition TreeBH.cpp:618
void RSort_Node3(TParticleCode *m, TParticleCode *m_temp, unsigned int n)
Поразрядная сортировка мортоновских кодов
Definition TreeBH.cpp:630
void FillMortonLowCells(int cell=0)
Функция заполнения списка нижних ячеек дерева (до которых производится обход)
Definition TreeBH.cpp:348
void CalcVeloTaylorExpansion(int lowCell)
Definition TreeBH.cpp:536
void ShareGabs(int numCell)
Definition TreeBH.cpp:854
void BuildMortonInternalTree()
Построение внутренних ячеек дерева
Definition TreeBH.cpp:314
std::pair< unsigned int, int > Prefix(int cell) const
Definition TreeBH.cpp:227
static std::pair< Point2D, Point2D > FindEnclosingRectangleOld(const std::vector< PointsCopy > &points)
Поиск габаритного прямоугольника системы точек (для OpenMP 2.0)
Definition TreeBH.cpp:82
double iQuadSideVar
Definition TreeBH.h:165
void MakeRootMortonTree()
Построение корня дерева и задание его общих параметров
Definition TreeBH.cpp:140
static unsigned int Morton2D(const Point2D &r)
Definition TreeBH.cpp:591
std::vector< int > mortonLowCells
Вектор индесков ячеек нижнего уровня в дереве
Definition TreeBH.h:247
static unsigned int ExpandBits(unsigned int v)
"Разрежение" двоичного представления беззнакового целого, вставляя по одному нулику между всеми битам...
Definition TreeBH.cpp:559
void FillMortonLowCellsA()
Definition TreeBH.cpp:362
Point2D lowLeft
Положение нижнего левого угла дерева
Definition TreeBH.h:168
void RSort_Parallel(TParticleCode *m, TParticleCode *m_temp, unsigned int n, unsigned int *s)
Параллельная версия поразрядной сортировки мортоновских кодов
Definition TreeBH.cpp:677
void getStatistics(int &numParticles, int &treshold, std::pair< int, int > &partLevel, std::pair< int, int > &leafsLevel, int &lowLevelCells) const
Definition TreeBH.cpp:749
const params & prm
Ссылка на параметры, загружаемые из файла
Definition TreeBH.h:150
std::vector< TParticleCode > mortonCodes
Список мортоновских кодов, вычисляемых по положениям частиц
Definition TreeBH.h:157
double & getBinom(int n, int k)
Definition TreeBH.h:238
int PrefixLength(int cell) const
Definition TreeBH.cpp:219
void AddVelo(int lowCell)
void GabsPyr(int numCell)
Definition TreeBH.cpp:870
void CalcVeloBiotSavart(int lowCell, std::unique_ptr< MortonTree > &treeInf)
Definition TreeBH.cpp:499
std::vector< double > binomCft
Definition TreeBH.h:223
const int maxTreeLevel
Максимальная глубина дерева при его обходе
Definition TreeBH.h:226
std::vector< double > iFact
Спикок обратных факториалов целых чисел
Definition TreeBH.h:219
void CalculateMortonTreeParams(int cell, int omplevel)
Вычисление параметров дерева (циркуляций и высших моментов)
Definition TreeBH.cpp:378
void SetCellGeometry(int cell)
Функция вычисления геометрических параметров внутренней ячейки дерева
Definition TreeBH.cpp:236
Класс, содержащий параметры метода Баонса - Хата для CPU.
Definition ParamsBH.h:65
int order
Точность расчета скоростей:
Definition ParamsBH.h:78
Шаблонный класс, определяющий вектор фиксированной длины Фактически представляет собой массив,...
Definition numvector.h:99
Вспомогательные функции для метода Барнса - Хата для CPU.
long long op
Глобальная переменная - счетчик количества операций
Definition TreeBH.cpp:46
Структура, соответствующая частице и ее мортоновскому коду
Definition TreeBH.h:65
double g
Циркуляция частицы (заносится в эту структуру для облегчения доступа)
Definition TreeBH.h:76
Point2D r
Положение частицы (заносится в эту структуру для облегчения доступа)
Definition TreeBH.h:73
unsigned int key
Мортоновский код частицы
Definition TreeBH.h:67
int originNumber
Индекс частицы в глобальном массиве частиц
Definition TreeBH.h:70
Структура, соответствующая ячейке мортоновского дерева
Definition TreeBH.h:92
bool leaf
Признак того, что при обходе дерева на более глубокие уровни спускаться не нужно
Definition TreeBH.h:107
numvector< int, 2 > range
Диапазон индексов частиц в ячейке: 0 - первая; 1 - последняя
Definition TreeBH.h:100
Point2D minx_miny
Definition TreeBH.h:130
int level
Уровень ячейки в дереве (не заполняется для частиц)
Definition TreeBH.h:116
Point2D center
Положение центра ячейки
Definition TreeBH.h:113
Point2D maxx_maxy
Definition TreeBH.h:130
Point2D size
Размер ячейки по x и по y (для листа (0;0)
Definition TreeBH.h:110
bool hasGabs
Definition TreeBH.h:133
std::vector< Point2D > E
Коэффициенты локального разложения в ячейке для вычисления скоростей
Definition TreeBH.h:122
std::vector< Point2D > mom
Мультипольные моменты ячейки
Definition TreeBH.h:119
std::vector< int > closeCells
Definition TreeBH.h:127
numvector< int, 2 > child
Индексы потомков (если это не particle): 0 - левый, 1 - правый
Definition TreeBH.h:97
int parent
Индекс родительской ячейки (0 — для корня дерева)
Definition TreeBH.h:94
bool particle
Definition TreeBH.h:104