xref: /petsc/src/tao/pde_constrained/tutorials/parabolic.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown #include <petsc/private/taoimpl.h>
2c4762a1bSJed Brown 
3c4762a1bSJed Brown typedef struct {
4c4762a1bSJed Brown   PetscInt  n;            /*  Number of variables */
5c4762a1bSJed Brown   PetscInt  m;            /*  Number of constraints per time step */
6c4762a1bSJed Brown   PetscInt  mx;           /*  grid points in each direction */
7c4762a1bSJed Brown   PetscInt  nt;           /*  Number of time steps; as of now, must be divisible by 8 */
8c4762a1bSJed Brown   PetscInt  ndata;        /*  Number of data points per sample */
9c4762a1bSJed Brown   PetscInt  ns;           /*  Number of samples */
10c4762a1bSJed Brown   PetscInt *sample_times; /*  Times of samples */
11c4762a1bSJed Brown   IS        s_is;
12c4762a1bSJed Brown   IS        d_is;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   VecScatter  state_scatter;
15c4762a1bSJed Brown   VecScatter  design_scatter;
16c4762a1bSJed Brown   VecScatter *yi_scatter;
17c4762a1bSJed Brown   VecScatter *di_scatter;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   Mat       Js, Jd, JsBlockPrec, JsInv, JsBlock;
20c4762a1bSJed Brown   PetscBool jformed, dsg_formed;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   PetscReal alpha; /*  Regularization parameter */
23c4762a1bSJed Brown   PetscReal beta;  /*  Weight attributed to ||u||^2 in regularization functional */
24c4762a1bSJed Brown   PetscReal noise; /*  Amount of noise to add to data */
25c4762a1bSJed Brown   PetscReal ht;    /*  Time step */
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   Mat Qblock, QblockT;
28c4762a1bSJed Brown   Mat L, LT;
29c4762a1bSJed Brown   Mat Div, Divwork;
30c4762a1bSJed Brown   Mat Grad;
31c4762a1bSJed Brown   Mat Av, Avwork, AvT;
32c4762a1bSJed Brown   Mat DSG;
33c4762a1bSJed Brown   Vec q;
34c4762a1bSJed Brown   Vec ur; /*  reference */
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   Vec  d;
37c4762a1bSJed Brown   Vec  dwork;
38c4762a1bSJed Brown   Vec *di;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   Vec y; /*  state variables */
41c4762a1bSJed Brown   Vec ywork;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   Vec  ytrue;
44c4762a1bSJed Brown   Vec *yi, *yiwork;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   Vec u; /*  design variables */
47c4762a1bSJed Brown   Vec uwork;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   Vec utrue;
50c4762a1bSJed Brown   Vec js_diag;
51c4762a1bSJed Brown   Vec c; /*  constraint vector */
52c4762a1bSJed Brown   Vec cwork;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   Vec lwork;
55c4762a1bSJed Brown   Vec S;
56c4762a1bSJed Brown   Vec Rwork, Swork, Twork;
57c4762a1bSJed Brown   Vec Av_u;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   KSP solver;
60c4762a1bSJed Brown   PC  prec;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscInt ksp_its;
63c4762a1bSJed Brown   PetscInt ksp_its_initial;
64c4762a1bSJed Brown } AppCtx;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
67c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
68c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
69c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void *);
70c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void *);
71c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
72c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
73c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
74c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
75c4762a1bSJed Brown PetscErrorCode ParabolicInitialize(AppCtx *user);
76c4762a1bSJed Brown PetscErrorCode ParabolicDestroy(AppCtx *user);
77c4762a1bSJed Brown PetscErrorCode ParabolicMonitor(Tao, void *);
78c4762a1bSJed Brown 
79c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat, Vec, Vec);
80c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat, Vec, Vec);
81c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat, Vec, Vec);
82c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat, Vec);
83c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat, MatDuplicateOption, Mat *);
84c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat, Vec, Vec);
85c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat, Vec, Vec);
86c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC, Vec, Vec);
87c4762a1bSJed Brown 
88c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat, Vec, Vec);
89c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat, Vec, Vec);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown PetscErrorCode Gather_i(Vec, Vec *, VecScatter *, PetscInt);
92c4762a1bSJed Brown PetscErrorCode Scatter_i(Vec, Vec *, VecScatter *, PetscInt);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown static char help[] = "";
95c4762a1bSJed Brown 
969371c9d4SSatish Balay int main(int argc, char **argv) {
97c4762a1bSJed Brown   Vec      x, x0;
98c4762a1bSJed Brown   Tao      tao;
99c4762a1bSJed Brown   AppCtx   user;
100c4762a1bSJed Brown   IS       is_allstate, is_alldesign;
101c4762a1bSJed Brown   PetscInt lo, hi, hi2, lo2, ksp_old;
102c4762a1bSJed Brown   PetscInt ntests = 1;
103c4762a1bSJed Brown   PetscInt i;
104c4762a1bSJed Brown #if defined(PETSC_USE_LOG)
105c4762a1bSJed Brown   PetscLogStage stages[1];
106c4762a1bSJed Brown #endif
107c4762a1bSJed Brown 
108327415f7SBarry Smith   PetscFunctionBeginUser;
1099566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
110c4762a1bSJed Brown   user.mx = 8;
111d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "parabolic example", NULL);
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL));
113c4762a1bSJed Brown   user.nt = 8;
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-nt", "Number of time steps", "", user.nt, &user.nt, NULL));
115c4762a1bSJed Brown   user.ndata = 64;
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL));
117c4762a1bSJed Brown   user.ns = 8;
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ns", "Number of samples", "", user.ns, &user.ns, NULL));
119c4762a1bSJed Brown   user.alpha = 1.0;
1209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL));
121c4762a1bSJed Brown   user.beta = 0.01;
1229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-beta", "Weight attributed to ||u||^2 in regularization functional", "", user.beta, &user.beta, NULL));
123c4762a1bSJed Brown   user.noise = 0.01;
1249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-noise", "Amount of noise to add to data", "", user.noise, &user.noise, NULL));
1259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL));
126d0609cedSBarry Smith   PetscOptionsEnd();
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   user.m  = user.mx * user.mx * user.mx; /*  number of constraints per time step */
129c4762a1bSJed Brown   user.n  = user.m * (user.nt + 1);      /*  number of variables */
130c4762a1bSJed Brown   user.ht = (PetscReal)1 / user.nt;
131c4762a1bSJed Brown 
1329566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.u));
1339566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.y));
1349566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.c));
1359566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.u, PETSC_DECIDE, user.n - user.m * user.nt));
1369566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.y, PETSC_DECIDE, user.m * user.nt));
1379566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.c, PETSC_DECIDE, user.m * user.nt));
1389566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.u));
1399566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.y));
1409566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user.c));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* Create scatters for reduced spaces.
143c4762a1bSJed Brown      If the state vector y and design vector u are partitioned as
144c4762a1bSJed Brown      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
145c4762a1bSJed Brown      then the solution vector x is organized as
146c4762a1bSJed Brown      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
147c4762a1bSJed Brown      The index sets user.s_is and user.d_is correspond to the indices of the
148c4762a1bSJed Brown      state and design variables owned by the current processor.
149c4762a1bSJed Brown   */
1509566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user.y, &lo, &hi));
1539566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user.u, &lo2, &hi2));
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate));
1569566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user.s_is));
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign));
1599566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user.d_is));
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, hi - lo + hi2 - lo2, user.n));
1629566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, user.s_is, user.y, is_allstate, &user.state_scatter));
1659566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, user.d_is, user.u, is_alldesign, &user.design_scatter));
1669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_alldesign));
1679566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_allstate));
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1709566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1719566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOLCL));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* Set up initial vectors and matrices */
1749566063dSJacob Faibussowitsch   PetscCall(ParabolicInitialize(&user));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(Gather(x, user.y, user.state_scatter, user.u, user.design_scatter));
1779566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &x0));
1789566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, x0));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /* Set solution vector with an initial guess */
1819566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
1829566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(tao, FormFunction, &user));
1839566063dSJacob Faibussowitsch   PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
1849566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
185c4762a1bSJed Brown 
1869566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, &user));
1879566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
188c4762a1bSJed Brown 
1899566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
1909566063dSJacob Faibussowitsch   PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1939566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Trials", &stages[0]));
1949566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stages[0]));
195c4762a1bSJed Brown   user.ksp_its_initial = user.ksp_its;
196c4762a1bSJed Brown   ksp_old              = user.ksp_its;
197c4762a1bSJed Brown   for (i = 0; i < ntests; i++) {
1989566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
19963a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its - ksp_old));
2009566063dSJacob Faibussowitsch     PetscCall(VecCopy(x0, x));
2019566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(tao, x));
202c4762a1bSJed Brown   }
2039566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
2049566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)x));
2059566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: "));
20663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial));
20763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total KSP iterations over %" PetscInt_FMT " trial(s): ", ntests));
20863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its));
2099566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations per trial: "));
21063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", (user.ksp_its - user.ksp_its_initial) / ntests));
211c4762a1bSJed Brown 
2129566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x0));
2159566063dSJacob Faibussowitsch   PetscCall(ParabolicDestroy(&user));
216c4762a1bSJed Brown 
2179566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
218b122ec5aSJacob Faibussowitsch   return 0;
219c4762a1bSJed Brown }
220c4762a1bSJed Brown /* ------------------------------------------------------------------- */
221c4762a1bSJed Brown /*
222c4762a1bSJed Brown    dwork = Qy - d
223c4762a1bSJed Brown    lwork = L*(u-ur)
224c4762a1bSJed Brown    f = 1/2 * (dwork.dork + alpha*lwork.lwork)
225c4762a1bSJed Brown */
2269371c9d4SSatish Balay PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) {
227c4762a1bSJed Brown   PetscReal d1 = 0, d2 = 0;
228c4762a1bSJed Brown   PetscInt  i, j;
229c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ptr;
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   PetscFunctionBegin;
2329566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2339566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
234c4762a1bSJed Brown   for (j = 0; j < user->ns; j++) {
235c4762a1bSJed Brown     i = user->sample_times[j];
2369566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j]));
237c4762a1bSJed Brown   }
2389566063dSJacob Faibussowitsch   PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns));
2399566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2409566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork, user->dwork, &d1));
241c4762a1bSJed Brown 
2429566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
2439566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L, user->uwork, user->lwork));
2449566063dSJacob Faibussowitsch   PetscCall(VecDot(user->lwork, user->lwork, &d2));
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha * d2);
247c4762a1bSJed Brown   PetscFunctionReturn(0);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown /* ------------------------------------------------------------------- */
251c4762a1bSJed Brown /*
252c4762a1bSJed Brown     state: g_s = Q' *(Qy - d)
253c4762a1bSJed Brown     design: g_d = alpha*L'*L*(u-ur)
254c4762a1bSJed Brown */
2559371c9d4SSatish Balay PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) {
256c4762a1bSJed Brown   PetscInt i, j;
257c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   PetscFunctionBegin;
2609566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2619566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
262c4762a1bSJed Brown   for (j = 0; j < user->ns; j++) {
263c4762a1bSJed Brown     i = user->sample_times[j];
2649566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j]));
265c4762a1bSJed Brown   }
2669566063dSJacob Faibussowitsch   PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns));
2679566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2689566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->dwork, user->di, user->di_scatter, user->ns));
2699566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ywork, 0.0));
2709566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
271c4762a1bSJed Brown   for (j = 0; j < user->ns; j++) {
272c4762a1bSJed Brown     i = user->sample_times[j];
2739566063dSJacob Faibussowitsch     PetscCall(MatMult(user->QblockT, user->di[j], user->yiwork[i]));
274c4762a1bSJed Brown   }
2759566063dSJacob Faibussowitsch   PetscCall(Gather_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
276c4762a1bSJed Brown 
2779566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
2789566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L, user->uwork, user->lwork));
2799566063dSJacob Faibussowitsch   PetscCall(MatMult(user->LT, user->lwork, user->uwork));
2809566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork, user->alpha));
2819566063dSJacob Faibussowitsch   PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
282c4762a1bSJed Brown   PetscFunctionReturn(0);
283c4762a1bSJed Brown }
284c4762a1bSJed Brown 
2859371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) {
286c4762a1bSJed Brown   PetscReal d1, d2;
287c4762a1bSJed Brown   PetscInt  i, j;
288c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ptr;
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   PetscFunctionBegin;
2919566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2929566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
293c4762a1bSJed Brown   for (j = 0; j < user->ns; j++) {
294c4762a1bSJed Brown     i = user->sample_times[j];
2959566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j]));
296c4762a1bSJed Brown   }
2979566063dSJacob Faibussowitsch   PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns));
2989566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2999566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork, user->dwork, &d1));
3009566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->dwork, user->di, user->di_scatter, user->ns));
3019566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ywork, 0.0));
3029566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
303c4762a1bSJed Brown   for (j = 0; j < user->ns; j++) {
304c4762a1bSJed Brown     i = user->sample_times[j];
3059566063dSJacob Faibussowitsch     PetscCall(MatMult(user->QblockT, user->di[j], user->yiwork[i]));
306c4762a1bSJed Brown   }
3079566063dSJacob Faibussowitsch   PetscCall(Gather_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
308c4762a1bSJed Brown 
3099566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
3109566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L, user->uwork, user->lwork));
3119566063dSJacob Faibussowitsch   PetscCall(VecDot(user->lwork, user->lwork, &d2));
3129566063dSJacob Faibussowitsch   PetscCall(MatMult(user->LT, user->lwork, user->uwork));
3139566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork, user->alpha));
314c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha * d2);
315c4762a1bSJed Brown 
3169566063dSJacob Faibussowitsch   PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
317c4762a1bSJed Brown   PetscFunctionReturn(0);
318c4762a1bSJed Brown }
319c4762a1bSJed Brown 
320c4762a1bSJed Brown /* ------------------------------------------------------------------- */
321c4762a1bSJed Brown /* A
322c4762a1bSJed Brown MatShell object
323c4762a1bSJed Brown */
3249371c9d4SSatish Balay PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr) {
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));
3299566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork, 0));
3309566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork, -1.0, user->u));
3319566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
3329566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
3339566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->Av_u, user->Swork));
3349566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Swork));
3359566063dSJacob Faibussowitsch   PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
3369566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork));
337c4762a1bSJed Brown   if (user->dsg_formed) {
3389566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(user->DSG));
339c4762a1bSJed Brown   } else {
3409566063dSJacob Faibussowitsch     PetscCall(MatMatMult(user->Divwork, user->Grad, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DSG));
341c4762a1bSJed Brown     user->dsg_formed = PETSC_TRUE;
342c4762a1bSJed Brown   }
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   /* B = speye(nx^3) + ht*DSG; */
3459566063dSJacob Faibussowitsch   PetscCall(MatScale(user->DSG, user->ht));
3469566063dSJacob Faibussowitsch   PetscCall(MatShift(user->DSG, 1.0));
347c4762a1bSJed Brown   PetscFunctionReturn(0);
348c4762a1bSJed Brown }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown /* ------------------------------------------------------------------- */
351c4762a1bSJed Brown /* B */
3529371c9d4SSatish Balay PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr) {
353c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   PetscFunctionBegin;
3569566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
357c4762a1bSJed Brown   PetscFunctionReturn(0);
358c4762a1bSJed Brown }
359c4762a1bSJed Brown 
3609371c9d4SSatish Balay PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y) {
361c4762a1bSJed Brown   PetscInt i;
362c4762a1bSJed Brown   AppCtx  *user;
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   PetscFunctionBegin;
3659566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
3669566063dSJacob Faibussowitsch   PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
3679566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
368c4762a1bSJed Brown   for (i = 1; i < user->nt; i++) {
3699566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3709566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i - 1]));
371c4762a1bSJed Brown   }
3729566063dSJacob Faibussowitsch   PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
373c4762a1bSJed Brown   PetscFunctionReturn(0);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown 
3769371c9d4SSatish Balay PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y) {
377c4762a1bSJed Brown   PetscInt i;
378c4762a1bSJed Brown   AppCtx  *user;
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   PetscFunctionBegin;
3819566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
3829566063dSJacob Faibussowitsch   PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
383c4762a1bSJed Brown   for (i = 0; i < user->nt - 1; i++) {
3849566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3859566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i + 1]));
386c4762a1bSJed Brown   }
387c4762a1bSJed Brown   i = user->nt - 1;
3889566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3899566063dSJacob Faibussowitsch   PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
390c4762a1bSJed Brown   PetscFunctionReturn(0);
391c4762a1bSJed Brown }
392c4762a1bSJed Brown 
3939371c9d4SSatish Balay PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y) {
394c4762a1bSJed Brown   AppCtx *user;
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   PetscFunctionBegin;
3979566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
3989566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Grad, X, user->Swork));
3999566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u));
4009566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Div, user->Swork, Y));
4019566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Y, user->ht, X));
402c4762a1bSJed Brown   PetscFunctionReturn(0);
403c4762a1bSJed Brown }
404c4762a1bSJed Brown 
4059371c9d4SSatish Balay PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y) {
406c4762a1bSJed Brown   PetscInt i;
407c4762a1bSJed Brown   AppCtx  *user;
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   PetscFunctionBegin;
4109566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   /* sdiag(1./v) */
4139566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork, 0));
4149566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork, -1.0, user->u));
4159566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
416c4762a1bSJed Brown 
417c4762a1bSJed Brown   /* sdiag(1./((Av*(1./v)).^2)) */
4189566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av, user->uwork, user->Swork));
4199566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->Swork, user->Swork, user->Swork));
4209566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Swork));
421c4762a1bSJed Brown 
422c4762a1bSJed Brown   /* (Av * (sdiag(1./v) * b)) */
4239566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork, user->uwork, X));
4249566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av, user->uwork, user->Twork));
425c4762a1bSJed Brown 
426c4762a1bSJed Brown   /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
4279566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork));
428c4762a1bSJed Brown 
4299566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
430c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
431c4762a1bSJed Brown     /* (sdiag(Grad*y(:,i)) */
4329566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad, user->yi[i], user->Twork));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown     /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
4359566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->Twork, user->Twork, user->Swork));
4369566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Div, user->Twork, user->yiwork[i]));
4379566063dSJacob Faibussowitsch     PetscCall(VecScale(user->yiwork[i], user->ht));
438c4762a1bSJed Brown   }
4399566063dSJacob Faibussowitsch   PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   PetscFunctionReturn(0);
442c4762a1bSJed Brown }
443c4762a1bSJed Brown 
4449371c9d4SSatish Balay PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y) {
445c4762a1bSJed Brown   PetscInt i;
446c4762a1bSJed Brown   AppCtx  *user;
447c4762a1bSJed Brown 
448c4762a1bSJed Brown   PetscFunctionBegin;
4499566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   /* sdiag(1./((Av*(1./v)).^2)) */
4529566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork, 0));
4539566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork, -1.0, user->u));
4549566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
4559566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av, user->uwork, user->Rwork));
4569566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->Rwork, user->Rwork, user->Rwork));
4579566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Rwork));
458c4762a1bSJed Brown 
4599566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, 0.0));
4609566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
4619566063dSJacob Faibussowitsch   PetscCall(Scatter_i(X, user->yiwork, user->yi_scatter, user->nt));
462c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
463c4762a1bSJed Brown     /* (Div' * b(:,i)) */
4649566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad, user->yiwork[i], user->Swork));
465c4762a1bSJed Brown 
466c4762a1bSJed Brown     /* sdiag(Grad*y(:,i)) */
4679566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad, user->yi[i], user->Twork));
468c4762a1bSJed Brown 
469c4762a1bSJed Brown     /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */
4709566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->Twork, user->Swork, user->Twork));
471c4762a1bSJed Brown 
472c4762a1bSJed Brown     /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */
4739566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->Twork, user->Rwork, user->Twork));
474c4762a1bSJed Brown 
475c4762a1bSJed Brown     /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
4769566063dSJacob Faibussowitsch     PetscCall(MatMult(user->AvT, user->Twork, user->yiwork[i]));
477c4762a1bSJed Brown 
478c4762a1bSJed Brown     /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
4799566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->yiwork[i], user->uwork, user->yiwork[i]));
4809566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, user->ht, user->yiwork[i]));
481c4762a1bSJed Brown   }
482c4762a1bSJed Brown   PetscFunctionReturn(0);
483c4762a1bSJed Brown }
484c4762a1bSJed Brown 
4859371c9d4SSatish Balay PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y) {
486c4762a1bSJed Brown   AppCtx *user;
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   PetscFunctionBegin;
4899566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(PC_shell, &user));
490c4762a1bSJed Brown 
491c4762a1bSJed Brown   if (user->dsg_formed) {
4929566063dSJacob Faibussowitsch     PetscCall(MatSOR(user->DSG, X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y));
493c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "DSG not formed");
494c4762a1bSJed Brown   PetscFunctionReturn(0);
495c4762a1bSJed Brown }
496c4762a1bSJed Brown 
4979371c9d4SSatish Balay PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y) {
498c4762a1bSJed Brown   AppCtx  *user;
499c4762a1bSJed Brown   PetscInt its, i;
500c4762a1bSJed Brown 
501c4762a1bSJed Brown   PetscFunctionBegin;
5029566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
503c4762a1bSJed Brown 
504c4762a1bSJed Brown   if (Y == user->ytrue) {
505c4762a1bSJed Brown     /* First solve is done with true solution to set up problem */
5069566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
507c4762a1bSJed Brown   } else {
5089566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
509c4762a1bSJed Brown   }
510c4762a1bSJed Brown 
5119566063dSJacob Faibussowitsch   PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
5129566063dSJacob Faibussowitsch   PetscCall(KSPSolve(user->solver, user->yi[0], user->yiwork[0]));
5139566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(user->solver, &its));
514c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
515c4762a1bSJed Brown 
516c4762a1bSJed Brown   for (i = 1; i < user->nt; i++) {
5179566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yi[i], 1.0, user->yiwork[i - 1]));
5189566063dSJacob Faibussowitsch     PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
5199566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(user->solver, &its));
520c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
521c4762a1bSJed Brown   }
5229566063dSJacob Faibussowitsch   PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
523c4762a1bSJed Brown   PetscFunctionReturn(0);
524c4762a1bSJed Brown }
525c4762a1bSJed Brown 
5269371c9d4SSatish Balay PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y) {
527c4762a1bSJed Brown   AppCtx  *user;
528c4762a1bSJed Brown   PetscInt its, i;
529c4762a1bSJed Brown 
530c4762a1bSJed Brown   PetscFunctionBegin;
5319566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
532c4762a1bSJed Brown 
5339566063dSJacob Faibussowitsch   PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
534c4762a1bSJed Brown 
535c4762a1bSJed Brown   i = user->nt - 1;
5369566063dSJacob Faibussowitsch   PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
537c4762a1bSJed Brown 
5389566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(user->solver, &its));
539c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
540c4762a1bSJed Brown 
541c4762a1bSJed Brown   for (i = user->nt - 2; i >= 0; i--) {
5429566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yi[i], 1.0, user->yiwork[i + 1]));
5439566063dSJacob Faibussowitsch     PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
544c4762a1bSJed Brown 
5459566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(user->solver, &its));
546c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
547c4762a1bSJed Brown   }
548c4762a1bSJed Brown 
5499566063dSJacob Faibussowitsch   PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
550c4762a1bSJed Brown   PetscFunctionReturn(0);
551c4762a1bSJed Brown }
552c4762a1bSJed Brown 
5539371c9d4SSatish Balay PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell) {
554c4762a1bSJed Brown   AppCtx *user;
555c4762a1bSJed Brown 
556c4762a1bSJed Brown   PetscFunctionBegin;
5579566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
558c4762a1bSJed Brown 
5599566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, new_shell));
5609566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT, (void (*)(void))StateMatMult));
5619566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate));
5629566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose));
5639566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*new_shell, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal));
564c4762a1bSJed Brown   PetscFunctionReturn(0);
565c4762a1bSJed Brown }
566c4762a1bSJed Brown 
5679371c9d4SSatish Balay PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X) {
568c4762a1bSJed Brown   AppCtx *user;
569c4762a1bSJed Brown 
570c4762a1bSJed Brown   PetscFunctionBegin;
5719566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell, &user));
5729566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->js_diag, X));
573c4762a1bSJed Brown   PetscFunctionReturn(0);
574c4762a1bSJed Brown }
575c4762a1bSJed Brown 
5769371c9d4SSatish Balay PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr) {
577c4762a1bSJed Brown   /* con = Ay - q, A = [B  0  0 ... 0;
578c4762a1bSJed Brown                        -I  B  0 ... 0;
579c4762a1bSJed Brown                         0 -I  B ... 0;
580c4762a1bSJed Brown                              ...     ;
581c4762a1bSJed Brown                         0    ... -I B]
582c4762a1bSJed Brown      B = ht * Div * Sigma * Grad + eye */
583c4762a1bSJed Brown   PetscInt i;
584c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
585c4762a1bSJed Brown 
586c4762a1bSJed Brown   PetscFunctionBegin;
5879566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
5889566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
5899566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
590c4762a1bSJed Brown   for (i = 1; i < user->nt; i++) {
5919566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
5929566063dSJacob Faibussowitsch     PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i - 1]));
593c4762a1bSJed Brown   }
5949566063dSJacob Faibussowitsch   PetscCall(Gather_i(C, user->yiwork, user->yi_scatter, user->nt));
5959566063dSJacob Faibussowitsch   PetscCall(VecAXPY(C, -1.0, user->q));
596c4762a1bSJed Brown   PetscFunctionReturn(0);
597c4762a1bSJed Brown }
598c4762a1bSJed Brown 
5999371c9d4SSatish Balay PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) {
600c4762a1bSJed Brown   PetscFunctionBegin;
6019566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6029566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6039566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
6049566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
605c4762a1bSJed Brown   PetscFunctionReturn(0);
606c4762a1bSJed Brown }
607c4762a1bSJed Brown 
6089371c9d4SSatish Balay PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) {
609c4762a1bSJed Brown   PetscInt i;
610c4762a1bSJed Brown 
611c4762a1bSJed Brown   PetscFunctionBegin;
612c4762a1bSJed Brown   for (i = 0; i < nt; i++) {
6139566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
6149566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
615c4762a1bSJed Brown   }
616c4762a1bSJed Brown   PetscFunctionReturn(0);
617c4762a1bSJed Brown }
618c4762a1bSJed Brown 
6199371c9d4SSatish Balay PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) {
620c4762a1bSJed Brown   PetscFunctionBegin;
6219566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6229566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6239566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
6249566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
625c4762a1bSJed Brown   PetscFunctionReturn(0);
626c4762a1bSJed Brown }
627c4762a1bSJed Brown 
6289371c9d4SSatish Balay PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) {
629c4762a1bSJed Brown   PetscInt i;
630c4762a1bSJed Brown 
631c4762a1bSJed Brown   PetscFunctionBegin;
632c4762a1bSJed Brown   for (i = 0; i < nt; i++) {
6339566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
6349566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
635c4762a1bSJed Brown   }
636c4762a1bSJed Brown   PetscFunctionReturn(0);
637c4762a1bSJed Brown }
638c4762a1bSJed Brown 
6399371c9d4SSatish Balay PetscErrorCode ParabolicInitialize(AppCtx *user) {
640c4762a1bSJed Brown   PetscInt    m, n, i, j, k, linear_index, istart, iend, iblock, lo, hi, lo2, hi2;
641c4762a1bSJed Brown   Vec         XX, YY, ZZ, XXwork, YYwork, ZZwork, UTwork, yi, di, bc;
642c4762a1bSJed Brown   PetscReal  *x, *y, *z;
643c4762a1bSJed Brown   PetscReal   h, stime;
644c4762a1bSJed Brown   PetscScalar hinv, neg_hinv, half = 0.5, sqrt_beta;
645c4762a1bSJed Brown   PetscInt    im, indx1, indx2, indy1, indy2, indz1, indz2, nx, ny, nz;
646c4762a1bSJed Brown   PetscReal   xri, yri, zri, xim, yim, zim, dx1, dx2, dy1, dy2, dz1, dz2, Dx, Dy, Dz;
647c4762a1bSJed Brown   PetscScalar v, vx, vy, vz;
648c4762a1bSJed Brown   IS          is_from_y, is_to_yi, is_from_d, is_to_di;
649c4762a1bSJed Brown   /* Data locations */
6509371c9d4SSatish Balay   PetscScalar xr[64] = {0.4970, 0.8498, 0.7814, 0.6268, 0.7782, 0.6402, 0.3617, 0.3160, 0.3610, 0.5298, 0.6987, 0.3331, 0.7962, 0.5596, 0.3866, 0.6774, 0.5407, 0.4518, 0.6702, 0.6061, 0.7580, 0.8997,
6519371c9d4SSatish Balay                         0.5198, 0.8326, 0.2138, 0.9198, 0.3000, 0.2833, 0.8288, 0.7076, 0.1820, 0.0728, 0.8447, 0.2367, 0.3239, 0.6413, 0.3114, 0.4731, 0.1192, 0.9273, 0.5724, 0.4331, 0.5136, 0.3547,
6529371c9d4SSatish Balay                         0.4413, 0.2602, 0.5698, 0.7278, 0.5261, 0.6230, 0.2454, 0.3948, 0.7479, 0.6582, 0.4660, 0.5594, 0.7574, 0.1143, 0.5900, 0.1065, 0.4260, 0.3294, 0.8276, 0.0756};
653c4762a1bSJed Brown 
6549371c9d4SSatish Balay   PetscScalar yr[64] = {0.7345, 0.9120, 0.9288, 0.7528, 0.4463, 0.4985, 0.2497, 0.6256, 0.3425, 0.9026, 0.6983, 0.4230, 0.7140, 0.2970, 0.4474, 0.8792, 0.6604, 0.2485, 0.7968, 0.6127, 0.1796, 0.2437,
6559371c9d4SSatish Balay                         0.5938, 0.6137, 0.3867, 0.5658, 0.4575, 0.1009, 0.0863, 0.3361, 0.0738, 0.3985, 0.6602, 0.1437, 0.0934, 0.5983, 0.5950, 0.0763, 0.0768, 0.2288, 0.5761, 0.1129, 0.3841, 0.6150,
6569371c9d4SSatish Balay                         0.6904, 0.6686, 0.1361, 0.4601, 0.4491, 0.3716, 0.1969, 0.6537, 0.6743, 0.6991, 0.4811, 0.5480, 0.1684, 0.4569, 0.6889, 0.8437, 0.3015, 0.2854, 0.8199, 0.2658};
657c4762a1bSJed Brown 
6589371c9d4SSatish Balay   PetscScalar zr[64] = {0.7668, 0.8573, 0.2654, 0.2719, 0.1060, 0.1311, 0.6232, 0.2295, 0.8009, 0.2147, 0.2119, 0.9325, 0.4473, 0.3600, 0.3374, 0.3819, 0.4066, 0.5801, 0.1673, 0.0959, 0.4638, 0.8236,
6599371c9d4SSatish Balay                         0.8800, 0.2939, 0.2028, 0.8262, 0.2706, 0.6276, 0.9085, 0.6443, 0.8241, 0.0712, 0.1824, 0.7789, 0.4389, 0.8415, 0.7055, 0.6639, 0.3653, 0.2078, 0.1987, 0.2297, 0.4321, 0.8115,
6609371c9d4SSatish Balay                         0.4915, 0.7764, 0.4657, 0.4627, 0.4569, 0.4232, 0.8514, 0.0674, 0.3227, 0.1055, 0.6690, 0.6313, 0.9226, 0.5461, 0.4126, 0.2364, 0.6096, 0.7042, 0.3914, 0.0711};
661c4762a1bSJed Brown 
662c4762a1bSJed Brown   PetscFunctionBegin;
6639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->mx, &x));
6649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->mx, &y));
6659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->mx, &z));
666c4762a1bSJed Brown   user->jformed    = PETSC_FALSE;
667c4762a1bSJed Brown   user->dsg_formed = PETSC_FALSE;
668c4762a1bSJed Brown 
669c4762a1bSJed Brown   n         = user->mx * user->mx * user->mx;
670c4762a1bSJed Brown   m         = 3 * user->mx * user->mx * (user->mx - 1);
671c4762a1bSJed Brown   sqrt_beta = PetscSqrtScalar(user->beta);
672c4762a1bSJed Brown 
673c4762a1bSJed Brown   user->ksp_its         = 0;
674c4762a1bSJed Brown   user->ksp_its_initial = 0;
675c4762a1bSJed Brown 
676c4762a1bSJed Brown   stime = (PetscReal)user->nt / user->ns;
6779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->ns, &user->sample_times));
6789371c9d4SSatish Balay   for (i = 0; i < user->ns; i++) { user->sample_times[i] = (PetscInt)(stime * ((PetscReal)i + 1.0) - 0.5); }
679c4762a1bSJed Brown 
6809566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &XX));
6819566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q));
6829566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(XX, PETSC_DECIDE, n));
6839566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->q, PETSC_DECIDE, n * user->nt));
6849566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(XX));
6859566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->q));
686c4762a1bSJed Brown 
6879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &YY));
6889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &ZZ));
6899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &XXwork));
6909566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &YYwork));
6919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &ZZwork));
6929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &UTwork));
6939566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &user->utrue));
6949566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX, &bc));
695c4762a1bSJed Brown 
696c4762a1bSJed Brown   /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
697c4762a1bSJed Brown   h        = 1.0 / user->mx;
698c4762a1bSJed Brown   hinv     = user->mx;
699c4762a1bSJed Brown   neg_hinv = -hinv;
700c4762a1bSJed Brown 
7019566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(XX, &istart, &iend));
702c4762a1bSJed Brown   for (linear_index = istart; linear_index < iend; linear_index++) {
703c4762a1bSJed Brown     i  = linear_index % user->mx;
704c4762a1bSJed Brown     j  = ((linear_index - i) / user->mx) % user->mx;
705c4762a1bSJed Brown     k  = ((linear_index - i) / user->mx - j) / user->mx;
706c4762a1bSJed Brown     vx = h * (i + 0.5);
707c4762a1bSJed Brown     vy = h * (j + 0.5);
708c4762a1bSJed Brown     vz = h * (k + 0.5);
7099566063dSJacob Faibussowitsch     PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES));
7109566063dSJacob Faibussowitsch     PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES));
7119566063dSJacob Faibussowitsch     PetscCall(VecSetValues(ZZ, 1, &linear_index, &vz, INSERT_VALUES));
712c4762a1bSJed Brown     if ((vx < 0.6) && (vx > 0.4) && (vy < 0.6) && (vy > 0.4) && (vy < 0.6) && (vz < 0.6) && (vz > 0.4)) {
713c4762a1bSJed Brown       v = 1000.0;
7149566063dSJacob Faibussowitsch       PetscCall(VecSetValues(bc, 1, &linear_index, &v, INSERT_VALUES));
715c4762a1bSJed Brown     }
716c4762a1bSJed Brown   }
717c4762a1bSJed Brown 
7189566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(XX));
7199566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(XX));
7209566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(YY));
7219566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(YY));
7229566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(ZZ));
7239566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(ZZ));
7249566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(bc));
7259566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(bc));
726c4762a1bSJed Brown 
727c4762a1bSJed Brown   /* Compute true parameter function
728c4762a1bSJed Brown      ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
7299566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX, XXwork));
7309566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY, YYwork));
7319566063dSJacob Faibussowitsch   PetscCall(VecCopy(ZZ, ZZwork));
732c4762a1bSJed Brown 
7339566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork, -0.5));
7349566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork, -0.5));
7359566063dSJacob Faibussowitsch   PetscCall(VecShift(ZZwork, -0.5));
736c4762a1bSJed Brown 
7379566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
7389566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
7399566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork));
740c4762a1bSJed Brown 
7419566063dSJacob Faibussowitsch   PetscCall(VecCopy(XXwork, user->utrue));
7429566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->utrue, 1.0, YYwork));
7439566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->utrue, 1.0, ZZwork));
7449566063dSJacob Faibussowitsch   PetscCall(VecScale(user->utrue, -10.0));
7459566063dSJacob Faibussowitsch   PetscCall(VecExp(user->utrue));
7469566063dSJacob Faibussowitsch   PetscCall(VecShift(user->utrue, 0.5));
747c4762a1bSJed Brown 
7489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XX));
7499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YY));
7509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ZZ));
7519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XXwork));
7529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YYwork));
7539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ZZwork));
7549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&UTwork));
755c4762a1bSJed Brown 
756c4762a1bSJed Brown   /* Initial guess and reference model */
7579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->utrue, &user->ur));
7589566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ur, 0.0));
759c4762a1bSJed Brown 
760c4762a1bSJed Brown   /* Generate Grad matrix */
7619566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad));
7629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, PETSC_DECIDE, m, n));
7639566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Grad));
7649566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Grad, 2, NULL, 2, NULL));
7659566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Grad, 2, NULL));
7669566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));
767c4762a1bSJed Brown 
768c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
769c4762a1bSJed Brown     if (i < m / 3) {
770c4762a1bSJed Brown       iblock = i / (user->mx - 1);
771c4762a1bSJed Brown       j      = iblock * user->mx + (i % (user->mx - 1));
7729566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
773c4762a1bSJed Brown       j = j + 1;
7749566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
775c4762a1bSJed Brown     }
776c4762a1bSJed Brown     if (i >= m / 3 && i < 2 * m / 3) {
777c4762a1bSJed Brown       iblock = (i - m / 3) / (user->mx * (user->mx - 1));
778c4762a1bSJed Brown       j      = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
7799566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
780c4762a1bSJed Brown       j = j + user->mx;
7819566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
782c4762a1bSJed Brown     }
783c4762a1bSJed Brown     if (i >= 2 * m / 3) {
784c4762a1bSJed Brown       j = i - 2 * m / 3;
7859566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
786c4762a1bSJed Brown       j = j + user->mx * user->mx;
7879566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
788c4762a1bSJed Brown     }
789c4762a1bSJed Brown   }
790c4762a1bSJed Brown 
7919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY));
7929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY));
793c4762a1bSJed Brown 
794c4762a1bSJed Brown   /* Generate arithmetic averaging matrix Av */
7959566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Av));
7969566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Av, PETSC_DECIDE, PETSC_DECIDE, m, n));
7979566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Av));
7989566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Av, 2, NULL, 2, NULL));
7999566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Av, 2, NULL));
8009566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Av, &istart, &iend));
801c4762a1bSJed Brown 
802c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
803c4762a1bSJed Brown     if (i < m / 3) {
804c4762a1bSJed Brown       iblock = i / (user->mx - 1);
805c4762a1bSJed Brown       j      = iblock * user->mx + (i % (user->mx - 1));
8069566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
807c4762a1bSJed Brown       j = j + 1;
8089566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
809c4762a1bSJed Brown     }
810c4762a1bSJed Brown     if (i >= m / 3 && i < 2 * m / 3) {
811c4762a1bSJed Brown       iblock = (i - m / 3) / (user->mx * (user->mx - 1));
812c4762a1bSJed Brown       j      = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
8139566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
814c4762a1bSJed Brown       j = j + user->mx;
8159566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
816c4762a1bSJed Brown     }
817c4762a1bSJed Brown     if (i >= 2 * m / 3) {
818c4762a1bSJed Brown       j = i - 2 * m / 3;
8199566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
820c4762a1bSJed Brown       j = j + user->mx * user->mx;
8219566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
822c4762a1bSJed Brown     }
823c4762a1bSJed Brown   }
824c4762a1bSJed Brown 
8259566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Av, MAT_FINAL_ASSEMBLY));
8269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Av, MAT_FINAL_ASSEMBLY));
827c4762a1bSJed Brown 
828c4762a1bSJed Brown   /* Generate transpose of averaging matrix Av */
8299566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Av, MAT_INITIAL_MATRIX, &user->AvT));
830c4762a1bSJed Brown 
8319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->L));
8329566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->L, PETSC_DECIDE, PETSC_DECIDE, m + n, n));
8339566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->L));
8349566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->L, 2, NULL, 2, NULL));
8359566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->L, 2, NULL));
8369566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->L, &istart, &iend));
837c4762a1bSJed Brown 
838c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
839c4762a1bSJed Brown     if (i < m / 3) {
840c4762a1bSJed Brown       iblock = i / (user->mx - 1);
841c4762a1bSJed Brown       j      = iblock * user->mx + (i % (user->mx - 1));
8429566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
843c4762a1bSJed Brown       j = j + 1;
8449566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
845c4762a1bSJed Brown     }
846c4762a1bSJed Brown     if (i >= m / 3 && i < 2 * m / 3) {
847c4762a1bSJed Brown       iblock = (i - m / 3) / (user->mx * (user->mx - 1));
848c4762a1bSJed Brown       j      = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
8499566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
850c4762a1bSJed Brown       j = j + user->mx;
8519566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
852c4762a1bSJed Brown     }
853c4762a1bSJed Brown     if (i >= 2 * m / 3 && i < m) {
854c4762a1bSJed Brown       j = i - 2 * m / 3;
8559566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
856c4762a1bSJed Brown       j = j + user->mx * user->mx;
8579566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
858c4762a1bSJed Brown     }
859c4762a1bSJed Brown     if (i >= m) {
860c4762a1bSJed Brown       j = i - m;
8619566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &sqrt_beta, INSERT_VALUES));
862c4762a1bSJed Brown     }
863c4762a1bSJed Brown   }
864c4762a1bSJed Brown 
8659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->L, MAT_FINAL_ASSEMBLY));
8669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->L, MAT_FINAL_ASSEMBLY));
867c4762a1bSJed Brown 
8689566063dSJacob Faibussowitsch   PetscCall(MatScale(user->L, PetscPowScalar(h, 1.5)));
869c4762a1bSJed Brown 
870c4762a1bSJed Brown   /* Generate Div matrix */
8719566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Grad, MAT_INITIAL_MATRIX, &user->Div));
872c4762a1bSJed Brown 
873c4762a1bSJed Brown   /* Build work vectors and matrices */
8749566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->S));
8759566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->S, PETSC_DECIDE, user->mx * user->mx * (user->mx - 1) * 3));
8769566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->S));
877c4762a1bSJed Brown 
8789566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork));
8799566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, m + user->mx * user->mx * user->mx));
8809566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->lwork));
881c4762a1bSJed Brown 
8829566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork));
8839566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(user->Av, MAT_SHARE_NONZERO_PATTERN, &user->Avwork));
884c4762a1bSJed Brown 
8859566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->d));
8869566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->d, PETSC_DECIDE, user->ndata * user->nt));
8879566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->d));
888c4762a1bSJed Brown 
8899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S, &user->Swork));
8909566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S, &user->Av_u));
8919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S, &user->Twork));
8929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S, &user->Rwork));
8939566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->y, &user->ywork));
8949566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u, &user->uwork));
8959566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u, &user->js_diag));
8969566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->c, &user->cwork));
8979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->d, &user->dwork));
898c4762a1bSJed Brown 
899c4762a1bSJed Brown   /* Create matrix-free shell user->Js for computing A*x */
9009566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m * user->nt, user, &user->Js));
9019566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (void (*)(void))StateMatMult));
9029566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate));
9039566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose));
9049566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal));
905c4762a1bSJed Brown 
906c4762a1bSJed Brown   /* Diagonal blocks of user->Js */
9079566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsBlock));
9089566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (void (*)(void))StateMatBlockMult));
909c4762a1bSJed Brown   /* Blocks are symmetric */
9109566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockMult));
911c4762a1bSJed Brown 
912c4762a1bSJed Brown   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
913c4762a1bSJed Brown      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
914c4762a1bSJed Brown      This is an SSOR preconditioner for user->JsBlock. */
9159566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsBlockPrec));
9169566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT, (void (*)(void))StateMatBlockPrecMult));
917c4762a1bSJed Brown   /* JsBlockPrec is symmetric */
9189566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockPrecMult));
919e9b2f0d4SBarry Smith   PetscCall(MatSetOption(user->JsBlockPrec, MAT_SYMMETRIC, PETSC_TRUE));
9209566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->JsBlockPrec, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
921c4762a1bSJed Brown 
922c4762a1bSJed Brown   /* Create a matrix-free shell user->Jd for computing B*x */
9239566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m, user, &user->Jd));
9249566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (void (*)(void))DesignMatMult));
9259566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (void (*)(void))DesignMatMultTranspose));
926c4762a1bSJed Brown 
927c4762a1bSJed Brown   /* User-defined routines for computing user->Js\x and user->Js^T\x*/
9289566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m * user->nt, user, &user->JsInv));
9299566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (void (*)(void))StateMatInvMult));
9309566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatInvTransposeMult));
931c4762a1bSJed Brown 
932c4762a1bSJed Brown   /* Solver options and tolerances */
9339566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver));
9349566063dSJacob Faibussowitsch   PetscCall(KSPSetType(user->solver, KSPCG));
9359566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->JsBlockPrec));
9369566063dSJacob Faibussowitsch   PetscCall(KSPSetInitialGuessNonzero(user->solver, PETSC_FALSE));
9379566063dSJacob Faibussowitsch   PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, 1e3, 500));
9389566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(user->solver));
9399566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(user->solver, &user->prec));
9409566063dSJacob Faibussowitsch   PetscCall(PCSetType(user->prec, PCSHELL));
941c4762a1bSJed Brown 
9429566063dSJacob Faibussowitsch   PetscCall(PCShellSetApply(user->prec, StateMatBlockPrecMult));
9439566063dSJacob Faibussowitsch   PetscCall(PCShellSetApplyTranspose(user->prec, StateMatBlockPrecMult));
9449566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(user->prec, user));
945c4762a1bSJed Brown 
946c4762a1bSJed Brown   /* Create scatter from y to y_1,y_2,...,y_nt */
9479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->nt * user->m, &user->yi_scatter));
9489566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &yi));
9499566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(yi, PETSC_DECIDE, user->mx * user->mx * user->mx));
9509566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(yi));
9519566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi, user->nt, &user->yi));
9529566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(yi, user->nt, &user->yiwork));
953c4762a1bSJed Brown 
9549566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->y, &lo2, &hi2));
955c4762a1bSJed Brown   istart = 0;
956c4762a1bSJed Brown   for (i = 0; i < user->nt; i++) {
9579566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->yi[i], &lo, &hi));
9589566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_yi));
9599566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_y));
9609566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->y, is_from_y, user->yi[i], is_to_yi, &user->yi_scatter[i]));
961c4762a1bSJed Brown     istart = istart + hi - lo;
9629566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_yi));
9639566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_y));
964c4762a1bSJed Brown   }
9659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&yi));
966c4762a1bSJed Brown 
967c4762a1bSJed Brown   /* Create scatter from d to d_1,d_2,...,d_ns */
9689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->ns * user->ndata, &user->di_scatter));
9699566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &di));
9709566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(di, PETSC_DECIDE, user->ndata));
9719566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(di));
9729566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(di, user->ns, &user->di));
9739566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->d, &lo2, &hi2));
974c4762a1bSJed Brown   istart = 0;
975c4762a1bSJed Brown   for (i = 0; i < user->ns; i++) {
9769566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->di[i], &lo, &hi));
9779566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_di));
9789566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_d));
9799566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->d, is_from_d, user->di[i], is_to_di, &user->di_scatter[i]));
980c4762a1bSJed Brown     istart = istart + hi - lo;
9819566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_to_di));
9829566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_d));
983c4762a1bSJed Brown   }
9849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&di));
985c4762a1bSJed Brown 
986c4762a1bSJed Brown   /* Assemble RHS of forward problem */
9879566063dSJacob Faibussowitsch   PetscCall(VecCopy(bc, user->yiwork[0]));
988*48a46eb9SPierre Jolivet   for (i = 1; i < user->nt; i++) PetscCall(VecSet(user->yiwork[i], 0.0));
9899566063dSJacob Faibussowitsch   PetscCall(Gather_i(user->q, user->yiwork, user->yi_scatter, user->nt));
9909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bc));
991c4762a1bSJed Brown 
992c4762a1bSJed Brown   /* Compute true state function yt given ut */
9939566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ytrue));
9949566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->ytrue, PETSC_DECIDE, n * user->nt));
9959566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->ytrue));
996c4762a1bSJed Brown 
997c4762a1bSJed Brown   /* First compute Av_u = Av*exp(-u) */
9989566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork, 0));
9999566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); /*  Note: user->utrue */
10009566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
10019566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
1002c4762a1bSJed Brown 
1003c20d7725SJed Brown   /* Symbolic DSG = Div * Grad */
10049566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(user->Div, user->Grad, NULL, &user->DSG));
10059566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(user->DSG, MATPRODUCT_AB));
10069566063dSJacob Faibussowitsch   PetscCall(MatProductSetAlgorithm(user->DSG, "default"));
10079566063dSJacob Faibussowitsch   PetscCall(MatProductSetFill(user->DSG, PETSC_DEFAULT));
10089566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(user->DSG));
10099566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(user->DSG));
1010c20d7725SJed Brown 
1011c4762a1bSJed Brown   user->dsg_formed = PETSC_TRUE;
1012c4762a1bSJed Brown 
1013c20d7725SJed Brown   /* Next form DSG = Div*Grad */
10149566063dSJacob Faibussowitsch   PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
10159566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u));
1016c4762a1bSJed Brown   if (user->dsg_formed) {
10179566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(user->DSG));
1018c4762a1bSJed Brown   } else {
10199566063dSJacob Faibussowitsch     PetscCall(MatMatMult(user->Div, user->Grad, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DSG));
1020c4762a1bSJed Brown     user->dsg_formed = PETSC_TRUE;
1021c4762a1bSJed Brown   }
1022c4762a1bSJed Brown   /* B = speye(nx^3) + ht*DSG; */
10239566063dSJacob Faibussowitsch   PetscCall(MatScale(user->DSG, user->ht));
10249566063dSJacob Faibussowitsch   PetscCall(MatShift(user->DSG, 1.0));
1025c4762a1bSJed Brown 
1026c4762a1bSJed Brown   /* Now solve for ytrue */
10279566063dSJacob Faibussowitsch   PetscCall(StateMatInvMult(user->Js, user->q, user->ytrue));
1028c4762a1bSJed Brown 
1029c4762a1bSJed Brown   /* Initial guess y0 for state given u0 */
1030c4762a1bSJed Brown 
1031c4762a1bSJed Brown   /* First compute Av_u = Av*exp(-u) */
10329566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork, 0));
10339566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork, -1.0, user->u)); /*  Note: user->u */
10349566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
10359566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
1036c4762a1bSJed Brown 
1037c4762a1bSJed Brown   /* Next form DSG = Div*S*Grad */
10389566063dSJacob Faibussowitsch   PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
10399566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u));
1040c4762a1bSJed Brown   if (user->dsg_formed) {
10419566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(user->DSG));
1042c4762a1bSJed Brown   } else {
10439566063dSJacob Faibussowitsch     PetscCall(MatMatMult(user->Div, user->Grad, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DSG));
1044c20d7725SJed Brown 
1045c4762a1bSJed Brown     user->dsg_formed = PETSC_TRUE;
1046c4762a1bSJed Brown   }
1047c4762a1bSJed Brown   /* B = speye(nx^3) + ht*DSG; */
10489566063dSJacob Faibussowitsch   PetscCall(MatScale(user->DSG, user->ht));
10499566063dSJacob Faibussowitsch   PetscCall(MatShift(user->DSG, 1.0));
1050c4762a1bSJed Brown 
1051c4762a1bSJed Brown   /* Now solve for y */
10529566063dSJacob Faibussowitsch   PetscCall(StateMatInvMult(user->Js, user->q, user->y));
1053c4762a1bSJed Brown 
1054c4762a1bSJed Brown   /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
10559566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Qblock));
10569566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Qblock, PETSC_DECIDE, PETSC_DECIDE, user->ndata, n));
10579566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Qblock));
10589566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Qblock, 8, NULL, 8, NULL));
10599566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Qblock, 8, NULL));
1060c4762a1bSJed Brown 
1061c4762a1bSJed Brown   for (i = 0; i < user->mx; i++) {
1062c4762a1bSJed Brown     x[i] = h * (i + 0.5);
1063c4762a1bSJed Brown     y[i] = h * (i + 0.5);
1064c4762a1bSJed Brown     z[i] = h * (i + 0.5);
1065c4762a1bSJed Brown   }
1066c4762a1bSJed Brown 
10679566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Qblock, &istart, &iend));
10689371c9d4SSatish Balay   nx = user->mx;
10699371c9d4SSatish Balay   ny = user->mx;
10709371c9d4SSatish Balay   nz = user->mx;
1071c4762a1bSJed Brown   for (i = istart; i < iend; i++) {
1072c4762a1bSJed Brown     xri = xr[i];
1073c4762a1bSJed Brown     im  = 0;
1074c4762a1bSJed Brown     xim = x[im];
1075c4762a1bSJed Brown     while (xri > xim && im < nx) {
1076c4762a1bSJed Brown       im  = im + 1;
1077c4762a1bSJed Brown       xim = x[im];
1078c4762a1bSJed Brown     }
1079c4762a1bSJed Brown     indx1 = im - 1;
1080c4762a1bSJed Brown     indx2 = im;
1081c4762a1bSJed Brown     dx1   = xri - x[indx1];
1082c4762a1bSJed Brown     dx2   = x[indx2] - xri;
1083c4762a1bSJed Brown 
1084c4762a1bSJed Brown     yri = yr[i];
1085c4762a1bSJed Brown     im  = 0;
1086c4762a1bSJed Brown     yim = y[im];
1087c4762a1bSJed Brown     while (yri > yim && im < ny) {
1088c4762a1bSJed Brown       im  = im + 1;
1089c4762a1bSJed Brown       yim = y[im];
1090c4762a1bSJed Brown     }
1091c4762a1bSJed Brown     indy1 = im - 1;
1092c4762a1bSJed Brown     indy2 = im;
1093c4762a1bSJed Brown     dy1   = yri - y[indy1];
1094c4762a1bSJed Brown     dy2   = y[indy2] - yri;
1095c4762a1bSJed Brown 
1096c4762a1bSJed Brown     zri = zr[i];
1097c4762a1bSJed Brown     im  = 0;
1098c4762a1bSJed Brown     zim = z[im];
1099c4762a1bSJed Brown     while (zri > zim && im < nz) {
1100c4762a1bSJed Brown       im  = im + 1;
1101c4762a1bSJed Brown       zim = z[im];
1102c4762a1bSJed Brown     }
1103c4762a1bSJed Brown     indz1 = im - 1;
1104c4762a1bSJed Brown     indz2 = im;
1105c4762a1bSJed Brown     dz1   = zri - z[indz1];
1106c4762a1bSJed Brown     dz2   = z[indz2] - zri;
1107c4762a1bSJed Brown 
1108c4762a1bSJed Brown     Dx = x[indx2] - x[indx1];
1109c4762a1bSJed Brown     Dy = y[indy2] - y[indy1];
1110c4762a1bSJed Brown     Dz = z[indz2] - z[indz1];
1111c4762a1bSJed Brown 
1112c4762a1bSJed Brown     j = indx1 + indy1 * nx + indz1 * nx * ny;
1113c4762a1bSJed Brown     v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
11149566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1115c4762a1bSJed Brown 
1116c4762a1bSJed Brown     j = indx1 + indy1 * nx + indz2 * nx * ny;
1117c4762a1bSJed Brown     v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
11189566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1119c4762a1bSJed Brown 
1120c4762a1bSJed Brown     j = indx1 + indy2 * nx + indz1 * nx * ny;
1121c4762a1bSJed Brown     v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
11229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1123c4762a1bSJed Brown 
1124c4762a1bSJed Brown     j = indx1 + indy2 * nx + indz2 * nx * ny;
1125c4762a1bSJed Brown     v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
11269566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1127c4762a1bSJed Brown 
1128c4762a1bSJed Brown     j = indx2 + indy1 * nx + indz1 * nx * ny;
1129c4762a1bSJed Brown     v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
11309566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1131c4762a1bSJed Brown 
1132c4762a1bSJed Brown     j = indx2 + indy1 * nx + indz2 * nx * ny;
1133c4762a1bSJed Brown     v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
11349566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1135c4762a1bSJed Brown 
1136c4762a1bSJed Brown     j = indx2 + indy2 * nx + indz1 * nx * ny;
1137c4762a1bSJed Brown     v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
11389566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1139c4762a1bSJed Brown 
1140c4762a1bSJed Brown     j = indx2 + indy2 * nx + indz2 * nx * ny;
1141c4762a1bSJed Brown     v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
11429566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1143c4762a1bSJed Brown   }
11449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Qblock, MAT_FINAL_ASSEMBLY));
11459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Qblock, MAT_FINAL_ASSEMBLY));
1146c4762a1bSJed Brown 
11479566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->Qblock, MAT_INITIAL_MATRIX, &user->QblockT));
11489566063dSJacob Faibussowitsch   PetscCall(MatTranspose(user->L, MAT_INITIAL_MATRIX, &user->LT));
1149c4762a1bSJed Brown 
1150c4762a1bSJed Brown   /* Add noise to the measurement data */
11519566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ywork, 1.0));
11529566063dSJacob Faibussowitsch   PetscCall(VecAYPX(user->ywork, user->noise, user->ytrue));
11539566063dSJacob Faibussowitsch   PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
1154c4762a1bSJed Brown   for (j = 0; j < user->ns; j++) {
1155c4762a1bSJed Brown     i = user->sample_times[j];
11569566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Qblock, user->yiwork[i], user->di[j]));
1157c4762a1bSJed Brown   }
11589566063dSJacob Faibussowitsch   PetscCall(Gather_i(user->d, user->di, user->di_scatter, user->ns));
1159c4762a1bSJed Brown 
1160c4762a1bSJed Brown   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11619566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(user->solver));
11629566063dSJacob Faibussowitsch   PetscCall(PetscFree(x));
11639566063dSJacob Faibussowitsch   PetscCall(PetscFree(y));
11649566063dSJacob Faibussowitsch   PetscCall(PetscFree(z));
1165c4762a1bSJed Brown   PetscFunctionReturn(0);
1166c4762a1bSJed Brown }
1167c4762a1bSJed Brown 
11689371c9d4SSatish Balay PetscErrorCode ParabolicDestroy(AppCtx *user) {
1169c4762a1bSJed Brown   PetscInt i;
1170c4762a1bSJed Brown 
1171c4762a1bSJed Brown   PetscFunctionBegin;
11729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Qblock));
11739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->QblockT));
11749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Div));
11759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Divwork));
11769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Grad));
11779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Av));
11789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Avwork));
11799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->AvT));
11809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->DSG));
11819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->L));
11829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->LT));
11839566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&user->solver));
11849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Js));
11859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Jd));
11869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsInv));
11879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsBlock));
11889566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsBlockPrec));
11899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->u));
11909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->uwork));
11919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->utrue));
11929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->y));
11939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ywork));
11949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ytrue));
11959566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->yi));
11969566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->nt, &user->yiwork));
11979566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(user->ns, &user->di));
11989566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->yi));
11999566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->yiwork));
12009566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->di));
12019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->c));
12029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->cwork));
12039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ur));
12049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->q));
12059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->d));
12069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->dwork));
12079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->lwork));
12089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->S));
12099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Swork));
12109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Av_u));
12119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Twork));
12129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Rwork));
12139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->js_diag));
12149566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->s_is));
12159566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->d_is));
12169566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->state_scatter));
12179566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->design_scatter));
1218*48a46eb9SPierre Jolivet   for (i = 0; i < user->nt; i++) PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
1219*48a46eb9SPierre Jolivet   for (i = 0; i < user->ns; i++) PetscCall(VecScatterDestroy(&user->di_scatter[i]));
12209566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->yi_scatter));
12219566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->di_scatter));
12229566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->sample_times));
1223c4762a1bSJed Brown   PetscFunctionReturn(0);
1224c4762a1bSJed Brown }
1225c4762a1bSJed Brown 
12269371c9d4SSatish Balay PetscErrorCode ParabolicMonitor(Tao tao, void *ptr) {
1227c4762a1bSJed Brown   Vec       X;
1228c4762a1bSJed Brown   PetscReal unorm, ynorm;
1229c4762a1bSJed Brown   AppCtx   *user = (AppCtx *)ptr;
1230c4762a1bSJed Brown 
1231c4762a1bSJed Brown   PetscFunctionBegin;
12329566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao, &X));
12339566063dSJacob Faibussowitsch   PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
12349566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue));
12359566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork, -1.0, user->utrue));
12369566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->uwork, NORM_2, &unorm));
12379566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->ywork, NORM_2, &ynorm));
12389566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm));
1239c4762a1bSJed Brown   PetscFunctionReturn(0);
1240c4762a1bSJed Brown }
1241c4762a1bSJed Brown 
1242c4762a1bSJed Brown /*TEST
1243c4762a1bSJed Brown 
1244c4762a1bSJed Brown    build:
1245c4762a1bSJed Brown       requires: !complex
1246c4762a1bSJed Brown 
1247c4762a1bSJed Brown    test:
1248c4762a1bSJed Brown       args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30
1249c4762a1bSJed Brown       requires: !single
1250c4762a1bSJed Brown 
1251c4762a1bSJed Brown TEST*/
1252