Actual source code: elliptic.c

petsc-master 2015-01-24
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 total variables */
 23:   PetscInt m; /* Number of constraints */
 24:   PetscInt nstate;
 25:   PetscInt ndesign;
 26:   PetscInt mx; /* grid points in each direction */
 27:   PetscInt ns; /* Number of data samples (1<=ns<=8)
 28:                   Currently only ns=1 is supported */
 29:   PetscInt ndata; /* Number of data points per sample */
 30:   IS       s_is;
 31:   IS       d_is;

 33:   VecScatter state_scatter;
 34:   VecScatter design_scatter;
 35:   VecScatter *yi_scatter, *di_scatter;
 36:   Vec        suby,subq,subd;
 37:   Mat        Js,Jd,JsPrec,JsInv,JsBlock;

 39:   PetscReal alpha; /* Regularization parameter */
 40:   PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */
 41:   PetscReal noise; /* Amount of noise to add to data */
 42:   PetscReal *ones;
 43:   Mat       Q;
 44:   Mat       MQ;
 45:   Mat       L;

 47:   Mat Grad;
 48:   Mat Av,Avwork;
 49:   Mat Div, Divwork;
 50:   Mat DSG;
 51:   Mat Diag,Ones;


 54:   Vec q;
 55:   Vec ur; /* reference */

 57:   Vec d;
 58:   Vec dwork;

 60:   Vec x; /* super vec of y,u */

 62:   Vec y; /* state variables */
 63:   Vec ywork;

 65:   Vec ytrue;

 67:   Vec u; /* design variables */
 68:   Vec uwork;

 70:   Vec utrue;

 72:   Vec js_diag;

 74:   Vec c; /* constraint vector */
 75:   Vec cwork;

 77:   Vec lwork;
 78:   Vec S;
 79:   Vec Swork,Twork,Sdiag,Ywork;
 80:   Vec Av_u;

 82:   KSP solver;
 83:   PC  prec;

 85:   PetscReal tola,tolb,tolc,told;
 86:   PetscInt  ksp_its;
 87:   PetscInt  ksp_its_initial;
 88:   int       stages[10];
 89:   PetscBool use_ptap;
 90:   PetscBool use_lrc;
 91: } AppCtx;

 93: PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*);
 94: PetscErrorCode FormGradient(Tao, Vec, Vec, void*);
 95: PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*);
 96: PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void*);
 97: PetscErrorCode FormJacobianDesign(Tao, Vec, Mat,void*);
 98: PetscErrorCode FormConstraints(Tao, Vec, Vec, void*);
 99: PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void*);
100: PetscErrorCode Gather(Vec, Vec, VecScatter, Vec, VecScatter);
101: PetscErrorCode Scatter(Vec, Vec, VecScatter, Vec, VecScatter);
102: PetscErrorCode EllipticInitialize(AppCtx*);
103: PetscErrorCode EllipticDestroy(AppCtx*);
104: PetscErrorCode EllipticMonitor(Tao, void*);

106: PetscErrorCode StateBlockMatMult(Mat,Vec,Vec);
107: PetscErrorCode StateMatMult(Mat,Vec,Vec);

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

113: PetscErrorCode QMatMult(Mat,Vec,Vec);
114: PetscErrorCode QMatMultTranspose(Mat,Vec,Vec);

116: static  char help[]="";

120: int main(int argc, char **argv)
121: {
122:   PetscErrorCode     ierr;
123:   Vec                x0;
124:   Tao                tao;
125:   TaoConvergedReason reason;
126:   AppCtx             user;
127:   PetscInt           ntests = 1;
128:   PetscInt           i;

130:   PetscInitialize(&argc, &argv, (char*)0,help);
131:   user.mx = 8;
132:   PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL);
133:   user.ns = 6;
134:   PetscOptionsInt("-ns","Number of data samples (1<=ns<=8)","",user.ns,&user.ns,NULL);
135:   user.ndata = 64;
136:   PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL);
137:   user.alpha = 0.1;
138:   PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL);
139:   user.beta = 0.00001;
140:   PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,NULL);
141:   user.noise = 0.01;
142:   PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,NULL);

144:   user.use_ptap = PETSC_FALSE;
145:   PetscOptionsBool("-use_ptap","Use ptap matrix for DSG","",user.use_ptap,&user.use_ptap,NULL);
146:   user.use_lrc = PETSC_FALSE;
147:   PetscOptionsBool("-use_lrc","Use lrc matrix for Js","",user.use_lrc,&user.use_lrc,NULL);

149:   user.m = user.ns*user.mx*user.mx*user.mx; /* number of constraints */
150:   user.nstate =  user.m;
151:   user.ndesign = user.mx*user.mx*user.mx;
152:   user.n = user.nstate + user.ndesign; /* number of variables */

154:   /* Create TAO solver and set desired solution method */
155:   TaoCreate(PETSC_COMM_WORLD,&tao);
156:   TaoSetType(tao,TAOLCL);

158:   /* Set up initial vectors and matrices */
159:   EllipticInitialize(&user);

161:   Gather(user.x,user.y,user.state_scatter,user.u,user.design_scatter);
162:   VecDuplicate(user.x,&x0);
163:   VecCopy(user.x,x0);

165:   /* Set solution vector with an initial guess */
166:   TaoSetInitialVector(tao,user.x);
167:   TaoSetObjectiveRoutine(tao, FormFunction, (void *)&user);
168:   TaoSetGradientRoutine(tao, FormGradient, (void *)&user);
169:   TaoSetConstraintsRoutine(tao, user.c, FormConstraints, (void *)&user);

171:   TaoSetJacobianStateRoutine(tao, user.Js, NULL, user.JsInv, FormJacobianState, (void *)&user);
172:   TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, (void *)&user);

174:   TaoSetStateDesignIS(tao,user.s_is,user.d_is);
175:   TaoSetFromOptions(tao);

177:   /* SOLVE THE APPLICATION */
178:   PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL);
179:   PetscLogStageRegister("Trials",&user.stages[1]);
180:   PetscLogStagePush(user.stages[1]);
181:   for (i=0; i<ntests; i++){
182:     TaoSolve(tao);
183:     PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its);
184:     VecCopy(x0,user.x);
185:   }
186:   PetscLogStagePop();
187:   PetscBarrier((PetscObject)user.x);
188:   PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ");
189:   PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial);

