Actual source code: parabolic.c

  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: TaoSetSolution();
  8:    Routines: TaoSetObjective();
  9:    Routines: TaoSetGradient();
 10:    Routines: TaoSetConstraintsRoutine();
 11:    Routines: TaoSetJacobianStateRoutine();
 12:    Routines: TaoSetJacobianDesignRoutine();
 13:    Routines: TaoSetStateDesignIS();
 14:    Routines: TaoSetFromOptions();
 15:    Routines: TaoSolve();
 16:    Routines: TaoDestroy();
 17:    Processors: n
 18: T*/

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

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

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

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

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

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

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

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

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

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

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

 76:   KSP solver;
 77:   PC prec;

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

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

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

105: PetscErrorCode DesignMatMult(Mat,Vec,Vec);
106: PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);

108: PetscErrorCode Gather_i(Vec,Vec*,VecScatter*,PetscInt);
109: PetscErrorCode Scatter_i(Vec,Vec*,VecScatter*,PetscInt);

111: static  char help[]="";

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

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

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:   TaoSetSolution(tao,x);
200:   TaoSetObjective(tao, FormFunction, &user);
201:   TaoSetGradient(tao, NULL, FormGradient, &user);
202:   TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user);

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

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

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

230:   TaoDestroy(&tao);
231:   VecDestroy(&x);
232:   VecDestroy(&x0);
233:   ParabolicDestroy(&user);

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

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

260:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
261:   MatMult(user->L,user->uwork,user->lwork);
262:   VecDot(user->lwork,user->lwork,&d2);

264:   *f = 0.5 * (d1 + user->alpha*d2);
265:   return 0;
266: }

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

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

295:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
296:   MatMult(user->L,user->uwork,user->lwork);
297:   MatMult(user->LT,user->lwork,user->uwork);
298:   VecScale(user->uwork, user->alpha);
299:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
300:   return 0;
301: }

303: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
304: {
305:   PetscReal      d1,d2;
306:   PetscInt       i,j;
307:   AppCtx         *user = (AppCtx*)ptr;

309:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
310:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
311:   for (j=0; j<user->ns; j++) {
312:     i = user->sample_times[j];
313:     MatMult(user->Qblock,user->yi[i],user->di[j]);
314:   }
315:   Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
316:   VecAXPY(user->dwork,-1.0,user->d);
317:   VecDot(user->dwork,user->dwork,&d1);
318:   Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
319:   VecSet(user->ywork,0.0);
320:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
321:   for (j=0; j<user->ns; j++) {
322:     i = user->sample_times[j];
323:     MatMult(user->QblockT,user->di[j],user->yiwork[i]);
324:   }
325:   Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);

327:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
328:   MatMult(user->L,user->uwork,user->lwork);
329:   VecDot(user->lwork,user->lwork,&d2);
330:   MatMult(user->LT,user->lwork,user->uwork);
331:   VecScale(user->uwork, user->alpha);
332:   *f = 0.5 * (d1 + user->alpha*d2);

334:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
335:   return 0;
336: }

338: /* ------------------------------------------------------------------- */
339: /* A
340: MatShell object
341: */
342: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
343: {
344:   AppCtx         *user = (AppCtx*)ptr;

346:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
347:   VecSet(user->uwork,0);
348:   VecAXPY(user->uwork,-1.0,user->u);
349:   VecExp(user->uwork);
350:   MatMult(user->Av,user->uwork,user->Av_u);
351:   VecCopy(user->Av_u,user->Swork);
352:   VecReciprocal(user->Swork);
353:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
354:   MatDiagonalScale(user->Divwork,NULL,user->Swork);
355:   if (user->dsg_formed) {
356:     MatProductNumeric(user->DSG);
357:   } else {
358:     MatMatMult(user->Divwork,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG);
359:     user->dsg_formed = PETSC_TRUE;
360:   }

362:   /* B = speye(nx^3) + ht*DSG; */
363:   MatScale(user->DSG,user->ht);
364:   MatShift(user->DSG,1.0);
365:   return 0;
366: }

