Actual source code: ex20.c
petsc-3.10.5 2019-03-28
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: */
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.0);
98: f[1] = 0.0;
99: VecRestoreArrayRead(X,&x);
100: VecRestoreArray(F,&f);
101: return(0);
102: }
104: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
105: {
106: PetscErrorCode ierr;
107: User user = (User)ctx;
108: const PetscScalar *x,*xdot;
109: PetscScalar *f;
112: VecGetArrayRead(X,&x);
113: VecGetArrayRead(Xdot,&xdot);
114: VecGetArray(F,&f);
115: f[0] = xdot[0] - (user->imex ? 0 : x[1]);
116: f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
117: VecRestoreArrayRead(X,&x);
118: VecRestoreArrayRead(Xdot,&xdot);
119: VecRestoreArray(F,&f);
120: return(0);
121: }
123: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
124: {
125: PetscErrorCode ierr;
126: User user = (User)ctx;
127: PetscInt rowcol[] = {0,1};
128: const PetscScalar *x;
129: PetscScalar J[2][2];
132: VecGetArrayRead(X,&x);
133: J[0][0] = a; J[0][1] = (user->imex ? 0 : -1.0);
134: 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]);
135: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
136: VecRestoreArrayRead(X,&x);
138: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
139: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
140: if (A != B) {
141: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
142: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
143: }
144: return(0);
145: }
147: /* This is an example of registering an user-provided ARKIMEX scheme */
148: static PetscErrorCode RegisterMyARK2(void)
149: {
153: {
154: const PetscReal
155: A[3][3] = {{0,0,0},
156: {0.41421356237309504880,0,0},
157: {0.75,0.25,0}},
158: At[3][3] = {{0,0,0},
159: {0.12132034355964257320,0.29289321881345247560,0},
160: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}};
161: TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);
162: }
163: return(0);
164: }
166: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
167: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
168: {
169: PetscErrorCode ierr;
170: const PetscScalar *x;
171: PetscReal tfinal, dt;
172: User user = (User)ctx;
173: Vec interpolatedX;
176: TSGetTimeStep(ts,&dt);
177: TSGetMaxTime(ts,&tfinal);
179: while (user->next_output <= t && user->next_output <= tfinal) {
180: VecDuplicate(X,&interpolatedX);
181: TSInterpolate(ts,user->next_output,interpolatedX);
182: VecGetArrayRead(interpolatedX,&x);
183: PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
184: user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
185: (double)PetscRealPart(x[1]));
186: VecRestoreArrayRead(interpolatedX,&x);
187: VecDestroy(&interpolatedX);
188: user->next_output += 0.1;
189: }
190: return(0);
191: }
193: int main(int argc,char **argv)
194: {
195: TS ts; /* nonlinear solver */
196: Vec x; /* solution, residual vectors */
197: Mat A; /* Jacobian matrix */
198: PetscInt steps;
199: PetscReal ftime = 0.5;
200: PetscBool monitor = PETSC_FALSE;
201: PetscScalar *x_ptr;
202: PetscMPIInt size;
203: struct _n_User user;
206: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207: Initialize program
208: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
210: MPI_Comm_size(PETSC_COMM_WORLD,&size);
211: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
213: /* Register user-specified ARKIMEX method */
214: RegisterMyARK2();
216: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217: Set runtime options
218: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219: user.imex = PETSC_TRUE;
220: user.next_output = 0.0;
221: user.mu = 1.0e6;
222: PetscOptionsGetBool(NULL,NULL,"-imex",&user.imex,NULL);
223: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
224: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
225: PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);
226: PetscOptionsEnd();
228: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229: Create necessary matrix and vectors, solve same ODE on every process
230: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231: MatCreate(PETSC_COMM_WORLD,&A);
232: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
233: MatSetFromOptions(A);
234: MatSetUp(A);
236: MatCreateVecs(A,&x,NULL);
238: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239: Create timestepping solver context
240: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
241: TSCreate(PETSC_COMM_WORLD,&ts);
242: TSSetType(ts,TSBEULER);
243: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
244: TSSetIFunction(ts,NULL,IFunction,&user);
245: TSSetIJacobian(ts,A,A,IJacobian,&user);
247: TSSetMaxTime(ts,ftime);
248: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
249: if (monitor) {
250: TSMonitorSet(ts,Monitor,&user,NULL);
251: }
253: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
254: Set initial conditions
255: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
256: VecGetArray(x,&x_ptr);
257: x_ptr[0] = 2.0; x_ptr[1] = -6.666665432100101e-01;
258: VecRestoreArray(x,&x_ptr);
259: TSSetTimeStep(ts,.001);
261: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262: Set runtime options
263: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
264: TSSetFromOptions(ts);
266: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267: Solve nonlinear system
268: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
269: TSSolve(ts,x);
270: TSGetSolveTime(ts,&ftime);
271: TSGetStepNumber(ts,&steps);
272: PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
273: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
275: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276: Free work space. All PETSc objects should be destroyed when they
277: are no longer needed.
278: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
279: MatDestroy(&A);
280: VecDestroy(&x);
281: TSDestroy(&ts);
283: PetscFinalize();
284: return(ierr);
285: }
287: /*TEST
289: test:
290: requires: !single
292: TEST*/