Salome HOME
Added a function to save all fields in signe phase flows
[tools/solverlab.git] / CoreFlows / Models / src / FiveEqsTwoFluid.cxx
1 /*
2  * FiveEqsTwoFluid.cxx
3  *
4  *  Created on: Jan 28, 2015
5  *      Author: M. Ndjinga
6  */
7 #include "FiveEqsTwoFluid.hxx"
8 #include <cstdlib>
9
10 using namespace std;
11
12 extern "C" int dgeev_(char *jobvl, char *jobvr, int *n, double *
13                 a, int *lda, double *wr, double *wi, double *vl,
14                 int *ldvl, double *vr, int *ldvr, double *work,
15                 int *lwork, int *info);
16
17
18 FiveEqsTwoFluid::FiveEqsTwoFluid(pressureEstimate pEstimate, int dim){
19         _Ndim=dim;
20         _nVar=2*(_Ndim+1)+1;
21         _nbPhases = 2;
22         _dragCoeffs=vector<double>(2,0);
23         _fluides.resize(2);
24         //Ne pas utiliser la loi de Stephane Dellacherie mais la stiffened gas standard
25         if (pEstimate==around1bar300K)//EOS at 1 bar and 373K
26         {
27                 cout<<"Fluid is water-Gas mixture around saturation point 1 bar and 373 K (100°C)"<<endl;
28                 *_runLogFile<<"Fluid is water-Gas mixture around saturation point 1 bar and 373 K (100°C)"<<endl;
29                 _fluides[0] = new StiffenedGas(1.34,1555,373,2.5e6);  //ideal gas law for Gas at pressure 1 bar and temperature 100°C: eref1=2.5e6
30                 _fluides[1] = new StiffenedGas(958,1e5,373,4.2e5,1543,3769);  //stiffened gas law for water at pressure 1 bar and temperature 100°C: eref2=5e5
31                 _Tsat=373;//saturation temperature at 1 bar
32                 _hsatl=4.2e5;//water enthalpy at saturation at 1 bar
33                 _hsatv=2.5e6;//Gas enthalpy at saturation at 1 bar
34         }
35         else//EOS at 155 bars and 618K
36         {
37                 cout<<"Fluid is water-Gas mixture around saturation point 155 bars and 618 K (345°C)"<<endl;
38                 *_runLogFile<<"Fluid is water-Gas mixture around saturation point 155 bars and 618 K (345°C)"<<endl;
39                 _fluides[0] = new StiffenedGas(102,1.55e7,618,2.44e6, 433,3633);  //stiffened gas law for Gas at pressure 155 bar and temperature 345°C: eref1=2.4e6
40                 _fluides[1] = new StiffenedGas(594,1.55e7,618,1.6e6, 621,3100);  //stiffened gas law for water at pressure 155 bar and temperature 345°C: eref2=1.6e6
41                 _Tsat=618;//saturation temperature at 155 bars
42                 _hsatl=1.63e6;//water enthalpy at saturation at 155 bars
43                 _hsatv=2.6e6;//Gas enthalpy at saturation at 155 bars
44         }
45         _latentHeat=_hsatv-_hsatl;
46         _intPressCoeff=1.5;
47 }
48
49 void FiveEqsTwoFluid::initialize()
50 {
51         cout<<"Initialising the five equation two fluid model"<<endl;
52         *_runLogFile<<"Initialising the five equation two fluid model"<<endl;
53
54         if(static_cast<StiffenedGas*>(_fluides[0])==NULL || static_cast<StiffenedGas*>(_fluides[1])==NULL)
55                 throw CdmathException("FiveEqsTwoFluid::initialize: both phase must have stiffened gas EOS");
56
57         _Uroe = new double[_nVar+1];
58
59         _lCon = new PetscScalar[_nVar];//should be deleted in ::terminate
60         _rCon = new PetscScalar[_nVar];//should be deleted in ::terminate
61         _JacoMat = new PetscScalar[_nVar*_nVar];//should be deleted in ::terminate
62
63         _gravite = vector<double>(_nVar,0);//Not to be confused with _GravityField3d (size _Ndim). _gravite (size _Nvar) is usefull for dealing with source term and implicitation of gravity vector
64         for(int i=0; i<_Ndim; i++)
65         {
66                 _gravite[i+1]=_GravityField3d[i];
67                 _gravite[i+1 +_Ndim+1]=_GravityField3d[i];
68         }
69         _GravityImplicitationMatrix = new PetscScalar[_nVar*_nVar];
70
71         if(_saveVelocity){
72                 _Vitesse1=Field("Gas velocity",CELLS,_mesh,3);//Forcement en dimension 3 pour le posttraitement des lignes de courant
73                 _Vitesse2=Field("Liquid velocity",CELLS,_mesh,3);//Forcement en dimension 3 pour le posttraitement des lignes de courant
74         }
75
76         if(_entropicCorrection)
77                 _entropicShift=vector<double>(_nVar);
78
79         ProblemFluid::initialize();
80 }
81
82 void FiveEqsTwoFluid::convectionState( const long &i, const long &j, const bool &IsBord){
83         //sortie: WRoe en (alpha, p, u1, u2, T, dm1,dm2,dalpha1,dp)
84         //entree: _conservativeVars en (rho1, rho1 u1, rho2, rho2 u2)
85
86         _idm[0] = _nVar*i;
87         for(int k=1; k<_nVar; k++)
88                 _idm[k] = _idm[k-1] + 1;
89         VecGetValues(_conservativeVars, _nVar, _idm, _Ui);
90
91         _idm[0] = _nVar*j;
92         for(int k=1; k<_nVar; k++)
93                 _idm[k] = _idm[k-1] + 1;
94         if(IsBord)
95                 VecGetValues(_Uext, _nVar, _idm, _Uj);
96         else
97                 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
98
99         if(_verbose && _nbTimeStep%_freqSave ==0)
100         {
101                 cout<<"Convection Left state cell " << i<< ": "<<endl;
102                 for(int k =0; k<_nVar; k++)
103                         cout<< _Ui[k]<<endl;
104                 cout<<"Convection Right state cell " << j<< ": "<<endl;
105                 for(int k =0; k<_nVar; k++)
106                         cout<< _Uj[k]<<endl;
107         }
108         /*
109         if(_Ui[0]<-(_precision) || _Uj[0]<-(_precision) || _Ui[_Ndim+1]<-(_precision) || _Uj[_Ndim+1]<-(_precision))
110         {
111                 cout<<"Warning: masse partielle negative!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!"<<endl;
112                 cout<< "valeurs a gauche: "<<_Ui[0]<<", "<<_Ui[_Ndim+1]<<", valeurs a droite: "<<_Uj[0]<<", "<<_Uj[_Ndim+1]<<endl;
113                 //        throw CdmathException(" Masse partielle negative, arret de calcul");
114         }
115         else
116                 _Ui[0]=max(0.,_Ui[0]);
117                 _Uj[0]=max(0.,_Uj[0]);
118                 _Ui[_Ndim+1]=max(0.,_Ui[_Ndim+1]);
119                 _Uj[_Ndim+1]=max(0.,_Uj[_Ndim+1]);
120          */
121
122         PetscScalar ri1, ri2, rj1, rj2, xi, xj;
123         //get _l and _r the primitive  states left and right of the interface
124         _idm[0] = _nVar*i;
125         for(int k=1; k<_nVar; k++)
126                 _idm[k] = _idm[k-1] + 1;
127         VecGetValues(_primitiveVars, _nVar, _idm, _l);
128
129         if(IsBord)
130         {
131                 //cout<<"_r is border"<<endl;
132                 //consToPrim(_Uj, _r);
133                 _idm[0] = 0;
134                 for(int k=1; k<_nVar; k++)
135                         _idm[k] = _idm[k-1] + 1;
136                 VecGetValues(_Vext, _nVar, _idm, _r);
137         }
138         else
139         {
140                 _idm[0] = _nVar*j;
141                 for(int k=1; k<_nVar; k++)
142                         _idm[k] = _idm[k-1] + 1;
143                 VecGetValues(_primitiveVars, _nVar, _idm, _r);
144         }
145         if(_verbose && _nbTimeStep%_freqSave ==0)
146         {
147                 cout<<"_l: "<<endl;
148                 for(int k=0;k<_nVar; k++)
149                         cout<< _l[k]<<endl;
150                 cout<<"_r: "<<endl;
151                 for(int k=0;k<_nVar; k++)
152                         cout<< _r[k]<<endl;
153         }
154         // _Uroe[0] = \tilde{\alpha_v} = 1 - \tilde{\alpha_l} (formula of Toumi)
155         if(2-_l[0]-_r[0] > _precision)
156                 _Uroe[0] = 1- 2*(1-_l[0])*(1-_r[0])/(2-_l[0]-_r[0]);
157         else
158                 _Uroe[0] = (_l[0]+_r[0])/2;
159
160         if(_l[0]+_r[0] > _precision)
161                 _Uroe[1] = (_l[1]*_l[0]+_r[1]*_r[0])/(_l[0]+_r[0]);
162         else
163                 _Uroe[1] = (_l[1]*(1-_l[0])+_r[1]*(1-_r[0]))/(2-_l[0]-_r[0]);
164
165         ri1 = sqrt(_Ui[0]); ri2 = sqrt(_Ui[_Ndim+1]);
166         rj1 = sqrt(_Uj[0]); rj2 = sqrt(_Uj[_Ndim+1]);
167         for(int k=0;k<_Ndim;k++)
168         {
169                 xi = _Ui[k+1];
170                 xj = _Uj[k+1];
171                 if(ri1>_precision && rj1>_precision)
172                         _Uroe[2+k] = (xi/ri1 + xj/rj1)/(ri1 + rj1);
173                 else if(ri1<_precision && rj1>_precision)
174                         _Uroe[2+k] =  xj/_Uj[0];
175                 else if(ri1>_precision && rj1<_precision)
176                         _Uroe[2+k] =  xi/_Ui[0];
177                 else
178                         _Uroe[2+k] =(_Ui[k+1+_Ndim+1]/ri2 + _Uj[k+1+_Ndim+1]/rj2)/(ri2 + rj2);
179
180                 xi = _Ui[k+1+_Ndim+1];
181                 xj = _Uj[k+1+_Ndim+1];
182                 if(ri2>_precision && rj2>_precision)
183                         _Uroe[1+k+_Ndim+1] = (xi/ri2 + xj/rj2)/(ri2 + rj2);
184                 else if(ri2<_precision && rj2>_precision)
185                         _Uroe[1+k+_Ndim+1] = xj/_Uj[_Ndim+1];
186                 else if(ri2>_precision && rj2<_precision)
187                         _Uroe[1+k+_Ndim+1] = xi/_Ui[_Ndim+1];
188                 else
189                         _Uroe[1+k+_Ndim+1] = (xi/ri1 + xj/rj1)/(ri1 + rj1);
190         }
191         _Uroe[_nVar-1]=.5*(_l[_nVar-1]+_r[_nVar-1]);
192
193         //Fin du remplissage dans la fonction convectionMatrices
194
195         if(_verbose && _nbTimeStep%_freqSave ==0)
196         {
197                 cout<<"Etat de Roe calcule: "<<endl;
198                 for(int k=0;k<_nVar; k++)
199                         cout<< _Uroe[k]<<endl;
200         }
201 }
202
203 void FiveEqsTwoFluid::diffusionStateAndMatrices(const long &i,const long &j, const bool &IsBord){
204         //sortie: matrices et etat Diffusion (alpha1 rho1, q1, alpha2 rho2, q2,T)
205         _idm[0] = _nVar*i;
206         for(int k=1; k<_nVar; k++)
207                 _idm[k] = _idm[k-1] + 1;
208
209         VecGetValues(_conservativeVars, _nVar, _idm, _Ui);
210         _idm[0] = _nVar*j;
211         for(int k=1; k<_nVar; k++)
212                 _idm[k] = _idm[k-1] + 1;
213
214         if(IsBord)
215                 VecGetValues(_Uextdiff, _nVar, _idm, _Uj);
216         else
217                 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
218
219         for(int k=0; k<_nVar; k++)
220                 _Udiff[k] = (_Ui[k]+_Uj[k])/2;
221         double q1_2,q2_2=0;
222         for (int i = 0; i<_Ndim;i++){
223                 q1_2+=_Udiff[        i+1]*_Udiff[        i+1];
224                 q2_2+=_Udiff[1+_Ndim+i+1]*_Udiff[1+_Ndim+i+1];
225         }
226         consToPrim(_Udiff,_phi);
227         _Udiff[_nVar-1]=_phi[_nVar-1];
228         double alpha=_phi[0];
229         double Tm=_phi[_nVar-1];
230         double mu1 = _fluides[0]->getViscosity(Tm);
231         double mu2 = _fluides[1]->getViscosity(Tm);
232         double lambda = alpha*_fluides[0]->getConductivity(Tm)+(1-alpha)*_fluides[1]->getConductivity(Tm);
233         double Cv1= _fluides[0]->constante("Cv");
234         double Cv2= _fluides[1]->constante("Cv");
235
236         if(_timeScheme==Implicit)
237         {
238                 for(int i=0; i<_nVar*_nVar;i++)
239                         _Diffusion[i] = 0;
240                 for(int idim=1;idim<_Ndim+1;idim++)
241                 {
242                         if(alpha>_precision){
243                                 _Diffusion[idim*_nVar] = alpha* mu1*_Udiff[idim]/(_Udiff[0]*_Udiff[0]);
244                                 _Diffusion[idim*_nVar+idim] = -alpha* mu1/_Udiff[0];
245                         }
246                         if(1-alpha>_precision){
247                                 _Diffusion[(idim+_Ndim+1)*_nVar] =  (1-alpha)* mu2*_Udiff[idim+_Ndim+1]/(_Udiff[_Ndim+1]*_Udiff[_Ndim+1]);
248                                 _Diffusion[(idim+_Ndim+1)*_nVar+idim+_Ndim+1] = -(1-alpha)* mu2/_Udiff[_Ndim+1];
249                         }
250                 }
251                 /*//Should correct the formula before using
252                 int i = (_nVar-1)*_nVar;
253                 _Diffusion[i]=lambda*(Tm/_Udiff[0]-q1_2/(2*Cv1*_Udiff[0]*_Udiff[0]*_Udiff[0]));
254                 _Diffusion[i+1+_Ndim]=lambda*(Tm/_Udiff[1+_Ndim]-q2_2/(2*Cv2*_Udiff[1+_Ndim]*_Udiff[1+_Ndim]*_Udiff[1+_Ndim]));
255                 for(int k=1;k<1+_Ndim;k++)
256                 {
257                         _Diffusion[i+k]= lambda*_Udiff[k]/(_Udiff[0]*_Udiff[0]*Cv1);
258                         _Diffusion[i+k+1+_Ndim]= lambda*_Udiff[k+1+_Ndim]/(_Udiff[1+_Ndim]*_Udiff[+1+_Ndim]*Cv2);
259                 }
260                 _Diffusion[_nVar*_nVar-1]=-lambda/(_Udiff[0]*Cv1+_Udiff[1+_Ndim]*Cv2);
261                  */
262         }
263 }
264
265 void FiveEqsTwoFluid::sourceVector(PetscScalar * Si,PetscScalar * Ui,PetscScalar * Vi, int i)
266 {
267         double m1=Ui[0],m2=Ui[1+_Ndim],rho=m1+m2, rhoE=Ui[_nVar-1], T=Vi[_nVar-1],alpha=Vi[0], P=Vi[1];
268         double norm_ur=0,norm_u1sq=0,norm_u2sq=0, Gamma;
269
270         for(int k=0; k<_Ndim; k++){
271                 norm_ur+=(Vi[2+k]-Vi[2+k+_Ndim])*(Vi[2+k]-Vi[2+k+_Ndim]);
272                 norm_u1sq+=Vi[2+k]*Vi[2+k];
273                 norm_u2sq+=Vi[2+k+_Ndim]*Vi[2+k+_Ndim];
274         }
275         norm_ur=sqrt(norm_ur);
276         double h=(rhoE-0.5*m1*norm_u1sq-0.5*m2*norm_u2sq+P)/rho;
277         double u_int[_Ndim];
278         for(int k=0; k<_Ndim; k++)
279                 u_int[k] = 0.5*(Vi[2+k]+Vi[2+k+_Ndim]);
280         //              u_int[k] = Vi[0]*Vi[2+k+_Ndim] + (1-Vi[0])*Vi[2+k];
281         if(i>=0 && T>_Tsat && alpha<1-_precision)//if(i>=0 && _hsatv>h  && h>_hsatl && alpha<1-_precision)
282                 Gamma=_heatPowerField(i)/_latentHeat;
283         else//boundary cell, no phase change
284                 Gamma=0;
285
286         for(int k=1; k<_Ndim+1; k++)
287         {
288                 Si[k] =_gravite[k]*m1-_dragCoeffs[0]*norm_ur*(Vi[1+k]-Vi[1+k+_Ndim]) + Gamma*u_int[k-1];//interfacial velocity= ul
289                 Si[k+_Ndim+1] =_gravite[k+_Ndim+1]*m2+ _dragCoeffs[0]*norm_ur*(Vi[1+k]-Vi[1+k+_Ndim]) - Gamma*u_int[k-1];
290         }
291         if(true){//heated boiling
292                 Si[0]=Gamma;
293                 Si[1+_Ndim]=-Gamma;
294         }
295         else if (P<_Psat && alpha<1-_precision){//flash boiling
296                 Si[0]=-_dHsatl_over_dp*_dp_over_dt(i)/_latentHeat;
297                 Si[1+_Ndim]=_dHsatl_over_dp*_dp_over_dt(i)/_latentHeat;
298         }
299         else{
300                 Si[0]=0;
301                 Si[1+_Ndim]=0;
302         }
303         if(i>=0)
304                 Si[_nVar-1]=_heatPowerField(i);
305         else//boundary cell, no heating
306                 Si[_nVar-1]=0;
307         for(int k=0; k<_Ndim; k++)
308                 Si[_nVar-1] +=_GravityField3d[k]*(Ui[1+k]+Ui[2+k+_Ndim])-_dragCoeffs[0]*norm_ur*(Vi[2+k]-Vi[2+k+_Ndim])*(Vi[2+k]-Vi[2+k+_Ndim]);
309
310         if(_timeScheme==Implicit)
311         {
312                 for(int i=0; i<_nVar*_nVar;i++)
313                         _GravityImplicitationMatrix[i] = 0;
314                 for(int i=0; i<_nVar/2;i++)
315                         _GravityImplicitationMatrix[i*_nVar]=-_gravite[i];
316                 for(int i=_nVar/2; i<_nVar;i++)
317                         _GravityImplicitationMatrix[i*_nVar+_nVar/2]=-_gravite[i];
318         }
319
320         if(_verbose && _nbTimeStep%_freqSave ==0)
321         {
322                 cout<<"FiveEqsTwoFluid::sourceVector"<<endl;
323                 cout<<"Ui="<<endl;
324                 for(int k=0;k<_nVar;k++)
325                         cout<<Ui[k]<<", ";
326                 cout<<endl;
327                 cout<<"Vi="<<endl;
328                 for(int k=0;k<_nVar;k++)
329                         cout<<Vi[k]<<", ";
330                 cout<<endl;
331                 cout<<"Si="<<endl;
332                 for(int k=0;k<_nVar;k++)
333                         cout<<Si[k]<<", ";
334                 cout<<endl;
335                 if(_timeScheme==Implicit)
336                         displayMatrix(_GravityImplicitationMatrix, _nVar, "Gravity implicitation matrix");
337         }
338 }
339
340 void FiveEqsTwoFluid::pressureLossVector(PetscScalar * pressureLoss, double K, PetscScalar * Ui, PetscScalar * Vi, PetscScalar * Uj, PetscScalar * Vj)
341 {
342         double norm_u1=0, u1_n=0, norm_u2=0, u2_n=0, m1, m2;
343         for(int i=0;i<_Ndim;i++){
344                 u1_n += _Uroe[1+i]      *_vec_normal[i];
345                 u2_n += _Uroe[1+i+_Ndim]*_vec_normal[i];
346         }
347         pressureLoss[0]=0;
348         pressureLoss[1+_Ndim]=0;
349         if(u1_n>0){
350                 for(int i=0;i<_Ndim;i++)
351                         norm_u1 += Vi[1+i]*Vi[1+i];
352                 norm_u1=sqrt(norm_u1);
353                 m1=Ui[0];
354                 for(int i=0;i<_Ndim;i++)
355                         pressureLoss[1+i]=-K*m1*norm_u1*Vi[1+i];
356         }
357         else{
358                 for(int i=0;i<_Ndim;i++)
359                         norm_u1 += Vj[1+i]*Vj[1+i];
360                 norm_u1=sqrt(norm_u1);
361                 m1=Uj[0];
362                 for(int i=0;i<_Ndim;i++)
363                         pressureLoss[1+i]=-K*m1*norm_u1*Vj[1+i];
364         }
365         if(u2_n>0){
366                 for(int i=0;i<_Ndim;i++)
367                         norm_u2 += Vi[2+i+_Ndim]*Vi[2+i+_Ndim];
368                 norm_u2=sqrt(norm_u2);
369                 m2=Ui[1+_Ndim];
370                 for(int i=0;i<_Ndim;i++)
371                         pressureLoss[2+i+_Ndim]=-K*m2*norm_u2*Vi[2+i+_Ndim];
372         }
373         else{
374                 for(int i=0;i<_Ndim;i++)
375                         norm_u2 += Vj[2+i+_Ndim]*Vj[2+i+_Ndim];
376                 norm_u2=sqrt(norm_u2);
377                 m2=Uj[1+_Ndim];
378                 for(int i=0;i<_Ndim;i++)
379                         pressureLoss[2+i+_Ndim]=-K*m2*norm_u2*Vj[2+i+_Ndim];
380         }
381         pressureLoss[_nVar-1]=-K*(m1*norm_u1*norm_u1*norm_u1+m2*norm_u2*norm_u2*norm_u2);
382
383         if(_verbose && _nbTimeStep%_freqSave ==0)
384         {
385                 cout<<"FiveEqsTwoFluid::pressureLossVector K= "<<K<<endl;
386                 cout<<"Ui="<<endl;
387                 for(int k=0;k<_nVar;k++)
388                         cout<<Ui[k]<<", ";
389                 cout<<endl;
390                 cout<<"Vi="<<endl;
391                 for(int k=0;k<_nVar;k++)
392                         cout<<Vi[k]<<", ";
393                 cout<<endl;
394                 cout<<"Uj="<<endl;
395                 for(int k=0;k<_nVar;k++)
396                         cout<<Uj[k]<<", ";
397                 cout<<endl;
398                 cout<<"Vj="<<endl;
399                 for(int k=0;k<_nVar;k++)
400                         cout<<Vj[k]<<", ";
401                 cout<<endl;
402                 cout<<"pressureLoss="<<endl;
403                 for(int k=0;k<_nVar;k++)
404                         cout<<pressureLoss[k]<<", ";
405                 cout<<endl;
406         }
407 }
408
409 void FiveEqsTwoFluid::porosityGradientSourceVector()
410 {
411         double u1_ni=0, u1_nj=0, u2_ni=0, u2_nj=0, rho1i, rho2i, rho1j, rho2j, pi=_Vi[1], pj=_Vj[1], Ti=_Vi[_nVar-1], Tj=_Vj[_nVar-1];
412         double pij1, pij2, alphaij=_Uroe[0];
413         for(int i=0;i<_Ndim;i++) {
414                 u1_ni += _Vi[2+i]*_vec_normal[i];
415                 u2_ni += _Vi[2+_Ndim+i]*_vec_normal[i];
416                 u1_nj += _Vj[2+i]*_vec_normal[i];
417                 u2_nj += _Vj[2+_Ndim+i]*_vec_normal[i];
418         }
419         _porosityGradientSourceVector[0]=0;
420         _porosityGradientSourceVector[1+_Ndim]=0;
421         rho1i = _fluides[0]->getDensity(pi, Ti);
422         rho2i = _fluides[1]->getDensity(pi, Ti);
423         rho1j = _fluides[0]->getDensity(pj, Tj);
424         rho2j = _fluides[1]->getDensity(pj, Tj);
425         pij1=(pi+pj)/2+rho1i*rho1j/2/(rho1i+rho1j)*(u1_ni-u1_nj)*(u1_ni-u1_nj);
426         pij2=(pi+pj)/2+rho2i*rho2j/2/(rho2i+rho2j)*(u2_ni-u2_nj)*(u2_ni-u2_nj);
427         for(int i=0;i<_Ndim;i++){
428                 _porosityGradientSourceVector[1+i]      =alphaij*pij1*(_porosityi-_porosityj)*2/(1/_inv_dxi+1/_inv_dxj);
429                 _porosityGradientSourceVector[2+_Ndim+i]=alphaij*pij2*(_porosityi-_porosityj)*2/(1/_inv_dxi+1/_inv_dxj);
430         }
431         _porosityGradientSourceVector[_nVar-1]=0;
432 }
433
434 double FiveEqsTwoFluid::intPressDef(double alpha, double u_r2, double rho1, double rho2, double Temperature)
435 {
436         return  _intPressCoeff*alpha*(1-alpha)*rho1*rho2*u_r2/( alpha*rho2+(1-alpha)*rho1);
437         +alpha*(1-alpha)*rho1*rho2*u_r2/((alpha*rho2+(1-alpha)*rho1)*(alpha*rho2+(1-alpha)*rho1)*(alpha*rho2+(1-alpha)*rho1)*(alpha*rho2+(1-alpha)*rho1))*u_r2
438                         *(alpha*alpha*rho2-(1-alpha)*(1-alpha)*rho1)
439                         *(alpha*alpha*rho2*rho2/(_fluides[0]->vitesseSonTemperature(Temperature,rho1)*_fluides[0]->vitesseSonTemperature(Temperature,rho1))
440                                         -(1-alpha)*(1-alpha)*rho1*rho1/(_fluides[1]->vitesseSonTemperature(Temperature,rho2)*_fluides[1]->vitesseSonTemperature(Temperature,rho2)));
441 }
442
443 Vector FiveEqsTwoFluid::convectionFlux(Vector U,Vector V, Vector normale, double porosity){
444         if(_verbose){
445                 cout<<"Ucons"<<endl;
446                 cout<<U<<endl;
447                 cout<<"Vprim"<<endl;
448                 cout<<V<<endl;
449         }
450
451         double phim1=U(0);//phi alpha1 rho1
452         double phim2=U(1+_Ndim);//phi alpha2 rho2
453         Vector phiq1(_Ndim),phiq2(_Ndim);//phi alpha1 rho1 u1, phi alpha2 rho2 u2
454         for(int i=0;i<_Ndim;i++){
455                 phiq1(i)=U(1+i);
456                 phiq2(i)=U(2+_Ndim+i);
457         }
458         double alpha=V(0);
459         double pression=V(1);
460         Vector vitesse1(_Ndim),vitesse2(_Ndim);
461         for(int i=0;i<_Ndim;i++){
462                 vitesse1(i)=V(2+i);
463                 vitesse2(i)=V(2+_Ndim+i);
464         }
465         double Temperature= V(_nVar-1);
466
467         double vitesse1n=vitesse1*normale;
468         double vitesse2n=vitesse2*normale;
469         double rho1=_fluides[0]->getDensity(pression,Temperature);
470         double rho2=_fluides[1]->getDensity(pression,Temperature);
471         double e1_int=_fluides[0]->getInternalEnergy(Temperature,rho1);
472         double e2_int=_fluides[1]->getInternalEnergy(Temperature,rho2);
473
474         double alpha_roe = _Uroe[0];//Toumi formula
475         // interfacial pressure term (hyperbolic correction)
476         double dpi = _Uroe[_nVar];
477
478         Vector F(_nVar);
479         F(0)=phim1*vitesse1n;
480         F(1+_Ndim)=phim2*vitesse2n;
481         for(int i=0;i<_Ndim;i++){
482                 F(1+i)=phim1*vitesse1n*vitesse1(i)+(alpha_roe*pression+dpi*alpha)*porosity*normale(i);
483                 F(2+_Ndim+i)=phim2*vitesse2n*vitesse2(i)+((1-alpha_roe)*pression+dpi*(1-alpha))*normale(i)*porosity;
484         }
485         F(_nVar-1)=phim1*(e1_int+0.5*vitesse1*vitesse1+pression/rho1)*vitesse1n+phim2*(e2_int+0.5*vitesse2*vitesse2+pression/rho2)*vitesse2n;
486
487         if(_verbose){
488                 cout<<"Flux F(U,V)"<<endl;
489                 cout<<F<<endl;
490         }
491
492         return F;
493 }
494
495 void FiveEqsTwoFluid::convectionJacobianMatrix(double *V, double *n)
496 {
497         complex< double > tmp;
498         // enter : V(nVar) : primitive variables
499         //            n    : normal vector
500         double alp = V[0];
501         double p = V[1];
502         double Tm = V[_nVar-1];
503         double rho1 = _fluides[0]->getDensity(p, Tm);
504         double rho2 = _fluides[1]->getDensity(p, Tm);
505         double ur_2 = 0;
506         for (int idim=0; idim<_Ndim; idim++){
507                 ur_2 += (V[2+idim]-V[2+idim+_Ndim])*(V[2+idim]-V[2+idim+_Ndim]);
508         }
509         // interfacial pressure term (hyperbolic correction)
510         double dpi1 = intPressDef(alp,ur_2, rho1, rho2,Tm);
511         double dpi2 = dpi1;
512
513         /********Prepare the parameters to compute the Jacobian Matrix********/
514         /**** coefficients a, b, c ****/
515         double inv_a1_2,inv_a2_2,b1,c1,a2,b2,c2;
516         double e1,e2;
517         e1 = _fluides[0]->getInternalEnergy(V[_nVar-1],rho1);// primitive variable _l[_nVar-1]=Tm
518         e2 = _fluides[1]->getInternalEnergy(V[_nVar-1],rho2);
519         inv_a1_2 = static_cast<StiffenedGas*>(_fluides[0])->getDiffDensPress(e1);
520         inv_a2_2 = static_cast<StiffenedGas*>(_fluides[1])->getDiffDensPress(e2);
521         //double getJumpDensInternalEnergy(const double p_l,const double p_r,const double e_l,const double e_r);
522         b1 = static_cast<StiffenedGas*>(_fluides[0])->getDiffDensInternalEnergy(p,e1);
523         b2 = static_cast<StiffenedGas*>(_fluides[1])->getDiffDensInternalEnergy(p,e2);
524         //double getJumpInternalEnergyTemperature();
525         c1 = static_cast<StiffenedGas*>(_fluides[0])->getDiffInternalEnergyTemperature();
526         c2 = static_cast<StiffenedGas*>(_fluides[1])->getDiffInternalEnergyTemperature();
527         /**** coefficients eta,  varrho_2 ****/
528         double eta[_Ndim], varrho_2;
529         // prefix m is arithmetic mean
530         double m1,m2, eta_n;
531         double u1[_Ndim],u2[_Ndim], alp_u1[_Ndim],alp_u2[_Ndim];
532         m1 = alp*rho1;
533         m2 = (1-alp)*rho2;
534         varrho_2 =1/((alp*rho2)*inv_a1_2+((1-alp)*rho1)*inv_a2_2);
535         eta_n = 0.;
536         for (int idim=0; idim<_Ndim; idim++){
537                 u1[idim] = V[idim+2];
538                 u2[idim] = V[_Ndim+idim+2];
539                 alp_u1[idim] = alp*u1[idim];
540                 alp_u2[idim] = (1-alp)*u2[idim];
541                 eta_n += (alp_u1[idim]*(1-p/rho1*inv_a1_2)+alp_u2[idim]*(1-p/rho2*inv_a2_2))*n[idim];
542         }
543         double eta_varrho_2n = eta_n*varrho_2;
544         /**** compute jump of Delta T, Delta e1, Delta e2 ****/
545         double inv_cm = 1/(c1*m1+c2*m2);
546         double DeltaT [_nVar], Delta_e1[_nVar], Delta_e2[_nVar];
547         // initialize DeltaT
548         DeltaT[0] =-e1;
549         DeltaT[1+_Ndim] =-e2;
550         DeltaT[_nVar-1] = 1 ;
551         for (int idim=0; idim<_Ndim; idim++){
552                 DeltaT[idim+1] = 0.;
553                 DeltaT[1+_Ndim+idim+1] = 0;
554         }
555         for (int idim=0; idim<_Ndim; idim++){
556                 // wrt mass gas
557                 DeltaT[0] += 0.5*u1[idim]*u1[idim];
558                 // wrt mass liquid
559                 DeltaT[_Ndim+1] += 0.5*u2[idim]*u2[idim];
560                 // wrt momentum gass
561                 DeltaT[idim+1] += - u1[idim];//*n[idim]
562                 // wrt momentum liquid
563                 DeltaT[_Ndim+idim+2] += - u2[idim];//*n[idim]
564         }
565         // finalize  DeltaT, Delta_e1 and Delta_e2
566         for (int i =0; i< _nVar; ++i){
567                 DeltaT[i] = inv_cm*DeltaT[i];
568                 Delta_e1[i] = c1*DeltaT[i];
569                 Delta_e2[i] = c2*DeltaT[i];
570         }
571         /**** compute differential flux (energy equation) A5 ****/
572
573         double dF5[_nVar];
574         // initialize A5
575         for (int i=0; i<_nVar; i++){
576                 dF5[i]=0;
577         }
578         dF5[0] = eta_varrho_2n*rho2; // mass gas
579         dF5[_Ndim+1] = eta_varrho_2n*rho1; // mass liquid
580         for (int idim=0; idim<_Ndim; idim++){
581                 // momentum gas
582                 dF5[idim+1]= (e1+p/rho1)*n[idim];
583                 // momentum liquid
584                 dF5[_Ndim+idim+2]=(e2+p/rho2)*n[idim];
585         }
586         // assign the value of A5 (last row of the Roe matrix)
587         for (int idim=0; idim<_Ndim; idim++){
588                 for (int jdim=0; jdim<_Ndim;jdim++){
589                         dF5[0] -= u1[idim]*u1[jdim]*u1[jdim]*n[idim];// -uin * ujn^2
590                         dF5[_Ndim+1] -= u2[idim]*u2[jdim]*u2[jdim]*n[idim];
591                         //momentum gas
592                         dF5[idim+1] += u1[idim]*u1[jdim]*n[jdim]+0.5*(u1[jdim]*u1[jdim])*n[idim];
593                         //momentum liquid
594                         dF5[_Ndim+idim+2] += u2[idim]*u2[jdim]*n[jdim]+0.5*(u2[jdim]*u2[jdim])*n[idim];
595                 }
596         }
597         // final dF5
598         double coef_e1, coef_e2;
599         coef_e1 = - eta_varrho_2n*alp*rho2*b1;
600         coef_e2 = - eta_varrho_2n*(1-alp)*rho1*b2;
601         for (int idim=0; idim<_Ndim; idim++){
602                 coef_e1 += (alp*rho1 - alp*p*b1/rho1)*u1[idim]*n[idim];
603                 coef_e2 += ((1-alp)*rho2 - (1-alp)*p*b2/rho2)*u2[idim]*n[idim];
604         }
605         for (int i =0; i< _nVar; i++){
606                 dF5[i] += coef_e1*Delta_e1[i] + coef_e2*Delta_e2[i];
607         }
608         /******** Construction de la matrice J *********/
609         //lignes de masse
610         for(int i=0; i<_nVar*_nVar;i++)
611                 _JacoMat[i]=0.;
612
613         for(int idim=0; idim<_Ndim;idim++)
614         {
615                 _JacoMat[1+idim]=n[idim];
616                 _JacoMat[1+idim+_Ndim+1]=0.;
617                 _JacoMat[(_Ndim+1)*_nVar+1+idim]=0.;
618                 _JacoMat[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=n[idim];
619         }
620         // Roe Matrix new version
621         for(int idim=0; idim<_Ndim;idim++)
622                 for (int jdim=0; jdim<_Ndim;jdim++){
623                         // momentum gas (neglect delta alpha and delta P)
624                         _JacoMat[ (1+idim)*_nVar] += -u1[idim]*u1[jdim]*n[jdim];
625                         _JacoMat[(1+idim)*_nVar+jdim+1] += u1[idim]*n[jdim];
626                         _JacoMat[(1+idim)*_nVar+idim+1] += u1[jdim]*n[jdim];
627                         // momentum liquid (neglect delta alpha and delta P)
628                         _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1] += -u2[idim]*u2[jdim]*n[jdim];
629                         _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1+jdim+1] += u2[idim]*n[jdim];
630                         _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1+idim+1] += u2[jdim]*n[jdim];
631                 }
632         // update \Delta alpha
633         /*
634          * (alpha *rho2*varrho_2+dpi1*(1-alpha)*inv_a2_2*varrho_2)*n[idim]
635          */
636         for (int idim=0; idim<_Ndim; idim++){
637                 _JacoMat[ (1+idim)*_nVar] += dpi1*varrho_2*(1-alp)*inv_a2_2*n[idim];
638                 _JacoMat[ (1+idim)*_nVar+_Ndim+1] += -dpi1*varrho_2*alp*inv_a1_2*n[idim];
639                 _JacoMat[(2+_Ndim+idim)*_nVar] += - dpi2*varrho_2*(1-alp)*inv_a2_2*n[idim];
640                 _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1] += dpi2*varrho_2*alp*inv_a1_2*n[idim];
641                 for  (int i=0; i<_nVar; i++){
642                         _JacoMat[ (1+idim)*_nVar+i]+=dpi1*varrho_2*(-alp*(1-alp)*inv_a2_2*b1*Delta_e1[i]+alp*(1-alp)*inv_a1_2*b2*Delta_e2[i])*n[idim];
643                         _JacoMat[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]+=-dpi2*varrho_2*(-alp*(1-alp)*inv_a2_2*b1*Delta_e1[i]+alp*(1-alp)*inv_a1_2*b2*Delta_e2[i])*n[idim];
644                 }
645         }
646         // update \Delta P
647         for (int idim=0; idim<_Ndim; idim++){
648                 _JacoMat[ (1+idim)*_nVar] += alp*varrho_2*rho2*n[idim];
649                 _JacoMat[ (1+idim)*_nVar+_Ndim+1] +=alp* varrho_2*rho1*n[idim];
650                 _JacoMat[(2+_Ndim+idim)*_nVar] +=  (1-alp)*varrho_2*rho2*n[idim];
651                 _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1] += (1-alp)* varrho_2*rho1*n[idim];
652                 for  (int i=0; i<_nVar; i++){
653                         _JacoMat[ (1+idim)*_nVar+i]+=alp*varrho_2*(-alp*rho2*b1*Delta_e1[i] -(1-alp)*rho1*b2*Delta_e2[i])*n[idim];
654                         _JacoMat[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]+=(1-alp)*varrho_2*(-alp*rho2*b1*Delta_e1[i] -(1-alp)*rho1*b2*Delta_e2[i])*n[idim];
655                 }
656         }
657         // last row (total energy)
658         for (int i=0; i<_nVar; i++){
659                 _JacoMat[(2*_Ndim+2)*_nVar +i] = dF5[i];
660         }
661 }
662
663 void FiveEqsTwoFluid::convectionMatrices()
664 {
665         //entree: URoe = alpha, p, u1, u2, T + ajout dpi
666         //sortie: matrices Roe+  et Roe- +Roe si schéma centre
667
668         if(_timeScheme==Implicit && _usePrimitiveVarsInNewton)
669                 throw CdmathException("Implicitation with primitive variables not yet available for FiveEqsTwoFluid model");
670
671         /*Definitions */
672         complex< double > tmp;
673         double  u_r2 = 0, u1_n=0, u2_n=0;
674         // u1_n = u1.n;u2_n = u2.n (scalar product)
675         for(int i=0;i<_Ndim;i++)
676         {
677                 //u_r2 += (_Uroe[2+i]-_Uroe[1+i+1+_Ndim])*(_Uroe[2+i]-_Uroe[1+i+1+_Ndim]);
678                 u_r2 += 0.5*((_l[2+i]-_l[2+i+_Ndim])*(_l[2+i]-_l[2+i+_Ndim])+(_r[2+i]-_r[2+i+_Ndim])*(_r[2+i]-_r[2+i+_Ndim])); //Kieu
679                 u1_n += _Uroe[2+i]*_vec_normal[i];
680                 u2_n += _Uroe[1+i+1+_Ndim]*_vec_normal[i];
681         }
682
683         //double alpha = (_l[0]+_r[0])*0.5;//Kieu formula
684         double alpha = _Uroe[0];//Toumi formula
685         //double p = _Uroe[1];//Toumi pressure
686
687         /***********Calcul des valeurs propres ********/
688
689         // ********Prepare the parameters to compute the Roe Matrix******** //
690         // **** coefficients eta,  varrho_2 **** //
691         double eta[_Ndim], varrho_2;
692         double rho1_l =  _fluides[0]->getDensity(_l[1], _l[_Ndim*2+2]);//(p,T)_Ndim*2+2
693         double rho2_l =  _fluides[1]->getDensity(_l[1], _l[_Ndim*2+2]);
694         double rho1_r =  _fluides[0]->getDensity(_r[1], _r[_Ndim*2+2]);
695         double rho2_r =  _fluides[1]->getDensity(_r[1], _r[_Ndim*2+2]);
696         // **** coefficients a, b, c **** //
697         double inv_a1_2,inv_a2_2,b1,c1,a2,b2,c2;
698         double e1_l,e1_r,e2_l,e2_r;
699         e1_l = _fluides[0]->getInternalEnergy(_l[_nVar-1],rho1_l);// primitive variable _l[_nVar-1]=Tm
700         e2_l = _fluides[1]->getInternalEnergy(_l[_nVar-1],rho2_l);
701         e1_r = _fluides[0]->getInternalEnergy(_r[_nVar-1],rho1_r);
702         e2_r = _fluides[1]->getInternalEnergy(_r[_nVar-1],rho2_r);
703         inv_a1_2 = static_cast<StiffenedGas*>(_fluides[0])->getJumpDensPress(e1_l,e1_r);
704         inv_a2_2 = static_cast<StiffenedGas*>(_fluides[1])->getJumpDensPress(e2_l,e2_r);
705         //double getJumpDensInternalEnergy(const double p_l,const double p_r,const double e_l,const double e_r);
706         b1 = static_cast<StiffenedGas*>(_fluides[0])->getJumpDensInternalEnergy(_l[1],_r[1],e1_l,e1_r);
707         b2 = static_cast<StiffenedGas*>(_fluides[1])->getJumpDensInternalEnergy(_l[1],_r[1],e2_l,e2_r);
708         //double getJumpInternalEnergyTemperature();
709         c1 = static_cast<StiffenedGas*>(_fluides[0])->getJumpInternalEnergyTemperature();
710         c2 = static_cast<StiffenedGas*>(_fluides[1])->getJumpInternalEnergyTemperature();
711
712         // prefix m is arithmetic mean
713         double m_alp1,m_rho1,m_rho2,m_P,m_e1,m_e2,m_m1,m_m2, eta_n;
714         double m_u1[_Ndim],m_u2[_Ndim], m_alp_u1[_Ndim],m_alp_u2[_Ndim];
715         double u1_l[_Ndim], u2_l[_Ndim],u1_r[_Ndim],u2_r[_Ndim];
716         double u1l_2 = 0, u1r_2 = 0, u2l_2 = 0, u2r_2 = 0;
717         m_alp1 =  0.5*(_l[0]+_r[0]);
718         m_rho1 = 0.5*(rho1_l+rho1_r);
719         m_rho2 = 0.5*(rho2_l+rho2_r);
720         m_P = 0.5*(_l[1]+_r[1]);
721         m_e1 = 0.5*(e1_l+e1_r);
722         m_e2 = 0.5*(e2_l+e2_r);
723         m_m1 = 0.5*(_l[0]*rho1_l+_r[0]*rho1_r);
724         m_m2 = 0.5*((1-_l[0])*rho2_l+(1-_r[0])*rho2_r);
725         varrho_2 =1/((m_alp1*m_rho2)*inv_a1_2+((1-m_alp1)*m_rho1)*inv_a2_2);
726         eta_n = 0.;
727
728         for (int idim=0; idim<_Ndim; idim++){
729                 u1_l[idim] = _l[idim+2];
730                 u1_r[idim] = _r[idim+2];
731                 u2_l[idim] = _l[_Ndim+idim+2];
732                 u2_r[idim] = _r[_Ndim+idim+2];
733                 m_u1[idim] = 0.5*(u1_l[idim] + u1_r[idim]);
734                 m_u2[idim] = 0.5*(u2_l[idim] + u2_r[idim]);
735                 m_alp_u1[idim] = 0.5*(_l[0]*u1_l[idim]+_r[0]*u1_r[idim]);
736                 m_alp_u2[idim] = 0.5*((1-_l[0])*u2_l[idim]+(1-_r[0])*u2_r[idim]);
737                 eta_n += (m_alp_u1[idim]*(1-m_P/m_rho1*inv_a1_2)+m_alp_u2[idim]*(1-m_P/m_rho2*inv_a2_2))*_vec_normal[idim];
738         }
739
740         double eta_varrho_2n = eta_n*varrho_2;
741         //  **** compute jump of Delta T, Delta e1, Delta e2 **** //
742         for  (int idim=0; idim<_Ndim; idim++){
743                 u1l_2 += u1_l[idim]*u1_l[idim];
744                 u1r_2 += u1_r[idim]*u1_r[idim];
745                 u2l_2 += u2_l[idim]*u2_l[idim];
746                 u2r_2 += u2_r[idim]*u2_r[idim];
747         }
748         double inv_m_cm = 1/(c1*m_m1+c2*m_m2);
749         double DeltaT [_nVar], Delta_e1[_nVar], Delta_e2[_nVar];
750         // initialize DeltaT
751         for (int i=0; i<_nVar; i++){
752                 DeltaT[i] = 0.;
753         }
754         DeltaT[0] += -m_e1;
755         DeltaT[1+_Ndim] += -m_e2;
756         DeltaT[_nVar-1] += 1 ;
757         for (int idim=0; idim<_Ndim; idim++){
758                 // wrt mass gas
759                 DeltaT[0] += 0.5*_l[idim+2]     *_r[idim+2];//0.5*\tilde{u_g}^2
760                 // wrt mass liquid
761                 DeltaT[_Ndim+1] += 0.5*_l[_Ndim+idim+2] *_r[_Ndim+idim+2];
762                 // wrt momentum gas
763                 DeltaT[idim+1] += - m_u1[idim];
764                 // wrt momentum liquid
765                 DeltaT[_Ndim+idim+2] += - m_u2[idim];
766         }
767
768         // finalize  DeltaT, Delta_e1 and Delta_e2
769         for (int i =0; i< _nVar; ++i){
770                 DeltaT[i] = inv_m_cm*DeltaT[i];
771                 Delta_e1[i] = c1*DeltaT[i];
772                 Delta_e2[i] = c2*DeltaT[i];
773         }
774
775         // *** compute jump flux (energy equation) A5 *** //
776
777         double A5[_nVar];
778         // initialize A5
779         for (int i=0; i<_nVar; i++){
780                 A5[i]=0;
781         }
782         A5[0] = eta_varrho_2n*m_rho2; // mass gas
783         A5[_Ndim+1] = eta_varrho_2n*m_rho1; // mass liquid
784         for (int idim=0; idim<_Ndim; idim++){
785                 // momentum gas
786                 A5[idim+1] = (m_e1+m_P/m_rho1)*_vec_normal[idim];
787                 // momentum liquid
788                 A5[_Ndim+idim+2] = (m_e2+m_P/m_rho2)*_vec_normal[idim];
789         }
790         // assign the value of A5 (last row of the Roe matrix)
791         for (int idim=0; idim<_Ndim; idim++){
792                 for (int jdim=0; jdim<_Ndim;jdim++){
793                         A5[0] += 0.5*(0.5*(_l[idim+2]*_l[jdim+2]*_l[jdim+2]+_r[idim+2]*_r[jdim+2]*_r[jdim+2])*_vec_normal[idim]);// m_(uin*uj^2)
794                         A5[0] -= m_u1[idim]*m_u1[jdim]*m_u1[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
795                         A5[0] -= 0.5*(m_u1[idim]*_vec_normal[idim]*0.5*(_l[jdim+2]*_l[jdim+2]+_r[jdim+2]*_r[jdim+2]));//-0.5*m_uin*m_uj^2
796                         A5[_Ndim+1]+= 0.5*(0.5*(_l[_Ndim+idim+2]*_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+idim+2]*_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2])*_vec_normal[idim]);// m_(uin*uj^2)
797                         A5[_Ndim+1] -= m_u2[idim]*m_u2[jdim]*m_u2[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
798                         A5[_Ndim+1] -= 0.5*(m_u2[idim]*0.5*(_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2]))*_vec_normal[idim];//-0.5*m_uin*m_uj^2
799                         //momentum gas
800                         A5[idim+1] += m_u1[jdim]*_vec_normal[jdim]*m_u1[idim]+0.5*_vec_normal[idim]*0.5*(_l[jdim+2]*_l[jdim+2]+_r[jdim+2]*_r[jdim+2]);
801                         //momentum liquid
802                         A5[_Ndim+idim+2] += m_u2[jdim]*_vec_normal[jdim]*m_u2[idim]+0.5*_vec_normal[idim]*0.5*(_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2]);
803                 }
804         }
805
806         // final A5
807         double coef_e1, coef_e2;
808         coef_e1 = - eta_varrho_2n*m_alp1*m_rho2*b1;
809         coef_e2 = - eta_varrho_2n*(1-m_alp1)*m_rho1*b2;
810         for (int idim=0; idim<_Ndim; idim++){
811                 coef_e1 += (0.5*(_l[0]*rho1_l*_l[idim+2]+_r[0]*rho1_r*_r[idim+2]) - m_alp_u1[idim]*m_P*b1/m_rho1)*_vec_normal[idim];
812                 coef_e2 += (0.5*((1-_l[0])*rho2_l*_l[_Ndim+idim+2]+(1-_r[0])*rho2_r*_r[_Ndim+idim+2])-m_alp_u2[idim]*m_P*b2/m_rho2)*_vec_normal[idim];
813         }
814         for (int i =0; i< _nVar; i++){
815                 A5[i] += coef_e1*Delta_e1[i] + coef_e2*Delta_e2[i];
816         }
817         // ******* Construction de la matrice de Roe ******** //
818         // interfacial pressure correction
819         double T=_Uroe[_nVar-1];
820         double dpi1 = intPressDef(alpha, u_r2, m_rho1,m_rho2,T);
821         double dpi2 = dpi1;
822         //saving dpi value for flux calculation later
823         _Uroe[_nVar]=dpi1 ;
824
825         //lignes de masse
826         for(int i=0; i<_nVar*_nVar;i++)
827                 _Aroe[i]=0;
828         //      alpha = 0.; dpi1 = 0.; dpi2 = 0.;
829         for(int idim=0; idim<_Ndim;idim++)
830         {
831                 _Aroe[1+idim]=_vec_normal[idim];
832                 _Aroe[1+idim+_Ndim+1]=0;
833                 _Aroe[(_Ndim+1)*_nVar+1+idim]=0;
834                 _Aroe[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=_vec_normal[idim];
835         }
836
837         // Take into account the convection term in the momentum eqts
838         for(int idim=0; idim<_Ndim;idim++)
839                 for (int jdim=0; jdim<_Ndim;jdim++){
840                         // momentum gas (neglect delta alpha and delta P)
841                         _Aroe[ (1+idim)*_nVar] += (0.5*(_l[2+idim]*_l[2+jdim]+_r[2+idim]*_r[2+jdim])-2*m_u1[idim]*m_u1[jdim])*_vec_normal[jdim];
842                         _Aroe[(1+idim)*_nVar+jdim+1] += m_u1[idim]*_vec_normal[jdim];
843                         _Aroe[(1+idim)*_nVar+idim+1] += m_u1[jdim]*_vec_normal[jdim];
844                         // momentum liquid (neglect delta alpha and delta P)
845                         _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += (0.5*(_l[_Ndim+2+idim]*_l[_Ndim+2+jdim]+_r[_Ndim+2+idim]*_r[_Ndim+2+jdim])-2*m_u2[idim]*m_u2[jdim])*_vec_normal[jdim];
846                         _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1+jdim+1] += m_u2[idim]*_vec_normal[jdim];
847                         _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1+idim+1] += m_u2[jdim]*_vec_normal[jdim];
848                 }
849
850         // update \Delta alpha
851         for (int idim=0; idim<_Ndim; idim++){
852                 _Aroe[ (1+idim)*_nVar] += dpi1*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
853                 _Aroe[ (1+idim)*_nVar+_Ndim+1] += -dpi1*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
854                 _Aroe[(2+_Ndim+idim)*_nVar] += - dpi2*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
855                 _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += dpi2*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
856                 for  (int i=0; i<_nVar; i++){
857                         _Aroe[ (1+idim)*_nVar+i]+=dpi1*varrho_2*(-m_alp1*(1-m_alp1)*inv_a2_2*b1*Delta_e1[i]+m_alp1*(1-m_alp1)*inv_a1_2*b2*Delta_e2[i])*_vec_normal[idim];
858                         _Aroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]+=-dpi2*varrho_2*(-m_alp1*(1-m_alp1)*inv_a2_2*b1*Delta_e1[i]+m_alp1*(1-m_alp1)*inv_a1_2*b2*Delta_e2[i])*_vec_normal[idim];
859                 }
860         }
861         // update \Delta P
862         for (int idim=0; idim<_Ndim; idim++){
863                 _Aroe[ (1+idim)*_nVar] += alpha*varrho_2*m_rho2*_vec_normal[idim];
864                 _Aroe[ (1+idim)*_nVar+_Ndim+1] += alpha* varrho_2*m_rho1*_vec_normal[idim];
865                 _Aroe[(2+_Ndim+idim)*_nVar] += (1-alpha)*varrho_2*m_rho2*_vec_normal[idim];
866                 _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += (1-alpha)* varrho_2*m_rho1*_vec_normal[idim];
867                 for  (int i=0; i<_nVar; i++){
868                         _Aroe[ (1+idim)*_nVar+i] += alpha*varrho_2*(-m_alp1*m_rho2*b1*Delta_e1[i] -(1-m_alp1)*m_rho1*b2*Delta_e2[i])*_vec_normal[idim];
869                         _Aroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i] += (1-alpha)*varrho_2*(-m_alp1*m_rho2*b1*Delta_e1[i] -(1-m_alp1)*m_rho1*b2*Delta_e2[i])*_vec_normal[idim];
870                 }
871         }
872
873         // last row (total energy)
874         for (int i=0; i<_nVar; i++){
875                 _Aroe[(2*_Ndim+2)*_nVar +i] += A5[i];
876         }
877
878         int LDA = _nVar;
879         int LDVL = _nVar;
880         int LDVR = _nVar;
881         int LWORK = 50*_nVar;
882         char jobvl[]="N", jobvr[]="N";
883         double WORK[LWORK], Aroe[_nVar*_nVar],egvaReal[_nVar],egvaImag[_nVar],
884         egVectorL[_nVar*_nVar],egVectorR[_nVar*_nVar];
885         int info=0;
886         double maxvploc=0;
887         std::vector< std::complex<double> > valeurs_propres_dist;
888         int taille_vp ;
889
890         for (int i=0; i<_nVar*_nVar; i++)
891                 Aroe[i] = _Aroe[i];
892
893         if (_verbose && _nbTimeStep%_freqSave ==0)
894         {
895                 cout<<endl<<"Matrice de Roe"<<endl;
896                 for(int i=0; i<_nVar;i++)
897                 {
898                         for(int j=0; j<_nVar;j++)
899                                 cout << _Aroe[i*_nVar+j]<< " , ";
900                         cout<<endl;
901                 }
902         }
903         /******** Compute the eigenvalues and eigenvectors of Roe Matrix (using lapack)*********/
904         Polynoms Poly;
905         dgeev_(jobvl, jobvr,  &_nVar,
906                         Aroe,&LDA,egvaReal,egvaImag, egVectorL,
907                         &LDVL,egVectorR,
908                         &LDVR, WORK,&LWORK,
909                         &info);
910         if(info <0){
911                 cout<<"FiveEqsTwoFluid::convectionMatrices: error dgeev_ : argument "<<-info<<" invalid"<<endl;
912                 *_runLogFile<<"FiveEqsTwoFluid::convectionMatrices: error dgeev_ : argument "<<-info<<" invalid"<<endl;
913                 throw CdmathException("FiveEqsTwoFluid::convectionMatrices: dgeev_ unsuccessful computation of the eigenvalues ");
914         }
915         else if(info >0)
916         {
917                 cout<<"Warning FiveEqsTwoFluid::convectionMatrices: dgeev_ did not compute all the eigenvalues, trying Rusanov scheme "<<endl;
918                 cout<<"Converged eigenvalues are ";
919                 for(int i=info; i<_nVar; i++)
920                         cout<<"("<< egvaReal[i]<<","<< egvaImag[i]<<"), ";
921                 cout<<endl;
922
923                 _maxvploc=0;
924                 for(int i =info; i<_nVar; i++)
925                         if (fabs(egvaReal[i])>_maxvploc)
926                                 _maxvploc=fabs(egvaReal[i]);
927                 if(_maxvploc>_maxvp)
928                         _maxvp=_maxvploc;
929
930                 valeurs_propres_dist=std::vector< std::complex<double> > (1,maxvploc);
931                 taille_vp =1;
932         }
933         else{
934                 if (_verbose && _nbTimeStep%_freqSave ==0)
935                 {
936                         for(int i=0; i<_nVar; i++)
937                                 cout<<" Vp real part " << egvaReal[i]<<", Imaginary part " << egvaImag[i]<<endl;
938                 }
939
940                 std::vector< std::complex<double> > valeurs_propres(_nVar);
941
942                 bool complexEigenvalues = false;
943                 for(int j=0; j<_nVar; j++){
944                         if (max(_l[0],_r[0])<_precision && abs(egvaImag[j])<_precision )// Kieu test Monophase
945                                 egvaImag[j] = 0.;
946                         else
947                                 if (egvaImag[j] >_precision){// Kieu
948                                         complexEigenvalues = true;
949                                 }
950                         if (abs(_l[0]-_r[0])<_precision*_precision && fabs(egvaImag[j])<_precision)// Kieu interfacial pressure
951                                 egvaImag[j] = 0.;
952                         valeurs_propres[j] = complex<double>(egvaReal[j],egvaImag[j]);
953                 }
954
955                 taille_vp =Poly.new_tri_selectif(valeurs_propres,valeurs_propres.size(),_precision);
956
957                 valeurs_propres_dist=vector< complex< double > >(taille_vp);
958                 for( int i=0 ; i<taille_vp ; i++)
959                         valeurs_propres_dist[i] = valeurs_propres[i];
960                 if(_verbose && _nbTimeStep%_freqSave ==0)
961                 {
962                         cout<<" Vp apres tri " << valeurs_propres_dist.size()<<endl;
963                         for(int ct =0; ct<taille_vp; ct++)
964                                 cout<< "("<<valeurs_propres_dist[ct].real()<< ", " <<valeurs_propres_dist[ct].imag() <<")  ";
965                         cout<< endl;
966                 }
967
968                 for(int i =0; i<taille_vp; i++){
969                         if (fabs(valeurs_propres_dist[i].real())>maxvploc)
970                                 maxvploc=fabs(valeurs_propres_dist[i].real());
971                         if(maxvploc>_maxvp)
972                                 _maxvp=maxvploc;
973                 }
974
975                 int existVpCplx = 0,pos_conj;
976                 double vp_imag_iter;
977                 for (int ct=0; ct<taille_vp; ct++) {
978                         vp_imag_iter = valeurs_propres_dist[ct].imag();
979                         if ( fabs(vp_imag_iter) > 100*_precision ) {
980                                 existVpCplx +=1;
981                                 if ( _part_imag_max < fabs(vp_imag_iter))
982                                         _part_imag_max = fabs(vp_imag_iter);
983                                 //On cherhe le conjugue
984                                 pos_conj = ct+1;
985                                 while(pos_conj<taille_vp && fabs(valeurs_propres_dist[pos_conj].imag()+vp_imag_iter)>_precision)
986                                         pos_conj++;
987                                 if(pos_conj!=ct+1 && pos_conj<taille_vp )
988                                 {
989                                         tmp=valeurs_propres_dist[ct+1];
990                                         valeurs_propres_dist[ct+1]=valeurs_propres_dist[pos_conj];
991                                         valeurs_propres_dist[pos_conj] = tmp;
992                                         ct++;
993                                 }
994                         }
995                 }
996                 if (existVpCplx >0)
997                         _nbVpCplx +=1;
998         }
999         /******* Construction des matrices de decentrement *******/
1000         if(_spaceScheme == centered )
1001         {
1002                 if(_entropicCorrection)
1003                         throw CdmathException("FiveEqsTwoFluid::convectionMatrices: entropic scheme not available for centered scheme");
1004                 for(int i=0; i<_nVar*_nVar;i++)
1005                         _absAroe[i]=0;
1006                 //  if(alpha<_precision || alpha>1-_precision)//rusanov
1007                 //          for(int i=0; i<_nVar;i++)
1008                 //            _absAroe[i*_nVar+i]=maxvploc;
1009         }
1010         if( _spaceScheme ==staggered){//To do: study entropic correction for staggered
1011                 if(_entropicCorrection)//To do: study entropic correction for staggered
1012                         throw CdmathException("FiveEqsTwoFluid::convectionMatrices: entropic scheme not yet available for staggered scheme");
1013                 /******** Construction de la matrice de decentrement staggered *********/
1014                 /***** Compute eta_n **************/
1015                 eta_n = 0.;
1016                 for (int idim=0; idim<_Ndim; idim++){
1017                         eta_n += (m_alp_u1[idim]*(-1-m_P/m_rho1*inv_a1_2)+m_alp_u2[idim]*(-1-m_P/m_rho2*inv_a2_2))*_vec_normal[idim];
1018                 }
1019                 double eta_varrho_2n = eta_n*varrho_2;
1020                 /**** compute jump flux (energy equation) A5 ****/
1021
1022                 double A5[_nVar];
1023                 // initialize A5
1024                 for (int i=0; i<_nVar; i++){
1025                         A5[i]=0;
1026                 }
1027                 A5[0] = eta_varrho_2n*m_rho2; // mass gas
1028                 A5[_Ndim+1] = eta_varrho_2n*m_rho1; // mass liquid
1029                 // assign the value of A5 (last row of the Roe matrix)
1030                 for (int idim=0; idim<_Ndim; idim++){
1031                         for (int jdim=0; jdim<_Ndim;jdim++){
1032                                 A5[0] += 0.5*(0.5*(_l[idim+2]*_l[jdim+2]*_l[jdim+2]+_r[idim+2]*_r[jdim+2]*_r[jdim+2])*_vec_normal[idim]);// m_(uin*uj^2)
1033                                 A5[0] -= m_u1[idim]*m_u1[jdim]*m_u1[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
1034                                 A5[0] -= 0.5*(m_u1[idim]*_vec_normal[idim]*0.5*(_l[jdim+2]*_l[jdim+2]+_r[jdim+2]*_r[jdim+2]));//-0.5*m_uin*m_uj^2
1035                                 A5[_Ndim+1]+= 0.5*(0.5*(_l[_Ndim+idim+2]*_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+idim+2]*_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2])*_vec_normal[idim]);// m_(uin*uj^2)
1036                                 A5[_Ndim+1] -= m_u2[idim]*m_u2[jdim]*m_u2[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
1037                                 A5[_Ndim+1] -= 0.5*(m_u2[idim]*_vec_normal[idim]*0.5*(_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2]));//-0.5*m_uin*m_uj^2
1038                         }
1039                 }
1040                 // final A5
1041                 double coef_e1, coef_e2;
1042                 coef_e1 = - eta_varrho_2n*m_alp1*m_rho2*b1;
1043                 coef_e2 = - eta_varrho_2n*(1-m_alp1)*m_rho1*b2;
1044                 for (int idim=0; idim<_Ndim; idim++){
1045                         coef_e1 += (0.5*(_l[0]*rho1_l*_l[idim+2]+_r[0]*rho1_r*_r[idim+2]) - m_alp_u1[idim]*m_P*b1/m_rho1)*_vec_normal[idim];
1046                         coef_e2 += (0.5*((1-_l[0])*rho2_l*_l[_Ndim+idim+2]+(1-_r[0])*rho2_r*_r[_Ndim+idim+2])-m_alp_u2[idim]*m_P*b2/m_rho2)*_vec_normal[idim];
1047                 }
1048                 for (int i =0; i< _nVar; i++){
1049                         A5[i] += coef_e1*Delta_e1[i] + coef_e2*Delta_e2[i];
1050                 }
1051                 /** Début remplissage matrice décentrement staggered **/
1052                 //lignes de masse
1053                 for(int i=0; i<_nVar*_nVar;i++)
1054                         _absAroe[i]=0;
1055
1056                 //Mass lines
1057                 for(int idim=0; idim<_Ndim;idim++)
1058                 {
1059                         _absAroe[1+idim]=_vec_normal[idim];
1060                         _absAroe[1+idim+_Ndim+1]=0;
1061                         _absAroe[(_Ndim+1)*_nVar+1+idim]=0;
1062                         _absAroe[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=_vec_normal[idim];
1063                 }
1064                 //Contribution of convection (rho u\times u) in the momentum equations
1065                 for(int idim=0; idim<_Ndim;idim++)
1066                         for (int jdim=0; jdim<_Ndim;jdim++){
1067                                 // momentum gas (neglect delta alpha and delta P)
1068                                 _absAroe[ (1+idim)*_nVar] += (0.5*(_l[2+idim]*_l[2+jdim]+_r[2+idim]*_r[2+jdim])-2*m_u1[idim]*m_u1[jdim])*_vec_normal[jdim];
1069                                 _absAroe[(1+idim)*_nVar+jdim+1] += m_u1[idim]*_vec_normal[jdim];
1070                                 _absAroe[(1+idim)*_nVar+idim+1] += m_u1[jdim]*_vec_normal[jdim];
1071                                 // momentum liquid (neglect delta alpha and delta P)
1072                                 _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += (0.5*(_l[_Ndim+2+idim]*_l[_Ndim+2+jdim]+_r[_Ndim+2+idim]*_r[_Ndim+2+jdim])-2*m_u2[idim]*m_u2[jdim])*_vec_normal[jdim];
1073                                 _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1+jdim+1] += m_u2[idim]*_vec_normal[jdim];
1074                                 _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1+idim+1] += m_u2[jdim]*_vec_normal[jdim];
1075                         }
1076                 // contribution of interfacial pressure in momentum equation
1077                 /*
1078                  * (alpha *rho2*varrho_2+dpi1*(1-alpha)*inv_a2_2*varrho_2)*_vec_normal[idim]
1079                  */
1080                 for (int idim=0; idim<_Ndim; idim++){
1081                         _absAroe[ (1+idim)*_nVar] += dpi1*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
1082                         _absAroe[ (1+idim)*_nVar+_Ndim+1] += -dpi1*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
1083                         _absAroe[(2+_Ndim+idim)*_nVar] += - dpi2*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
1084                         _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += dpi2*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
1085                         for  (int i=0; i<_nVar; i++){
1086                                 _absAroe[ (1+idim)*_nVar+i]+=dpi1*varrho_2*(-m_alp1*inv_a2_2*b1*Delta_e1[i]+(1-m_alp1)*inv_a1_2*b2*Delta_e2[i])*_vec_normal[idim];
1087                                 _absAroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]+=-dpi2*varrho_2*(-m_alp1*inv_a2_2*b1*Delta_e1[i]+(1-m_alp1)*inv_a1_2*b2*Delta_e2[i])*_vec_normal[idim];
1088                         }
1089                 }
1090                 // contribution of pressure gradient in momentum equation
1091                 for (int idim=0; idim<_Ndim; idim++){
1092                         _absAroe[ (1+idim)*_nVar] -= alpha*varrho_2*m_rho2*_vec_normal[idim];
1093                         _absAroe[ (1+idim)*_nVar+_Ndim+1] -=alpha* varrho_2*m_rho1*_vec_normal[idim];
1094                         _absAroe[(2+_Ndim+idim)*_nVar] -=  (1-alpha)*varrho_2*m_rho2*_vec_normal[idim];
1095                         _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1] -= (1-alpha)* varrho_2*m_rho1*_vec_normal[idim];
1096                         for  (int i=0; i<_nVar; i++){
1097                                 _absAroe[ (1+idim)*_nVar+i]-=alpha*varrho_2*(-m_alp1*m_rho2*b1*Delta_e1[i] -(1-m_alp1)*m_rho1*b2*Delta_e2[i])*_vec_normal[idim];
1098                                 _absAroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]-=(1-alpha)*varrho_2*(-m_alp1*m_rho2*b1*Delta_e1[i] -(1-m_alp1)*m_rho1*b2*Delta_e2[i])*_vec_normal[idim];
1099                         }
1100                 }
1101                 // last row (total energy) To be changed soon
1102                 for (int i=0; i<_nVar; i++){
1103                         _Aroe[(2*_Ndim+2)*_nVar +i] = A5[i];
1104                 }
1105                 double signu1,signu2;
1106                 if(u1_n>0)
1107                         signu1=1;
1108                 else if (u1_n<0)
1109                         signu1=-1;
1110                 else
1111                         signu1=0;
1112                 if(u2_n>0)
1113                         signu2=1;
1114                 else if (u2_n<0)
1115                         signu2=-1;
1116                 else
1117                         signu2=0;
1118
1119                 for(int i=0; i<(1+_Ndim)*_nVar;i++){
1120                         _absAroe[i] *= signu1;
1121                         _absAroe[i+(1+_Ndim)*_nVar] *= signu2;
1122                 }
1123
1124         }
1125         if(_spaceScheme ==upwind || _spaceScheme==pressureCorrection || _spaceScheme ==lowMach)//calcul de la valeur absolue
1126         {               //on ordonne les deux premieres valeurs
1127                 if(valeurs_propres_dist[1].real()<valeurs_propres_dist[0].real())
1128                 {
1129                         tmp=valeurs_propres_dist[0];
1130                         valeurs_propres_dist[0]=valeurs_propres_dist[1];
1131                         valeurs_propres_dist[1]=tmp;
1132                 }
1133                 vector< complex< double > > y (taille_vp,0);
1134                 for( int i=0 ; i<taille_vp ; i++)
1135                         y[i] = Poly.abs_generalise(valeurs_propres_dist[i]);
1136
1137                 if(_entropicCorrection)
1138                 {
1139                         double entShift0 = 0;
1140                         double entShift1 = 0;
1141                         entropicShift(_vec_normal,entShift0,entShift1);
1142                         //cout<<"entShift0= "<<entShift0<<endl;
1143                         for( int i=0 ; i<taille_vp ; i++)
1144                         {
1145                                 //cout<<"y["<<i<<"]="<<y[i].real()<<endl;
1146                                 y[i] += max(entShift0,entShift1);
1147                         }
1148                 }
1149
1150                 if(_verbose && _nbTimeStep%_freqSave ==0)
1151                 {
1152                         cout<<"valeurs propres"<<endl;
1153                         for( int i=0 ; i<taille_vp ; i++)
1154                                 cout<<valeurs_propres_dist[i] <<", "<<endl;
1155                         cout<<"valeurs à atteindre"<<endl;
1156                         for( int i=0 ; i<taille_vp ; i++)
1157                                 cout<<y[i] <<", "<<endl;
1158                 }
1159                 Poly.abs_par_interp_directe(taille_vp,valeurs_propres_dist, _Aroe, _nVar,_precision, _absAroe,y);
1160
1161                 if( _spaceScheme ==pressureCorrection){
1162                         for( int i=0 ; i<_Ndim ; i++)
1163                                 for( int j=0 ; j<_Ndim ; j++){
1164                                         _absAroe[(1+i)*_nVar+1+j]-=alpha*(valeurs_propres_dist[1].real()-valeurs_propres_dist[0].real())/2*_vec_normal[i]*_vec_normal[j];
1165                                         _absAroe[(2+_Ndim+i)*_nVar+2+_Ndim+j]-=(1-alpha)*(valeurs_propres_dist[1].real()-valeurs_propres_dist[0].real())/2*_vec_normal[i]*_vec_normal[j];
1166                                 }
1167                 }
1168                 else if( _spaceScheme ==lowMach){
1169                         double M=max(fabs(u1_n),fabs(u2_n))/maxvploc;
1170                         for( int i=0 ; i<_Ndim ; i++)
1171                                 for( int j=0 ; j<_Ndim ; j++){
1172                                         _absAroe[(1+i)*_nVar+1+j]-=(1-M)*alpha*(valeurs_propres_dist[1].real()-valeurs_propres_dist[0].real())/2*_vec_normal[i]*_vec_normal[j];
1173                                         _absAroe[(2+_Ndim+i)*_nVar+2+_Ndim+j]-=(1-M)*(1-alpha)*(valeurs_propres_dist[1].real()-valeurs_propres_dist[0].real())/2*_vec_normal[i]*_vec_normal[j];
1174                                 }
1175                 }
1176         }
1177
1178         //Calcul de la matrice signe pour VFFC, VFRoe et décentrement des termes source
1179
1180         if(_entropicCorrection || _spaceScheme ==pressureCorrection || _spaceScheme ==lowMach){
1181                 InvMatriceRoe( valeurs_propres_dist);
1182                 Poly.matrixProduct(_absAroe, _nVar, _nVar, _invAroe, _nVar, _nVar, _signAroe);
1183         }
1184         else if (_spaceScheme==upwind)//upwind sans entropic
1185                 SigneMatriceRoe( valeurs_propres_dist);
1186         else if (_spaceScheme==centered)//centre  sans entropic
1187         {
1188                 for(int i=0; i<_nVar*_nVar;i++)
1189                         _signAroe[i] = 0;
1190         }
1191         else  if(_spaceScheme ==staggered )
1192         {
1193                 double signu1,signu2;
1194                 if(u1_n>0)
1195                         signu1=1;
1196                 else if (u1_n<0)
1197                         signu1=-1;
1198                 else
1199                         signu1=0;
1200                 if(u2_n>0)
1201                         signu2=1;
1202                 else if (u2_n<0)
1203                         signu2=-1;
1204                 else
1205                         signu2=0;
1206                 for(int i=0; i<_nVar*_nVar;i++)
1207                         _signAroe[i] = 0;
1208                 _signAroe[0] = signu1;
1209                 _signAroe[(1+_Ndim)*_nVar +1+_Ndim] = signu2;
1210                 for(int i=2; i<_nVar-1;i++){
1211                         _signAroe[i*_nVar+i] = -signu1;
1212                         _signAroe[(i+1+_Ndim)*_nVar+i+1+_Ndim] = -signu2;
1213                 }
1214                 //_signAroe[_nVar*(_nVar-1)+_nVar-1] = signu;
1215         }
1216         else
1217                 throw CdmathException("FiveEqsTwoFluid::convectionMatrices: well balanced option not treated");
1218
1219         for(int i=0; i<_nVar*_nVar;i++)
1220         {
1221                 _AroeMinus[i] = (_Aroe[i]-_absAroe[i])/2;
1222                 _AroePlus[i]  = (_Aroe[i]+_absAroe[i])/2;
1223         }
1224         if(_timeScheme==Implicit && _usePrimitiveVarsInNewton)//Implicitation using primitive variables
1225                 for(int i=0; i<_nVar*_nVar;i++)
1226                         _AroeMinusImplicit[i] = (_AroeImplicit[i]-_absAroeImplicit[i])/2;
1227         else
1228                 for(int i=0; i<_nVar*_nVar;i++)
1229                         _AroeMinusImplicit[i] = _AroeMinus[i];
1230
1231         if(_verbose && _nbTimeStep%_freqSave ==0)
1232         {
1233                 cout<<"Matrice de Roe"<<endl;
1234                 for(int i=0; i<_nVar;i++){
1235                         for(int j=0; j<_nVar;j++)
1236                                 cout<<_Aroe[i*_nVar+j]<<" , ";
1237                         cout<<endl;
1238                 }
1239                 cout<<"Valeur absolue matrice de Roe"<<endl;
1240                 for(int i=0; i<_nVar;i++){
1241                         for(int j=0; j<_nVar;j++)
1242                                 cout<<_absAroe[i*_nVar+j]<<" , ";
1243                         cout<<endl;
1244                 }
1245                 cout<<"Signe matrice de Roe"<<endl;
1246                 for(int i=0; i<_nVar;i++){
1247                         for(int j=0; j<_nVar;j++)
1248                                 cout<<_signAroe[i*_nVar+j]<<" , ";
1249                         cout<<endl;
1250                 }
1251                 cout<<endl<<"Matrice _AroeMinus"<<endl;
1252                 for(int i=0; i<_nVar;i++)
1253                 {
1254                         for(int j=0; j<_nVar;j++)
1255                                 cout << _AroeMinus[i*_nVar+j]<< " , ";
1256                         cout<<endl;
1257                 }
1258                 cout<<endl<<"Matrice _AroePlus"<<endl;
1259                 for(int i=0; i<_nVar;i++)
1260                 {
1261                         for(int j=0; j<_nVar;j++)
1262                                 cout << _AroePlus[i*_nVar+j]<< " , ";
1263                         cout<<endl;
1264                 }
1265         }
1266 }
1267
1268 void FiveEqsTwoFluid::jacobianDiff(const int &j, string nameOfGroup)
1269 {
1270         int k;
1271         for(k=0; k<_nVar*_nVar;k++)
1272                 _JcbDiff[k] = 0;
1273         if (_limitField[nameOfGroup].bcType==Wall){
1274                 _idm[0] = _nVar*j;
1275                 for(k=1; k<_nVar; k++)
1276                         _idm[k] = _idm[k-1] + 1;
1277                 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
1278                 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
1279
1280                 double pression=_Vj[1];//pressure inside
1281                 double T=_Vj[_nVar-1];//temperature outside
1282                 double rho_v=_fluides[0]->getDensity(pression,T);
1283                 double rho_l=_fluides[1]->getDensity(pression,T);
1284                 _JcbDiff[0] = 1;
1285                 _JcbDiff[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1286                 _JcbDiff[_nVar]=_limitField[nameOfGroup].v_x[0];
1287                 _JcbDiff[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1288                 double v2_v=_limitField[nameOfGroup].v_x[0]*_limitField[nameOfGroup].v_x[0];
1289                 double v2_l=_limitField[nameOfGroup].v_x[1]*_limitField[nameOfGroup].v_x[1];
1290                 if(_Ndim>1)
1291                 {
1292                         _JcbDiff[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1293                         _JcbDiff[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1294                         v2_v+=_limitField[nameOfGroup].v_y[0]*_limitField[nameOfGroup].v_y[0];
1295                         v2_l+=_limitField[nameOfGroup].v_y[1]*_limitField[nameOfGroup].v_y[1];
1296                         if(_Ndim==3)
1297                         {
1298                                 _JcbDiff[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1299                                 _JcbDiff[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1300                                 v2_v+=_limitField[nameOfGroup].v_z[0]*_limitField[nameOfGroup].v_z[0];
1301                                 v2_l+=_limitField[nameOfGroup].v_z[1]*_limitField[nameOfGroup].v_z[1];
1302                         }
1303                 }
1304                 _JcbDiff[(_nVar-1)*_nVar]=         _fluides[0]->getInternalEnergy(_limitField[nameOfGroup].T,rho_v)+0.5*v2_v;
1305                 _JcbDiff[(_nVar-1)*_nVar +1+_Ndim]=_fluides[1]->getInternalEnergy(_limitField[nameOfGroup].T,rho_l)+0.5*v2_l;
1306         }
1307         else if (_limitField[nameOfGroup].bcType==Inlet){
1308                 /*
1309                 _JcbDiff[0] = 1;
1310                 _JcbDiff[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1311                 _JcbDiff[_nVar]=_limitField[nameOfGroup].v_x[0];
1312                 _JcbDiff[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1313                 double v2_v=_limitField[nameOfGroup].v_x[0]*_limitField[nameOfGroup].v_x[0];
1314                 double v2_l=_limitField[nameOfGroup].v_x[1]*_limitField[nameOfGroup].v_x[1];
1315                 if(_Ndim>1)
1316                 {
1317                         _JcbDiff[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1318                         _JcbDiff[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1319                         v2_v+=_limitField[nameOfGroup].v_y[0]*_limitField[nameOfGroup].v_y[0];
1320                         v2_l+=_limitField[nameOfGroup].v_y[1]*_limitField[nameOfGroup].v_y[1];
1321                         if(_Ndim==3)
1322                         {
1323                                 _JcbDiff[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1324                                 _JcbDiff[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1325                                 v2_v+=_limitField[nameOfGroup].v_z[0]*_limitField[nameOfGroup].v_z[0];
1326                                 v2_l+=_limitField[nameOfGroup].v_z[1]*_limitField[nameOfGroup].v_z[1];
1327                         }
1328                 }
1329                 _JcbDiff[(_nVar-1)*_nVar]=         _fluides[0]->getInternalEnergy(_limitField[nameOfGroup].T,rho_v)+0.5*v2_v;
1330                 _JcbDiff[(_nVar-1)*_nVar +1+_Ndim]=_fluides[1]->getInternalEnergy(_limitField[nameOfGroup].T,rho_l)+0.5*v2_l;
1331                  */
1332         } else if (_limitField[nameOfGroup].bcType==Outlet){
1333                 //extraction de l etat courant et primitives
1334                 /*
1335                 _idm[0] = j*_nVar;
1336                 for(k=1; k<_nVar;k++)
1337                 {_idm[k] = _idm[k-1] + 1;}
1338                 VecGetValues(_conservativeVars, _nVar, _idm, _phi);
1339                 VecGetValues(_primitiveVars, _nVar, _idm, _externalStates);
1340                  */
1341         }
1342         else if (_limitField[nameOfGroup].bcType!=Neumann && _limitField[nameOfGroup].bcType!=InletPressure){
1343                 cout<<"Condition  limite non traitee pour le bord "<<nameOfGroup<< endl;
1344                 throw CdmathException("FiveEqsTwoFluid::jacobianDiff: Condition  limite non traitee");
1345         }
1346 }
1347
1348
1349 void FiveEqsTwoFluid::setBoundaryState(string nameOfGroup, const int &j,double *normale){
1350         //To do controle signe des vitesses pour CL entree/sortie
1351         int k;
1352         double v1_2=0, v2_2=0, q1_n=0, q2_n=0, u1_n=0, u2_n=0;//quantités de mouvement et vitesses normales à la face limite;
1353         double v1[_Ndim], v2[_Ndim];
1354
1355         _idm[0] = _nVar*j;
1356         for(k=1; k<_nVar; k++)
1357                 _idm[k] = _idm[k-1] + 1;
1358
1359         VecGetValues(_conservativeVars, _nVar, _idm, _externalStates);//On initialise l'état fantôme avec l'état interne
1360         VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
1361
1362         for(k=0; k<_Ndim; k++){
1363                 q1_n+=_externalStates[(k+1)]*normale[k];
1364                 q2_n+=_externalStates[(k+1+1+_Ndim)]*normale[k];
1365                 u1_n+=_Vj[(k+2)]*normale[k];
1366                 u2_n+=_Vj[(k+2+_Ndim)]*normale[k];
1367         }
1368
1369         if(_verbose && _nbTimeStep%_freqSave ==0)
1370         {
1371                 cout << "Boundary conditions for group "<< nameOfGroup<< ", inner cell j= "<<j << " face unit normal vector "<<endl;
1372                 for(k=0; k<_Ndim; k++){
1373                         cout<<normale[k]<<", ";
1374                 }
1375                 cout<<endl;
1376         }
1377
1378         if (_limitField[nameOfGroup].bcType==Wall){
1379
1380                 //Pour la convection, inversion du sens des vitesses
1381                 for(k=0; k<_Ndim; k++){
1382                         _externalStates[(k+1)]-= 2*q1_n*normale[k];
1383                         _externalStates[(k+1+1+_Ndim)]-= 2*q2_n*normale[k];
1384                         _Vj[(k+2)]-= 2*u1_n*normale[k];
1385                         _Vj[(k+2+_Ndim)]-= 2*u2_n*normale[k];
1386                 }
1387                 _idm[0] = 0;
1388                 for(k=1; k<_nVar; k++)
1389                         _idm[k] = _idm[k-1] + 1;
1390
1391                 VecAssemblyBegin(_Uext);
1392                 VecAssemblyBegin(_Vext);
1393                 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1394                 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1395                 VecAssemblyEnd(_Uext);
1396                 VecAssemblyEnd(_Vext);
1397
1398                 //Pour la diffusion, paroi à vitesses et temperature imposees
1399                 double pression=_Vj[1];//pressure inside
1400                 double T=_Vj[_nVar-1];//temperature outside
1401                 double rho_v=_fluides[0]->getDensity(pression,T);
1402                 double rho_l=_fluides[1]->getDensity(pression,T);
1403                 _externalStates[1]=_externalStates[0]*_limitField[nameOfGroup].v_x[0];
1404                 _externalStates[2+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_x[1];
1405                 if(_Ndim>1)
1406                 {
1407                         _externalStates[2]=_externalStates[0]*_limitField[nameOfGroup].v_y[0];
1408                         _externalStates[3+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_y[1];
1409                         if(_Ndim==3)
1410                         {
1411                                 _externalStates[3]=_externalStates[0]*_limitField[nameOfGroup].v_z[0];
1412                                 _externalStates[4+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_z[1];
1413                         }
1414                 }
1415                 _externalStates[_nVar-1] = _externalStates[0]*(_fluides[0]->getInternalEnergy(_limitField[nameOfGroup].T,rho_v) + v1_2/2)
1416                                                                                                                                                                                                                                                                                                                                                                         +_externalStates[1+_Ndim]*(_fluides[1]->getInternalEnergy(_limitField[nameOfGroup].T,rho_l) + v2_2/2);
1417                 VecAssemblyBegin(_Uextdiff);
1418                 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1419                 VecAssemblyEnd(_Uextdiff);
1420         }
1421         else if (_limitField[nameOfGroup].bcType==Neumann){
1422                 _idm[0] = _nVar*j;
1423                 for(k=1; k<_nVar; k++)
1424                         _idm[k] = _idm[k-1] + 1;
1425
1426                 _idm[0] = 0;
1427                 for(k=1; k<_nVar; k++)
1428                         _idm[k] = _idm[k-1] + 1;
1429
1430                 VecAssemblyBegin(_Uext);
1431                 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1432                 VecAssemblyEnd(_Uext);
1433
1434                 VecAssemblyBegin(_Vext);
1435                 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1436                 VecAssemblyEnd(_Vext);
1437
1438                 VecAssemblyBegin(_Uextdiff);
1439                 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1440                 VecAssemblyEnd(_Uextdiff);
1441         }
1442         else if (_limitField[nameOfGroup].bcType==Inlet){
1443                 _idm[0] = _nVar*j;
1444                 for(k=1; k<_nVar; k++)
1445                         _idm[k] = _idm[k-1] + 1;
1446
1447                 double alpha=_limitField[nameOfGroup].alpha;//void fraction outside
1448                 double pression=_Vj[1];//pressure inside
1449                 double T=_limitField[nameOfGroup].T;//temperature outside
1450                 double rho_v=_fluides[0]->getDensity(pression,T);
1451                 double rho_l=_fluides[1]->getDensity(pression,T);
1452                 //cout<<"Inlet alpha= "<<alpha<<" pression= "<<pression<<" temperature= "<<T<<" velocity gas "<<_limitField[nameOfGroup].v_x[0]<<" velocity liq "<<_limitField[nameOfGroup].v_x[1]<<endl;
1453
1454                 _externalStates[0]=alpha*rho_v;
1455                 _externalStates[1 + _Ndim] = (1-alpha)*rho_l;
1456                 v1[0] = _limitField[nameOfGroup].v_x[0];
1457                 v2[0] = _limitField[nameOfGroup].v_x[1];
1458                 if (_Ndim >1){
1459                         v1[1] = _limitField[nameOfGroup].v_y[0];
1460                         v2[1] = _limitField[nameOfGroup].v_y[1];
1461                 }
1462                 if (_Ndim >2){
1463                         v1[2] = _limitField[nameOfGroup].v_z[0];
1464                         v2[2] = _limitField[nameOfGroup].v_z[1];
1465                 }
1466                 for (int idim=0;idim<_Ndim;idim++){
1467                         _externalStates[1 + idim] = v1[idim]* _externalStates[0]; // phase 1
1468                         _externalStates[2 + _Ndim + idim] = v2[idim]* _externalStates[1+_Ndim]; // phase 2
1469                         v1_2 += v1[idim]*v1[idim];
1470                         v2_2 += v2[idim]*v2[idim];
1471                         _Vj[2+idim] = v1[idim];
1472                         _Vj[2+_Ndim+idim] = v2[idim];
1473                 }
1474                 _externalStates[_nVar-1] = _externalStates[0]      *(_fluides[0]->getInternalEnergy(T,rho_v) + v1_2/2)
1475                                                                                                                                                                                                                                                                                                                                         +_externalStates[1+_Ndim]*(_fluides[1]->getInternalEnergy(T,rho_l) + v2_2/2);
1476                 // _Vj external primitives
1477                 _Vj[0] = alpha;
1478                 _Vj[_nVar-1] = T;
1479
1480                 _idm[0] = 0;
1481                 for(k=1; k<_nVar; k++)
1482                         _idm[k] = _idm[k-1] + 1;
1483
1484                 VecAssemblyBegin(_Uext);
1485                 VecAssemblyBegin(_Vext);
1486                 VecAssemblyBegin(_Uextdiff);
1487                 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1488                 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1489                 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1490                 VecAssemblyEnd(_Uext);
1491                 VecAssemblyEnd(_Vext);
1492                 VecAssemblyEnd(_Uextdiff);
1493         }
1494         else if (_limitField[nameOfGroup].bcType==InletPressure){
1495                 _idm[0] = _nVar*j;
1496                 for(k=1; k<_nVar; k++)
1497                         _idm[k] = _idm[k-1] + 1;
1498
1499                 //Computation of the hydrostatic contribution : scalar product between gravity vector and position vector
1500                 Cell Cj=_mesh.getCell(j);
1501                 double hydroPress=Cj.x()*_GravityField3d[0];
1502                 if(_Ndim>1){
1503                         hydroPress+=Cj.y()*_GravityField3d[1];
1504                         if(_Ndim>2)
1505                                 hydroPress+=Cj.z()*_GravityField3d[2];
1506                 }
1507                 hydroPress*=_externalStates[0]+_externalStates[_Ndim];//multiplication by rho the total density
1508
1509                 //Building the external state
1510                 double alpha=_limitField[nameOfGroup].alpha;
1511                 double pression=_limitField[nameOfGroup].p+hydroPress;
1512                 double T=_limitField[nameOfGroup].T;
1513                 double rho_v=_fluides[0]->getDensity(pression,T);
1514                 double rho_l=_fluides[1]->getDensity(pression,T);
1515                 _externalStates[0]=alpha*rho_v;
1516                 _externalStates[1+_Ndim]=(1-alpha)*rho_l;
1517
1518                 for(int idim=0;idim<_Ndim;idim++){
1519                         _externalStates[idim+1]=_externalStates[0]*_Vj[idim+2];
1520                         _externalStates[idim+2+_Ndim]=_externalStates[1+_Ndim]*_Vj[idim+2+_Ndim];
1521                         v1_2+=_Vj[2+idim]*_Vj[2+idim];
1522                         v2_2+=_Vj[2+_Ndim+idim]*_Vj[2+_Ndim+idim];
1523                 }
1524                 _externalStates[_nVar-1]=    alpha *rho_v*(_fluides[0]->getInternalEnergy(T,rho_v)+v1_2/2)
1525                                                                                                                                                                                                                                                                                                                                                 +(1-alpha)*rho_l*(_fluides[1]->getInternalEnergy(T,rho_l)+v2_2/2);
1526                 // _Vj external primitives
1527                 _Vj[0] = alpha;
1528                 _Vj[1] = pression;
1529                 _Vj[_nVar-1] = T;
1530
1531                 _idm[0] = 0;
1532                 for(k=1; k<_nVar; k++)
1533                         _idm[k] = _idm[k-1] + 1;
1534                 VecAssemblyBegin(_Uext);
1535                 VecAssemblyBegin(_Vext);
1536                 VecAssemblyBegin(_Uextdiff);
1537                 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1538                 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1539                 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1540                 VecAssemblyEnd(_Uext);
1541                 VecAssemblyEnd(_Vext);
1542                 VecAssemblyEnd(_Uextdiff);
1543         }
1544         else if (_limitField[nameOfGroup].bcType==Outlet){
1545                 _idm[0] = _nVar*j;
1546                 for(k=1; k<_nVar; k++)
1547                         _idm[k] = _idm[k-1] + 1;
1548
1549                 //Computation of the hydrostatic contribution : scalar product between gravity vector and position vector
1550                 Cell Cj=_mesh.getCell(j);
1551                 double hydroPress=Cj.x()*_GravityField3d[0];
1552                 if(_Ndim>1){
1553                         hydroPress+=Cj.y()*_GravityField3d[1];
1554                         if(_Ndim>2)
1555                                 hydroPress+=Cj.z()*_GravityField3d[2];
1556                 }
1557                 hydroPress*=_externalStates[0]+_externalStates[_Ndim];//multiplication by rho the total density
1558
1559                 //Building the external state
1560                 double pression_int=_Vj[1];
1561                 double pression_ext=_limitField[nameOfGroup].p+hydroPress;
1562                 double T=_Vj[_nVar-1];
1563                 double rho_v_int=_fluides[0]->getDensity(pression_int,T);
1564                 double rho_l_int=_fluides[1]->getDensity(pression_int,T);
1565                 double rho_v_ext=_fluides[0]->getDensity(pression_ext,T);
1566                 double rho_l_ext=_fluides[1]->getDensity(pression_ext,T);
1567
1568                 for(k=0;k<1+_Ndim;k++){
1569                         _externalStates[k]*=rho_v_ext/rho_v_int;
1570                         _externalStates[k+1+_Ndim]*=rho_l_ext/rho_l_int;
1571                 }
1572                 double alpha=_Vj[0];
1573                 //cout<<"Outlet alpha= "<<alpha<<" pression int= "<<pression_int<<" pression ext= "<<pression_ext<<" temperature= "<<T<<" velocity gas "<<_Uj[2]<<" velocity liq "<<_Uj[2+_Ndim]<<endl;
1574                 for(int idim=0;idim<_Ndim;idim++){
1575                         v1_2+=_Vj[2+idim]*_Vj[2+idim];
1576                         v2_2+=_Vj[2+_Ndim+idim]*_Vj[2+_Ndim+idim];
1577                 }
1578                 _externalStates[_nVar-1]=alpha*rho_v_ext*(_fluides[0]->getInternalEnergy(T,rho_v_int)+v1_2/2)+(1-alpha)*rho_l_ext*(_fluides[1]->getInternalEnergy(T,rho_l_int)+v2_2/2);
1579
1580                 // _Vj external primitives
1581                 _Vj[1] = pression_ext;
1582
1583                 _idm[0] = 0;
1584                 for(k=1; k<_nVar; k++)
1585                         _idm[k] = _idm[k-1] + 1;
1586                 VecAssemblyBegin(_Uext);
1587                 VecAssemblyBegin(_Vext);
1588                 VecAssemblyBegin(_Uextdiff);
1589                 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1590                 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1591                 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1592                 VecAssemblyEnd(_Uext);
1593                 VecAssemblyEnd(_Vext);
1594                 VecAssemblyEnd(_Uextdiff);
1595         }
1596         else{
1597                 cout<<"Boundary condition not set for boundary named "<<nameOfGroup<<endl;
1598                 cout<<"Accepted boundary condition are Neumann, Wall, Inlet, and Outlet"<<endl;
1599                 throw CdmathException("Unknown boundary condition");
1600         }
1601 }
1602
1603 void FiveEqsTwoFluid::addDiffusionToSecondMember
1604 (               const int &i,
1605                 const int &j,
1606                 bool isBord)
1607 {
1608         double Tm=_Udiff[_nVar-1];
1609         double lambdal=_fluides[1]->getConductivity(Tm);
1610         double lambdav=_fluides[0]->getConductivity(Tm);
1611         double mu1 =_fluides[0]->getViscosity(Tm);
1612         double mu2 = _fluides[1]->getViscosity(Tm);
1613
1614         if(mu1==0 && mu2 ==0 && lambdav==0 && lambdal==0 && _heatTransfertCoeff==0)
1615                 return;
1616
1617         //extraction des valeurs
1618         _idm[0] = _nVar*i;
1619         for(int k=1; k<_nVar; k++)
1620                 _idm[k] = _idm[k-1] + 1;
1621
1622         VecGetValues(_primitiveVars, _nVar, _idm, _Vi);
1623         if (_verbose && _nbTimeStep%_freqSave ==0)
1624         {
1625                 cout << "Contribution diffusion: variables primitives maille " << i<<endl;
1626                 for(int q=0; q<_nVar; q++)
1627                 {
1628                         cout << _Vi[q] << endl;
1629                 }
1630                 cout << endl;
1631         }
1632
1633         if(!isBord ){
1634                 for(int k=0; k<_nVar; k++)
1635                         _idm[k] = _nVar*j + k;
1636
1637                 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
1638         }
1639         else
1640         {
1641                 lambdal=max(lambdal,_heatTransfertCoeff);//wall nucleate boing -> larger heat transfer
1642
1643                 for(int k=0; k<_nVar; k++)
1644                         _idm[k] = k;
1645                 VecGetValues(_Uextdiff, _nVar, _idm, _phi);
1646                 consToPrim(_phi,_Vj);
1647         }
1648
1649         if (_verbose && _nbTimeStep%_freqSave ==0)
1650         {
1651                 cout << "Contribution diffusion: variables primitives maille " <<j <<endl;
1652                 for(int q=0; q<_nVar; q++)
1653                 {
1654                         cout << _Vj[q] << endl;
1655                 }
1656                 cout << endl;
1657         }
1658         double alpha=(_Vj[0]+_Vi[0])/2;
1659         double lambda = alpha*lambdav+(1-alpha)*lambdal;
1660         //on n'a pas de contribution sur la masse
1661         _phi[0]=0;
1662         _phi[_Ndim+1]=0;
1663         //contribution visqueuse sur la quantite de mouvement
1664         for(int k=1; k<_Ndim+1; k++)
1665         {
1666                 _phi[k]         = _inv_dxi*2/(1/_inv_dxi+1/_inv_dxj)*mu1*alpha    *(_Vj[2+k] - _Vi[2+k]);//attention car primitif=alpha p u1 u2 T
1667                 _phi[k+_Ndim+1] = _inv_dxi*2/(1/_inv_dxi+1/_inv_dxj)*mu2*(1-alpha)*(_Vj[2+k+_Ndim] - _Vi[2+k+_Ndim]);
1668         }
1669         _phi[_nVar-1] = _inv_dxi*2/(1/_inv_dxi+1/_inv_dxj)*lambda  *(_Vj[_nVar-1] - _Vi[_nVar-1]);
1670         _idm[0] = i;
1671         VecSetValuesBlocked(_b, 1, _idm, _phi, ADD_VALUES);
1672
1673         if(_verbose && _nbTimeStep%_freqSave ==0)
1674         {
1675                 cout << "Contribution diffusion au 2nd membre pour la maille " << i << ": "<<endl;
1676                 for(int q=0; q<_nVar; q++)
1677                 {
1678                         cout << _phi[q] << endl;
1679                 }
1680                 cout << endl;
1681         }
1682
1683         if(!isBord)
1684         {
1685                 //On change de signe pour l'autre contribution
1686                 for(int k=0; k<_nVar; k++)
1687                         _phi[k] *= -_inv_dxj/_inv_dxi;
1688                 _idm[0] = j;
1689
1690                 VecSetValuesBlocked(_b, 1, _idm, _phi, ADD_VALUES);
1691         }
1692
1693         if(_verbose && _nbTimeStep%_freqSave ==0)
1694         {
1695                 cout << "Contribution diffusion au 2nd membre pour la maille  " << j << ": "<<endl;
1696                 for(int q=0; q<_nVar; q++)
1697                 {
1698                         cout << _phi[q] << endl;
1699                 }
1700                 cout << endl;
1701
1702                 if(_timeScheme==Implicit)
1703                 {
1704                         cout << "Matrice de diffusion D, pour le couple (" << i << "," << j<< "):" << endl;
1705                         for(int i=0; i<_nVar; i++)
1706                         {
1707                                 for(int j=0; j<_nVar; j++)
1708                                         cout << _Diffusion[i*_nVar+j]<<", ";
1709                                 cout << endl;
1710                         }
1711                         cout << endl;
1712                 }
1713         }
1714 }
1715
1716 void FiveEqsTwoFluid::jacobian(const int &j, string nameOfGroup, double * normale)
1717 {
1718         int k;
1719         for(k=0; k<_nVar*_nVar;k++)
1720                 _Jcb[k] = 0;//No implicitation at this stage
1721
1722         // loop of boundary types
1723         if (_limitField[nameOfGroup].bcType==Wall)
1724         {
1725                 for(k=0; k<_nVar;k++)
1726                         _Jcb[k*_nVar + k] = 1;
1727                 for(k=1; k<1+_Ndim;k++)
1728                         for(int l=1; l<1+_Ndim;l++){
1729                                 _Jcb[k*_nVar + l] -= 2*normale[k-1]*normale[l-1];
1730                                 _Jcb[(k+1+_Ndim)*_nVar + l+1+_Ndim] -= 2*normale[k-1]*normale[l-1];
1731                         }
1732         }
1733         //not wall
1734         else if (_limitField[nameOfGroup].bcType==Inlet)
1735         {
1736                 _Jcb[0] = 1;
1737                 _Jcb[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1738                 _Jcb[_nVar]=_limitField[nameOfGroup].v_x[0];
1739                 _Jcb[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1740                 if(_Ndim>1)
1741                 {
1742                         _Jcb[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1743                         _Jcb[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1744                         if(_Ndim==3)
1745                         {
1746                                 _Jcb[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1747                                 _Jcb[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1748
1749                         }
1750                 }
1751         }
1752         // not wall, not inlet
1753         else if (_limitField[nameOfGroup].bcType==Outlet){
1754                 _idm[0] = j*_nVar;
1755                 for(k=1; k<_nVar;k++)
1756                 {_idm[k] = _idm[k-1] + 1;}
1757                 VecGetValues(_conservativeVars, _nVar, _idm, _phi);
1758                 VecGetValues(_primitiveVars, _nVar, _idm, _externalStates);
1759         }
1760         else  if (_limitField[nameOfGroup].bcType!=Neumann && _limitField[nameOfGroup].bcType!=InletPressure)// not wall, not inlet, not outlet
1761         {
1762                 cout << "group named "<<nameOfGroup << " : unknown boundary condition" << endl;
1763                 throw CdmathException("FiveEqs::jacobianDiff: This boundary condition is not treated");
1764         }
1765 }
1766
1767
1768 void FiveEqsTwoFluid::primToCons(const double *P, const int &i, double *W, const int &j){
1769         //P= alpha, p, u1, u2, T
1770         //W=m1,q1,m2,q2,rhoE =alpha1*rho1*(e1+u1^2/2)+alpha2*rho2*(e2+u2^2/2)
1771         double alpha=P[i*_nVar];
1772         double pression=P[i*_nVar+1];
1773         double temperature=P[i*_nVar+_nVar-1];
1774         double rho_v=_fluides[0]->getDensity(pression,temperature);
1775         double rho_l=_fluides[1]->getDensity(pression,temperature);
1776         double u1_sq=0, u2_sq=0;
1777
1778         W[j*_nVar] = alpha*rho_v;
1779         W[j*_nVar+1+_Ndim] = (1-alpha)*rho_l;
1780         // alpha*rho*u
1781         for(int k=0; k<_Ndim; k++)
1782         {
1783                 W[j*_nVar+(k+1)] = W[j*_nVar]*P[i*_nVar+(k+2)];//alpha1*rho1*u1
1784                 W[j*_nVar+(k+1)+1+_Ndim] = W[j*_nVar+1+_Ndim]*P[i*_nVar+(k+2)+_Ndim];//alpha2*rho2*u2
1785         }
1786         // total energy
1787         W[j*_nVar+_nVar-1] = W[j*(_nVar)]* _fluides[0]->getInternalEnergy(temperature,rho_v)+W[j*(_nVar)+1+_Ndim]* _fluides[1]->getInternalEnergy(temperature,rho_l);
1788         for(int k=0; k<_Ndim; k++){
1789                 u1_sq+=P[i*_nVar+(k+2)]*P[i*_nVar+(k+2)];
1790                 u2_sq+=P[i*_nVar+(k+2)+_Ndim]*P[i*_nVar+(k+2)+_Ndim];
1791         }
1792         W[j*_nVar+_nVar-1] += (W[j*_nVar]*u1_sq+W[j*_nVar+1+_Ndim]*u2_sq)*0.5;
1793 }
1794
1795 void FiveEqsTwoFluid::consToPrim(const double *Wcons, double* Wprim,double porosity)//To do: treat porosity
1796 {
1797         //Wprim= alpha, p, u1, u2, T
1798         //Wcons=m1,q1,m2,q2,rhoE
1799         double m_v=Wcons[0];
1800         double m_l=Wcons[1+_Ndim];
1801         double q1_sq = 0,q2_sq = 0;
1802         _minm1=min(m_v,_minm1);
1803         _minm2=min(m_l,_minm2);
1804
1805         if(m_v<-_precision || m_l<-_precision){
1806                 _nbMaillesNeg+=1;
1807         }
1808         for(int k=0;k<_Ndim;k++){
1809                 q1_sq += Wcons[k+1]*Wcons[k+1];
1810                 q2_sq += Wcons[k+1+1+_Ndim]*Wcons[k+1+1+_Ndim];
1811         }
1812         if(Wcons[0]>0)//_precision*_precision*_precision)
1813                 q1_sq /= Wcons[0];      //alpha1 rho1 u1²
1814         else
1815                 q1_sq = 0;
1816         if(Wcons[1+_Ndim]>0)//_precision*_precision*_precision)
1817                 q2_sq /= Wcons[1+_Ndim];        //alpha2 rho2 u1²
1818         else
1819                 q2_sq = 0;
1820         double rho_m_e_m=Wcons[_nVar-1] -0.5*(q1_sq+q2_sq);
1821         //calcul de la temperature et de la pression pour une loi stiffened gas
1822         double temperature= (rho_m_e_m-m_v*static_cast<StiffenedGas*>(_fluides[0])->getInternalEnergy(0)-m_l*static_cast<StiffenedGas*>(_fluides[1])->getInternalEnergy(0))/(m_v*_fluides[0]->constante("cv")+m_l*_fluides[1]->constante("cv"));
1823         double e_v=static_cast<StiffenedGas*>(_fluides[0])->getInternalEnergy(temperature);
1824         double e_l=static_cast<StiffenedGas*>(_fluides[1])->getInternalEnergy(temperature);
1825         double gamma_v=_fluides[0]->constante("gamma");
1826         double gamma_l=_fluides[1]->constante("gamma");
1827         double Pinf_v=- gamma_v*_fluides[0]->constante("p0");
1828         double Pinf_l=- gamma_l*_fluides[1]->constante("p0");
1829         double a=1;
1830         double b=-(Pinf_v+m_v*(gamma_v-1)*e_v+Pinf_l+m_l*(gamma_l-1)*e_l);
1831         double c=Pinf_v*Pinf_l+Pinf_v*m_l*(gamma_l-1)*e_l+ Pinf_l*m_v*(gamma_v-1)*e_v;
1832         double delta= b*b-4*a*c;
1833         if(delta<0){
1834                 cout<<"delta= "<<delta<<" <0"<<endl;
1835                 *_runLogFile<<"delta= "<<delta<<" <0"<<endl;
1836                 throw CdmathException("FiveEqsTwoFluid::consToPrim: Failed to compute pressure");
1837         }
1838         double pression=(-b+sqrt(delta))/(2*a);
1839         if (pression < 1){
1840                 cout << "pressure = "<< pression << " < 1 Pa " << endl;
1841                 cout << "Conservative state = ";
1842                 for(int k=0; k<_nVar; k++){
1843                         cout<<Wcons[k]<<", ";
1844                 }
1845                 cout<<endl;
1846                 *_runLogFile << "FiveEqsTwoFluid::consToPrim: Failed to compute pressure = "<< pression << " < 1 Pa " << endl;
1847                 throw CdmathException("FiveEqsTwoFluid::consToPrim: Failed to compute pressure");
1848         }
1849
1850         double rho_v=_fluides[0]->getDensity(pression,temperature);
1851         double alpha=m_v/rho_v;
1852         Wprim[0]= alpha;
1853         Wprim[1] =  pression;
1854         for(int k=0;k<_Ndim;k++){//vitesses
1855                 if(Wcons[0]>0)//_precision*_precision*_precision)
1856                         Wprim[k+2] = Wcons[k+1]/Wcons[0];
1857                 else
1858                         Wprim[k+2] = Wcons[k+2+_Ndim]/Wcons[1+_Ndim];
1859                 if(Wcons[1+_Ndim]>0)//_precision*_precision*_precision)
1860                         Wprim[k+2+_Ndim] = Wcons[k+2+_Ndim]/Wcons[1+_Ndim];
1861                 else
1862                         Wprim[k+2+_Ndim] = Wcons[k+1]/Wcons[0];
1863         }
1864         Wprim[_nVar-1] = temperature;
1865 }
1866
1867 void FiveEqsTwoFluid::entropicShift(double* n, double& vpcorr0, double& vpcorr1)
1868 {
1869
1870         // parameters of function dgeev_ (compute the eigenvalues)
1871         int LDA, LDVL,LWORK, SDIM,LDVR;
1872         LDA = _nVar;
1873         LDVL = _nVar;
1874         LDVR = _nVar;
1875         LWORK = 20*_nVar;
1876         char jobvl[]="N", jobvr[]="N";
1877         double WORK[LWORK], JacoMat[_nVar*_nVar],egvaReal[_nVar],egvaImag[_nVar],egVectorL[_nVar*_nVar],egVectorR[_nVar*_nVar];
1878         int info_l = 0, info_r = 0;
1879
1880         /******** Left: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1881         convectionJacobianMatrix(_l, n);
1882         for (int i=0; i<_nVar*_nVar; i++){
1883                 JacoMat[i] = _JacoMat[i];
1884         }
1885         dgeev_(jobvl, jobvl,  &_nVar,
1886                         JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1887                         &LDVL,egVectorR,
1888                         &LDVR, WORK,&LWORK,
1889                         &info_l);
1890
1891         //      /******** Right: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1892         convectionJacobianMatrix(_r, n);
1893         for (int i=0; i<_nVar*_nVar; i++){
1894                 JacoMat[i] = _JacoMat[i];
1895         }
1896         dgeev_(jobvl, jobvl,  &_nVar,
1897                         JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1898                         &LDVL,egVectorR,
1899                         &LDVR, WORK,&LWORK,
1900                         &info_r);
1901
1902         if (info_l < 0 || info_r < 0)
1903         {
1904                 cout<<"Warning FiveEqsTwoFluid::entropicShift: dgeev_ did not compute all the eigenvalues, trying heuristic entropy correction "<<endl;
1905                 double u1l_n=0, u2l_n=0, u1r_n=0, u2r_n=0;
1906                 for (int idim=0; idim<_Ndim; idim++){
1907                         u1l_n= _l[2+idim]      *n[idim];
1908                         u2l_n= _l[2+idim+_Ndim]*n[idim];
1909                         u1r_n= _r[2+idim]      *n[idim];
1910                         u2r_n= _r[2+idim+_Ndim]*n[idim];
1911                 }
1912
1913                 vpcorr0 =max(fabs(u1l_n-u1r_n),fabs(u2l_n-u2r_n));
1914                 vpcorr1 = vpcorr0;
1915         }
1916         else
1917         {
1918                 std::vector< std::complex<double> > eigValuesLeft(_nVar);
1919                 std::vector< std::complex<double> > eigValuesRight(_nVar);
1920                 for(int j=0; j<_nVar; j++){
1921                         eigValuesLeft[j] = complex<double>(egvaReal[j],egvaImag[j]);
1922                         eigValuesRight[j] = complex<double>(egvaReal[j],egvaImag[j]);
1923                 }
1924                 Polynoms Poly;
1925                 int sizeLeft =  Poly.new_tri_selectif(eigValuesLeft, eigValuesLeft.size(), _precision);
1926                 int sizeRight =  Poly.new_tri_selectif(eigValuesRight, eigValuesRight.size(), _precision);
1927                 if (_verbose && _nbTimeStep%_freqSave ==0)
1928                 {
1929                         cout<<" Eigenvalue of JacoMat Left: " << endl;
1930                         for(int i=0; i<sizeLeft; i++)
1931                                 cout<<eigValuesLeft[i] << ", "<<endl;
1932                 }
1933                 if (_verbose && _nbTimeStep%_freqSave ==0)
1934                 {
1935                         cout<<" Eigenvalue of JacoMat Right: " << endl;
1936                         for(int i=0; i<sizeRight; i++)
1937                                 cout<<eigValuesRight[i] << ", "<<endl;
1938                 }
1939                 vpcorr0 = 0;
1940                 for (int i=1; i<min(sizeLeft,sizeRight)-1; i++)
1941                         vpcorr0 = max(vpcorr0, abs(eigValuesRight[i]-eigValuesLeft[i]));// Kieu
1942                 vpcorr1 = vpcorr0;
1943         }
1944 }
1945
1946 void FiveEqsTwoFluid::entropicShift(double* n)
1947 {
1948
1949         // parameters of function dgeev_ (compute the eigenvalues)
1950         int LDA, LDVL,LWORK, SDIM,LDVR;
1951         LDA = _nVar;
1952         LDVL = _nVar;
1953         LDVR = _nVar;
1954         LWORK = 20*_nVar;
1955         char jobvl[]="N", jobvr[]="N";
1956         double WORK[LWORK], JacoMat[_nVar*_nVar],egvaReal[_nVar],egvaImag[_nVar],egVectorL[_nVar*_nVar],egVectorR[_nVar*_nVar];
1957         int info = 0;
1958         /******** Left: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1959         convectionJacobianMatrix(_l, n);
1960         for (int i=0; i<_nVar*_nVar; i++){
1961                 JacoMat[i] = _JacoMat[i];
1962         }
1963         dgeev_(jobvl, jobvl,  &_nVar,
1964                         JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1965                         &LDVL,egVectorR,
1966                         &LDVR, WORK,&LWORK,
1967                         &info);
1968         if (info != 0){
1969                 *_runLogFile<<"FiveEqsTwoFluid::JacoMat: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (left)"<<endl;
1970                 throw CdmathException(
1971                                 "FiveEqsTwoFluid::JacoMat: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (left)");
1972         }
1973
1974         std::vector< std::complex<double> > eigValuesLeft(_nVar,0.);
1975         for(int j=0; j<_nVar; j++){
1976                 eigValuesLeft[j] = complex<double>(egvaReal[j],egvaImag[j]);
1977         }
1978         //      /******** Right: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1979         convectionJacobianMatrix(_r, n);
1980         for (int i=0; i<_nVar*_nVar; i++){
1981                 JacoMat[i] = _JacoMat[i];
1982         }
1983         dgeev_(jobvl, jobvl,  &_nVar,
1984                         JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1985                         &LDVL,egVectorR,
1986                         &LDVR, WORK,&LWORK,
1987                         &info);
1988         if (info != 0)
1989                 {
1990                 *_runLogFile<<"FiveEqsTwoFluid::entropicShift: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (right)"<<endl;
1991                 throw CdmathException(
1992                                 "FiveEqsTwoFluid::entropicShift: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (right)");
1993                 }
1994
1995         std::vector< std::complex<double> > eigValuesRight(_nVar,0.);
1996         for(int j=0; j<_nVar; j++){
1997                 eigValuesRight[j] = complex<double>(egvaReal[j],egvaImag[j]);
1998         }
1999         Polynoms Poly;
2000         int sizeLeft =  Poly.new_tri_selectif(eigValuesLeft, eigValuesLeft.size(), _precision);
2001         int sizeRight =  Poly.new_tri_selectif(eigValuesRight, eigValuesRight.size(), _precision);
2002         if (_verbose && _nbTimeStep%_freqSave ==0)
2003         {
2004                 cout<<" Eigenvalue of JacoMat Left: " << endl;
2005                 for(int i=0; i<sizeLeft; i++)
2006                         cout<<eigValuesLeft[i] << ", "<<endl;
2007                 cout<<" Eigenvalue of JacoMat Right: " << endl;
2008                 for(int i=0; i<sizeRight; i++)
2009                         cout<<eigValuesRight[i] << ", "<<endl;
2010         }
2011         for (int i=0; i<min(sizeLeft,sizeRight)-1; i++)
2012                 _entropicShift[i]= abs(eigValuesRight[i]-eigValuesLeft[i]);
2013 }
2014
2015 Vector FiveEqsTwoFluid::staggeredVFFCFlux()
2016 {
2017         if(_spaceScheme!=staggered || _nonLinearFormulation!=VFFC)
2018                 throw CdmathException("IsothermalTwoFluid::staggeredVFFCFlux: staggeredVFFCFlux method should be called only for VFFC formulation and staggered upwinding");
2019         else//_spaceScheme==staggered
2020         {
2021                 Vector Fij(_nVar);
2022                 Fij(_nVar-1)=0;
2023                 double alpha_roe = _Uroe[0];//Toumi formula
2024                 // interfacial pressure term (hyperbolic correction)
2025                 double dpi = _Uroe[_nVar];
2026
2027                 double u1ijn=0, u2ijn=0, phialphaq1n=0, phialphaq2n=0,vitesse1n=0,vitesse2n=0;
2028                 for(int idim=0;idim<_Ndim;idim++){//URoe = alpha, p, u1, u2, T, dpi
2029                         u1ijn+=_vec_normal[idim]*_Uroe[2+idim];
2030                         u2ijn+=_vec_normal[idim]*_Uroe[2+_Ndim+idim];
2031                 }
2032                 if(u1ijn>=0)
2033                 {
2034                         for(int idim=0;idim<_Ndim;idim++){
2035                                 phialphaq1n+=_vec_normal[idim]*_Ui[1+idim];//phi alpha rho u n
2036                                 vitesse1n  +=_vec_normal[idim]*_Vi[2+idim];
2037                         }
2038                         Fij(0)=phialphaq1n;
2039                         for(int idim=0;idim<_Ndim;idim++)
2040                                 Fij(1+idim)=phialphaq1n*_Vi[2+idim]+(alpha_roe*_Vj[1]*_porosityj+dpi*_Vi[0]*_porosityi)*_vec_normal[idim];
2041
2042                         double pressioni=_Vi[1];
2043                         double Temperaturei= _Vi[_nVar-1];
2044                         double rho1=_fluides[0]->getDensity(pressioni,Temperaturei);
2045                         double e1_int=_fluides[0]->getInternalEnergy(Temperaturei,rho1);
2046                         Fij(_nVar-1)+=_Ui[0]*(e1_int+0.5*vitesse1n*vitesse1n+_Vj[1]/rho1)*vitesse1n;
2047                 }
2048                 else
2049                 {
2050                         for(int idim=0;idim<_Ndim;idim++){
2051                                 phialphaq2n+=_vec_normal[idim]*_Uj[1+idim];//phi alpha rho u n
2052                                 vitesse1n  +=_vec_normal[idim]*_Vj[2+idim];
2053                         }
2054                         Fij(0)=phialphaq2n;
2055                         for(int idim=0;idim<_Ndim;idim++)
2056                                 Fij(1+idim)=phialphaq2n*_Vj[2+idim]+(alpha_roe*_Vi[1]*_porosityi+dpi*_Vj[0]*_porosityj)*_vec_normal[idim];
2057
2058                         double pressionj=_Vj[1];
2059                         double Temperaturej= _Vj[_nVar-1];
2060                         double rho1=_fluides[0]->getDensity(pressionj,Temperaturej);
2061                         double e1_int=_fluides[0]->getInternalEnergy(Temperaturej,rho1);
2062                         Fij(_nVar-1)+=_Uj[0]*(e1_int+0.5*vitesse1n*vitesse1n+_Vi[1]/rho1)*vitesse1n;
2063                 }
2064
2065                 if(u2ijn>=0)
2066                 {
2067                         for(int idim=0;idim<_Ndim;idim++){
2068                                 phialphaq2n+=_vec_normal[idim]*_Ui[2+_Ndim+idim];//phi alpha rho u n
2069                                 vitesse2n  +=_vec_normal[idim]*_Vi[2+idim+_Ndim];
2070                         }
2071                         Fij(1+_Ndim)=phialphaq2n;
2072                         for(int idim=0;idim<_Ndim;idim++)
2073                                 Fij(2+_Ndim+idim)=phialphaq2n*_Vi[2+_Ndim+idim]+((1-alpha_roe)*_Vj[1]*_porosityj-dpi*_Vi[0]*_porosityi)*_vec_normal[idim];
2074
2075                         double pressioni=_Vi[1];
2076                         double Temperaturei= _Vi[_nVar-1];
2077                         double rho2=_fluides[1]->getDensity(pressioni,Temperaturei);
2078                         double e2_int=_fluides[1]->getInternalEnergy(Temperaturei,rho2);
2079                         Fij(_nVar-1)+=_Ui[1+_Ndim]*(e2_int+0.5*vitesse2n*vitesse2n+_Vj[1]/rho2)*vitesse2n;
2080                 }
2081                 else
2082                 {
2083                         for(int idim=0;idim<_Ndim;idim++){
2084                                 phialphaq2n+=_vec_normal[idim]*_Uj[2+_Ndim+idim];//phi alpha rho u n
2085                                 vitesse2n  +=_vec_normal[idim]*_Vj[2+idim+_Ndim];
2086                         }
2087                         Fij(1+_Ndim)=phialphaq2n;
2088                         for(int idim=0;idim<_Ndim;idim++)
2089                                 Fij(2+_Ndim+idim)=phialphaq2n*_Vj[2+_Ndim+idim]+((1-alpha_roe)*_Vi[1]*_porosityi-dpi*_Vj[0]*_porosityj)*_vec_normal[idim];
2090
2091                         double pressionj=_Vj[1];
2092                         double Temperaturej= _Vj[_nVar-1];
2093                         double rho2=_fluides[1]->getDensity(pressionj,Temperaturej);
2094                         double e2_int=_fluides[1]->getInternalEnergy(Temperaturej,rho2);
2095                         Fij(_nVar-1)+=_Uj[1+_Ndim]*(e2_int+0.5*vitesse2n*vitesse2n+_Vi[1]/rho2)*vitesse2n;
2096                 }
2097                 return Fij;
2098         }
2099 }
2100
2101 void FiveEqsTwoFluid::applyVFRoeLowMachCorrections(bool isBord, string groupname)
2102 {
2103         if(_nonLinearFormulation!=VFRoe)
2104                 throw CdmathException("FiveEqsTwoFluid::applyVFRoeLowMachCorrections: applyVFRoeLowMachCorrections method should be called only for VFRoe formulation");
2105         else//_nonLinearFormulation==VFRoe
2106         {
2107                 if(_spaceScheme==lowMach){
2108                         double u1_2=0, u2_2=0;
2109                         for(int i=0;i<_Ndim;i++){
2110                                 u1_2 += _Uroe[2+i]*_Uroe[2+i];
2111                                 u2_2 += _Uroe[2+i+_Ndim]*_Uroe[2+i+_Ndim];
2112                         }
2113
2114                         double  c = _maxvploc;//mixture sound speed
2115                         double M=max(sqrt(u1_2),sqrt(u2_2))/c;//Mach number
2116                         _Vij[1]=M*_Vij[1]+(1-M)*(_Vi[1]+_Vj[1])/2;
2117                         primToCons(_Vij,0,_Uij,0);
2118                 }
2119                 else if(_spaceScheme==pressureCorrection)
2120                 {
2121                         if(_pressureCorrectionOrder>2)
2122                                 throw CdmathException("FiveEqsTwoFluid::applyVFRoeLowMachCorrections pressure correction order can be only 1 or 2 for five equation two-fluid model");
2123
2124                         double norm_uij=0, uij_n=0, ui_n=0, uj_n=0;//mean velocities
2125                         double rho1 =  _fluides[0]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2126                         double rho2 =  _fluides[1]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2127                         double m1=_Uroe[0]*rho1,  m2=(1-_Uroe[0])*rho2;
2128                         double rhom=m1+m2;
2129                         for(int i=0;i<_Ndim;i++)
2130                         {
2131                                 norm_uij += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*(m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim]);
2132                                 uij_n += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*_vec_normal[i];
2133                                 ui_n += _Vi[2+i]*_vec_normal[i];
2134                                 uj_n += _Vj[2+i]*_vec_normal[i];
2135                         }
2136                         norm_uij=sqrt(norm_uij)/rhom;
2137                         uij_n/=rhom;
2138                         if(norm_uij>_precision)//avoid division by zero
2139                                 _Vij[1]=(_Vi[1]+_Vj[1])/2 + uij_n/norm_uij*(_Vj[1]-_Vi[1])/4 - rhom*norm_uij*(uj_n-ui_n)/4;
2140                         else
2141                                 _Vij[1]=(_Vi[1]+_Vj[1])/2                                    - rhom*norm_uij*(uj_n-ui_n)/4;
2142                 }
2143                 else if(_spaceScheme==staggered)
2144                 {
2145                         double qij_n=0;
2146                         double rho1 =  _fluides[0]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2147                         double rho2 =  _fluides[1]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2148                         double m1=_Uroe[0]*rho1,  m2=(1-_Uroe[0])*rho2;
2149                         double rhom=m1+m2;
2150                         for(int i=0;i<_Ndim;i++)
2151                                 qij_n += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*_vec_normal[i];
2152
2153                         if(qij_n>=0){
2154                                 _Vij[0]=_Vi[0];
2155                                 _Vij[1]=_Vj[1];
2156                                 for(int i=0;i<_Ndim;i++)
2157                                 {
2158                                         _Vij[2+i]               =_Vi[2+i];
2159                                         _Vij[2+i+_Ndim] =_Vi[2+i+_Ndim];
2160                                 }
2161                                 _Vij[_nVar-1]=_Vi[_nVar-1];
2162                         }
2163                         else{
2164                                 _Vij[0]=_Vj[0];
2165                                 _Vij[1]=_Vi[1];
2166                                 for(int i=0;i<_Ndim;i++)
2167                                 {
2168                                         _Vij[2+i]               =_Vj[2+i];
2169                                         _Vij[2+i+_Ndim] =_Vj[2+i+_Ndim];
2170                                 }
2171                                 _Vij[_nVar-1]=_Vj[_nVar-1];
2172                         }
2173                         primToCons(_Vij,0,_Uij,0);
2174                 }
2175         }
2176 }
2177
2178 void FiveEqsTwoFluid::computeScaling(double maxvp)
2179 {
2180         _blockDiag[0]=1;//alphaScaling;
2181         _invBlockDiag[0]=1;//_blockDiag[0];
2182         _blockDiag[1+_Ndim]=1;//-alphaScaling;
2183         _invBlockDiag[1+_Ndim]=1.0;//_blockDiag[1+_Ndim];
2184         for(int q=1; q<_Ndim+1; q++)
2185         {
2186                 _blockDiag[q]=1/maxvp;
2187                 _invBlockDiag[q]=1/_blockDiag[q];
2188                 _blockDiag[q+1+_Ndim]=1/maxvp;
2189                 _invBlockDiag[q+1+_Ndim]=1/_blockDiag[q+1+_Ndim];
2190         }
2191         _blockDiag[_nVar - 1]=1/(maxvp*maxvp);//1
2192         _invBlockDiag[_nVar - 1]=  1./_blockDiag[_nVar - 1] ;// 1.;//
2193 }
2194
2195 void FiveEqsTwoFluid::testConservation()
2196 {
2197         double SUM, DELTA, x;
2198         int I;
2199         for(int i=0; i<_nVar; i++)
2200         {
2201                 {
2202                         if(i == 0)
2203                                 cout << "Masse totale phase " << 0 <<" (kg): ";
2204                         else if( i == 1+_Ndim)
2205                                 cout << "Masse totale phase " << 1 <<" (kg): ";
2206                         else if( i == _nVar-1)
2207                                 cout << "Energie totale " <<" (J.m^-3): ";
2208                         else
2209                         {
2210                                 if(i < 1+_Ndim)
2211                                         cout << "Quantite de mouvement totale phase 0 (kg.m.s^-1): ";
2212                                 else
2213                                         cout << "Quantite de mouvement totale phase 1 (kg.m.s^-1): ";
2214                         }
2215                 }
2216                 SUM = 0;
2217                 DELTA = 0;
2218                 I =  i;
2219                 for(int j=0; j<_Nmailles; j++)
2220                 {
2221                         VecGetValues(_old, 1, &I, &x);//on recupere la valeur du champ
2222                         SUM += x*_mesh.getCell(j).getMeasure();
2223                         VecGetValues(_newtonVariation, 1, &I, &x);//on recupere la variation du champ
2224                         DELTA += x*_mesh.getCell(j).getMeasure();
2225                         I += _nVar;
2226                 }
2227                 if(fabs(SUM)>_precision)
2228                         cout << SUM << ", variation relative: " << fabs(DELTA /SUM)  << endl;
2229                 else
2230                         cout << " a une somme nulle,  variation absolue: " << fabs(DELTA) << endl;
2231         }
2232 }
2233
2234 void FiveEqsTwoFluid::save(){
2235         string prim(_path+"/FiveEqsTwoFluidPrim_");
2236         string cons(_path+"/FiveEqsTwoFluidCons_");
2237         prim+=_fileName;
2238         cons+=_fileName;
2239
2240         PetscInt Ii, lj;
2241         for (PetscInt i = 0; i < _Nmailles; i++){
2242                 /* j = 0 : void fraction
2243                            j = 1 : pressure
2244                            j = 2, 3, 4: velocity phase 1
2245                            j = 5, 6, 7: velocity phase 2
2246                            j = 8 : temperature */
2247                 for (int j = 0; j < _nVar; j++){
2248                         Ii = i*_nVar +j;
2249                         VecGetValues(_primitiveVars,1,&Ii,&_VV(i,j));
2250                 }
2251         }
2252         if(_saveConservativeField){
2253                 for (long i = 0; i < _Nmailles; i++){
2254                         for (int j = 0; j < _nVar; j++){
2255                                 Ii = i*_nVar +j;
2256                                 VecGetValues(_conservativeVars,1,&Ii,&_UU(i,j));
2257                         }
2258                 }
2259                 _UU.setTime(_time,_nbTimeStep+1);
2260         }
2261         _VV.setTime(_time,_nbTimeStep+1);
2262
2263         if (_nbTimeStep ==0 || _restartWithNewFileName){
2264                 string prim_suppress ="rm -rf "+prim+"_*";
2265                 string cons_suppress ="rm -rf "+cons+"_*";
2266                 system(prim_suppress.c_str());//Nettoyage des précédents calculs identiques
2267                 system(cons_suppress.c_str());//Nettoyage des précédents calculs identiques
2268                 _VV.setInfoOnComponent(0,"Void_fraction");
2269                 _VV.setInfoOnComponent(1,"Pressure_(Pa)");
2270                 _VV.setInfoOnComponent(2,"Velocity1_x_m/s");
2271
2272                 if (_Ndim>1)
2273                         _VV.setInfoOnComponent(3,"Velocity1_y_m/s");
2274                 if (_Ndim>2)
2275                         _VV.setInfoOnComponent(4,"Velocity1_z_m/s");
2276                 _VV.setInfoOnComponent(2+_Ndim,"Velocity2_x_m/s");
2277                 if (_Ndim>1)
2278                         _VV.setInfoOnComponent(3+_Ndim,"Velocity2_y_m/s");
2279                 if (_Ndim>2)
2280                         _VV.setInfoOnComponent(4+_Ndim,"Velocity2_z_m/s");
2281                 _VV.setInfoOnComponent(_nVar-1,"Temperature_(K)");
2282
2283                 switch(_saveFormat)
2284                 {
2285                 case VTK :
2286                         _VV.writeVTK(prim);
2287                         break;
2288                 case MED :
2289                         _VV.writeMED(prim);
2290                         break;
2291                 case CSV :
2292                         _VV.writeCSV(prim);
2293                         break;
2294                 }
2295                 if(_saveConservativeField){
2296                         _UU.setInfoOnComponent(0,"Partial_density1");// (kg/m^3)
2297                         _UU.setInfoOnComponent(1,"Momentum1_x");// phase1 (kg/m^2/s)
2298                         if (_Ndim>1)
2299                                 _UU.setInfoOnComponent(2,"Momentum1_y");// phase1 (kg/m^2/s)
2300                         if (_Ndim>2)
2301                                 _UU.setInfoOnComponent(3,"Momentum1_z");// phase1  (kg/m^2/s)
2302                         _UU.setInfoOnComponent(1+_Ndim,"Partial_density2");// phase2 (kg/m^3)
2303                         _UU.setInfoOnComponent(2+_Ndim,"Momentum2_x");// phase2 (kg/m^2/s)
2304
2305                         if (_Ndim>1)
2306                                 _UU.setInfoOnComponent(3+_Ndim,"Momentum2_y");// phase2 (kg/m^2/s)
2307                         if (_Ndim>2)
2308                                 _UU.setInfoOnComponent(4+_Ndim,"Momentum2_z");// phase2 (kg/m^2/s)
2309                         _UU.setInfoOnComponent(_nVar-1,"Total_energy");
2310
2311                         switch(_saveFormat)
2312                         {
2313                         case VTK :
2314                                 _UU.writeVTK(cons);
2315                                 break;
2316                         case MED :
2317                                 _UU.writeMED(cons);
2318                                 break;
2319                         case CSV :
2320                                 _UU.writeCSV(cons);
2321                                 break;
2322                         }
2323                 }
2324         }
2325         else{
2326                 switch(_saveFormat)
2327                 {
2328                 case VTK :
2329                         _VV.writeVTK(prim,false);
2330                         break;
2331                 case MED :
2332                         _VV.writeMED(prim,false);
2333                         break;
2334                 case CSV :
2335                         _VV.writeCSV(prim);
2336                         break;
2337                 }
2338
2339                 if(_saveConservativeField){
2340                         switch(_saveFormat)
2341                         {
2342                         case VTK :
2343                                 _UU.writeVTK(cons,false);
2344                                 break;
2345                         case MED :
2346                                 _UU.writeMED(cons,false);
2347                                 break;
2348                         case CSV :
2349                                 _UU.writeCSV(cons);
2350                                 break;
2351                         }
2352                 }
2353         }
2354         if(_saveVelocity){
2355                 for (long i = 0; i < _Nmailles; i++){
2356                         // j = 0 : concentration, j=1 : pressure; j = _nVar - 1: temperature; j = 2,..,_nVar-2: velocity
2357                         for (int j = 0; j < _Ndim; j++)//On récupère les composantes de vitesse
2358                         {
2359                                 Ii = i*_nVar +2+j;
2360                                 VecGetValues(_primitiveVars,1,&Ii,&_Vitesse1(i,j));
2361                                 Ii=i*_nVar +2+j+_Ndim;
2362                                 VecGetValues(_primitiveVars,1,&Ii,&_Vitesse2(i,j));
2363                         }
2364                         for (int j = _Ndim; j < 3; j++){//On met à zero les composantes de vitesse si la dimension est <3
2365                                 _Vitesse1(i,j)=0;
2366                                 _Vitesse2(i,j)=0;
2367                         }
2368                 }
2369                 _Vitesse1.setTime(_time,_nbTimeStep);
2370                 _Vitesse2.setTime(_time,_nbTimeStep);
2371                 if (_nbTimeStep ==0 || _restartWithNewFileName){                
2372                         _Vitesse1.setInfoOnComponent(0,"Velocity_x_(m/s)");
2373                         _Vitesse1.setInfoOnComponent(1,"Velocity_y_(m/s)");
2374                         _Vitesse1.setInfoOnComponent(2,"Velocity_z_(m/s)");
2375
2376                         _Vitesse2.setInfoOnComponent(0,"Velocity_x_(m/s)");
2377                         _Vitesse2.setInfoOnComponent(1,"Velocity_y_(m/s)");
2378                         _Vitesse2.setInfoOnComponent(2,"Velocity_z_(m/s)");
2379
2380                         switch(_saveFormat)
2381                         {
2382                         case VTK :
2383                                 _Vitesse1.writeVTK(prim+"_GasVelocity");
2384                                 _Vitesse2.writeVTK(prim+"_LiquidVelocity");
2385                                 break;
2386                         case MED :
2387                                 _Vitesse1.writeMED(prim+"_GasVelocity");
2388                                 _Vitesse2.writeMED(prim+"_LiquidVelocity");
2389                                 break;
2390                         case CSV :
2391                                 _Vitesse1.writeCSV(prim+"_GasVelocity");
2392                                 _Vitesse2.writeCSV(prim+"_LiquidVelocity");
2393                                 break;
2394                         }
2395                 }
2396                 else{
2397                         switch(_saveFormat)
2398                         {
2399                         case VTK :
2400                                 _Vitesse1.writeVTK(prim+"_GasVelocity",false);
2401                                 _Vitesse2.writeVTK(prim+"_LiquidVelocity",false);
2402                                 break;
2403                         case MED :
2404                                 _Vitesse1.writeMED(prim+"_GasVelocity",false);
2405                                 _Vitesse2.writeMED(prim+"_LiquidVelocity",false);
2406                                 break;
2407                         case CSV :
2408                                 _Vitesse1.writeCSV(prim+"_GasVelocity");
2409                                 _Vitesse2.writeCSV(prim+"_LiquidVelocity");
2410                                 break;
2411                         }
2412                 }
2413         }
2414
2415         if (_restartWithNewFileName)
2416                 _restartWithNewFileName=false;
2417 }