368: /* ------------------------------------------------------------------- */
369: /* B */
370: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
371: {
372:   AppCtx         *user = (AppCtx*)ptr;

374:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
375:   return 0;
376: }

378: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
379: {
380:   PetscInt       i;
381:   AppCtx         *user;

383:   MatShellGetContext(J_shell,&user);
384:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
385:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
386:   for (i=1; i<user->nt; i++) {
387:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
388:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
389:   }
390:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
391:   return 0;
392: }

394: PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
395: {
396:   PetscInt       i;
397:   AppCtx         *user;

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

411: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
412: {
413:   AppCtx         *user;

415:   MatShellGetContext(J_shell,&user);
416:   MatMult(user->Grad,X,user->Swork);
417:   VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
418:   MatMult(user->Div,user->Swork,Y);
419:   VecAYPX(Y,user->ht,X);
420:   return 0;
421: }

423: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
424: {
425:   PetscInt       i;
426:   AppCtx         *user;

428:   MatShellGetContext(J_shell,&user);

430:   /* sdiag(1./v) */
431:   VecSet(user->uwork,0);
432:   VecAXPY(user->uwork,-1.0,user->u);
433:   VecExp(user->uwork);

435:   /* sdiag(1./((Av*(1./v)).^2)) */
436:   MatMult(user->Av,user->uwork,user->Swork);
437:   VecPointwiseMult(user->Swork,user->Swork,user->Swork);
438:   VecReciprocal(user->Swork);

440:   /* (Av * (sdiag(1./v) * b)) */
441:   VecPointwiseMult(user->uwork,user->uwork,X);
442:   MatMult(user->Av,user->uwork,user->Twork);

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

447:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
448:   for (i=0; i<user->nt; i++) {
449:     /* (sdiag(Grad*y(:,i)) */
450:     MatMult(user->Grad,user->yi[i],user->Twork);

452:     /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
453:     VecPointwiseMult(user->Twork,user->Twork,user->Swork);
454:     MatMult(user->Div,user->Twork,user->yiwork[i]);
455:     VecScale(user->yiwork[i],user->ht);
456:   }
457:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);

459:   return 0;
460: }

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

467:   MatShellGetContext(J_shell,&user);

469:   /* sdiag(1./((Av*(1./v)).^2)) */
470:   VecSet(user->uwork,0);
471:   VecAXPY(user->uwork,-1.0,user->u);
472:   VecExp(user->uwork);
473:   MatMult(user->Av,user->uwork,user->Rwork);
474:   VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork);
475:   VecReciprocal(user->Rwork);

477:   VecSet(Y,0.0);
478:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
479:   Scatter_i(X,user->yiwork,user->yi_scatter,user->nt);
480:   for (i=0; i<user->nt; i++) {
481:     /* (Div' * b(:,i)) */
482:     MatMult(user->Grad,user->yiwork[i],user->Swork);

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

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

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

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

496:     /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
497:     VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i]);
498:     VecAXPY(Y,user->ht,user->yiwork[i]);
499:   }
500:   return 0;
501: }

