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

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

#include <World2D.h>

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

Public Member Functions

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

Public Attributes

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

Protected Attributes

LogStream info
 Поток для вывода логов и сообщений об ошибках More...
 
const PassportGen & passportGen
 Константная ссылка на паспорт конкретного расчета More...
 
std::unique_ptr< TimesGen > timestat
 Сведения о временах выполнения основных операций More...
 

Private Attributes

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

Detailed Description

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

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

Definition at line 68 of file World2D.h.

Constructor & Destructor Documentation

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

Конструктор

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

Definition at line 67 of file World2D.cpp.

67  :
68  WorldGen(passport_),
69  passport(dynamic_cast<const Passport&>(passport_)),
70  cuda(Gpu(*this))
71 {
72  std::stringstream ss;
73  ss << "#" << passport.problemNumber << " (" << passport.problemName << ")";
75 
77  currentStep = 0;
78 
79  timestat.reset(new Times(*this)),
80 
81  wake.reset(new Wake(*this));
82  // загрузка пелены из файла
84  wake->ReadFromFile(passport.wakesDir, passport.wakeDiscretizationProperties.fileWake); //Считываем из каталога с пеленой
85 
86  source.reset(new WakeDataBase(*this));
87  // загрузка положений источников из файла
89  source->ReadFromFile(passport.dir, passport.wakeDiscretizationProperties.fileSource); //Считываем из текущего каталога
90 
91  //считываем массив точек для подсчета и вывода поля скоростей и давлений
92  measureVP.reset(new MeasureVP(*this));
94  measureVP->ReadPointsFromFile(passport.dir);
95 
97  {
98  case 0:
99  velocity.reset(new VelocityBiotSavart(*this));
100  break;
101  //case 1:
102  // velocity.reset(new VelocityBarnesHut(*this));
103  // break;
104  }
105 
106  velocity->virtualVortexesParams.resize(passport.airfoilParams.size());
107 
108  for (size_t i = 0; i < passport.airfoilParams.size(); ++i)
109  {
110  switch (passport.numericalSchemes.panelsType.second)
111  {
112  case 0:
113  airfoil.emplace_back(new AirfoilRect(*this, i));
114  break;
115  //case 1:
116  // airfoil.emplace_back(new AirfoilCurv(*this, i));
117  // break;
118  }
119 
120  airfoil[i]->ReadFromFile(passport.airfoilsDir); //Считываем из каталога с коллекцией профилей
121 
123  {
124  case 0:
125  boundary.emplace_back(new BoundaryConstLayerAver(*this, i));
126  break;
127 
128  case 1:
129  boundary.emplace_back(new BoundaryLinLayerAver(*this, i));
130  break;
131 
132  case 10:
133  boundary.emplace_back(new BoundaryVortexCollocN(*this, i));
134  //info('e') << "BoundaryMDV is not implemented now! " << std::endl;
135  //exit(1);
136  break;
137 
138  default:
139  info('e') << "Unknown scheme!" << std::endl;
140  exit(1);
141  }
142 
143 
144  switch (passport.airfoilParams[i].mechanicalSystemType)
145  {
146  case 0:
147  mechanics.emplace_back(new MechanicsRigidImmovable(*this, i));
148  break;
149  case 1:
150  mechanics.emplace_back(new MechanicsRigidGivenLaw(*this, i));
151  break;
152  case 2:
153  mechanics.emplace_back(new MechanicsRigidOscillPart(*this, i));
154  break;
155 
156  case 3:
157  mechanics.emplace_back(new MechanicsRigidRotatePart(*this, i));
158  break;
159  /*
160  case 4:
161  mechanics.emplace_back(new MechanicsRigidOscillMon(*this, i));
162  break;
163  */
164  }
165 
166  }
167 
168  IQ.resize(passport.airfoilParams.size());
169  for (size_t i = 0; i < passport.airfoilParams.size(); ++i)
170  IQ[i].resize(passport.airfoilParams.size());
171 
172 
173  info.endl();
174  info('i') << "Start solving problem " << passport.problemName << std::endl;
175  info.endl();
176 
177  VMlib::PrintLogoToStream(info('_') << std::endl);
178 }//World2D(...)
std::unique_ptr< Velocity > velocity
Умный укзатель на объект, определяющий методику вычисления скоростей
Definition: World2D.h:87
static std::ostream * defaultWorld2DLogStream
Поток вывода логов и ошибок задачи
Definition: defs.h:213
std::pair< std::string, int > velocityComputation
Definition: Passport2D.h:190
Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле ...
std::string problemName
Название задачи
Definition: PassportGen.h:121
Класс для сбора статистики времени исполнения основных шагов алгоритма и вывода ее в файл ...
Definition: Times2D.h:59
std::pair< std::string, int > boundaryCondition
Метод аппроксимации граничных условий
Definition: Passport2D.h:199
void assignStream(std::ostream *pStr_, const std::string &label_)
Связывание потока логов с потоком вывода
Definition: LogStream.h:77
LogStream info
Поток для вывода логов и сообщений об ошибках
Definition: WorldGen.h:59
Класс, определяющий вид механической системы
Класс, определяющий вид механической системы
std::string dir
Рабочий каталог задачи
Definition: PassportGen.h:118
WorldGen(const VMlib::PassportGen &passport_)
Конструктор
Definition: WorldGen.cpp:44
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition: World2D.h:90
Класс, обеспечивающий возможность выполнения вычислений на GPU по технологии Nvidia CUDA...
Definition: Gpu2D.h:67
std::unique_ptr< WakeDataBase > source
Умный указатель на источники
Definition: World2D.h:93
Класс, определяющий вид механической системы
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
std::unique_ptr< MeasureVP > measureVP
Умный указатель на алгоритм вычисления полей скоростей и давления (для сохранения в файл) ...
Definition: World2D.h:96
Класс, определяющий тип обтекаемого профиля
Definition: Airfoil2DRect.h:64
Класс, опеделяющий вихревой след (пелену)
Definition: Wake2D.h:59
Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле ...
std::unique_ptr< TimesGen > timestat
Сведения о временах выполнения основных операций
Definition: WorldGen.h:65
double timeStart
Начальное время
Definition: PassportGen.h:59
std::string fileSource
Имя файла с положениями источников (без полного пути)
Definition: Passport2D.h:163
int saveVPstep
Шаг вычисления и сохранения скорости и давления
Definition: PassportGen.h:78
size_t problemNumber
Номер задачи
Definition: PassportGen.h:124
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition: Passport2D.h:302
void PrintLogoToStream(std::ostream &str)
Передача в поток вывода шапки программы VM2D/VM3D.
Definition: defs.cpp:105
Gpu cuda
Объект, управляющий графическим ускорителем
Definition: World2D.h:123
Класс, определяющий способ удовлетворения граничного условия на обтекаемом профиле ...
std::string fileWake
Имя файла с начальным состоянием вихревого следа (без полного пути)
Definition: Passport2D.h:160
Класс, отвечающий за вычисление поля скорости и давления в заданых точках для вывода ...
Definition: MeasureVP2D.h:64
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72
void setCurrTime(double t_) const
Установка текущего времени
Definition: Passport2D.h:108
std::pair< std::string, int > panelsType
Тип панелей (0 — прямые)
Definition: Passport2D.h:196
std::vector< std::vector< std::pair< Eigen::MatrixXd, Eigen::MatrixXd > > > IQ
Матрица, состоящая из пар матриц, в которых хранятся касательные и нормальные компоненты интегралов о...
Definition: World2D.h:105
Класс, определяющий вид механической системы
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
std::vector< AirfoilParams > airfoilParams
Список структур с параметрами профилей
Definition: Passport2D.h:280
Класс, опеделяющий набор вихрей
std::string wakesDir
Каталог с файлами вихревых следов
Definition: Passport2D.h:277
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
void endl()
Вывод в поток логов пустой строки
Definition: LogStream.h:100
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78
Класс, определяющий способ вычисления скоростей
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition: World2D.h:120
std::string airfoilsDir
Каталог с файлами профилей
Definition: Passport2D.h:274

