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

Класс, определяющий способ вычисления скоростей More...

#include <Velocity2DBiotSavart.h>

Inheritance diagram for VM2D::VelocityBiotSavart:
Collaboration diagram for VM2D::VelocityBiotSavart:

Public Member Functions

 VelocityBiotSavart (const World2D &W_)
 Конструктор More...
 
virtual ~VelocityBiotSavart ()
 Деструктор More...
 
void CalcConvVeloToSetOfPointsFromWake (const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius)
 Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа More...
 
virtual void CalcConvVelo () override
 Вычисление конвективных скоростей вихрей и виртуальных вихрей в вихревом следе, а также в точках wakeVP. More...
 
virtual void FillRhs (Eigen::VectorXd &rhs) const override
 Расчет вектора правой части (всего) More...
 
virtual void CalcVeloToWakeVP () override
 Вычисление скоростей в точках wakeVP. More...
 
void GetWakeInfluenceToRhs (const Airfoil &afl, std::vector< double > &wakeRhs) const
 Генерация вектора влияния вихревого следа на профиль More...
 
void CalcDiffVeloI1I2ToWakeFromSheets (const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2) override
 
void CalcDiffVeloI1I2ToWakeFromWake (const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2) override
 
void CalcDiffVeloI1I2ToSetOfPointsFromWake (const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)
 Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек More...
 
void CalcDiffVeloI1I2ToSetOfPointsFromSheets (const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)
 
void CalcDiffVeloI1I2 ()
 Вычисление диффузионных скоростей вихрей и виртуальных вихрей в вихревом следе More...
 
void CalcDiffVeloI0I3 ()
 
void LimitDiffVelo (std::vector< Point2D > &diffVel)
 Контроль больших значений диффузионных скоростей More...
 
void CalcDiffVelo ()
 Вычисление диффузионных скоростей More...
 
void ResizeAndZero ()
 Очистка старых массивов под хранение скоростей, выделение новой памяти и обнуление More...
 
void SaveVisStress ()
 Сохранение вязких напряжений More...
 

Public Attributes

VortexesParams wakeVortexesParams
 Струтура, определяющая параметры вихрей в следе More...
 
std::vector< VortexesParamsvirtualVortexesParams
 Вектор струтур, определяющий параметры виртуальных вихрей для профилей More...
 

Protected Attributes

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

Detailed Description

Класс, определяющий способ вычисления скоростей

Способ вычисления скоростей

  • напрямую по закону Био — Савара
Author
Марчевский Илья Константинович
Сокол Ксения Сергеевна
Рятина Евгения Павловна
Колганова Александра Олеговна
\Version 1.12
\date 14 января 2024 г.

Definition at line 64 of file Velocity2DBiotSavart.h.

Constructor & Destructor Documentation

VelocityBiotSavart::VelocityBiotSavart ( const World2D W_)

Конструктор

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

Definition at line 58 of file Velocity2DBiotSavart.cpp.

58  :
59  Velocity(W_)
60 {
61 };
Velocity(const World2D &W_)
Конструктор
Definition: Velocity2D.h:114
VelocityBiotSavart::~VelocityBiotSavart ( )
virtual

Деструктор

Definition at line 64 of file Velocity2DBiotSavart.cpp.

65 {
66 };

Member Function Documentation

void VelocityBiotSavart::CalcConvVelo ( )
overridevirtual

Вычисление конвективных скоростей вихрей и виртуальных вихрей в вихревом следе, а также в точках wakeVP.

Warning
скорости приплюсовываются к уже имеющимся

Implements VM2D::Velocity.

Definition at line 69 of file Velocity2DBiotSavart.cpp.

70 {
71  //const int checkStep = 3;
72 
73 
74  double tt1 = omp_get_wtime();
75 
76  W.getTimestat().timeCalcVortexConvVelo.first += omp_get_wtime();
77 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_CONV_TOWAKE))
78  //if (W.getCurrentStep() == checkStep)
79  GPUWakeToWakeFAST(W.getWake(), wakeVortexesParams.convVelo, wakeVortexesParams.epsastWake, true, true);
80  //else
81  //GPUCalcConvVeloToSetOfPointsFromWake(W.getWake(), wakeVortexesParams.convVelo, wakeVortexesParams.epsastWake, true, true);
82 #else
84 #endif
85 
86 
87 
88 
89  double tt2 = omp_get_wtime();
90 
91 
92 /*
93  if (W.getCurrentStep() == checkStep)
94  {
95  std::ofstream pos("checkPos.txt");
96  pos.precision(17);
97  pos << W.getWake().vtx.size() << std::endl;
98  for (int i = 0; i < W.getWake().vtx.size(); ++i)
99  pos << W.getWake().vtx[i].r()[0] << " " << W.getWake().vtx[i].r()[1] << " " << W.getWake().vtx[i].g() << std::endl;
100  pos.close();
101 
102  std::ofstream of("checkVelo.txt");
103  of.precision(16);
104  for (int i = 0; i < wakeVortexesParams.convVelo.size(); ++i)
105  of << i << " " << wakeVortexesParams.convVelo[i][0] << " " << wakeVortexesParams.convVelo[i][1] << " " \
106  << wakeVortexesParams.epsastWake[i] << std::endl;
107  of.close();
108  exit(100500);
109  }
110 */
111 
112  //W.getInfo('t') << "Convective velocities time = " << int(((tt2 - tt1) * 1000)*10) / 10.0 << " ms, bodies = " << W.getWake().vtx.size() << std::endl;
113 
114 
115  std::vector<Point2D> nullVector(0);
116 
117  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
118  {
119 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_CONV_TOBOU))
120  GPUCalcConvVeloToSetOfPointsFromWake(W.getBoundary(bou).virtualWake, nullVector, virtualVortexesParams[bou].epsastWake, false, true);
121 #else
122  CalcConvVeloToSetOfPointsFromWake(W.getBoundary(bou).virtualWake, nullVector, virtualVortexesParams[bou].epsastWake, false, true);
123 #endif
124  }
125 
126  //double tt3 = omp_get_wtime();
127 
128 //вычисление конвективных скоростей по закону Био-Савара от виртуальных вихрей
129  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
130  {
131  //Влияние на пелену от виртуальных вихрей
132 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_CONVVIRT))
133  W.getBoundary(bou).GPUCalcConvVelocityToSetOfPointsFromSheets(W.getWake(), wakeVortexesParams.convVelo);
134 #else
136 #endif
137  }
138 
139  //double tt4 = omp_get_wtime();
140 
141  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
142  {
143  //Скороcти самих виртуальных вихрей
145  }
146 
147  //double tt5 = omp_get_wtime();
148 
149  //std::cout << tt2 - tt1 << " " << tt3 - tt2 << " " << tt4 - tt3 << " " << tt5 - tt4 << std::endl;
150 
151  W.getTimestat().timeCalcVortexConvVelo.second += omp_get_wtime();
152 
153  W.getTimestat().timeVP.first += omp_get_wtime();
154 
155  //вычисление скоростей в заданных точках только на соответствующем шаге по времени
157  CalcVeloToWakeVP();
158  W.getTimestat().timeVP.second += omp_get_wtime();
159 }//CalcConvVelo()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
virtual void CalcConvVelocityAtVirtualVortexes(std::vector< Point2D > &velo) const =0
Вычисление конвективной скорости в точках расположения виртуальных вихрей
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
void CalcConvVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius)
Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа ...
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
timePeriod timeCalcVortexConvVelo
Начало и конец процесса вычисления конвективных скоростей вихрей
Definition: Times2D.h:86
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition: Velocity2D.h:108
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
virtual void CalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase &pointsDb, std::vector< Point2D > &velo) const =0
Вычисление конвективных скоростей в наборе точек, вызываемых наличием слоев вихрей и источников на пр...
int saveVPstep
Шаг вычисления и сохранения скорости и давления
Definition: PassportGen.h:78
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition: Velocity2D.h:105
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition: Velocity2D.h:85
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
timePeriod timeVP
Начало и конец процесса подсчета полей скоростей и давления и сохранения их в файл ...
Definition: Times2D.h:110
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition: WorldGen.h:91
std::vector< Point2D > convVelo
Вектор конвективных скоростей вихрей
Definition: Velocity2D.h:67
virtual void CalcVeloToWakeVP() override
Вычисление скоростей в точках wakeVP.

