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

Класс, опеделяющий текущую решаемую задачу More...

#include <World2D.h>

Inheritance diagram for VM2D::World2D:
Collaboration diagram for VM2D::World2D:

Public Member Functions

Point2D getV0 () const
 Возврат текущей скорости набегающего потока
 
bool isAnyMovable () const
 Возврат признака того, что хотя бы один из профилей подвижный
 
bool isAnyMovableOrDeformable () const
 Возврат признака того, что хотя бы один из профилей подвижный или деформируемый
 
const AirfoilgetAirfoil (size_t i) const
 Возврат константной ссылки на объект профиля
 
const AirfoilGeometrygetOldAirfoil (size_t i) const
 Возврат константной ссылки на объект старого профиля
 
AirfoilgetNonConstAirfoil (size_t i) const
 Возврат неконстантной ссылки на объект профиля
 
size_t getNumberOfAirfoil () const
 Возврат количества профилей в задаче
 
const BoundarygetBoundary (size_t i) const
 Возврат константной ссылки на объект граничного условия
 
BoundarygetNonConstBoundary (size_t i) const
 Возврат неконстантной ссылки на объект граничного условия
 
size_t getNumberOfBoundary () const
 Возврат количества граничных условий в задаче
 
size_t getDispBoundaryInSystem (size_t i) const
 Возврат смещения в системе dispBoundaryInSystem.
 
const MeasureVPgetMeasureVP () const
 Возврат константной ссылки на measureVP.
 
MeasureVPgetNonConstMeasureVP () const
 Возврат неконстантной ссылки на measureVP.
 
const MechanicsgetMechanics (size_t i) const
 Возврат константной ссылки на объект механики
 
MechanicsgetNonConstMechanics (size_t i) const
 Возврат неконстантной ссылки на объект механики
 
const std::vector< std::unique_ptr< Mechanics > > & getMechanicsVector () const
 
const WakegetWake () const
 Возврат константной ссылки на вихревой след
 
WakegetNonConstWake () const
 Возврат неконстантной ссылки на вихревой след
 
const WakeDataBasegetSource () const
 Возврат константной ссылки на источники в области течения
 
const VelocitygetVelocity () const
 Возврат константной ссылки на объект для вычисления скоростей
 
VelocitygetNonConstVelocity () const
 Возврат неконстантной ссылки на объект для вычисления скоростей
 
const PassportgetPassport () const
 Возврат константной ссылки на паспорт
 
PassportgetNonConstPassport () const
 Возврат неконстантной ссылки на паспорт
 
const GpugetCuda () const
 Возврат константной ссылки на объект, связанный с видеокартой (GPU)
 
GpugetNonConstCuda () const
 Возврат неконстантной ссылки на объект, связанный с видеокартой (GPU)
 
const std::pair< Eigen::MatrixXd, Eigen::MatrixXd > & getIQ (size_t i, size_t j) const
 Возврат константной ссылки на объект, связанный с матрицей интегралов от (r-xi)/|r-xi|^2.
 
VMlib::TimersGengetTimers () const
 Возврат ссылки на временную статистику выполнения шага расчета по времени
 
bool ifDivisible (int val) const
 
void SolveLinearSystem ()
 Решение системы линейных алгебраических уравнений
 
void FillIQ ()
 Заполнение матрицы, состоящей из интегралов от (r-xi) / |r-xi|^2.
 
void FillMatrixAndRhs ()
 Заполнение матрицы системы для всех профилей
 
void ReserveMemoryForMatrixAndRhs ()
 Вычисляем размер матрицы и резервируем память под нее и под правую часть
 
void CalcVortexVelo ()
 Вычисление скоростей (и конвективных, и диффузионных) вихрей (в пелене и виртуальных), а также в точках вычисления VP.
 
void CalcPanelsVeloAndAttachedSheets ()
 Вычисление скоростей панелей и интенсивностей присоединенных слоев вихрей и источников
 
void CalcAndSolveLinearSystem ()
 Набор матрицы, правой части и решение СЛАУ
 
void MoveVortexes (std::vector< Point2D > &newPos)
 Вычисляем новые положения вихрей (в пелене и виртуальных)
 
void CheckInside (std::vector< Point2D > &newPos, const std::vector< std::unique_ptr< AirfoilGeometry > > &oldAirfoil)
 Проверка проникновения вихрей внутрь профиля
 
void WakeAndAirfoilsMotion (bool dynamics)
 Перемещение вихрей и профилей на шаге
 
 World2D (const VMlib::PassportGen &passport_)
 Конструктор
 
 ~World2D ()
 Деструктор
 
void GenerateMechanicsHeader (size_t mechanicsNumber)
 
virtual void Step () override
 Функция выполнения предварительного шага
 
VMlib::LogStreamgetInfo () const
 Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода
 
std::ostream & getInfo (char x) const
 Возврат ссылки на поток вывода информации Необходимо для вывода телеметрической информации, информации об ошибках и т.п.
 
size_t getCurrentStep () const
 Возврат константной ссылки на параметры распараллеливания по MPI.
 
double getCurrentTime () const
 
const PassportGen & getPassportGen () const
 
bool isFinished () const
 Функция, возвращающая признак завершения счета в решаемой задаче
 

Public Attributes

int gabb
 
int check01
 
int check02
 
int checkPan
 
VMlib::vmTimer timerInitialBuild
 
VMlib::vmTimer timerRhs
 
VMlib::vmTimer timerFillMatrix
 
VMlib::vmTimer timerSlaeSolve
 
VMlib::vmTimer timerConvVelo
 
VMlib::vmTimer timerInside
 
VMlib::vmTimer timerMerging
 
size_t nVtxBeforeMerging
 

Protected Attributes

LogStream info
 Поток для вывода логов и сообщений об ошибках
 
const PassportGen & passportGen
 Константная ссылка на паспорт конкретного расчета
 
std::unique_ptr< TimersGen > timers
 Сведения о временах выполнения основных операций
 
size_t currentStep
 Текущий номер шага в решаемой задаче
 
double currentTime
 Текущее время в решаемой задаче
 

Private Attributes

std::vector< std::unique_ptr< Airfoil > > airfoil
 Список умных указателей на обтекаемые профили
 
std::vector< std::unique_ptr< AirfoilGeometry > > oldAirfoil
 Список умных указателей на обтекаемые профили для сохранения старого положения
 
std::vector< std::unique_ptr< Boundary > > boundary
 Список умных указателей на формирователи граничных условий на профилях
 
std::vector< size_t > dispBoundaryInSystem
 Список номеров, с которых начинаются элементы правой части (или матрицы) системы для профилей
 
std::vector< std::unique_ptr< Mechanics > > mechanics
 Список умных указателей на типы механической системы для каждого профиля
 
std::unique_ptr< Velocityvelocity
 Умный укзатель на объект, определяющий методику вычисления скоростей
 
std::unique_ptr< Wakewake
 Умный указатель на вихревой след
 
std::unique_ptr< WakeDataBasesource
 Умный указатель на источники
 
std::unique_ptr< MeasureVPmeasureVP
 Умный указатель на алгоритм вычисления полей скоростей и давления (для сохранения в файл)
 
Eigen::MatrixXd matrReord
 Матрица системы
 
std::vector< std::vector< std::pair< Eigen::MatrixXd, Eigen::MatrixXd > > > IQ
 Матрица, состоящая из пар матриц, в которых хранятся касательные и нормальные компоненты интегралов от ядра
 
Eigen::MatrixXd invMatr
 Обратная матрица
 
bool useInverseMatrix
 Признак использования обратной матрицы
 
Eigen::VectorXd rhsReord
 Правая часть системы
 
Eigen::VectorXd sol
 Решение системы
 
const Passportpassport
 Константная ссылка на паспорт конкретного расчета
 
Gpu cuda
 Объект, управляющий графическим ускорителем
 

Detailed Description

Класс, опеделяющий текущую решаемую задачу

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

\Version 1.14

Date
6 марта 2026 г.

Definition at line 73 of file World2D.h.

Constructor & Destructor Documentation

◆ World2D()

World2D::World2D ( const VMlib::PassportGen passport_)

Конструктор

Parameters
[in]passport_константная ссылка на паспорт расчета

Definition at line 75 of file World2D.cpp.

