Actual source code: hyperbolic.c

petsc-master 2014-10-23
Report Typos and Errors
  1: #include <petsctao.h>
  2: #include "../src/tao/pde_constrained/impls/lcl/lcl.h"

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

 22: typedef struct {
 23:   PetscInt n; /*  Number of variables */
 24:   PetscInt m; /*  Number of constraints */
 25:   PetscInt mx; /*  grid points in each direction */
 26:   PetscInt nt; /*  Number of time steps */
 27:   PetscInt ndata; /*  Number of data points per sample */
 28:   IS       s_is;
 29:   IS       d_is;
 30:   VecScatter state_scatter;
 31:   VecScatter design_scatter;
 32:   VecScatter *uxi_scatter,*uyi_scatter,*ux_scatter,*uy_scatter,*ui_scatter;
 33:   VecScatter *yi_scatter;

 35:   Mat       Js,Jd,JsBlockPrec,JsInv,JsBlock;
 36:   PetscBool jformed,c_formed;

 38:   PetscReal alpha; /*  Regularization parameter */
 39:   PetscReal gamma;
 40:   PetscReal ht; /*  Time step */
 41:   PetscReal T; /*  Final time */
 42:   Mat Q,QT;
 43:   Mat L,LT;
 44:   Mat Div,Divwork,Divxy[2];
 45:   Mat Grad,Gradxy[2];
 46:   Mat M;
 47:   Mat *C,*Cwork;
 48:   /* Mat Hs,Hd,Hsd; */
 49:   Vec q;
 50:   Vec ur; /*  reference */

 52:   Vec d;
 53:   Vec dwork;

 55:   Vec y; /*  state variables */
 56:   Vec ywork;
 57:   Vec ytrue;
 58:   Vec *yi,*yiwork,*ziwork;
 59:   Vec *uxi,*uyi,*uxiwork,*uyiwork,*ui,*uiwork;

 61:   Vec u; /*  design variables */
 62:   Vec uwork,vwork;
 63:   Vec utrue;

 65:   Vec js_diag;

 67:   Vec c; /*  constraint vector */
 68:   Vec cwork;

 70:   Vec lwork;

 72:   KSP      solver;
 73:   PC       prec;
 74:   PetscInt block_index;

 76:   TAO_LCL *lcl;
 77:   PetscInt ksp_its;
 78:   PetscInt ksp_its_initial;
 79: } AppCtx;


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

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

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

108: PetscErrorCode Scatter_yi(Vec,Vec*,VecScatter*,PetscInt); /*  y to y1,y2,...,y_nt */
109: PetscErrorCode Gather_yi(Vec,Vec*,VecScatter*,PetscInt);
110: PetscErrorCode Scatter_uxi_uyi(Vec,Vec*,VecScatter*,Vec*,VecScatter*,PetscInt); /*  u to ux_1,uy_1,ux_2,uy_2,...,u */
111: PetscErrorCode Gather_uxi_uyi(Vec,Vec*,VecScatter*,Vec*,VecScatter*,PetscInt);

113: static  char help[]="";

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

130:   PetscInitialize(&argc, &argv, (char*)0,help);

132:   user.mx = 32;
133:   PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL);
134:   user.nt = 16;
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.alpha = 10.0;
139:   PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL);
140:   user.T = 1.0/32.0;
141:   PetscOptionsReal("-Tfinal","Final time","",user.T,&user.T,NULL);

143:   user.m = user.mx*user.mx*user.nt; /*  number of constraints */
144:   user.n = user.mx*user.mx*3*user.nt; /*  number of variables */
145:   user.ht = user.T/user.nt; /*  Time step */
146:   user.gamma = user.T*user.ht / (user.mx*user.mx);

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

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

168:   VecGetOwnershipRange(user.y,&lo,&hi);
169:   VecGetOwnershipRange(user.u,&lo2,&hi2);

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

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

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

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

185:   /* Create TAO solver and set desired solution method */
186:   TaoCreate(PETSC_COMM_WORLD,&tao);
187:   TaoSetType(tao,TAOLCL);
188:   user.lcl = (TAO_LCL*)(tao->data);

190:   /* Set up initial vectors and matrices */
191:   HyperbolicInitialize(&user);

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

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

203:   TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, (void *)&user);

205:   TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, (void *)&user);

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

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

231:   TaoGetConvergedReason(tao,&reason);

233:   if (reason < 0) {
234:     PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge.\n");
235:   } else {
236:     PetscPrintf(MPI_COMM_WORLD, "Optimization terminated with status %D.\n", reason);
237:   }

239:   TaoDestroy(&tao);
240:   VecDestroy(&x);
241:   VecDestroy(&x0);
242:   HyperbolicDestroy(&user);
243:   PetscFinalize();
244:   return 0;
245: }
246: /* ------------------------------------------------------------------- */
249: /*
250:    dwork = Qy - d
251:    lwork = L*(u-ur).^2
252:    f = 1/2 * (dwork.dork + alpha*y.lwork)
253: */
254: PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
255: {
257:   PetscReal      d1=0,d2=0;
258:   AppCtx         *user = (AppCtx*)ptr;

261:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
262:   MatMult(user->Q,user->y,user->dwork);
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:   VecPointwiseMult(user->uwork,user->uwork,user->uwork);
268:   MatMult(user->L,user->uwork,user->lwork);
269:   VecDot(user->y,user->lwork,&d2);
270:   *f = 0.5 * (d1 + user->alpha*d2);
271:   return(0);
272: }

