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: TaoSetInitialVector();
  8:    Routines: TaoSetObjectiveRoutine();
  9:    Routines: TaoSetGradientRoutine();
 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);if (ierr) return ierr;
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:   TaoSetInitialVector(tao,x);
200:   TaoSetObjectiveRoutine(tao, FormFunction, &user);
201:   TaoSetGradientRoutine(tao, 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:     TaoSetInitialVector(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 ierr;
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: {
247:   PetscReal      d1=0,d2=0;
248:   PetscInt       i,j;
249:   AppCtx         *user = (AppCtx*)ptr;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

384:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
385:   return(0);
386: }

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

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

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

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

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

431:   MatShellGetContext(J_shell,&user);
432:   MatMult(user->Grad,X,user->Swork);
433:   VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
434:   MatMult(user->Div,user->Swork,Y);
435:   VecAYPX(Y,user->ht,X);
436:   return(0);
437: }

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

446:   MatShellGetContext(J_shell,&user);

448:   /* sdiag(1./v) */
449:   VecSet(user->uwork,0);
450:   VecAXPY(user->uwork,-1.0,user->u);
451:   VecExp(user->uwork);

453:   /* sdiag(1./((Av*(1./v)).^2)) */
454:   MatMult(user->Av,user->uwork,user->Swork);
455:   VecPointwiseMult(user->Swork,user->Swork,user->Swork);
456:   VecReciprocal(user->Swork);

458:   /* (Av * (sdiag(1./v) * b)) */
459:   VecPointwiseMult(user->uwork,user->uwork,X);
460:   MatMult(user->Av,user->uwork,user->Twork);

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

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

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

477:   return(0);
478: }

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

487:   MatShellGetContext(J_shell,&user);

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

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

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

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

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

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

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

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

529:   PCShellGetContext(PC_shell,&user);

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

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

544:   MatShellGetContext(J_shell,&user);

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

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

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

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

575:   MatShellGetContext(J_shell,&user);

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

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

582:   KSPGetIterationNumber(user->solver,&its);
583:   user->ksp_its = user->ksp_its + its;

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

589:     KSPGetIterationNumber(user->solver,&its);
590:     user->ksp_its = user->ksp_its + its;

592:   }

594:   Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
595:   return(0);
596: }

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

604:   MatShellGetContext(J_shell,&user);

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

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

620:   MatShellGetContext(J_shell,&user);
621:    VecCopy(user->js_diag,X);
622:   return(0);

624: }

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

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

651: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
652: {

656:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
657:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
658:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
659:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
660:   return(0);
661: }

663: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
664: {
666:   PetscInt       i;

669:   for (i=0; i<nt; i++) {
670:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
671:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
672:   }
673:   return(0);
674: }

676: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
677: {

681:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
682:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
683:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
684:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
685:   return(0);
686: }

688: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
689: {
691:   PetscInt       i;

694:   for (i=0; i<nt; i++) {
695:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
696:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
697:   }
698:   return(0);
699: }

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

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

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

742:   PetscMalloc1(user->mx,&x);
743:   PetscMalloc1(user->mx,&y);
744:   PetscMalloc1(user->mx,&z);
745:   user->jformed = PETSC_FALSE;
746:   user->dsg_formed = PETSC_FALSE;

748:   n = user->mx * user->mx * user->mx;
749:   m = 3 * user->mx * user->mx * (user->mx-1);
750:   sqrt_beta = PetscSqrtScalar(user->beta);

752:   user->ksp_its = 0;
753:   user->ksp_its_initial = 0;

755:   stime = (PetscReal)user->nt/user->ns;
756:   PetscMalloc1(user->ns,&user->sample_times);
757:   for (i=0; i<user->ns; i++) {
758:     user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
759:   }

761:   VecCreate(PETSC_COMM_WORLD,&XX);
762:   VecCreate(PETSC_COMM_WORLD,&user->q);
763:   VecSetSizes(XX,PETSC_DECIDE,n);
764:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
765:   VecSetFromOptions(XX);
766:   VecSetFromOptions(user->q);

768:   VecDuplicate(XX,&YY);
769:   VecDuplicate(XX,&ZZ);
770:   VecDuplicate(XX,&XXwork);
771:   VecDuplicate(XX,&YYwork);
772:   VecDuplicate(XX,&ZZwork);
773:   VecDuplicate(XX,&UTwork);
774:   VecDuplicate(XX,&user->utrue);
775:   VecDuplicate(XX,&bc);

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

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

799:   VecAssemblyBegin(XX);
800:   VecAssemblyEnd(XX);
801:   VecAssemblyBegin(YY);
802:   VecAssemblyEnd(YY);
803:   VecAssemblyBegin(ZZ);
804:   VecAssemblyEnd(ZZ);
805:   VecAssemblyBegin(bc);
806:   VecAssemblyEnd(bc);

808:   /* Compute true parameter function
809:      ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
810:   VecCopy(XX,XXwork);
811:   VecCopy(YY,YYwork);
812:   VecCopy(ZZ,ZZwork);

814:   VecShift(XXwork,-0.5);
815:   VecShift(YYwork,-0.5);
816:   VecShift(ZZwork,-0.5);

818:   VecPointwiseMult(XXwork,XXwork,XXwork);
819:   VecPointwiseMult(YYwork,YYwork,YYwork);
820:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

822:   VecCopy(XXwork,user->utrue);
823:   VecAXPY(user->utrue,1.0,YYwork);
824:   VecAXPY(user->utrue,1.0,ZZwork);
825:   VecScale(user->utrue,-10.0);
826:   VecExp(user->utrue);
827:   VecShift(user->utrue,0.5);

829:   VecDestroy(&XX);
830:   VecDestroy(&YY);
831:   VecDestroy(&ZZ);
832:   VecDestroy(&XXwork);
833:   VecDestroy(&YYwork);
834:   VecDestroy(&ZZwork);
835:   VecDestroy(&UTwork);

837:   /* Initial guess and reference model */
838:   VecDuplicate(user->utrue,&user->ur);
839:   VecSet(user->ur,0.0);

