]> SALOME platform Git repositories - tools/solverlab.git/blob - CoreFlows/Models/inc/TransportEquation.hxx
Salome HOME
Tested mesh fast equivalence when giving an input field
[tools/solverlab.git] / CoreFlows / Models / inc / TransportEquation.hxx
1 //============================================================================
2 /**
3  * \file TransportEquation.hxx
4  * \author Michael NDJINGA
5  * \version 1.0
6  * \date 24 March 2015
7  * \brief Fluid enthalpy transport equation
8  * */
9 //============================================================================
10
11 /*! \class TransportEquation TransportEquation.hxx "TransportEquation.hxx"
12  *  \brief Scalar advection equation for a fluid enthalpy
13  *  \details see \ref TransportEqPage for more details
14  */
15 #ifndef TransportEquation_HXX_
16 #define TransportEquation_HXX_
17
18 #include "ProblemCoreFlows.hxx"
19
20 using namespace std;
21
22
23 //! enumeration phase
24 /*! The fluid type can be LiquidPhase or water  */
25 enum phase
26 {
27         LiquidPhase,/**< Fluid considered is GasPhase */
28         GasPhase/**< Fluid considered is Gas */
29 };
30
31 //! enumeration pressureEstimate
32 /*! the pressure estimate needed to fit physical parameters  */
33 enum pressureMagnitude
34 {
35         around1bar300KTransport,/**< pressure is around 1 bar and temperature around 300K (for TransportEquation, SinglePhase and IsothermalTwoFluid) or 373 K (saturation for DriftModel and FiveEqsTwoFluid) */
36         around155bars600KTransport/**< pressure is around 155 bars  and temperature around 618 K (saturation) */
37 };
38
39 //! enumeration BoundaryType
40 /*! Boundary condition type  */
41 enum BoundaryTypeTransport      {InletTransport,  OutletTransport, NeumannTransport, DirichletTransport, NoneBCTransport};//Actually Inlet=Dirichlet and Outlet=Neumann
42
43 /** \struct LimitField
44  * \brief value of some fields on the boundary  */
45 struct LimitFieldTransport{
46         LimitFieldTransport(){bcType=NoneBCTransport; T=0; h=0; flux=0; }
47         LimitFieldTransport(BoundaryTypeTransport _bcType, double _T,   double _h,double _flux  ){
48                 bcType=_bcType; T=_T; h=_h; flux=_flux;
49         }
50
51         BoundaryTypeTransport bcType;
52         double T; //for inlet or Dirichlet
53         double h; //for inlet or Dirichlet
54         double flux; //for Neumann or outlet
55 };
56
57 class TransportEquation: public ProblemCoreFlows
58 {
59
60 public :
61         /** \fn TransportEquation
62                          * \brief Constructor for the enthalpy transport in a fluid
63                          * \param [in] phase : \ref Liquid or \ref Gas
64                          * \param [in] pressureMagnitude : \ref around1bar or \ref around155bars
65                          * \param [in] vector<double> : fluid velocity (assumed constant)
66                          *  */
67         TransportEquation(phase fluid, pressureMagnitude pEstimate,vector<double> vitesseTransport, MPI_Comm comm = MPI_COMM_WORLD);
68
69         //Gestion du calcul
70         virtual void initialize();
71         virtual void terminate();//vide la mémoire et enregistre le résultat final
72         bool initTimeStep(double dt);
73         double computeTimeStep(bool & stop);//propose un pas de temps pour le calcul. Celà nécessite de discrétiser les opérateur (convection, diffusion, sources) et pour chacun d'employer la condition cfl. En cas de problème durant ce calcul (exemple t=tmax), renvoie stop=true
74         void abortTimeStep();//efface les inconnues calculées par solveTimeStep() et reinitialise dt à 0
75         bool iterateTimeStep(bool &ok);
76         virtual void save();
77         virtual void validateTimeStep();
78
79         /* Boundary conditions */
80         /** \fn setIntletBoundaryCondition
81                          * \brief adds a new boundary condition of type Inlet
82                          * \details
83                          * \param [in] string : the name of the boundary
84                          * \param [in] double : the value of the temperature at the boundary
85                          * \param [out] void
86                          *  */
87         void setInletBoundaryCondition(string groupName,double enthalpy){
88                 _limitField[groupName]=LimitFieldTransport(InletTransport,-1,enthalpy,-1);
89         };
90
91         /** \fn setNeumannBoundaryCondition
92          * \brief adds a new boundary condition of type Neumann
93          * \details
94          * \param [in] string the name of the boundary
95          * \param [out] void
96          *  */
97         void setNeumannBoundaryCondition(string groupName, double flux=0){
98                 _limitField[groupName]=LimitFieldTransport(NeumannTransport,-1,flux,-1);
99         };
100
101         /** \fn setBoundaryFields
102          * \brief met à jour  _limitField  ( le type de condition limite )
103          * \details
104          * \param [in] string
105          * \param [out] void
106          *  */
107         void setBoundaryFields(map<string, LimitFieldTransport> boundaryFields){
108                 _limitField = boundaryFields;
109         };
110
111
112         /*Physical parameters*/
113         void setLiqSatEnthalpy(double hsatl){
114                 _hsatl=hsatl;
115         };
116         void setVapSatEnthalpy(double hsatv){
117                 _hsatv=hsatv;
118         };
119         void setLiqSatDensity(double rhosatl){
120                 _rhosatl=rhosatl;
121         };
122         void setVapSatDensity(double rhosatv){
123                 _rhosatv=rhosatv;
124         };
125         void setTransportVelocity(Vector v){
126                 _vitesseTransport=v;
127         };
128
129         /* set input fields to prepare the simulation */
130         vector<string> getInputFieldsNames();
131         void setInputField(const string& nameField, Field& inputField );//supply of a required input field
132         
133         /** \fn setRodTemperatureField
134          * \brief Set the rod temperature field
135          * \details
136          * \param [in] Field
137          * \param [out] void
138          *  */
139         void setRodTemperatureField(Field rodTemperature);
140
141         /** \fn setRodTemperature 
142          * \brief Set a constant rod temperature field
143          * \details
144          * \param [in] double
145          * \param [out] void
146          *  */
147         void setRodTemperature(double rodTemp){
148                 _rodTemperature=rodTemp;
149                 _isStationary=false;//Source term may be changed after previously reaching a stationary state
150         }
151
152         /** \fn getRodTemperatureField
153          * \brief
154          * \details
155          * \param [in] void
156          * \param [out] Field
157          *  */
158         Field& getRodTemperatureField(){ // ?? je ne retrouve pas cet attribut dans le file.cxx
159                 return _rodTemperatureField;
160         }
161
162         /* get output fields for postprocessing or coupling */
163         vector<string> getOutputFieldsNames() ;//liste tous les champs que peut fournir le code pour le postraitement
164         Field&         getOutputField(const string& nameField );//Renvoie un champs pour le postraitement
165
166         Field& getFluidTemperatureField(){
167                 return _TT;
168         }
169
170         Field& getEnthalpyField(){
171                 return _VV;
172         }
173
174         Field& getVoidFractionField(){
175                 return _Alpha;
176         }
177
178         Field& getDensityField(){
179                 return _Rho;
180         }
181
182 protected :
183         double computeTransportMatrix();
184         double computeRHS();
185         void updatePrimitives();
186
187         /* Postprocessing fields */
188         Field   _TT, _Alpha, _Rho;//Fields of temperature, void fraction, density. Unknown field is enthalpy (_VV)
189         double _rhosatv, _rhosatl;
190         double _Tref, _href, _cpref;
191
192         double temperature(double h){
193                 return _Tref+(h-_href)/_cpref;
194         };
195         double voidFraction(double h){
196                 double titre=(h-_href)/(_hsatv-_hsatl);
197                 if (titre<0)
198                         return 0;
199                 else if (titre>1)
200                         return 1;
201                 else return titre*_rhosatl/(titre*_rhosatl+(1-titre)*_rhosatv);
202         };
203         double density(double alpha){
204                 return alpha*_rhosatv+(1-alpha)*_rhosatl;
205         };
206
207         Vector _vitesseTransport, _normale;
208         bool _transportMatrixSet;
209         Vec _Hn, _deltaH, _Hk, _Hkm1, _b0;
210         Vec _Hn_seq; // Local sequential copy of the parallel vector _Hn, used for saving result files
211         double _dt_transport, _dt_src;
212
213         map<string, LimitFieldTransport> _limitField;
214         
215         /* source terms */
216         bool   _rodTemperatureFieldSet;
217         Field  _rodTemperatureField;
218         double _rodTemperature;
219 };
220
221 #endif /* TransportEquation_HXX_ */