VM2D 1.14
Vortex methods for 2D flows simulation
Loading...
Searching...
No Matches
Mechanics2DDeformable.cpp
Go to the documentation of this file.
1/*--------------------------------*- VM2D -*-----------------*---------------*\
2| ## ## ## ## #### ##### | | Version 1.14 |
3| ## ## ### ### ## ## ## ## | VM2D: Vortex Method | 2026/03/06 |
4| ## ## ## # ## ## ## ## | for 2D Flow Simulation *----------------*
5| #### ## ## ## ## ## | Open Source Code |
6| ## ## ## ###### ##### | https://www.github.com/vortexmethods/VM2D |
7| |
8| Copyright (C) 2017-2026 I. Marchevsky, K. Sokol, E. Ryatina, A. Kolganova |
9*-----------------------------------------------------------------------------*
10| File name: Mechanics2DDeformable.cpp |
11| Info: Source code of VM2D |
12| |
13| This file is part of VM2D. |
14| VM2D is free software: you can redistribute it and/or modify it |
15| under the terms of the GNU General Public License as published by |
16| the Free Software Foundation, either version 3 of the License, or |
17| (at your option) any later version. |
18| |
19| VM2D is distributed in the hope that it will be useful, but WITHOUT |
20| ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
21| FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
22| for more details. |
23| |
24| You should have received a copy of the GNU General Public License |
25| along with VM2D. If not, see <http://www.gnu.org/licenses/>. |
26\*---------------------------------------------------------------------------*/
27
28
42
43#include "Airfoil2D.h"
44#include "Airfoil2DDeformable.h"
45#include "Boundary2D.h"
46#include "MeasureVP2D.h"
47#include "StreamParser.h"
48#include "Velocity2D.h"
49#include "Wake2D.h"
50#include "World2D.h"
51
52using namespace VM2D;
53
54
55Beam::Beam(const World2D& W_, bool fsi_, double x0_, double L_, int R_) :
56 W(W_),
57 fsi(fsi_),
58 rho(10000.0),
59 F(0.02),
60 EJ(1.1703239289446188),
61 R(R_),
62 x0(x0_),
63 L(L_)
64{
65 qCoeff.resize(R);
66 currentPhi.resize(R);
67 currentDPhi.resize(R);
69};
70
71double Beam::phi(int n, double t) const
72{
73 return currentPhi[n];
74}
75
76void Beam::solveDU(int n, double dt)
77{
78 double cm = rho * F;
79 double ck = EJ * sqr(sqr(unitLambda[n] / L));
80 double cq = qCoeff[n];
81 //std::cout << "ck = " << ck << std::endl;
82 //std::cout << "lam = " << unitLambda[n] << std::endl;
83 //std::cout << "lam = " << unitLambda[n] / L << std::endl;
84
85 //double omega = sqrt(ck / cm);
86 //double cDamp = 0.0; //0.005;
87
88 double phiAst = currentPhi[n] + 0.5 * dt * currentDPhi[n];
89 double psiAst = currentDPhi[n] + 0.5 * dt * (-ck * currentPhi[n] /*- cDamp * omega * currentDPhi[n]*/ - cq) / cm;
90
91 currentPhi[n] += dt * psiAst;
92 currentDPhi[n] += dt * (-ck * phiAst /*- cDamp * omega * psiAst*/ - cq) / cm;
93}
94
95void Beam::solveDU_RK(int n, double dt) {
96 double cm = rho * F;
97 double ck = EJ * sqr(sqr(unitLambda[n] / L));
98 double cq = qCoeff[n];
99
100 double omega = sqrt(ck / cm);
101 double cDamp = 0.0;
102
103 auto acceleration = [&](double phi, double dphi) {
104 return (-ck * phi - cDamp * omega * dphi - cq) / cm;
105 };
106
107 double phi = currentPhi[n];
108 double dphi = currentDPhi[n];
109
110 double k1_phi = dphi;
111 double k1_dphi = acceleration(phi, dphi);
112
113 double k2_phi = dphi + 0.5 * dt * k1_dphi;
114 double k2_dphi = acceleration(phi + 0.5 * dt * k1_phi, dphi + 0.5 * dt * k1_dphi);
115
116 double k3_phi = dphi + 0.5 * dt * k2_dphi;
117 double k3_dphi = acceleration(phi + 0.5 * dt * k2_phi, dphi + 0.5 * dt * k2_dphi);
118
119 double k4_phi = dphi + dt * k3_dphi;
120 double k4_dphi = acceleration(phi + dt * k3_phi, dphi + dt * k3_dphi);
121
122 currentPhi[n] = phi + (dt / 6.0) * (k1_phi + 2 * k2_phi + 2 * k3_phi + k4_phi);
123 currentDPhi[n] = dphi + (dt / 6.0) * (k1_dphi + 2 * k2_dphi + 2 * k3_dphi + k4_dphi);
124
125 std::cout << "phi = " << currentPhi[n] << std::endl;
126}
127
128
129double Beam::getTotalDisp(double x, double t) const
130{
131 double result = 0.0;
132 for (int i = 0; i < R; ++i)
133 result += phi(i, t) * shape(i, x);
134 return result;
135}
136
137double Beam::getGivenLaw(double x, double t, double deformParam) const //имитатор деформации упругой линии
138{
139 const double c1 = -0.825;
140 const double c2 = 1.625;
141 double alpha;
142 if (fsi) //Turek
143 alpha = 0.0; //Turek
144 else
145 alpha = 0.1; //Fish
146
147
148 const double lambda = 1.0;
149 const double length = 1.0;
150 const double f = deformParam;
151
152 auto A = [c1, c2](double xi) {return 1.0 + (xi - 1.0) * c1 + (xi * xi - 1.0) * c2;};
153
154 //double result = alpha * A(x + 0.5) * sin(DPI * ((x + 0.5) / (length * lambda) - f * t));
155 double result = W.getPassport().physicalProperties.accelCft(t) * alpha * A(x + 0.5) * sin(DPI * ((x + 0.5) / (length * lambda) - f * t));
156
157 return result;
158}
159
160
161
162
163MechanicsDeformable::MechanicsDeformable(const World2D& W_, size_t numberInPassport_)
164 :
165 Mechanics(W_, numberInPassport_, true, true)
166{
167 const auto& airfoil = W_.getAirfoil(numberInPassport_);
168
169 Vcm0 = { 0.0, 0.0 };
170 Rcm0 = { airfoil.rcm[0], airfoil.rcm[1] };
171 Vcm = Vcm0;
172 Rcm = Rcm0;
173 VcmOld = Vcm0;
174 RcmOld = Rcm;
175
176
178
179 Point2D zero = { 0.0, 0.0 };
180
181 Initialize(zero, airfoil.rcm + zero, 0.0, airfoil.phiAfl + 0.0);
182
183 if (airfoil.phiAfl != 0)
184 {
185 W.getInfo('e') << "Airfoil rotation for Turek problem is not allowed" << std::endl;
186 exit(2345);
187 }
188
189 if (fsi)
190 {
191 //Turek
192 //Выделение упругой хорды
194 double x0 = airfoil.getR(indexOfUpperRightAngle)[0];
195 double x1 = x0;
196 while (fabs(x1 - x0) < 1e-12)
197 {
199 x0 = x1;
200 x1 = airfoil.getR(indexOfUpperRightAngle)[0];
201 }
203
205 for (size_t i = 0; i < upperShifts.size(); ++i)
206 upperShifts[i] = airfoil.getR(i + 1)[1] - airfoil.getR(0)[1];
207
208
210 double y0 = airfoil.getR(indexOfUpperLeftAngle)[1];
211 double y1 = y0;
212 while (fabs(y1 - y0) < 1e-12)
213 {
215 y0 = y1;
216 y1 = airfoil.getR(indexOfUpperLeftAngle)[1];
217 }
219
220
221 indexOfLowerRightAngle = airfoil.getNumberOfPanels() + 1;
222 x0 = airfoil.getR(indexOfLowerRightAngle)[0];
223 x1 = x0;
224 while (fabs(x1 - x0) < 1e-12)
225 {
227 x0 = x1;
228 x1 = airfoil.getR(indexOfLowerRightAngle)[0];
229 }
231
232
233 lowerShifts.resize(airfoil.getNumberOfPanels() - indexOfLowerRightAngle - 1);
234 for (size_t i = 0; i < lowerShifts.size(); ++i)
235 lowerShifts[i] = airfoil.getR(airfoil.getNumberOfPanels() - 1 - i)[1] - airfoil.getR(0)[1];
236
238 y0 = airfoil.getR(indexOfLowerLeftAngle)[1];
239 y1 = y0;
240 while (fabs(y1 - y0) < 1e-12)
241 {
243 y0 = y1;
244 y1 = airfoil.getR(indexOfLowerLeftAngle)[1];
245 }
247
248 //W.getInfo('i') << "UR: " << airfoil.getR(indexOfUpperRightAngle) << std::endl;
249 //W.getInfo('i') << "UL: " << airfoil.getR(indexOfUpperLeftAngle) << std::endl;
250 //W.getInfo('i') << "LR: " << airfoil.getR(indexOfLowerRightAngle) << std::endl;
251 //W.getInfo('i') << "LL: " << airfoil.getR(indexOfLowerLeftAngle) << std::endl;
252
254 {
255 W.getInfo('e') << "indexOfUpperLeftAngle - indexOfUpperRightAngle != indexOfLowerRightAngle - indexOfLowerLeftAngle" << std::endl;
256 exit(2346);
257 }
258
260
261 for (size_t i = 0; i < indexOfUpperLeftAngle - indexOfUpperRightAngle; ++i)
262 {
263 size_t idxUp = indexOfUpperLeftAngle - 1 - i;
264 size_t idxDn = indexOfLowerLeftAngle + i;
265 Point2D rUpLeft = airfoil.getR(idxUp + 1);
266 Point2D rUpRight = airfoil.getR(idxUp);
267 Point2D rDnLeft = airfoil.getR(idxDn);
268 Point2D rDnRight = airfoil.getR(idxDn + 1);
269
270 if (std::max(fabs(rUpLeft[0] - rDnLeft[0]), fabs(rUpRight[0] - rDnRight[0])) > \
271 0.01 * std::min((rUpRight[0] - rUpLeft[0]), (rDnRight[0] - rDnLeft[0])))
272 {
273 W.getInfo('e') << "x_up != x_dn" << std::endl;
274 exit(2347);
275 }
276
277 chord[i].beg = 0.5 * (rUpLeft + rDnLeft);
278 chord[i].end = 0.5 * (rUpRight + rDnRight);
279 chord[i].infPanels = { idxUp, idxDn };
280 chord[i].rightSemiWidth = 0.5 * (rUpRight - rDnRight)[1];
281 }
283 beam = std::make_unique<Beam>(W, fsi, chord[0].beg[0], initialChord.back().end[0] - initialChord[0].beg[0], 3); //Beam
284 }
285 else //Fish
286 {
287 int np = (int)airfoil.getNumberOfPanels();
288 chord.resize(np / 2);
289
290 chord[0].beg = airfoil.getR(np / 2);
291 chord[np / 2 - 1].end = airfoil.getR(0);
292 chord[np / 2 - 1].rightSemiWidth = 0.0;
293
294 for (size_t i = 0; i < np / 2; ++i)
295 {
296 if (i != 0)
297 chord[i].beg = 0.5 * (airfoil.getR(np / 2 - i) + airfoil.getR(np / 2 + i));
298
299 if (i != np / 2 - 1)
300 chord[i].end = 0.5 * (airfoil.getR(np / 2 - i - 1) + airfoil.getR(np / 2 + i + 1));
301
302 chord[i].infPanels = { np / 2 - i, np / 2 + i + 1 };
303
304 if (i != np / 2 - 1)
305 chord[i].rightSemiWidth = (airfoil.getR(np / 2 - i - 1) - airfoil.getR(np / 2 + i + 1)).length() * 0.5;
306 }
308 beam = std::make_unique<Beam>(W, fsi, chord[0].beg[0], initialChord.back().end[0] - initialChord[0].beg[0], 3); //Beam
309 }
310
311 //std::ofstream of("chord.txt");
312 //for (size_t i = 0; i < chord.size(); ++i)
313 // of << chord[i].beg[0] << " " << chord[i].beg[1] << " " << chord[i].end[0] << " " << chord[i].end[1] << " " << chord[i].infPanels.first << " " << chord[i].infPanels.second << std::endl;
314 //of.close();
315
316
317
318 initialPossibleWays = airfoil.possibleWays;
319};
320
321//Вычисление гидродинамической силы, действующей на профиль
323{
324 W.getTimers().start("Force");
325
326 const double& dt = W.getPassport().timeDiscretizationProperties.dt;
327
328 hydroDynamForce = { 0.0, 0.0 };
329 hydroDynamMoment = 0.0;
330
331 viscousForce = { 0.0, 0.0 };
332 viscousMoment = 0.0;
333
334 Point2D hDFGam = { 0.0, 0.0 }; //гидродинамические силы, обусловленные присоед.завихренностью
335 Point2D hDFdelta = { 0.0, 0.0 }; //гидродинамические силы, обусловленные приростом завихренности
336 Point2D hDFQ = { 0.0, 0.0 }; //гидродинамические силы, обусловленные присоед.источниками
337
338 double hDMGam = 0.0; //гидродинамический момент, обусловленный присоед.завихренностью
339 double hDMdelta = 0.0; //гидродинамический момент, обусловленный приростом завихренности
340 double hDMQ = 0.0; //гидродинамический момент, обусловленный присоед.источниками
341 for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
342 {
343 Point2D rK = 0.5 * (afl.getR(i + 1) + afl.getR(i)) - afl.rcm;
344
345 Point2D velK = 0.5 * (afl.getV(i) + afl.getV(i + 1));
346 double gAtt = (velK & afl.tau[i]);
347
348 double gAttOld = 0.0;
349 if (W.getCurrentStep() > 0)
350 {
351 auto oldAfl = W.getOldAirfoil(numberInPassport);
352 gAttOld = ((0.5 * (oldAfl.getV(i) + oldAfl.getV(i + 1))) & oldAfl.tau[i]);
353 }
354
355 double deltaGAtt = gAtt - gAttOld;
356
357 double qAtt = (velK & afl.nrm[i]);
358
360 double deltaK = boundary.sheets.freeVortexSheet(i, 0) * afl.len[i] - afl.gammaThrough[i] + deltaGAtt * afl.len[i];
361
362 /*1*/
363 hDFdelta += deltaK * Point2D({ -rK[1], rK[0] });
364 hDMdelta += 0.5 * deltaK * rK.length2();
365
366 /*2*/
367 hDFGam += 0.5 * velK.kcross() * gAtt * afl.len[i];
368 hDMGam += 0.5 * (rK ^ velK.kcross()) * gAtt * afl.len[i];
369
370 /*3*/
371 hDFQ -= 0.5 * velK * qAtt * afl.len[i];
372 hDMQ -= 0.5 * (rK ^ velK) * qAtt * afl.len[i];
373 }
374
375 const double rho = W.getPassport().physicalProperties.rho;
376
377 hydroDynamForce = rho * (hDFGam + hDFdelta * (1.0 / dt) + hDFQ);
378 hydroDynamMoment = rho * (hDMGam + hDMdelta / dt + hDMQ);
379
380 if ((W.getPassport().physicalProperties.nu > 0.0)/* && (W.currentStep > 0)*/)
381 for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
382 {
383 Point2D rK = 0.5 * (afl.getR(i + 1) + afl.getR(i)) - afl.rcm;
384 viscousForce += rho * afl.viscousStress[i] * afl.tau[i];
385 viscousMoment += rho * (afl.viscousStress[i] * afl.tau[i]) & rK;
386 }
387
388 W.getTimers().stop("Force");
389}// GetHydroDynamForce()
390
391// Вычисление скорости центра масс
393{
394 return Vcm;
395}//VeloOfAirfoilRcm(...)
396
397// Вычисление положения центра масс
399{
400 return Rcm;
401}//PositionOfAirfoilRcm(...)
402
404{
405 return Wcm;
406}//AngularVelocityOfAirfoil(...)
407
409{
410 if (afl.phiAfl != Phi)
411 {
412 std::cout << "afl.phiAfl != Phi" << std::endl;
413 exit(100600);
414 }
415
416 return afl.phiAfl;
417}//AngleOfAirfoil(...)
418
419// Вычисление скоростей начал панелей
421{
422 std::vector<Point2D> veloW(afl.getNumberOfPanels(), {0.0, 0.0});
423
424 //if (W.getCurrentStep() == 0)
425 // for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
426 // veloW[i] = { 0.0, 0.0 };//afl.getR(i).kcross();
427
428 if (W.getCurrentStep() > 0)
429 for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
430 veloW[i] = (1.0 / W.getPassport().timeDiscretizationProperties.dt) * (afl.getR(i) - W.getOldAirfoil(0).getR(i));
431
432 afl.setV(veloW);
433
434
435 //Циркуляция
437 circulation = 0.0;
438 for (size_t i = 0; i < afl.getNumberOfPanels(); ++i)
439 circulation += 0.5 * afl.len[i] * ((afl.getV(i) + afl.getV(i + 1)) & afl.tau[i]);
440
441}//VeloOfAirfoilPanels(...)
442
443
445{
446 double x0 = chord[0].beg[0];
447 double t = W.getCurrentTime();
448
449
450 if (fsi) //turek
451 {
452 for (int i = 0; i < beam->R; ++i)
454
455 for (size_t i = 0; i < chord.size(); ++i)
456 {
457 chord[i].beg[1] = beam->getTotalDisp(chord[i].beg[0], t);
458 chord[i].end[1] = beam->getTotalDisp(chord[i].end[0], t);
459 }
460
461 std::vector<Point2D> upperPoints(chord.size() + upperShifts.size());
462 std::vector<Point2D> lowerPoints(chord.size() + lowerShifts.size());
463
464 for (size_t i = 0; i < chord.size() - 1; ++i)
465 {
466 const Point2D& begIp1 = chord[i + 1].beg;
467 const Point2D& endI = chord[i].end;
468 Point2D normI = (endI - chord[i].beg).unit().kcross();
469 Point2D normIp1 = (chord[i + 1].end - begIp1).unit().kcross();
470 upperPoints[i] = endI + 0.5 * chord[i].rightSemiWidth * (normI + normIp1);
471 lowerPoints[i] = endI - 0.5 * chord[i].rightSemiWidth * (normI + normIp1);
472 }
473 Point2D normBack = (chord.back().end - chord.back().beg).unit().kcross();
474 upperPoints[chord.size() - 1] = chord.back().end + chord.back().rightSemiWidth * normBack;
475 lowerPoints[chord.size() - 1] = chord.back().end - chord.back().rightSemiWidth * normBack;
476
477 for (size_t i = 0; i < upperShifts.size(); ++i)
478 upperPoints[chord.size() + i] = chord.back().end + upperShifts[upperShifts.size() - 1 - i] * normBack;
479 for (size_t i = 0; i < lowerShifts.size(); ++i)
480 lowerPoints[chord.size() + i] = chord.back().end + lowerShifts[lowerShifts.size() - 1 - i] * normBack;
481
482 afl.setR(0) = { chord.back().end[0], beam->getTotalDisp(chord.back().end[0], t) };
483 for (size_t i = 0; i < upperPoints.size(); ++i)
484 afl.setR(indexOfUpperLeftAngle - 1 - i) = upperPoints[i];
485 for (size_t i = 0; i < lowerPoints.size(); ++i)
486 afl.setR(indexOfLowerLeftAngle + 1 + i) = lowerPoints[i];
487
488 /*
489 std::ofstream file;
490 file.open(W.getPassport().dir + "file" + std::to_string(W.getCurrentStep()) + ".txt");
491 for (size_t c = 0; c < chord.size(); ++c)
492 file << std::endl;
493 file.close();
494 //*/
495 }
496 else //Fish
497 {
498 for (size_t i = 0; i < chord.size(); ++i)
499 {
500 chord[i].beg[1] = beam->getGivenLaw(chord[i].beg[0], t, deformParam);
501 chord[i].end[1] = beam->getGivenLaw(chord[i].end[0], t, deformParam);
502 }
503
504 std::vector<Point2D> upperPoints(chord.size());
505 std::vector<Point2D> lowerPoints(chord.size());
506
507 for (size_t i = 0; i < chord.size() - 1; ++i)
508 {
509 const Point2D& begIp1 = chord[i + 1].beg;
510 const Point2D& endI = chord[i].end;
511 Point2D normI = (endI - chord[i].beg).unit().kcross();
512 Point2D normIp1 = (chord[i + 1].end - begIp1).unit().kcross();
513 upperPoints[i] = endI + 0.5 * chord[i].rightSemiWidth * (normI + normIp1);
514 lowerPoints[i] = endI - 0.5 * chord[i].rightSemiWidth * (normI + normIp1);
515 }
516
517 Point2D normFront = (chord[0].end - chord[0].beg).unit().kcross();
518 Point2D normBack = (chord.back().end - chord.back().beg).unit().kcross();
519
520 int nph = (int)chord.size();
521
522 afl.setR(0) = chord.back().end;
523 afl.setR(nph) = chord[0].beg;
524 for (size_t i = 0; i < nph - 1; ++i)
525 {
526 afl.setR(nph - 1 - i) = upperPoints[i];
527 afl.setR(nph + 1 + i) = lowerPoints[i];
528 }
529 }
530
531
534
535 afl.possibleWays.clear();
536 afl.possibleWays.resize(initialPossibleWays.size());
537
538 for (size_t w = 0; w < initialPossibleWays.size(); ++w)
539 {
540 std::vector<Point2D> initialWay = initialPossibleWays[w];
541 afl.possibleWays[w].resize(initialWay.size());
542 for (size_t p = 0; p < initialWay.size(); ++p)
543 {
544 double x = initialWay[p][0];
545 double y;
546
547 if (fsi) //turek
548 y = beam->getTotalDisp(x, t);
549 else //fish
550 y = beam->getGivenLaw(x, t, deformParam);
551
552 afl.possibleWays[w][p] = Point2D{ x,y };
553 }
554 }
555
557}//Move()
558
559
560
561
562#if defined(INITIAL) || defined(BRIDGE)
564{
565 mechParamsParser->get("deformParam", deformParam);
566 W.getInfo('i') << "mechDeformable: " << "deformParam = " << deformParam << std::endl;
567
568 mechParamsParser->get("fsi", fsi);
569 W.getInfo('i') << "mechDeformable: " << "fsi = " << fsi << std::endl;
570}//ReadSpecificParametersFromDictionary()
571#endif
Заголовочный файл с описанием класса Airfoil.
Заголовочный файл с описанием класса AirfoilDeformable.
Заголовочный файл с описанием класса Boundary.
Заголовочный файл с описанием класса MeasureVP.
Заголовочный файл с описанием класса MechanicsDeformable.
Заголовочный файл с описанием класса StreamParser.
Заголовочный файл с описанием класса Velocity.
Заголовочный файл с описанием класса Wake.
Заголовочный файл с описанием класса World2D.
double phiAfl
Поворот профиля
Definition Airfoil2D.h:100
std::vector< double > len
Длины панелей профиля
Definition Airfoil2D.h:94
void setV(const Point2D &vel)
Установка постоянной скорости всех вершин профиля
Definition Airfoil2D.h:145
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition Airfoil2D.h:113
const Point2D & getV(size_t q) const
Возврат константной ссылки на скорость вершины профиля
Definition Airfoil2D.h:137
std::vector< Point2D > nrm
Нормали к панелям профиля
Definition Airfoil2D.h:81
std::vector< Point2D > tau
Касательные к панелям профиля
Definition Airfoil2D.h:91
Point2D rcm
Положение центра масс профиля
Definition Airfoil2D.h:97
size_t getNumberOfPanels() const
Возврат количества панелей на профиле
Definition Airfoil2D.h:163
Point2D & setR(size_t q)
Возврат ссылки на вершину профиля
Definition Airfoil2D.h:125
std::vector< double > gammaThrough
Суммарные циркуляции вихрей, пересекших панели профиля на прошлом шаге
Definition Airfoil2D.h:276
virtual void GetGabarits(double gap=0.02)
Вычисляет габаритный прямоугольник профиля
void CalcNrmTauLen()
Вычисление нормалей, касательных и длин панелей по текущему положению вершин
void lightningTest()
Тест на "отвещенность".
std::vector< double > viscousStress
Нейросеть для коэффициентов I0 и I3 диффузионной скорости
Definition Airfoil2D.h:268
std::vector< std::vector< Point2D > > possibleWays
Возможные пути внутри профиля от точки (0, 0) к центрам всех панелей
Definition Airfoil2D.h:199
std::vector< double > qCoeff
double getGivenLaw(double x, double t, double deformParam) const
void solveDU_RK(int n, double dt)
Beam(const World2D &W_, bool fsi_, double x0_, double L_, int R_)
std::vector< std::vector< double > > presLastSteps
double shape(int n, double x) const
const std::vector< double > unitLambda
std::vector< double > currentPhi
void solveDU(int n, double dt)
std::vector< double > currentDPhi
double phi(int n, double t) const
double getTotalDisp(double x, double t) const
const size_t nLastSteps
Sheet sheets
Слои на профиле
Definition Boundary2D.h:96
virtual Point2D PositionOfAirfoilRcm(double currTime) override
Вычисление положения центра масс профиля
virtual Point2D VeloOfAirfoilRcm(double currTime) override
Вычисление скорости центра масс профиля
virtual void Move() override
Перемещение профиля в соответствии с законом
virtual void VeloOfAirfoilPanels(double currTime) override
Вычисление скоростей начал панелей
MechanicsDeformable(const World2D &W_, size_t numberInPassport_)
Конструктор
virtual void GetHydroDynamForce() override
Вычисление гидродинамической силы, действующей на профиль
virtual void ReadSpecificParametersFromDictionary() override
Чтение параметров конкретной механической системы
std::vector< ChordPanel > initialChord
virtual double AngularVelocityOfAirfoil(double currTime) override
Вычисление угловой скорости профиля
std::vector< std::vector< Point2D > > initialPossibleWays
std::vector< double > upperShifts
std::vector< double > lowerShifts
std::vector< ChordPanel > chord
std::unique_ptr< Beam > beam
virtual double AngleOfAirfoil(double currTime) override
Вычисление угла поворота профиля
Абстрактный класс, определяющий вид механической системы
Definition Mechanics2D.h:72
std::unique_ptr< VMlib::StreamParser > mechParamsParser
Умный указатель на парсер параметров механической системы
Definition Mechanics2D.h:98
Point2D hydroDynamForce
Вектор гидродинамической силы и момент, действующие на профиль
Point2D Vcm0
Начальная скорость центра и угловая скорость
const size_t numberInPassport
Номер профиля в паспорте
Definition Mechanics2D.h:82
Point2D RcmOld
Текущие положение профиля
Point2D VcmOld
Скорость и отклонение с предыдущего шага
const World2D & W
Константная ссылка на решаемую задачу
Definition Mechanics2D.h:79
void Initialize(Point2D Vcm0_, Point2D Rcm0_, double Wcm0_, double Phi0_)
Задание начального положения и начальной скорости
Point2D Rcm
Текущие положение профиля
Point2D Rcm0
Начальное положение профиля
Point2D viscousForce
Вектор силы и момент вязкого трения, действующие на профиль
double circulationOld
Циркуляция скорости по границе профиля с предыдущего шага
double hydroDynamMoment
Airfoil & afl
Definition Mechanics2D.h:87
Point2D Vcm
Текущие скорость центра и угловая скорость
double circulation
Текущая циркуляция скорости по границе профиля
const Boundary & boundary
Definition Mechanics2D.h:91
PhysicalProperties physicalProperties
Структура с физическими свойствами задачи
Definition Passport2D.h:289
const double & freeVortexSheet(size_t n, size_t moment) const
Definition Sheet2D.h:100
Класс, опеделяющий текущую решаемую задачу
Definition World2D.h:74
const Airfoil & getAirfoil(size_t i) const
Возврат константной ссылки на объект профиля
Definition World2D.h:157
const AirfoilGeometry & getOldAirfoil(size_t i) const
Возврат константной ссылки на объект старого профиля
Definition World2D.h:163
VMlib::TimersGen & getTimers() const
Возврат ссылки на временную статистику выполнения шага расчета по времени
Definition World2D.h:276
const Passport & getPassport() const
Возврат константной ссылки на паспорт
Definition World2D.h:251
TimeDiscretizationProperties timeDiscretizationProperties
Структура с параметрами процесса интегрирования по времени
void stop(const std::string &timerLabel)
Останов счетчика
Definition TimesGen.cpp:68
void start(const std::string &timerLabel)
Запуск счетчика
Definition TimesGen.cpp:55
VMlib::LogStream & getInfo() const
Возврат ссылки на объект LogStream Используется в техничеcких целях для организации вывода
Definition WorldGen.h:82
double getCurrentTime() const
Definition WorldGen.h:100
size_t getCurrentStep() const
Возврат константной ссылки на параметры распараллеливания по MPI.
Definition WorldGen.h:99
numvector< T, 2 > kcross() const
Геометрический поворот двумерного вектора на 90 градусов
Definition numvector.h:510
auto length2() const -> typename std::remove_const< typename std::remove_reference< decltype(this->data[0])>::type >::type
Вычисление квадрата нормы (длины) вектора
Definition numvector.h:386
size_t size() const
Definition numvector.h:114
const double DPI
Число .
Definition defs.h:88
double nu
Коэффициент кинематической вязкости среды
Definition Passport2D.h:99
double rho
Плотность потока
Definition Passport2D.h:75
double dt
Шаг по времени
Definition PassportGen.h:67