75 :
76 WorldGen(passport_),
77 passport(dynamic_cast<const Passport&>(passport_)),
78 cuda(Gpu(*this))
79{
80 std::stringstream ss;
81 ss << "#" << passport.problemNumber << " (" << passport.problemName << ")";
83
85 currentStep = 0;
86
87 std::vector<std::string> timerLabels = { "Step", "MatRhs", "Solve", "ConvVel", "DiffVel", "Force", "VelPres", "Inside", "Restr", "Save"};
88 timers = std::make_unique<VMlib::TimersGen>(*this, timerLabels);
89
90 wake.reset(new Wake(*this));
91 // загрузка пелены из файла
93 wake->ReadFromFile(passport.wakesDir, passport.wakeDiscretizationProperties.fileWake); //Считываем из каталога с пеленой
94
95 source.reset(new WakeDataBase(*this));
96 // загрузка положений источников из файла
98 source->ReadFromFile(passport.dir, passport.wakeDiscretizationProperties.fileSource); //Считываем из текущего каталога
99
100
101
103 {
104 case 0:
105 velocity.reset(new VelocityBiotSavart(*this));
106 break;
107 case 1:
108 velocity.reset(new VelocityBarnesHut(*this));
109 break;
110 }
111
112 velocity->virtualVortexesParams.resize(passport.airfoilParams.size());
113
114 auto CreateBoundary = [this](size_t i) {
116 {
117 case 0:
118 this->boundary.emplace_back(new BoundaryVortexCollocN(*this, i));
119 //info('e') << "BoundaryMDV is not implemented now! " << std::endl;
120 //exit(1);
121 break;
122
123 case 1:
124 this->boundary.emplace_back(new BoundaryConstLayerAver(*this, i));
125 break;
126
127 case 2:
128 this->boundary.emplace_back(new BoundaryLinLayerAver(*this, i));
129 break;
130
131 default:
132 info('e') << "Unknown scheme!" << std::endl;
133 exit(1);
134 }
135 };
136
137
138 for (size_t i = 0; i < passport.airfoilParams.size(); ++i)
139 {
140
141 switch (passport.airfoilParams[i].mechanicalSystemType)
142 {
143 case 0:
144 airfoil.emplace_back(new AirfoilRigid(*this, i));
145 airfoil[i]->ReadFromFile(passport.airfoilsDir);
146 CreateBoundary(i);
147 mechanics.emplace_back(new MechanicsRigidImmovable(*this, i));
148 break;
149
150 case 1:
151 airfoil.emplace_back(new AirfoilRigid(*this, i));
152 airfoil[i]->ReadFromFile(passport.airfoilsDir);
153 CreateBoundary(i);
154 mechanics.emplace_back(new MechanicsRigidGivenLaw(*this, i));
155 break;
156
157 case 2:
158 airfoil.emplace_back(new AirfoilRigid(*this, i));
159 airfoil[i]->ReadFromFile(passport.airfoilsDir);
160 CreateBoundary(i);
161 mechanics.emplace_back(new MechanicsRigidOscillPart(*this, i));
162 break;
163
164 case 3:
165 airfoil.emplace_back(new AirfoilRigid(*this, i));
166 airfoil[i]->ReadFromFile(passport.airfoilsDir);
167 CreateBoundary(i);
168 mechanics.emplace_back(new MechanicsRigidRotatePart(*this, i));
169 break;
170
171 case 4:
172 airfoil.emplace_back(new AirfoilDeformable(*this, i));
173 airfoil[i]->ReadFromFile(passport.airfoilsDir);
174 CreateBoundary(i);
175 mechanics.emplace_back(new MechanicsDeformable(*this, i));
176 break;
177 }
178 }
179
180 if (getPassport().wakeDiscretizationProperties.eps == 0)
181 {
182 double sumLength = 0.0;
183 size_t totNumPan = 0;
184 for (size_t bou = 0; bou < getNumberOfAirfoil(); ++bou)
185 {
186 for (size_t pnl = 0; pnl < getAirfoil(bou).getNumberOfPanels(); ++pnl)
187 sumLength += getAirfoil(bou).len[pnl];
188 totNumPan += getAirfoil(bou).getNumberOfPanels();
189 }
190 double eps = 0.5 * sumLength / totNumPan;
193 info('i') << "eps = " << getPassport().wakeDiscretizationProperties.eps << " is calculated automatically" << std::endl;
194 }
195
196 if (getPassport().wakeDiscretizationProperties.epscol == 0)
197 {
199
200 info('i') << "epscol = " << getPassport().wakeDiscretizationProperties.epscol << " is calculated automatically" << std::endl;
201 }
202
203 for (size_t afl = 0; afl < passport.airfoilParams.size(); ++afl)
204 if (passport.airfoilParams[afl].chord == 0)
205 {
206 auto& prm = getNonConstPassport().airfoilParams[afl];
207 prm.chord = (prm.initialGab.second[0] - prm.initialGab.first[0]) * prm.scale[0];
208 info('i') << "airfoil #" << afl << " chord = " << prm.chord << " is calculated automatically" << std::endl;
209 }
210
211 //2026-03-28
212 //считываем массив точек для подсчета и вывода поля скоростей и давлений
213 measureVP.reset(new MeasureVP(*this));
214
216 measureVP->ReadPointsFromFile(passport.dir);
217
218 //optimizer
219 //2026-03-28
220 //std::vector<Point2D> VPPoints, VPHistory;
221 //for (size_t a = 0; a < getNumberOfAirfoil(); ++a)
222 //{
223 // const auto& afl = *airfoil[a];
224 // for (size_t p = 0; p < afl.getNumberOfPanels(); ++p)
225 // VPHistory.push_back(0.5 * (afl.getR(p) + afl.getR(p + 1)) + afl.nrm[p] * passport.airfoilParams[a].chord * 0.01);
226 //}
227 //if (passport.timeDiscretizationProperties.saveVPstep != 0)
228 // measureVP->SetPoints(VPPoints, VPHistory);
229
230
231 IQ.resize(passport.airfoilParams.size());
232 for (size_t i = 0; i < passport.airfoilParams.size(); ++i)
233 IQ[i].resize(passport.airfoilParams.size());
234
235
236 info.endl();
237 info.endl();
238 info('i') << "Start solving problem " << passport.problemName << std::endl;
239 info.endl();
240
241 VMlib::PrintLogoToStream(info('_') << std::endl);
242}//World2D(...)
Класс, определяющий тип обтекаемого профиля
std::vector< double > len
Длины панелей профиля
Definition Airfoil2D.h:94
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition Airfoil2D.h:163
Класс, определяющий тип обтекаемого профиля
Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле
Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле
Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле
Класс, обеспечивающий возможность выполнения вычислений на GPU по технологии Nvidia CUDA.
Definition Gpu2D.h:69
Класс, отвечающий за вычисление поля скорости и давления в заданых точках для вывода
Definition MeasureVP2D.h:65
Класс, определяющий вид механической системы
Класс, определяющий вид механической системы
Класс, определяющий вид механической системы
Класс, определяющий вид механической системы
Класс, определяющий вид механической системы
std::string wakesDir
Каталог с файлами вихревых следов
Definition Passport2D.h:270
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition Passport2D.h:292
std::string airfoilsDir
Каталог с файлами профилей
Definition Passport2D.h:267
std::vector< AirfoilParams > airfoilParams
Список структур с параметрами профилей
Definition Passport2D.h:273
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition Passport2D.h:295
Класс, определяющий способ вычисления скоростей
Класс, определяющий способ вычисления скоростей
Класс, опеделяющий набор вихрей
Класс, опеделяющий вихревой след (пелену)
Definition Wake2D.h:63
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition World2D.h:174
Gpu cuda
Объект, управляющий графическим ускорителем
Definition World2D.h:137
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition World2D.h:77
std::unique_ptr< MeasureVP > measureVP
Умный указатель на алгоритм вычисления полей скоростей и давления (для сохранения в файл)
Definition World2D.h:101
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition World2D.h:157
std::unique_ptr< WakeDataBase > source
Умный указатель на источники
Definition World2D.h:98
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition World2D.h:89
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition World2D.h:83
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition World2D.h:251
Passport & getNonConstPassport() const
Возврат неконстантной ссылки на паспорт
Definition World2D.h:256
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition World2D.h:95
std::unique_ptr< Velocity > velocity
Умный укзатель на объект, определяющий методику вычисления скоростей
Definition World2D.h:92
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition World2D.h:134
std::vector< std::vector< std::pair< Eigen::MatrixXd, Eigen::MatrixXd > > > IQ
Матрица, состоящая из пар матриц, в которых хранятся касательные и нормальные компоненты интегралов о...
Definition World2D.h:117
void endl()
Вывод в поток логов пустой строки
Definition LogStream.h:103
void assignStream(std::ostream *pStr_, const std::string &label_)
Связывание потока логов с потоком вывода
Definition LogStream.h:80
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
std::string problemName
Название задачи
std::string dir
Рабочий каталог задачи
size_t problemNumber
Номер задачи
std::unique_ptr< TimersGen > timers
Сведения о временах выполнения основных операций
Definition WorldGen.h:66
double currentTime
Текущее время в решаемой задаче
Definition WorldGen.h:72
size_t currentStep
Текущий номер шага в решаемой задаче
Definition WorldGen.h:69
WorldGen(const VMlib::PassportGen &passport_)
Конструктор
Definition WorldGen.cpp:47
LogStream info
Поток для вывода логов и сообщений об ошибках
Definition WorldGen.h:60
void PrintLogoToStream(std::ostream &str)
Передача в поток вывода шапки программы VM2D/VM3D.
Definition defs.cpp:108
static std::ostream * defaultWorld2DLogStream
Поток вывода логов и ошибок задачи
Definition defs.h:224
std::pair< std::string, int > boundaryCondition
Метод аппроксимации граничных условий
Definition Passport2D.h:189
std::pair< std::string, int > velocityComputation
Definition Passport2D.h:180
std::string fileSource
Имя файла с положениями источников (без полного пути)
Definition Passport2D.h:150
double eps
Радиус вихря
Definition Passport2D.h:126
double epscol
Радиус коллапса
Definition Passport2D.h:132
std::string fileWake
Имя файла с начальным состоянием вихревого следа (без полного пути)
Definition Passport2D.h:147
double eps2
Квадрат радиуса вихря
Definition Passport2D.h:129
double timeStart
Начальное время
Definition PassportGen.h:61
int saveVPstep
Шаг вычисления и сохранения скорости и давления
Definition PassportGen.h:80
Here is the call graph for this function:

◆ ~World2D()

VM2D::World2D::~World2D ( )
inline

Деструктор

Definition at line 340 of file World2D.h.

340{};

Member Function Documentation

◆ CalcAndSolveLinearSystem()

void World2D::CalcAndSolveLinearSystem ( )

Набор матрицы, правой части и решение СЛАУ

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

Definition at line 1506 of file World2D.cpp.

1507{
1508 if (airfoil.size() > 0)
1509 {
1511
1512#if defined(__CUDACC__) || defined(USE_CUDA)
1515
1517
1518 if ((sch == 1) || (sch == 2) || (sch == 0))
1520 else
1521 {
1522 info('e') << "schemeSwitcher is not 0, or 1, or 2! " << std::endl;
1523 exit(1);
1524 }
1525
1526 cuda.RefreshWake(1);
1527 cuda.RefreshAfls(1);
1528 cuda.RefreshVirtualWakes(1);
1529
1530#endif
1531 const int linSystemScheme = passport.numericalSchemes.linearSystemSolver.second;
1532
1533 if (linSystemScheme == 0) //Gauss
1534 {
1535 if (currentStep == 0)
1536 {
1537#ifdef BRIDGE
1538 useInverseMatrix = true;
1539#endif //BRIDGE
1540
1541#ifdef INITIAL
1543 ((mechanics.size() == 1) && (!mechanics.front()->isDeform))
1544 ||
1545 (mechanics.size() > 1 && !isAnyMovableOrDeformable())
1546 );
1547#endif //INITIAL
1548 }
1549 }//if Gauss
1550
1551 if (linSystemScheme == 0 || linSystemScheme == 1 || (linSystemScheme == 2 && isAnyMovableOrDeformable()))
1552 FillIQ();
1553
1555
1556 //{
1557 // std::stringstream ss;
1558 // ss << "IQ1-" << currentStep;
1559 // std::ofstream of(passport.dir + "dbg/" + ss.str());
1560 // for (size_t i = 0; i < (IQ[0][0]).first.rows(); ++i)
1561 // {
1562 // for (size_t j = 0; j < (IQ[0][0]).first.cols(); ++j)
1563 // of << (IQ[0][0]).first(i, j) << " ";
1564 // of << std::endl;
1565 // }
1566 // of.close();
1567 //}
1568
1569 /*
1570 {
1571 std::stringstream ss;
1572 ss << "matr-" << currentStep;
1573 std::ofstream of(passport.dir + "dbg/" + ss.str());
1574 of.precision(16);
1575 for (size_t i = 0; i < matrReord.rows(); ++i)
1576 {
1577 for (size_t j = 0; j < matrReord.cols(); ++j)
1578 of << matrReord(i, j) << " ";
1579 of << std::endl;
1580 }
1581 of.close();
1582 }
1583
1584 {
1585 std::stringstream ss;
1586 ss << "rhs-" << currentStep;
1587 std::ofstream of(passport.dir + "dbg/" + ss.str());
1588 of.precision(16);
1589 for (size_t i = 0; i < rhsReord.rows(); ++i)
1590 {
1591 of << rhsReord(i) << std::endl;
1592 }
1593 of.close();
1594 }
1595
1596 //*/
1597
1599
1600 size_t currentRow = 0;
1601 for (size_t bou = 0; bou < boundary.size(); ++bou)
1602 {
1603 size_t nVars = boundary[bou]->GetUnknownsSize();
1604 Eigen::VectorXd locSol;
1605 locSol.resize(nVars);
1606 for (size_t i = 0; i < nVars; ++i)
1607 locSol(i) = sol(currentRow + i);
1608
1609 boundary[bou]->SolutionToFreeVortexSheetAndVirtualVortex(locSol);
1610 currentRow += nVars;// +1;
1611 }
1612
1613 //14-05-2024
1614 size_t nVirtVortices = 0;
1615 for (size_t bou = 0; bou < boundary.size(); ++bou)
1616 nVirtVortices += boundary[bou]->virtualWake.vtx.size();
1617
1618 wake->vtx.reserve(wake->vtx.size() + nVirtVortices);
1619 for (size_t bou = 0; bou < boundary.size(); ++bou)
1620 for (size_t v = 0; v < boundary[bou]->virtualWake.vtx.size(); ++v)
1621 wake->vtx.push_back(Vortex2D{ boundary[bou]->virtualWake.vtx[v].r(), 0.0 });
1622
1623 //for (size_t v = 0; v < wake->vtx.size(); ++v)
1624 //{
1625 // if (fabs(wake->vtx[v].g()) > 1.0)
1626 // std::cout << "Error gam! " << "v = " << v << ", gam[v] = " << wake->vtx[v].g() << std::endl;
1627 //}
1628
1629 //for (size_t bou = 0; bou < boundary.size(); ++bou)
1630 // for (size_t v = 0; v < airfoil[bou]->getNumberOfPanels(); ++v)
1631 // {
1632 // if (fabs(boundary[bou]->sheets.freeVortexSheet(v, 0)) > 1000.0)
1633 // std::cout << "Error gamma sheet! " << "bou = " << bou << ", v = " << v << ", gam[v] = " << boundary[bou]->sheets.freeVortexSheet(v, 0) << std::endl;
1634 // }
1635
1636 /*
1637 if (currentStep == 0)
1638 {
1639 Point2D addMass = { 0.0, 0.0 };
1640 for (size_t q = 0; q < airfoil[0]->getNumberOfPanels(); ++q)
1641 {
1642 addMass += (sol(q) + boundary[0]->sheets.attachedVortexSheet(q,0)) * 0.5 * (airfoil[0]->getR(q) + airfoil[0]->getR(q + 1)).kcross() * airfoil[0]->len[q];
1643 }
1644 addMass *= passport.physicalProperties.rho;
1645
1646 std::cout << "AddMass = " << addMass << std::endl;
1647 //exit(-42);
1648 }
1649 */
1650 }
1651}//CalcAndSolveLinearSystem()
void setAccelCoeff(double cft_)
Установка коэффициента разгона потока
Definition Gpu2D.h:272
void setMaxGamma(double gam_)
Установка максимально допустимой циркуляции вихря
Definition Gpu2D.h:295
void setSchemeSwitcher(int schemeSwitcher_)
Установка переключателя расчетных схем
Definition Gpu2D.h:309
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition Passport2D.h:289
void ReserveMemoryForMatrixAndRhs()
Вычисляем размер матрицы и резервируем память под нее и под правую часть
Definition World2D.cpp:1213
bool isAnyMovableOrDeformable() const
Возврат признака того, что хотя бы один из профилей подвижный или деформируемый
Definition World2D.cpp:1765
bool useInverseMatrix
Признак использования обратной матрицы
Definition World2D.h:123
void FillIQ()
Заполнение матрицы, состоящей из интегралов от (r-xi) / |r-xi|^2.
Definition World2D.cpp:1464
void SolveLinearSystem()
Решение системы линейных алгебраических уравнений
Definition World2D.cpp:826
void FillMatrixAndRhs()
Заполнение матрицы системы для всех профилей
Definition World2D.cpp:1115
Eigen::VectorXd sol
Решение системы
Definition World2D.h:131
Класс, опеделяющий двумерный вихревой элемент
Definition Vortex2D.h:59
double getCurrentTime() const
Definition WorldGen.h:100
std::pair< std::string, int > linearSystemSolver
Definition Passport2D.h:174
double accelCft(double currentTime) const
Функция-множитель, позволяющая моделировать разгон
double maxGamma
Максимально допустимая циркуляция вихря
Definition Passport2D.h:144
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalcPanelsVeloAndAttachedSheets()

void World2D::CalcPanelsVeloAndAttachedSheets ( )

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

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

Definition at line 1385 of file World2D.cpp.

