Actual source code: parabolic.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: #include <petsc/private/taoimpl.h>

  3: /*T
  4:    Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
  5:    Routines: TaoCreate();
  6:    Routines: TaoSetType();
  7:    Routines: TaoSetInitialVector();
  8:    Routines: TaoSetObjectiveRoutine();
  9:    Routines: TaoSetGradientRoutine();
 10:    Routines: TaoSetConstraintsRoutine();
 11:    Routines: TaoSetJacobianStateRoutine();
 12:    Routines: TaoSetJacobianDesignRoutine();
 13:    Routines: TaoSetStateDesignIS();
 14:    Routines: TaoSetFromOptions();
 15:    Routines: TaoSetHistory(); TaoGetHistory();
 16:    Routines: TaoSolve();
 17:    Routines: 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[]="";

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

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

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

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

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

170:   VecGetOwnershipRange(user.y,&lo,&hi);
171:   VecGetOwnershipRange(user.u,&lo2,&hi2);

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

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

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

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

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

191:   /* Set up initial vectors and matrices */
192:   ParabolicInitialize(&user);

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

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

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

207:   TaoSetFromOptions(tao);
208:   TaoSetStateDesignIS(tao,user.s_is,user.d_is);

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

232:   TaoDestroy(&tao);
233:   VecDestroy(&x);
234:   VecDestroy(&x0);
235:   ParabolicDestroy(&user);

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

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

264:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
265:   MatMult(user->L,user->uwork,user->lwork);
266:   VecDot(user->lwork,user->lwork,&d2);

268:   *f = 0.5 * (d1 + user->alpha*d2);
269:   return(0);
270: }

272: /* ------------------------------------------------------------------- */
273: /*
274:     state: g_s = Q' *(Qy - d)
275:     design: g_d = alpha*L'*L*(u-ur)
276: */
277: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
278: {
280:   PetscInt       i,j;
281:   AppCtx         *user = (AppCtx*)ptr;

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

301:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
302:   MatMult(user->L,user->uwork,user->lwork);
303:   MatMult(user->LT,user->lwork,user->uwork);
304:   VecScale(user->uwork, user->alpha);
305:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
306:   return(0);
307: }

309: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
310: {
312:   PetscReal      d1,d2;
313:   PetscInt       i,j;
314:   AppCtx         *user = (AppCtx*)ptr;

317:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
318:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
319:   for (j=0; j<user->ns; j++){
320:     i = user->sample_times[j];
321:     MatMult(user->Qblock,user->yi[i],user->di[j]);
322:   }
323:   Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
324:   VecAXPY(user->dwork,-1.0,user->d);
325:   VecDot(user->dwork,user->dwork,&d1);
326:   Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
327:   VecSet(user->ywork,0.0);
328:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
329:   for (j=0; j<user->ns; j++){
330:     i = user->sample_times[j];
331:     MatMult(user->QblockT,user->di[j],user->yiwork[i]);
332:   }
333:   Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);

335:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
336:   MatMult(user->L,user->uwork,user->lwork);
337:   VecDot(user->lwork,user->lwork,&d2);
338:   MatMult(user->LT,user->lwork,user->uwork);
339:   VecScale(user->uwork, user->alpha);
340:   *f = 0.5 * (d1 + user->alpha*d2);

342:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
343:   return(0);
344: }

346: /* ------------------------------------------------------------------- */
347: /* A
348: MatShell object
349: */
350: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
351: {
353:   AppCtx         *user = (AppCtx*)ptr;

356:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
357:   VecSet(user->uwork,0);
358:   VecAXPY(user->uwork,-1.0,user->u);
359:   VecExp(user->uwork);
360:   MatMult(user->Av,user->uwork,user->Av_u);
361:   VecCopy(user->Av_u,user->Swork);
362:   VecReciprocal(user->Swork);
363:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
364:   MatDiagonalScale(user->Divwork,NULL,user->Swork);
365:   if (user->dsg_formed) {
366:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
367:   } else {
368:     MatMatMultSymbolic(user->Divwork,user->Grad,PETSC_DECIDE,&user->DSG);
369:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
370:     user->dsg_formed = PETSC_TRUE;
371:   }

373:   /* B = speye(nx^3) + ht*DSG; */
374:   MatScale(user->DSG,user->ht);
375:   MatShift(user->DSG,1.0);
376:   return(0);
377: }

