Actual source code: ex16.c
petsc-3.6.4 2016-04-12
2: static char help[] = "Solves the van der Pol equation.\n\
3: Input parameters include:\n\
4: -mu : stiffness parameter\n\n";
6: /*
7: Concepts: TS^time-dependent nonlinear problems
8: Concepts: TS^van der Pol equation
9: Processors: 1
10: */
11: /* ------------------------------------------------------------------------
13: This program solves the van der Pol equation
14: y'' - \mu (1-y^2)*y' + y = 0 (1)
15: on the domain 0 <= x <= 1, with the boundary conditions
16: y(0) = 2, y'(0) = 0,
17: This is a nonlinear equation.
19: Notes:
20: This code demonstrates the TS solver interface to two variants of
21: linear problems, u_t = f(u,t), namely turning (1) into a system of
22: first order differential equations,
24: [ y' ] = [ z ]
25: [ z' ] [ \mu (1 - y^2) z - y ]
27: which then we can write as a vector equation
29: [ u_1' ] = [ u_2 ] (2)
30: [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ]
32: which is now in the desired form of u_t = f(u,t). One way that we
33: can split f(u,t) in (2) is to split by component,
35: [ u_1' ] = [ u_2 ] + [ 0 ]
36: [ u_2' ] [ 0 ] [ \mu (1 - u_1^2) u_2 - u_1 ]
38: where
40: [ F(u,t) ] = [ u_2 ]
41: [ 0 ]
43: and
45: [ G(u',u,t) ] = [ u_1' ] - [ 0 ]
46: [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ]
48: Using the definition of the Jacobian of G (from the PETSc user manual),
49: in the equation G(u',u,t) = F(u,t),
51: dG dG
52: J(G) = a * -- - --
53: du' du
55: where d is the partial derivative. In this example,
57: dG [ 1 ; 0 ]
58: -- = [ ]
59: du' [ 0 ; 1 ]
61: dG [ 0 ; 0 ]
62: -- = [ ]
63: du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ]
65: Hence,
67: [ a ; 0 ]
68: J(G) = [ ]
69: [ 2 \mu u_1*u_2 + 1; a - \mu (1 - u_1^2) ]
71: ------------------------------------------------------------------------- */
73: #include <petscts.h>
75: typedef struct _n_User *User;
76: struct _n_User {
77: PetscReal mu;
78: PetscBool imex;
79: PetscReal next_output;
80: };
82: /*
83: * User-defined routines
84: */
87: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
88: {
89: PetscErrorCode ierr;
90: User user = (User)ctx;
91: PetscScalar *f;
92: const PetscScalar *x;
95: VecGetArrayRead(X,&x);
96: VecGetArray(F,&f);
97: f[0] = (user->imex ? x[1] : 0);
98: f[1] = 0.0;
99: VecRestoreArrayRead(X,&x);
100: VecRestoreArray(F,&f);
101: return(0);
102: }
106: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
107: {
108: PetscErrorCode ierr;
109: User user = (User)ctx;
110: const PetscScalar *x,*xdot;
111: PetscScalar *f;
114: VecGetArrayRead(X,&x);
115: VecGetArrayRead(Xdot,&xdot);
116: VecGetArray(F,&f);
117: f[0] = xdot[0] + (user->imex ? 0 : x[1]);
118: f[1] = xdot[1] - user->mu*(1. - x[0]*x[0])*x[1] + x[0];
119: VecRestoreArrayRead(X,&x);
120: VecRestoreArrayRead(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: {
129: PetscErrorCode ierr;
130: User user = (User)ctx;
131: PetscReal mu = user->mu;
132: PetscInt rowcol[] = {0,1};
133: const PetscScalar *x;
134: PetscScalar J[2][2];
137: VecGetArrayRead(X,&x);
138: J[0][0] = a; J[0][1] = (user->imex ? 0 : 1.);
139: J[1][0] = 2.*mu*x[0]*x[1]+1.; J[1][1] = a - mu*(1. - x[0]*x[0]);
140: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
141: VecRestoreArrayRead(X,&x);
143: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
144: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
145: if (A != B) {
146: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
147: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
148: }
149: return(0);
150: }
154: static PetscErrorCode RegisterMyARK2(void)
155: {
159: {
160: const PetscReal
161: A[3][3] = {{0,0,0},
162: {0.41421356237309504880,0,0},
163: {0.75,0.25,0}},
164: At[3][3] = {{0,0,0},
165: {0.12132034355964257320,0.29289321881345247560,0},
166: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
167: *bembedt = NULL,*bembed = NULL;
168: TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL);
169: }
170: return(0);
171: }
175: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
176: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
177: {
178: PetscErrorCode ierr;
179: const PetscScalar *x;
180: PetscReal tfinal, dt;
181: User user = (User)ctx;
182: Vec interpolatedX;
185: TSGetTimeStep(ts,&dt);
186: TSGetDuration(ts,NULL,&tfinal);
188: while (user->next_output <= t && user->next_output <= tfinal) {
189: VecDuplicate(X,&interpolatedX);
190: TSInterpolate(ts,user->next_output,interpolatedX);
191: VecGetArrayRead(interpolatedX,&x);
192: PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",user->next_output,step,t,dt,(double)PetscRealPart(x[0]),(double)PetscRealPart(x[1]));
193: VecRestoreArrayRead(interpolatedX,&x);
194: VecDestroy(&interpolatedX);
196: user->next_output += 0.1;
197: }
198: return(0);
199: }
203: int main(int argc,char **argv)
204: {
205: TS ts; /* nonlinear solver */
206: Vec x; /* solution, residual vectors */
207: Mat A; /* Jacobian matrix */
208: PetscInt steps;
209: PetscReal ftime =0.5;
210: PetscBool monitor = PETSC_FALSE;
211: PetscScalar *x_ptr;
212: PetscMPIInt size;
213: struct _n_User user;
216: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217: Initialize program
218: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219: PetscInitialize(&argc,&argv,NULL,help);
221: MPI_Comm_size(PETSC_COMM_WORLD,&size);
222: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
224: RegisterMyARK2();
226: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227: Set runtime options
228: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229: user.mu = 1000;
230: user.imex = PETSC_TRUE;
231: user.next_output = 0.0;
233: PetscOptionsGetReal(NULL,"-mu",&user.mu,NULL);
234: PetscOptionsGetBool(NULL,"-imex",&user.imex,NULL);
235: PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
237: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238: Create necessary matrix and vectors, solve same ODE on every process
239: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240: MatCreate(PETSC_COMM_WORLD,&A);
241: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
242: MatSetFromOptions(A);
243: MatSetUp(A);
244: MatCreateVecs(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);
254: TSSetDuration(ts,PETSC_DEFAULT,ftime);
255: if (monitor) {
256: TSMonitorSet(ts,Monitor,&user,NULL);
257: }
259: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260: Set initial conditions
261: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262: VecGetArray(x,&x_ptr);
264: x_ptr[0] = 2; x_ptr[1] = 0.66666654321;
266: VecRestoreArray(x,&x_ptr);
267: TSSetInitialTimeStep(ts,0.0,.001);
269: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270: Set runtime options
271: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
272: TSSetFromOptions(ts);
274: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
275: Solve nonlinear system
276: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
277: TSSolve(ts,x);
278: TSGetSolveTime(ts,&ftime);
279: TSGetTimeStepNumber(ts,&steps);
280: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);
281: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
283: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
284: Free work space. All PETSc objects should be destroyed when they
285: are no longer needed.
286: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
287: MatDestroy(&A);
288: VecDestroy(&x);
289: TSDestroy(&ts);
291: PetscFinalize();
292: return(0);
293: }