503: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
504: {
505:   AppCtx         *user;

507:   PCShellGetContext(PC_shell,&user);

509:   if (user->dsg_formed) {
510:     MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
511:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed");
512:   return 0;
513: }

515: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
516: {
517:   AppCtx         *user;
518:   PetscInt       its,i;

520:   MatShellGetContext(J_shell,&user);

522:   if (Y == user->ytrue) {
523:     /* First solve is done with true solution to set up problem */
524:     KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
525:   } else {
526:     KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
527:   }

529:   Scatter_i(X,user->yi,user->yi_scatter,user->nt);
530:   KSPSolve(user->solver,user->yi[0],user->yiwork[0]);
531:   KSPGetIterationNumber(user->solver,&its);
532:   user->ksp_its = user->ksp_its + its;

534:   for (i=1; i<user->nt; i++) {
535:     VecAXPY(user->yi[i],1.0,user->yiwork[i-1]);
536:     KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
537:     KSPGetIterationNumber(user->solver,&its);
538:     user->ksp_its = user->ksp_its + its;
539:   }
540:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
541:   return 0;
542: }

544: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
545: {
546:   AppCtx         *user;
547:   PetscInt       its,i;

549:   MatShellGetContext(J_shell,&user);

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

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

556:   KSPGetIterationNumber(user->solver,&its);
557:   user->ksp_its = user->ksp_its + its;

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

563:     KSPGetIterationNumber(user->solver,&its);
564:     user->ksp_its = user->ksp_its + its;

566:   }

568:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
569:   return 0;
570: }

572: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
573: {
574:   AppCtx         *user;

576:   MatShellGetContext(J_shell,&user);

578:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
579:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
580:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
581:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
582:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
583:   return 0;
584: }

586: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
587: {
588:   AppCtx         *user;

590:   MatShellGetContext(J_shell,&user);
591:   VecCopy(user->js_diag,X);
592:   return 0;

594: }

596: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
597: {
598:   /* con = Ay - q, A = [B  0  0 ... 0;
599:                        -I  B  0 ... 0;
600:                         0 -I  B ... 0;
601:                              ...     ;
602:                         0    ... -I B]
603:      B = ht * Div * Sigma * Grad + eye */
604:   PetscInt       i;
605:   AppCtx         *user = (AppCtx*)ptr;

607:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
608:   Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
609:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
610:   for (i=1; i<user->nt; i++) {
611:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
612:     VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
613:   }
614:   Gather_i(C,user->yiwork,user->yi_scatter,user->nt);
615:   VecAXPY(C,-1.0,user->q);
616:   return 0;
617: }

619: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
620: {
621:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
622:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
623:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
624:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
625:   return 0;
626: }

628: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
629: {
630:   PetscInt       i;

632:   for (i=0; i<nt; i++) {
633:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
634:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
635:   }
636:   return 0;
637: }

639: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
640: {
641:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
642:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
643:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
644:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
645:   return 0;
646: }

648: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
649: {
650:   PetscInt       i;

652:   for (i=0; i<nt; i++) {
653:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
654:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
655:   }
656:   return 0;
657: }