Here is the call graph for this function:

void VelocityBiotSavart::CalcConvVeloToSetOfPointsFromWake ( const WakeDataBase pointsDb,
std::vector< Point2D > &  velo,
std::vector< double > &  domainRadius,
bool  calcVelo,
bool  calcRadius 
)

Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа

Parameters
[in]pointsDbконстантная ссылка на базу данных пелены из вихрей, в которых надо сосчитать конвективные скорости
[out]veloссылка на вектор скоростей в требуемых точках
[out]domainRadiusссылка на вектор радиусов вихревых доменов
[in]calcVeloпризнак вычисления скоростей в точках
[in]calcRadiusпризнак вычисления радиусов доменов

Definition at line 231 of file Velocity2DBiotSavart.cpp.

232 {
233  std::vector<Point2D> selfVelo(pointsDb.vtx.size());
234  domainRadius.resize(pointsDb.vtx.size());
235 
236  double cft = IDPI;
237 
238 #pragma warning (push)
239 #pragma warning (disable: 4101)
240  //Локальные переменные для цикла
241  Point2D velI;
242  Point2D tempVel;
243  double dst2eps, dst2;
244 #pragma warning (pop)
245 
247 
248  if (calcVelo)
249  {
250 #pragma omp parallel for default(none) shared(selfVelo, cft, calcVelo, calcRadius, eps2, pointsDb, domainRadius) private(tempVel, velI, dst2, dst2eps) schedule(dynamic, DYN_SCHEDULE)
251  for (int i = 0; i < pointsDb.vtx.size(); ++i)
252  {
253  double ee2[3] = { 10000.0, 10000.0, 10000.0 };
254 
255  velI.toZero();
256 
257  const Point2D& posI = pointsDb.vtx[i].r();
258 
259  for (size_t j = 0; j < W.getWake().vtx.size(); ++j)
260  {
261  const Point2D& posJ = W.getWake().vtx[j].r();
262 
263  dst2 = (posI-posJ).length2();
264 
265  //Модифицируем массив квадратов расстояний до ближайших вихрей из wake
266 #ifndef TESTONLYVELO
267  if (calcRadius)
268  VMlib::ModifyE2(ee2, dst2);
269 #endif // !TESTONLYVELO
270 
271  const double& gamJ = W.getWake().vtx[j].g();
272 
273  tempVel.toZero();
274  dst2eps = VMlib::boundDenom(dst2, eps2); //Сглаживать надо!!!
275 
276  tempVel = { -posI[1] + posJ[1], posI[0] - posJ[0] };
277  tempVel *= (gamJ / dst2eps);
278  velI += tempVel;
279  }
280 
281  for (size_t j = 0; j < W.getSource().vtx.size(); ++j)
282  {
283  const Point2D& posJ = W.getSource().vtx[j].r();
284  const double& gamJ = W.getPassport().physicalProperties.accelCft() * W.getSource().vtx[j].g();
285 
286  tempVel.toZero();
287 
288  dst2 = dist2(posI, posJ);
289  dst2eps = VMlib::boundDenom(dst2, W.getPassport().wakeDiscretizationProperties.eps2); //Сглаживать надо!!!
290 
291  tempVel = { posI[0] - posJ[0], posI[1] - posJ[1] };
292  tempVel *= (gamJ / dst2eps);
293  velI += tempVel;
294  }
295 #ifndef TESTONLYVELO
296  if (calcRadius)
297  {
298  for (size_t s = 0; s < W.getNumberOfBoundary(); ++s)
299  {
300  const auto& bou = W.getBoundary(s);
301  //Модифицируем массив квадратов расстояний до ближайших вихрей из virtualWake
302  for (size_t j = 0; j < bou.virtualWake.vtx.size(); ++j)
303  ModifyE2(ee2, dist2(posI, bou.virtualWake.vtx[j].r()));
304  }
305  }
306 #endif
307 
308  velI *= cft;
309  selfVelo[i] = velI;
310 
311 #ifndef TESTONLYVELO
312  if (calcRadius)
313  domainRadius[i] = 1.0 * sqrt((ee2[0] + ee2[1] + ee2[2]) / 3.0);
314 #endif
315  }
316  }
317  else if (calcRadius)
318  {
319 #pragma omp parallel for default(none) shared(selfVelo, cft, calcVelo, calcRadius, pointsDb, domainRadius) private(tempVel, velI, dst2) schedule(dynamic, DYN_SCHEDULE)
320  for (int i = 0; i < pointsDb.vtx.size(); ++i)
321  {
322  double ee2[3] = { 10000.0, 10000.0, 10000.0 };
323 
324  velI.toZero();
325 
326  const Point2D& posI = pointsDb.vtx[i].r();
327 
328  for (size_t j = 0; j < W.getWake().vtx.size(); ++j)
329  {
330  const Point2D& posJ = W.getWake().vtx[j].r();
331 
332  dst2 = dist2(posI, posJ);
333 
334  //Модифицируем массив квадратов расстояний до ближайших вихрей из wake
335  VMlib::ModifyE2(ee2, dst2);
336  }
337 
338  for (size_t s = 0; s < W.getNumberOfBoundary(); ++s)
339  {
340  for (size_t j = 0; j < W.getBoundary(s).virtualWake.vtx.size(); ++j)
341  {
342  const Point2D& posJ = W.getBoundary(s).virtualWake.vtx[j].r();
343  dst2 = dist2(posI, posJ);
344 
345  //Модифицируем массив квадратов расстояний до ближайших вихрей из virtualWake
346  ModifyE2(ee2, dst2);
347  /*
348  if (dst2 < ee2[0])
349  {
350  size_t pnl = W.getBoundary(s).virtualWake.aflPan[j].second;
351  ee2[0] = ee2[1] = ee2[2] = 0.5 * W.getBoundary(s).afl.len[pnl] / (W.getBoundary(s).vortexBeginEnd[pnl].second - W.getBoundary(s).vortexBeginEnd[pnl].first);
352  }
353  else
354  ModifyE2(ee2, dst2);
355  */
356 
357  }
358  }
359 
360  domainRadius[i] = 1.0 * sqrt((ee2[0] + ee2[1] + ee2[2]) / 3.0);
361  }
362  } //else
363 
364 
365  if (calcVelo)
366  for (size_t i = 0; i < velo.size(); ++i)
367  velo[i] += selfVelo[i];
368 }//CalcConvVeloToSetOfPointsFromWake(...)
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
const double IDPI
Число .
Definition: defs.h:76
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
void ModifyE2(double *ee2, double dst2)
Модифицирует массив квадратов расстояний до ближайших вихрей из wake.
Definition: defs.cpp:233
const WakeDataBase & getSource() const
Возврат константной ссылки на источники в области течения
Definition: World2D.h:207
double accelCft() const
Функция-множитель, позволяющая моделировать разгон
Definition: Passport2D.cpp:51
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
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
Список вихревых элементов
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
void ModifyE2(double *ee2, double dst2)
Модифицирует массив квадратов расстояний до ближайших вихрей из wake.
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
double boundDenom(double r2, double eps2)
Способ сглаживания скорости вихря (вихрь Рэнкина или вихрь Ламба)
Definition: defs.h:536
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
numvector< T, n > & toZero(P val=0)
Установка всех компонент вектора в константу (по умолчанию — нуль)
Definition: numvector.h:524
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
double eps2
Квадрат радиуса вихря
Definition: Passport2D.h:142

