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

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

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

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


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

 58:   Vec d;
 59:   Vec dwork;

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

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

 66:   Vec ytrue;

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

 71:   Vec utrue;

 73:   Vec js_diag;

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

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

 83:   KSP solver;
 84:   PC  prec;

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

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

108: PetscErrorCode StateBlockMatMult(Mat,Vec,Vec);
109: PetscErrorCode StateMatMult(Mat,Vec,Vec);

111: PetscErrorCode StateInvMatMult(Mat,Vec,Vec);
112: PetscErrorCode DesignMatMult(Mat,Vec,Vec);
113: PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);

115: PetscErrorCode QMatMult(Mat,Vec,Vec);
116: PetscErrorCode QMatMultTranspose(Mat,Vec,Vec);

118: static  char help[]="";

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

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

147:   PetscOptionsBool("-use_ptap","Use ptap matrix for DSG","",PETSC_FALSE,&user.use_ptap,&flag);
148:   PetscOptionsBool("-use_lrc","Use lrc matrix for Js","",PETSC_FALSE,&user.use_lrc,&flag);

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

155:   /* Create TAO solver and set desired solution method */
156:   TaoCreate(PETSC_COMM_WORLD,&tao);
157:   TaoSetType(tao,TAOLCL);
158:   user.lcl = (TAO_LCL*)(tao->data);

160:   /* Set up initial vectors and matrices */
161:   EllipticInitialize(&user);

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

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

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

176:   TaoSetStateDesignIS(tao,user.s_is,user.d_is);
177:   TaoSetFromOptions(tao);

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

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

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

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

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

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

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

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

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

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

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

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

