VM2D  1.12
Vortex methods for 2D flows simulation
VM2D::Wake Class Reference

Класс, опеделяющий вихревой след (пелену) More...

#include <Wake2D.h>

Inheritance diagram for VM2D::Wake:
Collaboration diagram for VM2D::Wake:

Public Member Functions

 Wake (const World2D &W_)
 Конструктор инициализации More...
 
 ~Wake ()
 Деструктор More...
 
void Inside (const std::vector< Point2D > &newPos, Airfoil &afl, bool isMoves, const Airfoil &oldAfl)
 Проверка пересечения вихрями следа профиля при перемещении More...
 
void Restruct ()
 Реструктуризация вихревого следа More...
 
int RemoveFar ()
 Зануление далеко улетевших вихрей More...
 
size_t RemoveZero ()
 Исключение нулевых и мелких вихрей More...
 
bool MoveInside (const Point2D &newPos, const Point2D &oldPos, const Airfoil &afl, size_t &panThrough)
 Проверка проникновения точки через границу профиля More...
 
bool MoveInsideMovingBoundary (const Point2D &newPos, const Point2D &oldPos, const Airfoil &oldAfl, const Airfoil &afl, size_t &panThrough)
 Проверка проникновения точки через границу профиля More...
 
void GetPairs (int type)
 Поиск ближайшего соседа More...
 
void GetPairsBS (int type)
 
void GetPairsBH (int type)
 
int Collaps (int type, int times)
 Коллапс вихрей More...
 
int CollapsNew (int type, int times)
 
void ReadFromFile (const std::string &dir, const std::string &fileName)
 Считывание вихревого следа из файла More...
 
void SaveKadrVtk (const std::string &filePrefix="Kadr") const
 Сохранение вихревого следа в файл .vtk. More...
 

Public Attributes

double collapseRightBorderParameter
 абсцисса, правее которой происходит линейный (вправо) рост радиуса коллапса More...
 
double collapseScaleParameter
 характерный масштаб, на котором происходит рост радиуса коллапса More...
 
const World2DW
 Константная ссылка на решаемую задачу More...
 
std::vector< Vortex2Dvtx
 Список вихревых элементов More...
 

Private Attributes

std::vector< int > neighb
 Вектор потенциальных соседей для будущего коллапса More...
 
std::vector< int > neighbNew
 

Static Private Attributes

static const int knb = 2
 

Detailed Description

Класс, опеделяющий вихревой след (пелену)

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

Definition at line 59 of file Wake2D.h.

Constructor & Destructor Documentation

VM2D::Wake::Wake ( const World2D W_)
inline

Конструктор инициализации

Parameters
[in]W_константная ссылка на решаемую задачу

Definition at line 73 of file Wake2D.h.

74  : WakeDataBase(W_)
75  {
76  };
WakeDataBase(const World2D &W_)
Конструктор инициализации
VM2D::Wake::~Wake ( )
inline

Деструктор

Definition at line 79 of file Wake2D.h.

79 { };

Here is the call graph for this function:

Member Function Documentation

int Wake::Collaps ( int  type,
int  times 
)

Коллапс вихрей

Parameters
[in]typeтип коллапса:
  • 0 — без приоритета знаков
  • 1 — коллапсировать только вихри разных знаков
  • 2 — коллапсировать только вихри одного знака
[in]timesчисло проходов алгоритма коллапса
Returns
число зануленных вихрей

Definition at line 412 of file Wake2D.cpp.

