Actual source code: ex16fwd.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  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  ]


 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[]={2};
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:   VecRestoreArrayRead(X,&x);
120:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

122:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
123:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
124:   return(0);
125: }

127: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
128: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
129: {
130:   PetscErrorCode    ierr;
131:   const PetscScalar *x;
132:   PetscReal         tfinal, dt, tprev;
133:   User              user = (User)ctx;

136:   TSGetTimeStep(ts,&dt);
137:   TSGetMaxTime(ts,&tfinal);
138:   TSGetPrevTime(ts,&tprev);
139:   VecGetArrayRead(X,&x);
140:   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]));
141:   PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
142:   VecRestoreArrayRead(X,&x);
143:   return(0);
144: }

146: int main(int argc,char **argv)
147: {
148:   TS             ts;            /* nonlinear solver */
149:   Vec            x;             /* solution, residual vectors */
150:   Mat            A;             /* Jacobian matrix */
151:   Mat            Jacp;          /* JacobianP matrix */
152:   PetscInt       steps;
153:   PetscReal      ftime   =0.5;
154:   PetscBool      monitor = PETSC_FALSE;
155:   PetscScalar    *x_ptr;
156:   PetscMPIInt    size;
157:   struct _n_User user;
159:   Mat            sp;

161:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162:      Initialize program
163:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164:   PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
165:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
166:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

168:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169:     Set runtime options
170:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
171:   user.mu          = 1;
172:   user.next_output = 0.0;


175:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
176:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);

178:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179:     Create necessary matrix and vectors, solve same ODE on every process
180:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181:   MatCreate(PETSC_COMM_WORLD,&A);
182:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
183:   MatSetFromOptions(A);
184:   MatSetUp(A);
185:   MatCreateVecs(A,&x,NULL);

187:   MatCreate(PETSC_COMM_WORLD,&Jacp);
188:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3);
189:   MatSetFromOptions(Jacp);
190:   MatSetUp(Jacp);

192:   MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp);
193:   MatZeroEntries(sp);
194:   MatShift(sp,1.0);

196:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197:      Create timestepping solver context
198:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
199:   TSCreate(PETSC_COMM_WORLD,&ts);
200:   TSSetType(ts,TSRK);
201:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
202:   /*   Set RHS Jacobian for the adjoint integration */
203:   TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);
204:   TSSetMaxTime(ts,ftime);
205:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
206:   if (monitor) {
207:     TSMonitorSet(ts,Monitor,&user,NULL);
208:   }
209:   TSForwardSetSensitivities(ts,3,sp);
210:   TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);

212:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213:      Set initial conditions
214:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
215:   VecGetArray(x,&x_ptr);

217:   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
218:   VecRestoreArray(x,&x_ptr);
219:   TSSetTimeStep(ts,.001);


222:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223:      Set runtime options
224:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225:   TSSetFromOptions(ts);

227:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228:      Solve nonlinear system
229:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230:   TSSolve(ts,x);
231:   TSGetSolveTime(ts,&ftime);
232:   TSGetStepNumber(ts,&steps);
233:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);
234:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

236:   PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n");
237:   MatView(sp,PETSC_VIEWER_STDOUT_WORLD);

239:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
240:      Free work space.  All PETSc objects should be destroyed when they
241:      are no longer needed.
242:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243:   MatDestroy(&A);
244:   MatDestroy(&Jacp);
245:   VecDestroy(&x);
246:   MatDestroy(&sp);
247:   TSDestroy(&ts);
248:   PetscFinalize();
249:   return ierr;
250: }

252: /*TEST

254:     test:
255:       args: -monitor 0 -ts_adapt_type none

257: TEST*/