Actual source code: hyperbolic.c

petsc-master 2016-02-06
Report Typos and Errors
  1: #include <petsctao.h>

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

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

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

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

 51:   Vec d;
 52:   Vec dwork;

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

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

 64:   Vec js_diag;

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

 69:   Vec lwork;

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

 75:   PetscInt ksp_its;
 76:   PetscInt ksp_its_initial;
 77: } AppCtx;


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

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

103: PetscErrorCode DesignMatMult(Mat,Vec,Vec);
104: PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);

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

111: static  char help[]="";

115: int main(int argc, char **argv)
116: {
117:   PetscErrorCode     ierr;
118:   Vec                x,x0;
119:   Tao                tao;
120:   TaoConvergedReason reason;
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
129: 
130:   PetscInitialize(&argc, &argv, (char*)0,help);

132:   user.mx = 32;
133:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,NULL,NULL);
134:   PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL);
135:   user.nt = 16;
136:   PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL);
137:   user.ndata = 64;
138:   PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL);
139:   user.alpha = 10.0;
140:   PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL);
141:   user.T = 1.0/32.0;
142:   PetscOptionsReal("-Tfinal","Final time","",user.T,&user.T,NULL);

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

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

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

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

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

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

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

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

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

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:   PetscOptionsEnd();
213:   PetscLogStageRegister("Trials",&stages[0]);
214:   PetscLogStagePush(stages[0]);
215:   user.ksp_its_initial = user.ksp_its;
216:   ksp_old = user.ksp_its;
217:   for (i=0; i<ntests; i++){
218:     TaoSolve(tao);
219:     PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its-ksp_old);
220:     VecCopy(x0,x);
221:     TaoSetInitialVector(tao,x);
222:   }
223:   PetscLogStagePop();
224:   PetscBarrier((PetscObject)x);
225:   PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ");
226:   PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial);
227:   PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %D trial(s): ",ntests);
228:   PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its);
229:   PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: ");
230:   PetscPrintf(PETSC_COMM_WORLD,"%D\n",(user.ksp_its-user.ksp_its_initial)/ntests);

232:   TaoGetConvergedReason(tao,&reason);

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

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

262:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
263:   MatMult(user->Q,user->y,user->dwork);
264:   VecAXPY(user->dwork,-1.0,user->d);
265:   VecDot(user->dwork,user->dwork,&d1);

267:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
268:   VecPointwiseMult(user->uwork,user->uwork,user->uwork);
269:   MatMult(user->L,user->uwork,user->lwork);
270:   VecDot(user->y,user->lwork,&d2);
271:   *f = 0.5 * (d1 + user->alpha*d2);
272:   return(0);
273: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

565:   if (Y == user->ytrue) {
566:     /* First solve is done using true solution to set up problem */
567:     KSPSetTolerances(user->solver,1e-4,1e-20,PETSC_DEFAULT,PETSC_DEFAULT);
568:   } else {
569:     KSPSetTolerances(user->solver,PETSC_DEFAULT,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:   PetscInt       its,i;

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

604:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);
605:   Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt);
606:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

608:   i = user->nt - 1;
609:   user->block_index = i;
610:   KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]);

612:   KSPGetIterationNumber(user->solver,&its);
613:   user->ksp_its = user->ksp_its + its;

615:   for (i=user->nt-2; i>=0; i--){
616:     MatMult(user->M,user->yiwork[i+1],user->ziwork[i+1]);
617:     VecAXPY(user->yi[i],1.0,user->ziwork[i+1]);
618:     user->block_index = i;
619:     KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]);

621:     KSPGetIterationNumber(user->solver,&its);
622:     user->ksp_its = user->ksp_its + its;
623:   }
624:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
625:   return(0);
626: }

630: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
631: {
633:   AppCtx         *user;

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

638:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
639:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
640:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
641:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
642:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
643:   return(0);
644: }

648: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
649: {
651:   AppCtx         *user;

654:   MatShellGetContext(J_shell,(void**)&user);
655:    VecCopy(user->js_diag,X);
656:   return(0);
657: }

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

674:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
675:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
676:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

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

681:   for (i=1; i<user->nt; i++){
682:     user->block_index = i;
683:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
684:     MatMult(user->M,user->yi[i-1],user->ziwork[i-1]);
685:     VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]);
686:   }

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

