Actual source code: parabolic.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: #include <petsc/private/taoimpl.h>

  3: /*T
  4:    Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
  5:    Routines: TaoCreate();
  6:    Routines: TaoSetType();
  7:    Routines: TaoSetInitialVector();
  8:    Routines: TaoSetObjectiveRoutine();
  9:    Routines: TaoSetGradientRoutine();
 10:    Routines: TaoSetConstraintsRoutine();
 11:    Routines: TaoSetJacobianStateRoutine();
 12:    Routines: TaoSetJacobianDesignRoutine();
 13:    Routines: TaoSetStateDesignIS();
 14:    Routines: TaoSetFromOptions();
 15:    Routines: TaoSetHistory(); TaoGetHistory();
 16:    Routines: TaoSolve();
 17:    Routines: TaoGetConvergedReason(); TaoDestroy();
 18:    Processors: n
 19: T*/

 21: typedef struct {
 22:   PetscInt n; /*  Number of variables */
 23:   PetscInt m; /*  Number of constraints per time step */
 24:   PetscInt mx; /*  grid points in each direction */
 25:   PetscInt nt; /*  Number of time steps; as of now, must be divisible by 8 */
 26:   PetscInt ndata; /*  Number of data points per sample */
 27:   PetscInt ns; /*  Number of samples */
 28:   PetscInt *sample_times; /*  Times of samples */
 29:   IS       s_is;
 30:   IS       d_is;

 32:   VecScatter state_scatter;
 33:   VecScatter design_scatter;
 34:   VecScatter *yi_scatter;
 35:   VecScatter *di_scatter;

 37:   Mat       Js,Jd,JsBlockPrec,JsInv,JsBlock;
 38:   PetscBool jformed,dsg_formed;

 40:   PetscReal alpha; /*  Regularization parameter */
 41:   PetscReal beta; /*  Weight attributed to ||u||^2 in regularization functional */
 42:   PetscReal noise; /*  Amount of noise to add to data */
 43:   PetscReal ht; /*  Time step */

 45:   Mat Qblock,QblockT;
 46:   Mat L,LT;
 47:   Mat Div,Divwork;
 48:   Mat Grad;
 49:   Mat Av,Avwork,AvT;
 50:   Mat DSG;
 51:   Vec q;
 52:   Vec ur; /*  reference */

 54:   Vec d;
 55:   Vec dwork;
 56:   Vec *di;

 58:   Vec y; /*  state variables */
 59:   Vec ywork;

 61:   Vec ytrue;
 62:   Vec *yi,*yiwork;

 64:   Vec u; /*  design variables */
 65:   Vec uwork;

 67:   Vec utrue;
 68:   Vec js_diag;
 69:   Vec c; /*  constraint vector */
 70:   Vec cwork;

 72:   Vec lwork;
 73:   Vec S;
 74:   Vec Rwork,Swork,Twork;
 75:   Vec Av_u;

 77:   KSP solver;
 78:   PC prec;

 80:   PetscInt ksp_its;
 81:   PetscInt ksp_its_initial;
 82: } AppCtx;


 85: PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*);
 86: PetscErrorCode FormGradient(Tao, Vec, Vec, void*);
 87: PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*);
 88: PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void*);
 89: PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void*);
 90: PetscErrorCode FormConstraints(Tao, Vec, Vec, void*);
 91: PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void*);
 92: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
 93: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
 94: PetscErrorCode ParabolicInitialize(AppCtx *user);
 95: PetscErrorCode ParabolicDestroy(AppCtx *user);
 96: PetscErrorCode ParabolicMonitor(Tao, void*);

 98: PetscErrorCode StateMatMult(Mat,Vec,Vec);
 99: PetscErrorCode StateMatBlockMult(Mat,Vec,Vec);
100: PetscErrorCode StateMatMultTranspose(Mat,Vec,Vec);
101: PetscErrorCode StateMatGetDiagonal(Mat,Vec);
102: PetscErrorCode StateMatDuplicate(Mat,MatDuplicateOption,Mat*);
103: PetscErrorCode StateMatInvMult(Mat,Vec,Vec);
104: PetscErrorCode StateMatInvTransposeMult(Mat,Vec,Vec);
105: PetscErrorCode StateMatBlockPrecMult(PC,Vec,Vec);

107: PetscErrorCode DesignMatMult(Mat,Vec,Vec);
108: PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);

110: PetscErrorCode Gather_i(Vec,Vec*,VecScatter*,PetscInt);
111: PetscErrorCode Scatter_i(Vec,Vec*,VecScatter*,PetscInt);

113: static  char help[]="";

117: int main(int argc, char **argv)
118: {
119:   PetscErrorCode     ierr;
120:   Vec                x,x0;
121:   Tao                tao;
122:   TaoConvergedReason reason;
123:   AppCtx             user;
124:   IS                 is_allstate,is_alldesign;
125:   PetscInt           lo,hi,hi2,lo2,ksp_old;
126:   PetscInt           ntests = 1;
127:   PetscInt           i;
128:   int                stages[1];

130:   PetscInitialize(&argc, &argv, (char*)0,help);

132:   user.mx = 8;
133:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,NULL,NULL);
134:   PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL);
135:   user.nt = 8;
136:   PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL);
137:   user.ndata = 64;
138:   PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL);
139:   user.ns = 8;
140:   PetscOptionsInt("-ns","Number of samples","",user.ns,&user.ns,NULL);
141:   user.alpha = 1.0;
142:   PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL);
143:   user.beta = 0.01;
144:   PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,NULL);
145:   user.noise = 0.01;
146:   PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,NULL);

148:   user.m = user.mx*user.mx*user.mx; /*  number of constraints per time step */
149:   user.n = user.m*(user.nt+1); /*  number of variables */
150:   user.ht = (PetscReal)1/user.nt;

152:   VecCreate(PETSC_COMM_WORLD,&user.u);
153:   VecCreate(PETSC_COMM_WORLD,&user.y);
154:   VecCreate(PETSC_COMM_WORLD,&user.c);
155:   VecSetSizes(user.u,PETSC_DECIDE,user.n-user.m*user.nt);
156:   VecSetSizes(user.y,PETSC_DECIDE,user.m*user.nt);
157:   VecSetSizes(user.c,PETSC_DECIDE,user.m*user.nt);
158:   VecSetFromOptions(user.u);
159:   VecSetFromOptions(user.y);
160:   VecSetFromOptions(user.c);

162:   /* Create scatters for reduced spaces.
163:      If the state vector y and design vector u are partitioned as
164:      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
165:      then the solution vector x is organized as
166:      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
167:      The index sets user.s_is and user.d_is correspond to the indices of the
168:      state and design variables owned by the current processor.
169:   */
170:   VecCreate(PETSC_COMM_WORLD,&x);

172:   VecGetOwnershipRange(user.y,&lo,&hi);
173:   VecGetOwnershipRange(user.u,&lo2,&hi2);

175:   ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate);
176:   ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user.s_is);

178:   ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign);
179:   ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user.d_is);

181:   VecSetSizes(x,hi-lo+hi2-lo2,user.n);
182:   VecSetFromOptions(x);

184:   VecScatterCreate(x,user.s_is,user.y,is_allstate,&user.state_scatter);
185:   VecScatterCreate(x,user.d_is,user.u,is_alldesign,&user.design_scatter);
186:   ISDestroy(&is_alldesign);
187:   ISDestroy(&is_allstate);

189:   /* Create TAO solver and set desired solution method */
190:   TaoCreate(PETSC_COMM_WORLD,&tao);
191:   TaoSetType(tao,TAOLCL);

193:   /* Set up initial vectors and matrices */
194:   ParabolicInitialize(&user);

196:   Gather(x,user.y,user.state_scatter,user.u,user.design_scatter);
197:   VecDuplicate(x,&x0);
198:   VecCopy(x,x0);

200:   /* Set solution vector with an initial guess */
201:   TaoSetInitialVector(tao,x);
202:   TaoSetObjectiveRoutine(tao, FormFunction, (void *)&user);
203:   TaoSetGradientRoutine(tao, FormGradient, (void *)&user);
204:   TaoSetConstraintsRoutine(tao, user.c, FormConstraints, (void *)&user);

206:   TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, (void *)&user);
207:   TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, (void *)&user);

209:   TaoSetFromOptions(tao);
210:   TaoSetStateDesignIS(tao,user.s_is,user.d_is);

212:  /* SOLVE THE APPLICATION */
213:   PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL);
214:   PetscOptionsEnd();
215:   PetscLogStageRegister("Trials",&stages[0]);
216:   PetscLogStagePush(stages[0]);
217:   user.ksp_its_initial = user.ksp_its;
218:   ksp_old = user.ksp_its;
219:   for (i=0; i<ntests; i++){
220:     TaoSolve(tao);
221:     PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its-ksp_old);
222:     VecCopy(x0,x);
223:     TaoSetInitialVector(tao,x);
224:   }
225:   PetscLogStagePop();
226:   PetscBarrier((PetscObject)x);
227:   PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ");
228:   PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial);
229:   PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %D trial(s): ",ntests);
230:   PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its);
231:   PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: ");
232:   PetscPrintf(PETSC_COMM_WORLD,"%D\n",(user.ksp_its-user.ksp_its_initial)/ntests);
233:   TaoGetConvergedReason(tao,&reason);

235:   if (reason < 0) {
236:     PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge.\n");
237:   } else {
238:     PetscPrintf(MPI_COMM_WORLD, "Optimization terminated with status %D.\n", reason);
239:   }
240:   TaoDestroy(&tao);
241:   VecDestroy(&x);
242:   VecDestroy(&x0);
243:   ParabolicDestroy(&user);

245:   PetscFinalize();
246:   return 0;
247: }
248: /* ------------------------------------------------------------------- */
251: /*
252:    dwork = Qy - d
253:    lwork = L*(u-ur)
254:    f = 1/2 * (dwork.dork + alpha*lwork.lwork)
255: */
256: PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
257: {
259:   PetscReal      d1=0,d2=0;
260:   PetscInt       i,j;
261:   AppCtx         *user = (AppCtx*)ptr;

264:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
265:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
266:   for (j=0; j<user->ns; j++){
267:     i = user->sample_times[j];
268:     MatMult(user->Qblock,user->yi[i],user->di[j]);
269:   }
270:   Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
271:   VecAXPY(user->dwork,-1.0,user->d);
272:   VecDot(user->dwork,user->dwork,&d1);

274:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
275:   MatMult(user->L,user->uwork,user->lwork);
276:   VecDot(user->lwork,user->lwork,&d2);

278:   *f = 0.5 * (d1 + user->alpha*d2);
279:   return(0);
280: }

282: /* ------------------------------------------------------------------- */
285: /*
286:     state: g_s = Q' *(Qy - d)
287:     design: g_d = alpha*L'*L*(u-ur)
288: */
289: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
290: {
292:   PetscInt       i,j;
293:   AppCtx         *user = (AppCtx*)ptr;

296:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
297:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
298:   for (j=0; j<user->ns; j++){
299:     i = user->sample_times[j];
300:     MatMult(user->Qblock,user->yi[i],user->di[j]);
301:   }
302:   Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
303:   VecAXPY(user->dwork,-1.0,user->d);
304:   Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
305:   VecSet(user->ywork,0.0);
306:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
307:   for (j=0; j<user->ns; j++){
308:     i = user->sample_times[j];
309:     MatMult(user->QblockT,user->di[j],user->yiwork[i]);
310:   }
311:   Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);

313:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
314:   MatMult(user->L,user->uwork,user->lwork);
315:   MatMult(user->LT,user->lwork,user->uwork);
316:   VecScale(user->uwork, user->alpha);
317:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
318:   return(0);
319: }

323: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
324: {
326:   PetscReal      d1,d2;
327:   PetscInt       i,j;
328:   AppCtx         *user = (AppCtx*)ptr;

331:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
332:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
333:   for (j=0; j<user->ns; j++){
334:     i = user->sample_times[j];
335:     MatMult(user->Qblock,user->yi[i],user->di[j]);
336:   }
337:   Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
338:   VecAXPY(user->dwork,-1.0,user->d);
339:   VecDot(user->dwork,user->dwork,&d1);
340:   Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
341:   VecSet(user->ywork,0.0);
342:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
343:   for (j=0; j<user->ns; j++){
344:     i = user->sample_times[j];
345:     MatMult(user->QblockT,user->di[j],user->yiwork[i]);
346:   }
347:   Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);

349:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
350:   MatMult(user->L,user->uwork,user->lwork);
351:   VecDot(user->lwork,user->lwork,&d2);
352:   MatMult(user->LT,user->lwork,user->uwork);
353:   VecScale(user->uwork, user->alpha);
354:   *f = 0.5 * (d1 + user->alpha*d2);

356:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
357:   return(0);
358: }

360: /* ------------------------------------------------------------------- */
363: /* A
364: MatShell object
365: */
366: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
367: {
369:   AppCtx         *user = (AppCtx*)ptr;

372:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
373:   VecSet(user->uwork,0);
374:   VecAXPY(user->uwork,-1.0,user->u);
375:   VecExp(user->uwork);
376:   MatMult(user->Av,user->uwork,user->Av_u);
377:   VecCopy(user->Av_u,user->Swork);
378:   VecReciprocal(user->Swork);
379:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
380:   MatDiagonalScale(user->Divwork,NULL,user->Swork);
381:   if (user->dsg_formed) {
382:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
383:   } else {
384:     MatMatMultSymbolic(user->Divwork,user->Grad,PETSC_DECIDE,&user->DSG);
385:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
386:     user->dsg_formed = PETSC_TRUE;
387:   }

389:   /* B = speye(nx^3) + ht*DSG; */
390:   MatScale(user->DSG,user->ht);
391:   MatShift(user->DSG,1.0);
392:   return(0);
393: }

395: /* ------------------------------------------------------------------- */
398: /* B */
399: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
400: {
402:   AppCtx         *user = (AppCtx*)ptr;

405:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
406:   return(0);
407: }

411: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
412: {
414:   PetscInt       i;
415:   AppCtx         *user;

418:   MatShellGetContext(J_shell,(void**)&user);
419:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
420:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
421:   for (i=1; i<user->nt; i++){
422:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
423:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
424:   }
425:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
426:   return(0);
427: }

431: PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
432: {
434:   PetscInt       i;
435:   AppCtx         *user;

438:   MatShellGetContext(J_shell,(void**)&user);
439:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
440:   for (i=0; i<user->nt-1; i++){
441:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
442:     VecAXPY(user->yiwork[i],-1.0,user->yi[i+1]);
443:   }
444:   i = user->nt-1;
445:   MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
446:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
447:   return(0);
448: }

452: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
453: {
455:   AppCtx         *user;

458:   MatShellGetContext(J_shell,(void**)&user);
459:   MatMult(user->Grad,X,user->Swork);
460:   VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
461:   MatMult(user->Div,user->Swork,Y);
462:   VecAYPX(Y,user->ht,X);
463:   return(0);
464: }

468: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
469: {
471:   PetscInt       i;
472:   AppCtx         *user;

475:   MatShellGetContext(J_shell,(void**)&user);

477:   /* sdiag(1./v) */
478:   VecSet(user->uwork,0);
479:   VecAXPY(user->uwork,-1.0,user->u);
480:   VecExp(user->uwork);

482:   /* sdiag(1./((Av*(1./v)).^2)) */
483:   MatMult(user->Av,user->uwork,user->Swork);
484:   VecPointwiseMult(user->Swork,user->Swork,user->Swork);
485:   VecReciprocal(user->Swork);

487:   /* (Av * (sdiag(1./v) * b)) */
488:   VecPointwiseMult(user->uwork,user->uwork,X);
489:   MatMult(user->Av,user->uwork,user->Twork);

491:   /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
492:   VecPointwiseMult(user->Swork,user->Twork,user->Swork);

494:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
495:   for (i=0; i<user->nt; i++){
496:     /* (sdiag(Grad*y(:,i)) */
497:     MatMult(user->Grad,user->yi[i],user->Twork);

499:     /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
500:     VecPointwiseMult(user->Twork,user->Twork,user->Swork);
501:     MatMult(user->Div,user->Twork,user->yiwork[i]);
502:     VecScale(user->yiwork[i],user->ht);
503:   }
504:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);

506:   return(0);
507: }

511: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
512: {
514:   PetscInt       i;
515:   AppCtx         *user;

518:   MatShellGetContext(J_shell,(void**)&user);

520:   /* sdiag(1./((Av*(1./v)).^2)) */
521:   VecSet(user->uwork,0);
522:   VecAXPY(user->uwork,-1.0,user->u);
523:   VecExp(user->uwork);
524:   MatMult(user->Av,user->uwork,user->Rwork);
525:   VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork);
526:   VecReciprocal(user->Rwork);

528:   VecSet(Y,0.0);
529:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
530:   Scatter_i(X,user->yiwork,user->yi_scatter,user->nt);
531:   for (i=0; i<user->nt; i++){
532:     /* (Div' * b(:,i)) */
533:     MatMult(user->Grad,user->yiwork[i],user->Swork);

535:     /* sdiag(Grad*y(:,i)) */
536:     MatMult(user->Grad,user->yi[i],user->Twork);

538:     /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */
539:     VecPointwiseMult(user->Twork,user->Swork,user->Twork);

541:     /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */
542:     VecPointwiseMult(user->Twork,user->Rwork,user->Twork);

544:     /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
545:     MatMult(user->AvT,user->Twork,user->yiwork[i]);

547:     /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
548:     VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i]);
549:     VecAXPY(Y,user->ht,user->yiwork[i]);
550:   }
551:   return(0);
552: }

556: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
557: {
559:   AppCtx         *user;

562:   PCShellGetContext(PC_shell,(void**)&user);

564:   if (user->dsg_formed) {
565:     MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
566:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed");
567:   return(0);
568: }

572: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
573: {
575:   AppCtx         *user;
576:   PetscInt       its,i;

579:   MatShellGetContext(J_shell,(void**)&user);

581:   if (Y == user->ytrue) {
582:     /* First solve is done with true solution to set up problem */
583:     KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
584:   } else {
585:     KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
586:   }

588:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
589:   KSPSolve(user->solver,user->yi[0],user->yiwork[0]);
590:   KSPGetIterationNumber(user->solver,&its);
591:   user->ksp_its = user->ksp_its + its;

593:   for (i=1; i<user->nt; i++){
594:     VecAXPY(user->yi[i],1.0,user->yiwork[i-1]);
595:     KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
596:     KSPGetIterationNumber(user->solver,&its);
597:     user->ksp_its = user->ksp_its + its;
598:   }
599:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
600:   return(0);
601: }

605: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
606: {
608:   AppCtx         *user;
609:   PetscInt       its,i;

612:   MatShellGetContext(J_shell,(void**)&user);

614:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);

616:   i = user->nt - 1;
617:   KSPSolve(user->solver,user->yi[i],user->yiwork[i]);

619:   KSPGetIterationNumber(user->solver,&its);
620:   user->ksp_its = user->ksp_its + its;

622:   for (i=user->nt-2; i>=0; i--){
623:     VecAXPY(user->yi[i],1.0,user->yiwork[i+1]);
624:     KSPSolve(user->solver,user->yi[i],user->yiwork[i]);

626:     KSPGetIterationNumber(user->solver,&its);
627:     user->ksp_its = user->ksp_its + its;

629:   }

631:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
632:   return(0);
633: }

637: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
638: {
640:   AppCtx         *user;

643:   MatShellGetContext(J_shell,(void**)&user);

645:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
646:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
647:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
648:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
649:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
650:   return(0);
651: }

655: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
656: {
658:   AppCtx         *user;

661:   MatShellGetContext(J_shell,(void**)&user);
662:    VecCopy(user->js_diag,X);
663:   return(0);

665: }

669: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
670: {
671:   /* con = Ay - q, A = [B  0  0 ... 0;
672:                        -I  B  0 ... 0;
673:                         0 -I  B ... 0;
674:                              ...     ;
675:                         0    ... -I B]
676:      B = ht * Div * Sigma * Grad + eye */
678:   PetscInt       i;
679:   AppCtx         *user = (AppCtx*)ptr;

682:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
683:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
684:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
685:   for (i=1; i<user->nt; i++){
686:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
687:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
688:   }
689:   Gather_i(C,user->yiwork,user->yi_scatter,user->nt);
690:   VecAXPY(C,-1.0,user->q);
691:   return(0);
692: }


697: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
698: {

702:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
703:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
704:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
705:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
706:   return(0);
707: }

711: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
712: {
714:   PetscInt       i;

717:   for (i=0; i<nt; i++){
718:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
719:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
720:   }
721:   return(0);
722: }


727: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
728: {

732:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
733:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
734:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
735:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
736:   return(0);
737: }

741: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
742: {
744:   PetscInt       i;

747:   for (i=0; i<nt; i++){
748:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
749:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
750:   }
751:   return(0);
752: }

756: PetscErrorCode ParabolicInitialize(AppCtx *user)
757: {
759:   PetscInt       m,n,i,j,k,linear_index,istart,iend,iblock,lo,hi,lo2,hi2;
760:   Vec            XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork,yi,di,bc;
761:   PetscReal      *x, *y, *z;
762:   PetscReal      h,stime;
763:   PetscScalar    hinv,neg_hinv,half = 0.5,sqrt_beta;
764:   PetscInt       im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
765:   PetscReal      xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
766:   PetscScalar    v,vx,vy,vz;
767:   IS             is_from_y,is_to_yi,is_from_d,is_to_di;
768:   /* Data locations */
769:   PetscScalar xr[64] = {0.4970,     0.8498,     0.7814,     0.6268,     0.7782,     0.6402,     0.3617,     0.3160,
770:                         0.3610,     0.5298,     0.6987,     0.3331,     0.7962,     0.5596,     0.3866,     0.6774,
771:                         0.5407,     0.4518,     0.6702,     0.6061,     0.7580,     0.8997,     0.5198,     0.8326,
772:                         0.2138,     0.9198,     0.3000,     0.2833,     0.8288,     0.7076,     0.1820,     0.0728,
773:                         0.8447,     0.2367,     0.3239,     0.6413,     0.3114,     0.4731,     0.1192,     0.9273,
774:                         0.5724,     0.4331,     0.5136,     0.3547,     0.4413,     0.2602,     0.5698,     0.7278,
775:                         0.5261,     0.6230,     0.2454,     0.3948,     0.7479,     0.6582,     0.4660,     0.5594,
776:                         0.7574,     0.1143,     0.5900,     0.1065,     0.4260,     0.3294,     0.8276,     0.0756};

778:   PetscScalar yr[64] = {0.7345,     0.9120,     0.9288,     0.7528,     0.4463,     0.4985,     0.2497,     0.6256,
779:                         0.3425,     0.9026,     0.6983,     0.4230,     0.7140,     0.2970,     0.4474,     0.8792,
780:                         0.6604,     0.2485,     0.7968,     0.6127,     0.1796,     0.2437,     0.5938,     0.6137,
781:                         0.3867,     0.5658,     0.4575,     0.1009,     0.0863,     0.3361,     0.0738,     0.3985,
782:                         0.6602,     0.1437,     0.0934,     0.5983,     0.5950,     0.0763,     0.0768,     0.2288,
783:                         0.5761,     0.1129,     0.3841,     0.6150,     0.6904,     0.6686,     0.1361,     0.4601,
784:                         0.4491,     0.3716,     0.1969,     0.6537,     0.6743,     0.6991,     0.4811,     0.5480,
785:                         0.1684,     0.4569,     0.6889,     0.8437,     0.3015,     0.2854,     0.8199,     0.2658};

787:   PetscScalar zr[64] = {0.7668,     0.8573,     0.2654,     0.2719,     0.1060,     0.1311,     0.6232,     0.2295,
788:                         0.8009,     0.2147,     0.2119,     0.9325,     0.4473,     0.3600,     0.3374,     0.3819,
789:                         0.4066,     0.5801,     0.1673,     0.0959,     0.4638,     0.8236,     0.8800,     0.2939,
790:                         0.2028,     0.8262,     0.2706,     0.6276,     0.9085,     0.6443,     0.8241,     0.0712,
791:                         0.1824,     0.7789,     0.4389,     0.8415,     0.7055,     0.6639,     0.3653,     0.2078,
792:                         0.1987,     0.2297,     0.4321,     0.8115,     0.4915,     0.7764,     0.4657,     0.4627,
793:                         0.4569,     0.4232,     0.8514,     0.0674,     0.3227,     0.1055,     0.6690,     0.6313,
794:                         0.9226,     0.5461,     0.4126,     0.2364,     0.6096,     0.7042,     0.3914,     0.0711};

797:   PetscMalloc1(user->mx,&x);
798:   PetscMalloc1(user->mx,&y);
799:   PetscMalloc1(user->mx,&z);
800:   user->jformed = PETSC_FALSE;
801:   user->dsg_formed = PETSC_FALSE;

803:   n = user->mx * user->mx * user->mx;
804:   m = 3 * user->mx * user->mx * (user->mx-1);
805:   sqrt_beta = PetscSqrtScalar(user->beta);

807:   user->ksp_its = 0;
808:   user->ksp_its_initial = 0;

810:   stime = (PetscReal)user->nt/user->ns;
811:   PetscMalloc1(user->ns,&user->sample_times);
812:   for (i=0; i<user->ns; i++){
813:     user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
814:   }

816:   VecCreate(PETSC_COMM_WORLD,&XX);
817:   VecCreate(PETSC_COMM_WORLD,&user->q);
818:   VecSetSizes(XX,PETSC_DECIDE,n);
819:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
820:   VecSetFromOptions(XX);
821:   VecSetFromOptions(user->q);

823:   VecDuplicate(XX,&YY);
824:   VecDuplicate(XX,&ZZ);
825:   VecDuplicate(XX,&XXwork);
826:   VecDuplicate(XX,&YYwork);
827:   VecDuplicate(XX,&ZZwork);
828:   VecDuplicate(XX,&UTwork);
829:   VecDuplicate(XX,&user->utrue);
830:   VecDuplicate(XX,&bc);

832:   /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
833:   h = 1.0/user->mx;
834:   hinv = user->mx;
835:   neg_hinv = -hinv;

837:   VecGetOwnershipRange(XX,&istart,&iend);
838:   for (linear_index=istart; linear_index<iend; linear_index++){
839:     i = linear_index % user->mx;
840:     j = ((linear_index-i)/user->mx) % user->mx;
841:     k = ((linear_index-i)/user->mx-j) / user->mx;
842:     vx = h*(i+0.5);
843:     vy = h*(j+0.5);
844:     vz = h*(k+0.5);
845:     VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
846:     VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
847:     VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES);
848:     if ((vx<0.6) && (vx>0.4) && (vy<0.6) && (vy>0.4) && (vy<0.6) && (vz<0.6) && (vz>0.4)){
849:       v = 1000.0;
850:       VecSetValues(bc,1,&linear_index,&v,INSERT_VALUES);
851:     }
852:   }

854:   VecAssemblyBegin(XX);
855:   VecAssemblyEnd(XX);
856:   VecAssemblyBegin(YY);
857:   VecAssemblyEnd(YY);
858:   VecAssemblyBegin(ZZ);
859:   VecAssemblyEnd(ZZ);
860:   VecAssemblyBegin(bc);
861:   VecAssemblyEnd(bc);

863:   /* Compute true parameter function
864:      ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
865:   VecCopy(XX,XXwork);
866:   VecCopy(YY,YYwork);
867:   VecCopy(ZZ,ZZwork);

869:   VecShift(XXwork,-0.5);
870:   VecShift(YYwork,-0.5);
871:   VecShift(ZZwork,-0.5);

873:   VecPointwiseMult(XXwork,XXwork,XXwork);
874:   VecPointwiseMult(YYwork,YYwork,YYwork);
875:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

877:   VecCopy(XXwork,user->utrue);
878:   VecAXPY(user->utrue,1.0,YYwork);
879:   VecAXPY(user->utrue,1.0,ZZwork);
880:   VecScale(user->utrue,-10.0);
881:   VecExp(user->utrue);
882:   VecShift(user->utrue,0.5);

884:   VecDestroy(&XX);
885:   VecDestroy(&YY);
886:   VecDestroy(&ZZ);
887:   VecDestroy(&XXwork);
888:   VecDestroy(&YYwork);
889:   VecDestroy(&ZZwork);
890:   VecDestroy(&UTwork);

892:   /* Initial guess and reference model */
893:   VecDuplicate(user->utrue,&user->ur);
894:   VecSet(user->ur,0.0);

896:   /* Generate Grad matrix */
897:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
898:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n);
899:   MatSetFromOptions(user->Grad);
900:   MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
901:   MatSeqAIJSetPreallocation(user->Grad,2,NULL);
902:   MatGetOwnershipRange(user->Grad,&istart,&iend);

904:   for (i=istart; i<iend; i++){
905:     if (i<m/3){
906:       iblock = i / (user->mx-1);
907:       j = iblock*user->mx + (i % (user->mx-1));
908:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
909:       j = j+1;
910:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
911:     }
912:     if (i>=m/3 && i<2*m/3){
913:       iblock = (i-m/3) / (user->mx*(user->mx-1));
914:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
915:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
916:       j = j + user->mx;
917:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
918:     }
919:     if (i>=2*m/3){
920:       j = i-2*m/3;
921:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
922:       j = j + user->mx*user->mx;
923:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
924:     }
925:   }

927:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
928:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);


931:   /* Generate arithmetic averaging matrix Av */
932:   MatCreate(PETSC_COMM_WORLD,&user->Av);
933:   MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n);
934:   MatSetFromOptions(user->Av);
935:   MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
936:   MatSeqAIJSetPreallocation(user->Av,2,NULL);
937:   MatGetOwnershipRange(user->Av,&istart,&iend);

939:   for (i=istart; i<iend; i++){
940:     if (i<m/3){
941:       iblock = i / (user->mx-1);
942:       j = iblock*user->mx + (i % (user->mx-1));
943:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
944:       j = j+1;
945:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
946:     }
947:     if (i>=m/3 && i<2*m/3){
948:       iblock = (i-m/3) / (user->mx*(user->mx-1));
949:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
950:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
951:       j = j + user->mx;
952:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
953:     }
954:     if (i>=2*m/3){
955:       j = i-2*m/3;
956:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
957:       j = j + user->mx*user->mx;
958:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
959:     }
960:   }

962:   MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
963:   MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);

965:   /* Generate transpose of averaging matrix Av */
966:   MatTranspose(user->Av,MAT_INITIAL_MATRIX,&user->AvT);

968:   MatCreate(PETSC_COMM_WORLD,&user->L);
969:   MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n);
970:   MatSetFromOptions(user->L);
971:   MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
972:   MatSeqAIJSetPreallocation(user->L,2,NULL);
973:   MatGetOwnershipRange(user->L,&istart,&iend);

975:   for (i=istart; i<iend; i++){
976:     if (i<m/3){
977:       iblock = i / (user->mx-1);
978:       j = iblock*user->mx + (i % (user->mx-1));
979:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
980:       j = j+1;
981:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
982:     }
983:     if (i>=m/3 && i<2*m/3){
984:       iblock = (i-m/3) / (user->mx*(user->mx-1));
985:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
986:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
987:       j = j + user->mx;
988:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
989:     }
990:     if (i>=2*m/3 && i<m){
991:       j = i-2*m/3;
992:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
993:       j = j + user->mx*user->mx;
994:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
995:     }
996:     if (i>=m){
997:       j = i - m;
998:       MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES);
999:     }
1000:   }

1002:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
1003:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);

1005:   MatScale(user->L,PetscPowScalar(h,1.5));

1007:   /* Generate Div matrix */
1008:   MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div);

1010:   /* Build work vectors and matrices */
1011:   VecCreate(PETSC_COMM_WORLD,&user->S);
1012:   VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
1013:   VecSetFromOptions(user->S);

1015:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1016:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
1017:   VecSetFromOptions(user->lwork);

1019:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
1020:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

1022:   VecCreate(PETSC_COMM_WORLD,&user->d);
1023:   VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
1024:   VecSetFromOptions(user->d);

1026:   VecDuplicate(user->S,&user->Swork);
1027:   VecDuplicate(user->S,&user->Av_u);
1028:   VecDuplicate(user->S,&user->Twork);
1029:   VecDuplicate(user->S,&user->Rwork);
1030:   VecDuplicate(user->y,&user->ywork);
1031:   VecDuplicate(user->u,&user->uwork);
1032:   VecDuplicate(user->u,&user->js_diag);
1033:   VecDuplicate(user->c,&user->cwork);
1034:   VecDuplicate(user->d,&user->dwork);

1036:   /* Create matrix-free shell user->Js for computing A*x */
1037:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
1038:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1039:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
1040:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
1041:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

1043:   /* Diagonal blocks of user->Js */
1044:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
1045:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
1046:   /* Blocks are symmetric */
1047:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);

1049:   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1050:      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1051:      This is an SSOR preconditioner for user->JsBlock. */
1052:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlockPrec);
1053:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult);
1054:   /* JsBlockPrec is symmetric */
1055:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMult);
1056:   MatSetOption(user->JsBlockPrec,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

1058:   /* Create a matrix-free shell user->Jd for computing B*x */
1059:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m,user,&user->Jd);
1060:   MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult);
1061:   MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose);

1063:   /* User-defined routines for computing user->Js\x and user->Js^T\x*/
1064:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->JsInv);
1065:   MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult);
1066:   MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult);

1068:   /* Solver options and tolerances */
1069:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1070:   KSPSetType(user->solver,KSPCG);
1071:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
1072:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
1073:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1074:   KSPGetPC(user->solver,&user->prec);
1075:   PCSetType(user->prec,PCSHELL);

1077:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1078:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
1079:   PCShellSetContext(user->prec,user);

1081:   /* Create scatter from y to y_1,y_2,...,y_nt */
1082:   PetscMalloc1(user->nt*user->m,&user->yi_scatter);
1083:   VecCreate(PETSC_COMM_WORLD,&yi);
1084:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
1085:   VecSetFromOptions(yi);
1086:   VecDuplicateVecs(yi,user->nt,&user->yi);
1087:   VecDuplicateVecs(yi,user->nt,&user->yiwork);

1089:   VecGetOwnershipRange(user->y,&lo2,&hi2);
1090:   istart = 0;
1091:   for (i=0; i<user->nt; i++){
1092:     VecGetOwnershipRange(user->yi[i],&lo,&hi);
1093:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
1094:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
1095:     VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
1096:     istart = istart + hi-lo;
1097:     ISDestroy(&is_to_yi);
1098:     ISDestroy(&is_from_y);
1099:   }
1100:   VecDestroy(&yi);

1102:   /* Create scatter from d to d_1,d_2,...,d_ns */
1103:   PetscMalloc1(user->ns*user->ndata,&user->di_scatter);
1104:   VecCreate(PETSC_COMM_WORLD,&di);
1105:   VecSetSizes(di,PETSC_DECIDE,user->ndata);
1106:   VecSetFromOptions(di);
1107:   VecDuplicateVecs(di,user->ns,&user->di);
1108:   VecGetOwnershipRange(user->d,&lo2,&hi2);
1109:   istart = 0;
1110:   for (i=0; i<user->ns; i++){
1111:     VecGetOwnershipRange(user->di[i],&lo,&hi);
1112:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_di);
1113:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d);
1114:     VecScatterCreate(user->d,is_from_d,user->di[i],is_to_di,&user->di_scatter[i]);
1115:     istart = istart + hi-lo;
1116:     ISDestroy(&is_to_di);
1117:     ISDestroy(&is_from_d);
1118:   }
1119:   VecDestroy(&di);

1121:   /* Assemble RHS of forward problem */
1122:   VecCopy(bc,user->yiwork[0]);
1123:   for (i=1; i<user->nt; i++){
1124:     VecSet(user->yiwork[i],0.0);
1125:   }
1126:   Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1127:   VecDestroy(&bc);

1129:   /* Compute true state function yt given ut */
1130:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1131:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1132:   VecSetFromOptions(user->ytrue);

1134:   /* First compute Av_u = Av*exp(-u) */
1135:   VecSet(user->uwork,0);
1136:   VecAXPY(user->uwork,-1.0,user->utrue); /*  Note: user->utrue */
1137:   VecExp(user->uwork);
1138:   MatMult(user->Av,user->uwork,user->Av_u);

1140:   MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1141:   user->dsg_formed = PETSC_TRUE;

1143:   /* Next form DSG = Div*S*Grad */
1144:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1145:   MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1146:   if (user->dsg_formed) {
1147:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1148:   } else {
1149:     MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1150:     MatMatMultNumeric(user->Div,user->Grad,user->DSG);
1151:     user->dsg_formed = PETSC_TRUE;
1152:   }
1153:   /* B = speye(nx^3) + ht*DSG; */
1154:   MatScale(user->DSG,user->ht);
1155:   MatShift(user->DSG,1.0);

1157:   /* Now solve for ytrue */
1158:   StateMatInvMult(user->Js,user->q,user->ytrue);

1160:   /* Initial guess y0 for state given u0 */

1162:   /* First compute Av_u = Av*exp(-u) */
1163:   VecSet(user->uwork,0);
1164:   VecAXPY(user->uwork,-1.0,user->u); /*  Note: user->u */
1165:   VecExp(user->uwork);
1166:   MatMult(user->Av,user->uwork,user->Av_u);

1168:   /* Next form DSG = Div*S*Grad */
1169:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1170:   MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1171:   if (user->dsg_formed) {
1172:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1173:   } else {
1174:     MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1175:     MatMatMultNumeric(user->Div,user->Grad,user->DSG);
1176:     user->dsg_formed = PETSC_TRUE;
1177:   }
1178:   /* B = speye(nx^3) + ht*DSG; */
1179:   MatScale(user->DSG,user->ht);
1180:   MatShift(user->DSG,1.0);

1182:   /* Now solve for y */
1183:   StateMatInvMult(user->Js,user->q,user->y);

1185:   /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
1186:   MatCreate(PETSC_COMM_WORLD,&user->Qblock);
1187:   MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n);
1188:   MatSetFromOptions(user->Qblock);
1189:   MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL);
1190:   MatSeqAIJSetPreallocation(user->Qblock,8,NULL);

1192:   for (i=0; i<user->mx; i++){
1193:     x[i] = h*(i+0.5);
1194:     y[i] = h*(i+0.5);
1195:     z[i] = h*(i+0.5);
1196:   }