191:   TaoGetConvergedReason(tao,&reason);
192:   if (reason < 0) {
193:     PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge.\n");
194:   } else {
195:     PetscPrintf(MPI_COMM_WORLD, "Optimization terminated with status %D.\n", reason);
196:   }

198:   TaoDestroy(&tao);
199:   VecDestroy(&x0);
200:   EllipticDestroy(&user);
201:   PetscFinalize();
202:   return 0;
203: }
204: /* ------------------------------------------------------------------- */
207: /*
208:    dwork = Qy - d
209:    lwork = L*(u-ur)
210:    f = 1/2 * (dwork.dwork + alpha*lwork.lwork)
211: */
212: PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
213: {
215:   PetscReal      d1=0,d2=0;
216:   AppCtx         *user = (AppCtx*)ptr;

219:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
220:   MatMult(user->MQ,user->y,user->dwork);
221:   VecAXPY(user->dwork,-1.0,user->d);
222:   VecDot(user->dwork,user->dwork,&d1);
223:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
224:   MatMult(user->L,user->uwork,user->lwork);
225:   VecDot(user->lwork,user->lwork,&d2);
226:   *f = 0.5 * (d1 + user->alpha*d2);
227:   return(0);
228: }

230: /* ------------------------------------------------------------------- */
233: /*
234:     state: g_s = Q' *(Qy - d)
235:     design: g_d = alpha*L'*L*(u-ur)
236: */
237: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
238: {
240:   AppCtx         *user = (AppCtx*)ptr;

243:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
244:   MatMult(user->MQ,user->y,user->dwork);
245:   VecAXPY(user->dwork,-1.0,user->d);
246:   MatMultTranspose(user->MQ,user->dwork,user->ywork);
247:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
248:   MatMult(user->L,user->uwork,user->lwork);
249:   MatMultTranspose(user->L,user->lwork,user->uwork);
250:   VecScale(user->uwork, user->alpha);
251:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
252:   return(0);
253: }

257: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
258: {
260:   PetscReal      d1,d2;
261:   AppCtx         *user = (AppCtx*)ptr;

264:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
265:   MatMult(user->MQ,user->y,user->dwork);
266:   VecAXPY(user->dwork,-1.0,user->d);
267:   VecDot(user->dwork,user->dwork,&d1);
268:   MatMultTranspose(user->MQ,user->dwork,user->ywork);

270:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
271:   MatMult(user->L,user->uwork,user->lwork);
272:   VecDot(user->lwork,user->lwork,&d2);
273:   MatMultTranspose(user->L,user->lwork,user->uwork);
274:   VecScale(user->uwork, user->alpha);
275:   *f = 0.5 * (d1 + user->alpha*d2);
276:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
277:   return(0);
278: }

280: /* ------------------------------------------------------------------- */
283: /* A
284: MatShell object
285: */
286: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
287: {
289:   AppCtx         *user = (AppCtx*)ptr;

292:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
293:   /* DSG = Div * (1/Av_u) * Grad */
294:   VecSet(user->uwork,0);
295:   VecAXPY(user->uwork,-1.0,user->u);
296:   VecExp(user->uwork);
297:   MatMult(user->Av,user->uwork,user->Av_u);
298:   VecCopy(user->Av_u,user->Swork);
299:   VecReciprocal(user->Swork);
300:   if (user->use_ptap) {
301:     MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);
302:     MatPtAP(user->Diag,user->Grad,MAT_REUSE_MATRIX,1.0,&user->DSG);
303:   } else {
304:     MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
305:     MatDiagonalScale(user->Divwork,NULL,user->Swork);
306:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
307:   }
308:   return(0);
309: }
310: /* ------------------------------------------------------------------- */
313: /* B */
314: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
315: {
317:   AppCtx         *user = (AppCtx*)ptr;

320:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
321:   return(0);
322: }

326: PetscErrorCode StateBlockMatMult(Mat J_shell, Vec X, Vec Y)
327: {
329:   PetscReal      sum;
330:   AppCtx         *user;

333:   MatShellGetContext(J_shell,(void**)&user);
334:   MatMult(user->DSG,X,Y);
335:   VecSum(X,&sum);
336:   sum /= user->ndesign;
337:   VecShift(Y,sum);
338:   return(0);
339: }

343: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
344: {
346:   PetscInt       i;
347:   AppCtx         *user;

350:   MatShellGetContext(J_shell,(void**)&user);
351:   if (user->ns == 1) {
352:     MatMult(user->JsBlock,X,Y);
353:   } else {
354:     for (i=0;i<user->ns;i++) {
355:       Scatter(X,user->subq,user->yi_scatter[i],0,0);
356:       Scatter(Y,user->suby,user->yi_scatter[i],0,0);
357:       MatMult(user->JsBlock,user->subq,user->suby);
358:       Gather(X,user->subq,user->yi_scatter[i],0,0);
359:       Gather(Y,user->suby,user->yi_scatter[i],0,0);
360:     }
361:   }
362:   return(0);
363: }

367: PetscErrorCode StateInvMatMult(Mat J_shell, Vec X, Vec Y)
368: {
370:   PetscInt       its,i;
371:   AppCtx         *user;

374:   MatShellGetContext(J_shell,(void**)&user);
375:   KSPSetOperators(user->solver,user->JsBlock,user->DSG);
376:   if (Y == user->ytrue) {
377:     /* First solve is done using true solution to set up problem */
378:     KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
379:   } else {
380:     KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
381:   }
382:   if (user->ns == 1) {
383:     KSPSolve(user->solver,X,Y);
384:     KSPGetIterationNumber(user->solver,&its);
385:     user->ksp_its+=its;
386:   } else {
387:     for (i=0;i<user->ns;i++) {
388:       Scatter(X,user->subq,user->yi_scatter[i],0,0);
389:       Scatter(Y,user->suby,user->yi_scatter[i],0,0);
390:       KSPSolve(user->solver,user->subq,user->suby);
391:       KSPGetIterationNumber(user->solver,&its);
392:       user->ksp_its+=its;
393:       Gather(X,user->subq,user->yi_scatter[i],0,0);
394:       Gather(Y,user->suby,user->yi_scatter[i],0,0);
395:     }
396:   }
397:   return(0);
398: }
401: PetscErrorCode QMatMult(Mat J_shell, Vec X, Vec Y)
402: {
404:   AppCtx         *user;
405:   PetscInt       i;

408:   MatShellGetContext(J_shell,(void**)&user);
409:   if (user->ns == 1) {
410:     MatMult(user->Q,X,Y);
411:   } else {
412:     for (i=0;i<user->ns;i++) {
413:       Scatter(X,user->subq,user->yi_scatter[i],0,0);
414:       Scatter(Y,user->subd,user->di_scatter[i],0,0);
415:       MatMult(user->Q,user->subq,user->subd);
416:       Gather(X,user->subq,user->yi_scatter[i],0,0);
417:       Gather(Y,user->subd,user->di_scatter[i],0,0);
418:     }
419:   }
420:   return(0);
421: }