Here is the call graph for this function:

Here is the caller graph for this function:

VM2D::World2D::~World2D ( )
inline

Деструктор

Definition at line 306 of file World2D.h.

306 { };

Here is the call graph for this function:

Member Function Documentation

void World2D::CalcAndSolveLinearSystem ( )

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

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

Definition at line 774 of file World2D.cpp.

775 {
776  if (airfoil.size() > 0)
777  {
779 
780 #if defined(__CUDACC__) || defined(USE_CUDA)
783 
785 
786  if ((sch == 0) || (sch == 1) || (sch == 10))
787  cuda.setSchemeSwitcher(sch + 1);
788  else
789  {
790  info('e') << "schemeSwitcher is not 0, or 1, or 10! " << std::endl;
791  exit(1);
792  }
793 
794  cuda.RefreshWake(1);
795  cuda.RefreshAfls(1);
796  cuda.RefreshVirtualWakes(1);
797 #endif
798 
799  if (currentStep == 0)
800  {
801  useInverseMatrix = (
802  (mechanics.size() == 1)
803  ||
804  (mechanics.size() > 1 && !std::any_of(mechanics.begin(), mechanics.end(), [](const std::unique_ptr<Mechanics>& m) { return m->isMoves; }))
805  );
806  }
807 
808  FillIQ();
810 
811  //{
812  // std::stringstream ss;
813  // ss << "IQ1-" << currentStep;
814  // std::ofstream of(passport.dir + "dbg/" + ss.str());
815  // for (size_t i = 0; i < (IQ[0][0]).first.rows(); ++i)
816  // {
817  // for (size_t j = 0; j < (IQ[0][0]).first.cols(); ++j)
818  // of << (IQ[0][0]).first(i, j) << " ";
819  // of << std::endl;
820  // }
821  // of.close();
822  //}
823 
824  /*
825  {
826  std::stringstream ss;
827  ss << "matr-" << currentStep;
828  std::ofstream of(passport.dir + "dbg/" + ss.str());
829  for (size_t i = 0; i < matr.rows(); ++i)
830  {
831  for (size_t j = 0; j < matr.cols(); ++j)
832  of << matr(i, j) << " ";
833  of << std::endl;
834  }
835  of.close();
836  }
837  */
838 
840 
841  getTimestat().timeOther.first += omp_get_wtime();
842 
843  size_t currentRow = 0;
844  for (size_t bou = 0; bou < boundary.size(); ++bou)
845  {
846  size_t nVars = boundary[bou]->GetUnknownsSize();
847  Eigen::VectorXd locSol;
848  locSol.resize(nVars);
849  for (size_t i = 0; i < nVars; ++i)
850  locSol(i) = sol(currentRow + i);
851 
852  boundary[bou]->SolutionToFreeVortexSheetAndVirtualVortex(locSol);
853  currentRow += nVars + 1;
854  }
855 
856  if (currentStep == 0)
857  {
858  Point2D addMass = { 0.0, 0.0 };
859  for (size_t q = 0; q < airfoil[0]->getNumberOfPanels(); ++q)
860  {
861  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];
862  }
863  addMass *= passport.physicalProperties.rho;
864 
865  std::cout << "AddMass = " << addMass << std::endl;
866  //exit(-42);
867  }
868 
869  getTimestat().timeOther.second += omp_get_wtime();
870  }
871 }//CalcAndSolveLinearSystem()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
void setSchemeSwitcher(int schemeSwitcher_)
Установка переключателя расчетных схем
Definition: Gpu2D.h:278
std::pair< std::string, int > boundaryCondition
Метод аппроксимации граничных условий
Definition: Passport2D.h:199
LogStream info
Поток для вывода логов и сообщений об ошибках
Definition: WorldGen.h:59
double maxGamma
Максимально допустимая циркуляция вихря
Definition: Passport2D.h:157
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
void setMaxGamma(double gam_)
Установка максимально допустимой циркуляции вихря
Definition: Gpu2D.h:264
double accelCft() const
Функция-множитель, позволяющая моделировать разгон
Definition: Passport2D.cpp:51
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
double rho
Плотность потока
Definition: Passport2D.h:75
void ReserveMemoryForMatrixAndRhs()
Вычисляем размер матрицы и резервируем память под нее и под правую часть
Definition: World2D.cpp:520
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition: Passport2D.h:302
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
void setAccelCoeff(double cft_)
Установка коэффициента разгона потока
Definition: Gpu2D.h:241
Gpu cuda
Объект, управляющий графическим ускорителем
Definition: World2D.h:123
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72
bool useInverseMatrix
Признак использования обратной матрицы
Definition: World2D.h:111
Eigen::VectorXd sol
Решение системы
Definition: World2D.h:117
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
void FillIQ()
Заполнение матрицы, состоящей из интегралов от (r-xi) / |r-xi|^2.
Definition: World2D.cpp:732
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
void FillMatrixAndRhs()
Заполнение матрицы системы для всех профилей
Definition: World2D.cpp:438
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78
void SolveLinearSystem()
Решение системы линейных алгебраических уравнений
Definition: World2D.cpp:333
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition: World2D.h:120
timePeriod timeOther
Все прочее
Definition: Times2D.h:113

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::CalcPanelsVeloAndAttachedSheets ( )

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

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

Definition at line 664 of file World2D.cpp.

