1c4762a1bSJed Brown #include <petsc/private/taoimpl.h>
2c4762a1bSJed Brown
3c4762a1bSJed Brown typedef struct {
4c4762a1bSJed Brown PetscInt n; /* Number of total variables */
5c4762a1bSJed Brown PetscInt m; /* Number of constraints */
6c4762a1bSJed Brown PetscInt nstate;
7c4762a1bSJed Brown PetscInt ndesign;
8c4762a1bSJed Brown PetscInt mx; /* grid points in each direction */
9c4762a1bSJed Brown PetscInt ns; /* Number of data samples (1<=ns<=8)
10c4762a1bSJed Brown Currently only ns=1 is supported */
11c4762a1bSJed Brown PetscInt ndata; /* Number of data points per sample */
12c4762a1bSJed Brown IS s_is;
13c4762a1bSJed Brown IS d_is;
14c4762a1bSJed Brown
15c4762a1bSJed Brown VecScatter state_scatter;
16c4762a1bSJed Brown VecScatter design_scatter;
17c4762a1bSJed Brown VecScatter *yi_scatter, *di_scatter;
18c4762a1bSJed Brown Vec suby, subq, subd;
19c4762a1bSJed Brown Mat Js, Jd, JsPrec, JsInv, JsBlock;
20c4762a1bSJed Brown
21c4762a1bSJed Brown PetscReal alpha; /* Regularization parameter */
22c4762a1bSJed Brown PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */
23c4762a1bSJed Brown PetscReal noise; /* Amount of noise to add to data */
24c4762a1bSJed Brown PetscReal *ones;
25c4762a1bSJed Brown Mat Q;
26c4762a1bSJed Brown Mat MQ;
27c4762a1bSJed Brown Mat L;
28c4762a1bSJed Brown
29c4762a1bSJed Brown Mat Grad;
30c4762a1bSJed Brown Mat Av, Avwork;
31c4762a1bSJed Brown Mat Div, Divwork;
32c4762a1bSJed Brown Mat DSG;
33c4762a1bSJed Brown Mat Diag, Ones;
34c4762a1bSJed Brown
35c4762a1bSJed Brown Vec q;
36c4762a1bSJed Brown Vec ur; /* reference */
37c4762a1bSJed Brown
38c4762a1bSJed Brown Vec d;
39c4762a1bSJed Brown Vec dwork;
40c4762a1bSJed Brown
41c4762a1bSJed Brown Vec x; /* super vec of y,u */
42c4762a1bSJed Brown
43c4762a1bSJed Brown Vec y; /* state variables */
44c4762a1bSJed Brown Vec ywork;
45c4762a1bSJed Brown
46c4762a1bSJed Brown Vec ytrue;
47c4762a1bSJed Brown
48c4762a1bSJed Brown Vec u; /* design variables */
49c4762a1bSJed Brown Vec uwork;
50c4762a1bSJed Brown
51c4762a1bSJed Brown Vec utrue;
52c4762a1bSJed Brown
53c4762a1bSJed Brown Vec js_diag;
54c4762a1bSJed Brown
55c4762a1bSJed Brown Vec c; /* constraint vector */
56c4762a1bSJed Brown Vec cwork;
57c4762a1bSJed Brown
58c4762a1bSJed Brown Vec lwork;
59c4762a1bSJed Brown Vec S;
60c4762a1bSJed Brown Vec Swork, Twork, Sdiag, Ywork;
61c4762a1bSJed Brown Vec Av_u;
62c4762a1bSJed Brown
63c4762a1bSJed Brown KSP solver;
64c4762a1bSJed Brown PC prec;
65c4762a1bSJed Brown
66c4762a1bSJed Brown PetscReal tola, tolb, tolc, told;
67c4762a1bSJed Brown PetscInt ksp_its;
68c4762a1bSJed Brown PetscInt ksp_its_initial;
69c4762a1bSJed Brown PetscLogStage stages[10];
70c4762a1bSJed Brown PetscBool use_ptap;
71c4762a1bSJed Brown PetscBool use_lrc;
72c4762a1bSJed Brown } AppCtx;
73c4762a1bSJed Brown
74c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
75c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
76c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
77c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void *);
78c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void *);
79c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
80c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
81c4762a1bSJed Brown PetscErrorCode Gather(Vec, Vec, VecScatter, Vec, VecScatter);
82c4762a1bSJed Brown PetscErrorCode Scatter(Vec, Vec, VecScatter, Vec, VecScatter);
83c4762a1bSJed Brown PetscErrorCode EllipticInitialize(AppCtx *);
84c4762a1bSJed Brown PetscErrorCode EllipticDestroy(AppCtx *);
85c4762a1bSJed Brown PetscErrorCode EllipticMonitor(Tao, void *);
86c4762a1bSJed Brown
87c4762a1bSJed Brown PetscErrorCode StateBlockMatMult(Mat, Vec, Vec);
88c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat, Vec, Vec);
89c4762a1bSJed Brown
90c4762a1bSJed Brown PetscErrorCode StateInvMatMult(Mat, Vec, Vec);
91c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat, Vec, Vec);
92c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat, Vec, Vec);
93c4762a1bSJed Brown
94c4762a1bSJed Brown PetscErrorCode QMatMult(Mat, Vec, Vec);
95c4762a1bSJed Brown PetscErrorCode QMatMultTranspose(Mat, Vec, Vec);
96c4762a1bSJed Brown
97c4762a1bSJed Brown static char help[] = "";
98c4762a1bSJed Brown
main(int argc,char ** argv)99d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown Vec x0;
102c4762a1bSJed Brown Tao tao;
103c4762a1bSJed Brown AppCtx user;
104c4762a1bSJed Brown PetscInt ntests = 1;
105c4762a1bSJed Brown PetscInt i;
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, "elliptic example", NULL);
1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL));
112c4762a1bSJed Brown user.ns = 6;
1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ns", "Number of data samples (1<=ns<=8)", "", user.ns, &user.ns, 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.alpha = 0.1;
1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL));
118c4762a1bSJed Brown user.beta = 0.00001;
1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-beta", "Weight attributed to ||u||^2 in regularization functional", "", user.beta, &user.beta, NULL));
120c4762a1bSJed Brown user.noise = 0.01;
1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-noise", "Amount of noise to add to data", "", user.noise, &user.noise, NULL));
122c4762a1bSJed Brown
123c4762a1bSJed Brown user.use_ptap = PETSC_FALSE;
1249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_ptap", "Use ptap matrix for DSG", "", user.use_ptap, &user.use_ptap, NULL));
125c4762a1bSJed Brown user.use_lrc = PETSC_FALSE;
1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_lrc", "Use lrc matrix for Js", "", user.use_lrc, &user.use_lrc, NULL));
1279566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL));
128d0609cedSBarry Smith PetscOptionsEnd();
12976280437SVaclav Hapla
130c4762a1bSJed Brown user.m = user.ns * user.mx * user.mx * user.mx; /* number of constraints */
131c4762a1bSJed Brown user.nstate = user.m;
132c4762a1bSJed Brown user.ndesign = user.mx * user.mx * user.mx;
133c4762a1bSJed Brown user.n = user.nstate + user.ndesign; /* number of variables */
134c4762a1bSJed Brown
135c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
1369566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1379566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLCL));
138c4762a1bSJed Brown
139c4762a1bSJed Brown /* Set up initial vectors and matrices */
1409566063dSJacob Faibussowitsch PetscCall(EllipticInitialize(&user));
141c4762a1bSJed Brown
1429566063dSJacob Faibussowitsch PetscCall(Gather(user.x, user.y, user.state_scatter, user.u, user.design_scatter));
1439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.x, &x0));
1449566063dSJacob Faibussowitsch PetscCall(VecCopy(user.x, x0));
145c4762a1bSJed Brown
146c4762a1bSJed Brown /* Set solution vector with an initial guess */
1479566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, user.x));
1489566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user));
1499566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
1509566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
151c4762a1bSJed Brown
1529566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, NULL, user.JsInv, FormJacobianState, &user));
1539566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
154c4762a1bSJed Brown
1559566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is));
1569566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
157c4762a1bSJed Brown
158c4762a1bSJed Brown /* SOLVE THE APPLICATION */
1599566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials", &user.stages[1]));
1609566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user.stages[1]));
161c4762a1bSJed Brown for (i = 0; i < ntests; i++) {
1629566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
16363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its));
1649566063dSJacob Faibussowitsch PetscCall(VecCopy(x0, user.x));
165c4762a1bSJed Brown }
1669566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop());
1679566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)user.x));
1689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: "));
16963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial));
170c4762a1bSJed Brown
1719566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
1729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0));
1739566063dSJacob Faibussowitsch PetscCall(EllipticDestroy(&user));
1749566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
175b122ec5aSJacob Faibussowitsch return 0;
176c4762a1bSJed Brown }
177c4762a1bSJed Brown /* ------------------------------------------------------------------- */
178c4762a1bSJed Brown /*
179c4762a1bSJed Brown dwork = Qy - d
180c4762a1bSJed Brown lwork = L*(u-ur)
181c4762a1bSJed Brown f = 1/2 * (dwork.dwork + alpha*lwork.lwork)
182c4762a1bSJed Brown */
FormFunction(Tao tao,Vec X,PetscReal * f,void * ptr)183d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr)
184d71ae5a4SJacob Faibussowitsch {
185c4762a1bSJed Brown PetscReal d1 = 0, d2 = 0;
186c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
187c4762a1bSJed Brown
188c4762a1bSJed Brown PetscFunctionBegin;
1899566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
1909566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->y, user->dwork));
1919566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
1929566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1));
1939566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
1949566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork));
1959566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork, user->lwork, &d2));
196c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2);
1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
198c4762a1bSJed Brown }
199c4762a1bSJed Brown
200c4762a1bSJed Brown /* ------------------------------------------------------------------- */
201c4762a1bSJed Brown /*
202c4762a1bSJed Brown state: g_s = Q' *(Qy - d)
203c4762a1bSJed Brown design: g_d = alpha*L'*L*(u-ur)
204c4762a1bSJed Brown */
FormGradient(Tao tao,Vec X,Vec G,void * ptr)205d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
208c4762a1bSJed Brown
209c4762a1bSJed Brown PetscFunctionBegin;
2109566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2119566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->y, user->dwork));
2129566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2139566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->MQ, user->dwork, user->ywork));
2149566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
2159566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork));
2169566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->L, user->lwork, user->uwork));
2179566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha));
2189566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
220c4762a1bSJed Brown }
221c4762a1bSJed Brown
FormFunctionGradient(Tao tao,Vec X,PetscReal * f,Vec G,void * ptr)222d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
223d71ae5a4SJacob Faibussowitsch {
224c4762a1bSJed Brown PetscReal d1, d2;
225c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
226c4762a1bSJed Brown
227c4762a1bSJed Brown PetscFunctionBegin;
2289566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2299566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->y, user->dwork));
2309566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d));
2319566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1));
2329566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->MQ, user->dwork, user->ywork));
233c4762a1bSJed Brown
2349566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
2359566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork));
2369566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork, user->lwork, &d2));
2379566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->L, user->lwork, user->uwork));
2389566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha));
239c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2);
2409566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown
244c4762a1bSJed Brown /* ------------------------------------------------------------------- */
245c4762a1bSJed Brown /* A
246c4762a1bSJed Brown MatShell object
247c4762a1bSJed Brown */
FormJacobianState(Tao tao,Vec X,Mat J,Mat JPre,Mat JInv,void * ptr)248d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
249d71ae5a4SJacob Faibussowitsch {
250c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
251c4762a1bSJed Brown
252c4762a1bSJed Brown PetscFunctionBegin;
2539566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
254c4762a1bSJed Brown /* DSG = Div * (1/Av_u) * Grad */
2559566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
2569566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u));
2579566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
2589566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
2599566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u, user->Swork));
2609566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork));
261c4762a1bSJed Brown if (user->use_ptap) {
2629566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES));
2639566063dSJacob Faibussowitsch PetscCall(MatPtAP(user->Diag, user->Grad, MAT_REUSE_MATRIX, 1.0, &user->DSG));
264c4762a1bSJed Brown } else {
2659566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
2669566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork));
2679566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG));
268c4762a1bSJed Brown }
2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
270c4762a1bSJed Brown }
271c4762a1bSJed Brown /* ------------------------------------------------------------------- */
272c4762a1bSJed Brown /* B */
FormJacobianDesign(Tao tao,Vec X,Mat J,void * ptr)273d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
274d71ae5a4SJacob Faibussowitsch {
275c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
276c4762a1bSJed Brown
277c4762a1bSJed Brown PetscFunctionBegin;
2789566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown
StateBlockMatMult(Mat J_shell,Vec X,Vec Y)282d71ae5a4SJacob Faibussowitsch PetscErrorCode StateBlockMatMult(Mat J_shell, Vec X, Vec Y)
283d71ae5a4SJacob Faibussowitsch {
284c4762a1bSJed Brown PetscReal sum;
285c4762a1bSJed Brown AppCtx *user;
286c4762a1bSJed Brown
287c4762a1bSJed Brown PetscFunctionBegin;
2889566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
2899566063dSJacob Faibussowitsch PetscCall(MatMult(user->DSG, X, Y));
2909566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum));
291c4762a1bSJed Brown sum /= user->ndesign;
2929566063dSJacob Faibussowitsch PetscCall(VecShift(Y, sum));
2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
294c4762a1bSJed Brown }
295c4762a1bSJed Brown
StateMatMult(Mat J_shell,Vec X,Vec Y)296d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
297d71ae5a4SJacob Faibussowitsch {
298c4762a1bSJed Brown PetscInt i;
299c4762a1bSJed Brown AppCtx *user;
300c4762a1bSJed Brown
301c4762a1bSJed Brown PetscFunctionBegin;
3029566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
303c4762a1bSJed Brown if (user->ns == 1) {
3049566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, X, Y));
305c4762a1bSJed Brown } else {
306c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
3079566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
3089566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0));
3099566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->subq, user->suby));
3109566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0));
311c4762a1bSJed Brown }
312c4762a1bSJed Brown }
3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
314c4762a1bSJed Brown }
315c4762a1bSJed Brown
StateInvMatMult(Mat J_shell,Vec X,Vec Y)316d71ae5a4SJacob Faibussowitsch PetscErrorCode StateInvMatMult(Mat J_shell, Vec X, Vec Y)
317d71ae5a4SJacob Faibussowitsch {
318c4762a1bSJed Brown PetscInt its, i;
319c4762a1bSJed Brown AppCtx *user;
320c4762a1bSJed Brown
321c4762a1bSJed Brown PetscFunctionBegin;
3229566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
3239566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->DSG));
324c4762a1bSJed Brown if (Y == user->ytrue) {
325c4762a1bSJed Brown /* First solve is done using true solution to set up problem */
326606f75f6SBarry Smith PetscCall(KSPSetTolerances(user->solver, 1e-8, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
327c4762a1bSJed Brown } else {
328606f75f6SBarry Smith PetscCall(KSPSetTolerances(user->solver, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
329c4762a1bSJed Brown }
330c4762a1bSJed Brown if (user->ns == 1) {
3319566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, X, Y));
3329566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
333c4762a1bSJed Brown user->ksp_its += its;
334c4762a1bSJed Brown } else {
335c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
3369566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
3379566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0));
3389566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->subq, user->suby));
3399566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its));
340c4762a1bSJed Brown user->ksp_its += its;
3419566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0));
342c4762a1bSJed Brown }
343c4762a1bSJed Brown }
3443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
345c4762a1bSJed Brown }
QMatMult(Mat J_shell,Vec X,Vec Y)346d71ae5a4SJacob Faibussowitsch PetscErrorCode QMatMult(Mat J_shell, Vec X, Vec Y)
347d71ae5a4SJacob Faibussowitsch {
348c4762a1bSJed Brown AppCtx *user;
349c4762a1bSJed Brown PetscInt i;
350c4762a1bSJed Brown
351c4762a1bSJed Brown PetscFunctionBegin;
3529566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
353c4762a1bSJed Brown if (user->ns == 1) {
3549566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, X, Y));
355c4762a1bSJed Brown } else {
356c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
3579566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
3589566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->subd, user->di_scatter[i], 0, 0));
3599566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->subq, user->subd));
3609566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->subd, user->di_scatter[i], 0, 0));
361c4762a1bSJed Brown }
362c4762a1bSJed Brown }
3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
364c4762a1bSJed Brown }
365c4762a1bSJed Brown
QMatMultTranspose(Mat J_shell,Vec X,Vec Y)366d71ae5a4SJacob Faibussowitsch PetscErrorCode QMatMultTranspose(Mat J_shell, Vec X, Vec Y)
367d71ae5a4SJacob Faibussowitsch {
368c4762a1bSJed Brown AppCtx *user;
369c4762a1bSJed Brown PetscInt i;
370c4762a1bSJed Brown
371c4762a1bSJed Brown PetscFunctionBegin;
3729566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
373c4762a1bSJed Brown if (user->ns == 1) {
3749566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Q, X, Y));
375c4762a1bSJed Brown } else {
376c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
3779566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subd, user->di_scatter[i], 0, 0));
3789566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0));
3799566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Q, user->subd, user->suby));
3809566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0));
381c4762a1bSJed Brown }
382c4762a1bSJed Brown }
3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
384c4762a1bSJed Brown }
385c4762a1bSJed Brown
DesignMatMult(Mat J_shell,Vec X,Vec Y)386d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
387d71ae5a4SJacob Faibussowitsch {
388c4762a1bSJed Brown PetscInt i;
389c4762a1bSJed Brown AppCtx *user;
390c4762a1bSJed Brown
391c4762a1bSJed Brown PetscFunctionBegin;
3929566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
393c4762a1bSJed Brown
394c4762a1bSJed Brown /* sdiag(1./v) */
3959566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
3969566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u));
3979566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
398c4762a1bSJed Brown
399c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */
4009566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Swork));
4019566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Swork, user->Swork));
4029566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork));
403c4762a1bSJed Brown
404c4762a1bSJed Brown /* (Av * (sdiag(1./v) * b)) */
4059566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->uwork, X));
4069566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Twork));
407c4762a1bSJed Brown
408c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
4099566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork));
410c4762a1bSJed Brown
411c4762a1bSJed Brown if (user->ns == 1) {
412c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) */
4139566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->y, user->Twork));
414c4762a1bSJed Brown
415c4762a1bSJed Brown /* Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
4169566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork));
4179566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Swork, Y));
418c4762a1bSJed Brown } else {
419c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
4209566063dSJacob Faibussowitsch PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0));
4219566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->subq, user->yi_scatter[i], 0, 0));
422c4762a1bSJed Brown
4239566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->suby, user->Twork));
4249566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork, user->Twork, user->Swork));
4259566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Twork, user->subq));
4269566063dSJacob Faibussowitsch PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0));
4279566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->subq, user->yi_scatter[i], 0, 0));
428c4762a1bSJed Brown }
429c4762a1bSJed Brown }
4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
431c4762a1bSJed Brown }
432c4762a1bSJed Brown
DesignMatMultTranspose(Mat J_shell,Vec X,Vec Y)433d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
434d71ae5a4SJacob Faibussowitsch {
435c4762a1bSJed Brown PetscInt i;
436c4762a1bSJed Brown AppCtx *user;
437c4762a1bSJed Brown
438c4762a1bSJed Brown PetscFunctionBegin;
4399566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user));
4409566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Y));
441c4762a1bSJed Brown
442c4762a1bSJed Brown /* Sdiag = 1./((Av*(1./v)).^2) */
4439566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
4449566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u));
4459566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
4469566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Swork));
4479566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Sdiag, user->Swork, user->Swork));
4489566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Sdiag));
449c4762a1bSJed Brown
450c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
4519566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
4529566063dSJacob Faibussowitsch PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0));
453c4762a1bSJed Brown
454c4762a1bSJed Brown /* Swork = (Div' * b(:,i)) */
4559566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->subq, user->Swork));
456c4762a1bSJed Brown
457c4762a1bSJed Brown /* Twork = Grad*y(:,i) */
4589566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->suby, user->Twork));
459c4762a1bSJed Brown
460c4762a1bSJed Brown /* Twork = sdiag(Twork) * Swork */
4619566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork, user->Swork, user->Twork));
462c4762a1bSJed Brown
463c4762a1bSJed Brown /* Swork = pointwisemult(Sdiag,Twork) */
4649566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Sdiag));
465c4762a1bSJed Brown
466c4762a1bSJed Brown /* Ywork = Av' * Swork */
4679566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Av, user->Swork, user->Ywork));
468c4762a1bSJed Brown
469c4762a1bSJed Brown /* Ywork = pointwisemult(uwork,Ywork) */
4709566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Ywork, user->uwork, user->Ywork));
4719566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, 1.0, user->Ywork));
4729566063dSJacob Faibussowitsch PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0));
473c4762a1bSJed Brown }
4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
475c4762a1bSJed Brown }
476c4762a1bSJed Brown
FormConstraints(Tao tao,Vec X,Vec C,void * ptr)477d71ae5a4SJacob Faibussowitsch PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
478d71ae5a4SJacob Faibussowitsch {
479c4762a1bSJed Brown /* C=Ay - q A = Div * Sigma * Grad + hx*hx*hx*ones(n,n) */
480c4762a1bSJed Brown PetscReal sum;
481c4762a1bSJed Brown PetscInt i;
482c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
483c4762a1bSJed Brown
484c4762a1bSJed Brown PetscFunctionBegin;
4859566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
486c4762a1bSJed Brown if (user->ns == 1) {
4879566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->y, user->Swork));
4889566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u));
4899566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Swork, C));
4909566063dSJacob Faibussowitsch PetscCall(VecSum(user->y, &sum));
491c4762a1bSJed Brown sum /= user->ndesign;
4929566063dSJacob Faibussowitsch PetscCall(VecShift(C, sum));
493c4762a1bSJed Brown } else {
494c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
4959566063dSJacob Faibussowitsch PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0));
4969566063dSJacob Faibussowitsch PetscCall(Scatter(C, user->subq, user->yi_scatter[i], 0, 0));
4979566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->suby, user->Swork));
4989566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u));
4999566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Swork, user->subq));
500c4762a1bSJed Brown
5019566063dSJacob Faibussowitsch PetscCall(VecSum(user->suby, &sum));
502c4762a1bSJed Brown sum /= user->ndesign;
5039566063dSJacob Faibussowitsch PetscCall(VecShift(user->subq, sum));
504c4762a1bSJed Brown
5059566063dSJacob Faibussowitsch PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0));
5069566063dSJacob Faibussowitsch PetscCall(Gather(C, user->subq, user->yi_scatter[i], 0, 0));
507c4762a1bSJed Brown }
508c4762a1bSJed Brown }
5099566063dSJacob Faibussowitsch PetscCall(VecAXPY(C, -1.0, user->q));
5103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
511c4762a1bSJed Brown }
512c4762a1bSJed Brown
Scatter(Vec x,Vec sub1,VecScatter scat1,Vec sub2,VecScatter scat2)513d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
514d71ae5a4SJacob Faibussowitsch {
515c4762a1bSJed Brown PetscFunctionBegin;
5169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat1, x, sub1, INSERT_VALUES, SCATTER_FORWARD));
5179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat1, x, sub1, INSERT_VALUES, SCATTER_FORWARD));
518c4762a1bSJed Brown if (sub2) {
5199566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat2, x, sub2, INSERT_VALUES, SCATTER_FORWARD));
5209566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat2, x, sub2, INSERT_VALUES, SCATTER_FORWARD));
521c4762a1bSJed Brown }
5223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
523c4762a1bSJed Brown }
524c4762a1bSJed Brown
Gather(Vec x,Vec sub1,VecScatter scat1,Vec sub2,VecScatter scat2)525d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
526d71ae5a4SJacob Faibussowitsch {
527c4762a1bSJed Brown PetscFunctionBegin;
5289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat1, sub1, x, INSERT_VALUES, SCATTER_REVERSE));
5299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat1, sub1, x, INSERT_VALUES, SCATTER_REVERSE));
530c4762a1bSJed Brown if (sub2) {
5319566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat2, sub2, x, INSERT_VALUES, SCATTER_REVERSE));
5329566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat2, sub2, x, INSERT_VALUES, SCATTER_REVERSE));
533c4762a1bSJed Brown }
5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
535c4762a1bSJed Brown }
536c4762a1bSJed Brown
EllipticInitialize(AppCtx * user)537d71ae5a4SJacob Faibussowitsch PetscErrorCode EllipticInitialize(AppCtx *user)
538d71ae5a4SJacob Faibussowitsch {
539c4762a1bSJed Brown PetscInt m, n, i, j, k, l, linear_index, is, js, ks, ls, istart, iend, iblock;
540c4762a1bSJed Brown Vec XX, YY, ZZ, XXwork, YYwork, ZZwork, UTwork;
541c4762a1bSJed Brown PetscReal *x, *y, *z;
542c4762a1bSJed Brown PetscReal h, meanut;
543c4762a1bSJed Brown PetscScalar hinv, neg_hinv, half = 0.5, sqrt_beta;
544c4762a1bSJed Brown PetscInt im, indx1, indx2, indy1, indy2, indz1, indz2, nx, ny, nz;
545c4762a1bSJed Brown IS is_alldesign, is_allstate;
546c4762a1bSJed Brown IS is_from_d;
547c4762a1bSJed Brown IS is_from_y;
548c4762a1bSJed Brown PetscInt lo, hi, hi2, lo2, ysubnlocal, dsubnlocal;
549c4762a1bSJed Brown const PetscInt *ranges, *subranges;
550c4762a1bSJed Brown PetscMPIInt size;
551c4762a1bSJed Brown PetscReal xri, yri, zri, xim, yim, zim, dx1, dx2, dy1, dy2, dz1, dz2, Dx, Dy, Dz;
552c4762a1bSJed Brown PetscScalar v, vx, vy, vz;
553c4762a1bSJed Brown PetscInt offset, subindex, subvec, nrank, kk;
554c4762a1bSJed Brown
5559371c9d4SSatish 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,
5569371c9d4SSatish 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,
5579371c9d4SSatish 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};
558c4762a1bSJed Brown
5599371c9d4SSatish 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,
5609371c9d4SSatish 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,
5619371c9d4SSatish 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};
562c4762a1bSJed Brown
5639371c9d4SSatish 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,
5649371c9d4SSatish 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,
5659371c9d4SSatish 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};
566c4762a1bSJed Brown
567c4762a1bSJed Brown PetscFunctionBegin;
5689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
5699566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Elliptic Setup", &user->stages[0]));
5709566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[0]));
571c4762a1bSJed Brown
572c4762a1bSJed Brown /* Create u,y,c,x */
5739566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->u));
5749566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->y));
5759566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->c));
5769566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->u, PETSC_DECIDE, user->ndesign));
5779566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->u));
5789566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(user->u, &ysubnlocal));
5799566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->y, ysubnlocal * user->ns, user->nstate));
5809566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->c, ysubnlocal * user->ns, user->m));
5819566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->y));
5829566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->c));
583c4762a1bSJed Brown
584c4762a1bSJed Brown /*
585c4762a1bSJed Brown *******************************
586c4762a1bSJed Brown Create scatters for x <-> y,u
587c4762a1bSJed Brown *******************************
588c4762a1bSJed Brown
589c4762a1bSJed Brown If the state vector y and design vector u are partitioned as
590c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
591c4762a1bSJed Brown then the solution vector x is organized as
592c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np].
593c4762a1bSJed Brown The index sets user->s_is and user->d_is correspond to the indices of the
594c4762a1bSJed Brown state and design variables owned by the current processor.
595c4762a1bSJed Brown */
5969566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->x));
597c4762a1bSJed Brown
5989566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->y, &lo, &hi));
5999566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->u, &lo2, &hi2));
600c4762a1bSJed Brown
6019566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate));
6029566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user->s_is));
6039566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign));
6049566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user->d_is));
605c4762a1bSJed Brown
6069566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->x, hi - lo + hi2 - lo2, user->n));
6079566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->x));
608c4762a1bSJed Brown
6099566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->x, user->s_is, user->y, is_allstate, &user->state_scatter));
6109566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->x, user->d_is, user->u, is_alldesign, &user->design_scatter));
6119566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign));
6129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate));
613c4762a1bSJed Brown /*
614c4762a1bSJed Brown *******************************
615c4762a1bSJed Brown Create scatter from y to y_1,y_2,...,y_ns
616c4762a1bSJed Brown *******************************
617c4762a1bSJed Brown */
6189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns, &user->yi_scatter));
6199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->suby));
6209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->subq));
621c4762a1bSJed Brown
6229566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->y, &lo2, &hi2));
623c4762a1bSJed Brown istart = 0;
624c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
6259566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->suby, &lo, &hi));
6269566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_y));
6279566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y, is_from_y, user->suby, NULL, &user->yi_scatter[i]));
628c4762a1bSJed Brown istart = istart + hi - lo;
6299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y));
630c4762a1bSJed Brown }
631c4762a1bSJed Brown /*
632c4762a1bSJed Brown *******************************
633c4762a1bSJed Brown Create scatter from d to d_1,d_2,...,d_ns
634c4762a1bSJed Brown *******************************
635c4762a1bSJed Brown */
6369566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->subd));
6379566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->subd, PETSC_DECIDE, user->ndata));
6389566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->subd));
6399566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->d));
6409566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(user->subd, &dsubnlocal));
6419566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->d, dsubnlocal * user->ns, user->ndata * user->ns));
6429566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->d));
6439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns, &user->di_scatter));
644c4762a1bSJed Brown
6459566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->d, &lo2, &hi2));
646c4762a1bSJed Brown istart = 0;
647c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
6489566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->subd, &lo, &hi));
6499566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_d));
6509566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->d, is_from_d, user->subd, NULL, &user->di_scatter[i]));
651c4762a1bSJed Brown istart = istart + hi - lo;
6529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_d));
653c4762a1bSJed Brown }
654c4762a1bSJed Brown
6559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &x));
6569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &y));
6579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &z));
658c4762a1bSJed Brown
659c4762a1bSJed Brown user->ksp_its = 0;
660c4762a1bSJed Brown user->ksp_its_initial = 0;
661c4762a1bSJed Brown
662c4762a1bSJed Brown n = user->mx * user->mx * user->mx;
663c4762a1bSJed Brown m = 3 * user->mx * user->mx * (user->mx - 1);
664c4762a1bSJed Brown sqrt_beta = PetscSqrtScalar(user->beta);
665c4762a1bSJed Brown
6669566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &XX));
6679566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q));
6689566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX, ysubnlocal, n));
6699566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q, ysubnlocal * user->ns, user->m));
6709566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX));
6719566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q));
672c4762a1bSJed Brown
6739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YY));
6749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &ZZ));
6759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &XXwork));
6769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YYwork));
6779566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &ZZwork));
6789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &UTwork));
6799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->utrue));
680c4762a1bSJed Brown
681c4762a1bSJed Brown /* map for striding q */
6829566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(user->q, &ranges));
6839566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(user->u, &subranges));
684c4762a1bSJed Brown
6859566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->q, &lo2, &hi2));
6869566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->u, &lo, &hi));
687c4762a1bSJed Brown /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
688c4762a1bSJed Brown h = 1.0 / user->mx;
689c4762a1bSJed Brown hinv = user->mx;
690c4762a1bSJed Brown neg_hinv = -hinv;
691c4762a1bSJed Brown
6929566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX, &istart, &iend));
693c4762a1bSJed Brown for (linear_index = istart; linear_index < iend; linear_index++) {
694c4762a1bSJed Brown i = linear_index % user->mx;
695c4762a1bSJed Brown j = ((linear_index - i) / user->mx) % user->mx;
696c4762a1bSJed Brown k = ((linear_index - i) / user->mx - j) / user->mx;
697c4762a1bSJed Brown vx = h * (i + 0.5);
698c4762a1bSJed Brown vy = h * (j + 0.5);
699c4762a1bSJed Brown vz = h * (k + 0.5);
7009566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES));
7019566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES));
7029566063dSJacob Faibussowitsch PetscCall(VecSetValues(ZZ, 1, &linear_index, &vz, INSERT_VALUES));
703c4762a1bSJed Brown for (is = 0; is < 2; is++) {
704c4762a1bSJed Brown for (js = 0; js < 2; js++) {
705c4762a1bSJed Brown for (ks = 0; ks < 2; ks++) {
706c4762a1bSJed Brown ls = is * 4 + js * 2 + ks;
707c4762a1bSJed Brown if (ls < user->ns) {
708c4762a1bSJed Brown l = ls * n + linear_index;
709c4762a1bSJed Brown /* remap */
710c4762a1bSJed Brown subindex = l % n;
711c4762a1bSJed Brown subvec = l / n;
712c4762a1bSJed Brown nrank = 0;
713c4762a1bSJed Brown while (subindex >= subranges[nrank + 1]) nrank++;
714c4762a1bSJed Brown offset = subindex - subranges[nrank];
715c4762a1bSJed Brown istart = 0;
716c4762a1bSJed Brown for (kk = 0; kk < nrank; kk++) istart += user->ns * (subranges[kk + 1] - subranges[kk]);
717c4762a1bSJed Brown istart += (subranges[nrank + 1] - subranges[nrank]) * subvec;
718c4762a1bSJed Brown l = istart + offset;
719c4762a1bSJed Brown v = 100 * PetscSinScalar(2 * PETSC_PI * (vx + 0.25 * is)) * PetscSinScalar(2 * PETSC_PI * (vy + 0.25 * js)) * PetscSinScalar(2 * PETSC_PI * (vz + 0.25 * ks));
7209566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->q, 1, &l, &v, INSERT_VALUES));
721c4762a1bSJed Brown }
722c4762a1bSJed Brown }
723c4762a1bSJed Brown }
724c4762a1bSJed Brown }
725c4762a1bSJed Brown }
726c4762a1bSJed Brown
7279566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX));
7289566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX));
7299566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY));
7309566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY));
7319566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(ZZ));
7329566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(ZZ));
7339566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->q));
7349566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->q));
735c4762a1bSJed Brown
736c4762a1bSJed Brown /* Compute true parameter function
737c4762a1bSJed Brown ut = exp(-((x-0.25)^2+(y-0.25)^2+(z-0.25)^2)/0.05) - exp((x-0.75)^2-(y-0.75)^2-(z-0.75))^2/0.05) */
7389566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork));
7399566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork));
7409566063dSJacob Faibussowitsch PetscCall(VecCopy(ZZ, ZZwork));
741c4762a1bSJed Brown
7429566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.25));
7439566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.25));
7449566063dSJacob Faibussowitsch PetscCall(VecShift(ZZwork, -0.25));
745c4762a1bSJed Brown
7469566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
7479566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
7489566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork));
749c4762a1bSJed Brown
7509566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, UTwork));
7519566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, YYwork));
7529566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, ZZwork));
7539566063dSJacob Faibussowitsch PetscCall(VecScale(UTwork, -20.0));
7549566063dSJacob Faibussowitsch PetscCall(VecExp(UTwork));
7559566063dSJacob Faibussowitsch PetscCall(VecCopy(UTwork, user->utrue));
756c4762a1bSJed Brown
7579566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork));
7589566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork));
7599566063dSJacob Faibussowitsch PetscCall(VecCopy(ZZ, ZZwork));
760c4762a1bSJed Brown
7619566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.75));
7629566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.75));
7639566063dSJacob Faibussowitsch PetscCall(VecShift(ZZwork, -0.75));
764c4762a1bSJed Brown
7659566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
7669566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
7679566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork));
768c4762a1bSJed Brown
7699566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, UTwork));
7709566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, YYwork));
7719566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, ZZwork));
7729566063dSJacob Faibussowitsch PetscCall(VecScale(UTwork, -20.0));
7739566063dSJacob Faibussowitsch PetscCall(VecExp(UTwork));
774c4762a1bSJed Brown
7759566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue, -1.0, UTwork));
776c4762a1bSJed Brown
7779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX));
7789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY));
7799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZ));
7809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork));
7819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork));
7829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZwork));
7839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&UTwork));
784c4762a1bSJed Brown
785c4762a1bSJed Brown /* Initial guess and reference model */
7869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue, &user->ur));
7879566063dSJacob Faibussowitsch PetscCall(VecSum(user->utrue, &meanut));
788c4762a1bSJed Brown meanut = meanut / n;
7899566063dSJacob Faibussowitsch PetscCall(VecSet(user->ur, meanut));
7909566063dSJacob Faibussowitsch PetscCall(VecCopy(user->ur, user->u));
791c4762a1bSJed Brown
792c4762a1bSJed Brown /* Generate Grad matrix */
7939566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad));
7949566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, ysubnlocal, m, n));
7959566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad));
7969566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad, 2, NULL, 2, NULL));
7979566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad, 2, NULL));
7989566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));
799c4762a1bSJed Brown
800c4762a1bSJed Brown for (i = istart; i < iend; i++) {
801c4762a1bSJed Brown if (i < m / 3) {
802c4762a1bSJed Brown iblock = i / (user->mx - 1);
803c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1));
8049566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
805c4762a1bSJed Brown j = j + 1;
8069566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
807c4762a1bSJed Brown }
808c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) {
809c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1));
810c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
8119566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
812c4762a1bSJed Brown j = j + user->mx;
8139566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
814c4762a1bSJed Brown }
815c4762a1bSJed Brown if (i >= 2 * m / 3) {
816c4762a1bSJed Brown j = i - 2 * m / 3;
8179566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
818c4762a1bSJed Brown j = j + user->mx * user->mx;
8199566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
820c4762a1bSJed Brown }
821c4762a1bSJed Brown }
822c4762a1bSJed Brown
8239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY));
8249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY));
825c4762a1bSJed Brown
826c4762a1bSJed Brown /* Generate arithmetic averaging matrix Av */
8279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Av));
8289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Av, PETSC_DECIDE, ysubnlocal, m, n));
8299566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Av));
8309566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Av, 2, NULL, 2, NULL));
8319566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Av, 2, NULL));
8329566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Av, &istart, &iend));
833c4762a1bSJed Brown
834c4762a1bSJed Brown for (i = istart; i < iend; i++) {
835c4762a1bSJed Brown if (i < m / 3) {
836c4762a1bSJed Brown iblock = i / (user->mx - 1);
837c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1));
8389566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
839c4762a1bSJed Brown j = j + 1;
8409566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
841c4762a1bSJed Brown }
842c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) {
843c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1));
844c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
8459566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
846c4762a1bSJed Brown j = j + user->mx;
8479566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
848c4762a1bSJed Brown }
849c4762a1bSJed Brown if (i >= 2 * m / 3) {
850c4762a1bSJed Brown j = i - 2 * m / 3;
8519566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
852c4762a1bSJed Brown j = j + user->mx * user->mx;
8539566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
854c4762a1bSJed Brown }
855c4762a1bSJed Brown }
856c4762a1bSJed Brown
8579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Av, MAT_FINAL_ASSEMBLY));
8589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Av, MAT_FINAL_ASSEMBLY));
859c4762a1bSJed Brown
8609566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->L));
8619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->L, PETSC_DECIDE, ysubnlocal, m + n, n));
8629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->L));
8639566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->L, 2, NULL, 2, NULL));
8649566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->L, 2, NULL));
8659566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->L, &istart, &iend));
866c4762a1bSJed Brown
867c4762a1bSJed Brown for (i = istart; i < iend; i++) {
868c4762a1bSJed Brown if (i < m / 3) {
869c4762a1bSJed Brown iblock = i / (user->mx - 1);
870c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1));
8719566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
872c4762a1bSJed Brown j = j + 1;
8739566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
874c4762a1bSJed Brown }
875c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) {
876c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1));
877c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
8789566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
879c4762a1bSJed Brown j = j + user->mx;
8809566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
881c4762a1bSJed Brown }
882c4762a1bSJed Brown if (i >= 2 * m / 3 && i < m) {
883c4762a1bSJed Brown j = i - 2 * m / 3;
8849566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
885c4762a1bSJed Brown j = j + user->mx * user->mx;
8869566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
887c4762a1bSJed Brown }
888c4762a1bSJed Brown if (i >= m) {
889c4762a1bSJed Brown j = i - m;
8909566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &sqrt_beta, INSERT_VALUES));
891c4762a1bSJed Brown }
892c4762a1bSJed Brown }
8939566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->L, MAT_FINAL_ASSEMBLY));
8949566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->L, MAT_FINAL_ASSEMBLY));
8959566063dSJacob Faibussowitsch PetscCall(MatScale(user->L, PetscPowScalar(h, 1.5)));
896c4762a1bSJed Brown
897c4762a1bSJed Brown /* Generate Div matrix */
898c4762a1bSJed Brown if (!user->use_ptap) {
899c4762a1bSJed Brown /* Generate Div matrix */
9009566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Div));
9019566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Div, ysubnlocal, PETSC_DECIDE, n, m));
9029566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Div));
9039566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Div, 4, NULL, 4, NULL));
9049566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Div, 6, NULL));
9059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));
906c4762a1bSJed Brown
907c4762a1bSJed Brown for (i = istart; i < iend; i++) {
908c4762a1bSJed Brown if (i < m / 3) {
909c4762a1bSJed Brown iblock = i / (user->mx - 1);
910c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1));
9119566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES));
912c4762a1bSJed Brown j = j + 1;
9139566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES));
914c4762a1bSJed Brown }
915c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) {
916c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1));
917c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
9189566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES));
919c4762a1bSJed Brown j = j + user->mx;
9209566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES));
921c4762a1bSJed Brown }
922c4762a1bSJed Brown if (i >= 2 * m / 3) {
923c4762a1bSJed Brown j = i - 2 * m / 3;
9249566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES));
925c4762a1bSJed Brown j = j + user->mx * user->mx;
9269566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES));
927c4762a1bSJed Brown }
928c4762a1bSJed Brown }
929c4762a1bSJed Brown
9309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Div, MAT_FINAL_ASSEMBLY));
9319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Div, MAT_FINAL_ASSEMBLY));
9329566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork));
933c4762a1bSJed Brown } else {
9349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Diag));
9359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Diag, PETSC_DECIDE, PETSC_DECIDE, m, m));
9369566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Diag));
9379566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Diag, 1, NULL, 0, NULL));
9389566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Diag, 1, NULL));
939c4762a1bSJed Brown }
940c4762a1bSJed Brown
941c4762a1bSJed Brown /* Build work vectors and matrices */
9429566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->S));
9439566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->S, PETSC_DECIDE, m));
9449566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->S));
945c4762a1bSJed Brown
9469566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork));
9479566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, m + user->mx * user->mx * user->mx));
9489566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork));
949c4762a1bSJed Brown
9509566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Av, MAT_SHARE_NONZERO_PATTERN, &user->Avwork));
951c4762a1bSJed Brown
9529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Swork));
9539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Sdiag));
9549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Av_u));
9559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Twork));
9569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ywork));
9579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->Ywork));
9589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->uwork));
9599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->js_diag));
9609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c, &user->cwork));
9619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->d, &user->dwork));
962c4762a1bSJed Brown
963c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */
9649566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal, user->nstate, user->ndesign, user, &user->Jd));
965*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (PetscErrorCodeFn *)DesignMatMult));
966*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)DesignMatMultTranspose));
967c4762a1bSJed Brown
968c4762a1bSJed Brown /* Compute true state function ytrue given utrue */
9699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ytrue));
970c4762a1bSJed Brown
971c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */
9729566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
9739566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); /* Note: user->utrue */
9749566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
9759566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
976c4762a1bSJed Brown
977c4762a1bSJed Brown /* Next form DSG = Div*S*Grad */
9789566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u, user->Swork));
9799566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork));
980c4762a1bSJed Brown if (user->use_ptap) {
9819566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES));
9829566063dSJacob Faibussowitsch PetscCall(MatPtAP(user->Diag, user->Grad, MAT_INITIAL_MATRIX, 1.0, &user->DSG));
983c4762a1bSJed Brown } else {
9849566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
9859566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork));
986c20d7725SJed Brown
9879566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Divwork, user->Grad, MAT_INITIAL_MATRIX, 1.0, &user->DSG));
988c4762a1bSJed Brown }
989c4762a1bSJed Brown
9909566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRIC, PETSC_TRUE));
9919566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
992c4762a1bSJed Brown
993c4762a1bSJed Brown if (user->use_lrc == PETSC_TRUE) {
994c4762a1bSJed Brown v = PetscSqrtReal(1.0 / user->ndesign);
9959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ndesign, &user->ones));
996c4762a1bSJed Brown
997ad540459SPierre Jolivet for (i = 0; i < user->ndesign; i++) user->ones[i] = v;
9989566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, ysubnlocal, PETSC_DECIDE, user->ndesign, 1, user->ones, &user->Ones));
9999566063dSJacob Faibussowitsch PetscCall(MatCreateLRC(user->DSG, user->Ones, NULL, user->Ones, &user->JsBlock));
10009566063dSJacob Faibussowitsch PetscCall(MatSetUp(user->JsBlock));
1001c4762a1bSJed Brown } else {
1002c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing (A + h^3*e*e^T)*x */
10039566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal, ysubnlocal, user->ndesign, user->ndesign, user, &user->JsBlock));
1004*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (PetscErrorCodeFn *)StateBlockMatMult));
1005*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateBlockMatMult));
1006c4762a1bSJed Brown }
10079566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsBlock, MAT_SYMMETRIC, PETSC_TRUE));
10089566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsBlock, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
10099566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal * user->ns, user->nstate, user->nstate, user, &user->Js));
1010*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (PetscErrorCodeFn *)StateMatMult));
1011*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatMult));
10129566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->Js, MAT_SYMMETRIC, PETSC_TRUE));
10139566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->Js, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
1014c4762a1bSJed Brown
10159566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal * user->ns, user->nstate, user->nstate, user, &user->JsInv));
1016*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (PetscErrorCodeFn *)StateInvMatMult));
1017*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateInvMatMult));
10189566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsInv, MAT_SYMMETRIC, PETSC_TRUE));
10199566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsInv, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
1020c4762a1bSJed Brown
10219566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRIC, PETSC_TRUE));
10229566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
1023c4762a1bSJed Brown /* Now solve for ytrue */
10249566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver));
10259566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver));
1026c4762a1bSJed Brown
10279566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->DSG));
1028c4762a1bSJed Brown
10299566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsInv, user->q, user->ytrue));
1030c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */
10319566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0));
10329566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u)); /* Note: user->u */
10339566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork));
10349566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
1035c4762a1bSJed Brown
1036c4762a1bSJed Brown /* Next update DSG = Div*S*Grad with user->u */
10379566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u, user->Swork));
10389566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork));
1039c4762a1bSJed Brown if (user->use_ptap) {
10409566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES));
10419566063dSJacob Faibussowitsch PetscCall(MatPtAP(user->Diag, user->Grad, MAT_REUSE_MATRIX, 1.0, &user->DSG));
1042c4762a1bSJed Brown } else {
10439566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
10449566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u));
10459566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG));
1046c4762a1bSJed Brown }
1047c4762a1bSJed Brown
1048c4762a1bSJed Brown /* Now solve for y */
1049c4762a1bSJed Brown
10509566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsInv, user->q, user->y));
1051c4762a1bSJed Brown
1052c4762a1bSJed Brown user->ksp_its_initial = user->ksp_its;
1053c4762a1bSJed Brown user->ksp_its = 0;
1054c4762a1bSJed Brown /* Construct projection matrix Q (blocks) */
10559566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Q));
10569566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Q, dsubnlocal, ysubnlocal, user->ndata, user->ndesign));
10579566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Q));
10589566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Q, 8, NULL, 8, NULL));
10599566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Q, 8, NULL));
1060c4762a1bSJed Brown
1061c4762a1bSJed Brown for (i = 0; i < user->mx; i++) {
1062c4762a1bSJed Brown x[i] = h * (i + 0.5);
1063c4762a1bSJed Brown y[i] = h * (i + 0.5);
1064c4762a1bSJed Brown z[i] = h * (i + 0.5);
1065c4762a1bSJed Brown }
10669566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Q, &istart, &iend));
1067c4762a1bSJed Brown
10689371c9d4SSatish Balay nx = user->mx;
10699371c9d4SSatish Balay ny = user->mx;
10709371c9d4SSatish Balay nz = user->mx;
1071c4762a1bSJed Brown for (i = istart; i < iend; i++) {
1072c4762a1bSJed Brown xri = xr[i];
1073c4762a1bSJed Brown im = 0;
1074c4762a1bSJed Brown xim = x[im];
1075c4762a1bSJed Brown while (xri > xim && im < nx) {
1076c4762a1bSJed Brown im = im + 1;
1077c4762a1bSJed Brown xim = x[im];
1078c4762a1bSJed Brown }
1079c4762a1bSJed Brown indx1 = im - 1;
1080c4762a1bSJed Brown indx2 = im;
1081c4762a1bSJed Brown dx1 = xri - x[indx1];
1082c4762a1bSJed Brown dx2 = x[indx2] - xri;
1083c4762a1bSJed Brown
1084c4762a1bSJed Brown yri = yr[i];
1085c4762a1bSJed Brown im = 0;
1086c4762a1bSJed Brown yim = y[im];
1087c4762a1bSJed Brown while (yri > yim && im < ny) {
1088c4762a1bSJed Brown im = im + 1;
1089c4762a1bSJed Brown yim = y[im];
1090c4762a1bSJed Brown }
1091c4762a1bSJed Brown indy1 = im - 1;
1092c4762a1bSJed Brown indy2 = im;
1093c4762a1bSJed Brown dy1 = yri - y[indy1];
1094c4762a1bSJed Brown dy2 = y[indy2] - yri;
1095c4762a1bSJed Brown
1096c4762a1bSJed Brown zri = zr[i];
1097c4762a1bSJed Brown im = 0;
1098c4762a1bSJed Brown zim = z[im];
1099c4762a1bSJed Brown while (zri > zim && im < nz) {
1100c4762a1bSJed Brown im = im + 1;
1101c4762a1bSJed Brown zim = z[im];
1102c4762a1bSJed Brown }
1103c4762a1bSJed Brown indz1 = im - 1;
1104c4762a1bSJed Brown indz2 = im;
1105c4762a1bSJed Brown dz1 = zri - z[indz1];
1106c4762a1bSJed Brown dz2 = z[indz2] - zri;
1107c4762a1bSJed Brown
1108c4762a1bSJed Brown Dx = x[indx2] - x[indx1];
1109c4762a1bSJed Brown Dy = y[indy2] - y[indy1];
1110c4762a1bSJed Brown Dz = z[indz2] - z[indz1];
1111c4762a1bSJed Brown
1112c4762a1bSJed Brown j = indx1 + indy1 * nx + indz1 * nx * ny;
1113c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
11149566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1115c4762a1bSJed Brown
1116c4762a1bSJed Brown j = indx1 + indy1 * nx + indz2 * nx * ny;
1117c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
11189566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1119c4762a1bSJed Brown
1120c4762a1bSJed Brown j = indx1 + indy2 * nx + indz1 * nx * ny;
1121c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
11229566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1123c4762a1bSJed Brown
1124c4762a1bSJed Brown j = indx1 + indy2 * nx + indz2 * nx * ny;
1125c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
11269566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1127c4762a1bSJed Brown
1128c4762a1bSJed Brown j = indx2 + indy1 * nx + indz1 * nx * ny;
1129c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
11309566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1131c4762a1bSJed Brown
1132c4762a1bSJed Brown j = indx2 + indy1 * nx + indz2 * nx * ny;
1133c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
11349566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1135c4762a1bSJed Brown
1136c4762a1bSJed Brown j = indx2 + indy2 * nx + indz1 * nx * ny;
1137c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
11389566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1139c4762a1bSJed Brown
1140c4762a1bSJed Brown j = indx2 + indy2 * nx + indz2 * nx * ny;
1141c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
11429566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1143c4762a1bSJed Brown }
1144c4762a1bSJed Brown
11459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Q, MAT_FINAL_ASSEMBLY));
11469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Q, MAT_FINAL_ASSEMBLY));
1147c4762a1bSJed Brown /* Create MQ (composed of blocks of Q */
11489566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, dsubnlocal * user->ns, PETSC_DECIDE, user->ndata * user->ns, user->nstate, user, &user->MQ));
1149*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->MQ, MATOP_MULT, (PetscErrorCodeFn *)QMatMult));
1150*57d50842SBarry Smith PetscCall(MatShellSetOperation(user->MQ, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)QMatMultTranspose));
1151c4762a1bSJed Brown
1152c4762a1bSJed Brown /* Add noise to the measurement data */
11539566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork, 1.0));
11549566063dSJacob Faibussowitsch PetscCall(VecAYPX(user->ywork, user->noise, user->ytrue));
11559566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->ywork, user->d));
1156c4762a1bSJed Brown
1157c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11589566063dSJacob Faibussowitsch PetscCall(PetscFree(x));
11599566063dSJacob Faibussowitsch PetscCall(PetscFree(y));
11609566063dSJacob Faibussowitsch PetscCall(PetscFree(z));
11619566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop());
11623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1163c4762a1bSJed Brown }
1164c4762a1bSJed Brown
EllipticDestroy(AppCtx * user)1165d71ae5a4SJacob Faibussowitsch PetscErrorCode EllipticDestroy(AppCtx *user)
1166d71ae5a4SJacob Faibussowitsch {
1167c4762a1bSJed Brown PetscInt i;
1168c4762a1bSJed Brown
1169c4762a1bSJed Brown PetscFunctionBegin;
11709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->DSG));
11719566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver));
11729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Q));
11739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->MQ));
1174c4762a1bSJed Brown if (!user->use_ptap) {
11759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div));
11769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork));
1177c4762a1bSJed Brown } else {
11789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Diag));
1179c4762a1bSJed Brown }
118048a46eb9SPierre Jolivet if (user->use_lrc) PetscCall(MatDestroy(&user->Ones));
1181c4762a1bSJed Brown
11829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad));
11839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Av));
11849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Avwork));
11859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L));
11869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js));
11879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd));
11889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock));
11899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv));
1190c4762a1bSJed Brown
11919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->x));
11929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u));
11939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork));
11949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue));
11959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y));
11969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork));
11979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue));
11989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c));
11999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork));
12009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur));
12019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q));
12029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d));
12039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork));
12049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork));
12059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->S));
12069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Swork));
12079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Sdiag));
12089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Ywork));
12099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Twork));
12109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Av_u));
12119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag));
12129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is));
12139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is));
12149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->suby));
12159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->subd));
12169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->subq));
12179566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter));
12189566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter));
1219c4762a1bSJed Brown for (i = 0; i < user->ns; i++) {
12209566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
12219566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->di_scatter[i]));
1222c4762a1bSJed Brown }
12239566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter));
12249566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di_scatter));
1225c4762a1bSJed Brown if (user->use_lrc) {
12269566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ones));
12279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Ones));
1228c4762a1bSJed Brown }
12293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1230c4762a1bSJed Brown }
1231c4762a1bSJed Brown
EllipticMonitor(Tao tao,void * ptr)1232d71ae5a4SJacob Faibussowitsch PetscErrorCode EllipticMonitor(Tao tao, void *ptr)
1233d71ae5a4SJacob Faibussowitsch {
1234c4762a1bSJed Brown Vec X;
1235c4762a1bSJed Brown PetscReal unorm, ynorm;
1236c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
1237c4762a1bSJed Brown
1238c4762a1bSJed Brown PetscFunctionBegin;
12399566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X));
12409566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
12419566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue));
12429566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue));
12439566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork, NORM_2, &unorm));
12449566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork, NORM_2, &ynorm));
12459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm));
12463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1247c4762a1bSJed Brown }
1248c4762a1bSJed Brown
1249c4762a1bSJed Brown /*TEST
1250c4762a1bSJed Brown
1251c4762a1bSJed Brown build:
1252c4762a1bSJed Brown requires: !complex
1253c4762a1bSJed Brown
1254c4762a1bSJed Brown test:
125510978b7dSBarry Smith args: -tao_monitor_constraint_norm -ns 1 -tao_type lcl -tao_gatol 1.e-3 -tao_max_it 11
1256c4762a1bSJed Brown requires: !single
1257c4762a1bSJed Brown
1258c4762a1bSJed Brown test:
1259c4762a1bSJed Brown suffix: 2
126010978b7dSBarry Smith args: -tao_monitor_constraint_norm -tao_type lcl -tao_max_it 11 -use_ptap -use_lrc -ns 1 -tao_gatol 1.e-3
1261c4762a1bSJed Brown requires: !single
1262c4762a1bSJed Brown
1263c4762a1bSJed Brown TEST*/
1264