Actual source code: parabolic.c

petsc-3.7.3 2016-08-01
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: TaoDestroy();
 18:    Processors: n
 19: T*/

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

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

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

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

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

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

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

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

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

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

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

 77:   KSP solver;
 78:   PC prec;

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


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

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

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

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

113: static  char help[]="";

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

234:   TaoDestroy(&tao);
235:   VecDestroy(&x);
236:   VecDestroy(&x0);
237:   ParabolicDestroy(&user);

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

258:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
259:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
260:   for (j=0; j<user->ns; j++){
261:     i = user->sample_times[j];
262:     MatMult(user->Qblock,user->yi[i],user->di[j]);
263:   }
264:   Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
265:   VecAXPY(user->dwork,-1.0,user->d);
266:   VecDot(user->dwork,user->dwork,&d1);

268:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
269:   MatMult(user->L,user->uwork,user->lwork);
270:   VecDot(user->lwork,user->lwork,&d2);

272:   *f = 0.5 * (d1 + user->alpha*d2);
273:   return(0);
274: }

276: /* ------------------------------------------------------------------- */
279: /*
280:     state: g_s = Q' *(Qy - d)
281:     design: g_d = alpha*L'*L*(u-ur)
282: */
283: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
284: {
286:   PetscInt       i,j;
287:   AppCtx         *user = (AppCtx*)ptr;

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

307:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
308:   MatMult(user->L,user->uwork,user->lwork);
309:   MatMult(user->LT,user->lwork,user->uwork);
310:   VecScale(user->uwork, user->alpha);
311:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
312:   return(0);
313: }

317: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
318: {
320:   PetscReal      d1,d2;
321:   PetscInt       i,j;
322:   AppCtx         *user = (AppCtx*)ptr;

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

343:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
344:   MatMult(user->L,user->uwork,user->lwork);
345:   VecDot(user->lwork,user->lwork,&d2);
346:   MatMult(user->LT,user->lwork,user->uwork);
347:   VecScale(user->uwork, user->alpha);
348:   *f = 0.5 * (d1 + user->alpha*d2);

350:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
351:   return(0);
352: }

354: /* ------------------------------------------------------------------- */
357: /* A
358: MatShell object
359: */
360: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
361: {
363:   AppCtx         *user = (AppCtx*)ptr;

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

383:   /* B = speye(nx^3) + ht*DSG; */
384:   MatScale(user->DSG,user->ht);
385:   MatShift(user->DSG,1.0);
386:   return(0);
387: }

389: /* ------------------------------------------------------------------- */
392: /* B */
393: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
394: {
396:   AppCtx         *user = (AppCtx*)ptr;

399:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
400:   return(0);
401: }

405: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
406: {
408:   PetscInt       i;
409:   AppCtx         *user;

412:   MatShellGetContext(J_shell,(void**)&user);
413:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
414:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
415:   for (i=1; i<user->nt; i++){
416:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
417:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
418:   }
419:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
420:   return(0);
421: }

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

432:   MatShellGetContext(J_shell,(void**)&user);
433:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
434:   for (i=0; i<user->nt-1; i++){
435:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
436:     VecAXPY(user->yiwork[i],-1.0,user->yi[i+1]);
437:   }
438:   i = user->nt-1;
439:   MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
440:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
441:   return(0);
442: }

446: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
447: {
449:   AppCtx         *user;

452:   MatShellGetContext(J_shell,(void**)&user);
453:   MatMult(user->Grad,X,user->Swork);
454:   VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
455:   MatMult(user->Div,user->Swork,Y);
456:   VecAYPX(Y,user->ht,X);
457:   return(0);
458: }

462: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
463: {
465:   PetscInt       i;
466:   AppCtx         *user;

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

471:   /* sdiag(1./v) */
472:   VecSet(user->uwork,0);
473:   VecAXPY(user->uwork,-1.0,user->u);
474:   VecExp(user->uwork);

476:   /* sdiag(1./((Av*(1./v)).^2)) */
477:   MatMult(user->Av,user->uwork,user->Swork);
478:   VecPointwiseMult(user->Swork,user->Swork,user->Swork);
479:   VecReciprocal(user->Swork);

481:   /* (Av * (sdiag(1./v) * b)) */
482:   VecPointwiseMult(user->uwork,user->uwork,X);
483:   MatMult(user->Av,user->uwork,user->Twork);

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

488:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
489:   for (i=0; i<user->nt; i++){
490:     /* (sdiag(Grad*y(:,i)) */
491:     MatMult(user->Grad,user->yi[i],user->Twork);

493:     /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
494:     VecPointwiseMult(user->Twork,user->Twork,user->Swork);
495:     MatMult(user->Div,user->Twork,user->yiwork[i]);
496:     VecScale(user->yiwork[i],user->ht);
497:   }
498:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);