Here is the call graph for this function:

Here is the caller graph for this function:

void Velocity::CalcDiffVelo ( )
inherited

Вычисление диффузионных скоростей

Вызывается в CalcVortexVelo()

omp

Definition at line 253 of file Velocity2D.cpp.

254 {
255  W.getTimestat().timeCalcVortexDiffVelo.first += omp_get_wtime();
256 
257  double t12start, t12finish;
258  double t03start, t03finish;
259  double tOtherstart, tOtherfinish;
260 
261  if (W.getPassport().physicalProperties.nu > 0.0)
262  {
263  t12start = omp_get_wtime();
265  t12finish = omp_get_wtime();
266 
267  t03start = omp_get_wtime();
269  t03finish = omp_get_wtime();
270 
271  //контроль застрелов диффузионной скорости
272  tOtherstart = omp_get_wtime();
274 
276  for (size_t bou = 0; bou < virtualVortexesParams.size(); ++bou)
277  LimitDiffVelo(virtualVortexesParams[bou].diffVelo);
278 
279 
280  for (size_t afl = 0; afl < W.getNumberOfAirfoil(); ++afl)
281  for (size_t i = 0; i < W.getAirfoil(afl).viscousStress.size(); ++i)
283 
284  SaveVisStress();
285  tOtherfinish = omp_get_wtime();
286 
287  }
288  W.getTimestat().timeCalcVortexDiffVelo.second += omp_get_wtime();
289 
290  //W.getInfo('t') << "diff_whole = " << W.getTimestat().timeCalcVortexDiffVelo.second - W.getTimestat().timeCalcVortexDiffVelo.first << std::endl;
291  //W.getInfo('t') << "diff_I12 = " << t12finish - t12start << std::endl;
292  //W.getInfo('t') << "diff_I03 = " << t03finish - t03start << std::endl;
293  //W.getInfo('t') << "diff_other = " << tOtherfinish - tOtherstart << std::endl;
294 
295 
296 }// CalcDiffVelo()
void CalcDiffVeloI1I2()
Вычисление диффузионных скоростей вихрей и виртуальных вихрей в вихревом следе
Definition: Velocity2D.cpp:56
std::vector< double > viscousStress
Касательные напряжения на панелях профиля
Definition: Airfoil2D.h:188
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
timePeriod timeCalcVortexDiffVelo
Начало и конец процесса вычисления диффузионных скоростей вихрей
Definition: Times2D.h:89
double nu
Коэффициент кинематической вязкости среды
Definition: Passport2D.h:96
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition: Velocity2D.h:108
void CalcDiffVeloI0I3()
Definition: Velocity2D.cpp:165
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition: Velocity2D.h:105
std::vector< Point2D > diffVelo
Вектор диффузионных скоростей вихрей
Definition: Velocity2D.h:70
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
void LimitDiffVelo(std::vector< Point2D > &diffVel)
Контроль больших значений диффузионных скоростей
Definition: Velocity2D.cpp:239
void SaveVisStress()
Сохранение вязких напряжений
Definition: Velocity2D.cpp:725
Airfoil & getNonConstAirfoil(size_t i) const
Возврат неконстантной ссылки на объект профиля
Definition: World2D.h:142

Here is the call graph for this function:

void Velocity::CalcDiffVeloI0I3 ( )
inherited

Definition at line 165 of file Velocity2D.cpp.

166 {
167  for (size_t afl = 0; afl < W.getNumberOfAirfoil(); ++afl)
168  {
169 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I0I3))
171  {
172  double tt1 = omp_get_wtime();
173  W.getNonConstAirfoil(afl).GPUGetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(W.getWake(), wakeVortexesParams.epsastWake, wakeVortexesParams.I0, wakeVortexesParams.I3);
174  double tt2 = omp_get_wtime();
175  //W.getInfo('t') << "I03_bou->wake = " << (tt2 - tt1) << std::endl;
176  }
177  else
179 #else
181 #endif
182  }
183 
184  //Порядок циклов именно такой, т.к. CUDA оптимизирует вызов ядер,
185  //и на один "слой" виртуальных вихрей считается влияние сразу от всех профилей
186  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
187  {
188  for (size_t afl = 0; afl < W.getNumberOfAirfoil(); ++afl)
189 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I0I3))
190  {
191  double tt1 = omp_get_wtime();
192  W.getNonConstAirfoil(afl).GPUGetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(W.getBoundary(bou).virtualWake, virtualVortexesParams[bou].epsastWake, virtualVortexesParams[bou].I0, virtualVortexesParams[bou].I3);
193  double tt2 = omp_get_wtime();
194  //W.getInfo('t') << "I03_bou->layer = " << (tt2 - tt1) << std::endl;
195  }
196 #else
198 #endif
199  }
200  //влияние поверхности
201  Point2D I3;
202  double I0;
203 
204  double domrad = 0.0;
205 
206 #pragma omp parallel for private(I0, I3, domrad)
207  for (int vt = 0; vt < (int)wakeVortexesParams.diffVelo.size(); ++vt)
208  {
210 
211  wakeVortexesParams.I0[vt] *= domrad;
212  wakeVortexesParams.I0[vt] += DPI * sqr(domrad);
213 
214  I3 = wakeVortexesParams.I3[vt];
215  I0 = wakeVortexesParams.I0[vt];
216 
217  if (fabs(I0) > 1.e-8)
218  wakeVortexesParams.diffVelo[vt] += I3 * (1.0 / I0);
219  }
220 
221  for (size_t targetBou = 0; targetBou < W.getNumberOfBoundary(); ++targetBou)
222 #pragma omp parallel for private(I0, I3, domrad)
223  for (int vt = 0; vt < (int)virtualVortexesParams[targetBou].diffVelo.size(); ++vt)
224  {
225  domrad = std::max(virtualVortexesParams[targetBou].epsastWake[vt], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
226 
227  virtualVortexesParams[targetBou].I0[vt] *= domrad;
228  virtualVortexesParams[targetBou].I0[vt] += DPI * sqr(domrad);
229 
230  I3 = virtualVortexesParams[targetBou].I3[vt];
231  I0 = virtualVortexesParams[targetBou].I0[vt];
232 
233  if (fabs(I0) > 1.e-8)
234  virtualVortexesParams[targetBou].diffVelo[vt] += I3 * (1.0 / I0);
235  }
236 }//CalcDiffVeloI0I3()
std::vector< double > I0
Вектор знаменателей (I0) диффузионных скоростей вихрей (обусловленных профилем)
Definition: Velocity2D.h:73
virtual void GetDiffVelocityI0I3ToWakeAndViscousStresses(const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)=0
std::vector< Point2D > I3
Вектор числителей (I3) диффузионных скоростей вихрей (обусловленных профилем)
Definition: Velocity2D.h:82
std::pair< std::string, int > velocityComputation
Definition: Passport2D.h:190
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
#define CU_I0I3
Definition: Gpudefs.h:78
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition: Velocity2D.h:108
T sqr(T x)
Возведение числа в квадрат
Definition: defs.h:430
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
virtual void GetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)=0
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек, обусловленных геометрией профиля, и вычисление вязкого трения
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition: Passport2D.h:302
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const double DPI
Число .
Definition: defs.h:79
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition: Velocity2D.h:105
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition: Velocity2D.h:85
std::vector< Point2D > diffVelo
Вектор диффузионных скоростей вихрей
Definition: Velocity2D.h:70
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
Airfoil & getNonConstAirfoil(size_t i) const
Возврат неконстантной ссылки на объект профиля
Definition: World2D.h:142
double getMinEpsAst() const
Функция минимально возможного значения для epsAst.
Definition: Passport2D.h:166

