Actual source code: ex16.c
petsc-3.5.4 2015-05-23
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 - 1; \mu (1 - u_1^2) ]
65: Hence,
67: [ a ; 0 ]
68: J(G) = [ ]
69: [ 2 \mu u_1 + 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: {
90: User user = (User)ctx;
91: PetscScalar *x,*f;
94: VecGetArray(X,&x);
95: VecGetArray(F,&f);
96: f[0] = (user->imex ? x[1] : 0);
97: f[1] = 0.0;
98: VecRestoreArray(X,&x);
99: VecRestoreArray(F,&f);
100: return(0);
101: }
105: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
106: {
108: User user = (User)ctx;
109: PetscScalar *x,*xdot,*f;
112: VecGetArray(X,&x);
113: VecGetArray(Xdot,&xdot);
114: VecGetArray(F,&f);
115: f[0] = xdot[0] + (user->imex ? 0 : x[1]);
116: f[1] = xdot[1] - user->mu*(1. - x[0]*x[0])*x[1] + x[0];
117: VecRestoreArray(X,&x);
118: VecRestoreArray(Xdot,&xdot);
119: VecRestoreArray(F,&f);
120: return(0);
121: }
125: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
126: {
128: User user = (User)ctx;
129: PetscReal mu = user->mu;
130: PetscInt rowcol[] = {0,1};
131: PetscScalar *x,J[2][2];
134: VecGetArray(X,&x);
135: J[0][0] = a; J[0][1] = (user->imex ? 0 : 1.);
136: J[1][0] = 2.*mu*x[0]*x[1]+1.; J[1][1] = a - mu*(1. - x[0]*x[0]);
137: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
138: VecRestoreArray(X,&x);
140: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
141: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
142: if (A != B) {
143: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
144: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
145: }
146: return(0);
147: }
151: static PetscErrorCode RegisterMyARK2(void)
152: {
156: {
157: const PetscReal
158: A[3][3] = {{0,0,0},
159: {0.41421356237309504880,0,0},
160: {0.75,0.25,0}},
161: At[3][3] = {{0,0,0},
162: {0.12132034355964257320,0.29289321881345247560,0},
163: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
164: *bembedt = NULL,*bembed = NULL;
165: TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL);
166: }
167: return(0);
168: }
172: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
173: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
174: {
175: PetscErrorCode ierr;
176: const PetscScalar *x;
177: PetscReal tfinal, dt;
178: User user = (User)ctx;
179: Vec interpolatedX;
182: TSGetTimeStep(ts,&dt);
183: TSGetDuration(ts,NULL,&tfinal);
185: while (user->next_output <= t && user->next_output <= tfinal) {
186: VecDuplicate(X,&interpolatedX);
187: TSInterpolate(ts,user->next_output,interpolatedX);
188: VecGetArrayRead(interpolatedX,&x);
189: 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]));
190: VecRestoreArrayRead(interpolatedX,&x);
191: VecDestroy(&interpolatedX);
193: user->next_output += 0.1;
194: }
195: return(0);
196: }
200: int main(int argc,char **argv)
201: {
202: TS ts; /* nonlinear solver */
203: Vec x; /* solution, residual vectors */
204: Mat A; /* Jacobian matrix */
205: PetscInt steps;
206: PetscReal ftime =0.5;
207: PetscBool monitor = PETSC_FALSE;
208: PetscScalar *x_ptr;
209: PetscMPIInt size;
210: struct _n_User user;
213: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: Initialize program
215: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216: PetscInitialize(&argc,&argv,NULL,help);
218: MPI_Comm_size(PETSC_COMM_WORLD,&size);
219: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
221: RegisterMyARK2();
223: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224: Set runtime options
225: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226: user.mu = 1000;
227: user.imex = PETSC_TRUE;
228: user.next_output = 0.0;
230: PetscOptionsGetReal(NULL,"-mu",&user.mu,NULL);
231: PetscOptionsGetBool(NULL,"-imex",&user.imex,NULL);
232: PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
234: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235: Create necessary matrix and vectors, solve same ODE on every process
236: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
237: MatCreate(PETSC_COMM_WORLD,&A);
238: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
239: MatSetFromOptions(A);
240: MatSetUp(A);
241: MatGetVecs(A,&x,NULL);
243: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244: Create timestepping solver context
245: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246: TSCreate(PETSC_COMM_WORLD,&ts);
247: TSSetType(ts,TSBEULER);
248: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
249: TSSetIFunction(ts,NULL,IFunction,&user);
250: TSSetIJacobian(ts,A,A,IJacobian,&user);
251: TSSetDuration(ts,PETSC_DEFAULT,ftime);
252: if (monitor) {
253: TSMonitorSet(ts,Monitor,&user,NULL);
254: }
256: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257: Set initial conditions
258: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
259: VecGetArray(x,&x_ptr);
261: x_ptr[0] = 2; x_ptr[1] = 0.66666654321;
263: VecRestoreArray(x,&x_ptr);
264: TSSetInitialTimeStep(ts,0.0,.001);
266: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267: Set runtime options
268: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
269: TSSetFromOptions(ts);
271: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272: Solve nonlinear system
273: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274: TSSolve(ts,x);
275: TSGetSolveTime(ts,&ftime);
276: TSGetTimeStepNumber(ts,&steps);
277: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);
278: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
280: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
281: Free work space. All PETSc objects should be destroyed when they
282: are no longer needed.
283: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
284: MatDestroy(&A);
285: VecDestroy(&x);
286: TSDestroy(&ts);
288: PetscFinalize();
289: return(0);
290: }