691:   return(0);
692: }


697: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
698: {

702:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
703:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
704:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
705:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
706:   return(0);
707: }

711: PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
712: {
714:   PetscInt       i;

717:   for (i=0; i<nt; i++){
718:     VecScatterBegin(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD);
719:     VecScatterEnd(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD);
720:     VecScatterBegin(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD);
721:     VecScatterEnd(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD);
722:   }
723:   return(0);
724: }

728: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
729: {

733:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
734:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
735:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
736:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
737:   return(0);
738: }

742: PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
743: {
745:   PetscInt       i;

748:   for (i=0; i<nt; i++){
749:     VecScatterBegin(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE);
750:     VecScatterEnd(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE);
751:     VecScatterBegin(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE);
752:     VecScatterEnd(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE);
753:   }
754:   return(0);
755: }

759: PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
760: {
762:   PetscInt       i;

765:   for (i=0; i<nt; i++){
766:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
767:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
768:   }
769:   return(0);
770: }

774: PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
775: {
777:   PetscInt       i;

780:   for (i=0; i<nt; i++){
781:     VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
782:     VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
783:   }
784:   return(0);
785: }

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

800:   user->jformed = PETSC_FALSE;
801:   user->c_formed = PETSC_FALSE;

803:   user->ksp_its = 0;
804:   user->ksp_its_initial = 0;

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

808:   h = 1.0/user->mx;
809:   hinv = user->mx;
810:   neg_hinv = -hinv;
811:   half_hinv = hinv / 2.0;
812:   neg_half_hinv = neg_hinv / 2.0;

814:   /* Generate Grad matrix */
815:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
816:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,2*n,n);
817:   MatSetFromOptions(user->Grad);
818:   MatMPIAIJSetPreallocation(user->Grad,3,NULL,3,NULL);
819:   MatSeqAIJSetPreallocation(user->Grad,3,NULL);
820:   MatGetOwnershipRange(user->Grad,&istart,&iend);

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

838:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
839:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);