Here is the call graph for this function:

Here is the caller graph for this function:

void Velocity::CalcDiffVeloI1I2 ( )
inherited

Вычисление диффузионных скоростей вихрей и виртуальных вихрей в вихревом следе

Вызывает 4 раза функцию CalcDiffVeloToSetOfPoints

Warning
скорости приплюсовываются к уже имеющимся

Definition at line 56 of file Velocity2D.cpp.

57 {
58  // !!! пелена на пелену
59 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I1I2))
61  {
62  double tt1 = omp_get_wtime();
63  //GPUCalcDiffVeloI1I2ToSetOfPointsFromWake(W.getWake(), wakeVortexesParams.epsastWake, W.getWake(), wakeVortexesParams.I1, wakeVortexesParams.I2);
65  double tt2 = omp_get_wtime();
66 
67  //W.getInfo('t') << "Diffusive velocities time = " << int(((tt2 - tt1) * 1000) * 10) / 10.0 << " ms, bodies = " << W.getWake().vtx.size() << std::endl;
68  //W.getInfo('t') << "I12_wake->wake = " << (tt2 - tt1) << std::endl;
69 
70  }
71  else
73 #else
75 #endif
76 
77 
78  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
79  {
80  //не нужно, т.к. сделано выше перед началом вычисления скоростей
81  //W.getNonConstBoundary(bou).virtualWake.WakeSynchronize();
82 
83 
84  //виртуальные на границе на след
85 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I1I2))
87  {
88  double tt1 = omp_get_wtime();
89  GPUCalcDiffVeloI1I2ToSetOfPointsFromSheets(W.getWake(), wakeVortexesParams.epsastWake, W.getBoundary(bou), wakeVortexesParams.I1, wakeVortexesParams.I2);
90  double tt2 = omp_get_wtime();
91  //W.getInfo('t') << "I12_layer->wake = " << (tt2 - tt1) << std::endl;
92  }
93  else
95 #else
97 #endif
98 
99  // !!! след на виртуальные
100 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I1I2))
101  {
102  double tt1 = omp_get_wtime();
103  GPUCalcDiffVeloI1I2ToSetOfPointsFromWake(W.getBoundary(bou).virtualWake, virtualVortexesParams[bou].epsastWake, W.getWake(), virtualVortexesParams[bou].I1, virtualVortexesParams[bou].I2);
104  double tt2 = omp_get_wtime();
105  //W.getInfo('t') << "I12_wake->layer = " << (tt2 - tt1) << std::endl;
106  }
107 #else
109 #endif
110  for (size_t targetBou = 0; targetBou < W.getNumberOfBoundary(); ++targetBou)
111  {
112  // виртуальные на виртуальные
113 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I1I2))
114  {
115  double tt1 = omp_get_wtime();
116  GPUCalcDiffVeloI1I2ToSetOfPointsFromSheets(W.getBoundary(targetBou).virtualWake, virtualVortexesParams[targetBou].epsastWake, W.getBoundary(bou), virtualVortexesParams[targetBou].I1, virtualVortexesParams[targetBou].I2);
117  double tt2 = omp_get_wtime();
118  //W.getInfo('t') << "I12_layer->layer = " << (tt2 - tt1) << std::endl;
119  }
120 #else
122 #endif
123  }
124 
125  } //for bou
126 
127 
128  double tDivstart, tDivfinish;
129 
130  Point2D I2;
131  double I1;
132 
133  tDivstart = omp_get_wtime();
134 
135 #pragma omp parallel for private(I1, I2)
136  for (int vt = 0; vt < (int)wakeVortexesParams.diffVelo.size(); ++vt)
137  {
138  I2 = wakeVortexesParams.I2[vt];
139  I1 = wakeVortexesParams.I1[vt];
140  if (fabs(I1) < 1.e-8)
141  wakeVortexesParams.diffVelo[vt] = { 0.0, 0.0 };
142  else
144  }
145 
146  for (size_t targetBou = 0; targetBou < W.getNumberOfBoundary(); ++targetBou)
147 #pragma omp parallel for private(I1, I2)
148  for (int vt = 0; vt < (int)virtualVortexesParams[targetBou].diffVelo.size(); ++vt)
149  {
150  I2 = virtualVortexesParams[targetBou].I2[vt];
151  I1 = virtualVortexesParams[targetBou].I1[vt];
152 
153  if (fabs(I1) < 1.e-8)
154  virtualVortexesParams[targetBou].diffVelo[vt] = { 0.0, 0.0 };
155  else
156  virtualVortexesParams[targetBou].diffVelo[vt] = I2 * (1.0 / (I1 * std::max(virtualVortexesParams[targetBou].epsastWake[vt], W.getPassport().wakeDiscretizationProperties.getMinEpsAst())));
157  }
158 
159  tDivfinish = omp_get_wtime();
160  //W.getInfo('t') << "I12_Div = " << (tDivfinish - tDivstart) << std::endl;
161 
162 }//CalcDiffVeloI1I2()
std::pair< std::string, int > velocityComputation
Definition: Passport2D.h:190
void CalcDiffVeloI1I2ToSetOfPointsFromWake(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек ...
Definition: Velocity2D.cpp:300
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
virtual void CalcDiffVeloI1I2ToWakeFromWake(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)=0
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition: Velocity2D.h:108
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
void CalcDiffVeloI1I2ToSetOfPointsFromSheets(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)
Definition: Velocity2D.cpp:365
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition: Passport2D.h:302
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< double > I1
Вектор знаменателей (I1) диффузионных скоростей вихрей (обусловленных завихренностью) ...
Definition: Velocity2D.h:76
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition: Velocity2D.h:105
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition: Velocity2D.h:85
std::vector< Point2D > diffVelo
Вектор диффузионных скоростей вихрей
Definition: Velocity2D.h:70
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
std::vector< Point2D > I2
Вектор числителей (I2) диффузионных скоростей вихрей (обусловленных завихренностью) ...
Definition: Velocity2D.h:79
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
virtual void CalcDiffVeloI1I2ToWakeFromSheets(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)=0
double getMinEpsAst() const
Функция минимально возможного значения для epsAst.
Definition: Passport2D.h:166

Here is the call graph for this function:

Here is the caller graph for this function:

void Velocity::CalcDiffVeloI1I2ToSetOfPointsFromSheets ( const WakeDataBase pointsDb,
const std::vector< double > &  domainRadius,
const Boundary bnd,
std::vector< double > &  I1,
std::vector< Point2D > &  I2 
)
inherited
Todo:
Понять природу магической константы 8.0 и синхронизировать с GPU
Todo:
Сделать переменной и синхронизировать с GPU
Todo:
Учитываем пока только нулевой момент решения

Definition at line 365 of file Velocity2D.cpp.

366 {
367  double tCPUSTART, tCPUEND;
368 
369  tCPUSTART = omp_get_wtime();
370 
371  std::vector<double> selfI1(pointsDb.vtx.size(), 0.0);
372  std::vector<Point2D> selfI2(pointsDb.vtx.size(), { 0.0, 0.0 });
373 
374 #pragma warning (push)
375 #pragma warning (disable: 4101)
376  //Локальные переменные для цикла
377  Point2D Rij;
378  double rij, expr;
379  double diffRadius;
380  double left;
381  double right;
382  double posJx;
383  double domRad;
384 #pragma warning (pop)
385 
386 
387 #pragma omp parallel for default(none) shared(selfI1, selfI2, domainRadius, bnd, pointsDb, std::cout) private(Rij, rij, expr, domRad, diffRadius, left, right, posJx)
388  for (int i = 0; i < pointsDb.vtx.size(); ++i)
389  {
390  const Vortex2D& vtxI = pointsDb.vtx[i];
391 
392  domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
393 
395  diffRadius = 8.0 * domRad;
396 
397  left = vtxI.r()[0] - diffRadius;
398  right = vtxI.r()[0] + diffRadius;
399 
400  for (size_t j = 0; j < bnd.afl.getNumberOfPanels(); ++j)
401  {
403  const int nQuadPt = 3;
404 
406  const double ptG = bnd.sheets.freeVortexSheet(j, 0) * bnd.afl.len[j] / nQuadPt;
407 
408  for (int q = 0; q < nQuadPt; ++q)
409  {
410  const Point2D& ptJ = bnd.afl.getR(j) + bnd.afl.tau[j] * (q + 0.5) * bnd.afl.len[j] * (1.0 / nQuadPt); // vorticesDb.vtx[j];
411  posJx = ptJ[0];
412 
413  if ((left < posJx) && (posJx < right))
414  {
415  Rij = vtxI.r() - ptJ;
416  rij = Rij.length();
417  if (rij < diffRadius && rij > 1.e-10)
418  {
419  expr = exp(-rij / domRad);
420  selfI2[i] += (ptG * expr / rij) * Rij;
421  selfI1[i] += ptG * expr;
422  }
423  }//if (rij>1e-6)
424  }
425  }//for j
426  } // for r
427 
428  for (size_t i = 0; i < I1.size(); ++i)
429  {
430  I1[i] += selfI1[i];
431  I2[i] += selfI2[i];
432  }
433 
434  tCPUEND = omp_get_wtime();
435  //W.getInfo('t') << "DIFF_CPU: " << tCPUEND - tCPUSTART << std::endl;
436 }//CalcDiffVeloI1I2ToSetOfPointsFromSheets(...)
const double & freeVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:99
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::vector< Vortex2D > vtx
Список вихревых элементов
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
double getMinEpsAst() const
Функция минимально возможного значения для epsAst.
Definition: Passport2D.h:166
const Airfoil & afl
Definition: Boundary2D.h:76

Here is the call graph for this function:

Here is the caller graph for this function:

void Velocity::CalcDiffVeloI1I2ToSetOfPointsFromWake ( const WakeDataBase pointsDb,
const std::vector< double > &  domainRadius,
const WakeDataBase vorticesDb,
std::vector< double > &  I1,
std::vector< Point2D > &  I2 
)
inherited

Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек

Parameters
[in]pointsDbконстантная ссылка на базу данных пелены из вихрей, в которых надо сосчитать диффузионные скорости
[in]domainRadiusконстантная ссылка на вектор радиусов вихревых доменов
[in]vorticesDbконстантная ссылка на на базу данных пелены из вихрей,от которых надо сосчитать влияния на points
[out]I1ссылка на вектор величин I1 (знаменателей в диффузионных скоростях) в требуемых точках
[out]I2ссылка на вектор величин I2 (числителей в диффузионных скоростях) в требуемых точках
Todo:
Понять природу магической константы 8.0 и синхронизировать с GPU

Definition at line 300 of file Velocity2D.cpp.

301 {
302  double tCPUSTART, tCPUEND;
303 
304  tCPUSTART = omp_get_wtime();
305 
306  std::vector<double> selfI1(pointsDb.vtx.size(), 0.0);
307  std::vector<Point2D> selfI2(pointsDb.vtx.size(), { 0.0, 0.0 });
308 
309 #pragma warning (push)
310 #pragma warning (disable: 4101)
311  //Локальные переменные для цикла
312  Point2D Rij;
313  double rij, expr;
314  double diffRadius, domRad;
315  double left;
316  double right;
317  double posJx;
318 #pragma warning (pop)
319 
320 #pragma omp parallel for default(none) shared(selfI1, selfI2, domainRadius, vorticesDb, pointsDb) private(Rij, rij, expr, diffRadius, domRad, left, right, posJx)
321  for (int i = 0; i < pointsDb.vtx.size(); ++i)
322  {
323  const Vortex2D& vtxI = pointsDb.vtx[i];
324 
325  domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
326 
328  diffRadius = 8.0 * domRad;
329 
330  left = vtxI.r()[0] - diffRadius;
331  right = vtxI.r()[0] + diffRadius;
332 
333  for (size_t j = 0; j < vorticesDb.vtx.size(); ++j)
334  {
335  const Vortex2D& vtxJ = vorticesDb.vtx[j];
336  posJx = vtxJ.r()[0];
337 
338  if ((left < posJx) && (posJx < right))
339  {
340  Rij = vtxI.r() - vtxJ.r();
341  rij = Rij.length();
342  if (rij < diffRadius && rij > 1.e-10)
343  {
344  expr = exp(-rij / domRad);
345  selfI2[i] += (vtxJ.g()* expr / rij) * Rij;
346  selfI1[i] += vtxJ.g()*expr;
347  }
348  }//if (rij>1e-6)
349  }//for j
350  } // for r
351 
352 
353  for (size_t i = 0; i < I1.size(); ++i)
354  {
355  I1[i] += selfI1[i];
356  I2[i] += selfI2[i];
357  }
358 
359  tCPUEND = omp_get_wtime();
360  //W.getInfo('t') << "DIFF_CPU: " << tCPUEND - tCPUSTART << std::endl;
361 }//CalcDiffVeloI1I2ToSetOfPointsFromWake(...)
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
HD double & g()
Функция для доступа к циркуляции вихря
Definition: Vortex2D.h:93
std::vector< Vortex2D > vtx
Список вихревых элементов
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
double getMinEpsAst() const
Функция минимально возможного значения для epsAst.
Definition: Passport2D.h:166

Here is the call graph for this function:

Here is the caller graph for this function:

void VelocityBiotSavart::CalcDiffVeloI1I2ToWakeFromSheets ( const WakeDataBase pointsDb,
const std::vector< double > &  domainRadius,
const Boundary bnd,
std::vector< double > &  I1,
std::vector< Point2D > &  I2 
)
overridevirtual

Implements VM2D::Velocity.

Definition at line 972 of file Velocity2DBiotSavart.cpp.

973 {
974  CalcDiffVeloI1I2ToSetOfPointsFromSheets(pointsDb, domainRadius, bnd, I1, I2);
975 }
void CalcDiffVeloI1I2ToSetOfPointsFromSheets(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)
Definition: Velocity2D.cpp:365

Here is the call graph for this function:

void VelocityBiotSavart::CalcDiffVeloI1I2ToWakeFromWake ( const WakeDataBase pointsDb,
const std::vector< double > &  domainRadius,
const WakeDataBase vorticesDb,
std::vector< double > &  I1,
std::vector< Point2D > &  I2 
)
overridevirtual

Implements VM2D::Velocity.

Definition at line 967 of file Velocity2DBiotSavart.cpp.

968 {
969  CalcDiffVeloI1I2ToSetOfPointsFromWake(pointsDb, domainRadius, vorticesDb, I1, I2);
970 }//CalcDiffVeloI1I2ToWakeFromWake(...)
void CalcDiffVeloI1I2ToSetOfPointsFromWake(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек ...
Definition: Velocity2D.cpp:300

Here is the call graph for this function:

void VelocityBiotSavart::CalcVeloToWakeVP ( )
overridevirtual

Вычисление скоростей в точках wakeVP.

Implements VM2D::Velocity.

Definition at line 163 of file Velocity2DBiotSavart.cpp.

164 {
165  std::vector<Point2D> velConvWake;
166  std::vector<std::vector<Point2D>> velConvBou;
167 
168  int addWSize = (int)W.getMeasureVP().getWakeVP().vtx.size();
169 
170  velConvWake.resize(addWSize, { 0.0, 0.0 });
171 
172  velConvBou.resize(W.getNumberOfBoundary());
173  for (size_t i = 0; i < W.getNumberOfBoundary(); ++i)
174  velConvBou[i].resize(addWSize, { 0.0, 0.0 });
175 
176 #if (defined(USE_CUDA))
177  W.getNonConstCuda().RefreshVP();
178 #endif
179 
180 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_VP))
181  GPUCalcConvVeloToSetOfPointsFromWake(W.getMeasureVP().getWakeVP(), velConvWake, W.getNonConstMeasureVP().getNonConstDomainRadius(), true, false);
182 #else
184 #endif
185 
186  for (size_t i = 0; i < W.getNumberOfBoundary(); ++i)
187 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_VP))
188  W.getNonConstBoundary(i).GPUCalcConvVelocityToSetOfPointsFromSheets(W.getMeasureVP().getWakeVP(), velConvBou[i]);
189 #else
191 #endif
192 
193  std::vector<Point2D>& velocityRef = W.getNonConstMeasureVP().getNonConstVelocity();
194  velocityRef.assign(addWSize, W.getPassport().physicalProperties.V0());
195 
196  for (int i = 0; i < addWSize; ++i)
197  velocityRef[i] += velConvWake[i];
198 
199  for (size_t bou = 0; bou < velConvBou.size(); ++bou)
200  for (int j = 0; j < addWSize; ++j)
201  velocityRef[j] += velConvBou[bou][j];
202 
203 }// CalcVeloToWakeVP()
Boundary & getNonConstBoundary(size_t i) const
Возврат неконстантной ссылки на объект граничного условия
Definition: World2D.h:159
Gpu & getNonConstCuda() const
Возврат неконстантной ссылки на объект, связанный с видеокартой (GPU)
Definition: World2D.h:232
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
#define CU_VP
Definition: Gpudefs.h:83
void CalcConvVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius)
Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа ...
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
virtual void CalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase &pointsDb, std::vector< Point2D > &velo) const =0
Вычисление конвективных скоростей в наборе точек, вызываемых наличием слоев вихрей и источников на пр...
const MeasureVP & getMeasureVP() const
Возврат константной ссылки на measureVP.
Definition: World2D.h:175
std::vector< Vortex2D > vtx
Список вихревых элементов
std::vector< Point2D > & getNonConstVelocity()
Возврат velocity.
Definition: MeasureVP2D.h:130
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
MeasureVP & getNonConstMeasureVP() const
Возврат неконстантной ссылки на measureVP.
Definition: World2D.h:180
Point2D V0() const
Функция скорости набегающего потока с учетом разгона
Definition: Passport2D.h:93
const WakeDataBase & getWakeVP() const
Возврат wakeVP.
Definition: MeasureVP2D.h:120
std::vector< double > & getNonConstDomainRadius()
Возврат domainRadius.
Definition: MeasureVP2D.h:140

