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