xref: /petsc/src/tao/pde_constrained/tutorials/elliptic.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
99c4762a1bSJed Brown int main(int argc, char **argv)
100c4762a1bSJed Brown {
101c4762a1bSJed Brown   Vec                x0;
102c4762a1bSJed Brown   Tao                tao;
103c4762a1bSJed Brown   AppCtx             user;
104c4762a1bSJed Brown   PetscInt           ntests = 1;
105c4762a1bSJed Brown   PetscInt           i;
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char*)0,help));
108c4762a1bSJed Brown   user.mx = 8;
109*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"elliptic example",NULL);
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL));
111c4762a1bSJed Brown   user.ns = 6;
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ns","Number of data samples (1<=ns<=8)","",user.ns,&user.ns,NULL));
113c4762a1bSJed Brown   user.ndata = 64;
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL));
115c4762a1bSJed Brown   user.alpha = 0.1;
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL));
117c4762a1bSJed Brown   user.beta = 0.00001;
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,NULL));
119c4762a1bSJed Brown   user.noise = 0.01;
1209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,NULL));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   user.use_ptap = PETSC_FALSE;
1239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_ptap","Use ptap matrix for DSG","",user.use_ptap,&user.use_ptap,NULL));
124c4762a1bSJed Brown   user.use_lrc = PETSC_FALSE;
1259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_lrc","Use lrc matrix for Js","",user.use_lrc,&user.use_lrc,NULL));
1269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL));
127*d0609cedSBarry Smith   PetscOptionsEnd();
12876280437SVaclav Hapla 
129c4762a1bSJed Brown   user.m = user.ns*user.mx*user.mx*user.mx; /* number of constraints */
130c4762a1bSJed Brown   user.nstate =  user.m;
131c4762a1bSJed Brown   user.ndesign = user.mx*user.mx*user.mx;
132c4762a1bSJed Brown   user.n = user.nstate + user.ndesign; /* number of variables */
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1359566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
1369566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOLCL));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* Set up initial vectors and matrices */
1399566063dSJacob Faibussowitsch   PetscCall(EllipticInitialize(&user));
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(Gather(user.x,user.y,user.state_scatter,user.u,user.design_scatter));
1429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.x,&x0));
1439566063dSJacob Faibussowitsch   PetscCall(VecCopy(user.x,x0));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Set solution vector with an initial guess */
1469566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,user.x));
1479566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(tao, FormFunction, &user));
1489566063dSJacob Faibussowitsch   PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
1499566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, NULL, user.JsInv, FormJacobianState, &user));
1529566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(TaoSetStateDesignIS(tao,user.s_is,user.d_is));
1559566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1589566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Trials",&user.stages[1]));
1599566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user.stages[1]));
160c4762a1bSJed Brown   for (i=0; i<ntests; i++) {
1619566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
1629566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its));
1639566063dSJacob Faibussowitsch     PetscCall(VecCopy(x0,user.x));
164c4762a1bSJed Brown   }
1659566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1669566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)user.x));
1679566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: "));
1689566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x0));
1729566063dSJacob Faibussowitsch   PetscCall(EllipticDestroy(&user));
1739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
174b122ec5aSJacob Faibussowitsch   return 0;
175c4762a1bSJed Brown }
176c4762a1bSJed Brown /* ------------------------------------------------------------------- */
177c4762a1bSJed Brown /*
178c4762a1bSJed Brown    dwork = Qy - d
179c4762a1bSJed Brown    lwork = L*(u-ur)
180c4762a1bSJed Brown    f = 1/2 * (dwork.dwork + alpha*lwork.lwork)
181c4762a1bSJed Brown */
182c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
183c4762a1bSJed Brown {
184c4762a1bSJed Brown   PetscReal      d1=0,d2=0;
185c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   PetscFunctionBegin;
1889566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
1899566063dSJacob Faibussowitsch   PetscCall(MatMult(user->MQ,user->y,user->dwork));
1909566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,-1.0,user->d));
1919566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork,user->dwork,&d1));
1929566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u));
1939566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L,user->uwork,user->lwork));
1949566063dSJacob Faibussowitsch   PetscCall(VecDot(user->lwork,user->lwork,&d2));
195c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
196c4762a1bSJed Brown   PetscFunctionReturn(0);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown /* ------------------------------------------------------------------- */
200c4762a1bSJed Brown /*
201c4762a1bSJed Brown     state: g_s = Q' *(Qy - d)
202c4762a1bSJed Brown     design: g_d = alpha*L'*L*(u-ur)
203c4762a1bSJed Brown */
204c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
205c4762a1bSJed Brown {
206c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   PetscFunctionBegin;
2099566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2109566063dSJacob Faibussowitsch   PetscCall(MatMult(user->MQ,user->y,user->dwork));
2119566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,-1.0,user->d));
2129566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(user->MQ,user->dwork,user->ywork));
2139566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u));
2149566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L,user->uwork,user->lwork));
2159566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(user->L,user->lwork,user->uwork));
2169566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork, user->alpha));
2179566063dSJacob Faibussowitsch   PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
218c4762a1bSJed Brown   PetscFunctionReturn(0);
219c4762a1bSJed Brown }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
222c4762a1bSJed Brown {
223c4762a1bSJed Brown   PetscReal      d1,d2;
224c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   PetscFunctionBegin;
2279566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2289566063dSJacob Faibussowitsch   PetscCall(MatMult(user->MQ,user->y,user->dwork));
2299566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->dwork,-1.0,user->d));
2309566063dSJacob Faibussowitsch   PetscCall(VecDot(user->dwork,user->dwork,&d1));
2319566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(user->MQ,user->dwork,user->ywork));
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u));
2349566063dSJacob Faibussowitsch   PetscCall(MatMult(user->L,user->uwork,user->lwork));
2359566063dSJacob Faibussowitsch   PetscCall(VecDot(user->lwork,user->lwork,&d2));
2369566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(user->L,user->lwork,user->uwork));
2379566063dSJacob Faibussowitsch   PetscCall(VecScale(user->uwork, user->alpha));
238c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
2399566063dSJacob Faibussowitsch   PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
240c4762a1bSJed Brown   PetscFunctionReturn(0);
241c4762a1bSJed Brown }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown /* ------------------------------------------------------------------- */
244c4762a1bSJed Brown /* A
245c4762a1bSJed Brown MatShell object
246c4762a1bSJed Brown */
247c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
248c4762a1bSJed Brown {
249c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
253c4762a1bSJed Brown   /* DSG = Div * (1/Av_u) * Grad */
2549566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork,0));
2559566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork,-1.0,user->u));
2569566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
2579566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av,user->uwork,user->Av_u));
2589566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->Av_u,user->Swork));
2599566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Swork));
260c4762a1bSJed Brown   if (user->use_ptap) {
2619566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES));
2629566063dSJacob Faibussowitsch     PetscCall(MatPtAP(user->Diag,user->Grad,MAT_REUSE_MATRIX,1.0,&user->DSG));
263c4762a1bSJed Brown   } else {
2649566063dSJacob Faibussowitsch     PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN));
2659566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Swork));
2669566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(user->DSG));
267c4762a1bSJed Brown   }
268c4762a1bSJed Brown   PetscFunctionReturn(0);
269c4762a1bSJed Brown }
270c4762a1bSJed Brown /* ------------------------------------------------------------------- */
271c4762a1bSJed Brown /* B */
272c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
273c4762a1bSJed Brown {
274c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   PetscFunctionBegin;
2779566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
278c4762a1bSJed Brown   PetscFunctionReturn(0);
279c4762a1bSJed Brown }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown PetscErrorCode StateBlockMatMult(Mat J_shell, Vec X, Vec Y)
282c4762a1bSJed Brown {
283c4762a1bSJed Brown   PetscReal      sum;
284c4762a1bSJed Brown   AppCtx         *user;
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   PetscFunctionBegin;
2879566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
2889566063dSJacob Faibussowitsch   PetscCall(MatMult(user->DSG,X,Y));
2899566063dSJacob Faibussowitsch   PetscCall(VecSum(X,&sum));
290c4762a1bSJed Brown   sum /= user->ndesign;
2919566063dSJacob Faibussowitsch   PetscCall(VecShift(Y,sum));
292c4762a1bSJed Brown   PetscFunctionReturn(0);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
296c4762a1bSJed Brown {
297c4762a1bSJed Brown   PetscInt       i;
298c4762a1bSJed Brown   AppCtx         *user;
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   PetscFunctionBegin;
3019566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
302c4762a1bSJed Brown   if (user->ns == 1) {
3039566063dSJacob Faibussowitsch     PetscCall(MatMult(user->JsBlock,X,Y));
304c4762a1bSJed Brown   } else {
305c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
3069566063dSJacob Faibussowitsch       PetscCall(Scatter(X,user->subq,user->yi_scatter[i],0,0));
3079566063dSJacob Faibussowitsch       PetscCall(Scatter(Y,user->suby,user->yi_scatter[i],0,0));
3089566063dSJacob Faibussowitsch       PetscCall(MatMult(user->JsBlock,user->subq,user->suby));
3099566063dSJacob Faibussowitsch       PetscCall(Gather(Y,user->suby,user->yi_scatter[i],0,0));
310c4762a1bSJed Brown     }
311c4762a1bSJed Brown   }
312c4762a1bSJed Brown   PetscFunctionReturn(0);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown PetscErrorCode StateInvMatMult(Mat J_shell, Vec X, Vec Y)
316c4762a1bSJed Brown {
317c4762a1bSJed Brown   PetscInt       its,i;
318c4762a1bSJed Brown   AppCtx         *user;
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   PetscFunctionBegin;
3219566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
3229566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(user->solver,user->JsBlock,user->DSG));
323c4762a1bSJed Brown   if (Y == user->ytrue) {
324c4762a1bSJed Brown     /* First solve is done using true solution to set up problem */
3259566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
326c4762a1bSJed Brown   } else {
3279566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
328c4762a1bSJed Brown   }
329c4762a1bSJed Brown   if (user->ns == 1) {
3309566063dSJacob Faibussowitsch     PetscCall(KSPSolve(user->solver,X,Y));
3319566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(user->solver,&its));
332c4762a1bSJed Brown     user->ksp_its+=its;
333c4762a1bSJed Brown   } else {
334c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
3359566063dSJacob Faibussowitsch       PetscCall(Scatter(X,user->subq,user->yi_scatter[i],0,0));
3369566063dSJacob Faibussowitsch       PetscCall(Scatter(Y,user->suby,user->yi_scatter[i],0,0));
3379566063dSJacob Faibussowitsch       PetscCall(KSPSolve(user->solver,user->subq,user->suby));
3389566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(user->solver,&its));
339c4762a1bSJed Brown       user->ksp_its+=its;
3409566063dSJacob Faibussowitsch       PetscCall(Gather(Y,user->suby,user->yi_scatter[i],0,0));
341c4762a1bSJed Brown     }
342c4762a1bSJed Brown   }
343c4762a1bSJed Brown   PetscFunctionReturn(0);
344c4762a1bSJed Brown }
345c4762a1bSJed Brown PetscErrorCode QMatMult(Mat J_shell, Vec X, Vec Y)
346c4762a1bSJed Brown {
347c4762a1bSJed Brown   AppCtx         *user;
348c4762a1bSJed Brown   PetscInt       i;
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
352c4762a1bSJed Brown   if (user->ns == 1) {
3539566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Q,X,Y));
354c4762a1bSJed Brown   } else {
355c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
3569566063dSJacob Faibussowitsch       PetscCall(Scatter(X,user->subq,user->yi_scatter[i],0,0));
3579566063dSJacob Faibussowitsch       PetscCall(Scatter(Y,user->subd,user->di_scatter[i],0,0));
3589566063dSJacob Faibussowitsch       PetscCall(MatMult(user->Q,user->subq,user->subd));
3599566063dSJacob Faibussowitsch       PetscCall(Gather(Y,user->subd,user->di_scatter[i],0,0));
360c4762a1bSJed Brown     }
361c4762a1bSJed Brown   }
362c4762a1bSJed Brown   PetscFunctionReturn(0);
363c4762a1bSJed Brown }
364c4762a1bSJed Brown 
365c4762a1bSJed Brown PetscErrorCode QMatMultTranspose(Mat J_shell, Vec X, Vec Y)
366c4762a1bSJed Brown {
367c4762a1bSJed Brown   AppCtx         *user;
368c4762a1bSJed Brown   PetscInt       i;
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   PetscFunctionBegin;
3719566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
372c4762a1bSJed Brown   if (user->ns == 1) {
3739566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(user->Q,X,Y));
374c4762a1bSJed Brown   } else {
375c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
3769566063dSJacob Faibussowitsch       PetscCall(Scatter(X,user->subd,user->di_scatter[i],0,0));
3779566063dSJacob Faibussowitsch       PetscCall(Scatter(Y,user->suby,user->yi_scatter[i],0,0));
3789566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(user->Q,user->subd,user->suby));
3799566063dSJacob Faibussowitsch       PetscCall(Gather(Y,user->suby,user->yi_scatter[i],0,0));
380c4762a1bSJed Brown     }
381c4762a1bSJed Brown   }
382c4762a1bSJed Brown   PetscFunctionReturn(0);
383c4762a1bSJed Brown }
384c4762a1bSJed Brown 
385c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
386c4762a1bSJed Brown {
387c4762a1bSJed Brown   PetscInt       i;
388c4762a1bSJed Brown   AppCtx         *user;
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   PetscFunctionBegin;
3919566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   /* sdiag(1./v) */
3949566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork,0));
3959566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork,-1.0,user->u));
3969566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   /* sdiag(1./((Av*(1./v)).^2)) */
3999566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av,user->uwork,user->Swork));
4009566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->Swork,user->Swork,user->Swork));
4019566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Swork));
402c4762a1bSJed Brown 
403c4762a1bSJed Brown   /* (Av * (sdiag(1./v) * b)) */
4049566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->uwork,user->uwork,X));
4059566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av,user->uwork,user->Twork));
406c4762a1bSJed Brown 
407c4762a1bSJed Brown   /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
4089566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->Swork,user->Twork,user->Swork));
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   if (user->ns == 1) {
411c4762a1bSJed Brown     /* (sdiag(Grad*y(:,i)) */
4129566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad,user->y,user->Twork));
413c4762a1bSJed Brown 
414c4762a1bSJed Brown     /* Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
4159566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->Swork,user->Twork,user->Swork));
4169566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(user->Grad,user->Swork,Y));
417c4762a1bSJed Brown   } else {
418c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
4199566063dSJacob Faibussowitsch       PetscCall(Scatter(user->y,user->suby,user->yi_scatter[i],0,0));
4209566063dSJacob Faibussowitsch       PetscCall(Scatter(Y,user->subq,user->yi_scatter[i],0,0));
421c4762a1bSJed Brown 
4229566063dSJacob Faibussowitsch       PetscCall(MatMult(user->Grad,user->suby,user->Twork));
4239566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(user->Twork,user->Twork,user->Swork));
4249566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(user->Grad,user->Twork,user->subq));
4259566063dSJacob Faibussowitsch       PetscCall(Gather(user->y,user->suby,user->yi_scatter[i],0,0));
4269566063dSJacob Faibussowitsch       PetscCall(Gather(Y,user->subq,user->yi_scatter[i],0,0));
427c4762a1bSJed Brown     }
428c4762a1bSJed Brown   }
429c4762a1bSJed Brown   PetscFunctionReturn(0);
430c4762a1bSJed Brown }
431c4762a1bSJed Brown 
432c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
433c4762a1bSJed Brown {
434c4762a1bSJed Brown   PetscInt       i;
435c4762a1bSJed Brown   AppCtx         *user;
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   PetscFunctionBegin;
4389566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(J_shell,&user));
4399566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Y));
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /* Sdiag = 1./((Av*(1./v)).^2) */
4429566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork,0));
4439566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork,-1.0,user->u));
4449566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
4459566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av,user->uwork,user->Swork));
4469566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(user->Sdiag,user->Swork,user->Swork));
4479566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Sdiag));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   for (i=0;i<user->ns;i++) {
4509566063dSJacob Faibussowitsch     PetscCall(Scatter(X,user->subq,user->yi_scatter[i],0,0));
4519566063dSJacob Faibussowitsch     PetscCall(Scatter(user->y,user->suby,user->yi_scatter[i],0,0));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown     /* Swork = (Div' * b(:,i)) */
4549566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad,user->subq,user->Swork));
455c4762a1bSJed Brown 
456c4762a1bSJed Brown     /* Twork = Grad*y(:,i) */
4579566063dSJacob Faibussowitsch     PetscCall(MatMult(user->Grad,user->suby,user->Twork));
458c4762a1bSJed Brown 
459c4762a1bSJed Brown     /* Twork = sdiag(Twork) * Swork */
4609566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->Twork,user->Swork,user->Twork));
461c4762a1bSJed Brown 
462c4762a1bSJed Brown     /* Swork = pointwisemult(Sdiag,Twork) */
4639566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->Swork,user->Twork,user->Sdiag));
464c4762a1bSJed Brown 
465c4762a1bSJed Brown     /* Ywork = Av' * Swork */
4669566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(user->Av,user->Swork,user->Ywork));
467c4762a1bSJed Brown 
468c4762a1bSJed Brown     /* Ywork = pointwisemult(uwork,Ywork) */
4699566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(user->Ywork,user->uwork,user->Ywork));
4709566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y,1.0,user->Ywork));
4719566063dSJacob Faibussowitsch     PetscCall(Gather(user->y,user->suby,user->yi_scatter[i],0,0));
472c4762a1bSJed Brown   }
473c4762a1bSJed Brown   PetscFunctionReturn(0);
474c4762a1bSJed Brown }
475c4762a1bSJed Brown 
476c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
477c4762a1bSJed Brown {
478c4762a1bSJed Brown    /* C=Ay - q      A = Div * Sigma * Grad + hx*hx*hx*ones(n,n) */
479c4762a1bSJed Brown    PetscReal      sum;
480c4762a1bSJed Brown    PetscInt       i;
481c4762a1bSJed Brown    AppCtx         *user = (AppCtx*)ptr;
482c4762a1bSJed Brown 
483c4762a1bSJed Brown    PetscFunctionBegin;
4849566063dSJacob Faibussowitsch    PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
485c4762a1bSJed Brown    if (user->ns == 1) {
4869566063dSJacob Faibussowitsch      PetscCall(MatMult(user->Grad,user->y,user->Swork));
4879566063dSJacob Faibussowitsch      PetscCall(VecPointwiseDivide(user->Swork,user->Swork,user->Av_u));
4889566063dSJacob Faibussowitsch      PetscCall(MatMultTranspose(user->Grad,user->Swork,C));
4899566063dSJacob Faibussowitsch      PetscCall(VecSum(user->y,&sum));
490c4762a1bSJed Brown      sum /= user->ndesign;
4919566063dSJacob Faibussowitsch      PetscCall(VecShift(C,sum));
492c4762a1bSJed Brown    } else {
493c4762a1bSJed Brown      for (i=0;i<user->ns;i++) {
4949566063dSJacob Faibussowitsch       PetscCall(Scatter(user->y,user->suby,user->yi_scatter[i],0,0));
4959566063dSJacob Faibussowitsch       PetscCall(Scatter(C,user->subq,user->yi_scatter[i],0,0));
4969566063dSJacob Faibussowitsch       PetscCall(MatMult(user->Grad,user->suby,user->Swork));
4979566063dSJacob Faibussowitsch       PetscCall(VecPointwiseDivide(user->Swork,user->Swork,user->Av_u));
4989566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(user->Grad,user->Swork,user->subq));
499c4762a1bSJed Brown 
5009566063dSJacob Faibussowitsch       PetscCall(VecSum(user->suby,&sum));
501c4762a1bSJed Brown       sum /= user->ndesign;
5029566063dSJacob Faibussowitsch       PetscCall(VecShift(user->subq,sum));
503c4762a1bSJed Brown 
5049566063dSJacob Faibussowitsch       PetscCall(Gather(user->y,user->suby,user->yi_scatter[i],0,0));
5059566063dSJacob Faibussowitsch       PetscCall(Gather(C,user->subq,user->yi_scatter[i],0,0));
506c4762a1bSJed Brown      }
507c4762a1bSJed Brown    }
5089566063dSJacob Faibussowitsch    PetscCall(VecAXPY(C,-1.0,user->q));
509c4762a1bSJed Brown    PetscFunctionReturn(0);
510c4762a1bSJed Brown }
511c4762a1bSJed Brown 
512c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
513c4762a1bSJed Brown {
514c4762a1bSJed Brown   PetscFunctionBegin;
5159566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat1,x,sub1,INSERT_VALUES,SCATTER_FORWARD));
5169566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat1,x,sub1,INSERT_VALUES,SCATTER_FORWARD));
517c4762a1bSJed Brown   if (sub2) {
5189566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat2,x,sub2,INSERT_VALUES,SCATTER_FORWARD));
5199566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat2,x,sub2,INSERT_VALUES,SCATTER_FORWARD));
520c4762a1bSJed Brown   }
521c4762a1bSJed Brown   PetscFunctionReturn(0);
522c4762a1bSJed Brown }
523c4762a1bSJed Brown 
524c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
525c4762a1bSJed Brown {
526c4762a1bSJed Brown   PetscFunctionBegin;
5279566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat1,sub1,x,INSERT_VALUES,SCATTER_REVERSE));
5289566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat1,sub1,x,INSERT_VALUES,SCATTER_REVERSE));
529c4762a1bSJed Brown   if (sub2) {
5309566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scat2,sub2,x,INSERT_VALUES,SCATTER_REVERSE));
5319566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scat2,sub2,x,INSERT_VALUES,SCATTER_REVERSE));
532c4762a1bSJed Brown   }
533c4762a1bSJed Brown   PetscFunctionReturn(0);
534c4762a1bSJed Brown }
535c4762a1bSJed Brown 
536c4762a1bSJed Brown PetscErrorCode EllipticInitialize(AppCtx *user)
537c4762a1bSJed Brown {
538c4762a1bSJed Brown   PetscInt       m,n,i,j,k,l,linear_index,is,js,ks,ls,istart,iend,iblock;
539c4762a1bSJed Brown   Vec            XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork;
540c4762a1bSJed Brown   PetscReal      *x,*y,*z;
541c4762a1bSJed Brown   PetscReal      h,meanut;
542c4762a1bSJed Brown   PetscScalar    hinv,neg_hinv,half = 0.5,sqrt_beta;
543c4762a1bSJed Brown   PetscInt       im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
544c4762a1bSJed Brown   IS             is_alldesign,is_allstate;
545c4762a1bSJed Brown   IS             is_from_d;
546c4762a1bSJed Brown   IS             is_from_y;
547c4762a1bSJed Brown   PetscInt       lo,hi,hi2,lo2,ysubnlocal,dsubnlocal;
548c4762a1bSJed Brown   const PetscInt *ranges, *subranges;
549c4762a1bSJed Brown   PetscMPIInt    size;
550c4762a1bSJed Brown   PetscReal      xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
551c4762a1bSJed Brown   PetscScalar    v,vx,vy,vz;
552c4762a1bSJed Brown   PetscInt       offset,subindex,subvec,nrank,kk;
553c4762a1bSJed Brown 
554c4762a1bSJed Brown   PetscScalar xr[64] = {0.4970,     0.8498,     0.7814,     0.6268,     0.7782,     0.6402,     0.3617,     0.3160,
555c4762a1bSJed Brown                         0.3610,     0.5298,     0.6987,     0.3331,     0.7962,     0.5596,     0.3866,     0.6774,
556c4762a1bSJed Brown                         0.5407,     0.4518,     0.6702,     0.6061,     0.7580,     0.8997,     0.5198,     0.8326,
557c4762a1bSJed Brown                         0.2138,     0.9198,     0.3000,     0.2833,     0.8288,     0.7076,     0.1820,     0.0728,
558c4762a1bSJed Brown                         0.8447,     0.2367,     0.3239,     0.6413,     0.3114,     0.4731,     0.1192,     0.9273,
559c4762a1bSJed Brown                         0.5724,     0.4331,     0.5136,     0.3547,     0.4413,     0.2602,     0.5698,     0.7278,
560c4762a1bSJed Brown                         0.5261,     0.6230,     0.2454,     0.3948,     0.7479,     0.6582,     0.4660,     0.5594,
561c4762a1bSJed Brown                         0.7574,     0.1143,     0.5900,     0.1065,     0.4260,     0.3294,     0.8276,     0.0756};
562c4762a1bSJed Brown 
563c4762a1bSJed Brown   PetscScalar yr[64] = {0.7345,     0.9120,     0.9288,     0.7528,     0.4463,     0.4985,     0.2497,     0.6256,
564c4762a1bSJed Brown                         0.3425,     0.9026,     0.6983,     0.4230,     0.7140,     0.2970,     0.4474,     0.8792,
565c4762a1bSJed Brown                         0.6604,     0.2485,     0.7968,     0.6127,     0.1796,     0.2437,     0.5938,     0.6137,
566c4762a1bSJed Brown                         0.3867,     0.5658,     0.4575,     0.1009,     0.0863,     0.3361,     0.0738,     0.3985,
567c4762a1bSJed Brown                         0.6602,     0.1437,     0.0934,     0.5983,     0.5950,     0.0763,     0.0768,     0.2288,
568c4762a1bSJed Brown                         0.5761,     0.1129,     0.3841,     0.6150,     0.6904,     0.6686,     0.1361,     0.4601,
569c4762a1bSJed Brown                         0.4491,     0.3716,     0.1969,     0.6537,     0.6743,     0.6991,     0.4811,     0.5480,
570c4762a1bSJed Brown                         0.1684,     0.4569,     0.6889,     0.8437,     0.3015,     0.2854,     0.8199,     0.2658};
571c4762a1bSJed Brown 
572c4762a1bSJed Brown   PetscScalar zr[64] = {0.7668,     0.8573,     0.2654,     0.2719,     0.1060,     0.1311,     0.6232,     0.2295,
573c4762a1bSJed Brown                         0.8009,     0.2147,     0.2119,     0.9325,     0.4473,     0.3600,     0.3374,     0.3819,
574c4762a1bSJed Brown                         0.4066,     0.5801,     0.1673,     0.0959,     0.4638,     0.8236,     0.8800,     0.2939,
575c4762a1bSJed Brown                         0.2028,     0.8262,     0.2706,     0.6276,     0.9085,     0.6443,     0.8241,     0.0712,
576c4762a1bSJed Brown                         0.1824,     0.7789,     0.4389,     0.8415,     0.7055,     0.6639,     0.3653,     0.2078,
577c4762a1bSJed Brown                         0.1987,     0.2297,     0.4321,     0.8115,     0.4915,     0.7764,     0.4657,     0.4627,
578c4762a1bSJed Brown                         0.4569,     0.4232,     0.8514,     0.0674,     0.3227,     0.1055,     0.6690,     0.6313,
579c4762a1bSJed Brown                         0.9226,     0.5461,     0.4126,     0.2364,     0.6096,     0.7042,     0.3914,     0.0711};
580c4762a1bSJed Brown 
581c4762a1bSJed Brown   PetscFunctionBegin;
5829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
5839566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Elliptic Setup",&user->stages[0]));
5849566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(user->stages[0]));
585c4762a1bSJed Brown 
586c4762a1bSJed Brown   /* Create u,y,c,x */
5879566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->u));
5889566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->y));
5899566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->c));
5909566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->u,PETSC_DECIDE,user->ndesign));
5919566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->u));
5929566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(user->u,&ysubnlocal));
5939566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->y,ysubnlocal*user->ns,user->nstate));
5949566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->c,ysubnlocal*user->ns,user->m));
5959566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->y));
5969566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->c));
597c4762a1bSJed Brown 
598c4762a1bSJed Brown   /*
599c4762a1bSJed Brown      *******************************
600c4762a1bSJed Brown      Create scatters for x <-> y,u
601c4762a1bSJed Brown      *******************************
602c4762a1bSJed Brown 
603c4762a1bSJed Brown      If the state vector y and design vector u are partitioned as
604c4762a1bSJed Brown      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
605c4762a1bSJed Brown      then the solution vector x is organized as
606c4762a1bSJed Brown      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
607c4762a1bSJed Brown      The index sets user->s_is and user->d_is correspond to the indices of the
608c4762a1bSJed Brown      state and design variables owned by the current processor.
609c4762a1bSJed Brown   */
6109566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->x));
611c4762a1bSJed Brown 
6129566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->y,&lo,&hi));
6139566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->u,&lo2,&hi2));
614c4762a1bSJed Brown 
6159566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate));
6169566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user->s_is));
6179566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign));
6189566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user->d_is));
619c4762a1bSJed Brown 
6209566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->x,hi-lo+hi2-lo2,user->n));
6219566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->x));
622c4762a1bSJed Brown 
6239566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(user->x,user->s_is,user->y,is_allstate,&user->state_scatter));
6249566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(user->x,user->d_is,user->u,is_alldesign,&user->design_scatter));
6259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_alldesign));
6269566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_allstate));
627c4762a1bSJed Brown   /*
628c4762a1bSJed Brown      *******************************
629c4762a1bSJed Brown      Create scatter from y to y_1,y_2,...,y_ns
630c4762a1bSJed Brown      *******************************
631c4762a1bSJed Brown   */
6329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->ns,&user->yi_scatter));
6339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->suby));
6349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->subq));
635c4762a1bSJed Brown 
6369566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->y,&lo2,&hi2));
637c4762a1bSJed Brown   istart = 0;
638c4762a1bSJed Brown   for (i=0; i<user->ns; i++) {
6399566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->suby,&lo,&hi));
6409566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y));
6419566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->y,is_from_y,user->suby,NULL,&user->yi_scatter[i]));
642c4762a1bSJed Brown     istart = istart + hi-lo;
6439566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_y));
644c4762a1bSJed Brown   }
645c4762a1bSJed Brown   /*
646c4762a1bSJed Brown      *******************************
647c4762a1bSJed Brown      Create scatter from d to d_1,d_2,...,d_ns
648c4762a1bSJed Brown      *******************************
649c4762a1bSJed Brown   */
6509566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->subd));
6519566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->subd,PETSC_DECIDE,user->ndata));
6529566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->subd));
6539566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->d));
6549566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(user->subd,&dsubnlocal));
6559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->d,dsubnlocal*user->ns,user->ndata*user->ns));
6569566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->d));
6579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->ns,&user->di_scatter));
658c4762a1bSJed Brown 
6599566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->d,&lo2,&hi2));
660c4762a1bSJed Brown   istart = 0;
661c4762a1bSJed Brown   for (i=0; i<user->ns; i++) {
6629566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(user->subd,&lo,&hi));
6639566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d));
6649566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(user->d,is_from_d,user->subd,NULL,&user->di_scatter[i]));
665c4762a1bSJed Brown     istart = istart + hi-lo;
6669566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is_from_d));
667c4762a1bSJed Brown   }
668c4762a1bSJed Brown 
6699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->mx,&x));
6709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->mx,&y));
6719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->mx,&z));
672c4762a1bSJed Brown 
673c4762a1bSJed Brown   user->ksp_its = 0;
674c4762a1bSJed Brown   user->ksp_its_initial = 0;
675c4762a1bSJed Brown 
676c4762a1bSJed Brown   n = user->mx * user->mx * user->mx;
677c4762a1bSJed Brown   m = 3 * user->mx * user->mx * (user->mx-1);
678c4762a1bSJed Brown   sqrt_beta = PetscSqrtScalar(user->beta);
679c4762a1bSJed Brown 
6809566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&XX));
6819566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->q));
6829566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(XX,ysubnlocal,n));
6839566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->q,ysubnlocal*user->ns,user->m));
6849566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(XX));
6859566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->q));
686c4762a1bSJed Brown 
6879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&YY));
6889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&ZZ));
6899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&XXwork));
6909566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&YYwork));
6919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&ZZwork));
6929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&UTwork));
6939566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(XX,&user->utrue));
694c4762a1bSJed Brown 
695c4762a1bSJed Brown   /* map for striding q */
6969566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRanges(user->q,&ranges));
6979566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRanges(user->u,&subranges));
698c4762a1bSJed Brown 
6999566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->q,&lo2,&hi2));
7009566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(user->u,&lo,&hi));
701c4762a1bSJed Brown   /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
702c4762a1bSJed Brown   h = 1.0/user->mx;
703c4762a1bSJed Brown   hinv = user->mx;
704c4762a1bSJed Brown   neg_hinv = -hinv;
705c4762a1bSJed Brown 
7069566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(XX,&istart,&iend));
707c4762a1bSJed Brown   for (linear_index=istart; linear_index<iend; linear_index++) {
708c4762a1bSJed Brown     i = linear_index % user->mx;
709c4762a1bSJed Brown     j = ((linear_index-i)/user->mx) % user->mx;
710c4762a1bSJed Brown     k = ((linear_index-i)/user->mx-j) / user->mx;
711c4762a1bSJed Brown     vx = h*(i+0.5);
712c4762a1bSJed Brown     vy = h*(j+0.5);
713c4762a1bSJed Brown     vz = h*(k+0.5);
7149566063dSJacob Faibussowitsch     PetscCall(VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES));
7159566063dSJacob Faibussowitsch     PetscCall(VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES));
7169566063dSJacob Faibussowitsch     PetscCall(VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES));
717c4762a1bSJed Brown     for (is=0; is<2; is++) {
718c4762a1bSJed Brown       for (js=0; js<2; js++) {
719c4762a1bSJed Brown         for (ks=0; ks<2; ks++) {
720c4762a1bSJed Brown           ls = is*4 + js*2 + ks;
721c4762a1bSJed Brown           if (ls<user->ns) {
722c4762a1bSJed Brown             l =ls*n + linear_index;
723c4762a1bSJed Brown             /* remap */
724c4762a1bSJed Brown             subindex = l%n;
725c4762a1bSJed Brown             subvec = l/n;
726c4762a1bSJed Brown             nrank=0;
727c4762a1bSJed Brown             while (subindex >= subranges[nrank+1]) nrank++;
728c4762a1bSJed Brown             offset = subindex - subranges[nrank];
729c4762a1bSJed Brown             istart=0;
730c4762a1bSJed Brown             for (kk=0;kk<nrank;kk++) istart+=user->ns*(subranges[kk+1]-subranges[kk]);
731c4762a1bSJed Brown             istart += (subranges[nrank+1]-subranges[nrank])*subvec;
732c4762a1bSJed Brown             l = istart+offset;
733c4762a1bSJed 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));
7349566063dSJacob Faibussowitsch             PetscCall(VecSetValues(user->q,1,&l,&v,INSERT_VALUES));
735c4762a1bSJed Brown           }
736c4762a1bSJed Brown         }
737c4762a1bSJed Brown       }
738c4762a1bSJed Brown     }
739c4762a1bSJed Brown   }
740c4762a1bSJed Brown 
7419566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(XX));
7429566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(XX));
7439566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(YY));
7449566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(YY));
7459566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(ZZ));
7469566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(ZZ));
7479566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(user->q));
7489566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(user->q));
749c4762a1bSJed Brown 
750c4762a1bSJed Brown   /* Compute true parameter function
751c4762a1bSJed 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) */
7529566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX,XXwork));
7539566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY,YYwork));
7549566063dSJacob Faibussowitsch   PetscCall(VecCopy(ZZ,ZZwork));
755c4762a1bSJed Brown 
7569566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork,-0.25));
7579566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork,-0.25));
7589566063dSJacob Faibussowitsch   PetscCall(VecShift(ZZwork,-0.25));
759c4762a1bSJed Brown 
7609566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork,XXwork,XXwork));
7619566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork,YYwork,YYwork));
7629566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(ZZwork,ZZwork,ZZwork));
763c4762a1bSJed Brown 
7649566063dSJacob Faibussowitsch   PetscCall(VecCopy(XXwork,UTwork));
7659566063dSJacob Faibussowitsch   PetscCall(VecAXPY(UTwork,1.0,YYwork));
7669566063dSJacob Faibussowitsch   PetscCall(VecAXPY(UTwork,1.0,ZZwork));
7679566063dSJacob Faibussowitsch   PetscCall(VecScale(UTwork,-20.0));
7689566063dSJacob Faibussowitsch   PetscCall(VecExp(UTwork));
7699566063dSJacob Faibussowitsch   PetscCall(VecCopy(UTwork,user->utrue));
770c4762a1bSJed Brown 
7719566063dSJacob Faibussowitsch   PetscCall(VecCopy(XX,XXwork));
7729566063dSJacob Faibussowitsch   PetscCall(VecCopy(YY,YYwork));
7739566063dSJacob Faibussowitsch   PetscCall(VecCopy(ZZ,ZZwork));
774c4762a1bSJed Brown 
7759566063dSJacob Faibussowitsch   PetscCall(VecShift(XXwork,-0.75));
7769566063dSJacob Faibussowitsch   PetscCall(VecShift(YYwork,-0.75));
7779566063dSJacob Faibussowitsch   PetscCall(VecShift(ZZwork,-0.75));
778c4762a1bSJed Brown 
7799566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(XXwork,XXwork,XXwork));
7809566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(YYwork,YYwork,YYwork));
7819566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(ZZwork,ZZwork,ZZwork));
782c4762a1bSJed Brown 
7839566063dSJacob Faibussowitsch   PetscCall(VecCopy(XXwork,UTwork));
7849566063dSJacob Faibussowitsch   PetscCall(VecAXPY(UTwork,1.0,YYwork));
7859566063dSJacob Faibussowitsch   PetscCall(VecAXPY(UTwork,1.0,ZZwork));
7869566063dSJacob Faibussowitsch   PetscCall(VecScale(UTwork,-20.0));
7879566063dSJacob Faibussowitsch   PetscCall(VecExp(UTwork));
788c4762a1bSJed Brown 
7899566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->utrue,-1.0,UTwork));
790c4762a1bSJed Brown 
7919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XX));
7929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YY));
7939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ZZ));
7949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&XXwork));
7959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&YYwork));
7969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ZZwork));
7979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&UTwork));
798c4762a1bSJed Brown 
799c4762a1bSJed Brown   /* Initial guess and reference model */
8009566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->utrue,&user->ur));
8019566063dSJacob Faibussowitsch   PetscCall(VecSum(user->utrue,&meanut));
802c4762a1bSJed Brown   meanut = meanut / n;
8039566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ur,meanut));
8049566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->ur,user->u));
805c4762a1bSJed Brown 
806c4762a1bSJed Brown   /* Generate Grad matrix */
8079566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Grad));
8089566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Grad,PETSC_DECIDE,ysubnlocal,m,n));
8099566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Grad));
8109566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL));
8119566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Grad,2,NULL));
8129566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Grad,&istart,&iend));
813c4762a1bSJed Brown 
814c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
815c4762a1bSJed Brown     if (i<m/3) {
816c4762a1bSJed Brown       iblock = i / (user->mx-1);
817c4762a1bSJed Brown       j = iblock*user->mx + (i % (user->mx-1));
8189566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES));
819c4762a1bSJed Brown       j = j+1;
8209566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES));
821c4762a1bSJed Brown     }
822c4762a1bSJed Brown     if (i>=m/3 && i<2*m/3) {
823c4762a1bSJed Brown       iblock = (i-m/3) / (user->mx*(user->mx-1));
824c4762a1bSJed Brown       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
8259566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES));
826c4762a1bSJed Brown       j = j + user->mx;
8279566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES));
828c4762a1bSJed Brown     }
829c4762a1bSJed Brown     if (i>=2*m/3) {
830c4762a1bSJed Brown       j = i-2*m/3;
8319566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES));
832c4762a1bSJed Brown       j = j + user->mx*user->mx;
8339566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES));
834c4762a1bSJed Brown     }
835c4762a1bSJed Brown   }
836c4762a1bSJed Brown 
8379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY));
8389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY));
839c4762a1bSJed Brown 
840c4762a1bSJed Brown   /* Generate arithmetic averaging matrix Av */
8419566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Av));
8429566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Av,PETSC_DECIDE,ysubnlocal,m,n));
8439566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Av));
8449566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL));
8459566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Av,2,NULL));
8469566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Av,&istart,&iend));
847c4762a1bSJed Brown 
848c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
849c4762a1bSJed Brown     if (i<m/3) {
850c4762a1bSJed Brown       iblock = i / (user->mx-1);
851c4762a1bSJed Brown       j = iblock*user->mx + (i % (user->mx-1));
8529566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES));
853c4762a1bSJed Brown       j = j+1;
8549566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES));
855c4762a1bSJed Brown     }
856c4762a1bSJed Brown     if (i>=m/3 && i<2*m/3) {
857c4762a1bSJed Brown       iblock = (i-m/3) / (user->mx*(user->mx-1));
858c4762a1bSJed Brown       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
8599566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES));
860c4762a1bSJed Brown       j = j + user->mx;
8619566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES));
862c4762a1bSJed Brown     }
863c4762a1bSJed Brown     if (i>=2*m/3) {
864c4762a1bSJed Brown       j = i-2*m/3;
8659566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES));
866c4762a1bSJed Brown       j = j + user->mx*user->mx;
8679566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES));
868c4762a1bSJed Brown     }
869c4762a1bSJed Brown   }
870c4762a1bSJed Brown 
8719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY));
8729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY));
873c4762a1bSJed Brown 
8749566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->L));
8759566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->L,PETSC_DECIDE,ysubnlocal,m+n,n));
8769566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->L));
8779566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL));
8789566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->L,2,NULL));
8799566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->L,&istart,&iend));
880c4762a1bSJed Brown 
881c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
882c4762a1bSJed Brown     if (i<m/3) {
883c4762a1bSJed Brown       iblock = i / (user->mx-1);
884c4762a1bSJed Brown       j = iblock*user->mx + (i % (user->mx-1));
8859566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES));
886c4762a1bSJed Brown       j = j+1;
8879566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES));
888c4762a1bSJed Brown     }
889c4762a1bSJed Brown     if (i>=m/3 && i<2*m/3) {
890c4762a1bSJed Brown       iblock = (i-m/3) / (user->mx*(user->mx-1));
891c4762a1bSJed Brown       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
8929566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES));
893c4762a1bSJed Brown       j = j + user->mx;
8949566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES));
895c4762a1bSJed Brown     }
896c4762a1bSJed Brown     if (i>=2*m/3 && i<m) {
897c4762a1bSJed Brown       j = i-2*m/3;
8989566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES));
899c4762a1bSJed Brown       j = j + user->mx*user->mx;
9009566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES));
901c4762a1bSJed Brown     }
902c4762a1bSJed Brown     if (i>=m) {
903c4762a1bSJed Brown       j = i - m;
9049566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES));
905c4762a1bSJed Brown     }
906c4762a1bSJed Brown   }
9079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY));
9089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY));
9099566063dSJacob Faibussowitsch   PetscCall(MatScale(user->L,PetscPowScalar(h,1.5)));
910c4762a1bSJed Brown 
911c4762a1bSJed Brown   /* Generate Div matrix */
912c4762a1bSJed Brown   if (!user->use_ptap) {
913c4762a1bSJed Brown     /* Generate Div matrix */
9149566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Div));
9159566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(user->Div,ysubnlocal,PETSC_DECIDE,n,m));
9169566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(user->Div));
9179566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(user->Div,4,NULL,4,NULL));
9189566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(user->Div,6,NULL));
9199566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(user->Grad,&istart,&iend));
920c4762a1bSJed Brown 
921c4762a1bSJed Brown     for (i=istart; i<iend; i++) {
922c4762a1bSJed Brown       if (i<m/3) {
923c4762a1bSJed Brown         iblock = i / (user->mx-1);
924c4762a1bSJed Brown         j = iblock*user->mx + (i % (user->mx-1));
9259566063dSJacob Faibussowitsch         PetscCall(MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES));
926c4762a1bSJed Brown         j = j+1;
9279566063dSJacob Faibussowitsch         PetscCall(MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES));
928c4762a1bSJed Brown       }
929c4762a1bSJed Brown       if (i>=m/3 && i<2*m/3) {
930c4762a1bSJed Brown         iblock = (i-m/3) / (user->mx*(user->mx-1));
931c4762a1bSJed Brown         j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
9329566063dSJacob Faibussowitsch         PetscCall(MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES));
933c4762a1bSJed Brown         j = j + user->mx;
9349566063dSJacob Faibussowitsch         PetscCall(MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES));
935c4762a1bSJed Brown       }
936c4762a1bSJed Brown       if (i>=2*m/3) {
937c4762a1bSJed Brown         j = i-2*m/3;
9389566063dSJacob Faibussowitsch         PetscCall(MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES));
939c4762a1bSJed Brown         j = j + user->mx*user->mx;
9409566063dSJacob Faibussowitsch         PetscCall(MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES));
941c4762a1bSJed Brown       }
942c4762a1bSJed Brown     }
943c4762a1bSJed Brown 
9449566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(user->Div,MAT_FINAL_ASSEMBLY));
9459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(user->Div,MAT_FINAL_ASSEMBLY));
9469566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork));
947c4762a1bSJed Brown   } else {
9489566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Diag));
9499566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(user->Diag,PETSC_DECIDE,PETSC_DECIDE,m,m));
9509566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(user->Diag));
9519566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(user->Diag,1,NULL,0,NULL));
9529566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(user->Diag,1,NULL));
953c4762a1bSJed Brown   }
954c4762a1bSJed Brown 
955c4762a1bSJed Brown   /* Build work vectors and matrices */
9569566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->S));
9579566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->S, PETSC_DECIDE, m));
9589566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->S));
959c4762a1bSJed Brown 
9609566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&user->lwork));
9619566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx));
9629566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->lwork));
963c4762a1bSJed Brown 
9649566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork));
965c4762a1bSJed Brown 
9669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S,&user->Swork));
9679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S,&user->Sdiag));
9689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S,&user->Av_u));
9699566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->S,&user->Twork));
9709566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->y,&user->ywork));
9719566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->Ywork));
9729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->uwork));
9739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->u,&user->js_diag));
9749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->c,&user->cwork));
9759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->d,&user->dwork));
976c4762a1bSJed Brown 
977c4762a1bSJed Brown   /* Create a matrix-free shell user->Jd for computing B*x */
9789566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal,user->nstate,user->ndesign,user,&user->Jd));
9799566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult));
9809566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose));
981c4762a1bSJed Brown 
982c4762a1bSJed Brown   /* Compute true state function ytrue given utrue */
9839566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->y,&user->ytrue));
984c4762a1bSJed Brown 
985c4762a1bSJed Brown   /* First compute Av_u = Av*exp(-u) */
9869566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork, 0));
9879566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork,-1.0,user->utrue)); /* Note: user->utrue */
9889566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
9899566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av,user->uwork,user->Av_u));
990c4762a1bSJed Brown 
991c4762a1bSJed Brown   /* Next form DSG = Div*S*Grad */
9929566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->Av_u,user->Swork));
9939566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Swork));
994c4762a1bSJed Brown   if (user->use_ptap) {
9959566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES));
9969566063dSJacob Faibussowitsch     PetscCall(MatPtAP(user->Diag,user->Grad,MAT_INITIAL_MATRIX,1.0,&user->DSG));
997c4762a1bSJed Brown   } else {
9989566063dSJacob Faibussowitsch     PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN));
9999566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Swork));
1000c20d7725SJed Brown 
10019566063dSJacob Faibussowitsch     PetscCall(MatMatMult(user->Divwork,user->Grad,MAT_INITIAL_MATRIX,1.0,&user->DSG));
1002c4762a1bSJed Brown   }
1003c4762a1bSJed Brown 
10049566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE));
10059566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
1006c4762a1bSJed Brown 
1007c4762a1bSJed Brown   if (user->use_lrc == PETSC_TRUE) {
1008c4762a1bSJed Brown     v=PetscSqrtReal(1.0 /user->ndesign);
10099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(user->ndesign,&user->ones));
1010c4762a1bSJed Brown 
1011c4762a1bSJed Brown     for (i=0;i<user->ndesign;i++) {
1012c4762a1bSJed Brown       user->ones[i]=v;
1013c4762a1bSJed Brown     }
10149566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PETSC_COMM_WORLD,ysubnlocal,PETSC_DECIDE,user->ndesign,1,user->ones,&user->Ones));
10159566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(user->Ones, MAT_FINAL_ASSEMBLY));
10169566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(user->Ones, MAT_FINAL_ASSEMBLY));
10179566063dSJacob Faibussowitsch     PetscCall(MatCreateLRC(user->DSG,user->Ones,NULL,user->Ones,&user->JsBlock));
10189566063dSJacob Faibussowitsch     PetscCall(MatSetUp(user->JsBlock));
1019c4762a1bSJed Brown   } else {
1020c4762a1bSJed Brown     /* Create matrix-free shell user->Js for computing (A + h^3*e*e^T)*x */
10219566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PETSC_COMM_WORLD,ysubnlocal,ysubnlocal,user->ndesign,user->ndesign,user,&user->JsBlock));
10229566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateBlockMatMult));
10239566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateBlockMatMult));
1024c4762a1bSJed Brown   }
10259566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->JsBlock,MAT_SYMMETRIC,PETSC_TRUE));
10269566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->JsBlock,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
10279566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal*user->ns,user->nstate,user->nstate,user,&user->Js));
10289566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult));
10299566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMult));
10309566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->Js,MAT_SYMMETRIC,PETSC_TRUE));
10319566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->Js,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
1032c4762a1bSJed Brown 
10339566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal*user->ns,user->nstate,user->nstate,user,&user->JsInv));
10349566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateInvMatMult));
10359566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateInvMatMult));
10369566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->JsInv,MAT_SYMMETRIC,PETSC_TRUE));
10379566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->JsInv,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
1038c4762a1bSJed Brown 
10399566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE));
10409566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
1041c4762a1bSJed Brown   /* Now solve for ytrue */
10429566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD,&user->solver));
10439566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(user->solver));
1044c4762a1bSJed Brown 
10459566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(user->solver,user->JsBlock,user->DSG));
1046c4762a1bSJed Brown 
10479566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsInv,user->q,user->ytrue));
1048c4762a1bSJed Brown   /* First compute Av_u = Av*exp(-u) */
10499566063dSJacob Faibussowitsch   PetscCall(VecSet(user->uwork,0));
10509566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork,-1.0,user->u)); /* Note: user->u */
10519566063dSJacob Faibussowitsch   PetscCall(VecExp(user->uwork));
10529566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Av,user->uwork,user->Av_u));
1053c4762a1bSJed Brown 
1054c4762a1bSJed Brown   /* Next update DSG = Div*S*Grad  with user->u */
10559566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->Av_u,user->Swork));
10569566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(user->Swork));
1057c4762a1bSJed Brown   if (user->use_ptap) {
10589566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES));
10599566063dSJacob Faibussowitsch     PetscCall(MatPtAP(user->Diag,user->Grad,MAT_REUSE_MATRIX,1.0,&user->DSG));
1060c4762a1bSJed Brown   } else {
10619566063dSJacob Faibussowitsch     PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN));
10629566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Av_u));
10639566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(user->DSG));
1064c4762a1bSJed Brown   }
1065c4762a1bSJed Brown 
1066c4762a1bSJed Brown   /* Now solve for y */
1067c4762a1bSJed Brown 
10689566063dSJacob Faibussowitsch   PetscCall(MatMult(user->JsInv,user->q,user->y));
1069c4762a1bSJed Brown 
1070c4762a1bSJed Brown   user->ksp_its_initial = user->ksp_its;
1071c4762a1bSJed Brown   user->ksp_its = 0;
1072c4762a1bSJed Brown   /* Construct projection matrix Q (blocks) */
10739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Q));
10749566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Q,dsubnlocal,ysubnlocal,user->ndata,user->ndesign));
10759566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Q));
10769566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(user->Q,8,NULL,8,NULL));
10779566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(user->Q,8,NULL));
1078c4762a1bSJed Brown 
1079c4762a1bSJed Brown   for (i=0; i<user->mx; i++) {
1080c4762a1bSJed Brown     x[i] = h*(i+0.5);
1081c4762a1bSJed Brown     y[i] = h*(i+0.5);
1082c4762a1bSJed Brown     z[i] = h*(i+0.5);
1083c4762a1bSJed Brown   }
10849566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(user->Q,&istart,&iend));
1085c4762a1bSJed Brown 
1086c4762a1bSJed Brown   nx = user->mx; ny = user->mx; nz = user->mx;
1087c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
1088c4762a1bSJed Brown 
1089c4762a1bSJed Brown     xri = xr[i];
1090c4762a1bSJed Brown     im = 0;
1091c4762a1bSJed Brown     xim = x[im];
1092c4762a1bSJed Brown     while (xri>xim && im<nx) {
1093c4762a1bSJed Brown       im = im+1;
1094c4762a1bSJed Brown       xim = x[im];
1095c4762a1bSJed Brown     }
1096c4762a1bSJed Brown     indx1 = im-1;
1097c4762a1bSJed Brown     indx2 = im;
1098c4762a1bSJed Brown     dx1 = xri - x[indx1];
1099c4762a1bSJed Brown     dx2 = x[indx2] - xri;
1100c4762a1bSJed Brown 
1101c4762a1bSJed Brown     yri = yr[i];
1102c4762a1bSJed Brown     im = 0;
1103c4762a1bSJed Brown     yim = y[im];
1104c4762a1bSJed Brown     while (yri>yim && im<ny) {
1105c4762a1bSJed Brown       im = im+1;
1106c4762a1bSJed Brown       yim = y[im];
1107c4762a1bSJed Brown     }
1108c4762a1bSJed Brown     indy1 = im-1;
1109c4762a1bSJed Brown     indy2 = im;
1110c4762a1bSJed Brown     dy1 = yri - y[indy1];
1111c4762a1bSJed Brown     dy2 = y[indy2] - yri;
1112c4762a1bSJed Brown 
1113c4762a1bSJed Brown     zri = zr[i];
1114c4762a1bSJed Brown     im = 0;
1115c4762a1bSJed Brown     zim = z[im];
1116c4762a1bSJed Brown     while (zri>zim && im<nz) {
1117c4762a1bSJed Brown       im = im+1;
1118c4762a1bSJed Brown       zim = z[im];
1119c4762a1bSJed Brown     }
1120c4762a1bSJed Brown     indz1 = im-1;
1121c4762a1bSJed Brown     indz2 = im;
1122c4762a1bSJed Brown     dz1 = zri - z[indz1];
1123c4762a1bSJed Brown     dz2 = z[indz2] - zri;
1124c4762a1bSJed Brown 
1125c4762a1bSJed Brown     Dx = x[indx2] - x[indx1];
1126c4762a1bSJed Brown     Dy = y[indy2] - y[indy1];
1127c4762a1bSJed Brown     Dz = z[indz2] - z[indz1];
1128c4762a1bSJed Brown 
1129c4762a1bSJed Brown     j = indx1 + indy1*nx + indz1*nx*ny;
1130c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
11319566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1132c4762a1bSJed Brown 
1133c4762a1bSJed Brown     j = indx1 + indy1*nx + indz2*nx*ny;
1134c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
11359566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1136c4762a1bSJed Brown 
1137c4762a1bSJed Brown     j = indx1 + indy2*nx + indz1*nx*ny;
1138c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
11399566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1140c4762a1bSJed Brown 
1141c4762a1bSJed Brown     j = indx1 + indy2*nx + indz2*nx*ny;
1142c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
11439566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1144c4762a1bSJed Brown 
1145c4762a1bSJed Brown     j = indx2 + indy1*nx + indz1*nx*ny;
1146c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
11479566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1148c4762a1bSJed Brown 
1149c4762a1bSJed Brown     j = indx2 + indy1*nx + indz2*nx*ny;
1150c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
11519566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1152c4762a1bSJed Brown 
1153c4762a1bSJed Brown     j = indx2 + indy2*nx + indz1*nx*ny;
1154c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
11559566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1156c4762a1bSJed Brown 
1157c4762a1bSJed Brown     j = indx2 + indy2*nx + indz2*nx*ny;
1158c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
11599566063dSJacob Faibussowitsch     PetscCall(MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES));
1160c4762a1bSJed Brown   }
1161c4762a1bSJed Brown 
11629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY));
11639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY));
1164c4762a1bSJed Brown   /* Create MQ (composed of blocks of Q */
11659566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,dsubnlocal*user->ns,PETSC_DECIDE,user->ndata*user->ns,user->nstate,user,&user->MQ));
11669566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->MQ,MATOP_MULT,(void(*)(void))QMatMult));
11679566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(user->MQ,MATOP_MULT_TRANSPOSE,(void(*)(void))QMatMultTranspose));
1168c4762a1bSJed Brown 
1169c4762a1bSJed Brown   /* Add noise to the measurement data */
11709566063dSJacob Faibussowitsch   PetscCall(VecSet(user->ywork,1.0));
11719566063dSJacob Faibussowitsch   PetscCall(VecAYPX(user->ywork,user->noise,user->ytrue));
11729566063dSJacob Faibussowitsch   PetscCall(MatMult(user->MQ,user->ywork,user->d));
1173c4762a1bSJed Brown 
1174c4762a1bSJed Brown   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11759566063dSJacob Faibussowitsch   PetscCall(PetscFree(x));
11769566063dSJacob Faibussowitsch   PetscCall(PetscFree(y));
11779566063dSJacob Faibussowitsch   PetscCall(PetscFree(z));
11789566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1179c4762a1bSJed Brown   PetscFunctionReturn(0);
1180c4762a1bSJed Brown }
1181c4762a1bSJed Brown 
1182c4762a1bSJed Brown PetscErrorCode EllipticDestroy(AppCtx *user)
1183c4762a1bSJed Brown {
1184c4762a1bSJed Brown   PetscInt       i;
1185c4762a1bSJed Brown 
1186c4762a1bSJed Brown   PetscFunctionBegin;
11879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->DSG));
11889566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&user->solver));
11899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Q));
11909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->MQ));
1191c4762a1bSJed Brown   if (!user->use_ptap) {
11929566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->Div));
11939566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->Divwork));
1194c4762a1bSJed Brown   } else {
11959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->Diag));
1196c4762a1bSJed Brown   }
1197c4762a1bSJed Brown   if (user->use_lrc) {
11989566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->Ones));
1199c4762a1bSJed Brown   }
1200c4762a1bSJed Brown 
12019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Grad));
12029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Av));
12039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Avwork));
12049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->L));
12059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Js));
12069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Jd));
12079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsBlock));
12089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->JsInv));
1209c4762a1bSJed Brown 
12109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->x));
12119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->u));
12129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->uwork));
12139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->utrue));
12149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->y));
12159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ywork));
12169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ytrue));
12179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->c));
12189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->cwork));
12199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ur));
12209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->q));
12219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->d));
12229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->dwork));
12239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->lwork));
12249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->S));
12259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Swork));
12269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Sdiag));
12279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Ywork));
12289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Twork));
12299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Av_u));
12309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->js_diag));
12319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->s_is));
12329566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user->d_is));
12339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->suby));
12349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->subd));
12359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->subq));
12369566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->state_scatter));
12379566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->design_scatter));
1238c4762a1bSJed Brown   for (i=0;i<user->ns;i++) {
12399566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
12409566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&user->di_scatter[i]));
1241c4762a1bSJed Brown   }
12429566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->yi_scatter));
12439566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->di_scatter));
1244c4762a1bSJed Brown   if (user->use_lrc) {
12459566063dSJacob Faibussowitsch     PetscCall(PetscFree(user->ones));
12469566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user->Ones));
1247c4762a1bSJed Brown   }
1248c4762a1bSJed Brown   PetscFunctionReturn(0);
1249c4762a1bSJed Brown }
1250c4762a1bSJed Brown 
1251c4762a1bSJed Brown PetscErrorCode EllipticMonitor(Tao tao, void *ptr)
1252c4762a1bSJed Brown {
1253c4762a1bSJed Brown   Vec            X;
1254c4762a1bSJed Brown   PetscReal      unorm,ynorm;
1255c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
1256c4762a1bSJed Brown 
1257c4762a1bSJed Brown   PetscFunctionBegin;
12589566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao,&X));
12599566063dSJacob Faibussowitsch   PetscCall(Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
12609566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->ywork,-1.0,user->ytrue));
12619566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->uwork,-1.0,user->utrue));
12629566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->uwork,NORM_2,&unorm));
12639566063dSJacob Faibussowitsch   PetscCall(VecNorm(user->ywork,NORM_2,&ynorm));
12649566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm));
1265c4762a1bSJed Brown   PetscFunctionReturn(0);
1266c4762a1bSJed Brown }
1267c4762a1bSJed Brown 
1268c4762a1bSJed Brown /*TEST
1269c4762a1bSJed Brown 
1270c4762a1bSJed Brown    build:
1271c4762a1bSJed Brown       requires: !complex
1272c4762a1bSJed Brown 
1273c4762a1bSJed Brown    test:
1274c4762a1bSJed Brown       args: -tao_cmonitor -ns 1 -tao_type lcl -tao_gatol 1.e-3 -tao_max_it 11
1275c4762a1bSJed Brown       requires: !single
1276c4762a1bSJed Brown 
1277c4762a1bSJed Brown    test:
1278c4762a1bSJed Brown       suffix: 2
1279c4762a1bSJed Brown       args: -tao_cmonitor -tao_type lcl -tao_max_it 11 -use_ptap -use_lrc -ns 1 -tao_gatol 1.e-3
1280c4762a1bSJed Brown       requires: !single
1281c4762a1bSJed Brown 
1282c4762a1bSJed Brown TEST*/
1283