659: PetscErrorCode ParabolicInitialize(AppCtx *user)
660: {
661:   PetscInt       m,n,i,j,k,linear_index,istart,iend,iblock,lo,hi,lo2,hi2;
662:   Vec            XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork,yi,di,bc;
663:   PetscReal      *x, *y, *z;
664:   PetscReal      h,stime;
665:   PetscScalar    hinv,neg_hinv,half = 0.5,sqrt_beta;
666:   PetscInt       im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
667:   PetscReal      xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
668:   PetscScalar    v,vx,vy,vz;
669:   IS             is_from_y,is_to_yi,is_from_d,is_to_di;
670:   /* Data locations */
671:   PetscScalar xr[64] = {0.4970,     0.8498,     0.7814,     0.6268,     0.7782,     0.6402,     0.3617,     0.3160,
672:                         0.3610,     0.5298,     0.6987,     0.3331,     0.7962,     0.5596,     0.3866,     0.6774,
673:                         0.5407,     0.4518,     0.6702,     0.6061,     0.7580,     0.8997,     0.5198,     0.8326,
674:                         0.2138,     0.9198,     0.3000,     0.2833,     0.8288,     0.7076,     0.1820,     0.0728,
675:                         0.8447,     0.2367,     0.3239,     0.6413,     0.3114,     0.4731,     0.1192,     0.9273,
676:                         0.5724,     0.4331,     0.5136,     0.3547,     0.4413,     0.2602,     0.5698,     0.7278,
677:                         0.5261,     0.6230,     0.2454,     0.3948,     0.7479,     0.6582,     0.4660,     0.5594,
678:                         0.7574,     0.1143,     0.5900,     0.1065,     0.4260,     0.3294,     0.8276,     0.0756};

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

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

698:   PetscMalloc1(user->mx,&x);
699:   PetscMalloc1(user->mx,&y);
700:   PetscMalloc1(user->mx,&z);
701:   user->jformed = PETSC_FALSE;
702:   user->dsg_formed = PETSC_FALSE;

704:   n = user->mx * user->mx * user->mx;
705:   m = 3 * user->mx * user->mx * (user->mx-1);
706:   sqrt_beta = PetscSqrtScalar(user->beta);

708:   user->ksp_its = 0;
709:   user->ksp_its_initial = 0;

711:   stime = (PetscReal)user->nt/user->ns;
712:   PetscMalloc1(user->ns,&user->sample_times);
713:   for (i=0; i<user->ns; i++) {
714:     user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
715:   }

717:   VecCreate(PETSC_COMM_WORLD,&XX);
718:   VecCreate(PETSC_COMM_WORLD,&user->q);
719:   VecSetSizes(XX,PETSC_DECIDE,n);
720:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
721:   VecSetFromOptions(XX);
722:   VecSetFromOptions(user->q);

724:   VecDuplicate(XX,&YY);
725:   VecDuplicate(XX,&ZZ);
726:   VecDuplicate(XX,&XXwork);
727:   VecDuplicate(XX,&YYwork);
728:   VecDuplicate(XX,&ZZwork);
729:   VecDuplicate(XX,&UTwork);
730:   VecDuplicate(XX,&user->utrue);
731:   VecDuplicate(XX,&bc);

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

738:   VecGetOwnershipRange(XX,&istart,&iend);
739:   for (linear_index=istart; linear_index<iend; linear_index++) {
740:     i = linear_index % user->mx;
741:     j = ((linear_index-i)/user->mx) % user->mx;
742:     k = ((linear_index-i)/user->mx-j) / user->mx;
743:     vx = h*(i+0.5);
744:     vy = h*(j+0.5);
745:     vz = h*(k+0.5);
746:     VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
747:     VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
748:     VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES);
749:     if ((vx<0.6) && (vx>0.4) && (vy<0.6) && (vy>0.4) && (vy<0.6) && (vz<0.6) && (vz>0.4)) {
750:       v = 1000.0;
751:       VecSetValues(bc,1,&linear_index,&v,INSERT_VALUES);
752:     }
753:   }

755:   VecAssemblyBegin(XX);
756:   VecAssemblyEnd(XX);
757:   VecAssemblyBegin(YY);
758:   VecAssemblyEnd(YY);
759:   VecAssemblyBegin(ZZ);
760:   VecAssemblyEnd(ZZ);
761:   VecAssemblyBegin(bc);
762:   VecAssemblyEnd(bc);

764:   /* Compute true parameter function
765:      ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
766:   VecCopy(XX,XXwork);
767:   VecCopy(YY,YYwork);
768:   VecCopy(ZZ,ZZwork);

770:   VecShift(XXwork,-0.5);
771:   VecShift(YYwork,-0.5);
772:   VecShift(ZZwork,-0.5);

774:   VecPointwiseMult(XXwork,XXwork,XXwork);
775:   VecPointwiseMult(YYwork,YYwork,YYwork);
776:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

778:   VecCopy(XXwork,user->utrue);
779:   VecAXPY(user->utrue,1.0,YYwork);
780:   VecAXPY(user->utrue,1.0,ZZwork);
781:   VecScale(user->utrue,-10.0);
782:   VecExp(user->utrue);
783:   VecShift(user->utrue,0.5);

785:   VecDestroy(&XX);
786:   VecDestroy(&YY);
787:   VecDestroy(&ZZ);
788:   VecDestroy(&XXwork);
789:   VecDestroy(&YYwork);
790:   VecDestroy(&ZZwork);
791:   VecDestroy(&UTwork);

793:   /* Initial guess and reference model */
794:   VecDuplicate(user->utrue,&user->ur);
795:   VecSet(user->ur,0.0);

