Droplet Lagrangian Transient One-dimensional Reacting Code Implementation of both liquid and gas phase governing equations.
Du kan inte välja fler än 25 ämnen Ämnen måste starta med en bokstav eller siffra, kan innehålla bindestreck ('-') och vara max 35 tecken långa.

1984 lines
60KB

  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. //Index here is 1-based
  124. int maxTemperatureIndex(const double* y,const size_t nt,const size_t nvar ,size_t nPts){
  125. double maxT = 0.0e0;
  126. double TempT = 0.0e0;
  127. int index = 0 ;
  128. for (size_t i=1;i<=nPts;i++){
  129. TempT = y[(i-1)*nvar+nt] ;
  130. if(TempT > maxT){
  131. maxT = TempT ;
  132. index = i;
  133. }
  134. }
  135. return(index);
  136. }
  137. double maxCurvPositionR(const double* y, const size_t nt,
  138. const size_t nvar, const size_t nr, size_t nPts){
  139. double maxCurvT=0.0e0;
  140. double gradTp=0.0e0;
  141. double gradTm=0.0e0;
  142. double curvT=0.0e0;
  143. double dx=0.0e0;
  144. double dr=0.0e0;
  145. double pos=0.0e0;
  146. size_t j,jm,jp;
  147. size_t r,rp,rm;
  148. for (size_t i = 1; i <nPts-1; i++) {
  149. j=i*nvar+nt;
  150. jm=(i-1)*nvar+nt;
  151. jp=(i+1)*nvar+nt;
  152. r = i*nvar+nr;
  153. rm = (i-1)*nvar+nr;
  154. rp = (i+1)*nvar+nr;
  155. //gradTp=fabs((y[jp]-y[j])/(x[i+1]-x[i]));
  156. //gradTm=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  157. //dx=0.5e0*((x[i]+x[i+1])-(x[i-1]+x[i]));
  158. //curvT=(gradTp-gradTm)/dx;
  159. gradTp = fabs((y[jp]-y[j])/(y[rp]-y[r]));
  160. gradTp = fabs((y[j]-y[jm])/(y[r]-y[rm]));
  161. dr = 0.5e0*(y[rp]-y[rm]);
  162. curvT=(gradTp-gradTm)/dr;
  163. if (curvT>=maxCurvT) {
  164. maxCurvT=curvT;
  165. pos=y[r];
  166. }
  167. }
  168. return(pos);
  169. }
  170. int maxCurvIndexR(const double* y, const size_t nt,
  171. const size_t nvar, const size_t nr, size_t nPts){
  172. double maxCurvT=0.0e0;
  173. double gradTp=0.0e0;
  174. double gradTm=0.0e0;
  175. double curvT=0.0e0;
  176. double dx=0.0e0;
  177. double dr=0.0e0;
  178. int pos=0;
  179. size_t j,jm,jp;
  180. size_t r,rm,rp;
  181. for (size_t i = 1; i <nPts-1; i++) {
  182. j=i*nvar+nt;
  183. jm=(i-1)*nvar+nt;
  184. jp=(i+1)*nvar+nt;
  185. r = i*nvar+nr;
  186. rm = (i-1)*nvar+nr;
  187. rp = (i+1)*nvar+nr;
  188. //gradTp=fabs((y[jp]-y[j])/(x[i+1]-x[i]));
  189. //gradTm=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  190. //dx=0.5e0*((x[i]+x[i+1])-(x[i-1]+x[i]));
  191. //curvT=(gradTp-gradTm)/dx;
  192. gradTp = fabs((y[jp]-y[j])/(y[rp]-y[r]));
  193. gradTp = fabs((y[j]-y[jm])/(y[r]-y[rm]));
  194. dr = 0.5e0*(y[rp]-y[rm]);
  195. curvT=(gradTp-gradTm)/dr;
  196. if (curvT>=maxCurvT) {
  197. maxCurvT=curvT;
  198. pos=i;
  199. }
  200. }
  201. return(pos);
  202. }
  203. double maxGradPosition(const double* y, const size_t nt,
  204. const size_t nvar, const double* x, size_t nPts){
  205. double maxGradT=0.0e0;
  206. double gradT=0.0e0;
  207. double pos=0.0e0;
  208. size_t j,jm;
  209. for (size_t i = 1; i <nPts; i++) {
  210. j=i*nvar+nt;
  211. jm=(i-1)*nvar+nt;
  212. gradT=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  213. if (gradT>=maxGradT) {
  214. maxGradT=gradT;
  215. pos=x[i];
  216. }
  217. }
  218. return(pos);
  219. }
  220. int maxGradIndex(const double* y, const size_t nt,
  221. const size_t nvar, const double* x, size_t nPts){
  222. double maxGradT=0.0e0;
  223. double gradT=0.0e0;
  224. int pos=0;
  225. size_t j,jm;
  226. for (size_t i = 1; i <nPts; i++) {
  227. j=i*nvar+nt;
  228. jm=(i-1)*nvar+nt;
  229. gradT=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  230. if (gradT>=maxGradT) {
  231. maxGradT=gradT;
  232. pos=i;
  233. }
  234. }
  235. return(pos);
  236. }
  237. double maxCurvPosition(const double* y, const size_t nt,
  238. const size_t nvar, const double* x, size_t nPts){
  239. double maxCurvT=0.0e0;
  240. double gradTp=0.0e0;
  241. double gradTm=0.0e0;
  242. double curvT=0.0e0;
  243. double dx=0.0e0;
  244. double pos=0.0e0;
  245. size_t j,jm,jp;
  246. for (size_t i = 1; i <nPts-1; i++) {
  247. j=i*nvar+nt;
  248. jm=(i-1)*nvar+nt;
  249. jp=(i+1)*nvar+nt;
  250. gradTp=fabs((y[jp]-y[j])/(x[i+1]-x[i]));
  251. gradTm=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  252. dx=0.5e0*((x[i]+x[i+1])-(x[i-1]+x[i]));
  253. curvT=(gradTp-gradTm)/dx;
  254. if (curvT>=maxCurvT) {
  255. maxCurvT=curvT;
  256. pos=x[i];
  257. }
  258. }
  259. return(pos);
  260. }
  261. int maxCurvIndex(const double* y, const size_t nt,
  262. const size_t nvar, const double* x, size_t nPts){
  263. double maxCurvT=0.0e0;
  264. double gradTp=0.0e0;
  265. double gradTm=0.0e0;
  266. double curvT=0.0e0;
  267. double dx=0.0e0;
  268. int pos=0;
  269. size_t j,jm,jp;
  270. for (size_t i = 1; i <nPts-1; i++) {
  271. j=i*nvar+nt;
  272. jm=(i-1)*nvar+nt;
  273. jp=(i+1)*nvar+nt;
  274. gradTp=fabs((y[jp]-y[j])/(x[i+1]-x[i]));
  275. gradTm=fabs((y[j]-y[jm])/(x[i]-x[i-1]));
  276. dx=0.5e0*((x[i]+x[i+1])-(x[i-1]+x[i]));
  277. curvT=(gradTp-gradTm)/dx;
  278. if (curvT>=maxCurvT) {
  279. maxCurvT=curvT;
  280. pos=i;
  281. }
  282. }
  283. return(pos);
  284. }
  285. //double isothermPosition(const double* y, const double T, const size_t nt,
  286. // const size_t nvar, const double* x, const size_t nPts){
  287. // double pos=x[nPts-1];
  288. // size_t j;
  289. // for (size_t i = 1; i <nPts; i++) {
  290. // j=i*nvar+nt;
  291. // if (y[j]<=T) {
  292. // pos=x[i];
  293. // break;
  294. // }
  295. // }
  296. // return(pos);
  297. //}
  298. //******************* scan the temperature from left to right ******************//
  299. double isothermPosition(const double* y, const double T, const size_t nt,
  300. const size_t nvar, const double* x, const size_t nPts){
  301. double pos=x[0];
  302. size_t j;
  303. for (size_t i = 0; i <nPts; i++) {
  304. j=i*nvar+nt;
  305. if (y[j]>=T) {
  306. pos=x[i];
  307. break;
  308. }
  309. }
  310. return(pos);
  311. }
  312. void updateSolution(double* y, double* ydot, const size_t nvar,
  313. const double xOld[],const double xNew[],const size_t nPts){
  314. double ytemp[nPts],ydottemp[nPts];
  315. gsl_interp_accel* acc;
  316. gsl_spline* spline;
  317. acc = gsl_interp_accel_alloc();
  318. spline = gsl_spline_alloc(gsl_interp_steffen, nPts);
  319. gsl_interp_accel* accdot;
  320. gsl_spline* splinedot;
  321. accdot = gsl_interp_accel_alloc();
  322. splinedot = gsl_spline_alloc(gsl_interp_steffen, nPts);
  323. for (size_t j = 0; j < nvar; j++) {
  324. for (size_t i = 0; i < nPts; i++) {
  325. ytemp[i]=y[j+i*nvar];
  326. ydottemp[i]=ydot[j+i*nvar];
  327. }
  328. gsl_spline_init(spline,xOld,ytemp,nPts);
  329. gsl_spline_init(splinedot,xOld,ydottemp,nPts);
  330. for (size_t i = 0; i < nPts; i++) {
  331. y[j+i*nvar]=gsl_spline_eval(spline,xNew[i],acc);
  332. ydot[j+i*nvar]=gsl_spline_eval(splinedot,xNew[i],accdot);
  333. }
  334. }
  335. //Exploring "fixing" boundary conditions:
  336. //for (size_t j = 1; j < nvar; j++) {
  337. // //printf("%15.6e\t%15.6e\n", y[j],y[j+nvar]);
  338. // y[j]=y[j+nvar];
  339. // //y[j+(nPts-1)*nvar]=y[j+(nPts-2)*nvar];
  340. // //ydot[j+nvar]=ydot[j];
  341. //}
  342. //y[0]=0.0e0;
  343. gsl_interp_accel_free(acc);
  344. gsl_spline_free(spline);
  345. gsl_interp_accel_free(accdot);
  346. gsl_spline_free(splinedot);
  347. }
  348. //Locate bath gas:
  349. size_t BathGasIndex(UserData data){
  350. size_t index=0;
  351. double max1;
  352. double max=data->gas->massFraction(0);
  353. for(size_t k=1;k<data->nsp;k++){
  354. max1=data->gas->massFraction(k);
  355. if(max1>=max){
  356. max=max1;
  357. index=k;
  358. }
  359. }
  360. return(index+1);
  361. }
  362. //Locate Oxidizer:
  363. size_t oxidizerIndex(UserData data){
  364. size_t index=0;
  365. for(size_t k=1;k<data->nsp;k++){
  366. if(data->gas->speciesName(k-1)=="O2"){
  367. index=k;
  368. }
  369. }
  370. return(index);
  371. }
  372. //Locate OH:
  373. size_t OHIndex(UserData data){
  374. size_t index=0;
  375. for(size_t k=1;k<data->nsp;k++){
  376. if(data->gas->speciesName(k-1)=="OH"){
  377. index=k;
  378. }
  379. }
  380. return(index);
  381. }
  382. //Locate HO2:
  383. size_t HO2Index(UserData data){
  384. size_t index=0;
  385. for(size_t k=1;k<data->nsp;k++){
  386. if(data->gas->speciesName(k-1)=="HO2"){
  387. index=k;
  388. }
  389. }
  390. return(index);
  391. }
  392. //Locate species index:
  393. size_t specIndex(UserData data,char *specName){
  394. size_t index=0;
  395. for(size_t k=1;k<data->nsp;k++){
  396. if(data->gas->speciesName(k-1)==specName){
  397. index=k;
  398. }
  399. }
  400. return(index);
  401. }
  402. int setAlgebraicVariables(N_Vector* id, UserData data){
  403. double *iddata;
  404. char* specPtr[2];
  405. N_VConst(ONE, *id);
  406. iddata = N_VGetArrayPointer_OpenMP(*id);
  407. data->k_bath=BathGasIndex(data);
  408. data->k_oxidizer=oxidizerIndex(data);
  409. data->k_OH=OHIndex(data);
  410. data->k_HO2=HO2Index(data);
  411. //use char* pointer to get the address
  412. for(int i=0;i<2;i++){
  413. specPtr[i] = data->dropSpec[i];
  414. }
  415. data->k_drop[0]=specIndex(data,specPtr[0]);
  416. data->k_drop[1]=specIndex(data,specPtr[1]);
  417. printf("Oxidizer index: %lu\n",data->k_oxidizer);
  418. printf("Bath gas index:%lu\n",data->k_bath);
  419. printf("Droplet species index:%lu,%lu\n",data->k_drop[0],data->k_drop[1]);
  420. for (size_t i = 1; i <=data->npts; i++) {
  421. /*Algebraic variables: indicated by ZERO.*/
  422. Rid(i)=ZERO;
  423. Pid(i)=ZERO;
  424. Yid(i,data->k_bath)=ZERO;
  425. Mdotid(i)=ZERO;
  426. }
  427. Mdotid(1)=ZERO;
  428. Rid(1)=ONE;
  429. Yid(1,data->k_drop[0])=ZERO;
  430. Yid(1,data->k_drop[1])=ZERO;
  431. Tid(data->npts)=ZERO;
  432. if(data->constantPressure){
  433. Pid(data->npts)=ONE;
  434. }else{
  435. Pid(data->npts)=ZERO;
  436. Rid(data->npts)=ONE;
  437. }
  438. if(data->dirichletInner){
  439. Tid(1)=ZERO;
  440. }else{
  441. Tid(1)=ONE;
  442. }
  443. if(data->dirichletOuter){
  444. for (size_t k = 1; k <=data->nsp; k++) {
  445. if(k!=data->k_bath){
  446. Yid(data->npts,k)=ONE;
  447. }
  448. Yid(1,k)=ZERO;
  449. Tid(data->npts)=ONE;
  450. }
  451. }else{
  452. for (size_t k = 1; k <=data->nsp; k++){
  453. Yid(1,k)=ZERO;
  454. Yid(data->npts,k)=ZERO;
  455. Tid(data->npts)=ZERO;
  456. }
  457. }
  458. return(0);
  459. }
  460. inline double calc_area(double x,int* i){
  461. switch (*i) {
  462. case 0:
  463. return(ONE);
  464. case 1:
  465. return(x);
  466. case 2:
  467. return(x*x);
  468. default:
  469. return(ONE);
  470. }
  471. }
  472. void readInitialCondition(FILE* input, double* ydata, const size_t nvar, const size_t nr, const size_t nPts, double Rg){
  473. FILE* output;output=fopen("test.dat","w");
  474. size_t bufLen=10000;
  475. size_t nRows=0;
  476. size_t nColumns=nvar;
  477. char buf[bufLen];
  478. char buf1[bufLen];
  479. char comment[1];
  480. char *ret;
  481. while (fgets(buf,bufLen, input)!=NULL){
  482. comment[0]=buf[0];
  483. if(strncmp(comment,"#",1)!=0){
  484. nRows++;
  485. }
  486. }
  487. rewind(input);
  488. printf("nRows: %ld\n", nRows);
  489. double y[nRows*nColumns];
  490. size_t i=0;
  491. while (fgets(buf,bufLen, input)!=NULL){
  492. comment[0]=buf[0];
  493. if(strncmp(comment,"#",1)==0){
  494. }
  495. else{
  496. ret=strtok(buf,"\t");
  497. size_t j=0;
  498. y[i*nColumns+j]=(double)(atof(ret));
  499. j++;
  500. while(ret!=NULL){
  501. ret=strtok(NULL,"\t");
  502. if(j<nColumns){
  503. y[i*nColumns+j]=(double)(atof(ret));
  504. }
  505. j++;
  506. }
  507. i++;
  508. }
  509. }
  510. for (i = 0; i < nRows; i++) {
  511. for (size_t j = 0; j < nColumns; j++) {
  512. fprintf(output, "%15.6e\t",y[i*nColumns+j]);
  513. }
  514. fprintf(output, "\n");
  515. }
  516. fclose(output);
  517. //double xOld[nRows],xNew[nPts],ytemp[nPts];
  518. double xOld[nRows],xNew[nPts],ytemp[nRows];
  519. for (size_t j = 0; j < nRows; j++) {
  520. xOld[j]=y[j*nColumns+nr];
  521. }
  522. //double dx=(xOld[nRows-1] - xOld[0])/((double)(nPts)-1.0e0);
  523. //for (size_t j = 0; j < nPts-1; j++) {
  524. // xNew[j]=(double)j*dx + xOld[0];
  525. //}
  526. //xNew[nPts-1]=xOld[nRows-1];
  527. double dX0;
  528. int NN = nPts-2;
  529. dX0 = (xOld[nRows-1]-xOld[0])*(pow(Rg,1.0/NN) - 1.0)/(pow(Rg,(NN+1.0)/NN) - 1.0);
  530. xNew[0] = xOld[0];
  531. for (size_t i = 1; i < nPts-1; i++) {
  532. xNew[i]=xNew[i-1]+dX0*pow(Rg,(i-1.0)/NN);
  533. }
  534. xNew[nPts-1] = xOld[nRows-1];
  535. gsl_interp_accel* acc;
  536. gsl_spline* spline;
  537. acc = gsl_interp_accel_alloc();
  538. spline = gsl_spline_alloc(gsl_interp_steffen, nRows);
  539. for (size_t j = 0; j < nColumns; j++) {
  540. for (size_t k = 0; k < nRows; k++) {
  541. ytemp[k]=y[j+k*nColumns];
  542. }
  543. gsl_spline_init(spline,xOld,ytemp,nRows);
  544. for (size_t k = 0; k < nPts; k++) {
  545. ydata[j+k*nColumns]=gsl_spline_eval(spline,xNew[k],acc);
  546. }
  547. }
  548. gsl_interp_accel_free(acc);
  549. gsl_spline_free(spline);
  550. }
  551. double systemMass(double* ydata,UserData data){
  552. double mass=0.0e0;
  553. double rho;
  554. for (size_t i = 2; i <=data->npts; i++) {
  555. data->gas->setState_TPY(T(i), P(i), &Y(i,1));
  556. rho=data->gas->density();
  557. //psi(i)=psi(i-1)+rho*(R(i)-R(i-1))*calc_area(HALF*(R(i)+R(i-1)),&m);
  558. mass+=rho*(R(i)-R(i-1))*calc_area(R(i),&data->metric);
  559. }
  560. return(mass);
  561. }
  562. int initializePsiGrid(double* ydata, double* psidata, UserData data){
  563. double rho,rhom;
  564. /*Create a psi grid that corresponds CONSISTENTLY to the spatial grid
  565. * "R" created above. Note that the Lagrangian variable psi has units
  566. * of kg. */
  567. psi(1)=ZERO;
  568. for (size_t i = 2; i <=data->npts; i++) {
  569. data->gas->setState_TPY(T(i), P(i), &Y(i,1));
  570. rho=data->gas->density();
  571. data->gas->setState_TPY(T(i-1), P(i-1), &Y(i-1,1));
  572. rhom=data->gas->density();
  573. //psi(i)=psi(i-1)+rho*(R(i)-R(i-1))*calc_area(HALF*(R(i)+R(i-1)),&data->metric);
  574. //psi(i)=psi(i-1)+rho*(R(i)-R(i-1))*calc_area(R(i),&data->metric);
  575. 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;
  576. }
  577. /*The mass of the entire system is the value of psi at the last grid
  578. * point. Normalize psi by this mass so that it varies from zero to
  579. * one. This makes psi dimensionless. So the mass needs to be
  580. * multiplied back in the approporiate places in the governing
  581. * equations so that units match.*/
  582. data->mass=psi(data->npts);
  583. for (size_t i = 1; i <=data->npts; i++) {
  584. psi(i)=psi(i)/data->mass;
  585. }
  586. return(0);
  587. }
  588. int setInitialCondition(N_Vector* y,
  589. N_Vector* ydot,
  590. UserData data){
  591. double* ydata;
  592. double* ydotdata;
  593. double* psidata;
  594. double* innerMassFractionsData, Rd, massDrop;
  595. double rhomhalf, lambdamhalf, YVmhalf[data->nsp];
  596. double f=ZERO;
  597. double g=ZERO;
  598. double perturb,rho;
  599. double epsilon=ZERO;
  600. int m,ier;
  601. ydata = N_VGetArrayPointer_OpenMP(*y);
  602. ydotdata = N_VGetArrayPointer_OpenMP(*ydot);
  603. innerMassFractionsData = data->innerMassFractions;
  604. Rd = data->Rd;
  605. massDrop = data->massDrop;
  606. //Mass of droplet
  607. //massDrop=1.0/3.0*Rd*Rd*Rd*997.0; //TODO: The density of the droplet should be a user input
  608. //massDrop=1.0/3.0*Rd*Rd*Rd*684.0; //TODO:The density of the droplet(n-heptane) should be a user input
  609. massDrop = 1.0/3.0*Rd*Rd*Rd*data->dropRho;
  610. if(data->adaptiveGrid){
  611. psidata = data->grid->xOld;
  612. }
  613. else{
  614. psidata = data->uniformGrid;
  615. }
  616. m=data->metric;
  617. data->innerTemperature=data->initialTemperature;
  618. for (size_t k = 1; k <=data->nsp; k++) {
  619. innerMassFractionsData[k-1]=data->gas->massFraction(k-1);
  620. }
  621. //Define Grid:
  622. double dR=(data->domainLength)/((double)(data->npts)-1.0e0);
  623. double dv=(pow(data->domainLength,1+data->metric)-pow(data->firstRadius*data->domainLength,1+data->metric))/((double)(data->npts)-1.0e0);
  624. //Initialize the R(i),T(i)(data->initialTemperature),Y(i,k),P(i)
  625. for (size_t i = 1; i <=data->npts; i++) {
  626. if(data->metric==0){
  627. R(i)=Rd+(double)((i-1)*dR);
  628. }else{
  629. if(i==1){
  630. R(i)=ZERO;
  631. }else if(i==2){
  632. R(i)=data->firstRadius*data->domainLength;
  633. }else{
  634. R(i)=pow(pow(R(i-1),1+data->metric)+dv,1.0/((double)(1+data->metric)));
  635. }
  636. }
  637. T(i)=data->initialTemperature;
  638. for (size_t k = 1; k <=data->nsp; k++) {
  639. Y(i,k)=data->gas->massFraction(k-1); //Indexing different in Cantera
  640. }
  641. P(i)=data->initialPressure*Cantera::OneAtm;
  642. }
  643. R(data->npts)=data->domainLength+data->Rd;
  644. // /********** test R(i) and the volumn between grids *****************/
  645. // printf("Print the first 4 R(i)s and volumn between them: \n") ;
  646. // printf("Grids: 1st:%2.6e,2nd:%2.6e,3rd:%2.6e,4th:%2.6e \n",R(1),R(2),R(3),R(4));
  647. // double v1,v2,v3,v4,v5;
  648. // v1 =pow(R(2),1+data->metric) - pow(R(1),1+data->metric);
  649. // v2 =pow(R(3),1+data->metric) - pow(R(2),1+data->metric);
  650. // v3 =pow(R(4),1+data->metric) - pow(R(3),1+data->metric);
  651. // v4 =pow(R(5),1+data->metric) - pow(R(4),1+data->metric);
  652. // v5 =pow(R(6),1+data->metric) - pow(R(5),1+data->metric);
  653. // printf("Volumn: 1st:%2.6e,2nd:%2.6e,3rd:%2.6e,4th:%2.6e,5th:%2.6e\n",v1,v2,v3,v4,v5) ;
  654. double Tmax;
  655. double Tmin=data->initialTemperature;
  656. double w=data->mixingWidth;
  657. double YN2=ZERO;
  658. double YO2=ZERO;
  659. double YFuel,YOxidizer,sum;
  660. //if(data->problemType==0){
  661. // data->gas->equilibrate("HP");
  662. // data->maxTemperature=data->gas->temperature();
  663. //}
  664. //else if(data->problemType==1){
  665. // /*Premixed Combustion: Equilibrium products comprise ignition
  666. // * kernel at t=0. The width of the kernel is "mixingWidth"
  667. // * shifted by "shift" from the center.*/
  668. // data->gas->equilibrate("HP");
  669. // Tmax=data->gas->temperature();
  670. // for (size_t i = 1; i <=data->npts; i++) {
  671. // g=HALF*(tanh((R(i)-data->shift)/w)+ONE); //increasing function of x
  672. // f=ONE-g; //decreasing function of x
  673. // T(i)=(Tmax-Tmin)*f+Tmin;
  674. // for (size_t k = 1; k <=data->nsp; k++) {
  675. // Y(i,k)=(data->gas->massFraction(k-1)-Y(i,k))*f+Y(i,k);
  676. // }
  677. // }
  678. // if(data->dirichletOuter){
  679. // T(data->npts)=data->wallTemperature;
  680. // }
  681. //}
  682. //else if(data->problemType==2){
  683. // FILE* input;
  684. // if(input=fopen("initialCondition.dat","r")){
  685. // readInitialCondition(input, ydata, data->nvar, data->nr, data->npts);
  686. // fclose(input);
  687. // }
  688. // else{
  689. // printf("file initialCondition.dat not found!\n");
  690. // return(-1);
  691. // }
  692. //}
  693. initializePsiGrid(ydata,psidata,data);
  694. if(data->adaptiveGrid){
  695. // if(data->problemType!=0){
  696. // data->grid->position=maxGradPosition(ydata, data->nt, data->nvar,
  697. // data->grid->xOld, data->npts);
  698. // //data->grid->position=maxCurvPosition(ydata, data->nt, data->nvar,
  699. // // data->grid->xOld, data->npts);
  700. // }
  701. // else{
  702. // }
  703. if(data->problemType!=3){
  704. data->grid->position=0.0e0;
  705. double x=data->grid->position+data->gridOffset*data->grid->leastMove;
  706. printf("New grid center:%15.6e\n",x);
  707. // /**************** TEST THE data->grid->xOld *******************/
  708. // double* ptr3 = data->grid->xOld ;
  709. // printf("SetInitialCondition function is called,before reGrid,Start print the first 5 elements of the xOld array : \n");
  710. // printf("1st:%.6f, 2nd:%.6f, 3rd:%.6f, 4th:%.6f, 5th:%.6f.\n",ptr3[0],ptr3[1],ptr3[2],ptr3[3],ptr3[4]);
  711. ier=reGrid(data->grid, x);
  712. if(ier==-1)return(-1);
  713. // /**************** TEST THE data->grid->xOld *******************/
  714. // double* ptr = data->grid->xOld ;
  715. // printf("SetInitialCondition function is called,after reGrid,Start print the first 5 elements of the xOld array : \n");
  716. // printf("1st:%.6f, 2nd:%.6f, 3rd:%.6f, 4th:%.6f, 5th:%.6f.\n",ptr[0],ptr[1],ptr[2],ptr[3],ptr[4]);
  717. updateSolution(ydata, ydotdata, data->nvar,
  718. data->grid->xOld,data->grid->x,data->npts);
  719. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  720. }
  721. }
  722. else{
  723. double Rg = data->Rg, dpsi0;
  724. int NN = data->npts-2;
  725. double psiNew[data->npts];
  726. dpsi0 = (pow(Rg,1.0/NN) - 1.0)/(pow(Rg,(NN+1.0)/NN) - 1.0);
  727. psiNew[0] = 0;
  728. for (size_t i = 1; i < data->npts-1; i++) {
  729. psiNew[i]=psiNew[i-1]+dpsi0*pow(Rg,(i-1.0)/NN);
  730. }
  731. psiNew[data->npts-1] = 1.0;
  732. printf("Last point:%15.6e\n",psiNew[data->npts-1]);
  733. updateSolution(ydata, ydotdata, data->nvar,
  734. data->uniformGrid,psiNew,data->npts);
  735. storeGrid(psiNew,data->uniformGrid,data->npts);
  736. //double psiNew[data->npts];
  737. //double dpsi=1.0e0/((double)(data->npts)-1.0e0);
  738. //for (size_t i = 0; i < data->npts; i++) {
  739. // psiNew[i]=(double)(i)*dpsi;
  740. //}
  741. //printf("Last point:%15.6e\n",psiNew[data->npts-1]);
  742. //updateSolution(ydata, ydotdata, data->nvar,
  743. // data->uniformGrid,psiNew,data->npts);
  744. //storeGrid(psiNew,data->uniformGrid,data->npts);
  745. }
  746. if(data->problemType==0){
  747. data->gas->equilibrate("HP");
  748. data->maxTemperature=data->gas->temperature();
  749. }
  750. else if(data->problemType==1){
  751. /*Premixed Combustion: Equilibrium products comprise ignition
  752. * kernel at t=0. The width of the kernel is "mixingWidth"
  753. * shifted by "shift" from the center.*/
  754. data->gas->equilibrate("HP");
  755. Tmax=data->gas->temperature();
  756. for (size_t i = 1; i <=data->npts; i++) {
  757. g=HALF*(tanh((R(i)-data->shift)/w)+ONE); //increasing function of x
  758. f=ONE-g; //decreasing function of x
  759. T(i)=(Tmax-Tmin)*f+Tmin;
  760. for (size_t k = 1; k <=data->nsp; k++) {
  761. Y(i,k)=(data->gas->massFraction(k-1)-Y(i,k))*f+Y(i,k);
  762. }
  763. }
  764. if(data->dirichletOuter){
  765. T(data->npts)=data->wallTemperature;
  766. }
  767. }
  768. else if(data->problemType==2){
  769. FILE* input;
  770. if(input=fopen("initialCondition.dat","r")){
  771. readInitialCondition(input, ydata, data->nvar, data->nr, data->npts, data->Rg);
  772. fclose(input);
  773. }
  774. else{
  775. printf("file initialCondition.dat not found!\n");
  776. return(-1);
  777. }
  778. initializePsiGrid(ydata,psidata,data);
  779. }
  780. else if(data->problemType==3){
  781. FILE* input;
  782. if(input=fopen("restart.bin","r")){
  783. readRestart(y, ydot, input, data);
  784. fclose(input);
  785. printf("Restart solution loaded!\n");
  786. printf("Problem starting at t=%15.6e\n",data->tNow);
  787. return(0);
  788. }
  789. else{
  790. printf("file restart.bin not found!\n");
  791. return(-1);
  792. }
  793. }
  794. if(data->reflectProblem){
  795. double temp;
  796. int j=1;
  797. while (data->npts+1-2*j>=0) {
  798. temp=T(j);
  799. T(j)=T(data->npts+1-j);
  800. T(data->npts+1-j)=temp;
  801. for (size_t k = 1; k <=data->nsp; k++) {
  802. temp=Y(j,k);
  803. Y(j,k)=Y(data->npts+1-j,k);
  804. Y(data->npts+1-j,k)=temp;
  805. }
  806. j=j+1;
  807. }
  808. }
  809. /*Floor small values to zero*/
  810. for (size_t i = 1; i <=data->npts; i++) {
  811. for (size_t k = 1; k <=data->nsp; k++) {
  812. if(fabs(Y(i,k))<=data->massFractionTolerance){
  813. Y(i,k)=0.0e0;
  814. }
  815. }
  816. }
  817. //Set grid to location of maximum curvature calculated by r instead of x:
  818. if(data->adaptiveGrid){
  819. //data->grid->position=maxCurvPosition(ydata, data->nt, data->nvar,
  820. // data->grid->x, data->npts);
  821. int maxCurvII = 0;
  822. maxCurvII = maxCurvIndexR(ydata,data->nt,data->nvar,data->nr,data->npts);
  823. data->grid->position = data->grid->x[maxCurvII];
  824. ier=reGrid(data->grid, data->grid->position);
  825. updateSolution(ydata, ydotdata, data->nvar,
  826. data->grid->xOld,data->grid->x,data->npts);
  827. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  828. /******** Test the maxCurvPosition and related variables *******/
  829. printf("The maxCurvPositition(based on r) is : %.6f,Temperature is : %.3f \n",data->grid->position,T(maxCurvII+1));
  830. }
  831. /*Ensure consistent boundary conditions*/
  832. //T(1)=T(2);
  833. //for (size_t k = 1; k <=data->nsp; k++) {
  834. // Y(1,k)=Y(2,k);
  835. // Y(data->npts,k)=Y(data->npts-1,k);
  836. //}
  837. return(0);
  838. }
  839. inline double Qdot(double* t,
  840. double* x,
  841. double* ignTime,
  842. double* kernelSize,
  843. double* maxQdot){
  844. double qdot;
  845. if(*x<=*kernelSize){
  846. if((*t)<=(*ignTime)){
  847. qdot=(*maxQdot);
  848. }
  849. else{
  850. qdot=0.0e0;
  851. }
  852. }else{
  853. qdot=0.0e0;
  854. }
  855. return(qdot);
  856. }
  857. inline void setGas(UserData data,
  858. double *ydata,
  859. size_t gridPoint){
  860. data->gas->setTemperature(T(gridPoint));
  861. data->gas->setMassFractions_NoNorm(&Y(gridPoint,1));
  862. data->gas->setPressure(P(gridPoint));
  863. }
  864. void getTransport(UserData data,
  865. double *ydata,
  866. size_t gridPoint,
  867. double *rho,
  868. double *lambda,
  869. double YV[]){
  870. double YAvg[data->nsp],
  871. XLeft[data->nsp],
  872. XRight[data->nsp],
  873. gradX[data->nsp];
  874. setGas(data,ydata,gridPoint);
  875. data->gas->getMoleFractions(XLeft);
  876. setGas(data,ydata,gridPoint+1);
  877. data->gas->getMoleFractions(XRight);
  878. for (size_t k = 1; k <=data->nsp; k++) {
  879. YAvg(k)=HALF*(Y(gridPoint,k)+
  880. Y(gridPoint+1,k));
  881. gradX(k)=(XRight(k)-XLeft(k))/
  882. (R(gridPoint+1)-R(gridPoint));
  883. }
  884. double TAvg = HALF*(T(gridPoint)+T(gridPoint+1));
  885. double gradT=(T(gridPoint+1)-T(gridPoint))/
  886. (R(gridPoint+1)-R(gridPoint));
  887. data->gas->setTemperature(TAvg);
  888. data->gas->setMassFractions_NoNorm(YAvg);
  889. data->gas->setPressure(P(gridPoint));
  890. *rho=data->gas->density();
  891. *lambda=data->trmix->thermalConductivity();
  892. data->trmix->getSpeciesFluxes(1,&gradT,data->nsp,
  893. gradX,data->nsp,YV);
  894. //setGas(data,ydata,gridPoint);
  895. }
  896. int residue(double t, N_Vector y, N_Vector ydot, N_Vector res, void *user_data){
  897. /*Declare and fetch nvectors and user data:*/
  898. double *ydata, *ydotdata, *resdata, *psidata, *innerMassFractionsData;
  899. UserData data;
  900. data = (UserData)user_data;
  901. size_t npts=data->npts;
  902. size_t nsp=data->nsp;
  903. size_t k_bath = data->k_bath;
  904. size_t k_drop[2] = {data->k_drop[0],data->k_drop[1]};
  905. ydata = N_VGetArrayPointer_OpenMP(y);
  906. ydotdata= N_VGetArrayPointer_OpenMP(ydot);
  907. resdata = N_VGetArrayPointer_OpenMP(res);
  908. if(data->adaptiveGrid==1){
  909. psidata = data->grid->x;
  910. }else{
  911. psidata = data->uniformGrid;
  912. }
  913. innerMassFractionsData = data->innerMassFractions;
  914. /* Grid stencil:*/
  915. /*-------|---------*---------|---------*---------|-------*/
  916. /*-------|---------*---------|---------*---------|-------*/
  917. /*-------|---------*---------|---------*---------|-------*/
  918. /*-------m-------mhalf-------j-------phalf-------p-------*/
  919. /*-------|---------*---------|---------*---------|-------*/
  920. /*-------|---------*---------|---------*---------|-------*/
  921. /*-------|<=======dxm=======>|<=======dxp=======>|-------*/
  922. /*-------|---------*<======dxav=======>*---------|-------*/
  923. /*-------|<================dxpm=================>|-------*/
  924. /* Various variables defined for book-keeping and storing previously
  925. * calculated values:
  926. * rho : densities at points m, mhalf, j, p, and phalf.
  927. * area : the matric at points m, mhalf, j, p, and phalf.
  928. * m : exponent that determines geometry;
  929. * lambda : thermal conductivities at mhalf and phalf.
  930. * mdot : mass flow rate at m, j, and p.
  931. * X : mole fractions at j and p.
  932. * YV : diffusion fluxes at mhalf and phalf.
  933. * Tgrad : temperature gradient at mhalf and phalf.
  934. * Tav : average temperature between two points.
  935. * Pav : average pressure between two points.
  936. * Yav : average mass fractions between two points.
  937. * Xgradhalf : mole fraction gradient at j.
  938. * Cpb : mass based bulk specific heat.
  939. * tranTerm : transient terms.
  940. * advTerm : advection terms.
  941. * diffTerm : diffusion terms.
  942. * srcTerm : source terms.
  943. */
  944. double rhomhalf, rhom, lambdamhalf, YVmhalf[nsp],
  945. rho,
  946. rhophalf, lambdaphalf, YVphalf[nsp],
  947. Cpb, Cvb, Cp[nsp], Cpl[2], dHvl[2], rhol, wdot[nsp], enthalpy[nsp], energy[nsp],
  948. tranTerm, diffTerm, srcTerm, advTerm,
  949. area,areamhalf,areaphalf,aream,areamhalfsq,areaphalfsq;
  950. /*Aliases for difference coefficients:*/
  951. double cendfm, cendfc, cendfp;
  952. /*Aliases for various grid spacings:*/
  953. double dpsip, dpsiav, dpsipm, dpsim, dpsimm;
  954. dpsip=dpsiav=dpsipm=dpsim=dpsimm=ONE;
  955. //double mass, mdotIn;
  956. double mass, massDrop;
  957. double sum, sum1, sum2, sum3;
  958. size_t j,k;
  959. int m;
  960. m=data->metric; //Unitless
  961. mass=data->mass; //Units: kg
  962. //massDrop=data->massDrop; //Units: kg
  963. //mdotIn=data->mdot*calc_area(R(npts),&m); //Units: kg/s
  964. // /*evaluate properties at j=1*************************/
  965. setGas(data,ydata,1);
  966. rhom=data->gas->density();
  967. Cpb=data->gas->cp_mass(); //J/kg/K
  968. Cvb=data->gas->cv_mass(); //J/kg/K
  969. //TODO: Create user input model for these. They should not be constant
  970. //Cpl=4.182e03; //J/kg/K for liquid water
  971. //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
  972. // - 9.62429538e+01*T(1) + 1.27131046e+04;
  973. /*Update specific heat:Cpl for liquid proprane&n-heptane&dodecane*/
  974. /*Based on the existiong data,a linear expression is used*/
  975. // Cpl= 12.10476*T(1) - 746.60143; // Unit:J/(kg*K), Need to be improved later
  976. //Cpl[0] = 2.423*T(1)+1661.074; //Unit:J/(kg*K),range:1atm,propane
  977. Cpl[0] = 3.336*T(1)+1289.5; //Unit:J/(kg*K),range:1atm,n-Heptane
  978. Cpl[1] = 3.815*T(1)+1081.8; //Unit:J/(kg*K),range:1atm,n-Dodecane
  979. double Cp_l = 0.0;
  980. for(size_t i=0;i<=1;i++){
  981. Cp_l=Cp_l+Cpl[i]*data->dropMassFrac[i];
  982. }
  983. //dHvl=2.260e6; //J/kg heat of vaporization of water
  984. //double Tr = T(1)/647.1; //Reduced Temperature: Droplet Temperature divided by critical temperature of water
  985. //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
  986. double Tr1 = T(1)/369.8; //Reduced Temperature;
  987. 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
  988. double Tr2 = T(1)/540.2; //Reduced Temperature:Droplet Temperature divided by critical temperature of liquid n-heptane, Source:NIST
  989. 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
  990. double dHv_l = 0.0;
  991. for(size_t i=0;i<=1;i++){
  992. dHv_l=dHv_l+dHvl[i]*data->dropMassFrac[i];
  993. }
  994. /*Following section is related to the property of water*/
  995. //rhol = 997.0;
  996. //massDrop=1.0/3.0*R(1)*R(1)*R(1)*997.0; //TODO: The density of the droplet should be a user input
  997. /*Following section is related to the property of liquid n-heptane and n-dodecane (mainly density rho)*/
  998. /*Density of these two species should be temperature dependent*/
  999. //data->dropDens[0] = -1.046*T(1)+823.794; //Unit:kg/m^3,density of propane @1atm
  1000. data->dropDens[0] = -0.858*T(1)+933.854; //Unit:kg/m^3,density of n-Heptane @1atm
  1001. data->dropDens[1] = -0.782*T(1)+979.643; //Unit:kg/m^3,density of n-Dodecane @1atm
  1002. //rhol = data->dropRho; //Unit:kg/m^3
  1003. rhol = 0.0;
  1004. for(size_t i=0;i<=1;i++){
  1005. rhol = rhol + data->dropMassFrac[i]*data->dropDens[i];
  1006. }
  1007. massDrop = 1.0/3.0*R(1)*R(1)*R(1)*rhol; //Unit:kg, this is the mass of liquid droplet
  1008. aream= calc_area(R(1),&m);
  1009. /*******************************************************************/
  1010. /*Calculate values at j=2's m and mhalf*****************************/
  1011. getTransport(data, ydata, 1, &rhomhalf,&lambdamhalf,YVmhalf);
  1012. areamhalf= calc_area(HALF*(R(1)+R(2)),&m);
  1013. aream = calc_area(R(1),&m);
  1014. areamhalfsq= areamhalf*areamhalf;
  1015. /*******************************************************************/
  1016. /*Fill up res with left side (center) boundary conditions:**********/
  1017. /*We impose zero fluxes at the center:*/
  1018. /*Mass:*/
  1019. if(data->quasiSteady){
  1020. Rres(1)=Rdot(1); // Should remain satisfied for quasi-steady droplet
  1021. }else{
  1022. Rres(1)=Rdot(1) + Mdot(1)/rhol/aream;
  1023. }
  1024. /*Species:*/
  1025. sum=ZERO;
  1026. //TODO: Antoine Parameters should be part of user input, not hardcoded
  1027. //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)
  1028. // / P(1) / data->gas->meanMolecularWeight();
  1029. //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)
  1030. // / P(1) / data->gas->meanMolecularWeight();
  1031. //Yres(1,k_drop)=Y(1,k_drop)-2.566785e-02;
  1032. /*Following using the Antoine Parameters to calculate the pressure of volatile component(n-heptane)*/
  1033. //Antoine parameter from NIST website
  1034. double p_i[2]={0.0,0.0} ;
  1035. p_i[0] = 1.0e5*pow(10,4.53678-(1149.36/(T(1)+24.906)))*data->dropMassFrac[0] ; //unit:Pa,For N-PROPANE
  1036. p_i[1] = 1.0e5*pow(10,4.02832-(1268.636/(T(1)-56.199)))*data->dropMassFrac[1] ; //unit:Pa,FOR N-HEPTANE
  1037. //p_i = 1.0e3*pow(2.71828,16.7-(4060.0/(T(1)-37.0))); //Unit:Pa,FOR WATER
  1038. //Yres(1,k_drop)=Y(1,k_drop) - p_i * data->gas->molecularWeight(k_drop-1)
  1039. // / P(1) / data->gas->meanMolecularWeight();
  1040. for(size_t i=0;i<=1;i++){
  1041. 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());
  1042. }
  1043. sum=sum+Y(1,k_drop[0])+Y(1,k_drop[1]);
  1044. //Mdotres(1)=Mdot(1) - (YVmhalf(k_drop)*areamhalf) / (1-Y(1,k_drop)); //Evaporating mass
  1045. 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
  1046. for (k = 1; k <=nsp; k++) {
  1047. if(k!=k_bath && k!=k_drop[0] && k!=k_drop[1]){
  1048. //TODO: May need to include chemical source term
  1049. Yres(1,k)=Y(1,k)*Mdot(1) + YVmhalf(k)*areamhalf;
  1050. //Yres(1,k)=YVmhalf(k)*areamhalf;
  1051. sum=sum+Y(1,k);
  1052. //if(fabs(Mdot(1))>1e-14){
  1053. ///Yres(1,k)=innerMassFractionsData[k-1]-
  1054. /// Y(1,k)-
  1055. /// (YVmhalf(k)*areamhalf)/Mdot(1);
  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. }
  1672. void getReactions(UserData data,N_Vector* y,FILE* output){
  1673. double Tmax;
  1674. double* ydata;
  1675. //double deltaT = 400.0;
  1676. //int i = 0;
  1677. size_t nRxns;
  1678. int index;
  1679. nRxns = data->gas->nReactions();
  1680. //DEBUG
  1681. printf("Total Number of Rxns:%zu \n",nRxns);
  1682. double fwdROP[nRxns],revROP[nRxns],netROP[nRxns];
  1683. std::string* rxnArr = new std::string[nRxns];
  1684. /*DEBUG*/
  1685. printf("Memory Allocation for Rxns Array Succeed!\n");
  1686. ydata = N_VGetArrayPointer_OpenMP(*y);
  1687. Tmax = maxTemperature(ydata,data->nt,data->nvar,data->npts);
  1688. //get the rxns' equation array
  1689. for(size_t ii=0;ii<nRxns;ii++){
  1690. rxnArr[ii]= data->gas->reactionString(ii);
  1691. }
  1692. if(Tmax>=(data->initialTemperature+data->deltaT)){
  1693. index = maxTemperatureIndex(ydata,data->nt,data->nvar,data->npts);
  1694. setGas(data,ydata,index);
  1695. /*Get forward/reverse/net rate of progress of rxns*/
  1696. data->gas->getRevRatesOfProgress(revROP);
  1697. data->gas->getFwdRatesOfProgress(fwdROP);
  1698. data->gas->getNetRatesOfProgress(netROP);
  1699. for(size_t j=0 ; j<nRxns; j++){
  1700. fprintf(output,"%30s\t",rxnArr[j]);
  1701. fprintf(output,"%15.9e\t%15.9e\t%15.9e\t\n",fwdROP[j],revROP[j],netROP[j]);
  1702. }
  1703. fprintf(output, "\n");
  1704. fclose(output);
  1705. }
  1706. delete[] rxnArr;
  1707. }
  1708. void getSpecies(UserData data,N_Vector* y,FILE* output){
  1709. double Tmax;
  1710. double* ydata;
  1711. //double deltaT = 400.0;
  1712. //int i = 0;
  1713. int index;
  1714. double fwdROP[data->nsp],revROP[data->nsp],netROP[data->nsp];
  1715. std::string* spArr = new std::string[data->nsp];
  1716. /*DEBUG*/
  1717. printf("Memory Allocation for Species Array Succeed!\n");
  1718. ydata = N_VGetArrayPointer_OpenMP(*y);
  1719. Tmax = maxTemperature(ydata,data->nt,data->nvar,data->npts);
  1720. //get the species name array
  1721. for(size_t ii=0;ii<data->nsp;ii++){
  1722. spArr[ii]= data->gas->speciesName(ii);
  1723. }
  1724. if(Tmax>=(data->initialTemperature+data->deltaT)){
  1725. //i = i + 1 ;
  1726. index = maxTemperatureIndex(ydata,data->nt,data->nvar,data->npts);
  1727. setGas(data,ydata,index);
  1728. /*Get forward/reverse/net rate of progress of rxns*/
  1729. data->gas->getDestructionRates(revROP);
  1730. data->gas->getCreationRates(fwdROP);
  1731. data->gas->getNetProductionRates(netROP);
  1732. /*Print data to the pertinent output file*/
  1733. for(size_t j=0 ; j<data->nsp; j++){
  1734. fprintf(output,"%15s\t",spArr[j]);
  1735. fprintf(output,"%15.9e\t%15.9e\t%15.9e\t\n",fwdROP[j],revROP[j],netROP[j]);
  1736. }
  1737. fprintf(output, "\n");
  1738. fclose(output);
  1739. }
  1740. delete[] spArr;
  1741. }