Salome HOME
Changed prototype of function setLinearSolver
[tools/solverlab.git] / CoreFlows / Models / inc / ProblemFluid.hxx
1 //============================================================================
2 // Name        : ProblemCoreFlows
3 // Author      : M. Ndjinga
4 // Version     :
5 // Copyright   : CEA Saclay 2014
6 // Description : Generic class for thermal hydraulics problems
7 //============================================================================
8
9 /*! \class ProblemFluid ProblemFluid.hxx "ProblemFluid.hxx"
10  *  \brief Factorises the methods that are common to the non scalar models (fluid models)
11  *  \details Common functions to fluid models
12  */
13 #ifndef PROBLEMFLUID_HXX_
14 #define PROBLEMFLUID_HXX_
15
16 #include "Fluide.h"
17 #include "ProblemCoreFlows.hxx"
18 #include "utilitaire_algebre.h"
19
20 using namespace std;
21
22 //! enumeration SpaceScheme
23 /*! Several numerical schemes are available */
24 enum SpaceScheme
25 {
26         upwind,/**<  classical full upwinding scheme (first order in space) */
27         centered,/**<  centered scheme (second order in space) */
28         pressureCorrection,/**<  include a pressure correction in the upwind scheme to increase precision at low Mach numbers */
29         lowMach,/**<  include an upwinding proportional to the Mach numer scheme to increase precision at low Mach numbers */
30         staggered,/**<  scheme inspired by staggered discretisations */
31 };
32
33 //! enumeration pressureEstimate
34 /*! the pressure estimate needed to fit physical parameters  */
35 enum pressureEstimate
36 {
37         around1bar300K,/**< pressure is around 1 bar and temperature around 300K (for TransportEquation, SinglePhase and IsothermalTwoFluid) or 373 K (saturation for DriftModel and FiveEqsTwoFluid) */
38         around155bars600K/**< pressure is around 155 bars  and temperature around 618 K (saturation) */
39 };
40
41 //! enumeration phaseType
42 /*! The fluid type can be Gas or water  */
43 enum phaseType
44 {
45         Liquid,/**< Fluid considered is water */
46         Gas/**< Fluid considered is Gas */
47 };
48
49 //! enumeration NonLinearFormulation
50 /*! the formulation used to compute the non viscous fluxes */
51 enum NonLinearFormulation
52 {
53         Roe,/**< Ph. Roe non linear formulation is used */
54         VFRoe,/**< Masella, Faille and Gallouet non linear formulation is used */
55         VFFC,/**< Ghidaglia, Kumbaro and Le Coq non linear formulation is used */
56         reducedRoe,/**< compacted formulation of Roe scheme without computation of the fluxes */
57 };
58
59 //! enumeration BoundaryType
60 /*! Boundary condition type  */
61 enum BoundaryType       {Wall, InnerWall, Inlet, InletPressure, InletRotationVelocity, InletEnthalpy, Outlet, Neumann, NoTypeSpecified};
62 /** \struct LimitField
63  * \brief value of some fields on the boundary  */
64 struct LimitField{
65         LimitField(){bcType=NoTypeSpecified; p=0; v_x=vector<double> (0,0); v_y=vector<double> (0,0); v_z=vector<double> (0,0); T=0; h=0; alpha=0;      conc=0;}
66         LimitField(BoundaryType _bcType,        double _p,      vector<double> _v_x, vector<double> _v_y, vector<double> _v_z,
67                         double _T,      double _h,      double _alpha,  double _conc){
68                 bcType=_bcType; p=_p; v_x=_v_x; v_y=_v_y; v_z=_v_z;     T=_T; h=_h; alpha=_alpha;       conc=_conc;
69         }
70
71         BoundaryType bcType;
72         double p;//For outlet (fluid models)
73         vector<double> v_x; vector<double> v_y; vector<double> v_z;//For wall and inlet (fluid models)
74         double T; //for wall and inlet (DriftModel and FiveEqsTwoFluid) and for Dirichlet (DiffusionEquation)
75         double h; //for inlet (TransportEquation)
76         double alpha; //For inlet (IsothermalTwoFluid and FiveEqsTwoFluid)
77         double conc;//For inlet (DriftModel)
78 };
79
80 class ProblemFluid: public ProblemCoreFlows
81 {
82
83 public :
84         /**\fn
85          * \brief constructeur de la classe ProblemFluid
86          */
87         ProblemFluid(void);
88
89         //Gestion du calcul (interface ICoCo)
90
91         /** \fn initialize
92          * \brief Alocates memory and checks that the mesh, the boundary and the intitial data are set
93          * \Details It is a pure virtual function overloaded y each model.
94          * @param  void
95          *  */
96         virtual void initialize();
97
98         /** \fn terminate
99          * \brief empties the memory
100          * @param void
101          *  */
102         virtual void terminate();
103
104         /** \fn initTimeStep
105          * \brief Sets a new time step dt to be solved later
106          *  @param  dt is  the value of the time step
107          *  \return false if dt <0 et True otherwise
108          *                         */
109         bool initTimeStep(double dt);
110
111         /** \fn computeTimeStep
112          * \brief Proposes a value for the next time step to be solved using mesh data and cfl coefficient
113          *  \return  double dt the proposed time step
114          *  \return  bool stop, true if the calculation should not be continued (stationary state, maximum time or time step numer reached)
115          *  */
116         double computeTimeStep(bool & stop);
117
118         /** \fn abortTimeStep
119          * \brief Reset the time step dt to 0
120          *  */
121         void abortTimeStep();
122
123         /** \fn iterateTimeStep
124          * \brief calls computeNewtonVariation to perform one Newton iteration and tests the convergence
125          * @param
126          * @return boolean ok is true is the newton iteration gave a physically acceptable result
127          * */
128         virtual bool iterateTimeStep(bool &ok);
129
130         /** \fn save
131          * \brief saves the current results in MED or VTK files
132          * \details It is a pure virtual function overloaded in each model
133          * @param  void
134          */
135         virtual void save()=0;
136
137         /** \fn validateTimeStep
138          * \brief Validates the solution computed y solveTimeStep
139          * \details updates the currens time t=t+dt, save unknown fields, resets the time step dt to 0, tests the stationnarity.
140          * c It is a pure virtual function overloaded in each model
141          * @param  void
142          *  */
143         virtual void validateTimeStep();
144
145         /* Boundary conditions */
146         /** \fn setNeumannBoundaryCondition
147          * \brief adds a new boundary condition of type Neumann
148          * \details
149          * \param [in] string the name of the boundary
150          * \param [out] void
151          *  */
152         void setNeumannBoundaryCondition(string groupName){
153                 _limitField[groupName]=LimitField(Neumann,-1,vector<double>(_Ndim,0),vector<double>(_Ndim,0),vector<double>(_Ndim,0),-1,-1,-1,-1);
154         };
155
156         /** \fn setOutletBoundaryCondition
157          * \brief Adds a new boundary condition of type Outlet
158          * \details
159          * \param [in] string : the name of the boundary
160          * \param [in] double : the value of the pressure at the boundary
161          * \param [out] void
162          *  */
163         void setOutletBoundaryCondition(string groupName,double Pressure){
164                 _limitField[groupName]=LimitField(Outlet,Pressure,vector<double>(_nbPhases,0),vector<double>(_nbPhases,0),vector<double>(_nbPhases,0),-1,-1,-1,-1);
165         };
166
167         /** \fn setOutletBoundaryCondition
168          * \brief Adds a new boundary condition of type Outlet taking into account the hydrostatic pressure variations
169          * \details The pressure is not constant on the boundary but varies linearly with a slope given by the gravity vector
170          * \param [in] string : the name of the boundary
171          * \param [in] double : the value of the pressure at the boundary
172          * \param [in] vector<double> : reference_point position on the boundary where the value Pressure will be imposed
173          * \param [out] void
174          *  */
175         void setOutletBoundaryCondition(string groupName,double referencePressure, vector<double> reference_point){
176                 /* On the boundary we have P-Pref=rho g\cdot(x-xref) */
177                 _gravityReferencePoint=reference_point;
178                 _limitField[groupName]=LimitField(Outlet,referencePressure,vector<double>(_nbPhases,0),vector<double>(_nbPhases,0),vector<double>(_nbPhases,0),-1,-1,-1,-1);
179         };
180
181         /** \fn setBoundaryFields
182          * \brief met à jour  _limitField  ( le type de condition limite )
183          * \details
184          * \param [in] string
185          * \param [out] void
186          *  */
187         void setBoundaryFields(map<string, LimitField> boundaryFields){
188                 _limitField = boundaryFields;
189         };
190
191         /** \fn setViscosity
192          * \brief sets the vector of viscosity coefficients
193          * @param viscosite is a vector of size equal to the number of phases and containing the viscosity of each phase
194          * @return throws an exception if the input vector size is not equal to the number of phases
195          *       * */
196         void setViscosity(vector<double> viscosite){
197                 if(_nbPhases!= viscosite.size())
198                         throw CdmathException("ProblemFluid::setViscosity: incorrect vector size vs number of phases");
199                 for(int i=0;i<_nbPhases;i++)
200                         _fluides[i]->setViscosity(viscosite[i]);
201         };
202
203         /** \fn setConductivity
204          * \brief sets the vector of conductivity coefficients
205          * @param conductivite is a vector of size equal to the number of phases and containing the conductivity of each phase
206          * @return throws an exception if the input vector size is not equal to the number of phases
207          * */
208         void setConductivity(vector<double> conductivite){
209                 if(_nbPhases!= conductivite.size())
210                         throw CdmathException("ProblemFluid::setConductivity: incorrect vector size vs number of phases");
211                 for(int i=0;i<_nbPhases;i++)
212                         _fluides[i]->setConductivity(conductivite[i]);
213         };
214
215         /** \fn setGravity
216          * \brief sets the gravity force in the model
217          * @param gravite is a vector of size equal to the space dimension
218          *                               * */
219         void setGravity(vector<double> gravite){
220                 _GravityField3d = gravite;
221         };
222
223         /** \fn setDragCoeffs
224          * \brief Sets the drag coefficients
225          * @param dragCoeffs is a  vector of size equal to the number of phases and containing the value of the friction coefficient of each phase
226          * @return throws an exception if the input vector size is not equal to the numer of phases
227          * */
228         void setDragCoeffs(vector<double> dragCoeffs){
229                 if(_nbPhases!= dragCoeffs.size())
230                         throw CdmathException("ProblemFluid::setDragCoeffs: incorrect vector size vs number of phases");
231                 for(int i=0;i<_nbPhases;i++)
232                 {
233                         _fluides[i]->setDragCoeffs(dragCoeffs[i]);
234                         _dragCoeffs[i]=dragCoeffs[i];
235                 }
236         };
237
238         /** \fn getNumberOfPhases
239          * \brief The numer of phase (one or two) depending on the model considered
240          * @param void
241          * @return the number of phases considered in the model
242          * */
243         int getNumberOfPhases() const {
244                 return _nbPhases;
245         };
246
247         /** \fn computeNewtonVariation
248          * \brief Builds and solves the linear system to obtain the variation Ukp1-Uk in a Newton scheme
249          * @param void
250          * */
251         virtual void computeNewtonVariation();
252
253         /** \fn testConservation
254          * \brief Teste et affiche la conservation de masse et de la quantité de mouvement
255          * \Details la fonction est virtuelle pure, on la surcharge dans chacun des modèles
256          * @param void
257          * */
258         virtual void testConservation()=0;
259
260         /** \fn saveVelocity
261          * \brief saves the velocity field in a separate 3D file so that paraview can display the streamlines
262          * @param bool
263          * */
264         void saveVelocity(bool save_v=true){
265                 _saveVelocity=save_v;
266         }
267
268         /** \fn saveConservativeField
269          * \brief saves the conservative fields (density, momentum etc...)
270          * */
271         void saveConservativeField(bool save=true){
272                 _saveConservativeField=save;
273         }
274         /** \fn setEntropicCorrection
275          * \brief include an entropy correction to avoid non entropic solutions
276          * @param boolean that is true if entropy correction should be applied
277          * */
278         void setEntropicCorrection(bool entropyCorr){
279                 _entropicCorrection=entropyCorr;
280         }
281
282         /** \fn setPressureCorrectionOrder
283          * \brief In case a pressure correction scheme is set by a call to setNumericalScheme(pressureCorrection) this function allows the setting of the type of pressure correction to be used
284          * @param int the order of the pressure correction
285         * \details The following treatment is applied depending on the value of the input parameter order
286         * \details 1 -> no pressure correction (pressure correction is applied nowhere in the domain), standard upwind scheme instead is used
287         * \details 2 -> standard pressure correction is applied everywhere in the domain, even at the boundary faces
288         * \details 3 -> standard pressure correction is applied only inside the domain (not at the boundary faces)
289         * \details 4 -> no pressure correction (pressure correction is applied nowhere in the domain), no Riemann problem at wall boundaries (boundary pressure = inner pressure)
290         * \details 5 -> standard pressure correction is applied everywhere in the domain, no Riemann problem at the boundary (boundary pressure = inner pressure)
291         * \details 6 -> standard pressure correction is applied inside the domain and a special pressure correction involving gravity is applied at the boundary, no Riemann problem at wall boundaries (boundary pressure = inner pressure+ source term)
292         * */
293         void setPressureCorrectionOrder(int order){
294                 if( order >6 || order <1)
295                         throw CdmathException("ProblemFluid::setPressureCorrectionOrder Pressure correction order must be an integer between 1 and 4");
296                 else
297                         _pressureCorrectionOrder=order;
298
299                 if(order==1)
300                         _spaceScheme=upwind;
301         }
302
303         // Petsc resolution
304
305         /** \fn setLinearSolver
306          * \brief sets the linear solver and preconditioner
307          * \details virtuelle function overloaded by intanciable classes
308          * @param kspType linear solver type (GMRES or BICGSTAB)
309          * @param pcType preconditioner (ILU,LU or NONE)
310          * @param scaling performs a bancing of the system matrix before calling th preconditioner
311          */
312         void setLinearSolver(linearSolver kspType, preconditioner pcType, double maxIts=50, bool scaling=false)
313         {
314                 ProblemCoreFlows::setLinearSolver(kspType, pcType, maxIts);
315                 _isScaling= scaling;
316         };
317
318         /** \fn setLatentHeat
319          * \brief Sets the value of the latent heat
320          * @param double L, the value of the latent heat
321          * */
322         void setLatentHeat(double L){
323                 _latentHeat=L;
324         }
325
326         /** \fn getLatentHeat
327          * \brief returns the value of the latent heat
328          * @param double L, the value of the latent heat
329          * */
330         double getLatentHeat() const{
331                 return _latentHeat;
332         }
333
334         /** \fn setSatTemp
335          * \brief sets the saturation temperature
336          * @param  Tsat double corresponds to saturation temperature
337          * */
338         void setSatTemp(double Tsat){
339                 _Tsat=Tsat;
340         }
341
342         /** \fn setSatTemp
343          * \brief sets the saturation temperature
344          * @param  Tsat double corresponds to saturation temperature
345          * */
346         double getSatTemp() const {
347                 return _Tsat;
348         }
349
350         /** \fn setSatPressure
351          * \brief sets the saturation pressure
352          * @param  Psat double corresponds to saturation pressure
353          * */
354         void setSatPressure(double Psat, double dHsatl_over_dp=0.05){
355                 _Psat=Psat;
356                 _dHsatl_over_dp=dHsatl_over_dp;
357         }
358
359         /** \fn setPorosityField
360          * \brief set the porosity field;
361          * @param [in] Field porosity field (field on CELLS)
362          * */
363         void setPorosityField(Field Porosity){
364                 _porosityField=Porosity;
365                 _porosityFieldSet=true;
366         }
367
368         /** \fn getPorosityField
369          * \brief returns the porosity field;
370          * @param
371          * */
372         Field getPorosityField() const {
373                 return _porosityField;
374         }
375
376         /** \fn setPorosityFieldFile
377          * \brief set the porosity field
378          * \details
379          * \param [in] string fileName (including file path)
380          * \param [in] string fieldName
381          * \param [out] void
382          *  */
383         void setPorosityField(string fileName, string fieldName){
384                 _porosityField=Field(fileName, CELLS,fieldName);
385                 _porosityFieldSet=true;
386         }
387
388         /** \fn setPressureLossField
389          * \brief set the pressure loss coefficients field;
390          * @param [in] Field pressure loss field (field on FACES)
391          * */
392         void setPressureLossField(Field PressureLoss){
393                 _pressureLossField=PressureLoss;
394                 _pressureLossFieldSet=true;
395         }
396         /** \fn setPressureLossField
397          * \brief set the pressure loss coefficient field
398          * \details localised friction force
399          * \param [in] string fileName (including file path)
400          * \param [in] string fieldName
401          * \param [out] void
402          *  */
403         void setPressureLossField(string fileName, string fieldName){
404                 _pressureLossField=Field(fileName, FACES,fieldName);
405                 _pressureLossFieldSet=true;
406         }
407
408         /** \fn setSectionField
409          * \brief set the cross section field;
410          * @param [in] Field cross section field (field on CELLS)
411          * */
412         void setSectionField(Field sectionField){
413                 _sectionField=sectionField;
414                 _sectionFieldSet=true;
415         }
416         /** \fn setSectionField
417          * \brief set the cross section field
418          * \details for variable cross section pipe network
419          * \param [in] string fileName (including file path)
420          * \param [in] string fieldName
421          * \param [out] void
422          *  */
423         void setSectionField(string fileName, string fieldName){
424                 _sectionField=Field(fileName, CELLS,fieldName);
425                 _sectionFieldSet=true;
426         }
427
428         /** \fn setNonLinearFormulation
429          * \brief sets the formulation used for the computation of non viscous fluxes
430          * \details Roe, VFRoe, VFFC
431          * \param [in] enum NonLinearFormulation
432          * \param [out] void
433          *  */
434         void setNonLinearFormulation(NonLinearFormulation nonLinearFormulation){
435                 //if(nonLinearFormulation != Roe && nonLinearFormulation != VFRoe && nonLinearFormulation != VFFC && nonLinearFormulation!=reducedRoe)
436                 //      throw CdmathException("nonLinearFormulation should be either Roe, VFRoe or VFFC");//extra security for swig binding
437                 _nonLinearFormulation=nonLinearFormulation;
438         }
439
440         /** \fn getNonLinearFormulation
441          * \brief returns the formulation used for the computation of non viscous fluxes
442          * \details Roe, VFRoe, VFFC
443          * \param [in] void
444          * \param [out] enum NonLinearFormulation
445          *  */
446         NonLinearFormulation getNonLinearFormulation() const{
447                 return _nonLinearFormulation;
448         }
449
450         /** \fn usePrimitiveVarsInNewton
451          * \brief use Primitive Vars instead of conservative vars In Newton scheme for implicit schemes
452          * \param [in] bool
453          *  */
454         void usePrimitiveVarsInNewton(bool usePrimitiveVarsInNewton){
455                 _usePrimitiveVarsInNewton=usePrimitiveVarsInNewton;
456         }
457
458         /** \fn getSpaceScheme
459          * \brief returns the  space scheme name
460          * \param [in] void
461          * \param [out] enum SpaceScheme(upwind, centred, pressureCorrection, pressureCorrection, staggered)
462          *  */
463         SpaceScheme getSpaceScheme() const;
464
465         /** \fn setNumericalScheme
466          * \brief sets the numerical method (upwind vs centered and explicit vs implicit)
467          * \details
468          * \param [in] SpaceScheme
469          * \param [in] TimeScheme
470          * \param [out] void
471          *  */
472         void setNumericalScheme(SpaceScheme scheme, TimeScheme method=Explicit);
473
474         //données initiales
475         /*
476         virtual vector<string> getInputFieldsNames()=0 ;//Renvoie les noms des champs dont le problème a besoin (données initiales)
477         virtual  Field& getInputFieldTemplate(const string& name)=0;//Renvoie le format de champs attendu (maillage, composantes etc)
478         virtual void setInputField(const string& name, const Field& afield)=0;//enregistre les valeurs d'une donnée initiale
479         virtual vector<string> getOutputFieldsNames()=0 ;//liste tous les champs que peut fournir le code pour le postraitement
480         virtual Field& getOutputField(const string& nameField )=0;//Renvoie un champs pour le postraitement
481          */
482
483         Field getConservativeField() const ;
484         Field getPrimitiveField() const;
485
486 protected :
487         /** Number of phases in the fluid **/
488         int _nbPhases;
489         /** Field of conservative variables (the primitive variables are defined in the mother class ProblemCoreFlows **/
490         Field  _UU;
491         /** Field of interfacial states of the VFRoe scheme **/
492         Field _UUstar, _VVstar;
493
494         SpaceScheme _spaceScheme;
495         /** the formulation used to compute the non viscous fluxes **/
496         NonLinearFormulation _nonLinearFormulation;
497
498         map<string, LimitField> _limitField;
499
500         /** boolean used to specify that an entropic correction should be used **/
501         bool _entropicCorrection;
502         /** Vector containing the eigenvalue jumps for the entropic correction **/
503         vector<double> _entropicShift;
504         /** In case a pressure correction scheme is used some more option regarding the type of pressure correction **/
505         int _pressureCorrectionOrder;
506
507         /** Fluid equation of state **/
508         vector< Fluide* > _fluides;
509         //!Viscosity coefficients 
510         vector<double> _viscosite;
511         //!Conductivity coefficients 
512         vector<double> _conductivite;
513
514         /** Source terms **/
515         vector<double> _gravite, _GravityField3d, _gravityReferencePoint, _dragCoeffs;//_GravityField3d has size _Ndim whereas _gravite has size _Nvar and is usefull for dealing with source term and implicitation of gravity vector
516         double _latentHeat, _Tsat,_Psat,_dHsatl_over_dp;
517         Field _porosityField, _pressureLossField, _dp_over_dt, _sectionField;
518         bool _porosityFieldSet, _pressureLossFieldSet, _sectionFieldSet;
519         double _porosityi, _porosityj;//porosity of the left and right states around an interface
520
521         /** User options **/
522         bool _isScaling;
523         bool _saveVelocity;/* In order to display streamlines with paraview */
524         bool _saveConservativeField;/* Save conservative fields such as density or momentum for instance */
525         bool _saveInterfacialField;/* Save interfacial fields of the VFRoe scheme */
526         bool _usePrimitiveVarsInNewton;
527
528         // Variables du schema numerique 
529         Vec _conservativeVars, _newtonVariation, _bScaling,_old, _primitiveVars, _Uext,_Uextdiff ,_vecScaling,_invVecScaling, _Vext;
530         //courant state vector, vector computed at next time step, second member of the equation
531         PetscScalar *_AroePlus, *_AroeMinus,*_Jcb,*_JcbDiff, *_a, *_blockDiag,  *_invBlockDiag,*_Diffusion, *_GravityImplicitationMatrix;
532         PetscScalar *_Aroe, *_absAroe, *_signAroe, *_invAroe;
533         PetscScalar *_AroeMinusImplicit,*_AroePlusImplicit,*_AroeImplicit,*_absAroeImplicit;//negative part of the roe matrix used in the implicit scheme matrix
534         PetscScalar * _primToConsJacoMat; //Jacobian matrix of the prim-->cons function
535         PetscScalar *_phi, *_Ui, *_Uj,  *_Vi, *_Vj,  *_Si, *_Sj, * _pressureLossVector, * _porosityGradientSourceVector, *_externalStates;
536         double *_Uroe, *_Udiff, *_temp, *_l, *_r,  *_vec_normal;
537         double * _Uij, *_Vij;//Conservative and primitive interfacial states of the VFRoe scheme
538         int *_idm, *_idn;
539         double _inv_dxi,_inv_dxj;//diametre des cellules i et j autour d'une face
540         double _err_press_max,_part_imag_max,_minm1,_minm2;
541         int _nbMaillesNeg, _nbVpCplx;
542         bool _isBoundary;// la face courante est elle une face de bord ?
543         double _maxvploc;
544
545         /** \fn convectionState
546          * \brief calcule l'etat de Roe entre deux cellules voisinnes
547          * \Details c'ets une fonction virtuelle, on la surcharge dans chacun des modèles
548          * @param i,j : entiers correspondant aux numéros des cellules à gauche et à droite de l'interface
549          * @param IsBord : est un booléen qui nous dit si la cellule voisine est sur le bord ou pas
550          * */
551         virtual void convectionState(const long &i, const long &j, const bool &IsBord)=0;
552
553         /** \fn convectionMatrices
554          * \brief calcule la matrice de convection de l'etat interfacial entre deux cellules voisinnes
555          * \Details convectionMatrices est une fonction virtuelle pure, on la surcharge dans chacun des modèles
556          * */
557         virtual void convectionMatrices()=0;
558
559         /** \fn diffusionStateAndMatrices
560          * \brief calcule la matrice de diffusion de l'etat interface pour la diffusion
561          * \Details est une fonction virtuelle pure, on la surcharge dans chacun eds modèles
562          * @param i,j : entier correspondent aux indices des cellules  à gauche et droite respectivement
563          * @param IsBord: bollean telling if (i,j) is a boundary face
564          * @return
565          * */
566         virtual void diffusionStateAndMatrices(const long &i,const long &j, const bool &IsBord)=0;
567
568         /** \fn addSourceTermToSecondMember
569          * \brief Adds the contribution of source terms to the right hand side of the system: gravity,
570          *  phase change, heat power and friction
571          * @param i,j : left and right cell number
572          * @param nbNeighboursi, integer giving the number of neighbours of cell i
573          * @param nbNeighboursj, integer giving the number of neighbours of cell j
574          * @param boolean isBoundary is true for a boundary face (i,j) and false otherwise
575          * @param double mesureFace the lenght or area of the face
576          * */
577         void addSourceTermToSecondMember(const int i, int nbNeighboursi,const int j, int nbNeighboursj,bool isBoundary, int ij, double mesureFace);
578
579         /** \fn sourceVector
580          * \brief Computes the source term (at the exclusion of pressure loss terms)
581          * \Details pure virtual function, overloaded by each model
582          * @param Si output vector containing the source  term
583          * @param Ui, Vi input conservative and primitive vectors
584          * @param i the cell number. Used to read the power field
585          * @return
586          * */
587         virtual void sourceVector(PetscScalar * Si,PetscScalar * Ui,PetscScalar * Vi, int i)=0;
588
589         /** \fn convectionFlux
590          * \brief Computes the convection flux F(U) projected on a vector n
591          * @param U is the conservative variable vector 
592          * @param V is the primitive variable vector
593          * @param normal is a unit vector giving the direction where the convection flux matrix F(U) is projected
594          * @param porosity is the ration of the volume occupied by the fluid in the cell (default value is 1)
595          * @return The convection flux projected in the direction given by the normal vector: F(U)*normal */
596         virtual Vector convectionFlux(Vector U,Vector V, Vector normale, double porosity)=0;
597
598         /** \fn pressureLossVector
599          * \brief Computes the contribution of pressure loss terms in the source term
600          * \Details pure virtual function, overloaded by each model
601          * @param pressureLoss output vector containing the pressure loss contributions
602          * @param K, input pressure loss coefficient
603          * @param Ui input primitive vectors
604          * @param Vi input conservative vectors
605          * @param Uj input primitive vectors
606          * @param Vj input conservative vectors
607          * @return
608          * */
609         virtual void pressureLossVector(PetscScalar * pressureLoss, double K, PetscScalar * Ui, PetscScalar * Vi, PetscScalar * Uj, PetscScalar * Vj)=0;
610
611         /** \fn porosityGradientSourceVector
612          * \brief Computes the contribution of the porosity variation in the source term
613          * \Details pure virtual function, overloaded by each model
614          * @param porosityGradientVector output vector containing the porosity variation contributions to the source term
615          * @param Ui input primitive vectors celli
616          * @param Vi input conservative vectors celli
617          * @param porosityi input porosity value celli
618          * @param deltaxi input diameter celli
619          * @param Uj input primitive vectors cellj
620          * @param Vj input conservative vectors cellj
621          * @param porosityj input porosity value cellj
622          * @param deltaxj input diameter cellj
623          * @return
624          * */
625         virtual void porosityGradientSourceVector()=0;
626
627         /** \fn jacobian
628          * \brief Calcule la jacobienne de la ConditionLimite convection
629          * \Details est une fonction virtuelle pure, qu'on surcharge dans chacun des modèles
630          * @param j entier , l'indice de la cellule sur le bord
631          * @param nameOfGroup  : chaine de caractères, correspond au type de la condition limite
632          * */
633         virtual void jacobian(const int &j, string nameOfGroup,double * normale)=0;
634
635         /** \fn jacobianDiff
636          * \brief Calcule la jacobienne de la CL de diffusion
637          * \Details est une fonction virtuelle pure, qu'on surcharge dans chacun des modèles
638          * @param j entier , l'indice de la cellule sur le bord
639          * @param nameOfGroup  : chaine de caractères, correspond au type de la condition limite
640          * */
641         virtual void jacobianDiff(const int &j, string nameOfGroup)=0;
642
643         /** \fn setBoundaryState
644          * \brief Calcule l'etat fictif a la frontiere
645          * \Details est une fonction virtuelle pure, qu'on surcharge dans chacun des modèles
646          * @param j entier , l'indice de la cellule sur le bord
647          * @param nameOfGroup  : chaine de caractères, correspond au type de la condition limite
648          * @param normale est un vecteur de double correspond au vecteur normal à l'interface
649          * @return
650          * */
651         virtual void setBoundaryState(string nameOfGroup, const int &j,double *normale)=0;
652
653         /** \fn addDiffusionToSecondMember
654          * \brief Compute the contribution of the diffusion operator to the right hand side of the system
655          * \Details this function is pure virtual, and overloaded in each physical model class
656          * @param i left cell number
657          * @param j right cell number
658          * @param boolean isBoundary is true for a boundary face (i,j) and false otherwise
659          * */
660         virtual void addDiffusionToSecondMember(const int &i,const int &j,bool isBoundary)=0;
661
662         //!Computes the interfacial flux for the VFFC formulation of the staggered upwinding
663         virtual Vector staggeredVFFCFlux()=0;
664         //!Compute the corrected interfacial state for lowMach, pressureCorrection and staggered versions of the VFRoe formulation
665         virtual void applyVFRoeLowMachCorrections(bool isBord, string groupname="")=0;
666
667         //remplit les vecteurs de scaling
668         /** \fn computeScaling
669          * \brief Special preconditioner based on a matrix scaling strategy
670          * \Details pure  virtual function overloaded in every model class
671          * @param offset double , correspond à la plus grande valeur propre
672          * @return
673          * */
674         virtual void computeScaling(double offset) =0;
675
676         // Fonctions utilisant la loi d'etat 
677
678         /** \fn consToPrim
679          * \brief computes the primitive vector state from a conservative vector state
680          * \Details ure virtual, implemented by each model
681          * @param Ucons : conservative variable vector
682          * @pram Vprim : primitive variable vector
683          * @param porosity is the porosity coefficient in case of a porous modeling
684          * */
685         virtual void consToPrim(const double *Ucons, double* Vprim,double porosity=1) = 0;
686
687         /** \fn primToCons
688          * \brief computes the conservative vector state from a primitive vector state
689          * \Details pure virtual, implemented by each model
690          * @param U : conservative variable vector, may contain several states
691          * @pram V : primitive variable vector, may contain several states
692          * @param i : index of the conservative state in the vector U
693          * @param j : index of the primitive state in the vector V
694          *       */
695         virtual void primToCons(const double *V, const int &i, double *U, const int &j)=0;
696
697         /** \fn primToConsJacobianMatrix
698          * \brief computes the jacobian matrix of the cons->prim function
699          * \Details pure virtual, implemented by each model
700          * @pram V : primitive vector state
701          *       */
702         //void primToConsJacobianMatrix(double *V)=0;
703
704         /** \fn getRacines
705          * \brief Computes the roots of a polynomial
706          * @param polynome is a vector containing the coefficients of the polynom
707          * @return vector containing the roots (complex numbers)
708          * */
709         vector< complex<double> > getRacines(vector< double > polynome);
710
711         // Some supplement functions ---------------------------------------------------------------------------------------------
712
713         /** \fn addConvectionToSecondMember
714          * \brief Adds the contribution of the convection to the system right hand side for a face (i,j) inside the domain
715          * @param i left cell number
716          * @param j right cell number
717          * @param isBord is a boolean that is true if the interface (i,j) is a boundary interface
718          * @param groupname : is a string that may be used when isBord is true to specify which boundary the face (i,j) belongs to
719          * */
720         virtual void addConvectionToSecondMember(const int &i,const int &j,bool isBord, string groupname="");
721
722         /** \fn updatePrimitives
723          * \brief updates the global primitive vector from the global conservative vector
724          * @param void
725          * */
726         void updatePrimitives();
727
728         /** \fn updateConservatives
729          * \brief updates the global conservative vector from the global primitive vector
730          * @param void
731          * */
732         void updateConservatives();
733
734         /** \fn AbsMatriceRoe
735          * \brief Computes the absolute value of the Roe matrix
736          * @param valeurs_propres_dist is the vector of distinct eigenvalues of the Roe matrix
737          * */
738         void AbsMatriceRoe(vector< complex<double> > valeurs_propres_dist);
739
740         /** \fn SigneMatriceRoe
741          * \brief Computes the sign of the Roe matrix
742          * @param valeurs_propres_dist is the vector of distinct eigenvalues of the Roe matrix
743          * */
744         void SigneMatriceRoe(vector< complex<double> > valeurs_propres_dist);
745
746         /** \fn InvMatriceRoe
747          * \brief Computes the inverse of the Roe matrix
748          * @param valeurs_propres_dist is the vector of distinct eigenvalues of the Roe matrix
749          * */
750         void InvMatriceRoe(vector< complex<double> > valeurs_propres_dist);
751
752         /** \fn entropicShift
753          * \brief computes the eigenvalue jumps for the entropy correction
754          * @param normal vector n to the interface between the two cells _l and _r
755          * */
756         virtual void entropicShift(double* n)=0;
757
758 };
759
760 #endif /* PROBLEMFLUID_HXX_ */