Actual source code: ex16adj.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: /*
  2:    Second step toward adjoint ERK solver
  3:    Features:
  4:      1. checkpointed data stored in the RK context;
  5:      2. adjoint solver is implemented inside TS;
  6:      3. support arrays of lambda and mu.
  7: */
  8: static char help[] = "Solves the van der Pol equation.\n\
  9: Input parameters include:\n\
 10:       -mu : stiffness parameter\n\n";

 12: /*
 13:    Concepts: TS^time-dependent nonlinear problems
 14:    Concepts: TS^van der Pol equation
 15:    Processors: 1
 16: */
 17: /* ------------------------------------------------------------------------

 19:    This program solves the van der Pol equation
 20:        y'' - \mu (1-y^2)*y' + y = 0        (1)
 21:    on the domain 0 <= x <= 1, with the boundary conditions
 22:        y(0) = 2, y'(0) = 0,
 23:    This is a nonlinear equation.

 25:    Notes:
 26:    This code demonstrates the TS solver interface to two variants of
 27:    linear problems, u_t = f(u,t), namely turning (1) into a system of
 28:    first order differential equations,

 30:    [ y' ] = [          z          ]
 31:    [ z' ]   [ \mu (1 - y^2) z - y ]

 33:    which then we can write as a vector equation

 35:    [ u_1' ] = [             u_2           ]  (2)
 36:    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]

 38:    which is now in the desired form of u_t = f(u,t). 
 39:   ------------------------------------------------------------------------- */

 41: #include <petscts.h>
 42: #include <petscmat.h>
 43: typedef struct _n_User *User;
 44: struct _n_User {
 45:   PetscReal mu;
 46:   PetscReal next_output;
 47:   PetscReal tprev;
 48: };

 50: /*
 51: *  User-defined routines
 52: */
 55: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 56: {
 57:   PetscErrorCode    ierr;
 58:   User              user = (User)ctx;
 59:   PetscScalar       *f;
 60:   const PetscScalar *x;

 63:   VecGetArrayRead(X,&x);
 64:   VecGetArray(F,&f);
 65:   f[0] = x[1];
 66:   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
 67:   VecRestoreArrayRead(X,&x);
 68:   VecRestoreArray(F,&f);
 69:   return(0);
 70: }

 74: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
 75: {
 76:   PetscErrorCode    ierr;
 77:   User              user = (User)ctx;
 78:   PetscReal         mu   = user->mu;
 79:   PetscInt          rowcol[] = {0,1};
 80:   PetscScalar       J[2][2];
 81:   const PetscScalar *x;

 84:   VecGetArrayRead(X,&x);
 85:   J[0][0] = 0;
 86:   J[1][0] = -2.*mu*x[1]*x[0]-1.;
 87:   J[0][1] = 1.0;
 88:   J[1][1] = mu*(1.0-x[0]*x[0]);
 89:   MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 90:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 91:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 92:   if (A != B) {
 93:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 94:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 95:   }
 96:   VecRestoreArrayRead(X,&x);
 97:   return(0);
 98: }

102: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
103: {
104:   PetscErrorCode    ierr;
105:   PetscInt          row[] = {0,1},col[]={0};
106:   PetscScalar       J[2][1];
107:   const PetscScalar *x;

110:   VecGetArrayRead(X,&x);
111:   J[0][0] = 0;
112:   J[1][0] = (1.-x[0]*x[0])*x[1];
113:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
114:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
115:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
116:   VecRestoreArrayRead(X,&x);
117:   return(0);
118: }

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

131:   TSGetTimeStep(ts,&dt);
132:   TSGetDuration(ts,NULL,&tfinal);
133:   TSGetPrevTime(ts,&tprev);
134:   VecGetArrayRead(X,&x);
135:   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]));
136:   PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
137:   VecRestoreArrayRead(X,&x);
138:   return(0);
139: }

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

158:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159:      Initialize program
160:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161:   PetscInitialize(&argc,&argv,NULL,help);

163:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
164:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

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


173:   PetscOptionsGetReal(NULL,"-mu",&user.mu,NULL);
174:   PetscOptionsGetBool(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,1);
187:   MatSetFromOptions(Jacp);
188:   MatSetUp(Jacp);
189: 
190:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191:      Create timestepping solver context
192:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193:   TSCreate(PETSC_COMM_WORLD,&ts);
194:   TSSetType(ts,TSRK);
195:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
196:   TSSetDuration(ts,PETSC_DEFAULT,ftime);
197:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
198:   if (monitor) {
199:     TSMonitorSet(ts,Monitor,&user,NULL);
200:   }
201: 
202:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203:      Set initial conditions
204:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205:   VecGetArray(x,&x_ptr);

207:   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
208:   VecRestoreArray(x,&x_ptr);
209:   TSSetInitialTimeStep(ts,0.0,.001);

211:   /*
212:     Have the TS save its trajectory so that TSAdjointSolve() may be used
213:   */
214:   TSSetSaveTrajectory(ts);

216:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217:      Set runtime options
218:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219:   TSSetFromOptions(ts);

221:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222:      Solve nonlinear system
223:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224:   TSSolve(ts,x);
225:   TSGetSolveTime(ts,&ftime);
226:   TSGetTimeStepNumber(ts,&steps);
227:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);
228:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

230:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231:      Start the Adjoint model
232:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233:   MatCreateVecs(A,&lambda[0],NULL);
234:   MatCreateVecs(A,&lambda[1],NULL);
235:   /*   Reset initial conditions for the adjoint integration */
236:   VecGetArray(lambda[0],&x_ptr);
237:   x_ptr[0] = 1.0;   x_ptr[1] = 0.0;
238:   VecRestoreArray(lambda[0],&x_ptr);
239:   VecGetArray(lambda[1],&x_ptr);
240:   x_ptr[0] = 0.0;   x_ptr[1] = 1.0;
241:   VecRestoreArray(lambda[1],&x_ptr);

243:   MatCreateVecs(Jacp,&mu[0],NULL);
244:   MatCreateVecs(Jacp,&mu[1],NULL);
245:   VecGetArray(mu[0],&x_ptr);
246:   x_ptr[0] = 0.0;
247:   VecRestoreArray(mu[0],&x_ptr);
248:   VecGetArray(mu[1],&x_ptr);
249:   x_ptr[0] = 0.0;
250:   VecRestoreArray(mu[1],&x_ptr);
251:   TSSetCostGradients(ts,2,lambda,mu);

253:   /*   Set RHS Jacobian for the adjoint integration */
254:   TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);

256:   /*   Set RHS JacobianP */
257:   TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,&user);

259:   TSAdjointSolve(ts);

261:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
262:   VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD);
263:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
264:   VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD);

266:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267:      Free work space.  All PETSc objects should be destroyed when they
268:      are no longer needed.
269:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
270:   MatDestroy(&A);
271:   MatDestroy(&Jacp);
272:   VecDestroy(&x);
273:   VecDestroy(&lambda[0]);
274:   VecDestroy(&lambda[1]);
275:   VecDestroy(&mu[0]);
276:   VecDestroy(&mu[1]);
277:   TSDestroy(&ts);

279:   PetscFinalize();
280:   return(0);
281: }