Salome HOME
update CoreFlows
[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         void save();
48
49         /** \fn setIntletBoundaryCondition
50          * \brief adds a new boundary condition of type Inlet
51          * \details
52          * \param [in] string : the name of the boundary
53          * \param [in] double : the value of the temperature at the boundary
54          * \param [in] double : the value of the x component of the velocity at the boundary
55          * \param [in] double : the value of the y component of the velocity at the boundary
56          * \param [in] double : the value of the z component of the velocity at the boundary
57          * \param [out] void
58          *  */
59         void setInletBoundaryCondition(string groupName,double Temperature,double v_x=0, double v_y=0, double v_z=0){
60                 _limitField[groupName]=LimitField(Inlet,-1,vector<double>(1,v_x),vector<double>(1,v_y),vector<double>(1,v_z),Temperature,-1,-1,-1);
61         };
62         /** \fn setIntletRotationBoundaryCondition
63          * \brief adds a new boundary condition of type InletRotationVelocity
64          * \details
65          * \param [in] string : the name of the boundary
66          * \param [in] double : the value of the temperature at the boundary
67          * \param [in] double : the components of the rotational of the inlet velocity
68          * \param [out] void
69          *  */
70         void setInletRotationBoundaryCondition(string groupName,double Temperature,double omega_x=0, double omega_y=0, double omega_z=0){
71                 _limitField[groupName]=LimitField(InletRotationVelocity,0,vector<double>(1,omega_x),vector<double>(1,omega_y),vector<double>(1,omega_z),Temperature,-1,-1,-1);
72         };
73         /** \fn setIntletPressureBoundaryCondition
74          * \brief adds a new boundary condition of type InletPressure
75          * \details
76          * \param [in] string : the name of the boundary
77          * \param [in] double : the value of the pressure at the boundary
78          * \param [in] double : the value of the temperature at the boundary
79          * \param [in] double : the components of the rotational of the inlet velocity
80          * \param [out] void
81          *  */
82         void setInletPressureBoundaryCondition(string groupName, double pressure,double Temperature,double omega_x=0, double omega_y=0, double omega_z=0){
83                 _limitField[groupName]=LimitField(InletPressure,pressure,vector<double>(1,omega_x),vector<double>(1,omega_y),vector<double>(1,omega_z),Temperature,-1,-1,-1);
84         };
85         /** \fn setIntletPressureBoundaryCondition
86          * \brief adds a new boundary condition of type InletPressure taking into account the hydrostatic pressure variations
87          * \details The pressure is not constant on the boundary but varies linearly with a slope given by the gravity vector
88          * \param [in] string : the name of the boundary
89          * \param [in] double : the value of the pressure at the boundary
90          * \param [in] double : the value of the temperature at the boundary
91          * \param [in] vector<double> : reference_point position on the boundary where the value Pressure will be imposed
92          * \param [out] void
93          *  */
94         void setInletPressureBoundaryCondition(string groupName, double pressure,double Temperature, vector<double> reference_point){
95                 /* On the boundary we have P-Pref=rho g\cdot(x-xref) hence P=Pref-g\cdot xref + g\cdot x */
96                 pressure-=reference_point[0]*_GravityField3d[0];
97                 if(_Ndim>1){
98                         pressure-=reference_point[1]*_GravityField3d[1];
99                         if(_Ndim>2)
100                                 pressure-=reference_point[2]*_GravityField3d[2];
101                 }
102
103                 _limitField[groupName]=LimitField(InletPressure,pressure,vector<double>(0,0),vector<double>(0,0),vector<double>(0,0),Temperature,-1,-1,-1);
104         };
105         /** \fn setWallBoundaryCondition
106          * \brief adds a new boundary condition of type Wall
107          * \details
108          * \param [in] string : the name of the boundary
109          * \param [in] double : the value of the temperature at the boundary
110          * \param [in] double : the value of the x component of the velocity at the boundary
111          * \param [in] double : the value of the y component of the velocity at the boundary
112          * \param [in] double : the value of the z component of the velocity at the boundary
113          * \param [out] void
114          *  */
115         void setWallBoundaryCondition(string groupName,double Temperature=-1, double v_x=0, double v_y=0, double v_z=0){
116                 if(Temperature!=-1)
117                         _limitField[groupName]=LimitField(Wall,-1,vector<double>(1,v_x),vector<double>(1,v_y),vector<double>(1,v_z),Temperature,-1,-1,-1);
118                 else
119                         _limitField[groupName]=LimitField(Wall,-1,vector<double>(1,v_x),vector<double>(1,v_y),vector<double>(1,v_z),getReferenceTemperature(),-1,-1,-1);
120         };
121
122         /** \fn computeNewtonVariation
123          * \brief Builds and solves the linear system to obtain the variation Vkp1-Vk in a Newton scheme using primitive variables
124          * @param
125          * */
126         void computeNewtonVariation();
127
128         /** \fn iterateTimeStep
129          * \brief calls computeNewtonVariation to perform one Newton iteration and tests the convergence
130          * @param
131          * @return boolean ok is true is the newton iteration gave a physically acceptable result
132          * */
133         bool iterateTimeStep(bool &ok);
134
135         double getReferencePressure()    { return _Pref; };
136         double getReferenceTemperature() { return _Tref; };
137         
138         //get output fields for postprocessing or coupling
139         vector<string> getOutputFieldsNames() ;//liste tous les champs que peut fournir le code pour le postraitement
140         Field&         getOutputField(const string& nameField );//Renvoie un champs pour le postraitement
141         Field& getPressureField();
142         Field& getVelocityField();
143         Field& getVelocityXField();
144         Field& getTemperatureField();
145         Field& getDensityField();
146         Field& getMomentumField();
147         Field& getTotalEnergyField();
148         Field& getEnthalpyField();
149
150 protected :
151         double  _drho_sur_dp,   _drho_sur_dT;//derivatives of the density rho wrt cv, p, T
152         double  _drhoE_sur_dp,  _drhoE_sur_dT;//derivatives of the total energy rho E wrt cv, p, T
153         bool _useDellacherieEOS;
154         double _Tref; //EOS reference temperature
155         double _Pref; //EOS reference pressure
156
157         //!calcule l'etat de Roe de deux etats
158         void convectionState( const long &i, const long &j, const bool &IsBord);
159         //!calcule la matrice de convection de l'etat interfacial entre deux cellules voisinnes
160         void convectionMatrices();
161         //!Calcule le flux pour un état et une porosité et une normale donnés
162         Vector convectionFlux(Vector U,Vector V, Vector normale, double porosity);
163         //!calcule la matrice de diffusion de l'etat interface pour la diffusion
164         void diffusionStateAndMatrices(const long &i,const long &j, const bool &IsBord);
165         //!Computes the source vector associated to the cell i
166         void sourceVector(PetscScalar * Si,PetscScalar * Ui,PetscScalar * Vi, int i);
167         //!Computes the pressure loss associated to the face ij
168         void pressureLossVector(PetscScalar * pressureLoss, double K, PetscScalar * Ui, PetscScalar * Vi, PetscScalar * Uj, PetscScalar * Vj);
169         //!Computes the contribution of the porosity gradient associated to the face ij to the source term
170         void porosityGradientSourceVector();
171         //!Calcule la jacobienne de la CL convection
172         void jacobian(const int &j, string nameOfGroup,double * normale);
173         //!Calcule la jacobienne de la CL de diffusion
174         void jacobianDiff(const int &j, string nameOfGroup);
175         //!Calcule l'etat fictif a la frontiere
176         void setBoundaryState(string nameOfGroup, const int &j,double *normale);// delete &nf Kieu
177         //!Adds the contribution of diffusion to the RHS
178         void addDiffusionToSecondMember(const int &i,const int &j,bool isBord);
179         //!Computation of the Roe matrix
180         void RoeMatrixConservativeVariables(double u_n, double total_enthalpy,Vector velocity, double k, double K);
181         void convectionMatrixPrimitiveVariables( double rho, double u_n, double H,Vector velocity);
182         //!Computation of the staggered Roe upwinding matrix in conservative variables
183         void staggeredRoeUpwindingMatrixConservativeVariables(  double u_n, double total_enthalpy, Vector velocity, double k, double K);
184         //!Computation of the staggered Roe upwinding matrix in primitive variables
185         void staggeredRoeUpwindingMatrixPrimitiveVariables(double density, double u_n,double total_enthalpy, Vector velocity);
186         /**** staggered VFFC scheme has some specific features (no need for Roe matrix), hence some specific functions ******/
187         //!Computes the interfacial flux for the VFFC formulation of the staggered upwinding
188         Vector staggeredVFFCFlux();
189         //!Computes the matrices A^+ and A^- for the VFFC formulation of the staggered upwinding
190         void staggeredVFFCMatricesConservativeVariables(double u_n);
191         //!Computes the matrices A^+ and A^- for the VFFC formulation of the staggered upwinding using Primitive Variables
192         void staggeredVFFCMatricesPrimitiveVariables(double u_n);
193         //!Compute the corrected interfacial state for lowMach, pressureCorrection and staggered versions of the VFRoe formulation
194         void applyVFRoeLowMachCorrections(bool isBord, string groupname="");
195         //!Special preconditioner based on a matrix scaling strategy
196         void computeScaling(double offset);
197         //!Calcule les saut de valeurs propres pour la correction entropique
198         void entropicShift(double* n);
199         // Fonctions utilisant la loi d'etat
200         void consToPrim(const double *Ucons, double* Vprim,double porosity=1);
201         void primToCons(const double *V, const int &i, double *U, const int &j);
202         void primToConsJacobianMatrix(double *V);
203         /** \fn getDensityDerivatives
204          * \brief Computes the partial derivatives of rho, and rho E with regard to the primitive variables  p and  T
205          * @param pressure
206          * @param temperature
207          * @param square of the velocity vector
208         */
209         void getDensityDerivatives( double pressure, double temperature, double v2);
210
211         bool _saveAllFields;
212         Field _Enthalpy, _Pressure, _Density, _Temperature, _Momentum, _TotalEnergy, _Vitesse, _VitesseX, _VitesseY, _VitesseZ;
213
214         };
215 #endif /* SINGLEPHASE_HXX_*/