You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

1523 lines
45KB

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