1386{
1387 //вычисляем скорости панелей
1388 for (size_t i = 0; i < airfoil.size(); ++i)
1390
1391 //вычисляем интенсивности присоединенных слоев
1392 for (size_t i = 0; i < airfoil.size(); ++i)
1393 boundary[i]->ComputeAttachedSheetsIntensity();
1394
1395}//CalcPanelsVeloAndAttachedSheets(...)
double dt
Шаг по времени
Definition PassportGen.h:67
Here is the caller graph for this function:

◆ CalcVortexVelo()

void World2D::CalcVortexVelo ( )

Вычисление скоростей (и конвективных, и диффузионных) вихрей (в пелене и виртуальных), а также в точках вычисления VP.

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

Definition at line 1273 of file World2D.cpp.

1274{
1275 velocity->ResizeAndZero();
1276
1277 getTimers().start("ConvVel");
1278
1279 //Конвективные скорости всех вихрей (и в пелене, и виртуальных), индуцируемые вихрями в пелене
1280 //Подготовка CUDA
1281#if defined(__CUDACC__) || defined(USE_CUDA)
1282 cuda.RefreshWake(2);
1283
1284 cuda.RefreshAfls(2);
1285 cuda.RefreshVirtualWakes(2);
1286#endif
1287 getTimers().stop("ConvVel");
1288
1289 velocity->CalcConvVelo();
1290
1291 //Расчет средних значений eps для каждой панели и их передача на видеокарту
1292
1293 getTimers().start("DiffVel");
1294
1295 for (size_t bou = 0; bou < getNumberOfBoundary(); ++bou)
1297
1298#if defined(__CUDACC__) || defined(USE_CUDA)
1299 for (size_t i = 0; i < airfoil.size(); ++i)
1300 cuda.CopyMemToDev<double, 1>(airfoil[i]->getNumberOfPanels(), airfoil[i]->meanEpsOverPanel.data(), airfoil[i]->devMeanEpsOverPanelPtr);
1301#endif
1302
1303 getTimers().stop("DiffVel");
1304
1305 //Вычисление диффузионных скоростей вихрей (и в пелене, и виртуальных)
1306 velocity->CalcDiffVelo();
1307
1308 //Обнуление вязких напряжений, если они не были вычислены
1309 for (auto& afl : airfoil)
1310 {
1311 if (afl->viscousStress.size() == 0)
1312 afl->viscousStress.resize(afl->getNumberOfPanels(), 0.0);
1313 }
1314
1315 getTimers().start("Save");
1316
1317 /*
1318 //Сохранение всех параметров для вихрей в пелене
1319 if(!(currentStep % 1))
1320 {
1321 VMlib::CreateDirectory(passport.dir, "dbg");
1322 std::ostringstream sss;
1323 sss << "prmWake";
1324 sss << currentStep;
1325 std::ofstream prmtFile(passport.dir + "dbg/" + sss.str());
1326 prmtFile.precision(11);
1327 prmtFile << "i x y g epsast convVeloX convVeloY diffVeloX diffVeloY I0 I1 I2X I2Y I3X I3Y" << std::endl;
1328 for (size_t i = 0; i < wake->vtx.size(); ++i)
1329 prmtFile << i << " " \
1330 << wake->vtx[i].r()[0] << " " << wake->vtx[i].r()[1] << " " \
1331 << wake->vtx[i].g() << " " \
1332 << velocity->wakeVortexesParams.epsastWake[i] << " " \
1333 << velocity->wakeVortexesParams.convVelo[i][0] << " " << velocity->wakeVortexesParams.convVelo[i][1] << " "\
1334 << velocity->wakeVortexesParams.diffVelo[i][0] << " " << velocity->wakeVortexesParams.diffVelo[i][1] << " "\
1335 //<< std::endl;
1336 << velocity->wakeVortexesParams.I0[i] << " " \
1337 << velocity->wakeVortexesParams.I1[i] << " " \
1338 << velocity->wakeVortexesParams.I2[i][0] << " " << velocity->wakeVortexesParams.I2[i][1] << " " \
1339 << velocity->wakeVortexesParams.I3[i][0] << " " << velocity->wakeVortexesParams.I3[i][1] << " " \
1340 << "\n";
1341
1342 prmtFile.close();
1343 }
1344//*/
1345
1346/*
1347 //Сохранение всех параметров для виртуальных вихрей
1348 if (!(currentStep % 1))
1349 {
1350 for (size_t b = 0; b < boundary.size(); ++b)
1351 {
1352 std::ostringstream sss;
1353 sss << "prmVirtual_";
1354 sss << b << "-";
1355 sss << currentStep;
1356 std::ofstream prmFileVirt(passport.dir + "dbg/" + sss.str());
1357 //prmFileVirt.precision(16);
1358 prmFileVirt << "i x y g epsast convVeloX convVeloY diffVeloX diffVeloY I0 I1 I2X I2Y I3X I3Y" << std::endl;
1359 for (size_t i = 0; i < boundary[b]->virtualWake.vtx.size(); ++i)
1360 prmFileVirt << i << " " \
1361 << boundary[b]->virtualWake.vtx[i].r()[0] << " " << boundary[b]->virtualWake.vtx[i].r()[1] << " " \
1362 << boundary[b]->virtualWake.vtx[i].g() << " " \
1363 << velocity->virtualVortexesParams[b].epsastWake[i] << " " \
1364 << velocity->virtualVortexesParams[b].convVelo[i][0] << " " << velocity->virtualVortexesParams[b].convVelo[i][1] << " "\
1365 << velocity->virtualVortexesParams[b].diffVelo[i][0] << " " << velocity->virtualVortexesParams[b].diffVelo[i][1] << " "\
1366 //<< std::endl;
1367 << velocity->virtualVortexesParams[b].I0[i] << " " \
1368 << velocity->virtualVortexesParams[b].I1[i] << " " \
1369 << velocity->virtualVortexesParams[b].I2[i][0] << " " << velocity->virtualVortexesParams[b].I2[i][1] << " " \
1370 << velocity->virtualVortexesParams[b].I3[i][0] << " " << velocity->virtualVortexesParams[b].I3[i][1] << " " \
1371 << std::endl;
1372 prmFileVirt.close();
1373 }
1374 //if (currentStep==3) exit(-123);
1375 }
1376//*/
1377
1378 getTimers().stop("Save");
1379
1380}//CalcVortexVelo()
void calcMeanEpsOverPanel()
Вычисление средних значений eps на панелях
Definition Airfoil2D.cpp:93
VMlib::TimersGen & getTimers() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition World2D.h:276
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition World2D.h:191
Airfoil & getNonConstAirfoil(size_t i) const
Возврат неконстантной ссылки на объект профиля
Definition World2D.h:169
void stop(const std::string &timerLabel)
Останов счетчика
Definition TimesGen.cpp:68
void start(const std::string &timerLabel)
Запуск счетчика
Definition TimesGen.cpp:55
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CheckInside()

void World2D::CheckInside ( std::vector< Point2D > &  newPos,
const std::vector< std::unique_ptr< AirfoilGeometry > > &  oldAirfoil 
)

Проверка проникновения вихрей внутрь профиля

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

Parameters
[in]newPosновые позиции вихрей
[in]oldAirfoilконстантная ссылка на вектор из умных указателей на старые положения профилей

Definition at line 713 of file World2D.cpp.

714{
715 getTimers().start("Inside");
716
717 gabb = 0;
718 check01 = 0;
719 check02 = 0;
720 checkPan = 0;
721
722#if (!defined(USE_CUDA))
723 nVtxBeforeMerging = newPos.size();
724 for (size_t afl = 0; afl < airfoil.size(); ++afl)
725 wake->Inside(newPos, *airfoil[afl], mechanics[afl]->isMoves, *oldAirfoil[afl]);
726
727 //std::cout << "gab: " << gabb << " check1: " << check01 << " check2: " << check02 << " checkPan: " << checkPan << std::endl;
728#else
729// //////////////////////// CUDA ///////////////////////
732
733 if ((newPos.size() > 0) && (getNumberOfAirfoil() > 0))
734 {
735 int nTotPanels = 0;
736 for (size_t afl = 0; afl < airfoil.size(); ++afl)
737 nTotPanels += (int)airfoil[afl]->getNumberOfPanels();
738
739 nVtxBeforeMerging = newPos.size();
740
741 auto& auxTree = *getNonConstCuda().auxTreePnl;
742
744 {
745 auxTree.MemoryAllocate((int)getCuda().n_CUDA_pnls);
746 auxTree.UpdatePanelGeometry(nTotPanels, (double4*)airfoil[0]->devRPtr);
747 auxTree.Build();
748 auxTree.UpwardTraversal(0);
749 }
750
751 std::vector<int> hit(newPos.size());
752
753 if (isAnyMovableOrDeformable()) //Если профили подвижны - метод псевдонормалей
754 {
755 double* devNewpos_ptr;
756 cudaMalloc(&devNewpos_ptr, newPos.size() * sizeof(double) * 2);
757 cudaMemcpy(devNewpos_ptr, newPos.data(), newPos.size() * sizeof(double) * 2, cudaMemcpyHostToDevice);
758
759 auto& cntrTreePnt = *getCuda().cntrTreePoint;
760 cntrTreePnt.MemoryAllocate((int)getCuda().n_CUDA_wake);
761 cntrTreePnt.Update((int)newPos.size(), devNewpos_ptr, 0.0);
762 cntrTreePnt.Build();
763
764 BHcu::treeClosestPanelToPointsCalculationWrapper(auxTree, cntrTreePnt, getWake().devNearestPanelPtr, true, getAirfoil(0).devPsnPtr);
765
766 cudaMemcpy(hit.data(), getWake().devNearestPanelPtr, newPos.size() * sizeof(int), cudaMemcpyDeviceToHost);
767 cudaFree(devNewpos_ptr);
768 }
769 else //иначе (для неподвижного профиля - метод трассировки лучей
770 {
771 std::vector<std::pair<Point2D, Point2D>> segments(newPos.size());
772 for (size_t i = 0; i < newPos.size(); ++i)
773 {
774 segments[i].first = wake->vtx[i].r();
775 segments[i].second = newPos[i];
776 }
777 double* devSegments_ptr;
778
779 cudaMalloc(&devSegments_ptr, newPos.size() * sizeof(double) * 4);
780 cudaMemcpy(devSegments_ptr, segments.data(), newPos.size() * sizeof(double) * 4, cudaMemcpyHostToDevice);
781
782 auto& cntrTreeSeg = *getCuda().cntrTreeSegment;
783 cntrTreeSeg.MemoryAllocate((int)getCuda().n_CUDA_wake);
784 cntrTreeSeg.UpdatePanelGeometry((int)newPos.size(), (double4*)devSegments_ptr);
785 cntrTreeSeg.Build();
786
787 BHcu::treePanelsSegmentsIntersectionCalculationWrapper(auxTree, cntrTreeSeg, getWake().devNearestPanelPtr);
788 cudaMemcpy(hit.data(), getWake().devNearestPanelPtr, newPos.size() * sizeof(int), cudaMemcpyDeviceToHost);
789 cudaFree(devSegments_ptr);
790 }
791
792 size_t panCounter = 0;
793 for (size_t afl = 0; afl < airfoil.size(); ++afl)
794 {
795 std::vector<double> gamma(airfoil[afl]->getNumberOfPanels(), 0.0);
796 for (int i = 0; i < newPos.size(); ++i)
797 {
798 if (hit[i] != -1)
799 {
800 if ((hit[i] >= panCounter) && (hit[i] < panCounter + airfoil[afl]->getNumberOfPanels()))
801 {
802 if (fabs(wake->vtx[i].g()) > 1.0)
803 std::cout << "Too large gamma is through: i = " << i << ", " << hit[i] << ", " << "gamma[hit[i] - panCounter] += " << wake->vtx[i].g() << std::endl;
804
805 gamma[hit[i] - panCounter] += wake->vtx[i].g();
806 wake->vtx[i].g() = 0.0;
807 }
808 }
809 }
810 airfoil[afl]->gammaThrough = gamma;
811 panCounter += airfoil[afl]->getNumberOfPanels();
812 }//for afl
813 }
815
816#endif
817
818 getTimers().stop("Inside");
819}
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition World2D.h:226
std::vector< std::unique_ptr< AirfoilGeometry > > oldAirfoil
Список умных указателей на обтекаемые профили для сохранения старого положения
Definition World2D.h:80
VMlib::vmTimer timerInside
Definition World2D.h:358
Gpu & getNonConstCuda() const
Возврат неконстантной ссылки на объект, связанный с видеокартой (GPU)
Definition World2D.h:266
const Gpu & getCuda() const
Возврат константной ссылки на объект, связанный с видеокартой (GPU)
Definition World2D.h:261
size_t nVtxBeforeMerging
Definition WorldGen.h:76
const vmTimer & stop() const
Останов работающего счетчика времени
Definition TimesGen.h:125
const vmTimer & start() const
Запуск (первый или повторный) счетчика времени
Definition TimesGen.h:109
const vmTimer & reset() const
Сброс счетчика времени
Definition TimesGen.h:101
Here is the call graph for this function:
Here is the caller graph for this function:

◆ FillIQ()

void World2D::FillIQ ( )

Заполнение матрицы, состоящей из интегралов от (r-xi) / |r-xi|^2.

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