413 {
414  int nHlop = 0; //общее число убитых вихрей
415 
416  //int loc_hlop = 0; //
417 
418  std::vector<bool> flag; //как только вихрь сколлапсировался flag = 1
419  for (int z = 0; z < times; ++z)
420  {
421 
422 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_PAIRS))
423 
425  {
426  const_cast<Gpu&>(W.getCuda()).RefreshWake(3);
427  GPUGetPairs(type);
428  }
429  else
430  GetPairs(type);
431 #else
432  GetPairs(type);
433 #endif
434 
435  //loc_hlop = 0;//число схлопнутых вихрей
436 
437  flag.clear();
438  flag.resize(vtx.size(), false);
439 
440  double sumAbsGam, iws;
441  Point2D newPos;
442 
443  for (size_t vt = 0; vt + 1 < vtx.size(); ++vt)
444  {
445  Vortex2D& vtxI = vtx[vt];
446 
447  if (!flag[vt])
448  {
449  //std::cout << "A" << std::endl;
450 
451  int ssd = neighb[vt];
452 
453  //std::cout << "B" << std::endl;
454 
455  if ((ssd < 0) || (ssd >= flag.size()))
456  std::cout << "ssd = " << ssd << ", flag.size() = " << flag.size() << ", vtx.size() = " << vtx.size() << std::endl;
457 
458  if ((ssd != 0) && (!flag[ssd]))
459  {
460  //std::cout << "C0" << std::endl;
461 
462  Vortex2D& vtxK = vtx[ssd];
463 
464  //std::cout << "C" << std::endl;
465 
466  flag[ssd] = true;
467 
468  Vortex2D sumVtx;
469  sumVtx.g() = vtxI.g() + vtxK.g();
470 
471  switch (type)
472  {
473  case 0:
474  case 2:
475  sumAbsGam = fabs(vtxI.g()) + fabs(vtxK.g());
476 
477  iws = sumAbsGam > 1e-10 ? 1.0 / sumAbsGam : 1.0;
478 
479  sumVtx.r() = (vtxI.r() * fabs(vtxI.g()) + vtxK.r() * fabs(vtxK.g())) * iws;
480  break;
481 
482  case 1:
483  sumVtx.r() = (fabs(vtxI.g()) > fabs(vtxK.g())) ? vtxI.r() : vtxK.r();
484  break;
485  }
486 
487  bool fl_hit = true;
488  //double Ch1[2];
489  size_t hitpan = -1;
490 
491  //std::cout << "D" << std::endl;
492 
493  for (size_t afl = 0; afl < W.getNumberOfAirfoil(); ++afl)
494  {
495  //проверим, не оказался ли новый вихрь внутри контура
496  if (MoveInside(sumVtx.r(), vtxI.r(), W.getAirfoil(afl), hitpan) || MoveInside(sumVtx.r(), vtxK.r(), W.getAirfoil(afl), hitpan))
497  fl_hit = false;
498  }//for
499 
500  //std::cout << "E" << std::endl;
501 
502  if (fl_hit)
503  {
504  vtxI = sumVtx;
505  vtxK.g() = 0.0;
506  nHlop++;
507  }//if (fl_hit)
508 
509  //std::cout << "F" << std::endl;
510 
511  }//if ((ssd!=0)&&(!flag[ssd]))
512  }//if !flag[vt]
513  }//for vt
514 
515  }//for z
516 
517  return nHlop;
518 }//Collaps(...)
std::pair< std::string, int > velocityComputation
Definition: Passport2D.h:190
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
HD double & g()
Функция для доступа к циркуляции вихря
Definition: Vortex2D.h:93
Класс, обеспечивающий возможность выполнения вычислений на GPU по технологии Nvidia CUDA...
Definition: Gpu2D.h:67
void GetPairs(int type)
Поиск ближайшего соседа
Definition: Wake2D.cpp:317
bool MoveInside(const Point2D &newPos, const Point2D &oldPos, const Airfoil &afl, size_t &panThrough)
Проверка проникновения точки через границу профиля
Definition: Wake2D.cpp:66
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
std::vector< Vortex2D > vtx
Список вихревых элементов
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition: Passport2D.h:302
const World2D & W
Константная ссылка на решаемую задачу
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
const Gpu & getCuda() const
Возврат константной ссылки на объект, связанный с видеокартой (GPU)
Definition: World2D.h:227
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
std::vector< int > neighb
Вектор потенциальных соседей для будущего коллапса
Definition: Wake2D.h:63
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130

Here is the call graph for this function:

Here is the caller graph for this function:

int Wake::CollapsNew ( int  type,
int  times 
)

Definition at line 521 of file Wake2D.cpp.

522 {
523  int nHlop = 0; //общее число убитых вихрей
524 
525  const double& cSP = collapseScaleParameter;
526  const double& cRBP = collapseRightBorderParameter;
527 
529 
530  //int loc_hlop = 0; // neigh
531 
532  std::vector<bool> flag; //как только вихрь сколлапсировался flag = 1
533  for (int z = 0; z < times; ++z)
534  {
535  //loc_hlop = 0;//число схлопнутых вихрей
536 
537  flag.clear();
538  flag.resize(vtx.size(), false);
539 
540  double sumAbsGam, iws, r2, r2test;
541  Point2D newPos;
542 
543  for (size_t vt = 0; vt + 1 < vtx.size(); ++vt)
544  {
545  Vortex2D& vtxI = vtx[vt];
546  if (!flag[vt])
547  {
548  for (int s = 0; s < knb - 1; ++s)
549  {
550  int ssd = neighbNew[vt * (knb - 1) + s];
551  Vortex2D& vtxK = vtx[ssd];
552 
553  r2 = dist2(vtxI.r(), vtxK.r());
554  //линейное увеличение радиуса коллапса
555  double mnog = std::max(1.0, /* 2.0 * */ (vtxI.r()[0] - cRBP) / cSP);
557  if (type == 1)
558  r2test *= 4.0; //Увеличение радиуса коллапса в 2 раза для коллапса вихрей разных знаков
559 
560  bool keyGamma = false;
561  switch (type)
562  {
563  case 0:
564  keyGamma = (fabs(vtxI.g() * vtxK.g()) != 0.0) && (fabs(vtxI.g() + vtxK.g()) < sqr(mnog) * maxG);
565  break;
566  case 1:
567  keyGamma = (vtxI.g() * vtxK.g() < 0.0);
568  break;
569  case 2:
570  keyGamma = (vtxI.g() * vtxK.g() > 0.0) && (fabs(vtxI.g() + vtxK.g()) < sqr(mnog) * maxG);
571  break;
572  }
573 
574  if ((r2 < r2test) && (keyGamma))
575  {
576  if ((ssd != 0) && (!flag[ssd]))
577  {
578  flag[ssd] = true;
579 
580  Vortex2D sumVtx;
581  sumVtx.g() = vtxI.g() + vtxK.g();
582 
583  switch (type)
584  {
585  case 0:
586  case 2:
587  sumAbsGam = fabs(vtxI.g()) + fabs(vtxK.g());
588 
589  if (sumAbsGam > 1e-10)
590  sumVtx.r() = (vtxI.r() * fabs(vtxI.g()) + vtxK.r() * fabs(vtxK.g())) * (1.0 / sumAbsGam);
591  else
592  sumVtx.r() = 0.5 * (vtxI.r() + vtxK.r());
593 
594  break;
595 
596  case 1:
597  sumVtx.r() = (fabs(vtxI.g()) > fabs(vtxK.g())) ? vtxI.r() : vtxK.r();
598  break;
599  }
600 
601  bool fl_hit = true;
602  size_t hitpan = -1;
603 
604  for (size_t afl = 0; afl < W.getNumberOfAirfoil(); ++afl)
605  {
606  //проверим, не оказался ли новый вихрь внутри контура
607  if (MoveInside(sumVtx.r(), vtxI.r(), W.getAirfoil(afl), hitpan) || MoveInside(sumVtx.r(), vtxK.r(), W.getAirfoil(afl), hitpan))
608  fl_hit = false;
609  }//for
610 
611  if (fl_hit)
612  {
613  vtxI = sumVtx;
614  vtxK.g() = 0.0;
615  nHlop++;
616  flag[vt] = true;
617  }//if (fl_hit)
618 
619  }//if ((ssd!=0)&&(!flag[ssd]))
620  }
621 
622  if (flag[vt])
623  break;
624  }// for s
625  }//if !flag[vt]
626  }//for vt
627 
628  }//for z
629 
630  return nHlop;
631 }//Collaps(...)
double maxGamma
Максимально допустимая циркуляция вихря
Definition: Passport2D.h:157
static const int knb
Definition: Wake2D.h:65
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
HD double & g()
Функция для доступа к циркуляции вихря
Definition: Vortex2D.h:93
double collapseRightBorderParameter
абсцисса, правее которой происходит линейный (вправо) рост радиуса коллапса
Definition: Wake2D.h:152
T sqr(T x)
Возведение числа в квадрат
Definition: defs.h:430
bool MoveInside(const Point2D &newPos, const Point2D &oldPos, const Airfoil &afl, size_t &panThrough)
Проверка проникновения точки через границу профиля
Definition: Wake2D.cpp:66
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
auto dist2(const numvector< T, n > &x, const numvector< P, n > &y) -> typename std::remove_const< decltype(x[0]-y[0])>::type
Вычисление квадрата расстояния между двумя точками
Definition: numvector.h:788
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
std::vector< int > neighbNew
Definition: Wake2D.h:66
double collapseScaleParameter
характерный масштаб, на котором происходит рост радиуса коллапса
Definition: Wake2D.h:155
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
double epscol
Радиус коллапса
Definition: Passport2D.h:145

