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

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

#include <MeasureVP2D.h>

Collaboration diagram for VM2D::MeasureVP:

Public Member Functions

 MeasureVP (const World2D &W_)
 Конструктор More...
 
 ~MeasureVP ()
 Деструктор More...
 
void ReadPointsFromFile (const std::string &dir)
 Чтение точек, в которых нужно посчитать давление и скорость More...
 
void Initialization ()
 Инициализация векторов для вычисления скоростей и давлений Вызывается только на тех шагах расчета, когда это необходимо More...
 
void CalcPressure ()
 Расчет поля давления More...
 
void SaveVP ()
 Сохранение в файл вычисленных скоростей и давлений More...
 
const WakeDataBasegetWakeVP () const
 Возврат wakeVP. More...
 
const std::vector< Point2D > & getVelocity () const
 Возврат velocity. More...
 
std::vector< Point2D > & getNonConstVelocity ()
 Возврат velocity. More...
 
const std::vector< double > & getDomainRadius () const
 Возврат domainRadius. More...
 
std::vector< double > & getNonConstDomainRadius ()
 Возврат domainRadius. More...
 

Protected Attributes

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

Private Attributes

std::vector< Point2DinitialPoints
 Точки, которые считываются из файла (давление пишется в vtk-файлы) More...
 
std::vector< Point2DhistoryPoints
 Точки, которые считываются из файла (давление пишется в vtk и csv-файлы) More...
 
std::unique_ptr< WakeDataBasewakeVP
 Умный указатель на точки, в которых нужно вычислять в данный момент времени Хранятся в виде "мнимой" системы вихрей, чтобы воспользоваться методом CalcConvVeloToSetOfPoints(...) More...
 
std::vector< Point2Dvelocity
 Скорости в нужных точках More...
 
std::vector< double > domainRadius
 Радиусы вихрей в нужных точках More...
 
std::vector< double > pressure
 Давление в нужных точках More...
 

Detailed Description

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

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

Definition at line 64 of file MeasureVP2D.h.

Constructor & Destructor Documentation

MeasureVP::MeasureVP ( const World2D W_)

Конструктор

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

Definition at line 65 of file MeasureVP2D.cpp.

66  : W(W_)
67 {
68  wakeVP.reset(new WakeDataBase(W_));
69 };//MeasureVP(...)
const World2D & W
Константная ссылка на решаемую задачу
Definition: MeasureVP2D.h:91
std::unique_ptr< WakeDataBase > wakeVP
Умный указатель на точки, в которых нужно вычислять в данный момент времени Хранятся в виде "мнимой" ...
Definition: MeasureVP2D.h:75
Класс, опеделяющий набор вихрей
VM2D::MeasureVP::~MeasureVP ( )
inline

Деструктор

Definition at line 100 of file MeasureVP2D.h.

100 { };

Here is the call graph for this function:

Member Function Documentation

void MeasureVP::CalcPressure ( )

Расчет поля давления

Todo:
Пока используем только средние значения свободного слоя на панелях

Definition at line 194 of file MeasureVP2D.cpp.