379: /* ------------------------------------------------------------------- */
380: /* B */
381: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
382: {
384:   AppCtx         *user = (AppCtx*)ptr;

387:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
388:   return(0);
389: }

391: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
392: {
394:   PetscInt       i;
395:   AppCtx         *user;

398:   MatShellGetContext(J_shell,(void**)&user);
399:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
400:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
401:   for (i=1; i<user->nt; i++){
402:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
403:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
404:   }
405:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
406:   return(0);
407: }

409: PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
410: {
412:   PetscInt       i;
413:   AppCtx         *user;

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

428: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
429: {
431:   AppCtx         *user;

434:   MatShellGetContext(J_shell,(void**)&user);
435:   MatMult(user->Grad,X,user->Swork);
436:   VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
437:   MatMult(user->Div,user->Swork,Y);
438:   VecAYPX(Y,user->ht,X);
439:   return(0);
440: }

442: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
443: {
445:   PetscInt       i;
446:   AppCtx         *user;

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

451:   /* sdiag(1./v) */
452:   VecSet(user->uwork,0);
453:   VecAXPY(user->uwork,-1.0,user->u);
454:   VecExp(user->uwork);

456:   /* sdiag(1./((Av*(1./v)).^2)) */
457:   MatMult(user->Av,user->uwork,user->Swork);
458:   VecPointwiseMult(user->Swork,user->Swork,user->Swork);
459:   VecReciprocal(user->Swork);

461:   /* (Av * (sdiag(1./v) * b)) */
462:   VecPointwiseMult(user->uwork,user->uwork,X);
463:   MatMult(user->Av,user->uwork,user->Twork);

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

468:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
469:   for (i=0; i<user->nt; i++){
470:     /* (sdiag(Grad*y(:,i)) */
471:     MatMult(user->Grad,user->yi[i],user->Twork);

473:     /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
474:     VecPointwiseMult(user->Twork,user->Twork,user->Swork);
475:     MatMult(user->Div,user->Twork,user->yiwork[i]);
476:     VecScale(user->yiwork[i],user->ht);
477:   }
478:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);

480:   return(0);
481: }

483: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
484: {
486:   PetscInt       i;
487:   AppCtx         *user;

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

492:   /* sdiag(1./((Av*(1./v)).^2)) */
493:   VecSet(user->uwork,0);
494:   VecAXPY(user->uwork,-1.0,user->u);
495:   VecExp(user->uwork);
496:   MatMult(user->Av,user->uwork,user->Rwork);
497:   VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork);
498:   VecReciprocal(user->Rwork);

500:   VecSet(Y,0.0);
501:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
502:   Scatter_i(X,user->yiwork,user->yi_scatter,user->nt);
503:   for (i=0; i<user->nt; i++){
504:     /* (Div' * b(:,i)) */
505:     MatMult(user->Grad,user->yiwork[i],user->Swork);

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

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

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

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

519:     /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
520:     VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i]);
521:     VecAXPY(Y,user->ht,user->yiwork[i]);
522:   }
523:   return(0);
524: }

526: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
527: {
529:   AppCtx         *user;

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

534:   if (user->dsg_formed) {
535:     MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
536:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed");
537:   return(0);
538: }

540: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
541: {
543:   AppCtx         *user;
544:   PetscInt       its,i;

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

549:   if (Y == user->ytrue) {
550:     /* First solve is done with true solution to set up problem */
551:     KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
552:   } else {
553:     KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
554:   }

556:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
557:   KSPSolve(user->solver,user->yi[0],user->yiwork[0]);
558:   KSPGetIterationNumber(user->solver,&its);
559:   user->ksp_its = user->ksp_its + its;

561:   for (i=1; i<user->nt; i++){
562:     VecAXPY(user->yi[i],1.0,user->yiwork[i-1]);
563:     KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
564:     KSPGetIterationNumber(user->solver,&its);
565:     user->ksp_its = user->ksp_its + its;
566:   }
567:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
568:   return(0);
569: }

571: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
572: {
574:   AppCtx         *user;
575:   PetscInt       its,i;

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

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

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

585:   KSPGetIterationNumber(user->solver,&its);
586:   user->ksp_its = user->ksp_its + its;

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

592:     KSPGetIterationNumber(user->solver,&its);
593:     user->ksp_its = user->ksp_its + its;

595:   }

597:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
598:   return(0);
599: }