Here is the call graph for this function:

Here is the caller graph for this function:

void Wake::GetPairs ( int  type)

Поиск ближайшего соседа

Parameters
[in]typeтип коллапса:
  • 0 — без приоритета знаков
  • 1 — коллапсировать только вихри разных знаков
  • 2 — коллапсировать только вихри одного знака

Definition at line 317 of file Wake2D.cpp.

318 {
319  GetPairsBS(type);
320 }
void GetPairsBS(int type)
Definition: Wake2D.cpp:323

Here is the call graph for this function:

Here is the caller graph for this function:

void VM2D::Wake::GetPairsBH ( int  type)

Here is the caller graph for this function:

void Wake::GetPairsBS ( int  type)
Todo:
Доделать

Definition at line 323 of file Wake2D.cpp.

324 {
325  neighb.resize(0);
326  neighb.resize(vtx.size(), 0);
327 
328  Point2D Ri, Rk;
329 
331 
332 #pragma omp parallel for default(none) shared(type, maxG) schedule(dynamic, DYN_SCHEDULE)
333  for (int i = 0; i < vtx.size(); ++i)
334  {
335  int s = i;
336  const Vortex2D& vtxI = vtx[i];
337 
338  bool found = false;
339 
340  double r2, r2test;
341 
342  const double& cSP = collapseScaleParameter;
343  const double& cRBP = collapseRightBorderParameter;
344 
345  while ( (!found) && ( s + 1 < (int)vtx.size() ) )
346  {
347  s++;
348  const Vortex2D& vtxK = vtx[s];
349 
350  r2 = dist2(vtxI.r(), vtxK.r());
351 
352  //линейное увеличение радиуса коллапса
353  double mnog = std::max(1.0, /* 2.0 * */ (vtxI.r()[0] - cRBP) / cSP);
354 
355  r2test = sqr( W.getPassport().wakeDiscretizationProperties.epscol * mnog );
356 
357  if (type == 1)
358  r2test *= 4.0; //Увеличение радиуса коллапса в 2 раза для коллапса вихрей разных знаков
359 
360  //if (mnog > 1.0)
361  // r2test = sqr(0.005*cSP);
362 
363  if (r2 < r2test)
364  {
365  switch (type)
366  {
367  case 0:
368  found = ( fabs(vtxI.g()*vtxK.g()) != 0.0) && (fabs(vtxI.g() + vtxK.g()) < sqr(mnog) * maxG);
369  break;
370  case 1:
371  found = (vtxI.g()*vtxK.g() < 0.0);
372  break;
373  case 2:
374  found = (vtxI.g()*vtxK.g() > 0.0) && (fabs(vtxI.g() + vtxK.g()) < sqr(mnog) * maxG);
375  break;
376  }
377  }//if r2 < r2_test
378  }//while
379 
380  if (found)
381  neighb[i] = s;
382  }//for locI
383 }//GetPairsBS(...)
double maxGamma
Максимально допустимая циркуляция вихря
Definition: Passport2D.h:157
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
HD double & g()
Функция для доступа к циркуляции вихря
Definition: Vortex2D.h:93
double collapseRightBorderParameter
абсцисса, правее которой происходит линейный (вправо) рост радиуса коллапса
Definition: Wake2D.h:152
T sqr(T x)
Возведение числа в квадрат
Definition: defs.h:430
auto dist2(const numvector< T, n > &x, const numvector< P, n > &y) -> typename std::remove_const< decltype(x[0]-y[0])>::type
Вычисление квадрата расстояния между двумя точками
Definition: numvector.h:788
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
std::vector< int > neighb
Вектор потенциальных соседей для будущего коллапса
Definition: Wake2D.h:63
double collapseScaleParameter
характерный масштаб, на котором происходит рост радиуса коллапса
Definition: Wake2D.h:155
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
double epscol
Радиус коллапса
Definition: Passport2D.h:145