500:   return(0);
501: }

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

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

514:   /* sdiag(1./((Av*(1./v)).^2)) */
515:   VecSet(user->uwork,0);
516:   VecAXPY(user->uwork,-1.0,user->u);
517:   VecExp(user->uwork);
518:   MatMult(user->Av,user->uwork,user->Rwork);
519:   VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork);
520:   VecReciprocal(user->Rwork);

522:   VecSet(Y,0.0);
523:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
524:   Scatter_i(X,user->yiwork,user->yi_scatter,user->nt);
525:   for (i=0; i<user->nt; i++){
526:     /* (Div' * b(:,i)) */
527:     MatMult(user->Grad,user->yiwork[i],user->Swork);

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

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

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

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

541:     /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
542:     VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i]);
543:     VecAXPY(Y,user->ht,user->yiwork[i]);
544:   }
545:   return(0);
546: }

550: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
551: {
553:   AppCtx         *user;

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

558:   if (user->dsg_formed) {
559:     MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
560:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed");
561:   return(0);
562: }

566: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
567: {
569:   AppCtx         *user;
570:   PetscInt       its,i;

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

575:   if (Y == user->ytrue) {
576:     /* First solve is done with true solution to set up problem */
577:     KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
578:   } else {
579:     KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
580:   }

582:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
583:   KSPSolve(user->solver,user->yi[0],user->yiwork[0]);
584:   KSPGetIterationNumber(user->solver,&its);
585:   user->ksp_its = user->ksp_its + its;

587:   for (i=1; i<user->nt; i++){
588:     VecAXPY(user->yi[i],1.0,user->yiwork[i-1]);
589:     KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
590:     KSPGetIterationNumber(user->solver,&its);
591:     user->ksp_its = user->ksp_its + its;
592:   }
593:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
594:   return(0);
595: }

599: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
600: {
602:   AppCtx         *user;
603:   PetscInt       its,i;

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

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

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

613:   KSPGetIterationNumber(user->solver,&its);
614:   user->ksp_its = user->ksp_its + its;

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

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

623:   }

625:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
626:   return(0);
627: }

631: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
632: {
634:   AppCtx         *user;

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

639:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
640:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
641:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
642:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
643:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
644:   return(0);
645: }

649: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
650: {
652:   AppCtx         *user;

655:   MatShellGetContext(J_shell,(void**)&user);
656:    VecCopy(user->js_diag,X);
657:   return(0);

659: }

663: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
664: {
665:   /* con = Ay - q, A = [B  0  0 ... 0;
666:                        -I  B  0 ... 0;
667:                         0 -I  B ... 0;
668:                              ...     ;
669:                         0    ... -I B]
670:      B = ht * Div * Sigma * Grad + eye */
672:   PetscInt       i;
673:   AppCtx         *user = (AppCtx*)ptr;

676:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
677:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
678:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
679:   for (i=1; i<user->nt; i++){
680:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
681:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
682:   }
683:   Gather_i(C,user->yiwork,user->yi_scatter,user->nt);
684:   VecAXPY(C,-1.0,user->q);
685:   return(0);
686: }


691: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
692: {

696:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
697:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
698:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
699:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
700:   return(0);
701: }

705: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
706: {
708:   PetscInt       i;

711:   for (i=0; i<nt; i++){
712:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
713:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
714:   }
715:   return(0);
716: }


721: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
722: {

726:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
727:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
728:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
729:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
730:   return(0);
731: }

735: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
736: {
738:   PetscInt       i;

741:   for (i=0; i<nt; i++){
742:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
743:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
744:   }
745:   return(0);
746: }

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

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

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