274: /* ------------------------------------------------------------------- */
277: /*
278:     state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2
279:     design: g_d = alpha*(L'y).*(u-ur)
280: */
281: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
282: {
284:   AppCtx         *user = (AppCtx*)ptr;

287:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
288:   MatMult(user->Q,user->y,user->dwork);
289:   VecAXPY(user->dwork,-1.0,user->d);

291:   MatMult(user->QT,user->dwork,user->ywork);

293:   MatMult(user->LT,user->y,user->uwork);
294:   VecWAXPY(user->vwork,-1.0,user->ur,user->u);
295:   VecPointwiseMult(user->uwork,user->vwork,user->uwork);
296:   VecScale(user->uwork,user->alpha);

298:   VecPointwiseMult(user->vwork,user->vwork,user->vwork);
299:   MatMult(user->L,user->vwork,user->lwork);
300:   VecAXPY(user->ywork,0.5*user->alpha,user->lwork);

302:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
303:   return(0);
304: }

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

315:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
316:   MatMult(user->Q,user->y,user->dwork);
317:   VecAXPY(user->dwork,-1.0,user->d);

319:   MatMult(user->QT,user->dwork,user->ywork);

321:   VecDot(user->dwork,user->dwork,&d1);

323:   MatMult(user->LT,user->y,user->uwork);
324:   VecWAXPY(user->vwork,-1.0,user->ur,user->u);
325:   VecPointwiseMult(user->uwork,user->vwork,user->uwork);
326:   VecScale(user->uwork,user->alpha);

328:   VecPointwiseMult(user->vwork,user->vwork,user->vwork);
329:   MatMult(user->L,user->vwork,user->lwork);
330:   VecAXPY(user->ywork,0.5*user->alpha,user->lwork);

332:   VecDot(user->y,user->lwork,&d2);

334:   *f = 0.5 * (d1 + user->alpha*d2);
335:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
336:   return(0);
337: }

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

352:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
353:   Scatter_yi(user->u,user->ui,user->ui_scatter,user->nt);
354:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);
355:   for (i=0; i<user->nt; i++){
356:     MatCopy(user->Divxy[0],user->C[i],SAME_NONZERO_PATTERN);
357:     MatCopy(user->Divxy[1],user->Cwork[i],SAME_NONZERO_PATTERN);

359:     MatDiagonalScale(user->C[i],NULL,user->uxi[i]);
360:     MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]);
361:     MatAXPY(user->C[i],1.0,user->Cwork[i],SUBSET_NONZERO_PATTERN);
362:     MatScale(user->C[i],user->ht);
363:     MatShift(user->C[i],1.0);
364:   }
365:   return(0);
366: }

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

378:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
379:   return(0);
380: }

384: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
385: {
387:   PetscInt       i;
388:   AppCtx         *user;

391:   MatShellGetContext(J_shell,(void**)&user);
392:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);
393:   user->block_index = 0;
394:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);

396:   for (i=1; i<user->nt; i++){
397:     user->block_index = i;
398:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
399:     MatMult(user->M,user->yi[i-1],user->ziwork[i-1]);
400:     VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]);
401:   }
402:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
403:   return(0);
404: }

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

415:   MatShellGetContext(J_shell,(void**)&user);
416:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);

418:   for (i=0; i<user->nt-1; i++){
419:     user->block_index = i;
420:     MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]);
421:     MatMult(user->M,user->yi[i+1],user->ziwork[i+1]);
422:     VecAXPY(user->yiwork[i],-1.0,user->ziwork[i+1]);
423:   }

425:   i = user->nt-1;
426:   user->block_index = i;
427:   MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]);
428:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
429:   return(0);
430: }

434: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
435: {
437:   PetscInt       i;
438:   AppCtx         *user;

441:   MatShellGetContext(J_shell,(void**)&user);
442:   i = user->block_index;
443:   VecPointwiseMult(user->uxiwork[i],X,user->uxi[i]);
444:   VecPointwiseMult(user->uyiwork[i],X,user->uyi[i]);
445:   Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
446:   MatMult(user->Div,user->uiwork[i],Y);
447:   VecAYPX(Y,user->ht,X);
448:   return(0);
449: }

453: PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y)
454: {
456:   PetscInt       i;
457:   AppCtx         *user;

460:   MatShellGetContext(J_shell,(void**)&user);
461:   i = user->block_index;
462:   MatMult(user->Grad,X,user->uiwork[i]);
463:   Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
464:   VecPointwiseMult(user->uxiwork[i],user->uxi[i],user->uxiwork[i]);
465:   VecPointwiseMult(user->uyiwork[i],user->uyi[i],user->uyiwork[i]);
466:   VecWAXPY(Y,1.0,user->uxiwork[i],user->uyiwork[i]);
467:   VecAYPX(Y,user->ht,X);
468:   return(0);
469: }

473: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
474: {
476:   PetscInt       i;
477:   AppCtx         *user;

480:   MatShellGetContext(J_shell,(void**)&user);
481:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
482:   Scatter_uxi_uyi(X,user->uxiwork,user->uxi_scatter,user->uyiwork,user->uyi_scatter,user->nt);
483:   for (i=0; i<user->nt; i++){
484:     VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]);
485:     VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]);
486:     Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
487:     MatMult(user->Div,user->uiwork[i],user->ziwork[i]);
488:     VecScale(user->ziwork[i],user->ht);
489:   }
490:   Gather_yi(Y,user->ziwork,user->yi_scatter,user->nt);
491:   return(0);
492: }

496: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
497: {
499:   PetscInt       i;
500:   AppCtx         *user;

503:   MatShellGetContext(J_shell,(void**)&user);
504:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
505:   Scatter_yi(X,user->yiwork,user->yi_scatter,user->nt);
506:   for (i=0; i<user->nt; i++){
507:     MatMult(user->Grad,user->yiwork[i],user->uiwork[i]);
508:     Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
509:     VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]);
510:     VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]);
511:     Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
512:     VecScale(user->uiwork[i],user->ht);
513:   }
514:   Gather_yi(Y,user->uiwork,user->ui_scatter,user->nt);
515:   return(0);
516: }

520: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
521: {
523:   PetscInt       i;
524:   AppCtx         *user;

527:   PCShellGetContext(PC_shell,(void**)&user);
528:   i = user->block_index;
529:   if (user->c_formed) {
530:     MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
531:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
532:   return(0);
533: }

537: PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y)
538: {
540:   PetscInt       i;
541:   AppCtx         *user;

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

546:   i = user->block_index;
547:   if (user->c_formed) {
548:     MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
549:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
550:   return(0);
551: }

555: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
556: {
558:   AppCtx         *user;
559:   PetscReal      tau;
560:   PetscInt       its,i;

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

565:   if (Y == user->ytrue) {
566:     KSPSetTolerances(user->solver,1e-4,1e-20,PETSC_DEFAULT,PETSC_DEFAULT);
567:   } else if (user->lcl) {
568:     tau = user->lcl->tau[user->lcl->solve_type];
569:     KSPSetTolerances(user->solver,tau,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
570:   }
571:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);
572:   Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt);
573:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

575:   user->block_index = 0;
576:   KSPSolve(user->solver,user->yi[0],user->yiwork[0]);

578:   KSPGetIterationNumber(user->solver,&its);
579:   user->ksp_its = user->ksp_its + its;
580:   for (i=1; i<user->nt; i++){
581:     MatMult(user->M,user->yiwork[i-1],user->ziwork[i-1]);
582:     VecAXPY(user->yi[i],1.0,user->ziwork[i-1]);
583:     user->block_index = i;
584:     KSPSolve(user->solver,user->yi[i],user->yiwork[i]);

586:     KSPGetIterationNumber(user->solver,&its);
587:     user->ksp_its = user->ksp_its + its;
588:   }
589:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
590:   return(0);
591: }

595: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
596: {
598:   AppCtx         *user;
599:   PetscReal      tau;
600:   PetscInt       its,i;

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

605:   if (user->lcl) {
606:     tau = user->lcl->tau[user->lcl->solve_type];
607:     KSPSetTolerances(user->solver,tau,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
608:   }
609:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);
610:   Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt);
611:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

613:   i = user->nt - 1;
614:   user->block_index = i;
615:   KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]);

617:   KSPGetIterationNumber(user->solver,&its);
618:   user->ksp_its = user->ksp_its + its;

620:   for (i=user->nt-2; i>=0; i--){
621:     MatMult(user->M,user->yiwork[i+1],user->ziwork[i+1]);
622:     VecAXPY(user->yi[i],1.0,user->ziwork[i+1]);
623:     user->block_index = i;
624:     KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]);

626:     KSPGetIterationNumber(user->solver,&its);
627:     user->ksp_its = user->ksp_its + its;
628:   }
629:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
630:   return(0);
631: }

635: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
636: {
638:   AppCtx         *user;

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

643:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
644:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
645:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
646:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
647:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
648:   return(0);
649: }

653: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
654: {
656:   AppCtx         *user;

659:   MatShellGetContext(J_shell,(void**)&user);
660:    VecCopy(user->js_diag,X);
661:   return(0);
662: }

666: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
667: {
668:   /* con = Ay - q, A = [C(u1)  0     0     ...   0;
669:                          -M  C(u2)   0     ...   0;
670:                           0   -M   C(u3)   ...   0;
671:                                       ...         ;
672:                           0    ...      -M C(u_nt)]
673:      C(u) = eye + ht*Div*[diag(u1); diag(u2)]       */
675:   PetscInt       i;
676:   AppCtx         *user = (AppCtx*)ptr;

679:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
680:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
681:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

683:   user->block_index = 0;
684:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);

686:   for (i=1; i<user->nt; i++){
687:     user->block_index = i;
688:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
689:     MatMult(user->M,user->yi[i-1],user->ziwork[i-1]);
690:     VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]);
691:   }

693:   Gather_yi(C,user->yiwork,user->yi_scatter,user->nt);
694:   VecAXPY(C,-1.0,user->q);

696:   return(0);
697: }


702: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
703: {

707:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
708:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
709:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
710:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
711:   return(0);
712: }

716: PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
717: {
719:   PetscInt       i;

722:   for (i=0; i<nt; i++){
723:     VecScatterBegin(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD);
724:     VecScatterEnd(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD);
725:     VecScatterBegin(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD);
726:     VecScatterEnd(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD);
727:   }
728:   return(0);
729: }

733: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
734: {

738:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
739:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
740:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
741:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
742:   return(0);
743: }

747: PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
748: {
750:   PetscInt       i;

753:   for (i=0; i<nt; i++){
754:     VecScatterBegin(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE);
755:     VecScatterEnd(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE);
756:     VecScatterBegin(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE);
757:     VecScatterEnd(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE);
758:   }
759:   return(0);
760: }

764: PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
765: {
767:   PetscInt       i;

770:   for (i=0; i<nt; i++){
771:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
772:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
773:   }
774:   return(0);
775: }

779: PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
780: {
782:   PetscInt       i;

785:   for (i=0; i<nt; i++){
786:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
787:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
788:   }
789:   return(0);
790: }

794: PetscErrorCode HyperbolicInitialize(AppCtx *user)
795: {
797:   PetscInt       n,i,j,linear_index,istart,iend,iblock,lo,hi;
798:   Vec            XX,YY,XXwork,YYwork,yi,uxi,ui,bc;
799:   PetscReal      h,sum;
800:   PetscScalar    hinv,neg_hinv,quarter=0.25,one=1.0,half_hinv,neg_half_hinv;
801:   PetscScalar    vx,vy,zero=0.0;
802:   IS             is_from_y,is_to_yi,is_from_u,is_to_uxi,is_to_uyi;

805:   user->jformed = PETSC_FALSE;
806:   user->c_formed = PETSC_FALSE;

808:   user->ksp_its = 0;
809:   user->ksp_its_initial = 0;

811:   n = user->mx * user->mx;

813:   h = 1.0/user->mx;
814:   hinv = user->mx;
815:   neg_hinv = -hinv;
816:   half_hinv = hinv / 2.0;
817:   neg_half_hinv = neg_hinv / 2.0;

819:   /* Generate Grad matrix */
820:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
821:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,2*n,n);
822:   MatSetFromOptions(user->Grad);
823:   MatMPIAIJSetPreallocation(user->Grad,3,NULL,3,NULL);
824:   MatSeqAIJSetPreallocation(user->Grad,3,NULL);
825:   MatGetOwnershipRange(user->Grad,&istart,&iend);

827:   for (i=istart; i<iend; i++){
828:     if (i<n){
829:       iblock = i / user->mx;
830:       j = iblock*user->mx + ((i+user->mx-1) % user->mx);
831:       MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES);
832:       j = iblock*user->mx + ((i+1) % user->mx);
833:       MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES);
834:     }
835:     if (i>=n){
836:       j = (i - user->mx) % n;
837:       MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES);
838:       j = (j + 2*user->mx) % n;
839:       MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES);
840:     }
841:   }

843:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
844:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);

846:   MatCreate(PETSC_COMM_WORLD,&user->Gradxy[0]);
847:   MatSetSizes(user->Gradxy[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
848:   MatSetFromOptions(user->Gradxy[0]);
849:   MatMPIAIJSetPreallocation(user->Gradxy[0],3,NULL,3,NULL);
850:   MatSeqAIJSetPreallocation(user->Gradxy[0],3,NULL);
851:   MatGetOwnershipRange(user->Gradxy[0],&istart,&iend);

853:   for (i=istart; i<iend; i++){
854:     iblock = i / user->mx;
855:     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
856:     MatSetValues(user->Gradxy[0],1,&i,1,&j,&half_hinv,INSERT_VALUES);
857:     j = iblock*user->mx + ((i+1) % user->mx);
858:     MatSetValues(user->Gradxy[0],1,&i,1,&j,&neg_half_hinv,INSERT_VALUES);
859:     MatSetValues(user->Gradxy[0],1,&i,1,&i,&zero,INSERT_VALUES);
860:   }
861:   MatAssemblyBegin(user->Gradxy[0],MAT_FINAL_ASSEMBLY);
862:   MatAssemblyEnd(user->Gradxy[0],MAT_FINAL_ASSEMBLY);

864:   MatCreate(PETSC_COMM_WORLD,&user->Gradxy[1]);
865:   MatSetSizes(user->Gradxy[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
866:   MatSetFromOptions(user->Gradxy[1]);
867:   MatMPIAIJSetPreallocation(user->Gradxy[1],3,NULL,3,NULL);
868:   MatSeqAIJSetPreallocation(user->Gradxy[1],3,NULL);
869:   MatGetOwnershipRange(user->Gradxy[1],&istart,&iend);

871:   for (i=istart; i<iend; i++){
872:     j = (i + n - user->mx) % n;
873:     MatSetValues(user->Gradxy[1],1,&i,1,&j,&half_hinv,INSERT_VALUES);
874:     j = (j + 2*user->mx) % n;
875:     MatSetValues(user->Gradxy[1],1,&i,1,&j,&neg_half_hinv,INSERT_VALUES);
876:     MatSetValues(user->Gradxy[1],1,&i,1,&i,&zero,INSERT_VALUES);
877:   }
878:   MatAssemblyBegin(user->Gradxy[1],MAT_FINAL_ASSEMBLY);
879:   MatAssemblyEnd(user->Gradxy[1],MAT_FINAL_ASSEMBLY);

881:   /* Generate Div matrix */
882:   MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div);
883:   MatTranspose(user->Gradxy[0],MAT_INITIAL_MATRIX,&user->Divxy[0]);
884:   MatTranspose(user->Gradxy[1],MAT_INITIAL_MATRIX,&user->Divxy[1]);

886:   /* Off-diagonal averaging matrix */
887:   MatCreate(PETSC_COMM_WORLD,&user->M);
888:   MatSetSizes(user->M,PETSC_DECIDE,PETSC_DECIDE,n,n);
889:   MatSetFromOptions(user->M);
890:   MatMPIAIJSetPreallocation(user->M,4,NULL,4,NULL);
891:   MatSeqAIJSetPreallocation(user->M,4,NULL);
892:   MatGetOwnershipRange(user->M,&istart,&iend);

894:   for (i=istart; i<iend; i++){
895:     /* kron(Id,Av) */
896:     iblock = i / user->mx;
897:     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
898:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
899:     j = iblock*user->mx + ((i+1) % user->mx);
900:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);

902:     /* kron(Av,Id) */
903:     j = (i + user->mx) % n;
904:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
905:     j = (i + n - user->mx) % n;
906:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
907:   }
908:   MatAssemblyBegin(user->M,MAT_FINAL_ASSEMBLY);
909:   MatAssemblyEnd(user->M,MAT_FINAL_ASSEMBLY);

911:   /* Generate 2D grid */
912:   VecCreate(PETSC_COMM_WORLD,&XX);
913:   VecCreate(PETSC_COMM_WORLD,&user->q);
914:   VecSetSizes(XX,PETSC_DECIDE,n);
915:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
916:   VecSetFromOptions(XX);
917:   VecSetFromOptions(user->q);

919:   VecDuplicate(XX,&YY);
920:   VecDuplicate(XX,&XXwork);
921:   VecDuplicate(XX,&YYwork);
922:   VecDuplicate(XX,&user->d);
923:   VecDuplicate(XX,&user->dwork);

925:   VecGetOwnershipRange(XX,&istart,&iend);
926:   for (linear_index=istart; linear_index<iend; linear_index++){
927:     i = linear_index % user->mx;
928:     j = (linear_index-i)/user->mx;
929:     vx = h*(i+0.5);
930:     vy = h*(j+0.5);
931:     VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
932:     VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
933:   }

935:   VecAssemblyBegin(XX);
936:   VecAssemblyEnd(XX);
937:   VecAssemblyBegin(YY);
938:   VecAssemblyEnd(YY);

940:   /* Compute final density function yT
941:      yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2))
942:      yT = yT / (h^2*sum(yT)) */
943:   VecCopy(XX,XXwork);
944:   VecCopy(YY,YYwork);

946:   VecShift(XXwork,-0.25);
947:   VecShift(YYwork,-0.25);

949:   VecPointwiseMult(XXwork,XXwork,XXwork);
950:   VecPointwiseMult(YYwork,YYwork,YYwork);

952:   VecCopy(XXwork,user->dwork);
953:   VecAXPY(user->dwork,1.0,YYwork);
954:   VecScale(user->dwork,-30.0);
955:   VecExp(user->dwork);
956:   VecCopy(user->dwork,user->d);

958:   VecCopy(XX,XXwork);
959:   VecCopy(YY,YYwork);

961:   VecShift(XXwork,-0.75);
962:   VecShift(YYwork,-0.75);

964:   VecPointwiseMult(XXwork,XXwork,XXwork);
965:   VecPointwiseMult(YYwork,YYwork,YYwork);

967:   VecCopy(XXwork,user->dwork);
968:   VecAXPY(user->dwork,1.0,YYwork);
969:   VecScale(user->dwork,-30.0);
970:   VecExp(user->dwork);

972:   VecAXPY(user->d,1.0,user->dwork);
973:   VecShift(user->d,1.0);
974:   VecSum(user->d,&sum);
975:   VecScale(user->d,1.0/(h*h*sum));

977:   /* Initial conditions of forward problem */
978:   VecDuplicate(XX,&bc);
979:   VecCopy(XX,XXwork);
980:   VecCopy(YY,YYwork);

982:   VecShift(XXwork,-0.5);
983:   VecShift(YYwork,-0.5);

985:   VecPointwiseMult(XXwork,XXwork,XXwork);
986:   VecPointwiseMult(YYwork,YYwork,YYwork);

988:   VecWAXPY(bc,1.0,XXwork,YYwork);
989:   VecScale(bc,-50.0);
990:   VecExp(bc);
991:   VecShift(bc,1.0);
992:   VecSum(bc,&sum);
993:   VecScale(bc,1.0/(h*h*sum));

995:   /* Create scatter from y to y_1,y_2,...,y_nt */
996:   /*  TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */
997:   PetscMalloc1(user->nt*user->mx*user->mx,&user->yi_scatter);
998:   VecCreate(PETSC_COMM_WORLD,&yi);
999:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx);
1000:   VecSetFromOptions(yi);
1001:   VecDuplicateVecs(yi,user->nt,&user->yi);
1002:   VecDuplicateVecs(yi,user->nt,&user->yiwork);
1003:   VecDuplicateVecs(yi,user->nt,&user->ziwork);
1004:   for (i=0; i<user->nt; i++){
1005:     VecGetOwnershipRange(user->yi[i],&lo,&hi);
1006:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
1007:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+i*user->mx*user->mx,1,&is_from_y);
1008:     VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
1009:     ISDestroy(&is_to_yi);
1010:     ISDestroy(&is_from_y);
1011:   }

1013:   /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */
1014:   /*  TODO: reorder for better parallelism */
1015:   PetscMalloc1(user->nt*user->mx*user->mx,&user->uxi_scatter);
1016:   PetscMalloc1(user->nt*user->mx*user->mx,&user->uyi_scatter);
1017:   PetscMalloc1(user->nt*user->mx*user->mx,&user->ux_scatter);
1018:   PetscMalloc1(user->nt*user->mx*user->mx,&user->uy_scatter);
1019:   PetscMalloc1(2*user->nt*user->mx*user->mx,&user->ui_scatter);
1020:   VecCreate(PETSC_COMM_WORLD,&uxi);
1021:   VecCreate(PETSC_COMM_WORLD,&ui);
1022:   VecSetSizes(uxi,PETSC_DECIDE,user->mx*user->mx);
1023:   VecSetSizes(ui,PETSC_DECIDE,2*user->mx*user->mx);
1024:   VecSetFromOptions(uxi);
1025:   VecSetFromOptions(ui);
1026:   VecDuplicateVecs(uxi,user->nt,&user->uxi);
1027:   VecDuplicateVecs(uxi,user->nt,&user->uyi);
1028:   VecDuplicateVecs(uxi,user->nt,&user->uxiwork);
1029:   VecDuplicateVecs(uxi,user->nt,&user->uyiwork);
1030:   VecDuplicateVecs(ui,user->nt,&user->ui);
1031:   VecDuplicateVecs(ui,user->nt,&user->uiwork);
1032:   for (i=0; i<user->nt; i++){
1033:     VecGetOwnershipRange(user->uxi[i],&lo,&hi);
1034:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1035:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u);
1036:     VecScatterCreate(user->u,is_from_u,user->uxi[i],is_to_uxi,&user->uxi_scatter[i]);

1038:     ISDestroy(&is_to_uxi);
1039:     ISDestroy(&is_from_u);

1041:     VecGetOwnershipRange(user->uyi[i],&lo,&hi);
1042:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi);
1043:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+(2*i+1)*user->mx*user->mx,1,&is_from_u);
1044:     VecScatterCreate(user->u,is_from_u,user->uyi[i],is_to_uyi,&user->uyi_scatter[i]);

1046:     ISDestroy(&is_to_uyi);
1047:     ISDestroy(&is_from_u);

1049:     VecGetOwnershipRange(user->uxi[i],&lo,&hi);
1050:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1051:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_from_u);
1052:     VecScatterCreate(user->ui[i],is_from_u,user->uxi[i],is_to_uxi,&user->ux_scatter[i]);

1054:     ISDestroy(&is_to_uxi);
1055:     ISDestroy(&is_from_u);

1057:     VecGetOwnershipRange(user->uyi[i],&lo,&hi);
1058:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi);
1059:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+user->mx*user->mx,1,&is_from_u);
1060:     VecScatterCreate(user->ui[i],is_from_u,user->uyi[i],is_to_uyi,&user->uy_scatter[i]);

1062:     ISDestroy(&is_to_uyi);
1063:     ISDestroy(&is_from_u);

1065:     VecGetOwnershipRange(user->ui[i],&lo,&hi);
1066:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1067:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u);
1068:     VecScatterCreate(user->u,is_from_u,user->ui[i],is_to_uxi,&user->ui_scatter[i]);

1070:     ISDestroy(&is_to_uxi);
1071:     ISDestroy(&is_from_u);
1072:   }

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

1081:   /* Compute true velocity field utrue */
1082:   VecDuplicate(user->u,&user->utrue);
1083:   for (i=0; i<user->nt; i++){
1084:     VecCopy(YY,user->uxi[i]);
1085:     VecScale(user->uxi[i],150.0*i*user->ht);
1086:     VecCopy(XX,user->uyi[i]);
1087:     VecShift(user->uyi[i],-10.0);
1088:     VecScale(user->uyi[i],15.0*i*user->ht);
1089:   }
1090:   Gather_uxi_uyi(user->utrue,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

1092:   /* Initial guess and reference model */
1093:   VecDuplicate(user->utrue,&user->ur);
1094:   for (i=0; i<user->nt; i++){
1095:     VecCopy(XX,user->uxi[i]);
1096:     VecShift(user->uxi[i],i*user->ht);
1097:     VecCopy(YY,user->uyi[i]);
1098:     VecShift(user->uyi[i],-i*user->ht);
1099:   }
1100:   Gather_uxi_uyi(user->ur,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

1102:   /* Generate regularization matrix L */
1103:   MatCreate(PETSC_COMM_WORLD,&user->LT);
1104:   MatSetSizes(user->LT,PETSC_DECIDE,PETSC_DECIDE,2*n*user->nt,n*user->nt);
1105:   MatSetFromOptions(user->LT);
1106:   MatMPIAIJSetPreallocation(user->LT,1,NULL,1,NULL);
1107:   MatSeqAIJSetPreallocation(user->LT,1,NULL);
1108:   MatGetOwnershipRange(user->LT,&istart,&iend);

1110:   for (i=istart; i<iend; i++){
1111:     iblock = (i+n) / (2*n);
1112:     j = i - iblock*n;
1113:     MatSetValues(user->LT,1,&i,1,&j,&user->gamma,INSERT_VALUES);
1114:   }

1116:   MatAssemblyBegin(user->LT,MAT_FINAL_ASSEMBLY);
1117:   MatAssemblyEnd(user->LT,MAT_FINAL_ASSEMBLY);

1119:   MatTranspose(user->LT,MAT_INITIAL_MATRIX,&user->L);

1121:   /* Build work vectors and matrices */
1122:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1123:   VecSetType(user->lwork,VECMPI);
1124:   VecSetSizes(user->lwork,PETSC_DECIDE,user->m);
1125:   VecSetFromOptions(user->lwork);

1127:   MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);

1129:   VecDuplicate(user->y,&user->ywork);
1130:   VecDuplicate(user->u,&user->uwork);
1131:   VecDuplicate(user->u,&user->vwork);
1132:   VecDuplicate(user->u,&user->js_diag);
1133:   VecDuplicate(user->c,&user->cwork);

1135:   /* Create matrix-free shell user->Js for computing A*x */
1136:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->Js);
1137:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1138:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
1139:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
1140:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

1142:   /* Diagonal blocks of user->Js */
1143:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlock);
1144:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
1145:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMultTranspose);

1147:   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1148:      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1149:      This is an SOR preconditioner for user->JsBlock. */
1150:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlockPrec);
1151:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult);
1152:   MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMultTranspose);

1154:   /* Create a matrix-free shell user->Jd for computing B*x */
1155:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->n-user->m,user,&user->Jd);
1156:   MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult);
1157:   MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose);

1159:   /* User-defined routines for computing user->Js\x and user->Js^T\x*/
1160:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsInv);
1161:   MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult);
1162:   MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult);

1164:   /* Build matrices for SOR preconditioner */
1165:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);
1166:   MatShift(user->Divxy[0],0.0); /*  Force C[i] and Divxy[0] to share same nonzero pattern */
1167:   MatAXPY(user->Divxy[0],0.0,user->Divxy[1],DIFFERENT_NONZERO_PATTERN);
1168:   PetscMalloc1(5*n,&user->C);
1169:   PetscMalloc1(2*n,&user->Cwork);
1170:   for (i=0; i<user->nt; i++){
1171:     MatDuplicate(user->Divxy[0],MAT_COPY_VALUES,&user->C[i]);
1172:     MatDuplicate(user->Divxy[1],MAT_COPY_VALUES,&user->Cwork[i]);

1174:     MatDiagonalScale(user->C[i],NULL,user->uxi[i]);
1175:     MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]);
1176:     MatAXPY(user->C[i],1.0,user->Cwork[i],SUBSET_NONZERO_PATTERN);
1177:     MatScale(user->C[i],user->ht);
1178:     MatShift(user->C[i],1.0);
1179:   }

1181:   /* Solver options and tolerances */
1182:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1183:   KSPSetType(user->solver,KSPGMRES);
1184:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
1185:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE); /*  TODO: why is true slower? */
1186:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1187:   /* KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500); */
1188:   KSPGetPC(user->solver,&user->prec);
1189:   PCSetType(user->prec,PCSHELL);

1191:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1192:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMultTranspose);
1193:   PCShellSetContext(user->prec,user);

1195:   /* Compute true state function yt given ut */
1196:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1197:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1198:   VecSetFromOptions(user->ytrue);
1199:   user->c_formed = PETSC_TRUE;
1200:   VecCopy(user->utrue,user->u); /*  Set u=utrue temporarily for StateMatInv */
1201:   VecSet(user->ytrue,0.0); /*  Initial guess */
1202:   StateMatInvMult(user->Js,user->q,user->ytrue);
1203:   VecCopy(user->ur,user->u); /*  Reset u=ur */

1205:   /* Initial guess y0 for state given u0 */
1206:   user->lcl->solve_type = LCL_FORWARD1;
1207:   StateMatInvMult(user->Js,user->q,user->y);