Here is the call graph for this function:

Here is the caller graph for this function:

void Wake::Inside ( const std::vector< Point2D > &  newPos,
Airfoil afl,
bool  isMoves,
const Airfoil oldAfl 
)

Проверка пересечения вихрями следа профиля при перемещении

Исполняется сразу для всех вихрей в пелене, осуществляет проверку для отдельного профиля
Вихри, попавшие внутрь профиля, получают нулевую циркуляцию, а их "бывшая" циркуляция передается в вектор gammaThrough в структуру данных профиля

Parameters
[in]newPosконстантная ссылка на вектор из новых положений вихрей в вихревом следе
[in]isMovesпризнак того, что профиль подвижный
[in]oldAflконстантная ссылка контролируемый профиль до перемещения (используется, если у профиля стоит признак того, что он движется)
[in,out]aflссылка на контролируемый профиль (происходит изменение afl->gammaThrough)

Definition at line 279 of file Wake2D.cpp.

280 {
281  std::vector<double> gamma;
282  gamma.resize(afl.getNumberOfPanels(), 0.0);
283 
284  std::vector<int> through;
285  through.resize(vtx.size(), -1);
286 
287 #pragma omp parallel for default(none) shared(afl, oldAfl, isMoves, through, newPos)
288  for (int i = 0; i < (int)vtx.size(); ++i)
289  {
290  size_t minN;
291 
292  bool crit = isMoves ? MoveInsideMovingBoundary(newPos[i], vtx[i].r(), oldAfl, afl, minN) : MoveInside(newPos[i], vtx[i].r(), afl, minN);
293 
294  if (crit)
295  through[i] = (int)minN;
296  }//for i
297 
298 
299  //std::stringstream sss;
300  //sss << "through_" << W.currentStep;
301  //std::ofstream of(W.getPassport().dir + "dbg/" + sss.str());
302  //for (size_t i = 0; i < gamma.size(); ++i)
303  // of << gamma[i] << std::endl;
304  //of.close();
305 
306  for (size_t q = 0; q < through.size(); ++q)
307  if (through[q] > -1)
308  {
309  gamma[through[q]] += vtx[q].g();
310  vtx[q].g() = 0.0;
311  }
312 
313  afl.gammaThrough = gamma;
314 }//Inside(...)
std::vector< double > gammaThrough
Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
Definition: Airfoil2D.h:196
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
bool MoveInside(const Point2D &newPos, const Point2D &oldPos, const Airfoil &afl, size_t &panThrough)
Проверка проникновения точки через границу профиля
Definition: Wake2D.cpp:66
bool MoveInsideMovingBoundary(const Point2D &newPos, const Point2D &oldPos, const Airfoil &oldAfl, const Airfoil &afl, size_t &panThrough)
Проверка проникновения точки через границу профиля
Definition: Wake2D.cpp:185
std::vector< Vortex2D > vtx
Список вихревых элементов

Here is the call graph for this function:

Here is the caller graph for this function:

bool Wake::MoveInside ( const Point2D newPos,
const Point2D oldPos,
const Airfoil afl,
size_t &  panThrough 
)

Проверка проникновения точки через границу профиля

Parameters
[in]newPosконстантная ссылка на смещенное (новое) положение
[in]oldPosконстантная ссылка на несмещенное (старое) положение
[in]aflконстантная ссылка на контролируемый профиль
[out]panThroughномер "протыкаемой" панели return признак пересечения профиля

Definition at line 66 of file Wake2D.cpp.

67 {
68  //const double porog_r = 1e-12;
69 
70  //double minDist = 1.0E+10; //расстояние до пробиваемой панели
71  panThrough = -1;
72 
73  //проверка габ. прямоугольника
74  if (afl.isOutsideGabarits(newPos) && afl.isOutsideGabarits(oldPos))
75  return false;
76 
77  //если внутри габ. прямоугольника - проводим контроль
78  bool hit = false;
79 
80  //std::pair<double, double> gabPointX = std::minmax(newPos[0], oldPos[0]);
81  //std::pair<double, double> gabPointY = std::minmax(newPos[1], oldPos[1]);
82 
83  //std::pair<double, double> gabPanelX, gabPanelY;
84 
85  //Проверка на пересечение
86  //double r0 = 0, r1 = 0, r2 = 0, r3 = 0;
87 
88  //int cnt = 0;
89  for (size_t j = 0; j < afl.getNumberOfPanels(); ++j)
90  {
91  const Point2D& aflRj = afl.getR(j);
92  const Point2D& aflRj1 = afl.getR(j + 1);
93 
94 // gabPanelX = std::minmax(aflRj[0], aflRj1[0]);
95 // gabPanelY = std::minmax(aflRj[1], aflRj1[1]);
96 
97 // if ((gabPointX.second >= gabPanelX.first) && (gabPanelX.second >= gabPointX.first) && (gabPointY.second >= gabPanelY.first) && (gabPanelY.second >= gabPointY.first))
98 // {
99  if ((((aflRj - oldPos) ^ (newPos - oldPos)) * ((aflRj1 - oldPos) ^ (newPos - oldPos)) <= 0) && \
100  (((oldPos - aflRj) ^ (aflRj1 - aflRj)) * ((newPos - aflRj) ^ (aflRj1 - aflRj)) <= 0))
101  {
102  hit = true;
103  panThrough = j;
104  break;
105  }
106 // }
107 // else { ++cnt; }
108  }//for j
109 
110 // std::cout << cnt << std::endl;
111 
112  return hit;
113 }//MoveInside(...)
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
bool isOutsideGabarits(const Point2D &r) const
Определяет, находится ли точка с радиус-вектором вне габаритного прямоугольника профиля ...
Definition: Airfoil2D.cpp:76

