2 * IsothermalTwoFluid.cxx
4 * Created on: Sep 16, 2014
8 #include "IsothermalTwoFluid.hxx"
12 IsothermalTwoFluid::IsothermalTwoFluid(pressureEstimate pEstimate, int dim){
16 _dragCoeffs=vector<double>(2,0);
18 if (pEstimate==around1bar300K)//EOS at 1 bar and 300K
20 cout<<"Fluid is air-water mixture around 1 bar and 300 K (27°C)"<<endl;
21 *_runLogFile<<"Fluid is air-water mixture around 1 bar and 300 K (27°C)"<<endl;
22 _Temperature=300;//Constant temperature of the model
23 _internalEnergy1=2.22e5;//nitrogen internal energy at 1bar, 300K
24 _internalEnergy2=1.12e5;//water internal energy at 1 bar, 300K
25 _fluides[0] = new StiffenedGas(1.4,743,_Temperature,_internalEnergy1); //ideal gas law for nitrogen at pressure 1 bar and temperature 27°C, c_v=743
26 _fluides[1] = new StiffenedGas(996,1e5,_Temperature,_internalEnergy2,1501,4130); //stiffened gas law for water at pressure 1 bar and temperature 27°C
28 else//EOS at 155 bars and 618K
30 cout<<"Fluid is water-Gas mixture around saturation point 155 bars and 618 K (345°C)"<<endl;
31 *_runLogFile<<"Fluid is water-Gas mixture around saturation point 155 bars and 618 K (345°C)"<<endl;
32 _Temperature=618;//Constant temperature of the model
33 _internalEnergy1=2.44e6;//Gas internal energy at saturation at 155 bar
34 _internalEnergy2=1.6e6;//water internal energy at saturation at 155 bar
35 _fluides[0] = new StiffenedGas(102,1.55e7,_Temperature,_internalEnergy1, 433,3633); //stiffened gas law for Gas at pressure 155 bar and temperature 345°C:
36 _fluides[1] = new StiffenedGas(594,1.55e7,_Temperature,_internalEnergy2, 621,3100); //stiffened gas law for water at pressure 155 bar and temperature 345°C:
41 void IsothermalTwoFluid::initialize(){
42 cout<<"Initialising the isothermal two-fluid model"<<endl;
43 *_runLogFile<<"Initialising the isothermal two-fluid model"<<endl;
45 _Uroe = new double[_nVar+1];
47 _guessalpha = _VV(0,0);
49 _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
50 for(int i=0; i<_Ndim; i++)
52 _gravite[i+1]=_GravityField3d[i];
53 _gravite[i+1 +_Ndim+1]=_GravityField3d[i];
55 _GravityImplicitationMatrix = new PetscScalar[_nVar*_nVar];
58 _Vitesse1=Field("Gas velocity",CELLS,_mesh,3);//Forcement en dimension 3 pour le posttraitement des lignes de courant
59 _Vitesse2=Field("Liquid velocity",CELLS,_mesh,3);//Forcement en dimension 3 pour le posttraitement des lignes de courant
61 if(_entropicCorrection)
62 _entropicShift=vector<double>(3,0);
64 ProblemFluid::initialize();
67 void IsothermalTwoFluid::convectionState( const long &i, const long &j, const bool &IsBord){
68 //sortie: WRoe en (alpha, p, u1, u2, dm1,dm2,dalpha1,dp)z
69 //entree: _conservativeVars en (alpha1 rho1, alpha1 rho1 u1, alpha2 rho2, alpha2 rho2 u2)
71 // _l always inside the domain (index i)
72 // _r is maybe the boundary cell (negative index)
73 // _l and _r are primative vectors (alp, P, u1, u2)
76 for(int k=1; k<_nVar; k++)
77 _idm[k] = _idm[k-1] + 1;
78 VecGetValues(_conservativeVars, _nVar, _idm, _Ui);
81 for(int k=1; k<_nVar; k++)
82 _idm[k] = _idm[k-1] + 1;
84 VecGetValues(_Uext, _nVar, _idm, _Uj);
86 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
87 if(_verbose && _nbTimeStep%_freqSave ==0)
89 cout<<"Convection Left state cell " << i<< ": "<<endl;
90 for(int k =0; k<_nVar; k++)
92 cout<<"Convection Right state cell " << j<< ": "<<endl;
93 for(int k =0; k<_nVar; k++)
96 // if(_Ui[0]<-(_precision) || _Uj[0]<-(_precision) || _Ui[_Ndim+1]<-(_precision) || _Uj[_Ndim+1]<-(_precision))
98 // cout<<"!!!!!!!!!!!!!!!!!!!!!!!! Masse partielle negative, arret de calcul!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!"<<endl;
99 // cout<< "valeurs a gauche: "<<_Ui[0]<<", "<<_Ui[_Ndim+1]<<"valeurs a droite: "<<_Uj[0]<<", "<<_Uj[_Ndim+1]<<endl;
100 // throw CdmathException(" Masse partielle negative, arret de calcul");
104 _Ui[0]=max(0.,_Ui[0]);// mass1 a gauche
105 _Uj[0]=max(0.,_Uj[0]);// mass1 a droite
106 _Ui[_Ndim+1]=max(0.,_Ui[_Ndim+1]);// mass2 a gauche
107 _Uj[_Ndim+1]=max(0.,_Uj[_Ndim+1]);// mass2 a droite
110 PetscScalar ri1, ri2, rj1, rj2, xi, xj;
112 for(int k=1; k<_nVar; k++)
113 _idm[k] = _idm[k-1] + 1;
114 VecGetValues(_primitiveVars, _nVar, _idm, _l);
116 if(_verbose && _nbTimeStep%_freqSave ==0)
118 cout<<"Etat de Roe, etat primitif gauche: "<<endl;
119 for(int i =0; i<_nVar; i++)
123 // boundary : compute the _r
129 // inside the domain : extract from the primative vector
133 for(int k=1; k<_nVar; k++)
134 _idm[k] = _idm[k-1] + 1;
135 VecGetValues(_primitiveVars, _nVar, _idm, _r);
138 if(_verbose && _nbTimeStep%_freqSave ==0)
140 cout<<"Etat de Roe, etat primitif droite: "<<endl;
141 for(int i =0; i<_nVar; i++)
145 // Using Toumi linearisation (read in the article of Toumi) : alpha^{Roe}_l
146 if(2-_l[0]-_r[0] > _precision)
147 _Uroe[0] = 1- 2*(1-_l[0])*(1-_r[0])/(2-_l[0]-_r[0]);
148 // Using an average : (alp_l + alp_r)/2 : suggestion of Michael (no theory)
150 _Uroe[0] = (_l[0]+_r[0])/2;
151 // Pressure is computed as function of alp and P_l, P_r (Toumi article)
152 if(_l[0]+_r[0] > _precision)
153 _Uroe[1] = (_l[1]*_l[0]+_r[1]*_r[0])/(_l[0]+_r[0]);
155 _Uroe[1] = (_l[1]*(1-_l[0])+_r[1]*(1-_r[0]))/(2-_l[0]-_r[0]);
156 // i :left, j : right (U is normally conservative variable)
157 ri1 = sqrt(_Ui[0]); ri2 = sqrt(_Ui[_Ndim+1]);
158 rj1 = sqrt(_Uj[0]); rj2 = sqrt(_Uj[_Ndim+1]);
159 // Roe average formula of the velocities
160 for(int k=0;k<_Ndim;k++)
164 // avoid dividing by zero, if mass is zero, do not consider the distribution of such a phase
165 if(ri1>_precision && rj1>_precision)
166 _Uroe[2+k] = (xi/ri1 + xj/rj1)/(ri1 + rj1);
167 else if(ri1<_precision && rj1>_precision)
168 _Uroe[2+k] = xj/_Uj[0];
169 else if(ri1>_precision && rj1<_precision)
170 _Uroe[2+k] = xi/_Ui[0];
172 _Uroe[2+k] =(_Ui[k+1+_Ndim+1]/ri2 + _Uj[k+1+_Ndim+1]/rj2)/(ri2 + rj2);
174 xi = _Ui[k+1+_Ndim+1];
175 xj = _Uj[k+1+_Ndim+1];
176 if(ri2>_precision && rj2>_precision)
177 _Uroe[2+k+_Ndim] = (xi/ri2 + xj/rj2)/(ri2 + rj2);
178 else if(ri2<_precision && rj2>_precision)
179 _Uroe[2+k+_Ndim] = xj/_Uj[_Ndim+1];
180 else if(ri2>_precision && rj2<_precision)
181 _Uroe[2+k+_Ndim] = xi/_Ui[_Ndim+1];
183 _Uroe[2+k+_Ndim] = (xi/ri1 + xj/rj1)/(ri1 + rj1);
185 if(_verbose && _nbTimeStep%_freqSave ==0)
187 cout<<"Etat de Roe calcule: "<<endl;
188 for(int k=0;k<_nVar; k++)//At this point _Uroe[_nVar] is not yet set
189 cout<< _Uroe[k]<<endl;
193 void IsothermalTwoFluid::diffusionStateAndMatrices(const long &i,const long &j, const bool &IsBord){
194 //sortie: matrices et etat de diffusion (m_v, q_v, m_l, q_l)
196 for(int k=1; k<_nVar; k++)
197 _idm[k] = _idm[k-1] + 1;
199 VecGetValues(_conservativeVars, _nVar, _idm, _Ui);
201 for(int k=1; k<_nVar; k++)
202 _idm[k] = _idm[k-1] + 1;
205 VecGetValues(_Uextdiff, _nVar, _idm, _Uj);
207 VecGetValues(_conservativeVars, _nVar, _idm, _Uj);
209 for(int k=0; k<_nVar; k++)
210 _Udiff[k] = (_Ui[k]+_Uj[k])/2;
212 for (int i = 0; i<_Ndim;i++)
213 q_2+=_Udiff[i+1]*_Udiff[i+1];
215 if(_timeScheme==Implicit)
217 for(int i=0; i<_nVar*_nVar;i++)
219 double mu1 = _fluides[0]->getViscosity(_Temperature);
220 double mu2 = _fluides[1]->getViscosity(_Temperature);
221 for(int idim=1;idim<_Ndim+1;idim++)
223 _Diffusion[idim*_nVar] = mu1*_Udiff[idim]/(_Udiff[0]*_Udiff[0]);
224 _Diffusion[idim*_nVar+idim] = -mu1/_Udiff[0];
225 _Diffusion[(idim+_Ndim+1)*_nVar] = mu2*_Udiff[idim+_Ndim+1]/(_Udiff[_Ndim+1]*_Udiff[_Ndim+1]);
226 _Diffusion[(idim+_Ndim+1)*_nVar+idim+_Ndim+1] = -mu2/_Udiff[_Ndim+1];
231 void IsothermalTwoFluid::convectionMatrices()
233 //entree: URoe = alpha, p, u1, u2 + ajout dpi pour calcul flux ultérieur
234 //sortie: matrices Roe+ et Roe- +Roe si scheme centre
236 if(_timeScheme==Implicit && _usePrimitiveVarsInNewton)
237 throw CdmathException("Implicitation with primitive variables not yet available for IsothermalTwoFluid model");
240 complex< double > tmp;
241 double u1_n, u1_2, u2_n, u2_2, u_r2;
246 for(int i=0;i<_Ndim;i++)
248 u1_2 += _Uroe[2+i]*_Uroe[2+i];
249 u1_n += _Uroe[2+i]*_vec_normal[i];
250 u2_2 += _Uroe[2+i+_Ndim]*_Uroe[2+i+_Ndim];
251 u2_n += _Uroe[2+i+_Ndim]*_vec_normal[i];
252 u_r2 += (_Uroe[2+i]-_Uroe[2+i+_Ndim])*(_Uroe[2+i]-_Uroe[2+i+_Ndim]);
254 //Ancienne construction Mat Roe (Dm1,Dm2,Dalp,Dp)
255 double alpha = _Uroe[0];
257 double rho1 = _fluides[0]->getDensity(p, _Temperature);
258 double rho2 = _fluides[1]->getDensity(p, _Temperature);
259 double dpi1 = intPressDef(alpha, u_r2, rho1, rho2);
261 double invcsq1 = 1/_fluides[0]->vitesseSonTemperature(_Temperature,rho1);
263 double invcsq2 = 1/_fluides[1]->vitesseSonTemperature(_Temperature,rho2);
265 double g2 = 1/(alpha*rho2*invcsq1+(1-alpha)*rho1*invcsq2);
266 double g2press=g2, g2alpha=g2;
267 //saving dpi value for flux calculation later
270 /***********Calcul des valeurs propres ********/
272 vector< double > pol_car= Poly.polynome_caracteristique(alpha, 1-alpha, u1_n, u2_n, rho1, rho2,invcsq1,invcsq2, dpi1, dpi2, g2press, g2alpha, g2,_precision);
273 for(int ct=0;ct<4;ct++){
274 if (abs(pol_car[ct])<_precision*_precision)
277 if(_verbose && _nbTimeStep%_freqSave ==0)
279 cout<<"pol caract= "<<endl;
280 for(int i =0; i<5; i++)
281 cout<<pol_car[i]<<" ";
283 cout<<"alpha= "<<alpha<<", p= " << p << ", rho1= " << rho1<< ", rho2= " << rho2<< ", c1= " <<sqrt(1/invcsq1)<<
284 ", c2= " <<sqrt(1/invcsq2)<<endl;
285 cout<< "u1_n= "<<u1_n<<", u2_n= "<<u2_n<< ", dpi1= "<<dpi1<< ", dpi2= "<<dpi2<< endl;
287 vector< complex<double> > valeurs_propres = getRacines(pol_car);
289 //On ajoute les valeurs propres triviales
292 if( !Poly.belongTo(u1_n,valeurs_propres, _precision) )
293 valeurs_propres.push_back(u1_n);//vp vapor energy
294 if( !Poly.belongTo(u2_n,valeurs_propres, _precision) )
295 valeurs_propres.push_back(u2_n);//vp liquid energy
297 bool doubleeigenval = norm(valeurs_propres[0] - valeurs_propres[1])<_precision;//norm= suqare of the magnitude
300 valeurs_propres[0] = valeurs_propres[valeurs_propres.size()-1];
301 valeurs_propres.pop_back();
304 int taille_vp = Poly.new_tri_selectif(valeurs_propres,valeurs_propres.size(),_precision);//valeurs_propres.size();//
307 for(int i =0; i<taille_vp; i++)
308 if(fabs(valeurs_propres[i].real())>_maxvploc)
309 _maxvploc=fabs(valeurs_propres[i].real());
313 int existVpCplx = 0,pos_conj;
315 for (int ct=0; ct<taille_vp; ct++) {
316 vp_imag_iter = valeurs_propres[ct].imag();
317 if ( fabs(vp_imag_iter) > _precision ) {
319 if ( _part_imag_max < fabs(vp_imag_iter))
320 _part_imag_max = fabs(vp_imag_iter);
321 //On cherhe le conjugue
323 while(pos_conj<taille_vp && fabs(valeurs_propres[pos_conj].imag()+vp_imag_iter)>_precision)
325 if(pos_conj!=ct+1 && pos_conj<taille_vp )
327 tmp=valeurs_propres[ct+1];
328 valeurs_propres[ct+1]=valeurs_propres[pos_conj];
329 valeurs_propres[pos_conj] = tmp;
337 //on ordonne les deux premieres valeurs
339 if(valeurs_propres[1].real()<valeurs_propres[0].real())
341 tmp=valeurs_propres[0];
342 valeurs_propres[0]=valeurs_propres[1];
343 valeurs_propres[1]=tmp;
346 if(_verbose && _nbTimeStep%_freqSave ==0)
348 cout<<" Vp apres tri " << valeurs_propres.size()<<endl;
349 for(int ct =0; ct<taille_vp; ct++)
350 cout<< "("<<valeurs_propres[ct].real()<< ", " <<valeurs_propres[ct].imag() <<") ";
354 /******** Construction de la matrice de Roe *********/
356 for(int i=0; i<_nVar*_nVar;i++)
359 for(int idim=0; idim<_Ndim;idim++)
361 _Aroe[1+idim]=_vec_normal[idim];
362 _Aroe[1+idim+_Ndim+1]=0;
363 _Aroe[(_Ndim+1)*_nVar+1+idim]=0;
364 _Aroe[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=_vec_normal[idim];
367 for(int idim=0; idim<_Ndim;idim++)
369 //premiere colonne (masse gaz)
370 _Aroe[ (1+idim)*_nVar]= (alpha *rho2*g2press+dpi1*(1-alpha)*invcsq2*g2alpha)*_vec_normal[idim] - u1_n*_Uroe[2+idim];
371 _Aroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar]=((1-alpha)*rho2*g2press-dpi2*(1-alpha)*invcsq2*g2alpha)*_vec_normal[idim];
372 //colonnes intermediaires
373 for(int jdim=0; jdim<_Ndim;jdim++)
375 _Aroe[ (1+idim)*_nVar + jdim + 1] = _Uroe[ 2+idim]*_vec_normal[jdim];
376 _Aroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar + jdim + 1+_Ndim+1] = _Uroe[_Ndim+2+idim]*_vec_normal[jdim];
379 _Aroe[ (1+idim)*_nVar + idim + 1] += u1_n;
380 _Aroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar + _Ndim+1+ idim + 1] += u2_n;
381 //troisieme colonne (masse liquide)
382 _Aroe[ (1+idim)*_nVar + _Ndim+1]= (alpha *rho1*g2press-dpi1*alpha*invcsq1*g2alpha)*_vec_normal[idim];
383 _Aroe[(_Ndim+1)*_nVar+(1+idim)*_nVar + _Ndim+1]=((1-alpha)*rho1*g2press+dpi2*alpha*invcsq1*g2alpha)*_vec_normal[idim] - u2_n*_Uroe[1+idim+_Ndim+1];
386 /******* Construction des matrices de decentrement *****/
387 if(_spaceScheme == centered){
388 if(_entropicCorrection)
389 throw CdmathException("IsothermalTwoFluid::convectionMatrices: entropic scheme not available for centered scheme");
390 for(int i=0; i<_nVar*_nVar;i++)
393 if( _spaceScheme ==staggered){
394 if(_entropicCorrection)//To do: study entropic correction for staggered
395 throw CdmathException("IsothermalTwoFluid::convectionMatrices: entropic scheme not yet available for staggered scheme");
396 /******** Construction du decentrement du type decale *********/
398 for(int i=0; i<_nVar*_nVar;i++)
401 for(int idim=0; idim<_Ndim;idim++)
403 _absAroe[1+idim]=_vec_normal[idim];
404 _absAroe[1+idim+_Ndim+1]=0;
405 _absAroe[(_Ndim+1)*_nVar+1+idim]=0;
406 _absAroe[(_Ndim+1)*_nVar+1+idim+_Ndim+1]=_vec_normal[idim];
409 for(int idim=0; idim<_Ndim;idim++)
411 //premiere colonne (masse gaz)
412 _absAroe[ (1+idim)*_nVar]= (-alpha *rho2*g2press+dpi1*(1-alpha)*invcsq2*g2alpha)*_vec_normal[idim] - u1_n*_Uroe[2+idim];
413 _absAroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar]=(-(1-alpha)*rho2*g2press-dpi2*(1-alpha)*invcsq2*g2alpha)*_vec_normal[idim];
414 //colonnes intermediaires
415 for(int jdim=0; jdim<_Ndim;jdim++)
417 _absAroe[ (1+idim)*_nVar + jdim + 1] = _Uroe[ 2+idim]*_vec_normal[jdim];
418 _absAroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar + jdim + 1+_Ndim+1] = _Uroe[_Ndim+2+idim]*_vec_normal[jdim];
421 _absAroe[ (1+idim)*_nVar + idim + 1] += u1_n;
422 _absAroe[(_Ndim+1)*_nVar+ (1+idim)*_nVar + _Ndim+1+ idim + 1] += u2_n;
423 //troisieme colonne (masse liquide)
424 _absAroe[ (1+idim)*_nVar + _Ndim+1]= (-alpha *rho1*g2press-dpi1*alpha*invcsq1*g2alpha)*_vec_normal[idim];
425 _absAroe[(_Ndim+1)*_nVar+(1+idim)*_nVar + _Ndim+1]=(-(1-alpha)*rho1*g2press+dpi2*alpha*invcsq1*g2alpha)*_vec_normal[idim] - u2_n*_Uroe[1+idim+_Ndim+1];
428 double signu1,signu2;
442 for(int i=0; i<(1+_Ndim)*_nVar;i++){
443 _absAroe[i] *= signu1;
444 _absAroe[i+(1+_Ndim)*_nVar] *= signu2;
447 if(_spaceScheme==upwind || _spaceScheme ==lowMach)
449 vector< complex< double > > y (taille_vp,0);
451 for( int i=0 ; i<taille_vp ; i++)
452 y[i] = Poly.abs_generalise(valeurs_propres[i]);
454 if(_entropicCorrection)
456 entropicShift(_vec_normal);
457 y[0] +=_entropicShift[0];
458 y[taille_vp-1] +=_entropicShift[2];
459 for( int i=1 ; i<taille_vp-1 ; i++)
460 y[i] +=_entropicShift[1];
463 Poly.abs_par_interp_directe(taille_vp,valeurs_propres, _Aroe, _nVar,_precision, _absAroe,y);
464 if( _spaceScheme ==pressureCorrection){
465 for( int i=0 ; i<_Ndim ; i++)
466 for( int j=0 ; j<_Ndim ; j++){
467 _absAroe[(1+i)*_nVar+1+j]-=alpha*(valeurs_propres[1].real()-valeurs_propres[0].real())/2*_vec_normal[i]*_vec_normal[j];
468 _absAroe[(2+_Ndim+i)*_nVar+2+_Ndim+j]-=(1-alpha)*(valeurs_propres[1].real()-valeurs_propres[0].real())/2*_vec_normal[i]*_vec_normal[j];
471 else if( _spaceScheme ==lowMach){
472 double M=max(fabs(u1_2),fabs(u2_2))/_maxvploc;
473 for( int i=0 ; i<_Ndim ; i++)
474 for( int j=0 ; j<_Ndim ; j++){
475 _absAroe[(1+i)*_nVar+1+j]-=(1-M)*alpha*(valeurs_propres[1].real()-valeurs_propres[0].real())/2*_vec_normal[i]*_vec_normal[j];
476 _absAroe[(2+_Ndim+i)*_nVar+2+_Ndim+j]-=(1-M)*(1-alpha)*(valeurs_propres[1].real()-valeurs_propres[0].real())/2*_vec_normal[i]*_vec_normal[j];
481 //Calcul de la matrice signe pour VFFC, VFRoe et décentrement des termes source
482 vector< complex< double > > valeurs_propres_dist(taille_vp,0);
483 for( int i=0 ; i<taille_vp ; i++)
484 valeurs_propres_dist[i] = valeurs_propres[i];
486 if(_entropicCorrection || _spaceScheme ==pressureCorrection || _spaceScheme ==lowMach){
487 InvMatriceRoe( valeurs_propres_dist);
488 Poly.matrixProduct(_absAroe, _nVar, _nVar, _invAroe, _nVar, _nVar, _signAroe);
490 else if (_spaceScheme==upwind)//upwind sans entropic
491 SigneMatriceRoe( valeurs_propres_dist);
492 else if (_spaceScheme==centered)//centre sans entropic
494 for(int i=0; i<_nVar*_nVar;i++)
497 else if(_spaceScheme ==staggered )
499 double signu1,signu2;
512 for(int i=0; i<_nVar*_nVar;i++)
514 _signAroe[0] = signu1;
515 _signAroe[(1+_Ndim)*_nVar +1+_Ndim] = signu2;
516 for(int i=2; i<_nVar-1;i++){
517 _signAroe[i*_nVar+i] = -signu1;
518 _signAroe[(i+1+_Ndim)*_nVar+i+1+_Ndim] = -signu2;
522 throw CdmathException("IsothermalTwoFluid::convectionMatrices: well balanced option not treated");
524 for(int i=0; i<_nVar*_nVar;i++)
526 _AroeMinus[i] = (_Aroe[i]-_absAroe[i])/2;
527 _AroePlus[i] = (_Aroe[i]+_absAroe[i])/2;
529 if(_timeScheme==Implicit && _usePrimitiveVarsInNewton)//Implicitation using primitive variables
530 for(int i=0; i<_nVar*_nVar;i++)
531 _AroeMinusImplicit[i] = (_AroeImplicit[i]-_absAroeImplicit[i])/2;
533 for(int i=0; i<_nVar*_nVar;i++)
534 _AroeMinusImplicit[i] = _AroeMinus[i];
536 if(_verbose && _nbTimeStep%_freqSave ==0)
538 cout<<endl<<"Matrice de Roe"<<endl;
539 for(int i=0; i<_nVar;i++)
541 for(int j=0; j<_nVar;j++)
542 cout << _Aroe[i*_nVar+j]<< " , ";
545 cout<<"Valeur absolue matrice de Roe"<<endl;
546 for(int i=0; i<_nVar;i++){
547 for(int j=0; j<_nVar;j++)
548 cout<<_absAroe[i*_nVar+j]<<" , ";
551 cout<<endl<<"Matrice _AroeMinus"<<endl;
552 for(int i=0; i<_nVar;i++)
554 for(int j=0; j<_nVar;j++)
555 cout << _AroeMinus[i*_nVar+j]<< " , ";
558 cout<<endl<<"Matrice _AroePlus"<<endl;
559 for(int i=0; i<_nVar;i++)
561 for(int j=0; j<_nVar;j++)
562 cout << _AroePlus[i*_nVar+j]<< " , ";
568 void IsothermalTwoFluid::setBoundaryState(string nameOfGroup, const int &j,double *normale){
569 //To do controle signe des vitesses pour CL entree/sortie
572 for(k=1; k<_nVar; k++)
573 _idm[k] = _idm[k-1] + 1;
575 VecGetValues(_conservativeVars, _nVar, _idm, _externalStates);//On initialise l'état fantôme avec l'état interne
576 double q1_n=0, q2_n=0;//quantité de mouvement normale à la face limite
577 for(k=0; k<_Ndim; k++){
578 q1_n+=_externalStates[(k+1)]*normale[k];
579 q2_n+=_externalStates[(k+1+1+_Ndim)]*normale[k];
582 if(_verbose && _nbTimeStep%_freqSave ==0)
584 cout << "Boundary conditions for group "<< nameOfGroup<< ", inner cell j= "<<j << " face unit normal vector "<<endl;
585 for(k=0; k<_Ndim; k++){
586 cout<<normale[k]<<", ";
590 if (_limitField[nameOfGroup].bcType==Wall){
592 //Pour la convection, inversion du sens des vitesses
593 for(k=0; k<_Ndim; k++){
594 _externalStates[(k+1)]-= 2*q1_n*normale[k];
595 _externalStates[(k+1+1+_Ndim)]-= 2*q2_n*normale[k];
599 for(k=1; k<_nVar; k++)
600 _idm[k] = _idm[k-1] + 1;
602 VecAssemblyBegin(_Uext);
603 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
604 VecAssemblyEnd(_Uext);
606 //Pour la diffusion, paroi à vitesses imposees
607 _externalStates[1]=_externalStates[0]*_limitField[nameOfGroup].v_x[0];
608 _externalStates[2+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_x[1];
611 _externalStates[2]=_externalStates[0]*_limitField[nameOfGroup].v_y[0];
612 _externalStates[3+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_y[1];
615 _externalStates[3]=_externalStates[0]*_limitField[nameOfGroup].v_z[0];
616 _externalStates[4+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_z[1];
619 VecAssemblyBegin(_Uextdiff);
620 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
621 VecAssemblyEnd(_Uextdiff);
623 else if (_limitField[nameOfGroup].bcType==Neumann){
625 for(k=1; k<_nVar; k++)
626 _idm[k] = _idm[k-1] + 1;
628 VecAssemblyBegin(_Uext);
629 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
630 VecAssemblyEnd(_Uext);
632 VecAssemblyBegin(_Uextdiff);
633 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
634 VecAssemblyEnd(_Uextdiff);
636 else if (_limitField[nameOfGroup].bcType==Inlet){
638 for(k=1; k<_nVar; k++)
639 _idm[k] = _idm[k-1] + 1;
641 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
642 double alpha=_limitField[nameOfGroup].alpha;
643 double pression=_Vj[1];
644 double T=_Temperature;
645 double rho_v=_fluides[0]->getDensity(pression,T);
646 double rho_l=_fluides[1]->getDensity(pression,T);
647 _externalStates[0]=alpha*rho_v;
648 _externalStates[1]=alpha*rho_v*(_limitField[nameOfGroup].v_x[0]);
649 _externalStates[1+_Ndim]=(1-alpha)*rho_l;
650 _externalStates[2+_Ndim]=(1-alpha)*rho_l*(_limitField[nameOfGroup].v_x[1]);
653 _externalStates[2]=_externalStates[0]*_limitField[nameOfGroup].v_y[0];
654 _externalStates[3+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_y[1];
657 _externalStates[3]=_externalStates[0]*_limitField[nameOfGroup].v_z[0];
658 _externalStates[4+_Ndim]=_externalStates[1+_Ndim]*_limitField[nameOfGroup].v_z[1];
662 for(k=1; k<_nVar; k++)
663 _idm[k] = _idm[k-1] + 1;
664 VecAssemblyBegin(_Uext);
665 VecAssemblyBegin(_Uextdiff);
666 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
667 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
668 VecAssemblyEnd(_Uext);
669 VecAssemblyEnd(_Uextdiff);
671 else if (_limitField[nameOfGroup].bcType==InletPressure){
672 //Computation of the hydrostatic contribution : scalar product between gravity vector and position vector
673 Cell Cj=_mesh.getCell(j);
674 double hydroPress=Cj.x()*_GravityField3d[0];
676 hydroPress+=Cj.y()*_GravityField3d[1];
678 hydroPress+=Cj.z()*_GravityField3d[2];
680 hydroPress*=_externalStates[0]+_externalStates[_Ndim];//multiplication by rho the total density
682 //Building the external state
684 for(k=1; k<_nVar; k++)
685 _idm[k] = _idm[k-1] + 1;
687 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
688 double pression=_limitField[nameOfGroup].p+hydroPress;
689 double alpha=_limitField[nameOfGroup].alpha;
690 double T=_Temperature;
691 double rho_v=_fluides[0]->getDensity(pression,T);
692 double rho_l=_fluides[1]->getDensity(pression,T);
693 _externalStates[0]=alpha*rho_v;
694 _externalStates[1+_Ndim]=(1-alpha)*rho_l;
696 for(k=1;k<1+_Ndim;k++){
697 _externalStates[k]=_externalStates[0]*_Vj[1+k];
698 _externalStates[k+1+_Ndim]=_externalStates[1+_Ndim]*_Vj[k+1+_Ndim];
701 for(k=1; k<_nVar; k++)
702 _idm[k] = _idm[k-1] + 1;
703 VecAssemblyBegin(_Uext);
704 VecAssemblyBegin(_Uextdiff);
705 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
706 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
707 VecAssemblyEnd(_Uext);
708 VecAssemblyEnd(_Uextdiff);
710 if(_verbose && _nbTimeStep%_freqSave ==0)
712 cout<<"Etat fantôme InletPressure"<<endl;
713 for(int k=0;k<_nVar;k++)
714 cout<<_externalStates[k]<<endl;
717 else if (_limitField[nameOfGroup].bcType==Outlet){
719 for(k=1; k<_nVar; k++)
720 _idm[k] = _idm[k-1] + 1;
722 //Computation of the hydrostatic contribution : scalar product between gravity vector and position vector
723 Cell Cj=_mesh.getCell(j);
724 double hydroPress=Cj.x()*_GravityField3d[0];
726 hydroPress+=Cj.y()*_GravityField3d[1];
728 hydroPress+=Cj.z()*_GravityField3d[2];
730 hydroPress*=_externalStates[0]+_externalStates[_Ndim];//multiplication by rho the total density
732 //Building the external state
733 VecGetValues(_primitiveVars, _nVar, _idm,_Vj);
734 double pression_int=_Vj[1];
735 double pression_ext=_limitField[nameOfGroup].p+hydroPress;
736 double T=_Vj[_nVar-1];
737 double rho_v_int=_fluides[0]->getDensity(pression_int,T);
738 double rho_l_int=_fluides[1]->getDensity(pression_int,T);
739 double rho_v_ext=_fluides[0]->getDensity(pression_ext,T);
740 double rho_l_ext=_fluides[1]->getDensity(pression_ext,T);
742 for(k=0;k<1+_Ndim;k++){
743 _externalStates[k]*=rho_v_ext/rho_v_int;
744 _externalStates[k+1+_Ndim]*=rho_l_ext/rho_l_int;
747 for(k=1; k<_nVar; k++)
748 _idm[k] = _idm[k-1] + 1;
749 VecAssemblyBegin(_Uext);
750 VecAssemblyBegin(_Uextdiff);
751 VecSetValues(_Uext, _nVar, _idm, _externalStates, INSERT_VALUES);
752 VecSetValues(_Uextdiff, _nVar, _idm, _externalStates, INSERT_VALUES);
753 VecAssemblyEnd(_Uext);
754 VecAssemblyEnd(_Uextdiff);
757 cout<<"!!!!!!!!!!!!!!!!! Error IsothermalTwoFluid::setBoundaryState !!!!!!!!!!"<<endl;
758 cout<<"!!!!!!!!!!!!! Boundary condition not set for boundary named"<<nameOfGroup<< ", _limitField[nameOfGroup].bcType= "<<_limitField[nameOfGroup].bcType<<endl;
759 cout<<"Accepted boundary condition are Neumann, Wall, Inlet, and Outlet"<<endl;
760 throw CdmathException("Unknown boundary condition");
764 void IsothermalTwoFluid::addDiffusionToSecondMember
770 double mu1 = _fluides[0]->getViscosity(_Temperature);
771 double mu2 = _fluides[1]->getViscosity(_Temperature);
776 //extraction des valeurs
778 for(int k=1; k<_nVar; k++)
779 _idm[k] = _idm[k-1] + 1;
781 VecGetValues(_primitiveVars, _nVar, _idm, _Vi);
782 if (_verbose && _nbTimeStep%_freqSave ==0)
784 cout << "Contribution diffusion: variables primitives maille " << i<<endl;
785 for(int q=0; q<_nVar; q++)
787 cout << _Vi[q] << endl;
793 for(int k=0; k<_nVar; k++)
794 _idm[k] = _nVar*j + k;
795 VecGetValues(_primitiveVars, _nVar, _idm, _Vj);
799 for(int k=0; k<_nVar; k++)
802 VecGetValues(_Uextdiff, _nVar, _idm, _phi);
803 consToPrim(_phi,_Vj);
806 if (_verbose && _nbTimeStep%_freqSave ==0)
808 cout << "Contribution diffusion: variables primitives maille " <<j <<endl;
809 for(int q=0; q<_nVar; q++)
811 cout << _Vj[q] << endl;
817 double alpha=(_Vj[0]+_Vi[0])/2;
818 //on n'a pas de contribution sur la masse
821 //contribution visqueuse sur la quantite de mouvement
822 for(int k=1; k<_Ndim+1; k++)
824 _phi[k] = _inv_dxi*2/(1/_inv_dxi+1/_inv_dxj)*alpha*mu1*(_Vj[k] - _Vi[k]);//attention car primitif=alpha p u1 u2
825 _phi[k+_Ndim+1] = _inv_dxi*2/(1/_inv_dxi+1/_inv_dxj)*(1-alpha)*mu2*(_Vj[1+k+_Ndim] - _Vi[1+k+_Ndim]);
829 VecSetValuesBlocked(_b, 1, _idm, _phi, ADD_VALUES);
831 if(_verbose && _nbTimeStep%_freqSave ==0)
833 cout << "Contribution diffusion au 2nd membre pour la maille " << i << ": "<<endl;
834 for(int q=0; q<_nVar; q++)
836 cout << _phi[q] << endl;
843 //On change de signe pour l'autre contribution
844 for(int k=0; k<_nVar; k++)
845 _phi[k] *= -_inv_dxj/_inv_dxi;
848 VecSetValuesBlocked(_b, 1, _idm, _phi, ADD_VALUES);
851 if(_verbose && _nbTimeStep%_freqSave ==0)
853 cout << "Contribution diffusion au 2nd membre pour la maille " << j << ": "<<endl;
854 for(int q=0; q<_nVar; q++)
856 cout << _phi[q] << endl;
860 if(_timeScheme==Implicit)
862 cout << "Matrice de diffusion D, pour le couple (" << i << "," << j<< "):" << endl;
863 for(int i=0; i<_nVar; i++)
865 for(int j=0; j<_nVar; j++)
866 cout << _Diffusion[i*_nVar+j]<<", ";
874 void IsothermalTwoFluid::sourceVector(PetscScalar * Si,PetscScalar * Ui,PetscScalar * Vi, int i)
876 double m1=Ui[0],m2=Ui[1+_Ndim],alpha=Vi[0], P=Vi[1];
877 double norm_ur=0, Gamma;
878 for(int k=0; k<_Ndim; k++)
879 norm_ur+=(Vi[2+k]-Vi[2+k+_Ndim])*(Vi[2+k]-Vi[2+k+_Ndim]);
880 norm_ur=sqrt(norm_ur);
882 if(i>=0 &&_Temperature>_Tsat && alpha<1-_precision)
883 Gamma=_heatPowerField(i)/_latentHeat;
886 for(int k=1; k<_Ndim+1; k++)
888 //cout<<"Vi[1+"<<k+_Ndim<<"]="<<Vi[1+k+_Ndim]<<endl;
889 Si[k] =_gravite[k]*m1 -_dragCoeffs[0]*norm_ur*(Vi[1+k]-Vi[1+k+_Ndim])+ Gamma*Vi[1+k+_Ndim];//interfacial velocity= ul
890 Si[k+_Ndim+1] =_gravite[k+_Ndim+1]*m2 + _dragCoeffs[0]*norm_ur*(Vi[1+k]-Vi[1+k+_Ndim])- Gamma*Vi[1+k+_Ndim];
892 if(true){//heated boiling
896 else if (P<_Psat && alpha<1-_precision){//flash boiling
897 Si[0]=-_dHsatl_over_dp*_dp_over_dt(i)/_latentHeat;
898 Si[1+_Ndim]=_dHsatl_over_dp*_dp_over_dt(i)/_latentHeat;
905 if(_timeScheme==Implicit)
907 for(int i=0; i<_nVar*_nVar;i++)
908 _GravityImplicitationMatrix[i] = 0;
909 for(int i=0; i<_nVar/2;i++)
910 _GravityImplicitationMatrix[i*_nVar]=-_gravite[i];
911 for(int i=_nVar/2; i<_nVar;i++)
912 _GravityImplicitationMatrix[i*_nVar+_nVar/2]=-_gravite[i];
915 if(_verbose && _nbTimeStep%_freqSave ==0)
917 cout<<"IsothermalTwoFluid::sourceVector"<<endl;
919 for(int k=0;k<_nVar;k++)
923 for(int k=0;k<_nVar;k++)
927 for(int k=0;k<_nVar;k++)
930 if(_timeScheme==Implicit)
931 displayMatrix(_GravityImplicitationMatrix, _nVar, "Gravity implicitation matrix");
934 void IsothermalTwoFluid::pressureLossVector(PetscScalar * pressureLoss, double K, PetscScalar * Ui, PetscScalar * Vi, PetscScalar * Uj, PetscScalar * Vj)
936 double norm_u1=0, u1_n=0, norm_u2=0, u2_n=0, m1, m2;
937 for(int i=0;i<_Ndim;i++){
938 u1_n += _Uroe[2+i] *_vec_normal[i];
939 u2_n += _Uroe[2+i+_Ndim]*_vec_normal[i];
942 pressureLoss[1+_Ndim]=0;
944 for(int i=0;i<_Ndim;i++)
945 norm_u1 += Vi[1+i]*Vi[1+i];
946 norm_u1=sqrt(norm_u1);
948 for(int i=0;i<_Ndim;i++)
949 pressureLoss[1+i]=-K*m1*norm_u1*Vi[1+i];
952 for(int i=0;i<_Ndim;i++)
953 norm_u1 += Vj[1+i]*Vj[1+i];
954 norm_u1=sqrt(norm_u1);
956 for(int i=0;i<_Ndim;i++)
957 pressureLoss[1+i]=-K*m1*norm_u1*Vj[1+i];
960 for(int i=0;i<_Ndim;i++)
961 norm_u2 += Vi[2+i+_Ndim]*Vi[2+i+_Ndim];
962 norm_u2=sqrt(norm_u2);
964 for(int i=0;i<_Ndim;i++)
965 pressureLoss[2+i+_Ndim]=-K*m2*norm_u2*Vi[2+i+_Ndim];
968 for(int i=0;i<_Ndim;i++)
969 norm_u2 += Vj[2+i+_Ndim]*Vj[2+i+_Ndim];
970 norm_u2=sqrt(norm_u2);
972 for(int i=0;i<_Ndim;i++)
973 pressureLoss[2+i+_Ndim]=-K*m2*norm_u2*Vj[2+i+_Ndim];
975 if(_verbose && _nbTimeStep%_freqSave ==0)
977 cout<<"IsothermalTwoFluid::pressureLossVector K= "<<K<<endl;
979 for(int k=0;k<_nVar;k++)
983 for(int k=0;k<_nVar;k++)
987 for(int k=0;k<_nVar;k++)
991 for(int k=0;k<_nVar;k++)
994 cout<<"pressureLoss="<<endl;
995 for(int k=0;k<_nVar;k++)
996 cout<<pressureLoss[k]<<", ";
1001 void IsothermalTwoFluid::porosityGradientSourceVector()
1003 double u1_ni=0, u1_nj=0, u2_ni=0, u2_nj=0, rho1i, rho2i, rho1j, rho2j, pi=_Vi[1], pj=_Vj[1], pij1, pij2, alphaij=_Uroe[0];
1004 for(int i=0;i<_Ndim;i++) {
1005 u1_ni += _Vi[2+i]*_vec_normal[i];
1006 u2_ni += _Vi[2+_Ndim+i]*_vec_normal[i];
1007 u1_nj += _Vj[2+i]*_vec_normal[i];
1008 u2_nj += _Vj[2+_Ndim+i]*_vec_normal[i];
1010 _porosityGradientSourceVector[0]=0;
1011 _porosityGradientSourceVector[1+_Ndim]=0;
1012 rho1i = _fluides[0]->getDensity(pi, _Temperature);
1013 rho2i = _fluides[1]->getDensity(pi, _Temperature);
1014 rho1j = _fluides[0]->getDensity(pj, _Temperature);
1015 rho2j = _fluides[1]->getDensity(pj, _Temperature);
1016 pij1=(pi+pj)/2+rho1i*rho1j/2/(rho1i+rho1j)*(u1_ni-u1_nj)*(u1_ni-u1_nj);
1017 pij2=(pi+pj)/2+rho2i*rho2j/2/(rho2i+rho2j)*(u2_ni-u2_nj)*(u2_ni-u2_nj);
1018 for(int i=0;i<_Ndim;i++){
1019 _porosityGradientSourceVector[1+i] =alphaij*pij1*(_porosityi-_porosityj)*2/(1/_inv_dxi+1/_inv_dxj);
1020 _porosityGradientSourceVector[2+_Ndim+i]=alphaij*pij2*(_porosityi-_porosityj)*2/(1/_inv_dxi+1/_inv_dxj);
1024 /* Funtion of equations of states */
1025 double IsothermalTwoFluid::ecartPression(double m1,double m2, double alpha, double e1, double e2){
1026 if(alpha>_precision*_precision&& alpha<1-_precision*_precision)
1027 return _fluides[0]->getPressure((m1/alpha)*e1,m1/alpha) - _fluides[1]->getPressure((m2/(1-alpha))*e2,m2/(1-alpha));
1028 else if(alpha<=_precision*_precision)
1030 //cout<<"Warning ecartPression, alpha close to 0"<<endl;
1031 // if(m1<_precision*_precision*_precision)
1034 return alpha/_precision*_precision;
1038 //cout<<"Warning ecartPression, alpha close to 1"<<endl;
1039 // if(m2<_precision*_precision*_precision)
1042 return -(1-alpha)/_precision*_precision;
1046 double IsothermalTwoFluid::ecartPressionDerivee(double m1,double m2, double alpha, double e1, double e2){
1047 if(alpha>_precision*_precision && alpha<1-_precision*_precision )
1048 return -(m1/alpha)/alpha*e1*_fluides[0]->getPressureDerivativeRhoE() - (m2/(1-alpha))/(1-alpha)*e2*_fluides[1]->getPressureDerivativeRhoE();
1049 else if (alpha<_precision*_precision)
1050 return -1/_precision*_precision;
1052 return 1/_precision*_precision;
1055 double IsothermalTwoFluid::intPressDef(double alpha, double u_r2, double rho1, double rho2)
1057 return _intPressCoeff*alpha*(1-alpha)*rho1*rho2*u_r2/( alpha*rho2+(1-alpha)*rho1)
1058 +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
1059 *(alpha*alpha*rho2-(1-alpha)*(1-alpha)*rho1)
1060 *(alpha*alpha*rho2*rho2/(_fluides[0]->vitesseSonTemperature(_Temperature,rho1)*_fluides[0]->vitesseSonTemperature(_Temperature,rho1))
1061 -(1-alpha)*(1-alpha)*rho1*rho1/(_fluides[1]->vitesseSonTemperature(_Temperature,rho2)*_fluides[1]->vitesseSonTemperature(_Temperature,rho2)));
1064 void IsothermalTwoFluid::entropicShift(double* n)
1066 vector< double > pol_car;
1069 double u1_n = 0, u1_2 = 0, u2_n = 0, u2_2 = 0, u_r2 = 0;
1070 for(int i=0;i<_Ndim;i++)
1072 u1_2 += _l[2+i]*_l[2+i];
1073 u1_n += _l[2+i]*n[i];
1074 u2_2 += _l[1+i+1+_Ndim]*_l[1+i+1+_Ndim];
1075 u2_n += _l[1+i+1+_Ndim]*n[i];
1076 u_r2 += (_l[2+i]-_l[1+i+1+_Ndim])*(_l[2+i]-_l[1+i+1+_Ndim]);
1078 double alpha = _l[0];
1080 double rho1 = _fluides[0]->getDensity(p, _Temperature);
1081 double rho2 = _fluides[1]->getDensity(p, _Temperature);
1082 double dpi1 = intPressDef(alpha, u_r2, rho1, rho2);
1084 double invcsq1 = 1/_fluides[0]->vitesseSonTemperature(_Temperature,rho1);
1086 double invcsq2 = 1/_fluides[1]->vitesseSonTemperature(_Temperature,rho2);
1090 pol_car= Poly.polynome_caracteristique(alpha, 1-alpha, u1_n, u2_n, rho1, rho2, invcsq1, invcsq2, dpi1, dpi2);
1091 for(int ct=0;ct<4;ct++){
1092 if (abs(pol_car[ct])<_precision*_precision)
1095 vector< complex<double> > vp_left = getRacines(pol_car);
1096 int taille_vp_left = Poly.new_tri_selectif(vp_left,vp_left.size(),_precision);
1098 if(_verbose && _nbTimeStep%_freqSave ==0)
1100 cout<<"Entropic shift left eigenvalues: "<<endl;
1101 for(unsigned int ct =0; ct<vp_left.size(); ct++)
1102 cout<<"("<< vp_left[ct].real()<< ", " <<vp_left[ct].imag() << ")";
1110 for(int i=0;i<_Ndim;i++)
1112 u1_2 += _r[2+i]*_r[2+i];
1113 u1_n += _r[2+i]*n[i];
1114 u2_2 += _r[1+i+1+_Ndim]*_r[1+i+1+_Ndim];
1115 u2_n += _r[1+i+1+_Ndim]*n[i];
1116 u_r2 += (_r[2+i]-_r[1+i+1+_Ndim])*(_r[2+i]-_r[1+i+1+_Ndim]);
1120 rho1 = _fluides[0]->getDensity(p, _Temperature);
1121 rho2 = _fluides[1]->getDensity(p, _Temperature);
1122 dpi1 = intPressDef(alpha, u_r2, rho1, rho2);
1124 invcsq1 = 1/_fluides[0]->vitesseSonTemperature(_Temperature,rho1);
1126 invcsq2 = 1/_fluides[1]->vitesseSonTemperature(_Temperature,rho2);
1129 pol_car= Poly.polynome_caracteristique(alpha, 1-alpha, u1_n, u2_n, rho1, rho2, invcsq1, invcsq2, dpi1, dpi2);
1130 for(int ct=0;ct<4;ct++){
1131 if (abs(pol_car[ct])<_precision*_precision)
1134 vector< complex<double> > vp_right = getRacines(pol_car);
1135 int taille_vp_right = Poly.new_tri_selectif(vp_right,vp_right.size(),_precision);
1137 if(_verbose && _nbTimeStep%_freqSave ==0)
1139 cout<<"Entropic shift right eigenvalues: "<<endl;
1140 for(unsigned int ct =0; ct<vp_right.size(); ct++)
1141 cout<<"("<<vp_right[ct].real()<< ", " <<vp_right[ct].imag() <<")";
1144 _entropicShift[0] = abs(vp_left[0]-vp_right[0]);
1145 _entropicShift[2] = abs(vp_left[taille_vp_left-1]-vp_right[taille_vp_right-1]);
1146 _entropicShift[1]=0;
1147 for(int i=1;i<min(taille_vp_right-1,taille_vp_left-1);i++)
1148 _entropicShift[1] = max(_entropicShift[1],abs(vp_left[i]-vp_right[i]));
1149 if(_verbose && _nbTimeStep%_freqSave ==0)
1151 cout<<"eigenvalue jumps "<<endl;
1152 cout<< _entropicShift[0] << ", " << _entropicShift[1] << ", "<< _entropicShift[2] <<endl;
1156 void IsothermalTwoFluid::computeScaling(double maxvp)
1158 // double alphaScaling;
1159 // if(_guessalpha>_precision && _guessalpha<1-_precision)
1160 // alphaScaling=_guessalpha;
1162 // alphaScaling=0.5;
1164 _blockDiag[0]=1;//alphaScaling;
1165 _invBlockDiag[0]=1;//_blockDiag[0];
1166 _blockDiag[1+_Ndim]=1;//-alphaScaling;
1167 _invBlockDiag[1+_Ndim]=1.0;//_blockDiag[1+_Ndim];
1168 for(int q=1; q<_Ndim+1; q++)
1170 _blockDiag[q]=1/(maxvp*maxvp);
1171 _invBlockDiag[q]=1;//_blockDiag[q];
1172 _blockDiag[q+1+_Ndim]=1/(maxvp*maxvp);
1173 _invBlockDiag[q+1+_Ndim]=1;//_blockDiag[q+1+_Ndim];
1177 void IsothermalTwoFluid::jacobian(const int &j, string nameOfGroup,double * normale)
1180 for(k=0; k<_nVar*_nVar;k++)
1181 _Jcb[k] = 0;//No implicitation at this stage
1183 // loop on boundaries
1184 if (_limitField[nameOfGroup].bcType==Wall)
1186 for(k=0; k<_nVar;k++)
1187 _Jcb[k*_nVar + k] = 1;
1188 for(k=1; k<1+_Ndim;k++)
1189 for(int l=1; l<1+_Ndim;l++){
1190 _Jcb[k*_nVar + l] -= 2*normale[k-1]*normale[l-1];
1191 _Jcb[(k+1+_Ndim)*_nVar + l+1+_Ndim] -= 2*normale[k-1]*normale[l-1];
1195 else if (_limitField[nameOfGroup].bcType==Inlet)
1198 _Jcb[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1199 _Jcb[_nVar]=_limitField[nameOfGroup].v_x[0];
1200 _Jcb[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1203 _Jcb[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1204 _Jcb[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1207 _Jcb[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1208 _Jcb[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1213 // not wall, not inlet
1214 else if (_limitField[nameOfGroup].bcType==Outlet){
1216 for(k=1; k<_nVar;k++)
1217 {_idm[k] = _idm[k-1] + 1;}
1218 VecGetValues(_conservativeVars, _nVar, _idm, _phi);
1219 VecGetValues(_primitiveVars, _nVar, _idm, _externalStates);
1221 else if (_limitField[nameOfGroup].bcType!=Neumann && _limitField[nameOfGroup].bcType!=InletPressure)// not wall, not inlet, not outlet
1223 cout << "group named "<<nameOfGroup << " : unknown boundary condition" << endl;
1224 throw CdmathException("IsothermalTwoFluid::jacobianDiff: This boundary condition is not treated");
1228 void IsothermalTwoFluid::jacobianDiff(const int &j, string nameOfGroup)
1232 for(k=0; k<_nVar*_nVar;k++)
1234 if (_limitField[nameOfGroup].bcType==Wall){
1236 _JcbDiff[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1237 _JcbDiff[_nVar]=_limitField[nameOfGroup].v_x[0];
1238 _JcbDiff[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1241 _JcbDiff[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1242 _JcbDiff[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1246 _JcbDiff[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1247 _JcbDiff[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1250 } else if (_limitField[nameOfGroup].bcType==Inlet){
1253 _JcbDiff[(1+_Ndim)*_nVar +1+_Ndim] = 1;
1254 _JcbDiff[_nVar]=_limitField[nameOfGroup].v_x[0];
1255 _JcbDiff[(2+_Ndim)*_nVar +1+_Ndim] =_limitField[nameOfGroup].v_x[1];
1258 _JcbDiff[2*_nVar]=_limitField[nameOfGroup].v_y[0];
1259 _JcbDiff[(3+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_y[1];
1263 _JcbDiff[3*_nVar]=_limitField[nameOfGroup].v_z[0];
1264 _JcbDiff[(4+_Ndim)*_nVar +1+_Ndim]= _limitField[nameOfGroup].v_z[1];
1268 } else if (_limitField[nameOfGroup].bcType==Outlet){
1270 //extraction de l etat courant et primitives
1272 for(k=1; k<_nVar;k++)
1273 {_idm[k] = _idm[k-1] + 1;}
1274 VecGetValues(_conservativeVars, _nVar, _idm, _phi);
1275 VecGetValues(_primitiveVars, _nVar, _idm, _externalStates);
1278 else if (_limitField[nameOfGroup].bcType!=Neumann && _limitField[nameOfGroup].bcType!=InletPressure){
1279 cout<<"Condition limite non traitee pour le bord "<<nameOfGroup<< endl;
1280 throw CdmathException("IsothermalTwoFluid::jacobianDiff: Condition limite non traitee");
1284 void IsothermalTwoFluid::primToCons(const double *P, const int &i, double *W, const int &j){
1287 W[j*_nVar] = P[i*_nVar]*_fluides[0]->getDensity(P[i*_nVar+1],_Temperature);
1288 W[j*_nVar+1+_Ndim] = (1-P[i*_nVar])*_fluides[1]->getDensity(P[i*_nVar+1],_Temperature);
1289 for(int k=0; k<_Ndim; k++)
1291 W[j*_nVar+(k+1)] = W[j*_nVar]*P[i*_nVar+(k+2)];
1292 W[j*_nVar+(k+1)+1+_Ndim] = W[j*_nVar+1+_Ndim]*P[i*_nVar+(k+2)+_Ndim];
1297 void IsothermalTwoFluid::consToPrim(const double *Wcons, double* Wprim,double porosity)//To do: treat porosity
1302 double m2=Wcons[_Ndim+1];
1303 double e1 = _internalEnergy1;
1304 double e2 = _internalEnergy2;
1306 _minm1=min(m1,_minm1);
1307 _minm2=min(m2,_minm2);
1308 if(m1<-_precision || m2<-_precision)
1310 if(fabs(m1)<_precision*_precision)
1313 Wprim[1]=_fluides[1]->getPressure(m2*e2,m2);
1315 for(int idim=0; idim<_Ndim; idim++)
1317 Wprim[2+_Ndim+idim] = Wcons[2+_Ndim+idim]/Wcons[1+_Ndim];
1318 Wprim[2+idim] = Wprim[2+_Ndim+idim];
1321 else if(fabs(m2)<_precision*_precision)
1324 Wprim[1]=_fluides[0]->getPressure(m1*e1,m1);
1325 for(int idim=0; idim<_Ndim; idim++)
1327 Wprim[2+idim] = Wcons[1+idim]/Wcons[0];
1328 Wprim[2+_Ndim+idim] = Wprim[2+idim];
1333 for(int idim=0; idim<_Ndim; idim++)
1335 Wprim[2+idim] = Wcons[1+idim]/Wcons[0];
1336 Wprim[2+_Ndim+idim] = Wcons[2+_Ndim+idim]/Wcons[1+_Ndim];
1338 double alphanewton, alphainf=0, alphasup=1;
1339 int iterMax=50, iter=0;
1341 double dp=ecartPression( m1, m2, _guessalpha, e1, e2);
1342 double dpprim=ecartPressionDerivee( m1, m2, _guessalpha, e1, e2);
1344 while(fabs(dp)>1e3*_precision && iter<iterMax)
1348 // alphainf = _guessalpha;
1350 // alphansup = _guessalpha;
1352 // _guessalpha = (alphainf+alphansup)/2;
1355 alphainf = _guessalpha;
1357 alphasup = _guessalpha;
1359 alphanewton = _guessalpha-dp/dpprim;
1361 if(alphanewton<=alphainf)
1363 _guessalpha = (9*alphainf+alphasup)/10;
1364 //cout<< "dichotomie"<<endl;
1366 else if(alphanewton>=alphasup)
1368 _guessalpha = (alphainf+9*alphasup)/10;
1369 //cout<< "dichotomie"<<endl;
1372 _guessalpha=alphanewton;
1375 if(_verbose && _nbTimeStep%_freqSave ==0)
1376 cout<<"consToPrim diphasique iter= " <<iter<<" dp= " << dp<<" dpprim= " << dpprim<< " _guessalpha= " << _guessalpha<<endl;
1377 dp=ecartPression( m1, m2, _guessalpha, e1, e2);
1378 dpprim=ecartPressionDerivee( m1, m2, _guessalpha, e1, e2);
1382 if(_verbose && _nbTimeStep%_freqSave ==0)
1387 _err_press_max = max(_err_press_max,fabs(dp/1.e5));
1388 // if(_guessalpha>0.5)
1394 Wprim[1]=_fluides[0]->getPressure((m1/_guessalpha)*e1,m1/_guessalpha);
1396 Wprim[1]=_fluides[1]->getPressure((m2/(1-_guessalpha))*e2,m2/(1-_guessalpha));
1397 Wprim[0]=_guessalpha;
1400 cout << "pressure = "<< Wprim[1] << " < 0 " << endl;
1401 cout << "Conservative state = ";
1402 for(int k=0; k<_nVar; k++){
1403 cout<<Wcons[k]<<", ";
1406 *_runLogFile<< "IsothermalTwoFluid::consToPrim: negative pressure = "<< Wprim[1] << " < 0 " << endl;
1407 throw CdmathException("IsothermalTwoFluid::consToPrim: negative pressure");
1411 Vector IsothermalTwoFluid::convectionFlux(Vector U,Vector V, Vector normale, double porosity){
1413 cout<<"Ucons"<<endl;
1415 cout<<"Vprim"<<endl;
1419 double phim1=U(0);//phi alpha1 rho1
1420 double phim2=U(1+_Ndim);//phi alpha2 rho2
1421 Vector phiq1(_Ndim),phiq2(_Ndim);//phi alpha1 rho1 u1, phi alpha2 rho2 u2
1422 for(int i=0;i<_Ndim;i++){
1424 phiq2(i)=U(2+_Ndim+i);
1427 double pression=V(1);
1428 Vector vitesse1(_Ndim),vitesse2(_Ndim);
1429 for(int i=0;i<_Ndim;i++){
1431 vitesse2(i)=V(2+_Ndim+i);
1434 double vitesse1n=vitesse1*normale;
1435 double vitesse2n=vitesse2*normale;
1437 double alpha_roe = _Uroe[0];//Toumi formula
1438 // interfacial pressure term (hyperbolic correction)
1439 double dpi = _Uroe[_nVar];
1442 F(0)=phim1*vitesse1n;
1443 F(1+_Ndim)=phim2*vitesse2n;
1444 for(int i=0;i<_Ndim;i++){
1445 F(1+i)=phim1*vitesse1n*vitesse1(i)+(alpha_roe*pression+dpi*alpha)*porosity*normale(i);
1446 F(2+_Ndim+i)=phim2*vitesse2n*vitesse2(i)+((1-alpha_roe)*pression-dpi*alpha)*normale(i)*porosity;
1450 cout<<"Flux F(U,V)"<<endl;
1457 Vector IsothermalTwoFluid::staggeredVFFCFlux()
1459 if(_spaceScheme!=staggered || _nonLinearFormulation!=VFFC)
1460 throw CdmathException("IsothermalTwoFluid::staggeredVFFCFlux: staggeredVFFCFlux method should be called only for VFFC formulation and staggered upwinding");
1461 else//_spaceScheme==staggered
1464 double alpha_roe = _Uroe[0];//Toumi formula
1465 // interfacial pressure term (hyperbolic correction)
1466 double dpi = _Uroe[_nVar];
1468 double u1ijn=0, u2ijn=0, phialphaq1n=0, phialphaq2n=0;
1469 for(int idim=0;idim<_Ndim;idim++){//URoe = alpha, p, u1, u2, dpi
1470 u1ijn+=_vec_normal[idim]*_Uroe[2+idim];
1471 u2ijn+=_vec_normal[idim]*_Uroe[2+_Ndim+idim];
1475 for(int idim=0;idim<_Ndim;idim++)
1476 phialphaq1n+=_vec_normal[idim]*_Ui[1+idim];//phi alpha rho u n
1478 for(int idim=0;idim<_Ndim;idim++)
1479 Fij(1+idim)=phialphaq1n*_Vi[2+idim]+(alpha_roe*_Vj[1]*_porosityj+dpi*_Vi[0]*_porosityi)*_vec_normal[idim];
1483 for(int idim=0;idim<_Ndim;idim++)
1484 phialphaq2n+=_vec_normal[idim]*_Uj[1+idim];//phi alpha rho u n
1486 for(int idim=0;idim<_Ndim;idim++)
1487 Fij(1+idim)=phialphaq2n*_Vj[2+idim]+(alpha_roe*_Vi[1]*_porosityi+dpi*_Vj[0]*_porosityj)*_vec_normal[idim];
1492 for(int idim=0;idim<_Ndim;idim++)
1493 phialphaq2n+=_vec_normal[idim]*_Ui[2+_Ndim+idim];//phi alpha rho u n
1494 Fij(1+_Ndim)=phialphaq2n;
1495 for(int idim=0;idim<_Ndim;idim++)
1496 Fij(2+_Ndim+idim)=phialphaq2n*_Vi[2+_Ndim+idim]+((1-alpha_roe)*_Vj[1]*_porosityj-dpi*_Vi[0]*_porosityi)*_vec_normal[idim];
1500 for(int idim=0;idim<_Ndim;idim++)
1501 phialphaq2n+=_vec_normal[idim]*_Uj[2+_Ndim+idim];//phi alpha rho u n
1502 Fij(1+_Ndim)=phialphaq2n;
1503 for(int idim=0;idim<_Ndim;idim++)
1504 Fij(2+_Ndim+idim)=phialphaq2n*_Vj[2+_Ndim+idim]+((1-alpha_roe)*_Vi[1]*_porosityi-dpi*_Vj[0]*_porosityj)*_vec_normal[idim];
1510 void IsothermalTwoFluid::applyVFRoeLowMachCorrections(bool isBord, string groupname)
1512 if(_nonLinearFormulation!=VFRoe)
1513 throw CdmathException("IsothermalTwoFluid::applyVFRoeLowMachCorrections: applyVFRoeLowMachCorrections method should be called only for VFRoe formulation");
1514 else//_nonLinearFormulation==VFRoe
1516 if(_spaceScheme==lowMach){
1517 double u1_2=0, u2_2=0;
1518 for(int i=0;i<_Ndim;i++){
1519 u1_2 += _Uroe[2+i]*_Uroe[2+i];
1520 u2_2 += _Uroe[2+i+_Ndim]*_Uroe[2+i+_Ndim];
1523 double c = _maxvploc;//mixture sound speed
1524 double M=max(sqrt(u1_2),sqrt(u2_2))/c;//Mach number
1525 _Vij[1]=M*_Vij[1]+(1-M)*(_Vi[1]+_Vj[1])/2;
1526 primToCons(_Vij,0,_Uij,0);
1528 else if(_spaceScheme==pressureCorrection)
1530 if(_pressureCorrectionOrder>2)
1531 throw CdmathException("IsothermalTwoFluid::applyVFRoeLowMachCorrections pressure correction order can be only 1 or 2 for Isothermal two-fluid model");
1533 double norm_uij=0, uij_n=0, ui_n=0, uj_n=0;//mean velocities
1534 double rho1 = _fluides[0]->getDensity(_Uroe[1],_Temperature);
1535 double rho2 = _fluides[1]->getDensity(_Uroe[1],_Temperature);
1536 double m1=_Uroe[0]*rho1, m2=(1-_Uroe[0])*rho2;
1538 for(int i=0;i<_Ndim;i++)
1540 norm_uij += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*(m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim]);
1541 uij_n += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*_vec_normal[i];
1542 ui_n += _Vi[2+i]*_vec_normal[i];
1543 uj_n += _Vj[2+i]*_vec_normal[i];
1545 norm_uij=sqrt(norm_uij)/rhom;
1547 if(norm_uij>_precision)//avoid division by zero
1548 _Vij[1]=(_Vi[1]+_Vj[1])/2 + uij_n/norm_uij*(_Vj[1]-_Vi[1])/4 - rhom*norm_uij*(uj_n-ui_n)/4;
1550 _Vij[1]=(_Vi[1]+_Vj[1])/2 - rhom*norm_uij*(uj_n-ui_n)/4;
1553 else if(_spaceScheme==staggered)
1556 double rho1 = _fluides[0]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
1557 double rho2 = _fluides[1]->getDensity(_Uroe[1],_Uroe[_nVar-1]);
1558 double m1=_Uroe[0]*rho1, m2=(1-_Uroe[0])*rho2;
1560 for(int i=0;i<_Ndim;i++)
1561 qij_n += (m1*_Uroe[2+i]+m2*_Uroe[2+i+_Ndim])*_vec_normal[i];
1566 for(int i=0;i<_Ndim;i++)
1568 _Vij[2+i] =_Vi[2+i];
1569 _Vij[2+i+_Ndim] =_Vi[2+i+_Ndim];
1575 for(int i=0;i<_Ndim;i++)
1577 _Vij[2+i] =_Vj[2+i];
1578 _Vij[2+i+_Ndim] =_Vj[2+i+_Ndim];
1581 primToCons(_Vij,0,_Uij,0);
1586 void IsothermalTwoFluid::testConservation()
1588 double SUM, DELTA, x;
1590 for(int i=0; i<_nVar; i++)
1594 cout << "Masse totale phase " << 0 <<" (kg): ";
1595 else if( i == 1+_Ndim)
1596 cout << "Masse totale phase " << 1 <<" (kg): ";
1600 cout << "Quantite de mouvement totale phase 0 (kg.m.s^-1): ";
1602 cout << "Quantite de mouvement totale phase 1 (kg.m.s^-1): ";
1608 for(int j=0; j<_Nmailles; j++)
1610 VecGetValues(_old, 1, &I, &x);//on recupere la valeur du champ
1611 SUM += x*_mesh.getCell(j).getMeasure();
1612 VecGetValues(_newtonVariation, 1, &I, &x);//on recupere la variation du champ
1613 DELTA += x*_mesh.getCell(j).getMeasure();
1616 if(fabs(SUM)>_precision)
1617 cout << SUM<< ", variation relative: " << fabs(DELTA /SUM) << endl;
1619 cout << " a une somme nulle, variation absolue: " << fabs(DELTA) << endl;
1623 void IsothermalTwoFluid::save(){
1624 string prim(_path+"/IsothermalTwoFluidPrim_");
1625 string cons(_path+"/IsothermalTwoFluidCons_");
1630 for (long i = 0; i < _Nmailles; i++){
1631 /* j = 0 : void fraction
1633 j = 2, 3, 4: velocity phase 1
1634 j = 5, 6, 7: velocity phase 2 */
1635 for (int j = 0; j < _nVar; j++){
1637 VecGetValues(_primitiveVars,1,&Ii,&_VV(i,j));
1640 if(_saveConservativeField){
1641 for (long i = 0; i < _Nmailles; i++){
1642 for (int j = 0; j < _nVar; j++){
1644 // cout<<"i= "<<i<< " j= "<<j<<" UU(i,j)= "<<_UU(i,j)<<endl;
1645 VecGetValues(_conservativeVars,1,&Ii,&_UU(i,j));
1648 _UU.setTime(_time,_nbTimeStep);
1650 _VV.setTime(_time,_nbTimeStep);
1651 if (_nbTimeStep ==0 || _restartWithNewFileName){
1652 string prim_suppress ="rm -rf "+prim+"_*";
1653 string cons_suppress ="rm -rf "+cons+"_*";
1654 system(prim_suppress.c_str());//Nettoyage des précédents calculs identiques
1655 system(cons_suppress.c_str());//Nettoyage des précédents calculs identiques
1656 _VV.setInfoOnComponent(0,"Void_fraction");
1657 _VV.setInfoOnComponent(1,"Pressure_(Pa)");
1658 _VV.setInfoOnComponent(2,"Velocity1_x_m/s");
1660 _VV.setInfoOnComponent(3,"Velocity1_y_m/s");
1662 _VV.setInfoOnComponent(4,"Velocity1_z_m/s");
1664 _VV.setInfoOnComponent(2+_Ndim,"Velocity2_x_m/s");
1666 _VV.setInfoOnComponent(3+_Ndim,"Velocity2_y_m/s");
1668 _VV.setInfoOnComponent(4+_Ndim,"Velocity2_z_m/s");
1683 if(_saveConservativeField){
1684 _UU.setInfoOnComponent(0,"Partial_density1");// (kg/m^3)
1685 _UU.setInfoOnComponent(1,"Momentum1_x");// phase1 (kg/m^2/s)
1687 _UU.setInfoOnComponent(2,"Momentum1_y");// phase1 (kg/m^2/s)
1689 _UU.setInfoOnComponent(3,"Momentum1_z");// phase1 (kg/m^2/s)
1691 _UU.setInfoOnComponent(1+_Ndim,"Partial_density2");// phase2 (kg/m^3)
1692 _UU.setInfoOnComponent(2+_Ndim,"Momentum2_x");// phase2 (kg/m^2/s)
1694 _UU.setInfoOnComponent(3+_Ndim,"Momentum2_y");// phase2 (kg/m^2/s)
1696 _UU.setInfoOnComponent(4+_Ndim,"Momentum2_z");// phase2 (kg/m^2/s)
1716 _VV.writeVTK(prim,false);
1719 _VV.writeMED(prim,false);
1726 if(_saveConservativeField){
1730 _UU.writeVTK(cons,false);
1733 _UU.writeMED(cons,false);
1742 for (long i = 0; i < _Nmailles; i++){
1743 // j = 0 : concentration, j=1 : pressure; j = _nVar - 1: temperature; j = 2,..,_nVar-2: velocity
1744 for (int j = 0; j < _Ndim; j++)//On récupère les composantes de vitesse
1747 VecGetValues(_primitiveVars,1,&Ii,&_Vitesse1(i,j));
1748 Ii=i*_nVar +2+j+_Ndim;
1749 VecGetValues(_primitiveVars,1,&Ii,&_Vitesse2(i,j));
1751 for (int j = _Ndim; j < 3; j++){//On met à zero les composantes de vitesse si la dimension est <3
1756 _Vitesse1.setTime(_time,_nbTimeStep);
1757 _Vitesse2.setTime(_time,_nbTimeStep);
1758 if (_nbTimeStep ==0 || _restartWithNewFileName){
1759 _Vitesse1.setInfoOnComponent(0,"Velocity_x_(m/s)");
1760 _Vitesse1.setInfoOnComponent(1,"Velocity_y_(m/s)");
1761 _Vitesse1.setInfoOnComponent(2,"Velocity_z_(m/s)");
1763 _Vitesse2.setInfoOnComponent(0,"Velocity_x_(m/s)");
1764 _Vitesse2.setInfoOnComponent(1,"Velocity_y_(m/s)");
1765 _Vitesse2.setInfoOnComponent(2,"Velocity_z_(m/s)");
1770 _Vitesse1.writeVTK(prim+"_GasVelocity");
1771 _Vitesse2.writeVTK(prim+"_LiquidVelocity");
1774 _Vitesse1.writeMED(prim+"_GasVelocity");
1775 _Vitesse2.writeMED(prim+"_LiquidVelocity");
1778 _Vitesse1.writeCSV(prim+"_GasVelocity");
1779 _Vitesse2.writeCSV(prim+"_LiquidVelocity");
1787 _Vitesse1.writeVTK(prim+"_GasVelocity",false);
1788 _Vitesse2.writeVTK(prim+"_LiquidVelocity",false);
1791 _Vitesse1.writeMED(prim+"_GasVelocity",false);
1792 _Vitesse2.writeMED(prim+"_LiquidVelocity",false);
1795 _Vitesse1.writeCSV(prim+"_GasVelocity");
1796 _Vitesse2.writeCSV(prim+"_LiquidVelocity");
1802 if (_restartWithNewFileName)
1803 _restartWithNewFileName=false;