xref: /petsc/src/ts/tutorials/eimex/allen_cahn.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 *);
23*2a8381b2SBarry Smith static PetscErrorCode FormIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscCtx ctx);
24c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS, Vec, void *);
25c4762a1bSJed Brown 
main(int argc,char ** argv)26d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
27d71ae5a4SJacob Faibussowitsch {
28c4762a1bSJed Brown   TS        ts;
29c4762a1bSJed Brown   Vec       x; /* solution vector */
30c4762a1bSJed Brown   Mat       A; /* Jacobian */
31c4762a1bSJed Brown   PetscInt  steps, mx;
32c4762a1bSJed Brown   PetscReal ftime;
33c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
34c4762a1bSJed Brown 
35327415f7SBarry Smith   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Initialize user application context */
39d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Allen-Cahn equation", "");
40c4762a1bSJed Brown   user.param  = 9e-4;
41c4762a1bSJed Brown   user.xleft  = -1.;
42c4762a1bSJed Brown   user.xright = 2.;
43c4762a1bSJed Brown   user.mx     = 400;
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-eps", "parameter", "", user.param, &user.param, NULL));
45d0609cedSBarry Smith   PetscOptionsEnd();
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48c4762a1bSJed Brown    Set runtime options
49c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50c4762a1bSJed Brown   /*
519566063dSJacob Faibussowitsch    * PetscCall(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    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
569566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
579566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, user.mx, user.mx));
589566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
599566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
609566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63c4762a1bSJed Brown    Create time stepping solver context
64c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
659566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
669566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSEIMEX));
679566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
689566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &user));
699566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, FormIJacobian, &user));
70c4762a1bSJed Brown   ftime = 22;
719566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
729566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75c4762a1bSJed Brown    Set initial conditions
76c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
779566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts, x, &user));
789566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
799566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &mx));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82c4762a1bSJed Brown    Set runtime options
83c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
849566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87c4762a1bSJed Brown    Solve nonlinear system
88c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
899566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
909566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &ftime));
919566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
9263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "eps %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.param, steps, (double)ftime));
939566063dSJacob Faibussowitsch   /*   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));*/
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96c4762a1bSJed Brown    Free work space.
97c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1009566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1019566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
102b122ec5aSJacob Faibussowitsch   return 0;
103c4762a1bSJed Brown }
104c4762a1bSJed Brown 
RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void * ptr)105d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr)
106d71ae5a4SJacob Faibussowitsch {
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);
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1189566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
119c4762a1bSJed Brown   f[0] = 2. * eps * (x[1] - x[0]) / (hx * hx); /*boundary*/
120ad540459SPierre Jolivet   for (i = 1; i < mx - 1; i++) f[i] = eps * (x[i + 1] - 2. * x[i] + x[i - 1]) / (hx * hx);
121c4762a1bSJed Brown   f[mx - 1] = 2. * eps * (x[mx - 2] - x[mx - 1]) / (hx * hx); /*boundary*/
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void * ptr)127d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ptr)
128d71ae5a4SJacob Faibussowitsch {
129c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
130c4762a1bSJed Brown   PetscScalar       *f;
131c4762a1bSJed Brown   const PetscScalar *x, *xdot;
132c4762a1bSJed Brown   PetscInt           i, mx;
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   PetscFunctionBegin;
135c4762a1bSJed Brown   mx = user->mx;
1369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
1389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
139c4762a1bSJed Brown 
140ad540459SPierre Jolivet   for (i = 0; i < mx; i++) f[i] = xdot[i] - x[i] * (1. - x[i] * x[i]);
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,PetscCtx ctx)148*2a8381b2SBarry Smith static PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, PetscCtx ctx)
149d71ae5a4SJacob Faibussowitsch {
150c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ctx;
151c4762a1bSJed Brown   PetscScalar        v;
152c4762a1bSJed Brown   const PetscScalar *x;
153c4762a1bSJed Brown   PetscInt           i, col;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &x));
157c4762a1bSJed Brown   for (i = 0; i < user->mx; i++) {
158c4762a1bSJed Brown     v   = a - 1. + 3. * x[i] * x[i];
159c4762a1bSJed Brown     col = i;
1609566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, &i, 1, &col, &v, INSERT_VALUES));
161c4762a1bSJed Brown   }
1629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &x));
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
166c4762a1bSJed Brown   if (J != Jpre) {
1679566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
1689566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
169c4762a1bSJed Brown   }
170c4762a1bSJed Brown   /*  MatView(J,PETSC_VIEWER_STDOUT_WORLD);*/
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
FormInitialSolution(TS ts,Vec U,PetscCtx ctx)174*2a8381b2SBarry Smith static PetscErrorCode FormInitialSolution(TS ts, Vec U, PetscCtx ctx)
175d71ae5a4SJacob Faibussowitsch {
176c4762a1bSJed Brown   AppCtx      *user = (AppCtx *)ctx;
177c4762a1bSJed Brown   PetscInt     i;
178c4762a1bSJed Brown   PetscScalar *x;
179c4762a1bSJed Brown   PetscReal    hx, x_map;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   PetscFunctionBegin;
182c4762a1bSJed Brown   hx = (user->xright - user->xleft) / (PetscReal)(user->mx - 1);
1839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &x));
184c4762a1bSJed Brown   for (i = 0; i < user->mx; i++) {
185c4762a1bSJed Brown     x_map = user->xleft + i * hx;
186c4762a1bSJed Brown     if (x_map >= 0.7065) {
187c4762a1bSJed Brown       x[i] = PetscTanhReal((x_map - 0.8) / (2. * PetscSqrtReal(user->param)));
188c4762a1bSJed Brown     } else if (x_map >= 0.4865) {
189c4762a1bSJed Brown       x[i] = PetscTanhReal((0.613 - x_map) / (2. * PetscSqrtReal(user->param)));
190c4762a1bSJed Brown     } else if (x_map >= 0.28) {
191c4762a1bSJed Brown       x[i] = PetscTanhReal((x_map - 0.36) / (2. * PetscSqrtReal(user->param)));
192c4762a1bSJed Brown     } else if (x_map >= -0.7) {
193c4762a1bSJed Brown       x[i] = PetscTanhReal((0.2 - x_map) / (2. * PetscSqrtReal(user->param)));
194c4762a1bSJed Brown     } else if (x_map >= -1) {
195c4762a1bSJed Brown       x[i] = PetscTanhReal((x_map + 0.9) / (2. * PetscSqrtReal(user->param)));
196c4762a1bSJed Brown     }
197c4762a1bSJed Brown   }
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &x));
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
200c4762a1bSJed Brown }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown /*TEST
203c4762a1bSJed Brown 
204c4762a1bSJed Brown      test:
205188af4bfSBarry Smith        args: -ts_rtol 1e-04 -ts_time_step 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
206c4762a1bSJed Brown        requires: x
207c4762a1bSJed Brown        timeoutfactor: 3
208c4762a1bSJed Brown 
209c4762a1bSJed Brown TEST*/
210