Droplet Lagrangian Transient One-dimensional Reacting Code Implementation of both liquid and gas phase governing equations.
No puede seleccionar más de 25 temas Los temas deben comenzar con una letra o número, pueden incluir guiones ('-') y pueden tener hasta 35 caracteres de largo.

1423 líneas
42KB

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