Definition at line 1464 of file World2D.cpp.

1465{
1466 getTimers().start("MatRhs");
1467
1468 for (size_t bou = 0; bou < boundary.size(); ++bou)
1469 {
1470 if (currentStep == 0 || mechanics[bou]->isDeform)
1471 {
1472 boundary[bou]->FillIQSelf(IQ[bou][bou]);
1473 }
1474
1475 for (size_t oth = 0; oth < boundary.size(); ++oth)
1476 {
1477
1478#ifdef INITIAL
1479 if (currentStep == 0 || !useInverseMatrix)
1480 {
1481 //size_t nVarsOther = boundary[oth]->GetUnknownsSize();
1482 if (bou != oth)
1483 {
1484 //std::cout << "matr!" << std::endl;
1485 boundary[bou]->FillIQFromOther(*boundary[oth], IQ[bou][oth]);
1486 }
1487 }// if (currentStep == 0 || !useInverseMatrix)
1488#endif
1489
1490#ifdef BRIDGE
1491 if (currentStep == 0)
1492 {
1493 //size_t nVarsOther = boundary[oth]->GetUnknownsSize();
1494 if (bou != oth)
1495 boundary[bou]->FillIQFromOther(*boundary[oth], IQ[bou][oth]);
1496 }// if (currentStep == 0)
1497#endif
1498
1499 }// for oth
1500
1501 }// for bou
1502
1503 getTimers().stop("MatRhs");
1504}//FillIQ()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ FillMatrixAndRhs()

void World2D::FillMatrixAndRhs ( )

Заполнение матрицы системы для всех профилей

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

Definition at line 1115 of file World2D.cpp.

1116{
1117 getTimers().start("MatRhs");
1118
1119 const int linSystemScheme = passport.numericalSchemes.linearSystemSolver.second;
1120
1121 if (linSystemScheme == 0 || linSystemScheme == 1)
1122 {
1125
1126 Eigen::MatrixXd locMatr;
1127 Eigen::MatrixXd otherMatr;
1128 Eigen::VectorXd locLastLine, locLastCol;
1129
1130 std::vector<std::vector<Point2D>> locIQ;
1131 std::vector<std::vector<Point2D>> locOtherIQ;
1132
1133 //обнуляем матрицу на первом шаге расчета
1134 if (currentStep == 0)
1135 {
1136 for (int i = 0; i < matrReord.rows(); ++i)
1137 for (int j = 0; j < matrReord.cols(); ++j)
1138 matrReord(i, j) = 0.0;
1139 }
1140
1141 size_t currentRow = 0;
1142 size_t currentSkosRow = 0;
1143
1144 size_t nAllVars = 0;
1145 for (size_t bou = 0; bou < boundary.size(); ++bou)
1146 nAllVars += boundary[bou]->GetUnknownsSize();
1147
1148 for (size_t bou = 0; bou < boundary.size(); ++bou)
1149 {
1150 size_t nVars = boundary[bou]->GetUnknownsSize();
1151 if (currentStep == 0 || mechanics[bou]->isDeform)
1152 {
1153 locMatr.resize(nVars, nVars);
1154 locLastLine.resize(nVars);
1155 locLastCol.resize(nVars);
1156 }
1157
1158 if (currentStep == 0 || mechanics[bou]->isDeform)
1159 boundary[bou]->FillMatrixSelf(locMatr, locLastLine, locLastCol);
1160
1161
1162 //размазываем матрицу
1163 for (size_t i = 0; i < nVars; ++i)
1164 {
1165 if (currentStep == 0 || mechanics[bou]->isDeform)
1166 {
1167 for (size_t j = 0; j < nVars; ++j)
1168 matrReord(i + currentSkosRow, j + currentSkosRow) = locMatr(i, j);
1169
1170 matrReord(nAllVars + bou, i + currentSkosRow) = locLastLine(i);
1171 matrReord(i + currentSkosRow, nAllVars + bou) = locLastCol(i);
1172 }
1173 }
1174
1175 if ((currentStep == 0) || (!useInverseMatrix))
1176 {
1177 size_t currentCol = 0;
1178 size_t currentSkosCol = 0;
1179 for (size_t oth = 0; oth < boundary.size(); ++oth)
1180 {
1181 size_t nVarsOther = boundary[oth]->GetUnknownsSize();
1182
1183 if (bou != oth)
1184 {
1185 otherMatr.resize(nVars, nVarsOther);
1186
1187 boundary[bou]->FillMatrixFromOther(*boundary[oth], otherMatr);
1188
1189 //размазываем матрицу
1190 for (size_t i = 0; i < nVars; ++i)
1191 {
1192 for (size_t j = 0; j < nVarsOther; ++j)
1193 matrReord(i + currentSkosRow, j + currentSkosCol) = otherMatr(i, j);
1194 }
1195 }// if (bou != oth)
1196 currentCol += nVarsOther + 1;
1197 currentSkosCol += nVarsOther;
1198 }// for oth
1199 }// if (currentStep == 0 || mechanics[oth]->isMoves)
1200
1201 currentRow += nVars + 1;
1202 currentSkosRow += nVars;
1203 }// for bou
1205 }
1206
1207 velocity->FillRhs(rhsReord);
1208
1209 getTimers().stop("MatRhs");
1210}//FillMatrixAndRhs()
Eigen::MatrixXd matrReord
Матрица системы
Definition World2D.h:113
Eigen::VectorXd rhsReord
Правая часть системы
Definition World2D.h:127
VMlib::vmTimer timerFillMatrix
Definition World2D.h:354
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GenerateMechanicsHeader()

void World2D::GenerateMechanicsHeader ( size_t  mechanicsNumber)

Метод-обертка для вызова метода генерации заголовка файла нагрузок и заголовка файла положения (последнее — если профиль движется)

Parameters
[in]mechanicsNumberномер профиля, для которого генерируется заголовок файла

Definition at line 1457 of file World2D.cpp.

1458{
1459 mechanics[mechanicsNumber]->GenerateForcesHeader();
1460 mechanics[mechanicsNumber]->GeneratePositionHeader();
1461}//GenerateMechanicsHeader(...)
Here is the caller graph for this function:

◆ getAirfoil()

const Airfoil & VM2D::World2D::getAirfoil ( size_t  i) const
inline

Возврат константной ссылки на объект профиля

Parameters
[in]iномер профиля, константная ссылка на который возвращается
Returns
константную ссылку на i-й профиль

Definition at line 157 of file World2D.h.

157{ return *airfoil[i]; };
Here is the caller graph for this function:

◆ getBoundary()

const Boundary & VM2D::World2D::getBoundary ( size_t  i) const
inline

Возврат константной ссылки на объект граничного условия

Parameters
[in]iномер граничного условия, константная ссылка на которое возвращается
Returns
константную ссылку на i-е граничное условие

Definition at line 180 of file World2D.h.

180{ return *boundary[i]; };
Here is the caller graph for this function:

◆ getCuda()

const Gpu & VM2D::World2D::getCuda ( ) const
inline

Возврат константной ссылки на объект, связанный с видеокартой (GPU)

Returns
константную ссылку на объект, связанный с видеокартой (GPU)

Definition at line 261 of file World2D.h.

261{ return cuda; };
Here is the caller graph for this function:

◆ getCurrentStep()

size_t VMlib::WorldGen::getCurrentStep ( ) const
inlineinherited

Возврат константной ссылки на параметры распараллеливания по MPI.

Returns
константную ссылку на параметры распараллеливания по MPI const Parallel& getParallel() const { return parallel; };

Возврат номера текущего временного шага

Returns
номера текущего временного шага

Definition at line 99 of file WorldGen.h.

99{ return currentStep; };
Here is the caller graph for this function:

◆ getCurrentTime()

double VMlib::WorldGen::getCurrentTime ( ) const
inlineinherited

Definition at line 100 of file WorldGen.h.

100{ return currentTime; };
Here is the caller graph for this function:

◆ getDispBoundaryInSystem()

size_t VM2D::World2D::getDispBoundaryInSystem ( size_t  i) const
inline

Возврат смещения в системе dispBoundaryInSystem.

Parameters
[in]iномер граничного условия, константная ссылка на которое возвращается
Returns
константную ссылку на i-е граничное условие

Definition at line 197 of file World2D.h.

197{ return dispBoundaryInSystem[i]; };
std::vector< size_t > dispBoundaryInSystem
Список номеров, с которых начинаются элементы правой части (или матрицы) системы для профилей
Definition World2D.h:86

◆ getInfo() [1/2]

VMlib::LogStream & VMlib::WorldGen::getInfo ( ) const
inlineinherited

Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода

Returns
ссылку на объект LogStream

Definition at line 82 of file WorldGen.h.

82{ return info; };
Here is the caller graph for this function:

◆ getInfo() [2/2]

std::ostream & VMlib::WorldGen::getInfo ( char  x) const
inlineinherited

Возврат ссылки на поток вывода информации Необходимо для вывода телеметрической информации, информации об ошибках и т.п.

Parameters
[in]xсимвол, определяющий стиль вывода сообщения
Returns
ссылку на поток вывода информации

Definition at line 89 of file WorldGen.h.

89{ return info(x); };

◆ getIQ()

const std::pair< Eigen::MatrixXd, Eigen::MatrixXd > & VM2D::World2D::getIQ ( size_t  i,
size_t  j 
) const
inline

Возврат константной ссылки на объект, связанный с матрицей интегралов от (r-xi)/|r-xi|^2.

Returns
константную ссылку на объект, связанный с матрицей интегралов от (r-xi)/|r-xi|^2

Definition at line 271 of file World2D.h.

271{ return IQ[i][j]; };
Here is the caller graph for this function:

◆ getMeasureVP()

const MeasureVP & VM2D::World2D::getMeasureVP ( ) const
inline

Возврат константной ссылки на measureVP.

Returns
константную ссылку на вихревой след

Definition at line 202 of file World2D.h.

202{ return *measureVP; };
Here is the caller graph for this function:

◆ getMechanics()

const Mechanics & VM2D::World2D::getMechanics ( size_t  i) const
inline

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

Parameters
[in]iномер механики, константная ссылка на который возвращается
Returns
константную ссылку на i-ю механику

Definition at line 213 of file World2D.h.

213{ return *mechanics[i]; };
Here is the caller graph for this function:

◆ getMechanicsVector()

const std::vector< std::unique_ptr< Mechanics > > & VM2D::World2D::getMechanicsVector ( ) const
inline

Definition at line 221 of file World2D.h.

221{ return mechanics; };
Here is the caller graph for this function:

◆ getNonConstAirfoil()

Airfoil & VM2D::World2D::getNonConstAirfoil ( size_t  i) const
inline

Возврат неконстантной ссылки на объект профиля

Parameters
[in]iномер профиля, неконстантная ссылка на который возвращается
Returns
неконстантную ссылку на i-й профиль

Definition at line 169 of file World2D.h.

169{ return *airfoil[i]; };
Here is the caller graph for this function:

◆ getNonConstBoundary()

Boundary & VM2D::World2D::getNonConstBoundary ( size_t  i) const
inline

Возврат неконстантной ссылки на объект граничного условия

Parameters
[in]iномер граничного условия, неконстантная ссылка на которое возвращается
Returns
неконстантную ссылку на i-е граничное условие

Definition at line 186 of file World2D.h.

186{ return *boundary[i]; };
Here is the caller graph for this function:

◆ getNonConstCuda()

Gpu & VM2D::World2D::getNonConstCuda ( ) const
inline

Возврат неконстантной ссылки на объект, связанный с видеокартой (GPU)

Returns
неконстантную ссылку на объект, связанный с видеокартой (GPU)

Definition at line 266 of file World2D.h.

266{ return cuda; };
Here is the caller graph for this function:

◆ getNonConstMeasureVP()

MeasureVP & VM2D::World2D::getNonConstMeasureVP ( ) const
inline

Возврат неконстантной ссылки на measureVP.

Returns
неконстантную ссылку на вихревой след

Definition at line 207 of file World2D.h.

207{ return *measureVP; };
Here is the caller graph for this function:

◆ getNonConstMechanics()

Mechanics & VM2D::World2D::getNonConstMechanics ( size_t  i) const
inline

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

Parameters
[in]iномер механики, константная ссылка на который возвращается
Returns
неконстантную ссылку на i-ю механику

Definition at line 219 of file World2D.h.

219{ return *mechanics[i]; };
Here is the caller graph for this function:

◆ getNonConstPassport()

Passport & VM2D::World2D::getNonConstPassport ( ) const
inline

Возврат неконстантной ссылки на паспорт

Returns
неконстантную ссылку на паспорт

Definition at line 256 of file World2D.h.

256{ return const_cast<Passport&>(passport); };
Here is the caller graph for this function:

◆ getNonConstVelocity()

Velocity & VM2D::World2D::getNonConstVelocity ( ) const
inline

Возврат неконстантной ссылки на объект для вычисления скоростей

Returns
неконстантную ссылку на объект для вычисления скоростей

Definition at line 246 of file World2D.h.

246{ return *velocity; };

◆ getNonConstWake()

Wake & VM2D::World2D::getNonConstWake ( ) const
inline

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

Returns
неконстантную ссылку на вихревой след

Definition at line 231 of file World2D.h.

