Actual source code: ex20.c
petsc-3.5.4 2015-05-23
2: static char help[] = "Solves the van der Pol equation.\n\
3: Input parameters include:\n";
5: /*
6: Concepts: TS^time-dependent nonlinear problems
7: Concepts: TS^van der Pol equation DAE equivalent
8: Processors: 1
9: */
10: /* ------------------------------------------------------------------------
12: This program solves the van der Pol DAE ODE equivalent
13: y' = z (1)
14: z' = mu[(1-y^2)z-y]
15: on the domain 0 <= x <= 1, with the boundary conditions
16: y(0) = 2, y'(0) = -6.666665432100101e-01,
17: and
18: mu = 10^6.
19: This is a nonlinear equation.
21: Notes:
22: This code demonstrates the TS solver interface to a variant of
23: linear problems, u_t = f(u,t), namely turning (1) into a system of
24: first order differential equations,
26: [ y' ] = [ z ]
27: [ z' ] [ mu[(1-y^2)z-y] ]
29: which then we can write as a vector equation
31: [ u_1' ] = [ u_2 ] (2)
32: [ u_2' ] [ mu[(1-u_1^2)u_2-u_1] ]
34: which is now in the desired form of u_t = f(u,t). One way that we
35: can split f(u,t) in (2) is to split by component,
37: [ u_1' ] = [ u_2 ] + [ 0 ]
38: [ u_2' ] [ 0 ] [ mu[(1-u_1^2)u_2-u_1] ]
40: where
42: [ F(u,t) ] = [ u_2 ]
43: [ 0 ]
45: and
47: [ G(u',u,t) ] = [ u_1' ] - [ 0 ]
48: [ u_2' ] [ mu[(1-u_1^2)u_2-u_1] ]
50: Using the definition of the Jacobian of G (from the PETSc user manual),
51: in the equation G(u',u,t) = F(u,t),
53: dG dG
54: J(G) = a * -- - --
55: du' du
57: where d is the partial derivative. In this example,
59: dG [ 1 ; 0 ]
60: -- = [ ]
61: du' [ 0 ; 1 ]
63: dG [ 0 ; 0 ]
64: -- = [ ]
65: du [ -mu*(1.0 + 2.0*u_1*u_2) ; mu*(1-u_1*u_1) ]
67: Hence,
69: [ a ; 0 ]
70: J(G) = [ ]
71: [ mu*(1.0 + 2.0*u_1*u_2) ; a - mu*(1-u_1*u_1) ]
73: ------------------------------------------------------------------------- */
75: #include <petscts.h>
77: typedef struct _n_User *User;
78: struct _n_User {
79: PetscReal mu;
80: PetscBool imex;
81: PetscReal next_output;
82: };
84: /*
85: * User-defined routines
86: */
89: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
90: {
92: User user = (User)ctx;
93: PetscScalar *x,*f;
96: VecGetArray(X,&x);
97: VecGetArray(F,&f);
98: f[0] = (user->imex ? x[1] : 0.0);
99: f[1] = 0.0;
100: VecRestoreArray(X,&x);
101: VecRestoreArray(F,&f);
102: return(0);
103: }
107: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
108: {
110: User user = (User)ctx;
111: PetscScalar *x,*xdot,*f;
114: VecGetArray(X,&x);
115: VecGetArray(Xdot,&xdot);
116: VecGetArray(F,&f);
117: f[0] = xdot[0] - (user->imex ? 0 : x[1]);
118: f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
119: VecRestoreArray(X,&x);
120: VecRestoreArray(Xdot,&xdot);
121: VecRestoreArray(F,&f);
122: return(0);
123: }
127: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
128: {
130: User user = (User)ctx;
131: PetscInt rowcol[] = {0,1};
132: PetscScalar *x,J[2][2];
135: VecGetArray(X,&x);
136: J[0][0] = a; J[0][1] = (user->imex ? 0 : -1.0);
137: J[1][0] = user->mu*(1.0 + 2.0*x[0]*x[1]); J[1][1] = a - user->mu*(1.0-x[0]*x[0]);
138: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
139: VecRestoreArray(X,&x);
141: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
142: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
143: if (A != B) {
144: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
145: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
146: }
147: return(0);
148: }
152: /* This is an example of registering an user-provided ARKIMEX scheme */
153: static PetscErrorCode RegisterMyARK2(void)
154: {
158: {
159: const PetscReal
160: A[3][3] = {{0,0,0},
161: {0.41421356237309504880,0,0},
162: {0.75,0.25,0}},
163: At[3][3] = {{0,0,0},
164: {0.12132034355964257320,0.29289321881345247560,0},
165: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}};
166: TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);
167: }
168: return(0);
169: }
173: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
174: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
175: {
176: PetscErrorCode ierr;
177: const PetscScalar *x;
178: PetscReal tfinal, dt;
179: User user = (User)ctx;
180: Vec interpolatedX;
183: TSGetTimeStep(ts,&dt);
184: TSGetDuration(ts,NULL,&tfinal);
186: while (user->next_output <= t && user->next_output <= tfinal) {
187: VecDuplicate(X,&interpolatedX);
188: TSInterpolate(ts,user->next_output,interpolatedX);
189: VecGetArrayRead(interpolatedX,&x);
190: PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
191: user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
192: (double)PetscRealPart(x[1]));
193: VecRestoreArrayRead(interpolatedX,&x);
194: VecDestroy(&interpolatedX);
195: user->next_output += 0.1;
196: }
197: return(0);
198: }
202: int main(int argc,char **argv)
203: {
204: TS ts; /* nonlinear solver */
205: Vec x; /* solution, residual vectors */
206: Mat A; /* Jacobian matrix */
207: PetscInt steps;
208: PetscReal ftime = 0.5;
209: PetscBool monitor = PETSC_FALSE;
210: PetscScalar *x_ptr;
211: PetscMPIInt size;
212: struct _n_User user;
215: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216: Initialize program
217: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
218: PetscInitialize(&argc,&argv,NULL,help);
220: MPI_Comm_size(PETSC_COMM_WORLD,&size);
221: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
223: /* Register user-specified ARKIMEX method */
224: RegisterMyARK2();
226: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227: Set runtime options
228: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229: user.imex = PETSC_TRUE;
230: user.next_output = 0.0;
231: user.mu = 1.0e6;
232: PetscOptionsGetBool(NULL,"-imex",&user.imex,NULL);
233: PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
234: PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);
236: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237: Create necessary matrix and vectors, solve same ODE on every process
238: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239: MatCreate(PETSC_COMM_WORLD,&A);
240: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
241: MatSetFromOptions(A);
242: MatSetUp(A);
244: MatGetVecs(A,&x,NULL);
246: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247: Create timestepping solver context
248: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
249: TSCreate(PETSC_COMM_WORLD,&ts);
250: TSSetType(ts,TSBEULER);
251: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
252: TSSetIFunction(ts,NULL,IFunction,&user);
253: TSSetIJacobian(ts,A,A,IJacobian,&user);
255: TSSetDuration(ts,PETSC_DEFAULT,ftime);
256: if (monitor) {
257: TSMonitorSet(ts,Monitor,&user,NULL);
258: }
260: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
261: Set initial conditions
262: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
263: VecGetArray(x,&x_ptr);
264: x_ptr[0] = 2.0; x_ptr[1] = -6.666665432100101e-01;
265: VecRestoreArray(x,&x_ptr);
266: TSSetInitialTimeStep(ts,0.0,.001);
268: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269: Set runtime options
270: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271: TSSetFromOptions(ts);
273: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
274: Solve nonlinear system
275: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
276: TSSolve(ts,x);
277: TSGetSolveTime(ts,&ftime);
278: TSGetTimeStepNumber(ts,&steps);
279: PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
280: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
282: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283: Free work space. All PETSc objects should be destroyed when they
284: are no longer needed.
285: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
286: MatDestroy(&A);
287: VecDestroy(&x);
288: TSDestroy(&ts);
290: PetscFinalize();
291: return(0);
292: }