xref: /petsc/src/ts/tutorials/eimex/allen_cahn.c (revision 2da392cc7c10228af19ad9843ce5155178acb644)
1 static char help[] ="Solves the time dependent Allen-Cahn equation with IMEX methods";
2 
3 /*
4  * allen_cahn.c
5  *
6  *  Created on: Jun 8, 2012
7  *      Author: Hong Zhang
8  */
9 
10 #include <petscts.h>
11 
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;
20 
21 
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*);
26 
27 
28 int main(int argc, char **argv)
29 {
30   TS                ts;
31   Vec               x; /*solution vector*/
32   Mat               A; /*Jacobian*/
33   PetscInt          steps,mx;
34   PetscErrorCode    ierr;
35   PetscReal         ftime;
36   AppCtx            user;       /* user-defined work context */
37 
38   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
39 
40   /* Initialize user application context */
41   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Allen-Cahn equation","");
42   user.param = 9e-4;
43   user.xleft = -1.;
44   user.xright = 2.;
45   user.mx = 400;
46   ierr = PetscOptionsReal("-eps","parameter","",user.param,&user.param,NULL);CHKERRQ(ierr);
47   ierr = PetscOptionsEnd();CHKERRQ(ierr);
48 
49   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50    Set runtime options
51    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52   /*
53    * ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
54    */
55   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56    Create necessary matrix and vectors, solve same ODE on every process
57    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
59   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,user.mx,user.mx);CHKERRQ(ierr);
60   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
61   ierr = MatSetUp(A);CHKERRQ(ierr);
62   ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr);
63 
64   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65    Create time stepping solver context
66    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
68   ierr = TSSetType(ts,TSEIMEX);CHKERRQ(ierr);
69   ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr);
70   ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr);
71   ierr = TSSetIJacobian(ts,A,A,FormIJacobian,&user);CHKERRQ(ierr);
72   ftime = 22;
73   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
74   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
75 
76   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77    Set initial conditions
78    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79   ierr = FormInitialSolution(ts,x,&user);CHKERRQ(ierr);
80   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
81   ierr = VecGetSize(x,&mx);CHKERRQ(ierr);
82 
83   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84    Set runtime options
85    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
87 
88   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89    Solve nonlinear system
90    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91   ierr = TSSolve(ts,x);CHKERRQ(ierr);
92   ierr = TSGetTime(ts,&ftime);CHKERRQ(ierr);
93   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
94   ierr = PetscPrintf(PETSC_COMM_WORLD,"eps %g, steps %D, ftime %g\n",(double)user.param,steps,(double)ftime);CHKERRQ(ierr);
95   /*   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/
96 
97   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98    Free work space.
99    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100   ierr = MatDestroy(&A);CHKERRQ(ierr);
101   ierr = VecDestroy(&x);CHKERRQ(ierr);
102   ierr = TSDestroy(&ts);CHKERRQ(ierr);
103   ierr = PetscFinalize();
104   return ierr;
105 }
106 
107 
108 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
109 {
110   PetscErrorCode    ierr;
111   AppCtx            *user = (AppCtx*)ptr;
112   PetscScalar       *f;
113   const PetscScalar *x;
114   PetscInt          i,mx;
115   PetscReal         hx,eps;
116 
117   PetscFunctionBegin;
118   mx = user->mx;
119   eps = user->param;
120   hx = (user->xright-user->xleft)/(mx-1);
121   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
122   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
123   f[0] = 2.*eps*(x[1]-x[0])/(hx*hx); /*boundary*/
124   for (i=1;i<mx-1;i++) {
125     f[i]= eps*(x[i+1]-2.*x[i]+x[i-1])/(hx*hx);
126   }
127   f[mx-1] = 2.*eps*(x[mx-2]- x[mx-1])/(hx*hx); /*boundary*/
128   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
129   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 
134 static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
135 {
136   PetscErrorCode    ierr;
137   AppCtx            *user = (AppCtx*)ptr;
138   PetscScalar       *f;
139   const PetscScalar *x,*xdot;
140   PetscInt          i,mx;
141 
142   PetscFunctionBegin;
143   mx = user->mx;
144   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
145   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
146   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
147 
148   for (i=0;i<mx;i++) {
149     f[i]= xdot[i] - x[i]*(1.-x[i]*x[i]);
150   }
151 
152   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
153   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
154   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
155   PetscFunctionReturn(0);
156 }
157 
158 
159 static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U, Vec Udot, PetscReal a, Mat J,Mat Jpre,void *ctx)
160 {
161   PetscErrorCode    ierr;
162   AppCtx            *user = (AppCtx *)ctx;
163   PetscScalar       v;
164   const PetscScalar *x;
165   PetscInt          i,col;
166 
167   PetscFunctionBegin;
168   ierr = VecGetArrayRead(U,&x);CHKERRQ(ierr);
169   for (i=0; i < user->mx; i++) {
170     v = a - 1. + 3.*x[i]*x[i];
171     col = i;
172     ierr = MatSetValues(J,1,&i,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr);
173   }
174   ierr = VecRestoreArrayRead(U,&x);CHKERRQ(ierr);
175 
176   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
178   if (J != Jpre) {
179     ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180     ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
181   }
182   /*  MatView(J,PETSC_VIEWER_STDOUT_WORLD);*/
183   PetscFunctionReturn(0);
184 }
185 
186 
187 static PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ctx)
188 {
189   AppCtx *user = (AppCtx*)ctx;
190   PetscInt       i;
191   PetscScalar    *x;
192   PetscReal      hx,x_map;
193   PetscErrorCode ierr;
194 
195   PetscFunctionBegin;
196   hx = (user->xright-user->xleft)/(PetscReal)(user->mx-1);
197   ierr = VecGetArray(U,&x);CHKERRQ(ierr);
198   for (i=0;i<user->mx;i++) {
199     x_map = user->xleft + i*hx;
200     if (x_map >= 0.7065) {
201       x[i] = PetscTanhReal((x_map-0.8)/(2.*PetscSqrtReal(user->param)));
202     } else if (x_map >= 0.4865) {
203       x[i] = PetscTanhReal((0.613-x_map)/(2.*PetscSqrtReal(user->param)));
204     } else if (x_map >= 0.28) {
205       x[i] = PetscTanhReal((x_map-0.36)/(2.*PetscSqrtReal(user->param)));
206     } else if (x_map >= -0.7) {
207       x[i] = PetscTanhReal((0.2-x_map)/(2.*PetscSqrtReal(user->param)));
208     } else if (x_map >= -1) {
209       x[i] = PetscTanhReal((x_map+0.9)/(2.*PetscSqrtReal(user->param)));
210     }
211   }
212   ierr = VecRestoreArray(U,&x);CHKERRQ(ierr);
213   PetscFunctionReturn(0);
214 }
215 
216 /*TEST
217 
218      test:
219        args:  -ts_rtol 1e-04 -ts_dt 0.025 -pc_type lu -ksp_error_if_not_converged TRUE  -ts_type eimex -ts_adapt_type none -ts_eimex_order_adapt -ts_eimex_max_rows 7 -ts_monitor_draw_solution
220        requires: x
221        timeoutfactor: 3
222 
223 TEST*/
224