Actual source code: parabolic.c

petsc-master 2014-12-21
Report Typos and Errors
  1: #include <petsc-private/taoimpl.h>
  2: #include "../src/tao/pde_constrained/impls/lcl/lcl.h"

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

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

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

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

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

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

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

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

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

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

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

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

 78:   KSP solver;
 79:   PC prec;

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


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

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

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

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

115: static  char help[]="";

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

506:   return(0);
507: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

634:   }

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

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

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

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

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

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

670: }

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

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


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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1305: PetscErrorCode ParabolicDestroy(AppCtx *user)
1306: {
1308:   PetscInt       i;

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

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

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