797:   /* Generate Grad matrix */
798:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
799:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n);
800:   MatSetFromOptions(user->Grad);
801:   MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
802:   MatSeqAIJSetPreallocation(user->Grad,2,NULL);
803:   MatGetOwnershipRange(user->Grad,&istart,&iend);

805:   for (i=istart; i<iend; i++) {
806:     if (i<m/3) {
807:       iblock = i / (user->mx-1);
808:       j = iblock*user->mx + (i % (user->mx-1));
809:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
810:       j = j+1;
811:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
812:     }
813:     if (i>=m/3 && i<2*m/3) {
814:       iblock = (i-m/3) / (user->mx*(user->mx-1));
815:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
816:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
817:       j = j + user->mx;
818:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
819:     }
820:     if (i>=2*m/3) {
821:       j = i-2*m/3;
822:       MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
823:       j = j + user->mx*user->mx;
824:       MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
825:     }
826:   }

828:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
829:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);

831:   /* Generate arithmetic averaging matrix Av */
832:   MatCreate(PETSC_COMM_WORLD,&user->Av);
833:   MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n);
834:   MatSetFromOptions(user->Av);
835:   MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
836:   MatSeqAIJSetPreallocation(user->Av,2,NULL);
837:   MatGetOwnershipRange(user->Av,&istart,&iend);

839:   for (i=istart; i<iend; i++) {
840:     if (i<m/3) {
841:       iblock = i / (user->mx-1);
842:       j = iblock*user->mx + (i % (user->mx-1));
843:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
844:       j = j+1;
845:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
846:     }
847:     if (i>=m/3 && i<2*m/3) {
848:       iblock = (i-m/3) / (user->mx*(user->mx-1));
849:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
850:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
851:       j = j + user->mx;
852:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
853:     }
854:     if (i>=2*m/3) {
855:       j = i-2*m/3;
856:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
857:       j = j + user->mx*user->mx;
858:       MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
859:     }
860:   }

862:   MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
863:   MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);

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

868:   MatCreate(PETSC_COMM_WORLD,&user->L);
869:   MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n);
870:   MatSetFromOptions(user->L);
871:   MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
872:   MatSeqAIJSetPreallocation(user->L,2,NULL);
873:   MatGetOwnershipRange(user->L,&istart,&iend);

875:   for (i=istart; i<iend; i++) {
876:     if (i<m/3) {
877:       iblock = i / (user->mx-1);
878:       j = iblock*user->mx + (i % (user->mx-1));
879:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
880:       j = j+1;
881:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
882:     }
883:     if (i>=m/3 && i<2*m/3) {
884:       iblock = (i-m/3) / (user->mx*(user->mx-1));
885:       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
886:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
887:       j = j + user->mx;
888:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
889:     }
890:     if (i>=2*m/3 && i<m) {
891:       j = i-2*m/3;
892:       MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
893:       j = j + user->mx*user->mx;
894:       MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
895:     }
896:     if (i>=m) {
897:       j = i - m;
898:       MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES);
899:     }
900:   }

902:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
903:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);

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

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

910:   /* Build work vectors and matrices */
911:   VecCreate(PETSC_COMM_WORLD,&user->S);
912:   VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
913:   VecSetFromOptions(user->S);

915:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
916:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
917:   VecSetFromOptions(user->lwork);

919:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
920:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

922:   VecCreate(PETSC_COMM_WORLD,&user->d);
923:   VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
924:   VecSetFromOptions(user->d);

