Actual source code: ex16adj.c
petsc-3.11.4 2019-09-28
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 adjoint.
19: Notes:
20: This code demonstrates the TSAdjoint 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 JacobianP (the Jacobian w.r.t. parameter) function
46: dF [ 0 ]
47: --- = [ ]
48: d\mu [ (1 - u_1^2) u_2 ]
51: ------------------------------------------------------------------------- */
53: #include <petscts.h>
54: #include <petscmat.h>
55: typedef struct _n_User *User;
56: struct _n_User {
57: PetscReal mu;
58: PetscReal next_output;
59: PetscReal tprev;
60: };
62: /*
63: * User-defined routines
64: */
65: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
66: {
67: PetscErrorCode ierr;
68: User user = (User)ctx;
69: PetscScalar *f;
70: const PetscScalar *x;
73: VecGetArrayRead(X,&x);
74: VecGetArray(F,&f);
75: f[0] = x[1];
76: f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
77: VecRestoreArrayRead(X,&x);
78: VecRestoreArray(F,&f);
79: return(0);
80: }
82: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
83: {
84: PetscErrorCode ierr;
85: User user = (User)ctx;
86: PetscReal mu = user->mu;
87: PetscInt rowcol[] = {0,1};
88: PetscScalar J[2][2];
89: const PetscScalar *x;
92: VecGetArrayRead(X,&x);
93: J[0][0] = 0;
94: J[1][0] = -2.*mu*x[1]*x[0]-1.;
95: J[0][1] = 1.0;
96: J[1][1] = mu*(1.0-x[0]*x[0]);
97: MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
98: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
99: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
100: if (A != B) {
101: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
102: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
103: }
104: VecRestoreArrayRead(X,&x);
105: return(0);
106: }
108: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
109: {
110: PetscErrorCode ierr;
111: PetscInt row[] = {0,1},col[]={0};
112: PetscScalar J[2][1];
113: const PetscScalar *x;
116: VecGetArrayRead(X,&x);
117: J[0][0] = 0;
118: J[1][0] = (1.-x[0]*x[0])*x[1];
119: MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
120: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
121: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
122: VecRestoreArrayRead(X,&x);
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: Vec lambda[2],mu[2];
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_SELF,1,"This is a uniprocessor example only!");
167: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168: Set runtime options
169: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: user.mu = 1;
171: user.next_output = 0.0;
174: PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
175: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
177: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: Create necessary matrix and vectors, solve same ODE on every process
179: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: MatCreate(PETSC_COMM_WORLD,&A);
181: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
182: MatSetFromOptions(A);
183: MatSetUp(A);
184: MatCreateVecs(A,&x,NULL);
186: MatCreate(PETSC_COMM_WORLD,&Jacp);
187: MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
188: MatSetFromOptions(Jacp);
189: MatSetUp(Jacp);
191: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192: Create timestepping solver context
193: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
194: TSCreate(PETSC_COMM_WORLD,&ts);
195: TSSetType(ts,TSRK);
196: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
197: /* Set RHS Jacobian for the adjoint integration */
198: TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);
199: TSSetMaxTime(ts,ftime);
200: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
201: if (monitor) {
202: TSMonitorSet(ts,Monitor,&user,NULL);
203: }
205: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: Set initial conditions
207: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208: VecGetArray(x,&x_ptr);
210: x_ptr[0] = 2; x_ptr[1] = 0.66666654321;
211: VecRestoreArray(x,&x_ptr);
212: TSSetTimeStep(ts,.001);
214: /*
215: Have the TS save its trajectory so that TSAdjointSolve() may be used
216: */
217: TSSetSaveTrajectory(ts);
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: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234: Start the Adjoint model
235: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236: MatCreateVecs(A,&lambda[0],NULL);
237: MatCreateVecs(A,&lambda[1],NULL);
238: /* Reset initial conditions for the adjoint integration */
239: VecGetArray(lambda[0],&x_ptr);
240: x_ptr[0] = 1.0; x_ptr[1] = 0.0;
241: VecRestoreArray(lambda[0],&x_ptr);
242: VecGetArray(lambda[1],&x_ptr);
243: x_ptr[0] = 0.0; x_ptr[1] = 1.0;
244: VecRestoreArray(lambda[1],&x_ptr);
246: MatCreateVecs(Jacp,&mu[0],NULL);
247: MatCreateVecs(Jacp,&mu[1],NULL);
248: VecGetArray(mu[0],&x_ptr);
249: x_ptr[0] = 0.0;
250: VecRestoreArray(mu[0],&x_ptr);
251: VecGetArray(mu[1],&x_ptr);
252: x_ptr[0] = 0.0;
253: VecRestoreArray(mu[1],&x_ptr);
254: TSSetCostGradients(ts,2,lambda,mu);
257: /* Set RHS JacobianP */
258: TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);
260: TSAdjointSolve(ts);
262: VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
263: VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD);
264: VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
265: VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD);
267: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268: Free work space. All PETSc objects should be destroyed when they
269: are no longer needed.
270: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271: MatDestroy(&A);
272: MatDestroy(&Jacp);
273: VecDestroy(&x);
274: VecDestroy(&lambda[0]);
275: VecDestroy(&lambda[1]);
276: VecDestroy(&mu[0]);
277: VecDestroy(&mu[1]);
278: TSDestroy(&ts);
279: PetscFinalize();
280: return ierr;
281: }
283: /*TEST
285: test:
286: args: -monitor 0 -viewer_binary_skip_info -ts_trajectory_dirname ex16adjdir
288: test:
289: suffix: 2
290: args: -monitor 0 -ts_trajectory_type memory
292: TEST*/