xref: /petsc/src/tao/pde_constrained/tutorials/hyperbolic.c (revision 10978b7d1d81984ed7edceb76c290f82705a4964)
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 
94d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
95d71ae5a4SJacob Faibussowitsch {
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   PetscLogStage stages[1];
104c4762a1bSJed Brown 
105327415f7SBarry Smith   PetscFunctionBeginUser;
1069566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
107c4762a1bSJed Brown   user.mx = 32;
108d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "hyperbolic example", NULL);
1099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL));
110c4762a1bSJed Brown   user.nt = 16;
1119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-nt", "Number of time steps", "", user.nt, &user.nt, NULL));
112c4762a1bSJed Brown   user.ndata = 64;
1139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL));
114c4762a1bSJed Brown   user.alpha = 10.0;
1159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL));
116c4762a1bSJed Brown   user.T = 1.0 / 32.0;
1179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-Tfinal", "Final time", "", user.T, &user.T, NULL));
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL));
119d0609cedSBarry Smith   PetscOptionsEnd();
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   user.m     = user.mx * user.mx * user.nt;     /*  number of constraints */
122c4762a1bSJed Brown   user.n     = user.mx * user.mx * 3 * user.nt; /*  number of variables */
123c4762a1bSJed Brown   user.ht    = user.T / user.nt;                /*  Time step */
124c4762a1bSJed Brown   user.gamma = user.T * user.ht / (user.mx * user.mx);
125c4762a1bSJed Brown 
1269566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.u));
1279566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.y));
1289566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.c));
1299566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.u, PETSC_DECIDE, user.n - user.m));
1309566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.y, PETSC_DECIDE, user.m));
1319566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.c, PETSC_DECIDE, user.m));
1329566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.u));
1339566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.y));
1349566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.c));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* Create scatters for reduced spaces.
137c4762a1bSJed Brown      If the state vector y and design vector u are partitioned as
138c4762a1bSJed Brown      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
139c4762a1bSJed Brown      then the solution vector x is organized as
140c4762a1bSJed Brown      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
141c4762a1bSJed Brown      The index sets user.s_is and user.d_is correspond to the indices of the
142c4762a1bSJed Brown      state and design variables owned by the current processor.
143c4762a1bSJed Brown   */
1449566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user.y, &lo, &hi));
1479566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user.u, &lo2, &hi2));
148c4762a1bSJed Brown 
1499566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate));
1509566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user.s_is));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign));
1539566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user.d_is));
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, hi - lo + hi2 - lo2, user.n));
1569566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, user.s_is, user.y, is_allstate, &user.state_scatter));
1599566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, user.d_is, user.u, is_alldesign, &user.design_scatter));
1609566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_alldesign));
1619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_allstate));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1649566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1659566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOLCL));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Set up initial vectors and matrices */
1689566063dSJacob Faibussowitsch   PetscCall(HyperbolicInitialize(&user));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch   PetscCall(Gather(x, user.y, user.state_scatter, user.u, user.design_scatter));
1719566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &x0));
1729566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, x0));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* Set solution vector with an initial guess */
1759566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
1769566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(tao, FormFunction, &user));
1779566063dSJacob Faibussowitsch   PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
1789566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
1799566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, &user));
1809566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
1819566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
1829566063dSJacob Faibussowitsch   PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1859566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Trials", &stages[0]));
1869566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stages[0]));
187c4762a1bSJed Brown   user.ksp_its_initial = user.ksp_its;
188c4762a1bSJed Brown   ksp_old              = user.ksp_its;
189c4762a1bSJed Brown   for (i = 0; i < ntests; i++) {
1909566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
19163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its - ksp_old));
1929566063dSJacob Faibussowitsch     PetscCall(VecCopy(x0, x));
1939566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(tao, x));
194c4762a1bSJed Brown   }
1959566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1969566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)x));
1979566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: "));
19863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial));
19963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total KSP iterations over %" PetscInt_FMT " trial(s): ", ntests));
20063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its));
2019566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations per trial: "));
20263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", (user.ksp_its - user.ksp_its_initial) / ntests));
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
2059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x0));
2079566063dSJacob Faibussowitsch   PetscCall(HyperbolicDestroy(&user));
2089566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
209b122ec5aSJacob Faibussowitsch   return 0;
210c4762a1bSJed Brown }
211c4762a1bSJed Brown /* ------------------------------------------------------------------- */
212c4762a1bSJed Brown /*
213c4762a1bSJed Brown    dwork = Qy - d
214c4762a1bSJed Brown    lwork = L*(u-ur).^2
215c4762a1bSJed Brown    f = 1/2 * (dwork.dork + alpha*y.lwork)
216c4762a1bSJed Brown */
217d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr)
218d71ae5a4SJacob Faibussowitsch {
219c4762a1bSJed Brown   PetscReal d1 = 0, d2 = 0;
220c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ptr;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   PetscFunctionBegin;
2239566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2249566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Q, user->y, user->dwork));
2259566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2269566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork, user->dwork, &d1));
227c4762a1bSJed Brown 
2289566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
2299566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork, user->uwork, user->uwork));
2309566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L, user->uwork, user->lwork));
2319566063dSJacob Faibussowitsch   PetscCall(VecDot(user->y, user->lwork, &d2));
232c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha * d2);
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234c4762a1bSJed Brown }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown /* ------------------------------------------------------------------- */
237c4762a1bSJed Brown /*
238c4762a1bSJed Brown     state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2
239c4762a1bSJed Brown     design: g_d = alpha*(L'y).*(u-ur)
240c4762a1bSJed Brown */
241d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr)
242d71ae5a4SJacob Faibussowitsch {
243c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   PetscFunctionBegin;
2469566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2479566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Q, user->y, user->dwork));
2489566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
249c4762a1bSJed Brown 
2509566063dSJacob Faibussowitsch   PetscCall(MatMult(user->QT, user->dwork, user->ywork));
251c4762a1bSJed Brown 
2529566063dSJacob Faibussowitsch   PetscCall(MatMult(user->LT, user->y, user->uwork));
2539566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u));
2549566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork));
2559566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork, user->alpha));
256c4762a1bSJed Brown 
2579566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork));
2589566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L, user->vwork, user->lwork));
2599566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork));
260c4762a1bSJed Brown 
2619566063dSJacob Faibussowitsch   PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263c4762a1bSJed Brown }
264c4762a1bSJed Brown 
265d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
266d71ae5a4SJacob Faibussowitsch {
267c4762a1bSJed Brown   PetscReal d1, d2;
268c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ptr;
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   PetscFunctionBegin;
2719566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2729566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Q, user->y, user->dwork));
2739566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
274c4762a1bSJed Brown 
2759566063dSJacob Faibussowitsch   PetscCall(MatMult(user->QT, user->dwork, user->ywork));
276c4762a1bSJed Brown 
2779566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork, user->dwork, &d1));
278c4762a1bSJed Brown 
2799566063dSJacob Faibussowitsch   PetscCall(MatMult(user->LT, user->y, user->uwork));
2809566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u));
2819566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork));
2829566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork, user->alpha));
283c4762a1bSJed Brown 
2849566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork));
2859566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L, user->vwork, user->lwork));
2869566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork));
287c4762a1bSJed Brown 
2889566063dSJacob Faibussowitsch   PetscCall(VecDot(user->y, user->lwork, &d2));
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha * d2);
2919566063dSJacob Faibussowitsch   PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown /* ------------------------------------------------------------------- */
296c4762a1bSJed Brown /* A
297c4762a1bSJed Brown MatShell object
298c4762a1bSJed Brown */
299d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
300d71ae5a4SJacob Faibussowitsch {
301c4762a1bSJed Brown   PetscInt i;
302c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
3069566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->u, user->ui, user->ui_scatter, user->nt));
3079566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
308c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
3099566063dSJacob Faibussowitsch     PetscCall(MatCopy(user->Divxy[0], user->C[i], SUBSET_NONZERO_PATTERN));
3109566063dSJacob Faibussowitsch     PetscCall(MatCopy(user->Divxy[1], user->Cwork[i], SAME_NONZERO_PATTERN));
311c4762a1bSJed Brown 
3129566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i]));
3139566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i]));
3149566063dSJacob Faibussowitsch     PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], SUBSET_NONZERO_PATTERN));
3159566063dSJacob Faibussowitsch     PetscCall(MatScale(user->C[i], user->ht));
3169566063dSJacob Faibussowitsch     PetscCall(MatShift(user->C[i], 1.0));
317c4762a1bSJed Brown   }
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
319c4762a1bSJed Brown }
320c4762a1bSJed Brown 
321c4762a1bSJed Brown /* ------------------------------------------------------------------- */
322c4762a1bSJed Brown /* B */
323d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
324d71ae5a4SJacob Faibussowitsch {
325c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
326c4762a1bSJed Brown 
327c4762a1bSJed Brown   PetscFunctionBegin;
3289566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
330c4762a1bSJed Brown }
331c4762a1bSJed Brown 
332d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
333d71ae5a4SJacob Faibussowitsch {
334c4762a1bSJed Brown   PetscInt i;
335c4762a1bSJed Brown   AppCtx  *user;
336c4762a1bSJed Brown 
337c4762a1bSJed Brown   PetscFunctionBegin;
3389566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
3399566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
340c4762a1bSJed Brown   user->block_index = 0;
3419566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   for (i = 1; i < user->nt; i++) {
344c4762a1bSJed Brown     user->block_index = i;
3459566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3469566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1]));
3479566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1]));
348c4762a1bSJed Brown   }
3499566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
3503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown 
353d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
354d71ae5a4SJacob Faibussowitsch {
355c4762a1bSJed Brown   PetscInt i;
356c4762a1bSJed Brown   AppCtx  *user;
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   PetscFunctionBegin;
3599566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
3609566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   for (i = 0; i < user->nt - 1; i++) {
363c4762a1bSJed Brown     user->block_index = i;
3649566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i]));
3659566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M, user->yi[i + 1], user->ziwork[i + 1]));
3669566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i + 1]));
367c4762a1bSJed Brown   }
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   i                 = user->nt - 1;
370c4762a1bSJed Brown   user->block_index = i;
3719566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i]));
3729566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
3733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown 
376d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
377d71ae5a4SJacob Faibussowitsch {
378c4762a1bSJed Brown   PetscInt i;
379c4762a1bSJed Brown   AppCtx  *user;
380c4762a1bSJed Brown 
381c4762a1bSJed Brown   PetscFunctionBegin;
3829566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
383c4762a1bSJed Brown   i = user->block_index;
3849566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uxiwork[i], X, user->uxi[i]));
3859566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uyiwork[i], X, user->uyi[i]));
3869566063dSJacob Faibussowitsch   PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
3879566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Div, user->uiwork[i], Y));
3889566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Y, user->ht, X));
3893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
390c4762a1bSJed Brown }
391c4762a1bSJed Brown 
392d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y)
393d71ae5a4SJacob Faibussowitsch {
394c4762a1bSJed Brown   PetscInt i;
395c4762a1bSJed Brown   AppCtx  *user;
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   PetscFunctionBegin;
3989566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
399c4762a1bSJed Brown   i = user->block_index;
4009566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Grad, X, user->uiwork[i]));
4019566063dSJacob Faibussowitsch   PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4029566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uxiwork[i], user->uxi[i], user->uxiwork[i]));
4039566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uyiwork[i], user->uyi[i], user->uyiwork[i]));
4049566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Y, 1.0, user->uxiwork[i], user->uyiwork[i]));
4059566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Y, user->ht, X));
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
407c4762a1bSJed Brown }
408c4762a1bSJed Brown 
409d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
410d71ae5a4SJacob Faibussowitsch {
411c4762a1bSJed Brown   PetscInt i;
412c4762a1bSJed Brown   AppCtx  *user;
413c4762a1bSJed Brown 
414c4762a1bSJed Brown   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
4169566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt));
4179566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(X, user->uxiwork, user->uxi_scatter, user->uyiwork, user->uyi_scatter, user->nt));
418c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
4199566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i]));
4209566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i]));
4219566063dSJacob Faibussowitsch     PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4229566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Div, user->uiwork[i], user->ziwork[i]));
4239566063dSJacob Faibussowitsch     PetscCall(VecScale(user->ziwork[i], user->ht));
424c4762a1bSJed Brown   }
4259566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y, user->ziwork, user->yi_scatter, user->nt));
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
427c4762a1bSJed Brown }
428c4762a1bSJed Brown 
429d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
430d71ae5a4SJacob Faibussowitsch {
431c4762a1bSJed Brown   PetscInt i;
432c4762a1bSJed Brown   AppCtx  *user;
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   PetscFunctionBegin;
4359566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
4369566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt));
4379566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X, user->yiwork, user->yi_scatter, user->nt));
438c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
4399566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad, user->yiwork[i], user->uiwork[i]));
4409566063dSJacob Faibussowitsch     PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4419566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i]));
4429566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i]));
4439566063dSJacob Faibussowitsch     PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4449566063dSJacob Faibussowitsch     PetscCall(VecScale(user->uiwork[i], user->ht));
445c4762a1bSJed Brown   }
4469566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y, user->uiwork, user->ui_scatter, user->nt));
4473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
448c4762a1bSJed Brown }
449c4762a1bSJed Brown 
450d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
451d71ae5a4SJacob Faibussowitsch {
452c4762a1bSJed Brown   PetscInt i;
453c4762a1bSJed Brown   AppCtx  *user;
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   PetscFunctionBegin;
4569566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(PC_shell, &user));
457c4762a1bSJed Brown   i = user->block_index;
458c4762a1bSJed Brown   if (user->c_formed) {
4599566063dSJacob Faibussowitsch     PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y));
460c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed");
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
462c4762a1bSJed Brown }
463c4762a1bSJed Brown 
464d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y)
465d71ae5a4SJacob Faibussowitsch {
466c4762a1bSJed Brown   PetscInt i;
467c4762a1bSJed Brown   AppCtx  *user;
468c4762a1bSJed Brown 
469c4762a1bSJed Brown   PetscFunctionBegin;
4709566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(PC_shell, &user));
471c4762a1bSJed Brown 
472c4762a1bSJed Brown   i = user->block_index;
473c4762a1bSJed Brown   if (user->c_formed) {
4749566063dSJacob Faibussowitsch     PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y));
475c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed");
4763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
477c4762a1bSJed Brown }
478c4762a1bSJed Brown 
479d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
480d71ae5a4SJacob Faibussowitsch {
481c4762a1bSJed Brown   AppCtx  *user;
482c4762a1bSJed Brown   PetscInt its, i;
483c4762a1bSJed Brown 
484c4762a1bSJed Brown   PetscFunctionBegin;
4859566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   if (Y == user->ytrue) {
488c4762a1bSJed Brown     /* First solve is done using true solution to set up problem */
4899566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, PETSC_DEFAULT, PETSC_DEFAULT));
490c4762a1bSJed Brown   } else {
4919566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
492c4762a1bSJed Brown   }
4939566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
4949566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt));
4959566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   user->block_index = 0;
4989566063dSJacob Faibussowitsch   PetscCall(KSPSolve(user->solver, user->yi[0], user->yiwork[0]));
499c4762a1bSJed Brown 
5009566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(user->solver, &its));
501c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
502c4762a1bSJed Brown   for (i = 1; i < user->nt; i++) {
5039566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M, user->yiwork[i - 1], user->ziwork[i - 1]));
5049566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i - 1]));
505c4762a1bSJed Brown     user->block_index = i;
5069566063dSJacob Faibussowitsch     PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
507c4762a1bSJed Brown 
5089566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(user->solver, &its));
509c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
510c4762a1bSJed Brown   }
5119566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
513c4762a1bSJed Brown }
514c4762a1bSJed Brown 
515d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
516d71ae5a4SJacob Faibussowitsch {
517c4762a1bSJed Brown   AppCtx  *user;
518c4762a1bSJed Brown   PetscInt its, i;
519c4762a1bSJed Brown 
520c4762a1bSJed Brown   PetscFunctionBegin;
5219566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
522c4762a1bSJed Brown 
5239566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
5249566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt));
5259566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
526c4762a1bSJed Brown 
527c4762a1bSJed Brown   i                 = user->nt - 1;
528c4762a1bSJed Brown   user->block_index = i;
5299566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(user->solver, user->yi[i], user->yiwork[i]));
530c4762a1bSJed Brown 
5319566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(user->solver, &its));
532c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
533c4762a1bSJed Brown 
534c4762a1bSJed Brown   for (i = user->nt - 2; i >= 0; i--) {
5359566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M, user->yiwork[i + 1], user->ziwork[i + 1]));
5369566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i + 1]));
537c4762a1bSJed Brown     user->block_index = i;
5389566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(user->solver, user->yi[i], user->yiwork[i]));
539c4762a1bSJed Brown 
5409566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(user->solver, &its));
541c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
542c4762a1bSJed Brown   }
5439566063dSJacob Faibussowitsch   PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
545c4762a1bSJed Brown }
546c4762a1bSJed Brown 
547d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
548d71ae5a4SJacob Faibussowitsch {
549c4762a1bSJed Brown   AppCtx *user;
550c4762a1bSJed Brown 
551c4762a1bSJed Brown   PetscFunctionBegin;
5529566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
553c4762a1bSJed Brown 
5549566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, new_shell));
5559566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT, (void (*)(void))StateMatMult));
5569566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate));
5579566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose));
5589566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal));
5593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
560c4762a1bSJed Brown }
561c4762a1bSJed Brown 
562d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
563d71ae5a4SJacob Faibussowitsch {
564c4762a1bSJed Brown   AppCtx *user;
565c4762a1bSJed Brown 
566c4762a1bSJed Brown   PetscFunctionBegin;
5679566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
5689566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->js_diag, X));
5693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
570c4762a1bSJed Brown }
571c4762a1bSJed Brown 
572d71ae5a4SJacob Faibussowitsch PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
573d71ae5a4SJacob Faibussowitsch {
574c4762a1bSJed Brown   /* con = Ay - q, A = [C(u1)  0     0     ...   0;
575c4762a1bSJed Brown                          -M  C(u2)   0     ...   0;
576c4762a1bSJed Brown                           0   -M   C(u3)   ...   0;
577c4762a1bSJed Brown                                       ...         ;
578c4762a1bSJed Brown                           0    ...      -M C(u_nt)]
579c4762a1bSJed Brown      C(u) = eye + ht*Div*[diag(u1); diag(u2)]       */
580c4762a1bSJed Brown   PetscInt i;
581c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
582c4762a1bSJed Brown 
583c4762a1bSJed Brown   PetscFunctionBegin;
5849566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
5859566063dSJacob Faibussowitsch   PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt));
5869566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   user->block_index = 0;
5899566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
590c4762a1bSJed Brown 
591c4762a1bSJed Brown   for (i = 1; i < user->nt; i++) {
592c4762a1bSJed Brown     user->block_index = i;
5939566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
5949566063dSJacob Faibussowitsch     PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1]));
5959566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1]));
596c4762a1bSJed Brown   }
597c4762a1bSJed Brown 
5989566063dSJacob Faibussowitsch   PetscCall(Gather_yi(C, user->yiwork, user->yi_scatter, user->nt));
5999566063dSJacob Faibussowitsch   PetscCall(VecAXPY(C, -1.0, user->q));
600c4762a1bSJed Brown 
6013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
602c4762a1bSJed Brown }
603c4762a1bSJed Brown 
604d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
605d71ae5a4SJacob Faibussowitsch {
606c4762a1bSJed Brown   PetscFunctionBegin;
6079566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6089566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6099566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
6109566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
612c4762a1bSJed Brown }
613c4762a1bSJed Brown 
614d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
615d71ae5a4SJacob Faibussowitsch {
616c4762a1bSJed Brown   PetscInt i;
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   PetscFunctionBegin;
619c4762a1bSJed Brown   for (i = 0; i < nt; i++) {
6209566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD));
6219566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD));
6229566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD));
6239566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD));
624c4762a1bSJed Brown   }
6253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
626c4762a1bSJed Brown }
627c4762a1bSJed Brown 
628d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
629d71ae5a4SJacob Faibussowitsch {
630c4762a1bSJed Brown   PetscFunctionBegin;
6319566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6329566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6339566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
6349566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
636c4762a1bSJed Brown }
637c4762a1bSJed Brown 
638d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
639d71ae5a4SJacob Faibussowitsch {
640c4762a1bSJed Brown   PetscInt i;
641c4762a1bSJed Brown 
642c4762a1bSJed Brown   PetscFunctionBegin;
643c4762a1bSJed Brown   for (i = 0; i < nt; i++) {
6449566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE));
6459566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE));
6469566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE));
6479566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE));
648c4762a1bSJed Brown   }
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
650c4762a1bSJed Brown }
651c4762a1bSJed Brown 
652d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
653d71ae5a4SJacob Faibussowitsch {
654c4762a1bSJed Brown   PetscInt i;
655c4762a1bSJed Brown 
656c4762a1bSJed Brown   PetscFunctionBegin;
657c4762a1bSJed Brown   for (i = 0; i < nt; i++) {
6589566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
6599566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
660c4762a1bSJed Brown   }
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
662c4762a1bSJed Brown }
663c4762a1bSJed Brown 
664d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
665d71ae5a4SJacob Faibussowitsch {
666c4762a1bSJed Brown   PetscInt i;
667c4762a1bSJed Brown 
668c4762a1bSJed Brown   PetscFunctionBegin;
669c4762a1bSJed Brown   for (i = 0; i < nt; i++) {
6709566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
6719566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
672c4762a1bSJed Brown   }
6733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
674c4762a1bSJed Brown }
675c4762a1bSJed Brown 
676d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicInitialize(AppCtx *user)
677d71ae5a4SJacob Faibussowitsch {
678c4762a1bSJed Brown   PetscInt    n, i, j, linear_index, istart, iend, iblock, lo, hi;
679c4762a1bSJed Brown   Vec         XX, YY, XXwork, YYwork, yi, uxi, ui, bc;
680c4762a1bSJed Brown   PetscReal   h, sum;
681c4762a1bSJed Brown   PetscScalar hinv, neg_hinv, quarter = 0.25, one = 1.0, half_hinv, neg_half_hinv;
682c4762a1bSJed Brown   PetscScalar vx, vy, zero = 0.0;
683c4762a1bSJed Brown   IS          is_from_y, is_to_yi, is_from_u, is_to_uxi, is_to_uyi;
684c4762a1bSJed Brown 
685c4762a1bSJed Brown   PetscFunctionBegin;
686c4762a1bSJed Brown   user->jformed  = PETSC_FALSE;
687c4762a1bSJed Brown   user->c_formed = PETSC_FALSE;
688c4762a1bSJed Brown 
689c4762a1bSJed Brown   user->ksp_its         = 0;
690c4762a1bSJed Brown   user->ksp_its_initial = 0;
691c4762a1bSJed Brown 
692c4762a1bSJed Brown   n = user->mx * user->mx;
693c4762a1bSJed Brown 
694c4762a1bSJed Brown   h             = 1.0 / user->mx;
695c4762a1bSJed Brown   hinv          = user->mx;
696c4762a1bSJed Brown   neg_hinv      = -hinv;
697c4762a1bSJed Brown   half_hinv     = hinv / 2.0;
698c4762a1bSJed Brown   neg_half_hinv = neg_hinv / 2.0;
699c4762a1bSJed Brown 
700c4762a1bSJed Brown   /* Generate Grad matrix */
7019566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad));
7029566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, PETSC_DECIDE, 2 * n, n));
7039566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Grad));
7049566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Grad, 3, NULL, 3, NULL));
7059566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Grad, 3, NULL));
7069566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));
707c4762a1bSJed Brown 
708c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
709c4762a1bSJed Brown     if (i < n) {
710c4762a1bSJed Brown       iblock = i / user->mx;
711c4762a1bSJed Brown       j      = iblock * user->mx + ((i + user->mx - 1) % user->mx);
7129566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
713c4762a1bSJed Brown       j = iblock * user->mx + ((i + 1) % user->mx);
7149566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
715c4762a1bSJed Brown     }
716c4762a1bSJed Brown     if (i >= n) {
717c4762a1bSJed Brown       j = (i - user->mx) % n;
7189566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
719c4762a1bSJed Brown       j = (j + 2 * user->mx) % n;
7209566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
721c4762a1bSJed Brown     }
722c4762a1bSJed Brown   }
723c4762a1bSJed Brown 
7249566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY));
7259566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY));
726c4762a1bSJed Brown 
7279566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[0]));
7289566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Gradxy[0], PETSC_DECIDE, PETSC_DECIDE, n, n));
7299566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Gradxy[0]));
7309566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[0], 3, NULL, 3, NULL));
7319566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[0], 3, NULL));
7329566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Gradxy[0], &istart, &iend));
733c4762a1bSJed Brown 
734c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
735c4762a1bSJed Brown     iblock = i / user->mx;
736c4762a1bSJed Brown     j      = iblock * user->mx + ((i + user->mx - 1) % user->mx);
7379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
738c4762a1bSJed Brown     j = iblock * user->mx + ((i + 1) % user->mx);
7399566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
7409566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &i, &zero, INSERT_VALUES));
741c4762a1bSJed Brown   }
7429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Gradxy[0], MAT_FINAL_ASSEMBLY));
7439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Gradxy[0], MAT_FINAL_ASSEMBLY));
744c4762a1bSJed Brown 
7459566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[1]));
7469566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Gradxy[1], PETSC_DECIDE, PETSC_DECIDE, n, n));
7479566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Gradxy[1]));
7489566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[1], 3, NULL, 3, NULL));
7499566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[1], 3, NULL));
7509566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Gradxy[1], &istart, &iend));
751c4762a1bSJed Brown 
752c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
753c4762a1bSJed Brown     j = (i + n - user->mx) % n;
7549566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
755c4762a1bSJed Brown     j = (j + 2 * user->mx) % n;
7569566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
7579566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &i, &zero, INSERT_VALUES));
758c4762a1bSJed Brown   }
7599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Gradxy[1], MAT_FINAL_ASSEMBLY));
7609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Gradxy[1], MAT_FINAL_ASSEMBLY));
761c4762a1bSJed Brown 
762c4762a1bSJed Brown   /* Generate Div matrix */
7639566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Grad, MAT_INITIAL_MATRIX, &user->Div));
7649566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Gradxy[0], MAT_INITIAL_MATRIX, &user->Divxy[0]));
7659566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Gradxy[1], MAT_INITIAL_MATRIX, &user->Divxy[1]));
766c4762a1bSJed Brown 
767c4762a1bSJed Brown   /* Off-diagonal averaging matrix */
7689566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->M));
7699566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->M, PETSC_DECIDE, PETSC_DECIDE, n, n));
7709566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->M));
7719566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->M, 4, NULL, 4, NULL));
7729566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->M, 4, NULL));
7739566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->M, &istart, &iend));
774c4762a1bSJed Brown 
775c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
776c4762a1bSJed Brown     /* kron(Id,Av) */
777c4762a1bSJed Brown     iblock = i / user->mx;
778c4762a1bSJed Brown     j      = iblock * user->mx + ((i + user->mx - 1) % user->mx);
7799566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
780c4762a1bSJed Brown     j = iblock * user->mx + ((i + 1) % user->mx);
7819566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
782c4762a1bSJed Brown 
783c4762a1bSJed Brown     /* kron(Av,Id) */
784c4762a1bSJed Brown     j = (i + user->mx) % n;
7859566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
786c4762a1bSJed Brown     j = (i + n - user->mx) % n;
7879566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
788c4762a1bSJed Brown   }
7899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->M, MAT_FINAL_ASSEMBLY));
7909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->M, MAT_FINAL_ASSEMBLY));
791c4762a1bSJed Brown 
792c4762a1bSJed Brown   /* Generate 2D grid */
7939566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &XX));
7949566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q));
7959566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(XX, PETSC_DECIDE, n));
7969566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->q, PETSC_DECIDE, n * user->nt));
7979566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(XX));
7989566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->q));
799c4762a1bSJed Brown 
8009566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &YY));
8019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &XXwork));
8029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &YYwork));
8039566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &user->d));
8049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &user->dwork));
805c4762a1bSJed Brown 
8069566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(XX, &istart, &iend));
807c4762a1bSJed Brown   for (linear_index = istart; linear_index < iend; linear_index++) {
808c4762a1bSJed Brown     i  = linear_index % user->mx;
809c4762a1bSJed Brown     j  = (linear_index - i) / user->mx;
810c4762a1bSJed Brown     vx = h * (i + 0.5);
811c4762a1bSJed Brown     vy = h * (j + 0.5);
8129566063dSJacob Faibussowitsch     PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES));
8139566063dSJacob Faibussowitsch     PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES));
814c4762a1bSJed Brown   }
815c4762a1bSJed Brown 
8169566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(XX));
8179566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(XX));
8189566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(YY));
8199566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(YY));
820c4762a1bSJed Brown 
821c4762a1bSJed Brown   /* Compute final density function yT
822c4762a1bSJed Brown      yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2))
823c4762a1bSJed Brown      yT = yT / (h^2*sum(yT)) */
8249566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX, XXwork));
8259566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY, YYwork));
826c4762a1bSJed Brown 
8279566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork, -0.25));
8289566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork, -0.25));
829c4762a1bSJed Brown 
8309566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
8319566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
832c4762a1bSJed Brown 
8339566063dSJacob Faibussowitsch   PetscCall(VecCopy(XXwork, user->dwork));
8349566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, 1.0, YYwork));
8359566063dSJacob Faibussowitsch   PetscCall(VecScale(user->dwork, -30.0));
8369566063dSJacob Faibussowitsch   PetscCall(VecExp(user->dwork));
8379566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->dwork, user->d));
838c4762a1bSJed Brown 
8399566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX, XXwork));
8409566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY, YYwork));
841c4762a1bSJed Brown 
8429566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork, -0.75));
8439566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork, -0.75));
844c4762a1bSJed Brown 
8459566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
8469566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
847c4762a1bSJed Brown 
8489566063dSJacob Faibussowitsch   PetscCall(VecCopy(XXwork, user->dwork));
8499566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, 1.0, YYwork));
8509566063dSJacob Faibussowitsch   PetscCall(VecScale(user->dwork, -30.0));
8519566063dSJacob Faibussowitsch   PetscCall(VecExp(user->dwork));
852c4762a1bSJed Brown 
8539566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->d, 1.0, user->dwork));
8549566063dSJacob Faibussowitsch   PetscCall(VecShift(user->d, 1.0));
8559566063dSJacob Faibussowitsch   PetscCall(VecSum(user->d, &sum));
8569566063dSJacob Faibussowitsch   PetscCall(VecScale(user->d, 1.0 / (h * h * sum)));
857c4762a1bSJed Brown 
858c4762a1bSJed Brown   /* Initial conditions of forward problem */
8599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &bc));
8609566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX, XXwork));
8619566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY, YYwork));
862c4762a1bSJed Brown 
8639566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork, -0.5));
8649566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork, -0.5));
865c4762a1bSJed Brown 
8669566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
8679566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
868c4762a1bSJed Brown 
8699566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(bc, 1.0, XXwork, YYwork));
8709566063dSJacob Faibussowitsch   PetscCall(VecScale(bc, -50.0));
8719566063dSJacob Faibussowitsch   PetscCall(VecExp(bc));
8729566063dSJacob Faibussowitsch   PetscCall(VecShift(bc, 1.0));
8739566063dSJacob Faibussowitsch   PetscCall(VecSum(bc, &sum));
8749566063dSJacob Faibussowitsch   PetscCall(VecScale(bc, 1.0 / (h * h * sum)));
875c4762a1bSJed Brown 
876c4762a1bSJed Brown   /* Create scatter from y to y_1,y_2,...,y_nt */
877c4762a1bSJed Brown   /*  TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */
8789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->yi_scatter));
8799566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &yi));
8809566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(yi, PETSC_DECIDE, user->mx * user->mx));
8819566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(yi));
8829566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi, user->nt, &user->yi));
8839566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi, user->nt, &user->yiwork));
8849566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi, user->nt, &user->ziwork));
885c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
8869566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->yi[i], &lo, &hi));
8879566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_yi));
8889566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + i * user->mx * user->mx, 1, &is_from_y));
8899566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->y, is_from_y, user->yi[i], is_to_yi, &user->yi_scatter[i]));
8909566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_yi));
8919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_y));
892c4762a1bSJed Brown   }
893c4762a1bSJed Brown 
894c4762a1bSJed Brown   /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */
895c4762a1bSJed Brown   /*  TODO: reorder for better parallelism */
8969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uxi_scatter));
8979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uyi_scatter));
8989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->ux_scatter));
8999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uy_scatter));
9009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * user->nt * user->mx * user->mx, &user->ui_scatter));
9019566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &uxi));
9029566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &ui));
9039566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(uxi, PETSC_DECIDE, user->mx * user->mx));
9049566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(ui, PETSC_DECIDE, 2 * user->mx * user->mx));
9059566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(uxi));
9069566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(ui));
9079566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxi));
9089566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyi));
9099566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxiwork));
9109566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyiwork));
9119566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ui, user->nt, &user->ui));
9129566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ui, user->nt, &user->uiwork));
913c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
9149566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi));
9159566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi));
9169566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u));
9179566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->u, is_from_u, user->uxi[i], is_to_uxi, &user->uxi_scatter[i]));
918c4762a1bSJed Brown 
9199566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uxi));
9209566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
921c4762a1bSJed Brown 
9229566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi));
9239566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi));
9249566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + (2 * i + 1) * user->mx * user->mx, 1, &is_from_u));
9259566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->u, is_from_u, user->uyi[i], is_to_uyi, &user->uyi_scatter[i]));
926c4762a1bSJed Brown 
9279566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uyi));
9289566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
929c4762a1bSJed Brown 
9309566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi));
9319566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi));
9329566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_from_u));
9339566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uxi[i], is_to_uxi, &user->ux_scatter[i]));
934c4762a1bSJed Brown 
9359566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uxi));
9369566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
937c4762a1bSJed Brown 
9389566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi));
9399566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi));
9409566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + user->mx * user->mx, 1, &is_from_u));
9419566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uyi[i], is_to_uyi, &user->uy_scatter[i]));
942c4762a1bSJed Brown 
9439566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uyi));
9449566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
945c4762a1bSJed Brown 
9469566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->ui[i], &lo, &hi));
9479566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi));
9489566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u));
9499566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->u, is_from_u, user->ui[i], is_to_uxi, &user->ui_scatter[i]));
950c4762a1bSJed Brown 
9519566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_uxi));
9529566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_u));
953c4762a1bSJed Brown   }
954c4762a1bSJed Brown 
955c4762a1bSJed Brown   /* RHS of forward problem */
9569566063dSJacob Faibussowitsch   PetscCall(MatMult(user->M, bc, user->yiwork[0]));
95748a46eb9SPierre Jolivet   for (i = 1; i < user->nt; i++) PetscCall(VecSet(user->yiwork[i], 0.0));
9589566063dSJacob Faibussowitsch   PetscCall(Gather_yi(user->q, user->yiwork, user->yi_scatter, user->nt));
959c4762a1bSJed Brown 
960c4762a1bSJed Brown   /* Compute true velocity field utrue */
9619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u, &user->utrue));
962c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
9639566063dSJacob Faibussowitsch     PetscCall(VecCopy(YY, user->uxi[i]));
9649566063dSJacob Faibussowitsch     PetscCall(VecScale(user->uxi[i], 150.0 * i * user->ht));
9659566063dSJacob Faibussowitsch     PetscCall(VecCopy(XX, user->uyi[i]));
9669566063dSJacob Faibussowitsch     PetscCall(VecShift(user->uyi[i], -10.0));
9679566063dSJacob Faibussowitsch     PetscCall(VecScale(user->uyi[i], 15.0 * i * user->ht));
968c4762a1bSJed Brown   }
9699566063dSJacob Faibussowitsch   PetscCall(Gather_uxi_uyi(user->utrue, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
970c4762a1bSJed Brown 
971c4762a1bSJed Brown   /* Initial guess and reference model */
9729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->utrue, &user->ur));
973c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
9749566063dSJacob Faibussowitsch     PetscCall(VecCopy(XX, user->uxi[i]));
9759566063dSJacob Faibussowitsch     PetscCall(VecShift(user->uxi[i], i * user->ht));
9769566063dSJacob Faibussowitsch     PetscCall(VecCopy(YY, user->uyi[i]));
9779566063dSJacob Faibussowitsch     PetscCall(VecShift(user->uyi[i], -i * user->ht));
978c4762a1bSJed Brown   }
9799566063dSJacob Faibussowitsch   PetscCall(Gather_uxi_uyi(user->ur, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
980c4762a1bSJed Brown 
981c4762a1bSJed Brown   /* Generate regularization matrix L */
9829566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->LT));
9839566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->LT, PETSC_DECIDE, PETSC_DECIDE, 2 * n * user->nt, n * user->nt));
9849566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->LT));
9859566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->LT, 1, NULL, 1, NULL));
9869566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->LT, 1, NULL));
9879566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->LT, &istart, &iend));
988c4762a1bSJed Brown 
989c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
990c4762a1bSJed Brown     iblock = (i + n) / (2 * n);
991c4762a1bSJed Brown     j      = i - iblock * n;
9929566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->LT, 1, &i, 1, &j, &user->gamma, INSERT_VALUES));
993c4762a1bSJed Brown   }
994c4762a1bSJed Brown 
9959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->LT, MAT_FINAL_ASSEMBLY));
9969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->LT, MAT_FINAL_ASSEMBLY));
997c4762a1bSJed Brown 
9989566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->LT, MAT_INITIAL_MATRIX, &user->L));
999c4762a1bSJed Brown 
1000c4762a1bSJed Brown   /* Build work vectors and matrices */
10019566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork));
10029566063dSJacob Faibussowitsch   PetscCall(VecSetType(user->lwork, VECMPI));
10039566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, user->m));
10049566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->lwork));
1005c4762a1bSJed Brown 
10069566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork));
1007c4762a1bSJed Brown 
10089566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->y, &user->ywork));
10099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u, &user->uwork));
10109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u, &user->vwork));
10119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u, &user->js_diag));
10129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->c, &user->cwork));
1013c4762a1bSJed Brown 
1014c4762a1bSJed Brown   /* Create matrix-free shell user->Js for computing A*x */
10159566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->Js));
10169566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (void (*)(void))StateMatMult));
10179566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate));
10189566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose));
10199566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal));
1020c4762a1bSJed Brown 
1021c4762a1bSJed Brown   /* Diagonal blocks of user->Js */
10229566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlock));
10239566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (void (*)(void))StateMatBlockMult));
10249566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockMultTranspose));
1025c4762a1bSJed Brown 
1026c4762a1bSJed Brown   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1027c4762a1bSJed Brown      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1028c4762a1bSJed Brown      This is an SOR preconditioner for user->JsBlock. */
10299566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlockPrec));
10309566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT, (void (*)(void))StateMatBlockPrecMult));
10319566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockPrecMultTranspose));
1032c4762a1bSJed Brown 
1033c4762a1bSJed Brown   /* Create a matrix-free shell user->Jd for computing B*x */
10349566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->n - user->m, user, &user->Jd));
10359566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (void (*)(void))DesignMatMult));
10369566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (void (*)(void))DesignMatMultTranspose));
1037c4762a1bSJed Brown 
1038c4762a1bSJed Brown   /* User-defined routines for computing user->Js\x and user->Js^T\x*/
10399566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsInv));
10409566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (void (*)(void))StateMatInvMult));
10419566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatInvTransposeMult));
1042c4762a1bSJed Brown 
1043c4762a1bSJed Brown   /* Build matrices for SOR preconditioner */
10449566063dSJacob Faibussowitsch   PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
10459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(5 * n, &user->C));
10469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * n, &user->Cwork));
1047c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
10489566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(user->Divxy[0], MAT_COPY_VALUES, &user->C[i]));
10499566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(user->Divxy[1], MAT_COPY_VALUES, &user->Cwork[i]));
1050c4762a1bSJed Brown 
10519566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i]));
10529566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i]));
10539566063dSJacob Faibussowitsch     PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], DIFFERENT_NONZERO_PATTERN));
10549566063dSJacob Faibussowitsch     PetscCall(MatScale(user->C[i], user->ht));
10559566063dSJacob Faibussowitsch     PetscCall(MatShift(user->C[i], 1.0));
1056c4762a1bSJed Brown   }
1057c4762a1bSJed Brown 
1058c4762a1bSJed Brown   /* Solver options and tolerances */
10599566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver));
10609566063dSJacob Faibussowitsch   PetscCall(KSPSetType(user->solver, KSPGMRES));
10619566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->JsBlockPrec));
10629566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, 1e3, 500));
10639566063dSJacob Faibussowitsch   /* PetscCall(KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500)); */
10649566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(user->solver, &user->prec));
10659566063dSJacob Faibussowitsch   PetscCall(PCSetType(user->prec, PCSHELL));
1066c4762a1bSJed Brown 
10679566063dSJacob Faibussowitsch   PetscCall(PCShellSetApply(user->prec, StateMatBlockPrecMult));
10689566063dSJacob Faibussowitsch   PetscCall(PCShellSetApplyTranspose(user->prec, StateMatBlockPrecMultTranspose));
10699566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(user->prec, user));
1070c4762a1bSJed Brown 
1071c4762a1bSJed Brown   /* Compute true state function yt given ut */
10729566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ytrue));
10739566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->ytrue, PETSC_DECIDE, n * user->nt));
10749566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->ytrue));
1075c4762a1bSJed Brown   user->c_formed = PETSC_TRUE;
10769566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->utrue, user->u)); /*  Set u=utrue temporarily for StateMatInv */
10779566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ytrue, 0.0));      /*  Initial guess */
10789566063dSJacob Faibussowitsch   PetscCall(StateMatInvMult(user->Js, user->q, user->ytrue));
10799566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->ur, user->u)); /*  Reset u=ur */
1080c4762a1bSJed Brown 
1081c4762a1bSJed Brown   /* Initial guess y0 for state given u0 */
10829566063dSJacob Faibussowitsch   PetscCall(StateMatInvMult(user->Js, user->q, user->y));
1083c4762a1bSJed Brown 
1084c4762a1bSJed Brown   /* Data discretization */
10859566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Q));
10869566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Q, PETSC_DECIDE, PETSC_DECIDE, user->mx * user->mx, user->m));
10879566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Q));
10889566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Q, 0, NULL, 1, NULL));
10899566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Q, 1, NULL));
1090c4762a1bSJed Brown 
10919566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Q, &istart, &iend));
1092c4762a1bSJed Brown 
1093c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
1094c4762a1bSJed Brown     j = i + user->m - user->mx * user->mx;
10959566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &one, INSERT_VALUES));
1096c4762a1bSJed Brown   }
1097c4762a1bSJed Brown 
10989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Q, MAT_FINAL_ASSEMBLY));
10999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Q, MAT_FINAL_ASSEMBLY));
1100c4762a1bSJed Brown 
11019566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Q, MAT_INITIAL_MATRIX, &user->QT));
1102c4762a1bSJed Brown 
11039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XX));
11049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YY));
11059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XXwork));
11069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YYwork));
11079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&yi));
11089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&uxi));
11099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ui));
11109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bc));
1111c4762a1bSJed Brown 
1112c4762a1bSJed Brown   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11139566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(user->solver));
11143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1115c4762a1bSJed Brown }
1116c4762a1bSJed Brown 
1117d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicDestroy(AppCtx *user)
1118d71ae5a4SJacob Faibussowitsch {
1119c4762a1bSJed Brown   PetscInt i;
1120c4762a1bSJed Brown 
1121c4762a1bSJed Brown   PetscFunctionBegin;
11229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Q));
11239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->QT));
11249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Div));
11259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Divwork));
11269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Grad));
11279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->L));
11289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->LT));
11299566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&user->solver));
11309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Js));
11319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Jd));
11329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsBlockPrec));
11339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsInv));
11349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsBlock));
11359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Divxy[0]));
11369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Divxy[1]));
11379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Gradxy[0]));
11389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Gradxy[1]));
11399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->M));
1140c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
11419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->C[i]));
11429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->Cwork[i]));
1143c4762a1bSJed Brown   }
11449566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->C));
11459566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->Cwork));
11469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->u));
11479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->uwork));
11489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->vwork));
11499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->utrue));
11509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->y));
11519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ywork));
11529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ytrue));
11539566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->yi));
11549566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->yiwork));
11559566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->ziwork));
11569566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->uxi));
11579566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->uyi));
11589566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->uxiwork));
11599566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->uyiwork));
11609566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->ui));
11619566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->uiwork));
11629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->c));
11639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->cwork));
11649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ur));
11659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->q));
11669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->d));
11679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->dwork));
11689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->lwork));
11699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->js_diag));
11709566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->s_is));
11719566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->d_is));
11729566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->state_scatter));
11739566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->design_scatter));
1174c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
11759566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->uxi_scatter[i]));
11769566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->uyi_scatter[i]));
11779566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->ux_scatter[i]));
11789566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->uy_scatter[i]));
11799566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->ui_scatter[i]));
11809566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
1181c4762a1bSJed Brown   }
11829566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->uxi_scatter));
11839566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->uyi_scatter));
11849566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->ux_scatter));
11859566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->uy_scatter));
11869566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->ui_scatter));
11879566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->yi_scatter));
11883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1189c4762a1bSJed Brown }
1190c4762a1bSJed Brown 
1191d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1192d71ae5a4SJacob Faibussowitsch {
1193c4762a1bSJed Brown   Vec       X;
1194c4762a1bSJed Brown   PetscReal unorm, ynorm;
1195c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ptr;
1196c4762a1bSJed Brown 
1197c4762a1bSJed Brown   PetscFunctionBegin;
11989566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao, &X));
11999566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
12009566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue));
12019566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork, -1.0, user->utrue));
12029566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->uwork, NORM_2, &unorm));
12039566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->ywork, NORM_2, &ynorm));
12049566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm));
12053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1206c4762a1bSJed Brown }
1207c4762a1bSJed Brown 
1208c4762a1bSJed Brown /*TEST
1209c4762a1bSJed Brown 
1210c4762a1bSJed Brown    build:
1211c4762a1bSJed Brown       requires: !complex
1212c4762a1bSJed Brown 
1213c4762a1bSJed Brown    test:
1214c4762a1bSJed Brown       requires: !single
1215*10978b7dSBarry Smith       args: -tao_monitor_constraint_norm -tao_max_funcs 10 -tao_type lcl -tao_gatol 1.e-5
1216c4762a1bSJed Brown 
1217c4762a1bSJed Brown    test:
1218c4762a1bSJed Brown       suffix: guess_pod
1219c4762a1bSJed Brown       requires: !single
1220*10978b7dSBarry Smith       args: -tao_monitor_constraint_norm -tao_max_funcs 10 -tao_type lcl -ksp_guess_type pod -tao_gatol 1.e-5
1221c4762a1bSJed Brown 
1222c4762a1bSJed Brown TEST*/
1223