195 {
196  int addWSize = (int)wakeVP->vtx.size();
197  pressure.resize(addWSize);
198 
199  Point2D dri;
200  Point2D Vi;
201  Point2D vi;
202 
203  const Point2D& V0 = W.getPassport().physicalProperties.V0();
204  const double& eps2 = W.getPassport().wakeDiscretizationProperties.eps2;
205  const double& dt = W.getPassport().timeDiscretizationProperties.dt;
206 
207  double P0 = /*1*/ 0.5 * V0.length2();
208 
209  Point2D cPan;
210 
211 #pragma warning (push)
212 #pragma warning (disable: 4101)
213  double alpha;
214  double dst2eps;
215 #pragma warning (pop)
216 
217 #pragma omp parallel for default(none) private(alpha, dri, Vi, vi, dst2eps, cPan) shared(P0, dt, eps2, V0, addWSize, std::cout, IDPI)
218  for (int i = 0; i < addWSize; ++i)
219  {
220  pressure[i] = P0;
221  pressure[i] -= 0.5 * velocity[i].length2(); //2
222 
223  const Point2D& pt = wakeVP->vtx[i].r();
224 
225  for (size_t j = 0; j < W.getWake().vtx.size(); ++j)
226  {
227  dri = pt - W.getWake().vtx[j].r();
228  Vi = W.getVelocity().wakeVortexesParams.convVelo[j] + V0;
229 
230  dst2eps = VMlib::boundDenom(dri.length2(), eps2);
231 
232  vi = IDPI * W.getWake().vtx[j].g() / dst2eps * dri.kcross();
233  pressure[i] += vi & Vi; //3
234  }
235 
236  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
237  for (size_t j = 0; j < W.getBoundary(bou).virtualWake.vtx.size(); ++j)
238  {
239  dri = pt - W.getBoundary(bou).virtualWake.vtx[j].r();
240 
241  dst2eps = VMlib::boundDenom(dri.length2(), eps2);
242 
243  Vi = W.getBoundary(bou).virtualWake.vecHalfGamma[j];
244  vi = IDPI * W.getBoundary(bou).virtualWake.vtx[j].g() / dst2eps * dri.kcross();
245 
246  pressure[i] += vi & Vi; //4
247  }
248 
249  for (size_t bou = 0; bou < W.getNumberOfBoundary(); ++bou)
250  {
251  const Point2D& rcm = W.getAirfoil(bou).rcm;
252 
253  for (size_t j = 0; j < W.getAirfoil(bou).getNumberOfPanels(); ++j)
254  {
255  cPan = 0.5 * (W.getAirfoil(bou).getR(j) + W.getAirfoil(bou).getR(j + 1));
256  alpha = atan2((cPan - pt) ^ (rcm - pt), (cPan - pt)&(rcm - pt));
257 
259  pressure[i] += IDPI * alpha * W.getBoundary(bou).sheets.freeVortexSheet(j, 0) * W.getAirfoil(bou).len[j] / dt;
260 
261  pressure[i] -= IDPI * alpha * W.getAirfoil(bou).gammaThrough[j] / dt;
262  }
263  }//for bou
264  }//for i
265 
266  for (size_t i = 0; i < pressure.size(); ++i)
268 
269 }//CalcPressure()
const double & freeVortexSheet(size_t n, size_t moment) const
Definition: Sheet2D.h:99
const World2D & W
Константная ссылка на решаемую задачу
Definition: MeasureVP2D.h:91
size_t getNumberOfBoundary() const
Возврат количества граничных условий в задаче
Definition: World2D.h:164
const Velocity & getVelocity() const
Возврат константной ссылки на объект для вычисления скоростей
Definition: World2D.h:212
const double IDPI
Число .
Definition: defs.h:76
std::vector< Point2D > vecHalfGamma
Скорость вихрей виртуального следа конкретного профиля (равна Gamma/2) используется для расчета давле...
Definition: VirtualWake2D.h:70
const Wake & getWake() const
Возврат константной ссылки на вихревой след
Definition: World2D.h:197
std::vector< double > gammaThrough
Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
Definition: Airfoil2D.h:196
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition: Airfoil2D.h:101
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition: Airfoil2D.h:139
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition: numvector.h:383
double dt
Шаг по времени
Definition: PassportGen.h:65
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
const Boundary & getBoundary(size_t i) const
Возврат константной ссылки на объект граничного условия
Definition: World2D.h:153
double rho
Плотность потока
Definition: Passport2D.h:75
std::vector< Vortex2D > vtx
Список вихревых элементов
Класс, опеделяющий двумерный вектор
Definition: Point2D.h:75
std::vector< Point2D > velocity
Скорости в нужных точках
Definition: MeasureVP2D.h:78
Sheet sheets
Слои на профиле
Definition: Boundary2D.h:95
std::vector< double > pressure
Давление в нужных точках
Definition: MeasureVP2D.h:84
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
VortexesParams wakeVortexesParams
Струтура, определяющая параметры вихрей в следе
Definition: Velocity2D.h:105
double boundDenom(double r2, double eps2)
Способ сглаживания скорости вихря (вихрь Рэнкина или вихрь Ламба)
Definition: defs.h:536
Point2D rcm
Положение центра масс профиля
Definition: Airfoil2D.h:77
std::unique_ptr< WakeDataBase > wakeVP
Умный указатель на точки, в которых нужно вычислять в данный момент времени Хранятся в виде "мнимой" ...
Definition: MeasureVP2D.h:75
VirtualWake virtualWake
Виртуальный вихревой след конкретного профиля
Definition: Boundary2D.h:85
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition: World2D.h:130
Point2D V0() const
Функция скорости набегающего потока с учетом разгона
Definition: Passport2D.h:93
WakeDiscretizationProperties wakeDiscretizationProperties
Структура с параметрами дискретизации вихревого следа
Definition: Passport2D.h:299
std::vector< double > len
Длины панелей профиля
Definition: Airfoil2D.h:185
std::vector< Point2D > convVelo
Вектор конвективных скоростей вихрей
Definition: Velocity2D.h:67
double eps2
Квадрат радиуса вихря
Definition: Passport2D.h:142
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition: numvector.h:507