665 {
666  getTimestat().timeOther.first += omp_get_wtime();
667 
668  //вычисляем скорости панелей
669  for (size_t i = 0; i < airfoil.size(); ++i)
670  mechanics[i]->VeloOfAirfoilPanels(currentStep * passport.timeDiscretizationProperties.dt);
671 
672  //вычисляем интенсивности присоединенных слоев
673  for (size_t i = 0; i < airfoil.size(); ++i)
674  boundary[i]->ComputeAttachedSheetsIntensity();
675 
676  getTimestat().timeOther.second += omp_get_wtime();
677 }//CalcPanelsVeloAndAttachedSheets(...)
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
double dt
Шаг по времени
Definition: PassportGen.h:65
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition: World2D.h:120
timePeriod timeOther
Все прочее
Definition: Times2D.h:113

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::CalcVortexVelo ( bool  shiftTime)

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

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

Definition at line 563 of file World2D.cpp.

564 {
565  getTimestat().timeCalcVortexConvVelo.first += omp_get_wtime();
566 
567  velocity->ResizeAndZero();
568 
569  //Конвективные скорости всех вихрей (и в пелене, и виртуальных), индуцируемые вихрями в пелене
570  //Подготовка CUDA
571 #if defined(__CUDACC__) || defined(USE_CUDA)
572  cuda.RefreshWake(2);
573  cuda.RefreshAfls(2);
574  cuda.RefreshVirtualWakes(2);
575 #endif
576  getTimestat().timeCalcVortexConvVelo.second += omp_get_wtime();
577 
578  //Вычисление конвективных скоростей вихрей (и в пелене, и виртуальных), а также в точках wakeVP
579  if (shiftTime)
580  --currentStep;
581 
582  velocity->CalcConvVelo();
583 
584  if (shiftTime)
585  ++currentStep;
586 
587  //Расчет средних значений eps для каждой панели и их передача на видеокарту
588  getTimestat().timeCalcVortexDiffVelo.first += omp_get_wtime();
589  for (size_t bou = 0; bou < getNumberOfBoundary(); ++bou)
591 
592 #if defined(__CUDACC__) || defined(USE_CUDA)
593  for (size_t i = 0; i < airfoil.size(); ++i)
594  cuda.CopyMemToDev<double, 1>(airfoil[i]->getNumberOfPanels(), airfoil[i]->meanEpsOverPanel.data(), airfoil[i]->devMeanEpsOverPanelPtr);
595 #endif
596  getTimestat().timeCalcVortexDiffVelo.second += omp_get_wtime();
597 
598  //Вычисление диффузионных скоростей вихрей (и в пелене, и виртуальных)
599  velocity->CalcDiffVelo();
600 
601  getTimestat().timeSaveKadr.first += omp_get_wtime();
602  /*//Сохранение всех параметров для вихрей в пелене
603 
604  {
605  VMlib::CreateDirectory(passport.dir, "dbg");
606  std::ostringstream sss;
607  sss << "prmWake";
608  sss << currentStep;
609  std::ofstream prmtFile(passport.dir + "dbg/" + sss.str());
610  prmtFile << "i x y g epsast convVeloX convVeloY diffVeloX diffVeloY I0 I1 I2X I2Y I3X I3Y" << std::endl;
611  for (size_t i = 0; i < wake->vtx.size(); ++i)
612  prmtFile << i << " " \
613  << wake->vtx[i].r()[0] << " " << wake->vtx[i].r()[1] << " " \
614  << wake->vtx[i].g() << " " \
615  << velocity->wakeVortexesParams.epsastWake[i] << " " \
616  << velocity->wakeVortexesParams.convVelo[i][0] << " " << velocity->wakeVortexesParams.convVelo[i][1] << " "\
617  << velocity->wakeVortexesParams.diffVelo[i][0] << " " << velocity->wakeVortexesParams.diffVelo[i][1] << " "\
618  //<< std::endl;
619  << velocity->wakeVortexesParams.I0[i] << " " \
620  << velocity->wakeVortexesParams.I1[i] << " " \
621  << velocity->wakeVortexesParams.I2[i][0] << " " << velocity->wakeVortexesParams.I2[i][1] << " " \
622  << velocity->wakeVortexesParams.I3[i][0] << " " << velocity->wakeVortexesParams.I3[i][1] << " " \
623  << std::endl;
624 
625  prmtFile.close();
626  }
627 //*/
628 
629 /*
630  //Сохранение всех параметров для виртуальных вихрей
631  {
632  for (size_t b = 0; b < boundary.size(); ++b)
633  {
634  std::ostringstream sss;
635  sss << "prmVirtual_";
636  sss << b << "-";
637  sss << currentStep;
638  std::ofstream prmFileVirt(passport.dir + "dbg/" + sss.str());
639  prmFileVirt << "i x y g epsast convVeloX convVeloY diffVeloX diffVeloY I0 I1 I2X I2Y I3X I3Y" << std::endl;
640  for (size_t i = 0; i < boundary[b]->virtualWake.vtx.size(); ++i)
641  prmFileVirt << i << " " \
642  << boundary[b]->virtualWake.vtx[i].r()[0] << " " << boundary[b]->virtualWake.vtx[i].r()[1] << " " \
643  << boundary[b]->virtualWake.vtx[i].g() << " " \
644  << velocity->virtualVortexesParams[b].epsastWake[i] << " " \
645  << velocity->virtualVortexesParams[b].convVelo[i][0] << " " << velocity->virtualVortexesParams[b].convVelo[i][1] << " "\
646  << velocity->virtualVortexesParams[b].diffVelo[i][0] << " " << velocity->virtualVortexesParams[b].diffVelo[i][1] << " "\
647  //<< std::endl;
648  << velocity->virtualVortexesParams[b].I0[i] << " " \
649  << velocity->virtualVortexesParams[b].I1[i] << " " \
650  << velocity->virtualVortexesParams[b].I2[i][0] << " " << velocity->virtualVortexesParams[b].I2[i][1] << " " \
651  << velocity->virtualVortexesParams[b].I3[i][0] << " " << velocity->virtualVortexesParams[b].I3[i][1] << " " \
652  << std::endl;
653  prmFileVirt.close();
654  }
655  //if (currentStep==2) exit(-123);
656  }
657 */
658  getTimestat().timeSaveKadr.second += omp_get_wtime();
659 }//CalcVortexVelo()
std::unique_ptr< Velocity > velocity
Умный укзатель на объект, определяющий методику вычисления скоростей
Definition: World2D.h:87
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
timePeriod timeCalcVortexDiffVelo
Начало и конец процесса вычисления диффузионных скоростей вихрей
Definition: Times2D.h:89
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
timePeriod timeCalcVortexConvVelo
Начало и конец процесса вычисления конвективных скоростей вихрей
Definition: Times2D.h:86
timePeriod timeSaveKadr
Начало и конец процесса сохранения кадра в файл
Definition: Times2D.h:107
Gpu cuda
Объект, управляющий графическим ускорителем
Definition: World2D.h:123
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72
void calcMeanEpsOverPanel()
Вычисление средних значений eps на панелях
Definition: Airfoil2D.cpp:82
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
Airfoil & getNonConstAirfoil(size_t i) const
Возврат неконстантной ссылки на объект профиля
Definition: World2D.h:142

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

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

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