Here is the call graph for this function:

Here is the caller graph for this function:

void VelocityBiotSavart::FillRhs ( Eigen::VectorXd &  rhs) const
overridevirtual

Расчет вектора правой части (всего)

Implements VM2D::Velocity.

Definition at line 832 of file Velocity2DBiotSavart.cpp.

833 {
834  Eigen::VectorXd locRhs;
835  std::vector<double> lastRhs(W.getNumberOfBoundary());
836 
837  size_t currentRow = 0;
838 
839  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
840  {
841  //double tt0 = omp_get_wtime();
842 
843  const Airfoil& afl = W.getAirfoil(bou);
844  size_t np = afl.getNumberOfPanels();
845 
846  size_t nVars;
847 
848  nVars = W.getBoundary(bou).GetUnknownsSize();
849  locRhs.resize(nVars);
850 
851 
852  std::vector<double> wakeRhs;
853 
854  double tt1 = omp_get_wtime();
855 
856 #if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_RHS))
857  //GPUGetWakeInfluenceToRhs(afl, wakeRhs);
858  GPUFASTGetWakeInfluenceToRhs(afl, wakeRhs);
859  //GetWakeInfluenceToRhs(afl, wakeRhs);
860 
861 #else
862  GetWakeInfluenceToRhs(afl, wakeRhs);
863 #endif
864 
865 
866  double tt2 = omp_get_wtime();
867  //W.getInfo('t') << "Rhs time = " << int(((tt2 - tt1) * 1000) * 10) / 10.0 << " ms, bodies = " << W.getWake().vtx.size() << std::endl;
868 
869 
870  std::vector<double> vInfRhs;
871  afl.GetInfluenceFromVInfToPanel(vInfRhs);
872 
873  /*
874  //if (W.currentStep == 200)
875  {
876  std::stringstream ss;
877  ss << "rhsWake-" << W.currentStep;
878  std::ofstream of(W.getPassport().dir + "dbg/" + ss.str());
879  for (size_t i = 0; i < wakeRhs.size(); ++i)
880  of << wakeRhs[i] << std::endl;
881  of.close();
882  }
883  */
884 
885 
886 
887 #pragma omp parallel for \
888  default(none) \
889  shared(locRhs, afl, bou, wakeRhs, vInfRhs, np) \
890  schedule(dynamic, DYN_SCHEDULE)
891  for (int i = 0; i < (int)afl.getNumberOfPanels(); ++i)
892  {
893  locRhs(i) = -vInfRhs[i] - wakeRhs[i] + 0.25 * ((afl.getV(i) + afl.getV(i + 1)) & afl.tau[i]); //0.25 * (afl.getV(i) + afl.getV(i + 1))*afl.tau[i] - прямолинейные
894  if (W.getBoundary(bou).sheetDim > 1)
895  locRhs(np + i) = -vInfRhs[np + i] - wakeRhs[np + i];
896 
897  //влияние присоединенных слоев от самого себя и от других профилей
898  for (size_t q = 0; q < W.getNumberOfBoundary(); ++q)
899  {
900  const auto& sht = W.getBoundary(q).sheets;
901  const auto& iq = W.getIQ(bou, q);
902 
903  const Airfoil& aflOther = W.getAirfoil(q);
905  {
906  for (size_t j = 0; j < aflOther.getNumberOfPanels(); ++j)
907  {
908  if ((i != j) || (bou != q))
909  {
910  locRhs(i) += -iq.first(i, j) * sht.attachedVortexSheet(j, 0);
911  locRhs(i) += -iq.second(i, j) * sht.attachedSourceSheet(j, 0); // getIQ(bou, q).second(i, j) пока забито нулем для криволинейных
912  }//if (i != j)
913  }//for j
914  }
915  }//for q
916  }//for i
917 
918  lastRhs[bou] = 0.0;
919 
920  for (size_t q = 0; q < afl.gammaThrough.size(); ++q)
921  lastRhs[bou] += afl.gammaThrough[q];
922 
923 
924 
925 
927 
928  double dwcm; //wcm0, wcm1, wcmOld0, wcmOld1, dwc0, dwc1;
929  const double currT = W.getPassport().timeDiscretizationProperties.currTime;
930  const double dt = W.getPassport().timeDiscretizationProperties.dt;
931 
932  //wcm0 = W.getNonConstMechanics(bou).AngularVelocityOfAirfoil(currT);
933  //wcmOld0 = W.getNonConstMechanics(bou).AngularVelocityOfAirfoil(currT - dt);
934  //dwc0 = wcm0 - wcmOld0;
935 
937  lastRhs[bou] += 2.0 * dwcm * W.getAirfoil(bou).area * (W.getAirfoil(bou).inverse ? 1.0 : -1.0);
938 
939  /*
940  if (bou == 0)
941  {
942  lastRhs[bou] -= (2.0 * PI * 0.5) * dwc0 * 0.5;
943  }
944 
945  if (bou == 1)
946  {
947  wcm1 = W.getNonConstMechanics(bou).AngularVelocityOfAirfoil(currT);
948  wcmOld1 = W.getNonConstMechanics(bou).AngularVelocityOfAirfoil(currT - dt);
949  dwc1 = wcm1 - wcmOld1;
950  //std::cout << "dwc0,1 = " << dwc0 << " " << dwc1 << std::endl;
951  lastRhs[bou] += (2.0 * PI * 1.0) * (dwc1) * 1.0;
952  }
953  */
955 
956  //размазываем правую часть
957  for (size_t i = 0; i < nVars; ++i)
958  rhs(i + currentRow) = locRhs(i);
959 
960  rhs(currentRow + nVars) = lastRhs[bou];
961 
962  currentRow += nVars + 1;
963  }// for bou
964 }
const std::pair< Eigen::MatrixXd, Eigen::MatrixXd > & getIQ(size_t i, size_t j) const
Возврат константной ссылки на объект, связанный с матрицей интегралов от (r-xi)/|r-xi|^2.
Definition: World2D.h:237
double area
Площадь профиля
Definition: Airfoil2D.h:83
bool inverse
Признак разворота нормалей (для расчета внутренних течений)
Definition: Airfoil2D.h:139
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
double currTime
Текущее время
Definition: PassportGen.h:56
std::vector< double > gammaThrough
Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
Definition: Airfoil2D.h:196
Mechanics & getNonConstMechanics(size_t i) const
Возврат неконстантной ссылки на объект механики
Definition: World2D.h:192
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
const Point2D & getV(size_t q) const
Возврат константной ссылки на скорость вершины профиля
Definition: Airfoil2D.h:113
const bool isMoves
Переменная, отвечающая за то, двигается профиль или нет
Definition: Mechanics2D.h:130
const Mechanics & getMechanics(size_t i) const
Возврат константной ссылки на объект механики
Definition: World2D.h:186
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
double dt
Шаг по времени
Definition: PassportGen.h:65
void GetWakeInfluenceToRhs(const Airfoil &afl, std::vector< double > &wakeRhs) const
Генерация вектора влияния вихревого следа на профиль
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
virtual double AngularVelocityOfAirfoil(double currTime)=0
Вычисление угловой скорости профиля
std::vector< Point2D > tau
Касательные к панелям профиля
Definition: Airfoil2D.h:182
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
virtual void GetInfluenceFromVInfToPanel(std::vector< double > &vInfRhs) const =0
Вычисление влияния набегающего потока на панель для правой части
const bool isDeform
Переменная, отвечающая за то, деформируется профиль или нет
Definition: Mechanics2D.h:133
Абстрактный класс, определяющий обтекаемый профиль
Definition: Airfoil2D.h:60
size_t sheetDim
Размерность параметров каждого из слоев на каждой из панелей
Definition: Boundary2D.h:92
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
size_t GetUnknownsSize() const
Возврат размерности вектора решения
Definition: Boundary2D.cpp:71

