Droplet Lagrangian Transient One-dimensional Reacting Code Implementation of both liquid and gas phase governing equations.
Nevar pievienot vairāk kā 25 tēmas Tēmai ir jāsākas ar burtu vai ciparu, tā var saturēt domu zīmes ('-') un var būt līdz 35 simboliem gara.

440 rindas
13KB

  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. if(data->adaptiveGrid){
  181. dxMin=data->grid->leastMove;
  182. xOld=maxCurvPosition(ydata, data->nt, data->nvar,
  183. data->grid->x, data->npts);
  184. //xOld=isothermPosition(ydata, data->isotherm, data->nt,
  185. // data->nvar, data->grid->x, data->npts);
  186. }
  187. while (tNow<=finalTime && R(1)>100e-9) {
  188. t1=tNow;
  189. /*Floor small value to zero*/
  190. floorSmallValue(data, &y);
  191. /*reset the tolerance after ignition*/
  192. resetTolerance(data,&y,&atolv);
  193. if(data->quasiSteady){
  194. ier = IDASolve(mem, finalTime, &tNow, y, ydot, IDA_ONE_STEP);
  195. }else{
  196. /*This prevents the solver from taking a step so large that
  197. *the droplet radius becomes negative.
  198. *TODO:Try adding the constraint that the droplet radius must
  199. * be a positive number*/
  200. ier = IDASolve(mem, tNow+dtMax, &tNow, y, ydot, IDA_ONE_STEP);
  201. //ier = IDASolve(mem, tNow+dtMax, &tNow, y, ydot, IDA_NORMAL);
  202. }
  203. if(check_flag(&ier, "IDASolve", 1)){
  204. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  205. return(-1);
  206. }
  207. dt=tNow-t1;
  208. ier = IDAGetCurrentOrder(mem, &kcur);
  209. /******** Print the max Temperature *********/
  210. double maxT = 0.00;
  211. maxT = maxTemperature(ydata,data->nt,data->nvar,data->npts);
  212. printf("Maximum temperature is : %.3f[K] \n",maxT);
  213. if(data->adaptiveGrid==1 && data->moveGrid==1){
  214. x=maxCurvPosition(ydata, data->nt, data->nvar,
  215. data->grid->x, data->npts);
  216. //x=isothermPosition(ydata, data->isotherm, data->nt,
  217. // data->nvar, data->grid->x, data->npts);
  218. //x=maxGradPosition(ydata, data->nt, data->nvar,
  219. // data->grid->x, data->npts);
  220. dx=x-xOld;
  221. if(dx*dxMin>0.0e0){
  222. move=1;
  223. }else{
  224. move=-1;
  225. }
  226. //if(fabs(dx)>=dxMin && x+(double)(move)*0.5e0*dxMin<=1.0e0){
  227. dxRatio=fabs(dx/dxMin);
  228. /************ Print xOld, x,dx,dxMin,dxRatio ******************/
  229. printf("xOld = %.6f, x = %.6f,",xOld,x);
  230. printf("dx = %.6f, dxMin = %.6f, dxRatio = %.3f\n",dx,dxMin,dxRatio);
  231. if(dxRatio>=1.0e0 && dxRatio<=2.0e0){
  232. printf("Regridding!\n");
  233. data->regrid=1;
  234. printSpaceTimeOutput(tNow, &y, data->gridOutput, data);
  235. ier=reGrid(data->grid, x+(double)(move)*0.5e0*dxMin);
  236. if(ier==-1){
  237. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  238. return(-1);
  239. }
  240. updateSolution(ydata, ydotdata, data->nvar,
  241. data->grid->xOld,data->grid->x,data->npts);
  242. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  243. xOld=x;
  244. printf("Regrid Complete! Restarting Problem at %15.6e s\n",tNow);
  245. ier = IDAReInit(mem, tNow, y, ydot);
  246. if(check_flag(&ier, "IDAReInit", 1)){
  247. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  248. return(-1);
  249. }
  250. ier = IDASetInitStep(mem,1e-01*dt);
  251. printf("Reinitialized! Calculating Initial Conditions:\n");
  252. printf("Cross your fingers...\n");
  253. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e-01*dt);
  254. if(check_flag(&ier, "IDACalcIC", 1)){
  255. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e+01*dt);
  256. }
  257. //Every once in a while, for reasons
  258. //that befuddle this mathematically
  259. //lesser author, IDACalcIC fails. Here,
  260. //I desperately try to make it converge
  261. //again by sampling increasingly larger
  262. //time-steps:
  263. for (int i = 0; i < 10; i++) {
  264. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+(1e-01+pow(10,i)*dt));
  265. if(ier==0){
  266. break;
  267. }
  268. }
  269. //Failure :( Back to the drawing board:
  270. if(check_flag(&ier, "IDACalcIC", 1)){
  271. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  272. return(-1);
  273. }
  274. printf("Initial (Consistent) Conditions Calculated!\n");
  275. printf("------------------------------------------\n\n");
  276. if(data->writeEveryRegrid){
  277. printSpaceTimeOutput(tNow, &y, data->output, data);
  278. FILE* fp;
  279. fp=fopen("restart.bin","w");
  280. writeRestart(tNow, &y, &ydot, fp, data);
  281. fclose(fp);
  282. }
  283. }
  284. }
  285. /*Floor small value to zero*/
  286. floorSmallValue(data, &y);
  287. if(count%data->nSaves==0 && !data->writeEveryRegrid){
  288. printSpaceTimeOutput(tNow, &y, data->output, data);
  289. //if(data->writeRates){
  290. // printSpaceTimeRates(tNow, ydot, data);
  291. //}
  292. }
  293. // getTimescale(data,&y);
  294. // if(count%data->nSaves==0){
  295. // printTimescaleOutput(tNow,&y, data->timescaleOutput,data);
  296. // //printSpaceTimeOutput(tNow, &y, data->output, data);
  297. // //if(data->writeRates){
  298. // // printSpaceTimeRates(tNow, ydot, data);
  299. // //}
  300. // }
  301. /*Print global variables only if time-step is of high order!*/
  302. if(data->nTimeSteps==0){
  303. data->flamePosition[0]=0.0e0;
  304. data->flamePosition[1]=0.0e0;
  305. data->flameTime[0]=tNow;
  306. data->flameTime[1]=tNow;
  307. }
  308. ier = IDAGetNumErrTestFails(mem, &netf);
  309. ier = IDAGetNumNonlinSolvConvFails(mem, &ncfn);
  310. ier = IDADlsGetNumJacEvals(mem, &njevals);
  311. ier = IDADlsGetNumResEvals(mem, &nrevals);
  312. printf("etf = %ld ,"
  313. "nlsf= %ld ,"
  314. "J=%ld ,"
  315. "R=%ld ,"
  316. "o=%d ,",netf, ncfn, njevals, nrevals, kcur);
  317. printf("Time=%15.6e s,",tNow);
  318. printf("dt=%15.6e s,",dt);
  319. printf("frac: %15.6e\n",dxRatio);
  320. count++;
  321. data->nTimeSteps=count;
  322. }
  323. }
  324. SUNLinSolFree(LS);
  325. SUNMatDestroy(A);
  326. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  327. return(0);
  328. }
  329. void freeAtLast(void* mem,
  330. N_Vector *y,
  331. N_Vector *ydot,
  332. N_Vector *res,
  333. N_Vector *id,
  334. N_Vector *atolv,
  335. N_Vector *constraints,UserData data){
  336. IDAFree(&mem);
  337. freeSolution(y,ydot,res,id,atolv,constraints);
  338. freeUserData(data);
  339. }
  340. static int check_flag(void *flagvalue, const char *funcname, int opt)
  341. {
  342. int *errflag;
  343. /* Check if SUNDIALS function returned NULL pointer - no memory allocated */
  344. if (opt == 0 && flagvalue == NULL) {
  345. fprintf(stderr,
  346. "\nSUNDIALS_ERROR: %s() failed - returned NULL pointer\n\n",
  347. funcname);
  348. return(1);
  349. }
  350. /* Check if flag < 0 */
  351. else if (opt == 1) {
  352. errflag = (int *) flagvalue;
  353. if (*errflag < 0) {
  354. fprintf(stderr,
  355. "\nSUNDIALS_ERROR: %s() failed with flag = %d\n\n",
  356. funcname, *errflag);
  357. return(1);
  358. }
  359. }
  360. /* Check if function returned NULL pointer - no memory allocated */
  361. else if (opt == 2 && flagvalue == NULL) {
  362. fprintf(stderr,
  363. "\nMEMORY_ERROR: %s() failed - returned NULL pointer\n\n",
  364. funcname);
  365. return(1);
  366. }
  367. return(0);
  368. }