Definition at line 294 of file World2D.cpp.

295 {
296  getTimestat().timeCheckInside.first += omp_get_wtime();
297 
298 
299 #if (!defined(USE_CUDA))
300  for (size_t afl = 0; afl < airfoil.size(); ++afl)
301  wake->Inside(newPos, *airfoil[afl], mechanics[afl]->isMoves, *oldAirfoil[afl]);
302 #else
303  if (newPos.size() > 0)
305  {
306  double* devNewpos_ptr;
307  cuReserveDevMem((void*&)devNewpos_ptr, newPos.size() * sizeof(double) * 2, 0);
308  cuCopyFixedArray(devNewpos_ptr, newPos.data(), newPos.size() * sizeof(double) * 2, 0);
309 
310  for (size_t afl = 0; afl < airfoil.size(); ++afl)
311  {
312  std::vector<double> gamma(airfoil[afl]->getNumberOfPanels(), 0.0);
313  std::vector<unsigned int> hit = lbvh_check_inside((int)currentStep, (int)newPos.size(), devNewpos_ptr, (int)airfoil[afl]->getNumberOfPanels(), airfoil[afl]->devRPtr);
314  for (int i = 0; i < newPos.size(); ++i)
315  {
316  if (hit[i] != (unsigned int)(-1))
317  {
318  gamma[hit[i]] += wake->vtx[i].g();
319  wake->vtx[i].g() = 0.0;
320  }
321  }
322  airfoil[afl]->gammaThrough = gamma;
323  }//for afl
324 
325  cuDeleteFromDev(devNewpos_ptr);
326  }
327 #endif
328 
329  getTimestat().timeCheckInside.second += omp_get_wtime();
330 }
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition: World2D.h:90
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
timePeriod timeCheckInside
Начало и конец процесса контроля протыкания
Definition: Times2D.h:98
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::FillIQ ( )

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

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

Definition at line 732 of file World2D.cpp.

733 {
734  getTimestat().timeFillMatrixAndRhs.first += omp_get_wtime();
735 
736  for (size_t bou = 0; bou < boundary.size(); ++bou)
737  {
738  if (currentStep == 0 || mechanics[bou]->isDeform)
739  {
740  boundary[bou]->FillIQSelf(IQ[bou][bou]);
741  }
742 
743  for (size_t oth = 0; oth < boundary.size(); ++oth)
744  {
745 
746 #ifdef INITIAL
747  if (currentStep == 0 || !useInverseMatrix)
748  {
749  //size_t nVarsOther = boundary[oth]->GetUnknownsSize();
750  if (bou != oth)
751  {
752  //std::cout << "matr!" << std::endl;
753  boundary[bou]->FillIQFromOther(*boundary[oth], IQ[bou][oth]);
754  }
755  }// if (currentStep == 0 || !useInverseMatrix)
756 #endif
757 
758 #ifdef BRIDGE
759  if (currentStep == 0)
760  {
761  size_t nVarsOther = boundary[oth]->GetUnknownsSize();
762  if (bou != oth)
763  boundary[bou]->FillIQFromOther(*boundary[oth], IQ[bou][oth]);
764  }// if (currentStep == 0)
765 #endif
766 
767  }// for oth
768 
769  }// for bou
770 
771  getTimestat().timeFillMatrixAndRhs.second += omp_get_wtime();
772 }//FillIQ()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
timePeriod timeFillMatrixAndRhs
Начало и конец процесса заполнения матрицы и формирования правой части
Definition: Times2D.h:80
bool useInverseMatrix
Признак использования обратной матрицы
Definition: World2D.h:111
std::vector< std::vector< std::pair< Eigen::MatrixXd, Eigen::MatrixXd > > > IQ
Матрица, состоящая из пар матриц, в которых хранятся касательные и нормальные компоненты интегралов о...
Definition: World2D.h:105
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::FillMatrixAndRhs ( )

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

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

Definition at line 438 of file World2D.cpp.

439 {
440  getTimestat().timeFillMatrixAndRhs.first += omp_get_wtime();
441 
442  Eigen::MatrixXd locMatr;
443  Eigen::MatrixXd otherMatr;
444  Eigen::VectorXd locLastLine, locLastCol;
445 
446  std::vector<std::vector<Point2D>> locIQ;
447  std::vector<std::vector<Point2D>> locOtherIQ;
448 
449  //обнуляем матрицу на первом шаге расчета
450  if (currentStep == 0)
451  {
452  for (int i = 0; i < matr.rows(); ++i)
453  for (int j = 0; j < matr.cols(); ++j)
454  matr(i, j) = 0.0;
455  }
456 
457  size_t currentRow = 0;
458 
459  for (size_t bou = 0; bou < boundary.size(); ++bou)
460  {
461  size_t nVars = boundary[bou]->GetUnknownsSize();
462  if (currentStep == 0)
463  {
464  locMatr.resize(nVars, nVars);
465  locLastLine.resize(nVars);
466  locLastCol.resize(nVars);
467  }
468 
469  if (currentStep == 0 || mechanics[bou]->isDeform)
470  {
471  boundary[bou]->FillMatrixSelf(locMatr, locLastLine, locLastCol);
472  }
473 
474  //размазываем матрицу
475 
476  for (size_t i = 0; i < nVars; ++i)
477  {
478  if (currentStep == 0 || mechanics[bou]->isDeform)
479  {
480  for (size_t j = 0; j < nVars; ++j)
481  matr(i + currentRow, j + currentRow) = locMatr(i, j);
482  matr(currentRow + nVars, i + currentRow) = locLastLine(i);
483  matr(i + currentRow, currentRow + nVars) = locLastCol(i);
484  }
485  }
486 
487  if ( (currentStep == 0) || (!useInverseMatrix) )
488  {
489  size_t currentCol = 0;
490  for (size_t oth = 0; oth < boundary.size(); ++oth)
491  {
492  size_t nVarsOther = boundary[oth]->GetUnknownsSize();
493 
494  if (bou != oth)
495  {
496  otherMatr.resize(nVars, nVarsOther);
497 
498  boundary[bou]->FillMatrixFromOther(*boundary[oth], otherMatr);
499 
500  //размазываем матрицу
501  for (size_t i = 0; i < nVars; ++i)
502  {
503  for (size_t j = 0; j < nVarsOther; ++j)
504  matr(i + currentRow, j + currentCol) = otherMatr(i, j);
505  }
506  }// if (bou != oth)
507  currentCol += nVarsOther + 1;
508  }// for oth
509  }// if (currentStep == 0 || mechanics[oth]->isMoves)
510 
511  currentRow += nVars + 1;
512  }// for bou
513 
514  velocity->FillRhs(rhs);
515 
516  getTimestat().timeFillMatrixAndRhs.second += omp_get_wtime();
517 }//FillMatrixAndRhs()
std::unique_ptr< Velocity > velocity
Умный укзатель на объект, определяющий методику вычисления скоростей
Definition: World2D.h:87
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
Eigen::VectorXd rhs
Правая часть системы
Definition: World2D.h:114
Eigen::MatrixXd matr
Матрица системы
Definition: World2D.h:102
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
timePeriod timeFillMatrixAndRhs
Начало и конец процесса заполнения матрицы и формирования правой части
Definition: Times2D.h:80
bool useInverseMatrix
Признак использования обратной матрицы
Definition: World2D.h:111
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::GenerateMechanicsHeader ( size_t  mechanicsNumber)

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

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

