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

Класс, опеделяющий вихревой след (пелену) More...

#include <VirtualWake2D.h>

Inheritance diagram for VM2D::VirtualWake:
Collaboration diagram for VM2D::VirtualWake:

Public Member Functions

 VirtualWake (const World2D &W_, const Boundary &bnd_)
 Конструктор инициализации More...
 
 ~VirtualWake ()
 Деструктор More...
 
void ReadFromFile (const std::string &dir, const std::string &fileName)
 Считывание вихревого следа из файла More...
 
void SaveKadrVtk (const std::string &filePrefix="Kadr") const
 Сохранение вихревого следа в файл .vtk. More...
 

Public Attributes

const Boundarybnd
 Константная ссылка на границу, которой принадлежит данный виртуальный след More...
 
std::vector< Point2DvecHalfGamma
 Скорость вихрей виртуального следа конкретного профиля (равна Gamma/2) используется для расчета давления More...
 
std::vector< std::pair< size_t, size_t > > aflPan
 Пара чисел: номер профиля и номер панели, на которой рожден виртуальный вихрь More...
 
const World2DW
 Константная ссылка на решаемую задачу More...
 
std::vector< Vortex2Dvtx
 Список вихревых элементов More...
 

Detailed Description

Класс, опеделяющий вихревой след (пелену)

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

Definition at line 60 of file VirtualWake2D.h.

Constructor & Destructor Documentation

VM2D::VirtualWake::VirtualWake ( const World2D W_,
const Boundary bnd_ 
)
inline

Конструктор инициализации

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

Definition at line 79 of file VirtualWake2D.h.

80  : WakeDataBase(W_), bnd(bnd_) { };
const Boundary & bnd
Константная ссылка на границу, которой принадлежит данный виртуальный след
Definition: VirtualWake2D.h:67
WakeDataBase(const World2D &W_)
Конструктор инициализации
VM2D::VirtualWake::~VirtualWake ( )
inline

Деструктор

Definition at line 83 of file VirtualWake2D.h.

83 { };

Member Function Documentation

void WakeDataBase::ReadFromFile ( const std::string &  dir,
const std::string &  fileName 
)
inherited

Считывание вихревого следа из файла

Parameters
[in]dirконстантная ссылка на строку, задающую каталог, где лежит файл с вихревым следом
[in]fileNameконстантная ссылка на строку, задающую имя файла с вихревым следом

Definition at line 56 of file WakeDataBase2D.cpp.

57 {
58  std::string filename = dir + fileName;
59  std::ifstream wakeFile;
60 
61  if (fileExistTest(filename, W.getInfo(), { "txt", "TXT" }))
62  {
63  std::stringstream wakeFile(VMlib::Preprocessor(filename).resultString);
64 
65  VMlib::LogStream XXX;
66  VMlib::StreamParser wakeParser(XXX, "vortex wake file parser", wakeFile);
67 
68  wakeParser.get("vtx", vtx);
69  }
70 }//ReadFromFile(...)
Класс, определяющий работу с потоком логов
Definition: LogStream.h:53
Класс, позволяющий выполнять предварительную обработку файлов
Definition: Preprocessor.h:59
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
bool fileExistTest(std::string &fileName, LogStream &info, const std::list< std::string > &extList={})
Проверка существования файла
Definition: defs.h:324
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода ...
Definition: WorldGen.h:74
Класс, позволяющий выполнять разбор файлов и строк с настройками и параметрами
Definition: StreamParser.h:151

Here is the call graph for this function:

Here is the caller graph for this function:

void WakeDataBase::SaveKadrVtk ( const std::string &  filePrefix = "Kadr") const
inherited

Сохранение вихревого следа в файл .vtk.

Definition at line 73 of file WakeDataBase2D.cpp.