Here is the call graph for this function:

Here is the caller graph for this function:

const std::vector<double>& VM2D::MeasureVP::getDomainRadius ( ) const
inline

Возврат domainRadius.

Returns
константную ссылку на начало domainRadius

Definition at line 135 of file MeasureVP2D.h.

135 { return domainRadius; };
std::vector< double > domainRadius
Радиусы вихрей в нужных точках
Definition: MeasureVP2D.h:81
std::vector<double>& VM2D::MeasureVP::getNonConstDomainRadius ( )
inline

Возврат domainRadius.

Returns
неконстантную ссылку на начало domainRadius

Definition at line 140 of file MeasureVP2D.h.

140 { return domainRadius; };
std::vector< double > domainRadius
Радиусы вихрей в нужных точках
Definition: MeasureVP2D.h:81

Here is the caller graph for this function:

std::vector<Point2D>& VM2D::MeasureVP::getNonConstVelocity ( )
inline

Возврат velocity.

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

Definition at line 130 of file MeasureVP2D.h.

130 { return velocity; };
std::vector< Point2D > velocity
Скорости в нужных точках
Definition: MeasureVP2D.h:78

Here is the caller graph for this function:

const std::vector<Point2D>& VM2D::MeasureVP::getVelocity ( ) const
inline

Возврат velocity.

Returns
константную ссылку на velocity

Definition at line 125 of file MeasureVP2D.h.

125 { return velocity; };
std::vector< Point2D > velocity
Скорости в нужных точках
Definition: MeasureVP2D.h:78
const WakeDataBase& VM2D::MeasureVP::getWakeVP ( ) const
inline

Возврат wakeVP.

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

Definition at line 120 of file MeasureVP2D.h.

120 { return *wakeVP; };
std::unique_ptr< WakeDataBase > wakeVP
Умный указатель на точки, в которых нужно вычислять в данный момент времени Хранятся в виде "мнимой" ...
Definition: MeasureVP2D.h:75

Here is the caller graph for this function:

void MeasureVP::Initialization ( )

Инициализация векторов для вычисления скоростей и давлений Вызывается только на тех шагах расчета, когда это необходимо

Definition at line 128 of file MeasureVP2D.cpp.

