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
main(int argc,char ** argv)96d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
97d71ae5a4SJacob Faibussowitsch {
98c4762a1bSJed Brown Vec x, x0;
99c4762a1bSJed Brown Tao tao;
100c4762a1bSJed Brown AppCtx user;
101c4762a1bSJed Brown IS is_allstate, is_alldesign;
102c4762a1bSJed Brown PetscInt lo, hi, hi2, lo2, ksp_old;
103c4762a1bSJed Brown PetscInt ntests = 1;
104c4762a1bSJed Brown PetscInt i;
105c4762a1bSJed Brown PetscLogStage stages[1];
106c4762a1bSJed Brown
107327415f7SBarry Smith PetscFunctionBeginUser;
108c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
109c4762a1bSJed Brown user.mx = 8;
110d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "parabolic example", NULL);
1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL));
112c4762a1bSJed Brown user.nt = 8;
1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-nt", "Number of time steps", "", user.nt, &user.nt, NULL));
114c4762a1bSJed Brown user.ndata = 64;
1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL));
116c4762a1bSJed Brown user.ns = 8;
1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ns", "Number of samples", "", user.ns, &user.ns, NULL));
118c4762a1bSJed Brown user.alpha = 1.0;
1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL));
120c4762a1bSJed Brown user.beta = 0.01;
1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-beta", "Weight attributed to ||u||^2 in regularization functional", "", user.beta, &user.beta, NULL));
122c4762a1bSJed Brown user.noise = 0.01;
1239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-noise", "Amount of noise to add to data", "", user.noise, &user.noise, NULL));
1249566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL));
125d0609cedSBarry Smith PetscOptionsEnd();
126c4762a1bSJed Brown
127c4762a1bSJed Brown user.m = user.mx * user.mx * user.mx; /* number of constraints per time step */
128c4762a1bSJed Brown user.n = user.m * (user.nt + 1); /* number of variables */
129c4762a1bSJed Brown user.ht = (PetscReal)1 / user.nt;
130c4762a1bSJed Brown
1319566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.u));
1329566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.y));
1339566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.c));
1349566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.u, PETSC_DECIDE, user.n - user.m * user.nt));
1359566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.y, PETSC_DECIDE, user.m * user.nt));
1369566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.c, PETSC_DECIDE, user.m * user.nt));
1379566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.u));
1389566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.y));
1399566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.c));
140c4762a1bSJed Brown
141c4762a1bSJed Brown /* Create scatters for reduced spaces.
142c4762a1bSJed Brown If the state vector y and design vector u are partitioned as
143c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
144c4762a1bSJed Brown then the solution vector x is organized as
145c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np].
146c4762a1bSJed Brown The index sets user.s_is and user.d_is correspond to the indices of the
147c4762a1bSJed Brown state and design variables owned by the current processor.
148c4762a1bSJed Brown */
1499566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
150c4762a1bSJed Brown
1519566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.y, &lo, &hi));
1529566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.u, &lo2, &hi2));
153c4762a1bSJed Brown
1549566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate));
1559566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user.s_is));
156c4762a1bSJed Brown
1579566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign));
1589566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user.d_is));
159c4762a1bSJed Brown
1609566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, hi - lo + hi2 - lo2, user.n));
1619566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x));
162c4762a1bSJed Brown
1639566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.s_is, user.y, is_allstate, &user.state_scatter));
1649566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.d_is, user.u, is_alldesign, &user.design_scatter));
1659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign));
1669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate));
167c4762a1bSJed Brown
168c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
1699566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1709566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLCL));
171c4762a1bSJed Brown
172c4762a1bSJed Brown /* Set up initial vectors and matrices */
1739566063dSJacob Faibussowitsch PetscCall(ParabolicInitialize(&user));
174c4762a1bSJed Brown
1759566063dSJacob Faibussowitsch PetscCall(Gather(x, user.y, user.state_scatter, user.u, user.design_scatter));
1769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &x0));
1779566063dSJacob Faibussowitsch PetscCall(VecCopy(x, x0));
178c4762a1bSJed Brown
179c4762a1bSJed Brown /* Set solution vector with an initial guess */
1809566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
1819566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user));
1829566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
1839566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
184c4762a1bSJed Brown
1859566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, &user));
1869566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
187c4762a1bSJed Brown
1889566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
1899566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is));
190c4762a1bSJed Brown
191c4762a1bSJed Brown /* SOLVE THE APPLICATION */
1929566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials", &stages[0]));
1939566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0]));
194c4762a1bSJed Brown user.ksp_its_initial = user.ksp_its;
195c4762a1bSJed Brown ksp_old = user.ksp_its;
196c4762a1bSJed Brown for (i = 0; i < ntests; i++) {
1979566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
19863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its - ksp_old));
1999566063dSJacob Faibussowitsch PetscCall(VecCopy(x0, x));
2009566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
201c4762a1bSJed Brown }
2029566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop());
2039566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)x));
2049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: "));
20563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial));
20663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total KSP iterations over %" PetscInt_FMT " trial(s): ", ntests));
20763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its));
2089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations per trial: "));
20963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", (user.ksp_its - user.ksp_its_initial) / ntests));
210c4762a1bSJed Brown
2119566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
2129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
2139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0));
2149566063dSJacob Faibussowitsch PetscCall(ParabolicDestroy(&user));
215c4762a1bSJed Brown
2169566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
217b122ec5aSJacob Faibussowitsch return 0;
218c4762a1bSJed Brown }
219c4762a1bSJed Brown /* ------------------------------------------------------------------- */
220c4762a1bSJed Brown /*
221c4762a1bSJed Brown dwork = Qy - d
222c4762a1bSJed Brown lwork = L*(u-ur)
223c4762a1bSJed Brown f = 1/2 * (dwork.dork + alpha*lwork.lwork)
224c4762a1bSJed Brown */
FormFunction(Tao tao,Vec X,PetscReal * f,void * ptr)225d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr)
226d71ae5a4SJacob Faibussowitsch {
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);
2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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 */
FormGradient(Tao tao,Vec X,Vec G,void * ptr)255d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr)
256d71ae5a4SJacob Faibussowitsch {
257c4762a1bSJed Brown PetscInt i, j;
258c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
259c4762a1bSJed Brown
260c4762a1bSJed Brown PetscFunctionBegin;
2619566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2629566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
263c4762a1bSJed Brown for (j = 0; j < user->ns; j++) {
264c4762a1bSJed Brown i = user->sample_times[j];
2659566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j]));
266c4762a1bSJed Brown }
2679566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns));
2689566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2699566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->dwork, user->di, user->di_scatter, user->ns));
2709566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork, 0.0));
2719566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
272c4762a1bSJed Brown for (j = 0; j < user->ns; j++) {
273c4762a1bSJed Brown i = user->sample_times[j];
2749566063dSJacob Faibussowitsch PetscCall(MatMult(user->QblockT, user->di[j], user->yiwork[i]));
275c4762a1bSJed Brown }
2769566063dSJacob Faibussowitsch PetscCall(Gather_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
277c4762a1bSJed Brown
2789566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
2799566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork));
2809566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->lwork, user->uwork));
2819566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha));
2829566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
284c4762a1bSJed Brown }
285c4762a1bSJed Brown
FormFunctionGradient(Tao tao,Vec X,PetscReal * f,Vec G,void * ptr)286d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
287d71ae5a4SJacob Faibussowitsch {
288c4762a1bSJed Brown PetscReal d1, d2;
289c4762a1bSJed Brown PetscInt i, j;
290c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
291c4762a1bSJed Brown
292c4762a1bSJed Brown PetscFunctionBegin;
2939566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2949566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
295c4762a1bSJed Brown for (j = 0; j < user->ns; j++) {
296c4762a1bSJed Brown i = user->sample_times[j];
2979566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j]));
298c4762a1bSJed Brown }
2999566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns));
3009566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
3019566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1));
3029566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->dwork, user->di, user->di_scatter, user->ns));
3039566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork, 0.0));
3049566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
305c4762a1bSJed Brown for (j = 0; j < user->ns; j++) {
306c4762a1bSJed Brown i = user->sample_times[j];
3079566063dSJacob Faibussowitsch PetscCall(MatMult(user->QblockT, user->di[j], user->yiwork[i]));
308c4762a1bSJed Brown }
3099566063dSJacob Faibussowitsch PetscCall(Gather_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
310c4762a1bSJed Brown
3119566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
3129566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork));
3139566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork, user->lwork, &d2));
3149566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->lwork, user->uwork));
3159566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha));
316c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2);
317c4762a1bSJed Brown
3189566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
320c4762a1bSJed Brown }
321c4762a1bSJed Brown
322c4762a1bSJed Brown /* ------------------------------------------------------------------- */
323c4762a1bSJed Brown /* A
324c4762a1bSJed Brown MatShell object
325c4762a1bSJed Brown */
FormJacobianState(Tao tao,Vec X,Mat J,Mat JPre,Mat JInv,void * ptr)326d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
327d71ae5a4SJacob Faibussowitsch {
328c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
329c4762a1bSJed Brown
330c4762a1bSJed Brown PetscFunctionBegin;
3319566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
3329566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
3339566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u));
3349566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
3359566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
3369566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u, user->Swork));
3379566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork));
3389566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
3399566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork));
340c4762a1bSJed Brown if (user->dsg_formed) {
3419566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG));
342c4762a1bSJed Brown } else {
343fb842aefSJose E. Roman PetscCall(MatMatMult(user->Divwork, user->Grad, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &user->DSG));
344c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE;
345c4762a1bSJed Brown }
346c4762a1bSJed Brown
347c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */
3489566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG, user->ht));
3499566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG, 1.0));
3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown
353c4762a1bSJed Brown /* ------------------------------------------------------------------- */
354c4762a1bSJed Brown /* B */
FormJacobianDesign(Tao tao,Vec X,Mat J,void * ptr)355d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
356d71ae5a4SJacob Faibussowitsch {
357c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
358c4762a1bSJed Brown
359c4762a1bSJed Brown PetscFunctionBegin;
3609566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
362c4762a1bSJed Brown }
363c4762a1bSJed Brown
StateMatMult(Mat J_shell,Vec X,Vec Y)364d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
365d71ae5a4SJacob Faibussowitsch {
366c4762a1bSJed Brown PetscInt i;
367c4762a1bSJed Brown AppCtx *user;
368c4762a1bSJed Brown
369c4762a1bSJed Brown PetscFunctionBegin;
3709566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
3719566063dSJacob Faibussowitsch PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
3729566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
373c4762a1bSJed Brown for (i = 1; i < user->nt; i++) {
3749566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3759566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i - 1]));
376c4762a1bSJed Brown }
3779566063dSJacob Faibussowitsch PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
3783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
379c4762a1bSJed Brown }
380c4762a1bSJed Brown
StateMatMultTranspose(Mat J_shell,Vec X,Vec Y)381d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
382d71ae5a4SJacob Faibussowitsch {
383c4762a1bSJed Brown PetscInt i;
384c4762a1bSJed Brown AppCtx *user;
385c4762a1bSJed Brown
386c4762a1bSJed Brown PetscFunctionBegin;
3879566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
3889566063dSJacob Faibussowitsch PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
389c4762a1bSJed Brown for (i = 0; i < user->nt - 1; i++) {
3909566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3919566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i + 1]));
392c4762a1bSJed Brown }
393c4762a1bSJed Brown i = user->nt - 1;
3949566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3959566063dSJacob Faibussowitsch PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
3963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
397c4762a1bSJed Brown }
398c4762a1bSJed Brown
StateMatBlockMult(Mat J_shell,Vec X,Vec Y)399d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
400d71ae5a4SJacob Faibussowitsch {
401c4762a1bSJed Brown AppCtx *user;
402c4762a1bSJed Brown
403c4762a1bSJed Brown PetscFunctionBegin;
4049566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
4059566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, X, user->Swork));
4069566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u));
4079566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->Swork, Y));
4089566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, user->ht, X));
4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
410c4762a1bSJed Brown }
411c4762a1bSJed Brown
DesignMatMult(Mat J_shell,Vec X,Vec Y)412d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
413d71ae5a4SJacob Faibussowitsch {
414c4762a1bSJed Brown PetscInt i;
415c4762a1bSJed Brown AppCtx *user;
416c4762a1bSJed Brown
417c4762a1bSJed Brown PetscFunctionBegin;
4189566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
419c4762a1bSJed Brown
420c4762a1bSJed Brown /* sdiag(1./v) */
4219566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
4229566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u));
4239566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
424c4762a1bSJed Brown
425c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */
4269566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Swork));
4279566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Swork, user->Swork));
4289566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork));
429c4762a1bSJed Brown
430c4762a1bSJed Brown /* (Av * (sdiag(1./v) * b)) */
4319566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->uwork, X));
4329566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Twork));
433c4762a1bSJed Brown
434c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
4359566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork));
436c4762a1bSJed Brown
4379566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
438c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
439c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) */
4409566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->yi[i], user->Twork));
441c4762a1bSJed Brown
442c4762a1bSJed Brown /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
4439566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork, user->Twork, user->Swork));
4449566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->Twork, user->yiwork[i]));
4459566063dSJacob Faibussowitsch PetscCall(VecScale(user->yiwork[i], user->ht));
446c4762a1bSJed Brown }
4479566063dSJacob Faibussowitsch PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
449c4762a1bSJed Brown }
450c4762a1bSJed Brown
DesignMatMultTranspose(Mat J_shell,Vec X,Vec Y)451d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
452d71ae5a4SJacob Faibussowitsch {
453c4762a1bSJed Brown PetscInt i;
454c4762a1bSJed Brown AppCtx *user;
455c4762a1bSJed Brown
456c4762a1bSJed Brown PetscFunctionBegin;
4579566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
458c4762a1bSJed Brown
459c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */
4609566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
4619566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u));
4629566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
4639566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Rwork));
4649566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Rwork, user->Rwork, user->Rwork));
4659566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Rwork));
466c4762a1bSJed Brown
4679566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 0.0));
4689566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
4699566063dSJacob Faibussowitsch PetscCall(Scatter_i(X, user->yiwork, user->yi_scatter, user->nt));
470c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
471c4762a1bSJed Brown /* (Div' * b(:,i)) */
4729566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->yiwork[i], user->Swork));
473c4762a1bSJed Brown
474c4762a1bSJed Brown /* sdiag(Grad*y(:,i)) */
4759566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->yi[i], user->Twork));
476c4762a1bSJed Brown
477c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */
4789566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork, user->Swork, user->Twork));
479c4762a1bSJed Brown
480c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */
4819566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork, user->Rwork, user->Twork));
482c4762a1bSJed Brown
483c4762a1bSJed Brown /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
4849566063dSJacob Faibussowitsch PetscCall(MatMult(user->AvT, user->Twork, user->yiwork[i]));
485c4762a1bSJed Brown
486c4762a1bSJed Brown /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
4879566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->yiwork[i], user->uwork, user->yiwork[i]));
4889566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, user->ht, user->yiwork[i]));
489c4762a1bSJed Brown }
4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
491c4762a1bSJed Brown }
492c4762a1bSJed Brown
StateMatBlockPrecMult(PC PC_shell,Vec X,Vec Y)493d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
494d71ae5a4SJacob Faibussowitsch {
495c4762a1bSJed Brown AppCtx *user;
496c4762a1bSJed Brown
497c4762a1bSJed Brown PetscFunctionBegin;
4989566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell, &user));
499c4762a1bSJed Brown
500966bd95aSPierre Jolivet PetscCheck(user->dsg_formed, PETSC_COMM_SELF, PETSC_ERR_SUP, "DSG not formed");
5019566063dSJacob Faibussowitsch PetscCall(MatSOR(user->DSG, X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y));
5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
503c4762a1bSJed Brown }
504c4762a1bSJed Brown
StateMatInvMult(Mat J_shell,Vec X,Vec Y)505d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
506d71ae5a4SJacob Faibussowitsch {
507c4762a1bSJed Brown AppCtx *user;
508c4762a1bSJed Brown PetscInt its, i;
509c4762a1bSJed Brown
510c4762a1bSJed Brown PetscFunctionBegin;
5119566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
512c4762a1bSJed Brown
513c4762a1bSJed Brown if (Y == user->ytrue) {
514c4762a1bSJed Brown /* First solve is done with true solution to set up problem */
515606f75f6SBarry Smith PetscCall(KSPSetTolerances(user->solver, 1e-8, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
516c4762a1bSJed Brown } else {
517606f75f6SBarry Smith PetscCall(KSPSetTolerances(user->solver, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
518c4762a1bSJed Brown }
519c4762a1bSJed Brown
5209566063dSJacob Faibussowitsch PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
5219566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[0], user->yiwork[0]));
5229566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
523c4762a1bSJed Brown user->ksp_its = user->ksp_its + its;
524c4762a1bSJed Brown
525c4762a1bSJed Brown for (i = 1; i < user->nt; i++) {
5269566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->yiwork[i - 1]));
5279566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
5289566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
529c4762a1bSJed Brown user->ksp_its = user->ksp_its + its;
530c4762a1bSJed Brown }
5319566063dSJacob Faibussowitsch PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
533c4762a1bSJed Brown }
534c4762a1bSJed Brown
StateMatInvTransposeMult(Mat J_shell,Vec X,Vec Y)535d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
536d71ae5a4SJacob Faibussowitsch {
537c4762a1bSJed Brown AppCtx *user;
538c4762a1bSJed Brown PetscInt its, i;
539c4762a1bSJed Brown
540c4762a1bSJed Brown PetscFunctionBegin;
5419566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
542c4762a1bSJed Brown
5439566063dSJacob Faibussowitsch PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt));
544c4762a1bSJed Brown
545c4762a1bSJed Brown i = user->nt - 1;
5469566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
547c4762a1bSJed Brown
5489566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
549c4762a1bSJed Brown user->ksp_its = user->ksp_its + its;
550c4762a1bSJed Brown
551c4762a1bSJed Brown for (i = user->nt - 2; i >= 0; i--) {
5529566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->yiwork[i + 1]));
5539566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
554c4762a1bSJed Brown
5559566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
556c4762a1bSJed Brown user->ksp_its = user->ksp_its + its;
557c4762a1bSJed Brown }
558c4762a1bSJed Brown
5599566063dSJacob Faibussowitsch PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt));
5603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
561c4762a1bSJed Brown }
562c4762a1bSJed Brown
StateMatDuplicate(Mat J_shell,MatDuplicateOption opt,Mat * new_shell)563d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
564d71ae5a4SJacob Faibussowitsch {
565c4762a1bSJed Brown AppCtx *user;
566c4762a1bSJed Brown
567c4762a1bSJed Brown PetscFunctionBegin;
5689566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
569c4762a1bSJed Brown
5709566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, new_shell));
571*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT, (PetscErrorCodeFn *)StateMatMult));
572*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_DUPLICATE, (PetscErrorCodeFn *)StateMatDuplicate));
573*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatMultTranspose));
574*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)StateMatGetDiagonal));
5753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
576c4762a1bSJed Brown }
577c4762a1bSJed Brown
StateMatGetDiagonal(Mat J_shell,Vec X)578d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
579d71ae5a4SJacob Faibussowitsch {
580c4762a1bSJed Brown AppCtx *user;
581c4762a1bSJed Brown
582c4762a1bSJed Brown PetscFunctionBegin;
5839566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
5849566063dSJacob Faibussowitsch PetscCall(VecCopy(user->js_diag, X));
5853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
586c4762a1bSJed Brown }
587c4762a1bSJed Brown
FormConstraints(Tao tao,Vec X,Vec C,void * ptr)588d71ae5a4SJacob Faibussowitsch PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
589d71ae5a4SJacob Faibussowitsch {
590c4762a1bSJed Brown /* con = Ay - q, A = [B 0 0 ... 0;
591c4762a1bSJed Brown -I B 0 ... 0;
592c4762a1bSJed Brown 0 -I B ... 0;
593c4762a1bSJed Brown ... ;
594c4762a1bSJed Brown 0 ... -I B]
595c4762a1bSJed Brown B = ht * Div * Sigma * Grad + eye */
596c4762a1bSJed Brown PetscInt i;
597c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
598c4762a1bSJed Brown
599c4762a1bSJed Brown PetscFunctionBegin;
6009566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
6019566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt));
6029566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
603c4762a1bSJed Brown for (i = 1; i < user->nt; i++) {
6049566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
6059566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i - 1]));
606c4762a1bSJed Brown }
6079566063dSJacob Faibussowitsch PetscCall(Gather_i(C, user->yiwork, user->yi_scatter, user->nt));
6089566063dSJacob Faibussowitsch PetscCall(VecAXPY(C, -1.0, user->q));
6093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
610c4762a1bSJed Brown }
611c4762a1bSJed Brown
Scatter(Vec x,Vec state,VecScatter s_scat,Vec design,VecScatter d_scat)612d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
613d71ae5a4SJacob Faibussowitsch {
614c4762a1bSJed Brown PetscFunctionBegin;
6159566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6169566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6179566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
6189566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
6193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
620c4762a1bSJed Brown }
621c4762a1bSJed Brown
Scatter_i(Vec y,Vec * yi,VecScatter * scat,PetscInt nt)622d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
623d71ae5a4SJacob Faibussowitsch {
624c4762a1bSJed Brown PetscInt i;
625c4762a1bSJed Brown
626c4762a1bSJed Brown PetscFunctionBegin;
627c4762a1bSJed Brown for (i = 0; i < nt; i++) {
6289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
6299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
630c4762a1bSJed Brown }
6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
632c4762a1bSJed Brown }
633c4762a1bSJed Brown
Gather(Vec x,Vec state,VecScatter s_scat,Vec design,VecScatter d_scat)634d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
635d71ae5a4SJacob Faibussowitsch {
636c4762a1bSJed Brown PetscFunctionBegin;
6379566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6389566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6399566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
6409566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
6413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
642c4762a1bSJed Brown }
643c4762a1bSJed Brown
Gather_i(Vec y,Vec * yi,VecScatter * scat,PetscInt nt)644d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
645d71ae5a4SJacob Faibussowitsch {
646c4762a1bSJed Brown PetscInt i;
647c4762a1bSJed Brown
648c4762a1bSJed Brown PetscFunctionBegin;
649c4762a1bSJed Brown for (i = 0; i < nt; i++) {
6509566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
6519566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
652c4762a1bSJed Brown }
6533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
654c4762a1bSJed Brown }
655c4762a1bSJed Brown
ParabolicInitialize(AppCtx * user)656d71ae5a4SJacob Faibussowitsch PetscErrorCode ParabolicInitialize(AppCtx *user)
657d71ae5a4SJacob Faibussowitsch {
658c4762a1bSJed Brown PetscInt m, n, i, j, k, linear_index, istart, iend, iblock, lo, hi, lo2, hi2;
659c4762a1bSJed Brown Vec XX, YY, ZZ, XXwork, YYwork, ZZwork, UTwork, yi, di, bc;
660c4762a1bSJed Brown PetscReal *x, *y, *z;
661c4762a1bSJed Brown PetscReal h, stime;
662c4762a1bSJed Brown PetscScalar hinv, neg_hinv, half = 0.5, sqrt_beta;
663c4762a1bSJed Brown PetscInt im, indx1, indx2, indy1, indy2, indz1, indz2, nx, ny, nz;
664c4762a1bSJed Brown PetscReal xri, yri, zri, xim, yim, zim, dx1, dx2, dy1, dy2, dz1, dz2, Dx, Dy, Dz;
665c4762a1bSJed Brown PetscScalar v, vx, vy, vz;
666c4762a1bSJed Brown IS is_from_y, is_to_yi, is_from_d, is_to_di;
667c4762a1bSJed Brown /* Data locations */
6689371c9d4SSatish 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,
6699371c9d4SSatish 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,
6709371c9d4SSatish 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};
671c4762a1bSJed Brown
6729371c9d4SSatish 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,
6739371c9d4SSatish 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,
6749371c9d4SSatish 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};
675c4762a1bSJed Brown
6769371c9d4SSatish 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,
6779371c9d4SSatish 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,
6789371c9d4SSatish 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};
679c4762a1bSJed Brown
680c4762a1bSJed Brown PetscFunctionBegin;
6819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &x));
6829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &y));
6839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &z));
684c4762a1bSJed Brown user->jformed = PETSC_FALSE;
685c4762a1bSJed Brown user->dsg_formed = PETSC_FALSE;
686c4762a1bSJed Brown
687c4762a1bSJed Brown n = user->mx * user->mx * user->mx;
688c4762a1bSJed Brown m = 3 * user->mx * user->mx * (user->mx - 1);
689c4762a1bSJed Brown sqrt_beta = PetscSqrtScalar(user->beta);
690c4762a1bSJed Brown
691c4762a1bSJed Brown user->ksp_its = 0;
692c4762a1bSJed Brown user->ksp_its_initial = 0;
693c4762a1bSJed Brown
694c4762a1bSJed Brown stime = (PetscReal)user->nt / user->ns;
6959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns, &user->sample_times));
696ad540459SPierre Jolivet for (i = 0; i < user->ns; i++) user->sample_times[i] = (PetscInt)(stime * ((PetscReal)i + 1.0) - 0.5);
697c4762a1bSJed Brown
6989566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &XX));
6999566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q));
7009566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX, PETSC_DECIDE, n));
7019566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q, PETSC_DECIDE, n * user->nt));
7029566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX));
7039566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q));
704c4762a1bSJed Brown
7059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YY));
7069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &ZZ));
7079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &XXwork));
7089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YYwork));
7099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &ZZwork));
7109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &UTwork));
7119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->utrue));
7129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &bc));
713c4762a1bSJed Brown
714c4762a1bSJed Brown /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
715c4762a1bSJed Brown h = 1.0 / user->mx;
716c4762a1bSJed Brown hinv = user->mx;
717c4762a1bSJed Brown neg_hinv = -hinv;
718c4762a1bSJed Brown
7199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX, &istart, &iend));
720c4762a1bSJed Brown for (linear_index = istart; linear_index < iend; linear_index++) {
721c4762a1bSJed Brown i = linear_index % user->mx;
722c4762a1bSJed Brown j = ((linear_index - i) / user->mx) % user->mx;
723c4762a1bSJed Brown k = ((linear_index - i) / user->mx - j) / user->mx;
724c4762a1bSJed Brown vx = h * (i + 0.5);
725c4762a1bSJed Brown vy = h * (j + 0.5);
726c4762a1bSJed Brown vz = h * (k + 0.5);
7279566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES));
7289566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES));
7299566063dSJacob Faibussowitsch PetscCall(VecSetValues(ZZ, 1, &linear_index, &vz, INSERT_VALUES));
730c4762a1bSJed Brown if ((vx < 0.6) && (vx > 0.4) && (vy < 0.6) && (vy > 0.4) && (vy < 0.6) && (vz < 0.6) && (vz > 0.4)) {
731c4762a1bSJed Brown v = 1000.0;
7329566063dSJacob Faibussowitsch PetscCall(VecSetValues(bc, 1, &linear_index, &v, INSERT_VALUES));
733c4762a1bSJed Brown }
734c4762a1bSJed Brown }
735c4762a1bSJed Brown
7369566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX));
7379566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX));
7389566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY));
7399566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY));
7409566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(ZZ));
7419566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(ZZ));
7429566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bc));
7439566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bc));
744c4762a1bSJed Brown
745c4762a1bSJed Brown /* Compute true parameter function
746c4762a1bSJed Brown ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
7479566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork));
7489566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork));
7499566063dSJacob Faibussowitsch PetscCall(VecCopy(ZZ, ZZwork));
750c4762a1bSJed Brown
7519566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.5));
7529566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.5));
7539566063dSJacob Faibussowitsch PetscCall(VecShift(ZZwork, -0.5));
754c4762a1bSJed Brown
7559566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
7569566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
7579566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork));
758c4762a1bSJed Brown
7599566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, user->utrue));
7609566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue, 1.0, YYwork));
7619566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue, 1.0, ZZwork));
7629566063dSJacob Faibussowitsch PetscCall(VecScale(user->utrue, -10.0));
7639566063dSJacob Faibussowitsch PetscCall(VecExp(user->utrue));
7649566063dSJacob Faibussowitsch PetscCall(VecShift(user->utrue, 0.5));
765c4762a1bSJed Brown
7669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX));
7679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY));
7689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZ));
7699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork));
7709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork));
7719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZwork));
7729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&UTwork));
773c4762a1bSJed Brown
774c4762a1bSJed Brown /* Initial guess and reference model */
7759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue, &user->ur));
7769566063dSJacob Faibussowitsch PetscCall(VecSet(user->ur, 0.0));
777c4762a1bSJed Brown
778c4762a1bSJed Brown /* Generate Grad matrix */
7799566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad));
7809566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, PETSC_DECIDE, m, n));
7819566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad));
7829566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad, 2, NULL, 2, NULL));
7839566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad, 2, NULL));
7849566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));
785c4762a1bSJed Brown
786c4762a1bSJed Brown for (i = istart; i < iend; i++) {
787c4762a1bSJed Brown if (i < m / 3) {
788c4762a1bSJed Brown iblock = i / (user->mx - 1);
789c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1));
7909566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
791c4762a1bSJed Brown j = j + 1;
7929566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
793c4762a1bSJed Brown }
794c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) {
795c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1));
796c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
7979566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
798c4762a1bSJed Brown j = j + user->mx;
7999566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
800c4762a1bSJed Brown }
801c4762a1bSJed Brown if (i >= 2 * m / 3) {
802c4762a1bSJed Brown j = i - 2 * m / 3;
8039566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
804c4762a1bSJed Brown j = j + user->mx * user->mx;
8059566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
806c4762a1bSJed Brown }
807c4762a1bSJed Brown }
808c4762a1bSJed Brown
8099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY));
8109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY));
811c4762a1bSJed Brown
812c4762a1bSJed Brown /* Generate arithmetic averaging matrix Av */
8139566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Av));
8149566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Av, PETSC_DECIDE, PETSC_DECIDE, m, n));
8159566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Av));
8169566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Av, 2, NULL, 2, NULL));
8179566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Av, 2, NULL));
8189566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Av, &istart, &iend));
819c4762a1bSJed Brown
820c4762a1bSJed Brown for (i = istart; i < iend; i++) {
821c4762a1bSJed Brown if (i < m / 3) {
822c4762a1bSJed Brown iblock = i / (user->mx - 1);
823c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1));
8249566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
825c4762a1bSJed Brown j = j + 1;
8269566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
827c4762a1bSJed Brown }
828c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) {
829c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1));
830c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
8319566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
832c4762a1bSJed Brown j = j + user->mx;
8339566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
834c4762a1bSJed Brown }
835c4762a1bSJed Brown if (i >= 2 * m / 3) {
836c4762a1bSJed Brown j = i - 2 * m / 3;
8379566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
838c4762a1bSJed Brown j = j + user->mx * user->mx;
8399566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
840c4762a1bSJed Brown }
841c4762a1bSJed Brown }
842c4762a1bSJed Brown
8439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Av, MAT_FINAL_ASSEMBLY));
8449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Av, MAT_FINAL_ASSEMBLY));
845c4762a1bSJed Brown
846c4762a1bSJed Brown /* Generate transpose of averaging matrix Av */
8479566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Av, MAT_INITIAL_MATRIX, &user->AvT));
848c4762a1bSJed Brown
8499566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->L));
8509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->L, PETSC_DECIDE, PETSC_DECIDE, m + n, n));
8519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->L));
8529566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->L, 2, NULL, 2, NULL));
8539566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->L, 2, NULL));
8549566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->L, &istart, &iend));
855c4762a1bSJed Brown
856c4762a1bSJed Brown for (i = istart; i < iend; i++) {
857c4762a1bSJed Brown if (i < m / 3) {
858c4762a1bSJed Brown iblock = i / (user->mx - 1);
859c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1));
8609566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
861c4762a1bSJed Brown j = j + 1;
8629566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
863c4762a1bSJed Brown }
864c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) {
865c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1));
866c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
8679566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
868c4762a1bSJed Brown j = j + user->mx;
8699566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
870c4762a1bSJed Brown }
871c4762a1bSJed Brown if (i >= 2 * m / 3 && i < m) {
872c4762a1bSJed Brown j = i - 2 * m / 3;
8739566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
874c4762a1bSJed Brown j = j + user->mx * user->mx;
8759566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
876c4762a1bSJed Brown }
877c4762a1bSJed Brown if (i >= m) {
878c4762a1bSJed Brown j = i - m;
8799566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &sqrt_beta, INSERT_VALUES));
880c4762a1bSJed Brown }
881c4762a1bSJed Brown }
882c4762a1bSJed Brown
8839566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->L, MAT_FINAL_ASSEMBLY));
8849566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->L, MAT_FINAL_ASSEMBLY));
885c4762a1bSJed Brown
8869566063dSJacob Faibussowitsch PetscCall(MatScale(user->L, PetscPowScalar(h, 1.5)));
887c4762a1bSJed Brown
888c4762a1bSJed Brown /* Generate Div matrix */
8899566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Grad, MAT_INITIAL_MATRIX, &user->Div));
890c4762a1bSJed Brown
891c4762a1bSJed Brown /* Build work vectors and matrices */
8929566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->S));
8939566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->S, PETSC_DECIDE, user->mx * user->mx * (user->mx - 1) * 3));
8949566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->S));
895c4762a1bSJed Brown
8969566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork));
8979566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, m + user->mx * user->mx * user->mx));
8989566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork));
899c4762a1bSJed Brown
9009566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork));
9019566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Av, MAT_SHARE_NONZERO_PATTERN, &user->Avwork));
902c4762a1bSJed Brown
9039566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->d));
9049566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->d, PETSC_DECIDE, user->ndata * user->nt));
9059566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->d));
906c4762a1bSJed Brown
9079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Swork));
9089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Av_u));
9099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Twork));
9109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Rwork));
9119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ywork));
9129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->uwork));
9139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->js_diag));
9149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c, &user->cwork));
9159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->d, &user->dwork));
916c4762a1bSJed Brown
917c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing A*x */
9189566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m * user->nt, user, &user->Js));
919*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (PetscErrorCodeFn *)StateMatMult));
920*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_DUPLICATE, (PetscErrorCodeFn *)StateMatDuplicate));
921*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatMultTranspose));
922*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)StateMatGetDiagonal));
923c4762a1bSJed Brown
924c4762a1bSJed Brown /* Diagonal blocks of user->Js */
9259566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsBlock));
926*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (PetscErrorCodeFn *)StateMatBlockMult));
927c4762a1bSJed Brown /* Blocks are symmetric */
928*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatBlockMult));
929c4762a1bSJed Brown
930c4762a1bSJed Brown /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
931c4762a1bSJed Brown D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
932c4762a1bSJed Brown This is an SSOR preconditioner for user->JsBlock. */
9339566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsBlockPrec));
934*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT, (PetscErrorCodeFn *)StateMatBlockPrecMult));
935c4762a1bSJed Brown /* JsBlockPrec is symmetric */
936*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatBlockPrecMult));
937e9b2f0d4SBarry Smith PetscCall(MatSetOption(user->JsBlockPrec, MAT_SYMMETRIC, PETSC_TRUE));
9389566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsBlockPrec, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
939c4762a1bSJed Brown
940c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */
9419566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m, user, &user->Jd));
942*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (PetscErrorCodeFn *)DesignMatMult));
943*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)DesignMatMultTranspose));
944c4762a1bSJed Brown
945c4762a1bSJed Brown /* User-defined routines for computing user->Js\x and user->Js^T\x*/
9469566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m * user->nt, user, &user->JsInv));
947*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (PetscErrorCodeFn *)StateMatInvMult));
948*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatInvTransposeMult));
949c4762a1bSJed Brown
950c4762a1bSJed Brown /* Solver options and tolerances */
9519566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver));
9529566063dSJacob Faibussowitsch PetscCall(KSPSetType(user->solver, KSPCG));
9539566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->JsBlockPrec));
9549566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(user->solver, PETSC_FALSE));
9559566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, 1e3, 500));
9569566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver));
9579566063dSJacob Faibussowitsch PetscCall(KSPGetPC(user->solver, &user->prec));
9589566063dSJacob Faibussowitsch PetscCall(PCSetType(user->prec, PCSHELL));
959c4762a1bSJed Brown
9609566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(user->prec, StateMatBlockPrecMult));
9619566063dSJacob Faibussowitsch PetscCall(PCShellSetApplyTranspose(user->prec, StateMatBlockPrecMult));
9629566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(user->prec, user));
963c4762a1bSJed Brown
964c4762a1bSJed Brown /* Create scatter from y to y_1,y_2,...,y_nt */
9659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->m, &user->yi_scatter));
9669566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &yi));
9679566063dSJacob Faibussowitsch PetscCall(VecSetSizes(yi, PETSC_DECIDE, user->mx * user->mx * user->mx));
9689566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(yi));
9699566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yi));
9709566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yiwork));
971c4762a1bSJed Brown
9729566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->y, &lo2, &hi2));
973c4762a1bSJed Brown istart = 0;
974c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
9759566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->yi[i], &lo, &hi));
9769566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_yi));
9779566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_y));
9789566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y, is_from_y, user->yi[i], is_to_yi, &user->yi_scatter[i]));
979c4762a1bSJed Brown istart = istart + hi - lo;
9809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_yi));
9819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y));
982c4762a1bSJed Brown }
9839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yi));
984c4762a1bSJed Brown
985c4762a1bSJed Brown /* Create scatter from d to d_1,d_2,...,d_ns */
9869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns * user->ndata, &user->di_scatter));
9879566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &di));
9889566063dSJacob Faibussowitsch PetscCall(VecSetSizes(di, PETSC_DECIDE, user->ndata));
9899566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(di));
9909566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(di, user->ns, &user->di));
9919566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->d, &lo2, &hi2));
992c4762a1bSJed Brown istart = 0;
993c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
9949566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->di[i], &lo, &hi));
9959566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_di));
9969566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_d));
9979566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->d, is_from_d, user->di[i], is_to_di, &user->di_scatter[i]));
998c4762a1bSJed Brown istart = istart + hi - lo;
9999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_di));
10009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_d));
1001c4762a1bSJed Brown }
10029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&di));
1003c4762a1bSJed Brown
1004c4762a1bSJed Brown /* Assemble RHS of forward problem */
10059566063dSJacob Faibussowitsch PetscCall(VecCopy(bc, user->yiwork[0]));
100648a46eb9SPierre Jolivet for (i = 1; i < user->nt; i++) PetscCall(VecSet(user->yiwork[i], 0.0));
10079566063dSJacob Faibussowitsch PetscCall(Gather_i(user->q, user->yiwork, user->yi_scatter, user->nt));
10089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bc));
1009c4762a1bSJed Brown
1010c4762a1bSJed Brown /* Compute true state function yt given ut */
10119566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ytrue));
10129566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->ytrue, PETSC_DECIDE, n * user->nt));
10139566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->ytrue));
1014c4762a1bSJed Brown
1015c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */
10169566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
10179566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); /* Note: user->utrue */
10189566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
10199566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
1020c4762a1bSJed Brown
1021c20d7725SJed Brown /* Symbolic DSG = Div * Grad */
10229566063dSJacob Faibussowitsch PetscCall(MatProductCreate(user->Div, user->Grad, NULL, &user->DSG));
10239566063dSJacob Faibussowitsch PetscCall(MatProductSetType(user->DSG, MATPRODUCT_AB));
10249566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(user->DSG, "default"));
1025fb842aefSJose E. Roman PetscCall(MatProductSetFill(user->DSG, PETSC_DETERMINE));
10269566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(user->DSG));
10279566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(user->DSG));
1028c20d7725SJed Brown
1029c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE;
1030c4762a1bSJed Brown
1031c20d7725SJed Brown /* Next form DSG = Div*Grad */
10329566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
10339566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u));
1034c4762a1bSJed Brown if (user->dsg_formed) {
10359566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG));
1036c4762a1bSJed Brown } else {
1037fb842aefSJose E. Roman PetscCall(MatMatMult(user->Div, user->Grad, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &user->DSG));
1038c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE;
1039c4762a1bSJed Brown }
1040c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */
10419566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG, user->ht));
10429566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG, 1.0));
1043c4762a1bSJed Brown
1044c4762a1bSJed Brown /* Now solve for ytrue */
10459566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->ytrue));
1046c4762a1bSJed Brown
1047c4762a1bSJed Brown /* Initial guess y0 for state given u0 */
1048c4762a1bSJed Brown
1049c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */
10509566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
10519566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u)); /* Note: user->u */
10529566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
10539566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
1054c4762a1bSJed Brown
1055c4762a1bSJed Brown /* Next form DSG = Div*S*Grad */
10569566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
10579566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u));
1058c4762a1bSJed Brown if (user->dsg_formed) {
10599566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG));
1060c4762a1bSJed Brown } else {
1061fb842aefSJose E. Roman PetscCall(MatMatMult(user->Div, user->Grad, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &user->DSG));
1062c20d7725SJed Brown
1063c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE;
1064c4762a1bSJed Brown }
1065c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */
10669566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG, user->ht));
10679566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG, 1.0));
1068c4762a1bSJed Brown
1069c4762a1bSJed Brown /* Now solve for y */
10709566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->y));
1071c4762a1bSJed Brown
1072c4762a1bSJed Brown /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
10739566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Qblock));
10749566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Qblock, PETSC_DECIDE, PETSC_DECIDE, user->ndata, n));
10759566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Qblock));
10769566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Qblock, 8, NULL, 8, NULL));
10779566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Qblock, 8, NULL));
1078c4762a1bSJed Brown
1079c4762a1bSJed Brown for (i = 0; i < user->mx; i++) {
1080c4762a1bSJed Brown x[i] = h * (i + 0.5);
1081c4762a1bSJed Brown y[i] = h * (i + 0.5);
1082c4762a1bSJed Brown z[i] = h * (i + 0.5);
1083c4762a1bSJed Brown }
1084c4762a1bSJed Brown
10859566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Qblock, &istart, &iend));
10869371c9d4SSatish Balay nx = user->mx;
10879371c9d4SSatish Balay ny = user->mx;
10889371c9d4SSatish Balay nz = user->mx;
1089c4762a1bSJed Brown for (i = istart; i < iend; i++) {
1090c4762a1bSJed Brown xri = xr[i];
1091c4762a1bSJed Brown im = 0;
1092c4762a1bSJed Brown xim = x[im];
1093c4762a1bSJed Brown while (xri > xim && im < nx) {
1094c4762a1bSJed Brown im = im + 1;
1095c4762a1bSJed Brown xim = x[im];
1096c4762a1bSJed Brown }
1097c4762a1bSJed Brown indx1 = im - 1;
1098c4762a1bSJed Brown indx2 = im;
1099c4762a1bSJed Brown dx1 = xri - x[indx1];
1100c4762a1bSJed Brown dx2 = x[indx2] - xri;
1101c4762a1bSJed Brown
1102c4762a1bSJed Brown yri = yr[i];
1103c4762a1bSJed Brown im = 0;
1104c4762a1bSJed Brown yim = y[im];
1105c4762a1bSJed Brown while (yri > yim && im < ny) {
1106c4762a1bSJed Brown im = im + 1;
1107c4762a1bSJed Brown yim = y[im];
1108c4762a1bSJed Brown }
1109c4762a1bSJed Brown indy1 = im - 1;
1110c4762a1bSJed Brown indy2 = im;
1111c4762a1bSJed Brown dy1 = yri - y[indy1];
1112c4762a1bSJed Brown dy2 = y[indy2] - yri;
1113c4762a1bSJed Brown
1114c4762a1bSJed Brown zri = zr[i];
1115c4762a1bSJed Brown im = 0;
1116c4762a1bSJed Brown zim = z[im];
1117c4762a1bSJed Brown while (zri > zim && im < nz) {
1118c4762a1bSJed Brown im = im + 1;
1119c4762a1bSJed Brown zim = z[im];
1120c4762a1bSJed Brown }
1121c4762a1bSJed Brown indz1 = im - 1;
1122c4762a1bSJed Brown indz2 = im;
1123c4762a1bSJed Brown dz1 = zri - z[indz1];
1124c4762a1bSJed Brown dz2 = z[indz2] - zri;
1125c4762a1bSJed Brown
1126c4762a1bSJed Brown Dx = x[indx2] - x[indx1];
1127c4762a1bSJed Brown Dy = y[indy2] - y[indy1];
1128c4762a1bSJed Brown Dz = z[indz2] - z[indz1];
1129c4762a1bSJed Brown
1130c4762a1bSJed Brown j = indx1 + indy1 * nx + indz1 * nx * ny;
1131c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
11329566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1133c4762a1bSJed Brown
1134c4762a1bSJed Brown j = indx1 + indy1 * nx + indz2 * nx * ny;
1135c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
11369566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1137c4762a1bSJed Brown
1138c4762a1bSJed Brown j = indx1 + indy2 * nx + indz1 * nx * ny;
1139c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
11409566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1141c4762a1bSJed Brown
1142c4762a1bSJed Brown j = indx1 + indy2 * nx + indz2 * nx * ny;
1143c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
11449566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1145c4762a1bSJed Brown
1146c4762a1bSJed Brown j = indx2 + indy1 * nx + indz1 * nx * ny;
1147c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
11489566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1149c4762a1bSJed Brown
1150c4762a1bSJed Brown j = indx2 + indy1 * nx + indz2 * nx * ny;
1151c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
11529566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1153c4762a1bSJed Brown
1154c4762a1bSJed Brown j = indx2 + indy2 * nx + indz1 * nx * ny;
1155c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
11569566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1157c4762a1bSJed Brown
1158c4762a1bSJed Brown j = indx2 + indy2 * nx + indz2 * nx * ny;
1159c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
11609566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES));
1161c4762a1bSJed Brown }
11629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Qblock, MAT_FINAL_ASSEMBLY));
11639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Qblock, MAT_FINAL_ASSEMBLY));
1164c4762a1bSJed Brown
11659566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Qblock, MAT_INITIAL_MATRIX, &user->QblockT));
11669566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->L, MAT_INITIAL_MATRIX, &user->LT));
1167c4762a1bSJed Brown
1168c4762a1bSJed Brown /* Add noise to the measurement data */
11699566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork, 1.0));
11709566063dSJacob Faibussowitsch PetscCall(VecAYPX(user->ywork, user->noise, user->ytrue));
11719566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt));
1172c4762a1bSJed Brown for (j = 0; j < user->ns; j++) {
1173c4762a1bSJed Brown i = user->sample_times[j];
11749566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock, user->yiwork[i], user->di[j]));
1175c4762a1bSJed Brown }
11769566063dSJacob Faibussowitsch PetscCall(Gather_i(user->d, user->di, user->di_scatter, user->ns));
1177c4762a1bSJed Brown
1178c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11799566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver));
11809566063dSJacob Faibussowitsch PetscCall(PetscFree(x));
11819566063dSJacob Faibussowitsch PetscCall(PetscFree(y));
11829566063dSJacob Faibussowitsch PetscCall(PetscFree(z));
11833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1184c4762a1bSJed Brown }
1185c4762a1bSJed Brown
ParabolicDestroy(AppCtx * user)1186d71ae5a4SJacob Faibussowitsch PetscErrorCode ParabolicDestroy(AppCtx *user)
1187d71ae5a4SJacob Faibussowitsch {
1188c4762a1bSJed Brown PetscInt i;
1189c4762a1bSJed Brown
1190c4762a1bSJed Brown PetscFunctionBegin;
11919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Qblock));
11929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->QblockT));
11939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div));
11949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork));
11959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad));
11969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Av));
11979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Avwork));
11989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->AvT));
11999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->DSG));
12009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L));
12019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->LT));
12029566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver));
12039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js));
12049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd));
12059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv));
12069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock));
12079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlockPrec));
12089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u));
12099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork));
12109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue));
12119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y));
12129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork));
12139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue));
12149566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yi));
12159566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yiwork));
12169566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->ns, &user->di));
12179566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi));
12189566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yiwork));
12199566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di));
12209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c));
12219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork));
12229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur));
12239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q));
12249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d));
12259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork));
12269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork));
12279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->S));
12289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Swork));
12299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Av_u));
12309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Twork));
12319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Rwork));
12329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag));
12339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is));
12349566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is));
12359566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter));
12369566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter));
123748a46eb9SPierre Jolivet for (i = 0; i < user->nt; i++) PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
123848a46eb9SPierre Jolivet for (i = 0; i < user->ns; i++) PetscCall(VecScatterDestroy(&user->di_scatter[i]));
12399566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter));
12409566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di_scatter));
12419566063dSJacob Faibussowitsch PetscCall(PetscFree(user->sample_times));
12423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1243c4762a1bSJed Brown }
1244c4762a1bSJed Brown
ParabolicMonitor(Tao tao,void * ptr)1245d71ae5a4SJacob Faibussowitsch PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1246d71ae5a4SJacob Faibussowitsch {
1247c4762a1bSJed Brown Vec X;
1248c4762a1bSJed Brown PetscReal unorm, ynorm;
1249c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
1250c4762a1bSJed Brown
1251c4762a1bSJed Brown PetscFunctionBegin;
12529566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X));
12539566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
12549566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue));
12559566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue));
12569566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork, NORM_2, &unorm));
12579566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork, NORM_2, &ynorm));
12589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm));
12593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1260c4762a1bSJed Brown }
1261c4762a1bSJed Brown
1262c4762a1bSJed Brown /*TEST
1263c4762a1bSJed Brown
1264c4762a1bSJed Brown build:
1265c4762a1bSJed Brown requires: !complex
1266c4762a1bSJed Brown
1267c4762a1bSJed Brown test:
126810978b7dSBarry Smith args: -tao_monitor_constraint_norm -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30
1269c4762a1bSJed Brown requires: !single
1270c4762a1bSJed Brown
1271c4762a1bSJed Brown TEST*/
1272