Actual source code: ex19.c
petsc-3.7.3 2016-08-01
2: static char help[] = "Solves the van der Pol DAE.\n\
3: Input parameters include:\n";
5: /*
6: Concepts: TS^time-dependent nonlinear problems
7: Concepts: TS^van der Pol DAE
8: Processors: 1
9: */
10: /* ------------------------------------------------------------------------
12: This program solves the van der Pol DAE
13: y' = -z = f(y,z) (1)
14: 0 = y-(z^3/3 - z) = g(y,z)
15: on the domain 0 <= x <= 1, with the boundary conditions
16: y(0) = -2, y'(0) = -2.355301397608119909925287735864250951918
17: This is a nonlinear equation.
19: Notes:
20: This code demonstrates the TS solver interface with the Van der Pol DAE,
21: namely it is the case when there is no RHS (meaning the RHS == 0), and the
22: equations are converted to two variants of linear problems, u_t = f(u,t),
23: namely turning (1) into a vector equation in terms of u,
25: [ y' + z ] = [ 0 ]
26: [ (z^3/3 - z) - y ] [ 0 ]
28: which then we can write as a vector equation
30: [ u_1' + u_2 ] = [ 0 ] (2)
31: [ (u_2^3/3 - u_2) - u_1 ] [ 0 ]
33: which is now in the desired form of u_t = f(u,t). As this is a DAE, and
34: there is no u_2', there is no need for a split,
36: so
38: [ G(u',u,t) ] = [ u_1' ] + [ u_2 ]
39: [ 0 ] [ (u_2^3/3 - u_2) - u_1 ]
41: Using the definition of the Jacobian of G (from the PETSc user manual),
42: in the equation G(u',u,t) = F(u,t),
44: dG dG
45: J(G) = a * -- - --
46: du' du
48: where d is the partial derivative. In this example,
50: dG [ 1 ; 0 ]
51: -- = [ ]
52: du' [ 0 ; 0 ]
54: dG [ 0 ; 1 ]
55: -- = [ ]
56: du [ -1 ; 1 - u_2^2 ]
58: Hence,
60: [ a ; -1 ]
61: J(G) = [ ]
62: [ 1 ; u_2^2 - 1 ]
64: ------------------------------------------------------------------------- */
66: #include <petscts.h>
68: typedef struct _n_User *User;
69: struct _n_User {
70: PetscReal next_output;
71: };
73: /*
74: * User-defined routines
75: */
79: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
80: {
81: PetscErrorCode ierr;
82: PetscScalar *f;
83: const PetscScalar *x,*xdot;
86: VecGetArrayRead(X,&x);
87: VecGetArrayRead(Xdot,&xdot);
88: VecGetArray(F,&f);
89: f[0] = xdot[0] + x[1];
90: f[1] = (x[1]*x[1]*x[1]/3.0 - x[1])-x[0];
91: VecRestoreArrayRead(X,&x);
92: VecRestoreArrayRead(Xdot,&xdot);
93: VecRestoreArray(F,&f);
94: return(0);
95: }
99: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
100: {
101: PetscErrorCode ierr;
102: PetscInt rowcol[] = {0,1};
103: PetscScalar J[2][2];
104: const PetscScalar *x;
107: VecGetArrayRead(X,&x);
108: J[0][0] = a; J[0][1] = -1.;
109: J[1][0] = 1.; J[1][1] = -1. + x[1]*x[1];
110: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
111: VecRestoreArrayRead(X,&x);
113: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
114: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
115: if (A != B) {
116: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
117: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
118: }
119: return(0);
120: }
124: static PetscErrorCode RegisterMyARK2(void)
125: {
129: {
130: const PetscReal
131: A[3][3] = {{0,0,0},
132: {0.41421356237309504880,0,0},
133: {0.75,0.25,0}},
134: At[3][3] = {{0,0,0},
135: {0.12132034355964257320,0.29289321881345247560,0},
136: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
137: *bembedt = NULL,*bembed = NULL;
138: TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL);
139: }
140: return(0);
141: }
145: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
146: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
147: {
148: PetscErrorCode ierr;
149: const PetscScalar *x;
150: PetscReal tfinal, dt;
151: User user = (User)ctx;
152: Vec interpolatedX;
155: TSGetTimeStep(ts,&dt);
156: TSGetDuration(ts,NULL,&tfinal);
158: while (user->next_output <= t && user->next_output <= tfinal) {
159: VecDuplicate(X,&interpolatedX);
160: TSInterpolate(ts,user->next_output,interpolatedX);
161: VecGetArrayRead(interpolatedX,&x);
162: 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]));
163: VecRestoreArrayRead(interpolatedX,&x);
164: VecDestroy(&interpolatedX);
165: user->next_output += 0.1;
166: }
167: return(0);
168: }
172: int main(int argc,char **argv)
173: {
174: TS ts; /* nonlinear solver */
175: Vec x; /* solution, residual vectors */
176: Mat A; /* Jacobian matrix */
177: PetscInt steps;
178: PetscReal ftime = 0.5;
179: PetscBool monitor = PETSC_FALSE;
180: PetscScalar *x_ptr;
181: PetscMPIInt size;
182: struct _n_User user;
185: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186: Initialize program
187: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188: PetscInitialize(&argc,&argv,NULL,help);
190: MPI_Comm_size(PETSC_COMM_WORLD,&size);
191: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
193: RegisterMyARK2();
195: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196: Set runtime options
197: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
199: user.next_output = 0.0;
200: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
202: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: Create necessary matrix and vectors, solve same ODE on every process
204: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205: MatCreate(PETSC_COMM_WORLD,&A);
206: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
207: MatSetFromOptions(A);
208: MatSetUp(A);
209: MatCreateVecs(A,&x,NULL);
211: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212: Create timestepping solver context
213: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214: TSCreate(PETSC_COMM_WORLD,&ts);
215: TSSetType(ts,TSBEULER);
216: TSSetIFunction(ts,NULL,IFunction,&user);
217: TSSetIJacobian(ts,A,A,IJacobian,&user);
218: TSSetDuration(ts,PETSC_DEFAULT,ftime);
219: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
220: if (monitor) {
221: TSMonitorSet(ts,Monitor,&user,NULL);
222: }
224: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: Set initial conditions
226: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227: VecGetArray(x,&x_ptr);
228: x_ptr[0] = -2; x_ptr[1] = -2.355301397608119909925287735864250951918;
229: VecRestoreArray(x,&x_ptr);
230: TSSetInitialTimeStep(ts,0.0,.001);
232: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233: Set runtime options
234: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
235: TSSetFromOptions(ts);
237: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238: Solve nonlinear system
239: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240: TSSolve(ts,x);
241: TSGetSolveTime(ts,&ftime);
242: TSGetTimeStepNumber(ts,&steps);
243: PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
244: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
246: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247: Free work space. All PETSc objects should be destroyed when they
248: are no longer needed.
249: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250: MatDestroy(&A);
251: VecDestroy(&x);
252: TSDestroy(&ts);
254: PetscFinalize();
255: return(0);
256: }