Actual source code: allen_cahn.c

petsc-3.7.2 2016-06-05
Report Typos and Errors
  1: static char help[] ="Solves the time dependent Allen-Cahn equation with IMEX methods";

  3: /*
  4:  * allen_cahn.c
  5:  *
  6:  *  Created on: Jun 8, 2012
  7:  *      Author: Hong Zhang
  8:  */

 10: #include <petscts.h>

 12: /*
 13:  * application context
 14:  */
 15: typedef struct {
 16:   PetscReal   param;        /* parameter */
 17:   PetscReal   xleft,xright;  /* range in x-direction */
 18:   PetscInt    mx;           /* Discretization in x-direction */
 19: }AppCtx;


 22: static PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
 23: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 24: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *ctx);
 25: static PetscErrorCode FormInitialSolution(TS,Vec,void*);


 30: int main(int argc, char **argv)
 31: {
 32:   TS                ts;
 33:   Vec               x; /*solution vector*/
 34:   Mat               A; /*Jacobian*/
 35:   PetscInt          steps,maxsteps,mx;
 36:   PetscErrorCode    ierr;
 37:   PetscReal         ftime;
 38:   AppCtx      user;       /* user-defined work context */

 40:   PetscInitialize(&argc,&argv,NULL,help);

 42:   /* Initialize user application context */
 43:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Allen-Cahn equation","");
 44:   user.param = 9e-4;
 45:   user.xleft = -1.;
 46:   user.xright = 2.;
 47:   user.mx = 400;
 48:   PetscOptionsReal("-eps","parameter","",user.param,&user.param,NULL);
 49:   PetscOptionsEnd();

 51:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 52:    Set runtime options
 53:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 54:   /*
 55:    * PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
 56:    */
 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:    Create necessary matrix and vectors, solve same ODE on every process
 59:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   MatCreate(PETSC_COMM_WORLD,&A);
 61:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,user.mx,user.mx);
 62:   MatSetFromOptions(A);
 63:   MatSetUp(A);
 64:   MatCreateVecs(A,&x,NULL);

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:    Create time stepping solver context
 68:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 69:   TSCreate(PETSC_COMM_WORLD,&ts);
 70:   TSSetType(ts,TSEIMEX);
 71:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
 72:   TSSetIFunction(ts,NULL,FormIFunction,&user);
 73:   TSSetIJacobian(ts,A,A,FormIJacobian,&user);
 74:   ftime = 142;
 75:   maxsteps = 100000;
 76:   TSSetDuration(ts,maxsteps,ftime);
 77:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:    Set initial conditions
 81:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 82:   FormInitialSolution(ts,x,&user);
 83:   TSSetSolution(ts,x);
 84:   VecGetSize(x,&mx);

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:    Set runtime options
 88:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 89:   TSSetFromOptions(ts);

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:    Solve nonlinear system
 93:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 94:   TSSolve(ts,x);
 95:   TSGetTime(ts,&ftime);
 96:   TSGetTimeStepNumber(ts,&steps);
 97:   PetscPrintf(PETSC_COMM_WORLD,"eps %g, steps %D, ftime %g\n",(double)user.param,steps,(double)ftime);
 98:   /*   VecView(x,PETSC_VIEWER_STDOUT_WORLD);*/

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:    Free work space.
102:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   MatDestroy(&A);
104:   VecDestroy(&x);
105:   TSDestroy(&ts);
106:   PetscFinalize();
107:   return(0);
108: }


113: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
114: {
115:   PetscErrorCode    ierr;
116:   AppCtx            *user = (AppCtx*)ptr;
117:   PetscScalar       *f;
118:   const PetscScalar *x;
119:   PetscInt          i,mx;
120:   PetscReal         hx,eps;

123:   mx = user->mx;
124:   eps = user->param;
125:   hx = (user->xright-user->xleft)/(mx-1);
126:   VecGetArrayRead(X,&x);
127:   VecGetArray(F,&f);
128:   f[0] = 2.*eps*(x[1]-x[0])/(hx*hx); /*boundary*/
129:   for(i=1;i<mx-1;i++) {
130:     f[i]= eps*(x[i+1]-2.*x[i]+x[i-1])/(hx*hx);
131:   }
132:   f[mx-1] = 2.*eps*(x[mx-2]- x[mx-1])/(hx*hx); /*boundary*/
133:   VecRestoreArrayRead(X,&x);
134:   VecRestoreArray(F,&f);
135:   return(0);
136: }


141: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
142: {
143:   PetscErrorCode    ierr;
144:   AppCtx            *user = (AppCtx*)ptr;
145:   PetscScalar       *f;
146:   const PetscScalar *x,*xdot;
147:   PetscInt          i,mx;

150:   mx = user->mx;
151:   VecGetArrayRead(X,&x);
152:   VecGetArrayRead(Xdot,&xdot);
153:   VecGetArray(F,&f);

155:   for(i=0;i<mx;i++) {
156:     f[i]= xdot[i] - x[i]*(1.-x[i]*x[i]);
157:   }

159:   VecRestoreArrayRead(X,&x);
160:   VecRestoreArrayRead(Xdot,&xdot);
161:   VecRestoreArray(F,&f);
162:   return(0);
163: }


168: static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U, Vec Udot, PetscReal a, Mat J,Mat Jpre,void *ctx)
169: {
170:   PetscErrorCode    ierr;
171:   AppCtx            *user = (AppCtx *)ctx;
172:   PetscScalar       v;
173:   const PetscScalar *x;
174:   PetscInt          i,col;

177:   VecGetArrayRead(U,&x);
178:   for(i=0; i < user->mx; i++) {
179:     v = a - 1. + 3.*x[i]*x[i];
180:     col = i;
181:     MatSetValues(J,1,&i,1,&col,&v,INSERT_VALUES);
182:   }
183:   VecRestoreArrayRead(U,&x);

185:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
186:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
187:   if (J != Jpre) {
188:     MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
189:     MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
190:   }
191:   /*  MatView(J,PETSC_VIEWER_STDOUT_WORLD);*/
192:   return(0);
193: }


198: static PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ctx)
199: {
200:   AppCtx *user = (AppCtx*)ctx;
201:   PetscInt       i;
202:   PetscScalar    *x;
203:   PetscReal      hx,x_map;

207:   hx = (user->xright-user->xleft)/(PetscReal)(user->mx-1);
208:   VecGetArray(U,&x);
209:   for(i=0;i<user->mx;i++) {
210:     x_map = user->xleft + i*hx;
211:     if(x_map >= 0.7065) {
212:       x[i] = tanh((x_map-0.8)/(2.*PetscSqrtReal(user->param)));
213:     } else if(x_map >= 0.4865) {
214:       x[i] = tanh((0.613-x_map)/(2.*PetscSqrtReal(user->param)));
215:     } else if(x_map >= 0.28) {
216:       x[i] = tanh((x_map-0.36)/(2.*PetscSqrtReal(user->param)));
217:     } else if(x_map >= -0.7) {
218:       x[i] = tanh((0.2-x_map)/(2.*PetscSqrtReal(user->param)));
219:     } else if(x_map >= -1) {
220:       x[i] = tanh((x_map+0.9)/(2.*PetscSqrtReal(user->param)));
221:     }
222:   }
223:   VecRestoreArray(U,&x);
224:   return(0);
225: }