601: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
602: {
604:   AppCtx         *user;

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

609:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
610:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
611:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
612:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
613:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
614:   return(0);
615: }

617: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
618: {
620:   AppCtx         *user;

623:   MatShellGetContext(J_shell,(void**)&user);
624:    VecCopy(user->js_diag,X);
625:   return(0);

627: }

629: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
630: {
631:   /* con = Ay - q, A = [B  0  0 ... 0;
632:                        -I  B  0 ... 0;
633:                         0 -I  B ... 0;
634:                              ...     ;
635:                         0    ... -I B]
636:      B = ht * Div * Sigma * Grad + eye */
638:   PetscInt       i;
639:   AppCtx         *user = (AppCtx*)ptr;

642:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
643:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
644:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
645:   for (i=1; i<user->nt; i++){
646:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
647:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
648:   }
649:   Gather_i(C,user->yiwork,user->yi_scatter,user->nt);
650:   VecAXPY(C,-1.0,user->q);
651:   return(0);
652: }


655: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
656: {

660:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
661:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
662:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
663:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
664:   return(0);
665: }

667: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
668: {
670:   PetscInt       i;

673:   for (i=0; i<nt; i++){
674:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
675:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
676:   }
677:   return(0);
678: }


681: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
682: {

686:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
687:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
688:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
689:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
690:   return(0);
691: }

693: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
694: {
696:   PetscInt       i;

699:   for (i=0; i<nt; i++){
700:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
701:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
702:   }
703:   return(0);
704: }

706: PetscErrorCode ParabolicInitialize(AppCtx *user)
707: {
709:   PetscInt       m,n,i,j,k,linear_index,istart,iend,iblock,lo,hi,lo2,hi2;
710:   Vec            XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork,yi,di,bc;
711:   PetscReal      *x, *y, *z;
712:   PetscReal      h,stime;
713:   PetscScalar    hinv,neg_hinv,half = 0.5,sqrt_beta;
714:   PetscInt       im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
715:   PetscReal      xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
716:   PetscScalar    v,vx,vy,vz;
717:   IS             is_from_y,is_to_yi,is_from_d,is_to_di;
718:   /* Data locations */
719:   PetscScalar xr[64] = {0.4970,     0.8498,     0.7814,     0.6268,     0.7782,     0.6402,     0.3617,     0.3160,
720:                         0.3610,     0.5298,     0.6987,     0.3331,     0.7962,     0.5596,     0.3866,     0.6774,
721:                         0.5407,     0.4518,     0.6702,     0.6061,     0.7580,     0.8997,     0.5198,     0.8326,
722:                         0.2138,     0.9198,     0.3000,     0.2833,     0.8288,     0.7076,     0.1820,     0.0728,
723:                         0.8447,     0.2367,     0.3239,     0.6413,     0.3114,     0.4731,     0.1192,     0.9273,
724:                         0.5724,     0.4331,     0.5136,     0.3547,     0.4413,     0.2602,     0.5698,     0.7278,
725:                         0.5261,     0.6230,     0.2454,     0.3948,     0.7479,     0.6582,     0.4660,     0.5594,
726:                         0.7574,     0.1143,     0.5900,     0.1065,     0.4260,     0.3294,     0.8276,     0.0756};

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

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

747:   PetscMalloc1(user->mx,&x);
748:   PetscMalloc1(user->mx,&y);
749:   PetscMalloc1(user->mx,&z);
750:   user->jformed = PETSC_FALSE;
751:   user->dsg_formed = PETSC_FALSE;

753:   n = user->mx * user->mx * user->mx;
754:   m = 3 * user->mx * user->mx * (user->mx-1);
755:   sqrt_beta = PetscSqrtScalar(user->beta);

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

760:   stime = (PetscReal)user->nt/user->ns;
761:   PetscMalloc1(user->ns,&user->sample_times);
762:   for (i=0; i<user->ns; i++){
763:     user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
764:   }

766:   VecCreate(PETSC_COMM_WORLD,&XX);
767:   VecCreate(PETSC_COMM_WORLD,&user->q);
768:   VecSetSizes(XX,PETSC_DECIDE,n);
769:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
770:   VecSetFromOptions(XX);
771:   VecSetFromOptions(user->q);

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

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

787:   VecGetOwnershipRange(XX,&istart,&iend);
788:   for (linear_index=istart; linear_index<iend; linear_index++){
789:     i = linear_index % user->mx;
790:     j = ((linear_index-i)/user->mx) % user->mx;
791:     k = ((linear_index-i)/user->mx-j) / user->mx;
792:     vx = h*(i+0.5);
793:     vy = h*(j+0.5);
794:     vz = h*(k+0.5);
795:     VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
796:     VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
797:     VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES);
798:     if ((vx<0.6) && (vx>0.4) && (vy<0.6) && (vy>0.4) && (vy<0.6) && (vz<0.6) && (vz>0.4)){
799:       v = 1000.0;
800:       VecSetValues(bc,1,&linear_index,&v,INSERT_VALUES);
801:     }
802:   }