231{ return *wake; };
Here is the caller graph for this function:

◆ getNumberOfAirfoil()

size_t VM2D::World2D::getNumberOfAirfoil ( ) const
inline

Возврат количества профилей в задаче

Returns
количество профилей в задаче

Definition at line 174 of file World2D.h.

174{ return airfoil.size(); };
Here is the caller graph for this function:

◆ getNumberOfBoundary()

size_t VM2D::World2D::getNumberOfBoundary ( ) const
inline

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

Returns
количество граничных условий в задаче

Definition at line 191 of file World2D.h.

191{ return boundary.size(); };
Here is the caller graph for this function:

◆ getOldAirfoil()

const AirfoilGeometry & VM2D::World2D::getOldAirfoil ( size_t  i) const
inline

Возврат константной ссылки на объект старого профиля

Parameters
[in]iномер старого профиля, константная ссылка на который возвращается
Returns
константную ссылку на i-й старый профиль

Definition at line 163 of file World2D.h.

163{ return *oldAirfoil[i]; };
Here is the caller graph for this function:

◆ getPassport()

const Passport & VM2D::World2D::getPassport ( ) const
inline

Возврат константной ссылки на паспорт

Returns
константную ссылку на паспорт

Definition at line 251 of file World2D.h.

251{ return passport; };
Here is the caller graph for this function:

◆ getPassportGen()

const PassportGen & VMlib::WorldGen::getPassportGen ( ) const
inlineinherited

Definition at line 102 of file WorldGen.h.

102{ return passportGen; };
const PassportGen & passportGen
Константная ссылка на паспорт конкретного расчета
Definition WorldGen.h:63
Here is the caller graph for this function:

◆ getSource()

const WakeDataBase & VM2D::World2D::getSource ( ) const
inline

Возврат константной ссылки на источники в области течения

Returns
константную ссылку на источники в области течения

Definition at line 236 of file World2D.h.

236{ return *source; };
Here is the caller graph for this function:

◆ getTimers()

VMlib::TimersGen & VM2D::World2D::getTimers ( ) const
inline

Возврат ссылки на временную статистику выполнения шага расчета по времени

Returns
ссылку на временную статистику выполнения шага расчета по времени

Definition at line 276 of file World2D.h.

276{ return *timers; };
Here is the caller graph for this function:

◆ getV0()

Point2D VM2D::World2D::getV0 ( ) const
inline

Возврат текущей скорости набегающего потока

Definition at line 142 of file World2D.h.

143 {
145 }
Point2D V0(double currentTime) const
Функция скорости набегающего потока с учетом разгона
Definition Passport2D.h:93
Here is the call graph for this function:
Here is the caller graph for this function:

◆ getVelocity()

const Velocity & VM2D::World2D::getVelocity ( ) const
inline

Возврат константной ссылки на объект для вычисления скоростей

Returns
константную ссылку на объект для вычисления скоростей

Definition at line 241 of file World2D.h.

241{ return *velocity; };
Here is the caller graph for this function:

◆ getWake()

const Wake & VM2D::World2D::getWake ( ) const
inline

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

Returns
константную ссылку на вихревой след

Definition at line 226 of file World2D.h.

226{ return *wake; };
Here is the caller graph for this function:

◆ ifDivisible()

bool VM2D::World2D::ifDivisible ( int  val) const
inline

Definition at line 278 of file World2D.h.

278{ return ((val > 0) && (!(currentStep % val))); };
Here is the caller graph for this function:

◆ isAnyMovable()

bool World2D::isAnyMovable ( ) const

Возврат признака того, что хотя бы один из профилей подвижный

Definition at line 1758 of file World2D.cpp.

1759{
1760 return std::any_of(getMechanicsVector().begin(), getMechanicsVector().end(),
1761 [](const std::unique_ptr<Mechanics>& m) {return (m->isMoves); });
1762}
const std::vector< std::unique_ptr< Mechanics > > & getMechanicsVector() const
Definition World2D.h:221
Here is the call graph for this function:

◆ isAnyMovableOrDeformable()

bool World2D::isAnyMovableOrDeformable ( ) const

Возврат признака того, что хотя бы один из профилей подвижный или деформируемый

Definition at line 1765 of file World2D.cpp.

1766{
1767 return std::any_of(getMechanicsVector().begin(), getMechanicsVector().end(),
1768 [](const std::unique_ptr<Mechanics>& m) {return (m->isMoves || m->isDeform); });
1769}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ isFinished()

bool WorldGen::isFinished ( ) const
inherited

Функция, возвращающая признак завершения счета в решаемой задаче

true, если задача решена и выполнен признак останова; false если требуется еще выполнять шаги по времени

Definition at line 52 of file WorldGen.cpp.

53{
55}
double timeStop
Конечное время
Definition PassportGen.h:64

◆ MoveVortexes()

void World2D::MoveVortexes ( std::vector< Point2D > &  newPos)

Вычисляем новые положения вихрей (в пелене и виртуальных)

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

Parameters
[out]newPosновые позиции вихрей

Definition at line 1401 of file World2D.cpp.