129 {
130  W.getTimestat().timeVP.first += omp_get_wtime();
131 
133  {
134  W.getInfo('i') << "Preparing VP points" << std::endl;
135 
136 
137  wakeVP->vtx.clear();
138  wakeVP->vtx.resize(0);
139 
140  Vortex2D addvtx;
141 
142  velocity.clear();
143  velocity.resize(0);
144 
145  pressure.clear();
146  pressure.resize(0);
147 
148  domainRadius.clear();
149  domainRadius.resize(0);
150 
151 
152  //отключаем проверку того, что точки расположены вне профиля
153  /*
154  //определяем точки, которые не находятся внутри профилей
155  bool inside;
156  for (size_t i = 0; i < initialPoints.size(); ++i)
157  {
158  inside = false;
159  for (size_t j = 0; j < W.getNumberOfAirfoil(); ++j)
160  {
161  if (W.getAirfoil(j).isInsideGabarits(initialPoints[i]) && W.getAirfoil(j).IsPointInAirfoil(initialPoints[i]))
162  {
163  inside = true;
164  break;
165  }
166  }
167  if (!inside)
168  {
169  addvtx.r() = initialPoints[i];
170  wakeVP.vtx.push_back(addvtx);
171  }
172  }
173  */
174 
175  //добавляем все точки, а не только те, которые вне профиля
176  for (size_t i = 0; i < initialPoints.size(); ++i)
177  {
178  addvtx.r() = initialPoints[i];
179  wakeVP->vtx.push_back(addvtx);
180  }//for i
181 
182  for (size_t i = 0; i < historyPoints.size(); ++i)
183  {
184  addvtx.r() = historyPoints[i];
185  wakeVP->vtx.push_back(addvtx);
186  }//for i
187 
188  }
189  W.getTimestat().timeVP.second += omp_get_wtime();
190 
191 }//Initialization()
std::vector< Point2D > initialPoints
Точки, которые считываются из файла (давление пишется в vtk-файлы)
Definition: MeasureVP2D.h:68
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
const World2D & W
Константная ссылка на решаемую задачу
Definition: MeasureVP2D.h:91
HD Point2D & r()
Функция для доступа к радиус-вектору вихря
Definition: Vortex2D.h:85
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
int saveVPstep
Шаг вычисления и сохранения скорости и давления
Definition: PassportGen.h:78
std::vector< Point2D > velocity
Скорости в нужных точках
Definition: MeasureVP2D.h:78
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода ...
Definition: WorldGen.h:74
std::vector< double > pressure
Давление в нужных точках
Definition: MeasureVP2D.h:84
std::vector< double > domainRadius
Радиусы вихрей в нужных точках
Definition: MeasureVP2D.h:81
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
Класс, опеделяющий двумерный вихревой элемент
Definition: Vortex2D.h:51
std::vector< Point2D > historyPoints
Точки, которые считываются из файла (давление пишется в vtk и csv-файлы)
Definition: MeasureVP2D.h:71
std::unique_ptr< WakeDataBase > wakeVP
Умный указатель на точки, в которых нужно вычислять в данный момент времени Хранятся в виде "мнимой" ...
Definition: MeasureVP2D.h:75
timePeriod timeVP
Начало и конец процесса подсчета полей скоростей и давления и сохранения их в файл ...
Definition: Times2D.h:110
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition: WorldGen.h:91

Here is the call graph for this function:

Here is the caller graph for this function:

void MeasureVP::ReadPointsFromFile ( const std::string &  dir)

Чтение точек, в которых нужно посчитать давление и скорость

Parameters
[in]dirконстантная ссылка на строку — имя каталога, где лежит cчитываемый файл

Definition at line 72 of file MeasureVP2D.cpp.