425: PetscErrorCode QMatMultTranspose(Mat J_shell, Vec X, Vec Y)
426: {
428:   AppCtx         *user;
429:   PetscInt       i;

432:   MatShellGetContext(J_shell,(void**)&user);
433:   if (user->ns == 1) {
434:     MatMultTranspose(user->Q,X,Y);
435:   } else {
436:     for (i=0;i<user->ns;i++) {
437:       Scatter(X,user->subd,user->di_scatter[i],0,0);
438:       Scatter(Y,user->suby,user->yi_scatter[i],0,0);
439:       MatMultTranspose(user->Q,user->subd,user->suby);
440:       Gather(X,user->subd,user->di_scatter[i],0,0);
441:       Gather(Y,user->suby,user->yi_scatter[i],0,0);
442:     }
443:   }
444:   return(0);
445: }

449: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
450: {
452:   PetscInt       i;
453:   AppCtx         *user;

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

458:   /* sdiag(1./v) */
459:   VecSet(user->uwork,0);
460:   VecAXPY(user->uwork,-1.0,user->u);
461:   VecExp(user->uwork);

463:   /* sdiag(1./((Av*(1./v)).^2)) */
464:   MatMult(user->Av,user->uwork,user->Swork);
465:   VecPointwiseMult(user->Swork,user->Swork,user->Swork);
466:   VecReciprocal(user->Swork);

468:   /* (Av * (sdiag(1./v) * b)) */
469:   VecPointwiseMult(user->uwork,user->uwork,X);
470:   MatMult(user->Av,user->uwork,user->Twork);

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

475:   if (user->ns == 1) {
476:     /* (sdiag(Grad*y(:,i)) */
477:     MatMult(user->Grad,user->y,user->Twork);

479:     /* Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
480:     VecPointwiseMult(user->Swork,user->Twork,user->Swork);
481:     MatMultTranspose(user->Grad,user->Swork,Y);
482:   } else {
483:     for (i=0;i<user->ns;i++) {
484:       Scatter(user->y,user->suby,user->yi_scatter[i],0,0);
485:       Scatter(Y,user->subq,user->yi_scatter[i],0,0);

487:       MatMult(user->Grad,user->suby,user->Twork);
488:       VecPointwiseMult(user->Twork,user->Twork,user->Swork);
489:       MatMultTranspose(user->Grad,user->Twork,user->subq);
490:       Gather(user->y,user->suby,user->yi_scatter[i],0,0);
491:       Gather(Y,user->subq,user->yi_scatter[i],0,0);
492:     }
493:   }
494:   return(0);
495: }

499: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
500: {
502:   PetscInt       i;
503:   AppCtx         *user;

506:   MatShellGetContext(J_shell,(void**)&user);
507:   VecZeroEntries(Y);

509:   /* Sdiag = 1./((Av*(1./v)).^2) */
510:   VecSet(user->uwork,0);
511:   VecAXPY(user->uwork,-1.0,user->u);
512:   VecExp(user->uwork);
513:   MatMult(user->Av,user->uwork,user->Swork);
514:   VecPointwiseMult(user->Sdiag,user->Swork,user->Swork);
515:   VecReciprocal(user->Sdiag);

517:   for (i=0;i<user->ns;i++) {
518:     Scatter(X,user->subq,user->yi_scatter[i],0,0);
519:     Scatter(user->y,user->suby,user->yi_scatter[i],0,0);

521:     /* Swork = (Div' * b(:,i)) */
522:     MatMult(user->Grad,user->subq,user->Swork);

524:     /* Twork = Grad*y(:,i) */
525:     MatMult(user->Grad,user->suby,user->Twork);

527:     /* Twork = sdiag(Twork) * Swork */
528:     VecPointwiseMult(user->Twork,user->Swork,user->Twork);


531:     /* Swork = pointwisemult(Sdiag,Twork) */
532:     VecPointwiseMult(user->Swork,user->Twork,user->Sdiag);

534:     /* Ywork = Av' * Swork */
535:     MatMultTranspose(user->Av,user->Swork,user->Ywork);

537:     /* Ywork = pointwisemult(uwork,Ywork) */
538:     VecPointwiseMult(user->Ywork,user->uwork,user->Ywork);
539:     VecAXPY(Y,1.0,user->Ywork);
540:     Gather(X,user->subq,user->yi_scatter[i],0,0);
541:     Gather(user->y,user->suby,user->yi_scatter[i],0,0);
542:   }
543:   return(0);
544: }

548: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
549: {
550:    /* C=Ay - q      A = Div * Sigma * Grad + hx*hx*hx*ones(n,n) */
552:    PetscReal      sum;
553:    PetscInt       i;
554:    AppCtx         *user = (AppCtx*)ptr;

557:    Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
558:    if (user->ns == 1) {
559:      MatMult(user->Grad,user->y,user->Swork);
560:      VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
561:      MatMultTranspose(user->Grad,user->Swork,C);
562:      VecSum(user->y,&sum);
563:      sum /= user->ndesign;
564:      VecShift(C,sum);
565:    } else {
566:      for (i=0;i<user->ns;i++) {
567:       Scatter(user->y,user->suby,user->yi_scatter[i],0,0);
568:       Scatter(C,user->subq,user->yi_scatter[i],0,0);
569:       MatMult(user->Grad,user->suby,user->Swork);
570:       VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
571:       MatMultTranspose(user->Grad,user->Swork,user->subq);

573:       VecSum(user->suby,&sum);
574:       sum /= user->ndesign;
575:       VecShift(user->subq,sum);

577:       Gather(user->y,user->suby,user->yi_scatter[i],0,0);
578:       Gather(C,user->subq,user->yi_scatter[i],0,0);
579:      }
580:    }
581:    VecAXPY(C,-1.0,user->q);
582:    return(0);
583: }

587: PetscErrorCode Scatter(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
588: {

592:   VecScatterBegin(scat1,x,sub1,INSERT_VALUES,SCATTER_FORWARD);
593:   VecScatterEnd(scat1,x,sub1,INSERT_VALUES,SCATTER_FORWARD);
594:   if (sub2) {
595:     VecScatterBegin(scat2,x,sub2,INSERT_VALUES,SCATTER_FORWARD);
596:     VecScatterEnd(scat2,x,sub2,INSERT_VALUES,SCATTER_FORWARD);
597:   }
598:   return(0);
599: }

603: PetscErrorCode Gather(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
604: {

608:   VecScatterBegin(scat1,sub1,x,INSERT_VALUES,SCATTER_REVERSE);
609:   VecScatterEnd(scat1,sub1,x,INSERT_VALUES,SCATTER_REVERSE);
610:   if (sub2) {
611:     VecScatterBegin(scat2,sub2,x,INSERT_VALUES,SCATTER_REVERSE);
612:     VecScatterEnd(scat2,sub2,x,INSERT_VALUES,SCATTER_REVERSE);
613:   }
614:   return(0);
615: }

619: PetscErrorCode EllipticInitialize(AppCtx *user)
620: {
622:   PetscInt       m,n,i,j,k,l,linear_index,is,js,ks,ls,istart,iend,iblock;
623:   Vec            XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork;
624:   PetscReal      *x,*y,*z;
625:   PetscReal      h,meanut;
626:   PetscScalar    hinv,neg_hinv,half = 0.5,sqrt_beta;
627:   PetscInt       im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
628:   IS             is_alldesign,is_allstate;
629:   IS             is_from_d;
630:   IS             is_from_y;
631:   PetscInt       lo,hi,hi2,lo2,ysubnlocal,dsubnlocal;
632:   const PetscInt *ranges, *subranges;
633:   PetscMPIInt    size;
634:   PetscReal      xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
635:   PetscScalar    v,vx,vy,vz;
636:   PetscInt       offset,subindex,subvec,nrank,kk;

638:   PetscScalar xr[64] = {0.4970,     0.8498,     0.7814,     0.6268,     0.7782,     0.6402,     0.3617,     0.3160,
639:                         0.3610,     0.5298,     0.6987,     0.3331,     0.7962,     0.5596,     0.3866,     0.6774,
640:                         0.5407,     0.4518,     0.6702,     0.6061,     0.7580,     0.8997,     0.5198,     0.8326,
641:                         0.2138,     0.9198,     0.3000,     0.2833,     0.8288,     0.7076,     0.1820,     0.0728,
642:                         0.8447,     0.2367,     0.3239,     0.6413,     0.3114,     0.4731,     0.1192,     0.9273,
643:                         0.5724,     0.4331,     0.5136,     0.3547,     0.4413,     0.2602,     0.5698,     0.7278,
644:                         0.5261,     0.6230,     0.2454,     0.3948,     0.7479,     0.6582,     0.4660,     0.5594,
645:                         0.7574,     0.1143,     0.5900,     0.1065,     0.4260,     0.3294,     0.8276,     0.0756};

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

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

666:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
667:   PetscLogStageRegister("Elliptic Setup",&user->stages[0]);
668:   PetscLogStagePush(user->stages[0]);

670:   /* Create u,y,c,x */
671:   VecCreate(PETSC_COMM_WORLD,&user->u);
672:   VecCreate(PETSC_COMM_WORLD,&user->y);
673:   VecCreate(PETSC_COMM_WORLD,&user->c);
674:   VecSetSizes(user->u,PETSC_DECIDE,user->ndesign);
675:   VecSetFromOptions(user->u);
676:   VecGetLocalSize(user->u,&ysubnlocal);
677:   VecSetSizes(user->y,ysubnlocal*user->ns,user->nstate);
678:   VecSetSizes(user->c,ysubnlocal*user->ns,user->m);
679:   VecSetFromOptions(user->y);
680:   VecSetFromOptions(user->c);

682:   /*
683:      *******************************
684:      Create scatters for x <-> y,u
685:      *******************************

687:      If the state vector y and design vector u are partitioned as
688:      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
689:      then the solution vector x is organized as
690:      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
691:      The index sets user->s_is and user->d_is correspond to the indices of the
692:      state and design variables owned by the current processor.
693:   */
694:   VecCreate(PETSC_COMM_WORLD,&user->x);

696:   VecGetOwnershipRange(user->y,&lo,&hi);
697:   VecGetOwnershipRange(user->u,&lo2,&hi2);

699:   ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate);
700:   ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user->s_is);
701:   ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign);
702:   ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user->d_is);

