Actual source code: ex16fwd.c
1: static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\
2: Input parameters include:\n\
3: -mu : stiffness parameter\n\n";
5: /* ------------------------------------------------------------------------
7: This program solves the van der Pol equation
8: y'' - \mu (1-y^2)*y' + y = 0 (1)
9: on the domain 0 <= x <= 1, with the boundary conditions
10: y(0) = 2, y'(0) = 0,
11: and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with an explicit Runge-Kutta method and its discrete tangent linear model.
13: Notes:
14: This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t).
16: (1) can be turned into a system of first order ODEs
17: [ y' ] = [ z ]
18: [ z' ] [ \mu (1 - y^2) z - y ]
20: which then we can write as a vector equation
22: [ u_1' ] = [ u_2 ] (2)
23: [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ]
25: which is now in the form of u_t = F(u,t).
27: The user provides the right-hand-side function
29: [ f(u,t) ] = [ u_2 ]
30: [ \mu (1 - u_1^2) u_2 - u_1 ]
32: the Jacobian function
34: df [ 0 ; 1 ]
35: -- = [ ]
36: du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ]
38: and the JacobainP (the Jacobian w.r.t. parameter) function
40: df [ 0; 0; 0 ]
41: --- = [ ]
42: d\mu [ 0; 0; (1 - u_1^2) u_2 ]
44: ------------------------------------------------------------------------- */
46: #include <petscts.h>
47: #include <petscmat.h>
48: typedef struct _n_User *User;
49: struct _n_User {
50: PetscReal mu;
51: PetscReal next_output;
52: PetscReal tprev;
53: };
55: /*
56: User-defined routines
57: */
58: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
59: {
60: User user = (User)ctx;
61: PetscScalar *f;
62: const PetscScalar *x;
65: VecGetArrayRead(X,&x);
66: VecGetArray(F,&f);
67: f[0] = x[1];
68: f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
69: VecRestoreArrayRead(X,&x);
70: VecRestoreArray(F,&f);
71: return 0;
72: }
74: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
75: {
76: User user = (User)ctx;
77: PetscReal mu = user->mu;
78: PetscInt rowcol[] = {0,1};
79: PetscScalar J[2][2];
80: const PetscScalar *x;
83: VecGetArrayRead(X,&x);
84: J[0][0] = 0;
85: J[1][0] = -2.*mu*x[1]*x[0]-1.;
86: J[0][1] = 1.0;
87: J[1][1] = mu*(1.0-x[0]*x[0]);
88: MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
89: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
91: if (A != B) {
92: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
93: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
94: }
95: VecRestoreArrayRead(X,&x);
96: return 0;
97: }
99: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
100: {
101: PetscInt row[] = {0,1},col[]={2};
102: PetscScalar J[2][1];
103: const PetscScalar *x;
106: VecGetArrayRead(X,&x);
107: J[0][0] = 0;
108: J[1][0] = (1.-x[0]*x[0])*x[1];
109: VecRestoreArrayRead(X,&x);
110: MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
112: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
113: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
114: return 0;
115: }
117: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
118: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
119: {
120: const PetscScalar *x;
121: PetscReal tfinal, dt, tprev;
122: User user = (User)ctx;
125: TSGetTimeStep(ts,&dt);
126: TSGetMaxTime(ts,&tfinal);
127: TSGetPrevTime(ts,&tprev);
128: VecGetArrayRead(X,&x);
129: PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",(double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),(double)PetscRealPart(x[1]));
130: PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
131: VecRestoreArrayRead(X,&x);
132: return 0;
133: }
135: int main(int argc,char **argv)
136: {
137: TS ts; /* nonlinear solver */
138: Vec x; /* solution, residual vectors */
139: Mat A; /* Jacobian matrix */
140: Mat Jacp; /* JacobianP matrix */
141: PetscInt steps;
142: PetscReal ftime =0.5;
143: PetscBool monitor = PETSC_FALSE;
144: PetscScalar *x_ptr;
145: PetscMPIInt size;
146: struct _n_User user;
147: Mat sp;
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Initialize program
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152: PetscInitialize(&argc,&argv,NULL,help);
153: MPI_Comm_size(PETSC_COMM_WORLD,&size);
156: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: Set runtime options
158: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159: user.mu = 1;
160: user.next_output = 0.0;
162: PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
163: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
165: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166: Create necessary matrix and vectors, solve same ODE on every process
167: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168: MatCreate(PETSC_COMM_WORLD,&A);
169: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
170: MatSetFromOptions(A);
171: MatSetUp(A);
172: MatCreateVecs(A,&x,NULL);
174: MatCreate(PETSC_COMM_WORLD,&Jacp);
175: MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3);
176: MatSetFromOptions(Jacp);
177: MatSetUp(Jacp);
179: MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp);
180: MatZeroEntries(sp);
181: MatShift(sp,1.0);
183: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184: Create timestepping solver context
185: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186: TSCreate(PETSC_COMM_WORLD,&ts);
187: TSSetType(ts,TSRK);
188: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
189: /* Set RHS Jacobian for the adjoint integration */
190: TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);
191: TSSetMaxTime(ts,ftime);
192: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
193: if (monitor) {
194: TSMonitorSet(ts,Monitor,&user,NULL);
195: }
196: TSForwardSetSensitivities(ts,3,sp);
197: TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);
199: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200: Set initial conditions
201: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202: VecGetArray(x,&x_ptr);
204: x_ptr[0] = 2; x_ptr[1] = 0.66666654321;
205: VecRestoreArray(x,&x_ptr);
206: TSSetTimeStep(ts,.001);
208: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209: Set runtime options
210: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211: TSSetFromOptions(ts);
213: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: Solve nonlinear system
215: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216: TSSolve(ts,x);
217: TSGetSolveTime(ts,&ftime);
218: TSGetStepNumber(ts,&steps);
219: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);
220: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
222: PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n");
223: MatView(sp,PETSC_VIEWER_STDOUT_WORLD);
225: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226: Free work space. All PETSc objects should be destroyed when they
227: are no longer needed.
228: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229: MatDestroy(&A);
230: MatDestroy(&Jacp);
231: VecDestroy(&x);
232: MatDestroy(&sp);
233: TSDestroy(&ts);
234: PetscFinalize();
235: return 0;
236: }
238: /*TEST
240: test:
241: args: -monitor 0 -ts_adapt_type none
243: TEST*/