841:   /* Generate Grad matrix */
842:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
843:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n);
844:   MatSetFromOptions(user->Grad);
845:   MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
846:   MatSeqAIJSetPreallocation(user->Grad,2,NULL);
847:   MatGetOwnershipRange(user->Grad,&istart,&iend);

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

872:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
873:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);

875:   /* Generate arithmetic averaging matrix Av */
876:   MatCreate(PETSC_COMM_WORLD,&user->Av);
877:   MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n);
878:   MatSetFromOptions(user->Av);
879:   MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
880:   MatSeqAIJSetPreallocation(user->Av,2,NULL);
881:   MatGetOwnershipRange(user->Av,&istart,&iend);

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

906:   MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
907:   MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);

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

912:   MatCreate(PETSC_COMM_WORLD,&user->L);
913:   MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n);
914:   MatSetFromOptions(user->L);
915:   MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
916:   MatSeqAIJSetPreallocation(user->L,2,NULL);
917:   MatGetOwnershipRange(user->L,&istart,&iend);

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

946:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
947:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);

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

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

954:   /* Build work vectors and matrices */
955:   VecCreate(PETSC_COMM_WORLD,&user->S);
956:   VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
957:   VecSetFromOptions(user->S);

959:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
960:   VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
961:   VecSetFromOptions(user->lwork);

963:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
964:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

966:   VecCreate(PETSC_COMM_WORLD,&user->d);
967:   VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
968:   VecSetFromOptions(user->d);

970:   VecDuplicate(user->S,&user->Swork);
971:   VecDuplicate(user->S,&user->Av_u);
972:   VecDuplicate(user->S,&user->Twork);
973:   VecDuplicate(user->S,&user->Rwork);
974:   VecDuplicate(user->y,&user->ywork);
975:   VecDuplicate(user->u,&user->uwork);
976:   VecDuplicate(user->u,&user->js_diag);
977:   VecDuplicate(user->c,&user->cwork);
978:   VecDuplicate(user->d,&user->dwork);

980:   /* Create matrix-free shell user->Js for computing A*x */
981:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
982:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
983:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
984:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
985:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

987:   /* Diagonal blocks of user->Js */
988:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
989:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
990:   /* Blocks are symmetric */
991:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);

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

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

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

1012:   /* Solver options and tolerances */
1013:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1014:   KSPSetType(user->solver,KSPCG);
1015:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
1016:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
1017:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1018:   KSPSetFromOptions(user->solver);
1019:   KSPGetPC(user->solver,&user->prec);
1020:   PCSetType(user->prec,PCSHELL);

1022:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1023:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
1024:   PCShellSetContext(user->prec,user);

1026:   /* Create scatter from y to y_1,y_2,...,y_nt */
1027:   PetscMalloc1(user->nt*user->m,&user->yi_scatter);
1028:   VecCreate(PETSC_COMM_WORLD,&yi);
1029:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
1030:   VecSetFromOptions(yi);
1031:   VecDuplicateVecs(yi,user->nt,&user->yi);
1032:   VecDuplicateVecs(yi,user->nt,&user->yiwork);

1034:   VecGetOwnershipRange(user->y,&lo2,&hi2);
1035:   istart = 0;
1036:   for (i=0; i<user->nt; i++) {
1037:     VecGetOwnershipRange(user->yi[i],&lo,&hi);
1038:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
1039:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
1040:     VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
1041:     istart = istart + hi-lo;
1042:     ISDestroy(&is_to_yi);
1043:     ISDestroy(&is_from_y);
1044:   }
1045:   VecDestroy(&yi);

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

1066:   /* Assemble RHS of forward problem */
1067:   VecCopy(bc,user->yiwork[0]);
1068:   for (i=1; i<user->nt; i++) {
1069:     VecSet(user->yiwork[i],0.0);
1070:   }
1071:   Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1072:   VecDestroy(&bc);

1074:   /* Compute true state function yt given ut */
1075:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1076:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1077:   VecSetFromOptions(user->ytrue);

1079:   /* First compute Av_u = Av*exp(-u) */
1080:   VecSet(user->uwork,0);
1081:   VecAXPY(user->uwork,-1.0,user->utrue); /*  Note: user->utrue */
1082:   VecExp(user->uwork);
1083:   MatMult(user->Av,user->uwork,user->Av_u);

1085:   /* Symbolic DSG = Div * Grad */
1086:   MatProductCreate(user->Div,user->Grad,NULL,&user->DSG);
1087:   MatProductSetType(user->DSG,MATPRODUCT_AB);
1088:   MatProductSetAlgorithm(user->DSG,"default");
1089:   MatProductSetFill(user->DSG,PETSC_DEFAULT);
1090:   MatProductSetFromOptions(user->DSG);
1091:   MatProductSymbolic(user->DSG);

1093:   user->dsg_formed = PETSC_TRUE;

1095:   /* Next form DSG = Div*Grad */
1096:   MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1097:   MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1098:   if (user->dsg_formed) {
1099:     MatProductNumeric(user->DSG);
1100:   } else {
1101:     MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&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:     MatProductNumeric(user->DSG);
1124:   } else {
1125:     MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&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: }

1330: /*TEST

1332:    build:
1333:       requires: !complex

1335:    test:
1336:       args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30
1337:       requires: !single

1339: TEST*/