704:   VecSetSizes(user->x,hi-lo+hi2-lo2,user->n);
705:   VecSetFromOptions(user->x);

707:   VecScatterCreate(user->x,user->s_is,user->y,is_allstate,&user->state_scatter);
708:   VecScatterCreate(user->x,user->d_is,user->u,is_alldesign,&user->design_scatter);
709:   ISDestroy(&is_alldesign);
710:   ISDestroy(&is_allstate);
711:   /*
712:      *******************************
713:      Create scatter from y to y_1,y_2,...,y_ns
714:      *******************************
715:   */
716:   PetscMalloc1(user->ns,&user->yi_scatter);
717:   VecDuplicate(user->u,&user->suby);
718:   VecDuplicate(user->u,&user->subq);

720:   VecGetOwnershipRange(user->y,&lo2,&hi2);
721:   istart = 0;
722:   for (i=0; i<user->ns; i++){
723:     VecGetOwnershipRange(user->suby,&lo,&hi);
724:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
725:     VecScatterCreate(user->y,is_from_y,user->suby,NULL,&user->yi_scatter[i]);
726:     istart = istart + hi-lo;
727:     ISDestroy(&is_from_y);
728:   }
729:   /*
730:      *******************************
731:      Create scatter from d to d_1,d_2,...,d_ns
732:      *******************************
733:   */
734:   VecCreate(PETSC_COMM_WORLD,&user->subd);
735:   VecSetSizes(user->subd,PETSC_DECIDE,user->ndata);
736:   VecSetFromOptions(user->subd);
737:   VecCreate(PETSC_COMM_WORLD,&user->d);
738:   VecGetLocalSize(user->subd,&dsubnlocal);
739:   VecSetSizes(user->d,dsubnlocal*user->ns,user->ndata*user->ns);
740:   VecSetFromOptions(user->d);
741:   PetscMalloc1(user->ns,&user->di_scatter);

