VM2D  1.12
Vortex methods for 2D flows simulation
VM2D::Velocity Class Referenceabstract

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

#include <Velocity2D.h>

Inheritance diagram for VM2D::Velocity:
Collaboration diagram for VM2D::Velocity:

Public Member Functions

 Velocity (const World2D &W_)
 Конструктор More...
 
virtual void CalcConvVelo ()=0
 Вычисление конвективных скоростей вихрей и виртуальных вихрей в вихревом следе, а также в точках wakeVP. More...
 
virtual void CalcVeloToWakeVP ()=0
 Вычисление скоростей в точках wakeVP. More...
 
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)
 
virtual void CalcDiffVeloI1I2ToWakeFromSheets (const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)=0
 
virtual void CalcDiffVeloI1I2ToWakeFromWake (const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)=0
 
void CalcDiffVeloI1I2 ()
 Вычисление диффузионных скоростей вихрей и виртуальных вихрей в вихревом следе More...
 
void CalcDiffVeloI0I3 ()
 
void LimitDiffVelo (std::vector< Point2D > &diffVel)
 Контроль больших значений диффузионных скоростей More...
 
void CalcDiffVelo ()
 Вычисление диффузионных скоростей More...
 
void ResizeAndZero ()
 Очистка старых массивов под хранение скоростей, выделение новой памяти и обнуление More...
 
void SaveVisStress ()
 Сохранение вязких напряжений More...
 
virtual void FillRhs (Eigen::VectorXd &rhs) const =0
 Расчет вектора правой части (всего) More...
 
virtual ~Velocity ()
 Деструктор More...
 

Public Attributes

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

Protected Attributes

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

Detailed Description

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

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

Definition at line 97 of file Velocity2D.h.

Constructor & Destructor Documentation

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

Конструктор

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

Definition at line 114 of file Velocity2D.h.

114  :
115  W(W_)
116  {
117  virtualVortexesParams.resize(0);
118  };
const World2D & W
Константная ссылка на решаемую задачу
Definition: Velocity2D.h:101
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition: Velocity2D.h:108
virtual VM2D::Velocity::~Velocity ( )
inlinevirtual

Деструктор

Definition at line 180 of file Velocity2D.h.

180 { };

Member Function Documentation

virtual void VM2D::Velocity::CalcConvVelo ( )
pure virtual

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

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

Implemented in VM2D::VelocityBiotSavart.

void Velocity::CalcDiffVelo ( )

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

Вызывается в 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 ( )

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 ( )

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

Вызывает 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 
)
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 
)

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

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:

virtual void VM2D::Velocity::CalcDiffVeloI1I2ToWakeFromSheets ( const WakeDataBase pointsDb,
const std::vector< double > &  domainRadius,
const Boundary bnd,
std::vector< double > &  I1,
std::vector< Point2D > &  I2 
)
pure virtual

Implemented in VM2D::VelocityBiotSavart.

Here is the caller graph for this function:

virtual void VM2D::Velocity::CalcDiffVeloI1I2ToWakeFromWake ( const WakeDataBase pointsDb,
const std::vector< double > &  domainRadius,
const WakeDataBase vorticesDb,
std::vector< double > &  I1,
std::vector< Point2D > &  I2 
)
pure virtual

Implemented in VM2D::VelocityBiotSavart.

Here is the caller graph for this function:

virtual void VM2D::Velocity::CalcVeloToWakeVP ( )
pure virtual

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

Implemented in VM2D::VelocityBiotSavart.

virtual void VM2D::Velocity::FillRhs ( Eigen::VectorXd &  rhs) const
pure virtual

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

Implemented in VM2D::VelocityBiotSavart.

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

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

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 ( )

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

Вызывается в 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 ( )

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

Вызывается в 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

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

Definition at line 108 of file Velocity2D.h.

const World2D& VM2D::Velocity::W
protected

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

Definition at line 101 of file Velocity2D.h.

VortexesParams VM2D::Velocity::wakeVortexesParams

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

Definition at line 105 of file Velocity2D.h.


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