Actual source code: parabolic.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: #include <petsc/private/taoimpl.h>

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



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

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

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

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

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

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

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

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

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

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

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

 79:   KSP solver;
 80:   PC prec;

 82:   PetscInt ksp_its;
 83:   PetscInt ksp_its_initial;
 84: } AppCtx;


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

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

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

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

115: static  char help[]="";

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

444: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
445: {
447:   PetscInt       i;
448:   AppCtx         *user;

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

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

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

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

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

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

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

482:   return(0);
483: }

485: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
486: {
488:   PetscInt       i;
489:   AppCtx         *user;

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

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

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

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

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

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

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

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

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

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

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

542: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
543: {
545:   AppCtx         *user;
546:   PetscInt       its,i;

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

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

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

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

573: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
574: {
576:   AppCtx         *user;
577:   PetscInt       its,i;

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

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

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

587:   KSPGetIterationNumber(user->solver,&its);
588:   user->ksp_its = user->ksp_its + its;

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

594:     KSPGetIterationNumber(user->solver,&its);
595:     user->ksp_its = user->ksp_its + its;

597:   }

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

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

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

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

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

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

629: }

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

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


657: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
658: {

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

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

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


683: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
684: {

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

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

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

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

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

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

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

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

759:   user->ksp_its = 0;
760:   user->ksp_its_initial = 0;

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

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

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

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

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

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

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

821:   VecShift(XXwork,-0.5);
822:   VecShift(YYwork,-0.5);
823:   VecShift(ZZwork,-0.5);

825:   VecPointwiseMult(XXwork,XXwork,XXwork);
826:   VecPointwiseMult(YYwork,YYwork,YYwork);
827:   VecPointwiseMult(ZZwork,ZZwork,ZZwork);

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

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

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

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

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

879:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
880:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);


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

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

914:   MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
915:   MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);

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

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

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

954:   MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
955:   MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);

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

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

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

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

971:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
972:   MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);

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

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

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

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

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

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

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

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

1030:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1031:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
1032:   PCShellSetContext(user->prec,user);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1229:   MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1230:   MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);

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

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

1250: PetscErrorCode ParabolicDestroy(AppCtx *user)
1251: {
1253:   PetscInt       i;

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

1314: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1315: {
1317:   Vec            X;
1318:   PetscReal      unorm,ynorm;
1319:   AppCtx         *user = (AppCtx*)ptr;

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


1333: /*TEST

1335:    build:
1336:       requires: !complex

1338:    test:
1339:       args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30
1340:       requires: !single

1342: TEST*/