743:   VecGetOwnershipRange(user->d,&lo2,&hi2);
744:   istart = 0;
745:   for (i=0; i<user->ns; i++){
746:     VecGetOwnershipRange(user->subd,&lo,&hi);
747:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d);
748:     VecScatterCreate(user->d,is_from_d,user->subd,NULL,&user->di_scatter[i]);
749:     istart = istart + hi-lo;
750:     ISDestroy(&is_from_d);
751:   }

753:   PetscMalloc1(user->mx,&x);
754:   PetscMalloc1(user->mx,&y);
755:   PetscMalloc1(user->mx,&z);

757:   user->ksp_its = 0;
758:   user->ksp_its_initial = 0;

760:   n = user->mx * user->mx * user->mx;
761:   m = 3 * user->mx * user->mx * (user->mx-1);
762:   sqrt_beta = PetscSqrtScalar(user->beta);

764:   VecCreate(PETSC_COMM_WORLD,&XX);
765:   VecCreate(PETSC_COMM_WORLD,&user->q);
766:   VecSetSizes(XX,ysubnlocal,n);
767:   VecSetSizes(user->q,ysubnlocal*user->ns,user->m);
768:   VecSetFromOptions(XX);
769:   VecSetFromOptions(user->q);

771:   VecDuplicate(XX,&YY);
772:   VecDuplicate(XX,&ZZ);
773:   VecDuplicate(XX,&XXwork);
774:   VecDuplicate(XX,&YYwork);
775:   VecDuplicate(XX,&ZZwork);
776:   VecDuplicate(XX,&UTwork);
777:   VecDuplicate(XX,&user->utrue);

779:   /* map for striding q */
780:   VecGetOwnershipRanges(user->q,&ranges);
781:   VecGetOwnershipRanges(user->u,&subranges);

783:   VecGetOwnershipRange(user->q,&lo2,&hi2);
784:   VecGetOwnershipRange(user->u,&lo,&hi);
785:   /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
786:   h = 1.0/user->mx;
787:   hinv = user->mx;
788:   neg_hinv = -hinv;

790:   VecGetOwnershipRange(XX,&istart,&iend);
791:   for (linear_index=istart; linear_index<iend; linear_index++){
792:     i = linear_index % user->mx;
793:     j = ((linear_index-i)/user->mx) % user->mx;
794:     k = ((linear_index-i)/user->mx-j) / user->mx;
795:     vx = h*(i+0.5);
796:     vy = h*(j+0.5);
797:     vz = h*(k+0.5);
798:     VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
799:     VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
800:     VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES);
801:     for (is=0; is<2; is++){
802:       for (js=0; js<2; js++){
803:         for (ks=0; ks<2; ks++){
804:           ls = is*4 + js*2 + ks;
805:           if (ls<user->ns){
806:             l =ls*n + linear_index;
807:             /* remap */
808:             subindex = l%n;
809:             subvec = l/n;
810:             nrank=0;
811:             while (subindex >= subranges[nrank+1]) nrank++;
812:             offset = subindex - subranges[nrank];
813:             istart=0;
814:             for (kk=0;kk<nrank;kk++) istart+=user->ns*(subranges[kk+1]-subranges[kk]);
815:             istart += (subranges[nrank+1]-subranges[nrank])*subvec;
816:             l = istart+offset;
817:             v = 100*PetscSinScalar(2*PETSC_PI*(vx+0.25*is))*sin(2*PETSC_PI*(vy+0.25*js))*sin(2*PETSC_PI*(vz+0.25*ks));
818:             VecSetValues(user->q,1,&l,&v,INSERT_VALUES);
819:           }
820:         }
821:       }
822:     }
823:   }

825:   VecAssemblyBegin(XX);
826:   VecAssemblyEnd(XX);
827:   VecAssemblyBegin(YY);
828:   VecAssemblyEnd(YY);
829:   VecAssemblyBegin(ZZ);
830:   VecAssemblyEnd(ZZ);
831:   VecAssemblyBegin(user->q);
832:   VecAssemblyEnd(user->q);

834:   /* Compute true parameter function
835:      ut = exp(-((x-0.25)^2+(y-0.25)^2+(z-0.25)^2)/0.05) - exp((x-0.75)^2-(y-0.75)^2-(z-0.75))^2/0.05) */
836:   VecCopy(XX,XXwork);
837:   VecCopy(YY,YYwork);
838:   VecCopy(ZZ,ZZwork);

840:   VecShift(XXwork,-0.25);
841:   VecShift(YYwork,-0.25);
842:   VecShift(ZZwork,-0.25);

844:   VecPointwiseMult(XXwork,XXwork,XXwork);
845:   VecPointwiseMult(YYwork,YYwork,YYwork);
846:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

848:   VecCopy(XXwork,UTwork);
849:   VecAXPY(UTwork,1.0,YYwork);
850:   VecAXPY(UTwork,1.0,ZZwork);
851:   VecScale(UTwork,-20.0);
852:   VecExp(UTwork);
853:   VecCopy(UTwork,user->utrue);

855:   VecCopy(XX,XXwork);
856:   VecCopy(YY,YYwork);
857:   VecCopy(ZZ,ZZwork);

859:   VecShift(XXwork,-0.75);
860:   VecShift(YYwork,-0.75);
861:   VecShift(ZZwork,-0.75);