74 {
75  W.getTimestat().timeSaveKadr.first += omp_get_wtime();
76 
78  {
79  std::ofstream outfile;
80  size_t numberNonZero = 0;
81 
82  if (vtx.size() > 0)
83  numberNonZero += vtx.size();
84  else
85  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
86  numberNonZero += W.getAirfoil(q).getNumberOfPanels();
87  VMlib::CreateDirectory(W.getPassport().dir, "snapshots");
88 
89  if (W.getPassport().timeDiscretizationProperties.fileTypeVtx.second == 0) //text format vtk
90  {
91  std::string fname = VMlib::fileNameStep(filePrefix, W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "vtk");
92  outfile.open(W.getPassport().dir + "snapshots/" + fname);
93 
94  outfile << "# vtk DataFile Version 2.0" << std::endl;
95  outfile << "VM2D VTK result: " << (W.getPassport().dir + "snapshots/" + fname).c_str() << " saved " << VMlib::CurrentDataTime() << std::endl;
96  outfile << "ASCII" << std::endl;
97  outfile << "DATASET UNSTRUCTURED_GRID" << std::endl;
98  outfile << "POINTS " << numberNonZero << " float" << std::endl;
99 
100  if (vtx.size() > 0)
101  for (auto& v : vtx)
102  {
103  const Point2D& r = v.r();
104  outfile << r[0] << " " << r[1] << " " << "0.0" << std::endl;
105  }//for v
106  else
107  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
108  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
109  {
110  const Point2D& r = W.getAirfoil(q).getR(s);
111  outfile << r[0] << " " << r[1] << " " << "0.0" << std::endl;
112  }
113 
114  outfile << "CELLS " << numberNonZero << " " << 2 * numberNonZero << std::endl;
115  for (size_t i = 0; i < numberNonZero; ++i)
116  outfile << "1 " << i << std::endl;
117 
118  outfile << "CELL_TYPES " << numberNonZero << std::endl;
119  for (size_t i = 0; i < numberNonZero; ++i)
120  outfile << "1" << std::endl;
121 
122  outfile << std::endl;
123  outfile << "POINT_DATA " << numberNonZero << std::endl;
124  outfile << "SCALARS Gamma float 1" << std::endl;
125  outfile << "LOOKUP_TABLE default" << std::endl;
126 
127  if (vtx.size() > 0)
128  for (auto& v : vtx)
129  outfile << v.g() << std::endl;
130  else
131  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
132  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
133  {
134  //const Point2D& r = W.getAirfoil(q).getR(s);
135  outfile << "0.0" << std::endl;
136  }
137 
138  outfile.close();
139  }//if fileType = text
140  else if (W.getPassport().timeDiscretizationProperties.fileTypeVtx.second == 1) //binary format vtk
141  {
142  //Тест способа хранения чисел
143  uint16_t x = 0x0001;
144  bool littleEndian = (*((uint8_t*)&x));
145  const char eolnBIN[] = "\n";
146 
147  std::string fname = VMlib::fileNameStep(filePrefix, W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "vtk");
148  outfile.open(W.getPassport().dir + "snapshots/" + fname, std::ios::out | std::ios::binary);
149 
150  outfile << "# vtk DataFile Version 3.0" << "\r\n" << "VM2D VTK result: " << (W.getPassport().dir + "snapshots/" + fname).c_str() << " saved " << VMlib::CurrentDataTime() << eolnBIN;
151  outfile << "BINARY" << eolnBIN;
152  outfile << "DATASET UNSTRUCTURED_GRID" << eolnBIN << "POINTS " << numberNonZero << " " << "float" << eolnBIN;
153 
154 
155  if (vtx.size() > 0)
156  {
157  Eigen::VectorXf rData = Eigen::VectorXf::Zero(vtx.size() * 3);
158  for (size_t i = 0; i < vtx.size(); ++i)
159  {
160  rData(3 * i) = (float)(vtx[i].r())[0];
161  rData(3 * i + 1) = (float)(vtx[i].r())[1];
162  }//for i
163 
164  if (littleEndian)
165  for (int i = 0; i < vtx.size() * 3; ++i)
166  VMlib::SwapEnd(rData(i));
167  outfile.write(reinterpret_cast<char*>(rData.data()), vtx.size() * 3 * sizeof(float));
168  }
169  else
170  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
171  {
172  Eigen::VectorXf rData = Eigen::VectorXf::Zero(W.getAirfoil(q).getNumberOfPanels() * 3);
173  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
174  {
175  rData(3 * s) = (float)(W.getAirfoil(q).getR(s))[0];
176  rData(3 * s + 1) = (float)(W.getAirfoil(q).getR(s))[1];
177  }//for i
178 
179  if (littleEndian)
180  for (int i = 0; i < W.getAirfoil(q).getNumberOfPanels() * 3; ++i)
181  VMlib::SwapEnd(rData(i));
182  outfile.write(reinterpret_cast<char*>(rData.data()), W.getAirfoil(q).getNumberOfPanels() * 3 * sizeof(float));
183  }
184 
185  // CELLS
186  std::vector<int> cells(2 * numberNonZero);
187  for (size_t i = 0; i < numberNonZero; ++i)
188  {
189  cells[2 * i] = 1;
190  cells[2 * i + 1] = (int)i;
191  }
192 
193  std::vector<int> cellsTypes;
194  cellsTypes.resize(numberNonZero, 1);
195 
196  if (littleEndian)
197  {
198  for (int i = 0; i < numberNonZero * 2; ++i)
199  VMlib::SwapEnd(cells[i]);
200 
201  for (int i = 0; i < numberNonZero; ++i)
202  VMlib::SwapEnd(cellsTypes[i]);
203  }
204 
205  outfile << eolnBIN << "CELLS " << numberNonZero << " " << numberNonZero * 2 << eolnBIN;
206  outfile.write(reinterpret_cast<char*>(cells.data()), numberNonZero * 2 * sizeof(int));
207  outfile << eolnBIN << "CELL_TYPES " << numberNonZero << eolnBIN;
208  outfile.write(reinterpret_cast<char*>(cellsTypes.data()), numberNonZero * sizeof(int));
209 
210  //gammas
211  outfile << eolnBIN << "POINT_DATA " << numberNonZero << eolnBIN;
212  outfile << eolnBIN << "SCALARS Gamma " << "float" << " 1" << eolnBIN;
213  outfile << "LOOKUP_TABLE default" << eolnBIN;
214 
215  if (vtx.size() > 0)
216  {
217  Eigen::VectorXf pData = Eigen::VectorXf::Zero(numberNonZero);
218  for (int s = 0; s < vtx.size(); ++s)
219  pData(s) = (float)vtx[s].g();
220 
221  if (littleEndian)
222  for (int i = 0; i < vtx.size(); ++i)
223  VMlib::SwapEnd(pData(i));
224  outfile.write(reinterpret_cast<char*>(pData.data()), vtx.size() * sizeof(float));
225  }
226  else
227  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
228  {
229  Eigen::VectorXf pData = Eigen::VectorXf::Zero(W.getAirfoil(q).getNumberOfPanels());
230  for (int s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
231  pData(s) = 0;
232 
233  if (littleEndian)
234  for (int i = 0; i < W.getAirfoil(q).getNumberOfPanels(); ++i)
235  VMlib::SwapEnd(pData(i));
236  outfile.write(reinterpret_cast<char*>(pData.data()), W.getAirfoil(q).getNumberOfPanels() * sizeof(float));
237  }
238 
239  outfile << eolnBIN;
240  outfile.close();
241  }//if binary
242  if (W.getPassport().timeDiscretizationProperties.fileTypeVtx.second == 2) //csv
243  {
244  std::string fname = VMlib::fileNameStep(filePrefix, W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "csv");
245  outfile.open(W.getPassport().dir + "snapshots/" + fname);
246 
247  outfile << "point,x,y,G " << std::endl;
248 
249  int counter = 0;
250  if (vtx.size() > 0)
251  for (auto& v : vtx)
252  outfile << counter++ << "," << v.r()[0] << "," << v.r()[1] << "," << v.g() << std::endl;
253  else
254  for (size_t q = 0; q < W.getNumberOfAirfoil(); ++q)
255  for (size_t s = 0; s < W.getAirfoil(q).getNumberOfPanels(); ++s)
256  {
257  const Point2D& r = W.getAirfoil(q).getR(s);
258  outfile << counter++ << "," << r[0] << "," << r[1] << "," << "0.0" << std::endl;
259  }
260  outfile.close();
261  }//if fileType = text
262 
263  }
264 
265  W.getTimestat().timeSaveKadr.second += omp_get_wtime();
266 }
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
int saveVtxStep
Шаг сохранения кадров в бинарные файлы
Definition: PassportGen.h:73
std::string CurrentDataTime()
Формирование строки с текущем временем и датой
Definition: defs.cpp:44
int nameLength
Число разрядов в имени файла
Definition: PassportGen.h:68
std::string dir
Рабочий каталог задачи
Definition: PassportGen.h:118
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
std::string fileNameStep(const std::string &name, int length, size_t number, const std::string &ext)
Формирование имени файла
Definition: defs.h:357
bool ifDivisible(int val) const
Definition: World2D.h:244
timePeriod timeSaveKadr
Начало и конец процесса сохранения кадра в файл
Definition: Times2D.h:107
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
size_t getNumberOfAirfoil() const
Возврат количества профилей в задаче
Definition: World2D.h:147
std::vector< Vortex2D > vtx
Список вихревых элементов
const World2D & W
Константная ссылка на решаемую задачу
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::pair< std::string, int > fileTypeVtx
Тип файлов для сохранения скорости и давления
Definition: PassportGen.h:71
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition: WorldGen.h:91
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
void SwapEnd(T &var)
Вспомогательная функция перестановки байт местами (нужно для сохранения бинарных VTK) ...
Definition: defs.h:528
void CreateDirectory(const std::string &dir, const std::string &name)
Создание каталога
Definition: defs.h:414

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

std::vector<std::pair<size_t, size_t> > VM2D::VirtualWake::aflPan

Пара чисел: номер профиля и номер панели, на которой рожден виртуальный вихрь

Definition at line 73 of file VirtualWake2D.h.

const Boundary& VM2D::VirtualWake::bnd

Константная ссылка на границу, которой принадлежит данный виртуальный след

Definition at line 67 of file VirtualWake2D.h.

std::vector<Point2D> VM2D::VirtualWake::vecHalfGamma

Скорость вихрей виртуального следа конкретного профиля (равна Gamma/2) используется для расчета давления

Definition at line 70 of file VirtualWake2D.h.

std::vector<Vortex2D> VM2D::WakeDataBase::vtx
inherited

Список вихревых элементов

Definition at line 76 of file WakeDataBase2D.h.

const World2D& VM2D::WakeDataBase::W
inherited

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

Definition at line 69 of file WakeDataBase2D.h.


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