791:   PetscMalloc1(user->mx,&x);
792:   PetscMalloc1(user->mx,&y);
793:   PetscMalloc1(user->mx,&z);
794:   user->jformed = PETSC_FALSE;
795:   user->dsg_formed = PETSC_FALSE;

797:   n = user->mx * user->mx * user->mx;
798:   m = 3 * user->mx * user->mx * (user->mx-1);
799:   sqrt_beta = PetscSqrtScalar(user->beta);

801:   user->ksp_its = 0;
802:   user->ksp_its_initial = 0;

804:   stime = (PetscReal)user->nt/user->ns;
805:   PetscMalloc1(user->ns,&user->sample_times);
806:   for (i=0; i<user->ns; i++){
807:     user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
808:   }

810:   VecCreate(PETSC_COMM_WORLD,&XX);
811:   VecCreate(PETSC_COMM_WORLD,&user->q);
812:   VecSetSizes(XX,PETSC_DECIDE,n);
813:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
814:   VecSetFromOptions(XX);
815:   VecSetFromOptions(user->q);

817:   VecDuplicate(XX,&YY);
818:   VecDuplicate(XX,&ZZ);
819:   VecDuplicate(XX,&XXwork);
820:   VecDuplicate(XX,&YYwork);
821:   VecDuplicate(XX,&ZZwork);
822:   VecDuplicate(XX,&UTwork);
823:   VecDuplicate(XX,&user->utrue);
824:   VecDuplicate(XX,&bc);

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

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

848:   VecAssemblyBegin(XX);
849:   VecAssemblyEnd(XX);
850:   VecAssemblyBegin(YY);
851:   VecAssemblyEnd(YY);
852:   VecAssemblyBegin(ZZ);
853:   VecAssemblyEnd(ZZ);
854:   VecAssemblyBegin(bc);
855:   VecAssemblyEnd(bc);

857:   /* Compute true parameter function
858:      ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
859:   VecCopy(XX,XXwork);
860:   VecCopy(YY,YYwork);
861:   VecCopy(ZZ,ZZwork);

863:   VecShift(XXwork,-0.5);
864:   VecShift(YYwork,-0.5);
865:   VecShift(ZZwork,-0.5);

867:   VecPointwiseMult(XXwork,XXwork,XXwork);
868:   VecPointwiseMult(YYwork,YYwork,YYwork);
869:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

871:   VecCopy(XXwork,user->utrue);
872:   VecAXPY(user->utrue,1.0,YYwork);
873:   VecAXPY(user->utrue,1.0,ZZwork);
874:   VecScale(user->utrue,-10.0);
875:   VecExp(user->utrue);
876:   VecShift(user->utrue,0.5);

878:   VecDestroy(&XX);
879:   VecDestroy(&YY);
880:   VecDestroy(&ZZ);
881:   VecDestroy(&XXwork);
882:   VecDestroy(&YYwork);
883:   VecDestroy(&ZZwork);
884:   VecDestroy(&UTwork);

886:   /* Initial guess and reference model */
887:   VecDuplicate(user->utrue,&user->ur);
888:   VecSet(user->ur,0.0);

890:   /* Generate Grad matrix */
891:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
892:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,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);


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

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

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

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

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

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

996:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
997:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);

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

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

1004:   /* Build work vectors and matrices */
1005:   VecCreate(PETSC_COMM_WORLD,&user->S);
1006:   VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
1007:   VecSetFromOptions(user->S);

1009:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1010:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
1011:   VecSetFromOptions(user->lwork);

1013:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
1014:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

1016:   VecCreate(PETSC_COMM_WORLD,&user->d);
1017:   VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
1018:   VecSetFromOptions(user->d);

1020:   VecDuplicate(user->S,&user->Swork);
1021:   VecDuplicate(user->S,&user->Av_u);
1022:   VecDuplicate(user->S,&user->Twork);
1023:   VecDuplicate(user->S,&user->Rwork);
1024:   VecDuplicate(user->y,&user->ywork);
1025:   VecDuplicate(user->u,&user->uwork);
1026:   VecDuplicate(user->u,&user->js_diag);
1027:   VecDuplicate(user->c,&user->cwork);
1028:   VecDuplicate(user->d,&user->dwork);

1030:   /* Create matrix-free shell user->Js for computing A*x */
1031:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
1032:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1033:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
1034:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
1035:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