926:   VecDuplicate(user->S,&user->Swork);
927:   VecDuplicate(user->S,&user->Av_u);
928:   VecDuplicate(user->S,&user->Twork);
929:   VecDuplicate(user->S,&user->Rwork);
930:   VecDuplicate(user->y,&user->ywork);
931:   VecDuplicate(user->u,&user->uwork);
932:   VecDuplicate(user->u,&user->js_diag);
933:   VecDuplicate(user->c,&user->cwork);
934:   VecDuplicate(user->d,&user->dwork);

936:   /* Create matrix-free shell user->Js for computing A*x */
937:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
938:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
939:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
940:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
941:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

943:   /* Diagonal blocks of user->Js */
944:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
945:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
946:   /* Blocks are symmetric */
947:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);

949:   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
950:      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
951:      This is an SSOR preconditioner for user->JsBlock. */
952:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlockPrec);
953:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult);
954:   /* JsBlockPrec is symmetric */
955:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMult);
956:   MatSetOption(user->JsBlockPrec,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);

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

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

968:   /* Solver options and tolerances */
969:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
970:   KSPSetType(user->solver,KSPCG);
971:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
972:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
973:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
974:   KSPSetFromOptions(user->solver);
975:   KSPGetPC(user->solver,&user->prec);
976:   PCSetType(user->prec,PCSHELL);

978:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
979:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
980:   PCShellSetContext(user->prec,user);

982:   /* Create scatter from y to y_1,y_2,...,y_nt */
983:   PetscMalloc1(user->nt*user->m,&user->yi_scatter);
984:   VecCreate(PETSC_COMM_WORLD,&yi);
985:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
986:   VecSetFromOptions(yi);
987:   VecDuplicateVecs(yi,user->nt,&user->yi);
988:   VecDuplicateVecs(yi,user->nt,&user->yiwork);

990:   VecGetOwnershipRange(user->y,&lo2,&hi2);
991:   istart = 0;
992:   for (i=0; i<user->nt; i++) {
993:     VecGetOwnershipRange(user->yi[i],&lo,&hi);
994:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
995:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
996:     VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
997:     istart = istart + hi-lo;
998:     ISDestroy(&is_to_yi);
999:     ISDestroy(&is_from_y);
1000:   }
1001:   VecDestroy(&yi);

1003:   /* Create scatter from d to d_1,d_2,...,d_ns */
1004:   PetscMalloc1(user->ns*user->ndata,&user->di_scatter);
1005:   VecCreate(PETSC_COMM_WORLD,&di);
1006:   VecSetSizes(di,PETSC_DECIDE,user->ndata);
1007:   VecSetFromOptions(di);
1008:   VecDuplicateVecs(di,user->ns,&user->di);
1009:   VecGetOwnershipRange(user->d,&lo2,&hi2);
1010:   istart = 0;
1011:   for (i=0; i<user->ns; i++) {
1012:     VecGetOwnershipRange(user->di[i],&lo,&hi);
1013:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_di);
1014:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d);
1015:     VecScatterCreate(user->d,is_from_d,user->di[i],is_to_di,&user->di_scatter[i]);
1016:     istart = istart + hi-lo;
1017:     ISDestroy(&is_to_di);
1018:     ISDestroy(&is_from_d);
1019:   }
1020:   VecDestroy(&di);

1022:   /* Assemble RHS of forward problem */
1023:   VecCopy(bc,user->yiwork[0]);
1024:   for (i=1; i<user->nt; i++) {
1025:     VecSet(user->yiwork[i],0.0);
1026:   }
1027:   Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1028:   VecDestroy(&bc);

1030:   /* Compute true state function yt given ut */
1031:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1032:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1033:   VecSetFromOptions(user->ytrue);

1035:   /* First compute Av_u = Av*exp(-u) */
1036:   VecSet(user->uwork,0);
1037:   VecAXPY(user->uwork,-1.0,user->utrue); /*  Note: user->utrue */
1038:   VecExp(user->uwork);
1039:   MatMult(user->Av,user->uwork,user->Av_u);

