Actual source code: allen_cahn.c
petsc-3.6.4 2016-04-12
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,"-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);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Set initial conditions
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: FormInitialSolution(ts,x,&user);
82: TSSetSolution(ts,x);
83: VecGetSize(x,&mx);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Set runtime options
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: TSSetFromOptions(ts);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Solve nonlinear system
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: TSSolve(ts,x);
94: TSGetTime(ts,&ftime);
95: TSGetTimeStepNumber(ts,&steps);
96: PetscPrintf(PETSC_COMM_WORLD,"eps %g, steps %D, ftime %g\n",(double)user.param,steps,(double)ftime);
97: /* VecView(x,PETSC_VIEWER_STDOUT_WORLD);*/
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Free work space.
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: MatDestroy(&A);
103: VecDestroy(&x);
104: TSDestroy(&ts);
105: PetscFinalize();
106: return(0);
107: }
112: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
113: {
114: PetscErrorCode ierr;
115: AppCtx *user = (AppCtx*)ptr;
116: PetscScalar *f;
117: const PetscScalar *x;
118: PetscInt i,mx;
119: PetscReal hx,eps;
122: mx = user->mx;
123: eps = user->param;
124: hx = (user->xright-user->xleft)/(mx-1);
125: VecGetArrayRead(X,&x);
126: VecGetArray(F,&f);
127: f[0] = 2.*eps*(x[1]-x[0])/(hx*hx); /*boundary*/
128: for(i=1;i<mx-1;i++) {
129: f[i]= eps*(x[i+1]-2.*x[i]+x[i-1])/(hx*hx);
130: }
131: f[mx-1] = 2.*eps*(x[mx-2]- x[mx-1])/(hx*hx); /*boundary*/
132: VecRestoreArrayRead(X,&x);
133: VecRestoreArray(F,&f);
134: return(0);
135: }
140: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
141: {
142: PetscErrorCode ierr;
143: AppCtx *user = (AppCtx*)ptr;
144: PetscScalar *f;
145: const PetscScalar *x,*xdot;
146: PetscInt i,mx;
149: mx = user->mx;
150: VecGetArrayRead(X,&x);
151: VecGetArrayRead(Xdot,&xdot);
152: VecGetArray(F,&f);
154: for(i=0;i<mx;i++) {
155: f[i]= xdot[i] - x[i]*(1.-x[i]*x[i]);
156: }
158: VecRestoreArrayRead(X,&x);
159: VecRestoreArrayRead(Xdot,&xdot);
160: VecRestoreArray(F,&f);
161: return(0);
162: }
167: static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U, Vec Udot, PetscReal a, Mat J,Mat Jpre,void *ctx)
168: {
169: PetscErrorCode ierr;
170: AppCtx *user = (AppCtx *)ctx;
171: PetscScalar v;
172: const PetscScalar *x;
173: PetscInt i,col;
176: VecGetArrayRead(U,&x);
177: for(i=0; i < user->mx; i++) {
178: v = a - 1. + 3.*x[i]*x[i];
179: col = i;
180: MatSetValues(J,1,&i,1,&col,&v,INSERT_VALUES);
181: }
182: VecRestoreArrayRead(U,&x);
184: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
185: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
186: if (J != Jpre) {
187: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
188: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
189: }
190: /* MatView(J,PETSC_VIEWER_STDOUT_WORLD);*/
191: return(0);
192: }
197: static PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ctx)
198: {
199: AppCtx *user = (AppCtx*)ctx;
200: PetscInt i;
201: PetscScalar *x;
202: PetscReal hx,x_map;
206: hx = (user->xright-user->xleft)/(PetscReal)(user->mx-1);
207: VecGetArray(U,&x);
208: for(i=0;i<user->mx;i++) {
209: x_map = user->xleft + i*hx;
210: if(x_map >= 0.7065) {
211: x[i] = tanh((x_map-0.8)/(2.*PetscSqrtReal(user->param)));
212: } else if(x_map >= 0.4865) {
213: x[i] = tanh((0.613-x_map)/(2.*PetscSqrtReal(user->param)));
214: } else if(x_map >= 0.28) {
215: x[i] = tanh((x_map-0.36)/(2.*PetscSqrtReal(user->param)));
216: } else if(x_map >= -0.7) {
217: x[i] = tanh((0.2-x_map)/(2.*PetscSqrtReal(user->param)));
218: } else if(x_map >= -1) {
219: x[i] = tanh((x_map+0.9)/(2.*PetscSqrtReal(user->param)));
220: }
221: }
222: VecRestoreArray(U,&x);
223: return(0);
224: }