863:   VecPointwiseMult(XXwork,XXwork,XXwork);
864:   VecPointwiseMult(YYwork,YYwork,YYwork);
865:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

867:   VecCopy(XXwork,UTwork);
868:   VecAXPY(UTwork,1.0,YYwork);
869:   VecAXPY(UTwork,1.0,ZZwork);
870:   VecScale(UTwork,-20.0);
871:   VecExp(UTwork);

873:   VecAXPY(user->utrue,-1.0,UTwork);

875:   VecDestroy(&XX);
876:   VecDestroy(&YY);
877:   VecDestroy(&ZZ);
878:   VecDestroy(&XXwork);
879:   VecDestroy(&YYwork);
880:   VecDestroy(&ZZwork);
881:   VecDestroy(&UTwork);

883:   /* Initial guess and reference model */
884:   VecDuplicate(user->utrue,&user->ur);
885:   VecSum(user->utrue,&meanut);
886:   meanut = meanut / n;
887:   VecSet(user->ur,meanut);
888:   VecCopy(user->ur,user->u);

890:   /* Generate Grad matrix */
891:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
892:   MatSetSizes(user->Grad,PETSC_DECIDE,ysubnlocal,m,n);
893:   MatSetFromOptions(user->Grad);
894:   MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
895:   MatSeqAIJSetPreallocation(user->Grad,2,NULL);
896:   MatGetOwnershipRange(user->Grad,&istart,&iend);

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

921:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
922:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);

924:   /* Generate arithmetic averaging matrix Av */
925:   MatCreate(PETSC_COMM_WORLD,&user->Av);
926:   MatSetSizes(user->Av,PETSC_DECIDE,ysubnlocal,m,n);
927:   MatSetFromOptions(user->Av);
928:   MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
929:   MatSeqAIJSetPreallocation(user->Av,2,NULL);
930:   MatGetOwnershipRange(user->Av,&istart,&iend);

932:   for (i=istart; i<iend; i++){
933:     if (i<m/3){
934:       iblock = i / (user->mx-1);
935:       j = iblock*user->mx + (i % (user->mx-1));
936:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
937:       j = j+1;
938:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
939:     }
940:     if (i>=m/3 && i<2*m/3){
941:       iblock = (i-m/3) / (user->mx*(user->mx-1));
942:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
943:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
944:       j = j + user->mx;
945:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
946:     }
947:     if (i>=2*m/3){
948:       j = i-2*m/3;
949:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
950:       j = j + user->mx*user->mx;
951:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
952:     }
953:   }

955:   MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
956:   MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);

958:   MatCreate(PETSC_COMM_WORLD,&user->L);
959:   MatSetSizes(user->L,PETSC_DECIDE,ysubnlocal,m+n,n);
960:   MatSetFromOptions(user->L);
961:   MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
962:   MatSeqAIJSetPreallocation(user->L,2,NULL);
963:   MatGetOwnershipRange(user->L,&istart,&iend);

965:   for (i=istart; i<iend; i++){
966:     if (i<m/3){
967:       iblock = i / (user->mx-1);
968:       j = iblock*user->mx + (i % (user->mx-1));
969:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
970:       j = j+1;
971:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
972:     }
973:     if (i>=m/3 && i<2*m/3){
974:       iblock = (i-m/3) / (user->mx*(user->mx-1));
975:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
976:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
977:       j = j + user->mx;
978:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
979:     }
980:     if (i>=2*m/3 && i<m){
981:       j = i-2*m/3;
982:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
983:       j = j + user->mx*user->mx;
984:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
985:     }
986:     if (i>=m){
987:       j = i - m;
988:       MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES);
989:     }
990:   }
991:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
992:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);
993:   MatScale(user->L,PetscPowScalar(h,1.5));

995:   /* Generate Div matrix */
996:   if (!user->use_ptap) {
997:     /* Generate Div matrix */
998:     MatCreate(PETSC_COMM_WORLD,&user->Div);
999:     MatSetSizes(user->Div,ysubnlocal,PETSC_DECIDE,n,m);
1000:     MatSetFromOptions(user->Div);
1001:     MatMPIAIJSetPreallocation(user->Div,4,NULL,4,NULL);
1002:     MatSeqAIJSetPreallocation(user->Div,6,NULL);
1003:     MatGetOwnershipRange(user->Grad,&istart,&iend);

1005:     for (i=istart; i<iend; i++){
1006:       if (i<m/3){
1007:         iblock = i / (user->mx-1);
1008:         j = iblock*user->mx + (i % (user->mx-1));
1009:         MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES);
1010:         j = j+1;
1011:         MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES);
1012:       }
1013:       if (i>=m/3 && i<2*m/3){
1014:         iblock = (i-m/3) / (user->mx*(user->mx-1));
1015:         j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
1016:         MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES);
1017:         j = j + user->mx;
1018:         MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES);
1019:       }
1020:       if (i>=2*m/3){
1021:         j = i-2*m/3;
1022:         MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES);
1023:         j = j + user->mx*user->mx;
1024:         MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES);
1025:       }
1026:     }

1028:     MatAssemblyBegin(user->Div,MAT_FINAL_ASSEMBLY);
1029:     MatAssemblyEnd(user->Div,MAT_FINAL_ASSEMBLY);
1030:     MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
1031:   } else {
1032:     MatCreate(PETSC_COMM_WORLD,&user->Diag);
1033:     MatSetSizes(user->Diag,PETSC_DECIDE,PETSC_DECIDE,m,m);
1034:     MatSetFromOptions(user->Diag);
1035:     MatMPIAIJSetPreallocation(user->Diag,1,NULL,0,NULL);
1036:     MatSeqAIJSetPreallocation(user->Diag,1,NULL);
1037:   }

1039:   /* Build work vectors and matrices */
1040:   VecCreate(PETSC_COMM_WORLD,&user->S);
1041:   VecSetSizes(user->S, PETSC_DECIDE, m);
1042:   VecSetFromOptions(user->S);

1044:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1045:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
1046:   VecSetFromOptions(user->lwork);

1048:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

