PageRenderTime 40ms CodeModel.GetById 12ms RepoModel.GetById 0ms app.codeStats 0ms

/src/ts/examples/tutorials/ex22.c

https://bitbucket.org/alexei-matveev/petsc-debian-pkg
C | 281 lines | 188 code | 36 blank | 57 comment | 20 complexity | 59ed5f1fc188aaf0354f63e9508b00b6 MD5 | raw file
  1. static const char help[] = "Time-dependent advection-reaction PDE in 1d, demonstrates IMEX methods.\n";
  2. /*
  3. u_t + a1*u_x = -k1*u + k2*v + s1
  4. v_t + a2*v_x = k1*u - k2*v + s2
  5. 0 < x < 1;
  6. a1 = 1, k1 = 10^6, s1 = 0,
  7. a2 = 0, k2 = 2*k1, s2 = 1
  8. Initial conditions:
  9. u(x,0) = 1 + s2*x
  10. v(x,0) = k0/k1*u(x,0) + s1/k1
  11. Upstream boundary conditions:
  12. u(0,t) = 1-sin(12*t)^4
  13. */
  14. #include <petscdmda.h>
  15. #include <petscts.h>
  16. typedef PetscScalar Field[2];
  17. typedef struct _User *User;
  18. struct _User {
  19. PetscReal a[2]; /* Advection speeds */
  20. PetscReal k[2]; /* Reaction coefficients */
  21. PetscReal s[2]; /* Source terms */
  22. };
  23. static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
  24. static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
  25. static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
  26. static PetscErrorCode FormInitialSolution(TS,Vec,void*);
  27. #undef __FUNCT__
  28. #define __FUNCT__ "main"
  29. int main(int argc,char **argv)
  30. {
  31. TS ts; /* nonlinear solver */
  32. Vec X; /* solution, residual vectors */
  33. Mat J; /* Jacobian matrix */
  34. PetscInt steps,maxsteps,mx;
  35. PetscErrorCode ierr;
  36. DM da;
  37. PetscReal ftime,dt;
  38. struct _User user; /* user-defined work context */
  39. TSConvergedReason reason;
  40. PetscInitialize(&argc,&argv,(char *)0,help);
  41. /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  42. Create distributed array (DMDA) to manage parallel grid and vectors
  43. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
  44. ierr = DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-11,2,2,PETSC_NULL,&da);CHKERRQ(ierr);
  45. /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  46. Extract global vectors from DMDA;
  47. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
  48. ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr);
  49. /* Initialize user application context */
  50. ierr = PetscOptionsBegin(PETSC_COMM_WORLD,PETSC_NULL,"Advection-reaction options",""); {
  51. user.a[0] = 1; ierr = PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],PETSC_NULL);CHKERRQ(ierr);
  52. user.a[1] = 0; ierr = PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],PETSC_NULL);CHKERRQ(ierr);
  53. user.k[0] = 1e6; ierr = PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],PETSC_NULL);CHKERRQ(ierr);
  54. user.k[1] = 2*user.k[0]; ierr = PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],PETSC_NULL);CHKERRQ(ierr);
  55. user.s[0] = 0; ierr = PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],PETSC_NULL);CHKERRQ(ierr);
  56. user.s[1] = 1; ierr = PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],PETSC_NULL);CHKERRQ(ierr);
  57. } ierr = PetscOptionsEnd();CHKERRQ(ierr);
  58. /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  59. Create timestepping solver context
  60. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
  61. ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
  62. ierr = TSSetDM(ts,da);CHKERRQ(ierr);
  63. ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
  64. ierr = TSSetRHSFunction(ts,PETSC_NULL,FormRHSFunction,&user);CHKERRQ(ierr);
  65. ierr = TSSetIFunction(ts,PETSC_NULL,FormIFunction,&user);CHKERRQ(ierr);
  66. ierr = DMCreateMatrix(da,MATAIJ,&J);CHKERRQ(ierr);
  67. ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr);
  68. ftime = 1.0;
  69. maxsteps = 10000;
  70. ierr = TSSetDuration(ts,maxsteps,ftime);CHKERRQ(ierr);
  71. /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  72. Set initial conditions
  73. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
  74. ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr);
  75. ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
  76. ierr = VecGetSize(X,&mx);CHKERRQ(ierr);
  77. dt = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
  78. ierr = TSSetInitialTimeStep(ts,0.0,dt);CHKERRQ(ierr);
  79. /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  80. Set runtime options
  81. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
  82. ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
  83. /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  84. Solve nonlinear system
  85. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
  86. ierr = TSSolve(ts,X,&ftime);CHKERRQ(ierr);
  87. ierr = TSGetTimeStepNumber(ts,&steps);CHKERRQ(ierr);
  88. ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
  89. ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %G after %D steps\n",TSConvergedReasons[reason],ftime,steps);CHKERRQ(ierr);
  90. /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  91. Free work space.
  92. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
  93. ierr = MatDestroy(&J);CHKERRQ(ierr);
  94. ierr = VecDestroy(&X);CHKERRQ(ierr);
  95. ierr = TSDestroy(&ts);CHKERRQ(ierr);
  96. ierr = DMDestroy(&da);CHKERRQ(ierr);
  97. ierr = PetscFinalize();
  98. return 0;
  99. }
  100. #undef __FUNCT__
  101. #define __FUNCT__ "FormIFunction"
  102. static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
  103. {
  104. User user = (User)ptr;
  105. DM da;
  106. DMDALocalInfo info;
  107. PetscInt i;
  108. Field *x,*xdot,*f;
  109. PetscErrorCode ierr;
  110. PetscFunctionBegin;
  111. ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
  112. ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
  113. /* Get pointers to vector data */
  114. ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
  115. ierr = DMDAVecGetArray(da,Xdot,&xdot);CHKERRQ(ierr);
  116. ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
  117. /* Compute function over the locally owned part of the grid */
  118. for (i=info.xs; i<info.xs+info.xm; i++) {
  119. f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0];
  120. f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1];
  121. }
  122. /* Restore vectors */
  123. ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
  124. ierr = DMDAVecRestoreArray(da,Xdot,&xdot);CHKERRQ(ierr);
  125. ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
  126. PetscFunctionReturn(0);
  127. }
  128. #undef __FUNCT__
  129. #define __FUNCT__ "FormRHSFunction"
  130. static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
  131. {
  132. User user = (User)ptr;
  133. DM da;
  134. Vec Xloc;
  135. DMDALocalInfo info;
  136. PetscInt i,j;
  137. PetscReal hx;
  138. Field *x,*f;
  139. PetscErrorCode ierr;
  140. PetscFunctionBegin;
  141. ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
  142. ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
  143. hx = 1.0/(PetscReal)info.mx;
  144. /*
  145. Scatter ghost points to local vector,using the 2-step process
  146. DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
  147. By placing code between these two statements, computations can be
  148. done while messages are in transition.
  149. */
  150. ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr);
  151. ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
  152. ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
  153. /* Get pointers to vector data */
  154. ierr = DMDAVecGetArray(da,Xloc,&x);CHKERRQ(ierr);
  155. ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
  156. /* Compute function over the locally owned part of the grid */
  157. for (i=info.xs; i<info.xs+info.xm; i++) {
  158. const PetscReal *a = user->a;
  159. PetscReal u0t[2] = {1. - PetscPowScalar(sin(12*t),4.),0};
  160. for (j=0; j<2; j++) {
  161. if (i == 0) {
  162. f[i][j] = a[j]/hx*(1./3*u0t[j] + 0.5*x[i][j] - x[i+1][j] + 1./6*x[i+2][j]);
  163. } else if (i == 1) {
  164. f[i][j] = a[j]/hx*(-1./12*u0t[j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
  165. } else if (i == info.mx-2) {
  166. f[i][j] = a[j]/hx*(-1./6*x[i-2][j] + x[i-1][j] - 0.5*x[i][j] - 1./3*x[i+1][j]);
  167. } else if (i == info.mx-1) {
  168. f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]);
  169. } else {
  170. f[i][j] = a[j]/hx*(-1./12*x[i-2][j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
  171. }
  172. }
  173. }
  174. /* Restore vectors */
  175. ierr = DMDAVecRestoreArray(da,Xloc,&x);CHKERRQ(ierr);
  176. ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
  177. ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr);
  178. PetscFunctionReturn(0);
  179. }
  180. /* --------------------------------------------------------------------- */
  181. /*
  182. IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
  183. */
  184. #undef __FUNCT__
  185. #define __FUNCT__ "FormIJacobian"
  186. PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ptr)
  187. {
  188. User user = (User)ptr;
  189. PetscErrorCode ierr;
  190. DMDALocalInfo info;
  191. PetscInt i;
  192. DM da;
  193. Field *x,*xdot;
  194. PetscFunctionBegin;
  195. ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
  196. ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
  197. /* Get pointers to vector data */
  198. ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
  199. ierr = DMDAVecGetArray(da,Xdot,&xdot);CHKERRQ(ierr);
  200. /* Compute function over the locally owned part of the grid */
  201. for (i=info.xs; i<info.xs+info.xm; i++) {
  202. const PetscReal *k = user->k;
  203. PetscScalar v[2][2] = {{a + k[0], -k[1]},
  204. {-k[0], a+k[1]}};
  205. ierr = MatSetValuesBlocked(*Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES);CHKERRQ(ierr);
  206. }
  207. /* Restore vectors */
  208. ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
  209. ierr = DMDAVecRestoreArray(da,Xdot,&xdot);CHKERRQ(ierr);
  210. ierr = MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  211. ierr = MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  212. if (*J != *Jpre) {
  213. ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  214. ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
  215. }
  216. PetscFunctionReturn(0);
  217. }
  218. #undef __FUNCT__
  219. #define __FUNCT__ "FormInitialSolution"
  220. PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
  221. {
  222. User user = (User)ctx;
  223. DM da;
  224. PetscInt i;
  225. DMDALocalInfo info;
  226. Field *x;
  227. PetscReal hx;
  228. PetscErrorCode ierr;
  229. PetscFunctionBegin;
  230. ierr = TSGetDM(ts,&da);
  231. ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
  232. hx = 1.0/(PetscReal)info.mx;
  233. /* Get pointers to vector data */
  234. ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
  235. /* Compute function over the locally owned part of the grid */
  236. for (i=info.xs; i<info.xs+info.xm; i++) {
  237. PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0;
  238. x[i][0] = 1 + user->s[1]*r;
  239. x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik;
  240. }
  241. ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
  242. PetscFunctionReturn(0);
  243. }