328: PetscErrorCode StateBlockMatMult(Mat J_shell, Vec X, Vec Y)
329: {
331:   PetscReal      sum;
332:   AppCtx         *user;

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

345: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
346: {
348:   PetscInt       i;
349:   AppCtx         *user;

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

369: PetscErrorCode StateInvMatMult(Mat J_shell, Vec X, Vec Y)
370: {
372:   PetscInt       its,i;
373:   PetscReal      tau;
374:   AppCtx         *user;

377:   MatShellGetContext(J_shell,(void**)&user);
378:   KSPSetOperators(user->solver,user->JsBlock,user->DSG);
379:   if (Y == user->ytrue) {
380:     KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
381:   } else if (user->lcl) {
382:     tau = user->lcl->tau[user->lcl->solve_type];
383:     KSPSetTolerances(user->solver,tau,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
384:   }
385:   if (user->ns == 1) {
386:     KSPSolve(user->solver,X,Y);
387:     KSPGetIterationNumber(user->solver,&its);
388:     user->ksp_its+=its;
389:   } else {
390:     for (i=0;i<user->ns;i++) {
391:       Scatter(X,user->subq,user->yi_scatter[i],0,0);
392:       Scatter(Y,user->suby,user->yi_scatter[i],0,0);
393:       KSPSolve(user->solver,user->subq,user->suby);
394:       KSPGetIterationNumber(user->solver,&its);
395:       user->ksp_its+=its;
396:       Gather(X,user->subq,user->yi_scatter[i],0,0);
397:       Gather(Y,user->suby,user->yi_scatter[i],0,0);
398:     }
399:   }
400:   return(0);
401: }
404: PetscErrorCode QMatMult(Mat J_shell, Vec X, Vec Y)
405: {
407:   AppCtx         *user;
408:   PetscInt       i;

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

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

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

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

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

461:   /* sdiag(1./v) */
462:   VecSet(user->uwork,0);
463:   VecAXPY(user->uwork,-1.0,user->u);
464:   VecExp(user->uwork);

466:   /* sdiag(1./((Av*(1./v)).^2)) */
467:   MatMult(user->Av,user->uwork,user->Swork);
468:   VecPointwiseMult(user->Swork,user->Swork,user->Swork);
469:   VecReciprocal(user->Swork);

471:   /* (Av * (sdiag(1./v) * b)) */
472:   VecPointwiseMult(user->uwork,user->uwork,X);
473:   MatMult(user->Av,user->uwork,user->Twork);

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

478:   if (user->ns == 1) {
479:     /* (sdiag(Grad*y(:,i)) */
480:     MatMult(user->Grad,user->y,user->Twork);

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

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

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

509:   MatShellGetContext(J_shell,(void**)&user);
510:   VecZeroEntries(Y);

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

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

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

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

530:     /* Twork = sdiag(Twork) * Swork */
531:     VecPointwiseMult(user->Twork,user->Swork,user->Twork);


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

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

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

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

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

576:       VecSum(user->suby,&sum);
577:       sum /= user->ndesign;
578:       VecShift(user->subq,sum);

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

590: PetscErrorCode Scatter(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
591: {

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

606: PetscErrorCode Gather(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
607: {

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

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

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

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

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

669:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
670:   PetscLogStageRegister("Elliptic Setup",&user->stages[0]);
671:   PetscLogStagePush(user->stages[0]);

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

685:   /*
686:      *******************************
687:      Create scatters for x <-> y,u
688:      *******************************

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

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

702:   ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate);
703:   ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user->s_is);
704:   ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign);
705:   ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user->d_is);

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

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

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

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

756:   PetscMalloc1(user->mx,&x);
757:   PetscMalloc1(user->mx,&y);
758:   PetscMalloc1(user->mx,&z);

760:   user->ksp_its = 0;
761:   user->ksp_its_initial = 0;

763:   n = user->mx * user->mx * user->mx;
764:   m = 3 * user->mx * user->mx * (user->mx-1);
765:   sqrt_beta = PetscSqrtScalar(user->beta);

767:   VecCreate(PETSC_COMM_WORLD,&XX);
768:   VecCreate(PETSC_COMM_WORLD,&user->q);
769:   VecSetSizes(XX,ysubnlocal,n);
770:   VecSetSizes(user->q,ysubnlocal*user->ns,user->m);
771:   VecSetFromOptions(XX);
772:   VecSetFromOptions(user->q);

774:   VecDuplicate(XX,&YY);
775:   VecDuplicate(XX,&ZZ);
776:   VecDuplicate(XX,&XXwork);
777:   VecDuplicate(XX,&YYwork);
778:   VecDuplicate(XX,&ZZwork);
779:   VecDuplicate(XX,&UTwork);
780:   VecDuplicate(XX,&user->utrue);

782:   /* map for striding q */
783:   /*
784:   PetscMalloc1((size+1),&ranges);
785:   PetscMalloc1((size+1),&subranges); */
786:   VecGetOwnershipRanges(user->q,&ranges);
787:   VecGetOwnershipRanges(user->u,&subranges);

789:   VecGetOwnershipRange(user->q,&lo2,&hi2);
790:   VecGetOwnershipRange(user->u,&lo,&hi);
791:   /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
792:   h = 1.0/user->mx;
793:   hinv = user->mx;
794:   neg_hinv = -hinv;

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

831:   VecAssemblyBegin(XX);
832:   VecAssemblyEnd(XX);
833:   VecAssemblyBegin(YY);
834:   VecAssemblyEnd(YY);
835:   VecAssemblyBegin(ZZ);
836:   VecAssemblyEnd(ZZ);
837:   VecAssemblyBegin(user->q);
838:   VecAssemblyEnd(user->q);

840:   /* Compute true parameter function
841:      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) */
842:   VecCopy(XX,XXwork);
843:   VecCopy(YY,YYwork);
844:   VecCopy(ZZ,ZZwork);

846:   VecShift(XXwork,-0.25);
847:   VecShift(YYwork,-0.25);
848:   VecShift(ZZwork,-0.25);

850:   VecPointwiseMult(XXwork,XXwork,XXwork);
851:   VecPointwiseMult(YYwork,YYwork,YYwork);
852:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

854:   VecCopy(XXwork,UTwork);
855:   VecAXPY(UTwork,1.0,YYwork);
856:   VecAXPY(UTwork,1.0,ZZwork);
857:   VecScale(UTwork,-20.0);
858:   VecExp(UTwork);
859:   VecCopy(UTwork,user->utrue);

861:   VecCopy(XX,XXwork);
862:   VecCopy(YY,YYwork);
863:   VecCopy(ZZ,ZZwork);

865:   VecShift(XXwork,-0.75);
866:   VecShift(YYwork,-0.75);
867:   VecShift(ZZwork,-0.75);

869:   VecPointwiseMult(XXwork,XXwork,XXwork);
870:   VecPointwiseMult(YYwork,YYwork,YYwork);
871:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

873:   VecCopy(XXwork,UTwork);
874:   VecAXPY(UTwork,1.0,YYwork);
875:   VecAXPY(UTwork,1.0,ZZwork);
876:   VecScale(UTwork,-20.0);
877:   VecExp(UTwork);

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

881:   VecDestroy(&XX);
882:   VecDestroy(&YY);
883:   VecDestroy(&ZZ);
884:   VecDestroy(&XXwork);
885:   VecDestroy(&YYwork);
886:   VecDestroy(&ZZwork);
887:   VecDestroy(&UTwork);

889:   /* Initial guess and reference model */
890:   VecDuplicate(user->utrue,&user->ur);
891:   VecSum(user->utrue,&meanut);
892:   meanut = meanut / n;
893:   VecSet(user->ur,meanut);
894:   VecCopy(user->ur,user->u);

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

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

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

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

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

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

964:   MatCreate(PETSC_COMM_WORLD,&user->L);
965:   MatSetSizes(user->L,PETSC_DECIDE,ysubnlocal,m+n,n);
966:   MatSetFromOptions(user->L);
967:   MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
968:   MatSeqAIJSetPreallocation(user->L,2,NULL);
969:   MatGetOwnershipRange(user->L,&istart,&iend);

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

1001:   /* Generate Div matrix */
1002:   if (!user->use_ptap) {
1003:     /* Generate Div matrix */
1004:     MatCreate(PETSC_COMM_WORLD,&user->Div);
1005:     MatSetSizes(user->Div,ysubnlocal,PETSC_DECIDE,n,m);
1006:     MatSetFromOptions(user->Div);
1007:     MatMPIAIJSetPreallocation(user->Div,4,NULL,4,NULL);
1008:     MatSeqAIJSetPreallocation(user->Div,6,NULL);
1009:     MatGetOwnershipRange(user->Grad,&istart,&iend);

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

1034:     MatAssemblyBegin(user->Div,MAT_FINAL_ASSEMBLY);
1035:     MatAssemblyEnd(user->Div,MAT_FINAL_ASSEMBLY);
1036:     MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
1037:   } else {
1038:     MatCreate(PETSC_COMM_WORLD,&user->Diag);
1039:     MatSetSizes(user->Diag,PETSC_DECIDE,PETSC_DECIDE,m,m);
1040:     MatSetFromOptions(user->Diag);
1041:     MatMPIAIJSetPreallocation(user->Diag,1,NULL,0,NULL);
1042:     MatSeqAIJSetPreallocation(user->Diag,1,NULL);
1043:   }

1045:   /* Build work vectors and matrices */
1046:   VecCreate(PETSC_COMM_WORLD,&user->S);
1047:   VecSetSizes(user->S, PETSC_DECIDE, m);
1048:   VecSetFromOptions(user->S);

1050:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1051:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
1052:   VecSetFromOptions(user->lwork);

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

1056:   VecDuplicate(user->S,&user->Swork);
1057:   VecDuplicate(user->S,&user->Sdiag);
1058:   VecDuplicate(user->S,&user->Av_u);
1059:   VecDuplicate(user->S,&user->Twork);
1060:   VecDuplicate(user->y,&user->ywork);
1061:   VecDuplicate(user->u,&user->Ywork);
1062:   VecDuplicate(user->u,&user->uwork);
1063:   VecDuplicate(user->u,&user->js_diag);
1064:   VecDuplicate(user->c,&user->cwork);
1065:   VecDuplicate(user->d,&user->dwork);

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

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

1075:   /* First compute Av_u = Av*exp(-u) */
1076:   VecSet(user->uwork, 0);
1077:   VecAXPY(user->uwork,-1.0,user->utrue); /* Note: user->utrue */
1078:   VecExp(user->uwork);
1079:   MatMult(user->Av,user->uwork,user->Av_u);

1081:   /* Next form DSG = Div*S*Grad */
1082:   VecCopy(user->Av_u,user->Swork);
1083:   VecReciprocal(user->Swork);
1084:   if (user->use_ptap) {
1085:     MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);
1086:     MatPtAP(user->Diag,user->Grad,MAT_INITIAL_MATRIX,1.0,&user->DSG);
1087:   } else {
1088:     MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1089:     MatDiagonalScale(user->Divwork,NULL,user->Swork);
1090:     MatMatMultSymbolic(user->Divwork,user->Grad,1.0,&user->DSG);
1091:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1092:   }

1094:   MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE);
1095:   MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

1097:   if (user->use_lrc == PETSC_TRUE) {
1098:     v=PetscSqrtReal(1.0 /user->ndesign);
1099:     PetscMalloc1(user->ndesign,&user->ones);

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

1123:   MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal*user->ns,user->nstate,user->nstate,user,&user->JsInv);
1124:   MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateInvMatMult);
1125:   MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateInvMatMult);
1126:   MatSetOption(user->JsInv,MAT_SYMMETRIC,PETSC_TRUE);
1127:   MatSetOption(user->JsInv,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

1129:   MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE);
1130:   MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
1131:   /* Now solve for ytrue */
1132:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1133:   KSPSetFromOptions(user->solver);

1135:   KSPSetOperators(user->solver,user->JsBlock,user->DSG);
1136:   user->lcl->solve_type = LCL_FORWARD1;
1137:   MatMult(user->JsInv,user->q,user->ytrue);
1138:   /* First compute Av_u = Av*exp(-u) */
1139:   VecSet(user->uwork,0);
1140:   VecAXPY(user->uwork,-1.0,user->u); /* Note: user->u */
1141:   VecExp(user->uwork);
1142:   MatMult(user->Av,user->uwork,user->Av_u);

1144:   /* Next update DSG = Div*S*Grad  with user->u */
1145:   VecCopy(user->Av_u,user->Swork);
1146:   VecReciprocal(user->Swork);
1147:   if (user->use_ptap) {
1148:     MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);
1149:     MatPtAP(user->Diag,user->Grad,MAT_REUSE_MATRIX,1.0,&user->DSG);
1150:   } else {
1151:     MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1152:     MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1153:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1154:   }

1156:   /* Now solve for y */
1157:   user->lcl->solve_type = LCL_FORWARD1;
1158:   MatMult(user->JsInv,user->q,user->y);

1160:   user->ksp_its_initial = user->ksp_its;
1161:   user->ksp_its = 0;
1162:   /* Construct projection matrix Q (blocks) */
1163:   MatCreate(PETSC_COMM_WORLD,&user->Q);
1164:   MatSetSizes(user->Q,dsubnlocal,ysubnlocal,user->ndata,user->ndesign);
1165:   MatSetFromOptions(user->Q);
1166:   MatMPIAIJSetPreallocation(user->Q,8,NULL,8,NULL);
1167:   MatSeqAIJSetPreallocation(user->Q,8,NULL);

1169:   for (i=0; i<user->mx; i++){
1170:     x[i] = h*(i+0.5);
1171:     y[i] = h*(i+0.5);
1172:     z[i] = h*(i+0.5);
1173:   }
1174:   MatGetOwnershipRange(user->Q,&istart,&iend);

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

1179:     xri = xr[i];
1180:     im = 0;
1181:     xim = x[im];
1182:     while (xri>xim && im<nx){
1183:       im = im+1;
1184:       xim = x[im];
1185:     }
1186:     indx1 = im-1;
1187:     indx2 = im;
1188:     dx1 = xri - x[indx1];
1189:     dx2 = x[indx2] - xri;

1191:     yri = yr[i];
1192:     im = 0;
1193:     yim = y[im];
1194:     while (yri>yim && im<ny){
1195:       im = im+1;
1196:       yim = y[im];
1197:     }
1198:     indy1 = im-1;
1199:     indy2 = im;
1200:     dy1 = yri - y[indy1];
1201:     dy2 = y[indy2] - yri;

1203:     zri = zr[i];
1204:     im = 0;
1205:     zim = z[im];
1206:     while (zri>zim && im<nz){
1207:       im = im+1;
1208:       zim = z[im];
1209:     }
1210:     indz1 = im-1;
1211:     indz2 = im;
1212:     dz1 = zri - z[indz1];
1213:     dz2 = z[indz2] - zri;

1215:     Dx = x[indx2] - x[indx1];
1216:     Dy = y[indy2] - y[indy1];
1217:     Dz = z[indz2] - z[indz1];

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

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

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

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

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

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

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

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

1252:   MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY);
1253:   MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY);
1254:   /* Create MQ (composed of blocks of Q */
1255:   MatCreateShell(PETSC_COMM_WORLD,dsubnlocal*user->ns,PETSC_DECIDE,user->ndata*user->ns,user->nstate,user,&user->MQ);
1256:   MatShellSetOperation(user->MQ,MATOP_MULT,(void(*)(void))QMatMult);
1257:   MatShellSetOperation(user->MQ,MATOP_MULT_TRANSPOSE,(void(*)(void))QMatMultTranspose);

1259:   /* Add noise to the measurement data */
1260:   VecSet(user->ywork,1.0);
1261:   VecAYPX(user->ywork,user->noise,user->ytrue);
1262:   MatMult(user->MQ,user->ywork,user->d);

1264:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1265:   PetscFree(x);
1266:   PetscFree(y);
1267:   PetscFree(z);
1268:   PetscLogStagePop();
1269:   return(0);
1270: }

1274: PetscErrorCode EllipticDestroy(AppCtx *user)
1275: {
1277:   PetscInt       i;

1280:   MatDestroy(&user->DSG);
1281:   KSPDestroy(&user->solver);
1282:   MatDestroy(&user->Q);
1283:   MatDestroy(&user->MQ);
1284:   if (!user->use_ptap) {
1285:     MatDestroy(&user->Div);
1286:     MatDestroy(&user->Divwork);
1287:   } else {
1288:     MatDestroy(&user->Diag);
1289:   }
1290:   if (user->use_lrc) {
1291:     MatDestroy(&user->Ones);
1292:   }

1294:   MatDestroy(&user->Grad);
1295:   MatDestroy(&user->Av);
1296:   MatDestroy(&user->Avwork);
1297:   MatDestroy(&user->L);
1298:   MatDestroy(&user->Js);
1299:   MatDestroy(&user->Jd);
1300:   MatDestroy(&user->JsBlock);
1301:   MatDestroy(&user->JsInv);

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

1346: PetscErrorCode EllipticMonitor(Tao tao, void *ptr)
1347: {
1349:   Vec            X;
1350:   PetscReal      unorm,ynorm;
1351:   AppCtx         *user = (AppCtx*)ptr;

1354:   TaoGetSolutionVector(tao,&X);
1355:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1356:   VecAXPY(user->ywork,-1.0,user->ytrue);
1357:   VecAXPY(user->uwork,-1.0,user->utrue);
1358:   VecNorm(user->uwork,NORM_2,&unorm);
1359:   VecNorm(user->ywork,NORM_2,&ynorm);
1360:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1361:   return(0);
1362: }