804:   VecAssemblyBegin(XX);
805:   VecAssemblyEnd(XX);
806:   VecAssemblyBegin(YY);
807:   VecAssemblyEnd(YY);
808:   VecAssemblyBegin(ZZ);
809:   VecAssemblyEnd(ZZ);
810:   VecAssemblyBegin(bc);
811:   VecAssemblyEnd(bc);

813:   /* Compute true parameter function
814:      ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
815:   VecCopy(XX,XXwork);
816:   VecCopy(YY,YYwork);
817:   VecCopy(ZZ,ZZwork);

819:   VecShift(XXwork,-0.5);
820:   VecShift(YYwork,-0.5);
821:   VecShift(ZZwork,-0.5);

823:   VecPointwiseMult(XXwork,XXwork,XXwork);
824:   VecPointwiseMult(YYwork,YYwork,YYwork);
825:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

827:   VecCopy(XXwork,user->utrue);
828:   VecAXPY(user->utrue,1.0,YYwork);
829:   VecAXPY(user->utrue,1.0,ZZwork);
830:   VecScale(user->utrue,-10.0);
831:   VecExp(user->utrue);
832:   VecShift(user->utrue,0.5);

834:   VecDestroy(&XX);
835:   VecDestroy(&YY);
836:   VecDestroy(&ZZ);
837:   VecDestroy(&XXwork);
838:   VecDestroy(&YYwork);
839:   VecDestroy(&ZZwork);
840:   VecDestroy(&UTwork);

842:   /* Initial guess and reference model */
843:   VecDuplicate(user->utrue,&user->ur);
844:   VecSet(user->ur,0.0);

846:   /* Generate Grad matrix */
847:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
848:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n);
849:   MatSetFromOptions(user->Grad);
850:   MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
851:   MatSeqAIJSetPreallocation(user->Grad,2,NULL);
852:   MatGetOwnershipRange(user->Grad,&istart,&iend);

854:   for (i=istart; i<iend; i++){
855:     if (i<m/3){
856:       iblock = i / (user->mx-1);
857:       j = iblock*user->mx + (i % (user->mx-1));
858:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
859:       j = j+1;
860:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
861:     }
862:     if (i>=m/3 && i<2*m/3){
863:       iblock = (i-m/3) / (user->mx*(user->mx-1));
864:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
865:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
866:       j = j + user->mx;
867:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
868:     }
869:     if (i>=2*m/3){
870:       j = i-2*m/3;
871:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
872:       j = j + user->mx*user->mx;
873:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
874:     }
875:   }

877:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
878:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);


881:   /* Generate arithmetic averaging matrix Av */
882:   MatCreate(PETSC_COMM_WORLD,&user->Av);
883:   MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n);
884:   MatSetFromOptions(user->Av);
885:   MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
886:   MatSeqAIJSetPreallocation(user->Av,2,NULL);
887:   MatGetOwnershipRange(user->Av,&istart,&iend);

889:   for (i=istart; i<iend; i++){
890:     if (i<m/3){
891:       iblock = i / (user->mx-1);
892:       j = iblock*user->mx + (i % (user->mx-1));
893:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
894:       j = j+1;
895:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
896:     }
897:     if (i>=m/3 && i<2*m/3){
898:       iblock = (i-m/3) / (user->mx*(user->mx-1));
899:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
900:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
901:       j = j + user->mx;
902:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
903:     }
904:     if (i>=2*m/3){
905:       j = i-2*m/3;
906:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
907:       j = j + user->mx*user->mx;
908:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
909:     }
910:   }

912:   MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
913:   MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);

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

