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: /*
6: Concepts: TS^time-dependent nonlinear problems
7: Concepts: TS^van der Pol equation
8: Concepts: TS^adjoint sensitivity analysis
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: 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.
19: Notes:
20: This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t).
22: (1) can be turned into a system of first order ODEs
23: [ y' ] = [ z ]
24: [ z' ] [ \mu (1 - y^2) z - y ]
26: which then we can write as a vector equation
28: [ u_1' ] = [ u_2 ] (2)
29: [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ]
31: which is now in the form of u_t = F(u,t).
33: The user provides the right-hand-side function
35: [ f(u,t) ] = [ u_2 ]
36: [ \mu (1 - u_1^2) u_2 - u_1 ]
38: the Jacobian function
40: df [ 0 ; 1 ]
41: -- = [ ]
42: du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ]
44: and the JacobainP (the Jacobian w.r.t. parameter) function
46: df [ 0; 0; 0 ]
47: --- = [ ]
48: d\mu [ 0; 0; (1 - u_1^2) u_2 ]
50: ------------------------------------------------------------------------- */
52: #include <petscts.h>
53: #include <petscmat.h>
54: typedef struct _n_User *User;
55: struct _n_User {
56: PetscReal mu;
57: PetscReal next_output;
58: PetscReal tprev;
59: };
61: /*
62: User-defined routines
63: */
64: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
65: {
66: PetscErrorCode ierr;
67: User user = (User)ctx;
68: PetscScalar *f;
69: const PetscScalar *x;
72: VecGetArrayRead(X,&x);
73: VecGetArray(F,&f);
74: f[0] = x[1];
75: f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
76: VecRestoreArrayRead(X,&x);
77: VecRestoreArray(F,&f);
78: return(0);
79: }
81: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
82: {
83: PetscErrorCode ierr;
84: User user = (User)ctx;
85: PetscReal mu = user->mu;
86: PetscInt rowcol[] = {0,1};
87: PetscScalar J[2][2];
88: const PetscScalar *x;
91: VecGetArrayRead(X,&x);
92: J[0][0] = 0;
93: J[1][0] = -2.*mu*x[1]*x[0]-1.;
94: J[0][1] = 1.0;
95: J[1][1] = mu*(1.0-x[0]*x[0]);
96: MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
97: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
98: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
99: if (A != B) {
100: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
101: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
102: }
103: VecRestoreArrayRead(X,&x);
104: return(0);
105: }
107: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
108: {
109: PetscErrorCode ierr;
110: PetscInt row[] = {0,1},col[]={2};
111: PetscScalar J[2][1];
112: const PetscScalar *x;
115: VecGetArrayRead(X,&x);
116: J[0][0] = 0;
117: J[1][0] = (1.-x[0]*x[0])*x[1];
118: VecRestoreArrayRead(X,&x);
119: MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
121: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
122: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
123: return(0);
124: }
126: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
127: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
128: {
129: PetscErrorCode ierr;
130: const PetscScalar *x;
131: PetscReal tfinal, dt, tprev;
132: User user = (User)ctx;
135: TSGetTimeStep(ts,&dt);
136: TSGetMaxTime(ts,&tfinal);
137: TSGetPrevTime(ts,&tprev);
138: VecGetArrayRead(X,&x);
139: 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]));
140: PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
141: VecRestoreArrayRead(X,&x);
142: return(0);
143: }
145: int main(int argc,char **argv)
146: {
147: TS ts; /* nonlinear solver */
148: Vec x; /* solution, residual vectors */
149: Mat A; /* Jacobian matrix */
150: Mat Jacp; /* JacobianP matrix */
151: PetscInt steps;
152: PetscReal ftime =0.5;
153: PetscBool monitor = PETSC_FALSE;
154: PetscScalar *x_ptr;
155: PetscMPIInt size;
156: struct _n_User user;
158: Mat sp;
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Initialize program
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
164: MPI_Comm_size(PETSC_COMM_WORLD,&size);
165: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
167: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168: Set runtime options
169: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: user.mu = 1;
171: user.next_output = 0.0;
173: PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
174: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
176: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177: Create necessary matrix and vectors, solve same ODE on every process
178: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
179: MatCreate(PETSC_COMM_WORLD,&A);
180: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
181: MatSetFromOptions(A);
182: MatSetUp(A);
183: MatCreateVecs(A,&x,NULL);
185: MatCreate(PETSC_COMM_WORLD,&Jacp);
186: MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3);
187: MatSetFromOptions(Jacp);
188: MatSetUp(Jacp);
190: MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp);
191: MatZeroEntries(sp);
192: MatShift(sp,1.0);
194: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195: Create timestepping solver context
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197: TSCreate(PETSC_COMM_WORLD,&ts);
198: TSSetType(ts,TSRK);
199: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
200: /* Set RHS Jacobian for the adjoint integration */
201: TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);
202: TSSetMaxTime(ts,ftime);
203: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
204: if (monitor) {
205: TSMonitorSet(ts,Monitor,&user,NULL);
206: }
207: TSForwardSetSensitivities(ts,3,sp);
208: TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);
210: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: Set initial conditions
212: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213: VecGetArray(x,&x_ptr);
215: x_ptr[0] = 2; x_ptr[1] = 0.66666654321;
216: VecRestoreArray(x,&x_ptr);
217: TSSetTimeStep(ts,.001);
219: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220: Set runtime options
221: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
222: TSSetFromOptions(ts);
224: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: Solve nonlinear system
226: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227: TSSolve(ts,x);
228: TSGetSolveTime(ts,&ftime);
229: TSGetStepNumber(ts,&steps);
230: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);
231: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
233: PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n");
234: MatView(sp,PETSC_VIEWER_STDOUT_WORLD);
236: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237: Free work space. All PETSc objects should be destroyed when they
238: are no longer needed.
239: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240: MatDestroy(&A);
241: MatDestroy(&Jacp);
242: VecDestroy(&x);
243: MatDestroy(&sp);
244: TSDestroy(&ts);
245: PetscFinalize();
246: return ierr;
247: }
249: /*TEST
251: test:
252: args: -monitor 0 -ts_adapt_type none
254: TEST*/