VM2D 1.14
Vortex methods for 2D flows simulation
Loading...
Searching...
No Matches
Gmres2D.cpp File Reference

GMRES. More...

#include "Gmres2D.h"
#include "Airfoil2D.h"
#include "Boundary2D.h"
#include "MeasureVP2D.h"
#include "Mechanics2D.h"
#include "Preprocessor.h"
#include "StreamParser.h"
#include "Velocity2D.h"
#include "Wake2D.h"
#include "World2D.h"
Include dependency graph for Gmres2D.cpp:

Go to the source code of this file.

Typedefs

typedef Eigen::Map< const Eigen::VectorXd > MapVec
 

Functions

void PreCalculateCoef (const Airfoil &aflP, size_t nPanelsP, std::vector< Point2D > &prea, std::vector< Point2D > &prec, std::vector< Point2D > &prea1, std::vector< Point2D > &prec1, bool linScheme)
 

Detailed Description

GMRES.

Author
\Version 1.14
Date
6 2026 .

Definition in file Gmres2D.cpp.

Typedef Documentation

◆ MapVec

typedef Eigen::Map<const Eigen::VectorXd> MapVec

Definition at line 57 of file Gmres2D.cpp.

Function Documentation

◆ PreCalculateCoef()

void PreCalculateCoef ( const Airfoil aflP,
size_t  nPanelsP,
std::vector< Point2D > &  prea,
std::vector< Point2D > &  prec,
std::vector< Point2D > &  prea1,
std::vector< Point2D > &  prec1,
bool  linScheme 
)

Definition at line 528 of file Gmres2D.cpp.

536{
537 Point2D p1, s1, p2, s2, di, dj;
538 numvector<double, 3> alpha, lambda;
539 numvector<Point2D, 3> v00, v11;
540 Point2D i00, i11;
541
542 for (size_t i = 0; i < nPanelsP; ++i)
543 {
544 std::array<int, 2> vecV = { ((int)i - 1 >= 0) ? int(i - 1) : int(nPanelsP - 1), (i + 1 < nPanelsP) ? int(i + 1) : 0 };
545 for (int j : vecV)
546 {
547 p1 = aflP.getR(i + 1) - aflP.getR(j + 1);
548 s1 = aflP.getR(i + 1) - aflP.getR(j);
549 p2 = aflP.getR(i) - aflP.getR(j + 1);
550 s2 = aflP.getR(i) - aflP.getR(j);
551 di = aflP.getR(i + 1) - aflP.getR(i);
552 dj = aflP.getR(j + 1) - aflP.getR(j);
553
554 alpha = { \
555 (aflP.isAfter(j, i)) ? 0.0 : VMlib::Alpha(s2, s1), \
556 VMlib::Alpha(s2, p1), \
557 (aflP.isAfter(i, j)) ? 0.0 : VMlib::Alpha(p1, p2) \
558 };
559
560 lambda = { \
561 (aflP.isAfter(j, i)) ? 0.0 : VMlib::Lambda(s2, s1), \
562 VMlib::Lambda(s2, p1), \
563 (aflP.isAfter(i, j)) ? 0.0 : VMlib::Lambda(p1, p2) \
564 };
565
566 v00 = {
567 VMlib::Omega(s1, di.unit(), dj.unit()),
568 -VMlib::Omega(di, di.unit(), dj.unit()),
569 VMlib::Omega(p2, di.unit(), dj.unit())
570 };
571
572 i00 = IDPI / di.length() * (-(alpha[0] * v00[0] + alpha[1] * v00[1] + alpha[2] * v00[2]).kcross() \
573 + (lambda[0] * v00[0] + lambda[1] * v00[1] + lambda[2] * v00[2]));
574
575 if (aflP.isAfter(j, i))
576 prec[i] = (i00.kcross()) * (1.0 / dj.length());
577 if (aflP.isAfter(i, j))
578 prea[i] = (i00.kcross()) * (1.0 / dj.length());
579
580 if (linScheme) {
581
582 v11 = { 1.0 / (12.0 * di.length() * dj.length()) * \
583 (2.0 * (s1 & Omega(s1 - 3.0 * p2, di.unit(), dj.unit())) * VMlib::Omega(s1, di.unit(), dj.unit()) - \
584 s1.length2() * (s1 - 3.0 * p2)) - 0.25 * VMlib::Omega(s1, di.unit(), dj.unit()),
585 -di.length() / (12.0 * dj.length()) * Omega(di, dj.unit(), dj.unit()),
586 -dj.length() / (12.0 * di.length()) * Omega(dj, di.unit(), di.unit()) };
587
588 i11 = (IDPI / di.length()) * ((alpha[0] + alpha[2]) * v11[0] + (alpha[1] + alpha[2]) * v11[1] + alpha[2] * v11[2]\
589 + ((lambda[0] + lambda[2]) * v11[0] + (lambda[1] + lambda[2]) * v11[1] + lambda[2] * v11[2] \
590 + 1.0 / 12.0 * (dj.length() * di.unit() + di.length() * dj.unit() - 2.0 * VMlib::Omega(s1, di.unit(), dj.unit()))).kcross());
591 if (aflP.isAfter(j, i))
592 prec1[i] = (i11) * (1.0 / dj.length());
593 if (aflP.isAfter(i, j))
594 prea1[i] = (i11) * (1.0 / dj.length());
595 }
596 }
597 }
598}
const Point2D & getR(size_t q) const
Возврат константной ссылки на вершину профиля
Definition Airfoil2D.h:113
bool isAfter(size_t i, size_t j) const
Проверка, идет ли вершина i следом за вершиной j.
Definition Airfoil2D.cpp:73
Шаблонный класс, определяющий вектор фиксированной длины Фактически представляет собой массив,...
Definition numvector.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
auto unit(P newlen=1) const -> numvector< typename std::remove_const< decltype(this->data[0] *newlen)>::type, n >
Вычисление орта вектора или вектора заданной длины, коллинеарного данному
Definition numvector.h:402
P length() const
Вычисление 2-нормы (длины) вектора
Definition numvector.h:374
const double IDPI
Число .
Definition defs.h:85
double Lambda(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления логарифма отношения норм векторов
Definition defs.cpp:268
double Alpha(const Point2D &p, const Point2D &s)
Вспомогательная функция вычисления угла между векторами
Definition defs.cpp:262
Point2D Omega(const Point2D &a, const Point2D &b, const Point2D &c)
Вспомогательная функция вычисления величины .
Definition defs.cpp:274
Here is the call graph for this function: