Actual source code: parabolic.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <petsc-private/taoimpl.h>
  2: #include "../src/tao/pde_constrained/impls/lcl/lcl.h"

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

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

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

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

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

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

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

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

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

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

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

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

 78:   KSP solver;
 79:   PC prec;

 81:   TAO_LCL *lcl;
 82:   PetscInt ksp_its;
 83:   PetscInt ksp_its_initial;
 84: } AppCtx;


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

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

109: PetscErrorCode DesignMatMult(Mat,Vec,Vec);
110: PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);

112: PetscErrorCode Gather_i(Vec,Vec*,VecScatter*,PetscInt);
113: PetscErrorCode Scatter_i(Vec,Vec*,VecScatter*,PetscInt);

115: static  char help[]="";

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

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

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

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

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

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

174:   VecGetOwnershipRange(user.y,&lo,&hi);
175:   VecGetOwnershipRange(user.u,&lo2,&hi2);

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

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

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

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

191:   /* Create TAO solver and set desired solution method */
192:   TaoCreate(PETSC_COMM_WORLD,&tao);
193:   TaoSetType(tao,TAOLCL);
194:   user.lcl = (TAO_LCL*)(tao->data);
195:   /* Set up initial vectors and matrices */
196:   ParabolicInitialize(&user);

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

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

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

211:   TaoSetFromOptions(tao);
212:   TaoSetStateDesignIS(tao,user.s_is,user.d_is);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

507:   return(0);
508: }

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

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

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

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

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

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

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

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

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

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

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

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

573: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
574: {
576:   AppCtx         *user;
577:   PetscReal      tau;
578:   PetscInt       its,i;

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

583:   if (Y == user->ytrue) {
584:     KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
585:   } else if (user->lcl) {
586:     tau = user->lcl->tau[user->lcl->solve_type];
587:     KSPSetTolerances(user->solver,tau,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
588:   }

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

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

607: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
608: {
610:   PetscReal      tau;
611:   AppCtx         *user;
612:   PetscInt       its,i;

615:   MatShellGetContext(J_shell,(void**)&user);
616:   if (user->lcl) {
617:     tau = user->lcl->tau[user->lcl->solve_type];
618:     KSPSetTolerances(user->solver,tau,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
619:   }
620:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);

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

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

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

632:     KSPGetIterationNumber(user->solver,&its);
633:     user->ksp_its = user->ksp_its + its;

635:   }

637:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
638:   return(0);
639: }

643: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
644: {
646:   AppCtx         *user;

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

651:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
652:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
653:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
654:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
655:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
656:   return(0);
657: }

661: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
662: {
664:   AppCtx         *user;

667:   MatShellGetContext(J_shell,(void**)&user);
668:    VecCopy(user->js_diag,X);
669:   return(0);

671: }

675: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
676: {
677:   /* con = Ay - q, A = [B  0  0 ... 0;
678:                        -I  B  0 ... 0;
679:                         0 -I  B ... 0;
680:                              ...     ;
681:                         0    ... -I B]
682:      B = ht * Div * Sigma * Grad + eye */
684:   PetscInt       i;
685:   AppCtx         *user = (AppCtx*)ptr;

688:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
689:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
690:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
691:   for (i=1; i<user->nt; i++){
692:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
693:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
694:   }
695:   Gather_i(C,user->yiwork,user->yi_scatter,user->nt);
696:   VecAXPY(C,-1.0,user->q);
697:   return(0);
698: }


703: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
704: {

708:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
709:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
710:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
711:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
712:   return(0);
713: }

717: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
718: {
720:   PetscInt       i;

723:   for (i=0; i<nt; i++){
724:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
725:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
726:   }
727:   return(0);
728: }


733: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
734: {

738:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
739:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
740:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
741:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
742:   return(0);
743: }

747: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
748: {
750:   PetscInt       i;

753:   for (i=0; i<nt; i++){
754:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
755:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
756:   }
757:   return(0);
758: }

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

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

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

803:   PetscMalloc1(user->mx,&x);
804:   PetscMalloc1(user->mx,&y);
805:   PetscMalloc1(user->mx,&z);
806:   user->jformed = PETSC_FALSE;
807:   user->dsg_formed = PETSC_FALSE;

809:   n = user->mx * user->mx * user->mx;
810:   m = 3 * user->mx * user->mx * (user->mx-1);
811:   sqrt_beta = PetscSqrtScalar(user->beta);

813:   user->ksp_its = 0;
814:   user->ksp_its_initial = 0;

816:   stime = (PetscReal)user->nt/user->ns;
817:   PetscMalloc1(user->ns,&user->sample_times);
818:   for (i=0; i<user->ns; i++){
819:     user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
820:   }

822:   VecCreate(PETSC_COMM_WORLD,&XX);
823:   VecCreate(PETSC_COMM_WORLD,&user->q);
824:   VecSetSizes(XX,PETSC_DECIDE,n);
825:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
826:   VecSetFromOptions(XX);
827:   VecSetFromOptions(user->q);

829:   VecDuplicate(XX,&YY);
830:   VecDuplicate(XX,&ZZ);
831:   VecDuplicate(XX,&XXwork);
832:   VecDuplicate(XX,&YYwork);
833:   VecDuplicate(XX,&ZZwork);
834:   VecDuplicate(XX,&UTwork);
835:   VecDuplicate(XX,&user->utrue);
836:   VecDuplicate(XX,&bc);

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

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

860:   VecAssemblyBegin(XX);
861:   VecAssemblyEnd(XX);
862:   VecAssemblyBegin(YY);
863:   VecAssemblyEnd(YY);
864:   VecAssemblyBegin(ZZ);
865:   VecAssemblyEnd(ZZ);
866:   VecAssemblyBegin(bc);
867:   VecAssemblyEnd(bc);

869:   /* Compute true parameter function
870:      ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
871:   VecCopy(XX,XXwork);
872:   VecCopy(YY,YYwork);
873:   VecCopy(ZZ,ZZwork);

875:   VecShift(XXwork,-0.5);
876:   VecShift(YYwork,-0.5);
877:   VecShift(ZZwork,-0.5);

879:   VecPointwiseMult(XXwork,XXwork,XXwork);
880:   VecPointwiseMult(YYwork,YYwork,YYwork);
881:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

883:   VecCopy(XXwork,user->utrue);
884:   VecAXPY(user->utrue,1.0,YYwork);
885:   VecAXPY(user->utrue,1.0,ZZwork);
886:   VecScale(user->utrue,-10.0);
887:   VecExp(user->utrue);
888:   VecShift(user->utrue,0.5);

890:   VecDestroy(&XX);
891:   VecDestroy(&YY);
892:   VecDestroy(&ZZ);
893:   VecDestroy(&XXwork);
894:   VecDestroy(&YYwork);
895:   VecDestroy(&ZZwork);
896:   VecDestroy(&UTwork);

898:   /* Initial guess and reference model */
899:   VecDuplicate(user->utrue,&user->ur);
900:   VecSet(user->ur,0.0);

902:   /* Generate Grad matrix */
903:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
904:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n);
905:   MatSetFromOptions(user->Grad);
906:   MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
907:   MatSeqAIJSetPreallocation(user->Grad,2,NULL);
908:   MatGetOwnershipRange(user->Grad,&istart,&iend);

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

933:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
934:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);


937:   /* Generate arithmetic averaging matrix Av */
938:   MatCreate(PETSC_COMM_WORLD,&user->Av);
939:   MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n);
940:   MatSetFromOptions(user->Av);
941:   MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
942:   MatSeqAIJSetPreallocation(user->Av,2,NULL);
943:   MatGetOwnershipRange(user->Av,&istart,&iend);

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

968:   MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
969:   MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);

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

974:   MatCreate(PETSC_COMM_WORLD,&user->L);
975:   MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n);
976:   MatSetFromOptions(user->L);
977:   MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
978:   MatSeqAIJSetPreallocation(user->L,2,NULL);
979:   MatGetOwnershipRange(user->L,&istart,&iend);

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

1008:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
1009:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);

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

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

1016:   /* Build work vectors and matrices */
1017:   VecCreate(PETSC_COMM_WORLD,&user->S);
1018:   VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
1019:   VecSetFromOptions(user->S);

1021:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1022:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
1023:   VecSetFromOptions(user->lwork);

1025:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
1026:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

1028:   VecCreate(PETSC_COMM_WORLD,&user->d);
1029:   VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
1030:   VecSetFromOptions(user->d);

1032:   VecDuplicate(user->S,&user->Swork);
1033:   VecDuplicate(user->S,&user->Av_u);
1034:   VecDuplicate(user->S,&user->Twork);
1035:   VecDuplicate(user->S,&user->Rwork);
1036:   VecDuplicate(user->y,&user->ywork);
1037:   VecDuplicate(user->u,&user->uwork);
1038:   VecDuplicate(user->u,&user->js_diag);
1039:   VecDuplicate(user->c,&user->cwork);
1040:   VecDuplicate(user->d,&user->dwork);

1042:   /* Create matrix-free shell user->Js for computing A*x */
1043:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
1044:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1045:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
1046:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
1047:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

1049:   /* Diagonal blocks of user->Js */
1050:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
1051:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
1052:   /* Blocks are symmetric */
1053:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);

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

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

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

1074:   /* Solver options and tolerances */
1075:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1076:   KSPSetType(user->solver,KSPCG);
1077:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
1078:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
1079:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1080:   KSPGetPC(user->solver,&user->prec);
1081:   PCSetType(user->prec,PCSHELL);

1083:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1084:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
1085:   PCShellSetContext(user->prec,user);

1087:   /* Create scatter from y to y_1,y_2,...,y_nt */
1088:   PetscMalloc1(user->nt*user->m,&user->yi_scatter);
1089:   VecCreate(PETSC_COMM_WORLD,&yi);
1090:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
1091:   VecSetFromOptions(yi);
1092:   VecDuplicateVecs(yi,user->nt,&user->yi);
1093:   VecDuplicateVecs(yi,user->nt,&user->yiwork);

1095:   VecGetOwnershipRange(user->y,&lo2,&hi2);
1096:   istart = 0;
1097:   for (i=0; i<user->nt; i++){
1098:     VecGetOwnershipRange(user->yi[i],&lo,&hi);
1099:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
1100:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
1101:     VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
1102:     istart = istart + hi-lo;
1103:     ISDestroy(&is_to_yi);
1104:     ISDestroy(&is_from_y);
1105:   }
1106:   VecDestroy(&yi);

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

1127:   /* Assemble RHS of forward problem */
1128:   VecCopy(bc,user->yiwork[0]);
1129:   for (i=1; i<user->nt; i++){
1130:     VecSet(user->yiwork[i],0.0);
1131:   }
1132:   Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1133:   VecDestroy(&bc);

1135:   /* Compute true state function yt given ut */
1136:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1137:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1138:   VecSetFromOptions(user->ytrue);

1140:   /* First compute Av_u = Av*exp(-u) */
1141:   VecSet(user->uwork,0);
1142:   VecAXPY(user->uwork,-1.0,user->utrue); /*  Note: user->utrue */
1143:   VecExp(user->uwork);
1144:   MatMult(user->Av,user->uwork,user->Av_u);

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

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

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

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

1168:   /* First compute Av_u = Av*exp(-u) */
1169:   VecSet(user->uwork,0);
1170:   VecAXPY(user->uwork,-1.0,user->u); /*  Note: user->u */
1171:   VecExp(user->uwork);
1172:   MatMult(user->Av,user->uwork,user->Av_u);

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

1188:   /* Now solve for y */
1189:   user->lcl->solve_type = LCL_FORWARD1;
1190:   StateMatInvMult(user->Js,user->q,user->y);

1192:   /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
1193:   MatCreate(PETSC_COMM_WORLD,&user->Qblock);
1194:   MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n);
1195:   MatSetFromOptions(user->Qblock);
1196:   MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL);
1197:   MatSeqAIJSetPreallocation(user->Qblock,8,NULL);

1199:   for (i=0; i<user->mx; i++){
1200:     x[i] = h*(i+0.5);
1201:     y[i] = h*(i+0.5);
1202:     z[i] = h*(i+0.5);
1203:   }

1205:   MatGetOwnershipRange(user->Qblock,&istart,&iend);
1206:   nx = user->mx; ny = user->mx; nz = user->mx;
1207:   for (i=istart; i<iend; i++){
1208:     xri = xr[i];
1209:     im = 0;
1210:     xim = x[im];
1211:     while (xri>xim && im<nx){
1212:       im = im+1;
1213:       xim = x[im];
1214:     }
1215:     indx1 = im-1;
1216:     indx2 = im;
1217:     dx1 = xri - x[indx1];
1218:     dx2 = x[indx2] - xri;

1220:     yri = yr[i];
1221:     im = 0;
1222:     yim = y[im];
1223:     while (yri>yim && im<ny){
1224:       im = im+1;
1225:       yim = y[im];
1226:     }
1227:     indy1 = im-1;
1228:     indy2 = im;
1229:     dy1 = yri - y[indy1];
1230:     dy2 = y[indy2] - yri;

1232:     zri = zr[i];
1233:     im = 0;
1234:     zim = z[im];
1235:     while (zri>zim && im<nz){
1236:       im = im+1;
1237:       zim = z[im];
1238:     }
1239:     indz1 = im-1;
1240:     indz2 = im;
1241:     dz1 = zri - z[indz1];
1242:     dz2 = z[indz2] - zri;

1244:     Dx = x[indx2] - x[indx1];
1245:     Dy = y[indy2] - y[indy1];
1246:     Dz = z[indz2] - z[indz1];

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

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

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

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

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

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

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

1276:     j = indx2 + indy2*nx + indz2*nx*ny;
1277:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1278:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1279:   }
1280:   MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY);
1281:   MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY);

1283:   MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1284:   MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);

1286:   /* Add noise to the measurement data */
1287:   VecSet(user->ywork,1.0);
1288:   VecAYPX(user->ywork,user->noise,user->ytrue);
1289:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
1290:   for (j=0; j<user->ns; j++){
1291:     i = user->sample_times[j];
1292:     MatMult(user->Qblock,user->yiwork[i],user->di[j]);
1293:   }
1294:   Gather_i(user->d,user->di,user->di_scatter,user->ns);

1296:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1297:   KSPSetFromOptions(user->solver);
1298:   PetscFree(x);
1299:   PetscFree(y);
1300:   PetscFree(z);
1301:   return(0);
1302: }

1306: PetscErrorCode ParabolicDestroy(AppCtx *user)
1307: {
1309:   PetscInt       i;

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

1372: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1373: {
1375:   Vec            X;
1376:   PetscReal      unorm,ynorm;
1377:   AppCtx         *user = (AppCtx*)ptr;

1380:   TaoGetSolutionVector(tao,&X);
1381:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1382:   VecAXPY(user->ywork,-1.0,user->ytrue);
1383:   VecAXPY(user->uwork,-1.0,user->utrue);
1384:   VecNorm(user->uwork,NORM_2,&unorm);
1385:   VecNorm(user->ywork,NORM_2,&ynorm);
1386:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1387:   return(0);
1388: }