Here is the call graph for this function:

void VelocityBiotSavart::GetWakeInfluenceToRhs ( const Airfoil afl,
std::vector< double > &  wakeRhs 
) const

Генерация вектора влияния вихревого следа на профиль

Генерирует вектор влияния вихревого следа на профиль, используемый затем для расчета вектора правой части.

Parameters
[in]aflконстантная ссылка на профиль, правая часть для которого вычисляется
[out]wakeRhsссылка на вектор влияния вихревого следа на ОДИН профиль

Definition at line 641 of file Velocity2DBiotSavart.cpp.

642 {
643  size_t np = afl.getNumberOfPanels();
644  size_t shDim = W.getBoundary(afl.numberInPassport).sheetDim;
645 
646  wakeRhs.resize(W.getBoundary(afl.numberInPassport).GetUnknownsSize());
647 
648  //локальные переменные для цикла
649  std::vector<double> velI(shDim, 0.0);
650 
651 #pragma omp parallel for default(none) shared(shDim, afl, np, wakeRhs, IDPI) private(velI)
652  for (int i = 0; i < np; ++i)
653  {
654  velI.assign(shDim, 0.0);
655 
656  if (W.getWake().vtx.size() > 0)
657  {
658  //Учет влияния следа
659  afl.GetInfluenceFromVorticesToPanel(i, W.getWake().vtx.data(), W.getWake().vtx.size(), velI);
660  }
661 
662  if (W.getSource().vtx.size() > 0)
663  {
664  //Учет влияния источников
665  afl.GetInfluenceFromSourcesToPanel(i, W.getSource().vtx.data(), W.getSource().vtx.size(), velI);
666  }
667 
668  for (size_t j = 0; j < shDim; ++j)
669  velI[j] *= IDPI / afl.len[i];
670 
671  wakeRhs[i] = velI[0];
672 
673  if (shDim != 1)
674  wakeRhs[np + i] = velI[1];
675  }//for i
676 }//GetWakeInfluenceToRhs(...)
const double IDPI
Число .
Definition: defs.h:76
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
if(currentStep%10==0)
Definition: gammaCirc.h:22
const size_t numberInPassport
Номер профиля в паспорте
Definition: Airfoil2D.h:74
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
const WakeDataBase & getSource() const
Возврат константной ссылки на источники в области течения
Definition: World2D.h:207
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
std::vector< Vortex2D > vtx
Список вихревых элементов
virtual void GetInfluenceFromVorticesToPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const =0
Вычисление влияния части подряд идущих вихрей из вихревого следа на панель для правой части ...
size_t sheetDim
Размерность параметров каждого из слоев на каждой из панелей
Definition: Boundary2D.h:92
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
virtual void GetInfluenceFromSourcesToPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const =0
Вычисление влияния части подряд идущих источников из области течения на панель для правой части ...
size_t GetUnknownsSize() const
Возврат размерности вектора решения
Definition: Boundary2D.cpp:71

