xref: /petsc/src/tao/pde_constrained/tutorials/hyperbolic.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown #include <petsctao.h>
2c4762a1bSJed Brown 
3c4762a1bSJed Brown typedef struct {
4c4762a1bSJed Brown   PetscInt n; /*  Number of variables */
5c4762a1bSJed Brown   PetscInt m; /*  Number of constraints */
6c4762a1bSJed Brown   PetscInt mx; /*  grid points in each direction */
7c4762a1bSJed Brown   PetscInt nt; /*  Number of time steps */
8c4762a1bSJed Brown   PetscInt ndata; /*  Number of data points per sample */
9c4762a1bSJed Brown   IS       s_is;
10c4762a1bSJed Brown   IS       d_is;
11c4762a1bSJed Brown   VecScatter state_scatter;
12c4762a1bSJed Brown   VecScatter design_scatter;
13c4762a1bSJed Brown   VecScatter *uxi_scatter,*uyi_scatter,*ux_scatter,*uy_scatter,*ui_scatter;
14c4762a1bSJed Brown   VecScatter *yi_scatter;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   Mat       Js,Jd,JsBlockPrec,JsInv,JsBlock;
17c4762a1bSJed Brown   PetscBool jformed,c_formed;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   PetscReal alpha; /*  Regularization parameter */
20c4762a1bSJed Brown   PetscReal gamma;
21c4762a1bSJed Brown   PetscReal ht; /*  Time step */
22c4762a1bSJed Brown   PetscReal T; /*  Final time */
23c4762a1bSJed Brown   Mat Q,QT;
24c4762a1bSJed Brown   Mat L,LT;
25c4762a1bSJed Brown   Mat Div,Divwork,Divxy[2];
26c4762a1bSJed Brown   Mat Grad,Gradxy[2];
27c4762a1bSJed Brown   Mat M;
28c4762a1bSJed Brown   Mat *C,*Cwork;
29c4762a1bSJed Brown   /* Mat Hs,Hd,Hsd; */
30c4762a1bSJed Brown   Vec q;
31c4762a1bSJed Brown   Vec ur; /*  reference */
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   Vec d;
34c4762a1bSJed Brown   Vec dwork;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   Vec y; /*  state variables */
37c4762a1bSJed Brown   Vec ywork;
38c4762a1bSJed Brown   Vec ytrue;
39c4762a1bSJed Brown   Vec *yi,*yiwork,*ziwork;
40c4762a1bSJed Brown   Vec *uxi,*uyi,*uxiwork,*uyiwork,*ui,*uiwork;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   Vec u; /*  design variables */
43c4762a1bSJed Brown   Vec uwork,vwork;
44c4762a1bSJed Brown   Vec utrue;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   Vec js_diag;
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   Vec c; /*  constraint vector */
49c4762a1bSJed Brown   Vec cwork;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   Vec lwork;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   KSP      solver;
54c4762a1bSJed Brown   PC       prec;
55c4762a1bSJed Brown   PetscInt block_index;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   PetscInt ksp_its;
58c4762a1bSJed Brown   PetscInt ksp_its_initial;
59c4762a1bSJed Brown } AppCtx;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*);
62c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void*);
63c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*);
64c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void*);
65c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat,void*);
66c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void*);
67c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void*);
68c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
69c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
70c4762a1bSJed Brown PetscErrorCode HyperbolicInitialize(AppCtx *user);
71c4762a1bSJed Brown PetscErrorCode HyperbolicDestroy(AppCtx *user);
72c4762a1bSJed Brown PetscErrorCode HyperbolicMonitor(Tao, void*);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat,Vec,Vec);
75c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat,Vec,Vec);
76c4762a1bSJed Brown PetscErrorCode StateMatBlockMultTranspose(Mat,Vec,Vec);
77c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat,Vec,Vec);
78c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat,Vec);
79c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat,MatDuplicateOption,Mat*);
80c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat,Vec,Vec);
81c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat,Vec,Vec);
82c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC,Vec,Vec);
83c4762a1bSJed Brown 
84c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat,Vec,Vec);
85c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown PetscErrorCode Scatter_yi(Vec,Vec*,VecScatter*,PetscInt); /*  y to y1,y2,...,y_nt */
88c4762a1bSJed Brown PetscErrorCode Gather_yi(Vec,Vec*,VecScatter*,PetscInt);
89c4762a1bSJed Brown PetscErrorCode Scatter_uxi_uyi(Vec,Vec*,VecScatter*,Vec*,VecScatter*,PetscInt); /*  u to ux_1,uy_1,ux_2,uy_2,...,u */
90c4762a1bSJed Brown PetscErrorCode Gather_uxi_uyi(Vec,Vec*,VecScatter*,Vec*,VecScatter*,PetscInt);
91c4762a1bSJed Brown 
92c4762a1bSJed Brown static  char help[]="";
93c4762a1bSJed Brown 
94c4762a1bSJed Brown int main(int argc, char **argv)
95c4762a1bSJed Brown {
96c4762a1bSJed Brown   Vec                x,x0;
97c4762a1bSJed Brown   Tao                tao;
98c4762a1bSJed Brown   AppCtx             user;
99c4762a1bSJed Brown   IS                 is_allstate,is_alldesign;
100c4762a1bSJed Brown   PetscInt           lo,hi,hi2,lo2,ksp_old;
101c4762a1bSJed Brown   PetscInt           ntests = 1;
102c4762a1bSJed Brown   PetscInt           i;
103c4762a1bSJed Brown #if defined(PETSC_USE_LOG)
104c4762a1bSJed Brown   PetscLogStage      stages[1];
105c4762a1bSJed Brown #endif
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char*)0,help));
108c4762a1bSJed Brown   user.mx = 32;
109*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"hyperbolic example",NULL);
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL));
111c4762a1bSJed Brown   user.nt = 16;
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL));
113c4762a1bSJed Brown   user.ndata = 64;
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL));
115c4762a1bSJed Brown   user.alpha = 10.0;
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL));
117c4762a1bSJed Brown   user.T = 1.0/32.0;
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-Tfinal","Final time","",user.T,&user.T,NULL));
1199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL));
120*d0609cedSBarry Smith   PetscOptionsEnd();
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   user.m = user.mx*user.mx*user.nt; /*  number of constraints */
123c4762a1bSJed Brown   user.n = user.mx*user.mx*3*user.nt; /*  number of variables */
124c4762a1bSJed Brown   user.ht = user.T/user.nt; /*  Time step */
125c4762a1bSJed Brown   user.gamma = user.T*user.ht / (user.mx*user.mx);
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user.u));
1289566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user.y));
1299566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user.c));
1309566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.u,PETSC_DECIDE,user.n-user.m));
1319566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.y,PETSC_DECIDE,user.m));
1329566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.c,PETSC_DECIDE,user.m));
1339566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.u));
1349566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.y));
1359566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.c));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* Create scatters for reduced spaces.
138c4762a1bSJed Brown      If the state vector y and design vector u are partitioned as
139c4762a1bSJed Brown      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
140c4762a1bSJed Brown      then the solution vector x is organized as
141c4762a1bSJed Brown      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
142c4762a1bSJed Brown      The index sets user.s_is and user.d_is correspond to the indices of the
143c4762a1bSJed Brown      state and design variables owned by the current processor.
144c4762a1bSJed Brown   */
1459566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user.y,&lo,&hi));
1489566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user.u,&lo2,&hi2));
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate));
1519566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user.s_is));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign));
1549566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user.d_is));
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,hi-lo+hi2-lo2,user.n));
1579566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
158c4762a1bSJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,user.s_is,user.y,is_allstate,&user.state_scatter));
1609566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,user.d_is,user.u,is_alldesign,&user.design_scatter));
1619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_alldesign));
1629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_allstate));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1659566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
1669566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOLCL));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* Set up initial vectors and matrices */
1699566063dSJacob Faibussowitsch   PetscCall(HyperbolicInitialize(&user));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch   PetscCall(Gather(x,user.y,user.state_scatter,user.u,user.design_scatter));
1729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&x0));
1739566063dSJacob Faibussowitsch   PetscCall(VecCopy(x,x0));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* Set solution vector with an initial guess */
1769566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,x));
1779566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(tao, FormFunction, &user));
1789566063dSJacob Faibussowitsch   PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
1799566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
1809566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, &user));
1819566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
1829566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
1839566063dSJacob Faibussowitsch   PetscCall(TaoSetStateDesignIS(tao,user.s_is,user.d_is));
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1869566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Trials",&stages[0]));
1879566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stages[0]));
188c4762a1bSJed Brown   user.ksp_its_initial = user.ksp_its;
189c4762a1bSJed Brown   ksp_old = user.ksp_its;
190c4762a1bSJed Brown   for (i=0; i<ntests; i++) {
1919566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
1929566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its-ksp_old));
1939566063dSJacob Faibussowitsch     PetscCall(VecCopy(x0,x));
1949566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(tao,x));
195c4762a1bSJed Brown   }
1969566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1979566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)x));
1989566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: "));
1999566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial));
2009566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %D trial(s): ",ntests));
2019566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its));
2029566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: "));
2039566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%D\n",(user.ksp_its-user.ksp_its_initial)/ntests));
204c4762a1bSJed Brown 
2059566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
2069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x0));
2089566063dSJacob Faibussowitsch   PetscCall(HyperbolicDestroy(&user));
2099566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
210b122ec5aSJacob Faibussowitsch   return 0;
211c4762a1bSJed Brown }
212c4762a1bSJed Brown /* ------------------------------------------------------------------- */
213c4762a1bSJed Brown /*
214c4762a1bSJed Brown    dwork = Qy - d
215c4762a1bSJed Brown    lwork = L*(u-ur).^2
216c4762a1bSJed Brown    f = 1/2 * (dwork.dork + alpha*y.lwork)
217c4762a1bSJed Brown */
218c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
219c4762a1bSJed Brown {
220c4762a1bSJed Brown   PetscReal      d1=0,d2=0;
221c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   PetscFunctionBegin;
2249566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2259566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Q,user->y,user->dwork));
2269566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,-1.0,user->d));
2279566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork,user->dwork,&d1));
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u));
2309566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork,user->uwork,user->uwork));
2319566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L,user->uwork,user->lwork));
2329566063dSJacob Faibussowitsch   PetscCall(VecDot(user->y,user->lwork,&d2));
233c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
234c4762a1bSJed Brown   PetscFunctionReturn(0);
235c4762a1bSJed Brown }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown /* ------------------------------------------------------------------- */
238c4762a1bSJed Brown /*
239c4762a1bSJed Brown     state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2
240c4762a1bSJed Brown     design: g_d = alpha*(L'y).*(u-ur)
241c4762a1bSJed Brown */
242c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
243c4762a1bSJed Brown {
244c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   PetscFunctionBegin;
2479566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2489566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Q,user->y,user->dwork));
2499566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,-1.0,user->d));
250c4762a1bSJed Brown 
2519566063dSJacob Faibussowitsch   PetscCall(MatMult(user->QT,user->dwork,user->ywork));
252c4762a1bSJed Brown 
2539566063dSJacob Faibussowitsch   PetscCall(MatMult(user->LT,user->y,user->uwork));
2549566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->vwork,-1.0,user->ur,user->u));
2559566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork,user->vwork,user->uwork));
2569566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork,user->alpha));
257c4762a1bSJed Brown 
2589566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->vwork,user->vwork,user->vwork));
2599566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L,user->vwork,user->lwork));
2609566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork,0.5*user->alpha,user->lwork));
261c4762a1bSJed Brown 
2629566063dSJacob Faibussowitsch   PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
263c4762a1bSJed Brown   PetscFunctionReturn(0);
264c4762a1bSJed Brown }
265c4762a1bSJed Brown 
266c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
267c4762a1bSJed Brown {
268c4762a1bSJed Brown   PetscReal      d1,d2;
269c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   PetscFunctionBegin;
2729566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2739566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Q,user->y,user->dwork));
2749566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,-1.0,user->d));
275c4762a1bSJed Brown 
2769566063dSJacob Faibussowitsch   PetscCall(MatMult(user->QT,user->dwork,user->ywork));
277c4762a1bSJed Brown 
2789566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork,user->dwork,&d1));
279c4762a1bSJed Brown 
2809566063dSJacob Faibussowitsch   PetscCall(MatMult(user->LT,user->y,user->uwork));
2819566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->vwork,-1.0,user->ur,user->u));
2829566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork,user->vwork,user->uwork));
2839566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork,user->alpha));
284c4762a1bSJed Brown 
2859566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->vwork,user->vwork,user->vwork));
2869566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L,user->vwork,user->lwork));
2879566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork,0.5*user->alpha,user->lwork));
288c4762a1bSJed Brown 
2899566063dSJacob Faibussowitsch   PetscCall(VecDot(user->y,user->lwork,&d2));
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
2929566063dSJacob Faibussowitsch   PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
293c4762a1bSJed Brown   PetscFunctionReturn(0);
294c4762a1bSJed Brown }
295c4762a1bSJed Brown 
296c4762a1bSJed Brown /* ------------------------------------------------------------------- */
297c4762a1bSJed Brown /* A
298c4762a1bSJed Brown MatShell object
299c4762a1bSJed Brown */
300c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
301c4762a1bSJed Brown {
302c4762a1bSJed Brown   PetscInt       i;
303c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
3079566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->u,user->ui,user->ui_scatter,user->nt));
3089566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
309c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
3109566063dSJacob Faibussowitsch     PetscCall(MatCopy(user->Divxy[0],user->C[i],SUBSET_NONZERO_PATTERN));
3119566063dSJacob Faibussowitsch     PetscCall(MatCopy(user->Divxy[1],user->Cwork[i],SAME_NONZERO_PATTERN));
312c4762a1bSJed Brown 
3139566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->C[i],NULL,user->uxi[i]));
3149566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]));
3159566063dSJacob Faibussowitsch     PetscCall(MatAXPY(user->C[i],1.0,user->Cwork[i],SUBSET_NONZERO_PATTERN));
3169566063dSJacob Faibussowitsch     PetscCall(MatScale(user->C[i],user->ht));
3179566063dSJacob Faibussowitsch     PetscCall(MatShift(user->C[i],1.0));
318c4762a1bSJed Brown   }
319c4762a1bSJed Brown   PetscFunctionReturn(0);
320c4762a1bSJed Brown }
321c4762a1bSJed Brown 
322c4762a1bSJed Brown /* ------------------------------------------------------------------- */
323c4762a1bSJed Brown /* B */
324c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
325c4762a1bSJed Brown {
326c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   PetscFunctionBegin;
3299566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
330c4762a1bSJed Brown   PetscFunctionReturn(0);
331c4762a1bSJed Brown }
332c4762a1bSJed Brown 
333c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
334c4762a1bSJed Brown {
335c4762a1bSJed Brown   PetscInt       i;
336c4762a1bSJed Brown   AppCtx         *user;
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   PetscFunctionBegin;
3399566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
3409566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
341c4762a1bSJed Brown   user->block_index = 0;
3429566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsBlock,user->yi[0],user->yiwork[0]));
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
345c4762a1bSJed Brown     user->block_index = i;
3469566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i]));
3479566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M,user->yi[i-1],user->ziwork[i-1]));
3489566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]));
349c4762a1bSJed Brown   }
3509566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
351c4762a1bSJed Brown   PetscFunctionReturn(0);
352c4762a1bSJed Brown }
353c4762a1bSJed Brown 
354c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
355c4762a1bSJed Brown {
356c4762a1bSJed Brown   PetscInt       i;
357c4762a1bSJed Brown   AppCtx         *user;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   PetscFunctionBegin;
3609566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
3619566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   for (i=0; i<user->nt-1; i++) {
364c4762a1bSJed Brown     user->block_index = i;
3659566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]));
3669566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M,user->yi[i+1],user->ziwork[i+1]));
3679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i],-1.0,user->ziwork[i+1]));
368c4762a1bSJed Brown   }
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   i = user->nt-1;
371c4762a1bSJed Brown   user->block_index = i;
3729566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]));
3739566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
374c4762a1bSJed Brown   PetscFunctionReturn(0);
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
378c4762a1bSJed Brown {
379c4762a1bSJed Brown   PetscInt       i;
380c4762a1bSJed Brown   AppCtx         *user;
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
384c4762a1bSJed Brown   i = user->block_index;
3859566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uxiwork[i],X,user->uxi[i]));
3869566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uyiwork[i],X,user->uyi[i]));
3879566063dSJacob Faibussowitsch   PetscCall(Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
3889566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Div,user->uiwork[i],Y));
3899566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Y,user->ht,X));
390c4762a1bSJed Brown   PetscFunctionReturn(0);
391c4762a1bSJed Brown }
392c4762a1bSJed Brown 
393c4762a1bSJed Brown PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y)
394c4762a1bSJed Brown {
395c4762a1bSJed Brown   PetscInt       i;
396c4762a1bSJed Brown   AppCtx         *user;
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   PetscFunctionBegin;
3999566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
400c4762a1bSJed Brown   i = user->block_index;
4019566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Grad,X,user->uiwork[i]));
4029566063dSJacob Faibussowitsch   PetscCall(Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4039566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uxiwork[i],user->uxi[i],user->uxiwork[i]));
4049566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uyiwork[i],user->uyi[i],user->uyiwork[i]));
4059566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Y,1.0,user->uxiwork[i],user->uyiwork[i]));
4069566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Y,user->ht,X));
407c4762a1bSJed Brown   PetscFunctionReturn(0);
408c4762a1bSJed Brown }
409c4762a1bSJed Brown 
410c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
411c4762a1bSJed Brown {
412c4762a1bSJed Brown   PetscInt       i;
413c4762a1bSJed Brown   AppCtx         *user;
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   PetscFunctionBegin;
4169566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
4179566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt));
4189566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(X,user->uxiwork,user->uxi_scatter,user->uyiwork,user->uyi_scatter,user->nt));
419c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
4209566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]));
4219566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]));
4229566063dSJacob Faibussowitsch     PetscCall(Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4239566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Div,user->uiwork[i],user->ziwork[i]));
4249566063dSJacob Faibussowitsch     PetscCall(VecScale(user->ziwork[i],user->ht));
425c4762a1bSJed Brown   }
4269566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y,user->ziwork,user->yi_scatter,user->nt));
427c4762a1bSJed Brown   PetscFunctionReturn(0);
428c4762a1bSJed Brown }
429c4762a1bSJed Brown 
430c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
431c4762a1bSJed Brown {
432c4762a1bSJed Brown   PetscInt       i;
433c4762a1bSJed Brown   AppCtx         *user;
434c4762a1bSJed Brown 
435c4762a1bSJed Brown   PetscFunctionBegin;
4369566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
4379566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt));
4389566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X,user->yiwork,user->yi_scatter,user->nt));
439c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
4409566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad,user->yiwork[i],user->uiwork[i]));
4419566063dSJacob Faibussowitsch     PetscCall(Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4429566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]));
4439566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]));
4449566063dSJacob Faibussowitsch     PetscCall(Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4459566063dSJacob Faibussowitsch     PetscCall(VecScale(user->uiwork[i],user->ht));
446c4762a1bSJed Brown   }
4479566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y,user->uiwork,user->ui_scatter,user->nt));
448c4762a1bSJed Brown   PetscFunctionReturn(0);
449c4762a1bSJed Brown }
450c4762a1bSJed Brown 
451c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
452c4762a1bSJed Brown {
453c4762a1bSJed Brown   PetscInt       i;
454c4762a1bSJed Brown   AppCtx         *user;
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   PetscFunctionBegin;
4579566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(PC_shell,&user));
458c4762a1bSJed Brown   i = user->block_index;
459c4762a1bSJed Brown   if (user->c_formed) {
4609566063dSJacob Faibussowitsch     PetscCall(MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y));
461c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
462c4762a1bSJed Brown   PetscFunctionReturn(0);
463c4762a1bSJed Brown }
464c4762a1bSJed Brown 
465c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y)
466c4762a1bSJed Brown {
467c4762a1bSJed Brown   PetscInt       i;
468c4762a1bSJed Brown   AppCtx         *user;
469c4762a1bSJed Brown 
470c4762a1bSJed Brown   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(PC_shell,&user));
472c4762a1bSJed Brown 
473c4762a1bSJed Brown   i = user->block_index;
474c4762a1bSJed Brown   if (user->c_formed) {
4759566063dSJacob Faibussowitsch     PetscCall(MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y));
476c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
477c4762a1bSJed Brown   PetscFunctionReturn(0);
478c4762a1bSJed Brown }
479c4762a1bSJed Brown 
480c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
481c4762a1bSJed Brown {
482c4762a1bSJed Brown   AppCtx         *user;
483c4762a1bSJed Brown   PetscInt       its,i;
484c4762a1bSJed Brown 
485c4762a1bSJed Brown   PetscFunctionBegin;
4869566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   if (Y == user->ytrue) {
489c4762a1bSJed Brown     /* First solve is done using true solution to set up problem */
4909566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver,1e-4,1e-20,PETSC_DEFAULT,PETSC_DEFAULT));
491c4762a1bSJed Brown   } else {
4929566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
493c4762a1bSJed Brown   }
4949566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
4959566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt));
4969566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
497c4762a1bSJed Brown 
498c4762a1bSJed Brown   user->block_index = 0;
4999566063dSJacob Faibussowitsch   PetscCall(KSPSolve(user->solver,user->yi[0],user->yiwork[0]));
500c4762a1bSJed Brown 
5019566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(user->solver,&its));
502c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
503c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
5049566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M,user->yiwork[i-1],user->ziwork[i-1]));
5059566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yi[i],1.0,user->ziwork[i-1]));
506c4762a1bSJed Brown     user->block_index = i;
5079566063dSJacob Faibussowitsch     PetscCall(KSPSolve(user->solver,user->yi[i],user->yiwork[i]));
508c4762a1bSJed Brown 
5099566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(user->solver,&its));
510c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
511c4762a1bSJed Brown   }
5129566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
513c4762a1bSJed Brown   PetscFunctionReturn(0);
514c4762a1bSJed Brown }
515c4762a1bSJed Brown 
516c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
517c4762a1bSJed Brown {
518c4762a1bSJed Brown   AppCtx         *user;
519c4762a1bSJed Brown   PetscInt       its,i;
520c4762a1bSJed Brown 
521c4762a1bSJed Brown   PetscFunctionBegin;
5229566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
523c4762a1bSJed Brown 
5249566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
5259566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt));
5269566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
527c4762a1bSJed Brown 
528c4762a1bSJed Brown   i = user->nt - 1;
529c4762a1bSJed Brown   user->block_index = i;
5309566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]));
531c4762a1bSJed Brown 
5329566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(user->solver,&its));
533c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
534c4762a1bSJed Brown 
535c4762a1bSJed Brown   for (i=user->nt-2; i>=0; i--) {
5369566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M,user->yiwork[i+1],user->ziwork[i+1]));
5379566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yi[i],1.0,user->ziwork[i+1]));
538c4762a1bSJed Brown     user->block_index = i;
5399566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]));
540c4762a1bSJed Brown 
5419566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(user->solver,&its));
542c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
543c4762a1bSJed Brown   }
5449566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
545c4762a1bSJed Brown   PetscFunctionReturn(0);
546c4762a1bSJed Brown }
547c4762a1bSJed Brown 
548c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
549c4762a1bSJed Brown {
550c4762a1bSJed Brown   AppCtx         *user;
551c4762a1bSJed Brown 
552c4762a1bSJed Brown   PetscFunctionBegin;
5539566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
554c4762a1bSJed Brown 
5559566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell));
5569566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult));
5579566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate));
5589566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose));
5599566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal));
560c4762a1bSJed Brown   PetscFunctionReturn(0);
561c4762a1bSJed Brown }
562c4762a1bSJed Brown 
563c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
564c4762a1bSJed Brown {
565c4762a1bSJed Brown   AppCtx         *user;
566c4762a1bSJed Brown 
567c4762a1bSJed Brown   PetscFunctionBegin;
5689566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
5699566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->js_diag,X));
570c4762a1bSJed Brown   PetscFunctionReturn(0);
571c4762a1bSJed Brown }
572c4762a1bSJed Brown 
573c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
574c4762a1bSJed Brown {
575c4762a1bSJed Brown   /* con = Ay - q, A = [C(u1)  0     0     ...   0;
576c4762a1bSJed Brown                          -M  C(u2)   0     ...   0;
577c4762a1bSJed Brown                           0   -M   C(u3)   ...   0;
578c4762a1bSJed Brown                                       ...         ;
579c4762a1bSJed Brown                           0    ...      -M C(u_nt)]
580c4762a1bSJed Brown      C(u) = eye + ht*Div*[diag(u1); diag(u2)]       */
581c4762a1bSJed Brown   PetscInt       i;
582c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
583c4762a1bSJed Brown 
584c4762a1bSJed Brown   PetscFunctionBegin;
5859566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
5869566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt));
5879566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
588c4762a1bSJed Brown 
589c4762a1bSJed Brown   user->block_index = 0;
5909566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsBlock,user->yi[0],user->yiwork[0]));
591c4762a1bSJed Brown 
592c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
593c4762a1bSJed Brown     user->block_index = i;
5949566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i]));
5959566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M,user->yi[i-1],user->ziwork[i-1]));
5969566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]));
597c4762a1bSJed Brown   }
598c4762a1bSJed Brown 
5999566063dSJacob Faibussowitsch   PetscCall(Gather_yi(C,user->yiwork,user->yi_scatter,user->nt));
6009566063dSJacob Faibussowitsch   PetscCall(VecAXPY(C,-1.0,user->q));
601c4762a1bSJed Brown 
602c4762a1bSJed Brown   PetscFunctionReturn(0);
603c4762a1bSJed Brown }
604c4762a1bSJed Brown 
605c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
606c4762a1bSJed Brown {
607c4762a1bSJed Brown   PetscFunctionBegin;
6089566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD));
6099566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD));
6109566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD));
6119566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD));
612c4762a1bSJed Brown   PetscFunctionReturn(0);
613c4762a1bSJed Brown }
614c4762a1bSJed Brown 
615c4762a1bSJed Brown PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
616c4762a1bSJed Brown {
617c4762a1bSJed Brown   PetscInt       i;
618c4762a1bSJed Brown 
619c4762a1bSJed Brown   PetscFunctionBegin;
620c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6219566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD));
6229566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD));
6239566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD));
6249566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD));
625c4762a1bSJed Brown   }
626c4762a1bSJed Brown   PetscFunctionReturn(0);
627c4762a1bSJed Brown }
628c4762a1bSJed Brown 
629c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
630c4762a1bSJed Brown {
631c4762a1bSJed Brown   PetscFunctionBegin;
6329566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE));
6339566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE));
6349566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE));
6359566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE));
636c4762a1bSJed Brown   PetscFunctionReturn(0);
637c4762a1bSJed Brown }
638c4762a1bSJed Brown 
639c4762a1bSJed Brown PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
640c4762a1bSJed Brown {
641c4762a1bSJed Brown   PetscInt       i;
642c4762a1bSJed Brown 
643c4762a1bSJed Brown   PetscFunctionBegin;
644c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6459566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE));
6469566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE));
6479566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE));
6489566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE));
649c4762a1bSJed Brown   }
650c4762a1bSJed Brown   PetscFunctionReturn(0);
651c4762a1bSJed Brown }
652c4762a1bSJed Brown 
653c4762a1bSJed Brown PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
654c4762a1bSJed Brown {
655c4762a1bSJed Brown   PetscInt       i;
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   PetscFunctionBegin;
658c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6599566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD));
6609566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD));
661c4762a1bSJed Brown   }
662c4762a1bSJed Brown   PetscFunctionReturn(0);
663c4762a1bSJed Brown }
664c4762a1bSJed Brown 
665c4762a1bSJed Brown PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
666c4762a1bSJed Brown {
667c4762a1bSJed Brown   PetscInt       i;
668c4762a1bSJed Brown 
669c4762a1bSJed Brown   PetscFunctionBegin;
670c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6719566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE));
6729566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE));
673c4762a1bSJed Brown   }
674c4762a1bSJed Brown   PetscFunctionReturn(0);
675c4762a1bSJed Brown }
676c4762a1bSJed Brown 
677c4762a1bSJed Brown PetscErrorCode HyperbolicInitialize(AppCtx *user)
678c4762a1bSJed Brown {
679c4762a1bSJed Brown   PetscInt       n,i,j,linear_index,istart,iend,iblock,lo,hi;
680c4762a1bSJed Brown   Vec            XX,YY,XXwork,YYwork,yi,uxi,ui,bc;
681c4762a1bSJed Brown   PetscReal      h,sum;
682c4762a1bSJed Brown   PetscScalar    hinv,neg_hinv,quarter=0.25,one=1.0,half_hinv,neg_half_hinv;
683c4762a1bSJed Brown   PetscScalar    vx,vy,zero=0.0;
684c4762a1bSJed Brown   IS             is_from_y,is_to_yi,is_from_u,is_to_uxi,is_to_uyi;
685c4762a1bSJed Brown 
686c4762a1bSJed Brown   PetscFunctionBegin;
687c4762a1bSJed Brown   user->jformed = PETSC_FALSE;
688c4762a1bSJed Brown   user->c_formed = PETSC_FALSE;
689c4762a1bSJed Brown 
690c4762a1bSJed Brown   user->ksp_its = 0;
691c4762a1bSJed Brown   user->ksp_its_initial = 0;
692c4762a1bSJed Brown 
693c4762a1bSJed Brown   n = user->mx * user->mx;
694c4762a1bSJed Brown 
695c4762a1bSJed Brown   h = 1.0/user->mx;
696c4762a1bSJed Brown   hinv = user->mx;
697c4762a1bSJed Brown   neg_hinv = -hinv;
698c4762a1bSJed Brown   half_hinv = hinv / 2.0;
699c4762a1bSJed Brown   neg_half_hinv = neg_hinv / 2.0;
700c4762a1bSJed Brown 
701c4762a1bSJed Brown   /* Generate Grad matrix */
7029566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Grad));
7039566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,2*n,n));
7049566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Grad));
7059566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Grad,3,NULL,3,NULL));
7069566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Grad,3,NULL));
7079566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Grad,&istart,&iend));
708c4762a1bSJed Brown 
709c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
710c4762a1bSJed Brown     if (i<n) {
711c4762a1bSJed Brown       iblock = i / user->mx;
712c4762a1bSJed Brown       j = iblock*user->mx + ((i+user->mx-1) % user->mx);
7139566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES));
714c4762a1bSJed Brown       j = iblock*user->mx + ((i+1) % user->mx);
7159566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
716c4762a1bSJed Brown     }
717c4762a1bSJed Brown     if (i>=n) {
718c4762a1bSJed Brown       j = (i - user->mx) % n;
7199566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES));
720c4762a1bSJed Brown       j = (j + 2*user->mx) % n;
7219566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
722c4762a1bSJed Brown     }
723c4762a1bSJed Brown   }
724c4762a1bSJed Brown 
7259566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY));
7269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY));
727c4762a1bSJed Brown 
7289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Gradxy[0]));
7299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Gradxy[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
7309566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Gradxy[0]));
7319566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[0],3,NULL,3,NULL));
7329566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[0],3,NULL));
7339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Gradxy[0],&istart,&iend));
734c4762a1bSJed Brown 
735c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
736c4762a1bSJed Brown     iblock = i / user->mx;
737c4762a1bSJed Brown     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
7389566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[0],1,&i,1,&j,&half_hinv,INSERT_VALUES));
739c4762a1bSJed Brown     j = iblock*user->mx + ((i+1) % user->mx);
7409566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[0],1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
7419566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[0],1,&i,1,&i,&zero,INSERT_VALUES));
742c4762a1bSJed Brown   }
7439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Gradxy[0],MAT_FINAL_ASSEMBLY));
7449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Gradxy[0],MAT_FINAL_ASSEMBLY));
745c4762a1bSJed Brown 
7469566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Gradxy[1]));
7479566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Gradxy[1],PETSC_DECIDE,PETSC_DECIDE,n,n));
7489566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Gradxy[1]));
7499566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[1],3,NULL,3,NULL));
7509566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[1],3,NULL));
7519566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Gradxy[1],&istart,&iend));
752c4762a1bSJed Brown 
753c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
754c4762a1bSJed Brown     j = (i + n - user->mx) % n;
7559566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[1],1,&i,1,&j,&half_hinv,INSERT_VALUES));
756c4762a1bSJed Brown     j = (j + 2*user->mx) % n;
7579566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[1],1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
7589566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[1],1,&i,1,&i,&zero,INSERT_VALUES));
759c4762a1bSJed Brown   }
7609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Gradxy[1],MAT_FINAL_ASSEMBLY));
7619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Gradxy[1],MAT_FINAL_ASSEMBLY));
762c4762a1bSJed Brown 
763c4762a1bSJed Brown   /* Generate Div matrix */
7649566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div));
7659566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Gradxy[0],MAT_INITIAL_MATRIX,&user->Divxy[0]));
7669566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Gradxy[1],MAT_INITIAL_MATRIX,&user->Divxy[1]));
767c4762a1bSJed Brown 
768c4762a1bSJed Brown   /* Off-diagonal averaging matrix */
7699566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->M));
7709566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->M,PETSC_DECIDE,PETSC_DECIDE,n,n));
7719566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->M));
7729566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->M,4,NULL,4,NULL));
7739566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->M,4,NULL));
7749566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->M,&istart,&iend));
775c4762a1bSJed Brown 
776c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
777c4762a1bSJed Brown     /* kron(Id,Av) */
778c4762a1bSJed Brown     iblock = i / user->mx;
779c4762a1bSJed Brown     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
7809566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
781c4762a1bSJed Brown     j = iblock*user->mx + ((i+1) % user->mx);
7829566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
783c4762a1bSJed Brown 
784c4762a1bSJed Brown     /* kron(Av,Id) */
785c4762a1bSJed Brown     j = (i + user->mx) % n;
7869566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
787c4762a1bSJed Brown     j = (i + n - user->mx) % n;
7889566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
789c4762a1bSJed Brown   }
7909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->M,MAT_FINAL_ASSEMBLY));
7919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->M,MAT_FINAL_ASSEMBLY));
792c4762a1bSJed Brown 
793c4762a1bSJed Brown   /* Generate 2D grid */
7949566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&XX));
7959566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->q));
7969566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(XX,PETSC_DECIDE,n));
7979566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->q,PETSC_DECIDE,n*user->nt));
7989566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(XX));
7999566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->q));
800c4762a1bSJed Brown 
8019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&YY));
8029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&XXwork));
8039566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&YYwork));
8049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&user->d));
8059566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&user->dwork));
806c4762a1bSJed Brown 
8079566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(XX,&istart,&iend));
808c4762a1bSJed Brown   for (linear_index=istart; linear_index<iend; linear_index++) {
809c4762a1bSJed Brown     i = linear_index % user->mx;
810c4762a1bSJed Brown     j = (linear_index-i)/user->mx;
811c4762a1bSJed Brown     vx = h*(i+0.5);
812c4762a1bSJed Brown     vy = h*(j+0.5);
8139566063dSJacob Faibussowitsch     PetscCall(VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES));
8149566063dSJacob Faibussowitsch     PetscCall(VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES));
815c4762a1bSJed Brown   }
816c4762a1bSJed Brown 
8179566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(XX));
8189566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(XX));
8199566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(YY));
8209566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(YY));
821c4762a1bSJed Brown 
822c4762a1bSJed Brown   /* Compute final density function yT
823c4762a1bSJed Brown      yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2))
824c4762a1bSJed Brown      yT = yT / (h^2*sum(yT)) */
8259566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX,XXwork));
8269566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY,YYwork));
827c4762a1bSJed Brown 
8289566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork,-0.25));
8299566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork,-0.25));
830c4762a1bSJed Brown 
8319566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork,XXwork,XXwork));
8329566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork,YYwork,YYwork));
833c4762a1bSJed Brown 
8349566063dSJacob Faibussowitsch   PetscCall(VecCopy(XXwork,user->dwork));
8359566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,1.0,YYwork));
8369566063dSJacob Faibussowitsch   PetscCall(VecScale(user->dwork,-30.0));
8379566063dSJacob Faibussowitsch   PetscCall(VecExp(user->dwork));
8389566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->dwork,user->d));
839c4762a1bSJed Brown 
8409566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX,XXwork));
8419566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY,YYwork));
842c4762a1bSJed Brown 
8439566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork,-0.75));
8449566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork,-0.75));
845c4762a1bSJed Brown 
8469566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork,XXwork,XXwork));
8479566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork,YYwork,YYwork));
848c4762a1bSJed Brown 
8499566063dSJacob Faibussowitsch   PetscCall(VecCopy(XXwork,user->dwork));
8509566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,1.0,YYwork));
8519566063dSJacob Faibussowitsch   PetscCall(VecScale(user->dwork,-30.0));
8529566063dSJacob Faibussowitsch   PetscCall(VecExp(user->dwork));
853c4762a1bSJed Brown 
8549566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->d,1.0,user->dwork));
8559566063dSJacob Faibussowitsch   PetscCall(VecShift(user->d,1.0));
8569566063dSJacob Faibussowitsch   PetscCall(VecSum(user->d,&sum));
8579566063dSJacob Faibussowitsch   PetscCall(VecScale(user->d,1.0/(h*h*sum)));
858c4762a1bSJed Brown 
859c4762a1bSJed Brown   /* Initial conditions of forward problem */
8609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&bc));
8619566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX,XXwork));
8629566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY,YYwork));
863c4762a1bSJed Brown 
8649566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork,-0.5));
8659566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork,-0.5));
866c4762a1bSJed Brown 
8679566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork,XXwork,XXwork));
8689566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork,YYwork,YYwork));
869c4762a1bSJed Brown 
8709566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(bc,1.0,XXwork,YYwork));
8719566063dSJacob Faibussowitsch   PetscCall(VecScale(bc,-50.0));
8729566063dSJacob Faibussowitsch   PetscCall(VecExp(bc));
8739566063dSJacob Faibussowitsch   PetscCall(VecShift(bc,1.0));
8749566063dSJacob Faibussowitsch   PetscCall(VecSum(bc,&sum));
8759566063dSJacob Faibussowitsch   PetscCall(VecScale(bc,1.0/(h*h*sum)));
876c4762a1bSJed Brown 
877c4762a1bSJed Brown   /* Create scatter from y to y_1,y_2,...,y_nt */
878c4762a1bSJed Brown   /*  TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */
8799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt*user->mx*user->mx,&user->yi_scatter));
8809566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&yi));
8819566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx));
8829566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(yi));
8839566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi,user->nt,&user->yi));
8849566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi,user->nt,&user->yiwork));
8859566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi,user->nt,&user->ziwork));
886c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
8879566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->yi[i],&lo,&hi));
8889566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi));
8899566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+i*user->mx*user->mx,1,&is_from_y));
8909566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]));
8919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_yi));
8929566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_y));
893c4762a1bSJed Brown   }
894c4762a1bSJed Brown 
895c4762a1bSJed Brown   /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */
896c4762a1bSJed Brown   /*  TODO: reorder for better parallelism */
8979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt*user->mx*user->mx,&user->uxi_scatter));
8989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt*user->mx*user->mx,&user->uyi_scatter));
8999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt*user->mx*user->mx,&user->ux_scatter));
9009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt*user->mx*user->mx,&user->uy_scatter));
9019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*user->nt*user->mx*user->mx,&user->ui_scatter));
9029566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&uxi));
9039566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&ui));
9049566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(uxi,PETSC_DECIDE,user->mx*user->mx));
9059566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(ui,PETSC_DECIDE,2*user->mx*user->mx));
9069566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(uxi));
9079566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(ui));
9089566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi,user->nt,&user->uxi));
9099566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi,user->nt,&user->uyi));
9109566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi,user->nt,&user->uxiwork));
9119566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi,user->nt,&user->uyiwork));
9129566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ui,user->nt,&user->ui));
9139566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ui,user->nt,&user->uiwork));
914c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
9159566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uxi[i],&lo,&hi));
9169566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi));
9179566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u));
9189566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->u,is_from_u,user->uxi[i],is_to_uxi,&user->uxi_scatter[i]));
919c4762a1bSJed Brown 
9209566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uxi));
9219566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
922c4762a1bSJed Brown 
9239566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uyi[i],&lo,&hi));
9249566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi));
9259566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+(2*i+1)*user->mx*user->mx,1,&is_from_u));
9269566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->u,is_from_u,user->uyi[i],is_to_uyi,&user->uyi_scatter[i]));
927c4762a1bSJed Brown 
9289566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uyi));
9299566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
930c4762a1bSJed Brown 
9319566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uxi[i],&lo,&hi));
9329566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi));
9339566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_from_u));
9349566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->ui[i],is_from_u,user->uxi[i],is_to_uxi,&user->ux_scatter[i]));
935c4762a1bSJed Brown 
9369566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uxi));
9379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
938c4762a1bSJed Brown 
9399566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uyi[i],&lo,&hi));
9409566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi));
9419566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+user->mx*user->mx,1,&is_from_u));
9429566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->ui[i],is_from_u,user->uyi[i],is_to_uyi,&user->uy_scatter[i]));
943c4762a1bSJed Brown 
9449566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uyi));
9459566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
946c4762a1bSJed Brown 
9479566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->ui[i],&lo,&hi));
9489566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi));
9499566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u));
9509566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->u,is_from_u,user->ui[i],is_to_uxi,&user->ui_scatter[i]));
951c4762a1bSJed Brown 
9529566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uxi));
9539566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
954c4762a1bSJed Brown   }
955c4762a1bSJed Brown 
956c4762a1bSJed Brown   /* RHS of forward problem */
9579566063dSJacob Faibussowitsch   PetscCall(MatMult(user->M,bc,user->yiwork[0]));
958c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
9599566063dSJacob Faibussowitsch     PetscCall(VecSet(user->yiwork[i],0.0));
960c4762a1bSJed Brown   }
9619566063dSJacob Faibussowitsch   PetscCall(Gather_yi(user->q,user->yiwork,user->yi_scatter,user->nt));
962c4762a1bSJed Brown 
963c4762a1bSJed Brown   /* Compute true velocity field utrue */
9649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->utrue));
965c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
9669566063dSJacob Faibussowitsch     PetscCall(VecCopy(YY,user->uxi[i]));
9679566063dSJacob Faibussowitsch     PetscCall(VecScale(user->uxi[i],150.0*i*user->ht));
9689566063dSJacob Faibussowitsch     PetscCall(VecCopy(XX,user->uyi[i]));
9699566063dSJacob Faibussowitsch     PetscCall(VecShift(user->uyi[i],-10.0));
9709566063dSJacob Faibussowitsch     PetscCall(VecScale(user->uyi[i],15.0*i*user->ht));
971c4762a1bSJed Brown   }
9729566063dSJacob Faibussowitsch   PetscCall(Gather_uxi_uyi(user->utrue,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
973c4762a1bSJed Brown 
974c4762a1bSJed Brown   /* Initial guess and reference model */
9759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->utrue,&user->ur));
976c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
9779566063dSJacob Faibussowitsch     PetscCall(VecCopy(XX,user->uxi[i]));
9789566063dSJacob Faibussowitsch     PetscCall(VecShift(user->uxi[i],i*user->ht));
9799566063dSJacob Faibussowitsch     PetscCall(VecCopy(YY,user->uyi[i]));
9809566063dSJacob Faibussowitsch     PetscCall(VecShift(user->uyi[i],-i*user->ht));
981c4762a1bSJed Brown   }
9829566063dSJacob Faibussowitsch   PetscCall(Gather_uxi_uyi(user->ur,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
983c4762a1bSJed Brown 
984c4762a1bSJed Brown   /* Generate regularization matrix L */
9859566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->LT));
9869566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->LT,PETSC_DECIDE,PETSC_DECIDE,2*n*user->nt,n*user->nt));
9879566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->LT));
9889566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->LT,1,NULL,1,NULL));
9899566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->LT,1,NULL));
9909566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->LT,&istart,&iend));
991c4762a1bSJed Brown 
992c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
993c4762a1bSJed Brown     iblock = (i+n) / (2*n);
994c4762a1bSJed Brown     j = i - iblock*n;
9959566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->LT,1,&i,1,&j,&user->gamma,INSERT_VALUES));
996c4762a1bSJed Brown   }
997c4762a1bSJed Brown 
9989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->LT,MAT_FINAL_ASSEMBLY));
9999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->LT,MAT_FINAL_ASSEMBLY));
1000c4762a1bSJed Brown 
10019566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->LT,MAT_INITIAL_MATRIX,&user->L));
1002c4762a1bSJed Brown 
1003c4762a1bSJed Brown   /* Build work vectors and matrices */
10049566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->lwork));
10059566063dSJacob Faibussowitsch   PetscCall(VecSetType(user->lwork,VECMPI));
10069566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->lwork,PETSC_DECIDE,user->m));
10079566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->lwork));
1008c4762a1bSJed Brown 
10099566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork));
1010c4762a1bSJed Brown 
10119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->y,&user->ywork));
10129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->uwork));
10139566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->vwork));
10149566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->js_diag));
10159566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->c,&user->cwork));
1016c4762a1bSJed Brown 
1017c4762a1bSJed Brown   /* Create matrix-free shell user->Js for computing A*x */
10189566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->Js));
10199566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult));
10209566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate));
10219566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose));
10229566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal));
1023c4762a1bSJed Brown 
1024c4762a1bSJed Brown   /* Diagonal blocks of user->Js */
10259566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlock));
10269566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult));
10279566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMultTranspose));
1028c4762a1bSJed Brown 
1029c4762a1bSJed Brown   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1030c4762a1bSJed Brown      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1031c4762a1bSJed Brown      This is an SOR preconditioner for user->JsBlock. */
10329566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlockPrec));
10339566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult));
10349566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMultTranspose));
1035c4762a1bSJed Brown 
1036c4762a1bSJed Brown   /* Create a matrix-free shell user->Jd for computing B*x */
10379566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->n-user->m,user,&user->Jd));
10389566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult));
10399566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose));
1040c4762a1bSJed Brown 
1041c4762a1bSJed Brown   /* User-defined routines for computing user->Js\x and user->Js^T\x*/
10429566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsInv));
10439566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult));
10449566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult));
1045c4762a1bSJed Brown 
1046c4762a1bSJed Brown   /* Build matrices for SOR preconditioner */
10479566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
10489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(5*n,&user->C));
10499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*n,&user->Cwork));
1050c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
10519566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(user->Divxy[0],MAT_COPY_VALUES,&user->C[i]));
10529566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(user->Divxy[1],MAT_COPY_VALUES,&user->Cwork[i]));
1053c4762a1bSJed Brown 
10549566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->C[i],NULL,user->uxi[i]));
10559566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]));
10569566063dSJacob Faibussowitsch     PetscCall(MatAXPY(user->C[i],1.0,user->Cwork[i],DIFFERENT_NONZERO_PATTERN));
10579566063dSJacob Faibussowitsch     PetscCall(MatScale(user->C[i],user->ht));
10589566063dSJacob Faibussowitsch     PetscCall(MatShift(user->C[i],1.0));
1059c4762a1bSJed Brown   }
1060c4762a1bSJed Brown 
1061c4762a1bSJed Brown   /* Solver options and tolerances */
10629566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD,&user->solver));
10639566063dSJacob Faibussowitsch   PetscCall(KSPSetType(user->solver,KSPGMRES));
10649566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec));
10659566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500));
10669566063dSJacob Faibussowitsch   /* PetscCall(KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500)); */
10679566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(user->solver,&user->prec));
10689566063dSJacob Faibussowitsch   PetscCall(PCSetType(user->prec,PCSHELL));
1069c4762a1bSJed Brown 
10709566063dSJacob Faibussowitsch   PetscCall(PCShellSetApply(user->prec,StateMatBlockPrecMult));
10719566063dSJacob Faibussowitsch   PetscCall(PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMultTranspose));
10729566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(user->prec,user));
1073c4762a1bSJed Brown 
1074c4762a1bSJed Brown   /* Compute true state function yt given ut */
10759566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->ytrue));
10769566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt));
10779566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->ytrue));
1078c4762a1bSJed Brown   user->c_formed = PETSC_TRUE;
10799566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->utrue,user->u)); /*  Set u=utrue temporarily for StateMatInv */
10809566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ytrue,0.0)); /*  Initial guess */
10819566063dSJacob Faibussowitsch   PetscCall(StateMatInvMult(user->Js,user->q,user->ytrue));
10829566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->ur,user->u)); /*  Reset u=ur */
1083c4762a1bSJed Brown 
1084c4762a1bSJed Brown   /* Initial guess y0 for state given u0 */
10859566063dSJacob Faibussowitsch   PetscCall(StateMatInvMult(user->Js,user->q,user->y));
1086c4762a1bSJed Brown 
1087c4762a1bSJed Brown   /* Data discretization */
10889566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Q));
10899566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Q,PETSC_DECIDE,PETSC_DECIDE,user->mx*user->mx,user->m));
10909566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Q));
10919566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Q,0,NULL,1,NULL));
10929566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Q,1,NULL));
1093c4762a1bSJed Brown 
10949566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Q,&istart,&iend));
1095c4762a1bSJed Brown 
1096c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
1097c4762a1bSJed Brown     j = i + user->m - user->mx*user->mx;
10989566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&one,INSERT_VALUES));
1099c4762a1bSJed Brown   }
1100c4762a1bSJed Brown 
11019566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY));
11029566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY));
1103c4762a1bSJed Brown 
11049566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Q,MAT_INITIAL_MATRIX,&user->QT));
1105c4762a1bSJed Brown 
11069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XX));
11079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YY));
11089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XXwork));
11099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YYwork));
11109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&yi));
11119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&uxi));
11129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ui));
11139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bc));
1114c4762a1bSJed Brown 
1115c4762a1bSJed Brown   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11169566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(user->solver));
1117c4762a1bSJed Brown   PetscFunctionReturn(0);
1118c4762a1bSJed Brown }
1119c4762a1bSJed Brown 
1120c4762a1bSJed Brown PetscErrorCode HyperbolicDestroy(AppCtx *user)
1121c4762a1bSJed Brown {
1122c4762a1bSJed Brown   PetscInt       i;
1123c4762a1bSJed Brown 
1124c4762a1bSJed Brown   PetscFunctionBegin;
11259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Q));
11269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->QT));
11279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Div));
11289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Divwork));
11299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Grad));
11309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->L));
11319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->LT));
11329566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&user->solver));
11339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Js));
11349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Jd));
11359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsBlockPrec));
11369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsInv));
11379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsBlock));
11389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Divxy[0]));
11399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Divxy[1]));
11409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Gradxy[0]));
11419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Gradxy[1]));
11429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->M));
1143c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
11449566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->C[i]));
11459566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->Cwork[i]));
1146c4762a1bSJed Brown   }
11479566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->C));
11489566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->Cwork));
11499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->u));
11509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->uwork));
11519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->vwork));
11529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->utrue));
11539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->y));
11549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ywork));
11559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ytrue));
11569566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->yi));
11579566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->yiwork));
11589566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->ziwork));
11599566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->uxi));
11609566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->uyi));
11619566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->uxiwork));
11629566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->uyiwork));
11639566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->ui));
11649566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt,&user->uiwork));
11659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->c));
11669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->cwork));
11679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ur));
11689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->q));
11699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->d));
11709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->dwork));
11719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->lwork));
11729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->js_diag));
11739566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->s_is));
11749566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->d_is));
11759566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->state_scatter));
11769566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->design_scatter));
1177c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
11789566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->uxi_scatter[i]));
11799566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->uyi_scatter[i]));
11809566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->ux_scatter[i]));
11819566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->uy_scatter[i]));
11829566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->ui_scatter[i]));
11839566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
1184c4762a1bSJed Brown   }
11859566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->uxi_scatter));
11869566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->uyi_scatter));
11879566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->ux_scatter));
11889566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->uy_scatter));
11899566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->ui_scatter));
11909566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->yi_scatter));
1191c4762a1bSJed Brown   PetscFunctionReturn(0);
1192c4762a1bSJed Brown }
1193c4762a1bSJed Brown 
1194c4762a1bSJed Brown PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1195c4762a1bSJed Brown {
1196c4762a1bSJed Brown   Vec            X;
1197c4762a1bSJed Brown   PetscReal      unorm,ynorm;
1198c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
1199c4762a1bSJed Brown 
1200c4762a1bSJed Brown   PetscFunctionBegin;
12019566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao,&X));
12029566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
12039566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork,-1.0,user->ytrue));
12049566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork,-1.0,user->utrue));
12059566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->uwork,NORM_2,&unorm));
12069566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->ywork,NORM_2,&ynorm));
12079566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm));
1208c4762a1bSJed Brown   PetscFunctionReturn(0);
1209c4762a1bSJed Brown }
1210c4762a1bSJed Brown 
1211c4762a1bSJed Brown /*TEST
1212c4762a1bSJed Brown 
1213c4762a1bSJed Brown    build:
1214c4762a1bSJed Brown       requires: !complex
1215c4762a1bSJed Brown 
1216c4762a1bSJed Brown    test:
1217c4762a1bSJed Brown       requires: !single
1218c4762a1bSJed Brown       args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -tao_gatol 1.e-5
1219c4762a1bSJed Brown 
1220c4762a1bSJed Brown    test:
1221c4762a1bSJed Brown       suffix: guess_pod
1222c4762a1bSJed Brown       requires: !single
1223c4762a1bSJed Brown       args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -ksp_guess_type pod -tao_gatol 1.e-5
1224c4762a1bSJed Brown 
1225c4762a1bSJed Brown TEST*/
1226