TRUST 1.9.8
HPC thermohydraulic platform
Loading...
Searching...
No Matches
Piso.h
1/****************************************************************************
2* Copyright (c) 2026, CEA
3* All rights reserved.
4*
5* Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
6* 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
7* 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.
8* 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission.
9*
10* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
11* IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
12* OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
13*
14*****************************************************************************/
15
16#ifndef Piso_included
17#define Piso_included
18
19#include <Navier_Stokes_std.h>
20#include <Simpler.h>
21
22//Description
23
24// Ref. Journal of Computational Physics 62 P. 40-65
25// R. I. Issa
26
27// A = (M/dt + C(Uk) + D)
28// Bt et -B designent respectivement les operateurs gradient et divergence
29// delta_x designe l operateur Laplacien
30// H designe l operateur convection+diffusion
31
32// L algorithme PISO (Pressure Implicit Splitting Operator) est non iteratif
33// et se decompose en trois etapes * ** et ***.
34// Les etapes de l algorithme sont :
35
36// -ETAPE DE PREDICTION
37//
38// Recherche du champ de vitesse U* staisfaisant l equation de q.d.m. (avec pression Pn)
39// (rho/dt)*(U*-Un) = H(U*) - BtP + S
40// U* ne satisfait pas l equation de continuite
41
42// L etape de prediction est realisee en traitant le systeme suivant :
43// AU* = -BtPn + Sv + Ss + (M/dt)Un -> U*
44
45
46// -PREMIERE ETAPE DE CORRECTION
47//
48// Recherche de U** et P* satisfaisant la q.d.m. et l equation de continuite (a partir de U* et Un)
49// (rho/dt)*(U**-Un) = H(U*) - BtP* + S
50// delta_x U** = 0
51
52// Afin de stabiliser le systeme a traiter on implicite la partie diagonale de convection-diffusion
53// (rho/dt-Ao)*U** - (rho/dt)*Un = H'(U*) -BtP* + S avec H(U) = H'(U) + AoU
54// En retranchant l equation de prediction :
55// (rho/dt-Ao)*(U**-U*) = -Bt(P*-Pn) et d autre part delta_x U** = 0
56//
57// L etape de premiere correction est realisee en traitant (Da = rho/dt-Ao):
58// (BDa-1Bt)P' = BU* -> P' -> P* = Pn + P'
59// puis Da[Un]U' = -BtP' -> U' -> U** = U* + U'
60
61
62//-SECONDE ETAPE DE CORRECTION
63//
64// Recherche de U*** et P** satisfaisant la q.d.m. et l equation de continuite (a partir de U** et Un)
65// (rho/dt)*(U***-Un) = H(U**) - BtP** + S
66// delta_x U*** = 0
67
68// Afin de stabiliser le systeme a traiter on implicite la partie diagonale de convection-diffusion
69// (rho/dt-Ao)*U*** - (rho/dt)*Un = H'(U**) -BtP** + S
70// En retranchant l equation de premiere correction reformulee a l equation de seconde correction reformulee :
71// (rho/dt-Ao)*(U***-U**) = H'(U**-U*) -Bt(P**-P*) et d autre part delta_x U*** = 0 et delta_x U** = 0
72//
73// L etape de seconde correction est realisee en traitant (E = H'):
74// (BDa-1Bt)P'' = (BDa-1E)U' -> P'' -> P** = P* + P''
75// DaU'' = EU' -BtP'' -> U'' -> U*** = U** + U''
76
77// Rq : Des etapes de correction supplementaires peuvent etre realisees.
78// La version codee ici poursuit les corrections (compt_max-1 maximum)
79// sauf si le residu augmente par rapport a la correction precedente.
80
81class Piso : public Simpler
82{
83 Declare_instanciable(Piso);
84public :
85 void iterer_NS(Equation_base&, DoubleTab& current, DoubleTab& pression, double, Matrice_Morse&, double, DoubleTrav&,int nb_iter,int& converge, int& ok) override;
86protected :
87
88 int nb_corrections_max_ = 21; //nombre de corrections maximum pour affinet la projection
89 int avancement_crank_ = 0; // on ne fait pas vraiment du piso mais plutot du CN
90 int with_sources_ = 0; //prise en compte des termes sources dans la matrice de pression -> plus de termes en implicite, moins de termes en PISO
91
92 Entree& lire(const Motcle&, Entree&) override;
93
94private:
95 virtual void iterer_NS_PolyMAC_CDO(Navier_Stokes_std& eqn,DoubleTab& current,DoubleTab& pression, double dt, Matrice_Morse& matrice, int& ok);
96 // IBM stuff
97 void add_penality_term(Navier_Stokes_std& , DoubleTrav& resu , DoubleTrav& gradP);
98 void correct_incr_pressure(Navier_Stokes_std& , DoubleTrav& secmem);
99 void correct_gradP(Navier_Stokes_std& , DoubleTrav& gradP);
100 void correct_pressure(Navier_Stokes_std& , DoubleTab& , DoubleTab& );
101};
102
103//Description
104// Ref. G. Fauchet
105
106// L algorithme IMPLICITE est non iteratif et se decompose en deux etapes * et **.
107// Les etapes de l algorithme sont :
108
109// -ETAPE DE PREDICTION
110//
111// Recherche du champ de vitesse U* staisfaisant l equation de q.d.m. (avec pression Pn)
112// De facon identique a l algorithme PISO le terme de diffusion est exprime sous forme implicite
113// et le terme de convection sous forme semi-implicite.
114
115// (U*-Un)/dt = H(U*) - BtPn + Sv + Ss
116// U* ne satisfait pas l equation de continuite
117
118// L etape de prediction est realisee en traitant le systeme suivant :
119// AU* = -BtPn + Sv + Ss + (M/dt)Un -> U*
120
121// -ETAPE DE CORRECTION
122//
123// L equation de q.d.m. a l etape n+1 peut s ecrire :
124// (Un+1-Un)/dt = H(Un+1) - BtPn+1 + Sv + Ss
125// En retranchant l 'equation de q.d.m. ecrite pour U* et en negligeant
126// les termes de convection et diffusion sur U' = Un+1 - U* on a la relation :
127// (Un+1-U*)/dt = -BtP' avec P' = Pn+1 - Pn
128//
129// L etape de correction est realisee en traitant (M matrice de masse) :
130// (BM-1Bt)dt*P' = BU* -> dt*P' -> P'
131// Un+1 = U* -dt*BtP' -> Un+1
132//
133// Rq : La matrice de masse est constante par consequent le systeme (BM-1Bt) n est assemble qu une seule fois.
134class Implicite : public Piso
135{
136 Declare_instanciable(Implicite);
137};
138
139#endif /* Piso_included */
Class defining operators and methods for all reading operation in an input flow (file,...
Definition Entree.h:42
classe Equation_base Le role d'une equation est le calcul d'un ou plusieurs champs....
Classe Matrice_Morse Represente une matrice M (creuse), non necessairement carree.
Une chaine de caractere (Nom) en majuscules.
Definition Motcle.h:26
classe Navier_Stokes_std Cette classe porte les termes de l'equation de la dynamique
Definition Piso.h:82
int avancement_crank_
Definition Piso.h:89
int nb_corrections_max_
Definition Piso.h:88
Entree & lire(const Motcle &, Entree &) override
Definition Piso.cpp:60
int with_sources_
Definition Piso.h:90
void iterer_NS(Equation_base &, DoubleTab &current, DoubleTab &pression, double, Matrice_Morse &, double, DoubleTrav &, int nb_iter, int &converge, int &ok) override
Definition Piso.cpp:118