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

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

#include <Velocity2DBarnesHut.h>

Inheritance diagram for VM2D::VelocityBarnesHut:
Collaboration diagram for VM2D::VelocityBarnesHut:

Public Member Functions

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

Public Attributes

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

Protected Attributes

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

Detailed Description

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

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

  • в соответствии с гибридным быстрым методом Барнса — Савара / мультиполей
Author
Марчевский Илья Константинович
Сокол Ксения Сергеевна
Рятина Евгения Павловна
Колганова Александра Олеговна

\Version 1.14

Date
6 марта 2026 г.

Definition at line 64 of file Velocity2DBarnesHut.h.

Constructor & Destructor Documentation

◆ VelocityBarnesHut()

VelocityBarnesHut::VelocityBarnesHut ( const World2D W_)

Конструктор

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

Definition at line 55 of file Velocity2DBarnesHut.cpp.

55 :
56 Velocity(W_)
57{
58};
Абстрактный класс, определяющий способ вычисления скоростей
Definition Velocity2D.h:105

◆ ~VelocityBarnesHut()

VelocityBarnesHut::~VelocityBarnesHut ( )
virtual

Деструктор

Definition at line 61 of file Velocity2DBarnesHut.cpp.

62{
63};

Member Function Documentation

◆ CalcConvVelo()

void Velocity::CalcConvVelo ( )
inherited

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

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

Definition at line 888 of file Velocity2D.cpp.

889{
890 W.getTimers().start("ConvVel");
891
892 //Влияние следа на след
893#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_CONV_TOWAKE))
894 //Вызов виртуальной функции
895
898
899 std::unique_ptr<BHcu::CudaTreeInfo>& cntrTree = W.getNonConstCuda().cntrTreeWake;
900 auto& treeWake = *W.getNonConstCuda().inflTreeWake;
902 {
903 cntrTree->MemoryAllocate((int)W.getCuda().n_CUDA_wake);
904 cntrTree->Update((int)W.getWake().vtx.size(), W.getWake().devVtxPtr, W.getPassport().wakeDiscretizationProperties.eps);
905 cntrTree->Build();
906
907 //Перестроение дерева для вихрей //todo //временно для корректного расчета epsast
908 treeWake.MemoryAllocate((int)W.getCuda().n_CUDA_wake);
909 treeWake.Update((int)W.getWake().vtx.size(), W.getWake().devVtxPtr, W.getPassport().wakeDiscretizationProperties.eps);
910 treeWake.Build();
911 treeWake.UpwardTraversal(W.getPassport().numericalSchemes.nbodyMultipoleOrder);
912 }
913 GPUCalcConvVeloToSetOfPointsFromWake(W.getNonConstCuda().cntrTreeWake, W.getWake(), wakeVortexesParams.convVelo, wakeVortexesParams.epsastWake, true, true);
914
916#else
917 //Вызов виртуальной функции
919#endif
920
921 std::vector<Point2D> nullVector(0);
922 //вычисление конвективных скоростей по закону Био-Савара от виртуальных вихрей (точнее, от слоев)
923 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
924 {
925 //Влияние на след от панелей
926#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_CONVVIRT))
927
929 {
930 case 0:
931 W.getBoundary(bou).GPUCalcConvVelocityToSetOfPointsFromSheets(W.getWake(), wakeVortexesParams.convVelo);
932 break;
933 case 1:
934 if (bou == 0)
935 GPUCalcConvVelocityToSetOfPointsFromSheets(W.getNonConstCuda().cntrTreeWake, W.getWake(), wakeVortexesParams.convVelo);
936 break;
937 }
938 //
939
940#else
942#endif
943 }
944
945 //Скорости только что рожденных вихрей
946 {
947 size_t nVirtVortices = 0;
948 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
949 nVirtVortices += W.getBoundary(bou).virtualWake.vtx.size();
950
951 size_t counter = wakeVortexesParams.convVelo.size() - nVirtVortices;
952
953 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
954 for (size_t v = 0; v < W.getBoundary(bou).virtualWake.vtx.size(); ++v)
955 {
958 - W.getV0();
959 ++counter;
960 }
961 }
962
963 //Заполнение структур данных виртуальных вихрей
964 {
965 size_t nVirtVortices = 0;
966 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
967 nVirtVortices += W.getBoundary(bou).virtualWake.vtx.size();
968
969 size_t counter = wakeVortexesParams.convVelo.size() - nVirtVortices;
970
971 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
972 for (size_t v = 0; v < W.getBoundary(bou).virtualWake.vtx.size(); ++v)
973 {
974 virtualVortexesParams[bou].convVelo[v] = wakeVortexesParams.convVelo[counter];
975 virtualVortexesParams[bou].epsastWake[v] = wakeVortexesParams.epsastWake[counter];
976 ++counter;
977 }
978 }
979
980
981 W.getTimers().stop("ConvVel");
982 W.getTimers().start("VelPres");
983
984 //вычисление скоростей в заданных точках только на соответствующем шаге по времени
986 //optimizer
987 // if ((W.getCurrentStep() >= 0 && W.getCurrentStep() <= 12000) || W.getCurrentStep() == 0)
989
991 /*if ((W.getPassport().timeDiscretizationProperties.saveVPstep > 0) && (!(W.getCurrentStep() % W.getPassport().timeDiscretizationProperties.saveVPstep)))
992 {
993 int currentStep = W.getCurrentStep();
994 if (currentStep >= 1100 && currentStep <= 1300)
995 {
996 CalcVeloToWakeVP();
997 optimizedVelocity.AddVelocities(W.getNonConstMeasureVP().getNonConstVelocity());
998 }
999
1000 if (currentStep == 1200)
1001 {
1002 optimizedVelocity.PerformAveragingVelo(W.getPassport().dir);
1003 }
1004 }*/
1006
1007 W.getTimers().stop("VelPres");
1008}//CalcConvVelo()
const Point2D & getV(size_t q) const
Возврат константной ссылки на скорость вершины профиля
Definition Airfoil2D.h:137
const Airfoil & afl
Definition Boundary2D.h:77
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition Boundary2D.h:86
virtual void CalcConvVelocityToSetOfPointsFromSheets(const WakeDataBase &pointsDb, std::vector< Point2D > &velo) const =0
Вычисление конвективных скоростей в наборе точек, вызываемых наличием слоев вихрей и источников на пр...
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition Passport2D.h:292
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition Passport2D.h:295
std::vector< VortexesParams > virtualVortexesParams
Вектор струтур, определяющий параметры виртуальных вихрей для профилей
Definition Velocity2D.h:115
virtual void CalcConvVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius)=0
Вычисление конвективных скоростей и радиусов вихревых доменов в заданном наборе точек от следа
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition Velocity2D.h:112
virtual void CalcVeloToWakeVP()
Вычисление скоростей в точках wakeVP.
const World2D & W
Константная ссылка на решаемую задачу
Definition Velocity2D.h:108
std::vector< Point2D > vecHalfGamma
Скорость вихрей виртуального следа конкретного профиля (равна Gamma/2) используется для расчета давле...
std::vector< std::pair< size_t, size_t > > aflPan
Пара чисел: номер профиля и номер панели, на которой рожден виртуальный вихрь
std::vector< Vortex2D > vtx
Список вихревых элементов
VMlib::vmTimer timerConvVelo
Definition World2D.h:357
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition World2D.h:226
Gpu & getNonConstCuda() const
Возврат неконстантной ссылки на объект, связанный с видеокартой (GPU)
Definition World2D.h:266
VMlib::TimersGen & getTimers() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition World2D.h:276
Point2D getV0() const
Возврат текущей скорости набегающего потока
Definition World2D.h:142
const Gpu & getCuda() const
Возврат константной ссылки на объект, связанный с видеокартой (GPU)
Definition World2D.h:261
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition World2D.h:251
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition World2D.h:180
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition World2D.h:191
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
void stop(const std::string &timerLabel)
Останов счетчика
Definition TimesGen.cpp:68
void start(const std::string &timerLabel)
Запуск счетчика
Definition TimesGen.cpp:55
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition WorldGen.h:99
const vmTimer & stop() const
Останов работающего счетчика времени
Definition TimesGen.h:125
const vmTimer & start() const
Запуск (первый или повторный) счетчика времени
Definition TimesGen.h:109
const vmTimer & reset() const
Сброс счетчика времени
Definition TimesGen.h:101
std::pair< std::string, int > velocityComputation
Definition Passport2D.h:180
std::vector< Point2D > convVelo
Вектор конвективных скоростей вихрей
Definition Velocity2D.h:72
std::vector< double > epsastWake
Вектор характерных радиусов вихревых доменов (eps*)
Definition Velocity2D.h:90
double eps
Радиус вихря
Definition Passport2D.h:126
int saveVPstep
Шаг вычисления и сохранения скорости и давления
Definition PassportGen.h:80
Here is the call graph for this function:

◆ CalcConvVeloToSetOfPointsFromWake()

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

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

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

Implements VM2D::Velocity.

Definition at line 66 of file Velocity2DBarnesHut.cpp.

67{
68 BH::params prm;
71
72 int nCores = omp_get_max_threads();
75
76 //Выбор неочевиден
77 //&&&
78 prm.NumOfLevelsVortex = (int)log2(pointsDb.vtx.size());
79
80 BH::BarnesHut BH(prm, pointsDb.vtx);
81 double timeBuildTree = 0.0, timeRect = 0.0, timeSummarization = 0.0, timeComputation = 0.0;
82 BH.BuildOneTree(BH.treeVrt, prm.NumOfLevelsVortex, BH.pointsCopyVrt, timeBuildTree);
83 BH.BuildEnclosingRectangle(timeRect);
84 BH.InfluenceComputation(velo, domainRadius, timeSummarization, timeComputation, calcRadius);
85 //std::cout << timeBuildTree << " " << timeRect << " " << timeSummarization << " " << timeComputation << std::endl;
86 //std::cout << "Tot_time = " << timeBuildTree + timeRect + timeSummarization + timeComputation << std::endl;
87
88#ifdef USE_CUDA
89 if (calcVelo)
90 W.getCuda().CopyMemToDev<double, 2>(velo.size(), (double*)velo.data(), pointsDb.devVelPtr);
91 if (calcRadius)
92 W.getCuda().CopyMemToDev<double, 1>(domainRadius.size(), domainRadius.data(), pointsDb.devRadPtr);
93#endif
94
95
96}
Класс, определяющий основной алгоритм модификации метода Барнса - Хата
Definition BarnesHut.h:59
Класс, содержащий параметры метода Баонса - Хата для CPU.
Definition ParamsBH.h:65
double eps
Радиус вихревого элемента
Definition ParamsBH.h:68
int NumOfLevelsVortex
Максимальное количество уровней дерева
Definition ParamsBH.h:81
double theta
Параметр точности
Definition ParamsBH.h:84
int order
Точность расчета скоростей:
Definition ParamsBH.h:78
double eps2
Квадрат радиуса вихревого элемента
Definition ParamsBH.h:71
double eps2
Квадрат радиуса вихря
Definition Passport2D.h:129
Here is the call graph for this function:

◆ CalcConvVPVeloToSetOfPointsFromWake()

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

Reimplemented from VM2D::Velocity.

Definition at line 101 of file Velocity2DBarnesHut.cpp.