1050:   VecDuplicate(user->S,&user->Swork);
1051:   VecDuplicate(user->S,&user->Sdiag);
1052:   VecDuplicate(user->S,&user->Av_u);
1053:   VecDuplicate(user->S,&user->Twork);
1054:   VecDuplicate(user->y,&user->ywork);
1055:   VecDuplicate(user->u,&user->Ywork);
1056:   VecDuplicate(user->u,&user->uwork);
1057:   VecDuplicate(user->u,&user->js_diag);
1058:   VecDuplicate(user->c,&user->cwork);
1059:   VecDuplicate(user->d,&user->dwork);

1061:   /* Create a matrix-free shell user->Jd for computing B*x */
1062:   MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal,user->nstate,user->ndesign,user,&user->Jd);
1063:   MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult);
1064:   MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose);

1066:   /* Compute true state function ytrue given utrue */
1067:   VecDuplicate(user->y,&user->ytrue);

1069:   /* First compute Av_u = Av*exp(-u) */
1070:   VecSet(user->uwork, 0);
1071:   VecAXPY(user->uwork,-1.0,user->utrue); /* Note: user->utrue */
1072:   VecExp(user->uwork);
1073:   MatMult(user->Av,user->uwork,user->Av_u);

1075:   /* Next form DSG = Div*S*Grad */
1076:   VecCopy(user->Av_u,user->Swork);
1077:   VecReciprocal(user->Swork);
1078:   if (user->use_ptap) {
1079:     MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);
1080:     MatPtAP(user->Diag,user->Grad,MAT_INITIAL_MATRIX,1.0,&user->DSG);
1081:   } else {
1082:     MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1083:     MatDiagonalScale(user->Divwork,NULL,user->Swork);
1084:     MatMatMultSymbolic(user->Divwork,user->Grad,1.0,&user->DSG);
1085:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1086:   }

1088:   MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE);
1089:   MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

1091:   if (user->use_lrc == PETSC_TRUE) {
1092:     v=PetscSqrtReal(1.0 /user->ndesign);
1093:     PetscMalloc1(user->ndesign,&user->ones);

1095:     for (i=0;i<user->ndesign;i++) {
1096:       user->ones[i]=v;
1097:     }
1098:     MatCreateDense(PETSC_COMM_WORLD,ysubnlocal,PETSC_DECIDE,user->ndesign,1,user->ones,&user->Ones);
1099:     MatAssemblyBegin(user->Ones, MAT_FINAL_ASSEMBLY);
1100:     MatAssemblyEnd(user->Ones, MAT_FINAL_ASSEMBLY);
1101:     MatCreateLRC(user->DSG,user->Ones,user->Ones,&user->JsBlock);
1102:     MatSetUp(user->JsBlock);
1103:   } else {
1104:     /* Create matrix-free shell user->Js for computing (A + h^3*e*e^T)*x */
1105:     MatCreateShell(PETSC_COMM_WORLD,ysubnlocal,ysubnlocal,user->ndesign,user->ndesign,user,&user->JsBlock);
1106:     MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateBlockMatMult);
1107:     MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateBlockMatMult);
1108:   }
1109:   MatSetOption(user->JsBlock,MAT_SYMMETRIC,PETSC_TRUE);
1110:   MatSetOption(user->JsBlock,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
1111:   MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal*user->ns,user->nstate,user->nstate,user,&user->Js);
1112:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1113:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMult);
1114:   MatSetOption(user->Js,MAT_SYMMETRIC,PETSC_TRUE);
1115:   MatSetOption(user->Js,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

1117:   MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal*user->ns,user->nstate,user->nstate,user,&user->JsInv);
1118:   MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateInvMatMult);
1119:   MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateInvMatMult);
1120:   MatSetOption(user->JsInv,MAT_SYMMETRIC,PETSC_TRUE);
1121:   MatSetOption(user->JsInv,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

1123:   MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE);
1124:   MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
1125:   /* Now solve for ytrue */
1126:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1127:   KSPSetFromOptions(user->solver);

1129:   KSPSetOperators(user->solver,user->JsBlock,user->DSG);

1131:   MatMult(user->JsInv,user->q,user->ytrue);
1132:   /* First compute Av_u = Av*exp(-u) */
1133:   VecSet(user->uwork,0);
1134:   VecAXPY(user->uwork,-1.0,user->u); /* Note: user->u */
1135:   VecExp(user->uwork);
1136:   MatMult(user->Av,user->uwork,user->Av_u);

1138:   /* Next update DSG = Div*S*Grad  with user->u */
1139:   VecCopy(user->Av_u,user->Swork);
1140:   VecReciprocal(user->Swork);
1141:   if (user->use_ptap) {
1142:     MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);
1143:     MatPtAP(user->Diag,user->Grad,MAT_REUSE_MATRIX,1.0,&user->DSG);
1144:   } else {
1145:     MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1146:     MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1147:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1148:   }

1150:   /* Now solve for y */

1152:   MatMult(user->JsInv,user->q,user->y);

1154:   user->ksp_its_initial = user->ksp_its;
1155:   user->ksp_its = 0;
1156:   /* Construct projection matrix Q (blocks) */
1157:   MatCreate(PETSC_COMM_WORLD,&user->Q);
1158:   MatSetSizes(user->Q,dsubnlocal,ysubnlocal,user->ndata,user->ndesign);
1159:   MatSetFromOptions(user->Q);
1160:   MatMPIAIJSetPreallocation(user->Q,8,NULL,8,NULL);
1161:   MatSeqAIJSetPreallocation(user->Q,8,NULL);

1163:   for (i=0; i<user->mx; i++){
1164:     x[i] = h*(i+0.5);
1165:     y[i] = h*(i+0.5);
1166:     z[i] = h*(i+0.5);
1167:   }
1168:   MatGetOwnershipRange(user->Q,&istart,&iend);

