Droplet Lagrangian Transient One-dimensional Reacting Code Implementation of both liquid and gas phase governing equations.
Você não pode selecionar mais de 25 tópicos Os tópicos devem começar com uma letra ou um número, podem incluir traços ('-') e podem ter até 35 caracteres.

1618 linhas
48KB

  1. #ifndef GSL_DEF
  2. #define GSL_DEF
  3. #include <gsl/gsl_math.h>
  4. #include <gsl/gsl_spline.h>
  5. #endif
  6. #include "residue.h"
  7. #include "macros.h"
  8. #include <cmath>
  9. #include <stdio.h>
  10. #include "timing.hpp"
  11. double maxGradPosition(const double* y, const size_t nt,
  12. const size_t nvar, const double* x, size_t nPts){
  13. double maxGradT=0.0e0;
  14. double gradT=0.0e0;
  15. double pos=0.0e0;
  16. size_t j,jm;
  17. for (size_t i = 1; i <nPts; i++) {
  18. j=i*nvar+nt;
  19. jm=(i-1)*nvar+nt;
  20. gradT=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  21. if (gradT>=maxGradT) {
  22. maxGradT=gradT;
  23. pos=x[i];
  24. }
  25. }
  26. return(pos);
  27. }
  28. int maxGradIndex(const double* y, const size_t nt,
  29. const size_t nvar, const double* x, size_t nPts){
  30. double maxGradT=0.0e0;
  31. double gradT=0.0e0;
  32. int pos=0.0e0;
  33. size_t j,jm;
  34. for (size_t i = 1; i <nPts; i++) {
  35. j=i*nvar+nt;
  36. jm=(i-1)*nvar+nt;
  37. gradT=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  38. if (gradT>=maxGradT) {
  39. maxGradT=gradT;
  40. pos=i;
  41. }
  42. }
  43. return(pos);
  44. }
  45. double maxCurvPosition(const double* y, const size_t nt,
  46. const size_t nvar, const double* x, size_t nPts){
  47. double maxCurvT=0.0e0;
  48. double gradTp=0.0e0;
  49. double gradTm=0.0e0;
  50. double curvT=0.0e0;
  51. double dx=0.0e0;
  52. double pos=0.0e0;
  53. size_t j,jm,jp;
  54. for (size_t i = 1; i <nPts-1; i++) {
  55. j=i*nvar+nt;
  56. jm=(i-1)*nvar+nt;
  57. jp=(i+1)*nvar+nt;
  58. gradTp=fabs((y[jp]-y[j])/(x[i+1]-x[i]));
  59. gradTm=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  60. dx=0.5e0*((x[i]+x[i+1])-(x[i-1]+x[i]));
  61. curvT=(gradTp-gradTm)/dx;
  62. if (curvT>=maxCurvT) {
  63. maxCurvT=curvT;
  64. pos=x[i];
  65. }
  66. }
  67. return(pos);
  68. }
  69. int maxCurvIndex(const double* y, const size_t nt,
  70. const size_t nvar, const double* x, size_t nPts){
  71. double maxCurvT=0.0e0;
  72. double gradTp=0.0e0;
  73. double gradTm=0.0e0;
  74. double curvT=0.0e0;
  75. double dx=0.0e0;
  76. int pos=0;
  77. size_t j,jm,jp;
  78. for (size_t i = 1; i <nPts-1; i++) {
  79. j=i*nvar+nt;
  80. jm=(i-1)*nvar+nt;
  81. jp=(i+1)*nvar+nt;
  82. gradTp=fabs((y[jp]-y[j])/(x[i+1]-x[i]));
  83. gradTm=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  84. dx=0.5e0*((x[i]+x[i+1])-(x[i-1]+x[i]));
  85. curvT=(gradTp-gradTm)/dx;
  86. if (curvT>=maxCurvT) {
  87. maxCurvT=curvT;
  88. pos=i;
  89. }
  90. }
  91. return(pos);
  92. }
  93. //double isothermPosition(const double* y, const double T, const size_t nt,
  94. // const size_t nvar, const double* x, const size_t nPts){
  95. // double pos=x[nPts-1];
  96. // size_t j;
  97. // for (size_t i = 1; i <nPts; i++) {
  98. // j=i*nvar+nt;
  99. // if (y[j]<=T) {
  100. // pos=x[i];
  101. // break;
  102. // }
  103. // }
  104. // return(pos);
  105. //}
  106. //******************* scan the temperature from left to right ******************//
  107. double isothermPosition(const double* y, const double T, const size_t nt,
  108. const size_t nvar, const double* x, const size_t nPts){
  109. double pos=x[0];
  110. size_t j;
  111. for (size_t i = 0; i <nPts; i++) {
  112. j=i*nvar+nt;
  113. if (y[j]>=T) {
  114. pos=x[i];
  115. break;
  116. }
  117. }
  118. return(pos);
  119. }
  120. void updateSolution(double* y, double* ydot, const size_t nvar,
  121. const double xOld[],const double xNew[],const size_t nPts){
  122. double ytemp[nPts],ydottemp[nPts];
  123. gsl_interp_accel* acc;
  124. gsl_spline* spline;
  125. acc = gsl_interp_accel_alloc();
  126. spline = gsl_spline_alloc(gsl_interp_steffen, nPts);
  127. gsl_interp_accel* accdot;
  128. gsl_spline* splinedot;
  129. accdot = gsl_interp_accel_alloc();
  130. splinedot = gsl_spline_alloc(gsl_interp_steffen, nPts);
  131. for (size_t j = 0; j < nvar; j++) {
  132. for (size_t i = 0; i < nPts; i++) {
  133. ytemp[i]=y[j+i*nvar];
  134. ydottemp[i]=ydot[j+i*nvar];
  135. }
  136. gsl_spline_init(spline,xOld,ytemp,nPts);
  137. gsl_spline_init(splinedot,xOld,ydottemp,nPts);
  138. for (size_t i = 0; i < nPts; i++) {
  139. y[j+i*nvar]=gsl_spline_eval(spline,xNew[i],acc);
  140. ydot[j+i*nvar]=gsl_spline_eval(splinedot,xNew[i],accdot);
  141. }
  142. }
  143. //Exploring "fixing" boundary conditions:
  144. //for (size_t j = 1; j < nvar; j++) {
  145. // //printf("%15.6e\t%15.6e\n", y[j],y[j+nvar]);
  146. // y[j]=y[j+nvar];
  147. // //y[j+(nPts-1)*nvar]=y[j+(nPts-2)*nvar];
  148. // //ydot[j+nvar]=ydot[j];
  149. //}
  150. //y[0]=0.0e0;
  151. gsl_interp_accel_free(acc);
  152. gsl_spline_free(spline);
  153. gsl_interp_accel_free(accdot);
  154. gsl_spline_free(splinedot);
  155. }
  156. //Locate bath gas:
  157. size_t BathGasIndex(UserData data){
  158. size_t index=0;
  159. double max1;
  160. double max=data->gas->massFraction(0);
  161. for(size_t k=1;k<data->nsp;k++){
  162. max1=data->gas->massFraction(k);
  163. if(max1>=max){
  164. max=max1;
  165. index=k;
  166. }
  167. }
  168. return(index+1);
  169. }
  170. //Locate Oxidizer:
  171. size_t oxidizerIndex(UserData data){
  172. size_t index=0;
  173. for(size_t k=1;k<data->nsp;k++){
  174. if(data->gas->speciesName(k-1)=="O2"){
  175. index=k;
  176. }
  177. }
  178. return(index);
  179. }
  180. //Locate OH:
  181. size_t OHIndex(UserData data){
  182. size_t index=0;
  183. for(size_t k=1;k<data->nsp;k++){
  184. if(data->gas->speciesName(k-1)=="OH"){
  185. index=k;
  186. }
  187. }
  188. return(index);
  189. }
  190. //Locate HO2:
  191. size_t HO2Index(UserData data){
  192. size_t index=0;
  193. for(size_t k=1;k<data->nsp;k++){
  194. if(data->gas->speciesName(k-1)=="HO2"){
  195. index=k;
  196. }
  197. }
  198. return(index);
  199. }
  200. //Locate species index:
  201. size_t specIndex(UserData data,char *specName){
  202. size_t index=0;
  203. for(size_t k=1;k<data->nsp;k++){
  204. if(data->gas->speciesName(k-1)==specName){
  205. index=k;
  206. }
  207. }
  208. return(index);
  209. }
  210. int setAlgebraicVariables(N_Vector* id, UserData data){
  211. double *iddata;
  212. N_VConst(ONE, *id);
  213. iddata = N_VGetArrayPointer_OpenMP(*id);
  214. data->k_bath=BathGasIndex(data);
  215. data->k_oxidizer=oxidizerIndex(data);
  216. data->k_OH=OHIndex(data);
  217. data->k_HO2=HO2Index(data);
  218. data->k_drop=specIndex(data,data->dropSpec);
  219. printf("Oxidizer index: %lu\n",data->k_oxidizer);
  220. printf("Bath gas index:%lu\n",data->k_bath);
  221. printf("Droplet species index:%lu\n",data->k_drop);
  222. for (size_t i = 1; i <=data->npts; i++) {
  223. /*Algebraic variables: indicated by ZERO.*/
  224. Rid(i)=ZERO;
  225. Pid(i)=ZERO;
  226. Yid(i,data->k_bath)=ZERO;
  227. Mdotid(i)=ZERO;
  228. }
  229. Mdotid(1)=ZERO;
  230. Rid(1)=ONE;
  231. Yid(1,data->k_drop)=ZERO;
  232. Tid(data->npts)=ZERO;
  233. if(data->constantPressure){
  234. Pid(data->npts)=ONE;
  235. }else{
  236. Pid(data->npts)=ZERO;
  237. Rid(data->npts)=ONE;
  238. }
  239. if(data->dirichletInner){
  240. Tid(1)=ZERO;
  241. }else{
  242. Tid(1)=ONE;
  243. }
  244. if(data->dirichletOuter){
  245. for (size_t k = 1; k <=data->nsp; k++) {
  246. if(k!=data->k_bath){
  247. Yid(data->npts,k)=ONE;
  248. }
  249. Yid(1,k)=ZERO;
  250. Tid(data->npts)=ONE;
  251. }
  252. }else{
  253. for (size_t k = 1; k <=data->nsp; k++){
  254. Yid(1,k)=ZERO;
  255. Yid(data->npts,k)=ZERO;
  256. Tid(data->npts)=ZERO;
  257. }
  258. }
  259. return(0);
  260. }
  261. inline double calc_area(double x,int* i){
  262. switch (*i) {
  263. case 0:
  264. return(ONE);
  265. case 1:
  266. return(x);
  267. case 2:
  268. return(x*x);
  269. default:
  270. return(ONE);
  271. }
  272. }
  273. void readInitialCondition(FILE* input, double* ydata, const size_t nvar, const size_t nr, const size_t nPts, double Rg){
  274. FILE* output;output=fopen("test.dat","w");
  275. size_t bufLen=10000;
  276. size_t nRows=0;
  277. size_t nColumns=nvar;
  278. char buf[bufLen];
  279. char buf1[bufLen];
  280. char comment[1];
  281. char *ret;
  282. while (fgets(buf,bufLen, input)!=NULL){
  283. comment[0]=buf[0];
  284. if(strncmp(comment,"#",1)!=0){
  285. nRows++;
  286. }
  287. }
  288. rewind(input);
  289. printf("nRows: %ld\n", nRows);
  290. double y[nRows*nColumns];
  291. size_t i=0;
  292. while (fgets(buf,bufLen, input)!=NULL){
  293. comment[0]=buf[0];
  294. if(strncmp(comment,"#",1)==0){
  295. }
  296. else{
  297. ret=strtok(buf,"\t");
  298. size_t j=0;
  299. y[i*nColumns+j]=(double)(atof(ret));
  300. j++;
  301. while(ret!=NULL){
  302. ret=strtok(NULL,"\t");
  303. if(j<nColumns){
  304. y[i*nColumns+j]=(double)(atof(ret));
  305. }
  306. j++;
  307. }
  308. i++;
  309. }
  310. }
  311. for (i = 0; i < nRows; i++) {
  312. for (size_t j = 0; j < nColumns; j++) {
  313. fprintf(output, "%15.6e\t",y[i*nColumns+j]);
  314. }
  315. fprintf(output, "\n");
  316. }
  317. fclose(output);
  318. //double xOld[nRows],xNew[nPts],ytemp[nPts];
  319. double xOld[nRows],xNew[nPts],ytemp[nRows];
  320. for (size_t j = 0; j < nRows; j++) {
  321. xOld[j]=y[j*nColumns+nr];
  322. }
  323. //double dx=(xOld[nRows-1] - xOld[0])/((double)(nPts)-1.0e0);
  324. //for (size_t j = 0; j < nPts-1; j++) {
  325. // xNew[j]=(double)j*dx + xOld[0];
  326. //}
  327. //xNew[nPts-1]=xOld[nRows-1];
  328. double dX0;
  329. int NN = nPts-2;
  330. dX0 = (xOld[nRows-1]-xOld[0])*(pow(Rg,1.0/NN) - 1.0)/(pow(Rg,(NN+1.0)/NN) - 1.0);
  331. xNew[0] = xOld[0];
  332. for (size_t i = 1; i < nPts-1; i++) {
  333. xNew[i]=xNew[i-1]+dX0*pow(Rg,(i-1.0)/NN);
  334. }
  335. xNew[nPts-1] = xOld[nRows-1];
  336. gsl_interp_accel* acc;
  337. gsl_spline* spline;
  338. acc = gsl_interp_accel_alloc();
  339. spline = gsl_spline_alloc(gsl_interp_steffen, nRows);
  340. for (size_t j = 0; j < nColumns; j++) {
  341. for (size_t k = 0; k < nRows; k++) {
  342. ytemp[k]=y[j+k*nColumns];
  343. }
  344. gsl_spline_init(spline,xOld,ytemp,nRows);
  345. for (size_t k = 0; k < nPts; k++) {
  346. ydata[j+k*nColumns]=gsl_spline_eval(spline,xNew[k],acc);
  347. }
  348. }
  349. gsl_interp_accel_free(acc);
  350. gsl_spline_free(spline);
  351. }
  352. double systemMass(double* ydata,UserData data){
  353. double mass=0.0e0;
  354. double rho;
  355. for (size_t i = 2; i <=data->npts; i++) {
  356. data->gas->setState_TPY(T(i), P(i), &Y(i,1));
  357. rho=data->gas->density();
  358. //psi(i)=psi(i-1)+rho*(R(i)-R(i-1))*calc_area(HALF*(R(i)+R(i-1)),&m);
  359. mass+=rho*(R(i)-R(i-1))*calc_area(R(i),&data->metric);
  360. }
  361. return(mass);
  362. }
  363. int initializePsiGrid(double* ydata, double* psidata, UserData data){
  364. double rho,rhom;
  365. /*Create a psi grid that corresponds CONSISTENTLY to the spatial grid
  366. * "R" created above. Note that the Lagrangian variable psi has units
  367. * of kg. */
  368. psi(1)=ZERO;
  369. for (size_t i = 2; i <=data->npts; i++) {
  370. data->gas->setState_TPY(T(i), P(i), &Y(i,1));
  371. rho=data->gas->density();
  372. data->gas->setState_TPY(T(i-1), P(i-1), &Y(i-1,1));
  373. rhom=data->gas->density();
  374. //psi(i)=psi(i-1)+rho*(R(i)-R(i-1))*calc_area(HALF*(R(i)+R(i-1)),&data->metric);
  375. //psi(i)=psi(i-1)+rho*(R(i)-R(i-1))*calc_area(R(i),&data->metric);
  376. psi(i)=psi(i-1)+(R(i)-R(i-1))*(rho*calc_area(R(i),&data->metric)+rhom*calc_area(R(i-1),&data->metric))/TWO;
  377. }
  378. /*The mass of the entire system is the value of psi at the last grid
  379. * point. Normalize psi by this mass so that it varies from zero to
  380. * one. This makes psi dimensionless. So the mass needs to be
  381. * multiplied back in the approporiate places in the governing
  382. * equations so that units match.*/
  383. data->mass=psi(data->npts);
  384. for (size_t i = 1; i <=data->npts; i++) {
  385. psi(i)=psi(i)/data->mass;
  386. }
  387. return(0);
  388. }
  389. int setInitialCondition(N_Vector* y,
  390. N_Vector* ydot,
  391. UserData data){
  392. double* ydata;
  393. double* ydotdata;
  394. double* psidata;
  395. double* innerMassFractionsData, Rd, massDrop;
  396. double rhomhalf, lambdamhalf, YVmhalf[data->nsp];
  397. double f=ZERO;
  398. double g=ZERO;
  399. double perturb,rho;
  400. double epsilon=ZERO;
  401. int m,ier;
  402. ydata = N_VGetArrayPointer_OpenMP(*y);
  403. ydotdata = N_VGetArrayPointer_OpenMP(*ydot);
  404. innerMassFractionsData = data->innerMassFractions;
  405. Rd = data->Rd;
  406. massDrop = data->massDrop;
  407. //Mass of droplet
  408. //massDrop=1.0/3.0*Rd*Rd*Rd*997.0; //TODO: The density of the droplet should be a user input
  409. massDrop=1.0/3.0*Rd*Rd*Rd*684.0; //TODO:The density of the droplet(n-heptane) should be a user input
  410. if(data->adaptiveGrid){
  411. psidata = data->grid->xOld;
  412. }
  413. else{
  414. psidata = data->uniformGrid;
  415. }
  416. m=data->metric;
  417. data->innerTemperature=data->initialTemperature;
  418. for (size_t k = 1; k <=data->nsp; k++) {
  419. innerMassFractionsData[k-1]=data->gas->massFraction(k-1);
  420. }
  421. //Define Grid:
  422. double dR=(data->domainLength)/((double)(data->npts)-1.0e0);
  423. double dv=(pow(data->domainLength,1+data->metric)-pow(data->firstRadius*data->domainLength,1+data->metric))/((double)(data->npts)-1.0e0);
  424. for (size_t i = 1; i <=data->npts; i++) {
  425. if(data->metric==0){
  426. R(i)=Rd+(double)((i-1)*dR);
  427. }else{
  428. if(i==1){
  429. R(i)=ZERO;
  430. }else if(i==2){
  431. R(i)=data->firstRadius*data->domainLength;
  432. }else{
  433. R(i)=pow(pow(R(i-1),1+data->metric)+dv,1.0/((double)(1+data->metric)));
  434. }
  435. }
  436. T(i)=data->initialTemperature;
  437. for (size_t k = 1; k <=data->nsp; k++) {
  438. Y(i,k)=data->gas->massFraction(k-1); //Indexing different in Cantera
  439. }
  440. P(i)=data->initialPressure*Cantera::OneAtm;
  441. }
  442. R(data->npts)=data->domainLength+data->Rd;
  443. double Tmax;
  444. double Tmin=data->initialTemperature;
  445. double w=data->mixingWidth;
  446. double YN2=ZERO;
  447. double YO2=ZERO;
  448. double YFuel,YOxidizer,sum;
  449. //if(data->problemType==0){
  450. // data->gas->equilibrate("HP");
  451. // data->maxTemperature=data->gas->temperature();
  452. //}
  453. //else if(data->problemType==1){
  454. // /*Premixed Combustion: Equilibrium products comprise ignition
  455. // * kernel at t=0. The width of the kernel is "mixingWidth"
  456. // * shifted by "shift" from the center.*/
  457. // data->gas->equilibrate("HP");
  458. // Tmax=data->gas->temperature();
  459. // for (size_t i = 1; i <=data->npts; i++) {
  460. // g=HALF*(tanh((R(i)-data->shift)/w)+ONE); //increasing function of x
  461. // f=ONE-g; //decreasing function of x
  462. // T(i)=(Tmax-Tmin)*f+Tmin;
  463. // for (size_t k = 1; k <=data->nsp; k++) {
  464. // Y(i,k)=(data->gas->massFraction(k-1)-Y(i,k))*f+Y(i,k);
  465. // }
  466. // }
  467. // if(data->dirichletOuter){
  468. // T(data->npts)=data->wallTemperature;
  469. // }
  470. //}
  471. //else if(data->problemType==2){
  472. // FILE* input;
  473. // if(input=fopen("initialCondition.dat","r")){
  474. // readInitialCondition(input, ydata, data->nvar, data->nr, data->npts);
  475. // fclose(input);
  476. // }
  477. // else{
  478. // printf("file initialCondition.dat not found!\n");
  479. // return(-1);
  480. // }
  481. //}
  482. initializePsiGrid(ydata,psidata,data);
  483. if(data->adaptiveGrid){
  484. // if(data->problemType!=0){
  485. // data->grid->position=maxGradPosition(ydata, data->nt, data->nvar,
  486. // data->grid->xOld, data->npts);
  487. // //data->grid->position=maxCurvPosition(ydata, data->nt, data->nvar,
  488. // // data->grid->xOld, data->npts);
  489. // }
  490. // else{
  491. // }
  492. if(data->problemType!=3){
  493. data->grid->position=0.0e0;
  494. double x=data->grid->position+data->gridOffset*data->grid->leastMove;
  495. printf("New grid center:%15.6e\n",x);
  496. ier=reGrid(data->grid, x);
  497. if(ier==-1)return(-1);
  498. updateSolution(ydata, ydotdata, data->nvar,
  499. data->grid->xOld,data->grid->x,data->npts);
  500. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  501. }
  502. }
  503. else{
  504. double Rg = data->Rg, dpsi0;
  505. int NN = data->npts-2;
  506. double psiNew[data->npts];
  507. dpsi0 = (pow(Rg,1.0/NN) - 1.0)/(pow(Rg,(NN+1.0)/NN) - 1.0);
  508. psiNew[0] = 0;
  509. for (size_t i = 1; i < data->npts-1; i++) {
  510. psiNew[i]=psiNew[i-1]+dpsi0*pow(Rg,(i-1.0)/NN);
  511. }
  512. psiNew[data->npts-1] = 1.0;
  513. printf("Last point:%15.6e\n",psiNew[data->npts-1]);
  514. updateSolution(ydata, ydotdata, data->nvar,
  515. data->uniformGrid,psiNew,data->npts);
  516. storeGrid(psiNew,data->uniformGrid,data->npts);
  517. //double psiNew[data->npts];
  518. //double dpsi=1.0e0/((double)(data->npts)-1.0e0);
  519. //for (size_t i = 0; i < data->npts; i++) {
  520. // psiNew[i]=(double)(i)*dpsi;
  521. //}
  522. //printf("Last point:%15.6e\n",psiNew[data->npts-1]);
  523. //updateSolution(ydata, ydotdata, data->nvar,
  524. // data->uniformGrid,psiNew,data->npts);
  525. //storeGrid(psiNew,data->uniformGrid,data->npts);
  526. }
  527. if(data->problemType==0){
  528. data->gas->equilibrate("HP");
  529. data->maxTemperature=data->gas->temperature();
  530. }
  531. else if(data->problemType==1){
  532. /*Premixed Combustion: Equilibrium products comprise ignition
  533. * kernel at t=0. The width of the kernel is "mixingWidth"
  534. * shifted by "shift" from the center.*/
  535. data->gas->equilibrate("HP");
  536. Tmax=data->gas->temperature();
  537. for (size_t i = 1; i <=data->npts; i++) {
  538. g=HALF*(tanh((R(i)-data->shift)/w)+ONE); //increasing function of x
  539. f=ONE-g; //decreasing function of x
  540. T(i)=(Tmax-Tmin)*f+Tmin;
  541. for (size_t k = 1; k <=data->nsp; k++) {
  542. Y(i,k)=(data->gas->massFraction(k-1)-Y(i,k))*f+Y(i,k);
  543. }
  544. }
  545. if(data->dirichletOuter){
  546. T(data->npts)=data->wallTemperature;
  547. }
  548. }
  549. else if(data->problemType==2){
  550. FILE* input;
  551. if(input=fopen("initialCondition.dat","r")){
  552. readInitialCondition(input, ydata, data->nvar, data->nr, data->npts, data->Rg);
  553. fclose(input);
  554. }
  555. else{
  556. printf("file initialCondition.dat not found!\n");
  557. return(-1);
  558. }
  559. initializePsiGrid(ydata,psidata,data);
  560. }
  561. else if(data->problemType==3){
  562. FILE* input;
  563. if(input=fopen("restart.bin","r")){
  564. readRestart(y, ydot, input, data);
  565. fclose(input);
  566. printf("Restart solution loaded!\n");
  567. printf("Problem starting at t=%15.6e\n",data->tNow);
  568. return(0);
  569. }
  570. else{
  571. printf("file restart.bin not found!\n");
  572. return(-1);
  573. }
  574. }
  575. if(data->reflectProblem){
  576. double temp;
  577. int j=1;
  578. while (data->npts+1-2*j>=0) {
  579. temp=T(j);
  580. T(j)=T(data->npts+1-j);
  581. T(data->npts+1-j)=temp;
  582. for (size_t k = 1; k <=data->nsp; k++) {
  583. temp=Y(j,k);
  584. Y(j,k)=Y(data->npts+1-j,k);
  585. Y(data->npts+1-j,k)=temp;
  586. }
  587. j=j+1;
  588. }
  589. }
  590. /*Floor small values to zero*/
  591. for (size_t i = 1; i <=data->npts; i++) {
  592. for (size_t k = 1; k <=data->nsp; k++) {
  593. if(fabs(Y(i,k))<=data->massFractionTolerance){
  594. Y(i,k)=0.0e0;
  595. }
  596. }
  597. }
  598. //Set grid to location of maximum curvature:
  599. if(data->adaptiveGrid){
  600. data->grid->position=maxCurvPosition(ydata, data->nt, data->nvar,
  601. data->grid->x, data->npts);
  602. ier=reGrid(data->grid, data->grid->position);
  603. updateSolution(ydata, ydotdata, data->nvar,
  604. data->grid->xOld,data->grid->x,data->npts);
  605. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  606. }
  607. /*Ensure consistent boundary conditions*/
  608. //T(1)=T(2);
  609. //for (size_t k = 1; k <=data->nsp; k++) {
  610. // Y(1,k)=Y(2,k);
  611. // Y(data->npts,k)=Y(data->npts-1,k);
  612. //}
  613. return(0);
  614. }
  615. inline double Qdot(double* t,
  616. double* x,
  617. double* ignTime,
  618. double* kernelSize,
  619. double* maxQdot){
  620. double qdot;
  621. if(*x<=*kernelSize){
  622. if((*t)<=(*ignTime)){
  623. qdot=(*maxQdot);
  624. }
  625. else{
  626. qdot=0.0e0;
  627. }
  628. }else{
  629. qdot=0.0e0;
  630. }
  631. return(qdot);
  632. }
  633. inline void setGas(UserData data,
  634. double *ydata,
  635. size_t gridPoint){
  636. data->gas->setTemperature(T(gridPoint));
  637. data->gas->setMassFractions_NoNorm(&Y(gridPoint,1));
  638. data->gas->setPressure(P(gridPoint));
  639. }
  640. void getTransport(UserData data,
  641. double *ydata,
  642. size_t gridPoint,
  643. double *rho,
  644. double *lambda,
  645. double YV[]){
  646. double YAvg[data->nsp],
  647. XLeft[data->nsp],
  648. XRight[data->nsp],
  649. gradX[data->nsp];
  650. setGas(data,ydata,gridPoint);
  651. data->gas->getMoleFractions(XLeft);
  652. setGas(data,ydata,gridPoint+1);
  653. data->gas->getMoleFractions(XRight);
  654. for (size_t k = 1; k <=data->nsp; k++) {
  655. YAvg(k)=HALF*(Y(gridPoint,k)+
  656. Y(gridPoint+1,k));
  657. gradX(k)=(XRight(k)-XLeft(k))/
  658. (R(gridPoint+1)-R(gridPoint));
  659. }
  660. double TAvg = HALF*(T(gridPoint)+T(gridPoint+1));
  661. double gradT=(T(gridPoint+1)-T(gridPoint))/
  662. (R(gridPoint+1)-R(gridPoint));
  663. data->gas->setTemperature(TAvg);
  664. data->gas->setMassFractions_NoNorm(YAvg);
  665. data->gas->setPressure(P(gridPoint));
  666. *rho=data->gas->density();
  667. *lambda=data->trmix->thermalConductivity();
  668. data->trmix->getSpeciesFluxes(1,&gradT,data->nsp,
  669. gradX,data->nsp,YV);
  670. //setGas(data,ydata,gridPoint);
  671. }
  672. int residue(double t, N_Vector y, N_Vector ydot, N_Vector res, void *user_data){
  673. /*Declare and fetch nvectors and user data:*/
  674. double *ydata, *ydotdata, *resdata, *psidata, *innerMassFractionsData;
  675. UserData data;
  676. data = (UserData)user_data;
  677. size_t npts=data->npts;
  678. size_t nsp=data->nsp;
  679. size_t k_bath = data->k_bath;
  680. size_t k_drop = data->k_drop;
  681. ydata = N_VGetArrayPointer_OpenMP(y);
  682. ydotdata= N_VGetArrayPointer_OpenMP(ydot);
  683. resdata = N_VGetArrayPointer_OpenMP(res);
  684. if(data->adaptiveGrid==1){
  685. psidata = data->grid->x;
  686. }else{
  687. psidata = data->uniformGrid;
  688. }
  689. innerMassFractionsData = data->innerMassFractions;
  690. /* Grid stencil:*/
  691. /*-------|---------*---------|---------*---------|-------*/
  692. /*-------|---------*---------|---------*---------|-------*/
  693. /*-------|---------*---------|---------*---------|-------*/
  694. /*-------m-------mhalf-------j-------phalf-------p-------*/
  695. /*-------|---------*---------|---------*---------|-------*/
  696. /*-------|---------*---------|---------*---------|-------*/
  697. /*-------|<=======dxm=======>|<=======dxp=======>|-------*/
  698. /*-------|---------*<======dxav=======>*---------|-------*/
  699. /*-------|<================dxpm=================>|-------*/
  700. /* Various variables defined for book-keeping and storing previously
  701. * calculated values:
  702. * rho : densities at points m, mhalf, j, p, and phalf.
  703. * area : the matric at points m, mhalf, j, p, and phalf.
  704. * m : exponent that determines geometry;
  705. * lambda : thermal conductivities at mhalf and phalf.
  706. * mdot : mass flow rate at m, j, and p.
  707. * X : mole fractions at j and p.
  708. * YV : diffusion fluxes at mhalf and phalf.
  709. * Tgrad : temperature gradient at mhalf and phalf.
  710. * Tav : average temperature between two points.
  711. * Pav : average pressure between two points.
  712. * Yav : average mass fractions between two points.
  713. * Xgradhalf : mole fraction gradient at j.
  714. * Cpb : mass based bulk specific heat.
  715. * tranTerm : transient terms.
  716. * advTerm : advection terms.
  717. * diffTerm : diffusion terms.
  718. * srcTerm : source terms.
  719. */
  720. double rhomhalf, rhom, lambdamhalf, YVmhalf[nsp],
  721. rho,
  722. rhophalf, lambdaphalf, YVphalf[nsp],
  723. Cpb, Cvb, Cp[nsp], Cpl, dHvl, rhol, wdot[nsp], enthalpy[nsp], energy[nsp],
  724. tranTerm, diffTerm, srcTerm, advTerm,
  725. area,areamhalf,areaphalf,aream,areamhalfsq,areaphalfsq;
  726. /*Aliases for difference coefficients:*/
  727. double cendfm, cendfc, cendfp;
  728. /*Aliases for various grid spacings:*/
  729. double dpsip, dpsiav, dpsipm, dpsim, dpsimm;
  730. dpsip=dpsiav=dpsipm=dpsim=dpsimm=ONE;
  731. //double mass, mdotIn;
  732. double mass, massDrop;
  733. double sum, sum1, sum2, sum3;
  734. size_t j,k;
  735. int m;
  736. m=data->metric; //Unitless
  737. mass=data->mass; //Units: kg
  738. //massDrop=data->massDrop; //Units: kg
  739. //mdotIn=data->mdot*calc_area(R(npts),&m); //Units: kg/s
  740. // /*evaluate properties at j=1*************************/
  741. setGas(data,ydata,1);
  742. rhom=data->gas->density();
  743. Cpb=data->gas->cp_mass(); //J/kg/K
  744. Cvb=data->gas->cv_mass(); //J/kg/K
  745. //TODO: Create user input model for these. They should not be constant
  746. //Cpl=4.182e03; //J/kg/K for liquid water
  747. //Cpl= 5.67508633e-07*pow(T(1),4) - 7.78060597e-04*pow(T(1),3) + 4.08310544e-01*pow(T(1),2) //J/kg/K for liquid water
  748. // - 9.62429538e+01*T(1) + 1.27131046e+04;
  749. /*Update specific heat:Cpl for liquid n-heptane*/
  750. /*Based on the existiong data,a linear expression is used*/
  751. // Cpl= 12.10476*T(1) - 746.60143; // Unit:J/(kg*K), Need to be improved later
  752. Cpl = 2242.871642; //Unit:J/(kg*K),range:25~520K,Ginnings&Furukawa,NIST
  753. //dHvl=2.260e6; //J/kg heat of vaporization of water
  754. //double Tr = T(1)/647.1; //Reduced Temperature: Droplet Temperature divided by critical temperature of water
  755. //dHvl= 5.2053e07*pow(1-Tr,0.3199 - 0.212*Tr + 0.25795*Tr*Tr)/18.01; //J/kg latent heat of vaporization of water
  756. double Tr = T(1)/540.2; //Reduced Temperature:Droplet Temperature divided by critical temperature of liquid n-heptane, Source:NIST
  757. dHvl = 5.366e5*exp(-0.2831*Tr) * pow((1-Tr),0.2831); //Unit:J/kg, latent heat of vaporization of liquid n-heptane, Source: NIST
  758. /*Following section is related to the property of water*/
  759. //rhol = 997.0;
  760. //massDrop=1.0/3.0*R(1)*R(1)*R(1)*997.0; //TODO: The density of the droplet should be a user input
  761. /*Following section is related to the property of liquid n-heptane (mainly density rho)*/
  762. rhol = 684.0; //Unit:kg/m^3
  763. massDrop = 1.0/3.0*R(1)*R(1)*R(1)*rhol; //Unit:kg, this is the mass of liquid droplet
  764. aream= calc_area(R(1),&m);
  765. /*******************************************************************/
  766. /*Calculate values at j=2's m and mhalf*****************************/
  767. getTransport(data, ydata, 1, &rhomhalf,&lambdamhalf,YVmhalf);
  768. areamhalf= calc_area(HALF*(R(1)+R(2)),&m);
  769. aream = calc_area(R(1),&m);
  770. areamhalfsq= areamhalf*areamhalf;
  771. /*******************************************************************/
  772. /*Fill up res with left side (center) boundary conditions:**********/
  773. /*We impose zero fluxes at the center:*/
  774. /*Mass:*/
  775. if(data->quasiSteady){
  776. Rres(1)=Rdot(1); // Should remain satisfied for quasi-steady droplet
  777. }else{
  778. Rres(1)=Rdot(1) + Mdot(1)/rhol/aream;
  779. }
  780. /*Species:*/
  781. sum=ZERO;
  782. //TODO: Antoine Parameters should be part of user input, not hardcoded
  783. //Yres(1,k_drop)=Y(1,k_drop) - 1.0e5*pow(10,4.6543-(1435.264/(T(1)-64.848)))*data->gas->molecularWeight(k_drop-1)
  784. // / P(1) / data->gas->meanMolecularWeight();
  785. //Yres(1,k_drop)=Y(1,k_drop) - 1.0e3*pow(2.71828,16.7-(4060.0/(T(1)-37.0)))*data->gas->molecularWeight(k_drop-1)
  786. // / P(1) / data->gas->meanMolecularWeight();
  787. //Yres(1,k_drop)=Y(1,k_drop)-2.566785e-02;
  788. /*Following using the Antoine Parameters to calculate the pressure of volatile component(n-heptane)*/
  789. //Antoine parameter from NIST website
  790. double p_i=0.0 ;
  791. p_i = 1.0e5*pow(10,4.02832-(1268.636/(T(1)-56.199))) ; //unit:Pa,FOR N-HEPTANE
  792. //p_i = 1.0e3*pow(2.71828,16.7-(4060.0/(T(1)-37.0))); //Unit:Pa,FOR WATER
  793. Yres(1,k_drop)=Y(1,k_drop) - p_i * data->gas->molecularWeight(k_drop-1)
  794. / P(1) / data->gas->meanMolecularWeight();
  795. sum=sum+Y(1,k_drop);
  796. Mdotres(1)=Mdot(1) - (YVmhalf(k_drop)*areamhalf) / (1-Y(1,k_drop)); //Evaporating mass
  797. for (k = 1; k <=nsp; k++) {
  798. if(k!=k_bath && k!=k_drop){
  799. //TODO: May need to include chemical source term
  800. Yres(1,k)=Y(1,k)*Mdot(1) + YVmhalf(k)*areamhalf;
  801. //Yres(1,k)=YVmhalf(k)*areamhalf;
  802. sum=sum+Y(1,k);
  803. //if(fabs(Mdot(1))>1e-14){
  804. ///Yres(1,k)=innerMassFractionsData[k-1]-
  805. /// Y(1,k)-
  806. /// (YVmhalf(k)*areamhalf)/Mdot(1);
  807. //Yres(1,k)=Y(1,k)*Mdot(1) + YVmhalf(k)*areamhalf;
  808. //}
  809. //else{
  810. // //Yres(1,k)=Y(1,k)-innerMassFractionsData[k-1];
  811. // Yres(1,k)=Y(2,k)-Y(1,k);
  812. // /*The below flux boundary condition makes the
  813. // * problem more prone to diverge. How does one
  814. // * fix this?*/
  815. // //Yres(1,k)=YVmhalf(k);
  816. //}
  817. }
  818. }
  819. Yres(1,k_bath)=ONE-sum-Y(1,k_bath);
  820. /*Energy:*/
  821. if(data->dirichletInner){
  822. //Tres(1)=Tdot(1);
  823. //Tres(1)=T(1)-data->innerTemperature;
  824. Tres(1)=lambdamhalf*rhomhalf*areamhalfsq*(T(2)-T(1))/((psi(2)-psi(1))*mass)/dHvl - YVmhalf(k_drop)*areamhalf/(1-Y(1,k_drop));
  825. }else{
  826. //Tres(1)=Tdot(1) -
  827. // (rhomhalf*areamhalfsq*lambdamhalf*(T(2)-T(1))/((psi(2)-psi(1))*mass) - Mdot(1)*dHvl)
  828. // / (massDrop * Cpl);
  829. Tres(1)=Tdot(1) -
  830. (areamhalf*lambdamhalf*(T(2)-T(1))/(R(2)-R(1)) - Mdot(1)*dHvl)
  831. / (massDrop * Cpl);
  832. //Tres(1)=T(2)-T(1);
  833. //Tres(1)=Tdot(1) - (Pdot(1)/(rhom*Cpb))
  834. // +(double)(data->metric+1)*(rhomhalf*lambdamhalf*areamhalfsq*(T(2)-T(1))/psi(2)-psi(1));
  835. }
  836. /*Pressure:*/
  837. Pres(1)=P(2)-P(1);
  838. /*Fill up res with governing equations at inner points:*************/
  839. for (j = 2; j < npts; j++) {
  840. /*evaluate various mesh differences*///
  841. dpsip = (psi(j+1) - psi(j) )*mass;
  842. dpsim = (psi(j) - psi(j-1))*mass;
  843. dpsiav = HALF*(psi(j+1) - psi(j-1))*mass;
  844. dpsipm = (psi(j+1) - psi(j-1))*mass;
  845. /***********************************///
  846. /*evaluate various central difference coefficients*/
  847. cendfm = - dpsip / (dpsim*dpsipm);
  848. cendfc = (dpsip-dpsim) / (dpsip*dpsim);
  849. cendfp = dpsim / (dpsip*dpsipm);
  850. /**************************************************/
  851. /*evaluate properties at j*************************/
  852. setGas(data,ydata,j);
  853. rho=data->gas->density(); //kg/m^3
  854. Cpb=data->gas->cp_mass(); //J/kg/K
  855. Cvb=data->gas->cv_mass(); //J/kg/K
  856. data->gas->getNetProductionRates(wdot); //kmol/m^3
  857. data->gas->getEnthalpy_RT(enthalpy); //unitless
  858. data->gas->getCp_R(Cp); //unitless
  859. area = calc_area(R(j),&m); //m^2
  860. /*evaluate properties at p*************************/
  861. getTransport(data, ydata, j, &rhophalf,&lambdaphalf,YVphalf);
  862. areaphalf= calc_area(HALF*(R(j)+R(j+1)),&m);
  863. areaphalfsq= areaphalf*areaphalf;
  864. /**************************************************///
  865. /*Evaporating Mass*/
  866. Mdotres(j)=Mdot(j) - Mdot(j-1);
  867. /*Mass:*/
  868. /* ∂r/∂ψ = 1/ρA */
  869. Rres(j)=((R(j)-R(j-1))/dpsim)-(TWO/(rhom*aream+rho*area));
  870. /*Energy:*/
  871. /* ∂T/∂t = - ṁ(∂T/∂ψ)
  872. * + (∂/∂ψ)(λρA²∂T/∂ψ)
  873. * - (A/cₚ) ∑ YᵢVᵢcₚᵢ(∂T/∂ψ)
  874. * - (1/ρcₚ)∑ ώᵢhᵢ
  875. * + (1/ρcₚ)(∂P/∂t) */
  876. /*Notes:
  877. * λ has units J/m/s/K.
  878. * YᵢVᵢ has units kg/m^2/s.
  879. * hᵢ has units J/kmol, so we must multiply the enthalpy
  880. * defined above (getEnthalpy_RT) by T (K) and the gas constant
  881. * (J/kmol/K) to get the right units.
  882. * cₚᵢ has units J/kg/K, so we must multiply the specific heat
  883. * defined above (getCp_R) by the gas constant (J/kmol/K) and
  884. * divide by the molecular weight (kg/kmol) to get the right
  885. * units.
  886. * */
  887. //enthalpy formulation:
  888. //tranTerm = Tdot(j) - (Pdot(j)/(rho*Cpb));
  889. tranTerm = Tdot(j);
  890. sum=ZERO;
  891. sum1=ZERO;
  892. for (k = 1; k <=nsp; k++) {
  893. sum=sum+wdot(k)*enthalpy(k);
  894. sum1=sum1+(Cp(k)/data->gas->molecularWeight(k-1))*rho
  895. *HALF*(YVmhalf(k)+YVphalf(k));
  896. }
  897. sum=sum*Cantera::GasConstant*T(j);
  898. sum1=sum1*Cantera::GasConstant;
  899. diffTerm =(( (rhophalf*areaphalfsq*lambdaphalf*(T(j+1)-T(j))/dpsip)
  900. -(rhomhalf*areamhalfsq*lambdamhalf*(T(j)-T(j-1))/dpsim) )
  901. /(dpsiav*Cpb) )
  902. -(sum1*area*(cendfp*T(j+1)
  903. +cendfc*T(j)
  904. +cendfm*T(j-1))/Cpb);
  905. srcTerm = (sum-Qdot(&t,&R(j),&data->ignTime,&data->kernelSize,&data->maxQDot))/(rho*Cpb); // Qdot is forced heating to cause ignition, sum is the conversion of chemical enthalpy to sensible enthalpy
  906. advTerm = (Mdot(j)*(T(j)-T(j-1))/dpsim);
  907. Tres(j)= tranTerm - (Pdot(j)/(rho*Cpb))
  908. +advTerm
  909. -diffTerm
  910. +srcTerm;
  911. // //energy formulation:
  912. // tranTerm = Tdot(j);
  913. // sum=ZERO;
  914. // sum1=ZERO;
  915. // sum2=ZERO;
  916. // sum3=ZERO;
  917. // for (k = 1; k <=nsp; k++) {
  918. // energy(k)=enthalpy(k)-ONE;
  919. // sum=sum+wdot(k)*energy(k);
  920. // sum1=sum1+(Cp(k)/data->gas->molecularWeight(k-1))*rho
  921. // *HALF*(YVmhalf(k)+YVphalf(k));
  922. // sum2=sum2+(YVmhalf(k)/data->gas->molecularWeight(k-1));
  923. // sum3=sum3+(YVphalf(k)/data->gas->molecularWeight(k-1));
  924. // }
  925. // sum=sum*Cantera::GasConstant*T(j);
  926. // sum1=sum1*Cantera::GasConstant;
  927. // diffTerm =(( (rhophalf*areaphalfsq*lambdaphalf*(T(j+1)-T(j))/dpsip)
  928. // -(rhomhalf*areamhalfsq*lambdamhalf*(T(j)-T(j-1))/dpsim) )
  929. // /(dpsiav*Cvb) )
  930. // -(sum1*area*(cendfp*T(j+1)
  931. // +cendfc*T(j)
  932. // +cendfm*T(j-1))/Cvb);
  933. // srcTerm = (sum-Qdot(&t,&R(j),&data->ignTime,&data->kernelSize,&data->maxQDot))/(rho*Cvb);
  934. // advTerm = (mdotIn*(T(j)-T(j-1))/dpsim);
  935. // advTerm = advTerm + (Cantera::GasConstant*T(j)*area/Cvb)*((sum3-sum2)/dpsiav);
  936. // Tres(j)= tranTerm
  937. // +advTerm
  938. // -diffTerm
  939. // +srcTerm;
  940. /*Species:*/
  941. /* ∂Yᵢ/∂t = - ṁ(∂Yᵢ/∂ψ)
  942. * - (∂/∂ψ)(AYᵢVᵢ)
  943. * + (ώᵢWᵢ/ρ) */
  944. sum=ZERO;
  945. for (k = 1; k <=nsp; k++) {
  946. if(k!=k_bath){
  947. tranTerm = Ydot(j,k);
  948. diffTerm = (YVphalf(k)*areaphalf
  949. -YVmhalf(k)*areamhalf)/dpsiav;
  950. srcTerm = wdot(k)
  951. *(data->gas->molecularWeight(k-1))/rho;
  952. advTerm = (Mdot(j)*(Y(j,k)-Y(j-1,k))/dpsim);
  953. Yres(j,k)= tranTerm
  954. +advTerm
  955. +diffTerm
  956. -srcTerm;
  957. sum=sum+Y(j,k);
  958. }
  959. }
  960. Yres(j,k_bath)=ONE-sum-Y(j,k_bath);
  961. /*Pressure:*/
  962. Pres(j) = P(j+1)-P(j);
  963. /*Assign values evaluated at p and phalf to m
  964. * and mhalf to save some cpu cost:****************/
  965. areamhalf=areaphalf;
  966. areamhalfsq=areaphalfsq;
  967. aream=area;
  968. rhom=rho;
  969. rhomhalf=rhophalf;
  970. lambdamhalf=lambdaphalf;
  971. for (k = 1; k <=nsp; k++) {
  972. YVmhalf(k)=YVphalf(k);
  973. }
  974. /**************************************************/
  975. }
  976. /*******************************************************************///
  977. /*Fill up res with right side (wall) boundary conditions:***********/
  978. /*We impose zero fluxes at the wall:*/
  979. setGas(data,ydata,npts);
  980. rho=data->gas->density();
  981. area = calc_area(R(npts),&m);
  982. /*Mass:*/
  983. dpsim=(psi(npts)-psi(npts-1))*mass;
  984. Rres(npts)=((R(npts)-R(npts-1))/dpsim)-(TWO/(rhom*aream+rho*area));
  985. /*Energy:*/
  986. if(data->dirichletOuter){
  987. //Tres(npts)=T(npts)-data->wallTemperature;
  988. Tres(npts)=Tdot(npts);
  989. }
  990. else{
  991. Tres(npts)=T(npts)-T(npts-1);
  992. }
  993. /*Species:*/
  994. sum=ZERO;
  995. if(data->dirichletOuter){
  996. for (k = 1; k <=nsp; k++) {
  997. if(k!=k_bath){
  998. Yres(npts,k)=Ydot(npts,k);
  999. sum=sum+Y(npts,k);
  1000. }
  1001. }
  1002. }
  1003. else{
  1004. for (k = 1; k <=nsp; k++) {
  1005. if(k!=k_bath){
  1006. Yres(npts,k)=Y(npts,k)-Y(npts-1,k);
  1007. //Yres(npts,k)=YVmhalf(k);
  1008. sum=sum+Y(npts,k);
  1009. }
  1010. }
  1011. }
  1012. Yres(npts,k_bath)=ONE-sum-Y(npts,k_bath);
  1013. /*Pressure:*/
  1014. if(data->constantPressure){
  1015. Pres(npts)=Pdot(npts)-data->dPdt;
  1016. }
  1017. else{
  1018. Pres(npts)=R(npts)-data->domainLength;
  1019. //Pres(npts)=Rdot(npts);
  1020. }
  1021. /*Evaporating Mass*/
  1022. Mdotres(npts)=Mdot(npts) - Mdot(npts-1);
  1023. //for (j = 1; j <=npts; j++) {
  1024. // //for (k = 1; k <=nsp; k++) {
  1025. // // Yres(j,k)=Ydot(j,k);
  1026. // //}
  1027. // //Tres(j)=Tdot(j);
  1028. //}
  1029. return(0);
  1030. }
  1031. void printSpaceTimeHeader(UserData data)
  1032. {
  1033. fprintf((data->output), "%15s\t","#1");
  1034. for (size_t k = 1; k <=data->nvar+1; k++) {
  1035. fprintf((data->output), "%15lu\t",k+1);
  1036. }
  1037. fprintf((data->output), "%15lu\n",data->nvar+3);
  1038. fprintf((data->output), "%15s\t%15s\t%15s\t","#psi","time(s)","dpsi");
  1039. fprintf((data->output), "%15s\t%15s\t","radius(m)","Temp(K)");
  1040. for (size_t k = 1; k <=data->nsp; k++) {
  1041. fprintf((data->output), "%15s\t",data->gas->speciesName(k-1).c_str());
  1042. }
  1043. fprintf((data->output), "%15s\t","Pressure(Pa)");
  1044. fprintf((data->output), "%15s\n","Mdot (kg/s)");
  1045. }
  1046. void printSpaceTimeOutput(double t, N_Vector* y, FILE* output, UserData data)
  1047. {
  1048. double *ydata,*psidata;
  1049. ydata = N_VGetArrayPointer_OpenMP(*y);
  1050. if(data->adaptiveGrid){
  1051. psidata = data->grid->x;
  1052. }else{
  1053. psidata = data->uniformGrid;
  1054. }
  1055. for (size_t i = 0; i < data->npts; i++) {
  1056. fprintf(output, "%15.9e\t%15.9e\t",psi(i+1),t);
  1057. if(i==0){
  1058. fprintf(output, "%15.9e\t",psi(2)-psi(1));
  1059. }
  1060. else{
  1061. fprintf(output, "%15.6e\t",psi(i+1)-psi(i));
  1062. }
  1063. for (size_t j = 0; j < data->nvar; j++) {
  1064. if(j!=data->nvar-1){
  1065. fprintf(output, "%15.9e\t",ydata[j+i*data->nvar]);
  1066. }else{
  1067. fprintf(output, "%15.9e",ydata[j+i*data->nvar]);
  1068. }
  1069. }
  1070. fprintf(output, "\n");
  1071. }
  1072. fprintf(output, "\n");
  1073. }
  1074. void writeRestart(double t, N_Vector* y, N_Vector* ydot, FILE* output, UserData data){
  1075. double *ydata,*psidata, *ydotdata;
  1076. ydata = N_VGetArrayPointer_OpenMP(*y);
  1077. ydotdata = N_VGetArrayPointer_OpenMP(*ydot);
  1078. if(data->adaptiveGrid){
  1079. psidata = data->grid->x;
  1080. }else{
  1081. psidata = data->uniformGrid;
  1082. }
  1083. fwrite(&t, sizeof(t), 1, output); //write time
  1084. fwrite(psidata, data->npts*sizeof(psidata), 1, output); //write grid
  1085. fwrite(ydata, data->neq*sizeof(ydata), 1, output); //write solution
  1086. fwrite(ydotdata, data->neq*sizeof(ydotdata), 1, output); //write solutiondot
  1087. }
  1088. void readRestart(N_Vector* y, N_Vector* ydot, FILE* input, UserData data){
  1089. double *ydata,*psidata, *ydotdata;
  1090. double t;
  1091. if(data->adaptiveGrid){
  1092. psidata = data->grid->x;
  1093. }else{
  1094. psidata = data->uniformGrid;
  1095. }
  1096. ydata = N_VGetArrayPointer_OpenMP(*y);
  1097. ydotdata = N_VGetArrayPointer_OpenMP(*ydot);
  1098. fread(&t, sizeof(t), 1, input);
  1099. data->tNow=t;
  1100. fread(psidata, data->npts*sizeof(psidata), 1, input);
  1101. fread(ydata, data->neq*sizeof(ydata), 1, input);
  1102. fread(ydotdata, data->neq*sizeof(ydotdata), 1, input);
  1103. if(data->adaptiveGrid){
  1104. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  1105. }
  1106. }
  1107. void printGlobalHeader(UserData data)
  1108. {
  1109. fprintf((data->globalOutput), "%8s\t","#Time(s)");
  1110. //fprintf((data->globalOutput), "%15s","S_u(exp*)(m/s)");
  1111. fprintf((data->globalOutput), "%15s","bFlux(kg/m^2/s)");
  1112. fprintf((data->globalOutput), "%16s"," IsothermPos(m)");
  1113. fprintf((data->globalOutput), "%15s\t","Pressure(Pa)");
  1114. fprintf((data->globalOutput), "%15s\t","Pdot(Pa/s)");
  1115. fprintf((data->globalOutput), "%15s\t","gamma");
  1116. fprintf((data->globalOutput), "%15s\t","S_u(m/s)");
  1117. fprintf((data->globalOutput), "%15s\t","Tu(K)");
  1118. fprintf((data->globalOutput), "\n");
  1119. }
  1120. //
  1121. //void printSpaceTimeRates(double t, N_Vector ydot, UserData data)
  1122. //{
  1123. // double *ydotdata,*psidata;
  1124. // ydotdata = N_VGetArrayPointer_OpenMP(ydot);
  1125. // psidata = N_VGetArrayPointer_OpenMP(data->grid);
  1126. // for (int i = 0; i < data->npts; i++) {
  1127. // fprintf((data->ratesOutput), "%15.6e\t%15.6e\t",psi(i+1),t);
  1128. // for (int j = 0; j < data->nvar; j++) {
  1129. // fprintf((data->ratesOutput), "%15.6e\t",ydotdata[j+i*data->nvar]);
  1130. // }
  1131. // fprintf((data->ratesOutput), "\n");
  1132. // }
  1133. // fprintf((data->ratesOutput), "\n\n");
  1134. //}
  1135. //
  1136. void printGlobalVariables(double t, N_Vector* y, N_Vector* ydot, UserData data)
  1137. {
  1138. double *ydata,*ydotdata, *innerMassFractionsData, *psidata;
  1139. innerMassFractionsData = data->innerMassFractions;
  1140. if(data->adaptiveGrid){
  1141. psidata = data->grid->x;
  1142. }else{
  1143. psidata = data->uniformGrid;
  1144. }
  1145. double TAvg, RAvg, YAvg, psiAvg;
  1146. ydata = N_VGetArrayPointer_OpenMP(*y);
  1147. ydotdata = N_VGetArrayPointer_OpenMP(*ydot);
  1148. TAvg=data->isotherm;
  1149. double sum=ZERO;
  1150. double dpsim,area,aream,drdt;
  1151. double Cpb,Cvb,gamma,rho,flameArea,Tu;
  1152. /*Find the isotherm chosen by the user*/
  1153. size_t j=1;
  1154. size_t jj=1;
  1155. size_t jjj=1;
  1156. double wdot[data->nsp];
  1157. double wdotMax=0.0e0;
  1158. double advTerm=0.0e0;
  1159. psiAvg=0.0e0;
  1160. if(T(data->npts)>T(1)){
  1161. while (T(j)<TAvg) {
  1162. j=j+1;
  1163. }
  1164. YAvg=innerMassFractionsData[data->k_oxidizer-1]-Y(data->npts,data->k_oxidizer);
  1165. while (fabs((T(jj+1)-T(jj))/T(jj))>1e-08) {
  1166. jj=jj+1;
  1167. }
  1168. setGas(data,ydata,jj);
  1169. Tu=T(jj);
  1170. rho=data->gas->density();
  1171. Cpb=data->gas->cp_mass(); //J/kg/K
  1172. Cvb=data->gas->cv_mass(); //J/kg/K
  1173. gamma=Cpb/Cvb;
  1174. }
  1175. else{
  1176. while (T(j)>TAvg) {
  1177. j=j+1;
  1178. }
  1179. YAvg=innerMassFractionsData[data->k_oxidizer-1]-Y(1,data->k_oxidizer);
  1180. while (fabs((T(data->npts-jj-1)-T(data->npts-jj))/T(data->npts-jj))>1e-08) {
  1181. jj=jj+1;
  1182. }
  1183. setGas(data,ydata,data->npts-jj);
  1184. Tu=T(data->npts-jj);
  1185. rho=data->gas->density();
  1186. Cpb=data->gas->cp_mass(); //J/kg/K
  1187. Cvb=data->gas->cv_mass(); //J/kg/K
  1188. gamma=Cpb/Cvb;
  1189. }
  1190. if(T(j)<TAvg){
  1191. RAvg=((R(j+1)-R(j))/(T(j+1)-T(j)))*(TAvg-T(j))+R(j);
  1192. }
  1193. else{
  1194. RAvg=((R(j)-R(j-1))/(T(j)-T(j-1)))*(TAvg-T(j-1))+R(j-1);
  1195. }
  1196. ////Experimental burning speed calculation:
  1197. //int nMax=0;
  1198. ////nMax=maxCurvIndex(ydata, data->nt, data->nvar,
  1199. //// data->grid->x, data->npts);
  1200. //nMax=maxGradIndex(ydata, data->nt, data->nvar,
  1201. // data->grid->x, data->npts);
  1202. //advTerm=(T(nMax)-T(nMax-1))/(data->mass*(psi(nMax)-psi(nMax-1)));
  1203. //aream=calc_area(R(nMax),&data->metric);
  1204. ////setGas(data,ydata,nMax);
  1205. ////rho=data->gas->density();
  1206. //psiAvg=-Tdot(nMax)/(rho*aream*advTerm);
  1207. ////if(t>data->ignTime){
  1208. //// for(size_t n=2;n<data->npts;n++){
  1209. //// setGas(data,ydata,n);
  1210. //// data->gas->getNetProductionRates(wdot); //kmol/m^3
  1211. //// advTerm=(T(n)-T(n-1))/(data->mass*(psi(n)-psi(n-1)));
  1212. //// if(fabs(wdot[data->k_oxidizer-1])>=wdotMax){
  1213. //// aream=calc_area(R(n),&data->metric);
  1214. //// psiAvg=-Tdot(n)/(rho*aream*advTerm);
  1215. //// wdotMax=fabs(wdot[data->k_oxidizer-1]);
  1216. //// }
  1217. //// }
  1218. ////}
  1219. ////else{
  1220. //// psiAvg=0.0e0;
  1221. ////}
  1222. //drdt=(RAvg-data->flamePosition[1])/(t-data->flameTime[1]);
  1223. //data->flamePosition[0]=data->flamePosition[1];
  1224. //data->flamePosition[1]=RAvg;
  1225. //data->flameTime[0]=data->flameTime[1];
  1226. //data->flameTime[1]=t;
  1227. //flameArea=calc_area(RAvg,&data->metric);
  1228. /*Use the Trapezoidal rule to calculate the mass burning rate based on
  1229. * the consumption of O2*/
  1230. aream= calc_area(R(1)+1e-03*data->domainLength,&data->metric);
  1231. for (j = 2; j <data->npts; j++) {
  1232. dpsim=(psi(j)-psi(j-1))*data->mass;
  1233. area= calc_area(R(j),&data->metric);
  1234. sum=sum+HALF*dpsim*((Ydot(j-1,data->k_oxidizer)/aream)
  1235. +(Ydot(j,data->k_oxidizer)/area));
  1236. aream=area;
  1237. }
  1238. //double maxOH,maxHO2;
  1239. //maxOH=0.0e0;
  1240. //maxHO2=0.0e0;
  1241. //for(j=1;j<data->npts;j++){
  1242. // if(Y(j,data->k_OH)>maxOH){
  1243. // maxOH=Y(j,data->k_OH);
  1244. // }
  1245. //}
  1246. //for(j=1;j<data->npts;j++){
  1247. // if(Y(j,data->k_HO2)>maxHO2){
  1248. // maxHO2=Y(j,data->k_HO2);
  1249. // }
  1250. //}
  1251. fprintf((data->globalOutput), "%15.6e\t",t);
  1252. //fprintf((data->globalOutput), "%15.6e\t",psiAvg);
  1253. fprintf((data->globalOutput), "%15.6e\t",fabs(sum)/YAvg);
  1254. fprintf((data->globalOutput), "%15.6e\t",RAvg);
  1255. fprintf((data->globalOutput), "%15.6e\t",P(data->npts));
  1256. fprintf((data->globalOutput), "%15.6e\t",Pdot(data->npts));
  1257. fprintf((data->globalOutput), "%15.6e\t",gamma);
  1258. fprintf((data->globalOutput), "%15.6e\t",fabs(sum)/(YAvg*rho));
  1259. fprintf((data->globalOutput), "%15.6e\t",Tu);
  1260. fprintf((data->globalOutput), "\n");
  1261. }
  1262. //
  1263. //void printSpaceTimeOutputInterpolated(double t, N_Vector y, UserData data)
  1264. //{
  1265. // double *ydata,*psidata;
  1266. // ydata = N_VGetArrayPointer_OpenMP(y);
  1267. // psidata = N_VGetArrayPointer_OpenMP(data->grid);
  1268. // for (int i = 0; i < data->npts; i++) {
  1269. // fprintf((data->gridOutput), "%15.6e\t%15.6e\t",psi(i+1),t);
  1270. // for (int j = 0; j < data->nvar; j++) {
  1271. // fprintf((data->gridOutput), "%15.6e\t",ydata[j+i*data->nvar]);
  1272. // }
  1273. // fprintf((data->gridOutput), "\n");
  1274. // }
  1275. // fprintf((data->gridOutput), "\n\n");
  1276. //}
  1277. //
  1278. //
  1279. ////void repairSolution(N_Vector y, N_Vector ydot, UserData data){
  1280. //// int npts=data->npts;
  1281. //// double *ydata;
  1282. //// double *ydotdata;
  1283. //// ydata = N_VGetArrayPointer_OpenMP(y);
  1284. //// ydotdata = N_VGetArrayPointer_OpenMP(ydot);
  1285. ////
  1286. //// T(2)=T(1);
  1287. //// T(npts-1)=T(npts);
  1288. //// Tdot(2)=Tdot(1);
  1289. //// Tdot(npts-1)=Tdot(npts);
  1290. //// for (int k = 1; k <=data->nsp; k++) {
  1291. //// Y(2,k)=Y(1,k);
  1292. //// Y(npts-1,k)=Y(npts,k);
  1293. ////
  1294. //// Ydot(2,k)=Ydot(1,k);
  1295. //// Ydot(npts-1,k)=Ydot(npts,k);
  1296. //// }
  1297. ////}
  1298. /*Following functions are added to derive the characteristic time scale of species*/
  1299. void getTimescale(UserData data, N_Vector* y){
  1300. size_t i, k, nsp,npts ;
  1301. nsp = data->nsp ;
  1302. npts = data->npts ;
  1303. double rho, wdot_mole[nsp], wdot_mass[nsp],MW[nsp],concentra[nsp];
  1304. //double time_scale[npts*nsp] ;
  1305. double* ydata ;
  1306. ydata = N_VGetArrayPointer_OpenMP(*y); //access the data stored in N_Vector* y
  1307. for(i=1; i<= npts;i++){
  1308. setGas(data,ydata,i); //set the gas state at each grid point
  1309. rho = data->gas->density(); //get the averaged density at each grid point, Unit:kg/m^3
  1310. data->gas->getNetProductionRates(wdot_mole) ; //Unit:kmol/(m^3 * s)
  1311. for(k=1; k<= nsp; k++){
  1312. MW(k) = data->gas->molecularWeight(k-1) ; //Unit:kg/kmol
  1313. }
  1314. for(k=1; k<= nsp; k++){
  1315. wdot_mass(k) = wdot_mole(k) * MW(k) ; //Unit:kg/(m^3 * s)
  1316. }
  1317. for(k=1; k<= nsp; k++){
  1318. concentra(k) = Y(i,k) * rho ; //Unit:kg/m^3
  1319. }
  1320. for(k=1;k<= nsp;k++){
  1321. data->time_scale(i,k) = concentra(k)/(wdot_mass(k)+ 1.00e-16) ;
  1322. }
  1323. }
  1324. }
  1325. void printTimescaleHeader(UserData data)
  1326. {
  1327. fprintf((data->timescaleOutput), "%15s\t","#1");
  1328. for (size_t k = 1; k <=data->nsp+1; k++) {
  1329. fprintf((data->timescaleOutput), "%15lu\t",k+1);
  1330. }
  1331. fprintf((data->timescaleOutput), "%15lu\n",data->nsp+3);
  1332. fprintf((data->timescaleOutput), "%15s\t%15s\t%15s\t","#time","radius","Temp(K)");
  1333. //fprintf((data->output), "%15s\t%15s\t","radius(m)","Temp(K)");
  1334. for (size_t k = 1; k <=(data->nsp); k++) {
  1335. fprintf((data->timescaleOutput), "%15s\t",data->gas->speciesName(k-1).c_str());
  1336. }
  1337. //fprintf((data->output), "%15s\t","Pressure(Pa)");
  1338. //fprintf((data->timescaleOutput), "%15s\n",data->gas->speciesName(data->nsp-1).c_str());
  1339. //fprintf((data->output), "%15s\n","Mdot (kg/s)");
  1340. fprintf((data->timescaleOutput), "\n");
  1341. }
  1342. void printTimescaleOutput(double t,N_Vector* y,FILE* output,UserData data)
  1343. {
  1344. double* ydata;
  1345. ydata = N_VGetArrayPointer_OpenMP(*y);
  1346. for(size_t i=1 ; i<= data->npts; i++){
  1347. fprintf(output, "%15.9e\t%15.9e\t",t,R(i));
  1348. fprintf(output, "%15.9e\t",T(i));
  1349. for(size_t k=1; k<=data->nsp;k++){
  1350. fprintf(output, "%15.9e\t",data->time_scale(i,k));
  1351. }
  1352. fprintf(output, "\n");
  1353. }
  1354. fprintf(output, "\n");
  1355. }
  1356. void floorSmallValue(UserData data, N_Vector* y)
  1357. {
  1358. double* ydata;
  1359. ydata = N_VGetArrayPointer_OpenMP(*y);
  1360. //double sum = 0.00;
  1361. size_t k_bath = data->k_bath;
  1362. /*Floor small values to zero*/
  1363. for (size_t i = 1; i <=data->npts; i++) {
  1364. for (size_t k = 1; k <=data->nsp; k++) {
  1365. if(fabs(Y(i,k))<=data->massFractionTolerance){
  1366. Y(i,k)=0.0e0;
  1367. }
  1368. }
  1369. }
  1370. /*Dump the error to the bath gas*/
  1371. for (size_t i = 1; i <=data->npts; i++) {
  1372. double sum = 0.00;
  1373. for (size_t k = 1; k <=data->nsp; k++) {
  1374. if(k!=k_bath){
  1375. sum = sum + Y(i,k);
  1376. }
  1377. }
  1378. Y(i,k_bath) = ONE-sum;
  1379. }
  1380. }
  1381. void resetTolerance(UserData data, N_Vector* y,N_Vector* atolv)
  1382. {
  1383. double* ydata;
  1384. realtype* atolvdata;
  1385. ydata = N_VGetArrayPointer_OpenMP(*y);
  1386. atolvdata = N_VGetArrayPointer_OpenMP(*atolv);
  1387. double relTol,radTol,tempTol,presTol,massFracTol,bathGasTol,mdotTol;
  1388. double maxT=0.00, ignT=1800;
  1389. relTol = 1.0e-03 ;
  1390. radTol = 1.0e-05 ;
  1391. tempTol = 1.0e-01 ;
  1392. presTol = 1.0e00 ;
  1393. massFracTol = 1.0e-08;
  1394. bathGasTol = 1.0e-04 ;
  1395. mdotTol = 1.0e-09 ;
  1396. /*Get the maximum Tempture*/
  1397. for (size_t i =1;i <= data->npts;i++){
  1398. if(T(i) > maxT){
  1399. maxT = T(i);
  1400. }
  1401. }
  1402. /*reset the tolerance when maxT > ignT*/
  1403. if(maxT >= ignT){
  1404. data->relativeTolerance = relTol;
  1405. for(size_t i =1; i<=data->npts; i++){
  1406. atolT(i) = tempTol;
  1407. atolR(i) = radTol ;
  1408. atolP(i) = presTol ;
  1409. atolMdot(i) = mdotTol ;
  1410. for(size_t k =1; k<= data->nsp; k++){
  1411. if(k!=data->k_bath){
  1412. atolY(i,k) = massFracTol;
  1413. }else{
  1414. atolY(i,k) = bathGasTol;
  1415. }
  1416. }
  1417. }
  1418. }
  1419. }