Here is the call graph for this function:

Here is the caller graph for this function:

bool Wake::MoveInsideMovingBoundary ( const Point2D newPos,
const Point2D oldPos,
const Airfoil oldAfl,
const Airfoil afl,
size_t &  panThrough 
)

Проверка проникновения точки через границу профиля

Parameters
[in]newPosконстантная ссылка на смещенное (новое) положение вихря
[in]oldPosконстантная ссылка на несмещенное (старое) положение вихря
[in]oldAflконстантная ссылка на состояние контролируемого профиля до перемещения
[in]aflконстантная ссылка на контролируемый профиль
[out]panThroughномер "протыкаемой" панели return признак пересечения профиля
Todo:
сравнить производительности двух inside-ов

Definition at line 185 of file Wake2D.cpp.

186 {
187  panThrough = -1;
188 
190 
191  //проверка габ. прямоугольника
192  if (!afl.inverse && afl.isOutsideGabarits(newPos) && afl.isOutsideGabarits(oldPos))
193  return false;
194 
195  bool hit = false;
196 
197  double angle = 0;
198  double cs, sn;
199 
200  double dist2 = 1000000000.0;
201 
202  for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
203  {
204  Point2D v1, v2, vv;
205  v1 = afl.getR(i) - newPos;
206 
207  v2 = afl.getR(i + 1) - newPos;
208 
209  vv = afl.getR(i + 1) - afl.getR(i);
210 
211  double dst = v1.length2() + v2.length2();
212  if (dst < dist2)
213  {
214  dist2 = dst;
215  panThrough = i;
216  }
217 
218  cs = v1 & v2;
219  sn = v1 ^ v2;
220 
221  angle += atan2(sn, cs);
222  }//for i
223 
224  hit = ((angle > 3.14) || (angle < -3.14));
225 
226  if (afl.inverse)
227  hit = !hit;
228 
229  return hit;
230 }//MoveInsideMovingBoundary(...)
bool inverse
Признак разворота нормалей (для расчета внутренних течений)
Definition: Airfoil2D.h:139
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition: numvector.h:383
auto dist2(const numvector< T, n > &x, const numvector< P, n > &y) -> typename std::remove_const< decltype(x[0]-y[0])>::type
Вычисление квадрата расстояния между двумя точками
Definition: numvector.h:788
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
bool isOutsideGabarits(const Point2D &r) const
Определяет, находится ли точка с радиус-вектором вне габаритного прямоугольника профиля ...
Definition: Airfoil2D.cpp:76

Here is the call graph for this function:

Here is the caller graph for this function:

void WakeDataBase::ReadFromFile ( const std::string &  dir,
const std::string &  fileName 
)
inherited

Считывание вихревого следа из файла

Parameters
[in]dirконстантная ссылка на строку, задающую каталог, где лежит файл с вихревым следом
[in]fileNameконстантная ссылка на строку, задающую имя файла с вихревым следом

Definition at line 56 of file WakeDataBase2D.cpp.

57 {
58  std::string filename = dir + fileName;
59  std::ifstream wakeFile;
60 
61  if (fileExistTest(filename, W.getInfo(), { "txt", "TXT" }))
62  {
63  std::stringstream wakeFile(VMlib::Preprocessor(filename).resultString);
64 
65  VMlib::LogStream XXX;
66  VMlib::StreamParser wakeParser(XXX, "vortex wake file parser", wakeFile);
67 
68  wakeParser.get("vtx", vtx);
69  }
70 }//ReadFromFile(...)
Класс, определяющий работу с потоком логов
Definition: LogStream.h:53
Класс, позволяющий выполнять предварительную обработку файлов
Definition: Preprocessor.h:59
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
bool fileExistTest(std::string &fileName, LogStream &info, const std::list< std::string > &extList={})
Проверка существования файла
Definition: defs.h:324
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода ...
Definition: WorldGen.h:74
Класс, позволяющий выполнять разбор файлов и строк с настройками и параметрами
Definition: StreamParser.h:151

Here is the call graph for this function:

Here is the caller graph for this function:

int Wake::RemoveFar ( )

Зануление далеко улетевших вихрей

Returns
число вихрей в дальнем следе, которые занулены
Todo:
Пока профиль 1, расстояние от его центра; сделать от самого правого профиля

Definition at line 635 of file Wake2D.cpp.

636 {
637  int nFar = 0;
640  Point2D zerovec = { 0.0, 0.0 };
641 #pragma omp parallel for default(none) shared(distFar2, zerovec) reduction(+:nFar)
642  for (int i = 0; i <static_cast<int>(vtx.size()); ++i)
643  {
644  if (dist2(vtx[i].r(), zerovec) > distFar2)
645  {
646  vtx[i].g() = 0.0;
647  nFar++;
648  }
649  }
650 
651  return nFar;
652 }//RemoveFar()
T sqr(T x)
Возведение числа в квадрат
Definition: defs.h:430
auto dist2(const numvector< T, n > &x, const numvector< P, n > &y) -> typename std::remove_const< decltype(x[0]-y[0])>::type
Вычисление квадрата расстояния между двумя точками
Definition: numvector.h:788
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
double distFar
Расстояние от центра самого подветренного (правого) профиля, на котором вихри уничтожаются ...
Definition: Passport2D.h:148