Definition at line 725 of file World2D.cpp.

726 {
727  mechanics[mechanicsNumber]->GenerateForcesHeader();
728  mechanics[mechanicsNumber]->GeneratePositionHeader();
729 }//GenerateMechanicsHeader(...)
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84

Here is the caller graph for this function:

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

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

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

Definition at line 130 of file World2D.h.

130 { return *airfoil[i]; };
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72

Here is the caller graph for this function:

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

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

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

Definition at line 153 of file World2D.h.

153 { return *boundary[i]; };
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78

Here is the caller graph for this function:

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

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

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

Definition at line 227 of file World2D.h.

227 { return cuda; };
Gpu cuda
Объект, управляющий графическим ускорителем
Definition: World2D.h:123

Here is the caller graph for this function:

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

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

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

Definition at line 91 of file WorldGen.h.

91 { return currentStep; };
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69

Here is the caller graph for this function:

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

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

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

Definition at line 170 of file World2D.h.

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

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

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

Definition at line 74 of file WorldGen.h.

74 { return info; };
LogStream info
Поток для вывода логов и сообщений об ошибках
Definition: WorldGen.h:59

Here is the caller graph for this function:

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

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

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

Definition at line 81 of file WorldGen.h.

81 { return info(x); };
LogStream info
Поток для вывода логов и сообщений об ошибках
Definition: WorldGen.h:59
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 237 of file World2D.h.

237 { return IQ[i][j]; };
std::vector< std::vector< std::pair< Eigen::MatrixXd, Eigen::MatrixXd > > > IQ
Матрица, состоящая из пар матриц, в которых хранятся касательные и нормальные компоненты интегралов о...
Definition: World2D.h:105

Here is the caller graph for this function:

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

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

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

Definition at line 175 of file World2D.h.

175 { return *measureVP; };
std::unique_ptr< MeasureVP > measureVP
Умный указатель на алгоритм вычисления полей скоростей и давления (для сохранения в файл) ...
Definition: World2D.h:96

Here is the caller graph for this function:

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

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

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

Definition at line 186 of file World2D.h.

186 { return *mechanics[i]; };
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84

Here is the caller graph for this function:

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

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

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

Definition at line 142 of file World2D.h.

142 { return *airfoil[i]; };
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72

Here is the caller graph for this function:

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

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

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

Definition at line 159 of file World2D.h.

159 { return *boundary[i]; };
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78

Here is the caller graph for this function:

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

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

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

Definition at line 232 of file World2D.h.

232 { return cuda; };
Gpu cuda
Объект, управляющий графическим ускорителем
Definition: World2D.h:123

Here is the caller graph for this function:

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

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

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

Definition at line 180 of file World2D.h.

180 { return *measureVP; };
std::unique_ptr< MeasureVP > measureVP
Умный указатель на алгоритм вычисления полей скоростей и давления (для сохранения в файл) ...
Definition: World2D.h:96

Here is the caller graph for this function:

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

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

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

Definition at line 192 of file World2D.h.

192 { return *mechanics[i]; };
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84

Here is the caller graph for this function:

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

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

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

Definition at line 217 of file World2D.h.

217 { return *velocity; };
std::unique_ptr< Velocity > velocity
Умный укзатель на объект, определяющий методику вычисления скоростей
Definition: World2D.h:87

Here is the caller graph for this function:

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

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

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

Definition at line 202 of file World2D.h.

202 { return *wake; };
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition: World2D.h:90

Here is the caller graph for this function:

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

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

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

Definition at line 147 of file World2D.h.

147 { return airfoil.size(); };
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72

Here is the caller graph for this function:

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

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

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

Definition at line 164 of file World2D.h.

164 { return boundary.size(); };
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78

Here is the caller graph for this function:

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

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

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

Definition at line 136 of file World2D.h.

136 { return *oldAirfoil[i]; };
std::vector< std::unique_ptr< Airfoil > > oldAirfoil
Список умных указателей на обтекаемые профили для сохранения старого положения
Definition: World2D.h:75

Here is the caller graph for this function:

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

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

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

Definition at line 222 of file World2D.h.

222 { return passport; };
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition: World2D.h:120

Here is the caller graph for this function:

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

Definition at line 93 of file WorldGen.h.

93 { return passportGen; };
const PassportGen & passportGen
Константная ссылка на паспорт конкретного расчета
Definition: WorldGen.h:62

Here is the call graph for this function:

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

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

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

Definition at line 207 of file World2D.h.

207 { return *source; };
std::unique_ptr< WakeDataBase > source
Умный указатель на источники
Definition: World2D.h:93

Here is the caller graph for this function:

Times& VM2D::World2D::getTimestat ( ) const
inline

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

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

Definition at line 242 of file World2D.h.

242 { return dynamic_cast<Times&>(*timestat); };

Here is the caller graph for this function:

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

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

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

Definition at line 212 of file World2D.h.

212 { return *velocity; };
std::unique_ptr< Velocity > velocity
Умный укзатель на объект, определяющий методику вычисления скоростей
Definition: World2D.h:87

Here is the caller graph for this function:

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

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

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

Definition at line 197 of file World2D.h.

197 { return *wake; };
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition: World2D.h:90

Here is the caller graph for this function:

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

Definition at line 244 of file World2D.h.

244 { return ((val > 0) && (!(currentStep % val))); };
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69

Here is the call graph for this function:

Here is the caller graph for this function:

bool WorldGen::isFinished ( ) const
inherited

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

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

Definition at line 49 of file WorldGen.cpp.

50 {
52 }
double timeStop
Конечное время
Definition: PassportGen.h:62
double currTime
Текущее время
Definition: PassportGen.h:56
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
const PassportGen & passportGen
Константная ссылка на паспорт конкретного расчета
Definition: WorldGen.h:62