1041:   /* Symbolic DSG = Div * Grad */
1042:   MatProductCreate(user->Div,user->Grad,NULL,&user->DSG);
1043:   MatProductSetType(user->DSG,MATPRODUCT_AB);
1044:   MatProductSetAlgorithm(user->DSG,"default");
1045:   MatProductSetFill(user->DSG,PETSC_DEFAULT);
1046:   MatProductSetFromOptions(user->DSG);
1047:   MatProductSymbolic(user->DSG);

1049:   user->dsg_formed = PETSC_TRUE;

1051:   /* Next form DSG = Div*Grad */
1052:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1053:   MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1054:   if (user->dsg_formed) {
1055:     MatProductNumeric(user->DSG);
1056:   } else {
1057:     MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG);
1058:     user->dsg_formed = PETSC_TRUE;
1059:   }
1060:   /* B = speye(nx^3) + ht*DSG; */
1061:   MatScale(user->DSG,user->ht);
1062:   MatShift(user->DSG,1.0);

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

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

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

1075:   /* Next form DSG = Div*S*Grad */
1076:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1077:   MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1078:   if (user->dsg_formed) {
1079:     MatProductNumeric(user->DSG);
1080:   } else {
1081:     MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG);

1083:     user->dsg_formed = PETSC_TRUE;
1084:   }
1085:   /* B = speye(nx^3) + ht*DSG; */
1086:   MatScale(user->DSG,user->ht);
1087:   MatShift(user->DSG,1.0);

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

1092:   /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
1093:   MatCreate(PETSC_COMM_WORLD,&user->Qblock);
1094:   MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n);
1095:   MatSetFromOptions(user->Qblock);
1096:   MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL);
1097:   MatSeqAIJSetPreallocation(user->Qblock,8,NULL);

1099:   for (i=0; i<user->mx; i++) {
1100:     x[i] = h*(i+0.5);
1101:     y[i] = h*(i+0.5);
1102:     z[i] = h*(i+0.5);
1103:   }

1105:   MatGetOwnershipRange(user->Qblock,&istart,&iend);
1106:   nx = user->mx; ny = user->mx; nz = user->mx;
1107:   for (i=istart; i<iend; i++) {
1108:     xri = xr[i];
1109:     im = 0;
1110:     xim = x[im];
1111:     while (xri>xim && im<nx) {
1112:       im = im+1;
1113:       xim = x[im];
1114:     }
1115:     indx1 = im-1;
1116:     indx2 = im;
1117:     dx1 = xri - x[indx1];
1118:     dx2 = x[indx2] - xri;

1120:     yri = yr[i];
1121:     im = 0;
1122:     yim = y[im];
1123:     while (yri>yim && im<ny) {
1124:       im = im+1;
1125:       yim = y[im];
1126:     }
1127:     indy1 = im-1;
1128:     indy2 = im;
1129:     dy1 = yri - y[indy1];
1130:     dy2 = y[indy2] - yri;

1132:     zri = zr[i];
1133:     im = 0;
1134:     zim = z[im];
1135:     while (zri>zim && im<nz) {
1136:       im = im+1;
1137:       zim = z[im];
1138:     }
1139:     indz1 = im-1;
1140:     indz2 = im;
1141:     dz1 = zri - z[indz1];
1142:     dz2 = z[indz2] - zri;

1144:     Dx = x[indx2] - x[indx1];
1145:     Dy = y[indy2] - y[indy1];
1146:     Dz = z[indz2] - z[indz1];

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

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

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

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

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

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

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

1176:     j = indx2 + indy2*nx + indz2*nx*ny;
1177:     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1178:     MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1179:   }
1180:   MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY);
1181:   MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY);

1183:   MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1184:   MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);

1186:   /* Add noise to the measurement data */
1187:   VecSet(user->ywork,1.0);
1188:   VecAYPX(user->ywork,user->noise,user->ytrue);
1189:   Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
1190:   for (j=0; j<user->ns; j++) {
1191:     i = user->sample_times[j];
1192:     MatMult(user->Qblock,user->yiwork[i],user->di[j]);
1193:   }
1194:   Gather_i(user->d,user->di,user->di_scatter,user->ns);

