Actual source code: hyperbolic.c

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

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

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

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

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

 52:   Vec d;
 53:   Vec dwork;

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

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

 65:   Vec js_diag;

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

 70:   Vec lwork;

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

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


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

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

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

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

113: static  char help[]="";

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

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

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

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,"tao_lcl");
189:   user.lcl = (TAO_LCL*)(tao->data);

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

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

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

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

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

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

211:   /* SOLVE THE APPLICATION */
212:   PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,&flag);
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:   TaoGetTerminationReason(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, MatStructure* flag, 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:   *JInv = user->JsInv;
367:   return(0);
368: }

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

380:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
381:   *J = user->Jd;
382:   return(0);
383: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

558: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
559: {
561:   AppCtx         *user;
562:   PetscReal      tau;
563:   PetscInt       its,i;

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

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

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

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

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

598: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
599: {
601:   AppCtx         *user;
602:   PetscReal      tau;
603:   PetscInt       its,i;

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

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

616:   i = user->nt - 1;
617:   user->block_index = i;
618:   KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]);

620:   KSPGetIterationNumber(user->solver,&its);
621:   user->ksp_its = user->ksp_its + its;

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

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

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

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

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

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

662:   MatShellGetContext(J_shell,(void**)&user);
663:    VecCopy(user->js_diag,X);
664:   return(0);
665: }

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

682:   Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
683:   Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt);
684:   Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt);

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

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

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

699:   return(0);
700: }


705: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
706: {

710:   VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
711:   VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
712:   VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
713:   VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
714:   return(0);
715: }

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

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

736: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
737: {

741:   VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
742:   VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
743:   VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
744:   VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
745:   return(0);
746: }

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

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

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

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

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

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

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

808:   user->jformed = PETSC_FALSE;
809:   user->c_formed = PETSC_FALSE;

811:   user->ksp_its = 0;
812:   user->ksp_its_initial = 0;

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

816:   h = 1.0/user->mx;
817:   hinv = user->mx;
818:   neg_hinv = -hinv;
819:   half_hinv = hinv / 2.0;
820:   neg_half_hinv = neg_hinv / 2.0;

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

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

846:   MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
847:   MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);

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

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

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

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

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

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

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

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

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

922:   VecDuplicate(XX,&YY);
923:   VecDuplicate(XX,&XXwork);
924:   VecDuplicate(XX,&YYwork);
925:   VecDuplicate(XX,&user->d);
926:   VecDuplicate(XX,&user->dwork);

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

938:   VecAssemblyBegin(XX);
939:   VecAssemblyEnd(XX);
940:   VecAssemblyBegin(YY);
941:   VecAssemblyEnd(YY);

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

949:   VecShift(XXwork,-0.25);
950:   VecShift(YYwork,-0.25);

952:   VecPointwiseMult(XXwork,XXwork,XXwork);
953:   VecPointwiseMult(YYwork,YYwork,YYwork);

955:   VecCopy(XXwork,user->dwork);
956:   VecAXPY(user->dwork,1.0,YYwork);
957:   VecScale(user->dwork,-30.0);
958:   VecExp(user->dwork);
959:   VecCopy(user->dwork,user->d);

961:   VecCopy(XX,XXwork);
962:   VecCopy(YY,YYwork);

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

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

970:   VecCopy(XXwork,user->dwork);
971:   VecAXPY(user->dwork,1.0,YYwork);
972:   VecScale(user->dwork,-30.0);
973:   VecExp(user->dwork);

975:   VecAXPY(user->d,1.0,user->dwork);
976:   VecShift(user->d,1.0);
977:   VecSum(user->d,&sum);
978:   VecScale(user->d,1.0/(h*h*sum));

980:   /* Initial conditions of forward problem */
981:   VecDuplicate(XX,&bc);
982:   VecCopy(XX,XXwork);
983:   VecCopy(YY,YYwork);

985:   VecShift(XXwork,-0.5);
986:   VecShift(YYwork,-0.5);

988:   VecPointwiseMult(XXwork,XXwork,XXwork);
989:   VecPointwiseMult(YYwork,YYwork,YYwork);

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

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

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

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

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

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

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

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

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

1065:     ISDestroy(&is_to_uyi);
1066:     ISDestroy(&is_from_u);

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

1073:     ISDestroy(&is_to_uxi);
1074:     ISDestroy(&is_from_u);
1075:   }

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

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

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

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

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

1119:   MatAssemblyBegin(user->LT,MAT_FINAL_ASSEMBLY);
1120:   MatAssemblyEnd(user->LT,MAT_FINAL_ASSEMBLY);

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

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

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

1132:   VecDuplicate(user->y,&user->ywork);
1133:   VecDuplicate(user->u,&user->uwork);
1134:   VecDuplicate(user->u,&user->vwork);
1135:   VecDuplicate(user->u,&user->js_diag);
1136:   VecDuplicate(user->c,&user->cwork);

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

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

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

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

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

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

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

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

1194:   PCShellSetApply(user->prec,StateMatBlockPrecMult);
1195:   PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMultTranspose);
1196:   PCShellSetContext(user->prec,user);

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

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

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

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

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

1226:   MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY);
1227:   MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY);

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

1231:   VecDestroy(&XX);
1232:   VecDestroy(&YY);
1233:   VecDestroy(&XXwork);
1234:   VecDestroy(&YYwork);
1235:   VecDestroy(&yi);
1236:   VecDestroy(&uxi);
1237:   VecDestroy(&ui);
1238:   VecDestroy(&bc);

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

1247: PetscErrorCode HyperbolicDestroy(AppCtx *user)
1248: {
1250:   PetscInt       i;

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

1324: PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1325: {
1327:   Vec            X;
1328:   PetscReal      unorm,ynorm;
1329:   AppCtx         *user = (AppCtx*)ptr;

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