841:   MatCreate(PETSC_COMM_WORLD,&user->Gradxy[0]);
842:   MatSetSizes(user->Gradxy[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
843:   MatSetFromOptions(user->Gradxy[0]);
844:   MatMPIAIJSetPreallocation(user->Gradxy[0],3,NULL,3,NULL);
845:   MatSeqAIJSetPreallocation(user->Gradxy[0],3,NULL);
846:   MatGetOwnershipRange(user->Gradxy[0],&istart,&iend);

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

859:   MatCreate(PETSC_COMM_WORLD,&user->Gradxy[1]);
860:   MatSetSizes(user->Gradxy[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
861:   MatSetFromOptions(user->Gradxy[1]);
862:   MatMPIAIJSetPreallocation(user->Gradxy[1],3,NULL,3,NULL);
863:   MatSeqAIJSetPreallocation(user->Gradxy[1],3,NULL);
864:   MatGetOwnershipRange(user->Gradxy[1],&istart,&iend);

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

876:   /* Generate Div matrix */
877:   MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div);
878:   MatTranspose(user->Gradxy[0],MAT_INITIAL_MATRIX,&user->Divxy[0]);
879:   MatTranspose(user->Gradxy[1],MAT_INITIAL_MATRIX,&user->Divxy[1]);

881:   /* Off-diagonal averaging matrix */
882:   MatCreate(PETSC_COMM_WORLD,&user->M);
883:   MatSetSizes(user->M,PETSC_DECIDE,PETSC_DECIDE,n,n);
884:   MatSetFromOptions(user->M);
885:   MatMPIAIJSetPreallocation(user->M,4,NULL,4,NULL);
886:   MatSeqAIJSetPreallocation(user->M,4,NULL);
887:   MatGetOwnershipRange(user->M,&istart,&iend);

889:   for (i=istart; i<iend; i++){
890:     /* kron(Id,Av) */
891:     iblock = i / user->mx;
892:     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
893:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
894:     j = iblock*user->mx + ((i+1) % user->mx);
895:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);

897:     /* kron(Av,Id) */
898:     j = (i + user->mx) % n;
899:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
900:     j = (i + n - user->mx) % n;
901:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
902:   }
903:   MatAssemblyBegin(user->M,MAT_FINAL_ASSEMBLY);
904:   MatAssemblyEnd(user->M,MAT_FINAL_ASSEMBLY);

906:   /* Generate 2D grid */
907:   VecCreate(PETSC_COMM_WORLD,&XX);
908:   VecCreate(PETSC_COMM_WORLD,&user->q);
909:   VecSetSizes(XX,PETSC_DECIDE,n);
910:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
911:   VecSetFromOptions(XX);
912:   VecSetFromOptions(user->q);

914:   VecDuplicate(XX,&YY);
915:   VecDuplicate(XX,&XXwork);
916:   VecDuplicate(XX,&YYwork);
917:   VecDuplicate(XX,&user->d);
918:   VecDuplicate(XX,&user->dwork);

920:   VecGetOwnershipRange(XX,&istart,&iend);
921:   for (linear_index=istart; linear_index<iend; linear_index++){
922:     i = linear_index % user->mx;
923:     j = (linear_index-i)/user->mx;
924:     vx = h*(i+0.5);
925:     vy = h*(j+0.5);
926:     VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
927:     VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
928:   }

930:   VecAssemblyBegin(XX);
931:   VecAssemblyEnd(XX);
932:   VecAssemblyBegin(YY);
933:   VecAssemblyEnd(YY);

935:   /* Compute final density function yT
936:      yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2))
937:      yT = yT / (h^2*sum(yT)) */
938:   VecCopy(XX,XXwork);
939:   VecCopy(YY,YYwork);

941:   VecShift(XXwork,-0.25);
942:   VecShift(YYwork,-0.25);

944:   VecPointwiseMult(XXwork,XXwork,XXwork);
945:   VecPointwiseMult(YYwork,YYwork,YYwork);

947:   VecCopy(XXwork,user->dwork);
948:   VecAXPY(user->dwork,1.0,YYwork);
949:   VecScale(user->dwork,-30.0);
950:   VecExp(user->dwork);
951:   VecCopy(user->dwork,user->d);

953:   VecCopy(XX,XXwork);
954:   VecCopy(YY,YYwork);

956:   VecShift(XXwork,-0.75);
957:   VecShift(YYwork,-0.75);

959:   VecPointwiseMult(XXwork,XXwork,XXwork);
960:   VecPointwiseMult(YYwork,YYwork,YYwork);

962:   VecCopy(XXwork,user->dwork);
963:   VecAXPY(user->dwork,1.0,YYwork);
964:   VecScale(user->dwork,-30.0);
965:   VecExp(user->dwork);

967:   VecAXPY(user->d,1.0,user->dwork);
968:   VecShift(user->d,1.0);
969:   VecSum(user->d,&sum);
970:   VecScale(user->d,1.0/(h*h*sum));

972:   /* Initial conditions of forward problem */
973:   VecDuplicate(XX,&bc);
974:   VecCopy(XX,XXwork);
975:   VecCopy(YY,YYwork);

977:   VecShift(XXwork,-0.5);
978:   VecShift(YYwork,-0.5);

980:   VecPointwiseMult(XXwork,XXwork,XXwork);
981:   VecPointwiseMult(YYwork,YYwork,YYwork);

983:   VecWAXPY(bc,1.0,XXwork,YYwork);
984:   VecScale(bc,-50.0);
985:   VecExp(bc);
986:   VecShift(bc,1.0);
987:   VecSum(bc,&sum);
988:   VecScale(bc,1.0/(h*h*sum));

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

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

1033:     ISDestroy(&is_to_uxi);
1034:     ISDestroy(&is_from_u);

1036:     VecGetOwnershipRange(user->uyi[i],&lo,&hi);
1037:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi);
1038:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+(2*i+1)*user->mx*user->mx,1,&is_from_u);
1039:     VecScatterCreate(user->u,is_from_u,user->uyi[i],is_to_uyi,&user->uyi_scatter[i]);

1041:     ISDestroy(&is_to_uyi);
1042:     ISDestroy(&is_from_u);

1044:     VecGetOwnershipRange(user->uxi[i],&lo,&hi);
1045:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1046:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_from_u);
1047:     VecScatterCreate(user->ui[i],is_from_u,user->uxi[i],is_to_uxi,&user->ux_scatter[i]);

