4 * Created on: Jan 28, 2015
7 #include "FiveEqsTwoFluid.hxx"
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);
18 FiveEqsTwoFluid::FiveEqsTwoFluid(pressureEstimate pEstimate, int dim){
22 _dragCoeffs=vector<double>(2,0);
24 //Ne pas utiliser la loi de Stephane Dellacherie mais la stiffened gas standard
25 if (pEstimate==around1bar300K)//EOS at 1 bar and 373K
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
35 else//EOS at 155 bars and 618K
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
45 _latentHeat=_hsatv-_hsatl;
49 void FiveEqsTwoFluid::initialize()
51 cout<<"Initialising the five equation two fluid model"<<endl;
52 *_runLogFile<<"Initialising the five equation two fluid model"<<endl;
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");
57 _Uroe = new double[_nVar+1];
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
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++)
66 _gravite[i+1]=_GravityField3d[i];
67 _gravite[i+1 +_Ndim+1]=_GravityField3d[i];
69 _GravityImplicitationMatrix = new PetscScalar[_nVar*_nVar];
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
76 if(_entropicCorrection)
77 _entropicShift=vector<double>(_nVar);
79 ProblemFluid::initialize();
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)
87 for(int k=1; k<_nVar; k++)
88 _idm[k] = _idm[k-1] + 1;
89 VecGetValues(_conservativeVars, _nVar, _idm, _Ui);
92 for(int k=1; k<_nVar; k++)
93 _idm[k] = _idm[k-1] + 1;
95 VecGetValues(_Uext, _nVar, _idm, _Uj);
97 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
99 if(_verbose && _nbTimeStep%_freqSave ==0)
101 cout<<"Convection Left state cell " << i<< ": "<<endl;
102 for(int k =0; k<_nVar; k++)
104 cout<<"Convection Right state cell " << j<< ": "<<endl;
105 for(int k =0; k<_nVar; k++)
109 if(_Ui[0]<-(_precision) || _Uj[0]<-(_precision) || _Ui[_Ndim+1]<-(_precision) || _Uj[_Ndim+1]<-(_precision))
111 cout<<"Warning: masse partielle negative!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!"<<endl;
112 cout<< "valeurs a gauche: "<<_Ui[0]<<", "<<_Ui[_Ndim+1]<<", valeurs a droite: "<<_Uj[0]<<", "<<_Uj[_Ndim+1]<<endl;
113 // throw CdmathException(" Masse partielle negative, arret de calcul");
116 _Ui[0]=max(0.,_Ui[0]);
117 _Uj[0]=max(0.,_Uj[0]);
118 _Ui[_Ndim+1]=max(0.,_Ui[_Ndim+1]);
119 _Uj[_Ndim+1]=max(0.,_Uj[_Ndim+1]);
122 PetscScalar ri1, ri2, rj1, rj2, xi, xj;
123 //get _l and _r the primitive states left and right of the interface
125 for(int k=1; k<_nVar; k++)
126 _idm[k] = _idm[k-1] + 1;
127 VecGetValues(_primitiveVars, _nVar, _idm, _l);
131 //cout<<"_r is border"<<endl;
132 //consToPrim(_Uj, _r);
134 for(int k=1; k<_nVar; k++)
135 _idm[k] = _idm[k-1] + 1;
136 VecGetValues(_Vext, _nVar, _idm, _r);
141 for(int k=1; k<_nVar; k++)
142 _idm[k] = _idm[k-1] + 1;
143 VecGetValues(_primitiveVars, _nVar, _idm, _r);
145 if(_verbose && _nbTimeStep%_freqSave ==0)
148 for(int k=0;k<_nVar; k++)
151 for(int k=0;k<_nVar; k++)
154 // _Uroe[0] = \tilde{\alpha_v} = 1 - \tilde{\alpha_l} (formula of Toumi)
155 if(2-_l[0]-_r[0] > _precision)
156 _Uroe[0] = 1- 2*(1-_l[0])*(1-_r[0])/(2-_l[0]-_r[0]);
158 _Uroe[0] = (_l[0]+_r[0])/2;
160 if(_l[0]+_r[0] > _precision)
161 _Uroe[1] = (_l[1]*_l[0]+_r[1]*_r[0])/(_l[0]+_r[0]);
163 _Uroe[1] = (_l[1]*(1-_l[0])+_r[1]*(1-_r[0]))/(2-_l[0]-_r[0]);
165 ri1 = sqrt(_Ui[0]); ri2 = sqrt(_Ui[_Ndim+1]);
166 rj1 = sqrt(_Uj[0]); rj2 = sqrt(_Uj[_Ndim+1]);
167 for(int k=0;k<_Ndim;k++)
171 if(ri1>_precision && rj1>_precision)
172 _Uroe[2+k] = (xi/ri1 + xj/rj1)/(ri1 + rj1);
173 else if(ri1<_precision && rj1>_precision)
174 _Uroe[2+k] = xj/_Uj[0];
175 else if(ri1>_precision && rj1<_precision)
176 _Uroe[2+k] = xi/_Ui[0];
178 _Uroe[2+k] =(_Ui[k+1+_Ndim+1]/ri2 + _Uj[k+1+_Ndim+1]/rj2)/(ri2 + rj2);
180 xi = _Ui[k+1+_Ndim+1];
181 xj = _Uj[k+1+_Ndim+1];
182 if(ri2>_precision && rj2>_precision)
183 _Uroe[1+k+_Ndim+1] = (xi/ri2 + xj/rj2)/(ri2 + rj2);
184 else if(ri2<_precision && rj2>_precision)
185 _Uroe[1+k+_Ndim+1] = xj/_Uj[_Ndim+1];
186 else if(ri2>_precision && rj2<_precision)
187 _Uroe[1+k+_Ndim+1] = xi/_Ui[_Ndim+1];
189 _Uroe[1+k+_Ndim+1] = (xi/ri1 + xj/rj1)/(ri1 + rj1);
191 _Uroe[_nVar-1]=.5*(_l[_nVar-1]+_r[_nVar-1]);
193 //Fin du remplissage dans la fonction convectionMatrices
195 if(_verbose && _nbTimeStep%_freqSave ==0)
197 cout<<"Etat de Roe calcule: "<<endl;
198 for(int k=0;k<_nVar; k++)
199 cout<< _Uroe[k]<<endl;
203 void FiveEqsTwoFluid::diffusionStateAndMatrices(const long &i,const long &j, const bool &IsBord){
204 //sortie: matrices et etat Diffusion (alpha1 rho1, q1, alpha2 rho2, q2,T)
206 for(int k=1; k<_nVar; k++)
207 _idm[k] = _idm[k-1] + 1;
209 VecGetValues(_conservativeVars, _nVar, _idm, _Ui);
211 for(int k=1; k<_nVar; k++)
212 _idm[k] = _idm[k-1] + 1;
215 VecGetValues(_Uextdiff, _nVar, _idm, _Uj);
217 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
219 for(int k=0; k<_nVar; k++)
220 _Udiff[k] = (_Ui[k]+_Uj[k])/2;
222 for (int i = 0; i<_Ndim;i++){
223 q1_2+=_Udiff[ i+1]*_Udiff[ i+1];
224 q2_2+=_Udiff[1+_Ndim+i+1]*_Udiff[1+_Ndim+i+1];
226 consToPrim(_Udiff,_phi);
227 _Udiff[_nVar-1]=_phi[_nVar-1];
228 double alpha=_phi[0];
229 double Tm=_phi[_nVar-1];
230 double mu1 = _fluides[0]->getViscosity(Tm);
231 double mu2 = _fluides[1]->getViscosity(Tm);
232 double lambda = alpha*_fluides[0]->getConductivity(Tm)+(1-alpha)*_fluides[1]->getConductivity(Tm);
233 double Cv1= _fluides[0]->constante("Cv");
234 double Cv2= _fluides[1]->constante("Cv");
236 if(_timeScheme==Implicit)
238 for(int i=0; i<_nVar*_nVar;i++)
240 for(int idim=1;idim<_Ndim+1;idim++)
242 if(alpha>_precision){
243 _Diffusion[idim*_nVar] = alpha* mu1*_Udiff[idim]/(_Udiff[0]*_Udiff[0]);
244 _Diffusion[idim*_nVar+idim] = -alpha* mu1/_Udiff[0];
246 if(1-alpha>_precision){
247 _Diffusion[(idim+_Ndim+1)*_nVar] = (1-alpha)* mu2*_Udiff[idim+_Ndim+1]/(_Udiff[_Ndim+1]*_Udiff[_Ndim+1]);
248 _Diffusion[(idim+_Ndim+1)*_nVar+idim+_Ndim+1] = -(1-alpha)* mu2/_Udiff[_Ndim+1];
251 /*//Should correct the formula before using
252 int i = (_nVar-1)*_nVar;
253 _Diffusion[i]=lambda*(Tm/_Udiff[0]-q1_2/(2*Cv1*_Udiff[0]*_Udiff[0]*_Udiff[0]));
254 _Diffusion[i+1+_Ndim]=lambda*(Tm/_Udiff[1+_Ndim]-q2_2/(2*Cv2*_Udiff[1+_Ndim]*_Udiff[1+_Ndim]*_Udiff[1+_Ndim]));
255 for(int k=1;k<1+_Ndim;k++)
257 _Diffusion[i+k]= lambda*_Udiff[k]/(_Udiff[0]*_Udiff[0]*Cv1);
258 _Diffusion[i+k+1+_Ndim]= lambda*_Udiff[k+1+_Ndim]/(_Udiff[1+_Ndim]*_Udiff[+1+_Ndim]*Cv2);
260 _Diffusion[_nVar*_nVar-1]=-lambda/(_Udiff[0]*Cv1+_Udiff[1+_Ndim]*Cv2);
265 void FiveEqsTwoFluid::sourceVector(PetscScalar * Si,PetscScalar * Ui,PetscScalar * Vi, int i)
267 double m1=Ui[0],m2=Ui[1+_Ndim],rho=m1+m2, rhoE=Ui[_nVar-1], T=Vi[_nVar-1],alpha=Vi[0], P=Vi[1];
268 double norm_ur=0,norm_u1sq=0,norm_u2sq=0, Gamma;
270 for(int k=0; k<_Ndim; k++){
271 norm_ur+=(Vi[2+k]-Vi[2+k+_Ndim])*(Vi[2+k]-Vi[2+k+_Ndim]);
272 norm_u1sq+=Vi[2+k]*Vi[2+k];
273 norm_u2sq+=Vi[2+k+_Ndim]*Vi[2+k+_Ndim];
275 norm_ur=sqrt(norm_ur);
276 double h=(rhoE-0.5*m1*norm_u1sq-0.5*m2*norm_u2sq+P)/rho;
278 for(int k=0; k<_Ndim; k++)
279 u_int[k] = 0.5*(Vi[2+k]+Vi[2+k+_Ndim]);
280 // u_int[k] = Vi[0]*Vi[2+k+_Ndim] + (1-Vi[0])*Vi[2+k];
281 if(i>=0 && T>_Tsat && alpha<1-_precision)//if(i>=0 && _hsatv>h && h>_hsatl && alpha<1-_precision)
282 Gamma=_heatPowerField(i)/_latentHeat;
283 else//boundary cell, no phase change
286 for(int k=1; k<_Ndim+1; k++)
288 Si[k] =_gravite[k]*m1-_dragCoeffs[0]*norm_ur*(Vi[1+k]-Vi[1+k+_Ndim]) + Gamma*u_int[k-1];//interfacial velocity= ul
289 Si[k+_Ndim+1] =_gravite[k+_Ndim+1]*m2+ _dragCoeffs[0]*norm_ur*(Vi[1+k]-Vi[1+k+_Ndim]) - Gamma*u_int[k-1];
291 if(true){//heated boiling
295 else if (P<_Psat && alpha<1-_precision){//flash boiling
296 Si[0]=-_dHsatl_over_dp*_dp_over_dt(i)/_latentHeat;
297 Si[1+_Ndim]=_dHsatl_over_dp*_dp_over_dt(i)/_latentHeat;
304 Si[_nVar-1]=_heatPowerField(i);
305 else//boundary cell, no heating
307 for(int k=0; k<_Ndim; k++)
308 Si[_nVar-1] +=_GravityField3d[k]*(Ui[1+k]+Ui[2+k+_Ndim])-_dragCoeffs[0]*norm_ur*(Vi[2+k]-Vi[2+k+_Ndim])*(Vi[2+k]-Vi[2+k+_Ndim]);
310 if(_timeScheme==Implicit)
312 for(int i=0; i<_nVar*_nVar;i++)
313 _GravityImplicitationMatrix[i] = 0;
314 for(int i=0; i<_nVar/2;i++)
315 _GravityImplicitationMatrix[i*_nVar]=-_gravite[i];
316 for(int i=_nVar/2; i<_nVar;i++)
317 _GravityImplicitationMatrix[i*_nVar+_nVar/2]=-_gravite[i];
320 if(_verbose && _nbTimeStep%_freqSave ==0)
322 cout<<"FiveEqsTwoFluid::sourceVector"<<endl;
324 for(int k=0;k<_nVar;k++)
328 for(int k=0;k<_nVar;k++)
332 for(int k=0;k<_nVar;k++)
335 if(_timeScheme==Implicit)
336 displayMatrix(_GravityImplicitationMatrix, _nVar, "Gravity implicitation matrix");
340 void FiveEqsTwoFluid::pressureLossVector(PetscScalar * pressureLoss, double K, PetscScalar * Ui, PetscScalar * Vi, PetscScalar * Uj, PetscScalar * Vj)
342 double norm_u1=0, u1_n=0, norm_u2=0, u2_n=0, m1, m2;
343 for(int i=0;i<_Ndim;i++){
344 u1_n += _Uroe[1+i] *_vec_normal[i];
345 u2_n += _Uroe[1+i+_Ndim]*_vec_normal[i];
348 pressureLoss[1+_Ndim]=0;
350 for(int i=0;i<_Ndim;i++)
351 norm_u1 += Vi[1+i]*Vi[1+i];
352 norm_u1=sqrt(norm_u1);
354 for(int i=0;i<_Ndim;i++)
355 pressureLoss[1+i]=-K*m1*norm_u1*Vi[1+i];
358 for(int i=0;i<_Ndim;i++)
359 norm_u1 += Vj[1+i]*Vj[1+i];
360 norm_u1=sqrt(norm_u1);
362 for(int i=0;i<_Ndim;i++)
363 pressureLoss[1+i]=-K*m1*norm_u1*Vj[1+i];
366 for(int i=0;i<_Ndim;i++)
367 norm_u2 += Vi[2+i+_Ndim]*Vi[2+i+_Ndim];
368 norm_u2=sqrt(norm_u2);
370 for(int i=0;i<_Ndim;i++)
371 pressureLoss[2+i+_Ndim]=-K*m2*norm_u2*Vi[2+i+_Ndim];
374 for(int i=0;i<_Ndim;i++)
375 norm_u2 += Vj[2+i+_Ndim]*Vj[2+i+_Ndim];
376 norm_u2=sqrt(norm_u2);
378 for(int i=0;i<_Ndim;i++)
379 pressureLoss[2+i+_Ndim]=-K*m2*norm_u2*Vj[2+i+_Ndim];
381 pressureLoss[_nVar-1]=-K*(m1*norm_u1*norm_u1*norm_u1+m2*norm_u2*norm_u2*norm_u2);
383 if(_verbose && _nbTimeStep%_freqSave ==0)
385 cout<<"FiveEqsTwoFluid::pressureLossVector K= "<<K<<endl;
387 for(int k=0;k<_nVar;k++)
391 for(int k=0;k<_nVar;k++)
395 for(int k=0;k<_nVar;k++)
399 for(int k=0;k<_nVar;k++)
402 cout<<"pressureLoss="<<endl;
403 for(int k=0;k<_nVar;k++)
404 cout<<pressureLoss[k]<<", ";
409 void FiveEqsTwoFluid::porosityGradientSourceVector()
411 double u1_ni=0, u1_nj=0, u2_ni=0, u2_nj=0, rho1i, rho2i, rho1j, rho2j, pi=_Vi[1], pj=_Vj[1], Ti=_Vi[_nVar-1], Tj=_Vj[_nVar-1];
412 double pij1, pij2, alphaij=_Uroe[0];
413 for(int i=0;i<_Ndim;i++) {
414 u1_ni += _Vi[2+i]*_vec_normal[i];
415 u2_ni += _Vi[2+_Ndim+i]*_vec_normal[i];
416 u1_nj += _Vj[2+i]*_vec_normal[i];
417 u2_nj += _Vj[2+_Ndim+i]*_vec_normal[i];
419 _porosityGradientSourceVector[0]=0;
420 _porosityGradientSourceVector[1+_Ndim]=0;
421 rho1i = _fluides[0]->getDensity(pi, Ti);
422 rho2i = _fluides[1]->getDensity(pi, Ti);
423 rho1j = _fluides[0]->getDensity(pj, Tj);
424 rho2j = _fluides[1]->getDensity(pj, Tj);
425 pij1=(pi+pj)/2+rho1i*rho1j/2/(rho1i+rho1j)*(u1_ni-u1_nj)*(u1_ni-u1_nj);
426 pij2=(pi+pj)/2+rho2i*rho2j/2/(rho2i+rho2j)*(u2_ni-u2_nj)*(u2_ni-u2_nj);
427 for(int i=0;i<_Ndim;i++){
428 _porosityGradientSourceVector[1+i] =alphaij*pij1*(_porosityi-_porosityj)*2/(1/_inv_dxi+1/_inv_dxj);
429 _porosityGradientSourceVector[2+_Ndim+i]=alphaij*pij2*(_porosityi-_porosityj)*2/(1/_inv_dxi+1/_inv_dxj);
431 _porosityGradientSourceVector[_nVar-1]=0;
434 double FiveEqsTwoFluid::intPressDef(double alpha, double u_r2, double rho1, double rho2, double Temperature)
436 return _intPressCoeff*alpha*(1-alpha)*rho1*rho2*u_r2/( alpha*rho2+(1-alpha)*rho1);
437 +alpha*(1-alpha)*rho1*rho2*u_r2/((alpha*rho2+(1-alpha)*rho1)*(alpha*rho2+(1-alpha)*rho1)*(alpha*rho2+(1-alpha)*rho1)*(alpha*rho2+(1-alpha)*rho1))*u_r2
438 *(alpha*alpha*rho2-(1-alpha)*(1-alpha)*rho1)
439 *(alpha*alpha*rho2*rho2/(_fluides[0]->vitesseSonTemperature(Temperature,rho1)*_fluides[0]->vitesseSonTemperature(Temperature,rho1))
440 -(1-alpha)*(1-alpha)*rho1*rho1/(_fluides[1]->vitesseSonTemperature(Temperature,rho2)*_fluides[1]->vitesseSonTemperature(Temperature,rho2)));
443 Vector FiveEqsTwoFluid::convectionFlux(Vector U,Vector V, Vector normale, double porosity){
451 double phim1=U(0);//phi alpha1 rho1
452 double phim2=U(1+_Ndim);//phi alpha2 rho2
453 Vector phiq1(_Ndim),phiq2(_Ndim);//phi alpha1 rho1 u1, phi alpha2 rho2 u2
454 for(int i=0;i<_Ndim;i++){
456 phiq2(i)=U(2+_Ndim+i);
459 double pression=V(1);
460 Vector vitesse1(_Ndim),vitesse2(_Ndim);
461 for(int i=0;i<_Ndim;i++){
463 vitesse2(i)=V(2+_Ndim+i);
465 double Temperature= V(_nVar-1);
467 double vitesse1n=vitesse1*normale;
468 double vitesse2n=vitesse2*normale;
469 double rho1=_fluides[0]->getDensity(pression,Temperature);
470 double rho2=_fluides[1]->getDensity(pression,Temperature);
471 double e1_int=_fluides[0]->getInternalEnergy(Temperature,rho1);
472 double e2_int=_fluides[1]->getInternalEnergy(Temperature,rho2);
474 double alpha_roe = _Uroe[0];//Toumi formula
475 // interfacial pressure term (hyperbolic correction)
476 double dpi = _Uroe[_nVar];
479 F(0)=phim1*vitesse1n;
480 F(1+_Ndim)=phim2*vitesse2n;
481 for(int i=0;i<_Ndim;i++){
482 F(1+i)=phim1*vitesse1n*vitesse1(i)+(alpha_roe*pression+dpi*alpha)*porosity*normale(i);
483 F(2+_Ndim+i)=phim2*vitesse2n*vitesse2(i)+((1-alpha_roe)*pression+dpi*(1-alpha))*normale(i)*porosity;
485 F(_nVar-1)=phim1*(e1_int+0.5*vitesse1*vitesse1+pression/rho1)*vitesse1n+phim2*(e2_int+0.5*vitesse2*vitesse2+pression/rho2)*vitesse2n;
488 cout<<"Flux F(U,V)"<<endl;
495 void FiveEqsTwoFluid::convectionJacobianMatrix(double *V, double *n)
497 complex< double > tmp;
498 // enter : V(nVar) : primitive variables
502 double Tm = V[_nVar-1];
503 double rho1 = _fluides[0]->getDensity(p, Tm);
504 double rho2 = _fluides[1]->getDensity(p, Tm);
506 for (int idim=0; idim<_Ndim; idim++){
507 ur_2 += (V[2+idim]-V[2+idim+_Ndim])*(V[2+idim]-V[2+idim+_Ndim]);
509 // interfacial pressure term (hyperbolic correction)
510 double dpi1 = intPressDef(alp,ur_2, rho1, rho2,Tm);
513 /********Prepare the parameters to compute the Jacobian Matrix********/
514 /**** coefficients a, b, c ****/
515 double inv_a1_2,inv_a2_2,b1,c1,a2,b2,c2;
517 e1 = _fluides[0]->getInternalEnergy(V[_nVar-1],rho1);// primitive variable _l[_nVar-1]=Tm
518 e2 = _fluides[1]->getInternalEnergy(V[_nVar-1],rho2);
519 inv_a1_2 = static_cast<StiffenedGas*>(_fluides[0])->getDiffDensPress(e1);
520 inv_a2_2 = static_cast<StiffenedGas*>(_fluides[1])->getDiffDensPress(e2);
521 //double getJumpDensInternalEnergy(const double p_l,const double p_r,const double e_l,const double e_r);
522 b1 = static_cast<StiffenedGas*>(_fluides[0])->getDiffDensInternalEnergy(p,e1);
523 b2 = static_cast<StiffenedGas*>(_fluides[1])->getDiffDensInternalEnergy(p,e2);
524 //double getJumpInternalEnergyTemperature();
525 c1 = static_cast<StiffenedGas*>(_fluides[0])->getDiffInternalEnergyTemperature();
526 c2 = static_cast<StiffenedGas*>(_fluides[1])->getDiffInternalEnergyTemperature();
527 /**** coefficients eta, varrho_2 ****/
528 double eta[_Ndim], varrho_2;
529 // prefix m is arithmetic mean
531 double u1[_Ndim],u2[_Ndim], alp_u1[_Ndim],alp_u2[_Ndim];
534 varrho_2 =1/((alp*rho2)*inv_a1_2+((1-alp)*rho1)*inv_a2_2);
536 for (int idim=0; idim<_Ndim; idim++){
537 u1[idim] = V[idim+2];
538 u2[idim] = V[_Ndim+idim+2];
539 alp_u1[idim] = alp*u1[idim];
540 alp_u2[idim] = (1-alp)*u2[idim];
541 eta_n += (alp_u1[idim]*(1-p/rho1*inv_a1_2)+alp_u2[idim]*(1-p/rho2*inv_a2_2))*n[idim];
543 double eta_varrho_2n = eta_n*varrho_2;
544 /**** compute jump of Delta T, Delta e1, Delta e2 ****/
545 double inv_cm = 1/(c1*m1+c2*m2);
546 double DeltaT [_nVar], Delta_e1[_nVar], Delta_e2[_nVar];
549 DeltaT[1+_Ndim] =-e2;
550 DeltaT[_nVar-1] = 1 ;
551 for (int idim=0; idim<_Ndim; idim++){
553 DeltaT[1+_Ndim+idim+1] = 0;
555 for (int idim=0; idim<_Ndim; idim++){
557 DeltaT[0] += 0.5*u1[idim]*u1[idim];
559 DeltaT[_Ndim+1] += 0.5*u2[idim]*u2[idim];
561 DeltaT[idim+1] += - u1[idim];//*n[idim]
562 // wrt momentum liquid
563 DeltaT[_Ndim+idim+2] += - u2[idim];//*n[idim]
565 // finalize DeltaT, Delta_e1 and Delta_e2
566 for (int i =0; i< _nVar; ++i){
567 DeltaT[i] = inv_cm*DeltaT[i];
568 Delta_e1[i] = c1*DeltaT[i];
569 Delta_e2[i] = c2*DeltaT[i];
571 /**** compute differential flux (energy equation) A5 ****/
575 for (int i=0; i<_nVar; i++){
578 dF5[0] = eta_varrho_2n*rho2; // mass gas
579 dF5[_Ndim+1] = eta_varrho_2n*rho1; // mass liquid
580 for (int idim=0; idim<_Ndim; idim++){
582 dF5[idim+1]= (e1+p/rho1)*n[idim];
584 dF5[_Ndim+idim+2]=(e2+p/rho2)*n[idim];
586 // assign the value of A5 (last row of the Roe matrix)
587 for (int idim=0; idim<_Ndim; idim++){
588 for (int jdim=0; jdim<_Ndim;jdim++){
589 dF5[0] -= u1[idim]*u1[jdim]*u1[jdim]*n[idim];// -uin * ujn^2
590 dF5[_Ndim+1] -= u2[idim]*u2[jdim]*u2[jdim]*n[idim];
592 dF5[idim+1] += u1[idim]*u1[jdim]*n[jdim]+0.5*(u1[jdim]*u1[jdim])*n[idim];
594 dF5[_Ndim+idim+2] += u2[idim]*u2[jdim]*n[jdim]+0.5*(u2[jdim]*u2[jdim])*n[idim];
598 double coef_e1, coef_e2;
599 coef_e1 = - eta_varrho_2n*alp*rho2*b1;
600 coef_e2 = - eta_varrho_2n*(1-alp)*rho1*b2;
601 for (int idim=0; idim<_Ndim; idim++){
602 coef_e1 += (alp*rho1 - alp*p*b1/rho1)*u1[idim]*n[idim];
603 coef_e2 += ((1-alp)*rho2 - (1-alp)*p*b2/rho2)*u2[idim]*n[idim];
605 for (int i =0; i< _nVar; i++){
606 dF5[i] += coef_e1*Delta_e1[i] + coef_e2*Delta_e2[i];
608 /******** Construction de la matrice J *********/
610 for(int i=0; i<_nVar*_nVar;i++)
613 for(int idim=0; idim<_Ndim;idim++)
615 _JacoMat[1+idim]=n[idim];
616 _JacoMat[1+idim+_Ndim+1]=0.;
617 _JacoMat[(_Ndim+1)*_nVar+1+idim]=0.;
618 _JacoMat[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=n[idim];
620 // Roe Matrix new version
621 for(int idim=0; idim<_Ndim;idim++)
622 for (int jdim=0; jdim<_Ndim;jdim++){
623 // momentum gas (neglect delta alpha and delta P)
624 _JacoMat[ (1+idim)*_nVar] += -u1[idim]*u1[jdim]*n[jdim];
625 _JacoMat[(1+idim)*_nVar+jdim+1] += u1[idim]*n[jdim];
626 _JacoMat[(1+idim)*_nVar+idim+1] += u1[jdim]*n[jdim];
627 // momentum liquid (neglect delta alpha and delta P)
628 _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1] += -u2[idim]*u2[jdim]*n[jdim];
629 _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1+jdim+1] += u2[idim]*n[jdim];
630 _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1+idim+1] += u2[jdim]*n[jdim];
632 // update \Delta alpha
634 * (alpha *rho2*varrho_2+dpi1*(1-alpha)*inv_a2_2*varrho_2)*n[idim]
636 for (int idim=0; idim<_Ndim; idim++){
637 _JacoMat[ (1+idim)*_nVar] += dpi1*varrho_2*(1-alp)*inv_a2_2*n[idim];
638 _JacoMat[ (1+idim)*_nVar+_Ndim+1] += -dpi1*varrho_2*alp*inv_a1_2*n[idim];
639 _JacoMat[(2+_Ndim+idim)*_nVar] += - dpi2*varrho_2*(1-alp)*inv_a2_2*n[idim];
640 _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1] += dpi2*varrho_2*alp*inv_a1_2*n[idim];
641 for (int i=0; i<_nVar; i++){
642 _JacoMat[ (1+idim)*_nVar+i]+=dpi1*varrho_2*(-alp*(1-alp)*inv_a2_2*b1*Delta_e1[i]+alp*(1-alp)*inv_a1_2*b2*Delta_e2[i])*n[idim];
643 _JacoMat[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]+=-dpi2*varrho_2*(-alp*(1-alp)*inv_a2_2*b1*Delta_e1[i]+alp*(1-alp)*inv_a1_2*b2*Delta_e2[i])*n[idim];
647 for (int idim=0; idim<_Ndim; idim++){
648 _JacoMat[ (1+idim)*_nVar] += alp*varrho_2*rho2*n[idim];
649 _JacoMat[ (1+idim)*_nVar+_Ndim+1] +=alp* varrho_2*rho1*n[idim];
650 _JacoMat[(2+_Ndim+idim)*_nVar] += (1-alp)*varrho_2*rho2*n[idim];
651 _JacoMat[(2+_Ndim+idim)*_nVar+_Ndim+1] += (1-alp)* varrho_2*rho1*n[idim];
652 for (int i=0; i<_nVar; i++){
653 _JacoMat[ (1+idim)*_nVar+i]+=alp*varrho_2*(-alp*rho2*b1*Delta_e1[i] -(1-alp)*rho1*b2*Delta_e2[i])*n[idim];
654 _JacoMat[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]+=(1-alp)*varrho_2*(-alp*rho2*b1*Delta_e1[i] -(1-alp)*rho1*b2*Delta_e2[i])*n[idim];
657 // last row (total energy)
658 for (int i=0; i<_nVar; i++){
659 _JacoMat[(2*_Ndim+2)*_nVar +i] = dF5[i];
663 void FiveEqsTwoFluid::convectionMatrices()
665 //entree: URoe = alpha, p, u1, u2, T + ajout dpi
666 //sortie: matrices Roe+ et Roe- +Roe si schéma centre
668 if(_timeScheme==Implicit && _usePrimitiveVarsInNewton)
669 throw CdmathException("Implicitation with primitive variables not yet available for FiveEqsTwoFluid model");
672 complex< double > tmp;
673 double u_r2 = 0, u1_n=0, u2_n=0;
674 // u1_n = u1.n;u2_n = u2.n (scalar product)
675 for(int i=0;i<_Ndim;i++)
677 //u_r2 += (_Uroe[2+i]-_Uroe[1+i+1+_Ndim])*(_Uroe[2+i]-_Uroe[1+i+1+_Ndim]);
678 u_r2 += 0.5*((_l[2+i]-_l[2+i+_Ndim])*(_l[2+i]-_l[2+i+_Ndim])+(_r[2+i]-_r[2+i+_Ndim])*(_r[2+i]-_r[2+i+_Ndim])); //Kieu
679 u1_n += _Uroe[2+i]*_vec_normal[i];
680 u2_n += _Uroe[1+i+1+_Ndim]*_vec_normal[i];
683 //double alpha = (_l[0]+_r[0])*0.5;//Kieu formula
684 double alpha = _Uroe[0];//Toumi formula
685 //double p = _Uroe[1];//Toumi pressure
687 /***********Calcul des valeurs propres ********/
689 // ********Prepare the parameters to compute the Roe Matrix******** //
690 // **** coefficients eta, varrho_2 **** //
691 double eta[_Ndim], varrho_2;
692 double rho1_l = _fluides[0]->getDensity(_l[1], _l[_Ndim*2+2]);//(p,T)_Ndim*2+2
693 double rho2_l = _fluides[1]->getDensity(_l[1], _l[_Ndim*2+2]);
694 double rho1_r = _fluides[0]->getDensity(_r[1], _r[_Ndim*2+2]);
695 double rho2_r = _fluides[1]->getDensity(_r[1], _r[_Ndim*2+2]);
696 // **** coefficients a, b, c **** //
697 double inv_a1_2,inv_a2_2,b1,c1,a2,b2,c2;
698 double e1_l,e1_r,e2_l,e2_r;
699 e1_l = _fluides[0]->getInternalEnergy(_l[_nVar-1],rho1_l);// primitive variable _l[_nVar-1]=Tm
700 e2_l = _fluides[1]->getInternalEnergy(_l[_nVar-1],rho2_l);
701 e1_r = _fluides[0]->getInternalEnergy(_r[_nVar-1],rho1_r);
702 e2_r = _fluides[1]->getInternalEnergy(_r[_nVar-1],rho2_r);
703 inv_a1_2 = static_cast<StiffenedGas*>(_fluides[0])->getJumpDensPress(e1_l,e1_r);
704 inv_a2_2 = static_cast<StiffenedGas*>(_fluides[1])->getJumpDensPress(e2_l,e2_r);
705 //double getJumpDensInternalEnergy(const double p_l,const double p_r,const double e_l,const double e_r);
706 b1 = static_cast<StiffenedGas*>(_fluides[0])->getJumpDensInternalEnergy(_l[1],_r[1],e1_l,e1_r);
707 b2 = static_cast<StiffenedGas*>(_fluides[1])->getJumpDensInternalEnergy(_l[1],_r[1],e2_l,e2_r);
708 //double getJumpInternalEnergyTemperature();
709 c1 = static_cast<StiffenedGas*>(_fluides[0])->getJumpInternalEnergyTemperature();
710 c2 = static_cast<StiffenedGas*>(_fluides[1])->getJumpInternalEnergyTemperature();
712 // prefix m is arithmetic mean
713 double m_alp1,m_rho1,m_rho2,m_P,m_e1,m_e2,m_m1,m_m2, eta_n;
714 double m_u1[_Ndim],m_u2[_Ndim], m_alp_u1[_Ndim],m_alp_u2[_Ndim];
715 double u1_l[_Ndim], u2_l[_Ndim],u1_r[_Ndim],u2_r[_Ndim];
716 double u1l_2 = 0, u1r_2 = 0, u2l_2 = 0, u2r_2 = 0;
717 m_alp1 = 0.5*(_l[0]+_r[0]);
718 m_rho1 = 0.5*(rho1_l+rho1_r);
719 m_rho2 = 0.5*(rho2_l+rho2_r);
720 m_P = 0.5*(_l[1]+_r[1]);
721 m_e1 = 0.5*(e1_l+e1_r);
722 m_e2 = 0.5*(e2_l+e2_r);
723 m_m1 = 0.5*(_l[0]*rho1_l+_r[0]*rho1_r);
724 m_m2 = 0.5*((1-_l[0])*rho2_l+(1-_r[0])*rho2_r);
725 varrho_2 =1/((m_alp1*m_rho2)*inv_a1_2+((1-m_alp1)*m_rho1)*inv_a2_2);
728 for (int idim=0; idim<_Ndim; idim++){
729 u1_l[idim] = _l[idim+2];
730 u1_r[idim] = _r[idim+2];
731 u2_l[idim] = _l[_Ndim+idim+2];
732 u2_r[idim] = _r[_Ndim+idim+2];
733 m_u1[idim] = 0.5*(u1_l[idim] + u1_r[idim]);
734 m_u2[idim] = 0.5*(u2_l[idim] + u2_r[idim]);
735 m_alp_u1[idim] = 0.5*(_l[0]*u1_l[idim]+_r[0]*u1_r[idim]);
736 m_alp_u2[idim] = 0.5*((1-_l[0])*u2_l[idim]+(1-_r[0])*u2_r[idim]);
737 eta_n += (m_alp_u1[idim]*(1-m_P/m_rho1*inv_a1_2)+m_alp_u2[idim]*(1-m_P/m_rho2*inv_a2_2))*_vec_normal[idim];
740 double eta_varrho_2n = eta_n*varrho_2;
741 // **** compute jump of Delta T, Delta e1, Delta e2 **** //
742 for (int idim=0; idim<_Ndim; idim++){
743 u1l_2 += u1_l[idim]*u1_l[idim];
744 u1r_2 += u1_r[idim]*u1_r[idim];
745 u2l_2 += u2_l[idim]*u2_l[idim];
746 u2r_2 += u2_r[idim]*u2_r[idim];
748 double inv_m_cm = 1/(c1*m_m1+c2*m_m2);
749 double DeltaT [_nVar], Delta_e1[_nVar], Delta_e2[_nVar];
751 for (int i=0; i<_nVar; i++){
755 DeltaT[1+_Ndim] += -m_e2;
756 DeltaT[_nVar-1] += 1 ;
757 for (int idim=0; idim<_Ndim; idim++){
759 DeltaT[0] += 0.5*_l[idim+2] *_r[idim+2];//0.5*\tilde{u_g}^2
761 DeltaT[_Ndim+1] += 0.5*_l[_Ndim+idim+2] *_r[_Ndim+idim+2];
763 DeltaT[idim+1] += - m_u1[idim];
764 // wrt momentum liquid
765 DeltaT[_Ndim+idim+2] += - m_u2[idim];
768 // finalize DeltaT, Delta_e1 and Delta_e2
769 for (int i =0; i< _nVar; ++i){
770 DeltaT[i] = inv_m_cm*DeltaT[i];
771 Delta_e1[i] = c1*DeltaT[i];
772 Delta_e2[i] = c2*DeltaT[i];
775 // *** compute jump flux (energy equation) A5 *** //
779 for (int i=0; i<_nVar; i++){
782 A5[0] = eta_varrho_2n*m_rho2; // mass gas
783 A5[_Ndim+1] = eta_varrho_2n*m_rho1; // mass liquid
784 for (int idim=0; idim<_Ndim; idim++){
786 A5[idim+1] = (m_e1+m_P/m_rho1)*_vec_normal[idim];
788 A5[_Ndim+idim+2] = (m_e2+m_P/m_rho2)*_vec_normal[idim];
790 // assign the value of A5 (last row of the Roe matrix)
791 for (int idim=0; idim<_Ndim; idim++){
792 for (int jdim=0; jdim<_Ndim;jdim++){
793 A5[0] += 0.5*(0.5*(_l[idim+2]*_l[jdim+2]*_l[jdim+2]+_r[idim+2]*_r[jdim+2]*_r[jdim+2])*_vec_normal[idim]);// m_(uin*uj^2)
794 A5[0] -= m_u1[idim]*m_u1[jdim]*m_u1[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
795 A5[0] -= 0.5*(m_u1[idim]*_vec_normal[idim]*0.5*(_l[jdim+2]*_l[jdim+2]+_r[jdim+2]*_r[jdim+2]));//-0.5*m_uin*m_uj^2
796 A5[_Ndim+1]+= 0.5*(0.5*(_l[_Ndim+idim+2]*_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+idim+2]*_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2])*_vec_normal[idim]);// m_(uin*uj^2)
797 A5[_Ndim+1] -= m_u2[idim]*m_u2[jdim]*m_u2[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
798 A5[_Ndim+1] -= 0.5*(m_u2[idim]*0.5*(_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2]))*_vec_normal[idim];//-0.5*m_uin*m_uj^2
800 A5[idim+1] += m_u1[jdim]*_vec_normal[jdim]*m_u1[idim]+0.5*_vec_normal[idim]*0.5*(_l[jdim+2]*_l[jdim+2]+_r[jdim+2]*_r[jdim+2]);
802 A5[_Ndim+idim+2] += m_u2[jdim]*_vec_normal[jdim]*m_u2[idim]+0.5*_vec_normal[idim]*0.5*(_l[_Ndim+jdim+2]*_l[_Ndim+jdim+2]+_r[_Ndim+jdim+2]*_r[_Ndim+jdim+2]);
807 double coef_e1, coef_e2;
808 coef_e1 = - eta_varrho_2n*m_alp1*m_rho2*b1;
809 coef_e2 = - eta_varrho_2n*(1-m_alp1)*m_rho1*b2;
810 for (int idim=0; idim<_Ndim; idim++){
811 coef_e1 += (0.5*(_l[0]*rho1_l*_l[idim+2]+_r[0]*rho1_r*_r[idim+2]) - m_alp_u1[idim]*m_P*b1/m_rho1)*_vec_normal[idim];
812 coef_e2 += (0.5*((1-_l[0])*rho2_l*_l[_Ndim+idim+2]+(1-_r[0])*rho2_r*_r[_Ndim+idim+2])-m_alp_u2[idim]*m_P*b2/m_rho2)*_vec_normal[idim];
814 for (int i =0; i< _nVar; i++){
815 A5[i] += coef_e1*Delta_e1[i] + coef_e2*Delta_e2[i];
817 // ******* Construction de la matrice de Roe ******** //
818 // interfacial pressure correction
819 double T=_Uroe[_nVar-1];
820 double dpi1 = intPressDef(alpha, u_r2, m_rho1,m_rho2,T);
822 //saving dpi value for flux calculation later
826 for(int i=0; i<_nVar*_nVar;i++)
828 // alpha = 0.; dpi1 = 0.; dpi2 = 0.;
829 for(int idim=0; idim<_Ndim;idim++)
831 _Aroe[1+idim]=_vec_normal[idim];
832 _Aroe[1+idim+_Ndim+1]=0;
833 _Aroe[(_Ndim+1)*_nVar+1+idim]=0;
834 _Aroe[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=_vec_normal[idim];
837 // Take into account the convection term in the momentum eqts
838 for(int idim=0; idim<_Ndim;idim++)
839 for (int jdim=0; jdim<_Ndim;jdim++){
840 // momentum gas (neglect delta alpha and delta P)
841 _Aroe[ (1+idim)*_nVar] += (0.5*(_l[2+idim]*_l[2+jdim]+_r[2+idim]*_r[2+jdim])-2*m_u1[idim]*m_u1[jdim])*_vec_normal[jdim];
842 _Aroe[(1+idim)*_nVar+jdim+1] += m_u1[idim]*_vec_normal[jdim];
843 _Aroe[(1+idim)*_nVar+idim+1] += m_u1[jdim]*_vec_normal[jdim];
844 // momentum liquid (neglect delta alpha and delta P)
845 _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += (0.5*(_l[_Ndim+2+idim]*_l[_Ndim+2+jdim]+_r[_Ndim+2+idim]*_r[_Ndim+2+jdim])-2*m_u2[idim]*m_u2[jdim])*_vec_normal[jdim];
846 _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1+jdim+1] += m_u2[idim]*_vec_normal[jdim];
847 _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1+idim+1] += m_u2[jdim]*_vec_normal[jdim];
850 // update \Delta alpha
851 for (int idim=0; idim<_Ndim; idim++){
852 _Aroe[ (1+idim)*_nVar] += dpi1*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
853 _Aroe[ (1+idim)*_nVar+_Ndim+1] += -dpi1*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
854 _Aroe[(2+_Ndim+idim)*_nVar] += - dpi2*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
855 _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += dpi2*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
856 for (int i=0; i<_nVar; i++){
857 _Aroe[ (1+idim)*_nVar+i]+=dpi1*varrho_2*(-m_alp1*(1-m_alp1)*inv_a2_2*b1*Delta_e1[i]+m_alp1*(1-m_alp1)*inv_a1_2*b2*Delta_e2[i])*_vec_normal[idim];
858 _Aroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i]+=-dpi2*varrho_2*(-m_alp1*(1-m_alp1)*inv_a2_2*b1*Delta_e1[i]+m_alp1*(1-m_alp1)*inv_a1_2*b2*Delta_e2[i])*_vec_normal[idim];
862 for (int idim=0; idim<_Ndim; idim++){
863 _Aroe[ (1+idim)*_nVar] += alpha*varrho_2*m_rho2*_vec_normal[idim];
864 _Aroe[ (1+idim)*_nVar+_Ndim+1] += alpha* varrho_2*m_rho1*_vec_normal[idim];
865 _Aroe[(2+_Ndim+idim)*_nVar] += (1-alpha)*varrho_2*m_rho2*_vec_normal[idim];
866 _Aroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += (1-alpha)* varrho_2*m_rho1*_vec_normal[idim];
867 for (int i=0; i<_nVar; i++){
868 _Aroe[ (1+idim)*_nVar+i] += alpha*varrho_2*(-m_alp1*m_rho2*b1*Delta_e1[i] -(1-m_alp1)*m_rho1*b2*Delta_e2[i])*_vec_normal[idim];
869 _Aroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar+i] += (1-alpha)*varrho_2*(-m_alp1*m_rho2*b1*Delta_e1[i] -(1-m_alp1)*m_rho1*b2*Delta_e2[i])*_vec_normal[idim];
873 // last row (total energy)
874 for (int i=0; i<_nVar; i++){
875 _Aroe[(2*_Ndim+2)*_nVar +i] += A5[i];
881 int LWORK = 50*_nVar;
882 char jobvl[]="N", jobvr[]="N";
883 double WORK[LWORK], Aroe[_nVar*_nVar],egvaReal[_nVar],egvaImag[_nVar],
884 egVectorL[_nVar*_nVar],egVectorR[_nVar*_nVar];
887 std::vector< std::complex<double> > valeurs_propres_dist;
890 for (int i=0; i<_nVar*_nVar; i++)
893 if (_verbose && _nbTimeStep%_freqSave ==0)
895 cout<<endl<<"Matrice de Roe"<<endl;
896 for(int i=0; i<_nVar;i++)
898 for(int j=0; j<_nVar;j++)
899 cout << _Aroe[i*_nVar+j]<< " , ";
903 /******** Compute the eigenvalues and eigenvectors of Roe Matrix (using lapack)*********/
904 dgeev_(jobvl, jobvr, &_nVar,
905 Aroe,&LDA,egvaReal,egvaImag, egVectorL,
910 cout<<"FiveEqsTwoFluid::convectionMatrices: error dgeev_ : argument "<<-info<<" invalid"<<endl;
911 *_runLogFile<<"FiveEqsTwoFluid::convectionMatrices: error dgeev_ : argument "<<-info<<" invalid"<<endl;
912 throw CdmathException("FiveEqsTwoFluid::convectionMatrices: dgeev_ unsuccessful computation of the eigenvalues ");
916 cout<<"Warning FiveEqsTwoFluid::convectionMatrices: dgeev_ did not compute all the eigenvalues, trying Rusanov scheme "<<endl;
917 cout<<"Converged eigenvalues are ";
918 for(int i=info; i<_nVar; i++)
919 cout<<"("<< egvaReal[i]<<","<< egvaImag[i]<<"), ";
923 for(int i =info; i<_nVar; i++)
924 if (fabs(egvaReal[i])>_maxvploc)
925 _maxvploc=fabs(egvaReal[i]);
929 valeurs_propres_dist=std::vector< std::complex<double> > (1,maxvploc);
933 if (_verbose && _nbTimeStep%_freqSave ==0)
935 for(int i=0; i<_nVar; i++)
936 cout<<" Vp real part " << egvaReal[i]<<", Imaginary part " << egvaImag[i]<<endl;
939 std::vector< std::complex<double> > valeurs_propres(_nVar);
941 bool complexEigenvalues = false;
942 for(int j=0; j<_nVar; j++){
943 if (max(_l[0],_r[0])<_precision && abs(egvaImag[j])<_precision )// Kieu test Monophase
946 if (egvaImag[j] >_precision){// Kieu
947 complexEigenvalues = true;
949 if (abs(_l[0]-_r[0])<_precision*_precision && fabs(egvaImag[j])<_precision)// Kieu interfacial pressure
951 valeurs_propres[j] = complex<double>(egvaReal[j],egvaImag[j]);
954 taille_vp =Polynoms::new_tri_selectif(valeurs_propres,valeurs_propres.size(),_precision);
956 valeurs_propres_dist=vector< complex< double > >(taille_vp);
957 for( int i=0 ; i<taille_vp ; i++)
958 valeurs_propres_dist[i] = valeurs_propres[i];
959 if(_verbose && _nbTimeStep%_freqSave ==0)
961 cout<<" Vp apres tri " << valeurs_propres_dist.size()<<endl;
962 for(int ct =0; ct<taille_vp; ct++)
963 cout<< "("<<valeurs_propres_dist[ct].real()<< ", " <<valeurs_propres_dist[ct].imag() <<") ";
967 for(int i =0; i<taille_vp; i++){
968 if (fabs(valeurs_propres_dist[i].real())>maxvploc)
969 maxvploc=fabs(valeurs_propres_dist[i].real());
974 int existVpCplx = 0,pos_conj;
976 for (int ct=0; ct<taille_vp; ct++) {
977 vp_imag_iter = valeurs_propres_dist[ct].imag();
978 if ( fabs(vp_imag_iter) > 100*_precision ) {
980 if ( _part_imag_max < fabs(vp_imag_iter))
981 _part_imag_max = fabs(vp_imag_iter);
982 //On cherhe le conjugue
984 while(pos_conj<taille_vp && fabs(valeurs_propres_dist[pos_conj].imag()+vp_imag_iter)>_precision)
986 if(pos_conj!=ct+1 && pos_conj<taille_vp )
988 tmp=valeurs_propres_dist[ct+1];
989 valeurs_propres_dist[ct+1]=valeurs_propres_dist[pos_conj];
990 valeurs_propres_dist[pos_conj] = tmp;
998 /******* Construction des matrices de decentrement *******/
999 if(_spaceScheme == centered )
1001 if(_entropicCorrection)
1002 throw CdmathException("FiveEqsTwoFluid::convectionMatrices: entropic scheme not available for centered scheme");
1003 for(int i=0; i<_nVar*_nVar;i++)
1005 // if(alpha<_precision || alpha>1-_precision)//rusanov
1006 // for(int i=0; i<_nVar;i++)
1007 // _absAroe[i*_nVar+i]=maxvploc;
1009 if( _spaceScheme ==staggered){//To do: study entropic correction for staggered
1010 if(_entropicCorrection)//To do: study entropic correction for staggered
1011 throw CdmathException("FiveEqsTwoFluid::convectionMatrices: entropic scheme not yet available for staggered scheme");
1012 /******** Construction de la matrice de decentrement staggered *********/
1013 /***** Compute eta_n **************/
1015 for (int idim=0; idim<_Ndim; idim++){
1016 eta_n += (m_alp_u1[idim]*(-1-m_P/m_rho1*inv_a1_2)+m_alp_u2[idim]*(-1-m_P/m_rho2*inv_a2_2))*_vec_normal[idim];
1018 double eta_varrho_2n = eta_n*varrho_2;
1019 /**** compute jump flux (energy equation) A5 ****/
1023 for (int i=0; i<_nVar; i++){
1026 A5[0] = eta_varrho_2n*m_rho2; // mass gas
1027 A5[_Ndim+1] = eta_varrho_2n*m_rho1; // mass liquid
1028 // assign the value of A5 (last row of the Roe matrix)
1029 for (int idim=0; idim<_Ndim; idim++){
1030 for (int jdim=0; jdim<_Ndim;jdim++){
1031 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)
1032 A5[0] -= m_u1[idim]*m_u1[jdim]*m_u1[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
1033 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
1034 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)
1035 A5[_Ndim+1] -= m_u2[idim]*m_u2[jdim]*m_u2[jdim]*_vec_normal[idim];// -m_uin * m_ujn^2
1036 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
1040 double coef_e1, coef_e2;
1041 coef_e1 = - eta_varrho_2n*m_alp1*m_rho2*b1;
1042 coef_e2 = - eta_varrho_2n*(1-m_alp1)*m_rho1*b2;
1043 for (int idim=0; idim<_Ndim; idim++){
1044 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];
1045 coef_e2 += (0.5*((1-_l[0])*rho2_l*_l[_Ndim+idim+2]+(1-_r[0])*rho2_r*_r[_Ndim+idim+2])-m_alp_u2[idim]*m_P*b2/m_rho2)*_vec_normal[idim];
1047 for (int i =0; i< _nVar; i++){
1048 A5[i] += coef_e1*Delta_e1[i] + coef_e2*Delta_e2[i];
1050 /** Début remplissage matrice décentrement staggered **/
1052 for(int i=0; i<_nVar*_nVar;i++)
1056 for(int idim=0; idim<_Ndim;idim++)
1058 _absAroe[1+idim]=_vec_normal[idim];
1059 _absAroe[1+idim+_Ndim+1]=0;
1060 _absAroe[(_Ndim+1)*_nVar+1+idim]=0;
1061 _absAroe[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=_vec_normal[idim];
1063 //Contribution of convection (rho u\times u) in the momentum equations
1064 for(int idim=0; idim<_Ndim;idim++)
1065 for (int jdim=0; jdim<_Ndim;jdim++){
1066 // momentum gas (neglect delta alpha and delta P)
1067 _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];
1068 _absAroe[(1+idim)*_nVar+jdim+1] += m_u1[idim]*_vec_normal[jdim];
1069 _absAroe[(1+idim)*_nVar+idim+1] += m_u1[jdim]*_vec_normal[jdim];
1070 // momentum liquid (neglect delta alpha and delta P)
1071 _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];
1072 _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1+jdim+1] += m_u2[idim]*_vec_normal[jdim];
1073 _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1+idim+1] += m_u2[jdim]*_vec_normal[jdim];
1075 // contribution of interfacial pressure in momentum equation
1077 * (alpha *rho2*varrho_2+dpi1*(1-alpha)*inv_a2_2*varrho_2)*_vec_normal[idim]
1079 for (int idim=0; idim<_Ndim; idim++){
1080 _absAroe[ (1+idim)*_nVar] += dpi1*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
1081 _absAroe[ (1+idim)*_nVar+_Ndim+1] += -dpi1*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
1082 _absAroe[(2+_Ndim+idim)*_nVar] += - dpi2*varrho_2*(1-m_alp1)*inv_a2_2*_vec_normal[idim];
1083 _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1] += dpi2*varrho_2*(m_alp1)*inv_a1_2*_vec_normal[idim];
1084 for (int i=0; i<_nVar; i++){
1085 _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];
1086 _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];
1089 // contribution of pressure gradient in momentum equation
1090 for (int idim=0; idim<_Ndim; idim++){
1091 _absAroe[ (1+idim)*_nVar] -= alpha*varrho_2*m_rho2*_vec_normal[idim];
1092 _absAroe[ (1+idim)*_nVar+_Ndim+1] -=alpha* varrho_2*m_rho1*_vec_normal[idim];
1093 _absAroe[(2+_Ndim+idim)*_nVar] -= (1-alpha)*varrho_2*m_rho2*_vec_normal[idim];
1094 _absAroe[(2+_Ndim+idim)*_nVar+_Ndim+1] -= (1-alpha)* varrho_2*m_rho1*_vec_normal[idim];
1095 for (int i=0; i<_nVar; i++){
1096 _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];
1097 _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];
1100 // last row (total energy) To be changed soon
1101 for (int i=0; i<_nVar; i++){
1102 _Aroe[(2*_Ndim+2)*_nVar +i] = A5[i];
1104 double signu1,signu2;
1118 for(int i=0; i<(1+_Ndim)*_nVar;i++){
1119 _absAroe[i] *= signu1;
1120 _absAroe[i+(1+_Ndim)*_nVar] *= signu2;
1124 if(_spaceScheme ==upwind || _spaceScheme==pressureCorrection || _spaceScheme ==lowMach)//calcul de la valeur absolue
1125 { //on ordonne les deux premieres valeurs
1126 if(valeurs_propres_dist[1].real()<valeurs_propres_dist[0].real())
1128 tmp=valeurs_propres_dist[0];
1129 valeurs_propres_dist[0]=valeurs_propres_dist[1];
1130 valeurs_propres_dist[1]=tmp;
1132 vector< complex< double > > y (taille_vp,0);
1133 for( int i=0 ; i<taille_vp ; i++)
1134 y[i] = Polynoms::abs_generalise(valeurs_propres_dist[i]);
1136 if(_entropicCorrection)
1138 double entShift0 = 0;
1139 double entShift1 = 0;
1140 entropicShift(_vec_normal,entShift0,entShift1);
1141 //cout<<"entShift0= "<<entShift0<<endl;
1142 for( int i=0 ; i<taille_vp ; i++)
1144 //cout<<"y["<<i<<"]="<<y[i].real()<<endl;
1145 y[i] += max(entShift0,entShift1);
1149 if(_verbose && _nbTimeStep%_freqSave ==0)
1151 cout<<"valeurs propres"<<endl;
1152 for( int i=0 ; i<taille_vp ; i++)
1153 cout<<valeurs_propres_dist[i] <<", "<<endl;
1154 cout<<"valeurs à atteindre"<<endl;
1155 for( int i=0 ; i<taille_vp ; i++)
1156 cout<<y[i] <<", "<<endl;
1158 Polynoms::abs_par_interp_directe(taille_vp,valeurs_propres_dist, _Aroe, _nVar,_precision, _absAroe,y);
1160 if( _spaceScheme ==pressureCorrection){
1161 for( int i=0 ; i<_Ndim ; i++)
1162 for( int j=0 ; j<_Ndim ; j++){
1163 _absAroe[(1+i)*_nVar+1+j]-=alpha*(valeurs_propres_dist[1].real()-valeurs_propres_dist[0].real())/2*_vec_normal[i]*_vec_normal[j];
1164 _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];
1167 else if( _spaceScheme ==lowMach){
1168 double M=max(fabs(u1_n),fabs(u2_n))/maxvploc;
1169 for( int i=0 ; i<_Ndim ; i++)
1170 for( int j=0 ; j<_Ndim ; j++){
1171 _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];
1172 _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];
1177 //Calcul de la matrice signe pour VFFC, VFRoe et décentrement des termes source
1179 if(_entropicCorrection || _spaceScheme ==pressureCorrection || _spaceScheme ==lowMach){
1180 InvMatriceRoe( valeurs_propres_dist);
1181 Polynoms::matrixProduct(_absAroe, _nVar, _nVar, _invAroe, _nVar, _nVar, _signAroe);
1183 else if (_spaceScheme==upwind)//upwind sans entropic
1184 SigneMatriceRoe( valeurs_propres_dist);
1185 else if (_spaceScheme==centered)//centre sans entropic
1187 for(int i=0; i<_nVar*_nVar;i++)
1190 else if(_spaceScheme ==staggered )
1192 double signu1,signu2;
1205 for(int i=0; i<_nVar*_nVar;i++)
1207 _signAroe[0] = signu1;
1208 _signAroe[(1+_Ndim)*_nVar +1+_Ndim] = signu2;
1209 for(int i=2; i<_nVar-1;i++){
1210 _signAroe[i*_nVar+i] = -signu1;
1211 _signAroe[(i+1+_Ndim)*_nVar+i+1+_Ndim] = -signu2;
1213 //_signAroe[_nVar*(_nVar-1)+_nVar-1] = signu;
1216 throw CdmathException("FiveEqsTwoFluid::convectionMatrices: well balanced option not treated");
1218 for(int i=0; i<_nVar*_nVar;i++)
1220 _AroeMinus[i] = (_Aroe[i]-_absAroe[i])/2;
1221 _AroePlus[i] = (_Aroe[i]+_absAroe[i])/2;
1223 if(_timeScheme==Implicit && _usePrimitiveVarsInNewton)//Implicitation using primitive variables
1224 for(int i=0; i<_nVar*_nVar;i++)
1225 _AroeMinusImplicit[i] = (_AroeImplicit[i]-_absAroeImplicit[i])/2;
1227 for(int i=0; i<_nVar*_nVar;i++)
1228 _AroeMinusImplicit[i] = _AroeMinus[i];
1230 if(_verbose && _nbTimeStep%_freqSave ==0)
1232 cout<<"Matrice de Roe"<<endl;
1233 for(int i=0; i<_nVar;i++){
1234 for(int j=0; j<_nVar;j++)
1235 cout<<_Aroe[i*_nVar+j]<<" , ";
1238 cout<<"Valeur absolue matrice de Roe"<<endl;
1239 for(int i=0; i<_nVar;i++){
1240 for(int j=0; j<_nVar;j++)
1241 cout<<_absAroe[i*_nVar+j]<<" , ";
1244 cout<<"Signe matrice de Roe"<<endl;
1245 for(int i=0; i<_nVar;i++){
1246 for(int j=0; j<_nVar;j++)
1247 cout<<_signAroe[i*_nVar+j]<<" , ";
1250 cout<<endl<<"Matrice _AroeMinus"<<endl;
1251 for(int i=0; i<_nVar;i++)
1253 for(int j=0; j<_nVar;j++)
1254 cout << _AroeMinus[i*_nVar+j]<< " , ";
1257 cout<<endl<<"Matrice _AroePlus"<<endl;
1258 for(int i=0; i<_nVar;i++)
1260 for(int j=0; j<_nVar;j++)
1261 cout << _AroePlus[i*_nVar+j]<< " , ";
1267 void FiveEqsTwoFluid::jacobianDiff(const int &j, string nameOfGroup)
1270 for(k=0; k<_nVar*_nVar;k++)
1272 if (_limitField[nameOfGroup].bcType==Wall){
1274 for(k=1; k<_nVar; k++)
1275 _idm[k] = _idm[k-1] + 1;
1276 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
1277 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
1279 double pression=_Vj[1];//pressure inside
1280 double T=_Vj[_nVar-1];//temperature outside
1281 double rho_v=_fluides[0]->getDensity(pression,T);
1282 double rho_l=_fluides[1]->getDensity(pression,T);
1284 _JcbDiff[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1285 _JcbDiff[_nVar]=_limitField[nameOfGroup].v_x[0];
1286 _JcbDiff[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1287 double v2_v=_limitField[nameOfGroup].v_x[0]*_limitField[nameOfGroup].v_x[0];
1288 double v2_l=_limitField[nameOfGroup].v_x[1]*_limitField[nameOfGroup].v_x[1];
1291 _JcbDiff[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1292 _JcbDiff[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1293 v2_v+=_limitField[nameOfGroup].v_y[0]*_limitField[nameOfGroup].v_y[0];
1294 v2_l+=_limitField[nameOfGroup].v_y[1]*_limitField[nameOfGroup].v_y[1];
1297 _JcbDiff[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1298 _JcbDiff[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1299 v2_v+=_limitField[nameOfGroup].v_z[0]*_limitField[nameOfGroup].v_z[0];
1300 v2_l+=_limitField[nameOfGroup].v_z[1]*_limitField[nameOfGroup].v_z[1];
1303 _JcbDiff[(_nVar-1)*_nVar]= _fluides[0]->getInternalEnergy(_limitField[nameOfGroup].T,rho_v)+0.5*v2_v;
1304 _JcbDiff[(_nVar-1)*_nVar +1+_Ndim]=_fluides[1]->getInternalEnergy(_limitField[nameOfGroup].T,rho_l)+0.5*v2_l;
1306 else if (_limitField[nameOfGroup].bcType==Inlet){
1309 _JcbDiff[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1310 _JcbDiff[_nVar]=_limitField[nameOfGroup].v_x[0];
1311 _JcbDiff[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1312 double v2_v=_limitField[nameOfGroup].v_x[0]*_limitField[nameOfGroup].v_x[0];
1313 double v2_l=_limitField[nameOfGroup].v_x[1]*_limitField[nameOfGroup].v_x[1];
1316 _JcbDiff[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1317 _JcbDiff[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1318 v2_v+=_limitField[nameOfGroup].v_y[0]*_limitField[nameOfGroup].v_y[0];
1319 v2_l+=_limitField[nameOfGroup].v_y[1]*_limitField[nameOfGroup].v_y[1];
1322 _JcbDiff[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1323 _JcbDiff[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1324 v2_v+=_limitField[nameOfGroup].v_z[0]*_limitField[nameOfGroup].v_z[0];
1325 v2_l+=_limitField[nameOfGroup].v_z[1]*_limitField[nameOfGroup].v_z[1];
1328 _JcbDiff[(_nVar-1)*_nVar]= _fluides[0]->getInternalEnergy(_limitField[nameOfGroup].T,rho_v)+0.5*v2_v;
1329 _JcbDiff[(_nVar-1)*_nVar +1+_Ndim]=_fluides[1]->getInternalEnergy(_limitField[nameOfGroup].T,rho_l)+0.5*v2_l;
1331 } else if (_limitField[nameOfGroup].bcType==Outlet){
1332 //extraction de l etat courant et primitives
1335 for(k=1; k<_nVar;k++)
1336 {_idm[k] = _idm[k-1] + 1;}
1337 VecGetValues(_conservativeVars, _nVar, _idm, _phi);
1338 VecGetValues(_primitiveVars, _nVar, _idm, _externalStates);
1341 else if (_limitField[nameOfGroup].bcType!=Neumann && _limitField[nameOfGroup].bcType!=InletPressure){
1342 cout<<"Condition limite non traitee pour le bord "<<nameOfGroup<< endl;
1343 throw CdmathException("FiveEqsTwoFluid::jacobianDiff: Condition limite non traitee");
1348 void FiveEqsTwoFluid::setBoundaryState(string nameOfGroup, const int &j,double *normale){
1349 //To do controle signe des vitesses pour CL entree/sortie
1351 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;
1352 double v1[_Ndim], v2[_Ndim];
1355 for(k=1; k<_nVar; k++)
1356 _idm[k] = _idm[k-1] + 1;
1358 VecGetValues(_conservativeVars, _nVar, _idm, _externalStates);//On initialise l'état fantôme avec l'état interne
1359 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
1361 for(k=0; k<_Ndim; k++){
1362 q1_n+=_externalStates[(k+1)]*normale[k];
1363 q2_n+=_externalStates[(k+1+1+_Ndim)]*normale[k];
1364 u1_n+=_Vj[(k+2)]*normale[k];
1365 u2_n+=_Vj[(k+2+_Ndim)]*normale[k];
1368 if(_verbose && _nbTimeStep%_freqSave ==0)
1370 cout << "Boundary conditions for group "<< nameOfGroup<< ", inner cell j= "<<j << " face unit normal vector "<<endl;
1371 for(k=0; k<_Ndim; k++){
1372 cout<<normale[k]<<", ";
1377 if (_limitField[nameOfGroup].bcType==Wall){
1379 //Pour la convection, inversion du sens des vitesses
1380 for(k=0; k<_Ndim; k++){
1381 _externalStates[(k+1)]-= 2*q1_n*normale[k];
1382 _externalStates[(k+1+1+_Ndim)]-= 2*q2_n*normale[k];
1383 _Vj[(k+2)]-= 2*u1_n*normale[k];
1384 _Vj[(k+2+_Ndim)]-= 2*u2_n*normale[k];
1387 for(k=1; k<_nVar; k++)
1388 _idm[k] = _idm[k-1] + 1;
1390 VecAssemblyBegin(_Uext);
1391 VecAssemblyBegin(_Vext);
1392 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1393 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1394 VecAssemblyEnd(_Uext);
1395 VecAssemblyEnd(_Vext);
1397 //Pour la diffusion, paroi à vitesses et temperature imposees
1398 double pression=_Vj[1];//pressure inside
1399 double T=_Vj[_nVar-1];//temperature outside
1400 double rho_v=_fluides[0]->getDensity(pression,T);
1401 double rho_l=_fluides[1]->getDensity(pression,T);
1402 _externalStates[1]=_externalStates[0]*_limitField[nameOfGroup].v_x[0];
1403 _externalStates[2+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_x[1];
1406 _externalStates[2]=_externalStates[0]*_limitField[nameOfGroup].v_y[0];
1407 _externalStates[3+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_y[1];
1410 _externalStates[3]=_externalStates[0]*_limitField[nameOfGroup].v_z[0];
1411 _externalStates[4+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_z[1];
1414 _externalStates[_nVar-1] = _externalStates[0]*(_fluides[0]->getInternalEnergy(_limitField[nameOfGroup].T,rho_v) + v1_2/2)
1415 +_externalStates[1+_Ndim]*(_fluides[1]->getInternalEnergy(_limitField[nameOfGroup].T,rho_l) + v2_2/2);
1416 VecAssemblyBegin(_Uextdiff);
1417 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1418 VecAssemblyEnd(_Uextdiff);
1420 else if (_limitField[nameOfGroup].bcType==Neumann){
1422 for(k=1; k<_nVar; k++)
1423 _idm[k] = _idm[k-1] + 1;
1426 for(k=1; k<_nVar; k++)
1427 _idm[k] = _idm[k-1] + 1;
1429 VecAssemblyBegin(_Uext);
1430 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1431 VecAssemblyEnd(_Uext);
1433 VecAssemblyBegin(_Vext);
1434 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1435 VecAssemblyEnd(_Vext);
1437 VecAssemblyBegin(_Uextdiff);
1438 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1439 VecAssemblyEnd(_Uextdiff);
1441 else if (_limitField[nameOfGroup].bcType==Inlet){
1443 for(k=1; k<_nVar; k++)
1444 _idm[k] = _idm[k-1] + 1;
1446 double alpha=_limitField[nameOfGroup].alpha;//void fraction outside
1447 double pression=_Vj[1];//pressure inside
1448 double T=_limitField[nameOfGroup].T;//temperature outside
1449 double rho_v=_fluides[0]->getDensity(pression,T);
1450 double rho_l=_fluides[1]->getDensity(pression,T);
1451 //cout<<"Inlet alpha= "<<alpha<<" pression= "<<pression<<" temperature= "<<T<<" velocity gas "<<_limitField[nameOfGroup].v_x[0]<<" velocity liq "<<_limitField[nameOfGroup].v_x[1]<<endl;
1453 _externalStates[0]=alpha*rho_v;
1454 _externalStates[1 + _Ndim] = (1-alpha)*rho_l;
1455 v1[0] = _limitField[nameOfGroup].v_x[0];
1456 v2[0] = _limitField[nameOfGroup].v_x[1];
1458 v1[1] = _limitField[nameOfGroup].v_y[0];
1459 v2[1] = _limitField[nameOfGroup].v_y[1];
1462 v1[2] = _limitField[nameOfGroup].v_z[0];
1463 v2[2] = _limitField[nameOfGroup].v_z[1];
1465 for (int idim=0;idim<_Ndim;idim++){
1466 _externalStates[1 + idim] = v1[idim]* _externalStates[0]; // phase 1
1467 _externalStates[2 + _Ndim + idim] = v2[idim]* _externalStates[1+_Ndim]; // phase 2
1468 v1_2 += v1[idim]*v1[idim];
1469 v2_2 += v2[idim]*v2[idim];
1470 _Vj[2+idim] = v1[idim];
1471 _Vj[2+_Ndim+idim] = v2[idim];
1473 _externalStates[_nVar-1] = _externalStates[0] *(_fluides[0]->getInternalEnergy(T,rho_v) + v1_2/2)
1474 +_externalStates[1+_Ndim]*(_fluides[1]->getInternalEnergy(T,rho_l) + v2_2/2);
1475 // _Vj external primitives
1480 for(k=1; k<_nVar; k++)
1481 _idm[k] = _idm[k-1] + 1;
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);
1493 else if (_limitField[nameOfGroup].bcType==InletPressure){
1495 for(k=1; k<_nVar; k++)
1496 _idm[k] = _idm[k-1] + 1;
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];
1502 hydroPress+=Cj.y()*_GravityField3d[1];
1504 hydroPress+=Cj.z()*_GravityField3d[2];
1506 hydroPress*=_externalStates[0]+_externalStates[_Ndim];//multiplication by rho the total density
1508 //Building the external state
1509 double alpha=_limitField[nameOfGroup].alpha;
1510 double pression=_limitField[nameOfGroup].p+hydroPress;
1511 double T=_limitField[nameOfGroup].T;
1512 double rho_v=_fluides[0]->getDensity(pression,T);
1513 double rho_l=_fluides[1]->getDensity(pression,T);
1514 _externalStates[0]=alpha*rho_v;
1515 _externalStates[1+_Ndim]=(1-alpha)*rho_l;
1517 for(int idim=0;idim<_Ndim;idim++){
1518 _externalStates[idim+1]=_externalStates[0]*_Vj[idim+2];
1519 _externalStates[idim+2+_Ndim]=_externalStates[1+_Ndim]*_Vj[idim+2+_Ndim];
1520 v1_2+=_Vj[2+idim]*_Vj[2+idim];
1521 v2_2+=_Vj[2+_Ndim+idim]*_Vj[2+_Ndim+idim];
1523 _externalStates[_nVar-1]= alpha *rho_v*(_fluides[0]->getInternalEnergy(T,rho_v)+v1_2/2)
1524 +(1-alpha)*rho_l*(_fluides[1]->getInternalEnergy(T,rho_l)+v2_2/2);
1525 // _Vj external primitives
1531 for(k=1; k<_nVar; k++)
1532 _idm[k] = _idm[k-1] + 1;
1533 VecAssemblyBegin(_Uext);
1534 VecAssemblyBegin(_Vext);
1535 VecAssemblyBegin(_Uextdiff);
1536 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1537 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1538 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1539 VecAssemblyEnd(_Uext);
1540 VecAssemblyEnd(_Vext);
1541 VecAssemblyEnd(_Uextdiff);
1543 else if (_limitField[nameOfGroup].bcType==Outlet){
1545 for(k=1; k<_nVar; k++)
1546 _idm[k] = _idm[k-1] + 1;
1548 //Computation of the hydrostatic contribution : scalar product between gravity vector and position vector
1549 Cell Cj=_mesh.getCell(j);
1550 double hydroPress=Cj.x()*_GravityField3d[0];
1552 hydroPress+=Cj.y()*_GravityField3d[1];
1554 hydroPress+=Cj.z()*_GravityField3d[2];
1556 hydroPress*=_externalStates[0]+_externalStates[_Ndim];//multiplication by rho the total density
1558 //Building the external state
1559 double pression_int=_Vj[1];
1560 double pression_ext=_limitField[nameOfGroup].p+hydroPress;
1561 double T=_Vj[_nVar-1];
1562 double rho_v_int=_fluides[0]->getDensity(pression_int,T);
1563 double rho_l_int=_fluides[1]->getDensity(pression_int,T);
1564 double rho_v_ext=_fluides[0]->getDensity(pression_ext,T);
1565 double rho_l_ext=_fluides[1]->getDensity(pression_ext,T);
1567 for(k=0;k<1+_Ndim;k++){
1568 _externalStates[k]*=rho_v_ext/rho_v_int;
1569 _externalStates[k+1+_Ndim]*=rho_l_ext/rho_l_int;
1571 double alpha=_Vj[0];
1572 //cout<<"Outlet alpha= "<<alpha<<" pression int= "<<pression_int<<" pression ext= "<<pression_ext<<" temperature= "<<T<<" velocity gas "<<_Uj[2]<<" velocity liq "<<_Uj[2+_Ndim]<<endl;
1573 for(int idim=0;idim<_Ndim;idim++){
1574 v1_2+=_Vj[2+idim]*_Vj[2+idim];
1575 v2_2+=_Vj[2+_Ndim+idim]*_Vj[2+_Ndim+idim];
1577 _externalStates[_nVar-1]=alpha*rho_v_ext*(_fluides[0]->getInternalEnergy(T,rho_v_int)+v1_2/2)+(1-alpha)*rho_l_ext*(_fluides[1]->getInternalEnergy(T,rho_l_int)+v2_2/2);
1579 // _Vj external primitives
1580 _Vj[1] = pression_ext;
1583 for(k=1; k<_nVar; k++)
1584 _idm[k] = _idm[k-1] + 1;
1585 VecAssemblyBegin(_Uext);
1586 VecAssemblyBegin(_Vext);
1587 VecAssemblyBegin(_Uextdiff);
1588 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
1589 VecSetValues(_Vext, _nVar, _idm, _Vj, INSERT_VALUES);
1590 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
1591 VecAssemblyEnd(_Uext);
1592 VecAssemblyEnd(_Vext);
1593 VecAssemblyEnd(_Uextdiff);
1596 cout<<"Boundary condition not set for boundary named "<<nameOfGroup<<endl;
1597 cout<<"Accepted boundary condition are Neumann, Wall, Inlet, and Outlet"<<endl;
1598 throw CdmathException("Unknown boundary condition");
1602 void FiveEqsTwoFluid::addDiffusionToSecondMember
1607 double Tm=_Udiff[_nVar-1];
1608 double lambdal=_fluides[1]->getConductivity(Tm);
1609 double lambdav=_fluides[0]->getConductivity(Tm);
1610 double mu1 =_fluides[0]->getViscosity(Tm);
1611 double mu2 = _fluides[1]->getViscosity(Tm);
1613 if(mu1==0 && mu2 ==0 && lambdav==0 && lambdal==0 && _heatTransfertCoeff==0)
1616 //extraction des valeurs
1618 for(int k=1; k<_nVar; k++)
1619 _idm[k] = _idm[k-1] + 1;
1621 VecGetValues(_primitiveVars, _nVar, _idm, _Vi);
1622 if (_verbose && _nbTimeStep%_freqSave ==0)
1624 cout << "Contribution diffusion: variables primitives maille " << i<<endl;
1625 for(int q=0; q<_nVar; q++)
1627 cout << _Vi[q] << endl;
1633 for(int k=0; k<_nVar; k++)
1634 _idm[k] = _nVar*j + k;
1636 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
1640 lambdal=max(lambdal,_heatTransfertCoeff);//wall nucleate boing -> larger heat transfer
1642 for(int k=0; k<_nVar; k++)
1644 VecGetValues(_Uextdiff, _nVar, _idm, _phi);
1645 consToPrim(_phi,_Vj);
1648 if (_verbose && _nbTimeStep%_freqSave ==0)
1650 cout << "Contribution diffusion: variables primitives maille " <<j <<endl;
1651 for(int q=0; q<_nVar; q++)
1653 cout << _Vj[q] << endl;
1657 double alpha=(_Vj[0]+_Vi[0])/2;
1658 double lambda = alpha*lambdav+(1-alpha)*lambdal;
1659 //on n'a pas de contribution sur la masse
1662 //contribution visqueuse sur la quantite de mouvement
1663 for(int k=1; k<_Ndim+1; k++)
1665 _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
1666 _phi[k+_Ndim+1] = _inv_dxi*2/(1/_inv_dxi+1/_inv_dxj)*mu2*(1-alpha)*(_Vj[2+k+_Ndim] - _Vi[2+k+_Ndim]);
1668 _phi[_nVar-1] = _inv_dxi*2/(1/_inv_dxi+1/_inv_dxj)*lambda *(_Vj[_nVar-1] - _Vi[_nVar-1]);
1670 VecSetValuesBlocked(_b, 1, _idm, _phi, ADD_VALUES);
1672 if(_verbose && _nbTimeStep%_freqSave ==0)
1674 cout << "Contribution diffusion au 2nd membre pour la maille " << i << ": "<<endl;
1675 for(int q=0; q<_nVar; q++)
1677 cout << _phi[q] << endl;
1684 //On change de signe pour l'autre contribution
1685 for(int k=0; k<_nVar; k++)
1686 _phi[k] *= -_inv_dxj/_inv_dxi;
1689 VecSetValuesBlocked(_b, 1, _idm, _phi, ADD_VALUES);
1692 if(_verbose && _nbTimeStep%_freqSave ==0)
1694 cout << "Contribution diffusion au 2nd membre pour la maille " << j << ": "<<endl;
1695 for(int q=0; q<_nVar; q++)
1697 cout << _phi[q] << endl;
1701 if(_timeScheme==Implicit)
1703 cout << "Matrice de diffusion D, pour le couple (" << i << "," << j<< "):" << endl;
1704 for(int i=0; i<_nVar; i++)
1706 for(int j=0; j<_nVar; j++)
1707 cout << _Diffusion[i*_nVar+j]<<", ";
1715 void FiveEqsTwoFluid::jacobian(const int &j, string nameOfGroup, double * normale)
1718 for(k=0; k<_nVar*_nVar;k++)
1719 _Jcb[k] = 0;//No implicitation at this stage
1721 // loop of boundary types
1722 if (_limitField[nameOfGroup].bcType==Wall)
1724 for(k=0; k<_nVar;k++)
1725 _Jcb[k*_nVar + k] = 1;
1726 for(k=1; k<1+_Ndim;k++)
1727 for(int l=1; l<1+_Ndim;l++){
1728 _Jcb[k*_nVar + l] -= 2*normale[k-1]*normale[l-1];
1729 _Jcb[(k+1+_Ndim)*_nVar + l+1+_Ndim] -= 2*normale[k-1]*normale[l-1];
1733 else if (_limitField[nameOfGroup].bcType==Inlet)
1736 _Jcb[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1737 _Jcb[_nVar]=_limitField[nameOfGroup].v_x[0];
1738 _Jcb[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1741 _Jcb[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1742 _Jcb[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1745 _Jcb[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1746 _Jcb[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1751 // not wall, not inlet
1752 else if (_limitField[nameOfGroup].bcType==Outlet){
1754 for(k=1; k<_nVar;k++)
1755 {_idm[k] = _idm[k-1] + 1;}
1756 VecGetValues(_conservativeVars, _nVar, _idm, _phi);
1757 VecGetValues(_primitiveVars, _nVar, _idm, _externalStates);
1759 else if (_limitField[nameOfGroup].bcType!=Neumann && _limitField[nameOfGroup].bcType!=InletPressure)// not wall, not inlet, not outlet
1761 cout << "group named "<<nameOfGroup << " : unknown boundary condition" << endl;
1762 throw CdmathException("FiveEqs::jacobianDiff: This boundary condition is not treated");
1767 void FiveEqsTwoFluid::primToCons(const double *P, const int &i, double *W, const int &j){
1768 //P= alpha, p, u1, u2, T
1769 //W=m1,q1,m2,q2,rhoE =alpha1*rho1*(e1+u1^2/2)+alpha2*rho2*(e2+u2^2/2)
1770 double alpha=P[i*_nVar];
1771 double pression=P[i*_nVar+1];
1772 double temperature=P[i*_nVar+_nVar-1];
1773 double rho_v=_fluides[0]->getDensity(pression,temperature);
1774 double rho_l=_fluides[1]->getDensity(pression,temperature);
1775 double u1_sq=0, u2_sq=0;
1777 W[j*_nVar] = alpha*rho_v;
1778 W[j*_nVar+1+_Ndim] = (1-alpha)*rho_l;
1780 for(int k=0; k<_Ndim; k++)
1782 W[j*_nVar+(k+1)] = W[j*_nVar]*P[i*_nVar+(k+2)];//alpha1*rho1*u1
1783 W[j*_nVar+(k+1)+1+_Ndim] = W[j*_nVar+1+_Ndim]*P[i*_nVar+(k+2)+_Ndim];//alpha2*rho2*u2
1786 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);
1787 for(int k=0; k<_Ndim; k++){
1788 u1_sq+=P[i*_nVar+(k+2)]*P[i*_nVar+(k+2)];
1789 u2_sq+=P[i*_nVar+(k+2)+_Ndim]*P[i*_nVar+(k+2)+_Ndim];
1791 W[j*_nVar+_nVar-1] += (W[j*_nVar]*u1_sq+W[j*_nVar+1+_Ndim]*u2_sq)*0.5;
1794 void FiveEqsTwoFluid::consToPrim(const double *Wcons, double* Wprim,double porosity)//To do: treat porosity
1796 //Wprim= alpha, p, u1, u2, T
1797 //Wcons=m1,q1,m2,q2,rhoE
1798 double m_v=Wcons[0];
1799 double m_l=Wcons[1+_Ndim];
1800 double q1_sq = 0,q2_sq = 0;
1801 _minm1=min(m_v,_minm1);
1802 _minm2=min(m_l,_minm2);
1804 if(m_v<-_precision || m_l<-_precision){
1807 for(int k=0;k<_Ndim;k++){
1808 q1_sq += Wcons[k+1]*Wcons[k+1];
1809 q2_sq += Wcons[k+1+1+_Ndim]*Wcons[k+1+1+_Ndim];
1811 if(Wcons[0]>0)//_precision*_precision*_precision)
1812 q1_sq /= Wcons[0]; //alpha1 rho1 u1²
1815 if(Wcons[1+_Ndim]>0)//_precision*_precision*_precision)
1816 q2_sq /= Wcons[1+_Ndim]; //alpha2 rho2 u1²
1819 double rho_m_e_m=Wcons[_nVar-1] -0.5*(q1_sq+q2_sq);
1820 //calcul de la temperature et de la pression pour une loi stiffened gas
1821 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"));
1822 double e_v=static_cast<StiffenedGas*>(_fluides[0])->getInternalEnergy(temperature);
1823 double e_l=static_cast<StiffenedGas*>(_fluides[1])->getInternalEnergy(temperature);
1824 double gamma_v=_fluides[0]->constante("gamma");
1825 double gamma_l=_fluides[1]->constante("gamma");
1826 double Pinf_v=- gamma_v*_fluides[0]->constante("p0");
1827 double Pinf_l=- gamma_l*_fluides[1]->constante("p0");
1829 double b=-(Pinf_v+m_v*(gamma_v-1)*e_v+Pinf_l+m_l*(gamma_l-1)*e_l);
1830 double c=Pinf_v*Pinf_l+Pinf_v*m_l*(gamma_l-1)*e_l+ Pinf_l*m_v*(gamma_v-1)*e_v;
1831 double delta= b*b-4*a*c;
1833 cout<<"delta= "<<delta<<" <0"<<endl;
1834 *_runLogFile<<"delta= "<<delta<<" <0"<<endl;
1835 throw CdmathException("FiveEqsTwoFluid::consToPrim: Failed to compute pressure");
1837 double pression=(-b+sqrt(delta))/(2*a);
1839 cout << "pressure = "<< pression << " < 1 Pa " << endl;
1840 cout << "Conservative state = ";
1841 for(int k=0; k<_nVar; k++){
1842 cout<<Wcons[k]<<", ";
1845 *_runLogFile << "FiveEqsTwoFluid::consToPrim: Failed to compute pressure = "<< pression << " < 1 Pa " << endl;
1846 throw CdmathException("FiveEqsTwoFluid::consToPrim: Failed to compute pressure");
1849 double rho_v=_fluides[0]->getDensity(pression,temperature);
1850 double alpha=m_v/rho_v;
1852 Wprim[1] = pression;
1853 for(int k=0;k<_Ndim;k++){//vitesses
1854 if(Wcons[0]>0)//_precision*_precision*_precision)
1855 Wprim[k+2] = Wcons[k+1]/Wcons[0];
1857 Wprim[k+2] = Wcons[k+2+_Ndim]/Wcons[1+_Ndim];
1858 if(Wcons[1+_Ndim]>0)//_precision*_precision*_precision)
1859 Wprim[k+2+_Ndim] = Wcons[k+2+_Ndim]/Wcons[1+_Ndim];
1861 Wprim[k+2+_Ndim] = Wcons[k+1]/Wcons[0];
1863 Wprim[_nVar-1] = temperature;
1866 void FiveEqsTwoFluid::entropicShift(double* n, double& vpcorr0, double& vpcorr1)
1869 // parameters of function dgeev_ (compute the eigenvalues)
1870 int LDA, LDVL,LWORK, SDIM,LDVR;
1875 char jobvl[]="N", jobvr[]="N";
1876 double WORK[LWORK], JacoMat[_nVar*_nVar],egvaReal[_nVar],egvaImag[_nVar],egVectorL[_nVar*_nVar],egVectorR[_nVar*_nVar];
1877 int info_l = 0, info_r = 0;
1879 /******** Left: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1880 convectionJacobianMatrix(_l, n);
1881 for (int i=0; i<_nVar*_nVar; i++){
1882 JacoMat[i] = _JacoMat[i];
1884 dgeev_(jobvl, jobvl, &_nVar,
1885 JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1890 // /******** Right: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1891 convectionJacobianMatrix(_r, n);
1892 for (int i=0; i<_nVar*_nVar; i++){
1893 JacoMat[i] = _JacoMat[i];
1895 dgeev_(jobvl, jobvl, &_nVar,
1896 JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1901 if (info_l < 0 || info_r < 0)
1903 cout<<"Warning FiveEqsTwoFluid::entropicShift: dgeev_ did not compute all the eigenvalues, trying heuristic entropy correction "<<endl;
1904 double u1l_n=0, u2l_n=0, u1r_n=0, u2r_n=0;
1905 for (int idim=0; idim<_Ndim; idim++){
1906 u1l_n= _l[2+idim] *n[idim];
1907 u2l_n= _l[2+idim+_Ndim]*n[idim];
1908 u1r_n= _r[2+idim] *n[idim];
1909 u2r_n= _r[2+idim+_Ndim]*n[idim];
1912 vpcorr0 =max(fabs(u1l_n-u1r_n),fabs(u2l_n-u2r_n));
1917 std::vector< std::complex<double> > eigValuesLeft(_nVar);
1918 std::vector< std::complex<double> > eigValuesRight(_nVar);
1919 for(int j=0; j<_nVar; j++){
1920 eigValuesLeft[j] = complex<double>(egvaReal[j],egvaImag[j]);
1921 eigValuesRight[j] = complex<double>(egvaReal[j],egvaImag[j]);
1923 int sizeLeft = Polynoms::new_tri_selectif(eigValuesLeft, eigValuesLeft.size(), _precision);
1924 int sizeRight = Polynoms::new_tri_selectif(eigValuesRight, eigValuesRight.size(), _precision);
1925 if (_verbose && (_nbTimeStep-1)%_freqSave ==0)
1927 cout<<" Eigenvalue of JacoMat Left: " << endl;
1928 for(int i=0; i<sizeLeft; i++)
1929 cout<<eigValuesLeft[i] << ", "<<endl;
1931 if (_verbose && _nbTimeStep%_freqSave ==0)
1933 cout<<" Eigenvalue of JacoMat Right: " << endl;
1934 for(int i=0; i<sizeRight; i++)
1935 cout<<eigValuesRight[i] << ", "<<endl;
1938 for (int i=1; i<min(sizeLeft,sizeRight)-1; i++)
1939 vpcorr0 = max(vpcorr0, abs(eigValuesRight[i]-eigValuesLeft[i]));// Kieu
1944 void FiveEqsTwoFluid::entropicShift(double* n)
1947 // parameters of function dgeev_ (compute the eigenvalues)
1948 int LDA, LDVL,LWORK, SDIM,LDVR;
1953 char jobvl[]="N", jobvr[]="N";
1954 double WORK[LWORK], JacoMat[_nVar*_nVar],egvaReal[_nVar],egvaImag[_nVar],egVectorL[_nVar*_nVar],egVectorR[_nVar*_nVar];
1956 /******** Left: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1957 convectionJacobianMatrix(_l, n);
1958 for (int i=0; i<_nVar*_nVar; i++){
1959 JacoMat[i] = _JacoMat[i];
1961 dgeev_(jobvl, jobvl, &_nVar,
1962 JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1967 *_runLogFile<<"FiveEqsTwoFluid::JacoMat: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (left)"<<endl;
1968 throw CdmathException(
1969 "FiveEqsTwoFluid::JacoMat: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (left)");
1972 std::vector< std::complex<double> > eigValuesLeft(_nVar,0.);
1973 for(int j=0; j<_nVar; j++){
1974 eigValuesLeft[j] = complex<double>(egvaReal[j],egvaImag[j]);
1976 // /******** Right: Compute the eigenvalues and eigenvectors of Jacobian Matrix (using lapack)********/
1977 convectionJacobianMatrix(_r, n);
1978 for (int i=0; i<_nVar*_nVar; i++){
1979 JacoMat[i] = _JacoMat[i];
1981 dgeev_(jobvl, jobvl, &_nVar,
1982 JacoMat,&LDA,egvaReal,egvaImag, egVectorL,
1988 *_runLogFile<<"FiveEqsTwoFluid::entropicShift: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (right)"<<endl;
1989 throw CdmathException(
1990 "FiveEqsTwoFluid::entropicShift: dgeev_ unsuccessful computation of the eigenvalues of Jacobian Matrix (right)");
1993 std::vector< std::complex<double> > eigValuesRight(_nVar,0.);
1994 for(int j=0; j<_nVar; j++){
1995 eigValuesRight[j] = complex<double>(egvaReal[j],egvaImag[j]);
1997 int sizeLeft = Polynoms::new_tri_selectif(eigValuesLeft, eigValuesLeft.size(), _precision);
1998 int sizeRight = Polynoms::new_tri_selectif(eigValuesRight, eigValuesRight.size(), _precision);
1999 if (_verbose && (_nbTimeStep-1)%_freqSave ==0)
2001 cout<<" Eigenvalue of JacoMat Left: " << endl;
2002 for(int i=0; i<sizeLeft; i++)
2003 cout<<eigValuesLeft[i] << ", "<<endl;
2004 cout<<" Eigenvalue of JacoMat Right: " << endl;
2005 for(int i=0; i<sizeRight; i++)
2006 cout<<eigValuesRight[i] << ", "<<endl;
2008 for (int i=0; i<min(sizeLeft,sizeRight)-1; i++)
2009 _entropicShift[i]= abs(eigValuesRight[i]-eigValuesLeft[i]);
2012 Vector FiveEqsTwoFluid::staggeredVFFCFlux()
2014 if(_spaceScheme!=staggered || _nonLinearFormulation!=VFFC)
2015 throw CdmathException("IsothermalTwoFluid::staggeredVFFCFlux: staggeredVFFCFlux method should be called only for VFFC formulation and staggered upwinding");
2016 else//_spaceScheme==staggered
2020 double alpha_roe = _Uroe[0];//Toumi formula
2021 // interfacial pressure term (hyperbolic correction)
2022 double dpi = _Uroe[_nVar];
2024 double u1ijn=0, u2ijn=0, phialphaq1n=0, phialphaq2n=0,vitesse1n=0,vitesse2n=0;
2025 for(int idim=0;idim<_Ndim;idim++){//URoe = alpha, p, u1, u2, T, dpi
2026 u1ijn+=_vec_normal[idim]*_Uroe[2+idim];
2027 u2ijn+=_vec_normal[idim]*_Uroe[2+_Ndim+idim];
2031 for(int idim=0;idim<_Ndim;idim++){
2032 phialphaq1n+=_vec_normal[idim]*_Ui[1+idim];//phi alpha rho u n
2033 vitesse1n +=_vec_normal[idim]*_Vi[2+idim];
2036 for(int idim=0;idim<_Ndim;idim++)
2037 Fij(1+idim)=phialphaq1n*_Vi[2+idim]+(alpha_roe*_Vj[1]*_porosityj+dpi*_Vi[0]*_porosityi)*_vec_normal[idim];
2039 double pressioni=_Vi[1];
2040 double Temperaturei= _Vi[_nVar-1];
2041 double rho1=_fluides[0]->getDensity(pressioni,Temperaturei);
2042 double e1_int=_fluides[0]->getInternalEnergy(Temperaturei,rho1);
2043 Fij(_nVar-1)+=_Ui[0]*(e1_int+0.5*vitesse1n*vitesse1n+_Vj[1]/rho1)*vitesse1n;
2047 for(int idim=0;idim<_Ndim;idim++){
2048 phialphaq2n+=_vec_normal[idim]*_Uj[1+idim];//phi alpha rho u n
2049 vitesse1n +=_vec_normal[idim]*_Vj[2+idim];
2052 for(int idim=0;idim<_Ndim;idim++)
2053 Fij(1+idim)=phialphaq2n*_Vj[2+idim]+(alpha_roe*_Vi[1]*_porosityi+dpi*_Vj[0]*_porosityj)*_vec_normal[idim];
2055 double pressionj=_Vj[1];
2056 double Temperaturej= _Vj[_nVar-1];
2057 double rho1=_fluides[0]->getDensity(pressionj,Temperaturej);
2058 double e1_int=_fluides[0]->getInternalEnergy(Temperaturej,rho1);
2059 Fij(_nVar-1)+=_Uj[0]*(e1_int+0.5*vitesse1n*vitesse1n+_Vi[1]/rho1)*vitesse1n;
2064 for(int idim=0;idim<_Ndim;idim++){
2065 phialphaq2n+=_vec_normal[idim]*_Ui[2+_Ndim+idim];//phi alpha rho u n
2066 vitesse2n +=_vec_normal[idim]*_Vi[2+idim+_Ndim];
2068 Fij(1+_Ndim)=phialphaq2n;
2069 for(int idim=0;idim<_Ndim;idim++)
2070 Fij(2+_Ndim+idim)=phialphaq2n*_Vi[2+_Ndim+idim]+((1-alpha_roe)*_Vj[1]*_porosityj-dpi*_Vi[0]*_porosityi)*_vec_normal[idim];
2072 double pressioni=_Vi[1];
2073 double Temperaturei= _Vi[_nVar-1];
2074 double rho2=_fluides[1]->getDensity(pressioni,Temperaturei);
2075 double e2_int=_fluides[1]->getInternalEnergy(Temperaturei,rho2);
2076 Fij(_nVar-1)+=_Ui[1+_Ndim]*(e2_int+0.5*vitesse2n*vitesse2n+_Vj[1]/rho2)*vitesse2n;
2080 for(int idim=0;idim<_Ndim;idim++){
2081 phialphaq2n+=_vec_normal[idim]*_Uj[2+_Ndim+idim];//phi alpha rho u n
2082 vitesse2n +=_vec_normal[idim]*_Vj[2+idim+_Ndim];
2084 Fij(1+_Ndim)=phialphaq2n;
2085 for(int idim=0;idim<_Ndim;idim++)
2086 Fij(2+_Ndim+idim)=phialphaq2n*_Vj[2+_Ndim+idim]+((1-alpha_roe)*_Vi[1]*_porosityi-dpi*_Vj[0]*_porosityj)*_vec_normal[idim];
2088 double pressionj=_Vj[1];
2089 double Temperaturej= _Vj[_nVar-1];
2090 double rho2=_fluides[1]->getDensity(pressionj,Temperaturej);
2091 double e2_int=_fluides[1]->getInternalEnergy(Temperaturej,rho2);
2092 Fij(_nVar-1)+=_Uj[1+_Ndim]*(e2_int+0.5*vitesse2n*vitesse2n+_Vi[1]/rho2)*vitesse2n;
2098 void FiveEqsTwoFluid::applyVFRoeLowMachCorrections(bool isBord, string groupname)
2100 if(_nonLinearFormulation!=VFRoe)
2101 throw CdmathException("FiveEqsTwoFluid::applyVFRoeLowMachCorrections: applyVFRoeLowMachCorrections method should be called only for VFRoe formulation");
2102 else//_nonLinearFormulation==VFRoe
2104 if(_spaceScheme==lowMach){
2105 double u1_2=0, u2_2=0;
2106 for(int i=0;i<_Ndim;i++){
2107 u1_2 += _Uroe[2+i]*_Uroe[2+i];
2108 u2_2 += _Uroe[2+i+_Ndim]*_Uroe[2+i+_Ndim];
2111 double c = _maxvploc;//mixture sound speed
2112 double M=max(sqrt(u1_2),sqrt(u2_2))/c;//Mach number
2113 _Vij[1]=M*_Vij[1]+(1-M)*(_Vi[1]+_Vj[1])/2;
2114 primToCons(_Vij,0,_Uij,0);
2116 else if(_spaceScheme==pressureCorrection)
2118 if(_pressureCorrectionOrder>2)
2119 throw CdmathException("FiveEqsTwoFluid::applyVFRoeLowMachCorrections pressure correction order can be only 1 or 2 for five equation two-fluid model");
2121 double norm_uij=0, uij_n=0, ui_n=0, uj_n=0;//mean velocities
2122 double rho1 = _fluides[0]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2123 double rho2 = _fluides[1]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2124 double m1=_Uroe[0]*rho1, m2=(1-_Uroe[0])*rho2;
2126 for(int i=0;i<_Ndim;i++)
2128 norm_uij += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*(m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim]);
2129 uij_n += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*_vec_normal[i];
2130 ui_n += _Vi[2+i]*_vec_normal[i];
2131 uj_n += _Vj[2+i]*_vec_normal[i];
2133 norm_uij=sqrt(norm_uij)/rhom;
2135 if(norm_uij>_precision)//avoid division by zero
2136 _Vij[1]=(_Vi[1]+_Vj[1])/2 + uij_n/norm_uij*(_Vj[1]-_Vi[1])/4 - rhom*norm_uij*(uj_n-ui_n)/4;
2138 _Vij[1]=(_Vi[1]+_Vj[1])/2 - rhom*norm_uij*(uj_n-ui_n)/4;
2140 else if(_spaceScheme==staggered)
2143 double rho1 = _fluides[0]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2144 double rho2 = _fluides[1]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
2145 double m1=_Uroe[0]*rho1, m2=(1-_Uroe[0])*rho2;
2147 for(int i=0;i<_Ndim;i++)
2148 qij_n += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*_vec_normal[i];
2153 for(int i=0;i<_Ndim;i++)
2155 _Vij[2+i] =_Vi[2+i];
2156 _Vij[2+i+_Ndim] =_Vi[2+i+_Ndim];
2158 _Vij[_nVar-1]=_Vi[_nVar-1];
2163 for(int i=0;i<_Ndim;i++)
2165 _Vij[2+i] =_Vj[2+i];
2166 _Vij[2+i+_Ndim] =_Vj[2+i+_Ndim];
2168 _Vij[_nVar-1]=_Vj[_nVar-1];
2170 primToCons(_Vij,0,_Uij,0);
2175 void FiveEqsTwoFluid::computeScaling(double maxvp)
2177 _blockDiag[0]=1;//alphaScaling;
2178 _invBlockDiag[0]=1;//_blockDiag[0];
2179 _blockDiag[1+_Ndim]=1;//-alphaScaling;
2180 _invBlockDiag[1+_Ndim]=1.0;//_blockDiag[1+_Ndim];
2181 for(int q=1; q<_Ndim+1; q++)
2183 _blockDiag[q]=1/maxvp;
2184 _invBlockDiag[q]=1/_blockDiag[q];
2185 _blockDiag[q+1+_Ndim]=1/maxvp;
2186 _invBlockDiag[q+1+_Ndim]=1/_blockDiag[q+1+_Ndim];
2188 _blockDiag[_nVar - 1]=1/(maxvp*maxvp);//1
2189 _invBlockDiag[_nVar - 1]= 1./_blockDiag[_nVar - 1] ;// 1.;//
2192 void FiveEqsTwoFluid::testConservation()
2194 double SUM, DELTA, x;
2196 for(int i=0; i<_nVar; i++)
2200 cout << "Masse totale phase " << 0 <<" (kg): ";
2201 else if( i == 1+_Ndim)
2202 cout << "Masse totale phase " << 1 <<" (kg): ";
2203 else if( i == _nVar-1)
2204 cout << "Energie totale " <<" (J.m^-3): ";
2208 cout << "Quantite de mouvement totale phase 0 (kg.m.s^-1): ";
2210 cout << "Quantite de mouvement totale phase 1 (kg.m.s^-1): ";
2216 for(int j=0; j<_Nmailles; j++)
2218 VecGetValues(_old, 1, &I, &x);//on recupere la valeur du champ
2219 SUM += x*_mesh.getCell(j).getMeasure();
2220 VecGetValues(_newtonVariation, 1, &I, &x);//on recupere la variation du champ
2221 DELTA += x*_mesh.getCell(j).getMeasure();
2224 if(fabs(SUM)>_precision)
2225 cout << SUM << ", variation relative: " << fabs(DELTA /SUM) << endl;
2227 cout << " a une somme nulle, variation absolue: " << fabs(DELTA) << endl;
2231 void FiveEqsTwoFluid::save(){
2232 string prim(_path+"/FiveEqsTwoFluidPrim_");
2233 string cons(_path+"/FiveEqsTwoFluidCons_");
2238 for (PetscInt i = 0; i < _Nmailles; i++){
2239 /* j = 0 : void fraction
2241 j = 2, 3, 4: velocity phase 1
2242 j = 5, 6, 7: velocity phase 2
2243 j = 8 : temperature */
2244 for (int j = 0; j < _nVar; j++){
2246 VecGetValues(_primitiveVars,1,&Ii,&_VV(i,j));
2249 if(_saveConservativeField){
2250 for (long i = 0; i < _Nmailles; i++){
2251 for (int j = 0; j < _nVar; j++){
2253 VecGetValues(_conservativeVars,1,&Ii,&_UU(i,j));
2256 _UU.setTime(_time,_nbTimeStep+1);
2258 _VV.setTime(_time,_nbTimeStep+1);
2260 if (_nbTimeStep ==0 || _restartWithNewFileName){
2261 string prim_suppress ="rm -rf "+prim+"_*";
2262 string cons_suppress ="rm -rf "+cons+"_*";
2263 system(prim_suppress.c_str());//Nettoyage des précédents calculs identiques
2264 system(cons_suppress.c_str());//Nettoyage des précédents calculs identiques
2265 _VV.setInfoOnComponent(0,"Void_fraction");
2266 _VV.setInfoOnComponent(1,"Pressure_(Pa)");
2267 _VV.setInfoOnComponent(2,"Velocity1_x_m/s");
2270 _VV.setInfoOnComponent(3,"Velocity1_y_m/s");
2272 _VV.setInfoOnComponent(4,"Velocity1_z_m/s");
2273 _VV.setInfoOnComponent(2+_Ndim,"Velocity2_x_m/s");
2275 _VV.setInfoOnComponent(3+_Ndim,"Velocity2_y_m/s");
2277 _VV.setInfoOnComponent(4+_Ndim,"Velocity2_z_m/s");
2278 _VV.setInfoOnComponent(_nVar-1,"Temperature_(K)");
2292 if(_saveConservativeField){
2293 _UU.setInfoOnComponent(0,"Partial_density1");// (kg/m^3)
2294 _UU.setInfoOnComponent(1,"Momentum1_x");// phase1 (kg/m^2/s)
2296 _UU.setInfoOnComponent(2,"Momentum1_y");// phase1 (kg/m^2/s)
2298 _UU.setInfoOnComponent(3,"Momentum1_z");// phase1 (kg/m^2/s)
2299 _UU.setInfoOnComponent(1+_Ndim,"Partial_density2");// phase2 (kg/m^3)
2300 _UU.setInfoOnComponent(2+_Ndim,"Momentum2_x");// phase2 (kg/m^2/s)
2303 _UU.setInfoOnComponent(3+_Ndim,"Momentum2_y");// phase2 (kg/m^2/s)
2305 _UU.setInfoOnComponent(4+_Ndim,"Momentum2_z");// phase2 (kg/m^2/s)
2306 _UU.setInfoOnComponent(_nVar-1,"Total_energy");
2326 _VV.writeVTK(prim,false);
2329 _VV.writeMED(prim,false);
2336 if(_saveConservativeField){
2340 _UU.writeVTK(cons,false);
2343 _UU.writeMED(cons,false);
2352 for (long i = 0; i < _Nmailles; i++){
2353 // j = 0 : concentration, j=1 : pressure; j = _nVar - 1: temperature; j = 2,..,_nVar-2: velocity
2354 for (int j = 0; j < _Ndim; j++)//On récupère les composantes de vitesse
2357 VecGetValues(_primitiveVars,1,&Ii,&_Vitesse1(i,j));
2358 Ii=i*_nVar +2+j+_Ndim;
2359 VecGetValues(_primitiveVars,1,&Ii,&_Vitesse2(i,j));
2361 for (int j = _Ndim; j < 3; j++){//On met à zero les composantes de vitesse si la dimension est <3
2366 _Vitesse1.setTime(_time,_nbTimeStep);
2367 _Vitesse2.setTime(_time,_nbTimeStep);
2368 if (_nbTimeStep ==0 || _restartWithNewFileName){
2369 _Vitesse1.setInfoOnComponent(0,"Velocity_x_(m/s)");
2370 _Vitesse1.setInfoOnComponent(1,"Velocity_y_(m/s)");
2371 _Vitesse1.setInfoOnComponent(2,"Velocity_z_(m/s)");
2373 _Vitesse2.setInfoOnComponent(0,"Velocity_x_(m/s)");
2374 _Vitesse2.setInfoOnComponent(1,"Velocity_y_(m/s)");
2375 _Vitesse2.setInfoOnComponent(2,"Velocity_z_(m/s)");
2380 _Vitesse1.writeVTK(prim+"_GasVelocity");
2381 _Vitesse2.writeVTK(prim+"_LiquidVelocity");
2384 _Vitesse1.writeMED(prim+"_GasVelocity");
2385 _Vitesse2.writeMED(prim+"_LiquidVelocity");
2388 _Vitesse1.writeCSV(prim+"_GasVelocity");
2389 _Vitesse2.writeCSV(prim+"_LiquidVelocity");
2397 _Vitesse1.writeVTK(prim+"_GasVelocity",false);
2398 _Vitesse2.writeVTK(prim+"_LiquidVelocity",false);
2401 _Vitesse1.writeMED(prim+"_GasVelocity",false);
2402 _Vitesse2.writeMED(prim+"_LiquidVelocity",false);
2405 _Vitesse1.writeCSV(prim+"_GasVelocity");
2406 _Vitesse2.writeCSV(prim+"_LiquidVelocity");
2412 if (_restartWithNewFileName)
2413 _restartWithNewFileName=false;