1 //============================================================================
2 // Name : ProblemCoreFlows
5 // Copyright : CEA Saclay 2014
6 // Description : Generic class for thermal hydraulics problems
7 //============================================================================
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
13 #ifndef PROBLEMFLUID_HXX_
14 #define PROBLEMFLUID_HXX_
17 #include "ProblemCoreFlows.hxx"
18 #include "utilitaire_algebre.h"
22 //! enumeration SpaceScheme
23 /*! Several numerical schemes are available */
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 */
33 //! enumeration pressureEstimate
34 /*! the pressure estimate needed to fit physical parameters */
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) */
41 //! enumeration phaseType
42 /*! The fluid type can be Gas or water */
45 Liquid,/**< Fluid considered is water */
46 Gas/**< Fluid considered is Gas */
49 //! enumeration NonLinearFormulation
50 /*! the formulation used to compute the non viscous fluxes */
51 enum NonLinearFormulation
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 */
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 */
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;
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)
80 class ProblemFluid: public ProblemCoreFlows
85 * \brief constructeur de la classe ProblemFluid
89 //Gestion du calcul (interface ICoCo)
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.
96 virtual void initialize();
99 * \brief empties the memory
102 virtual void terminate();
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
109 bool initTimeStep(double dt);
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)
116 double computeTimeStep(bool & stop);
118 /** \fn abortTimeStep
119 * \brief Reset the time step dt to 0
121 void abortTimeStep();
123 /** \fn iterateTimeStep
124 * \brief calls computeNewtonVariation to perform one Newton iteration and tests the convergence
126 * @return boolean ok is true is the newton iteration gave a physically acceptable result
128 virtual bool iterateTimeStep(bool &ok);
131 * \brief saves the current results in MED or VTK files
132 * \details It is a pure virtual function overloaded in each model
135 virtual void save()=0;
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
143 virtual void validateTimeStep();
145 /* Boundary conditions */
146 /** \fn setNeumannBoundaryCondition
147 * \brief adds a new boundary condition of type Neumann
149 * \param [in] string the name of the boundary
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);
156 /** \fn setOutletBoundaryCondition
157 * \brief Adds a new boundary condition of type Outlet
159 * \param [in] string : the name of the boundary
160 * \param [in] double : the value of the pressure at the boundary
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);
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
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);
181 /** \fn setBoundaryFields
182 * \brief met à jour _limitField ( le type de condition limite )
187 void setBoundaryFields(map<string, LimitField> boundaryFields){
188 _limitField = boundaryFields;
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
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]);
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
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]);
216 * \brief sets the gravity force in the model
217 * @param gravite is a vector of size equal to the space dimension
219 void setGravity(vector<double> gravite){
220 _GravityField3d = gravite;
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
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++)
233 _fluides[i]->setDragCoeffs(dragCoeffs[i]);
234 _dragCoeffs[i]=dragCoeffs[i];
238 /** \fn getNumberOfPhases
239 * \brief The numer of phase (one or two) depending on the model considered
241 * @return the number of phases considered in the model
243 int getNumberOfPhases() const {
247 /** \fn computeNewtonVariation
248 * \brief Builds and solves the linear system to obtain the variation Ukp1-Uk in a Newton scheme
251 virtual void computeNewtonVariation();
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
258 virtual void testConservation()=0;
261 * \brief saves the velocity field in a separate 3D file so that paraview can display the streamlines
264 void saveVelocity(bool save_v=true){
265 _saveVelocity=save_v;
268 /** \fn saveConservativeField
269 * \brief saves the conservative fields (density, momentum etc...)
271 void saveConservativeField(bool save=true){
272 _saveConservativeField=save;
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
278 void setEntropicCorrection(bool entropyCorr){
279 _entropicCorrection=entropyCorr;
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)
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");
297 _pressureCorrectionOrder=order;
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
312 void setLinearSolver(linearSolver kspType, preconditioner pcType, bool scaling=false)
314 ProblemCoreFlows::setLinearSolver(kspType, pcType);
318 /** \fn setLatentHeat
319 * \brief Sets the value of the latent heat
320 * @param double L, the value of the latent heat
322 void setLatentHeat(double L){
326 /** \fn getLatentHeat
327 * \brief returns the value of the latent heat
328 * @param double L, the value of the latent heat
330 double getLatentHeat() const{
335 * \brief sets the saturation temperature
336 * @param Tsat double corresponds to saturation temperature
338 void setSatTemp(double Tsat){
343 * \brief sets the saturation temperature
344 * @param Tsat double corresponds to saturation temperature
346 double getSatTemp() const {
350 /** \fn setSatPressure
351 * \brief sets the saturation pressure
352 * @param Psat double corresponds to saturation pressure
354 void setSatPressure(double Psat, double dHsatl_over_dp=0.05){
356 _dHsatl_over_dp=dHsatl_over_dp;
359 /** \fn setPorosityField
360 * \brief set the porosity field;
361 * @param [in] Field porosity field (field on CELLS)
363 void setPorosityField(Field Porosity){
364 _porosityField=Porosity;
365 _porosityFieldSet=true;
368 /** \fn getPorosityField
369 * \brief returns the porosity field;
372 Field getPorosityField() const {
373 return _porosityField;
376 /** \fn setPorosityFieldFile
377 * \brief set the porosity field
379 * \param [in] string fileName (including file path)
380 * \param [in] string fieldName
383 void setPorosityField(string fileName, string fieldName){
384 _porosityField=Field(fileName, CELLS,fieldName);
385 _porosityFieldSet=true;
388 /** \fn setPressureLossField
389 * \brief set the pressure loss coefficients field;
390 * @param [in] Field pressure loss field (field on FACES)
392 void setPressureLossField(Field PressureLoss){
393 _pressureLossField=PressureLoss;
394 _pressureLossFieldSet=true;
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
403 void setPressureLossField(string fileName, string fieldName){
404 _pressureLossField=Field(fileName, FACES,fieldName);
405 _pressureLossFieldSet=true;
408 /** \fn setSectionField
409 * \brief set the cross section field;
410 * @param [in] Field cross section field (field on CELLS)
412 void setSectionField(Field sectionField){
413 _sectionField=sectionField;
414 _sectionFieldSet=true;
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
423 void setSectionField(string fileName, string fieldName){
424 _sectionField=Field(fileName, CELLS,fieldName);
425 _sectionFieldSet=true;
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
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;
440 /** \fn getNonLinearFormulation
441 * \brief returns the formulation used for the computation of non viscous fluxes
442 * \details Roe, VFRoe, VFFC
444 * \param [out] enum NonLinearFormulation
446 NonLinearFormulation getNonLinearFormulation() const{
447 return _nonLinearFormulation;
450 /** \fn usePrimitiveVarsInNewton
451 * \brief use Primitive Vars instead of conservative vars In Newton scheme for implicit schemes
454 void usePrimitiveVarsInNewton(bool usePrimitiveVarsInNewton){
455 _usePrimitiveVarsInNewton=usePrimitiveVarsInNewton;
458 /** \fn getSpaceScheme
459 * \brief returns the space scheme name
461 * \param [out] enum SpaceScheme(upwind, centred, pressureCorrection, pressureCorrection, staggered)
463 SpaceScheme getSpaceScheme() const;
465 /** \fn setNumericalScheme
466 * \brief sets the numerical method (upwind vs centered and explicit vs implicit)
468 * \param [in] SpaceScheme
469 * \param [in] TimeScheme
472 void setNumericalScheme(SpaceScheme scheme, TimeScheme method=Explicit);
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
483 Field getConservativeField() const ;
484 Field getPrimitiveField() const;
487 /** Number of phases in the fluid **/
489 /** Field of conservative variables (the primitive variables are defined in the mother class ProblemCoreFlows **/
491 /** Field of interfacial states of the VFRoe scheme **/
492 Field _UUstar, _VVstar;
494 SpaceScheme _spaceScheme;
495 /** the formulation used to compute the non viscous fluxes **/
496 NonLinearFormulation _nonLinearFormulation;
498 map<string, LimitField> _limitField;
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;
507 /** Fluid equation of state **/
508 vector< Fluide* > _fluides;
509 //!Viscosity coefficients
510 vector<double> _viscosite;
511 //!Conductivity coefficients
512 vector<double> _conductivite;
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
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;
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
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 ?
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
551 virtual void convectionState(const long &i, const long &j, const bool &IsBord)=0;
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
557 virtual void convectionMatrices()=0;
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
566 virtual void diffusionStateAndMatrices(const long &i,const long &j, const bool &IsBord)=0;
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
577 void addSourceTermToSecondMember(const int i, int nbNeighboursi,const int j, int nbNeighboursj,bool isBoundary, int ij, double mesureFace);
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
587 virtual void sourceVector(PetscScalar * Si,PetscScalar * Ui,PetscScalar * Vi, int i)=0;
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;
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
609 virtual void pressureLossVector(PetscScalar * pressureLoss, double K, PetscScalar * Ui, PetscScalar * Vi, PetscScalar * Uj, PetscScalar * Vj)=0;
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
625 virtual void porosityGradientSourceVector()=0;
630 * \brief Calcule la jacobienne de la ConditionLimite convection
631 * \Details est une fonction virtuelle pure, qu'on surcharge dans chacun des modèles
632 * @param j entier , l'indice de la cellule sur le bord
633 * @param nameOfGroup : chaine de caractères, correspond au type de la condition limite
635 virtual void jacobian(const int &j, string nameOfGroup,double * normale)=0;
638 * \brief Calcule la jacobienne de la CL de diffusion
639 * \Details est une fonction virtuelle pure, qu'on surcharge dans chacun des modèles
640 * @param j entier , l'indice de la cellule sur le bord
641 * @param nameOfGroup : chaine de caractères, correspond au type de la condition limite
643 virtual void jacobianDiff(const int &j, string nameOfGroup)=0;
645 /** \fn setBoundaryState
646 * \brief Calcule l'etat fictif a la frontiere
647 * \Details est une fonction virtuelle pure, qu'on surcharge dans chacun des modèles
648 * @param j entier , l'indice de la cellule sur le bord
649 * @param nameOfGroup : chaine de caractères, correspond au type de la condition limite
650 * @param normale est un vecteur de double correspond au vecteur normal à l'interface
653 virtual void setBoundaryState(string nameOfGroup, const int &j,double *normale)=0;
655 /** \fn addDiffusionToSecondMember
656 * \brief Compute the contribution of the diffusion operator to the right hand side of the system
657 * \Details this function is pure virtual, and overloaded in each physical model class
658 * @param i left cell number
659 * @param j right cell number
660 * @param boolean isBoundary is true for a boundary face (i,j) and false otherwise
662 virtual void addDiffusionToSecondMember(const int &i,const int &j,bool isBoundary)=0;
664 //!Computes the interfacial flux for the VFFC formulation of the staggered upwinding
665 virtual Vector staggeredVFFCFlux()=0;
666 //!Compute the corrected interfacial state for lowMach, pressureCorrection and staggered versions of the VFRoe formulation
667 virtual void applyVFRoeLowMachCorrections(bool isBord, string groupname="")=0;
669 //remplit les vecteurs de scaling
670 /** \fn computeScaling
671 * \brief Special preconditioner based on a matrix scaling strategy
672 * \Details pure virtual function overloaded in every model class
673 * @param offset double , correspond à la plus grande valeur propre
676 virtual void computeScaling(double offset) =0;
678 // Fonctions utilisant la loi d'etat
681 * \brief computes the primitive vector state from a conservative vector state
682 * \Details ure virtual, implemented by each model
683 * @param Ucons : conservative variable vector
684 * @pram Vprim : primitive variable vector
685 * @param porosity is the porosity coefficient in case of a porous modeling
687 virtual void consToPrim(const double *Ucons, double* Vprim,double porosity=1) = 0;
690 * \brief computes the conservative vector state from a primitive vector state
691 * \Details pure virtual, implemented by each model
692 * @param U : conservative variable vector, may contain several states
693 * @pram V : primitive variable vector, may contain several states
694 * @param i : index of the conservative state in the vector U
695 * @param j : index of the primitive state in the vector V
697 virtual void primToCons(const double *V, const int &i, double *U, const int &j)=0;
699 /** \fn primToConsJacobianMatrix
700 * \brief computes the jacobian matrix of the cons->prim function
701 * \Details pure virtual, implemented by each model
702 * @pram V : primitive vector state
704 //void primToConsJacobianMatrix(double *V)=0;
707 * \brief Computes the roots of a polynomial
708 * @param polynome is a vector containing the coefficients of the polynom
709 * @return vector containing the roots (complex numbers)
711 vector< complex<double> > getRacines(vector< double > polynome);
713 // Some supplement functions ---------------------------------------------------------------------------------------------
715 /** \fn addConvectionToSecondMember
716 * \brief Adds the contribution of the convection to the system right hand side for a face (i,j) inside the domain
717 * @param i left cell number
718 * @param j right cell number
719 * @param isBord is a boolean that is true if the interface (i,j) is a boundary interface
720 * @param groupname : is a string that may be used when isBord is true to specify which boundary the face (i,j) belongs to
722 virtual void addConvectionToSecondMember(const int &i,const int &j,bool isBord, string groupname="");
724 /** \fn updatePrimitives
725 * \brief updates the global primitive vector from the global conservative vector
728 void updatePrimitives();
730 /** \fn updateConservatives
731 * \brief updates the global conservative vector from the global primitive vector
734 void updateConservatives();
736 /** \fn AbsMatriceRoe
737 * \brief Computes the absolute value of the Roe matrix
738 * @param valeurs_propres_dist is the vector of distinct eigenvalues of the Roe matrix
740 void AbsMatriceRoe(vector< complex<double> > valeurs_propres_dist);
742 /** \fn SigneMatriceRoe
743 * \brief Computes the sign of the Roe matrix
744 * @param valeurs_propres_dist is the vector of distinct eigenvalues of the Roe matrix
746 void SigneMatriceRoe(vector< complex<double> > valeurs_propres_dist);
748 /** \fn InvMatriceRoe
749 * \brief Computes the inverse of the Roe matrix
750 * @param valeurs_propres_dist is the vector of distinct eigenvalues of the Roe matrix
752 void InvMatriceRoe(vector< complex<double> > valeurs_propres_dist);
754 /** \fn entropicShift
755 * \brief computes the eigenvalue jumps for the entropy correction
756 * @param normal vector n to the interface between the two cells _l and _r
758 virtual void entropicShift(double* n)=0;
762 #endif /* PROBLEMFLUID_HXX_ */