102{
103 std::vector<Point2D> selfVelo(pointsDb.vtx.size());
104 domainRadius.resize(pointsDb.vtx.size());
105
106 double cft = IDPI;
107
108#pragma warning (push)
109#pragma warning (disable: 4101)
110 //Локальные переменные для цикла
111 Point2D velI;
112 Point2D tempVel;
113 double dst2eps, dst2;
114#pragma warning (pop)
115
117
118 if (calcVelo)
119 {
120#pragma omp parallel for default(none) shared(selfVelo, cft, calcVelo, calcRadius, eps2, pointsDb, domainRadius) private(tempVel, velI, dst2, dst2eps) schedule(dynamic, DYN_SCHEDULE)
121 for (int i = 0; i < pointsDb.vtx.size(); ++i)
122 {
123 double ee2[3] = { 10000.0, 10000.0, 10000.0 };
124
125 velI.toZero();
126
127 const Point2D& posI = pointsDb.vtx[i].r();
128
129 for (size_t j = 0; j < W.getWake().vtx.size(); ++j)
130 {
131 const Point2D& posJ = W.getWake().vtx[j].r();
132
133 dst2 = (posI - posJ).length2();
134
135// //Модифицируем массив квадратов расстояний до ближайших вихрей из wake
136//#ifndef TESTONLYVELO
137// if (calcRadius)
138// VMlib::ModifyE2(ee2, dst2);
139//#endif //!TESTONLYVELO
140
141 const double& gamJ = W.getWake().vtx[j].g();
142
143 tempVel.toZero();
144 dst2eps = VMlib::boundDenom(dst2, eps2); //Сглаживать надо!!!
145
146 tempVel = { -posI[1] + posJ[1], posI[0] - posJ[0] };
147 tempVel *= (gamJ / dst2eps);
148 velI += tempVel;
149 }
150
151 for (size_t j = 0; j < W.getSource().vtx.size(); ++j)
152 {
153 const Point2D& posJ = W.getSource().vtx[j].r();
154 const double& gamJ = W.getPassport().physicalProperties.accelCft(W.getCurrentTime()) * W.getSource().vtx[j].g();
155
156 tempVel.toZero();
157
158 dst2 = dist2(posI, posJ);
159 dst2eps = VMlib::boundDenom(dst2, W.getPassport().wakeDiscretizationProperties.eps2); //Сглаживать надо!!!
160
161 tempVel = { posI[0] - posJ[0], posI[1] - posJ[1] };
162 tempVel *= (gamJ / dst2eps);
163 velI += tempVel;
164 }
165//#ifndef TESTONLYVELO
166// if (calcRadius)
167// {
168// for (size_t s = 0; s < W.getNumberOfBoundary(); ++s)
169// {
170// const auto& bou = W.getBoundary(s);
171// //Модифицируем массив квадратов расстояний до ближайших вихрей из virtualWake
172// for (size_t j = 0; j < bou.virtualWake.vtx.size(); ++j)
173// ModifyE2(ee2, dist2(posI, bou.virtualWake.vtx[j].r()));
174// }
175// }
176//#endif
177
178 velI *= cft;
179 selfVelo[i] = velI;
180
181//#ifndef TESTONLYVELO
182// if (calcRadius)
183// domainRadius[i] = 1.0 * sqrt((ee2[0] + ee2[1] + ee2[2]) / 3.0);
184//#endif
185 }
186 }
187// else if (calcRadius)
188// {
189//#pragma omp parallel for default(none) shared(selfVelo, cft, calcVelo, calcRadius, pointsDb, domainRadius) private(tempVel, velI, dst2) schedule(dynamic, DYN_SCHEDULE)
190// for (int i = 0; i < pointsDb.vtx.size(); ++i)
191// {
192// double ee2[3] = { 10000.0, 10000.0, 10000.0 };
193//
194// velI.toZero();
195//
196// const Point2D& posI = pointsDb.vtx[i].r();
197//
198// for (size_t j = 0; j < W.getWake().vtx.size(); ++j)
199// {
200// const Point2D& posJ = W.getWake().vtx[j].r();
201//
202// dst2 = dist2(posI, posJ);
203//
204// //Модифицируем массив квадратов расстояний до ближайших вихрей из wake
205// VMlib::ModifyE2(ee2, dst2);
206// }
207//
208// for (size_t s = 0; s < W.getNumberOfBoundary(); ++s)
209// {
210// for (size_t j = 0; j < W.getBoundary(s).virtualWake.vtx.size(); ++j)
211// {
212// const Point2D& posJ = W.getBoundary(s).virtualWake.vtx[j].r();
213// dst2 = dist2(posI, posJ);
214//
215// //Модифицируем массив квадратов расстояний до ближайших вихрей из virtualWake
216// ModifyE2(ee2, dst2);
217// }
218// }
219//
220// domainRadius[i] = 1.0 * sqrt((ee2[0] + ee2[1] + ee2[2]) / 3.0);
221// }
222// } //else
223
224
225 if (calcVelo)
226 for (size_t i = 0; i < velo.size(); ++i)
227 velo[i] += selfVelo[i];
228
229#ifdef USE_CUDA
230 if (calcVelo)
231 W.getCuda().CopyMemToDev<double, 2>(velo.size(), (double*)velo.data(), pointsDb.devVelPtr);
232// if (calcRadius)
233// W.getCuda().CopyMemToDev<double, 1>(domainRadius.size(), domainRadius.data(), pointsDb.devRadPtr);
234#endif
235
236}//CalcConvVPVeloToSetOfPointsFromWake(...)
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition Passport2D.h:289
const WakeDataBase & getSource() const
Возврат константной ссылки на источники в области течения
Definition World2D.h:236
double getCurrentTime() const
Definition WorldGen.h:100
numvector< T, n > & toZero(P val=0)
Установка всех компонент вектора в константу (по умолчанию — нуль)
Definition numvector.h:527
const double IDPI
Число .
Definition defs.h:85
double boundDenom(double r2, double eps2)
Способ сглаживания скорости вихря (вихрь Рэнкина или вихрь Ламба)
Definition defs.h:563
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:835
double accelCft(double currentTime) const
Функция-множитель, позволяющая моделировать разгон
Here is the call graph for this function:

◆ CalcDiffVelo()

void Velocity::CalcDiffVelo ( )
inherited

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

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

Definition at line 148 of file Velocity2D.cpp.

