VM2D 1.14
Vortex methods for 2D flows simulation
Loading...
Searching...
No Matches
BH::BarnesHut Class Reference

Класс, определяющий основной алгоритм модификации метода Барнса - Хата More...

#include <BarnesHut.h>

Collaboration diagram for BH::BarnesHut:

Public Member Functions

 BarnesHut (const params &prm_, const std::vector< Vortex2D > &pointsVrt)
 Конструктор для решения задачи NBODY о вычислении скоростей вихревых частиц
 
 ~BarnesHut ()
 Деструктор
 
void BuildOneTree (std::unique_ptr< MortonTree > &tree, int maxTreeLevel, std::vector< PointsCopy > &pointsCopy, double &time)
 Построение одного дерева tree на основе заданных точек pointsCopy

 
void BuildNecessaryTrees (double &time)
 Построение всех нужных деревьев на основе заданных точек pointsCopy

 
void BuildEnclosingRectangle (double &time)
 
void InfluenceComputation (std::vector< Point2D > &result, std::vector< double > &epsast, double &timeParams, double &timeInfl, bool calcRadius)
 Расчет влияния в точках дерева, характерных для решаемой задачи (определяется внутри функции)
 

Public Attributes

const paramsprm
 Ссылка на параметры, считываемые из файла
 
std::vector< PointsCopypointsCopyVrt
 Список оберток положений вихрей
 
std::unique_ptr< MortonTreetreeVrt
 Умный yказатель на дерево вихрей
 

Detailed Description

Класс, определяющий основной алгоритм модификации метода Барнса - Хата

Author
Марчевский Илья Константинович
Сокол Ксения Сергеевна
Рятина Евгения Павловна
Колганова Александра Олеговна

\Version 1.14

Date
6 марта 2026 г.

Definition at line 58 of file BarnesHut.h.

Constructor & Destructor Documentation

◆ BarnesHut()

BH::BarnesHut::BarnesHut ( const params prm_,
const std::vector< Vortex2D > &  pointsVrt 
)

Конструктор для решения задачи NBODY о вычислении скоростей вихревых частиц

Parameters
[in]prmконстантная ссылка на параметры, считываемые из файла
[in]pointsVrtконстантная ссылка на список частиц

Definition at line 47 of file BarnesHut.cpp.

48 : prm(prm_)
49 {
50 pointsCopyVrt.insert(pointsCopyVrt.end(), pointsVrt.begin(), pointsVrt.end());
51 }//BarnesHut(...)
const params & prm
Ссылка на параметры, считываемые из файла
Definition BarnesHut.h:62
std::vector< PointsCopy > pointsCopyVrt
Список оберток положений вихрей
Definition BarnesHut.h:65

◆ ~BarnesHut()

BH::BarnesHut::~BarnesHut ( )
inline

Деструктор

Definition at line 77 of file BarnesHut.h.

77{};

Member Function Documentation

◆ BuildEnclosingRectangle()

void BH::BarnesHut::BuildEnclosingRectangle ( double &  time)

Definition at line 84 of file BarnesHut.cpp.

85 {
86 double t1 = omp_get_wtime();
87
88#ifdef needTreeVrt
89 //treeVrt->MinMaxCellsSearch(0); //габаритные прямоугольники для листовых ячеек
90
91 for (auto& cell : treeVrt->mortonLowCells)
92 treeVrt->MinMaxCellsSearch(cell);
93
94 treeVrt->GabsPyr(0); //расчет габаритных прямоугольников при подъеме от листов к корню (но на самом деле наоборот)
95#endif
96
97 double t2 = omp_get_wtime();
98 time += t2 - t1;
99
100
101
102 }
std::unique_ptr< MortonTree > treeVrt
Умный yказатель на дерево вихрей
Definition BarnesHut.h:68

◆ BuildNecessaryTrees()

void BH::BarnesHut::BuildNecessaryTrees ( double &  time)

Построение всех нужных деревьев на основе заданных точек pointsCopy

Parameters
[in,out]timeвремя, затрачиваемое на построение дерева (накопительный итог)

Definition at line 76 of file BarnesHut.cpp.

77 {
78#ifdef needTreeVrt
79 BuildOneTree(treeVrt, prm.NumOfLevelsVortex, pointsCopyVrt, time); // дерево вихрей
80#endif
81 }
void BuildOneTree(std::unique_ptr< MortonTree > &tree, int maxTreeLevel, std::vector< PointsCopy > &pointsCopy, double &time)
Построение одного дерева tree на основе заданных точек pointsCopy
Definition BarnesHut.cpp:56
int NumOfLevelsVortex
Максимальное количество уровней дерева
Definition ParamsBH.h:81
Here is the call graph for this function:

◆ BuildOneTree()

void BH::BarnesHut::BuildOneTree ( std::unique_ptr< MortonTree > &  tree,
int  maxTreeLevel,
std::vector< PointsCopy > &  pointsCopy,
double &  time 
)

Построение одного дерева tree на основе заданных точек pointsCopy