918:   MatCreate(PETSC_COMM_WORLD,&user->L);
919:   MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n);
920:   MatSetFromOptions(user->L);
921:   MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
922:   MatSeqAIJSetPreallocation(user->L,2,NULL);
923:   MatGetOwnershipRange(user->L,&istart,&iend);

925:   for (i=istart; i<iend; i++){
926:     if (i<m/3){
927:       iblock = i / (user->mx-1);
928:       j = iblock*user->mx + (i % (user->mx-1));
929:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
930:       j = j+1;
931:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
932:     }
933:     if (i>=m/3 && i<2*m/3){
934:       iblock = (i-m/3) / (user->mx*(user->mx-1));
935:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
936:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
937:       j = j + user->mx;
938:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
939:     }
940:     if (i>=2*m/3 && i<m){
941:       j = i-2*m/3;
942:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
943:       j = j + user->mx*user->mx;
944:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
945:     }
946:     if (i>=m){
947:       j = i - m;
948:       MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES);
949:     }
950:   }

952:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
953:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);

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

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

960:   /* Build work vectors and matrices */
961:   VecCreate(PETSC_COMM_WORLD,&user->S);
962:   VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
963:   VecSetFromOptions(user->S);

965:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
966:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
967:   VecSetFromOptions(user->lwork);

969:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
970:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

972:   VecCreate(PETSC_COMM_WORLD,&user->d);
973:   VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
974:   VecSetFromOptions(user->d);

976:   VecDuplicate(user->S,&user->Swork);
977:   VecDuplicate(user->S,&user->Av_u);
978:   VecDuplicate(user->S,&user->Twork);
979:   VecDuplicate(user->S,&user->Rwork);
980:   VecDuplicate(user->y,&user->ywork);
981:   VecDuplicate(user->u,&user->uwork);
982:   VecDuplicate(user->u,&user->js_diag);
983:   VecDuplicate(user->c,&user->cwork);
984:   VecDuplicate(user->d,&user->dwork);

986:   /* Create matrix-free shell user->Js for computing A*x */
987:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
988:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
989:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
990:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
991:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

993:   /* Diagonal blocks of user->Js */
994:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
995:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
996:   /* Blocks are symmetric */
997:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);

999:   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1000:      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1001:      This is an SSOR preconditioner for user->JsBlock. */
1002:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlockPrec);
1003:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult);
1004:   /* JsBlockPrec is symmetric */
1005:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMult);
1006:   MatSetOption(user->JsBlockPrec,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

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

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

1018:   /* Solver options and tolerances */
1019:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1020:   KSPSetType(user->solver,KSPCG);
1021:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
1022:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
1023:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1024:   KSPSetFromOptions(user->solver);
1025:   KSPGetPC(user->solver,&user->prec);
1026:   PCSetType(user->prec,PCSHELL);

1028:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1029:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
1030:   PCShellSetContext(user->prec,user);

1032:   /* Create scatter from y to y_1,y_2,...,y_nt */
1033:   PetscMalloc1(user->nt*user->m,&user->yi_scatter);
1034:   VecCreate(PETSC_COMM_WORLD,&yi);
1035:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
1036:   VecSetFromOptions(yi);
1037:   VecDuplicateVecs(yi,user->nt,&user->yi);
1038:   VecDuplicateVecs(yi,user->nt,&user->yiwork);

1040:   VecGetOwnershipRange(user->y,&lo2,&hi2);
1041:   istart = 0;
1042:   for (i=0; i<user->nt; i++){
1043:     VecGetOwnershipRange(user->yi[i],&lo,&hi);
1044:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
1045:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
1046:     VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
1047:     istart = istart + hi-lo;
1048:     ISDestroy(&is_to_yi);
1049:     ISDestroy(&is_from_y);
1050:   }
1051:   VecDestroy(&yi);

1053:   /* Create scatter from d to d_1,d_2,...,d_ns */
1054:   PetscMalloc1(user->ns*user->ndata,&user->di_scatter);
1055:   VecCreate(PETSC_COMM_WORLD,&di);
1056:   VecSetSizes(di,PETSC_DECIDE,user->ndata);
1057:   VecSetFromOptions(di);
1058:   VecDuplicateVecs(di,user->ns,&user->di);
1059:   VecGetOwnershipRange(user->d,&lo2,&hi2);
1060:   istart = 0;
1061:   for (i=0; i<user->ns; i++){
1062:     VecGetOwnershipRange(user->di[i],&lo,&hi);
1063:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_di);
1064:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d);
1065:     VecScatterCreate(user->d,is_from_d,user->di[i],is_to_di,&user->di_scatter[i]);
1066:     istart = istart + hi-lo;
1067:     ISDestroy(&is_to_di);
1068:     ISDestroy(&is_from_d);
1069:   }
1070:   VecDestroy(&di);

