Actual source code: ex19.c
petsc-dev 2014-02-02
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: {
82: PetscScalar *x,*xdot,*f;
85: VecGetArray(X,&x);
86: VecGetArray(Xdot,&xdot);
87: VecGetArray(F,&f);
88: f[0] = xdot[0] + x[1];
89: f[1] = (x[1]*x[1]*x[1]/3 - x[1])-x[0];
90: VecRestoreArray(X,&x);
91: VecRestoreArray(Xdot,&xdot);
92: VecRestoreArray(F,&f);
93: return(0);
94: }
98: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx)
99: {
101: PetscInt rowcol[] = {0,1};
102: PetscScalar *x,J[2][2];
105: VecGetArray(X,&x);
106: J[0][0] = a; J[0][1] = -1.;
107: J[1][0] = 1.; J[1][1] = -1. + x[1]*x[1];
108: MatSetValues(*B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
109: VecRestoreArray(X,&x);
111: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
112: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
113: if (*A != *B) {
114: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
115: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
116: }
117: *flag = SAME_NONZERO_PATTERN;
118: return(0);
119: }
123: static PetscErrorCode RegisterMyARK2(void)
124: {
128: {
129: const PetscReal
130: A[3][3] = {{0,0,0},
131: {0.41421356237309504880,0,0},
132: {0.75,0.25,0}},
133: At[3][3] = {{0,0,0},
134: {0.12132034355964257320,0.29289321881345247560,0},
135: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
136: *bembedt = NULL,*bembed = NULL;
137: TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL);
138: }
139: return(0);
140: }
144: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
145: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
146: {
147: PetscErrorCode ierr;
148: const PetscScalar *x;
149: PetscReal tfinal, dt;
150: User user = (User)ctx;
151: Vec interpolatedX;
154: TSGetTimeStep(ts,&dt);
155: TSGetDuration(ts,NULL,&tfinal);
157: while (user->next_output <= t && user->next_output <= tfinal) {
158: VecDuplicate(X,&interpolatedX);
159: TSInterpolate(ts,user->next_output,interpolatedX);
160: VecGetArrayRead(interpolatedX,&x);
161: 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]));
162: VecRestoreArrayRead(interpolatedX,&x);
163: VecDestroy(&interpolatedX);
164: user->next_output += 0.1;
165: }
166: return(0);
167: }
171: int main(int argc,char **argv)
172: {
173: TS ts; /* nonlinear solver */
174: Vec x; /* solution, residual vectors */
175: Mat A; /* Jacobian matrix */
176: PetscInt steps;
177: PetscReal ftime = 0.5;
178: PetscBool monitor = PETSC_FALSE;
179: PetscScalar *x_ptr;
180: PetscMPIInt size;
181: struct _n_User user;
184: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: Initialize program
186: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187: PetscInitialize(&argc,&argv,NULL,help);
189: MPI_Comm_size(PETSC_COMM_WORLD,&size);
190: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
192: RegisterMyARK2();
194: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195: Set runtime options
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198: user.next_output = 0.0;
199: PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
201: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202: Create necessary matrix and vectors, solve same ODE on every process
203: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204: MatCreate(PETSC_COMM_WORLD,&A);
205: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
206: MatSetFromOptions(A);
207: MatSetUp(A);
208: MatGetVecs(A,&x,NULL);
210: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: Create timestepping solver context
212: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213: TSCreate(PETSC_COMM_WORLD,&ts);
214: TSSetType(ts,TSBEULER);
215: TSSetIFunction(ts,NULL,IFunction,&user);
216: TSSetIJacobian(ts,A,A,IJacobian,&user);
217: TSSetDuration(ts,PETSC_DEFAULT,ftime);
218: if (monitor) {
219: TSMonitorSet(ts,Monitor,&user,NULL);
220: }
222: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223: Set initial conditions
224: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225: VecGetArray(x,&x_ptr);
226: x_ptr[0] = -2; x_ptr[1] = -2.355301397608119909925287735864250951918;
227: VecRestoreArray(x,&x_ptr);
228: TSSetInitialTimeStep(ts,0.0,.001);
230: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231: Set runtime options
232: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233: TSSetFromOptions(ts);
235: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236: Solve nonlinear system
237: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238: TSSolve(ts,x);
239: TSGetSolveTime(ts,&ftime);
240: TSGetTimeStepNumber(ts,&steps);
241: PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
242: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
244: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245: Free work space. All PETSc objects should be destroyed when they
246: are no longer needed.
247: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
248: MatDestroy(&A);
249: VecDestroy(&x);
250: TSDestroy(&ts);
252: PetscFinalize();
253: return(0);
254: }