149{
150 W.getTimers().start("DiffVel");
151
152 double t12start, t12finish;
153 double t03start, t03finish;
154 double tOtherstart, tOtherfinish;
155
157 {
158 t12start = omp_get_wtime();
160 t12finish = omp_get_wtime();
161
162 t03start = omp_get_wtime();
164 t03finish = omp_get_wtime();
165
166 tOtherstart = omp_get_wtime();
167 for (Point2D& diffV : wakeVortexesParams.diffVelo)
168 diffV *= W.getPassport().physicalProperties.nu;
169
170 //контроль застрелов диффузионной скорости
172
173 for (size_t afl = 0; afl < W.getNumberOfAirfoil(); ++afl)
174 for (size_t i = 0; i < W.getAirfoil(afl).viscousStress.size(); ++i)
176
177 //SaveVisStress();
178
179
180 //Заполнение структуры данных виртуальных вихрей
181 size_t nVirtVortices = 0;
182 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
183 nVirtVortices += W.getBoundary(bou).virtualWake.vtx.size();
184
185 size_t counter = wakeVortexesParams.convVelo.size() - nVirtVortices;
186
187 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
188 for (size_t v = 0; v < W.getBoundary(bou).virtualWake.vtx.size(); ++v)
189 {
190 virtualVortexesParams[bou].diffVelo[v] = wakeVortexesParams.diffVelo[counter];
191 virtualVortexesParams[bou].I0[v] = wakeVortexesParams.I0[counter];
192 virtualVortexesParams[bou].I1[v] = wakeVortexesParams.I1[counter];
193 virtualVortexesParams[bou].I2[v] = wakeVortexesParams.I2[counter];
194 virtualVortexesParams[bou].I3[v] = wakeVortexesParams.I3[counter];
195 ++counter;
196 }
197
198 tOtherfinish = omp_get_wtime();
199 }
200 W.getTimers().stop("DiffVel");
201
202}// CalcDiffVelo()
std::vector< double > viscousStress
Нейросеть для коэффициентов I0 и I3 диффузионной скорости
Definition Airfoil2D.h:268
void LimitDiffVelo(std::vector< Point2D > &diffVel)
Контроль больших значений диффузионных скоростей
void CalcDiffVeloI0I3()
void CalcDiffVeloI1I2()
Вычисление диффузионных скоростей вихрей и виртуальных вихрей в вихревом следе
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition World2D.h:174
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition World2D.h:157
Airfoil & getNonConstAirfoil(size_t i) const
Возврат неконстантной ссылки на объект профиля
Definition World2D.h:169
double nu
Коэффициент кинематической вязкости среды
Definition Passport2D.h:99
std::vector< Point2D > I2
Вектор числителей (I2) диффузионных скоростей вихрей (обусловленных завихренностью)
Definition Velocity2D.h:84
std::vector< double > I0
Вектор знаменателей (I0) диффузионных скоростей вихрей (обусловленных профилем)
Definition Velocity2D.h:78
std::vector< Point2D > diffVelo
Вектор диффузионных скоростей вихрей
Definition Velocity2D.h:75
std::vector< Point2D > I3
Вектор числителей (I3) диффузионных скоростей вихрей (обусловленных профилем)
Definition Velocity2D.h:87
std::vector< double > I1
Вектор знаменателей (I1) диффузионных скоростей вихрей (обусловленных завихренностью)
Definition Velocity2D.h:81
Here is the call graph for this function:

◆ CalcDiffVeloI0I3()

void Velocity::CalcDiffVeloI0I3 ( )
inherited

Definition at line 99 of file Velocity2D.cpp.

100{
101 for (size_t afl = 0; afl < W.getNumberOfAirfoil(); ++afl)
102 {
103#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I0I3))
104 W.getNonConstAirfoil(afl).GPUGetDiffVelocityI0I3ToSetOfPointsAndViscousStresses(wakeVortexesParams.epsastWake, wakeVortexesParams.I0, wakeVortexesParams.I3);
105#else
107#endif
108 }
109
110
111 //влияние поверхности
112 Point2D I3;
113 double I0;
114
115 double domrad = 0.0;
116
117#pragma omp parallel for private(I0, I3, domrad)
118 for (int vt = 0; vt < (int)wakeVortexesParams.diffVelo.size(); ++vt)
119 {
121
122 wakeVortexesParams.I0[vt] *= domrad;
123 wakeVortexesParams.I0[vt] += DPI * sqr(domrad);
124
125
126 I3 = wakeVortexesParams.I3[vt];
127 I0 = wakeVortexesParams.I0[vt];
128
129
130 if (fabs(I0) > 1.e-8)
131 wakeVortexesParams.diffVelo[vt] += I3 * (1.0 / I0);
132 }
133}//CalcDiffVeloI0I3()
virtual void GetDiffVelocityI0I3ToWakeAndViscousStresses(const WakeDataBase &pointsDb, std::vector< double > &domainRadius, std::vector< double > &I0, std::vector< Point2D > &I3)
const double DPI
Число .
Definition defs.h:88
T sqr(T x)
Умножение a на комплексно сопряженноe к b.
Definition defsBH.h:101
double getMinEpsAst() const
Функция минимально возможного значения для epsAst.
Definition Passport2D.h:153
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcDiffVeloI1I2()

void Velocity::CalcDiffVeloI1I2 ( )
inherited

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

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

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

Definition at line 55 of file Velocity2D.cpp.