Here is the call graph for this function:

Here is the caller graph for this function:

void Velocity::LimitDiffVelo ( std::vector< Point2D > &  diffVel)
inherited

Контроль больших значений диффузионных скоростей

Parameters
[in,out]diffVelссылка на вектор диффузионных скоростей

Definition at line 239 of file Velocity2D.cpp.

240 {
241  for (size_t i = 0; i < diffVel.size(); ++i)
242  {
243  Point2D& diffV = diffVel[i];
244 
245  diffV *= W.getPassport().physicalProperties.nu;
246 
247  if (diffV.length() > 1.5 * W.getPassport().physicalProperties.vRef)
249  }
250 }
double vRef
Референсная скорость
Definition: Passport2D.h:81
double nu
Коэффициент кинематической вязкости среды
Definition: Passport2D.h:96
P length() const
Вычисление 2-нормы (длины) вектора
Definition: numvector.h:371
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
void normalize(P newlen=1.0)
Нормирование вектора на заданную длину
Definition: numvector.h:416
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222

Here is the call graph for this function:

Here is the caller graph for this function:

void Velocity::ResizeAndZero ( )
inherited

Очистка старых массивов под хранение скоростей, выделение новой памяти и обнуление

Вызывается в CalcVortexVelo() на каждом шаге расчета перед непосредственно расчетом скоростей