Here is the caller graph for this function:

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

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

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

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

Definition at line 681 of file World2D.cpp.

682 {
683  getTimestat().timeMoveVortexes.first += omp_get_wtime();
684 
685  size_t nvt = wake->vtx.size();
686  for (size_t i = 0; i < getNumberOfBoundary(); ++i)
687  nvt += boundary[i]->virtualWake.vtx.size();
688 
689  //newPos.clear();
690  newPos.resize(nvt);
691 
692 
693 #pragma omp parallel for
694  for (int i = 0; i < (int)wake->vtx.size(); ++i)
695  {
696  newPos[i] = wake->vtx[i].r() \
697  + (velocity->wakeVortexesParams.convVelo[i] +
698  velocity->wakeVortexesParams.diffVelo[i] +
700  }
701 
702  int curCounterNewPos = (int)wake->vtx.size();
703 
704  wake->vtx.resize(nvt);
705 
706 
707  for (size_t bou = 0; bou < boundary.size(); ++bou)
708  {
709 #pragma omp parallel for
710  for (int i = 0; i < (int)boundary[bou]->virtualWake.vtx.size(); ++i)
711  {
712  wake->vtx[curCounterNewPos + i] = (boundary[bou]->virtualWake.vtx[i]);
713  newPos[curCounterNewPos + i] = (boundary[bou]->virtualWake.vtx[i].r() \
714  + (velocity->virtualVortexesParams[bou].convVelo[i] +
715  velocity->virtualVortexesParams[bou].diffVelo[i] +
717  }
718  curCounterNewPos += (int)boundary[bou]->virtualWake.vtx.size();
719  }
720 
721  getTimestat().timeMoveVortexes.second += omp_get_wtime();
722 }//MoveVortexes(...)
std::unique_ptr< Velocity > velocity
Умный укзатель на объект, определяющий методику вычисления скоростей
Definition: World2D.h:87
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition: World2D.h:90
double dt
Шаг по времени
Definition: PassportGen.h:65
timePeriod timeMoveVortexes
Начало и конец процесса перемещения вихрей
Definition: Times2D.h:95
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
Point2D V0() const
Функция скорости набегающего потока с учетом разгона
Definition: Passport2D.h:93
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition: World2D.h:120

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::ReserveMemoryForMatrixAndRhs ( )

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

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

Definition at line 520 of file World2D.cpp.

521 {
522  getTimestat().timeReserveMemoryForMatrixAndRhs.first += omp_get_wtime();
523 
524 
525  if (currentStep == 0)
526  {
527  dispBoundaryInSystem.resize(boundary.size());
528  dispBoundaryInSystem[0] = 0;
529 
530  for (size_t i = 1; i < boundary.size(); ++i)
531  {
532  dispBoundaryInSystem[i] = dispBoundaryInSystem[i - 1] + boundary[i - 1]->GetUnknownsSize() + 1;
533  }
534 
535  size_t matrSize = boundary.size();
536  for (auto it = boundary.begin(); it != boundary.end(); ++it)
537  matrSize += (*it)->GetUnknownsSize();
538 
539  matr.resize(matrSize, matrSize);
540  matr.setZero();
541 
542  for (size_t i = 0; i < getNumberOfAirfoil(); ++i)
543  {
544  size_t nVari = boundary[i]->GetUnknownsSize();
545  for (size_t j = 0; j < getNumberOfAirfoil(); ++j)
546  {
547  size_t nVarj = boundary[j]->GetUnknownsSize();
548  IQ[i][j].first.resize(nVari, nVarj);
549  IQ[i][j].second.resize(nVari, nVarj);
550  }
551  }
552 
553  rhs.resize(matrSize);
554  }
555  rhs.setZero();
556 
557  getTimestat().timeReserveMemoryForMatrixAndRhs.second += omp_get_wtime();
558 }//ReserveMemoryForMatrixAndRhs()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
Eigen::VectorXd rhs
Правая часть системы
Definition: World2D.h:114
Eigen::MatrixXd matr
Матрица системы
Definition: World2D.h:102
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
std::vector< std::vector< std::pair< Eigen::MatrixXd, Eigen::MatrixXd > > > IQ
Матрица, состоящая из пар матриц, в которых хранятся касательные и нормальные компоненты интегралов о...
Definition: World2D.h:105
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
timePeriod timeReserveMemoryForMatrixAndRhs
Начало и конец процесса выделения памяти под матрицу и правую часть
Definition: Times2D.h:77
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78
std::vector< size_t > dispBoundaryInSystem
Список номеров, с которых начинаются элементы правой части (или матрицы) системы для профилей ...
Definition: World2D.h:81

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::SolveLinearSystem ( )

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

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

Definition at line 333 of file World2D.cpp.

334 {
335  getTimestat().timeSolveLinearSystem.first += omp_get_wtime();
336 
337 /*
338  if (currentStep == 0)
339  {
340  invMatr = matr;
341  std::ifstream fileMatrix("matrix.txt");
342  int nx, ny;
343  fileMatrix >> nx;
344  fileMatrix >> ny;
345  for (int i = 0; i < matr.rows(); ++i)
346  {
347  for (int j = 0; j < matr.cols(); ++j)
348  {
349  fileMatrix >> matr(i, j);
350  }
351  }
352  fileMatrix.close();
353  }
354 */
355 
356 
357  if (useInverseMatrix && (currentStep == 0))
358  {
359  info('t') << "Inverting matrix" << std::endl;
360 
361 #if (defined(USE_CUDA))
362  invMatr.resize(matr.rows(), matr.cols());
363  for (int i = 0; i < (int)matr.rows(); ++i)
364  for (int j = 0; j < (int)matr.cols(); ++j)
365  invMatr(i, j) = (i == j) ? 1.0 : 0.0;
366  cuInverseMatrix((int)matr.rows(), matr.data(), invMatr.data());
367 #else
368  invMatr = matr.inverse();
369 #endif
370  }
371 
372 
373 
374 /*
375  if (currentStep == 0)
376  {
377  invMatr = matr;
378  std::ifstream fileMatrix("invMatrix.txt");
379  int nx, ny;
380  fileMatrix >> nx;
381  fileMatrix >> ny;
382  for (int i = 0; i < invMatr.rows(); ++i)
383  {
384  for (int j = 0; j < invMatr.cols(); ++j)
385  {
386  fileMatrix >> invMatr(i, j);
387  }
388  }
389  fileMatrix.close();
390  }
391 */
392 
393 /*
394  if (currentStep == 0)
395  {
396  Eigen::MatrixXd mmul = matr * invMatr;
397  for (int i = 0; i < invMatr.rows(); ++i)
398  mmul(i,i) -= 1.0;
399 
400  double maxElem = 0.0;
401  for (int i = 0; i < mmul.rows(); ++i)
402  for (int j = 0; j < mmul.cols(); ++j)
403  if (fabs(mmul(i,j)) > maxElem )
404  maxElem = mmul(i,j);
405 
406  std::cout << "A*A^-1: MAX ELEMENT = " << maxElem << std::endl;
407  }
408 */
409 
410 /*
411  if (currentStep == 0)
412  {
413  std::ofstream fileMatrix("invMatrix4x3200.txt");
414  fileMatrix << invMatr.rows() << " " << invMatr.cols() << std::endl;
415  fileMatrix.precision(18);
416  for (int i = 0; i < invMatr.rows(); ++i)
417  {
418  for (int j = 0; j < invMatr.cols(); ++j)
419  {
420  fileMatrix << invMatr(i, j) << " ";
421  }
422  fileMatrix << std::endl;
423  }
424  fileMatrix.close();
425  }
426 */
427 
428 
429  if (useInverseMatrix)
430  sol = invMatr * rhs;
431  else
432  sol = matr.partialPivLu().solve(rhs);
433 
434  getTimestat().timeSolveLinearSystem.second += omp_get_wtime();
435 }//SolveLinearSystem()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
Eigen::MatrixXd invMatr
Обратная матрица
Definition: World2D.h:108
Eigen::VectorXd rhs
Правая часть системы
Definition: World2D.h:114
LogStream info
Поток для вывода логов и сообщений об ошибках
Definition: WorldGen.h:59
Eigen::MatrixXd matr
Матрица системы
Definition: World2D.h:102
timePeriod timeSolveLinearSystem
Начало и конец процесса решения системы линейных алгебраических уравнений
Definition: Times2D.h:83
bool useInverseMatrix
Признак использования обратной матрицы
Definition: World2D.h:111
Eigen::VectorXd sol
Решение системы
Definition: World2D.h:117
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::Step ( )
overridevirtual

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

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

Implements VMlib::WorldGen.

Definition at line 196 of file World2D.cpp.

197 {
198  try {
199  //Очистка статистики
200  getTimestat().ToZero();
201 
202  //Засечка времени в начале шага
203  getTimestat().timeWholeStep.first += omp_get_wtime();
204 
206  measureVP->Initialization();
207 
208  //for (int i = 0; i < wake->vtx.size(); ++i)
209  //{
210  // if (std::isnan(wake->vtx[i].r()[0]) || std::isnan(wake->vtx[i].r()[1]) || std::isnan(wake->vtx[i].g()))
211  // {
212  // std::cout << wake->vtx[i].r() << " " << wake->vtx[i].g() << std::endl;
213  // exit(-4);
214  // }
215  //}
216 
217 
219 
220  //Вычисление скоростей вихрей: для тех, которые в следе, и виртуальных, а также в точках wakeVP
221  CalcVortexVelo(false);
222 
223  //Расчет и сохранение поля давления
225  {
226  measureVP->CalcPressure();
227  measureVP->SaveVP();
228  }
229 
230  std::vector<MechanicsRigidOscillPart*> mechOscilPart;
231  if (mechanics.size() > 0)
232  {
233  mechOscilPart.resize(mechanics.size(), nullptr);
234 
235  for (size_t s = 0; s < mechanics.size(); ++s)
236  mechOscilPart[s] = dynamic_cast<MechanicsRigidOscillPart*>(mechanics[0].get());
237 
238  if (getPassport().airfoilParams[0].addedMass.length2() > 0)
239  {
240  for (size_t s = 0; s < mechanics.size(); ++s)
241  mechanics[s]->hydroDynamForce = { 0.0, 0.0 }; //12-01
242 
243  WakeAndAirfoilsMotion(); //12-01
245  CalcAndSolveLinearSystem(); //12-01
246  }
247 
248  }
249  //Вычисление сил, действующих на профиль и сохранение в файл
250 
251  for (auto& mech : mechanics)
252  {
253  mech->GetHydroDynamForce();
254  mech->GenerateForcesString();
255  mech->GeneratePositionString();
256  }
257 
258  for (size_t s = 0; s < mechanics.size(); ++s)
259  if (mechOscilPart[s])
260  {
261  if (getPassport().airfoilParams[0].addedMass.length2() > 0)
262  mechOscilPart[s]->UpdateU(); //12-01
263  }
264 
265  //Движение вихрей (сброс вихрей + передвижение пелены)
266  WakeAndAirfoilsMotion(); //12-01
267 
268  wake->Restruct();
269 
270  wake->SaveKadrVtk();
271 
272  //Засечка времени в конце шага
273  getTimestat().timeWholeStep.second += omp_get_wtime();
274 
275 
276  info('i') << "Step = " << currentStep \
277  << " PhysTime = " << passport.physicalProperties.getCurrTime() \
278  << " StepTime = " << Times::dT(getTimestat().timeWholeStep) << std::endl;
280 
281 
283  ++currentStep;
284  }
285  catch (...)
286  {
287  info('e') << "!!! Exception from unknown source !!!" << std::endl;
288  exit(-1);
289  }
290 }//Step()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
virtual void ToZero() override
Обнуление состояния временной статистики
Definition: Times2D.cpp:100
virtual void GenerateStatString() const override
Сохранение строки со статистикой в файл временной статистики
Definition: Times2D.cpp:72
double getCurrTime() const
Возвращает текуще время
Definition: Passport2D.h:102
LogStream info
Поток для вывода логов и сообщений об ошибках
Definition: WorldGen.h:59
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition: World2D.h:90
double dt
Шаг по времени
Definition: PassportGen.h:65
bool ifDivisible(int val) const
Definition: World2D.h:244
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
std::unique_ptr< MeasureVP > measureVP
Умный указатель на алгоритм вычисления полей скоростей и давления (для сохранения в файл) ...
Definition: World2D.h:96
static double dT(const timePeriod &t)
Definition: TimesGen.h:82
int saveVPstep
Шаг вычисления и сохранения скорости и давления
Definition: PassportGen.h:78
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
void CalcAndSolveLinearSystem()
Набор матрицы, правой части и решение СЛАУ
Definition: World2D.cpp:774
size_t currentStep
Текущий номер шага в решаемой задаче
Definition: WorldGen.h:69
void CalcPanelsVeloAndAttachedSheets()
Вычисление скоростей панелей и интенсивностей присоединенных слоев вихрей и источников ...
Definition: World2D.cpp:664
void addCurrTime(double dt_) const
Добавление шага к текущему времени
Definition: Passport2D.h:114
std::vector< AirfoilParams > airfoilParams
Список структур с параметрами профилей
Definition: Passport2D.h:280
void CalcVortexVelo(bool shiftTime)
Вычисление скоростей (и конвективных, и диффузионных) вихрей (в пелене и виртуальных), а также в точках вычисления VP.
Definition: World2D.cpp:563
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition: World2D.h:120
timePeriod timeWholeStep
Начало и конец процесса выполнения шага целиком
Definition: Times2D.h:71
void WakeAndAirfoilsMotion()
Перемещение вихрей и профилей на шаге
Definition: World2D.cpp:873

Here is the call graph for this function:

Here is the caller graph for this function:

void World2D::WakeAndAirfoilsMotion ( )

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

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

Definition at line 873 of file World2D.cpp.

874 {
875  std::vector<Point2D> newPos;
876 
877  MoveVortexes(newPos);
878 
879 #ifdef BRIDGE
880  double totalForce = 0;
881  for (size_t afl = 0; afl < airfoil.size(); ++afl)
882  {
883  totalForce += mechanics[afl]->hydroDynamForce[1];
884  }
885  totalForce *= 1.2;
886 
887  mechanics[0]->hydroDynamForce[1] = totalForce;
888 #endif
889 
890  getTimestat().timeOther.first += omp_get_wtime();
891 
892  oldAirfoil.resize(0);
893  for (auto& afl : airfoil)
894  {
895  switch (passport.numericalSchemes.panelsType.second)
896  {
897  case 0:
898  oldAirfoil.emplace_back(new AirfoilRect(*afl));
899  break;
900  //case 1:
901  // oldAirfoil.emplace_back(new AirfoilCurv(*afl));
902  // break;
903  }
904 
905 
906 #ifdef BRIDGE
907  if (afl->numberInPassport == 0)
908  {
909  mechanics[afl->numberInPassport]->Move();
910  }
911  else
912  {
913  Mechanics& mechGen = *mechanics[0];
914  MechanicsRigidOscillPart& mech0 = dynamic_cast<MechanicsRigidOscillPart&>(mechGen);
915  double dy = mech0.getY() - mech0.getYOld();
916  double du = mech0.getU() - mech0.getUOld();
917 
918  Mechanics& mechGenI = *mechanics[afl->numberInPassport];
919  MechanicsRigidOscillPart& mechI = dynamic_cast<MechanicsRigidOscillPart&>(mechGenI);
920 
921  mechI.getUOld() = mechI.getU();
922  mechI.getYOld() = mechI.getY();
923 
924  //std::cout << "afl = " << afl << ", dy = " << dy << std::endl;
925 
926  airfoil[afl->numberInPassport]->Move({ 0.0, dy });
927  mechI.Vcm[1] += du;
928 
929  mechI.getY() += dy;
930  mechI.getU() += du;
931  }
932 #endif
933 
934 #ifdef INITIAL
935  mechanics[afl->numberInPassport]->Move();
936 #endif
937  }//for
938 
939  for (auto& bou : boundary)
940  bou->virtualWake.vtx.clear();
941 
942  getTimestat().timeOther.second += omp_get_wtime();
943 
944  CheckInside(newPos, oldAirfoil);
945 
946  getTimestat().timeOther.first += omp_get_wtime();
947 
948  //передача новых положений вихрей в пелену
949  for (size_t i = 0; i < wake->vtx.size(); ++i)
950  wake->vtx[i].r() = newPos[i];
951 
952 // getWake().SaveKadrVtk();
953  getTimestat().timeOther.second += omp_get_wtime();
954 }//WakeAndAirfoilsMotion()
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
void MoveVortexes(std::vector< Point2D > &newPos)
Вычисляем новые положения вихрей (в пелене и виртуальных)
Definition: World2D.cpp:681
Класс, определяющий вид механической системы
std::vector< std::unique_ptr< Airfoil > > oldAirfoil
Список умных указателей на обтекаемые профили для сохранения старого положения
Definition: World2D.h:75
std::unique_ptr< Wake > wake
Умный указатель на вихревой след
Definition: World2D.h:90
std::vector< std::unique_ptr< Mechanics > > mechanics
Список умных указателей на типы механической системы для каждого профиля
Definition: World2D.h:84
Класс, определяющий тип обтекаемого профиля
Definition: Airfoil2DRect.h:64
void CheckInside(std::vector< Point2D > &newPos, const std::vector< std::unique_ptr< Airfoil >> &oldAirfoil)
Проверка проникновения вихрей внутрь профиля
Definition: World2D.cpp:294
NumericalSchemes numericalSchemes
Структура с используемыми численными схемами
Definition: Passport2D.h:302
double & getY()
текущее отклонение профиля
std::vector< std::unique_ptr< Airfoil > > airfoil
Список умных указателей на обтекаемые профили
Definition: World2D.h:72
std::pair< std::string, int > panelsType
Тип панелей (0 — прямые)
Definition: Passport2D.h:196
std::vector< std::unique_ptr< Boundary > > boundary
Список умных указателей на формирователи граничных условий на профилях
Definition: World2D.h:78
const Passport & passport
Константная ссылка на паспорт конкретного расчета
Definition: World2D.h:120
timePeriod timeOther
Все прочее
Definition: Times2D.h:113
Абстрактный класс, определяющий вид механической системы
Definition: Mechanics2D.h:71

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

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

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

Definition at line 72 of file World2D.h.

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

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

Definition at line 78 of file World2D.h.

Gpu VM2D::World2D::cuda
mutableprivate

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

Definition at line 123 of file World2D.h.

size_t VMlib::WorldGen::currentStep
inherited

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

Definition at line 69 of file WorldGen.h.

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

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

Definition at line 81 of file World2D.h.

LogStream VMlib::WorldGen::info
mutableprotectedinherited

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

Definition at line 59 of file WorldGen.h.

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

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

Definition at line 108 of file World2D.h.

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

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

Definition at line 105 of file World2D.h.

Eigen::MatrixXd VM2D::World2D::matr
private

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

Definition at line 102 of file World2D.h.

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

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

Definition at line 96 of file World2D.h.

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

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

Definition at line 84 of file World2D.h.

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

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

Definition at line 75 of file World2D.h.

const Passport& VM2D::World2D::passport
private

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

Definition at line 120 of file World2D.h.

const PassportGen& VMlib::WorldGen::passportGen
protectedinherited

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

Definition at line 62 of file WorldGen.h.

Eigen::VectorXd VM2D::World2D::rhs
private

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

Definition at line 114 of file World2D.h.

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

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

Definition at line 117 of file World2D.h.

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

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

Definition at line 93 of file World2D.h.

std::unique_ptr<TimesGen> VMlib::WorldGen::timestat
protectedinherited

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

Definition at line 65 of file WorldGen.h.

bool VM2D::World2D::useInverseMatrix
private

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

Definition at line 111 of file World2D.h.

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

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

Definition at line 87 of file World2D.h.

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

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

Definition at line 90 of file World2D.h.


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