56{
57 // !!! пелена на пелену
58#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I1I2))
60 GPUCalcDiffVeloI1I2ToSetOfPointsFromWake(W.getWake(), wakeVortexesParams.epsastWake, W.getWake(), wakeVortexesParams.I1, wakeVortexesParams.I2);
61 else
63#else
65#endif
66
67
68 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
69 {
70 //не нужно, т.к. сделано выше перед началом вычисления скоростей
71 //W.getNonConstBoundary(bou).virtualWake.WakeSynchronize();
72
73
74 //виртуальные на границе на след
75#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_I1I2))
76 GPUCalcDiffVeloI1I2ToSetOfPointsFromSheets(W.getWake(), wakeVortexesParams.epsastWake, W.getBoundary(bou), wakeVortexesParams.I1, wakeVortexesParams.I2);
77#else
79#endif
80 } //for bou
81
82 Point2D I2;
83 double I1;
84
85#pragma omp parallel for private(I1, I2)
86 for (int vt = 0; vt < (int)wakeVortexesParams.diffVelo.size(); ++vt)
87 {
88 I2 = wakeVortexesParams.I2[vt];
89 I1 = wakeVortexesParams.I1[vt];
90
91 if (fabs(I1) < 1.e-8)
92 wakeVortexesParams.diffVelo[vt] = { 0.0, 0.0 };
93 else
95 }
96}//CalcDiffVeloI1I2()
void CalcDiffVeloI1I2ToWakeFromSheets(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)
void CalcDiffVeloI1I2ToWakeFromWake(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcDiffVeloI1I2ToSetOfPointsFromSheets()

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 264 of file Velocity2D.cpp.

265{
266 double tCPUSTART, tCPUEND;
267
268 tCPUSTART = omp_get_wtime();
269
270 std::vector<double> selfI1(pointsDb.vtx.size(), 0.0);
271 std::vector<Point2D> selfI2(pointsDb.vtx.size(), { 0.0, 0.0 });
272
273#pragma warning (push)
274#pragma warning (disable: 4101)
275 //Локальные переменные для цикла
276 Point2D Rij;
277 double rij, expr;
278 double diffRadius;
279 double left;
280 double right;
281 double posJx;
282 double domRad;
283#pragma warning (pop)
284
285
286#pragma omp parallel for default(none) shared(selfI1, selfI2, domainRadius, bnd, pointsDb, std::cout) private(Rij, rij, expr, domRad, diffRadius, left, right, posJx)
287 for (int i = 0; i < pointsDb.vtx.size(); ++i)
288 {
289 const Vortex2D& vtxI = pointsDb.vtx[i];
290
291 domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
292
294 diffRadius = 8.0 * domRad;
295
296 left = vtxI.r()[0] - diffRadius;
297 right = vtxI.r()[0] + diffRadius;
298
299 for (size_t j = 0; j < bnd.afl.getNumberOfPanels(); ++j)
300 {
302 const int nQuadPt = 3;
303
305 const double ptG = bnd.sheets.freeVortexSheet(j, 0) * bnd.afl.len[j] / nQuadPt;
306
307 for (int q = 0; q < nQuadPt; ++q)
308 {
309 const Point2D& ptJ = bnd.afl.getR(j) + bnd.afl.tau[j] * (q + 0.5) * bnd.afl.len[j] * (1.0 / nQuadPt); // vorticesDb.vtx[j];
310 posJx = ptJ[0];
311
312 if ((left < posJx) && (posJx < right))
313 {
314 Rij = vtxI.r() - ptJ;
315 rij = Rij.length();
316 if (rij < diffRadius && rij > 1.e-10)
317 {
318 expr = exp(-rij / domRad);
319 selfI2[i] += (ptG * expr / rij) * Rij;
320 selfI1[i] += ptG * expr;
321 }
322 }//if (rij>1e-6)
323 }
324 }//for j
325 } // for r
326
327 for (size_t i = 0; i < I1.size(); ++i)
328 {
329 I1[i] += selfI1[i];
330 I2[i] += selfI2[i];
331 }
332
333 tCPUEND = omp_get_wtime();
334 //W.getInfo('t') << "DIFF_CPU: " << tCPUEND - tCPUSTART << std::endl;
335}//CalcDiffVeloI1I2ToSetOfPointsFromSheets(...)
std::vector< double > len
Длины панелей профиля
Definition Airfoil2D.h:94
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition Airfoil2D.h:113
std::vector< Point2D > tau
Касательные к панелям профиля
Definition Airfoil2D.h:91
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition Airfoil2D.h:163
Sheet sheets
Слои на профиле
Definition Boundary2D.h:96
const double & freeVortexSheet(size_t n, size_t moment) const
Definition Sheet2D.h:100
Класс, опеделяющий двумерный вихревой элемент
Definition Vortex2D.h:59
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition Vortex2D.h:87
P length() const
Вычисление 2-нормы (длины) вектора
Definition numvector.h:374
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcDiffVeloI1I2ToSetOfPointsFromWake()

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 206 of file Velocity2D.cpp.

207{
208 std::vector<double> selfI1(pointsDb.vtx.size(), 0.0);
209 std::vector<Point2D> selfI2(pointsDb.vtx.size(), { 0.0, 0.0 });
210
211#pragma warning (push)
212#pragma warning (disable: 4101)
213 //Локальные переменные для цикла
214 Point2D Rij;
215 double rij, expr;
216 double diffRadius, domRad;
217 double left;
218 double right;
219 double posJx;
220#pragma warning (pop)
221
222#pragma omp parallel for default(none) shared(selfI1, selfI2, domainRadius, vorticesDb, pointsDb) private(Rij, rij, expr, diffRadius, domRad, left, right, posJx)
223 for (int i = 0; i < pointsDb.vtx.size(); ++i)
224 {
225 const Vortex2D& vtxI = pointsDb.vtx[i];
226
227 domRad = std::max(domainRadius[i], W.getPassport().wakeDiscretizationProperties.getMinEpsAst());
228
230 diffRadius = 8.0 * domRad;
231
232 left = vtxI.r()[0] - diffRadius;
233 right = vtxI.r()[0] + diffRadius;
234
235 for (size_t j = 0; j < vorticesDb.vtx.size(); ++j)
236 {
237 const Vortex2D& vtxJ = vorticesDb.vtx[j];
238 posJx = vtxJ.r()[0];
239
240 if ((left < posJx) && (posJx < right))
241 {
242 Rij = vtxI.r() - vtxJ.r();
243 rij = Rij.length();
244 if (rij < diffRadius && rij > 1.e-10)
245 {
246 expr = exp(-rij / domRad);
247 selfI2[i] += (vtxJ.g()* expr / rij) * Rij;
248 selfI1[i] += vtxJ.g()*expr;
249 }
250 }//if (rij>1e-6)
251 }//for j
252 } // for r
253
254
255 for (size_t i = 0; i < I1.size(); ++i)
256 {
257 I1[i] += selfI1[i];
258 I2[i] += selfI2[i];
259 }
260}//CalcDiffVeloI1I2ToSetOfPointsFromWake(...)
HD double & g()
Функция для доступа к циркуляции вихря
Definition Vortex2D.h:95
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcDiffVeloI1I2ToWakeFromSheets()

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

Definition at line 1093 of file Velocity2D.cpp.

1094{
1095 CalcDiffVeloI1I2ToSetOfPointsFromSheets(pointsDb, domainRadius, bnd, I1, I2);
1096}
void CalcDiffVeloI1I2ToSetOfPointsFromSheets(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const Boundary &bnd, std::vector< double > &I1, std::vector< Point2D > &I2)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcDiffVeloI1I2ToWakeFromWake()

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

Definition at line 1088 of file Velocity2D.cpp.

1089{
1090 CalcDiffVeloI1I2ToSetOfPointsFromWake(pointsDb, domainRadius, vorticesDb, I1, I2);
1091}//CalcDiffVeloI1I2ToWakeFromWake(...)
void CalcDiffVeloI1I2ToSetOfPointsFromWake(const WakeDataBase &pointsDb, const std::vector< double > &domainRadius, const WakeDataBase &vorticesDb, std::vector< double > &I1, std::vector< Point2D > &I2)
Вычисление числителей и знаменателей диффузионных скоростей в заданном наборе точек
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcVeloToWakeVP()

void Velocity::CalcVeloToWakeVP ( )
virtualinherited

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

Definition at line 1013 of file Velocity2D.cpp.

1014{
1015 std::vector<Point2D> velConvWake;
1016 std::vector<std::vector<Point2D>> velConvBou;
1017
1018 int addWSize = (int)W.getMeasureVP().getWakeVP().vtx.size();
1019
1020 velConvWake.resize(addWSize, { 0.0, 0.0 });
1021
1022 velConvBou.resize(W.getNumberOfBoundary());
1023 for (size_t i = 0; i < W.getNumberOfBoundary(); ++i)
1024 velConvBou[i].resize(addWSize, { 0.0, 0.0 });
1025
1026 std::vector<Point2D>& velocityRef = W.getNonConstMeasureVP().getNonConstVelocity();
1027 velocityRef.assign(addWSize, W.getV0());
1028
1029
1030#if (defined(USE_CUDA))
1031 W.getNonConstCuda().RefreshVP();
1032#endif
1033
1034#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_VP))
1035 std::unique_ptr<BHcu::CudaTreeInfo>& cntrTree = W.getNonConstCuda().cntrTreeVP;
1037 {
1039 {
1040 cntrTree->MemoryAllocate((int)W.getCuda().n_CUDA_velVP);
1041 cntrTree->Update((int)W.getMeasureVP().getWakeVP().vtx.size(), W.getMeasureVP().getWakeVP().devVtxPtr, 0.0);
1042 cntrTree->Build();
1043 }
1044 }
1045 GPUCalcConvVeloToSetOfPointsFromWake(cntrTree, W.getMeasureVP().getWakeVP(), velocityRef, W.getNonConstMeasureVP().getNonConstDomainRadius(), true, false);
1046#else
1048 {
1049 case 0:
1051 break;
1052 case 1:
1054 }
1055#endif
1056
1057 for (size_t i = 0; i < W.getNumberOfBoundary(); ++i)
1058#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_VP))
1060 {
1061 case 0:
1062 W.getNonConstBoundary(i).GPUCalcConvVelocityToSetOfPointsFromSheets(W.getMeasureVP().getWakeVP(), velConvBou[i]);
1063
1064 for (int s = 0; s < addWSize; ++s)
1065 velocityRef[s] += velConvWake[s];
1066
1067 for (size_t bou = 0; bou < velConvBou.size(); ++bou)
1068 for (int j = 0; j < addWSize; ++j)
1069 velocityRef[j] += velConvBou[bou][j];
1070 break;
1071 case 1:
1072 if (i == 0)
1073 GPUCalcConvVelocityToSetOfPointsFromSheets(cntrTree, W.getMeasureVP().getWakeVP(), velocityRef);
1074 break;
1075 }
1076#else
1078 for (int i = 0; i < addWSize; ++i)
1079 velocityRef[i] += velConvWake[i];
1080
1081 for (size_t bou = 0; bou < velConvBou.size(); ++bou)
1082 for (int j = 0; j < addWSize; ++j)
1083 velocityRef[j] += velConvBou[bou][j];
1084#endif
1085}// CalcVeloToWakeVP()
#define CU_VP
Definition Gpudefs.h:86
const WakeDataBase & getWakeVP() const
Возврат wakeVP.
std::vector< Point2D > & getNonConstVelocity()
Возврат velocity.
std::vector< double > & getNonConstDomainRadius()
Возврат domainRadius.
virtual void CalcConvVPVeloToSetOfPointsFromWake(const WakeDataBase &pointsDb, std::vector< Point2D > &velo, std::vector< double > &domainRadius, bool calcVelo, bool calcRadius)
Definition Velocity2D.h:141
bool isAnyMovableOrDeformable() const
Возврат признака того, что хотя бы один из профилей подвижный или деформируемый
Definition World2D.cpp:1765
MeasureVP & getNonConstMeasureVP() const
Возврат неконстантной ссылки на measureVP.
Definition World2D.h:207
Boundary & getNonConstBoundary(size_t i) const
Возврат неконстантной ссылки на объект граничного условия
Definition World2D.h:186
const MeasureVP & getMeasureVP() const
Возврат константной ссылки на measureVP.
Definition World2D.h:202
Here is the call graph for this function:
Here is the caller graph for this function:

◆ FillRhs()

void Velocity::FillRhs ( Eigen::VectorXd &  rhsReord) const
inherited

!!!!!!!!!!влияние присоединенных слоев от самого себя и от других профилей

Definition at line 781 of file Velocity2D.cpp.

782{
783 Eigen::VectorXd locRhs;
784 std::vector<double> lastRhs(W.getNumberOfBoundary());
785
786 size_t currentRow = 0;
787 size_t currentSkosRow = 0;
788
789 size_t nAllVars = 0;
790 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
791 nAllVars += W.getBoundary(bou).GetUnknownsSize();
792
793 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
794 {
795 const Airfoil& afl = W.getAirfoil(bou);
796 size_t np = afl.getNumberOfPanels();
797
798 size_t nVars;
799
800 nVars = W.getBoundary(bou).GetUnknownsSize();
801 locRhs.resize(nVars);
802
803
804 std::vector<double> wakeRhs;
805
806 double tt1 = omp_get_wtime();
807
808#if (defined(__CUDACC__) || defined(USE_CUDA)) && (defined(CU_RHS))
809 W.timerRhs.reset();
810 W.timerRhs.start();
812 GPUFASTGetWakeInfluenceToRhs(afl, wakeRhs);
813 else
814 GPUGetWakeInfluenceToRhs(afl, wakeRhs);
815 W.timerRhs.stop();
816#else
817 GetWakeInfluenceToRhs(afl, wakeRhs);
818#endif
819 std::vector<double> vInfRhs;
820 afl.GetInfluenceFromVInfToPanel(vInfRhs);
821
822#pragma omp parallel for \
823 default(none) \
824 shared(locRhs, afl, bou, wakeRhs, vInfRhs, np) \
825 schedule(dynamic, DYN_SCHEDULE)
826 for (int i = 0; i < (int)afl.getNumberOfPanels(); ++i)
827 {
828 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] - прямолинейные
829 if (W.getBoundary(bou).sheetDim > 1)
830 locRhs(np + i) = -vInfRhs[np + i] - wakeRhs[np + i];
831
833 //Если включен быстрый метод и есть подвижные тела, то пока что все равно вычисляем IQ
835 {
836 for (size_t q = 0; q < W.getNumberOfBoundary(); ++q)
837 {
838 const auto& sht = W.getBoundary(q).sheets;
839 const auto& iq = W.getIQ(bou, q);
840
841 const Airfoil& aflOther = W.getAirfoil(q);
842 if (W.getMechanics(q).isMoves)
843 {
844 for (size_t j = 0; j < aflOther.getNumberOfPanels(); ++j)
845 {
846 if ((i != j) || (bou != q))
847 {
848 locRhs(i) += -iq.first(i, j) * sht.attachedVortexSheet(j, 0);
849 locRhs(i) += -iq.second(i, j) * sht.attachedSourceSheet(j, 0);
850 }//if (i != j)
851
852 }//for j
853 }
854 }//for q
855 }
856 }//for i
857
858 lastRhs[bou] = 0.0;
859
860#pragma omp for
861 for (int q = 0; q < (int)afl.gammaThrough.size(); ++q)
862 lastRhs[bou] += afl.gammaThrough[q];
863
865
866 const double currT = W.getCurrentTime();
867 const double dt = W.getPassport().timeDiscretizationProperties.dt;
868
869 lastRhs[bou] += (W.getMechanics(bou).circulation - W.getMechanics(bou).circulationOld) * (W.getAirfoil(bou).inverse ? 1.0 : -1.0);
870
872
873 //размазываем правую часть
874 for (size_t i = 0; i < nVars; ++i)
875 rhsReord(i + currentSkosRow) = locRhs(i);
876
877 rhsReord(nAllVars + bou) = lastRhs[bou];
878
879 currentRow += nVars + 1;
880 currentSkosRow += nVars;
881 }// for bou
882}
bool inverse
Признак разворота нормалей (для расчета внутренних течений)
Definition Airfoil2D.h:76
Абстрактный класс, определяющий обтекаемый профиль
Definition Airfoil2D.h:182
std::vector< double > gammaThrough
Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
Definition Airfoil2D.h:276
virtual void GetInfluenceFromVInfToPanel(std::vector< double > &vInfRhs) const
Вычисление влияния набегающего потока на панель для правой части
size_t GetUnknownsSize() const
Возврат размерности вектора решения
size_t sheetDim
Размерность параметров каждого из слоев на каждой из панелей
Definition Boundary2D.h:93
const bool isMoves
Переменная, отвечающая за то, двигается профиль или нет
double circulationOld
Циркуляция скорости по границе профиля с предыдущего шага
double circulation
Текущая циркуляция скорости по границе профиля
void GetWakeInfluenceToRhs(const Airfoil &afl, std::vector< double > &wakeRhs) const
Генерация вектора влияния вихревого следа на профиль
const Mechanics & getMechanics(size_t i) const
Возврат константной ссылки на объект механики
Definition World2D.h:213
const std::pair< Eigen::MatrixXd, Eigen::MatrixXd > & getIQ(size_t i, size_t j) const
Возврат константной ссылки на объект, связанный с матрицей интегралов от (r-xi)/|r-xi|^2.
Definition World2D.h:271
VMlib::vmTimer timerRhs
Definition World2D.h:353
double dt
Шаг по времени
Definition PassportGen.h:67
Here is the call graph for this function:

◆ GetWakeInfluenceToRhs()

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

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

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

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

Definition at line 612 of file Velocity2D.cpp.

613{
614 size_t np = afl.getNumberOfPanels();
615 size_t shDim = W.getBoundary(afl.numberInPassport).sheetDim;
616
617 wakeRhs.resize(W.getBoundary(afl.numberInPassport).GetUnknownsSize());
618
619 //локальные переменные для цикла
620 std::vector<double> velI(shDim, 0.0);
621
622#pragma omp parallel for default(none) shared(shDim, afl, np, wakeRhs, IDPI) private(velI)
623 for (int i = 0; i < np; ++i)
624 {
625 velI.assign(shDim, 0.0);
626
627 if (W.getWake().vtx.size() > 0)
628 {
629 //Учет влияния следа
630 afl.GetInfluenceFromVorticesToPanel(i, W.getWake().vtx.data(), W.getWake().vtx.size(), velI);
631 }
632
633 if (W.getSource().vtx.size() > 0)
634 {
635 //Учет влияния источников
636 afl.GetInfluenceFromSourcesToPanel(i, W.getSource().vtx.data(), W.getSource().vtx.size(), velI);
637 }
638
639 for (size_t j = 0; j < shDim; ++j)
640 velI[j] *= IDPI / afl.len[i];
641
642 wakeRhs[i] = velI[0];
643
644 if (shDim != 1)
645 wakeRhs[np + i] = velI[1];
646 }//for i
647}//GetWakeInfluenceToRhs(...)
virtual void GetInfluenceFromSourcesToPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const
Вычисление влияния части подряд идущих источников из области течения на панель для правой части
virtual void GetInfluenceFromVorticesToPanel(size_t panel, const Vortex2D *ptr, ptrdiff_t count, std::vector< double > &panelRhs) const
Вычисление влияния части подряд идущих вихрей из вихревого следа на панель для правой части
const size_t numberInPassport
Номер профиля в паспорте
Definition Airfoil2D.h:188
Here is the call graph for this function:
Here is the caller graph for this function:

◆ LimitDiffVelo()

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

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

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

Definition at line 136 of file Velocity2D.cpp.

137{
138 for (size_t i = 0; i < diffVel.size(); ++i)
139 {
140 Point2D& diffV = diffVel[i];
141
142 if (diffV.length() > 1.5 * W.getPassport().physicalProperties.vRef)
144 }
145}
void normalize(P newlen=1.0)
Нормирование вектора на заданную длину
Definition numvector.h:419
double vRef
Референсная скорость
Definition Passport2D.h:81
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ResizeAndZero()

void Velocity::ResizeAndZero ( )
inherited

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

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

Definition at line 516 of file Velocity2D.cpp.

517{
518 W.getTimers().start("ConvVel");
520 wakeVortexesParams.convVelo.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
521 W.getTimers().stop("ConvVel");
522
523 W.getTimers().start("DiffVel");
524 wakeVortexesParams.I0.clear();
525 wakeVortexesParams.I0.resize(W.getWake().vtx.size(), 0.0);
526
527 wakeVortexesParams.I1.clear();
528 wakeVortexesParams.I1.resize(W.getWake().vtx.size(), 0.0);
529
530 wakeVortexesParams.I2.clear();
531 wakeVortexesParams.I2.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
532
533 wakeVortexesParams.I3.clear();
534 wakeVortexesParams.I3.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
535
537 wakeVortexesParams.diffVelo.resize(W.getWake().vtx.size(), { 0.0, 0.0 });
538
540 wakeVortexesParams.epsastWake.resize(W.getWake().vtx.size(), 0.0);
541 W.getTimers().stop("DiffVel");
542
543 //Создаем массивы под виртуальные вихри
544 W.getTimers().start("ConvVel");
545 virtualVortexesParams.clear();
547 W.getTimers().stop("ConvVel");
548
549 for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
550 {
551 W.getTimers().start("ConvVel");
552 virtualVortexesParams[bou].convVelo.clear();
553 virtualVortexesParams[bou].convVelo.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
554 W.getTimers().stop("ConvVel");
555
556 W.getTimers().start("DiffVel");
557 virtualVortexesParams[bou].I0.clear();
558 virtualVortexesParams[bou].I0.resize(W.getBoundary(bou).virtualWake.vtx.size(), 0.0);
559
560 virtualVortexesParams[bou].I1.clear();
561 virtualVortexesParams[bou].I1.resize(W.getBoundary(bou).virtualWake.vtx.size(), 0.0);
562
563 virtualVortexesParams[bou].I2.clear();
564 virtualVortexesParams[bou].I2.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
565
566 virtualVortexesParams[bou].I3.clear();
567 virtualVortexesParams[bou].I3.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
568
569 virtualVortexesParams[bou].diffVelo.clear();
570 virtualVortexesParams[bou].diffVelo.resize(W.getBoundary(bou).virtualWake.vtx.size(), { 0.0, 0.0 });
571
572 virtualVortexesParams[bou].epsastWake.clear();
573 virtualVortexesParams[bou].epsastWake.resize(W.getBoundary(bou).virtualWake.vtx.size(), 0.0);
574
575 W.getTimers().stop("DiffVel");
576 }
577}//ResizeAndZero()
Here is the call graph for this function:

◆ SaveVisStress()

void Velocity::SaveVisStress ( )
inherited

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

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

Definition at line 580 of file Velocity2D.cpp.

581{
583 // if ((W.getPassport().timeDiscretizationProperties.saveVTK > 0) && (W.ifDivisible(10)) && (W.getNumberOfAirfoil() > 0))
584 {
585 for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
586 {
587 if (q == 0)
588 VMlib::CreateDirectory(W.getPassport().dir, "visStress");
589
590 std::stringstream ss;
591 ss << "VisStress_" << q << "-";
593 std::ofstream outfile;
594 outfile.open(W.getPassport().dir + "visStress/" + fname);
595
596 outfile << W.getAirfoil(q).viscousStress.size() << std::endl; //Сохранение числа вихрей в пелене
597
598 for (size_t i = 0; i < W.getAirfoil(q).viscousStress.size(); ++i)
599 {
600 const Point2D& r = 0.5 * (W.getAirfoil(q).getR(i + 1) + W.getAirfoil(q).getR(i));
601 double gi = W.getAirfoil(q).viscousStress[i];
602 outfile << static_cast<int>(i) << " " << r[0] << " " << r[1] << " " << gi << std::endl;
603 }//for i
604 outfile.close();
605 }
606 }
607}//SaveVisStress()
bool ifDivisible(int val) const
Definition World2D.h:278
std::string dir
Рабочий каталог задачи
std::string fileNameStep(const std::string &name, int length, size_t number, const std::string &ext)
Формирование имени файла
Definition defs.h:381
void CreateDirectory(const std::string &dir, const std::string &name)
Создание каталога
Definition defs.h:441
int saveVtxStep
Шаг сохранения кадров в бинарные файлы
Definition PassportGen.h:75
int saveVisStress
Шаг вычисления и сохранения скорости и давления
Definition PassportGen.h:83
int nameLength
Число разрядов в имени файла
Definition PassportGen.h:70
Here is the call graph for this function:

Member Data Documentation

◆ optimizedVelocity

OptimizedVelocity VM2D::Velocity::optimizedVelocity
inherited

Definition at line 117 of file Velocity2D.h.

◆ virtualVortexesParams

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

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

Definition at line 115 of file Velocity2D.h.

◆ W

const World2D& VM2D::Velocity::W
protectedinherited

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

Definition at line 108 of file Velocity2D.h.

◆ wakeVortexesParams

VortexesParams VM2D::Velocity::wakeVortexesParams
inherited

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

Definition at line 112 of file Velocity2D.h.


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