TrioCFD 1.9.8
TrioCFD documentation
Loading...
Searching...
No Matches
Schema_RK_Rationnel.cpp
1/****************************************************************************
2* Copyright (c) 2025, 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#include <Schema_RK_Rationnel.h>
17#include <Equation_base.h>
18
19// XD runge_kutta_rationnel_ordre_2 schema_temps_base runge_kutta_rationnel_ordre_2 INHERITS_BRACE This is the
20// XD_CONT Runge-Kutta rational scheme of second order. The method is described in the note: Wambeck - Rational
21// XD_CONT Runge-Kutta methods for solving systems of ordinary differential equations, at the link:
22// XD_CONT https://link.springer.com/article/10.1007/BF02252381. Although rational methods require more computational
23// XD_CONT work than linear ones, they can have some other properties, such as a stable behaviour with explicitness,
24// XD_CONT which make them preferable. The CFD application of this RRK2 scheme is described in the note:
25// XD_CONT https://link.springer.com/content/pdf/10.1007\%2F3-540-13917-6_112.pdf.
26Implemente_instanciable(RRK2,"Runge_Kutta_Rationnel_ordre_2",TRUSTSchema_RK<Ordre_RK::RATIO_DEUX>);
27
29
31
32/*! @brief Effectue un pas de temps de Runge Kutta rationnel d'ordre 2, sur l'equation passee en parametre.
33 *
34 * Le schema de Runge Kutta rationel d'ordre 2:
35 * g1=hf(y0)
36 * g2=hf(y0+c2g1)
37 * y1=y0+(g1g1)/(b1g1+b2g2)
38 * ou ab/d = (a(b,d)+b(d,a)-d(a,b)/(d,d)
39 * y1=y0+(2g1(g1,b1g1+b2g2)-(b1g1+b2g2)(g1,g1)/(b1g1+b2g2,b1g1+b2g2)
40 * y1=y0+(2g1(g1,"g2")-("g2")(g1,g1)/("g2","g2")
41 * ordre2 si b2c2=-1/2
42 * b2c2<=-1/2 A0 stabilite et I stabilite
43 * b2c2<= 1/(2cos(alpha)(2-cos(alpha))) O<=alpha<pi/2 Aalpha stabilite
44 *
45 */
47{
48 // Warning sur les 100 premiers pas de temps si facsec est egal a 1 pour faire reflechir l'utilisateur
49 if (nb_pas_dt() >= 0 && nb_pas_dt() <= NW && facsec_ == 1) print_warning(NW);
50
51 const double b1 = 2.0, b2 = -1, c2 = 0.5;
52
53 DoubleTab& present = eqn.inconnue().valeurs(), &futur = eqn.inconnue().futur();
54
55 // g1=futur=f(y0)
56 DoubleTrav g1(present), g2(present); // just for initializing the array structure ...
57
58 // sauv=y0
59 DoubleTrav sauv(present);
60 sauv = present;
61
63
64 // g1=hf(y0)
65 g1 *= dt_;
66
67 // present=y0+c2g1
68 present.ajoute(c2, g1);
69
70 // g2=futur=f(y0+c2g1)
72
73 // g2=b2"g2"
74 g2 *= (b2 * dt_);
75
76 // g2=b2"g2" + b1g1
77 g2.ajoute(b1, g1);
78 // g2.axpby(b2*dt_, g2, b1, g1); to fuse 2 kernels
79
80 // normeg2=("g2","g2")
81 double normeg2 = mp_carre_norme_vect(g2) + DMINFLOAT;
82 // normeg1=-(g1,g1)
83 double normeg1 = -mp_carre_norme_vect(g1);
84 // psc1=2(g1,"g2")
85 double psc1 = 2.0 * mp_prodscal(g1, g2);
86
87 // y1=y0+(2g1(g1,"g2")-("g2")(g1,g1)/("g2","g2")
88 // ToDo: implement axpby(a, x, b, y, result);
89 //futur.axpby(psc1/(normeg2*dt), g1, normeg1, g2/(normeg2*dt)); to fuse 4 kernels
90 futur = g1;
91 futur *= psc1;
92 futur.ajoute(normeg1, g2);
93 futur /= (normeg2 * dt_);
94 // futur = (g1 * psc1 + g2 * normeg1)/(normeg2 * dt_)
95 present = sauv;
96 update_critere_statio(futur, eqn);
97 futur *= dt_;
98 futur += sauv;
99 // futur = sauv + (g1 * psc1 + g2 * normeg1)/normeg2
101 futur.echange_espace_virtuel();
102
103 return 1;
104}
DoubleTab & futur(int i=1) override
Renvoie les valeurs du champs a l'instant t+i.
DoubleTab & valeurs() override
Renvoie le tableau des valeurs du champ au temps courant.
virtual void imposer_cond_lim(Champ_Inc_base &, double)=0
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....
virtual const Champ_Inc_base & inconnue() const =0
virtual DoubleTab & derivee_en_temps_inco(DoubleTab &)
Returns the time derivative of the unknown I of the equation: dI/dt = M-1*(sum(operators(I) + sources...
virtual Domaine_Cl_dis_base & domaine_Cl_dis()
Renvoie le domaine des conditions aux limite discretisee associee a l'equation.
virtual Entree & readOn(Entree &)
Lecture d'un Objet_U sur un flot d'entree Methode a surcharger.
Definition Objet_U.cpp:293
virtual Sortie & printOn(Sortie &) const
Ecriture de l'objet sur un flot de sortie Methode a surcharger.
Definition Objet_U.cpp:282
: classe RRK2 Cette classe represente un schema de Runge Kutta
int faire_un_pas_de_temps_eqn_base(Equation_base &) override
Effectue un pas de temps de Runge Kutta rationnel d'ordre 2, sur l'equation passee en parametre.
double temps_courant() const
Renvoie le temps courant.
double dt_
Pas de temps de calcul.
double pas_de_temps() const
Renvoie le pas de temps (delta_t) courant.
int nb_pas_dt() const
Renvoie le nombre de pas de temps effectues.
void update_critere_statio(const DoubleTab &tab_critere, Equation_base &equation)
//Actualisation de stationnaire_atteint_ et residu_ (critere residu_<seuil_statio_)
Classe de base des flux de sortie.
Definition Sortie.h:52
void ajoute(_SCALAR_TYPE_ alpha, const TRUSTVect &y, Mp_vect_options opt=VECT_ALL_ITEMS)
Definition TRUSTVect.tpp:52