1196:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1197:   KSPSetFromOptions(user->solver);
1198:   PetscFree(x);
1199:   PetscFree(y);
1200:   PetscFree(z);
1201:   return 0;
1202: }

1204: PetscErrorCode ParabolicDestroy(AppCtx *user)
1205: {
1206:   PetscInt       i;

1208:   MatDestroy(&user->Qblock);
1209:   MatDestroy(&user->QblockT);
1210:   MatDestroy(&user->Div);
1211:   MatDestroy(&user->Divwork);
1212:   MatDestroy(&user->Grad);
1213:   MatDestroy(&user->Av);
1214:   MatDestroy(&user->Avwork);
1215:   MatDestroy(&user->AvT);
1216:   MatDestroy(&user->DSG);
1217:   MatDestroy(&user->L);
1218:   MatDestroy(&user->LT);
1219:   KSPDestroy(&user->solver);
1220:   MatDestroy(&user->Js);
1221:   MatDestroy(&user->Jd);
1222:   MatDestroy(&user->JsInv);
1223:   MatDestroy(&user->JsBlock);
1224:   MatDestroy(&user->JsBlockPrec);
1225:   VecDestroy(&user->u);
1226:   VecDestroy(&user->uwork);
1227:   VecDestroy(&user->utrue);
1228:   VecDestroy(&user->y);
1229:   VecDestroy(&user->ywork);
1230:   VecDestroy(&user->ytrue);
1231:   VecDestroyVecs(user->nt,&user->yi);
1232:   VecDestroyVecs(user->nt,&user->yiwork);
1233:   VecDestroyVecs(user->ns,&user->di);
1234:   PetscFree(user->yi);
1235:   PetscFree(user->yiwork);
1236:   PetscFree(user->di);
1237:   VecDestroy(&user->c);
1238:   VecDestroy(&user->cwork);
1239:   VecDestroy(&user->ur);
1240:   VecDestroy(&user->q);
1241:   VecDestroy(&user->d);
1242:   VecDestroy(&user->dwork);
1243:   VecDestroy(&user->lwork);
1244:   VecDestroy(&user->S);
1245:   VecDestroy(&user->Swork);
1246:   VecDestroy(&user->Av_u);
1247:   VecDestroy(&user->Twork);
1248:   VecDestroy(&user->Rwork);
1249:   VecDestroy(&user->js_diag);
1250:   ISDestroy(&user->s_is);
1251:   ISDestroy(&user->d_is);
1252:   VecScatterDestroy(&user->state_scatter);
1253:   VecScatterDestroy(&user->design_scatter);
1254:   for (i=0; i<user->nt; i++) {
1255:     VecScatterDestroy(&user->yi_scatter[i]);
1256:   }
1257:   for (i=0; i<user->ns; i++) {
1258:     VecScatterDestroy(&user->di_scatter[i]);
1259:   }
1260:   PetscFree(user->yi_scatter);
1261:   PetscFree(user->di_scatter);
1262:   PetscFree(user->sample_times);
1263:   return 0;
1264: }

1266: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1267: {
1268:   Vec            X;
1269:   PetscReal      unorm,ynorm;
1270:   AppCtx         *user = (AppCtx*)ptr;

1272:   TaoGetSolution(tao,&X);
1273:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1274:   VecAXPY(user->ywork,-1.0,user->ytrue);
1275:   VecAXPY(user->uwork,-1.0,user->utrue);
1276:   VecNorm(user->uwork,NORM_2,&unorm);
1277:   VecNorm(user->ywork,NORM_2,&ynorm);
1278:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1279:   return 0;
1280: }

1282: /*TEST

1284:    build:
1285:       requires: !complex

1287:    test:
1288:       args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30
1289:       requires: !single

1291: TEST*/