Actual source code: hyperbolic.c

petsc-master 2016-06-28
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: 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:   AppCtx             user;
121:   IS                 is_allstate,is_alldesign;
122:   PetscInt           lo,hi,hi2,lo2,ksp_old;
123:   PetscInt           ntests = 1;
124:   PetscInt           i;
125: #if defined(PETSC_USE_LOG)
126:   PetscLogStage      stages[1];
127: #endif

129:   PetscInitialize(&argc, &argv, (char*)0,help);if (ierr) return ierr;
130:   user.mx = 32;
131:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,NULL,NULL);
132:   PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL);
133:   user.nt = 16;
134:   PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL);
135:   user.ndata = 64;
136:   PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL);
137:   user.alpha = 10.0;
138:   PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL);
139:   user.T = 1.0/32.0;
140:   PetscOptionsReal("-Tfinal","Final time","",user.T,&user.T,NULL);

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

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

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

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

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

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

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

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

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

188:   /* Set up initial vectors and matrices */
189:   HyperbolicInitialize(&user);

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

195:   /* Set solution vector with an initial guess */
196:   TaoSetInitialVector(tao,x);
197:   TaoSetObjectiveRoutine(tao, FormFunction, (void *)&user);
198:   TaoSetGradientRoutine(tao, FormGradient, (void *)&user);
199:   TaoSetConstraintsRoutine(tao, user.c, FormConstraints, (void *)&user);
200:   TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, (void *)&user);
201:   TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, (void *)&user);
202:   TaoSetFromOptions(tao);
203:   TaoSetStateDesignIS(tao,user.s_is,user.d_is);

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

227:   TaoDestroy(&tao);
228:   VecDestroy(&x);
229:   VecDestroy(&x0);
230:   HyperbolicDestroy(&user);
231:   PetscFinalize();
232:   return ierr;
233: }
234: /* ------------------------------------------------------------------- */
237: /*
238:    dwork = Qy - d
239:    lwork = L*(u-ur).^2
240:    f = 1/2 * (dwork.dork + alpha*y.lwork)
241: */
242: PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
243: {
245:   PetscReal      d1=0,d2=0;
246:   AppCtx         *user = (AppCtx*)ptr;

249:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
250:   MatMult(user->Q,user->y,user->dwork);
251:   VecAXPY(user->dwork,-1.0,user->d);
252:   VecDot(user->dwork,user->dwork,&d1);

254:   VecWAXPY(user->uwork,-1.0,user->ur,user->u);
255:   VecPointwiseMult(user->uwork,user->uwork,user->uwork);
256:   MatMult(user->L,user->uwork,user->lwork);
257:   VecDot(user->y,user->lwork,&d2);
258:   *f = 0.5 * (d1 + user->alpha*d2);
259:   return(0);
260: }

262: /* ------------------------------------------------------------------- */
265: /*
266:     state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2
267:     design: g_d = alpha*(L'y).*(u-ur)
268: */
269: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
270: {
272:   AppCtx         *user = (AppCtx*)ptr;

275:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
276:   MatMult(user->Q,user->y,user->dwork);
277:   VecAXPY(user->dwork,-1.0,user->d);

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

281:   MatMult(user->LT,user->y,user->uwork);
282:   VecWAXPY(user->vwork,-1.0,user->ur,user->u);
283:   VecPointwiseMult(user->uwork,user->vwork,user->uwork);
284:   VecScale(user->uwork,user->alpha);

286:   VecPointwiseMult(user->vwork,user->vwork,user->vwork);
287:   MatMult(user->L,user->vwork,user->lwork);
288:   VecAXPY(user->ywork,0.5*user->alpha,user->lwork);

290:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
291:   return(0);
292: }

296: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
297: {
299:   PetscReal      d1,d2;
300:   AppCtx         *user = (AppCtx*)ptr;

303:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
304:   MatMult(user->Q,user->y,user->dwork);
305:   VecAXPY(user->dwork,-1.0,user->d);

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

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

311:   MatMult(user->LT,user->y,user->uwork);
312:   VecWAXPY(user->vwork,-1.0,user->ur,user->u);
313:   VecPointwiseMult(user->uwork,user->vwork,user->uwork);
314:   VecScale(user->uwork,user->alpha);

316:   VecPointwiseMult(user->vwork,user->vwork,user->vwork);
317:   MatMult(user->L,user->vwork,user->lwork);
318:   VecAXPY(user->ywork,0.5*user->alpha,user->lwork);

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

322:   *f = 0.5 * (d1 + user->alpha*d2);
323:   Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
324:   return(0);
325: }

327: /* ------------------------------------------------------------------- */
330: /* A
331: MatShell object
332: */
333: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
334: {
336:   PetscInt       i;
337:   AppCtx         *user = (AppCtx*)ptr;

340:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
341:   Scatter_yi(user->u,user->ui,user->ui_scatter,user->nt);
342:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);
343:   for (i=0; i<user->nt; i++){
344:     MatCopy(user->Divxy[0],user->C[i],SAME_NONZERO_PATTERN);
345:     MatCopy(user->Divxy[1],user->Cwork[i],SAME_NONZERO_PATTERN);

347:     MatDiagonalScale(user->C[i],NULL,user->uxi[i]);
348:     MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]);
349:     MatAXPY(user->C[i],1.0,user->Cwork[i],SUBSET_NONZERO_PATTERN);
350:     MatScale(user->C[i],user->ht);
351:     MatShift(user->C[i],1.0);
352:   }
353:   return(0);
354: }

356: /* ------------------------------------------------------------------- */
359: /* B */
360: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
361: {
363:   AppCtx         *user = (AppCtx*)ptr;

366:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
367:   return(0);
368: }

372: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
373: {
375:   PetscInt       i;
376:   AppCtx         *user;

379:   MatShellGetContext(J_shell,(void**)&user);
380:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);
381:   user->block_index = 0;
382:   MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);

384:   for (i=1; i<user->nt; i++){
385:     user->block_index = i;
386:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
387:     MatMult(user->M,user->yi[i-1],user->ziwork[i-1]);
388:     VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]);
389:   }
390:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
391:   return(0);
392: }

396: PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
397: {
399:   PetscInt       i;
400:   AppCtx         *user;

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

406:   for (i=0; i<user->nt-1; i++){
407:     user->block_index = i;
408:     MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]);
409:     MatMult(user->M,user->yi[i+1],user->ziwork[i+1]);
410:     VecAXPY(user->yiwork[i],-1.0,user->ziwork[i+1]);
411:   }

413:   i = user->nt-1;
414:   user->block_index = i;
415:   MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]);
416:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
417:   return(0);
418: }

422: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
423: {
425:   PetscInt       i;
426:   AppCtx         *user;

429:   MatShellGetContext(J_shell,(void**)&user);
430:   i = user->block_index;
431:   VecPointwiseMult(user->uxiwork[i],X,user->uxi[i]);
432:   VecPointwiseMult(user->uyiwork[i],X,user->uyi[i]);
433:   Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
434:   MatMult(user->Div,user->uiwork[i],Y);
435:   VecAYPX(Y,user->ht,X);
436:   return(0);
437: }

441: PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y)
442: {
444:   PetscInt       i;
445:   AppCtx         *user;

448:   MatShellGetContext(J_shell,(void**)&user);
449:   i = user->block_index;
450:   MatMult(user->Grad,X,user->uiwork[i]);
451:   Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
452:   VecPointwiseMult(user->uxiwork[i],user->uxi[i],user->uxiwork[i]);
453:   VecPointwiseMult(user->uyiwork[i],user->uyi[i],user->uyiwork[i]);
454:   VecWAXPY(Y,1.0,user->uxiwork[i],user->uyiwork[i]);
455:   VecAYPX(Y,user->ht,X);
456:   return(0);
457: }

461: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
462: {
464:   PetscInt       i;
465:   AppCtx         *user;

468:   MatShellGetContext(J_shell,(void**)&user);
469:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
470:   Scatter_uxi_uyi(X,user->uxiwork,user->uxi_scatter,user->uyiwork,user->uyi_scatter,user->nt);
471:   for (i=0; i<user->nt; i++){
472:     VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]);
473:     VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]);
474:     Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
475:     MatMult(user->Div,user->uiwork[i],user->ziwork[i]);
476:     VecScale(user->ziwork[i],user->ht);
477:   }
478:   Gather_yi(Y,user->ziwork,user->yi_scatter,user->nt);
479:   return(0);
480: }

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

491:   MatShellGetContext(J_shell,(void**)&user);
492:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
493:   Scatter_yi(X,user->yiwork,user->yi_scatter,user->nt);
494:   for (i=0; i<user->nt; i++){
495:     MatMult(user->Grad,user->yiwork[i],user->uiwork[i]);
496:     Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
497:     VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]);
498:     VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]);
499:     Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]);
500:     VecScale(user->uiwork[i],user->ht);
501:   }
502:   Gather_yi(Y,user->uiwork,user->ui_scatter,user->nt);
503:   return(0);
504: }

508: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
509: {
511:   PetscInt       i;
512:   AppCtx         *user;

515:   PCShellGetContext(PC_shell,(void**)&user);
516:   i = user->block_index;
517:   if (user->c_formed) {
518:     MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
519:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
520:   return(0);
521: }

525: PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y)
526: {
528:   PetscInt       i;
529:   AppCtx         *user;

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

534:   i = user->block_index;
535:   if (user->c_formed) {
536:     MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
537:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
538:   return(0);
539: }

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

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

552:   if (Y == user->ytrue) {
553:     /* First solve is done using true solution to set up problem */
554:     KSPSetTolerances(user->solver,1e-4,1e-20,PETSC_DEFAULT,PETSC_DEFAULT);
555:   } else {
556:     KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
557:   }
558:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);
559:   Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt);
560:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

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

565:   KSPGetIterationNumber(user->solver,&its);
566:   user->ksp_its = user->ksp_its + its;
567:   for (i=1; i<user->nt; i++){
568:     MatMult(user->M,user->yiwork[i-1],user->ziwork[i-1]);
569:     VecAXPY(user->yi[i],1.0,user->ziwork[i-1]);
570:     user->block_index = i;
571:     KSPSolve(user->solver,user->yi[i],user->yiwork[i]);

573:     KSPGetIterationNumber(user->solver,&its);
574:     user->ksp_its = user->ksp_its + its;
575:   }
576:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
577:   return(0);
578: }

582: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
583: {
585:   AppCtx         *user;
586:   PetscInt       its,i;

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

591:   Scatter_yi(X,user->yi,user->yi_scatter,user->nt);
592:   Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt);
593:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

595:   i = user->nt - 1;
596:   user->block_index = i;
597:   KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]);

599:   KSPGetIterationNumber(user->solver,&its);
600:   user->ksp_its = user->ksp_its + its;

602:   for (i=user->nt-2; i>=0; i--){
603:     MatMult(user->M,user->yiwork[i+1],user->ziwork[i+1]);
604:     VecAXPY(user->yi[i],1.0,user->ziwork[i+1]);
605:     user->block_index = i;
606:     KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]);

608:     KSPGetIterationNumber(user->solver,&its);
609:     user->ksp_its = user->ksp_its + its;
610:   }
611:   Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt);
612:   return(0);
613: }

617: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
618: {
620:   AppCtx         *user;

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

625:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
626:   MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
627:   MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
628:   MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
629:   MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
630:   return(0);
631: }

635: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
636: {
638:   AppCtx         *user;

641:   MatShellGetContext(J_shell,(void**)&user);
642:    VecCopy(user->js_diag,X);
643:   return(0);
644: }

648: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
649: {
650:   /* con = Ay - q, A = [C(u1)  0     0     ...   0;
651:                          -M  C(u2)   0     ...   0;
652:                           0   -M   C(u3)   ...   0;
653:                                       ...         ;
654:                           0    ...      -M C(u_nt)]
655:      C(u) = eye + ht*Div*[diag(u1); diag(u2)]       */
657:   PetscInt       i;
658:   AppCtx         *user = (AppCtx*)ptr;

661:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
662:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
663:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

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

668:   for (i=1; i<user->nt; i++){
669:     user->block_index = i;
670:     MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
671:     MatMult(user->M,user->yi[i-1],user->ziwork[i-1]);
672:     VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]);
673:   }

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

678:   return(0);
679: }


