Actual source code: parabolic.c

petsc-3.12.5 2020-03-29
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: TaoSolve();
 16:    Routines: TaoDestroy();
 17:    Processors: n
 18: T*/



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

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

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

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

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

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

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

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

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

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

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

 78:   KSP solver;
 79:   PC prec;

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


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

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

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

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

114: static  char help[]="";

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

481:   return(0);
482: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

596:   }

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

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

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

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

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

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

628: }

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

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


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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1249: PetscErrorCode ParabolicDestroy(AppCtx *user)
1250: {
1252:   PetscInt       i;

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

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

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


1332: /*TEST

1334:    build:
1335:       requires: !complex

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

1341: TEST*/