Here is the call graph for this function:

Here is the caller graph for this function:

size_t Wake::RemoveZero ( )

Исключение нулевых и мелких вихрей

Returns
число исключенных вихрей

Definition at line 655 of file Wake2D.cpp.

656 {
657  const double porog_g = 1e-15;
658 
659  std::vector<Vortex2D/*, VM2D::MyAlloc<VMlib::Vortex2D>*/> newWake;
660 
661  newWake.reserve(vtx.size());
662 
663  for (size_t q = 0; q < vtx.size(); ++q)
664  if (fabs(vtx[q].g()) > porog_g)
665  newWake.push_back(vtx[q]);
666 
667  size_t delta = vtx.size() - newWake.size();
668 
669  newWake.swap(vtx);
670 
671  return delta;
672 }//RemoveZero()
std::vector< Vortex2D > vtx
Список вихревых элементов
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51

Here is the caller graph for this function:

void Wake::Restruct ( )

Реструктуризация вихревого следа

Исполняется сразу для всех вихрей в пелене
Вихри, находящиеся далеко от профилей, удаляются
Вихри, которые сильно сблизились, коллапсируются

Definition at line 685 of file Wake2D.cpp.

686 {
687  W.getTimestat().timeRestruct.first += omp_get_wtime();
688 
689  // Определение параметров, отвечающих за увеличение радиуса коллапса
690  std::vector<double> rightBorder, horizSpan;
691  rightBorder.reserve(W.getNumberOfAirfoil());
692  horizSpan.reserve(W.getNumberOfAirfoil());
693 
694  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
695  {
696 
697  rightBorder.emplace_back(W.getAirfoil(q).upRight[0]);
698  horizSpan.emplace_back(W.getAirfoil(q).upRight[0] - W.getAirfoil(q).lowLeft[0]);
699  }
700 
701  if (W.getNumberOfAirfoil() > 0)
702  {
703  W.getNonConstWake().collapseRightBorderParameter = *std::max_element(rightBorder.begin(), rightBorder.end());
704  W.getNonConstWake().collapseScaleParameter = *std::max_element(horizSpan.begin(), horizSpan.end());
705  }
706  else
707  {
710  }
711 #if defined(__CUDACC__) || defined(USE_CUDA)
713 #endif
714 
715 
716 
718 
719 
720  //CPU/GPU - прямой алгоритм
721  //Collaps(1, 1);
722  //Collaps(2, 1);
723 
724 
725  //*
726  //быстрый алгоритм
727  for (int collapsStep = 1; collapsStep <= 2; ++collapsStep)
728  {
729 #if defined(USE_CUDA)
730 #define knnGPU
731 #endif
732 
733  neighbNew.resize(vtx.size() * (knb - 1));
734 
735  if (vtx.size() > 2 * knb)
736  {
737 #ifndef knnGPU
738  //CPU
739  std::vector<std::vector<std::pair<double, size_t>>> initdist(vtx.size());
740  for (auto& d : initdist)
741  d.resize(2 * knb, { -1.0, -1 });
742 
743  WakekNN(vtx, knb, initdist);//CPU
744 #else
745  std::vector<std::pair<double, size_t>> initdistcuda(knb * vtx.size(), { -1.0, -1 }); //CUDA
746  kNNcuda(vtx, knb, initdistcuda, vecForKnn); //CUDA
747 #endif
748 
749  for (int i = 0; i < vtx.size(); ++i)
750  {
751 #ifndef knnGPU
752  neighbNew[i] = (int)initdist[i][1].second;
753 #else
754  for (int j = 0; j < knb - 1; ++j)
755  neighbNew[i * (knb - 1) + j] = (int)initdistcuda[(i * knb) + (j + 1)].second;
756 #endif
757  }
758  }
759 
760  CollapsNew(collapsStep, 1);
761 
762  }
763 
765 //*/
766 
767  RemoveFar();
768  RemoveZero();
769 
770  W.getTimestat().timeRestruct.second += omp_get_wtime();
771 }//Restruct()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
Gpu & getNonConstCuda() const
Возврат неконстантной ссылки на объект, связанный с видеокартой (GPU)
Definition: World2D.h:232
timePeriod timeRestruct
Начало и конец процесса реструктуризации пелены
Definition: Times2D.h:101
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
Point2D lowLeft
Левый нижний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:190
static const int knb
Definition: Wake2D.h:65
int RemoveFar()
Зануление далеко улетевших вихрей
Definition: Wake2D.cpp:635
double collapseRightBorderParameter
абсцисса, правее которой происходит линейный (вправо) рост радиуса коллапса
Definition: Wake2D.h:152
size_t RemoveZero()
Исключение нулевых и мелких вихрей
Definition: Wake2D.cpp:655
Point2D upRight
Правый верхний угол габаритного прямоугольника профиля
Definition: Airfoil2D.h:191
Wake & getNonConstWake() const
Возврат неконстантной ссылки на вихревой след
Definition: World2D.h:202
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
void setCollapseCoeff(double pos_, double refLength_)
Установка правой границы самого правого профиля (для организации увеличения радиуса коллапса) ...
Definition: Gpu2D.h:253
std::vector< int > neighbNew
Definition: Wake2D.h:66
double collapseScaleParameter
характерный масштаб, на котором происходит рост радиуса коллапса
Definition: Wake2D.h:155
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
void WakekNN(const std::vector< Vortex2D > &vtx, const size_t k, std::vector< std::vector< std::pair< double, size_t >>> &initdist)
Definition: knnCPU.cpp:447
int CollapsNew(int type, int times)
Definition: Wake2D.cpp:521