684: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
685: {

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

698: PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
699: {
701:   PetscInt       i;

704:   for (i=0; i<nt; i++){
705:     VecScatterBegin(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD);
706:     VecScatterEnd(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD);
707:     VecScatterBegin(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD);
708:     VecScatterEnd(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD);
709:   }
710:   return(0);
711: }

715: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
716: {

720:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
721:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
722:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
723:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
724:   return(0);
725: }

729: PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
730: {
732:   PetscInt       i;

735:   for (i=0; i<nt; i++){
736:     VecScatterBegin(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE);
737:     VecScatterEnd(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE);
738:     VecScatterBegin(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE);
739:     VecScatterEnd(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE);
740:   }
741:   return(0);
742: }

746: PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
747: {
749:   PetscInt       i;

752:   for (i=0; i<nt; i++){
753:     VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
754:     VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
755:   }
756:   return(0);
757: }

761: PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
762: {
764:   PetscInt       i;

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

776: PetscErrorCode HyperbolicInitialize(AppCtx *user)
777: {
779:   PetscInt       n,i,j,linear_index,istart,iend,iblock,lo,hi;
780:   Vec            XX,YY,XXwork,YYwork,yi,uxi,ui,bc;
781:   PetscReal      h,sum;
782:   PetscScalar    hinv,neg_hinv,quarter=0.25,one=1.0,half_hinv,neg_half_hinv;
783:   PetscScalar    vx,vy,zero=0.0;
784:   IS             is_from_y,is_to_yi,is_from_u,is_to_uxi,is_to_uyi;

787:   user->jformed = PETSC_FALSE;
788:   user->c_formed = PETSC_FALSE;

790:   user->ksp_its = 0;
791:   user->ksp_its_initial = 0;

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

795:   h = 1.0/user->mx;
796:   hinv = user->mx;
797:   neg_hinv = -hinv;
798:   half_hinv = hinv / 2.0;
799:   neg_half_hinv = neg_hinv / 2.0;

801:   /* Generate Grad matrix */
802:   MatCreate(PETSC_COMM_WORLD,&user->Grad);
803:   MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,2*n,n);
804:   MatSetFromOptions(user->Grad);
805:   MatMPIAIJSetPreallocation(user->Grad,3,NULL,3,NULL);
806:   MatSeqAIJSetPreallocation(user->Grad,3,NULL);
807:   MatGetOwnershipRange(user->Grad,&istart,&iend);

809:   for (i=istart; i<iend; i++){
810:     if (i<n){
811:       iblock = i / user->mx;
812:       j = iblock*user->mx + ((i+user->mx-1) % user->mx);
813:       MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES);
814:       j = iblock*user->mx + ((i+1) % user->mx);
815:       MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES);
816:     }
817:     if (i>=n){
818:       j = (i - user->mx) % n;
819:       MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES);
820:       j = (j + 2*user->mx) % n;
821:       MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES);
822:     }
823:   }

825:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
826:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);

828:   MatCreate(PETSC_COMM_WORLD,&user->Gradxy[0]);
829:   MatSetSizes(user->Gradxy[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
830:   MatSetFromOptions(user->Gradxy[0]);
831:   MatMPIAIJSetPreallocation(user->Gradxy[0],3,NULL,3,NULL);
832:   MatSeqAIJSetPreallocation(user->Gradxy[0],3,NULL);
833:   MatGetOwnershipRange(user->Gradxy[0],&istart,&iend);

835:   for (i=istart; i<iend; i++){
836:     iblock = i / user->mx;
837:     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
838:     MatSetValues(user->Gradxy[0],1,&i,1,&j,&half_hinv,INSERT_VALUES);
839:     j = iblock*user->mx + ((i+1) % user->mx);
840:     MatSetValues(user->Gradxy[0],1,&i,1,&j,&neg_half_hinv,INSERT_VALUES);
841:     MatSetValues(user->Gradxy[0],1,&i,1,&i,&zero,INSERT_VALUES);
842:   }
843:   MatAssemblyBegin(user->Gradxy[0],MAT_FINAL_ASSEMBLY);
844:   MatAssemblyEnd(user->Gradxy[0],MAT_FINAL_ASSEMBLY);

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

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

863:   /* Generate Div matrix */
864:   MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div);
865:   MatTranspose(user->Gradxy[0],MAT_INITIAL_MATRIX,&user->Divxy[0]);
866:   MatTranspose(user->Gradxy[1],MAT_INITIAL_MATRIX,&user->Divxy[1]);

868:   /* Off-diagonal averaging matrix */
869:   MatCreate(PETSC_COMM_WORLD,&user->M);
870:   MatSetSizes(user->M,PETSC_DECIDE,PETSC_DECIDE,n,n);
871:   MatSetFromOptions(user->M);
872:   MatMPIAIJSetPreallocation(user->M,4,NULL,4,NULL);
873:   MatSeqAIJSetPreallocation(user->M,4,NULL);
874:   MatGetOwnershipRange(user->M,&istart,&iend);

876:   for (i=istart; i<iend; i++){
877:     /* kron(Id,Av) */
878:     iblock = i / user->mx;
879:     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
880:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
881:     j = iblock*user->mx + ((i+1) % user->mx);
882:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);

884:     /* kron(Av,Id) */
885:     j = (i + user->mx) % n;
886:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
887:     j = (i + n - user->mx) % n;
888:     MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES);
889:   }
890:   MatAssemblyBegin(user->M,MAT_FINAL_ASSEMBLY);
891:   MatAssemblyEnd(user->M,MAT_FINAL_ASSEMBLY);