1072:   /* Assemble RHS of forward problem */
1073:   VecCopy(bc,user->yiwork[0]);
1074:   for (i=1; i<user->nt; i++){
1075:     VecSet(user->yiwork[i],0.0);
1076:   }
1077:   Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1078:   VecDestroy(&bc);

1080:   /* Compute true state function yt given ut */
1081:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1082:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1083:   VecSetFromOptions(user->ytrue);

1085:   /* First compute Av_u = Av*exp(-u) */
1086:   VecSet(user->uwork,0);
1087:   VecAXPY(user->uwork,-1.0,user->utrue); /*  Note: user->utrue */
1088:   VecExp(user->uwork);
1089:   MatMult(user->Av,user->uwork,user->Av_u);

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

1094:   /* Next form DSG = Div*S*Grad */
1095:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1096:   MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1097:   if (user->dsg_formed) {
1098:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1099:   } else {
1100:     MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1101:     MatMatMultNumeric(user->Div,user->Grad,user->DSG);
1102:     user->dsg_formed = PETSC_TRUE;
1103:   }
1104:   /* B = speye(nx^3) + ht*DSG; */
1105:   MatScale(user->DSG,user->ht);
1106:   MatShift(user->DSG,1.0);

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

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

1113:   /* First compute Av_u = Av*exp(-u) */
1114:   VecSet(user->uwork,0);
1115:   VecAXPY(user->uwork,-1.0,user->u); /*  Note: user->u */
1116:   VecExp(user->uwork);
1117:   MatMult(user->Av,user->uwork,user->Av_u);

1119:   /* Next form DSG = Div*S*Grad */
1120:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1121:   MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1122:   if (user->dsg_formed) {
1123:     MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1124:   } else {
1125:     MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1126:     MatMatMultNumeric(user->Div,user->Grad,user->DSG);
1127:     user->dsg_formed = PETSC_TRUE;
1128:   }
1129:   /* B = speye(nx^3) + ht*DSG; */
1130:   MatScale(user->DSG,user->ht);
1131:   MatShift(user->DSG,1.0);

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

1136:   /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
1137:   MatCreate(PETSC_COMM_WORLD,&user->Qblock);
1138:   MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n);
1139:   MatSetFromOptions(user->Qblock);
1140:   MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL);
1141:   MatSeqAIJSetPreallocation(user->Qblock,8,NULL);

1143:   for (i=0; i<user->mx; i++){
1144:     x[i] = h*(i+0.5);
1145:     y[i] = h*(i+0.5);
1146:     z[i] = h*(i+0.5);
1147:   }

1149:   MatGetOwnershipRange(user->Qblock,&istart,&iend);
1150:   nx = user->mx; ny = user->mx; nz = user->mx;
1151:   for (i=istart; i<iend; i++){
1152:     xri = xr[i];
1153:     im = 0;
1154:     xim = x[im];
1155:     while (xri>xim && im<nx){
1156:       im = im+1;
1157:       xim = x[im];
1158:     }
1159:     indx1 = im-1;
1160:     indx2 = im;
1161:     dx1 = xri - x[indx1];
1162:     dx2 = x[indx2] - xri;

1164:     yri = yr[i];
1165:     im = 0;
1166:     yim = y[im];
1167:     while (yri>yim && im<ny){
1168:       im = im+1;
1169:       yim = y[im];
1170:     }
1171:     indy1 = im-1;
1172:     indy2 = im;
1173:     dy1 = yri - y[indy1];
1174:     dy2 = y[indy2] - yri;

1176:     zri = zr[i];
1177:     im = 0;
1178:     zim = z[im];
1179:     while (zri>zim && im<nz){
1180:       im = im+1;
1181:       zim = z[im];
1182:     }
1183:     indz1 = im-1;
1184:     indz2 = im;
1185:     dz1 = zri - z[indz1];
1186:     dz2 = z[indz2] - zri;

1188:     Dx = x[indx2] - x[indx1];
1189:     Dy = y[indy2] - y[indy1];
1190:     Dz = z[indz2] - z[indz1];

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

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

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

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

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

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

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

1220:     j = indx2 + indy2*nx + indz2*nx*ny;
1221:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1222:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1223:   }
1224:   MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY);
1225:   MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY);