73 {
74  std::string filename = dir + "pointsVP";
75  std::ifstream VPFile;
76 
77  if (fileExistTest(filename, W.getInfo(), { "txt", "TXT" }))
78  {
79  std::stringstream VPFile(VMlib::Preprocessor(filename).resultString);
80 
81  VMlib::StreamParser VPParser(W.getInfo(), "velocity & pressure parser", VPFile);
82 
83  VPParser.get("points", initialPoints);
84  VPParser.get("history", historyPoints);
85 
86  for (auto& pt : initialPoints)
87  pt = pt.rotated(-W.getPassport().rotateAngleVpPoints * PI / 180.0);
88 
89  for (auto& pt : historyPoints)
90  pt = pt.rotated(-W.getPassport().rotateAngleVpPoints * PI / 180.0);
91 
92 
93  if (historyPoints.size() > 0)
94  {
95  VMlib::CreateDirectory(dir, "velPres");
96  std::string VPFileNameList = W.getPassport().dir + "velPres/listPoints";
97  std::ofstream VPFileList(VPFileNameList.c_str());
98 
99  VMlib::PrintLogoToTextFile(VPFileList, VPFileNameList.c_str(), "List of points, where velocity and pressure are measured in csv-files");
100 
101  VMlib::PrintHeaderToTextFile(VPFileList, "No. FileName X Y");
102 
103  for (size_t q = 0; q < historyPoints.size(); ++q)
104  {
105  VPFileList << '\n' << q << '\t' << VMlib::fileNameStep("VP-atPoint-", 2, q, "csv") << '\t' << historyPoints[q][0] << '\t' << historyPoints[q][1];
106 
107  std::string VPFileNameCsv;
108  VPFileNameCsv = W.getPassport().dir + "velPres/" + VMlib::fileNameStep("VP-atPoint-", 2, q, "csv");
109 
110  std::ofstream VPFileCsv(VPFileNameCsv.c_str());
111 
113  VPFileCsv << "t,Vx,Vy,p" << std::endl;
114  else
115  VPFileCsv << "t,CVx,CVy,Cp" << std::endl;
116 
117  VPFileCsv.close();
118  VPFileCsv.clear();
119  }
120 
121  VPFileList.close();
122  VPFileList.clear();
123  }
124  }
125 }//ReadPointsFromFile(...)
std::vector< Point2D > initialPoints
Точки, которые считываются из файла (давление пишется в vtk-файлы)
Definition: MeasureVP2D.h:68
const World2D & W
Константная ссылка на решаемую задачу
Definition: MeasureVP2D.h:91
const double PI
Число .
Definition: defs.h:73
std::string dir
Рабочий каталог задачи
Definition: PassportGen.h:118
Класс, позволяющий выполнять предварительную обработку файлов
Definition: Preprocessor.h:59
bool calcCoefficients
Признак вычисления коэффициентов вместо сил
Definition: Passport2D.h:289
void PrintHeaderToTextFile(std::ofstream &str, const std::string &header)
Формирование подзаголовка в текстовом файле вывода программы VM2D/VM3D.
Definition: defs.cpp:172
std::string fileNameStep(const std::string &name, int length, size_t number, const std::string &ext)
Формирование имени файла
Definition: defs.h:357
double rotateAngleVpPoints
Угол поворота точек VP.
Definition: Passport2D.h:292
bool fileExistTest(std::string &fileName, LogStream &info, const std::list< std::string > &extList={})
Проверка существования файла
Definition: defs.h:324
void PrintLogoToTextFile(std::ofstream &str, const std::string &fileName, const std::string &descr)
Формирование заголовка файла программы VM2D/VM3D.
Definition: defs.cpp:136
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода ...
Definition: WorldGen.h:74
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
std::vector< Point2D > historyPoints
Точки, которые считываются из файла (давление пишется в vtk и csv-файлы)
Definition: MeasureVP2D.h:71
void CreateDirectory(const std::string &dir, const std::string &name)
Создание каталога
Definition: defs.h:414
Класс, позволяющий выполнять разбор файлов и строк с настройками и параметрами
Definition: StreamParser.h:151
bool get(const std::string &name, std::vector< Point2D > &res, const std::vector< Point2D > *defValue=nullptr, bool echoDefault=true) const
Считывание вектора из двумерных точек из базы данных
Definition: StreamParser.h:314

Here is the call graph for this function:

Here is the caller graph for this function:

void MeasureVP::SaveVP ( )

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

Definition at line 273 of file MeasureVP2D.cpp.