893:   /* Generate 2D grid */
894:   VecCreate(PETSC_COMM_WORLD,&XX);
895:   VecCreate(PETSC_COMM_WORLD,&user->q);
896:   VecSetSizes(XX,PETSC_DECIDE,n);
897:   VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
898:   VecSetFromOptions(XX);
899:   VecSetFromOptions(user->q);

901:   VecDuplicate(XX,&YY);
902:   VecDuplicate(XX,&XXwork);
903:   VecDuplicate(XX,&YYwork);
904:   VecDuplicate(XX,&user->d);
905:   VecDuplicate(XX,&user->dwork);

907:   VecGetOwnershipRange(XX,&istart,&iend);
908:   for (linear_index=istart; linear_index<iend; linear_index++){
909:     i = linear_index % user->mx;
910:     j = (linear_index-i)/user->mx;
911:     vx = h*(i+0.5);
912:     vy = h*(j+0.5);
913:     VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
914:     VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
915:   }

917:   VecAssemblyBegin(XX);
918:   VecAssemblyEnd(XX);
919:   VecAssemblyBegin(YY);
920:   VecAssemblyEnd(YY);

922:   /* Compute final density function yT
923:      yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2))
924:      yT = yT / (h^2*sum(yT)) */
925:   VecCopy(XX,XXwork);
926:   VecCopy(YY,YYwork);

928:   VecShift(XXwork,-0.25);
929:   VecShift(YYwork,-0.25);

931:   VecPointwiseMult(XXwork,XXwork,XXwork);
932:   VecPointwiseMult(YYwork,YYwork,YYwork);

934:   VecCopy(XXwork,user->dwork);
935:   VecAXPY(user->dwork,1.0,YYwork);
936:   VecScale(user->dwork,-30.0);
937:   VecExp(user->dwork);
938:   VecCopy(user->dwork,user->d);

940:   VecCopy(XX,XXwork);
941:   VecCopy(YY,YYwork);

943:   VecShift(XXwork,-0.75);
944:   VecShift(YYwork,-0.75);

946:   VecPointwiseMult(XXwork,XXwork,XXwork);
947:   VecPointwiseMult(YYwork,YYwork,YYwork);

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

954:   VecAXPY(user->d,1.0,user->dwork);
955:   VecShift(user->d,1.0);
956:   VecSum(user->d,&sum);
957:   VecScale(user->d,1.0/(h*h*sum));

959:   /* Initial conditions of forward problem */
960:   VecDuplicate(XX,&bc);
961:   VecCopy(XX,XXwork);
962:   VecCopy(YY,YYwork);

964:   VecShift(XXwork,-0.5);
965:   VecShift(YYwork,-0.5);

967:   VecPointwiseMult(XXwork,XXwork,XXwork);
968:   VecPointwiseMult(YYwork,YYwork,YYwork);

970:   VecWAXPY(bc,1.0,XXwork,YYwork);
971:   VecScale(bc,-50.0);
972:   VecExp(bc);
973:   VecShift(bc,1.0);
974:   VecSum(bc,&sum);
975:   VecScale(bc,1.0/(h*h*sum));

977:   /* Create scatter from y to y_1,y_2,...,y_nt */
978:   /*  TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */
979:   PetscMalloc1(user->nt*user->mx*user->mx,&user->yi_scatter);
980:   VecCreate(PETSC_COMM_WORLD,&yi);
981:   VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx);
982:   VecSetFromOptions(yi);
983:   VecDuplicateVecs(yi,user->nt,&user->yi);
984:   VecDuplicateVecs(yi,user->nt,&user->yiwork);
985:   VecDuplicateVecs(yi,user->nt,&user->ziwork);
986:   for (i=0; i<user->nt; i++){
987:     VecGetOwnershipRange(user->yi[i],&lo,&hi);
988:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
989:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+i*user->mx*user->mx,1,&is_from_y);
990:     VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
991:     ISDestroy(&is_to_yi);
992:     ISDestroy(&is_from_y);
993:   }