1049:     ISDestroy(&is_to_uxi);
1050:     ISDestroy(&is_from_u);

1052:     VecGetOwnershipRange(user->uyi[i],&lo,&hi);
1053:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi);
1054:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+user->mx*user->mx,1,&is_from_u);
1055:     VecScatterCreate(user->ui[i],is_from_u,user->uyi[i],is_to_uyi,&user->uy_scatter[i]);

1057:     ISDestroy(&is_to_uyi);
1058:     ISDestroy(&is_from_u);

1060:     VecGetOwnershipRange(user->ui[i],&lo,&hi);
1061:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1062:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u);
1063:     VecScatterCreate(user->u,is_from_u,user->ui[i],is_to_uxi,&user->ui_scatter[i]);

1065:     ISDestroy(&is_to_uxi);
1066:     ISDestroy(&is_from_u);
1067:   }

1069:   /* RHS of forward problem */
1070:   MatMult(user->M,bc,user->yiwork[0]);
1071:   for (i=1; i<user->nt; i++){
1072:     VecSet(user->yiwork[i],0.0);
1073:   }
1074:   Gather_yi(user->q,user->yiwork,user->yi_scatter,user->nt);

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

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

1097:   /* Generate regularization matrix L */
1098:   MatCreate(PETSC_COMM_WORLD,&user->LT);
1099:   MatSetSizes(user->LT,PETSC_DECIDE,PETSC_DECIDE,2*n*user->nt,n*user->nt);
1100:   MatSetFromOptions(user->LT);
1101:   MatMPIAIJSetPreallocation(user->LT,1,NULL,1,NULL);
1102:   MatSeqAIJSetPreallocation(user->LT,1,NULL);
1103:   MatGetOwnershipRange(user->LT,&istart,&iend);

1105:   for (i=istart; i<iend; i++){
1106:     iblock = (i+n) / (2*n);
1107:     j = i - iblock*n;
1108:     MatSetValues(user->LT,1,&i,1,&j,&user->gamma,INSERT_VALUES);
1109:   }

1111:   MatAssemblyBegin(user->LT,MAT_FINAL_ASSEMBLY);
1112:   MatAssemblyEnd(user->LT,MAT_FINAL_ASSEMBLY);

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

1116:   /* Build work vectors and matrices */
1117:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1118:   VecSetType(user->lwork,VECMPI);
1119:   VecSetSizes(user->lwork,PETSC_DECIDE,user->m);
1120:   VecSetFromOptions(user->lwork);

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

1124:   VecDuplicate(user->y,&user->ywork);
1125:   VecDuplicate(user->u,&user->uwork);
1126:   VecDuplicate(user->u,&user->vwork);
1127:   VecDuplicate(user->u,&user->js_diag);
1128:   VecDuplicate(user->c,&user->cwork);

1130:   /* Create matrix-free shell user->Js for computing A*x */
1131:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->Js);
1132:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1133:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
1134:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
1135:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

1137:   /* Diagonal blocks of user->Js */
1138:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlock);
1139:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
1140:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMultTranspose);

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

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

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

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

1169:     MatDiagonalScale(user->C[i],NULL,user->uxi[i]);
1170:     MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]);
1171:     MatAXPY(user->C[i],1.0,user->Cwork[i],SUBSET_NONZERO_PATTERN);
1172:     MatScale(user->C[i],user->ht);
1173:     MatShift(user->C[i],1.0);
1174:   }

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

1186:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1187:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMultTranspose);
1188:   PCShellSetContext(user->prec,user);

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

1200:   /* Initial guess y0 for state given u0 */
1201:   StateMatInvMult(user->Js,user->q,user->y);

1203:   /* Data discretization */
1204:   MatCreate(PETSC_COMM_WORLD,&user->Q);
1205:   MatSetSizes(user->Q,PETSC_DECIDE,PETSC_DECIDE,user->mx*user->mx,user->m);
1206:   MatSetFromOptions(user->Q);
1207:   MatMPIAIJSetPreallocation(user->Q,0,NULL,1,NULL);
1208:   MatSeqAIJSetPreallocation(user->Q,1,NULL);

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