1170:   nx = user->mx; ny = user->mx; nz = user->mx;
1171:   for (i=istart; i<iend; i++){

1173:     xri = xr[i];
1174:     im = 0;
1175:     xim = x[im];
1176:     while (xri>xim && im<nx){
1177:       im = im+1;
1178:       xim = x[im];
1179:     }
1180:     indx1 = im-1;
1181:     indx2 = im;
1182:     dx1 = xri - x[indx1];
1183:     dx2 = x[indx2] - xri;

1185:     yri = yr[i];
1186:     im = 0;
1187:     yim = y[im];
1188:     while (yri>yim && im<ny){
1189:       im = im+1;
1190:       yim = y[im];
1191:     }
1192:     indy1 = im-1;
1193:     indy2 = im;
1194:     dy1 = yri - y[indy1];
1195:     dy2 = y[indy2] - yri;

1197:     zri = zr[i];
1198:     im = 0;
1199:     zim = z[im];
1200:     while (zri>zim && im<nz){
1201:       im = im+1;
1202:       zim = z[im];
1203:     }
1204:     indz1 = im-1;
1205:     indz2 = im;
1206:     dz1 = zri - z[indz1];
1207:     dz2 = z[indz2] - zri;

1209:     Dx = x[indx2] - x[indx1];
1210:     Dy = y[indy2] - y[indy1];
1211:     Dz = z[indz2] - z[indz1];

1213:     j = indx1 + indy1*nx + indz1*nx*ny;
1214:     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1215:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);

1217:     j = indx1 + indy1*nx + indz2*nx*ny;
1218:     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1219:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);

1221:     j = indx1 + indy2*nx + indz1*nx*ny;
1222:     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1223:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);

1225:     j = indx1 + indy2*nx + indz2*nx*ny;
1226:     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1227:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);

1229:     j = indx2 + indy1*nx + indz1*nx*ny;
1230:     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1231:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);

1233:     j = indx2 + indy1*nx + indz2*nx*ny;
1234:     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1235:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);

1237:     j = indx2 + indy2*nx + indz1*nx*ny;
1238:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1239:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);

1241:     j = indx2 + indy2*nx + indz2*nx*ny;
1242:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1243:     MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);
1244:   }

1246:   MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY);
1247:   MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY);
1248:   /* Create MQ (composed of blocks of Q */
1249:   MatCreateShell(PETSC_COMM_WORLD,dsubnlocal*user->ns,PETSC_DECIDE,user->ndata*user->ns,user->nstate,user,&user->MQ);
1250:   MatShellSetOperation(user->MQ,MATOP_MULT,(void(*)(void))QMatMult);
1251:   MatShellSetOperation(user->MQ,MATOP_MULT_TRANSPOSE,(void(*)(void))QMatMultTranspose);

1253:   /* Add noise to the measurement data */
1254:   VecSet(user->ywork,1.0);
1255:   VecAYPX(user->ywork,user->noise,user->ytrue);
1256:   MatMult(user->MQ,user->ywork,user->d);

1258:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1259:   PetscFree(x);
1260:   PetscFree(y);
1261:   PetscFree(z);
1262:   PetscLogStagePop();
1263:   return(0);
1264: }

1268: PetscErrorCode EllipticDestroy(AppCtx *user)
1269: {
1271:   PetscInt       i;

1274:   MatDestroy(&user->DSG);
1275:   KSPDestroy(&user->solver);
1276:   MatDestroy(&user->Q);
1277:   MatDestroy(&user->MQ);
1278:   if (!user->use_ptap) {
1279:     MatDestroy(&user->Div);
1280:     MatDestroy(&user->Divwork);
1281:   } else {
1282:     MatDestroy(&user->Diag);
1283:   }
1284:   if (user->use_lrc) {
1285:     MatDestroy(&user->Ones);
1286:   }

1288:   MatDestroy(&user->Grad);
1289:   MatDestroy(&user->Av);
1290:   MatDestroy(&user->Avwork);
1291:   MatDestroy(&user->L);
1292:   MatDestroy(&user->Js);
1293:   MatDestroy(&user->Jd);
1294:   MatDestroy(&user->JsBlock);
1295:   MatDestroy(&user->JsInv);

1297:   VecDestroy(&user->x);
1298:   VecDestroy(&user->u);
1299:   VecDestroy(&user->uwork);
1300:   VecDestroy(&user->utrue);
1301:   VecDestroy(&user->y);
1302:   VecDestroy(&user->ywork);
1303:   VecDestroy(&user->ytrue);
1304:   VecDestroy(&user->c);
1305:   VecDestroy(&user->cwork);
1306:   VecDestroy(&user->ur);
1307:   VecDestroy(&user->q);
1308:   VecDestroy(&user->d);
1309:   VecDestroy(&user->dwork);
1310:   VecDestroy(&user->lwork);
1311:   VecDestroy(&user->S);
1312:   VecDestroy(&user->Swork);
1313:   VecDestroy(&user->Sdiag);
1314:   VecDestroy(&user->Ywork);
1315:   VecDestroy(&user->Twork);
1316:   VecDestroy(&user->Av_u);
1317:   VecDestroy(&user->js_diag);
1318:   ISDestroy(&user->s_is);
1319:   ISDestroy(&user->d_is);
1320:   VecDestroy(&user->suby);
1321:   VecDestroy(&user->subd);
1322:   VecDestroy(&user->subq);
1323:   VecScatterDestroy(&user->state_scatter);
1324:   VecScatterDestroy(&user->design_scatter);
1325:   for (i=0;i<user->ns;i++) {
1326:     VecScatterDestroy(&user->yi_scatter[i]);
1327:     VecScatterDestroy(&user->di_scatter[i]);
1328:   }
1329:   PetscFree(user->yi_scatter);
1330:   PetscFree(user->di_scatter);
1331:   if (user->use_lrc) {
1332:     PetscFree(user->ones);
1333:     MatDestroy(&user->Ones);
1334:   }
1335:   return(0);
1336: }

1340: PetscErrorCode EllipticMonitor(Tao tao, void *ptr)
1341: {
1343:   Vec            X;
1344:   PetscReal      unorm,ynorm;
1345:   AppCtx         *user = (AppCtx*)ptr;

1348:   TaoGetSolutionVector(tao,&X);
1349:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1350:   VecAXPY(user->ywork,-1.0,user->ytrue);
1351:   VecAXPY(user->uwork,-1.0,user->utrue);
1352:   VecNorm(user->uwork,NORM_2,&unorm);
1353:   VecNorm(user->ywork,NORM_2,&ynorm);
1354:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1355:   return(0);
1356: }