Parameters
[in,out]treeумный указатель на дерево
[in]maxTreeLevelмаксимальная глубина при обходе дерева
[in]pointsCopyнеконстантная ссылка на список данных (оберток) точек, по которым строится дерево
[in,out]timeвремя, затрачиваемое на построение дерева (накопительный итог)

Definition at line 56 of file BarnesHut.cpp.

57 {
58 tree = std::make_unique<MortonTree>(prm, maxTreeLevel, pointsCopy);
59
60 double t1 = omp_get_wtime();
61 tree->MakeRootMortonTree();
62 tree->BuildMortonInternalTree();
63 tree->BuildMortonParticlesTree();
64
65 // Один из двух вариантов заполнения нижних вершин
66 //tree->FillMortonLowCellsA(); // Рекурсивный
67 tree->FillMortonLowCells(); // Нерекурсивный
68
69 double t2 = omp_get_wtime();
70
71 time += t2-t1;
72 }
Here is the caller graph for this function:

◆ InfluenceComputation()

void BH::BarnesHut::InfluenceComputation ( std::vector< Point2D > &  result,
std::vector< double > &  epsast,
double &  timeParams,
double &  timeInfl,
bool  calcRadius 
)

Расчет влияния в точках дерева, характерных для решаемой задачи (определяется внутри функции)

Parameters
[out]resultссылка на вектор, в который сохраняются вычисленные скорости
[in,out]timeParamsвремя расчета параметров деревьев (накопительный итог)
[in,out]timeInflвремя расчета влияния (накопительный итог)

Definition at line 106 of file BarnesHut.cpp.

107 {
108 double tTreeParamsStart = omp_get_wtime();
109
110#ifdef OLD_OMP
111 omp_set_nested(1);
112#else
113 omp_set_max_active_levels(3);
114 std::cout << "OMP_LEVEL = 2" << std::endl;
115 //omp_set_max_active_levels(prm.maxLevelOmp + 1);
116#endif
117
118
119 auto& treeContr = treeVrt;
120 auto& pointsCopy = pointsCopyVrt;
121
122
123 treeVrt->CalculateMortonTreeParams(0, 0);
124
125 double tTreeParamsFinish = omp_get_wtime();
126 timeParams += tTreeParamsFinish - tTreeParamsStart;
127
128 double tInflStart = omp_get_wtime();
129
130#pragma omp parallel for schedule(dynamic, 10) //reduction(+:t1, t2, t3)
131 for (int i = 0; i < (int)treeContr->mortonLowCells.size(); ++i)
132 {
133 auto& lci = treeContr->mortonLowCells[i];
134 auto& lowCell = treeContr->mortonTree[lci];
135
136 for (auto& e : lowCell.E)
137 e.toZero();
138
139 treeContr->CalcLocalCoeffToLowLevel(lci, treeVrt, 0, true);
140 treeContr->CalcVeloBiotSavart(lci, treeVrt);
141 treeContr->CalcVeloTaylorExpansion(lci);
142 }
143
144
145 //CPU - neib
146 const size_t knb = 3;
147 std::vector<std::vector<std::pair<double, size_t>>> initdist(pointsCopy.size());
148 for (auto& d : initdist)
149 d.resize(2 * knb, { -1.0, -1 });
150 double timeKnn = -omp_get_wtime();
151 WakekNNnewForEpsast(pointsCopy, knb, initdist);//CPU
152
153#pragma omp parallel for
154 for (int j = 0; j < pointsCopy.size(); ++j)
155 {
156 double sd2 = (initdist[j][0].first + initdist[j][1].first + initdist[j][2].first) / 3.0;
157 pointsCopy[j].epsast = (sd2 > 0) ? sqrt(sd2) : 1000.0;
158 }
159 timeKnn += omp_get_wtime();
160
161 double tInflStop = omp_get_wtime();
162 timeInfl += tInflStop - tInflStart;
163
164 int n = (int)pointsCopy.size();
165
166#pragma omp parallel for
167 for (int i = 0; i < n; ++i)
168 result[i] = IDPI * pointsCopy[i].veloCopy;
169
170 if (calcRadius)
171 {
172#pragma omp parallel for
173 for (int i = 0; i < n; ++i)
174 epsast[i] = pointsCopy[i].epsast;
175 }
176 }//InfluenceComputation(...)
void WakekNNnewForEpsast(const std::vector< BH::PointsCopy > &vtx, const size_t k, std::vector< std::vector< std::pair< double, size_t > > > &initdist)
Definition knnCPU.cpp:458
static const double IDPI
Definition defsBH.h:60

Member Data Documentation

◆ pointsCopyVrt

std::vector<PointsCopy> BH::BarnesHut::pointsCopyVrt

Список оберток положений вихрей

Definition at line 65 of file BarnesHut.h.

◆ prm

const params& BH::BarnesHut::prm

Ссылка на параметры, считываемые из файла

Definition at line 62 of file BarnesHut.h.

◆ treeVrt

std::unique_ptr<MortonTree> BH::BarnesHut::treeVrt
mutable

Умный yказатель на дерево вихрей

Definition at line 68 of file BarnesHut.h.


The documentation for this class was generated from the following files: