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.

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