995:   /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */
996:   /*  TODO: reorder for better parallelism */
997:   PetscMalloc1(user->nt*user->mx*user->mx,&user->uxi_scatter);
998:   PetscMalloc1(user->nt*user->mx*user->mx,&user->uyi_scatter);
999:   PetscMalloc1(user->nt*user->mx*user->mx,&user->ux_scatter);
1000:   PetscMalloc1(user->nt*user->mx*user->mx,&user->uy_scatter);
1001:   PetscMalloc1(2*user->nt*user->mx*user->mx,&user->ui_scatter);
1002:   VecCreate(PETSC_COMM_WORLD,&uxi);
1003:   VecCreate(PETSC_COMM_WORLD,&ui);
1004:   VecSetSizes(uxi,PETSC_DECIDE,user->mx*user->mx);
1005:   VecSetSizes(ui,PETSC_DECIDE,2*user->mx*user->mx);
1006:   VecSetFromOptions(uxi);
1007:   VecSetFromOptions(ui);
1008:   VecDuplicateVecs(uxi,user->nt,&user->uxi);
1009:   VecDuplicateVecs(uxi,user->nt,&user->uyi);
1010:   VecDuplicateVecs(uxi,user->nt,&user->uxiwork);
1011:   VecDuplicateVecs(uxi,user->nt,&user->uyiwork);
1012:   VecDuplicateVecs(ui,user->nt,&user->ui);
1013:   VecDuplicateVecs(ui,user->nt,&user->uiwork);
1014:   for (i=0; i<user->nt; i++){
1015:     VecGetOwnershipRange(user->uxi[i],&lo,&hi);
1016:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1017:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u);
1018:     VecScatterCreate(user->u,is_from_u,user->uxi[i],is_to_uxi,&user->uxi_scatter[i]);

1020:     ISDestroy(&is_to_uxi);
1021:     ISDestroy(&is_from_u);

1023:     VecGetOwnershipRange(user->uyi[i],&lo,&hi);
1024:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi);
1025:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+(2*i+1)*user->mx*user->mx,1,&is_from_u);
1026:     VecScatterCreate(user->u,is_from_u,user->uyi[i],is_to_uyi,&user->uyi_scatter[i]);

1028:     ISDestroy(&is_to_uyi);
1029:     ISDestroy(&is_from_u);

1031:     VecGetOwnershipRange(user->uxi[i],&lo,&hi);
1032:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1033:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_from_u);
1034:     VecScatterCreate(user->ui[i],is_from_u,user->uxi[i],is_to_uxi,&user->ux_scatter[i]);

1036:     ISDestroy(&is_to_uxi);
1037:     ISDestroy(&is_from_u);

1039:     VecGetOwnershipRange(user->uyi[i],&lo,&hi);
1040:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi);
1041:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+user->mx*user->mx,1,&is_from_u);
1042:     VecScatterCreate(user->ui[i],is_from_u,user->uyi[i],is_to_uyi,&user->uy_scatter[i]);

1044:     ISDestroy(&is_to_uyi);
1045:     ISDestroy(&is_from_u);

1047:     VecGetOwnershipRange(user->ui[i],&lo,&hi);
1048:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi);
1049:     ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u);
1050:     VecScatterCreate(user->u,is_from_u,user->ui[i],is_to_uxi,&user->ui_scatter[i]);

1052:     ISDestroy(&is_to_uxi);
1053:     ISDestroy(&is_from_u);
1054:   }

1056:   /* RHS of forward problem */
1057:   MatMult(user->M,bc,user->yiwork[0]);
1058:   for (i=1; i<user->nt; i++){
1059:     VecSet(user->yiwork[i],0.0);
1060:   }
1061:   Gather_yi(user->q,user->yiwork,user->yi_scatter,user->nt);

1063:   /* Compute true velocity field utrue */
1064:   VecDuplicate(user->u,&user->utrue);
1065:   for (i=0; i<user->nt; i++){
1066:     VecCopy(YY,user->uxi[i]);
1067:     VecScale(user->uxi[i],150.0*i*user->ht);
1068:     VecCopy(XX,user->uyi[i]);
1069:     VecShift(user->uyi[i],-10.0);
1070:     VecScale(user->uyi[i],15.0*i*user->ht);
1071:   }
1072:   Gather_uxi_uyi(user->utrue,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

1074:   /* Initial guess and reference model */
1075:   VecDuplicate(user->utrue,&user->ur);
1076:   for (i=0; i<user->nt; i++){
1077:     VecCopy(XX,user->uxi[i]);
1078:     VecShift(user->uxi[i],i*user->ht);
1079:     VecCopy(YY,user->uyi[i]);
1080:     VecShift(user->uyi[i],-i*user->ht);
1081:   }
1082:   Gather_uxi_uyi(user->ur,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

1084:   /* Generate regularization matrix L */
1085:   MatCreate(PETSC_COMM_WORLD,&user->LT);
1086:   MatSetSizes(user->LT,PETSC_DECIDE,PETSC_DECIDE,2*n*user->nt,n*user->nt);
1087:   MatSetFromOptions(user->LT);
1088:   MatMPIAIJSetPreallocation(user->LT,1,NULL,1,NULL);
1089:   MatSeqAIJSetPreallocation(user->LT,1,NULL);
1090:   MatGetOwnershipRange(user->LT,&istart,&iend);

1092:   for (i=istart; i<iend; i++){
1093:     iblock = (i+n) / (2*n);
1094:     j = i - iblock*n;
1095:     MatSetValues(user->LT,1,&i,1,&j,&user->gamma,INSERT_VALUES);
1096:   }

1098:   MatAssemblyBegin(user->LT,MAT_FINAL_ASSEMBLY);
1099:   MatAssemblyEnd(user->LT,MAT_FINAL_ASSEMBLY);

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

1103:   /* Build work vectors and matrices */
1104:   VecCreate(PETSC_COMM_WORLD,&user->lwork);
1105:   VecSetType(user->lwork,VECMPI);
1106:   VecSetSizes(user->lwork,PETSC_DECIDE,user->m);
1107:   VecSetFromOptions(user->lwork);

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

1111:   VecDuplicate(user->y,&user->ywork);
1112:   VecDuplicate(user->u,&user->uwork);
1113:   VecDuplicate(user->u,&user->vwork);
1114:   VecDuplicate(user->u,&user->js_diag);
1115:   VecDuplicate(user->c,&user->cwork);

1117:   /* Create matrix-free shell user->Js for computing A*x */
1118:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->Js);
1119:   MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1120:   MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
1121:   MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
1122:   MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);