274 {
275  W.getTimestat().timeSaveKadr.first += omp_get_wtime();
276 
277  double scaleV = 1.0, scaleP = 1.0;
279  {
280  const double& vRef = W.getPassport().physicalProperties.vRef;
281  scaleV = 1.0 / vRef;
282  scaleP = 1.0 / (0.5 * W.getPassport().physicalProperties.rho * sqr(vRef));
283  }
284 
285 
286  std::ofstream outfile;
287 
288  VMlib::CreateDirectory(W.getPassport().dir, "velPres");
289 
290  if (W.getPassport().timeDiscretizationProperties.fileTypeVP.second == 0) //text format VTK
291  {
292  std::string fname = VMlib::fileNameStep("VelPres", W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "vtk");
293  outfile.open(W.getPassport().dir + "velPres/" + fname);
294 
295  outfile << "# vtk DataFile Version 2.0" << std::endl;
296  outfile << "VM2D VTK result: " << (W.getPassport().dir + "velPres/" + fname).c_str() << " saved " << VMlib::CurrentDataTime() << std::endl;
297  outfile << "ASCII" << std::endl;
298  outfile << "DATASET UNSTRUCTURED_GRID" << std::endl;
299  outfile << "POINTS " << wakeVP->vtx.size() << " float" << std::endl;
300 
301  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
302  {
303  double xi = (wakeVP->vtx[i].r())[0];
304  double yi = (wakeVP->vtx[i].r())[1];
305  outfile << xi << " " << yi << " " << "0.0" << std::endl;
306  }//for i
307 
308  outfile << "CELLS " << wakeVP->vtx.size() << " " << 2 * wakeVP->vtx.size() << std::endl;
309  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
310  outfile << "1 " << i << std::endl;
311 
312  outfile << "CELL_TYPES " << wakeVP->vtx.size() << std::endl;
313  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
314  outfile << "1" << std::endl;
315 
316  outfile << std::endl;
317  outfile << "POINT_DATA " << wakeVP->vtx.size() << std::endl;
318 
320  outfile << "VECTORS V float" << std::endl;
321  else
322  outfile << "VECTORS CV float" << std::endl;
323 
324  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
325  {
326  outfile << velocity[i][0] * scaleV << " " << velocity[i][1] * scaleV << " 0.0" << std::endl;
327  }//for i
328 
329  outfile << std::endl;
330 
332  outfile << "SCALARS P float 1" << std::endl;
333  else
334  outfile << "SCALARS CP float 1" << std::endl;
335 
336  outfile << "LOOKUP_TABLE default" << std::endl;
337 
338  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
339  {
340  outfile << pressure[i] * scaleP << std::endl;
341  }//for i
342 
343  outfile.close();
344  }
345  else if (W.getPassport().timeDiscretizationProperties.fileTypeVP.second == 1) //binary format VTK
346  {
347  //Тест способа хранения чисел
348  uint16_t x = 0x0001;
349  bool littleEndian = (*((uint8_t*)&x));
350  const char eolnBIN[] = "\n";
351 
352  std::string fname = VMlib::fileNameStep("VelPres", W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "vtk");
353  outfile.open(W.getPassport().dir + "velPres/" + fname, std::ios::out | std::ios::binary);
354 
355  outfile << "# vtk DataFile Version 3.0" << "\r\n" << "VM2D VTK result: " << (W.getPassport().dir + "velPres/" + fname).c_str() << " saved " << VMlib::CurrentDataTime() << eolnBIN;
356  outfile << "BINARY" << eolnBIN;
357  outfile << "DATASET UNSTRUCTURED_GRID" << eolnBIN << "POINTS " << wakeVP->vtx.size() << " " << "float" << eolnBIN;
358 
359 
360  Eigen::VectorXf pData = Eigen::VectorXf::Zero(wakeVP->vtx.size());
361  Eigen::VectorXf vData = Eigen::VectorXf::Zero(wakeVP->vtx.size() * 3);
362  Eigen::VectorXf rData = Eigen::VectorXf::Zero(wakeVP->vtx.size() * 3);
363 
364 
365  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
366  {
367  rData(3 * i + 0) = (float)(wakeVP->vtx[i].r())[0];
368  rData(3 * i + 1) = (float)(wakeVP->vtx[i].r())[1];
369  }//for i
370 
371  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
372  {
373  vData(3 * i + 0) = (float)(velocity[i][0] * scaleV);
374  vData(3 * i + 1) = (float)(velocity[i][1] * scaleV);
375  }//for i
376 
377  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
378  {
379  pData(i) = (float)(pressure[i] * scaleP);
380  }//for i
381 
382  //POINTS
383  if (littleEndian)
384  for (int i = 0; i < wakeVP->vtx.size() * 3; ++i)
385  VMlib::SwapEnd(rData(i));
386  outfile.write(reinterpret_cast<char*>(rData.data()), wakeVP->vtx.size() * 3 * sizeof(float));
387 
388  // CELLS
389  std::vector<int> cells(2 * wakeVP->vtx.size());
390  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
391  {
392  cells[2 * i + 0] = 1;
393  cells[2 * i + 1] = (int)i;
394  }
395 
396  std::vector<int> cellsTypes;
397  cellsTypes.resize(wakeVP->vtx.size(), 1);
398 
399  if (littleEndian)
400  {
401  for (int i = 0; i < wakeVP->vtx.size() * 2; ++i)
402  VMlib::SwapEnd(cells[i]);
403 
404  for (int i = 0; i < wakeVP->vtx.size(); ++i)
405  VMlib::SwapEnd(cellsTypes[i]);
406  }
407 
408  outfile << eolnBIN << "CELLS " << wakeVP->vtx.size() << " " << wakeVP->vtx.size() * 2 << eolnBIN;
409  outfile.write(reinterpret_cast<char*>(cells.data()), wakeVP->vtx.size() * 2 * sizeof(int));
410  outfile << eolnBIN << "CELL_TYPES " << wakeVP->vtx.size() << eolnBIN;
411  outfile.write(reinterpret_cast<char*>(cellsTypes.data()), wakeVP->vtx.size() * sizeof(int));
412 
413  //VECTORS V
414  if (littleEndian)
415  for (int i = 0; i < wakeVP->vtx.size() * 3; ++i)
416  VMlib::SwapEnd(vData(i));
417 
418  outfile << eolnBIN << "POINT_DATA " << wakeVP->vtx.size() << eolnBIN;
419 
421  outfile << "VECTORS V " << "float" << eolnBIN;
422  else
423  outfile << "VECTORS CV " << "float" << eolnBIN;
424 
425  outfile.write(reinterpret_cast<char*>(vData.data()), wakeVP->vtx.size() * 3 * sizeof(float));
426 
427  //SCALARS P
428  if (littleEndian)
429  for (int i = 0; i < wakeVP->vtx.size(); ++i)
430  VMlib::SwapEnd(pData(i));
431 
433  outfile << eolnBIN << "SCALARS P " << "float" << " 1" << eolnBIN;
434  else
435  outfile << eolnBIN << "SCALARS CP " << "float" << " 1" << eolnBIN;
436 
437  outfile << "LOOKUP_TABLE default" << eolnBIN;
438  outfile.write(reinterpret_cast<char*>(pData.data()), wakeVP->vtx.size() * sizeof(float));
439 
440  outfile << eolnBIN;
441 
442  outfile.close();
443  }
444  else if (W.getPassport().timeDiscretizationProperties.fileTypeVP.second == 2) //csv
445  {
446  std::string fname = VMlib::fileNameStep("VelPres", W.getPassport().timeDiscretizationProperties.nameLength, W.getCurrentStep(), "csv");
447  outfile.open(W.getPassport().dir + "velPres/" + fname);
448 
450  outfile << "point,x,y,Vx,Vy,P" << std::endl;
451  else
452  outfile << "point,x,y,CVx,CVy,CP" << std::endl;
453 
454  for (size_t i = 0; i < wakeVP->vtx.size(); ++i)
455  {
456  double xi = (wakeVP->vtx[i].r())[0];
457  double yi = (wakeVP->vtx[i].r())[1];
458  outfile << i << "," << xi << "," << yi \
459  << "," << velocity[i][0] * scaleV << "," << velocity[i][1] * scaleV \
460  << "," << pressure[i] * scaleP << std::endl;
461  }//for i
462  outfile.close();
463  }
464  else if (W.getPassport().timeDiscretizationProperties.fileTypeVP.second == 3) //csvBundle
465  {
466  std::string fnameBunCsv = W.getPassport().dir + "velPres/" + "velPresBundle.csv";
467  std::ofstream VPFileBunCsv(fnameBunCsv.c_str(), W.getCurrentStep() ? std::ios::app : std::ios::out);
468  {
469  if (W.getCurrentStep() == 0)
470  VPFileBunCsv << wakeVP->vtx.size() << std::endl;
471 
472 
473  for (size_t q = 0; q < wakeVP->vtx.size(); ++q)
474  {
475  VPFileBunCsv << W.getCurrentStep() << "," \
477  << q << "," \
478  << wakeVP->vtx[q].r()[0] << "," << wakeVP->vtx[q].r()[1] << "," \
479  << velocity[q][0] * scaleV << "," << velocity[q][1] * scaleV << "," \
480  << pressure[q] * scaleP << std::endl;
481  }
482  }
483 
484  }
485 
486 
487  //Вывод в csv-файлы точек "historyPoints"
488  for (size_t q = 0; q < historyPoints.size(); ++q)
489  {
490  std::string VPFileNameCsv;
491  VPFileNameCsv = W.getPassport().dir + "velPres/" + VMlib::fileNameStep("VP-atPoint-", 2, q, "csv");
492 
493  std::ofstream VPFileCsv(VPFileNameCsv.c_str(), W.getCurrentStep() ? std::ios::app : std::ios::out);
494 
495  if (W.getCurrentStep() == 0)
496  {
498  VPFileCsv << "point,time,Vx,Vy,P" << std::endl;
499  else
500  VPFileCsv << "point,time,CVx,CVy,CP" << std::endl;
501  }
502 
503  VPFileCsv << q << "," << W.getPassport().physicalProperties.getCurrTime() << ","
504  << velocity[initialPoints.size() + q][0] * scaleV << ","
505  << velocity[initialPoints.size() + q][1] * scaleV << ","
506  << pressure[initialPoints.size() + q] * scaleP << std::endl;
507  VPFileCsv.close();
508  }
509 
510  W.getTimestat().timeSaveKadr.second += omp_get_wtime();
511 }
std::vector< Point2D > initialPoints
Точки, которые считываются из файла (давление пишется в vtk-файлы)
Definition: MeasureVP2D.h:68
std::pair< std::string, int > fileTypeVP
Тип файлов для сохранения скорости и давления
Definition: PassportGen.h:76
Times & getTimestat() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition: World2D.h:242
const World2D & W
Константная ссылка на решаемую задачу
Definition: MeasureVP2D.h:91
std::string CurrentDataTime()
Формирование строки с текущем временем и датой
Definition: defs.cpp:44
double vRef
Референсная скорость
Definition: Passport2D.h:81
double getCurrTime() const
Возвращает текуще время
Definition: Passport2D.h:102
int nameLength
Число разрядов в имени файла
Definition: PassportGen.h:68
std::string dir
Рабочий каталог задачи
Definition: PassportGen.h:118
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition: Passport2D.h:296
bool calcCoefficients
Признак вычисления коэффициентов вместо сил
Definition: Passport2D.h:289
std::string fileNameStep(const std::string &name, int length, size_t number, const std::string &ext)
Формирование имени файла
Definition: defs.h:357
T sqr(T x)
Возведение числа в квадрат
Definition: defs.h:430
timePeriod timeSaveKadr
Начало и конец процесса сохранения кадра в файл
Definition: Times2D.h:107
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
Definition: PassportGen.h:127
double rho
Плотность потока
Definition: Passport2D.h:75
std::vector< Point2D > velocity
Скорости в нужных точках
Definition: MeasureVP2D.h:78
std::vector< double > pressure
Давление в нужных точках
Definition: MeasureVP2D.h:84
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition: World2D.h:222
std::vector< Point2D > historyPoints
Точки, которые считываются из файла (давление пишется в vtk и csv-файлы)
Definition: MeasureVP2D.h:71
std::unique_ptr< WakeDataBase > wakeVP
Умный указатель на точки, в которых нужно вычислять в данный момент времени Хранятся в виде "мнимой" ...
Definition: MeasureVP2D.h:75
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition: WorldGen.h:91
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<double> VM2D::MeasureVP::domainRadius
private