1037:   /* Diagonal blocks of user->Js */
1038:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
1039:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
1040:   /* Blocks are symmetric */
1041:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);

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

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

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

1062:   /* Solver options and tolerances */
1063:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1064:   KSPSetType(user->solver,KSPCG);
1065:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
1066:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
1067:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1068:   KSPSetFromOptions(user->solver);
1069:   KSPGetPC(user->solver,&user->prec);
1070:   PCSetType(user->prec,PCSHELL);

1072:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1073:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
1074:   PCShellSetContext(user->prec,user);

1076:   /* Create scatter from y to y_1,y_2,...,y_nt */
1077:   PetscMalloc1(user->nt*user->m,&user->yi_scatter);
1078:   VecCreate(PETSC_COMM_WORLD,&yi);
1079:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
1080:   VecSetFromOptions(yi);
1081:   VecDuplicateVecs(yi,user->nt,&user->yi);
1082:   VecDuplicateVecs(yi,user->nt,&user->yiwork);

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

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

1116:   /* Assemble RHS of forward problem */
1117:   VecCopy(bc,user->yiwork[0]);
1118:   for (i=1; i<user->nt; i++){
1119:     VecSet(user->yiwork[i],0.0);
1120:   }
1121:   Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1122:   VecDestroy(&bc);

1124:   /* Compute true state function yt given ut */
1125:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1126:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1127:   VecSetFromOptions(user->ytrue);

1129:   /* First compute Av_u = Av*exp(-u) */
1130:   VecSet(user->uwork,0);
1131:   VecAXPY(user->uwork,-1.0,user->utrue); /*  Note: user->utrue */
1132:   VecExp(user->uwork);
1133:   MatMult(user->Av,user->uwork,user->Av_u);

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

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

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

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

1157:   /* First compute Av_u = Av*exp(-u) */
1158:   VecSet(user->uwork,0);
1159:   VecAXPY(user->uwork,-1.0,user->u); /*  Note: user->u */
1160:   VecExp(user->uwork);
1161:   MatMult(user->Av,user->uwork,user->Av_u);

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

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

1180:   /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
1181:   MatCreate(PETSC_COMM_WORLD,&user->Qblock);
1182:   MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n);
1183:   MatSetFromOptions(user->Qblock);
1184:   MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL);
1185:   MatSeqAIJSetPreallocation(user->Qblock,8,NULL);

1187:   for (i=0; i<user->mx; i++){
1188:     x[i] = h*(i+0.5);
1189:     y[i] = h*(i+0.5);
1190:     z[i] = h*(i+0.5);
1191:   }

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

1208:     yri = yr[i];
1209:     im = 0;
1210:     yim = y[im];
1211:     while (yri>yim && im<ny){
1212:       im = im+1;
1213:       yim = y[im];
1214:     }
1215:     indy1 = im-1;
1216:     indy2 = im;
1217:     dy1 = yri - y[indy1];
1218:     dy2 = y[indy2] - yri;

1220:     zri = zr[i];
1221:     im = 0;
1222:     zim = z[im];
1223:     while (zri>zim && im<nz){
1224:       im = im+1;
1225:       zim = z[im];
1226:     }
1227:     indz1 = im-1;
1228:     indz2 = im;
1229:     dz1 = zri - z[indz1];
1230:     dz2 = z[indz2] - zri;

1232:     Dx = x[indx2] - x[indx1];
1233:     Dy = y[indy2] - y[indy1];
1234:     Dz = z[indz2] - z[indz1];

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

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

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

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

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

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

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

1264:     j = indx2 + indy2*nx + indz2*nx*ny;
1265:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1266:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1267:   }
1268:   MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY);
1269:   MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY);

1271:   MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1272:   MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);

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

1284:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1285:   KSPSetFromOptions(user->solver);
1286:   PetscFree(x);
1287:   PetscFree(y);
1288:   PetscFree(z);
1289:   return(0);
1290: }

1294: PetscErrorCode ParabolicDestroy(AppCtx *user)
1295: {
1297:   PetscInt       i;

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

1360: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1361: {
1363:   Vec            X;
1364:   PetscReal      unorm,ynorm;
1365:   AppCtx         *user = (AppCtx*)ptr;

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