Droplet Lagrangian Transient One-dimensional Reacting Code Implementation of both liquid and gas phase governing equations.
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

489 lines
15KB

  1. /*
  2. _____ ___ ____ ____
  3. |_ _/ _ \| _ \ / ___|
  4. | || | | | |_) | |
  5. | || |_| | _ <| |___
  6. |_| \___/|_| \_\\____|
  7. */
  8. #include "UserData.h"
  9. #include "solution.h"
  10. #include "residue.h"
  11. #include "macros.h"
  12. #include "timing.h"
  13. #include <ida/ida.h>
  14. #include <ida/ida_direct.h>
  15. #include <sunmatrix/sunmatrix_band.h>
  16. #include <sunlinsol/sunlinsol_lapackband.h>
  17. //#include <ida/ida_band.h>
  18. static int check_flag(void *flagvalue,
  19. const char *funcname,
  20. int opt);
  21. void freeAtLast(void* mem, N_Vector *y,
  22. N_Vector *ydot,
  23. N_Vector *res,
  24. N_Vector *id,
  25. N_Vector *atolv,
  26. N_Vector *constraints,UserData data);
  27. int main(){
  28. // Read input file specifying the details of the case and store them
  29. FILE *input;input=fopen("input.dat","r");
  30. UserData data;data=NULL;data=allocateUserData(input);
  31. fclose(input);
  32. data->clockStart=get_wall_time();
  33. // /**************** TEST THE xOld *******************/
  34. // double* ptr1 = data->grid->xOld ;
  35. // printf("After allocateUserData in main.cpp,Start print the first 5 elements of the xOld array : \n");
  36. // printf("1st:%.6f, 2nd:%.6f, 3rd:%.6f, 4th:%.6f, 5th:%.6f.\n",ptr1[0],ptr1[1],ptr1[2],ptr1[3],ptr1[4]);
  37. if(data==NULL){
  38. printf("check input file!\n");
  39. freeUserData(data);
  40. return(-1);
  41. }
  42. // Allocate solution variables
  43. long int ier,mu,ml,count,netf,ncfn,njevals,nrevals;
  44. realtype tNow,*atolvdata,*constraintsdata,finalTime,dtMax,tolsfac;
  45. N_Vector y,ydot,id,res,atolv,constraints;
  46. y=ydot=id=res=atolv=constraints=NULL;
  47. ier=allocateSolution(data->neq,data->nThreads,&y,&ydot,&id,&res,&atolv,&constraints);
  48. ier=setAlgebraicVariables(&id,data);
  49. // /**************** TEST THE xOld *******************/
  50. // double* ptr2 = data->grid->xOld ;
  51. // printf("Before setInitialCondition in main.cpp,Start print the first 5 elements of the xOld array : \n");
  52. // printf("1st:%.6f, 2nd:%.6f, 3rd:%.6f, 4th:%.6f, 5th:%.6f.\n",ptr2[0],ptr2[1],ptr2[2],ptr2[3],ptr2[4]);
  53. ier=setInitialCondition(&y,&ydot,data);
  54. if(ier==-1)return(-1);
  55. tNow=data->tNow;
  56. finalTime=data->finalTime;
  57. //TODO: dtMax should be a user input
  58. dtMax = 1e-4;
  59. double* ydata;
  60. double* ydotdata;
  61. ydata = N_VGetArrayPointer_OpenMP(y);
  62. ydotdata = N_VGetArrayPointer_OpenMP(ydot);
  63. ////////// DEBUG ///////////////////
  64. //double* resdata;
  65. //double* iddata;
  66. //resdata = N_VGetArrayPointer_OpenMP(res);
  67. //iddata = N_VGetArrayPointer_OpenMP(id);
  68. ///////////////////////////////////
  69. void *mem;mem=NULL;mem = IDACreate();
  70. ier = IDASetUserData(mem, data);
  71. ier = IDASetId(mem, id);
  72. ier = IDAInit(mem, residue, tNow, y, ydot);
  73. // Atol array
  74. atolvdata = N_VGetArrayPointer_OpenMP(atolv);
  75. for (size_t i = 1; i <=data->npts; i++) {
  76. atolT(i) = data->temperatureTolerance;
  77. for (size_t k = 1; k <=data->nsp; k++) {
  78. if(k!=data->k_bath){
  79. atolY(i,k) = data->massFractionTolerance;
  80. }
  81. else{
  82. atolY(i,k) = data->bathGasTolerance;
  83. }
  84. }
  85. atolR(i) = data->radiusTolerance;
  86. atolP(i) = data->pressureTolerance;
  87. atolMdot(i) = data->MdotTolerance;
  88. }
  89. ier = IDASVtolerances(mem, data->relativeTolerance, atolv);
  90. mu = 2*data->nvar; ml = mu;
  91. SUNMatrix A; A=NULL;
  92. A=SUNBandMatrix(data->neq,mu,ml,mu+ml);
  93. SUNLinearSolver LS; LS=NULL;
  94. //LS=SUNBandLinearSolver(y,A);
  95. LS=SUNLapackBand(y,A);
  96. ier=IDADlsSetLinearSolver(mem,LS,A);
  97. //ier = IDABand(mem, data->neq, mu, ml);
  98. constraintsdata = N_VGetArrayPointer_OpenMP(constraints);
  99. if(data->setConstraints){
  100. for (size_t i = 1; i <=data->npts; i++) {
  101. for (size_t k = 1; k <=data->nsp; k++) {
  102. constraintsY(i,k) = ONE;
  103. }
  104. }
  105. ier=IDASetConstraints(mem, constraints);
  106. }
  107. if(!data->quasiSteady){
  108. constraintsR(1) = ONE;
  109. ier=IDASetConstraints(mem, constraints);
  110. }
  111. ier = IDASetSuppressAlg(mem, data->suppressAlg);
  112. if(check_flag(&ier, "IDASetSuppressAlg", 1)) return(1);
  113. //ier= IDASetMaxNumStepsIC(mem, 1);
  114. //ier= IDASetMaxNumJacsIC(mem,8);
  115. //ier= IDASetMaxNumItersIC(mem,100);
  116. //ier= IDASetMaxBacksIC(mem,2000);
  117. //ier = IDASetLineSearchOffIC(mem,SUNTRUE);
  118. //////// DEBUG ///////////
  119. //if(data->dryRun){
  120. // ier = residue(tNow,y, ydot, res, data);
  121. // for(size_t k=0; k < data->nvar*data->npts; k++){
  122. // printf("%i: %15.6e\n",k,resdata[k]);
  123. // }
  124. // for(size_t k=0; k < data->nvar*data->npts; k++){
  125. // if(iddata[k] == 1){
  126. // ydotdata[k] = -resdata[k];
  127. // }
  128. // }
  129. // ier = residue(tNow,y, ydot, res, data);
  130. // for(size_t k=0; k < data->nvar*data->npts; k++){
  131. // printf("%i: %15.6e\n",k,resdata[k]);
  132. // }
  133. //}
  134. //for(size_t k=0; k < data->neq; k++){
  135. // if(iddata[k] == 1){
  136. // ydotdata[k] = -resdata[k];
  137. // }
  138. //}
  139. //////////////////////////
  140. if(!data->dryRun){
  141. printf("Calculating Initial Conditions:\n");
  142. printf("Cross your fingers...\n");
  143. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, 1e-5*finalTime);
  144. //If at first IDACalcIC doesn't succeed, try, try, try again:
  145. for (int i = 0; i < 10; i++) {
  146. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, (1e-01+pow(10,i)*finalTime));
  147. /************* Print the #of iterations ************/
  148. //printf("This the %dth try of calculating initial conditions. \n",i);
  149. if(ier==0){
  150. break;
  151. }
  152. }
  153. //...or else cry again :(
  154. if(check_flag(&ier, "IDACalcIC", 1)){
  155. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  156. return(-1);
  157. }else{
  158. printf("Initial (Consistent) Conditions Calculated!\n");
  159. }
  160. ier = IDASetInitStep(mem,1e-12);
  161. }
  162. printSpaceTimeHeader(data);
  163. printGlobalHeader(data);
  164. printTimescaleHeader(data);
  165. printSpaceTimeOutput(tNow, &y, data->output, data);
  166. printSpaceTimeOutput(tNow, &y, data->gridOutput, data);
  167. // getTimescale(data,&y) ;
  168. // printTimescaleOutput(tNow, &y, data->timescaleOutput,data);
  169. if(!data->dryRun){
  170. count=0;
  171. double dt=1e-08;
  172. double t1=0.0e0;
  173. double xOld=0.0e0;
  174. double x=0.0e0;
  175. double dx=0.0e0;
  176. double dxMin=1.0e0;
  177. double dxRatio=dx/dxMin;
  178. int move=0;
  179. int kcur=0;
  180. int RGCOUNT=0;
  181. if(data->adaptiveGrid){
  182. dxMin=data->grid->leastMove;
  183. xOld=maxCurvPosition(ydata, data->nt, data->nvar,
  184. data->grid->x, data->npts);
  185. //xOld=isothermPosition(ydata, data->isotherm, data->nt,
  186. // data->nvar, data->grid->x, data->npts);
  187. }
  188. while (tNow<=finalTime && R(1)>100e-9) {
  189. t1=tNow;
  190. /*Floor small value to zero*/
  191. floorSmallValue(data, &y);
  192. if(data->quasiSteady){
  193. ier = IDASolve(mem, finalTime, &tNow, y, ydot, IDA_ONE_STEP);
  194. }else{
  195. /*This prevents the solver from taking a step so large that
  196. *the droplet radius becomes negative.
  197. *TODO:Try adding the constraint that the droplet radius must
  198. * be a positive number*/
  199. ier = IDASolve(mem, tNow+dtMax, &tNow, y, ydot, IDA_ONE_STEP);
  200. //ier = IDASolve(mem, tNow+dtMax, &tNow, y, ydot, IDA_NORMAL);
  201. }
  202. if(check_flag(&ier, "IDASolve", 1)){
  203. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  204. return(-1);
  205. }
  206. dt=tNow-t1;
  207. ier = IDAGetCurrentOrder(mem, &kcur);
  208. /******** Print the max Temperature *********/
  209. double maxT = 0.00;
  210. maxT = maxTemperature(ydata,data->nt,data->nvar,data->npts);
  211. printf("Maximum temperature is : %.3f[K] \n",maxT);
  212. if(data->adaptiveGrid==1 && data->moveGrid==1){
  213. x=maxCurvPosition(ydata, data->nt, data->nvar,
  214. data->grid->x, data->npts);
  215. //x=isothermPosition(ydata, data->isotherm, data->nt,
  216. // data->nvar, data->grid->x, data->npts);
  217. //x=maxGradPosition(ydata, data->nt, data->nvar,
  218. // data->grid->x, data->npts);
  219. dx=x-xOld;
  220. if(dx*dxMin>0.0e0){
  221. move=1;
  222. }else{
  223. move=-1;
  224. }
  225. //if(fabs(dx)>=dxMin && x+(double)(move)*0.5e0*dxMin<=1.0e0){
  226. dxRatio=fabs(dx/dxMin);
  227. /************ Print xOld, x,dx,dxMin,dxRatio ******************/
  228. printf("xOld = %.6f, x = %.6f,",xOld,x);
  229. printf("dx = %.6f, dxMin = %.6f, dxRatio = %.3f\n",dx,dxMin,dxRatio);
  230. if(dxRatio>=1.0e0 && dxRatio<=2.0e0){
  231. printf("Regridding!\n");
  232. data->regrid=1;
  233. printSpaceTimeOutput(tNow, &y, data->gridOutput, data);
  234. ier=reGrid(data->grid, x+(double)(move)*0.5e0*dxMin);
  235. if(ier==-1){
  236. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  237. return(-1);
  238. }
  239. updateSolution(ydata, ydotdata, data->nvar,
  240. data->grid->xOld,data->grid->x,data->npts);
  241. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  242. xOld=x;
  243. printf("Regrid Complete! Restarting Problem at %15.6e s\n",tNow);
  244. ier = IDAReInit(mem, tNow, y, ydot);
  245. if(check_flag(&ier, "IDAReInit", 1)){
  246. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  247. return(-1);
  248. }
  249. ier = IDASetInitStep(mem,1e-01*dt);
  250. printf("Reinitialized! Calculating Initial Conditions:\n");
  251. printf("Cross your fingers...\n");
  252. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e-01*dt);
  253. if(check_flag(&ier, "IDACalcIC", 1)){
  254. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e+01*dt);
  255. }
  256. //Every once in a while, for reasons
  257. //that befuddle this mathematically
  258. //lesser author, IDACalcIC fails. Here,
  259. //I desperately try to make it converge
  260. //again by sampling increasingly larger
  261. //time-steps:
  262. for (int i = 0; i < 10; i++) {
  263. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+(1e-01+pow(10,i)*dt));
  264. if(ier==0){
  265. break;
  266. }
  267. }
  268. //Failure :( Back to the drawing board:
  269. if(check_flag(&ier, "IDACalcIC", 1)){
  270. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  271. return(-1);
  272. }
  273. printf("Initial (Consistent) Conditions Calculated!\n");
  274. printf("------------------------------------------\n\n");
  275. if(data->writeEveryRegrid){
  276. printSpaceTimeOutput(tNow, &y, data->output, data);
  277. FILE* fp;
  278. fp=fopen("restart.bin","w");
  279. writeRestart(tNow, &y, &ydot, fp, data);
  280. fclose(fp);
  281. }
  282. }
  283. }
  284. /*reset the tolerance after ignition*/
  285. //resetTolerance(data,&y,&atolv);
  286. /*regrid and update the solution based on R,re-initialize the problem*/
  287. /*For the time being,we only allow TORC to REGRID once for each run*/
  288. if(data->JJRG ==1 && (maxT >= data->initialTemperature+data->deltaT)){
  289. if(RGCOUNT<1){
  290. RGCOUNT = RGCOUNT +1;
  291. REGRID(ydata,ydotdata,data);
  292. initializePsiGrid(ydata,data->uniformGrid,data);
  293. printf("REGRID Complete!Restarting Problem at %15.6e s\n",tNow);
  294. ier = IDAReInit(mem,tNow,y,ydot);
  295. if(check_flag(&ier,"IDAReInit",1)){
  296. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  297. return(-1);
  298. }
  299. ier= IDASetInitStep(mem,1e-01*dt);
  300. //ier= IDASetInitStep(mem,0.0);
  301. printf("Reinitialized!Calculating Initial Conditions:\n");
  302. printf("Cross your fingers...\n");
  303. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e-01*dt);
  304. if(check_flag(&ier, "IDACalcIC", 1)){
  305. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e+01*dt);
  306. }
  307. //Every once in a while, for reasons
  308. //that befuddle this mathematically
  309. //lesser author, IDACalcIC fails. Here,
  310. //I desperately try to make it converge
  311. //again by sampling increasingly larger
  312. //time-steps:
  313. for (int i = 0; i < 10; i++) {
  314. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+(1e-01+pow(10,i)*dt));
  315. if(ier==0){
  316. break;
  317. }
  318. }
  319. //Failure :( Back to the drawing board:
  320. if(check_flag(&ier, "IDACalcIC", 1)){
  321. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  322. return(-1);
  323. }
  324. printf("Initial (Consistent) Conditions Calculated!\n");
  325. printf("------------------------------------------\n\n");
  326. }
  327. }
  328. /*Floor small value to zero*/
  329. floorSmallValue(data, &y);
  330. if(count%data->nSaves==0 && !data->writeEveryRegrid){
  331. printSpaceTimeOutput(tNow, &y, data->output, data);
  332. //if(data->writeRates){
  333. // printSpaceTimeRates(tNow, ydot, data);
  334. //}
  335. }
  336. // getTimescale(data,&y);
  337. // if(count%data->nSaves==0){
  338. // printTimescaleOutput(tNow,&y, data->timescaleOutput,data);
  339. // //printSpaceTimeOutput(tNow, &y, data->output, data);
  340. // //if(data->writeRates){
  341. // // printSpaceTimeRates(tNow, ydot, data);
  342. // //}
  343. // }
  344. /*Print global variables only if time-step is of high order!*/
  345. if(data->nTimeSteps==0){
  346. data->flamePosition[0]=0.0e0;
  347. data->flamePosition[1]=0.0e0;
  348. data->flameTime[0]=tNow;
  349. data->flameTime[1]=tNow;
  350. }
  351. ier = IDAGetNumErrTestFails(mem, &netf);
  352. ier = IDAGetNumNonlinSolvConvFails(mem, &ncfn);
  353. ier = IDADlsGetNumJacEvals(mem, &njevals);
  354. ier = IDADlsGetNumResEvals(mem, &nrevals);
  355. printf("etf = %ld ,"
  356. "nlsf= %ld ,"
  357. "J=%ld ,"
  358. "R=%ld ,"
  359. "o=%d ,",netf, ncfn, njevals, nrevals, kcur);
  360. printf("Time=%15.6e s,",tNow);
  361. printf("dt=%15.6e s,",dt);
  362. printf("frac: %15.6e\n",dxRatio);
  363. count++;
  364. data->nTimeSteps=count;
  365. }
  366. }
  367. SUNLinSolFree(LS);
  368. SUNMatDestroy(A);
  369. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  370. return(0);
  371. }
  372. void freeAtLast(void* mem,
  373. N_Vector *y,
  374. N_Vector *ydot,
  375. N_Vector *res,
  376. N_Vector *id,
  377. N_Vector *atolv,
  378. N_Vector *constraints,UserData data){
  379. IDAFree(&mem);
  380. freeSolution(y,ydot,res,id,atolv,constraints);
  381. freeUserData(data);
  382. }
  383. static int check_flag(void *flagvalue, const char *funcname, int opt)
  384. {
  385. int *errflag;
  386. /* Check if SUNDIALS function returned NULL pointer - no memory allocated */
  387. if (opt == 0 && flagvalue == NULL) {
  388. fprintf(stderr,
  389. "\nSUNDIALS_ERROR: %s() failed - returned NULL pointer\n\n",
  390. funcname);
  391. return(1);
  392. }
  393. /* Check if flag < 0 */
  394. else if (opt == 1) {
  395. errflag = (int *) flagvalue;
  396. if (*errflag < 0) {
  397. fprintf(stderr,
  398. "\nSUNDIALS_ERROR: %s() failed with flag = %d\n\n",
  399. funcname, *errflag);
  400. return(1);
  401. }
  402. }
  403. /* Check if function returned NULL pointer - no memory allocated */
  404. else if (opt == 2 && flagvalue == NULL) {
  405. fprintf(stderr,
  406. "\nMEMORY_ERROR: %s() failed - returned NULL pointer\n\n",
  407. funcname);
  408. return(1);
  409. }
  410. return(0);
  411. }