Droplet Lagrangian Transient One-dimensional Reacting Code Implementation of both liquid and gas phase governing equations.
Nie możesz wybrać więcej, niż 25 tematów Tematy muszą się zaczynać od litery lub cyfry, mogą zawierać myślniki ('-') i mogą mieć do 35 znaków.

384 wiersze
11KB

  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. if(data==NULL){
  34. printf("check input file!\n");
  35. freeUserData(data);
  36. return(-1);
  37. }
  38. // Allocate solution variables
  39. long int ier,mu,ml,count,netf,ncfn,njevals,nrevals;
  40. realtype tNow,*atolvdata,*constraintsdata,finalTime,dtMax,tolsfac;
  41. N_Vector y,ydot,id,res,atolv,constraints;
  42. y=ydot=id=res=atolv=constraints=NULL;
  43. ier=allocateSolution(data->neq,data->nThreads,&y,&ydot,&id,&res,&atolv,&constraints);
  44. ier=setAlgebraicVariables(&id,data);
  45. ier=setInitialCondition(&y,&ydot,data);
  46. if(ier==-1)return(-1);
  47. tNow=data->tNow;
  48. finalTime=data->finalTime;
  49. //TODO: dtMax should be a user input
  50. dtMax = 1e-4;
  51. double* ydata;
  52. double* ydotdata;
  53. ydata = N_VGetArrayPointer_OpenMP(y);
  54. ydotdata = N_VGetArrayPointer_OpenMP(ydot);
  55. ////////// DEBUG ///////////////////
  56. //double* resdata;
  57. //double* iddata;
  58. //resdata = N_VGetArrayPointer_OpenMP(res);
  59. //iddata = N_VGetArrayPointer_OpenMP(id);
  60. ///////////////////////////////////
  61. void *mem;mem=NULL;mem = IDACreate();
  62. ier = IDASetUserData(mem, data);
  63. ier = IDASetId(mem, id);
  64. ier = IDAInit(mem, residue, tNow, y, ydot);
  65. // Atol array
  66. atolvdata = N_VGetArrayPointer_OpenMP(atolv);
  67. for (size_t i = 1; i <=data->npts; i++) {
  68. atolT(i) = data->temperatureTolerance;
  69. for (size_t k = 1; k <=data->nsp; k++) {
  70. if(k!=data->k_bath){
  71. atolY(i,k) = data->massFractionTolerance;
  72. }
  73. else{
  74. atolY(i,k) = data->bathGasTolerance;
  75. }
  76. }
  77. atolR(i) = data->radiusTolerance;
  78. atolP(i) = data->pressureTolerance;
  79. atolMdot(i) = data->MdotTolerance;
  80. }
  81. ier = IDASVtolerances(mem, data->relativeTolerance, atolv);
  82. mu = 2*data->nvar; ml = mu;
  83. SUNMatrix A; A=NULL;
  84. A=SUNBandMatrix(data->neq,mu,ml,mu+ml);
  85. SUNLinearSolver LS; LS=NULL;
  86. //LS=SUNBandLinearSolver(y,A);
  87. LS=SUNLapackBand(y,A);
  88. ier=IDADlsSetLinearSolver(mem,LS,A);
  89. //ier = IDABand(mem, data->neq, mu, ml);
  90. constraintsdata = N_VGetArrayPointer_OpenMP(constraints);
  91. if(data->setConstraints){
  92. for (size_t i = 1; i <=data->npts; i++) {
  93. for (size_t k = 1; k <=data->nsp; k++) {
  94. constraintsY(i,k) = ONE;
  95. }
  96. }
  97. ier=IDASetConstraints(mem, constraints);
  98. }
  99. if(!data->quasiSteady){
  100. constraintsR(1) = ONE;
  101. ier=IDASetConstraints(mem, constraints);
  102. }
  103. ier = IDASetSuppressAlg(mem, data->suppressAlg);
  104. if(check_flag(&ier, "IDASetSuppressAlg", 1)) return(1);
  105. //ier= IDASetMaxNumStepsIC(mem, 1);
  106. //ier= IDASetMaxNumJacsIC(mem,8);
  107. //ier= IDASetMaxNumItersIC(mem,100);
  108. //ier= IDASetMaxBacksIC(mem,2000);
  109. //ier = IDASetLineSearchOffIC(mem,SUNTRUE);
  110. //////// DEBUG ///////////
  111. //if(data->dryRun){
  112. // ier = residue(tNow,y, ydot, res, data);
  113. // for(size_t k=0; k < data->nvar*data->npts; k++){
  114. // printf("%i: %15.6e\n",k,resdata[k]);
  115. // }
  116. // for(size_t k=0; k < data->nvar*data->npts; k++){
  117. // if(iddata[k] == 1){
  118. // ydotdata[k] = -resdata[k];
  119. // }
  120. // }
  121. // ier = residue(tNow,y, ydot, res, data);
  122. // for(size_t k=0; k < data->nvar*data->npts; k++){
  123. // printf("%i: %15.6e\n",k,resdata[k]);
  124. // }
  125. //}
  126. //for(size_t k=0; k < data->neq; k++){
  127. // if(iddata[k] == 1){
  128. // ydotdata[k] = -resdata[k];
  129. // }
  130. //}
  131. //////////////////////////
  132. if(!data->dryRun){
  133. printf("Calculating Initial Conditions:\n");
  134. printf("Cross your fingers...\n");
  135. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, 1e-05*finalTime);
  136. //If at first IDACalcIC doesn't succeed, try, try, try again:
  137. for (int i = 0; i < 10; i++) {
  138. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, (1e-01+pow(10,i)*finalTime));
  139. if(ier==0){
  140. break;
  141. }
  142. }
  143. //...or else cry again :(
  144. if(check_flag(&ier, "IDACalcIC", 1)){
  145. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  146. return(-1);
  147. }else{
  148. printf("Initial (Consistent) Conditions Calculated!\n");
  149. }
  150. ier = IDASetInitStep(mem,1e-12);
  151. }
  152. printSpaceTimeHeader(data);
  153. printGlobalHeader(data);
  154. printSpaceTimeOutput(tNow, &y, data->output, data);
  155. printSpaceTimeOutput(tNow, &y, data->gridOutput, data);
  156. if(!data->dryRun){
  157. count=0;
  158. double dt=1e-08;
  159. double t1=0.0e0;
  160. double xOld=0.0e0;
  161. double x=0.0e0;
  162. double dx=0.0e0;
  163. double dxMin=1.0e0;
  164. double dxRatio=dx/dxMin;
  165. int move=0;
  166. int kcur=0;
  167. if(data->adaptiveGrid){
  168. dxMin=data->grid->leastMove;
  169. //xOld=maxCurvPosition(ydata, data->nt, data->nvar,
  170. // data->grid->x, data->npts);
  171. xOld=isothermPosition(ydata, data->isotherm, data->nt,
  172. data->nvar, data->grid->x, data->npts);
  173. }
  174. while (tNow<=finalTime && R(1)>100e-9) {
  175. t1=tNow;
  176. if(data->quasiSteady){
  177. ier = IDASolve(mem, finalTime, &tNow, y, ydot, IDA_ONE_STEP);
  178. }else{
  179. /*This prevents the solver from taking a step so large that
  180. *the droplet radius becomes negative.
  181. *TODO:Try adding the constraint that the droplet radius must
  182. * be a positive number*/
  183. ier = IDASolve(mem, tNow+dtMax, &tNow, y, ydot, IDA_ONE_STEP);
  184. //ier = IDASolve(mem, tNow+dtMax, &tNow, y, ydot, IDA_NORMAL);
  185. }
  186. if(check_flag(&ier, "IDASolve", 1)){
  187. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  188. return(-1);
  189. }
  190. dt=tNow-t1;
  191. ier = IDAGetCurrentOrder(mem, &kcur);
  192. if(data->adaptiveGrid==1 && data->moveGrid==1){
  193. x=maxCurvPosition(ydata, data->nt, data->nvar,
  194. data->grid->x, data->npts);
  195. //x=isothermPosition(ydata, data->isotherm, data->nt,
  196. // data->nvar, data->grid->x, data->npts);
  197. //x=maxGradPosition(ydata, data->nt, data->nvar,
  198. // data->grid->x, data->npts);
  199. dx=x-xOld;
  200. if(dx*dxMin>0.0e0){
  201. move=1;
  202. }else{
  203. move=-1;
  204. }
  205. //if(fabs(dx)>=dxMin && x+(double)(move)*0.5e0*dxMin<=1.0e0){
  206. dxRatio=fabs(dx/dxMin);
  207. if(dxRatio>=1.0e0 && dxRatio<=2.0e0){
  208. printf("Regridding!\n");
  209. data->regrid=1;
  210. printSpaceTimeOutput(tNow, &y, data->gridOutput, data);
  211. ier=reGrid(data->grid, x+(double)(move)*0.5e0*dxMin);
  212. if(ier==-1){
  213. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  214. return(-1);
  215. }
  216. updateSolution(ydata, ydotdata, data->nvar,
  217. data->grid->xOld,data->grid->x,data->npts);
  218. storeGrid(data->grid->x,data->grid->xOld,data->npts);
  219. xOld=x;
  220. printf("Regrid Complete! Restarting Problem at %15.6e s\n",tNow);
  221. ier = IDAReInit(mem, tNow, y, ydot);
  222. if(check_flag(&ier, "IDAReInit", 1)){
  223. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  224. return(-1);
  225. }
  226. ier = IDASetInitStep(mem,1e-01*dt);
  227. printf("Reinitialized! Calculating Initial Conditions:\n");
  228. printf("Cross your fingers...\n");
  229. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e-01*dt);
  230. if(check_flag(&ier, "IDACalcIC", 1)){
  231. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+1e+01*dt);
  232. }
  233. //Every once in a while, for reasons
  234. //that befuddle this mathematically
  235. //lesser author, IDACalcIC fails. Here,
  236. //I desperately try to make it converge
  237. //again by sampling increasingly larger
  238. //time-steps:
  239. for (int i = 0; i < 10; i++) {
  240. ier = IDACalcIC(mem, IDA_YA_YDP_INIT, tNow+(1e-01+pow(10,i)*dt));
  241. if(ier==0){
  242. break;
  243. }
  244. }
  245. //Failure :( Back to the drawing board:
  246. if(check_flag(&ier, "IDACalcIC", 1)){
  247. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  248. return(-1);
  249. }
  250. printf("Initial (Consistent) Conditions Calculated!\n");
  251. printf("------------------------------------------\n\n");
  252. if(data->writeEveryRegrid){
  253. printSpaceTimeOutput(tNow, &y, data->output, data);
  254. FILE* fp;
  255. fp=fopen("restart.bin","w");
  256. writeRestart(tNow, &y, &ydot, fp, data);
  257. fclose(fp);
  258. }
  259. }
  260. }
  261. if(count%data->nSaves==0 && !data->writeEveryRegrid){
  262. printSpaceTimeOutput(tNow, &y, data->output, data);
  263. //if(data->writeRates){
  264. // printSpaceTimeRates(tNow, ydot, data);
  265. //}
  266. }
  267. /*Print global variables only if time-step is of high order!*/
  268. if(data->nTimeSteps==0){
  269. data->flamePosition[0]=0.0e0;
  270. data->flamePosition[1]=0.0e0;
  271. data->flameTime[0]=tNow;
  272. data->flameTime[1]=tNow;
  273. }
  274. ier = IDAGetNumErrTestFails(mem, &netf);
  275. ier = IDAGetNumNonlinSolvConvFails(mem, &ncfn);
  276. ier = IDADlsGetNumJacEvals(mem, &njevals);
  277. ier = IDADlsGetNumResEvals(mem, &nrevals);
  278. printf("etf = %ld ,"
  279. "nlsf= %ld ,"
  280. "J=%ld ,"
  281. "R=%ld ,"
  282. "o=%d ,",netf, ncfn, njevals, nrevals, kcur);
  283. printf("Time=%15.6e s,",tNow);
  284. printf("frac: %15.6e\n",dxRatio);
  285. count++;
  286. data->nTimeSteps=count;
  287. }
  288. }
  289. SUNLinSolFree(LS);
  290. SUNMatDestroy(A);
  291. freeAtLast(mem,&y,&ydot,&res,&id,&atolv,&constraints,data);
  292. return(0);
  293. }
  294. void freeAtLast(void* mem,
  295. N_Vector *y,
  296. N_Vector *ydot,
  297. N_Vector *res,
  298. N_Vector *id,
  299. N_Vector *atolv,
  300. N_Vector *constraints,UserData data){
  301. IDAFree(&mem);
  302. freeSolution(y,ydot,res,id,atolv,constraints);
  303. freeUserData(data);
  304. }
  305. static int check_flag(void *flagvalue, const char *funcname, int opt)
  306. {
  307. int *errflag;
  308. /* Check if SUNDIALS function returned NULL pointer - no memory allocated */
  309. if (opt == 0 && flagvalue == NULL) {
  310. fprintf(stderr,
  311. "\nSUNDIALS_ERROR: %s() failed - returned NULL pointer\n\n",
  312. funcname);
  313. return(1);
  314. }
  315. /* Check if flag < 0 */
  316. else if (opt == 1) {
  317. errflag = (int *) flagvalue;
  318. if (*errflag < 0) {
  319. fprintf(stderr,
  320. "\nSUNDIALS_ERROR: %s() failed with flag = %d\n\n",
  321. funcname, *errflag);
  322. return(1);
  323. }
  324. }
  325. /* Check if function returned NULL pointer - no memory allocated */
  326. else if (opt == 2 && flagvalue == NULL) {
  327. fprintf(stderr,
  328. "\nMEMORY_ERROR: %s() failed - returned NULL pointer\n\n",
  329. funcname);
  330. return(1);
  331. }
  332. return(0);
  333. }