TRUST 1.9.8
HPC thermohydraulic platform
Loading...
Searching...
No Matches
Pb_Dilatable_base.cpp
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#include <Schema_Euler_explicite.h>
17#include <Schema_Euler_Implicite.h>
18#include <Fluide_Dilatable_base.h>
19#include <Schema_RK_Williamson.h>
20#include <Schema_RK_Rationnel.h>
21#include <Loi_Fermeture_base.h>
22#include <Pb_Dilatable_base.h>
23#include <Probleme_Couple.h>
24#include <Equation_base.h>
25#include <Perf_counters.h>
26#include <Domaine.h>
27#include <Debog.h>
28
29Implemente_base(Pb_Dilatable_base,"Pb_Dilatable_base",Pb_Fluide_base);
30
32
34
40
42{
43 if (!sub_type(Schema_Euler_explicite,sch) && !sub_type(Schema_Euler_Implicite,sch)
44 && !sub_type(RRK2,sch) && !sub_type(RK2,sch) && !sub_type(RK3,sch) && !sub_type(RK4,sch))
45 {
46 Cerr << "TRUST can't solve a " << que_suis_je() << " with a " << sch.que_suis_je() << " time scheme." << finl;
48 }
49 if ( sub_type(Schema_Euler_Implicite,sch) && is_coupled())
50 {
51 Cerr << finl;
52 Cerr << "Coupled problem with unique Euler implicit time scheme with " << que_suis_je() << "fluid are not allowed!" << finl;
53 Cerr << "Choose another time scheme or set a time scheme for each of your problems." << finl;
55 }
57}
58
60{
61 Fluide_Dilatable_base& le_fluide = ref_cast(Fluide_Dilatable_base,milieu());
62 if (le_fluide.type_fluide()=="Gaz_Reel") equation(1).inconnue().nommer("temperature");
63 if (le_fluide.type_fluide()=="Melange_Binaire") equation(1).inconnue().nommer("fraction_massique");
64
65 le_fluide.completer(*this);
66 le_fluide.preparer_calcul();
68 le_fluide.calculer_masse_volumique();
69 le_fluide.preparer_calcul();
70}
71
73{
74 bool ok = Pb_Fluide_base::initTimeStep(dt);
75 le_fluide_->preparer_pas_temps();
76 return ok;
77}
78
80{
82 equation(1).mettre_a_jour(temps); // thermique
83 equation(0).mettre_a_jour(temps); // NS
84
85 for(int i=2; i<nombre_d_equations(); i++) // if species ...
86 equation(i).mettre_a_jour(temps);
87
88 les_postraitements_.mettre_a_jour(temps);
89 domaine().mettre_a_jour(temps,domaine_dis(),*this);
90 for (auto& itr : liste_loi_fermeture_)
91 {
92 Loi_Fermeture_base& loi=itr.valeur();
93 loi.mettre_a_jour(temps);
94 }
95}
96
98{
101 double temps_present=sch.temps_courant();
102 double temps_futur=temps_present+sch.pas_de_temps();
103
104 //1. Solve all the equations except the first one (Navier Stokes, solved later at //6)
105 for (int i=1; i<nombre_d_equations(); i++)
106 {
108 statistics().begin_count(STD_COUNTERS::update_variables,statistics().get_last_opened_counter_level()+1);
109 equation(i).milieu().mettre_a_jour(temps_futur);;
110 equation(i).inconnue().mettre_a_jour(temps_futur);
111 statistics().end_count(STD_COUNTERS::update_variables);
112 }
113
114 statistics().begin_count(STD_COUNTERS::update_variables,statistics().get_last_opened_counter_level()+1);
115
116 //2. Compute temperature-dependent coefficients
117 le_fluide_->calculer_coeff_T();
118
119 //3. Compute volumic mass (update Cp)
120 le_fluide_->calculer_masse_volumique();
121
122 //4. Solve EDO equation (pressure) if needed (ie. QC)
123 le_fluide_->Resoudre_EDO_PT();
124
125 //5. Compute volumic mass
126 le_fluide_->calculer_masse_volumique();
127 statistics().end_count(STD_COUNTERS::update_variables);
128 //6. Solve Navier Stokes equation
130 statistics().begin_count(STD_COUNTERS::update_variables,statistics().get_last_opened_counter_level()+1);
131 equation(0).milieu().mettre_a_jour(temps_futur);
132 equation(0).inconnue().mettre_a_jour(temps_futur);
133 statistics().end_count(STD_COUNTERS::update_variables);
134 // Update pressure fields (total/thermo/hydro) if necessary
135 update_pressure_fields(temps_futur);
136
137 // on recule les inconnues (le pb mettra a jour les equations)
138 for (int i=0; i<nombre_d_equations(); i++)
140
141 // Calculs coeffs echange sur l'instant sur lequel doivent agir les operateurs.
142 double tps=schema_temps().temps_defaut();
143 for(int i=0; i<nombre_d_equations(); i++)
145
146 converged=true;
147 return true;
148}
149
151{
152 le_fluide_->update_pressure_fields(t);
153}
void mettre_a_jour(double temps) override
Effectue une mise a jour en temps du champ inconnue.
Champ_Inc_base & reculer(int i=1)
Recule le pointeur courant de i pas de temps, dans la liste des valeurs temporelles conservees.
static void set_nom_pb_actuel(const Nom &nom)
Definition Debog.cpp:44
virtual int calculer_coeffs_echange(double temps)
Calcul des coefficients d'echange pour les problemes couples thermiques.
virtual void mettre_a_jour(double temps, Domaine_dis_base &, Probleme_base &)
Class defining operators and methods for all reading operation in an input flow (file,...
Definition Entree.h:42
virtual const Milieu_base & milieu() const =0
virtual const Champ_Inc_base & inconnue() const =0
virtual void mettre_a_jour(double temps)
La valeur de l'inconnue sur le pas de temps a ete calculee.
virtual Domaine_Cl_dis_base & domaine_Cl_dis()
Renvoie le domaine des conditions aux limite discretisee associee a l'equation.
void nommer(const Nom &) override
Donne un nom au champ.
classe Fluide_Dilatable_base Cette classe represente un d'un fluide dilatable,
void preparer_calcul() override
Prepare le fluide au calcul.
const Nom type_fluide() const
virtual void completer(const Probleme_base &)
Complete le fluide avec les champs inconnus associes au probleme.
: Classe de base des lois de fermetures.
virtual void mettre_a_jour(double temps)
Cette methode est appelee par le probleme apres mettre_a_jour() des equations et du milieu.
classe Milieu_base Cette classe est la base de la hierarchie des milieux (physiques)
Definition Milieu_base.h:50
virtual void mettre_a_jour(double temps)
const Nom & que_suis_je() const
renvoie la chaine identifiant la classe.
Definition Objet_U.cpp:104
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 Pb_Dilatable_base Cette classe est censee factoriser ce qui est commun a l'ensemble
void mettre_a_jour(double temps) override
Effectue une mise a jour en temps du probleme.
void associer_sch_tps_base(const Schema_Temps_base &) override
Associe un schema en temps au probleme.
bool iterateTimeStep(bool &converged) override
In the case solveTimeStep uses an iterative process, this method executes a single iteration.
bool initTimeStep(double dt) override
This method allocates and initializes the unknown and given fields for the future time step.
virtual void update_pressure_fields(double)
void preparer_calcul() override
Prepare le calcul: initialise les parametres du milieu et prepare le calcul de chacune des equations.
void associer_milieu_base(const Milieu_base &) override
Associe un milieu physique aux equations du probleme.
classe Pb_Fluide_base Cette classe a pour but de disposer d une classe amont pour
const Nom & le_nom() const override
Donne le nom de l'Objet_U Methode a surcharger : renvoie "neant" dans cette implementation.
Definition Probleme_U.h:109
virtual void associer_milieu_base(const Milieu_base &)
Associe un milieu physique aux equations du probleme.
const Domaine & domaine() const
Renvoie le domaine associe au probleme.
virtual void preparer_calcul()
Prepare le calcul: initialise les parametres du milieu et prepare le calcul de chacune des equations.
Postraitements les_postraitements_
virtual void associer_sch_tps_base(const Schema_Temps_base &)
Associe un schema en temps au probleme.
bool is_coupled() const
virtual const Milieu_base & milieu() const
Renvoie le milieu physique associe au probleme.
const Schema_Temps_base & schema_temps() const
Renvoie le schema en temps associe au probleme.
virtual int nombre_d_equations() const =0
virtual const Equation_base & equation(int) const =0
bool initTimeStep(double dt) override
This method allocates and initializes the unknown and given fields for the future time step.
const Domaine_dis_base & domaine_dis() const
Renvoie le domaine discretise associe au probleme.
static void exit(int exit_code=-1)
Routine de sortie de TRUST dans une region Kokkos.
Definition Process.cpp:455
: classe RK2 Cette classe represente un schema en temps de Runge Kutta d'ordre 2, cas 1 de Williamson...
: classe RK3 Cette classe represente un schema en temps de Runge Kutta d'ordre 3, cas 7 de Williamson...
: classe RK4 Cette classe represente un schema en temps de Runge Kutta d'ordre 4 degnere (schema a tr...
: classe RRK2 Cette classe represente un schema de Runge Kutta
: classe Schema_Euler_explicite Cette classe represente un schema en temps d'Euler explicite: U(n+1) ...
class Schema_Temps_base
double temps_courant() const
Renvoie le temps courant.
double pas_de_temps() const
Renvoie le pas de temps (delta_t) courant.
virtual int faire_un_pas_de_temps_eqn_base(Equation_base &)=0
virtual double temps_defaut() const =0
Classe de base des flux de sortie.
Definition Sortie.h:52