1227:   MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1228:   MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);

1230:   /* Add noise to the measurement data */
1231:   VecSet(user->ywork,1.0);
1232:   VecAYPX(user->ywork,user->noise,user->ytrue);
1233:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
1234:   for (j=0; j<user->ns; j++){
1235:     i = user->sample_times[j];
1236:     MatMult(user->Qblock,user->yiwork[i],user->di[j]);
1237:   }
1238:   Gather_i(user->d,user->di,user->di_scatter,user->ns);

1240:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1241:   KSPSetFromOptions(user->solver);
1242:   PetscFree(x);
1243:   PetscFree(y);
1244:   PetscFree(z);
1245:   return(0);
1246: }

1248: PetscErrorCode ParabolicDestroy(AppCtx *user)
1249: {
1251:   PetscInt       i;

1254:   MatDestroy(&user->Qblock);
1255:   MatDestroy(&user->QblockT);
1256:   MatDestroy(&user->Div);
1257:   MatDestroy(&user->Divwork);
1258:   MatDestroy(&user->Grad);
1259:   MatDestroy(&user->Av);
1260:   MatDestroy(&user->Avwork);
1261:   MatDestroy(&user->AvT);
1262:   MatDestroy(&user->DSG);
1263:   MatDestroy(&user->L);
1264:   MatDestroy(&user->LT);
1265:   KSPDestroy(&user->solver);
1266:   MatDestroy(&user->Js);
1267:   MatDestroy(&user->Jd);
1268:   MatDestroy(&user->JsInv);
1269:   MatDestroy(&user->JsBlock);
1270:   MatDestroy(&user->JsBlockPrec);
1271:   VecDestroy(&user->u);
1272:   VecDestroy(&user->uwork);
1273:   VecDestroy(&user->utrue);
1274:   VecDestroy(&user->y);
1275:   VecDestroy(&user->ywork);
1276:   VecDestroy(&user->ytrue);
1277:   VecDestroyVecs(user->nt,&user->yi);
1278:   VecDestroyVecs(user->nt,&user->yiwork);
1279:   VecDestroyVecs(user->ns,&user->di);
1280:   PetscFree(user->yi);
1281:   PetscFree(user->yiwork);
1282:   PetscFree(user->di);
1283:   VecDestroy(&user->c);
1284:   VecDestroy(&user->cwork);
1285:   VecDestroy(&user->ur);
1286:   VecDestroy(&user->q);
1287:   VecDestroy(&user->d);
1288:   VecDestroy(&user->dwork);
1289:   VecDestroy(&user->lwork);
1290:   VecDestroy(&user->S);
1291:   VecDestroy(&user->Swork);
1292:   VecDestroy(&user->Av_u);
1293:   VecDestroy(&user->Twork);
1294:   VecDestroy(&user->Rwork);
1295:   VecDestroy(&user->js_diag);
1296:   ISDestroy(&user->s_is);
1297:   ISDestroy(&user->d_is);
1298:   VecScatterDestroy(&user->state_scatter);
1299:   VecScatterDestroy(&user->design_scatter);
1300:   for (i=0; i<user->nt; i++){
1301:     VecScatterDestroy(&user->yi_scatter[i]);
1302:   }
1303:   for (i=0; i<user->ns; i++){
1304:     VecScatterDestroy(&user->di_scatter[i]);
1305:   }
1306:   PetscFree(user->yi_scatter);
1307:   PetscFree(user->di_scatter);
1308:   PetscFree(user->sample_times);
1309:   return(0);
1310: }

1312: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1313: {
1315:   Vec            X;
1316:   PetscReal      unorm,ynorm;
1317:   AppCtx         *user = (AppCtx*)ptr;

1320:   TaoGetSolutionVector(tao,&X);
1321:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1322:   VecAXPY(user->ywork,-1.0,user->ytrue);
1323:   VecAXPY(user->uwork,-1.0,user->utrue);
1324:   VecNorm(user->uwork,NORM_2,&unorm);
1325:   VecNorm(user->ywork,NORM_2,&ynorm);
1326:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1327:   return(0);
1328: }