xref: /petsc/src/ts/tutorials/eimex/allen_cahn.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] ="Solves the time dependent Allen-Cahn equation with IMEX methods";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown  * allen_cahn.c
5c4762a1bSJed Brown  *
6c4762a1bSJed Brown  *  Created on: Jun 8, 2012
7c4762a1bSJed Brown  *      Author: Hong Zhang
8c4762a1bSJed Brown  */
9c4762a1bSJed Brown 
10c4762a1bSJed Brown #include <petscts.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown  * application context
14c4762a1bSJed Brown  */
15c4762a1bSJed Brown typedef struct {
16c4762a1bSJed Brown   PetscReal   param;        /* parameter */
17c4762a1bSJed Brown   PetscReal   xleft,xright;  /* range in x-direction */
18c4762a1bSJed Brown   PetscInt    mx;           /* Discretization in x-direction */
19c4762a1bSJed Brown }AppCtx;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
22c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
23c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *ctx);
24c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown int main(int argc, char **argv)
27c4762a1bSJed Brown {
28c4762a1bSJed Brown   TS                ts;
29c4762a1bSJed Brown   Vec               x; /*solution vector*/
30c4762a1bSJed Brown   Mat               A; /*Jacobian*/
31c4762a1bSJed Brown   PetscInt          steps,mx;
32c4762a1bSJed Brown   PetscErrorCode    ierr;
33c4762a1bSJed Brown   PetscReal         ftime;
34c4762a1bSJed Brown   AppCtx            user;       /* user-defined work context */
35c4762a1bSJed Brown 
36*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Initialize user application context */
391e1ea65dSPierre Jolivet   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Allen-Cahn equation","");CHKERRQ(ierr);
40c4762a1bSJed Brown   user.param = 9e-4;
41c4762a1bSJed Brown   user.xleft = -1.;
42c4762a1bSJed Brown   user.xright = 2.;
43c4762a1bSJed Brown   user.mx = 400;
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-eps","parameter","",user.param,&user.param,NULL));
45c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48c4762a1bSJed Brown    Set runtime options
49c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50c4762a1bSJed Brown   /*
515f80ce2aSJacob Faibussowitsch    * CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
52c4762a1bSJed Brown    */
53c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54c4762a1bSJed Brown    Create necessary matrix and vectors, solve same ODE on every process
55c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,user.mx,user.mx));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,NULL));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63c4762a1bSJed Brown    Create time stepping solver context
64c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
655f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSEIMEX));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&user));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,FormIFunction,&user));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,A,A,FormIJacobian,&user));
70c4762a1bSJed Brown   ftime = 22;
715f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,ftime));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75c4762a1bSJed Brown    Set initial conditions
76c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
775f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialSolution(ts,x,&user));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(x,&mx));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82c4762a1bSJed Brown    Set runtime options
83c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
845f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87c4762a1bSJed Brown    Solve nonlinear system
88c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
895f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts,&ftime));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"eps %g, steps %D, ftime %g\n",(double)user.param,steps,(double)ftime));
935f80ce2aSJacob Faibussowitsch   /*   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));*/
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96c4762a1bSJed Brown    Free work space.
97c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
101*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
102*b122ec5aSJacob Faibussowitsch   return 0;
103c4762a1bSJed Brown }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
106c4762a1bSJed Brown {
107c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
108c4762a1bSJed Brown   PetscScalar       *f;
109c4762a1bSJed Brown   const PetscScalar *x;
110c4762a1bSJed Brown   PetscInt          i,mx;
111c4762a1bSJed Brown   PetscReal         hx,eps;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   PetscFunctionBegin;
114c4762a1bSJed Brown   mx = user->mx;
115c4762a1bSJed Brown   eps = user->param;
116c4762a1bSJed Brown   hx = (user->xright-user->xleft)/(mx-1);
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
119c4762a1bSJed Brown   f[0] = 2.*eps*(x[1]-x[0])/(hx*hx); /*boundary*/
120c4762a1bSJed Brown   for (i=1;i<mx-1;i++) {
121c4762a1bSJed Brown     f[i]= eps*(x[i+1]-2.*x[i]+x[i-1])/(hx*hx);
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown   f[mx-1] = 2.*eps*(x[mx-2]- x[mx-1])/(hx*hx); /*boundary*/
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
126c4762a1bSJed Brown   PetscFunctionReturn(0);
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
130c4762a1bSJed Brown {
131c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
132c4762a1bSJed Brown   PetscScalar       *f;
133c4762a1bSJed Brown   const PetscScalar *x,*xdot;
134c4762a1bSJed Brown   PetscInt          i,mx;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   PetscFunctionBegin;
137c4762a1bSJed Brown   mx = user->mx;
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Xdot,&xdot));
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   for (i=0;i<mx;i++) {
143c4762a1bSJed Brown     f[i]= xdot[i] - x[i]*(1.-x[i]*x[i]);
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Xdot,&xdot));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U, Vec Udot, PetscReal a, Mat J,Mat Jpre,void *ctx)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ctx;
155c4762a1bSJed Brown   PetscScalar       v;
156c4762a1bSJed Brown   const PetscScalar *x;
157c4762a1bSJed Brown   PetscInt          i,col;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBegin;
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&x));
161c4762a1bSJed Brown   for (i=0; i < user->mx; i++) {
162c4762a1bSJed Brown     v = a - 1. + 3.*x[i]*x[i];
163c4762a1bSJed Brown     col = i;
1645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(J,1,&i,1,&col,&v,INSERT_VALUES));
165c4762a1bSJed Brown   }
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&x));
167c4762a1bSJed Brown 
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
170c4762a1bSJed Brown   if (J != Jpre) {
1715f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown   /*  MatView(J,PETSC_VIEWER_STDOUT_WORLD);*/
175c4762a1bSJed Brown   PetscFunctionReturn(0);
176c4762a1bSJed Brown }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ctx)
179c4762a1bSJed Brown {
180c4762a1bSJed Brown   AppCtx *user = (AppCtx*)ctx;
181c4762a1bSJed Brown   PetscInt       i;
182c4762a1bSJed Brown   PetscScalar    *x;
183c4762a1bSJed Brown   PetscReal      hx,x_map;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   PetscFunctionBegin;
186c4762a1bSJed Brown   hx = (user->xright-user->xleft)/(PetscReal)(user->mx-1);
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(U,&x));
188c4762a1bSJed Brown   for (i=0;i<user->mx;i++) {
189c4762a1bSJed Brown     x_map = user->xleft + i*hx;
190c4762a1bSJed Brown     if (x_map >= 0.7065) {
191c4762a1bSJed Brown       x[i] = PetscTanhReal((x_map-0.8)/(2.*PetscSqrtReal(user->param)));
192c4762a1bSJed Brown     } else if (x_map >= 0.4865) {
193c4762a1bSJed Brown       x[i] = PetscTanhReal((0.613-x_map)/(2.*PetscSqrtReal(user->param)));
194c4762a1bSJed Brown     } else if (x_map >= 0.28) {
195c4762a1bSJed Brown       x[i] = PetscTanhReal((x_map-0.36)/(2.*PetscSqrtReal(user->param)));
196c4762a1bSJed Brown     } else if (x_map >= -0.7) {
197c4762a1bSJed Brown       x[i] = PetscTanhReal((0.2-x_map)/(2.*PetscSqrtReal(user->param)));
198c4762a1bSJed Brown     } else if (x_map >= -1) {
199c4762a1bSJed Brown       x[i] = PetscTanhReal((x_map+0.9)/(2.*PetscSqrtReal(user->param)));
200c4762a1bSJed Brown     }
201c4762a1bSJed Brown   }
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(U,&x));
203c4762a1bSJed Brown   PetscFunctionReturn(0);
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /*TEST
207c4762a1bSJed Brown 
208c4762a1bSJed Brown      test:
209c4762a1bSJed Brown        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
210c4762a1bSJed Brown        requires: x
211c4762a1bSJed Brown        timeoutfactor: 3
212c4762a1bSJed Brown 
213c4762a1bSJed Brown TEST*/
214