Радиусы вихрей в нужных точках

Definition at line 81 of file MeasureVP2D.h.

std::vector<Point2D> VM2D::MeasureVP::historyPoints
private

Точки, которые считываются из файла (давление пишется в vtk и csv-файлы)

Definition at line 71 of file MeasureVP2D.h.

std::vector<Point2D> VM2D::MeasureVP::initialPoints
private

Точки, которые считываются из файла (давление пишется в vtk-файлы)

Definition at line 68 of file MeasureVP2D.h.

std::vector<double> VM2D::MeasureVP::pressure
private

Давление в нужных точках

Definition at line 84 of file MeasureVP2D.h.

std::vector<Point2D> VM2D::MeasureVP::velocity
private

Скорости в нужных точках

Definition at line 78 of file MeasureVP2D.h.

const World2D& VM2D::MeasureVP::W
protected

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

Todo:
Сделать учет давления на бесконечности

Definition at line 91 of file MeasureVP2D.h.

std::unique_ptr<WakeDataBase> VM2D::MeasureVP::wakeVP
private

Умный указатель на точки, в которых нужно вычислять в данный момент времени Хранятся в виде "мнимой" системы вихрей, чтобы воспользоваться методом CalcConvVeloToSetOfPoints(...)

Definition at line 75 of file MeasureVP2D.h.


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