1212:   for (i=istart; i<iend; i++){
1213:     j = i + user->m - user->mx*user->mx;
1214:     MatSetValues(user->Q,1,&i,1,&j,&one,INSERT_VALUES);
1215:   }

1217:   MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY);
1218:   MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY);

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

1222:   VecDestroy(&XX);
1223:   VecDestroy(&YY);
1224:   VecDestroy(&XXwork);
1225:   VecDestroy(&YYwork);
1226:   VecDestroy(&yi);
1227:   VecDestroy(&uxi);
1228:   VecDestroy(&ui);
1229:   VecDestroy(&bc);

1231:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1232:   KSPSetFromOptions(user->solver);
1233:   return(0);
1234: }

1238: PetscErrorCode HyperbolicDestroy(AppCtx *user)
1239: {
1241:   PetscInt       i;

1244:   MatDestroy(&user->Q);
1245:   MatDestroy(&user->QT);
1246:   MatDestroy(&user->Div);
1247:   MatDestroy(&user->Divwork);
1248:   MatDestroy(&user->Grad);
1249:   MatDestroy(&user->L);
1250:   MatDestroy(&user->LT);
1251:   KSPDestroy(&user->solver);
1252:   MatDestroy(&user->Js);
1253:   MatDestroy(&user->Jd);
1254:   MatDestroy(&user->JsBlockPrec);
1255:   MatDestroy(&user->JsInv);
1256:   MatDestroy(&user->JsBlock);
1257:   MatDestroy(&user->Divxy[0]);
1258:   MatDestroy(&user->Divxy[1]);
1259:   MatDestroy(&user->Gradxy[0]);
1260:   MatDestroy(&user->Gradxy[1]);
1261:   MatDestroy(&user->M);
1262:   for (i=0; i<user->nt; i++){
1263:     MatDestroy(&user->C[i]);
1264:     MatDestroy(&user->Cwork[i]);
1265:   }
1266:   PetscFree(user->C);
1267:   PetscFree(user->Cwork);
1268:   VecDestroy(&user->u);
1269:   VecDestroy(&user->uwork);
1270:   VecDestroy(&user->vwork);
1271:   VecDestroy(&user->utrue);
1272:   VecDestroy(&user->y);
1273:   VecDestroy(&user->ywork);
1274:   VecDestroy(&user->ytrue);
1275:   VecDestroyVecs(user->nt,&user->yi);
1276:   VecDestroyVecs(user->nt,&user->yiwork);
1277:   VecDestroyVecs(user->nt,&user->ziwork);
1278:   VecDestroyVecs(user->nt,&user->uxi);
1279:   VecDestroyVecs(user->nt,&user->uyi);
1280:   VecDestroyVecs(user->nt,&user->uxiwork);
1281:   VecDestroyVecs(user->nt,&user->uyiwork);
1282:   VecDestroyVecs(user->nt,&user->ui);
1283:   VecDestroyVecs(user->nt,&user->uiwork);
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->js_diag);
1292:   ISDestroy(&user->s_is);
1293:   ISDestroy(&user->d_is);
1294:   VecScatterDestroy(&user->state_scatter);
1295:   VecScatterDestroy(&user->design_scatter);
1296:   for (i=0; i<user->nt; i++){
1297:     VecScatterDestroy(&user->uxi_scatter[i]);
1298:     VecScatterDestroy(&user->uyi_scatter[i]);
1299:     VecScatterDestroy(&user->ux_scatter[i]);
1300:     VecScatterDestroy(&user->uy_scatter[i]);
1301:     VecScatterDestroy(&user->ui_scatter[i]);
1302:     VecScatterDestroy(&user->yi_scatter[i]);
1303:   }
1304:   PetscFree(user->uxi_scatter);
1305:   PetscFree(user->uyi_scatter);
1306:   PetscFree(user->ux_scatter);
1307:   PetscFree(user->uy_scatter);
1308:   PetscFree(user->ui_scatter);
1309:   PetscFree(user->yi_scatter);
1310:   return(0);
1311: }

1315: PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1316: {
1318:   Vec            X;
1319:   PetscReal      unorm,ynorm;
1320:   AppCtx         *user = (AppCtx*)ptr;

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