1c4762a1bSJed Brown #include <petsctao.h>
2c4762a1bSJed Brown
3c4762a1bSJed Brown typedef struct {
4c4762a1bSJed Brown PetscInt n; /* Number of variables */
5c4762a1bSJed Brown PetscInt m; /* Number of constraints */
6c4762a1bSJed Brown PetscInt mx; /* grid points in each direction */
7c4762a1bSJed Brown PetscInt nt; /* Number of time steps */
8c4762a1bSJed Brown PetscInt ndata; /* Number of data points per sample */
9c4762a1bSJed Brown IS s_is;
10c4762a1bSJed Brown IS d_is;
11c4762a1bSJed Brown VecScatter state_scatter;
12c4762a1bSJed Brown VecScatter design_scatter;
13c4762a1bSJed Brown VecScatter *uxi_scatter, *uyi_scatter, *ux_scatter, *uy_scatter, *ui_scatter;
14c4762a1bSJed Brown VecScatter *yi_scatter;
15c4762a1bSJed Brown
16c4762a1bSJed Brown Mat Js, Jd, JsBlockPrec, JsInv, JsBlock;
17c4762a1bSJed Brown PetscBool jformed, c_formed;
18c4762a1bSJed Brown
19c4762a1bSJed Brown PetscReal alpha; /* Regularization parameter */
20c4762a1bSJed Brown PetscReal gamma;
21c4762a1bSJed Brown PetscReal ht; /* Time step */
22c4762a1bSJed Brown PetscReal T; /* Final time */
23c4762a1bSJed Brown Mat Q, QT;
24c4762a1bSJed Brown Mat L, LT;
25c4762a1bSJed Brown Mat Div, Divwork, Divxy[2];
26c4762a1bSJed Brown Mat Grad, Gradxy[2];
27c4762a1bSJed Brown Mat M;
28c4762a1bSJed Brown Mat *C, *Cwork;
29c4762a1bSJed Brown /* Mat Hs,Hd,Hsd; */
30c4762a1bSJed Brown Vec q;
31c4762a1bSJed Brown Vec ur; /* reference */
32c4762a1bSJed Brown
33c4762a1bSJed Brown Vec d;
34c4762a1bSJed Brown Vec dwork;
35c4762a1bSJed Brown
36c4762a1bSJed Brown Vec y; /* state variables */
37c4762a1bSJed Brown Vec ywork;
38c4762a1bSJed Brown Vec ytrue;
39c4762a1bSJed Brown Vec *yi, *yiwork, *ziwork;
40c4762a1bSJed Brown Vec *uxi, *uyi, *uxiwork, *uyiwork, *ui, *uiwork;
41c4762a1bSJed Brown
42c4762a1bSJed Brown Vec u; /* design variables */
43c4762a1bSJed Brown Vec uwork, vwork;
44c4762a1bSJed Brown Vec utrue;
45c4762a1bSJed Brown
46c4762a1bSJed Brown Vec js_diag;
47c4762a1bSJed Brown
48c4762a1bSJed Brown Vec c; /* constraint vector */
49c4762a1bSJed Brown Vec cwork;
50c4762a1bSJed Brown
51c4762a1bSJed Brown Vec lwork;
52c4762a1bSJed Brown
53c4762a1bSJed Brown KSP solver;
54c4762a1bSJed Brown PC prec;
55c4762a1bSJed Brown PetscInt block_index;
56c4762a1bSJed Brown
57c4762a1bSJed Brown PetscInt ksp_its;
58c4762a1bSJed Brown PetscInt ksp_its_initial;
59c4762a1bSJed Brown } AppCtx;
60c4762a1bSJed Brown
61c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
62c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
63c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
64c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void *);
65c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void *);
66c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
67c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
68c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
69c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
70c4762a1bSJed Brown PetscErrorCode HyperbolicInitialize(AppCtx *user);
71c4762a1bSJed Brown PetscErrorCode HyperbolicDestroy(AppCtx *user);
72c4762a1bSJed Brown PetscErrorCode HyperbolicMonitor(Tao, void *);
73c4762a1bSJed Brown
74c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat, Vec, Vec);
75c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat, Vec, Vec);
76c4762a1bSJed Brown PetscErrorCode StateMatBlockMultTranspose(Mat, Vec, Vec);
77c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat, Vec, Vec);
78c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat, Vec);
79c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat, MatDuplicateOption, Mat *);
80c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat, Vec, Vec);
81c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat, Vec, Vec);
82c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC, Vec, Vec);
83c4762a1bSJed Brown
84c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat, Vec, Vec);
85c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat, Vec, Vec);
86c4762a1bSJed Brown
87c4762a1bSJed Brown PetscErrorCode Scatter_yi(Vec, Vec *, VecScatter *, PetscInt); /* y to y1,y2,...,y_nt */
88c4762a1bSJed Brown PetscErrorCode Gather_yi(Vec, Vec *, VecScatter *, PetscInt);
89c4762a1bSJed Brown PetscErrorCode Scatter_uxi_uyi(Vec, Vec *, VecScatter *, Vec *, VecScatter *, PetscInt); /* u to ux_1,uy_1,ux_2,uy_2,...,u */
90c4762a1bSJed Brown PetscErrorCode Gather_uxi_uyi(Vec, Vec *, VecScatter *, Vec *, VecScatter *, PetscInt);
91c4762a1bSJed Brown
92c4762a1bSJed Brown static char help[] = "";
93c4762a1bSJed Brown
main(int argc,char ** argv)94d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
95d71ae5a4SJacob Faibussowitsch {
96c4762a1bSJed Brown Vec x, x0;
97c4762a1bSJed Brown Tao tao;
98c4762a1bSJed Brown AppCtx user;
99c4762a1bSJed Brown IS is_allstate, is_alldesign;
100c4762a1bSJed Brown PetscInt lo, hi, hi2, lo2, ksp_old;
101c4762a1bSJed Brown PetscInt ntests = 1;
102c4762a1bSJed Brown PetscInt i;
103c4762a1bSJed Brown PetscLogStage stages[1];
104c4762a1bSJed Brown
105327415f7SBarry Smith PetscFunctionBeginUser;
106c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
107c4762a1bSJed Brown user.mx = 32;
108d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "hyperbolic example", NULL);
1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL));
110c4762a1bSJed Brown user.nt = 16;
1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-nt", "Number of time steps", "", user.nt, &user.nt, NULL));
112c4762a1bSJed Brown user.ndata = 64;
1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL));
114c4762a1bSJed Brown user.alpha = 10.0;
1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL));
116c4762a1bSJed Brown user.T = 1.0 / 32.0;
1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Tfinal", "Final time", "", user.T, &user.T, NULL));
1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL));
119d0609cedSBarry Smith PetscOptionsEnd();
120c4762a1bSJed Brown
121c4762a1bSJed Brown user.m = user.mx * user.mx * user.nt; /* number of constraints */
122c4762a1bSJed Brown user.n = user.mx * user.mx * 3 * user.nt; /* number of variables */
123c4762a1bSJed Brown user.ht = user.T / user.nt; /* Time step */
124c4762a1bSJed Brown user.gamma = user.T * user.ht / (user.mx * user.mx);
125c4762a1bSJed Brown
1269566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.u));
1279566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.y));
1289566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.c));
1299566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.u, PETSC_DECIDE, user.n - user.m));
1309566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.y, PETSC_DECIDE, user.m));
1319566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.c, PETSC_DECIDE, user.m));
1329566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.u));
1339566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.y));
1349566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.c));
135c4762a1bSJed Brown
136c4762a1bSJed Brown /* Create scatters for reduced spaces.
137c4762a1bSJed Brown If the state vector y and design vector u are partitioned as
138c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
139c4762a1bSJed Brown then the solution vector x is organized as
140c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np].
141c4762a1bSJed Brown The index sets user.s_is and user.d_is correspond to the indices of the
142c4762a1bSJed Brown state and design variables owned by the current processor.
143c4762a1bSJed Brown */
1449566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
145c4762a1bSJed Brown
1469566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.y, &lo, &hi));
1479566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.u, &lo2, &hi2));
148c4762a1bSJed Brown
1499566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate));
1509566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user.s_is));
151c4762a1bSJed Brown
1529566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign));
1539566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user.d_is));
154c4762a1bSJed Brown
1559566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, hi - lo + hi2 - lo2, user.n));
1569566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x));
157c4762a1bSJed Brown
1589566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.s_is, user.y, is_allstate, &user.state_scatter));
1599566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.d_is, user.u, is_alldesign, &user.design_scatter));
1609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign));
1619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate));
162c4762a1bSJed Brown
163c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
1649566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1659566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLCL));
166c4762a1bSJed Brown
167c4762a1bSJed Brown /* Set up initial vectors and matrices */
1689566063dSJacob Faibussowitsch PetscCall(HyperbolicInitialize(&user));
169c4762a1bSJed Brown
1709566063dSJacob Faibussowitsch PetscCall(Gather(x, user.y, user.state_scatter, user.u, user.design_scatter));
1719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &x0));
1729566063dSJacob Faibussowitsch PetscCall(VecCopy(x, x0));
173c4762a1bSJed Brown
174c4762a1bSJed Brown /* Set solution vector with an initial guess */
1759566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
1769566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user));
1779566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
1789566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
1799566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, &user));
1809566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
1819566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
1829566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is));
183c4762a1bSJed Brown
184c4762a1bSJed Brown /* SOLVE THE APPLICATION */
1859566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials", &stages[0]));
1869566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0]));
187c4762a1bSJed Brown user.ksp_its_initial = user.ksp_its;
188c4762a1bSJed Brown ksp_old = user.ksp_its;
189c4762a1bSJed Brown for (i = 0; i < ntests; i++) {
1909566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
19163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its - ksp_old));
1929566063dSJacob Faibussowitsch PetscCall(VecCopy(x0, x));
1939566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
194c4762a1bSJed Brown }
1959566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop());
1969566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)x));
1979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: "));
19863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial));
19963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total KSP iterations over %" PetscInt_FMT " trial(s): ", ntests));
20063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its));
2019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations per trial: "));
20263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", (user.ksp_its - user.ksp_its_initial) / ntests));
203c4762a1bSJed Brown
2049566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
2059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
2069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0));
2079566063dSJacob Faibussowitsch PetscCall(HyperbolicDestroy(&user));
2089566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
209b122ec5aSJacob Faibussowitsch return 0;
210c4762a1bSJed Brown }
211c4762a1bSJed Brown /* ------------------------------------------------------------------- */
212c4762a1bSJed Brown /*
213c4762a1bSJed Brown dwork = Qy - d
214c4762a1bSJed Brown lwork = L*(u-ur).^2
215c4762a1bSJed Brown f = 1/2 * (dwork.dork + alpha*y.lwork)
216c4762a1bSJed Brown */
FormFunction(Tao tao,Vec X,PetscReal * f,void * ptr)217d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr)
218d71ae5a4SJacob Faibussowitsch {
219c4762a1bSJed Brown PetscReal d1 = 0, d2 = 0;
220c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
221c4762a1bSJed Brown
222c4762a1bSJed Brown PetscFunctionBegin;
2239566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2249566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork));
2259566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2269566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1));
227c4762a1bSJed Brown
2289566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
2299566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->uwork, user->uwork));
2309566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork));
2319566063dSJacob Faibussowitsch PetscCall(VecDot(user->y, user->lwork, &d2));
232c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2);
2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
234c4762a1bSJed Brown }
235c4762a1bSJed Brown
236c4762a1bSJed Brown /* ------------------------------------------------------------------- */
237c4762a1bSJed Brown /*
238c4762a1bSJed Brown state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2
239c4762a1bSJed Brown design: g_d = alpha*(L'y).*(u-ur)
240c4762a1bSJed Brown */
FormGradient(Tao tao,Vec X,Vec G,void * ptr)241d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr)
242d71ae5a4SJacob Faibussowitsch {
243c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
244c4762a1bSJed Brown
245c4762a1bSJed Brown PetscFunctionBegin;
2469566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2479566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork));
2489566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
249c4762a1bSJed Brown
2509566063dSJacob Faibussowitsch PetscCall(MatMult(user->QT, user->dwork, user->ywork));
251c4762a1bSJed Brown
2529566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->y, user->uwork));
2539566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u));
2549566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork));
2559566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha));
256c4762a1bSJed Brown
2579566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork));
2589566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->vwork, user->lwork));
2599566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork));
260c4762a1bSJed Brown
2619566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
263c4762a1bSJed Brown }
264c4762a1bSJed Brown
FormFunctionGradient(Tao tao,Vec X,PetscReal * f,Vec G,void * ptr)265d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
266d71ae5a4SJacob Faibussowitsch {
267c4762a1bSJed Brown PetscReal d1, d2;
268c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
269c4762a1bSJed Brown
270c4762a1bSJed Brown PetscFunctionBegin;
2719566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2729566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork));
2739566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
274c4762a1bSJed Brown
2759566063dSJacob Faibussowitsch PetscCall(MatMult(user->QT, user->dwork, user->ywork));
276c4762a1bSJed Brown
2779566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1));
278c4762a1bSJed Brown
2799566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->y, user->uwork));
2809566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u));
2819566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork));
2829566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha));
283c4762a1bSJed Brown
2849566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork));
2859566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->vwork, user->lwork));
2869566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork));
287c4762a1bSJed Brown
2889566063dSJacob Faibussowitsch PetscCall(VecDot(user->y, user->lwork, &d2));
289c4762a1bSJed Brown
290c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2);
2919566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown
295c4762a1bSJed Brown /* ------------------------------------------------------------------- */
296c4762a1bSJed Brown /* A
297c4762a1bSJed Brown MatShell object
298c4762a1bSJed Brown */
FormJacobianState(Tao tao,Vec X,Mat J,Mat JPre,Mat JInv,void * ptr)299d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
300d71ae5a4SJacob Faibussowitsch {
301c4762a1bSJed Brown PetscInt i;
302c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
303c4762a1bSJed Brown
304c4762a1bSJed Brown PetscFunctionBegin;
3059566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
3069566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->u, user->ui, user->ui_scatter, user->nt));
3079566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
308c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
3099566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Divxy[0], user->C[i], SUBSET_NONZERO_PATTERN));
3109566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Divxy[1], user->Cwork[i], SAME_NONZERO_PATTERN));
311c4762a1bSJed Brown
3129566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i]));
3139566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i]));
3149566063dSJacob Faibussowitsch PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], SUBSET_NONZERO_PATTERN));
3159566063dSJacob Faibussowitsch PetscCall(MatScale(user->C[i], user->ht));
3169566063dSJacob Faibussowitsch PetscCall(MatShift(user->C[i], 1.0));
317c4762a1bSJed Brown }
3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
319c4762a1bSJed Brown }
320c4762a1bSJed Brown
321c4762a1bSJed Brown /* ------------------------------------------------------------------- */
322c4762a1bSJed Brown /* B */
FormJacobianDesign(Tao tao,Vec X,Mat J,void * ptr)323d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
324d71ae5a4SJacob Faibussowitsch {
325c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
326c4762a1bSJed Brown
327c4762a1bSJed Brown PetscFunctionBegin;
3289566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
330c4762a1bSJed Brown }
331c4762a1bSJed Brown
StateMatMult(Mat J_shell,Vec X,Vec Y)332d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
333d71ae5a4SJacob Faibussowitsch {
334c4762a1bSJed Brown PetscInt i;
335c4762a1bSJed Brown AppCtx *user;
336c4762a1bSJed Brown
337c4762a1bSJed Brown PetscFunctionBegin;
3389566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
3399566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
340c4762a1bSJed Brown user->block_index = 0;
3419566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
342c4762a1bSJed Brown
343c4762a1bSJed Brown for (i = 1; i < user->nt; i++) {
344c4762a1bSJed Brown user->block_index = i;
3459566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
3469566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1]));
3479566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1]));
348c4762a1bSJed Brown }
3499566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown
StateMatMultTranspose(Mat J_shell,Vec X,Vec Y)353d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
354d71ae5a4SJacob Faibussowitsch {
355c4762a1bSJed Brown PetscInt i;
356c4762a1bSJed Brown AppCtx *user;
357c4762a1bSJed Brown
358c4762a1bSJed Brown PetscFunctionBegin;
3599566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
3609566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
361c4762a1bSJed Brown
362c4762a1bSJed Brown for (i = 0; i < user->nt - 1; i++) {
363c4762a1bSJed Brown user->block_index = i;
3649566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i]));
3659566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i + 1], user->ziwork[i + 1]));
3669566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i + 1]));
367c4762a1bSJed Brown }
368c4762a1bSJed Brown
369c4762a1bSJed Brown i = user->nt - 1;
370c4762a1bSJed Brown user->block_index = i;
3719566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i]));
3729566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
3733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown
StateMatBlockMult(Mat J_shell,Vec X,Vec Y)376d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
377d71ae5a4SJacob Faibussowitsch {
378c4762a1bSJed Brown PetscInt i;
379c4762a1bSJed Brown AppCtx *user;
380c4762a1bSJed Brown
381c4762a1bSJed Brown PetscFunctionBegin;
3829566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
383c4762a1bSJed Brown i = user->block_index;
3849566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], X, user->uxi[i]));
3859566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], X, user->uyi[i]));
3869566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
3879566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->uiwork[i], Y));
3889566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, user->ht, X));
3893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
390c4762a1bSJed Brown }
391c4762a1bSJed Brown
StateMatBlockMultTranspose(Mat J_shell,Vec X,Vec Y)392d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y)
393d71ae5a4SJacob Faibussowitsch {
394c4762a1bSJed Brown PetscInt i;
395c4762a1bSJed Brown AppCtx *user;
396c4762a1bSJed Brown
397c4762a1bSJed Brown PetscFunctionBegin;
3989566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
399c4762a1bSJed Brown i = user->block_index;
4009566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, X, user->uiwork[i]));
4019566063dSJacob Faibussowitsch PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4029566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->uxi[i], user->uxiwork[i]));
4039566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->uyi[i], user->uyiwork[i]));
4049566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Y, 1.0, user->uxiwork[i], user->uyiwork[i]));
4059566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, user->ht, X));
4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
407c4762a1bSJed Brown }
408c4762a1bSJed Brown
DesignMatMult(Mat J_shell,Vec X,Vec Y)409d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
410d71ae5a4SJacob Faibussowitsch {
411c4762a1bSJed Brown PetscInt i;
412c4762a1bSJed Brown AppCtx *user;
413c4762a1bSJed Brown
414c4762a1bSJed Brown PetscFunctionBegin;
4159566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
4169566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt));
4179566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(X, user->uxiwork, user->uxi_scatter, user->uyiwork, user->uyi_scatter, user->nt));
418c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
4199566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i]));
4209566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i]));
4219566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4229566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->uiwork[i], user->ziwork[i]));
4239566063dSJacob Faibussowitsch PetscCall(VecScale(user->ziwork[i], user->ht));
424c4762a1bSJed Brown }
4259566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->ziwork, user->yi_scatter, user->nt));
4263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
427c4762a1bSJed Brown }
428c4762a1bSJed Brown
DesignMatMultTranspose(Mat J_shell,Vec X,Vec Y)429d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
430d71ae5a4SJacob Faibussowitsch {
431c4762a1bSJed Brown PetscInt i;
432c4762a1bSJed Brown AppCtx *user;
433c4762a1bSJed Brown
434c4762a1bSJed Brown PetscFunctionBegin;
4359566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
4369566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt));
4379566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yiwork, user->yi_scatter, user->nt));
438c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
4399566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->yiwork[i], user->uiwork[i]));
4409566063dSJacob Faibussowitsch PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4419566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i]));
4429566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i]));
4439566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i]));
4449566063dSJacob Faibussowitsch PetscCall(VecScale(user->uiwork[i], user->ht));
445c4762a1bSJed Brown }
4469566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->uiwork, user->ui_scatter, user->nt));
4473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
448c4762a1bSJed Brown }
449c4762a1bSJed Brown
StateMatBlockPrecMult(PC PC_shell,Vec X,Vec Y)450d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
451d71ae5a4SJacob Faibussowitsch {
452c4762a1bSJed Brown PetscInt i;
453c4762a1bSJed Brown AppCtx *user;
454c4762a1bSJed Brown
455c4762a1bSJed Brown PetscFunctionBegin;
4569566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell, &user));
457c4762a1bSJed Brown i = user->block_index;
458966bd95aSPierre Jolivet PetscCheck(user->c_formed, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed");
4599566063dSJacob Faibussowitsch PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y));
4603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
461c4762a1bSJed Brown }
462c4762a1bSJed Brown
StateMatBlockPrecMultTranspose(PC PC_shell,Vec X,Vec Y)463d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y)
464d71ae5a4SJacob Faibussowitsch {
465c4762a1bSJed Brown PetscInt i;
466c4762a1bSJed Brown AppCtx *user;
467c4762a1bSJed Brown
468c4762a1bSJed Brown PetscFunctionBegin;
4699566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell, &user));
470c4762a1bSJed Brown
471c4762a1bSJed Brown i = user->block_index;
472966bd95aSPierre Jolivet PetscCheck(user->c_formed, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed");
4739566063dSJacob Faibussowitsch PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y));
4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
475c4762a1bSJed Brown }
476c4762a1bSJed Brown
StateMatInvMult(Mat J_shell,Vec X,Vec Y)477d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
478d71ae5a4SJacob Faibussowitsch {
479c4762a1bSJed Brown AppCtx *user;
480c4762a1bSJed Brown PetscInt its, i;
481c4762a1bSJed Brown
482c4762a1bSJed Brown PetscFunctionBegin;
4839566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
484c4762a1bSJed Brown
485c4762a1bSJed Brown if (Y == user->ytrue) {
486c4762a1bSJed Brown /* First solve is done using true solution to set up problem */
487606f75f6SBarry Smith PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, PETSC_CURRENT, PETSC_CURRENT));
488c4762a1bSJed Brown } else {
489606f75f6SBarry Smith PetscCall(KSPSetTolerances(user->solver, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
490c4762a1bSJed Brown }
4919566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
4929566063dSJacob Faibussowitsch PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt));
4939566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
494c4762a1bSJed Brown
495c4762a1bSJed Brown user->block_index = 0;
4969566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[0], user->yiwork[0]));
497c4762a1bSJed Brown
4989566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
499c4762a1bSJed Brown user->ksp_its = user->ksp_its + its;
500c4762a1bSJed Brown for (i = 1; i < user->nt; i++) {
5019566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yiwork[i - 1], user->ziwork[i - 1]));
5029566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i - 1]));
503c4762a1bSJed Brown user->block_index = i;
5049566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i]));
505c4762a1bSJed Brown
5069566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
507c4762a1bSJed Brown user->ksp_its = user->ksp_its + its;
508c4762a1bSJed Brown }
5099566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
5103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
511c4762a1bSJed Brown }
512c4762a1bSJed Brown
StateMatInvTransposeMult(Mat J_shell,Vec X,Vec Y)513d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
514d71ae5a4SJacob Faibussowitsch {
515c4762a1bSJed Brown AppCtx *user;
516c4762a1bSJed Brown PetscInt its, i;
517c4762a1bSJed Brown
518c4762a1bSJed Brown PetscFunctionBegin;
5199566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
520c4762a1bSJed Brown
5219566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt));
5229566063dSJacob Faibussowitsch PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt));
5239566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
524c4762a1bSJed Brown
525c4762a1bSJed Brown i = user->nt - 1;
526c4762a1bSJed Brown user->block_index = i;
5279566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(user->solver, user->yi[i], user->yiwork[i]));
528c4762a1bSJed Brown
5299566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
530c4762a1bSJed Brown user->ksp_its = user->ksp_its + its;
531c4762a1bSJed Brown
532c4762a1bSJed Brown for (i = user->nt - 2; i >= 0; i--) {
5339566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yiwork[i + 1], user->ziwork[i + 1]));
5349566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i + 1]));
535c4762a1bSJed Brown user->block_index = i;
5369566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(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 }
5419566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt));
5423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
543c4762a1bSJed Brown }
544c4762a1bSJed Brown
StateMatDuplicate(Mat J_shell,MatDuplicateOption opt,Mat * new_shell)545d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
546d71ae5a4SJacob Faibussowitsch {
547c4762a1bSJed Brown AppCtx *user;
548c4762a1bSJed Brown
549c4762a1bSJed Brown PetscFunctionBegin;
5509566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
551c4762a1bSJed Brown
5529566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, new_shell));
553*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT, (PetscErrorCodeFn *)StateMatMult));
554*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_DUPLICATE, (PetscErrorCodeFn *)StateMatDuplicate));
555*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatMultTranspose));
556*57d50842SBarry Smith PetscCall(MatShellSetOperation(*new_shell, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)StateMatGetDiagonal));
5573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
558c4762a1bSJed Brown }
559c4762a1bSJed Brown
StateMatGetDiagonal(Mat J_shell,Vec X)560d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
561d71ae5a4SJacob Faibussowitsch {
562c4762a1bSJed Brown AppCtx *user;
563c4762a1bSJed Brown
564c4762a1bSJed Brown PetscFunctionBegin;
5659566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
5669566063dSJacob Faibussowitsch PetscCall(VecCopy(user->js_diag, X));
5673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
568c4762a1bSJed Brown }
569c4762a1bSJed Brown
FormConstraints(Tao tao,Vec X,Vec C,void * ptr)570d71ae5a4SJacob Faibussowitsch PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
571d71ae5a4SJacob Faibussowitsch {
572c4762a1bSJed Brown /* con = Ay - q, A = [C(u1) 0 0 ... 0;
573c4762a1bSJed Brown -M C(u2) 0 ... 0;
574c4762a1bSJed Brown 0 -M C(u3) ... 0;
575c4762a1bSJed Brown ... ;
576c4762a1bSJed Brown 0 ... -M C(u_nt)]
577c4762a1bSJed Brown C(u) = eye + ht*Div*[diag(u1); diag(u2)] */
578c4762a1bSJed Brown PetscInt i;
579c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
580c4762a1bSJed Brown
581c4762a1bSJed Brown PetscFunctionBegin;
5829566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
5839566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt));
5849566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
585c4762a1bSJed Brown
586c4762a1bSJed Brown user->block_index = 0;
5879566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0]));
588c4762a1bSJed Brown
589c4762a1bSJed Brown for (i = 1; i < user->nt; i++) {
590c4762a1bSJed Brown user->block_index = i;
5919566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i]));
5929566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1]));
5939566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1]));
594c4762a1bSJed Brown }
595c4762a1bSJed Brown
5969566063dSJacob Faibussowitsch PetscCall(Gather_yi(C, user->yiwork, user->yi_scatter, user->nt));
5979566063dSJacob Faibussowitsch PetscCall(VecAXPY(C, -1.0, user->q));
5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
599c4762a1bSJed Brown }
600c4762a1bSJed Brown
Scatter(Vec x,Vec state,VecScatter s_scat,Vec design,VecScatter d_scat)601d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
602d71ae5a4SJacob Faibussowitsch {
603c4762a1bSJed Brown PetscFunctionBegin;
6049566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6059566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD));
6069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
6079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD));
6083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
609c4762a1bSJed Brown }
610c4762a1bSJed Brown
Scatter_uxi_uyi(Vec u,Vec * uxi,VecScatter * scatx,Vec * uyi,VecScatter * scaty,PetscInt nt)611d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
612d71ae5a4SJacob Faibussowitsch {
613c4762a1bSJed Brown PetscInt i;
614c4762a1bSJed Brown
615c4762a1bSJed Brown PetscFunctionBegin;
616c4762a1bSJed Brown for (i = 0; i < nt; i++) {
6179566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD));
6189566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD));
6199566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD));
6209566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD));
621c4762a1bSJed Brown }
6223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
623c4762a1bSJed Brown }
624c4762a1bSJed Brown
Gather(Vec x,Vec state,VecScatter s_scat,Vec design,VecScatter d_scat)625d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
626d71ae5a4SJacob Faibussowitsch {
627c4762a1bSJed Brown PetscFunctionBegin;
6289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE));
6309566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
6319566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE));
6323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
633c4762a1bSJed Brown }
634c4762a1bSJed Brown
Gather_uxi_uyi(Vec u,Vec * uxi,VecScatter * scatx,Vec * uyi,VecScatter * scaty,PetscInt nt)635d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
636d71ae5a4SJacob Faibussowitsch {
637c4762a1bSJed Brown PetscInt i;
638c4762a1bSJed Brown
639c4762a1bSJed Brown PetscFunctionBegin;
640c4762a1bSJed Brown for (i = 0; i < nt; i++) {
6419566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE));
6429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE));
6439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE));
6449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE));
645c4762a1bSJed Brown }
6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
647c4762a1bSJed Brown }
648c4762a1bSJed Brown
Scatter_yi(Vec y,Vec * yi,VecScatter * scat,PetscInt nt)649d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
650d71ae5a4SJacob Faibussowitsch {
651c4762a1bSJed Brown PetscInt i;
652c4762a1bSJed Brown
653c4762a1bSJed Brown PetscFunctionBegin;
654c4762a1bSJed Brown for (i = 0; i < nt; i++) {
6559566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
6569566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD));
657c4762a1bSJed Brown }
6583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
659c4762a1bSJed Brown }
660c4762a1bSJed Brown
Gather_yi(Vec y,Vec * yi,VecScatter * scat,PetscInt nt)661d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
662d71ae5a4SJacob Faibussowitsch {
663c4762a1bSJed Brown PetscInt i;
664c4762a1bSJed Brown
665c4762a1bSJed Brown PetscFunctionBegin;
666c4762a1bSJed Brown for (i = 0; i < nt; i++) {
6679566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
6689566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE));
669c4762a1bSJed Brown }
6703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
671c4762a1bSJed Brown }
672c4762a1bSJed Brown
HyperbolicInitialize(AppCtx * user)673d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicInitialize(AppCtx *user)
674d71ae5a4SJacob Faibussowitsch {
675c4762a1bSJed Brown PetscInt n, i, j, linear_index, istart, iend, iblock, lo, hi;
676c4762a1bSJed Brown Vec XX, YY, XXwork, YYwork, yi, uxi, ui, bc;
677c4762a1bSJed Brown PetscReal h, sum;
678c4762a1bSJed Brown PetscScalar hinv, neg_hinv, quarter = 0.25, one = 1.0, half_hinv, neg_half_hinv;
679c4762a1bSJed Brown PetscScalar vx, vy, zero = 0.0;
680c4762a1bSJed Brown IS is_from_y, is_to_yi, is_from_u, is_to_uxi, is_to_uyi;
681c4762a1bSJed Brown
682c4762a1bSJed Brown PetscFunctionBegin;
683c4762a1bSJed Brown user->jformed = PETSC_FALSE;
684c4762a1bSJed Brown user->c_formed = PETSC_FALSE;
685c4762a1bSJed Brown
686c4762a1bSJed Brown user->ksp_its = 0;
687c4762a1bSJed Brown user->ksp_its_initial = 0;
688c4762a1bSJed Brown
689c4762a1bSJed Brown n = user->mx * user->mx;
690c4762a1bSJed Brown
691c4762a1bSJed Brown h = 1.0 / user->mx;
692c4762a1bSJed Brown hinv = user->mx;
693c4762a1bSJed Brown neg_hinv = -hinv;
694c4762a1bSJed Brown half_hinv = hinv / 2.0;
695c4762a1bSJed Brown neg_half_hinv = neg_hinv / 2.0;
696c4762a1bSJed Brown
697c4762a1bSJed Brown /* Generate Grad matrix */
6989566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad));
6999566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, PETSC_DECIDE, 2 * n, n));
7009566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad));
7019566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad, 3, NULL, 3, NULL));
7029566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad, 3, NULL));
7039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));
704c4762a1bSJed Brown
705c4762a1bSJed Brown for (i = istart; i < iend; i++) {
706c4762a1bSJed Brown if (i < n) {
707c4762a1bSJed Brown iblock = i / user->mx;
708c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx);
7099566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
710c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx);
7119566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
712c4762a1bSJed Brown }
713c4762a1bSJed Brown if (i >= n) {
714c4762a1bSJed Brown j = (i - user->mx) % n;
7159566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
716c4762a1bSJed Brown j = (j + 2 * user->mx) % n;
7179566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
718c4762a1bSJed Brown }
719c4762a1bSJed Brown }
720c4762a1bSJed Brown
7219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY));
7229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY));
723c4762a1bSJed Brown
7249566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[0]));
7259566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Gradxy[0], PETSC_DECIDE, PETSC_DECIDE, n, n));
7269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Gradxy[0]));
7279566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[0], 3, NULL, 3, NULL));
7289566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[0], 3, NULL));
7299566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Gradxy[0], &istart, &iend));
730c4762a1bSJed Brown
731c4762a1bSJed Brown for (i = istart; i < iend; i++) {
732c4762a1bSJed Brown iblock = i / user->mx;
733c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx);
7349566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
735c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx);
7369566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
7379566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &i, &zero, INSERT_VALUES));
738c4762a1bSJed Brown }
7399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Gradxy[0], MAT_FINAL_ASSEMBLY));
7409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Gradxy[0], MAT_FINAL_ASSEMBLY));
741c4762a1bSJed Brown
7429566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[1]));
7439566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Gradxy[1], PETSC_DECIDE, PETSC_DECIDE, n, n));
7449566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Gradxy[1]));
7459566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[1], 3, NULL, 3, NULL));
7469566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[1], 3, NULL));
7479566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Gradxy[1], &istart, &iend));
748c4762a1bSJed Brown
749c4762a1bSJed Brown for (i = istart; i < iend; i++) {
750c4762a1bSJed Brown j = (i + n - user->mx) % n;
7519566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &half_hinv, INSERT_VALUES));
752c4762a1bSJed Brown j = (j + 2 * user->mx) % n;
7539566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES));
7549566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &i, &zero, INSERT_VALUES));
755c4762a1bSJed Brown }
7569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Gradxy[1], MAT_FINAL_ASSEMBLY));
7579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Gradxy[1], MAT_FINAL_ASSEMBLY));
758c4762a1bSJed Brown
759c4762a1bSJed Brown /* Generate Div matrix */
7609566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Grad, MAT_INITIAL_MATRIX, &user->Div));
7619566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Gradxy[0], MAT_INITIAL_MATRIX, &user->Divxy[0]));
7629566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Gradxy[1], MAT_INITIAL_MATRIX, &user->Divxy[1]));
763c4762a1bSJed Brown
764c4762a1bSJed Brown /* Off-diagonal averaging matrix */
7659566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->M));
7669566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->M, PETSC_DECIDE, PETSC_DECIDE, n, n));
7679566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->M));
7689566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->M, 4, NULL, 4, NULL));
7699566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->M, 4, NULL));
7709566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->M, &istart, &iend));
771c4762a1bSJed Brown
772c4762a1bSJed Brown for (i = istart; i < iend; i++) {
773c4762a1bSJed Brown /* kron(Id,Av) */
774c4762a1bSJed Brown iblock = i / user->mx;
775c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx);
7769566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
777c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx);
7789566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
779c4762a1bSJed Brown
780c4762a1bSJed Brown /* kron(Av,Id) */
781c4762a1bSJed Brown j = (i + user->mx) % n;
7829566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
783c4762a1bSJed Brown j = (i + n - user->mx) % n;
7849566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES));
785c4762a1bSJed Brown }
7869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->M, MAT_FINAL_ASSEMBLY));
7879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->M, MAT_FINAL_ASSEMBLY));
788c4762a1bSJed Brown
789c4762a1bSJed Brown /* Generate 2D grid */
7909566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &XX));
7919566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q));
7929566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX, PETSC_DECIDE, n));
7939566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q, PETSC_DECIDE, n * user->nt));
7949566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX));
7959566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q));
796c4762a1bSJed Brown
7979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YY));
7989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &XXwork));
7999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YYwork));
8009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->d));
8019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->dwork));
802c4762a1bSJed Brown
8039566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX, &istart, &iend));
804c4762a1bSJed Brown for (linear_index = istart; linear_index < iend; linear_index++) {
805c4762a1bSJed Brown i = linear_index % user->mx;
806c4762a1bSJed Brown j = (linear_index - i) / user->mx;
807c4762a1bSJed Brown vx = h * (i + 0.5);
808c4762a1bSJed Brown vy = h * (j + 0.5);
8099566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES));
8109566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES));
811c4762a1bSJed Brown }
812c4762a1bSJed Brown
8139566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX));
8149566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX));
8159566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY));
8169566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY));
817c4762a1bSJed Brown
818c4762a1bSJed Brown /* Compute final density function yT
819c4762a1bSJed Brown yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2))
820c4762a1bSJed Brown yT = yT / (h^2*sum(yT)) */
8219566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork));
8229566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork));
823c4762a1bSJed Brown
8249566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.25));
8259566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.25));
826c4762a1bSJed Brown
8279566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
8289566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
829c4762a1bSJed Brown
8309566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, user->dwork));
8319566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, 1.0, YYwork));
8329566063dSJacob Faibussowitsch PetscCall(VecScale(user->dwork, -30.0));
8339566063dSJacob Faibussowitsch PetscCall(VecExp(user->dwork));
8349566063dSJacob Faibussowitsch PetscCall(VecCopy(user->dwork, user->d));
835c4762a1bSJed Brown
8369566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork));
8379566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork));
838c4762a1bSJed Brown
8399566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.75));
8409566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.75));
841c4762a1bSJed Brown
8429566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
8439566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
844c4762a1bSJed Brown
8459566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, user->dwork));
8469566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, 1.0, YYwork));
8479566063dSJacob Faibussowitsch PetscCall(VecScale(user->dwork, -30.0));
8489566063dSJacob Faibussowitsch PetscCall(VecExp(user->dwork));
849c4762a1bSJed Brown
8509566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->d, 1.0, user->dwork));
8519566063dSJacob Faibussowitsch PetscCall(VecShift(user->d, 1.0));
8529566063dSJacob Faibussowitsch PetscCall(VecSum(user->d, &sum));
8539566063dSJacob Faibussowitsch PetscCall(VecScale(user->d, 1.0 / (h * h * sum)));
854c4762a1bSJed Brown
855c4762a1bSJed Brown /* Initial conditions of forward problem */
8569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &bc));
8579566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork));
8589566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork));
859c4762a1bSJed Brown
8609566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.5));
8619566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.5));
862c4762a1bSJed Brown
8639566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
8649566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
865c4762a1bSJed Brown
8669566063dSJacob Faibussowitsch PetscCall(VecWAXPY(bc, 1.0, XXwork, YYwork));
8679566063dSJacob Faibussowitsch PetscCall(VecScale(bc, -50.0));
8689566063dSJacob Faibussowitsch PetscCall(VecExp(bc));
8699566063dSJacob Faibussowitsch PetscCall(VecShift(bc, 1.0));
8709566063dSJacob Faibussowitsch PetscCall(VecSum(bc, &sum));
8719566063dSJacob Faibussowitsch PetscCall(VecScale(bc, 1.0 / (h * h * sum)));
872c4762a1bSJed Brown
873c4762a1bSJed Brown /* Create scatter from y to y_1,y_2,...,y_nt */
874c4762a1bSJed Brown /* TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */
8759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->yi_scatter));
8769566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &yi));
8779566063dSJacob Faibussowitsch PetscCall(VecSetSizes(yi, PETSC_DECIDE, user->mx * user->mx));
8789566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(yi));
8799566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yi));
8809566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yiwork));
8819566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->ziwork));
882c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
8839566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->yi[i], &lo, &hi));
8849566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_yi));
8859566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + i * user->mx * user->mx, 1, &is_from_y));
8869566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y, is_from_y, user->yi[i], is_to_yi, &user->yi_scatter[i]));
8879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_yi));
8889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y));
889c4762a1bSJed Brown }
890c4762a1bSJed Brown
891c4762a1bSJed Brown /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */
892c4762a1bSJed Brown /* TODO: reorder for better parallelism */
8939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uxi_scatter));
8949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uyi_scatter));
8959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->ux_scatter));
8969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uy_scatter));
8979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * user->nt * user->mx * user->mx, &user->ui_scatter));
8989566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &uxi));
8999566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &ui));
9009566063dSJacob Faibussowitsch PetscCall(VecSetSizes(uxi, PETSC_DECIDE, user->mx * user->mx));
9019566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ui, PETSC_DECIDE, 2 * user->mx * user->mx));
9029566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(uxi));
9039566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ui));
9049566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxi));
9059566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyi));
9069566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxiwork));
9079566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyiwork));
9089566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ui, user->nt, &user->ui));
9099566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ui, user->nt, &user->uiwork));
910c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
9119566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi));
9129566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi));
9139566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u));
9149566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->uxi[i], is_to_uxi, &user->uxi_scatter[i]));
915c4762a1bSJed Brown
9169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi));
9179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u));
918c4762a1bSJed Brown
9199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi));
9209566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi));
9219566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + (2 * i + 1) * user->mx * user->mx, 1, &is_from_u));
9229566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->uyi[i], is_to_uyi, &user->uyi_scatter[i]));
923c4762a1bSJed Brown
9249566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uyi));
9259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u));
926c4762a1bSJed Brown
9279566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi));
9289566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi));
9299566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_from_u));
9309566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uxi[i], is_to_uxi, &user->ux_scatter[i]));
931c4762a1bSJed Brown
9329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi));
9339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u));
934c4762a1bSJed Brown
9359566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi));
9369566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi));
9379566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + user->mx * user->mx, 1, &is_from_u));
9389566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uyi[i], is_to_uyi, &user->uy_scatter[i]));
939c4762a1bSJed Brown
9409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uyi));
9419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u));
942c4762a1bSJed Brown
9439566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->ui[i], &lo, &hi));
9449566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi));
9459566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u));
9469566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->ui[i], is_to_uxi, &user->ui_scatter[i]));
947c4762a1bSJed Brown
9489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi));
9499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u));
950c4762a1bSJed Brown }
951c4762a1bSJed Brown
952c4762a1bSJed Brown /* RHS of forward problem */
9539566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, bc, user->yiwork[0]));
95448a46eb9SPierre Jolivet for (i = 1; i < user->nt; i++) PetscCall(VecSet(user->yiwork[i], 0.0));
9559566063dSJacob Faibussowitsch PetscCall(Gather_yi(user->q, user->yiwork, user->yi_scatter, user->nt));
956c4762a1bSJed Brown
957c4762a1bSJed Brown /* Compute true velocity field utrue */
9589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->utrue));
959c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
9609566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, user->uxi[i]));
9619566063dSJacob Faibussowitsch PetscCall(VecScale(user->uxi[i], 150.0 * i * user->ht));
9629566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, user->uyi[i]));
9639566063dSJacob Faibussowitsch PetscCall(VecShift(user->uyi[i], -10.0));
9649566063dSJacob Faibussowitsch PetscCall(VecScale(user->uyi[i], 15.0 * i * user->ht));
965c4762a1bSJed Brown }
9669566063dSJacob Faibussowitsch PetscCall(Gather_uxi_uyi(user->utrue, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
967c4762a1bSJed Brown
968c4762a1bSJed Brown /* Initial guess and reference model */
9699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue, &user->ur));
970c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
9719566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, user->uxi[i]));
9729566063dSJacob Faibussowitsch PetscCall(VecShift(user->uxi[i], i * user->ht));
9739566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, user->uyi[i]));
9749566063dSJacob Faibussowitsch PetscCall(VecShift(user->uyi[i], -i * user->ht));
975c4762a1bSJed Brown }
9769566063dSJacob Faibussowitsch PetscCall(Gather_uxi_uyi(user->ur, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
977c4762a1bSJed Brown
978c4762a1bSJed Brown /* Generate regularization matrix L */
9799566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->LT));
9809566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->LT, PETSC_DECIDE, PETSC_DECIDE, 2 * n * user->nt, n * user->nt));
9819566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->LT));
9829566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->LT, 1, NULL, 1, NULL));
9839566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->LT, 1, NULL));
9849566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->LT, &istart, &iend));
985c4762a1bSJed Brown
986c4762a1bSJed Brown for (i = istart; i < iend; i++) {
987c4762a1bSJed Brown iblock = (i + n) / (2 * n);
988c4762a1bSJed Brown j = i - iblock * n;
9899566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->LT, 1, &i, 1, &j, &user->gamma, INSERT_VALUES));
990c4762a1bSJed Brown }
991c4762a1bSJed Brown
9929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->LT, MAT_FINAL_ASSEMBLY));
9939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->LT, MAT_FINAL_ASSEMBLY));
994c4762a1bSJed Brown
9959566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->LT, MAT_INITIAL_MATRIX, &user->L));
996c4762a1bSJed Brown
997c4762a1bSJed Brown /* Build work vectors and matrices */
9989566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork));
9999566063dSJacob Faibussowitsch PetscCall(VecSetType(user->lwork, VECMPI));
10009566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, user->m));
10019566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork));
1002c4762a1bSJed Brown
10039566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork));
1004c4762a1bSJed Brown
10059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ywork));
10069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->uwork));
10079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->vwork));
10089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->js_diag));
10099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c, &user->cwork));
1010c4762a1bSJed Brown
1011c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing A*x */
10129566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->Js));
1013*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (PetscErrorCodeFn *)StateMatMult));
1014*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_DUPLICATE, (PetscErrorCodeFn *)StateMatDuplicate));
1015*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatMultTranspose));
1016*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)StateMatGetDiagonal));
1017c4762a1bSJed Brown
1018c4762a1bSJed Brown /* Diagonal blocks of user->Js */
10199566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlock));
1020*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (PetscErrorCodeFn *)StateMatBlockMult));
1021*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatBlockMultTranspose));
1022c4762a1bSJed Brown
1023c4762a1bSJed Brown /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1024c4762a1bSJed Brown D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1025c4762a1bSJed Brown This is an SOR preconditioner for user->JsBlock. */
10269566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlockPrec));
1027*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT, (PetscErrorCodeFn *)StateMatBlockPrecMult));
1028*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatBlockPrecMultTranspose));
1029c4762a1bSJed Brown
1030c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */
10319566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->n - user->m, user, &user->Jd));
1032*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (PetscErrorCodeFn *)DesignMatMult));
1033*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)DesignMatMultTranspose));
1034c4762a1bSJed Brown
1035c4762a1bSJed Brown /* User-defined routines for computing user->Js\x and user->Js^T\x*/
10369566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsInv));
1037*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (PetscErrorCodeFn *)StateMatInvMult));
1038*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatInvTransposeMult));
1039c4762a1bSJed Brown
1040c4762a1bSJed Brown /* Build matrices for SOR preconditioner */
10419566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt));
10429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * n, &user->C));
10439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * n, &user->Cwork));
1044c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
10459566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Divxy[0], MAT_COPY_VALUES, &user->C[i]));
10469566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Divxy[1], MAT_COPY_VALUES, &user->Cwork[i]));
1047c4762a1bSJed Brown
10489566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i]));
10499566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i]));
10509566063dSJacob Faibussowitsch PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], DIFFERENT_NONZERO_PATTERN));
10519566063dSJacob Faibussowitsch PetscCall(MatScale(user->C[i], user->ht));
10529566063dSJacob Faibussowitsch PetscCall(MatShift(user->C[i], 1.0));
1053c4762a1bSJed Brown }
1054c4762a1bSJed Brown
1055c4762a1bSJed Brown /* Solver options and tolerances */
10569566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver));
10579566063dSJacob Faibussowitsch PetscCall(KSPSetType(user->solver, KSPGMRES));
10589566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->JsBlockPrec));
10599566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, 1e3, 500));
10609566063dSJacob Faibussowitsch /* PetscCall(KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500)); */
10619566063dSJacob Faibussowitsch PetscCall(KSPGetPC(user->solver, &user->prec));
10629566063dSJacob Faibussowitsch PetscCall(PCSetType(user->prec, PCSHELL));
1063c4762a1bSJed Brown
10649566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(user->prec, StateMatBlockPrecMult));
10659566063dSJacob Faibussowitsch PetscCall(PCShellSetApplyTranspose(user->prec, StateMatBlockPrecMultTranspose));
10669566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(user->prec, user));
1067c4762a1bSJed Brown
1068c4762a1bSJed Brown /* Compute true state function yt given ut */
10699566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ytrue));
10709566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->ytrue, PETSC_DECIDE, n * user->nt));
10719566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->ytrue));
1072c4762a1bSJed Brown user->c_formed = PETSC_TRUE;
10739566063dSJacob Faibussowitsch PetscCall(VecCopy(user->utrue, user->u)); /* Set u=utrue temporarily for StateMatInv */
10749566063dSJacob Faibussowitsch PetscCall(VecSet(user->ytrue, 0.0)); /* Initial guess */
10759566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->ytrue));
10769566063dSJacob Faibussowitsch PetscCall(VecCopy(user->ur, user->u)); /* Reset u=ur */
1077c4762a1bSJed Brown
1078c4762a1bSJed Brown /* Initial guess y0 for state given u0 */
10799566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->y));
1080c4762a1bSJed Brown
1081c4762a1bSJed Brown /* Data discretization */
10829566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Q));
10839566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Q, PETSC_DECIDE, PETSC_DECIDE, user->mx * user->mx, user->m));
10849566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Q));
10859566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Q, 0, NULL, 1, NULL));
10869566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Q, 1, NULL));
1087c4762a1bSJed Brown
10889566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Q, &istart, &iend));
1089c4762a1bSJed Brown
1090c4762a1bSJed Brown for (i = istart; i < iend; i++) {
1091c4762a1bSJed Brown j = i + user->m - user->mx * user->mx;
10929566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &one, INSERT_VALUES));
1093c4762a1bSJed Brown }
1094c4762a1bSJed Brown
10959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Q, MAT_FINAL_ASSEMBLY));
10969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Q, MAT_FINAL_ASSEMBLY));
1097c4762a1bSJed Brown
10989566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Q, MAT_INITIAL_MATRIX, &user->QT));
1099c4762a1bSJed Brown
11009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX));
11019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY));
11029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork));
11039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork));
11049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yi));
11059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uxi));
11069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ui));
11079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bc));
1108c4762a1bSJed Brown
1109c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11109566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver));
11113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1112c4762a1bSJed Brown }
1113c4762a1bSJed Brown
HyperbolicDestroy(AppCtx * user)1114d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicDestroy(AppCtx *user)
1115d71ae5a4SJacob Faibussowitsch {
1116c4762a1bSJed Brown PetscInt i;
1117c4762a1bSJed Brown
1118c4762a1bSJed Brown PetscFunctionBegin;
11199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Q));
11209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->QT));
11219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div));
11229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork));
11239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad));
11249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L));
11259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->LT));
11269566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver));
11279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js));
11289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd));
11299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlockPrec));
11309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv));
11319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock));
11329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divxy[0]));
11339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divxy[1]));
11349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Gradxy[0]));
11359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Gradxy[1]));
11369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->M));
1137c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
11389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->C[i]));
11399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Cwork[i]));
1140c4762a1bSJed Brown }
11419566063dSJacob Faibussowitsch PetscCall(PetscFree(user->C));
11429566063dSJacob Faibussowitsch PetscCall(PetscFree(user->Cwork));
11439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u));
11449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork));
11459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->vwork));
11469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue));
11479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y));
11489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork));
11499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue));
11509566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yi));
11519566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yiwork));
11529566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->ziwork));
11539566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uxi));
11549566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uyi));
11559566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uxiwork));
11569566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uyiwork));
11579566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->ui));
11589566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uiwork));
11599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c));
11609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork));
11619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur));
11629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q));
11639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d));
11649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork));
11659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork));
11669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag));
11679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is));
11689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is));
11699566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter));
11709566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter));
1171c4762a1bSJed Brown for (i = 0; i < user->nt; i++) {
11729566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uxi_scatter[i]));
11739566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uyi_scatter[i]));
11749566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->ux_scatter[i]));
11759566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uy_scatter[i]));
11769566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->ui_scatter[i]));
11779566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
1178c4762a1bSJed Brown }
11799566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uxi_scatter));
11809566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uyi_scatter));
11819566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ux_scatter));
11829566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uy_scatter));
11839566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ui_scatter));
11849566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter));
11853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1186c4762a1bSJed Brown }
1187c4762a1bSJed Brown
HyperbolicMonitor(Tao tao,void * ptr)1188d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1189d71ae5a4SJacob Faibussowitsch {
1190c4762a1bSJed Brown Vec X;
1191c4762a1bSJed Brown PetscReal unorm, ynorm;
1192c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
1193c4762a1bSJed Brown
1194c4762a1bSJed Brown PetscFunctionBegin;
11959566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X));
11969566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
11979566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue));
11989566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue));
11999566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork, NORM_2, &unorm));
12009566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork, NORM_2, &ynorm));
12019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm));
12023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1203c4762a1bSJed Brown }
1204c4762a1bSJed Brown
1205c4762a1bSJed Brown /*TEST
1206c4762a1bSJed Brown
1207c4762a1bSJed Brown build:
1208c4762a1bSJed Brown requires: !complex
1209c4762a1bSJed Brown
1210c4762a1bSJed Brown test:
1211c4762a1bSJed Brown requires: !single
121210978b7dSBarry Smith args: -tao_monitor_constraint_norm -tao_max_funcs 10 -tao_type lcl -tao_gatol 1.e-5
1213c4762a1bSJed Brown
1214c4762a1bSJed Brown test:
1215c4762a1bSJed Brown suffix: guess_pod
1216c4762a1bSJed Brown requires: !single
121710978b7dSBarry Smith args: -tao_monitor_constraint_norm -tao_max_funcs 10 -tao_type lcl -ksp_guess_type pod -tao_gatol 1.e-5
1218c4762a1bSJed Brown
1219c4762a1bSJed Brown TEST*/
1220