Salome HOME
Added a function to save all fields in signe phase flows
[tools/solverlab.git] / CoreFlows / Models / inc / SinglePhase.hxx
1 //============================================================================
2 /**
3  * \file DriftModel.hxx
4  * \author Michael NDJINGA, Kieu Nguyen
5  * \version 1.0
6  * \date 01 janv. 2016
7  * \brief The compressible Navier-Stokes equations
8  * */
9 //============================================================================
10
11 /*! \class SinglePhase SinglePhase.hxx "SinglePhase.hxx"
12  *  \brief The compressible Navier-Stokes equations
13  *  \details The model consists in one mass, one momentum and one energy equation, see \ref NSModelsPage for more details
14  */
15 #ifndef SINGLEPHASE_HXX_
16 #define SINGLEPHASE_HXX_
17
18 #include "ProblemFluid.hxx"
19
20 class SinglePhase : public ProblemFluid{
21 public :
22         /** \fn SinglePhase
23          * \brief Constructor for the Navier-Stokes system
24          * \param [in] phaseType : \ref Liquid or \ref Gas
25          * \param [in] pressureEstimate : \ref around1bar or \ref around155bars
26          * \param [in] int : mesh dimension
27          * \param [in] bool : There are two possible equations of state for the fluid
28          *  */
29         SinglePhase(phaseType fluid, pressureEstimate pEstimate,int dim,bool useDellacherieEOS=false);
30
31         /** \fn setViscosity
32          * \brief sets the viscosity
33          * @param viscosite : value of the dynamic viscosity
34          *       * */
35         void setViscosityConstant( double viscosite ){
36                 _fluides[0]->setViscosity(viscosite);
37         };
38
39         //! system initialisation
40         void initialize();
41
42         //fonctions d'echange de flux
43         //      void getOutputField(const Vec &Flux, const string Champ, const int numBord)=0;//, PetscInt *indices_Flux, PetscInt *indices_Bord, const long range)=0;
44         //      double trace(const int &numBord, Vec &out)=0;
45         void testConservation();
46
47         /** \fn saveAllFields
48          * \brief saves every interesting field in a separate file
49          * @param boolean saveAllFields
50          * */
51         void saveAllFields(bool saveAllFields=true){
52                 _saveAllFields=saveAllFields;
53         }
54
55         void save();
56
57         /** \fn setIntletBoundaryCondition
58          * \brief adds a new boundary condition of type Inlet
59          * \details
60          * \param [in] string : the name of the boundary
61          * \param [in] double : the value of the temperature at the boundary
62          * \param [in] double : the value of the x component of the velocity at the boundary
63          * \param [in] double : the value of the y component of the velocity at the boundary
64          * \param [in] double : the value of the z component of the velocity at the boundary
65          * \param [out] void
66          *  */
67         void setInletBoundaryCondition(string groupName,double Temperature,double v_x=0, double v_y=0, double v_z=0){
68                 _limitField[groupName]=LimitField(Inlet,-1,vector<double>(1,v_x),vector<double>(1,v_y),vector<double>(1,v_z),Temperature,-1,-1,-1);
69         };
70         /** \fn setIntletRotationBoundaryCondition
71          * \brief adds a new boundary condition of type InletRotationVelocity
72          * \details
73          * \param [in] string : the name of the boundary
74          * \param [in] double : the value of the temperature at the boundary
75          * \param [in] double : the components of the rotational of the inlet velocity
76          * \param [out] void
77          *  */
78         void setInletRotationBoundaryCondition(string groupName,double Temperature,double omega_x=0, double omega_y=0, double omega_z=0){
79                 _limitField[groupName]=LimitField(InletRotationVelocity,0,vector<double>(1,omega_x),vector<double>(1,omega_y),vector<double>(1,omega_z),Temperature,-1,-1,-1);
80         };
81         /** \fn setIntletPressureBoundaryCondition
82          * \brief adds a new boundary condition of type InletPressure
83          * \details
84          * \param [in] string : the name of the boundary
85          * \param [in] double : the value of the pressure at the boundary
86          * \param [in] double : the value of the temperature at the boundary
87          * \param [in] double : the components of the rotational of the inlet velocity
88          * \param [out] void
89          *  */
90         void setInletPressureBoundaryCondition(string groupName, double pressure,double Temperature,double omega_x=0, double omega_y=0, double omega_z=0){
91                 _limitField[groupName]=LimitField(InletPressure,pressure,vector<double>(1,omega_x),vector<double>(1,omega_y),vector<double>(1,omega_z),Temperature,-1,-1,-1);
92         };
93         /** \fn setIntletPressureBoundaryCondition
94          * \brief adds a new boundary condition of type InletPressure taking into account the hydrostatic pressure variations
95          * \details The pressure is not constant on the boundary but varies linearly with a slope given by the gravity vector
96          * \param [in] string : the name of the boundary
97          * \param [in] double : the value of the pressure at the boundary
98          * \param [in] double : the value of the temperature at the boundary
99          * \param [in] vector<double> : reference_point position on the boundary where the value Pressure will be imposed
100          * \param [out] void
101          *  */
102         void setInletPressureBoundaryCondition(string groupName, double pressure,double Temperature, vector<double> reference_point){
103                 /* On the boundary we have P-Pref=rho g\cdot(x-xref) hence P=Pref-g\cdot xref + g\cdot x */
104                 pressure-=reference_point[0]*_GravityField3d[0];
105                 if(_Ndim>1){
106                         pressure-=reference_point[1]*_GravityField3d[1];
107                         if(_Ndim>2)
108                                 pressure-=reference_point[2]*_GravityField3d[2];
109                 }
110
111                 _limitField[groupName]=LimitField(InletPressure,pressure,vector<double>(0,0),vector<double>(0,0),vector<double>(0,0),Temperature,-1,-1,-1);
112         };
113         /** \fn setWallBoundaryCondition
114          * \brief adds a new boundary condition of type Wall
115          * \details
116          * \param [in] string : the name of the boundary
117          * \param [in] double : the value of the temperature at the boundary
118          * \param [in] double : the value of the x component of the velocity at the boundary
119          * \param [in] double : the value of the y component of the velocity at the boundary
120          * \param [in] double : the value of the z component of the velocity at the boundary
121          * \param [out] void
122          *  */
123         void setWallBoundaryCondition(string groupName,double Temperature=-1, double v_x=0, double v_y=0, double v_z=0){
124                 if(Temperature!=-1)
125                         _limitField[groupName]=LimitField(Wall,-1,vector<double>(1,v_x),vector<double>(1,v_y),vector<double>(1,v_z),Temperature,-1,-1,-1);
126                 else
127                         _limitField[groupName]=LimitField(Wall,-1,vector<double>(1,v_x),vector<double>(1,v_y),vector<double>(1,v_z),getReferenceTemperature(),-1,-1,-1);
128         };
129
130         /** \fn computeNewtonVariation
131          * \brief Builds and solves the linear system to obtain the variation Vkp1-Vk in a Newton scheme using primitive variables
132          * @param
133          * */
134         void computeNewtonVariation();
135
136         /** \fn iterateTimeStep
137          * \brief calls computeNewtonVariation to perform one Newton iteration and tests the convergence
138          * @param
139          * @return boolean ok is true is the newton iteration gave a physically acceptable result
140          * */
141         bool iterateTimeStep(bool &ok);
142
143         double getReferencePressure()    { return _Pref; };
144         double getReferenceTemperature() { return _Tref; };
145         
146         //get output fields for postprocessing or coupling
147         vector<string> getOutputFieldsNames() ;//liste tous les champs que peut fournir le code pour le postraitement
148         Field&         getOutputField(const string& nameField );//Renvoie un champs pour le postraitement
149         Field& getPressureField();
150         Field& getVelocityField();
151         Field& getVelocityXField();
152         Field& getTemperatureField();
153         Field& getDensityField();
154         Field& getMomentumField();
155         Field& getTotalEnergyField();
156         Field& getEnthalpyField();
157
158 protected :
159         double  _drho_sur_dp,   _drho_sur_dT;//derivatives of the density rho wrt cv, p, T
160         double  _drhoE_sur_dp,  _drhoE_sur_dT;//derivatives of the total energy rho E wrt cv, p, T
161         bool _useDellacherieEOS;
162         double _Tref; //EOS reference temperature
163         double _Pref; //EOS reference pressure
164
165         //!calcule l'etat de Roe de deux etats
166         void convectionState( const long &i, const long &j, const bool &IsBord);
167         //!calcule la matrice de convection de l'etat interfacial entre deux cellules voisinnes
168         void convectionMatrices();
169         //!Calcule le flux pour un état et une porosité et une normale donnés
170         Vector convectionFlux(Vector U,Vector V, Vector normale, double porosity);
171         //!calcule la matrice de diffusion de l'etat interface pour la diffusion
172         void diffusionStateAndMatrices(const long &i,const long &j, const bool &IsBord);
173         //!Computes the source vector associated to the cell i
174         void sourceVector(PetscScalar * Si,PetscScalar * Ui,PetscScalar * Vi, int i);
175         //!Computes the pressure loss associated to the face ij
176         void pressureLossVector(PetscScalar * pressureLoss, double K, PetscScalar * Ui, PetscScalar * Vi, PetscScalar * Uj, PetscScalar * Vj);
177         //!Computes the contribution of the porosity gradient associated to the face ij to the source term
178         void porosityGradientSourceVector();
179         //!Calcule la jacobienne de la CL convection
180         void jacobian(const int &j, string nameOfGroup,double * normale);
181         //!Calcule la jacobienne de la CL de diffusion
182         void jacobianDiff(const int &j, string nameOfGroup);
183         //!Calcule l'etat fictif a la frontiere
184         void setBoundaryState(string nameOfGroup, const int &j,double *normale);// delete &nf Kieu
185         //!Adds the contribution of diffusion to the RHS
186         void addDiffusionToSecondMember(const int &i,const int &j,bool isBord);
187         //!Computation of the Roe matrix
188         void RoeMatrixConservativeVariables(double u_n, double total_enthalpy,Vector velocity, double k, double K);
189         void convectionMatrixPrimitiveVariables( double rho, double u_n, double H,Vector velocity);
190         //!Computation of the staggered Roe upwinding matrix in conservative variables
191         void staggeredRoeUpwindingMatrixConservativeVariables(  double u_n, double total_enthalpy, Vector velocity, double k, double K);
192         //!Computation of the staggered Roe upwinding matrix in primitive variables
193         void staggeredRoeUpwindingMatrixPrimitiveVariables(double density, double u_n,double total_enthalpy, Vector velocity);
194         /**** staggered VFFC scheme has some specific features (no need for Roe matrix), hence some specific functions ******/
195         //!Computes the interfacial flux for the VFFC formulation of the staggered upwinding
196         Vector staggeredVFFCFlux();
197         //!Computes the matrices A^+ and A^- for the VFFC formulation of the staggered upwinding
198         void staggeredVFFCMatricesConservativeVariables(double u_n);
199         //!Computes the matrices A^+ and A^- for the VFFC formulation of the staggered upwinding using Primitive Variables
200         void staggeredVFFCMatricesPrimitiveVariables(double u_n);
201         //!Compute the corrected interfacial state for lowMach, pressureCorrection and staggered versions of the VFRoe formulation
202         void applyVFRoeLowMachCorrections(bool isBord, string groupname="");
203         //!Special preconditioner based on a matrix scaling strategy
204         void computeScaling(double offset);
205         //!Calcule les saut de valeurs propres pour la correction entropique
206         void entropicShift(double* n);
207         // Fonctions utilisant la loi d'etat
208         void consToPrim(const double *Ucons, double* Vprim,double porosity=1);
209         void primToCons(const double *V, const int &i, double *U, const int &j);
210         void primToConsJacobianMatrix(double *V);
211         /** \fn getDensityDerivatives
212          * \brief Computes the partial derivatives of rho, and rho E with regard to the primitive variables  p and  T
213          * @param pressure
214          * @param temperature
215          * @param square of the velocity vector
216         */
217         void getDensityDerivatives( double pressure, double temperature, double v2);
218
219         bool _saveAllFields;
220         Field _Enthalpy, _Pressure, _Density, _Temperature, _Momentum, _TotalEnergy, _Vitesse, _VitesseX, _VitesseY, _VitesseZ;
221
222         };
223 #endif /* SINGLEPHASE_HXX_*/