1209:   /* Data discretization */
1210:   MatCreate(PETSC_COMM_WORLD,&user->Q);
1211:   MatSetSizes(user->Q,PETSC_DECIDE,PETSC_DECIDE,user->mx*user->mx,user->m);
1212:   MatSetFromOptions(user->Q);
1213:   MatMPIAIJSetPreallocation(user->Q,0,NULL,1,NULL);
1214:   MatSeqAIJSetPreallocation(user->Q,1,NULL);

1216:   MatGetOwnershipRange(user->Q,&istart,&iend);

1218:   for (i=istart; i<iend; i++){
1219:     j = i + user->m - user->mx*user->mx;
1220:     MatSetValues(user->Q,1,&i,1,&j,&one,INSERT_VALUES);
1221:   }

1223:   MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY);
1224:   MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY);

1226:   MatTranspose(user->Q,MAT_INITIAL_MATRIX,&user->QT);

1228:   VecDestroy(&XX);
1229:   VecDestroy(&YY);
1230:   VecDestroy(&XXwork);
1231:   VecDestroy(&YYwork);
1232:   VecDestroy(&yi);
1233:   VecDestroy(&uxi);
1234:   VecDestroy(&ui);
1235:   VecDestroy(&bc);

1237:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1238:   KSPSetFromOptions(user->solver);
1239:   return(0);
1240: }

1244: PetscErrorCode HyperbolicDestroy(AppCtx *user)
1245: {
1247:   PetscInt       i;

1250:   MatDestroy(&user->Q);
1251:   MatDestroy(&user->QT);
1252:   MatDestroy(&user->Div);
1253:   MatDestroy(&user->Divwork);
1254:   MatDestroy(&user->Grad);
1255:   MatDestroy(&user->L);
1256:   MatDestroy(&user->LT);
1257:   KSPDestroy(&user->solver);
1258:   MatDestroy(&user->Js);
1259:   MatDestroy(&user->Jd);
1260:   MatDestroy(&user->JsBlockPrec);
1261:   MatDestroy(&user->JsInv);
1262:   MatDestroy(&user->JsBlock);
1263:   MatDestroy(&user->Divxy[0]);
1264:   MatDestroy(&user->Divxy[1]);
1265:   MatDestroy(&user->Gradxy[0]);
1266:   MatDestroy(&user->Gradxy[1]);
1267:   MatDestroy(&user->M);
1268:   for (i=0; i<user->nt; i++){
1269:     MatDestroy(&user->C[i]);
1270:     MatDestroy(&user->Cwork[i]);
1271:   }
1272:   PetscFree(user->C);
1273:   PetscFree(user->Cwork);
1274:   VecDestroy(&user->u);
1275:   VecDestroy(&user->uwork);
1276:   VecDestroy(&user->vwork);
1277:   VecDestroy(&user->utrue);
1278:   VecDestroy(&user->y);
1279:   VecDestroy(&user->ywork);
1280:   VecDestroy(&user->ytrue);
1281:   VecDestroyVecs(user->nt,&user->yi);
1282:   VecDestroyVecs(user->nt,&user->yiwork);
1283:   VecDestroyVecs(user->nt,&user->ziwork);
1284:   VecDestroyVecs(user->nt,&user->uxi);
1285:   VecDestroyVecs(user->nt,&user->uyi);
1286:   VecDestroyVecs(user->nt,&user->uxiwork);
1287:   VecDestroyVecs(user->nt,&user->uyiwork);
1288:   VecDestroyVecs(user->nt,&user->ui);
1289:   VecDestroyVecs(user->nt,&user->uiwork);
1290:   VecDestroy(&user->c);
1291:   VecDestroy(&user->cwork);
1292:   VecDestroy(&user->ur);
1293:   VecDestroy(&user->q);
1294:   VecDestroy(&user->d);
1295:   VecDestroy(&user->dwork);
1296:   VecDestroy(&user->lwork);
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->uxi_scatter[i]);
1304:     VecScatterDestroy(&user->uyi_scatter[i]);
1305:     VecScatterDestroy(&user->ux_scatter[i]);
1306:     VecScatterDestroy(&user->uy_scatter[i]);
1307:     VecScatterDestroy(&user->ui_scatter[i]);
1308:     VecScatterDestroy(&user->yi_scatter[i]);
1309:   }
1310:   PetscFree(user->uxi_scatter);
1311:   PetscFree(user->uyi_scatter);
1312:   PetscFree(user->ux_scatter);
1313:   PetscFree(user->uy_scatter);
1314:   PetscFree(user->ui_scatter);
1315:   PetscFree(user->yi_scatter);
1316:   return(0);
1317: }

1321: PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1322: {
1324:   Vec            X;
1325:   PetscReal      unorm,ynorm;
1326:   AppCtx         *user = (AppCtx*)ptr;

1329:   TaoGetSolutionVector(tao,&X);
1330:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1331:   VecAXPY(user->ywork,-1.0,user->ytrue);
1332:   VecAXPY(user->uwork,-1.0,user->utrue);
1333:   VecNorm(user->uwork,NORM_2,&unorm);
1334:   VecNorm(user->ywork,NORM_2,&ynorm);
1335:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1336:   return(0);
1337: }