Here is the call graph for this function:

Here is the caller graph for this function:

void WakeDataBase::SaveKadrVtk ( const std::string &  filePrefix = "Kadr") const
inherited

Сохранение вихревого следа в файл .vtk.

Definition at line 73 of file WakeDataBase2D.cpp.

74 {
75  W.getTimestat().timeSaveKadr.first += omp_get_wtime();
76 
78  {
79  std::ofstream outfile;
80  size_t numberNonZero = 0;
81 
82  if (vtx.size() > 0)
83  numberNonZero += vtx.size();
84  else
85  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
86  numberNonZero += W.getAirfoil(q).getNumberOfPanels();
87  VMlib::CreateDirectory(W.getPassport().dir, "snapshots");
88 
89  if (W.getPassport().timeDiscretizationProperties.fileTypeVtx.second == 0) //text format vtk
90  {
91  std::string fname = VMlib::fileNameStep(filePrefix, W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "vtk");
92  outfile.open(W.getPassport().dir + "snapshots/" + fname);
93 
94  outfile << "# vtk DataFile Version 2.0" << std::endl;
95  outfile << "VM2D VTK result: " << (W.getPassport().dir + "snapshots/" + fname).c_str() << " saved " << VMlib::CurrentDataTime() << std::endl;
96  outfile << "ASCII" << std::endl;
97  outfile << "DATASET UNSTRUCTURED_GRID" << std::endl;
98  outfile << "POINTS " << numberNonZero << " float" << std::endl;
99 
100  if (vtx.size() > 0)
101  for (auto& v : vtx)
102  {
103  const Point2D& r = v.r();
104  outfile << r[0] << " " << r[1] << " " << "0.0" << std::endl;
105  }//for v
106  else
107  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
108  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
109  {
110  const Point2D& r = W.getAirfoil(q).getR(s);
111  outfile << r[0] << " " << r[1] << " " << "0.0" << std::endl;
112  }
113 
114  outfile << "CELLS " << numberNonZero << " " << 2 * numberNonZero << std::endl;
115  for (size_t i = 0; i < numberNonZero; ++i)
116  outfile << "1 " << i << std::endl;
117 
118  outfile << "CELL_TYPES " << numberNonZero << std::endl;
119  for (size_t i = 0; i < numberNonZero; ++i)
120  outfile << "1" << std::endl;
121 
122  outfile << std::endl;
123  outfile << "POINT_DATA " << numberNonZero << std::endl;
124  outfile << "SCALARS Gamma float 1" << std::endl;
125  outfile << "LOOKUP_TABLE default" << std::endl;
126 
127  if (vtx.size() > 0)
128  for (auto& v : vtx)
129  outfile << v.g() << std::endl;
130  else
131  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
132  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
133  {
134  //const Point2D& r = W.getAirfoil(q).getR(s);
135  outfile << "0.0" << std::endl;
136  }
137 
138  outfile.close();
139  }//if fileType = text
140  else if (W.getPassport().timeDiscretizationProperties.fileTypeVtx.second == 1) //binary format vtk
141  {
142  //Тест способа хранения чисел
143  uint16_t x = 0x0001;
144  bool littleEndian = (*((uint8_t*)&x));
145  const char eolnBIN[] = "\n";
146 
147  std::string fname = VMlib::fileNameStep(filePrefix, W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "vtk");
148  outfile.open(W.getPassport().dir + "snapshots/" + fname, std::ios::out | std::ios::binary);
149 
150  outfile << "# vtk DataFile Version 3.0" << "\r\n" << "VM2D VTK result: " << (W.getPassport().dir + "snapshots/" + fname).c_str() << " saved " << VMlib::CurrentDataTime() << eolnBIN;
151  outfile << "BINARY" << eolnBIN;
152  outfile << "DATASET UNSTRUCTURED_GRID" << eolnBIN << "POINTS " << numberNonZero << " " << "float" << eolnBIN;
153 
154 
155  if (vtx.size() > 0)
156  {
157  Eigen::VectorXf rData = Eigen::VectorXf::Zero(vtx.size() * 3);
158  for (size_t i = 0; i < vtx.size(); ++i)
159  {
160  rData(3 * i) = (float)(vtx[i].r())[0];
161  rData(3 * i + 1) = (float)(vtx[i].r())[1];
162  }//for i
163 
164  if (littleEndian)
165  for (int i = 0; i < vtx.size() * 3; ++i)
166  VMlib::SwapEnd(rData(i));
167  outfile.write(reinterpret_cast<char*>(rData.data()), vtx.size() * 3 * sizeof(float));
168  }
169  else
170  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
171  {
172  Eigen::VectorXf rData = Eigen::VectorXf::Zero(W.getAirfoil(q).getNumberOfPanels() * 3);
173  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
174  {
175  rData(3 * s) = (float)(W.getAirfoil(q).getR(s))[0];
176  rData(3 * s + 1) = (float)(W.getAirfoil(q).getR(s))[1];
177  }//for i
178 
179  if (littleEndian)
180  for (int i = 0; i < W.getAirfoil(q).getNumberOfPanels() * 3; ++i)
181  VMlib::SwapEnd(rData(i));
182  outfile.write(reinterpret_cast<char*>(rData.data()), W.getAirfoil(q).getNumberOfPanels() * 3 * sizeof(float));
183  }
184 
185  // CELLS
186  std::vector<int> cells(2 * numberNonZero);
187  for (size_t i = 0; i < numberNonZero; ++i)
188  {
189  cells[2 * i] = 1;
190  cells[2 * i + 1] = (int)i;
191  }
192 
193  std::vector<int> cellsTypes;
194  cellsTypes.resize(numberNonZero, 1);
195 
196  if (littleEndian)
197  {
198  for (int i = 0; i < numberNonZero * 2; ++i)
199  VMlib::SwapEnd(cells[i]);
200 
201  for (int i = 0; i < numberNonZero; ++i)
202  VMlib::SwapEnd(cellsTypes[i]);
203  }
204 
205  outfile << eolnBIN << "CELLS " << numberNonZero << " " << numberNonZero * 2 << eolnBIN;
206  outfile.write(reinterpret_cast<char*>(cells.data()), numberNonZero * 2 * sizeof(int));
207  outfile << eolnBIN << "CELL_TYPES " << numberNonZero << eolnBIN;
208  outfile.write(reinterpret_cast<char*>(cellsTypes.data()), numberNonZero * sizeof(int));
209 
210  //gammas
211  outfile << eolnBIN << "POINT_DATA " << numberNonZero << eolnBIN;
212  outfile << eolnBIN << "SCALARS Gamma " << "float" << " 1" << eolnBIN;
213  outfile << "LOOKUP_TABLE default" << eolnBIN;
214 
215  if (vtx.size() > 0)
216  {
217  Eigen::VectorXf pData = Eigen::VectorXf::Zero(numberNonZero);
218  for (int s = 0; s < vtx.size(); ++s)
219  pData(s) = (float)vtx[s].g();
220 
221  if (littleEndian)
222  for (int i = 0; i < vtx.size(); ++i)
223  VMlib::SwapEnd(pData(i));
224  outfile.write(reinterpret_cast<char*>(pData.data()), vtx.size() * sizeof(float));
225  }
226  else
227  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
228  {
229  Eigen::VectorXf pData = Eigen::VectorXf::Zero(W.getAirfoil(q).getNumberOfPanels());
230  for (int s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
231  pData(s) = 0;
232 
233  if (littleEndian)
234  for (int i = 0; i < W.getAirfoil(q).getNumberOfPanels(); ++i)
235  VMlib::SwapEnd(pData(i));
236  outfile.write(reinterpret_cast<char*>(pData.data()), W.getAirfoil(q).getNumberOfPanels() * sizeof(float));
237  }
238 
239  outfile << eolnBIN;
240  outfile.close();
241  }//if binary
242  if (W.getPassport().timeDiscretizationProperties.fileTypeVtx.second == 2) //csv
243  {
244  std::string fname = VMlib::fileNameStep(filePrefix, W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "csv");
245  outfile.open(W.getPassport().dir + "snapshots/" + fname);
246 
247  outfile << "point,x,y,G " << std::endl;
248 
249  int counter = 0;
250  if (vtx.size() > 0)
251  for (auto& v : vtx)
252  outfile << counter++ << "," << v.r()[0] << "," << v.r()[1] << "," << v.g() << std::endl;
253  else
254  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
255  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
256  {
257  const Point2D& r = W.getAirfoil(q).getR(s);
258  outfile << counter++ << "," << r[0] << "," << r[1] << "," << "0.0" << std::endl;
259  }
260  outfile.close();
261  }//if fileType = text
262 
263  }
264 
265  W.getTimestat().timeSaveKadr.second += omp_get_wtime();
266 }
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
int saveVtxStep
Шаг сохранения кадров в бинарные файлы
Definition: PassportGen.h:73
std::string CurrentDataTime()
Формирование строки с текущем временем и датой
Definition: defs.cpp:44
int nameLength
Число разрядов в имени файла
Definition: PassportGen.h:68
std::string dir
Рабочий каталог задачи
Definition: PassportGen.h:118
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::string fileNameStep(const std::string &name, int length, size_t number, const std::string &ext)
Формирование имени файла
Definition: defs.h:357
bool ifDivisible(int val) const
Definition: World2D.h:244
timePeriod timeSaveKadr
Начало и конец процесса сохранения кадра в файл
Definition: Times2D.h:107
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::pair< std::string, int > fileTypeVtx
Тип файлов для сохранения скорости и давления
Definition: PassportGen.h:71
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition: WorldGen.h:91
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
void SwapEnd(T &var)
Вспомогательная функция перестановки байт местами (нужно для сохранения бинарных VTK) ...
Definition: defs.h:528
void CreateDirectory(const std::string &dir, const std::string &name)
Создание каталога
Definition: defs.h:414

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

double VM2D::Wake::collapseRightBorderParameter

абсцисса, правее которой происходит линейный (вправо) рост радиуса коллапса

Definition at line 152 of file Wake2D.h.

double VM2D::Wake::collapseScaleParameter

характерный масштаб, на котором происходит рост радиуса коллапса

Definition at line 155 of file Wake2D.h.

const int VM2D::Wake::knb = 2
staticprivate

Definition at line 65 of file Wake2D.h.

std::vector<int> VM2D::Wake::neighb
private

Вектор потенциальных соседей для будущего коллапса

Definition at line 63 of file Wake2D.h.

std::vector<int> VM2D::Wake::neighbNew
private

Definition at line 66 of file Wake2D.h.

std::vector<Vortex2D> VM2D::WakeDataBase::vtx
inherited

Список вихревых элементов

Definition at line 76 of file WakeDataBase2D.h.

const World2D& VM2D::WakeDataBase::W
inherited

Константная ссылка на решаемую задачу

Definition at line 69 of file WakeDataBase2D.h.


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