1198:   MatGetOwnershipRange(user->Qblock,&istart,&iend);
1199:   nx = user->mx; ny = user->mx; nz = user->mx;
1200:   for (i=istart; i<iend; i++){
1201:     xri = xr[i];
1202:     im = 0;
1203:     xim = x[im];
1204:     while (xri>xim && im<nx){
1205:       im = im+1;
1206:       xim = x[im];
1207:     }
1208:     indx1 = im-1;
1209:     indx2 = im;
1210:     dx1 = xri - x[indx1];
1211:     dx2 = x[indx2] - xri;

1213:     yri = yr[i];
1214:     im = 0;
1215:     yim = y[im];
1216:     while (yri>yim && im<ny){
1217:       im = im+1;
1218:       yim = y[im];
1219:     }
1220:     indy1 = im-1;
1221:     indy2 = im;
1222:     dy1 = yri - y[indy1];
1223:     dy2 = y[indy2] - yri;

1225:     zri = zr[i];
1226:     im = 0;
1227:     zim = z[im];
1228:     while (zri>zim && im<nz){
1229:       im = im+1;
1230:       zim = z[im];
1231:     }
1232:     indz1 = im-1;
1233:     indz2 = im;
1234:     dz1 = zri - z[indz1];
1235:     dz2 = z[indz2] - zri;

1237:     Dx = x[indx2] - x[indx1];
1238:     Dy = y[indy2] - y[indy1];
1239:     Dz = z[indz2] - z[indz1];

1241:     j = indx1 + indy1*nx + indz1*nx*ny;
1242:     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1243:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);

1245:     j = indx1 + indy1*nx + indz2*nx*ny;
1246:     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1247:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);

1249:     j = indx1 + indy2*nx + indz1*nx*ny;
1250:     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1251:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);

1253:     j = indx1 + indy2*nx + indz2*nx*ny;
1254:     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1255:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);

1257:     j = indx2 + indy1*nx + indz1*nx*ny;
1258:     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1259:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);

1261:     j = indx2 + indy1*nx + indz2*nx*ny;
1262:     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1263:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);

1265:     j = indx2 + indy2*nx + indz1*nx*ny;
1266:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1267:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);

1269:     j = indx2 + indy2*nx + indz2*nx*ny;
1270:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1271:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1272:   }
1273:   MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY);
1274:   MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY);

1276:   MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1277:   MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);

1279:   /* Add noise to the measurement data */
1280:   VecSet(user->ywork,1.0);
1281:   VecAYPX(user->ywork,user->noise,user->ytrue);
1282:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
1283:   for (j=0; j<user->ns; j++){
1284:     i = user->sample_times[j];
1285:     MatMult(user->Qblock,user->yiwork[i],user->di[j]);
1286:   }
1287:   Gather_i(user->d,user->di,user->di_scatter,user->ns);

1289:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1290:   KSPSetFromOptions(user->solver);
1291:   PetscFree(x);
1292:   PetscFree(y);
1293:   PetscFree(z);
1294:   return(0);
1295: }

1299: PetscErrorCode ParabolicDestroy(AppCtx *user)
1300: {
1302:   PetscInt       i;

1305:   MatDestroy(&user->Qblock);
1306:   MatDestroy(&user->QblockT);
1307:   MatDestroy(&user->Div);
1308:   MatDestroy(&user->Divwork);
1309:   MatDestroy(&user->Grad);
1310:   MatDestroy(&user->Av);
1311:   MatDestroy(&user->Avwork);
1312:   MatDestroy(&user->AvT);
1313:   MatDestroy(&user->DSG);
1314:   MatDestroy(&user->L);
1315:   MatDestroy(&user->LT);
1316:   KSPDestroy(&user->solver);
1317:   MatDestroy(&user->Js);
1318:   MatDestroy(&user->Jd);
1319:   MatDestroy(&user->JsInv);
1320:   MatDestroy(&user->JsBlock);
1321:   MatDestroy(&user->JsBlockPrec);
1322:   VecDestroy(&user->u);
1323:   VecDestroy(&user->uwork);
1324:   VecDestroy(&user->utrue);
1325:   VecDestroy(&user->y);
1326:   VecDestroy(&user->ywork);
1327:   VecDestroy(&user->ytrue);
1328:   VecDestroyVecs(user->nt,&user->yi);
1329:   VecDestroyVecs(user->nt,&user->yiwork);
1330:   VecDestroyVecs(user->ns,&user->di);
1331:   PetscFree(user->yi);
1332:   PetscFree(user->yiwork);
1333:   PetscFree(user->di);
1334:   VecDestroy(&user->c);
1335:   VecDestroy(&user->cwork);
1336:   VecDestroy(&user->ur);
1337:   VecDestroy(&user->q);
1338:   VecDestroy(&user->d);
1339:   VecDestroy(&user->dwork);
1340:   VecDestroy(&user->lwork);
1341:   VecDestroy(&user->S);
1342:   VecDestroy(&user->Swork);
1343:   VecDestroy(&user->Av_u);
1344:   VecDestroy(&user->Twork);
1345:   VecDestroy(&user->Rwork);
1346:   VecDestroy(&user->js_diag);
1347:   ISDestroy(&user->s_is);
1348:   ISDestroy(&user->d_is);
1349:   VecScatterDestroy(&user->state_scatter);
1350:   VecScatterDestroy(&user->design_scatter);
1351:   for (i=0; i<user->nt; i++){
1352:     VecScatterDestroy(&user->yi_scatter[i]);
1353:   }
1354:   for (i=0; i<user->ns; i++){
1355:     VecScatterDestroy(&user->di_scatter[i]);
1356:   }
1357:   PetscFree(user->yi_scatter);
1358:   PetscFree(user->di_scatter);
1359:   PetscFree(user->sample_times);
1360:   return(0);
1361: }

1365: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1366: {
1368:   Vec            X;
1369:   PetscReal      unorm,ynorm;
1370:   AppCtx         *user = (AppCtx*)ptr;

1373:   TaoGetSolutionVector(tao,&X);
1374:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1375:   VecAXPY(user->ywork,-1.0,user->ytrue);
1376:   VecAXPY(user->uwork,-1.0,user->utrue);
1377:   VecNorm(user->uwork,NORM_2,&unorm);
1378:   VecNorm(user->ywork,NORM_2,&ynorm);
1379:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1380:   return(0);
1381: }