Droplet Lagrangian Transient One-dimensional Reacting Code Implementation of both liquid and gas phase governing equations.
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.

1918 lines
57KB

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