1124:   /* Diagonal blocks of user->Js */
1125:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlock);
1126:   MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
1127:   MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMultTranspose);

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

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

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

1146:   /* Build matrices for SOR preconditioner */
1147:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);
1148:   MatShift(user->Divxy[0],0.0); /*  Force C[i] and Divxy[0] to share same nonzero pattern */
1149:   MatAXPY(user->Divxy[0],0.0,user->Divxy[1],DIFFERENT_NONZERO_PATTERN);
1150:   PetscMalloc1(5*n,&user->C);
1151:   PetscMalloc1(2*n,&user->Cwork);
1152:   for (i=0; i<user->nt; i++){
1153:     MatDuplicate(user->Divxy[0],MAT_COPY_VALUES,&user->C[i]);
1154:     MatDuplicate(user->Divxy[1],MAT_COPY_VALUES,&user->Cwork[i]);

1156:     MatDiagonalScale(user->C[i],NULL,user->uxi[i]);
1157:     MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]);
1158:     MatAXPY(user->C[i],1.0,user->Cwork[i],SUBSET_NONZERO_PATTERN);
1159:     MatScale(user->C[i],user->ht);
1160:     MatShift(user->C[i],1.0);
1161:   }

1163:   /* Solver options and tolerances */
1164:   KSPCreate(PETSC_COMM_WORLD,&user->solver);
1165:   KSPSetType(user->solver,KSPGMRES);
1166:   KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
1167:   KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE); /*  TODO: why is true slower? */
1168:   KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1169:   /* KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500); */
1170:   KSPGetPC(user->solver,&user->prec);
1171:   PCSetType(user->prec,PCSHELL);

1173:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1174:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMultTranspose);
1175:   PCShellSetContext(user->prec,user);

1177:   /* Compute true state function yt given ut */
1178:   VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1179:   VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1180:   VecSetFromOptions(user->ytrue);
1181:   user->c_formed = PETSC_TRUE;
1182:   VecCopy(user->utrue,user->u); /*  Set u=utrue temporarily for StateMatInv */
1183:   VecSet(user->ytrue,0.0); /*  Initial guess */
1184:   StateMatInvMult(user->Js,user->q,user->ytrue);
1185:   VecCopy(user->ur,user->u); /*  Reset u=ur */

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

1190:   /* Data discretization */
1191:   MatCreate(PETSC_COMM_WORLD,&user->Q);
1192:   MatSetSizes(user->Q,PETSC_DECIDE,PETSC_DECIDE,user->mx*user->mx,user->m);
1193:   MatSetFromOptions(user->Q);
1194:   MatMPIAIJSetPreallocation(user->Q,0,NULL,1,NULL);
1195:   MatSeqAIJSetPreallocation(user->Q,1,NULL);

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

1199:   for (i=istart; i<iend; i++){
1200:     j = i + user->m - user->mx*user->mx;
1201:     MatSetValues(user->Q,1,&i,1,&j,&one,INSERT_VALUES);
1202:   }

1204:   MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY);
1205:   MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY);

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