Definition at line 661 of file Velocity2D.cpp.

662 {
663  W.getTimestat().timeCalcVortexConvVelo.first += omp_get_wtime();
665  wakeVortexesParams.convVelo.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
666  W.getTimestat().timeCalcVortexConvVelo.second += omp_get_wtime();
667 
668  W.getTimestat().timeCalcVortexDiffVelo.first += omp_get_wtime();
669  wakeVortexesParams.I0.clear();
670  wakeVortexesParams.I0.resize(W.getWake().vtx.size(), 0.0);
671 
672  wakeVortexesParams.I1.clear();
673  wakeVortexesParams.I1.resize(W.getWake().vtx.size(), 0.0);
674 
675  wakeVortexesParams.I2.clear();
676  wakeVortexesParams.I2.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
677 
678  wakeVortexesParams.I3.clear();
679  wakeVortexesParams.I3.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
680 
682  wakeVortexesParams.diffVelo.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
683 
685  wakeVortexesParams.epsastWake.resize(W.getWake().vtx.size(), 0.0);
686  W.getTimestat().timeCalcVortexDiffVelo.second += omp_get_wtime();
687 
688  //Создаем массивы под виртуальные вихри
689  W.getTimestat().timeCalcVortexConvVelo.first += omp_get_wtime();
690  virtualVortexesParams.clear();
692  W.getTimestat().timeCalcVortexConvVelo.second += omp_get_wtime();
693 
694  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
695  {
696  W.getTimestat().timeCalcVortexConvVelo.first += omp_get_wtime();
697  virtualVortexesParams[bou].convVelo.clear();
698  virtualVortexesParams[bou].convVelo.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
699  W.getTimestat().timeCalcVortexConvVelo.second += omp_get_wtime();
700 
701  W.getTimestat().timeCalcVortexDiffVelo.first += omp_get_wtime();
702  virtualVortexesParams[bou].I0.clear();
703  virtualVortexesParams[bou].I0.resize(W.getBoundary(bou).virtualWake.vtx.size(), 0.0);
704 
705  virtualVortexesParams[bou].I1.clear();
706  virtualVortexesParams[bou].I1.resize(W.getBoundary(bou).virtualWake.vtx.size(), 0.0);
707 
708  virtualVortexesParams[bou].I2.clear();
709  virtualVortexesParams[bou].I2.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
710 
711  virtualVortexesParams[bou].I3.clear();
712  virtualVortexesParams[bou].I3.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
713 
714  virtualVortexesParams[bou].diffVelo.clear();
715  virtualVortexesParams[bou].diffVelo.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
716 
717  virtualVortexesParams[bou].epsastWake.clear();
718  virtualVortexesParams[bou].epsastWake.resize(W.getBoundary(bou).virtualWake.vtx.size(), 0.0);
719 
720  W.getTimestat().timeCalcVortexDiffVelo.second += omp_get_wtime();
721  }
722 }//ResizeAndZero()
std::vector< double > I0
Вектор знаменателей (I0) диффузионных скоростей вихрей (обусловленных профилем)
Definition: Velocity2D.h:73
std::vector< Point2D > I3
Вектор числителей (I3) диффузионных скоростей вихрей (обусловленных профилем)
Definition: Velocity2D.h:82
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
timePeriod timeCalcVortexDiffVelo
Начало и конец процесса вычисления диффузионных скоростей вихрей
Definition: Times2D.h:89
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
timePeriod timeCalcVortexConvVelo
Начало и конец процесса вычисления конвективных скоростей вихрей
Definition: Times2D.h:86
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition: Velocity2D.h:108
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
std::vector< Vortex2D > vtx
Список вихревых элементов
std::vector< double > I1
Вектор знаменателей (I1) диффузионных скоростей вихрей (обусловленных завихренностью) ...
Definition: Velocity2D.h:76
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition: Velocity2D.h:105
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition: Velocity2D.h:85
std::vector< Point2D > diffVelo
Вектор диффузионных скоростей вихрей
Definition: Velocity2D.h:70
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
std::vector< Point2D > I2
Вектор числителей (I2) диффузионных скоростей вихрей (обусловленных завихренностью) ...
Definition: Velocity2D.h:79
std::vector< Point2D > convVelo
Вектор конвективных скоростей вихрей
Definition: Velocity2D.h:67

Here is the call graph for this function:

void Velocity::SaveVisStress ( )
inherited

Сохранение вязких напряжений

Вызывается в CalcDiffVelo()

Definition at line 725 of file Velocity2D.cpp.

726 {
728  // if ((W.getPassport().timeDiscretizationProperties.saveVTK > 0) && (W.ifDivisible(10)) && (W.getNumberOfAirfoil() > 0))
729  {
730  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
731  {
732  if (q == 0)
733  VMlib::CreateDirectory(W.getPassport().dir, "visStress");
734 
735  std::stringstream ss;
736  ss << "VisStress_" << q << "-";
737  std::string fname = VMlib::fileNameStep(ss.str(), W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "txt");
738  std::ofstream outfile;
739  outfile.open(W.getPassport().dir + "visStress/" + fname);
740 
741  outfile << W.getAirfoil(q).viscousStress.size() << std::endl; //Сохранение числа вихрей в пелене
742 
743  for (size_t i = 0; i < W.getAirfoil(q).viscousStress.size(); ++i)
744  {
745  const Point2D& r = 0.5 * (W.getAirfoil(q).getR(i + 1) + W.getAirfoil(q).getR(i));
746  double gi = W.getAirfoil(q).viscousStress[i];
747  outfile << static_cast<int>(i) << " " << r[0] << " " << r[1] << " " << gi << std::endl;
748  }//for i
749  outfile.close();
750  }
751  }
752 }
std::vector< double > viscousStress
Касательные напряжения на панелях профиля
Definition: Airfoil2D.h:188
int saveVtxStep
Шаг сохранения кадров в бинарные файлы
Definition: PassportGen.h:73
int nameLength
Число разрядов в имени файла
Definition: PassportGen.h:68
std::string dir
Рабочий каталог задачи
Definition: PassportGen.h:118
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
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
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
int saveVisStress
Шаг вычисления и сохранения скорости и давления
Definition: PassportGen.h:81
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition: WorldGen.h:91
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
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

std::vector<VortexesParams> VM2D::Velocity::virtualVortexesParams
inherited

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

Definition at line 108 of file Velocity2D.h.

const World2D& VM2D::Velocity::W
protectedinherited

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

Definition at line 101 of file Velocity2D.h.

VortexesParams VM2D::Velocity::wakeVortexesParams
inherited

Струтура, определяющая параметры вихрей в следе

Definition at line 105 of file Velocity2D.h.


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