1402{
1403 //getTimers().start("Move");
1404
1405 size_t nvt = wake->vtx.size();
1406 size_t nVirtVortex = 0;
1407 for (size_t i = 0; i < getNumberOfBoundary(); ++i)
1408 nVirtVortex += boundary[i]->virtualWake.vtx.size();
1409
1410 //newPos.clear();
1411 newPos.resize(nvt);
1412
1413
1414 // 1. Инициализация генератора случайных чисел (ПСЧГ)
1415 //std::random_device rd; // Источник энтропии для инициализации
1416 //std::mt19937 gen(rd()); // Mersenne Twister, инициализированный случайным значением
1417
1418 // 2. Определение параметров нормального распределения
1419 // μ (среднее) = 0.0, σ (стандартное отклонение) = 1.0
1420 //std::normal_distribution<> d(0.0, sqrt(2.0 * getPassport().physicalProperties.nu * getPassport().timeDiscretizationProperties.dt));
1421
1422
1423#pragma omp parallel for
1424 for (int i = 0; i < (int)wake->vtx.size(); ++i)
1425 {
1426 newPos[i] = wake->vtx[i].r()
1427 + (velocity->wakeVortexesParams.convVelo[i] +
1428 velocity->wakeVortexesParams.diffVelo[i] +
1430
1431 //newPos[i] = wake->vtx[i].r() + (velocity->wakeVortexesParams.convVelo[i] + getV0()) * passport.timeDiscretizationProperties.dt + Point2D{ d(gen), d(gen) };
1432 }
1433
1434 //for (int i = 0; i < (int)wake->vtx.size(); ++i)
1435 //{
1436 // std::cout << i << " " << wake->vtx[i].r() << " " << \
1437 // velocity->wakeVortexesParams.convVelo[i] << " " << \
1438 // velocity->wakeVortexesParams.diffVelo[i] << " " << \
1439 // getV0() << "\n";
1440 //}
1441
1442
1443 size_t counter = wake->vtx.size() - nVirtVortex;
1444 for (size_t bou = 0; bou < boundary.size(); ++bou)
1445 {
1446 for (int i = 0; i < (int)boundary[bou]->virtualWake.vtx.size(); ++i)
1447 {
1448 wake->vtx[counter].g() = boundary[bou]->virtualWake.vtx[i].g();
1449 ++counter;
1450 }
1451 }
1452
1453 //getTimers().stop("Move");
1454}//MoveVortexes(...)
Point2D getV0() const
Возврат текущей скорости набегающего потока
Definition World2D.h:142
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReserveMemoryForMatrixAndRhs()

void World2D::ReserveMemoryForMatrixAndRhs ( )

Вычисляем размер матрицы и резервируем память под нее и под правую часть

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

Definition at line 1213 of file World2D.cpp.

1214{
1215 //getTimers().start("Mem");
1216
1217 if (currentStep == 0)
1218 {
1219 dispBoundaryInSystem.resize(boundary.size());
1220 dispBoundaryInSystem[0] = 0;
1221
1222 for (size_t i = 1; i < boundary.size(); ++i)
1223 {
1224 dispBoundaryInSystem[i] = dispBoundaryInSystem[i - 1] + boundary[i - 1]->GetUnknownsSize() + 1;
1225 }
1226
1227 size_t matrSize = boundary.size();
1228 size_t matrSkosSize = 0;
1229
1230 for (auto it = boundary.begin(); it != boundary.end(); ++it)
1231 {
1232 matrSize += (*it)->GetUnknownsSize();
1233 matrSkosSize += (*it)->GetUnknownsSize();
1234 }
1235
1236 if ((getPassport().numericalSchemes.linearSystemSolver.second == 0 || getPassport().numericalSchemes.linearSystemSolver.second == 1 || (getPassport().numericalSchemes.linearSystemSolver.second == 2 && isAnyMovableOrDeformable())))
1237 {
1238 //matr.resize(matrSize, matrSize);
1239 matrReord.resize(matrSize, matrSize);
1240 //matrSkos.resize(matrSkosSize, matrSkosSize);
1241
1242 //matr.setZero();
1243 matrReord.setZero();
1244 //matrSkos.setZero();
1245
1246
1247 for (size_t i = 0; i < getNumberOfAirfoil(); ++i)
1248 {
1249 size_t nVari = boundary[i]->GetUnknownsSize();
1250 for (size_t j = 0; j < getNumberOfAirfoil(); ++j)
1251 {
1252 size_t nVarj = boundary[j]->GetUnknownsSize();
1253 IQ[i][j].first.resize(nVari, nVarj);
1254 IQ[i][j].second.resize(nVari, nVarj);
1255 }
1256 }
1257 }
1258
1259 //rhs.resize(matrSize);
1260 rhsReord.resize(matrSize);
1261 //rhsSkos.resize(matrSkosSize);
1262 }
1263 //rhs.setZero();
1264 rhsReord.setZero();
1265 //rhsSkos.setZero();
1266
1267 //getTimers().stop("Mem");
1268}//ReserveMemoryForMatrixAndRhs()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SolveLinearSystem()

void World2D::SolveLinearSystem ( )

Решение системы линейных алгебраических уравнений

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

Definition at line 826 of file World2D.cpp.

827{
828 getTimers().start("Solve");
829/*
830 if (currentStep == 0)
831 {
832 invMatr = matr;
833 std::ifstream fileMatrix("matrix.txt");
834 int nx, ny;
835 fileMatrix >> nx;
836 fileMatrix >> ny;
837 for (int i = 0; i < matr.rows(); ++i)
838 {
839 for (int j = 0; j < matr.cols(); ++j)
840 {
841 fileMatrix >> matr(i, j);
842 }
843 }
844 fileMatrix.close();
845 }
846*/
847
848/*
849 if (currentStep == 0)
850 {
851 std::ofstream fileMatrix(passport.dir + "/dbg/matrix.txt");
852 //int nx, ny;
853 //fileMatrix >> nx;
854 //fileMatrix >> ny;
855 fileMatrix.precision(16);
856 for (int i = 0; i < matrReord.rows(); ++i)
857 {
858 for (int j = 0; j < matrReord.cols(); ++j)
859 {
860 fileMatrix << matrReord(i, j) << " ";
861 }
862 fileMatrix << std::endl;
863 }
864 fileMatrix.close();
865 }
866
867
868
869 {
870 std::ofstream fileMatrix(passport.dir + "/dbg/rhs"+std::to_string(currentStep)+".txt");
871 //int nx, ny;
872 //fileMatrix >> nx;
873 //fileMatrix >> ny;
874 fileMatrix.precision(16);
875 for (int i = 0; i < rhsReord.size(); ++i)
876 {
877 fileMatrix << rhsReord(i) << std::endl;
878 }
879 fileMatrix.close();
880 }//*/
881
882 //exit(-100500);
883
885 const int linSystemScheme = passport.numericalSchemes.linearSystemSolver.second;
886
887 if (linSystemScheme == 0) // Gaussian elimination
888 {
889 double t1 = -omp_get_wtime();
890 if (useInverseMatrix && (currentStep == 0))
891 {
892 info('t') << "Inverting matrix... ";
893
894#if (defined(USE_CUDA))
895 invMatr.resize(matrReord.rows(), matrReord.cols());
896 for (int i = 0; i < (int)matrReord.rows(); ++i)
897 for (int j = 0; j < (int)matrReord.cols(); ++j)
898 invMatr(i, j) = (i == j) ? 1.0 : 0.0;
899 cuInverseMatrix((int)matrReord.rows(), matrReord.data(), invMatr.data());
900#else
901 invMatr = matrReord.inverse();
902#endif
903 info('t') << "done" << std::endl;
904 }
905
906 if (currentStep == 0)
907 info('t') << "Solving system at step #0... ";
908
910 {
911 /*
912 std::cout << "Solution via invMatrix" << std::endl;
913 for (int i = 0; i < invMatr.rows(); ++i)
914 for (int j = 0; j < invMatr.cols(); ++j)
915 if (fabs(invMatr(i, j)) > 1e+5)
916 std::cout << "invMatrixError: " << "invA(" << i << ", " << j << ") = " << invMatr(i, j) << std::endl;
917
918 for (int i = 0; i < rhsReord.size(); ++i)
919 if (fabs(rhsReord(i)) > 1e+5)
920 std::cout << "rhsReordError: " << "rhsReord(" << i << ") = " << rhsReord(i) << std::endl;
921 */
924 sol = invMatr * rhsReord;
926 }
927 else
928 {
931 sol = matrReord.partialPivLu().solve(rhsReord);
933 }
934
935 if (currentStep == 0)
936 info('t') << "done" << std::endl;
937
938 t1 += omp_get_wtime();
939
940 //std::cout << " Time in Gauss = " << t1 << std::endl;
941/*
942 std::ofstream solFile(getPassport().dir + "/sol" + std::to_string(currentStep) + "-Gauss.txt");
943 solFile.precision(16);
944 for (int i = 0; i < sol.size(); ++i)
945 solFile << sol(i) << std::endl;
946 solFile.close();
947 exit(10101);
948//*/
949 }//Gauss
950
951/*
952 if (currentStep == 0)
953 { N
954 invMatr = matr;
955 std::ifstream fileMatrix("invMatrix.txt");
956 int nx, ny;
957 fileMatrix >> nx;
958 fileMatrix >> ny;
959 for (int i = 0; i < invMatr.rows(); ++i)
960 {
961 for (int j = 0; j < invMatr.cols(); ++j)
962 {
963 fileMatrix >> invMatr(i, j);
964 }
965 }
966 fileMatrix.close();
967 }
968*/
969
970/*
971 if (currentStep == 0)
972 {
973 Eigen::MatrixXd mmul = matr * invMatr;
974 for (int i = 0; i < invMatr.rows(); ++i)
975 mmul(i,i) -= 1.0;
976
977 double maxElem = 0.0;
978 for (int i = 0; i < mmul.rows(); ++i)
979 for (int j = 0; j < mmul.cols(); ++j)
980 if (fabs(mmul(i,j)) > maxElem )
981 maxElem = mmul(i,j);
982
983 std::cout << "A*A^-1: MAX ELEMENT = " << maxElem << std::endl;
984 }
985*/
986
987/*
988 if (currentStep == 0)
989 {
990 std::ofstream fileMatrix("invMatrix4x3200.txt");
991 fileMatrix << invMatr.rows() << " " << invMatr.cols() << std::endl;
992 fileMatrix.precision(18);
993 for (int i = 0; i < invMatr.rows(); ++i)
994 {
995 for (int j = 0; j < invMatr.cols(); ++j)
996 {
997 fileMatrix << invMatr(i, j) << " ";
998 }
999 fileMatrix << std::endl;
1000 }
1001 fileMatrix.close();
1002 }
1003*/
1004
1005
1006
1007 if ((linSystemScheme == 1 || linSystemScheme == 2)) // GMRES
1008 {
1009 int nFullVars = (int)getNumberOfBoundary();
1010 for (int i = 0; i < getNumberOfBoundary(); ++i)
1011 nFullVars += (int)(boundary[i]->GetUnknownsSize());
1012
1013 std::vector<std::vector<double>> Ggam(getNumberOfBoundary());
1014 for (int i = 0; i < getNumberOfBoundary(); ++i)
1015 Ggam[i].resize(boundary[i]->GetUnknownsSize());
1016
1017 std::vector<double> GR(getNumberOfBoundary());
1018
1019#ifndef USE_CUDA
1020 std::vector<double> Grhs(rhsReord.size());
1021 for (int i = 0; i < rhsReord.size(); ++i)
1022 Grhs[i] = rhsReord(i);
1023
1024 std::vector<int> Gpos(getNumberOfBoundary(), 0);
1025 for (int i = 1; i < getNumberOfBoundary(); ++i)
1026 Gpos[i] = Gpos[i - 1] + (int)(boundary[i - 1]->GetUnknownsSize());
1027
1028 std::vector<int> Gvsize(getNumberOfBoundary());
1029 for (int i = 0; i < getNumberOfBoundary(); ++i)
1030 Gvsize[i] = (int)(boundary[i]->GetUnknownsSize());
1031
1033 {
1034 //for direct GMRES
1035 std::vector<double> Gmatr(nFullVars * nFullVars);
1036 for (size_t i = 0; i < nFullVars; ++i)
1037 for (size_t j = 0; j < nFullVars; ++j)
1038 Gmatr[i * nFullVars + j] = matrReord(i, j);
1039
1040 for (int i = 0; i < getNumberOfBoundary(); ++i)
1041 for (int j = 0; j < boundary[i]->GetUnknownsSize(); ++j)
1042 auto y_test = airfoil[i]->len[j % airfoil[i]->getNumberOfPanels()];
1043
1044 GMRES_Direct(*this, nFullVars, (int)getNumberOfBoundary(), Gmatr, Grhs, Gpos, Gvsize, Ggam, GR);
1045
1046 sol.resize(nFullVars);
1047 int cntr = 0;
1048 for (int i = 0; i < getNumberOfBoundary(); ++i)
1049 for (int j = 0; j < boundary[i]->GetUnknownsSize(); ++j)
1050 sol(cntr++) = Ggam[i][j] /*/ airfoil[i]->len[j % airfoil[i]->getNumberOfPanels()]*/;
1051 for (int i = 0; i < getNumberOfBoundary(); ++i)
1052 sol(cntr++) = GR[i];
1053
1054 //std::ofstream solFile(getPassport().dir + "/dbg/sol-direct-gmres" + std::to_string(currentStep) + ".txt");
1055 //solFile.precision(16);
1056 //for (int i = 0; i < sol.size(); ++i)
1057 // solFile << sol(i) << std::endl;
1058 //solFile.close();
1059 }
1060#else
1062 {
1063 getInfo('e') << "Direct GMRES for CUDA is not implemented" << std::endl;
1064 exit(-2);
1065 }
1066
1067 if ( (passport.numericalSchemes.linearSystemSolver.second == 2))
1068 {
1069 // for fast GMRES
1070 std::vector<std::vector<double>> GGrhs(getNumberOfBoundary());
1071 int cntrRhs = 0;
1072 for (int i = 0; i < getNumberOfBoundary(); ++i)
1073 GGrhs[i].resize(boundary[i]->GetUnknownsSize());
1074 for (int i = 0; i < getNumberOfBoundary(); ++i)
1075 for (int j = 0; j < boundary[i]->GetUnknownsSize(); ++j)
1076 GGrhs[i][j] = rhsReord(cntrRhs++);
1077
1078 std::vector<double> GrhsReg(getNumberOfBoundary());
1079 for (int i = 0; i < getNumberOfBoundary(); ++i)
1080 GrhsReg[i] = rhsReord[nFullVars - (getNumberOfBoundary() - i)];
1081
1082
1083 int niter;
1085
1086 double time_GMRES = -omp_get_wtime();
1087 GMRES(*this, Ggam, GR, GGrhs, GrhsReg, niter, linScheme);
1088 time_GMRES += omp_get_wtime();
1089 std::cout << "Time_GMRES = " << time_GMRES << std::endl;
1090
1091 sol.resize(nFullVars);
1092 int cntr = 0;
1093 for (int i = 0; i < getNumberOfBoundary(); ++i)
1094 for (int j = 0; j < boundary[i]->GetUnknownsSize(); ++j)
1095 sol(cntr++) = Ggam[i][j] / airfoil[i]->len[j % airfoil[i]->getNumberOfPanels()];
1096 for (int i = 0; i < getNumberOfBoundary(); ++i)
1097 sol(cntr++) = GR[i];
1098
1099 /*
1100 std::ofstream solFile(getPassport().dir + "/sol-fast-new-gmres" + std::to_string(currentStep) + ".txt");
1101 solFile.precision(16);
1102 for (int i = 0; i < sol.size(); ++i)
1103 solFile << sol(i) << std::endl;
1104 solFile.close();
1105 exit(1010110);
1106 //*/
1107 }
1108#endif
1109 }
1110
1111 getTimers().stop("Solve");
1112}//SolveLinearSystem()
Eigen::MatrixXd invMatr
Обратная матрица
Definition World2D.h:120
VMlib::vmTimer timerSlaeSolve
Definition World2D.h:355
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода
Definition WorldGen.h:82
void GMRES_Direct(const World2D &W, int nAllVars, int nafl, const std::vector< double > &mtrDir, const std::vector< double > &rhsDir, const std::vector< int > &pos, const std::vector< int > &vsize, std::vector< std::vector< double > > &gam, std::vector< double > &R)
Definition Gmres2D.cpp:206
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Step()

void World2D::Step ( )
overridevirtual

Функция выполнения предварительного шага

Основная функция выполнения одного шага по времени

Implements VMlib::WorldGen.

Definition at line 260 of file World2D.cpp.

261{
262 //try {
263 //Очистка статистики
265
266#ifdef USE_CUDA
267 cuSetCurrentStep((int)currentStep, 1);
268#endif // USE_CUDA
269
270
271 //Засечка времени в начале шага
272 getTimers().start("Step");
273
275
276 size_t countStrongCoupling = 0;
277 for (size_t m = 0; m < mechanics.size(); ++m)
278 {
279 MechanicsRigidOscillPart* mechVar = dynamic_cast<MechanicsRigidOscillPart*>(mechanics[m].get());
280 if (mechVar && getPassport().airfoilParams[m].addedMass.length2() > 0)
281 {
282 mechVar->getStrongCoupling() = true;
283 ++countStrongCoupling;
284 }
285 }
286
287 bool semiImplicitStrategy = ((countStrongCoupling == mechanics.size()) && (mechanics.size() > 0));
288 if ((currentStep == 0) && (semiImplicitStrategy))
289 info('i') << "Strong (semi-implicit) coupling strategy" << std::endl;
290
291//НЕПОДВИЖНЫЕ ТЕЛА
292//*
293 int nTotPan = 0;
294 for (size_t s = 0; s < getNumberOfAirfoil(); ++s)
295 nTotPan += (int)getAirfoil(s).getNumberOfPanels();
296
297 if (!semiImplicitStrategy)
298 {
300
303 if (getPassport().numericalSchemes.velocityComputation.second == 0 && getPassport().numericalSchemes.linearSystemSolver.second == 2)
304 {
305#ifdef USE_CUDA
306 cuda.RefreshAfls(3);
307 if (nTotPan > 0)
308 {
309 auto& afl = getAirfoil(0);
310 auto& treePnlVrt = *getCuda().inflTreePnlVortex;
311 if (getCurrentStep() == 0)
312 treePnlVrt.MemoryAllocate((int)getCuda().n_CUDA_pnls);
313
315 {
316 //Построение дерева для влияющих панелей (вихревых слоев)
317 treePnlVrt.UpdatePanelGeometry(nTotPan, (double4*)afl.devRPtr);
318 treePnlVrt.Build();
319 }
320
321 treePnlVrt.UpdatePanelAttachedVortexIntensity(afl.devAttachedVortexSheetPtr, afl.devAttachedVortexSheetLinPtr);
322 treePnlVrt.UpwardTraversal(multipoleOrder);
323 }
324#endif
325 }
326
327 if (getPassport().numericalSchemes.velocityComputation.second == 1)
328 {
329#ifdef USE_CUDA
330 cuda.RefreshWake(3);
331 cuda.RefreshAfls(3);
332
333 //Построение дерева для вихрей
334 auto& treeWake = *getCuda().inflTreeWake;
335 treeWake.MemoryAllocate((int)getCuda().n_CUDA_wake);
336 treeWake.Update((int)getWake().vtx.size(), getWake().devVtxPtr, getPassport().wakeDiscretizationProperties.eps);
337 treeWake.Build();
338 treeWake.UpwardTraversal(getPassport().numericalSchemes.nbodyMultipoleOrder);
339
340 if (nTotPan > 0)
341 {
342 auto& afl = getAirfoil(0);
343 auto& treePnl = *getCuda().cntrTreePnl;
344 auto& treePnlVrt = *getCuda().inflTreePnlVortex;
345 auto& treePnlSrc = *getCuda().inflTreePnlSource;
346 auto& treePnlAux = *getCuda().auxTreePnl;
347
348 if (getCurrentStep() == 0)
349 {
350 treePnl.MemoryAllocate((int)getCuda().n_CUDA_pnls);
351 treePnlAux.MemoryAllocate((int)getCuda().n_CUDA_pnls);
352 treePnlVrt.MemoryAllocate((int)getCuda().n_CUDA_pnls);
354 treePnlSrc.MemoryAllocate((int)getCuda().n_CUDA_pnls);
355 }
356
358 {
359 //Построение дерева для влияющих панелей (вихревых слоев)
360 treePnlVrt.UpdatePanelGeometry(nTotPan, (double4*)afl.devRPtr);
361 treePnlVrt.Build();
362
363 //Построение контрольного дерева панелей
364 treePnl.UpdatePanelGeometry((int)nTotPan, (double4*)afl.devRPtr);
365 treePnl.Build();
366
367 //Построение дерева для влияющих панелей (источников)
369 {
370 treePnlSrc.UpdatePanelGeometry(nTotPan, (double4*)afl.devRPtr);
371 treePnlSrc.UpdatePanelAttachedSourceIntensity(afl.devAttachedSourceSheetPtr, afl.devAttachedSourceSheetLinPtr);
372 treePnlSrc.Build();
373 treePnlSrc.UpwardTraversal(multipoleOrder);
374 }
375 }
376
377 treePnlVrt.UpdatePanelAttachedVortexIntensity(afl.devAttachedVortexSheetPtr, afl.devAttachedVortexSheetLinPtr);
378 treePnlVrt.UpwardTraversal(multipoleOrder);
379 }
380#endif
381 }
383
384 measureVP->Initialization();
386
387 //added masses
388 if (getPassport().physicalProperties.typeAccel.second == 3)
389 {
390 std::vector<Point2D> lambdaAdd(getNumberOfAirfoil(), {0.0, 0.0});
391 std::vector<double> muAdd(getNumberOfAirfoil());
392 for (size_t bou = 0; bou < getNumberOfAirfoil(); ++bou)
393 {
394 if (dynamic_cast<MechanicsRigidGivenLaw*>(&getNonConstMechanics(bou)))
395 {
396 for (size_t pnl = 0; pnl < getAirfoil(bou).getNumberOfPanels(); ++pnl)
397 {
398 double sumGam = getBoundary(bou).sheets.freeVortexSheet(pnl, 0) + getBoundary(bou).sheets.attachedVortexSheet(pnl, 0);
399 Point2D rpnl = 0.5 * (getAirfoil(bou).getR(pnl) + getAirfoil(bou).getR(pnl + 1));
400 lambdaAdd[bou] += rpnl.kcross() * sumGam * getAirfoil(bou).len[pnl];
401 muAdd[bou] += (rpnl - getAirfoil(bou).rcm).length2() * sumGam * getAirfoil(bou).len[pnl];
402 }
403 lambdaAdd[bou] *= getPassport().physicalProperties.rho;
404 muAdd[bou] *= 0.5 * getPassport().physicalProperties.rho;
405 info('i') << "Added Masses for airfoil #" << bou << " = { " << lambdaAdd[bou][0] << ", " << lambdaAdd[bou][1] << ", " << muAdd[bou] << " }" << std::endl;
406
407 char direction;
408 switch ((int)(getPassport().physicalProperties.timeAccel))
409 {
410 case 0:
411 direction = 'x';
412 break;
413 case 1:
414 direction = 'y';
415 break;
416 case 2:
417 direction = 'w';
418 break;
419 default:
420 direction = '?';
421 info('e') << "Wrong dirfection is specified!" << std::endl;
422 exit(1);
423 }
424
425 std::string addMassFileName = getPassport().dir + "addMass-" + std::to_string(getAirfoil(bou).numberInPassport);
426 std::ofstream addMassFile;
427 if (!VMlib::fileExistTest(addMassFileName, info) || getPassport().physicalProperties.timeAccel == 0)
428 {
429 addMassFile.open(addMassFileName);
430 addMassFile << "Added Masses for airfoil:" << std::endl;
431 }
432 else
433 addMassFile.open(addMassFileName, std::ios_base::app);
434
435 addMassFile << direction << "-direction: " << lambdaAdd[bou][0] << " " << lambdaAdd[bou][1] << " " << muAdd[bou] << std::endl;
436
437 addMassFile.close();
438 }
439 }
440 //info('i') << "Goodbye!" << std::endl;
441 //exit(0);
443 }
444
446 if (nTotPan > 0 && getPassport().numericalSchemes.velocityComputation.second == 1)
447 {
448 auto& afl = getAirfoil(0);
449#if USE_CUDA
450 getCuda().inflTreePnlVortex->UpdatePanelFreeAndAttachedVortexIntensity(afl.devFreeVortexSheetPtr, afl.devFreeVortexSheetLinPtr, afl.devAttachedVortexSheetPtr, afl.devAttachedVortexSheetLinPtr);
451 getCuda().inflTreePnlVortex->UpwardTraversal(multipoleOrder);
452#endif
453 }
455
456 //Вычисление скоростей вихрей: для тех, которые в следе, и виртуальных, а также в точках wakeVP
458
459//#include "gammaCirc.h"
460
461 //Расчет и сохранение поля давления
463 //optimizer
464 //if ((getCurrentStep() >= 0 && getCurrentStep() <= 12000) || getCurrentStep()==0) //2026-03-28
465 {
466 const double& vRef = getPassport().physicalProperties.vRef;
467 //scaleV = 1.0 / vRef;
468 double scaleP = 1.0 / (0.5 * getPassport().physicalProperties.rho * sqr(vRef));
469
470 measureVP->CalcPressure();
471
472 {
473 MechanicsDeformable* ptr = dynamic_cast<MechanicsDeformable*>(mechanics[0].get());
474 if (measureVP->getTotalNumberOfRealPoints() > 0)
475 measureVP->SaveVP();
476
477 if (ptr && !ptr->beam->fsi)
478
479 //сохранение главного вектора силы, вычисленного как интеграл от давления
480 //*
481 if (measureVP->elasticPoints.size() > 0)
482 {
483 std::ofstream presForcesFile;
484 if (currentStep == 0)
485 {
486 presForcesFile.open(getPassport().dir + "presForcesFile.csv");
487 presForcesFile << "time,Fx,Fy,Py" << std::endl;
488 }
489 else
490 presForcesFile.open(getPassport().dir + "presForcesFile.csv", std::ios_base::app);
491
492 auto r = measureVP->GetVPinElasticPoints();
493 Point2D presForce = { 0.0, 0.0 };
494 double yPower = 0.0;
495
496 for (int q = 0; q < r.size(); ++q)
497 {
498 presForce -= r[q].second * getAirfoil(0).len[q] * getAirfoil(0).nrm[q];
499 yPower -= (r[q].second * getAirfoil(0).len[q] * getAirfoil(0).nrm[q])[1] * (0.5 * (getAirfoil(0).getV(q) + getAirfoil(0).getV(q + 1))[1]);
500 }
501
502 presForcesFile << currentTime << "," << scaleP * presForce[0] << "," << scaleP * presForce[1] << "," << scaleP * yPower << std::endl;
503 presForcesFile.close();
504 }
505 }
506 //*/
507
508 //Коэффициенты разложения силы по балочным функциям
509 MechanicsDeformable* ptr = dynamic_cast<MechanicsDeformable*>(mechanics[0].get());
510 if (ptr && ptr->beam->fsi)
511 {
512 auto r = measureVP->GetVPinElasticPoints();
513
514 /*
515 std::ofstream elastPresFile;
516 elastPresFile.open(getPassport().dir + "elastPresFile-" + std::to_string(currentStep) + ".txt");
517 //if (currentStep == 0)
518 // elastPresFile.open(getPassport().dir + "elastPresFile.txt");
519 //else
520 // elastPresFile.open(getPassport().dir + "elastPresFile.txt", std::ios_base::app);
521 //elastPresFile << currentStep << " ";
522
523 for (size_t q = 0; q < r.size(); ++q)
524 elastPresFile << measureVP->elasticPoints[q][0] << " " << measureVP->elasticPoints[q][1] << " " << r[q].second << std::endl;
525 //elastPresFile << std::endl;
526 elastPresFile.close();
527 */
528
529 //Сохраняем давление на последних nLastSteps шагах по дисциплине очереди
530 std::vector<double> currentPres(ptr->chord.size());
531 for (size_t j = 0; j < ptr->chord.size(); ++j)
532 {
533 //currentPres[j] = -(ptr->beam->rho * 2 * ptr->beam->F); // gravity force for g = 2.0;
534 currentPres[j] = -(r[2 * j + 0].second - r[2 * j + 1].second);
535 }
536
537 if (ptr->beam->presLastSteps.size() < ptr->beam->nLastSteps)
538 ptr->beam->presLastSteps.push_back(currentPres);
539 else
540 {
541 for (int w = 1; w < ptr->beam->nLastSteps; ++w)
542 ptr->beam->presLastSteps[w - 1] = std::move(ptr->beam->presLastSteps[w]);
543 ptr->beam->presLastSteps.back() = currentPres;
544 }
545
546 for (int q = 0; q < ptr->beam->R; ++q)
547 {
548 ptr->beam->qCoeff[q] = 0;
549
550 if (ptr->beam->presLastSteps.size() == ptr->beam->nLastSteps)
551 {
552 for (size_t j = 0; j < ptr->chord.size(); ++j)
553 {
554 double averpres = 0.0;
555
556 //осредняем давление
557 for (int i = 0; i < ptr->beam->nLastSteps; ++i)
558 averpres += ptr->beam->presLastSteps[i][j];
559 averpres /= ptr->beam->nLastSteps;
560
561 ptr->beam->qCoeff[q] += -averpres * (ptr->initialChord[j].beg - ptr->initialChord[j].end).length() * ptr->beam->shape(q, 0.5 * (ptr->initialChord[j].beg + ptr->initialChord[j].end)[0]);
562 }
563 ptr->beam->qCoeff[q] /= (ptr->beam->intSqUnitShape * ptr->beam->L);
564 }
565
566 //std::cout << "q[" << q << "] = " << ptr->beam->qCoeff[q] << std::endl;
567 }
568
569 std::ofstream phiFile;
570 if (currentStep == 0)
571 {
572 phiFile.open(getPassport().dir + "phiFile.csv");
573 phiFile << "time";
574 for (int p = 0; p < ptr->beam->R; ++p)
575 phiFile << ",phi-" << std::to_string(p + 1);
576 phiFile << std::endl;
577 }
578 else
579 phiFile.open(getPassport().dir + "phiFile.csv", std::ios_base::app);
580
581 phiFile << currentTime;
582 for (int p = 0; p < ptr->beam->R; ++p)
583 phiFile << "," << ptr->beam->phi(p, 0);
584 phiFile << std::endl;
585
586 phiFile.close();
587 }
588//*/
589 //getInfo('e') << "Deformable airfoils are not supported now!";
590 }
591
592 //Вычисление сил, действующих на профиль и сохранение в файл
593 for (auto& mech : mechanics)
594 {
595 mech->GetHydroDynamForce();
596 mech->GenerateForcesString();
597 mech->GeneratePositionString();
598 }
599
600 //Движение вихрей (сброс вихрей + передвижение пелены)
602 wake->Restruct();
603 wake->SaveKadrVtk();
604
605 //std::string fname = getPassport().dir + "/airfoil-" + std::to_string(getAirfoil(0).numberInPassport) + "-" + std::to_string(currentStep) + ".pfl";
607 //{
608 // std::ofstream of(fname);
609 // for (size_t i = 0; i < getAirfoil(0).getNumberOfPanels(); ++i)
610 // of << getAirfoil(0).getR(i)[0] << " " << getAirfoil(0).getR(i)[1] << std::endl;
611 // of.close();
612 //}
613
614 }//if (!semiImplicitStrategy)
615
616
617//СОПРЯЖЕННАЯ ПОСТАНОВКА
618 if (semiImplicitStrategy)
619 {
621 measureVP->Initialization();
622
624
625 //Вычисление скоростей вихрей: для тех, которые в следе, и виртуальных, а также в точках wakeVP
627
628 //for (size_t s = 0; s < mechanics.size(); ++s) {// проверка: вычисление полной силы по циркуляциям новых вихрей
629 // mechanics[s]->GetHydroDynamForce();
630 // mechanics[s]->GenerateForcesString();
631 //}
632
633 //Расчет и сохранение поля давления
635 {
636 measureVP->CalcPressure();
637 measureVP->SaveVP();
638 }
639
640 WakeAndAirfoilsMotion(false); //профиль - кинематически, здесь Vold := V;
641 wake->Restruct();
642
645
646 //Вычисление сил, действующих на профиль и сохранение в файл
647
648 for (size_t m = 0; m < mechanics.size(); ++m)
649 {
650 mechanics[m]->GetHydroDynamForce();
651
652 MechanicsRigidOscillPart* mechVar = dynamic_cast<MechanicsRigidOscillPart*>(mechanics[m].get());
653 if (mechVar)
654 {
655 if (getPassport().airfoilParams[m].addedMass.length2() == 0)
656 {
657 getInfo('e') << "Added mass of the airfoil should be non-zero!" << std::endl;
658 exit(1);
659 }
660 mechVar->MoveOnlyVelo();
661 Point2D accel = (mechVar->getV() - mechVar->getVOld()) * (1.0 / getPassport().timeDiscretizationProperties.dt);
662 mechVar->hydroDynamForce -=
663 Point2D{ accel[0] * getPassport().airfoilParams[m].addedMass[0],
664 accel[1] * getPassport().airfoilParams[m].addedMass[1] };
665 mechVar->GenerateForcesString();
666 mechVar->GeneratePositionString();
667 }
668 else
669 exit(3333);
670 }
671 }
672
673
674 wake->SaveKadrVtk();
675//*/
677
678
679 //Засечка времени в конце шага
680 getTimers().stop("Step");
681 /*
682 std::cout << "nvt = " << nVtxBeforeMerging << std::endl;
683 std::cout << "tIBT = " << timerInitialBuild.duration() << std::endl;
684 std::cout << "tRHS = " << timerRhs.duration() << std::endl;
685 std::cout << "tFIL = " << timerFillMatrix.duration() << std::endl;
686 std::cout << "tLIN = " << timerSlaeSolve.duration() << std::endl;
687 std::cout << "tVEL = " << timerConvVelo.duration() << std::endl;
688 std::cout << "tINS = " << timerInside.duration() << std::endl;
689 std::cout << "tMRG = " << timerMerging.duration() << std::endl;
690 */
691
692 info('i') << "Step = " << getCurrentStep() \
693 << " PhysTime = " << getCurrentTime() \
694 << std::setprecision(3) \
695 << " StepTime = " << getTimers().durationStep() \
696 << std::setprecision(6) \
697 << std::endl;
698
700
702 currentStep += 1;
703 //}
704 //catch (...)
705 //{
706 // info('e') << "!!! Exception from unknown source !!!" << std::endl;
707 // exit(-1);
708 //}
709}//Step()
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition Airfoil2D.h:113
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition Airfoil2D.h:81
Point2D rcm
Положение центра масс профиля
Definition Airfoil2D.h:97
Sheet sheets
Слои на профиле
Definition Boundary2D.h:96
std::vector< ChordPanel > initialChord
std::vector< ChordPanel > chord
std::unique_ptr< Beam > beam
Point2D hydroDynamForce
Вектор гидродинамической силы и момент, действующие на профиль
void GeneratePositionString()
Сохранение строки со статистикой в файл нагрузок
void GenerateForcesString()
Сохранение строки со статистикой в файл нагрузок
Point2D & getV()
текущая скорость профиля
const double & attachedVortexSheet(size_t n, size_t moment) const
Definition Sheet2D.h:105
const double & freeVortexSheet(size_t n, size_t moment) const
Definition Sheet2D.h:100
VMlib::vmTimer timerInitialBuild
Definition World2D.h:352
void CalcPanelsVeloAndAttachedSheets()
Вычисление скоростей панелей и интенсивностей присоединенных слоев вихрей и источников
Definition World2D.cpp:1385
void CalcVortexVelo()
Вычисление скоростей (и конвективных, и диффузионных) вихрей (в пелене и виртуальных),...
Definition World2D.cpp:1273
void CalcAndSolveLinearSystem()
Набор матрицы, правой части и решение СЛАУ
Definition World2D.cpp:1506
void WakeAndAirfoilsMotion(bool dynamics)
Перемещение вихрей и профилей на шаге
Definition World2D.cpp:1653
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition World2D.h:180
bool ifDivisible(int val) const
Definition World2D.h:278
Mechanics & getNonConstMechanics(size_t i) const
Возврат неконстантной ссылки на объект механики
Definition World2D.h:219
void GenerateStatString(size_t stepNo, double curTime, size_t N)
Формирование очередной строки файла временной статистики
Definition TimesGen.cpp:111
void resetAll()
Сброс всех счетчиков
Definition TimesGen.cpp:81
double durationStep() const
Вывод счетчика всего шага в секундах
Definition TimesGen.h:184
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition WorldGen.h:99
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition numvector.h:510
const int multipoleOrder
Definition defs.h:65
T sqr(T x)
Умножение a на комплексно сопряженноe к b.
Definition defsBH.h:101
bool fileExistTest(std::string &fileName, LogStream &info, bool exitKey=false, const std::list< std::string > &extList={})
Проверка существования файла
Definition defs.h:342
double vRef
Референсная скорость
Definition Passport2D.h:81
double rho
Плотность потока
Definition Passport2D.h:75
Here is the call graph for this function:

◆ WakeAndAirfoilsMotion()

void World2D::WakeAndAirfoilsMotion ( bool  dynamics)

Перемещение вихрей и профилей на шаге

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

Definition at line 1653 of file World2D.cpp.

1654{
1655 std::vector<Point2D> newPos;
1656
1657 MoveVortexes(newPos);
1658
1659//#ifdef BRIDGE
1660// double totalForce = 0;
1661// for (size_t afl = 0; afl < airfoil.size(); ++afl)
1662// {
1663// totalForce += mechanics[afl]->hydroDynamForce[1];
1664// }
1665// totalForce *= 1.2;
1666//
1667// mechanics[0]->hydroDynamForce[1] = totalForce;
1668//#endif
1669
1670 oldAirfoil.resize(0);
1671 for (auto& afl : airfoil)
1672 {
1673 if (dynamic_cast<AirfoilRigid*>(afl.get()))
1674 oldAirfoil.emplace_back(new AirfoilGeometry(*afl));
1675
1676 if (dynamic_cast<AirfoilDeformable*>(afl.get()))
1677 oldAirfoil.emplace_back(new AirfoilGeometry(*afl));
1678
1679
1680#ifdef BRIDGE
1681 if (dynamics)
1682 {
1683 if (afl->numberInPassport == 0)
1684 {
1685 mechanics[afl->numberInPassport]->Move();
1686 }
1687 else
1688 {
1689 MechanicsRigidOscillPart* mechTest = dynamic_cast<MechanicsRigidOscillPart*>(mechanics[afl->numberInPassport].get());
1690 if (mechTest != nullptr)
1691 {
1692 Mechanics& mechGen = *mechanics[0];
1693 MechanicsRigidOscillPart& mech0 = dynamic_cast<MechanicsRigidOscillPart&>(mechGen);
1694
1695 Point2D dr = mech0.getR() - mech0.getROld();
1696 Point2D dv = mech0.getV() - mech0.getVOld();
1697
1698 double dphi = mech0.getPhi() - mech0.getPhiOld();
1699 double dw = mech0.getW() - mech0.getWOld();
1700
1701 Mechanics& mechGenI = *mechanics[afl->numberInPassport];
1702 MechanicsRigidOscillPart& mechI = dynamic_cast<MechanicsRigidOscillPart&>(mechGenI);
1703
1704 mechI.getROld() = mechI.getR();
1705 mechI.getVOld() = mechI.getV();
1706 mechI.getPhiOld() = mechI.getPhi();
1707 mechI.getWOld() = mechI.getW();
1708
1709 //std::cout << "afl = " << afl << ", dy = " << dy << std::endl;
1710
1711 airfoil[afl->numberInPassport]->Move(dr);
1712 airfoil[afl->numberInPassport]->Rotate(dphi);
1713
1714 mechI.getR() += dr;
1715 mechI.getV() += dv;
1716
1717 mechI.getPhi() += dphi;
1718 mechI.getW() += dw;
1719 }
1720 }
1721 }
1722#endif
1723
1724#ifdef INITIAL
1725 if (dynamics)
1726 mechanics[afl->numberInPassport]->Move();
1727 else
1728 {
1729 MechanicsRigidOscillPart* mech = dynamic_cast<MechanicsRigidOscillPart*>(mechanics[afl->numberInPassport].get());
1730 if (mech)
1731 mech->MoveKinematic();
1732 else
1733 exit(2222);
1734 }
1735#endif
1736 }//for
1737
1738#if defined(__CUDACC__) || defined(USE_CUDA)
1739 cuda.RefreshAfls(2);
1740#endif
1741
1742 for (auto& bou : boundary)
1743 bou->virtualWake.vtx.clear();
1744
1745
1746 CheckInside(newPos, oldAirfoil);
1747
1748
1749 //передача новых положений вихрей в пелену
1750 for (size_t i = 0; i < wake->vtx.size(); ++i)
1751 wake->vtx[i].r() = newPos[i];
1752
1753// getWake().SaveKadrVtk();
1754
1755}//WakeAndAirfoilsMotion()
Класс, определяющий форму профиля
Definition Airfoil2D.h:66
Абстрактный класс, определяющий вид механической системы
Definition Mechanics2D.h:72
Point2D & getR()
текущее отклонение профиля
double & getPhi()
текущий угол поворота профиля
double & getW()
текущая угловая скорость профиля
void MoveVortexes(std::vector< Point2D > &newPos)
Вычисляем новые положения вихрей (в пелене и виртуальных)
Definition World2D.cpp:1401
void CheckInside(std::vector< Point2D > &newPos, const std::vector< std::unique_ptr< AirfoilGeometry > > &oldAirfoil)
Проверка проникновения вихрей внутрь профиля
Definition World2D.cpp:713
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ airfoil

std::vector<std::unique_ptr<Airfoil> > VM2D::World2D::airfoil
private

Список умных указателей на обтекаемые профили

Definition at line 77 of file World2D.h.

◆ boundary

std::vector<std::unique_ptr<Boundary> > VM2D::World2D::boundary
private

Список умных указателей на формирователи граничных условий на профилях

Definition at line 83 of file World2D.h.

◆ check01

int VM2D::World2D::check01
mutable

Definition at line 105 of file World2D.h.

◆ check02

int VM2D::World2D::check02
mutable

Definition at line 106 of file World2D.h.

◆ checkPan

int VM2D::World2D::checkPan
mutable

Definition at line 107 of file World2D.h.

◆ cuda

Gpu VM2D::World2D::cuda
mutableprivate

Объект, управляющий графическим ускорителем

Definition at line 137 of file World2D.h.

◆ currentStep

size_t VMlib::WorldGen::currentStep
protectedinherited

Текущий номер шага в решаемой задаче

Definition at line 69 of file WorldGen.h.

◆ currentTime

double VMlib::WorldGen::currentTime
protectedinherited

Текущее время в решаемой задаче

Definition at line 72 of file WorldGen.h.

◆ dispBoundaryInSystem

std::vector<size_t> VM2D::World2D::dispBoundaryInSystem
private

Список номеров, с которых начинаются элементы правой части (или матрицы) системы для профилей

Definition at line 86 of file World2D.h.

◆ gabb

int VM2D::World2D::gabb
mutable

Definition at line 104 of file World2D.h.

◆ info

LogStream VMlib::WorldGen::info
mutableprotectedinherited

Поток для вывода логов и сообщений об ошибках

Definition at line 60 of file WorldGen.h.

◆ invMatr

Eigen::MatrixXd VM2D::World2D::invMatr
private

Обратная матрица

Definition at line 120 of file World2D.h.

◆ IQ

std::vector<std::vector<std::pair<Eigen::MatrixXd, Eigen::MatrixXd> > > VM2D::World2D::IQ
private

Матрица, состоящая из пар матриц, в которых хранятся касательные и нормальные компоненты интегралов от ядра

Definition at line 117 of file World2D.h.

◆ matrReord

Eigen::MatrixXd VM2D::World2D::matrReord
private

Матрица системы

Definition at line 113 of file World2D.h.

◆ measureVP

std::unique_ptr<MeasureVP> VM2D::World2D::measureVP
private

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

Definition at line 101 of file World2D.h.

◆ mechanics

std::vector<std::unique_ptr<Mechanics> > VM2D::World2D::mechanics
private

Список умных указателей на типы механической системы для каждого профиля

Definition at line 89 of file World2D.h.

◆ nVtxBeforeMerging

size_t VMlib::WorldGen::nVtxBeforeMerging
inherited

Definition at line 76 of file WorldGen.h.

◆ oldAirfoil

std::vector<std::unique_ptr<AirfoilGeometry> > VM2D::World2D::oldAirfoil
private

Список умных указателей на обтекаемые профили для сохранения старого положения

Definition at line 80 of file World2D.h.

◆ passport

const Passport& VM2D::World2D::passport
private

Константная ссылка на паспорт конкретного расчета

Definition at line 134 of file World2D.h.

◆ passportGen

const PassportGen& VMlib::WorldGen::passportGen
protectedinherited

Константная ссылка на паспорт конкретного расчета

Definition at line 63 of file WorldGen.h.

◆ rhsReord

Eigen::VectorXd VM2D::World2D::rhsReord
private

Правая часть системы

Definition at line 127 of file World2D.h.

◆ sol

Eigen::VectorXd VM2D::World2D::sol
private

Решение системы

Definition at line 131 of file World2D.h.

◆ source

std::unique_ptr<WakeDataBase> VM2D::World2D::source
private

Умный указатель на источники

Definition at line 98 of file World2D.h.

◆ timerConvVelo

VMlib::vmTimer VM2D::World2D::timerConvVelo

Definition at line 357 of file World2D.h.

◆ timerFillMatrix

VMlib::vmTimer VM2D::World2D::timerFillMatrix

Definition at line 354 of file World2D.h.

◆ timerInitialBuild

VMlib::vmTimer VM2D::World2D::timerInitialBuild

Definition at line 352 of file World2D.h.

◆ timerInside

VMlib::vmTimer VM2D::World2D::timerInside

Definition at line 358 of file World2D.h.

◆ timerMerging

VMlib::vmTimer VM2D::World2D::timerMerging

Definition at line 359 of file World2D.h.

◆ timerRhs

VMlib::vmTimer VM2D::World2D::timerRhs

Definition at line 353 of file World2D.h.

◆ timers

std::unique_ptr<TimersGen> VMlib::WorldGen::timers
protectedinherited

Сведения о временах выполнения основных операций

Definition at line 66 of file WorldGen.h.

◆ timerSlaeSolve

VMlib::vmTimer VM2D::World2D::timerSlaeSolve

Definition at line 355 of file World2D.h.

◆ useInverseMatrix

bool VM2D::World2D::useInverseMatrix
private

Признак использования обратной матрицы

Definition at line 123 of file World2D.h.

◆ velocity

std::unique_ptr<Velocity> VM2D::World2D::velocity
private

Умный укзатель на объект, определяющий методику вычисления скоростей

Definition at line 92 of file World2D.h.

◆ wake

std::unique_ptr<Wake> VM2D::World2D::wake
private

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

Definition at line 95 of file World2D.h.


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