1209:   VecDestroy(&XX);
1210:   VecDestroy(&YY);
1211:   VecDestroy(&XXwork);
1212:   VecDestroy(&YYwork);
1213:   VecDestroy(&yi);
1214:   VecDestroy(&uxi);
1215:   VecDestroy(&ui);
1216:   VecDestroy(&bc);

1218:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1219:   KSPSetFromOptions(user->solver);
1220:   return(0);
1221: }

1225: PetscErrorCode HyperbolicDestroy(AppCtx *user)
1226: {
1228:   PetscInt       i;

1231:   MatDestroy(&user->Q);
1232:   MatDestroy(&user->QT);
1233:   MatDestroy(&user->Div);
1234:   MatDestroy(&user->Divwork);
1235:   MatDestroy(&user->Grad);
1236:   MatDestroy(&user->L);
1237:   MatDestroy(&user->LT);
1238:   KSPDestroy(&user->solver);
1239:   MatDestroy(&user->Js);
1240:   MatDestroy(&user->Jd);
1241:   MatDestroy(&user->JsBlockPrec);
1242:   MatDestroy(&user->JsInv);
1243:   MatDestroy(&user->JsBlock);
1244:   MatDestroy(&user->Divxy[0]);
1245:   MatDestroy(&user->Divxy[1]);
1246:   MatDestroy(&user->Gradxy[0]);
1247:   MatDestroy(&user->Gradxy[1]);
1248:   MatDestroy(&user->M);
1249:   for (i=0; i<user->nt; i++){
1250:     MatDestroy(&user->C[i]);
1251:     MatDestroy(&user->Cwork[i]);
1252:   }
1253:   PetscFree(user->C);
1254:   PetscFree(user->Cwork);
1255:   VecDestroy(&user->u);
1256:   VecDestroy(&user->uwork);
1257:   VecDestroy(&user->vwork);
1258:   VecDestroy(&user->utrue);
1259:   VecDestroy(&user->y);
1260:   VecDestroy(&user->ywork);
1261:   VecDestroy(&user->ytrue);
1262:   VecDestroyVecs(user->nt,&user->yi);
1263:   VecDestroyVecs(user->nt,&user->yiwork);
1264:   VecDestroyVecs(user->nt,&user->ziwork);
1265:   VecDestroyVecs(user->nt,&user->uxi);
1266:   VecDestroyVecs(user->nt,&user->uyi);
1267:   VecDestroyVecs(user->nt,&user->uxiwork);
1268:   VecDestroyVecs(user->nt,&user->uyiwork);
1269:   VecDestroyVecs(user->nt,&user->ui);
1270:   VecDestroyVecs(user->nt,&user->uiwork);
1271:   VecDestroy(&user->c);
1272:   VecDestroy(&user->cwork);
1273:   VecDestroy(&user->ur);
1274:   VecDestroy(&user->q);
1275:   VecDestroy(&user->d);
1276:   VecDestroy(&user->dwork);
1277:   VecDestroy(&user->lwork);
1278:   VecDestroy(&user->js_diag);
1279:   ISDestroy(&user->s_is);
1280:   ISDestroy(&user->d_is);
1281:   VecScatterDestroy(&user->state_scatter);
1282:   VecScatterDestroy(&user->design_scatter);
1283:   for (i=0; i<user->nt; i++){
1284:     VecScatterDestroy(&user->uxi_scatter[i]);
1285:     VecScatterDestroy(&user->uyi_scatter[i]);
1286:     VecScatterDestroy(&user->ux_scatter[i]);
1287:     VecScatterDestroy(&user->uy_scatter[i]);
1288:     VecScatterDestroy(&user->ui_scatter[i]);
1289:     VecScatterDestroy(&user->yi_scatter[i]);
1290:   }
1291:   PetscFree(user->uxi_scatter);
1292:   PetscFree(user->uyi_scatter);
1293:   PetscFree(user->ux_scatter);
1294:   PetscFree(user->uy_scatter);
1295:   PetscFree(user->ui_scatter);
1296:   PetscFree(user->yi_scatter);
1297:   return(0);
1298: }

1302: PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1303: {
1305:   Vec            X;
1306:   PetscReal      unorm,ynorm;
1307:   AppCtx         *user = (AppCtx*)ptr;

1310:   TaoGetSolutionVector(tao,&X);
1311:   Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1312:   VecAXPY(user->ywork,-1.0,user->ytrue);
1313:   VecAXPY(user->uwork,-1.0,user->utrue);
1314:   VecNorm(user->uwork,NORM_2,&unorm);
1315:   VecNorm(user->ywork,NORM_2,&ynorm);
1316:   PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1317:   return(0);
1318: }