xref: /petsc/src/tao/pde_constrained/tutorials/elliptic.c (revision a82e8c82ed9474375a7f877f23dfa96948657643)
1c4762a1bSJed Brown #include <petsc/private/taoimpl.h>
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*T
4c4762a1bSJed Brown    Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
5c4762a1bSJed Brown    Routines: TaoCreate();
6c4762a1bSJed Brown    Routines: TaoSetType();
7*a82e8c82SStefano Zampini    Routines: TaoSetSolution();
8*a82e8c82SStefano Zampini    Routines: TaoSetObjective();
9*a82e8c82SStefano Zampini    Routines: TaoSetGradient();
10c4762a1bSJed Brown    Routines: TaoSetConstraintsRoutine();
11c4762a1bSJed Brown    Routines: TaoSetJacobianStateRoutine();
12c4762a1bSJed Brown    Routines: TaoSetJacobianDesignRoutine();
13c4762a1bSJed Brown    Routines: TaoSetStateDesignIS();
14c4762a1bSJed Brown    Routines: TaoSetFromOptions();
15c4762a1bSJed Brown    Routines: TaoSolve();
16c4762a1bSJed Brown    Routines: TaoDestroy();
17c4762a1bSJed Brown    Processors: n
18c4762a1bSJed Brown T*/
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct {
21c4762a1bSJed Brown   PetscInt n; /* Number of total variables */
22c4762a1bSJed Brown   PetscInt m; /* Number of constraints */
23c4762a1bSJed Brown   PetscInt nstate;
24c4762a1bSJed Brown   PetscInt ndesign;
25c4762a1bSJed Brown   PetscInt mx; /* grid points in each direction */
26c4762a1bSJed Brown   PetscInt ns; /* Number of data samples (1<=ns<=8)
27c4762a1bSJed Brown                   Currently only ns=1 is supported */
28c4762a1bSJed Brown   PetscInt ndata; /* Number of data points per sample */
29c4762a1bSJed Brown   IS       s_is;
30c4762a1bSJed Brown   IS       d_is;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   VecScatter state_scatter;
33c4762a1bSJed Brown   VecScatter design_scatter;
34c4762a1bSJed Brown   VecScatter *yi_scatter, *di_scatter;
35c4762a1bSJed Brown   Vec        suby,subq,subd;
36c4762a1bSJed Brown   Mat        Js,Jd,JsPrec,JsInv,JsBlock;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   PetscReal alpha; /* Regularization parameter */
39c4762a1bSJed Brown   PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */
40c4762a1bSJed Brown   PetscReal noise; /* Amount of noise to add to data */
41c4762a1bSJed Brown   PetscReal *ones;
42c4762a1bSJed Brown   Mat       Q;
43c4762a1bSJed Brown   Mat       MQ;
44c4762a1bSJed Brown   Mat       L;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   Mat Grad;
47c4762a1bSJed Brown   Mat Av,Avwork;
48c4762a1bSJed Brown   Mat Div, Divwork;
49c4762a1bSJed Brown   Mat DSG;
50c4762a1bSJed Brown   Mat Diag,Ones;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   Vec q;
53c4762a1bSJed Brown   Vec ur; /* reference */
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   Vec d;
56c4762a1bSJed Brown   Vec dwork;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   Vec x; /* super vec of y,u */
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   Vec y; /* state variables */
61c4762a1bSJed Brown   Vec ywork;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   Vec ytrue;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   Vec u; /* design variables */
66c4762a1bSJed Brown   Vec uwork;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   Vec utrue;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   Vec js_diag;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   Vec c; /* constraint vector */
73c4762a1bSJed Brown   Vec cwork;
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   Vec lwork;
76c4762a1bSJed Brown   Vec S;
77c4762a1bSJed Brown   Vec Swork,Twork,Sdiag,Ywork;
78c4762a1bSJed Brown   Vec Av_u;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   KSP solver;
81c4762a1bSJed Brown   PC  prec;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   PetscReal tola,tolb,tolc,told;
84c4762a1bSJed Brown   PetscInt  ksp_its;
85c4762a1bSJed Brown   PetscInt  ksp_its_initial;
86c4762a1bSJed Brown   PetscLogStage stages[10];
87c4762a1bSJed Brown   PetscBool use_ptap;
88c4762a1bSJed Brown   PetscBool use_lrc;
89c4762a1bSJed Brown } AppCtx;
90c4762a1bSJed Brown 
91c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*);
92c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void*);
93c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*);
94c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void*);
95c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat,void*);
96c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void*);
97c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void*);
98c4762a1bSJed Brown PetscErrorCode Gather(Vec, Vec, VecScatter, Vec, VecScatter);
99c4762a1bSJed Brown PetscErrorCode Scatter(Vec, Vec, VecScatter, Vec, VecScatter);
100c4762a1bSJed Brown PetscErrorCode EllipticInitialize(AppCtx*);
101c4762a1bSJed Brown PetscErrorCode EllipticDestroy(AppCtx*);
102c4762a1bSJed Brown PetscErrorCode EllipticMonitor(Tao, void*);
103c4762a1bSJed Brown 
104c4762a1bSJed Brown PetscErrorCode StateBlockMatMult(Mat,Vec,Vec);
105c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat,Vec,Vec);
106c4762a1bSJed Brown 
107c4762a1bSJed Brown PetscErrorCode StateInvMatMult(Mat,Vec,Vec);
108c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat,Vec,Vec);
109c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);
110c4762a1bSJed Brown 
111c4762a1bSJed Brown PetscErrorCode QMatMult(Mat,Vec,Vec);
112c4762a1bSJed Brown PetscErrorCode QMatMultTranspose(Mat,Vec,Vec);
113c4762a1bSJed Brown 
114c4762a1bSJed Brown static  char help[]="";
115c4762a1bSJed Brown 
116c4762a1bSJed Brown int main(int argc, char **argv)
117c4762a1bSJed Brown {
118c4762a1bSJed Brown   PetscErrorCode     ierr;
119c4762a1bSJed Brown   Vec                x0;
120c4762a1bSJed Brown   Tao                tao;
121c4762a1bSJed Brown   AppCtx             user;
122c4762a1bSJed Brown   PetscInt           ntests = 1;
123c4762a1bSJed Brown   PetscInt           i;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, (char*)0,help);if (ierr) return ierr;
126c4762a1bSJed Brown   user.mx = 8;
12776280437SVaclav Hapla   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"elliptic example",NULL);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL);CHKERRQ(ierr);
129c4762a1bSJed Brown   user.ns = 6;
130c4762a1bSJed Brown   ierr = PetscOptionsInt("-ns","Number of data samples (1<=ns<=8)","",user.ns,&user.ns,NULL);CHKERRQ(ierr);
131c4762a1bSJed Brown   user.ndata = 64;
132c4762a1bSJed Brown   ierr = PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL);CHKERRQ(ierr);
133c4762a1bSJed Brown   user.alpha = 0.1;
134c4762a1bSJed Brown   ierr = PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL);CHKERRQ(ierr);
135c4762a1bSJed Brown   user.beta = 0.00001;
136c4762a1bSJed Brown   ierr = PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,NULL);CHKERRQ(ierr);
137c4762a1bSJed Brown   user.noise = 0.01;
138c4762a1bSJed Brown   ierr = PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,NULL);CHKERRQ(ierr);
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   user.use_ptap = PETSC_FALSE;
141c4762a1bSJed Brown   ierr = PetscOptionsBool("-use_ptap","Use ptap matrix for DSG","",user.use_ptap,&user.use_ptap,NULL);CHKERRQ(ierr);
142c4762a1bSJed Brown   user.use_lrc = PETSC_FALSE;
143c4762a1bSJed Brown   ierr = PetscOptionsBool("-use_lrc","Use lrc matrix for Js","",user.use_lrc,&user.use_lrc,NULL);CHKERRQ(ierr);
14476280437SVaclav Hapla   ierr = PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL);CHKERRQ(ierr);
14576280437SVaclav Hapla   ierr = PetscOptionsEnd();CHKERRQ(ierr);
14676280437SVaclav Hapla 
147c4762a1bSJed Brown   user.m = user.ns*user.mx*user.mx*user.mx; /* number of constraints */
148c4762a1bSJed Brown   user.nstate =  user.m;
149c4762a1bSJed Brown   user.ndesign = user.mx*user.mx*user.mx;
150c4762a1bSJed Brown   user.n = user.nstate + user.ndesign; /* number of variables */
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
153c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
154c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOLCL);CHKERRQ(ierr);
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Set up initial vectors and matrices */
157c4762a1bSJed Brown   ierr = EllipticInitialize(&user);CHKERRQ(ierr);
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   ierr = Gather(user.x,user.y,user.state_scatter,user.u,user.design_scatter);CHKERRQ(ierr);
160c4762a1bSJed Brown   ierr = VecDuplicate(user.x,&x0);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = VecCopy(user.x,x0);CHKERRQ(ierr);
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Set solution vector with an initial guess */
164*a82e8c82SStefano Zampini   ierr = TaoSetSolution(tao,user.x);CHKERRQ(ierr);
165*a82e8c82SStefano Zampini   ierr = TaoSetObjective(tao, FormFunction, &user);CHKERRQ(ierr);
166*a82e8c82SStefano Zampini   ierr = TaoSetGradient(tao, NULL, FormGradient, &user);CHKERRQ(ierr);
1673ec1f749SStefano Zampini   ierr = TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user);CHKERRQ(ierr);
168c4762a1bSJed Brown 
1693ec1f749SStefano Zampini   ierr = TaoSetJacobianStateRoutine(tao, user.Js, NULL, user.JsInv, FormJacobianState, &user);CHKERRQ(ierr);
1703ec1f749SStefano Zampini   ierr = TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user);CHKERRQ(ierr);
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   ierr = TaoSetStateDesignIS(tao,user.s_is,user.d_is);CHKERRQ(ierr);
173c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
176c4762a1bSJed Brown   ierr = PetscLogStageRegister("Trials",&user.stages[1]);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = PetscLogStagePush(user.stages[1]);CHKERRQ(ierr);
178c4762a1bSJed Brown   for (i=0; i<ntests; i++) {
179c4762a1bSJed Brown     ierr = TaoSolve(tao);CHKERRQ(ierr);
180c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its);CHKERRQ(ierr);
181c4762a1bSJed Brown     ierr = VecCopy(x0,user.x);CHKERRQ(ierr);
182c4762a1bSJed Brown   }
183c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = PetscBarrier((PetscObject)user.x);CHKERRQ(ierr);
185c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ");CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial);CHKERRQ(ierr);
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
189c4762a1bSJed Brown   ierr = VecDestroy(&x0);CHKERRQ(ierr);
190c4762a1bSJed Brown   ierr = EllipticDestroy(&user);CHKERRQ(ierr);
191c4762a1bSJed Brown   ierr = PetscFinalize();
192c4762a1bSJed Brown   return ierr;
193c4762a1bSJed Brown }
194c4762a1bSJed Brown /* ------------------------------------------------------------------- */
195c4762a1bSJed Brown /*
196c4762a1bSJed Brown    dwork = Qy - d
197c4762a1bSJed Brown    lwork = L*(u-ur)
198c4762a1bSJed Brown    f = 1/2 * (dwork.dwork + alpha*lwork.lwork)
199c4762a1bSJed Brown */
200c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   PetscErrorCode ierr;
203c4762a1bSJed Brown   PetscReal      d1=0,d2=0;
204c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   PetscFunctionBegin;
207c4762a1bSJed Brown   ierr = Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);CHKERRQ(ierr);
208c4762a1bSJed Brown   ierr = MatMult(user->MQ,user->y,user->dwork);CHKERRQ(ierr);
209c4762a1bSJed Brown   ierr = VecAXPY(user->dwork,-1.0,user->d);CHKERRQ(ierr);
210c4762a1bSJed Brown   ierr = VecDot(user->dwork,user->dwork,&d1);CHKERRQ(ierr);
211c4762a1bSJed Brown   ierr = VecWAXPY(user->uwork,-1.0,user->ur,user->u);CHKERRQ(ierr);
212c4762a1bSJed Brown   ierr = MatMult(user->L,user->uwork,user->lwork);CHKERRQ(ierr);
213c4762a1bSJed Brown   ierr = VecDot(user->lwork,user->lwork,&d2);CHKERRQ(ierr);
214c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
215c4762a1bSJed Brown   PetscFunctionReturn(0);
216c4762a1bSJed Brown }
217c4762a1bSJed Brown 
218c4762a1bSJed Brown /* ------------------------------------------------------------------- */
219c4762a1bSJed Brown /*
220c4762a1bSJed Brown     state: g_s = Q' *(Qy - d)
221c4762a1bSJed Brown     design: g_d = alpha*L'*L*(u-ur)
222c4762a1bSJed Brown */
223c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
224c4762a1bSJed Brown {
225c4762a1bSJed Brown   PetscErrorCode ierr;
226c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   PetscFunctionBegin;
229c4762a1bSJed Brown   ierr = Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);CHKERRQ(ierr);
230c4762a1bSJed Brown   ierr = MatMult(user->MQ,user->y,user->dwork);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = VecAXPY(user->dwork,-1.0,user->d);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = MatMultTranspose(user->MQ,user->dwork,user->ywork);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = VecWAXPY(user->uwork,-1.0,user->ur,user->u);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = MatMult(user->L,user->uwork,user->lwork);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = MatMultTranspose(user->L,user->lwork,user->uwork);CHKERRQ(ierr);
236c4762a1bSJed Brown   ierr = VecScale(user->uwork, user->alpha);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);CHKERRQ(ierr);
238c4762a1bSJed Brown   PetscFunctionReturn(0);
239c4762a1bSJed Brown }
240c4762a1bSJed Brown 
241c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
242c4762a1bSJed Brown {
243c4762a1bSJed Brown   PetscErrorCode ierr;
244c4762a1bSJed Brown   PetscReal      d1,d2;
245c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   PetscFunctionBegin;
248c4762a1bSJed Brown   ierr = Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = MatMult(user->MQ,user->y,user->dwork);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = VecAXPY(user->dwork,-1.0,user->d);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = VecDot(user->dwork,user->dwork,&d1);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = MatMultTranspose(user->MQ,user->dwork,user->ywork);CHKERRQ(ierr);
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   ierr = VecWAXPY(user->uwork,-1.0,user->ur,user->u);CHKERRQ(ierr);
255c4762a1bSJed Brown   ierr = MatMult(user->L,user->uwork,user->lwork);CHKERRQ(ierr);
256c4762a1bSJed Brown   ierr = VecDot(user->lwork,user->lwork,&d2);CHKERRQ(ierr);
257c4762a1bSJed Brown   ierr = MatMultTranspose(user->L,user->lwork,user->uwork);CHKERRQ(ierr);
258c4762a1bSJed Brown   ierr = VecScale(user->uwork, user->alpha);CHKERRQ(ierr);
259c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
260c4762a1bSJed Brown   ierr = Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);CHKERRQ(ierr);
261c4762a1bSJed Brown   PetscFunctionReturn(0);
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
264c4762a1bSJed Brown /* ------------------------------------------------------------------- */
265c4762a1bSJed Brown /* A
266c4762a1bSJed Brown MatShell object
267c4762a1bSJed Brown */
268c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
269c4762a1bSJed Brown {
270c4762a1bSJed Brown   PetscErrorCode ierr;
271c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   PetscFunctionBegin;
274c4762a1bSJed Brown   ierr = Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);CHKERRQ(ierr);
275c4762a1bSJed Brown   /* DSG = Div * (1/Av_u) * Grad */
276c4762a1bSJed Brown   ierr = VecSet(user->uwork,0);CHKERRQ(ierr);
277c4762a1bSJed Brown   ierr = VecAXPY(user->uwork,-1.0,user->u);CHKERRQ(ierr);
278c4762a1bSJed Brown   ierr = VecExp(user->uwork);CHKERRQ(ierr);
279c4762a1bSJed Brown   ierr = MatMult(user->Av,user->uwork,user->Av_u);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr = VecCopy(user->Av_u,user->Swork);CHKERRQ(ierr);
281c4762a1bSJed Brown   ierr = VecReciprocal(user->Swork);CHKERRQ(ierr);
282c4762a1bSJed Brown   if (user->use_ptap) {
283c4762a1bSJed Brown     ierr = MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);CHKERRQ(ierr);
284c4762a1bSJed Brown     ierr = MatPtAP(user->Diag,user->Grad,MAT_REUSE_MATRIX,1.0,&user->DSG);CHKERRQ(ierr);
285c4762a1bSJed Brown   } else {
286c4762a1bSJed Brown     ierr = MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
287c4762a1bSJed Brown     ierr = MatDiagonalScale(user->Divwork,NULL,user->Swork);CHKERRQ(ierr);
288c20d7725SJed Brown     ierr = MatProductNumeric(user->DSG);CHKERRQ(ierr);
289c4762a1bSJed Brown   }
290c4762a1bSJed Brown   PetscFunctionReturn(0);
291c4762a1bSJed Brown }
292c4762a1bSJed Brown /* ------------------------------------------------------------------- */
293c4762a1bSJed Brown /* B */
294c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
295c4762a1bSJed Brown {
296c4762a1bSJed Brown   PetscErrorCode ierr;
297c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   PetscFunctionBegin;
300c4762a1bSJed Brown   ierr = Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);CHKERRQ(ierr);
301c4762a1bSJed Brown   PetscFunctionReturn(0);
302c4762a1bSJed Brown }
303c4762a1bSJed Brown 
304c4762a1bSJed Brown PetscErrorCode StateBlockMatMult(Mat J_shell, Vec X, Vec Y)
305c4762a1bSJed Brown {
306c4762a1bSJed Brown   PetscErrorCode ierr;
307c4762a1bSJed Brown   PetscReal      sum;
308c4762a1bSJed Brown   AppCtx         *user;
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   PetscFunctionBegin;
3113ec1f749SStefano Zampini   ierr = MatShellGetContext(J_shell,&user);CHKERRQ(ierr);
312c4762a1bSJed Brown   ierr = MatMult(user->DSG,X,Y);CHKERRQ(ierr);
313c4762a1bSJed Brown   ierr = VecSum(X,&sum);CHKERRQ(ierr);
314c4762a1bSJed Brown   sum /= user->ndesign;
315c4762a1bSJed Brown   ierr = VecShift(Y,sum);CHKERRQ(ierr);
316c4762a1bSJed Brown   PetscFunctionReturn(0);
317c4762a1bSJed Brown }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
320c4762a1bSJed Brown {
321c4762a1bSJed Brown   PetscErrorCode ierr;
322c4762a1bSJed Brown   PetscInt       i;
323c4762a1bSJed Brown   AppCtx         *user;
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   PetscFunctionBegin;
3263ec1f749SStefano Zampini   ierr = MatShellGetContext(J_shell,&user);CHKERRQ(ierr);
327c4762a1bSJed Brown   if (user->ns == 1) {
328c4762a1bSJed Brown     ierr = MatMult(user->JsBlock,X,Y);CHKERRQ(ierr);
329c4762a1bSJed Brown   } else {
330c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
331c4762a1bSJed Brown       ierr = Scatter(X,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
332c4762a1bSJed Brown       ierr = Scatter(Y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
333c4762a1bSJed Brown       ierr = MatMult(user->JsBlock,user->subq,user->suby);CHKERRQ(ierr);
334c4762a1bSJed Brown       ierr = Gather(Y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
335c4762a1bSJed Brown     }
336c4762a1bSJed Brown   }
337c4762a1bSJed Brown   PetscFunctionReturn(0);
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown PetscErrorCode StateInvMatMult(Mat J_shell, Vec X, Vec Y)
341c4762a1bSJed Brown {
342c4762a1bSJed Brown   PetscErrorCode ierr;
343c4762a1bSJed Brown   PetscInt       its,i;
344c4762a1bSJed Brown   AppCtx         *user;
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   PetscFunctionBegin;
3473ec1f749SStefano Zampini   ierr = MatShellGetContext(J_shell,&user);CHKERRQ(ierr);
348c4762a1bSJed Brown   ierr = KSPSetOperators(user->solver,user->JsBlock,user->DSG);CHKERRQ(ierr);
349c4762a1bSJed Brown   if (Y == user->ytrue) {
350c4762a1bSJed Brown     /* First solve is done using true solution to set up problem */
351c4762a1bSJed Brown     ierr = KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
352c4762a1bSJed Brown   } else {
353c4762a1bSJed Brown     ierr = KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
354c4762a1bSJed Brown   }
355c4762a1bSJed Brown   if (user->ns == 1) {
356c4762a1bSJed Brown     ierr = KSPSolve(user->solver,X,Y);CHKERRQ(ierr);
357c4762a1bSJed Brown     ierr = KSPGetIterationNumber(user->solver,&its);CHKERRQ(ierr);
358c4762a1bSJed Brown     user->ksp_its+=its;
359c4762a1bSJed Brown   } else {
360c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
361c4762a1bSJed Brown       ierr = Scatter(X,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
362c4762a1bSJed Brown       ierr = Scatter(Y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
363c4762a1bSJed Brown       ierr = KSPSolve(user->solver,user->subq,user->suby);CHKERRQ(ierr);
364c4762a1bSJed Brown       ierr = KSPGetIterationNumber(user->solver,&its);CHKERRQ(ierr);
365c4762a1bSJed Brown       user->ksp_its+=its;
366c4762a1bSJed Brown       ierr = Gather(Y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
367c4762a1bSJed Brown     }
368c4762a1bSJed Brown   }
369c4762a1bSJed Brown   PetscFunctionReturn(0);
370c4762a1bSJed Brown }
371c4762a1bSJed Brown PetscErrorCode QMatMult(Mat J_shell, Vec X, Vec Y)
372c4762a1bSJed Brown {
373c4762a1bSJed Brown   PetscErrorCode ierr;
374c4762a1bSJed Brown   AppCtx         *user;
375c4762a1bSJed Brown   PetscInt       i;
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   PetscFunctionBegin;
3783ec1f749SStefano Zampini   ierr = MatShellGetContext(J_shell,&user);CHKERRQ(ierr);
379c4762a1bSJed Brown   if (user->ns == 1) {
380c4762a1bSJed Brown     ierr = MatMult(user->Q,X,Y);CHKERRQ(ierr);
381c4762a1bSJed Brown   } else {
382c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
383c4762a1bSJed Brown       ierr = Scatter(X,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
384c4762a1bSJed Brown       ierr = Scatter(Y,user->subd,user->di_scatter[i],0,0);CHKERRQ(ierr);
385c4762a1bSJed Brown       ierr = MatMult(user->Q,user->subq,user->subd);CHKERRQ(ierr);
386c4762a1bSJed Brown       ierr = Gather(Y,user->subd,user->di_scatter[i],0,0);CHKERRQ(ierr);
387c4762a1bSJed Brown     }
388c4762a1bSJed Brown   }
389c4762a1bSJed Brown   PetscFunctionReturn(0);
390c4762a1bSJed Brown }
391c4762a1bSJed Brown 
392c4762a1bSJed Brown PetscErrorCode QMatMultTranspose(Mat J_shell, Vec X, Vec Y)
393c4762a1bSJed Brown {
394c4762a1bSJed Brown   PetscErrorCode ierr;
395c4762a1bSJed Brown   AppCtx         *user;
396c4762a1bSJed Brown   PetscInt       i;
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   PetscFunctionBegin;
3993ec1f749SStefano Zampini   ierr = MatShellGetContext(J_shell,&user);CHKERRQ(ierr);
400c4762a1bSJed Brown   if (user->ns == 1) {
401c4762a1bSJed Brown     ierr = MatMultTranspose(user->Q,X,Y);CHKERRQ(ierr);
402c4762a1bSJed Brown   } else {
403c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
404c4762a1bSJed Brown       ierr = Scatter(X,user->subd,user->di_scatter[i],0,0);CHKERRQ(ierr);
405c4762a1bSJed Brown       ierr = Scatter(Y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
406c4762a1bSJed Brown       ierr = MatMultTranspose(user->Q,user->subd,user->suby);CHKERRQ(ierr);
407c4762a1bSJed Brown       ierr = Gather(Y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
408c4762a1bSJed Brown     }
409c4762a1bSJed Brown   }
410c4762a1bSJed Brown   PetscFunctionReturn(0);
411c4762a1bSJed Brown }
412c4762a1bSJed Brown 
413c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
414c4762a1bSJed Brown {
415c4762a1bSJed Brown   PetscErrorCode ierr;
416c4762a1bSJed Brown   PetscInt       i;
417c4762a1bSJed Brown   AppCtx         *user;
418c4762a1bSJed Brown 
419c4762a1bSJed Brown   PetscFunctionBegin;
4203ec1f749SStefano Zampini   ierr = MatShellGetContext(J_shell,&user);CHKERRQ(ierr);
421c4762a1bSJed Brown 
422c4762a1bSJed Brown   /* sdiag(1./v) */
423c4762a1bSJed Brown   ierr = VecSet(user->uwork,0);CHKERRQ(ierr);
424c4762a1bSJed Brown   ierr = VecAXPY(user->uwork,-1.0,user->u);CHKERRQ(ierr);
425c4762a1bSJed Brown   ierr = VecExp(user->uwork);CHKERRQ(ierr);
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   /* sdiag(1./((Av*(1./v)).^2)) */
428c4762a1bSJed Brown   ierr = MatMult(user->Av,user->uwork,user->Swork);CHKERRQ(ierr);
429c4762a1bSJed Brown   ierr = VecPointwiseMult(user->Swork,user->Swork,user->Swork);CHKERRQ(ierr);
430c4762a1bSJed Brown   ierr = VecReciprocal(user->Swork);CHKERRQ(ierr);
431c4762a1bSJed Brown 
432c4762a1bSJed Brown   /* (Av * (sdiag(1./v) * b)) */
433c4762a1bSJed Brown   ierr = VecPointwiseMult(user->uwork,user->uwork,X);CHKERRQ(ierr);
434c4762a1bSJed Brown   ierr = MatMult(user->Av,user->uwork,user->Twork);CHKERRQ(ierr);
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
437c4762a1bSJed Brown   ierr = VecPointwiseMult(user->Swork,user->Twork,user->Swork);CHKERRQ(ierr);
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   if (user->ns == 1) {
440c4762a1bSJed Brown     /* (sdiag(Grad*y(:,i)) */
441c4762a1bSJed Brown     ierr = MatMult(user->Grad,user->y,user->Twork);CHKERRQ(ierr);
442c4762a1bSJed Brown 
443c4762a1bSJed Brown     /* Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
444c4762a1bSJed Brown     ierr = VecPointwiseMult(user->Swork,user->Twork,user->Swork);CHKERRQ(ierr);
445c4762a1bSJed Brown     ierr = MatMultTranspose(user->Grad,user->Swork,Y);CHKERRQ(ierr);
446c4762a1bSJed Brown   } else {
447c4762a1bSJed Brown     for (i=0;i<user->ns;i++) {
448c4762a1bSJed Brown       ierr = Scatter(user->y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
449c4762a1bSJed Brown       ierr = Scatter(Y,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
450c4762a1bSJed Brown 
451c4762a1bSJed Brown       ierr = MatMult(user->Grad,user->suby,user->Twork);CHKERRQ(ierr);
452c4762a1bSJed Brown       ierr = VecPointwiseMult(user->Twork,user->Twork,user->Swork);CHKERRQ(ierr);
453c4762a1bSJed Brown       ierr = MatMultTranspose(user->Grad,user->Twork,user->subq);CHKERRQ(ierr);
454c4762a1bSJed Brown       ierr = Gather(user->y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
455c4762a1bSJed Brown       ierr = Gather(Y,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
456c4762a1bSJed Brown     }
457c4762a1bSJed Brown   }
458c4762a1bSJed Brown   PetscFunctionReturn(0);
459c4762a1bSJed Brown }
460c4762a1bSJed Brown 
461c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
462c4762a1bSJed Brown {
463c4762a1bSJed Brown   PetscErrorCode ierr;
464c4762a1bSJed Brown   PetscInt       i;
465c4762a1bSJed Brown   AppCtx         *user;
466c4762a1bSJed Brown 
467c4762a1bSJed Brown   PetscFunctionBegin;
4683ec1f749SStefano Zampini   ierr = MatShellGetContext(J_shell,&user);CHKERRQ(ierr);
469c4762a1bSJed Brown   ierr = VecZeroEntries(Y);CHKERRQ(ierr);
470c4762a1bSJed Brown 
471c4762a1bSJed Brown   /* Sdiag = 1./((Av*(1./v)).^2) */
472c4762a1bSJed Brown   ierr = VecSet(user->uwork,0);CHKERRQ(ierr);
473c4762a1bSJed Brown   ierr = VecAXPY(user->uwork,-1.0,user->u);CHKERRQ(ierr);
474c4762a1bSJed Brown   ierr = VecExp(user->uwork);CHKERRQ(ierr);
475c4762a1bSJed Brown   ierr = MatMult(user->Av,user->uwork,user->Swork);CHKERRQ(ierr);
476c4762a1bSJed Brown   ierr = VecPointwiseMult(user->Sdiag,user->Swork,user->Swork);CHKERRQ(ierr);
477c4762a1bSJed Brown   ierr = VecReciprocal(user->Sdiag);CHKERRQ(ierr);
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   for (i=0;i<user->ns;i++) {
480c4762a1bSJed Brown     ierr = Scatter(X,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
481c4762a1bSJed Brown     ierr = Scatter(user->y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
482c4762a1bSJed Brown 
483c4762a1bSJed Brown     /* Swork = (Div' * b(:,i)) */
484c4762a1bSJed Brown     ierr = MatMult(user->Grad,user->subq,user->Swork);CHKERRQ(ierr);
485c4762a1bSJed Brown 
486c4762a1bSJed Brown     /* Twork = Grad*y(:,i) */
487c4762a1bSJed Brown     ierr = MatMult(user->Grad,user->suby,user->Twork);CHKERRQ(ierr);
488c4762a1bSJed Brown 
489c4762a1bSJed Brown     /* Twork = sdiag(Twork) * Swork */
490c4762a1bSJed Brown     ierr = VecPointwiseMult(user->Twork,user->Swork,user->Twork);CHKERRQ(ierr);
491c4762a1bSJed Brown 
492c4762a1bSJed Brown     /* Swork = pointwisemult(Sdiag,Twork) */
493c4762a1bSJed Brown     ierr = VecPointwiseMult(user->Swork,user->Twork,user->Sdiag);CHKERRQ(ierr);
494c4762a1bSJed Brown 
495c4762a1bSJed Brown     /* Ywork = Av' * Swork */
496c4762a1bSJed Brown     ierr = MatMultTranspose(user->Av,user->Swork,user->Ywork);CHKERRQ(ierr);
497c4762a1bSJed Brown 
498c4762a1bSJed Brown     /* Ywork = pointwisemult(uwork,Ywork) */
499c4762a1bSJed Brown     ierr = VecPointwiseMult(user->Ywork,user->uwork,user->Ywork);CHKERRQ(ierr);
500c4762a1bSJed Brown     ierr = VecAXPY(Y,1.0,user->Ywork);CHKERRQ(ierr);
501c4762a1bSJed Brown     ierr = Gather(user->y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
502c4762a1bSJed Brown   }
503c4762a1bSJed Brown   PetscFunctionReturn(0);
504c4762a1bSJed Brown }
505c4762a1bSJed Brown 
506c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
507c4762a1bSJed Brown {
508c4762a1bSJed Brown    /* C=Ay - q      A = Div * Sigma * Grad + hx*hx*hx*ones(n,n) */
509c4762a1bSJed Brown    PetscErrorCode ierr;
510c4762a1bSJed Brown    PetscReal      sum;
511c4762a1bSJed Brown    PetscInt       i;
512c4762a1bSJed Brown    AppCtx         *user = (AppCtx*)ptr;
513c4762a1bSJed Brown 
514c4762a1bSJed Brown    PetscFunctionBegin;
515c4762a1bSJed Brown    ierr = Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);CHKERRQ(ierr);
516c4762a1bSJed Brown    if (user->ns == 1) {
517c4762a1bSJed Brown      ierr = MatMult(user->Grad,user->y,user->Swork);CHKERRQ(ierr);
518c4762a1bSJed Brown      ierr = VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);CHKERRQ(ierr);
519c4762a1bSJed Brown      ierr = MatMultTranspose(user->Grad,user->Swork,C);CHKERRQ(ierr);
520c4762a1bSJed Brown      ierr = VecSum(user->y,&sum);CHKERRQ(ierr);
521c4762a1bSJed Brown      sum /= user->ndesign;
522c4762a1bSJed Brown      ierr = VecShift(C,sum);CHKERRQ(ierr);
523c4762a1bSJed Brown    } else {
524c4762a1bSJed Brown      for (i=0;i<user->ns;i++) {
525c4762a1bSJed Brown       ierr = Scatter(user->y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
526c4762a1bSJed Brown       ierr = Scatter(C,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
527c4762a1bSJed Brown       ierr = MatMult(user->Grad,user->suby,user->Swork);CHKERRQ(ierr);
528c4762a1bSJed Brown       ierr = VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);CHKERRQ(ierr);
529c4762a1bSJed Brown       ierr = MatMultTranspose(user->Grad,user->Swork,user->subq);CHKERRQ(ierr);
530c4762a1bSJed Brown 
531c4762a1bSJed Brown       ierr = VecSum(user->suby,&sum);CHKERRQ(ierr);
532c4762a1bSJed Brown       sum /= user->ndesign;
533c4762a1bSJed Brown       ierr = VecShift(user->subq,sum);CHKERRQ(ierr);
534c4762a1bSJed Brown 
535c4762a1bSJed Brown       ierr = Gather(user->y,user->suby,user->yi_scatter[i],0,0);CHKERRQ(ierr);
536c4762a1bSJed Brown       ierr = Gather(C,user->subq,user->yi_scatter[i],0,0);CHKERRQ(ierr);
537c4762a1bSJed Brown      }
538c4762a1bSJed Brown    }
539c4762a1bSJed Brown    ierr = VecAXPY(C,-1.0,user->q);CHKERRQ(ierr);
540c4762a1bSJed Brown    PetscFunctionReturn(0);
541c4762a1bSJed Brown }
542c4762a1bSJed Brown 
543c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
544c4762a1bSJed Brown {
545c4762a1bSJed Brown   PetscErrorCode ierr;
546c4762a1bSJed Brown 
547c4762a1bSJed Brown   PetscFunctionBegin;
548c4762a1bSJed Brown   ierr = VecScatterBegin(scat1,x,sub1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
549c4762a1bSJed Brown   ierr = VecScatterEnd(scat1,x,sub1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
550c4762a1bSJed Brown   if (sub2) {
551c4762a1bSJed Brown     ierr = VecScatterBegin(scat2,x,sub2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
552c4762a1bSJed Brown     ierr = VecScatterEnd(scat2,x,sub2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
553c4762a1bSJed Brown   }
554c4762a1bSJed Brown   PetscFunctionReturn(0);
555c4762a1bSJed Brown }
556c4762a1bSJed Brown 
557c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
558c4762a1bSJed Brown {
559c4762a1bSJed Brown   PetscErrorCode ierr;
560c4762a1bSJed Brown 
561c4762a1bSJed Brown   PetscFunctionBegin;
562c4762a1bSJed Brown   ierr = VecScatterBegin(scat1,sub1,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
563c4762a1bSJed Brown   ierr = VecScatterEnd(scat1,sub1,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
564c4762a1bSJed Brown   if (sub2) {
565c4762a1bSJed Brown     ierr = VecScatterBegin(scat2,sub2,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
566c4762a1bSJed Brown     ierr = VecScatterEnd(scat2,sub2,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
567c4762a1bSJed Brown   }
568c4762a1bSJed Brown   PetscFunctionReturn(0);
569c4762a1bSJed Brown }
570c4762a1bSJed Brown 
571c4762a1bSJed Brown PetscErrorCode EllipticInitialize(AppCtx *user)
572c4762a1bSJed Brown {
573c4762a1bSJed Brown   PetscErrorCode ierr;
574c4762a1bSJed Brown   PetscInt       m,n,i,j,k,l,linear_index,is,js,ks,ls,istart,iend,iblock;
575c4762a1bSJed Brown   Vec            XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork;
576c4762a1bSJed Brown   PetscReal      *x,*y,*z;
577c4762a1bSJed Brown   PetscReal      h,meanut;
578c4762a1bSJed Brown   PetscScalar    hinv,neg_hinv,half = 0.5,sqrt_beta;
579c4762a1bSJed Brown   PetscInt       im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
580c4762a1bSJed Brown   IS             is_alldesign,is_allstate;
581c4762a1bSJed Brown   IS             is_from_d;
582c4762a1bSJed Brown   IS             is_from_y;
583c4762a1bSJed Brown   PetscInt       lo,hi,hi2,lo2,ysubnlocal,dsubnlocal;
584c4762a1bSJed Brown   const PetscInt *ranges, *subranges;
585c4762a1bSJed Brown   PetscMPIInt    size;
586c4762a1bSJed Brown   PetscReal      xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
587c4762a1bSJed Brown   PetscScalar    v,vx,vy,vz;
588c4762a1bSJed Brown   PetscInt       offset,subindex,subvec,nrank,kk;
589c4762a1bSJed Brown 
590c4762a1bSJed Brown   PetscScalar xr[64] = {0.4970,     0.8498,     0.7814,     0.6268,     0.7782,     0.6402,     0.3617,     0.3160,
591c4762a1bSJed Brown                         0.3610,     0.5298,     0.6987,     0.3331,     0.7962,     0.5596,     0.3866,     0.6774,
592c4762a1bSJed Brown                         0.5407,     0.4518,     0.6702,     0.6061,     0.7580,     0.8997,     0.5198,     0.8326,
593c4762a1bSJed Brown                         0.2138,     0.9198,     0.3000,     0.2833,     0.8288,     0.7076,     0.1820,     0.0728,
594c4762a1bSJed Brown                         0.8447,     0.2367,     0.3239,     0.6413,     0.3114,     0.4731,     0.1192,     0.9273,
595c4762a1bSJed Brown                         0.5724,     0.4331,     0.5136,     0.3547,     0.4413,     0.2602,     0.5698,     0.7278,
596c4762a1bSJed Brown                         0.5261,     0.6230,     0.2454,     0.3948,     0.7479,     0.6582,     0.4660,     0.5594,
597c4762a1bSJed Brown                         0.7574,     0.1143,     0.5900,     0.1065,     0.4260,     0.3294,     0.8276,     0.0756};
598c4762a1bSJed Brown 
599c4762a1bSJed Brown   PetscScalar yr[64] = {0.7345,     0.9120,     0.9288,     0.7528,     0.4463,     0.4985,     0.2497,     0.6256,
600c4762a1bSJed Brown                         0.3425,     0.9026,     0.6983,     0.4230,     0.7140,     0.2970,     0.4474,     0.8792,
601c4762a1bSJed Brown                         0.6604,     0.2485,     0.7968,     0.6127,     0.1796,     0.2437,     0.5938,     0.6137,
602c4762a1bSJed Brown                         0.3867,     0.5658,     0.4575,     0.1009,     0.0863,     0.3361,     0.0738,     0.3985,
603c4762a1bSJed Brown                         0.6602,     0.1437,     0.0934,     0.5983,     0.5950,     0.0763,     0.0768,     0.2288,
604c4762a1bSJed Brown                         0.5761,     0.1129,     0.3841,     0.6150,     0.6904,     0.6686,     0.1361,     0.4601,
605c4762a1bSJed Brown                         0.4491,     0.3716,     0.1969,     0.6537,     0.6743,     0.6991,     0.4811,     0.5480,
606c4762a1bSJed Brown                         0.1684,     0.4569,     0.6889,     0.8437,     0.3015,     0.2854,     0.8199,     0.2658};
607c4762a1bSJed Brown 
608c4762a1bSJed Brown   PetscScalar zr[64] = {0.7668,     0.8573,     0.2654,     0.2719,     0.1060,     0.1311,     0.6232,     0.2295,
609c4762a1bSJed Brown                         0.8009,     0.2147,     0.2119,     0.9325,     0.4473,     0.3600,     0.3374,     0.3819,
610c4762a1bSJed Brown                         0.4066,     0.5801,     0.1673,     0.0959,     0.4638,     0.8236,     0.8800,     0.2939,
611c4762a1bSJed Brown                         0.2028,     0.8262,     0.2706,     0.6276,     0.9085,     0.6443,     0.8241,     0.0712,
612c4762a1bSJed Brown                         0.1824,     0.7789,     0.4389,     0.8415,     0.7055,     0.6639,     0.3653,     0.2078,
613c4762a1bSJed Brown                         0.1987,     0.2297,     0.4321,     0.8115,     0.4915,     0.7764,     0.4657,     0.4627,
614c4762a1bSJed Brown                         0.4569,     0.4232,     0.8514,     0.0674,     0.3227,     0.1055,     0.6690,     0.6313,
615c4762a1bSJed Brown                         0.9226,     0.5461,     0.4126,     0.2364,     0.6096,     0.7042,     0.3914,     0.0711};
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   PetscFunctionBegin;
618ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
619c4762a1bSJed Brown   ierr = PetscLogStageRegister("Elliptic Setup",&user->stages[0]);CHKERRQ(ierr);
620c4762a1bSJed Brown   ierr = PetscLogStagePush(user->stages[0]);CHKERRQ(ierr);
621c4762a1bSJed Brown 
622c4762a1bSJed Brown   /* Create u,y,c,x */
623c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->u);CHKERRQ(ierr);
624c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->y);CHKERRQ(ierr);
625c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->c);CHKERRQ(ierr);
626c4762a1bSJed Brown   ierr = VecSetSizes(user->u,PETSC_DECIDE,user->ndesign);CHKERRQ(ierr);
627c4762a1bSJed Brown   ierr = VecSetFromOptions(user->u);CHKERRQ(ierr);
628c4762a1bSJed Brown   ierr = VecGetLocalSize(user->u,&ysubnlocal);CHKERRQ(ierr);
629c4762a1bSJed Brown   ierr = VecSetSizes(user->y,ysubnlocal*user->ns,user->nstate);CHKERRQ(ierr);
630c4762a1bSJed Brown   ierr = VecSetSizes(user->c,ysubnlocal*user->ns,user->m);CHKERRQ(ierr);
631c4762a1bSJed Brown   ierr = VecSetFromOptions(user->y);CHKERRQ(ierr);
632c4762a1bSJed Brown   ierr = VecSetFromOptions(user->c);CHKERRQ(ierr);
633c4762a1bSJed Brown 
634c4762a1bSJed Brown   /*
635c4762a1bSJed Brown      *******************************
636c4762a1bSJed Brown      Create scatters for x <-> y,u
637c4762a1bSJed Brown      *******************************
638c4762a1bSJed Brown 
639c4762a1bSJed Brown      If the state vector y and design vector u are partitioned as
640c4762a1bSJed Brown      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
641c4762a1bSJed Brown      then the solution vector x is organized as
642c4762a1bSJed Brown      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
643c4762a1bSJed Brown      The index sets user->s_is and user->d_is correspond to the indices of the
644c4762a1bSJed Brown      state and design variables owned by the current processor.
645c4762a1bSJed Brown   */
646c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->x);CHKERRQ(ierr);
647c4762a1bSJed Brown 
648c4762a1bSJed Brown   ierr = VecGetOwnershipRange(user->y,&lo,&hi);CHKERRQ(ierr);
649c4762a1bSJed Brown   ierr = VecGetOwnershipRange(user->u,&lo2,&hi2);CHKERRQ(ierr);
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate);CHKERRQ(ierr);
652c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user->s_is);CHKERRQ(ierr);
653c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign);CHKERRQ(ierr);
654c4762a1bSJed Brown   ierr = ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user->d_is);CHKERRQ(ierr);
655c4762a1bSJed Brown 
656c4762a1bSJed Brown   ierr = VecSetSizes(user->x,hi-lo+hi2-lo2,user->n);CHKERRQ(ierr);
657c4762a1bSJed Brown   ierr = VecSetFromOptions(user->x);CHKERRQ(ierr);
658c4762a1bSJed Brown 
659c4762a1bSJed Brown   ierr = VecScatterCreate(user->x,user->s_is,user->y,is_allstate,&user->state_scatter);CHKERRQ(ierr);
660c4762a1bSJed Brown   ierr = VecScatterCreate(user->x,user->d_is,user->u,is_alldesign,&user->design_scatter);CHKERRQ(ierr);
661c4762a1bSJed Brown   ierr = ISDestroy(&is_alldesign);CHKERRQ(ierr);
662c4762a1bSJed Brown   ierr = ISDestroy(&is_allstate);CHKERRQ(ierr);
663c4762a1bSJed Brown   /*
664c4762a1bSJed Brown      *******************************
665c4762a1bSJed Brown      Create scatter from y to y_1,y_2,...,y_ns
666c4762a1bSJed Brown      *******************************
667c4762a1bSJed Brown   */
668c4762a1bSJed Brown   ierr = PetscMalloc1(user->ns,&user->yi_scatter);CHKERRQ(ierr);
669c4762a1bSJed Brown   ierr = VecDuplicate(user->u,&user->suby);CHKERRQ(ierr);
670c4762a1bSJed Brown   ierr = VecDuplicate(user->u,&user->subq);CHKERRQ(ierr);
671c4762a1bSJed Brown 
672c4762a1bSJed Brown   ierr = VecGetOwnershipRange(user->y,&lo2,&hi2);CHKERRQ(ierr);
673c4762a1bSJed Brown   istart = 0;
674c4762a1bSJed Brown   for (i=0; i<user->ns; i++) {
675c4762a1bSJed Brown     ierr = VecGetOwnershipRange(user->suby,&lo,&hi);CHKERRQ(ierr);
676c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);CHKERRQ(ierr);
677c4762a1bSJed Brown     ierr = VecScatterCreate(user->y,is_from_y,user->suby,NULL,&user->yi_scatter[i]);CHKERRQ(ierr);
678c4762a1bSJed Brown     istart = istart + hi-lo;
679c4762a1bSJed Brown     ierr = ISDestroy(&is_from_y);CHKERRQ(ierr);
680c4762a1bSJed Brown   }
681c4762a1bSJed Brown   /*
682c4762a1bSJed Brown      *******************************
683c4762a1bSJed Brown      Create scatter from d to d_1,d_2,...,d_ns
684c4762a1bSJed Brown      *******************************
685c4762a1bSJed Brown   */
686c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->subd);CHKERRQ(ierr);
687c4762a1bSJed Brown   ierr = VecSetSizes(user->subd,PETSC_DECIDE,user->ndata);CHKERRQ(ierr);
688c4762a1bSJed Brown   ierr = VecSetFromOptions(user->subd);CHKERRQ(ierr);
689c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->d);CHKERRQ(ierr);
690c4762a1bSJed Brown   ierr = VecGetLocalSize(user->subd,&dsubnlocal);CHKERRQ(ierr);
691c4762a1bSJed Brown   ierr = VecSetSizes(user->d,dsubnlocal*user->ns,user->ndata*user->ns);CHKERRQ(ierr);
692c4762a1bSJed Brown   ierr = VecSetFromOptions(user->d);CHKERRQ(ierr);
693c4762a1bSJed Brown   ierr = PetscMalloc1(user->ns,&user->di_scatter);CHKERRQ(ierr);
694c4762a1bSJed Brown 
695c4762a1bSJed Brown   ierr = VecGetOwnershipRange(user->d,&lo2,&hi2);CHKERRQ(ierr);
696c4762a1bSJed Brown   istart = 0;
697c4762a1bSJed Brown   for (i=0; i<user->ns; i++) {
698c4762a1bSJed Brown     ierr = VecGetOwnershipRange(user->subd,&lo,&hi);CHKERRQ(ierr);
699c4762a1bSJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d);CHKERRQ(ierr);
700c4762a1bSJed Brown     ierr = VecScatterCreate(user->d,is_from_d,user->subd,NULL,&user->di_scatter[i]);CHKERRQ(ierr);
701c4762a1bSJed Brown     istart = istart + hi-lo;
702c4762a1bSJed Brown     ierr = ISDestroy(&is_from_d);CHKERRQ(ierr);
703c4762a1bSJed Brown   }
704c4762a1bSJed Brown 
705c4762a1bSJed Brown   ierr = PetscMalloc1(user->mx,&x);CHKERRQ(ierr);
706c4762a1bSJed Brown   ierr = PetscMalloc1(user->mx,&y);CHKERRQ(ierr);
707c4762a1bSJed Brown   ierr = PetscMalloc1(user->mx,&z);CHKERRQ(ierr);
708c4762a1bSJed Brown 
709c4762a1bSJed Brown   user->ksp_its = 0;
710c4762a1bSJed Brown   user->ksp_its_initial = 0;
711c4762a1bSJed Brown 
712c4762a1bSJed Brown   n = user->mx * user->mx * user->mx;
713c4762a1bSJed Brown   m = 3 * user->mx * user->mx * (user->mx-1);
714c4762a1bSJed Brown   sqrt_beta = PetscSqrtScalar(user->beta);
715c4762a1bSJed Brown 
716c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&XX);CHKERRQ(ierr);
717c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->q);CHKERRQ(ierr);
718c4762a1bSJed Brown   ierr = VecSetSizes(XX,ysubnlocal,n);CHKERRQ(ierr);
719c4762a1bSJed Brown   ierr = VecSetSizes(user->q,ysubnlocal*user->ns,user->m);CHKERRQ(ierr);
720c4762a1bSJed Brown   ierr = VecSetFromOptions(XX);CHKERRQ(ierr);
721c4762a1bSJed Brown   ierr = VecSetFromOptions(user->q);CHKERRQ(ierr);
722c4762a1bSJed Brown 
723c4762a1bSJed Brown   ierr = VecDuplicate(XX,&YY);CHKERRQ(ierr);
724c4762a1bSJed Brown   ierr = VecDuplicate(XX,&ZZ);CHKERRQ(ierr);
725c4762a1bSJed Brown   ierr = VecDuplicate(XX,&XXwork);CHKERRQ(ierr);
726c4762a1bSJed Brown   ierr = VecDuplicate(XX,&YYwork);CHKERRQ(ierr);
727c4762a1bSJed Brown   ierr = VecDuplicate(XX,&ZZwork);CHKERRQ(ierr);
728c4762a1bSJed Brown   ierr = VecDuplicate(XX,&UTwork);CHKERRQ(ierr);
729c4762a1bSJed Brown   ierr = VecDuplicate(XX,&user->utrue);CHKERRQ(ierr);
730c4762a1bSJed Brown 
731c4762a1bSJed Brown   /* map for striding q */
732c4762a1bSJed Brown   ierr = VecGetOwnershipRanges(user->q,&ranges);CHKERRQ(ierr);
733c4762a1bSJed Brown   ierr = VecGetOwnershipRanges(user->u,&subranges);CHKERRQ(ierr);
734c4762a1bSJed Brown 
735c4762a1bSJed Brown   ierr = VecGetOwnershipRange(user->q,&lo2,&hi2);CHKERRQ(ierr);
736c4762a1bSJed Brown   ierr = VecGetOwnershipRange(user->u,&lo,&hi);CHKERRQ(ierr);
737c4762a1bSJed Brown   /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
738c4762a1bSJed Brown   h = 1.0/user->mx;
739c4762a1bSJed Brown   hinv = user->mx;
740c4762a1bSJed Brown   neg_hinv = -hinv;
741c4762a1bSJed Brown 
742c4762a1bSJed Brown   ierr = VecGetOwnershipRange(XX,&istart,&iend);CHKERRQ(ierr);
743c4762a1bSJed Brown   for (linear_index=istart; linear_index<iend; linear_index++) {
744c4762a1bSJed Brown     i = linear_index % user->mx;
745c4762a1bSJed Brown     j = ((linear_index-i)/user->mx) % user->mx;
746c4762a1bSJed Brown     k = ((linear_index-i)/user->mx-j) / user->mx;
747c4762a1bSJed Brown     vx = h*(i+0.5);
748c4762a1bSJed Brown     vy = h*(j+0.5);
749c4762a1bSJed Brown     vz = h*(k+0.5);
750c4762a1bSJed Brown     ierr = VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);CHKERRQ(ierr);
751c4762a1bSJed Brown     ierr = VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);CHKERRQ(ierr);
752c4762a1bSJed Brown     ierr = VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES);CHKERRQ(ierr);
753c4762a1bSJed Brown     for (is=0; is<2; is++) {
754c4762a1bSJed Brown       for (js=0; js<2; js++) {
755c4762a1bSJed Brown         for (ks=0; ks<2; ks++) {
756c4762a1bSJed Brown           ls = is*4 + js*2 + ks;
757c4762a1bSJed Brown           if (ls<user->ns) {
758c4762a1bSJed Brown             l =ls*n + linear_index;
759c4762a1bSJed Brown             /* remap */
760c4762a1bSJed Brown             subindex = l%n;
761c4762a1bSJed Brown             subvec = l/n;
762c4762a1bSJed Brown             nrank=0;
763c4762a1bSJed Brown             while (subindex >= subranges[nrank+1]) nrank++;
764c4762a1bSJed Brown             offset = subindex - subranges[nrank];
765c4762a1bSJed Brown             istart=0;
766c4762a1bSJed Brown             for (kk=0;kk<nrank;kk++) istart+=user->ns*(subranges[kk+1]-subranges[kk]);
767c4762a1bSJed Brown             istart += (subranges[nrank+1]-subranges[nrank])*subvec;
768c4762a1bSJed Brown             l = istart+offset;
769c4762a1bSJed 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));
770c4762a1bSJed Brown             ierr = VecSetValues(user->q,1,&l,&v,INSERT_VALUES);CHKERRQ(ierr);
771c4762a1bSJed Brown           }
772c4762a1bSJed Brown         }
773c4762a1bSJed Brown       }
774c4762a1bSJed Brown     }
775c4762a1bSJed Brown   }
776c4762a1bSJed Brown 
777c4762a1bSJed Brown   ierr = VecAssemblyBegin(XX);CHKERRQ(ierr);
778c4762a1bSJed Brown   ierr = VecAssemblyEnd(XX);CHKERRQ(ierr);
779c4762a1bSJed Brown   ierr = VecAssemblyBegin(YY);CHKERRQ(ierr);
780c4762a1bSJed Brown   ierr = VecAssemblyEnd(YY);CHKERRQ(ierr);
781c4762a1bSJed Brown   ierr = VecAssemblyBegin(ZZ);CHKERRQ(ierr);
782c4762a1bSJed Brown   ierr = VecAssemblyEnd(ZZ);CHKERRQ(ierr);
783c4762a1bSJed Brown   ierr = VecAssemblyBegin(user->q);CHKERRQ(ierr);
784c4762a1bSJed Brown   ierr = VecAssemblyEnd(user->q);CHKERRQ(ierr);
785c4762a1bSJed Brown 
786c4762a1bSJed Brown   /* Compute true parameter function
787c4762a1bSJed 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) */
788c4762a1bSJed Brown   ierr = VecCopy(XX,XXwork);CHKERRQ(ierr);
789c4762a1bSJed Brown   ierr = VecCopy(YY,YYwork);CHKERRQ(ierr);
790c4762a1bSJed Brown   ierr = VecCopy(ZZ,ZZwork);CHKERRQ(ierr);
791c4762a1bSJed Brown 
792c4762a1bSJed Brown   ierr = VecShift(XXwork,-0.25);CHKERRQ(ierr);
793c4762a1bSJed Brown   ierr = VecShift(YYwork,-0.25);CHKERRQ(ierr);
794c4762a1bSJed Brown   ierr = VecShift(ZZwork,-0.25);CHKERRQ(ierr);
795c4762a1bSJed Brown 
796c4762a1bSJed Brown   ierr = VecPointwiseMult(XXwork,XXwork,XXwork);CHKERRQ(ierr);
797c4762a1bSJed Brown   ierr = VecPointwiseMult(YYwork,YYwork,YYwork);CHKERRQ(ierr);
798c4762a1bSJed Brown   ierr = VecPointwiseMult(ZZwork,ZZwork,ZZwork);CHKERRQ(ierr);
799c4762a1bSJed Brown 
800c4762a1bSJed Brown   ierr = VecCopy(XXwork,UTwork);CHKERRQ(ierr);
801c4762a1bSJed Brown   ierr = VecAXPY(UTwork,1.0,YYwork);CHKERRQ(ierr);
802c4762a1bSJed Brown   ierr = VecAXPY(UTwork,1.0,ZZwork);CHKERRQ(ierr);
803c4762a1bSJed Brown   ierr = VecScale(UTwork,-20.0);CHKERRQ(ierr);
804c4762a1bSJed Brown   ierr = VecExp(UTwork);CHKERRQ(ierr);
805c4762a1bSJed Brown   ierr = VecCopy(UTwork,user->utrue);CHKERRQ(ierr);
806c4762a1bSJed Brown 
807c4762a1bSJed Brown   ierr = VecCopy(XX,XXwork);CHKERRQ(ierr);
808c4762a1bSJed Brown   ierr = VecCopy(YY,YYwork);CHKERRQ(ierr);
809c4762a1bSJed Brown   ierr = VecCopy(ZZ,ZZwork);CHKERRQ(ierr);
810c4762a1bSJed Brown 
811c4762a1bSJed Brown   ierr = VecShift(XXwork,-0.75);CHKERRQ(ierr);
812c4762a1bSJed Brown   ierr = VecShift(YYwork,-0.75);CHKERRQ(ierr);
813c4762a1bSJed Brown   ierr = VecShift(ZZwork,-0.75);CHKERRQ(ierr);
814c4762a1bSJed Brown 
815c4762a1bSJed Brown   ierr = VecPointwiseMult(XXwork,XXwork,XXwork);CHKERRQ(ierr);
816c4762a1bSJed Brown   ierr = VecPointwiseMult(YYwork,YYwork,YYwork);CHKERRQ(ierr);
817c4762a1bSJed Brown   ierr = VecPointwiseMult(ZZwork,ZZwork,ZZwork);CHKERRQ(ierr);
818c4762a1bSJed Brown 
819c4762a1bSJed Brown   ierr = VecCopy(XXwork,UTwork);CHKERRQ(ierr);
820c4762a1bSJed Brown   ierr = VecAXPY(UTwork,1.0,YYwork);CHKERRQ(ierr);
821c4762a1bSJed Brown   ierr = VecAXPY(UTwork,1.0,ZZwork);CHKERRQ(ierr);
822c4762a1bSJed Brown   ierr = VecScale(UTwork,-20.0);CHKERRQ(ierr);
823c4762a1bSJed Brown   ierr = VecExp(UTwork);CHKERRQ(ierr);
824c4762a1bSJed Brown 
825c4762a1bSJed Brown   ierr = VecAXPY(user->utrue,-1.0,UTwork);CHKERRQ(ierr);
826c4762a1bSJed Brown 
827c4762a1bSJed Brown   ierr = VecDestroy(&XX);CHKERRQ(ierr);
828c4762a1bSJed Brown   ierr = VecDestroy(&YY);CHKERRQ(ierr);
829c4762a1bSJed Brown   ierr = VecDestroy(&ZZ);CHKERRQ(ierr);
830c4762a1bSJed Brown   ierr = VecDestroy(&XXwork);CHKERRQ(ierr);
831c4762a1bSJed Brown   ierr = VecDestroy(&YYwork);CHKERRQ(ierr);
832c4762a1bSJed Brown   ierr = VecDestroy(&ZZwork);CHKERRQ(ierr);
833c4762a1bSJed Brown   ierr = VecDestroy(&UTwork);CHKERRQ(ierr);
834c4762a1bSJed Brown 
835c4762a1bSJed Brown   /* Initial guess and reference model */
836c4762a1bSJed Brown   ierr = VecDuplicate(user->utrue,&user->ur);CHKERRQ(ierr);
837c4762a1bSJed Brown   ierr = VecSum(user->utrue,&meanut);CHKERRQ(ierr);
838c4762a1bSJed Brown   meanut = meanut / n;
839c4762a1bSJed Brown   ierr = VecSet(user->ur,meanut);CHKERRQ(ierr);
840c4762a1bSJed Brown   ierr = VecCopy(user->ur,user->u);CHKERRQ(ierr);
841c4762a1bSJed Brown 
842c4762a1bSJed Brown   /* Generate Grad matrix */
843c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user->Grad);CHKERRQ(ierr);
844c4762a1bSJed Brown   ierr = MatSetSizes(user->Grad,PETSC_DECIDE,ysubnlocal,m,n);CHKERRQ(ierr);
845c4762a1bSJed Brown   ierr = MatSetFromOptions(user->Grad);CHKERRQ(ierr);
846c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);CHKERRQ(ierr);
847c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(user->Grad,2,NULL);CHKERRQ(ierr);
848c4762a1bSJed Brown   ierr = MatGetOwnershipRange(user->Grad,&istart,&iend);CHKERRQ(ierr);
849c4762a1bSJed Brown 
850c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
851c4762a1bSJed Brown     if (i<m/3) {
852c4762a1bSJed Brown       iblock = i / (user->mx-1);
853c4762a1bSJed Brown       j = iblock*user->mx + (i % (user->mx-1));
854c4762a1bSJed Brown       ierr = MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
855c4762a1bSJed Brown       j = j+1;
856c4762a1bSJed Brown       ierr = MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);CHKERRQ(ierr);
857c4762a1bSJed Brown     }
858c4762a1bSJed Brown     if (i>=m/3 && i<2*m/3) {
859c4762a1bSJed Brown       iblock = (i-m/3) / (user->mx*(user->mx-1));
860c4762a1bSJed Brown       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
861c4762a1bSJed Brown       ierr = MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
862c4762a1bSJed Brown       j = j + user->mx;
863c4762a1bSJed Brown       ierr = MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);CHKERRQ(ierr);
864c4762a1bSJed Brown     }
865c4762a1bSJed Brown     if (i>=2*m/3) {
866c4762a1bSJed Brown       j = i-2*m/3;
867c4762a1bSJed Brown       ierr = MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
868c4762a1bSJed Brown       j = j + user->mx*user->mx;
869c4762a1bSJed Brown       ierr = MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);CHKERRQ(ierr);
870c4762a1bSJed Brown     }
871c4762a1bSJed Brown   }
872c4762a1bSJed Brown 
873c4762a1bSJed Brown   ierr = MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
874c4762a1bSJed Brown   ierr = MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
875c4762a1bSJed Brown 
876c4762a1bSJed Brown   /* Generate arithmetic averaging matrix Av */
877c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user->Av);CHKERRQ(ierr);
878c4762a1bSJed Brown   ierr = MatSetSizes(user->Av,PETSC_DECIDE,ysubnlocal,m,n);CHKERRQ(ierr);
879c4762a1bSJed Brown   ierr = MatSetFromOptions(user->Av);CHKERRQ(ierr);
880c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);CHKERRQ(ierr);
881c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(user->Av,2,NULL);CHKERRQ(ierr);
882c4762a1bSJed Brown   ierr = MatGetOwnershipRange(user->Av,&istart,&iend);CHKERRQ(ierr);
883c4762a1bSJed Brown 
884c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
885c4762a1bSJed Brown     if (i<m/3) {
886c4762a1bSJed Brown       iblock = i / (user->mx-1);
887c4762a1bSJed Brown       j = iblock*user->mx + (i % (user->mx-1));
888c4762a1bSJed Brown       ierr = MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);CHKERRQ(ierr);
889c4762a1bSJed Brown       j = j+1;
890c4762a1bSJed Brown       ierr = MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);CHKERRQ(ierr);
891c4762a1bSJed Brown     }
892c4762a1bSJed Brown     if (i>=m/3 && i<2*m/3) {
893c4762a1bSJed Brown       iblock = (i-m/3) / (user->mx*(user->mx-1));
894c4762a1bSJed Brown       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
895c4762a1bSJed Brown       ierr = MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);CHKERRQ(ierr);
896c4762a1bSJed Brown       j = j + user->mx;
897c4762a1bSJed Brown       ierr = MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);CHKERRQ(ierr);
898c4762a1bSJed Brown     }
899c4762a1bSJed Brown     if (i>=2*m/3) {
900c4762a1bSJed Brown       j = i-2*m/3;
901c4762a1bSJed Brown       ierr = MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);CHKERRQ(ierr);
902c4762a1bSJed Brown       j = j + user->mx*user->mx;
903c4762a1bSJed Brown       ierr = MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);CHKERRQ(ierr);
904c4762a1bSJed Brown     }
905c4762a1bSJed Brown   }
906c4762a1bSJed Brown 
907c4762a1bSJed Brown   ierr = MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
908c4762a1bSJed Brown   ierr = MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
909c4762a1bSJed Brown 
910c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user->L);CHKERRQ(ierr);
911c4762a1bSJed Brown   ierr = MatSetSizes(user->L,PETSC_DECIDE,ysubnlocal,m+n,n);CHKERRQ(ierr);
912c4762a1bSJed Brown   ierr = MatSetFromOptions(user->L);CHKERRQ(ierr);
913c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);CHKERRQ(ierr);
914c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(user->L,2,NULL);CHKERRQ(ierr);
915c4762a1bSJed Brown   ierr = MatGetOwnershipRange(user->L,&istart,&iend);CHKERRQ(ierr);
916c4762a1bSJed Brown 
917c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
918c4762a1bSJed Brown     if (i<m/3) {
919c4762a1bSJed Brown       iblock = i / (user->mx-1);
920c4762a1bSJed Brown       j = iblock*user->mx + (i % (user->mx-1));
921c4762a1bSJed Brown       ierr = MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
922c4762a1bSJed Brown       j = j+1;
923c4762a1bSJed Brown       ierr = MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);CHKERRQ(ierr);
924c4762a1bSJed Brown     }
925c4762a1bSJed Brown     if (i>=m/3 && i<2*m/3) {
926c4762a1bSJed Brown       iblock = (i-m/3) / (user->mx*(user->mx-1));
927c4762a1bSJed Brown       j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
928c4762a1bSJed Brown       ierr = MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
929c4762a1bSJed Brown       j = j + user->mx;
930c4762a1bSJed Brown       ierr = MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);CHKERRQ(ierr);
931c4762a1bSJed Brown     }
932c4762a1bSJed Brown     if (i>=2*m/3 && i<m) {
933c4762a1bSJed Brown       j = i-2*m/3;
934c4762a1bSJed Brown       ierr = MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
935c4762a1bSJed Brown       j = j + user->mx*user->mx;
936c4762a1bSJed Brown       ierr = MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);CHKERRQ(ierr);
937c4762a1bSJed Brown     }
938c4762a1bSJed Brown     if (i>=m) {
939c4762a1bSJed Brown       j = i - m;
940c4762a1bSJed Brown       ierr = MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES);CHKERRQ(ierr);
941c4762a1bSJed Brown     }
942c4762a1bSJed Brown   }
943c4762a1bSJed Brown   ierr = MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
944c4762a1bSJed Brown   ierr = MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
945c4762a1bSJed Brown   ierr = MatScale(user->L,PetscPowScalar(h,1.5));CHKERRQ(ierr);
946c4762a1bSJed Brown 
947c4762a1bSJed Brown   /* Generate Div matrix */
948c4762a1bSJed Brown   if (!user->use_ptap) {
949c4762a1bSJed Brown     /* Generate Div matrix */
950c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&user->Div);CHKERRQ(ierr);
951c4762a1bSJed Brown     ierr = MatSetSizes(user->Div,ysubnlocal,PETSC_DECIDE,n,m);CHKERRQ(ierr);
952c4762a1bSJed Brown     ierr = MatSetFromOptions(user->Div);CHKERRQ(ierr);
953c4762a1bSJed Brown     ierr = MatMPIAIJSetPreallocation(user->Div,4,NULL,4,NULL);CHKERRQ(ierr);
954c4762a1bSJed Brown     ierr = MatSeqAIJSetPreallocation(user->Div,6,NULL);CHKERRQ(ierr);
955c4762a1bSJed Brown     ierr = MatGetOwnershipRange(user->Grad,&istart,&iend);CHKERRQ(ierr);
956c4762a1bSJed Brown 
957c4762a1bSJed Brown     for (i=istart; i<iend; i++) {
958c4762a1bSJed Brown       if (i<m/3) {
959c4762a1bSJed Brown         iblock = i / (user->mx-1);
960c4762a1bSJed Brown         j = iblock*user->mx + (i % (user->mx-1));
961c4762a1bSJed Brown         ierr = MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
962c4762a1bSJed Brown         j = j+1;
963c4762a1bSJed Brown         ierr = MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES);CHKERRQ(ierr);
964c4762a1bSJed Brown       }
965c4762a1bSJed Brown       if (i>=m/3 && i<2*m/3) {
966c4762a1bSJed Brown         iblock = (i-m/3) / (user->mx*(user->mx-1));
967c4762a1bSJed Brown         j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
968c4762a1bSJed Brown         ierr = MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
969c4762a1bSJed Brown         j = j + user->mx;
970c4762a1bSJed Brown         ierr = MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES);CHKERRQ(ierr);
971c4762a1bSJed Brown       }
972c4762a1bSJed Brown       if (i>=2*m/3) {
973c4762a1bSJed Brown         j = i-2*m/3;
974c4762a1bSJed Brown         ierr = MatSetValues(user->Div,1,&j,1,&i,&neg_hinv,INSERT_VALUES);CHKERRQ(ierr);
975c4762a1bSJed Brown         j = j + user->mx*user->mx;
976c4762a1bSJed Brown         ierr = MatSetValues(user->Div,1,&j,1,&i,&hinv,INSERT_VALUES);CHKERRQ(ierr);
977c4762a1bSJed Brown       }
978c4762a1bSJed Brown     }
979c4762a1bSJed Brown 
980c4762a1bSJed Brown     ierr = MatAssemblyBegin(user->Div,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
981c4762a1bSJed Brown     ierr = MatAssemblyEnd(user->Div,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
982c4762a1bSJed Brown     ierr = MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);CHKERRQ(ierr);
983c4762a1bSJed Brown   } else {
984c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&user->Diag);CHKERRQ(ierr);
985c4762a1bSJed Brown     ierr = MatSetSizes(user->Diag,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr);
986c4762a1bSJed Brown     ierr = MatSetFromOptions(user->Diag);CHKERRQ(ierr);
987c4762a1bSJed Brown     ierr = MatMPIAIJSetPreallocation(user->Diag,1,NULL,0,NULL);CHKERRQ(ierr);
988c4762a1bSJed Brown     ierr = MatSeqAIJSetPreallocation(user->Diag,1,NULL);CHKERRQ(ierr);
989c4762a1bSJed Brown   }
990c4762a1bSJed Brown 
991c4762a1bSJed Brown   /* Build work vectors and matrices */
992c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->S);CHKERRQ(ierr);
993c4762a1bSJed Brown   ierr = VecSetSizes(user->S, PETSC_DECIDE, m);CHKERRQ(ierr);
994c4762a1bSJed Brown   ierr = VecSetFromOptions(user->S);CHKERRQ(ierr);
995c4762a1bSJed Brown 
996c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&user->lwork);CHKERRQ(ierr);
997c4762a1bSJed Brown   ierr = VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);CHKERRQ(ierr);
998c4762a1bSJed Brown   ierr = VecSetFromOptions(user->lwork);CHKERRQ(ierr);
999c4762a1bSJed Brown 
1000c4762a1bSJed Brown   ierr = MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);CHKERRQ(ierr);
1001c4762a1bSJed Brown 
1002c4762a1bSJed Brown   ierr = VecDuplicate(user->S,&user->Swork);CHKERRQ(ierr);
1003c4762a1bSJed Brown   ierr = VecDuplicate(user->S,&user->Sdiag);CHKERRQ(ierr);
1004c4762a1bSJed Brown   ierr = VecDuplicate(user->S,&user->Av_u);CHKERRQ(ierr);
1005c4762a1bSJed Brown   ierr = VecDuplicate(user->S,&user->Twork);CHKERRQ(ierr);
1006c4762a1bSJed Brown   ierr = VecDuplicate(user->y,&user->ywork);CHKERRQ(ierr);
1007c4762a1bSJed Brown   ierr = VecDuplicate(user->u,&user->Ywork);CHKERRQ(ierr);
1008c4762a1bSJed Brown   ierr = VecDuplicate(user->u,&user->uwork);CHKERRQ(ierr);
1009c4762a1bSJed Brown   ierr = VecDuplicate(user->u,&user->js_diag);CHKERRQ(ierr);
1010c4762a1bSJed Brown   ierr = VecDuplicate(user->c,&user->cwork);CHKERRQ(ierr);
1011c4762a1bSJed Brown   ierr = VecDuplicate(user->d,&user->dwork);CHKERRQ(ierr);
1012c4762a1bSJed Brown 
1013c4762a1bSJed Brown   /* Create a matrix-free shell user->Jd for computing B*x */
1014c4762a1bSJed Brown   ierr = MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal,user->nstate,user->ndesign,user,&user->Jd);CHKERRQ(ierr);
1015c4762a1bSJed Brown   ierr = MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult);CHKERRQ(ierr);
1016c4762a1bSJed Brown   ierr = MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose);CHKERRQ(ierr);
1017c4762a1bSJed Brown 
1018c4762a1bSJed Brown   /* Compute true state function ytrue given utrue */
1019c4762a1bSJed Brown   ierr = VecDuplicate(user->y,&user->ytrue);CHKERRQ(ierr);
1020c4762a1bSJed Brown 
1021c4762a1bSJed Brown   /* First compute Av_u = Av*exp(-u) */
1022c4762a1bSJed Brown   ierr = VecSet(user->uwork, 0);CHKERRQ(ierr);
1023c4762a1bSJed Brown   ierr = VecAXPY(user->uwork,-1.0,user->utrue);CHKERRQ(ierr); /* Note: user->utrue */
1024c4762a1bSJed Brown   ierr = VecExp(user->uwork);CHKERRQ(ierr);
1025c4762a1bSJed Brown   ierr = MatMult(user->Av,user->uwork,user->Av_u);CHKERRQ(ierr);
1026c4762a1bSJed Brown 
1027c4762a1bSJed Brown   /* Next form DSG = Div*S*Grad */
1028c4762a1bSJed Brown   ierr = VecCopy(user->Av_u,user->Swork);CHKERRQ(ierr);
1029c4762a1bSJed Brown   ierr = VecReciprocal(user->Swork);CHKERRQ(ierr);
1030c4762a1bSJed Brown   if (user->use_ptap) {
1031c4762a1bSJed Brown     ierr = MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);CHKERRQ(ierr);
1032c4762a1bSJed Brown     ierr = MatPtAP(user->Diag,user->Grad,MAT_INITIAL_MATRIX,1.0,&user->DSG);CHKERRQ(ierr);
1033c4762a1bSJed Brown   } else {
1034c4762a1bSJed Brown     ierr = MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1035c4762a1bSJed Brown     ierr = MatDiagonalScale(user->Divwork,NULL,user->Swork);CHKERRQ(ierr);
1036c20d7725SJed Brown 
1037c20d7725SJed Brown     ierr = MatMatMult(user->Divwork,user->Grad,MAT_INITIAL_MATRIX,1.0,&user->DSG);CHKERRQ(ierr);
1038c4762a1bSJed Brown   }
1039c4762a1bSJed Brown 
1040c4762a1bSJed Brown   ierr = MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
1041c4762a1bSJed Brown   ierr = MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr);
1042c4762a1bSJed Brown 
1043c4762a1bSJed Brown   if (user->use_lrc == PETSC_TRUE) {
1044c4762a1bSJed Brown     v=PetscSqrtReal(1.0 /user->ndesign);
1045c4762a1bSJed Brown     ierr = PetscMalloc1(user->ndesign,&user->ones);CHKERRQ(ierr);
1046c4762a1bSJed Brown 
1047c4762a1bSJed Brown     for (i=0;i<user->ndesign;i++) {
1048c4762a1bSJed Brown       user->ones[i]=v;
1049c4762a1bSJed Brown     }
1050c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,ysubnlocal,PETSC_DECIDE,user->ndesign,1,user->ones,&user->Ones);CHKERRQ(ierr);
1051c4762a1bSJed Brown     ierr = MatAssemblyBegin(user->Ones, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1052c4762a1bSJed Brown     ierr = MatAssemblyEnd(user->Ones, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1053c4762a1bSJed Brown     ierr = MatCreateLRC(user->DSG,user->Ones,NULL,user->Ones,&user->JsBlock);CHKERRQ(ierr);
1054c4762a1bSJed Brown     ierr = MatSetUp(user->JsBlock);CHKERRQ(ierr);
1055c4762a1bSJed Brown   } else {
1056c4762a1bSJed Brown     /* Create matrix-free shell user->Js for computing (A + h^3*e*e^T)*x */
1057c4762a1bSJed Brown     ierr = MatCreateShell(PETSC_COMM_WORLD,ysubnlocal,ysubnlocal,user->ndesign,user->ndesign,user,&user->JsBlock);CHKERRQ(ierr);
1058c4762a1bSJed Brown     ierr = MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateBlockMatMult);CHKERRQ(ierr);
1059c4762a1bSJed Brown     ierr = MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateBlockMatMult);CHKERRQ(ierr);
1060c4762a1bSJed Brown   }
1061c4762a1bSJed Brown   ierr = MatSetOption(user->JsBlock,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
1062c4762a1bSJed Brown   ierr = MatSetOption(user->JsBlock,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr);
1063c4762a1bSJed Brown   ierr = MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal*user->ns,user->nstate,user->nstate,user,&user->Js);CHKERRQ(ierr);
1064c4762a1bSJed Brown   ierr = MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);CHKERRQ(ierr);
1065c4762a1bSJed Brown   ierr = MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMult);CHKERRQ(ierr);
1066c4762a1bSJed Brown   ierr = MatSetOption(user->Js,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
1067c4762a1bSJed Brown   ierr = MatSetOption(user->Js,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr);
1068c4762a1bSJed Brown 
1069c4762a1bSJed Brown   ierr = MatCreateShell(PETSC_COMM_WORLD,ysubnlocal*user->ns,ysubnlocal*user->ns,user->nstate,user->nstate,user,&user->JsInv);CHKERRQ(ierr);
1070c4762a1bSJed Brown   ierr = MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateInvMatMult);CHKERRQ(ierr);
1071c4762a1bSJed Brown   ierr = MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateInvMatMult);CHKERRQ(ierr);
1072c4762a1bSJed Brown   ierr = MatSetOption(user->JsInv,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
1073c4762a1bSJed Brown   ierr = MatSetOption(user->JsInv,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr);
1074c4762a1bSJed Brown 
1075c4762a1bSJed Brown   ierr = MatSetOption(user->DSG,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
1076c4762a1bSJed Brown   ierr = MatSetOption(user->DSG,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr);
1077c4762a1bSJed Brown   /* Now solve for ytrue */
1078c4762a1bSJed Brown   ierr = KSPCreate(PETSC_COMM_WORLD,&user->solver);CHKERRQ(ierr);
1079c4762a1bSJed Brown   ierr = KSPSetFromOptions(user->solver);CHKERRQ(ierr);
1080c4762a1bSJed Brown 
1081c4762a1bSJed Brown   ierr = KSPSetOperators(user->solver,user->JsBlock,user->DSG);CHKERRQ(ierr);
1082c4762a1bSJed Brown 
1083c4762a1bSJed Brown   ierr = MatMult(user->JsInv,user->q,user->ytrue);CHKERRQ(ierr);
1084c4762a1bSJed Brown   /* First compute Av_u = Av*exp(-u) */
1085c4762a1bSJed Brown   ierr = VecSet(user->uwork,0);CHKERRQ(ierr);
1086c4762a1bSJed Brown   ierr = VecAXPY(user->uwork,-1.0,user->u);CHKERRQ(ierr); /* Note: user->u */
1087c4762a1bSJed Brown   ierr = VecExp(user->uwork);CHKERRQ(ierr);
1088c4762a1bSJed Brown   ierr = MatMult(user->Av,user->uwork,user->Av_u);CHKERRQ(ierr);
1089c4762a1bSJed Brown 
1090c4762a1bSJed Brown   /* Next update DSG = Div*S*Grad  with user->u */
1091c4762a1bSJed Brown   ierr = VecCopy(user->Av_u,user->Swork);CHKERRQ(ierr);
1092c4762a1bSJed Brown   ierr = VecReciprocal(user->Swork);CHKERRQ(ierr);
1093c4762a1bSJed Brown   if (user->use_ptap) {
1094c4762a1bSJed Brown     ierr = MatDiagonalSet(user->Diag,user->Swork,INSERT_VALUES);CHKERRQ(ierr);
1095c4762a1bSJed Brown     ierr = MatPtAP(user->Diag,user->Grad,MAT_REUSE_MATRIX,1.0,&user->DSG);CHKERRQ(ierr);
1096c4762a1bSJed Brown   } else {
1097c4762a1bSJed Brown     ierr = MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1098c4762a1bSJed Brown     ierr = MatDiagonalScale(user->Divwork,NULL,user->Av_u);CHKERRQ(ierr);
1099c20d7725SJed Brown     ierr = MatProductNumeric(user->DSG);CHKERRQ(ierr);
1100c4762a1bSJed Brown   }
1101c4762a1bSJed Brown 
1102c4762a1bSJed Brown   /* Now solve for y */
1103c4762a1bSJed Brown 
1104c4762a1bSJed Brown   ierr = MatMult(user->JsInv,user->q,user->y);CHKERRQ(ierr);
1105c4762a1bSJed Brown 
1106c4762a1bSJed Brown   user->ksp_its_initial = user->ksp_its;
1107c4762a1bSJed Brown   user->ksp_its = 0;
1108c4762a1bSJed Brown   /* Construct projection matrix Q (blocks) */
1109c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user->Q);CHKERRQ(ierr);
1110c4762a1bSJed Brown   ierr = MatSetSizes(user->Q,dsubnlocal,ysubnlocal,user->ndata,user->ndesign);CHKERRQ(ierr);
1111c4762a1bSJed Brown   ierr = MatSetFromOptions(user->Q);CHKERRQ(ierr);
1112c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(user->Q,8,NULL,8,NULL);CHKERRQ(ierr);
1113c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(user->Q,8,NULL);CHKERRQ(ierr);
1114c4762a1bSJed Brown 
1115c4762a1bSJed Brown   for (i=0; i<user->mx; i++) {
1116c4762a1bSJed Brown     x[i] = h*(i+0.5);
1117c4762a1bSJed Brown     y[i] = h*(i+0.5);
1118c4762a1bSJed Brown     z[i] = h*(i+0.5);
1119c4762a1bSJed Brown   }
1120c4762a1bSJed Brown   ierr = MatGetOwnershipRange(user->Q,&istart,&iend);CHKERRQ(ierr);
1121c4762a1bSJed Brown 
1122c4762a1bSJed Brown   nx = user->mx; ny = user->mx; nz = user->mx;
1123c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
1124c4762a1bSJed Brown 
1125c4762a1bSJed Brown     xri = xr[i];
1126c4762a1bSJed Brown     im = 0;
1127c4762a1bSJed Brown     xim = x[im];
1128c4762a1bSJed Brown     while (xri>xim && im<nx) {
1129c4762a1bSJed Brown       im = im+1;
1130c4762a1bSJed Brown       xim = x[im];
1131c4762a1bSJed Brown     }
1132c4762a1bSJed Brown     indx1 = im-1;
1133c4762a1bSJed Brown     indx2 = im;
1134c4762a1bSJed Brown     dx1 = xri - x[indx1];
1135c4762a1bSJed Brown     dx2 = x[indx2] - xri;
1136c4762a1bSJed Brown 
1137c4762a1bSJed Brown     yri = yr[i];
1138c4762a1bSJed Brown     im = 0;
1139c4762a1bSJed Brown     yim = y[im];
1140c4762a1bSJed Brown     while (yri>yim && im<ny) {
1141c4762a1bSJed Brown       im = im+1;
1142c4762a1bSJed Brown       yim = y[im];
1143c4762a1bSJed Brown     }
1144c4762a1bSJed Brown     indy1 = im-1;
1145c4762a1bSJed Brown     indy2 = im;
1146c4762a1bSJed Brown     dy1 = yri - y[indy1];
1147c4762a1bSJed Brown     dy2 = y[indy2] - yri;
1148c4762a1bSJed Brown 
1149c4762a1bSJed Brown     zri = zr[i];
1150c4762a1bSJed Brown     im = 0;
1151c4762a1bSJed Brown     zim = z[im];
1152c4762a1bSJed Brown     while (zri>zim && im<nz) {
1153c4762a1bSJed Brown       im = im+1;
1154c4762a1bSJed Brown       zim = z[im];
1155c4762a1bSJed Brown     }
1156c4762a1bSJed Brown     indz1 = im-1;
1157c4762a1bSJed Brown     indz2 = im;
1158c4762a1bSJed Brown     dz1 = zri - z[indz1];
1159c4762a1bSJed Brown     dz2 = z[indz2] - zri;
1160c4762a1bSJed Brown 
1161c4762a1bSJed Brown     Dx = x[indx2] - x[indx1];
1162c4762a1bSJed Brown     Dy = y[indy2] - y[indy1];
1163c4762a1bSJed Brown     Dz = z[indz2] - z[indz1];
1164c4762a1bSJed Brown 
1165c4762a1bSJed Brown     j = indx1 + indy1*nx + indz1*nx*ny;
1166c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1167c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1168c4762a1bSJed Brown 
1169c4762a1bSJed Brown     j = indx1 + indy1*nx + indz2*nx*ny;
1170c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1171c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1172c4762a1bSJed Brown 
1173c4762a1bSJed Brown     j = indx1 + indy2*nx + indz1*nx*ny;
1174c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1175c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1176c4762a1bSJed Brown 
1177c4762a1bSJed Brown     j = indx1 + indy2*nx + indz2*nx*ny;
1178c4762a1bSJed Brown     v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1179c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1180c4762a1bSJed Brown 
1181c4762a1bSJed Brown     j = indx2 + indy1*nx + indz1*nx*ny;
1182c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1183c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1184c4762a1bSJed Brown 
1185c4762a1bSJed Brown     j = indx2 + indy1*nx + indz2*nx*ny;
1186c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1187c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1188c4762a1bSJed Brown 
1189c4762a1bSJed Brown     j = indx2 + indy2*nx + indz1*nx*ny;
1190c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1191c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1192c4762a1bSJed Brown 
1193c4762a1bSJed Brown     j = indx2 + indy2*nx + indz2*nx*ny;
1194c4762a1bSJed Brown     v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1195c4762a1bSJed Brown     ierr = MatSetValues(user->Q,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
1196c4762a1bSJed Brown   }
1197c4762a1bSJed Brown 
1198c4762a1bSJed Brown   ierr = MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1199c4762a1bSJed Brown   ierr = MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1200c4762a1bSJed Brown   /* Create MQ (composed of blocks of Q */
1201c4762a1bSJed Brown   ierr = MatCreateShell(PETSC_COMM_WORLD,dsubnlocal*user->ns,PETSC_DECIDE,user->ndata*user->ns,user->nstate,user,&user->MQ);CHKERRQ(ierr);
1202c4762a1bSJed Brown   ierr = MatShellSetOperation(user->MQ,MATOP_MULT,(void(*)(void))QMatMult);CHKERRQ(ierr);
1203c4762a1bSJed Brown   ierr = MatShellSetOperation(user->MQ,MATOP_MULT_TRANSPOSE,(void(*)(void))QMatMultTranspose);CHKERRQ(ierr);
1204c4762a1bSJed Brown 
1205c4762a1bSJed Brown   /* Add noise to the measurement data */
1206c4762a1bSJed Brown   ierr = VecSet(user->ywork,1.0);CHKERRQ(ierr);
1207c4762a1bSJed Brown   ierr = VecAYPX(user->ywork,user->noise,user->ytrue);CHKERRQ(ierr);
1208c4762a1bSJed Brown   ierr = MatMult(user->MQ,user->ywork,user->d);CHKERRQ(ierr);
1209c4762a1bSJed Brown 
1210c4762a1bSJed Brown   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1211c4762a1bSJed Brown   ierr = PetscFree(x);CHKERRQ(ierr);
1212c4762a1bSJed Brown   ierr = PetscFree(y);CHKERRQ(ierr);
1213c4762a1bSJed Brown   ierr = PetscFree(z);CHKERRQ(ierr);
1214c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
1215c4762a1bSJed Brown   PetscFunctionReturn(0);
1216c4762a1bSJed Brown }
1217c4762a1bSJed Brown 
1218c4762a1bSJed Brown PetscErrorCode EllipticDestroy(AppCtx *user)
1219c4762a1bSJed Brown {
1220c4762a1bSJed Brown   PetscErrorCode ierr;
1221c4762a1bSJed Brown   PetscInt       i;
1222c4762a1bSJed Brown 
1223c4762a1bSJed Brown   PetscFunctionBegin;
1224c4762a1bSJed Brown   ierr = MatDestroy(&user->DSG);CHKERRQ(ierr);
1225c4762a1bSJed Brown   ierr = KSPDestroy(&user->solver);CHKERRQ(ierr);
1226c4762a1bSJed Brown   ierr = MatDestroy(&user->Q);CHKERRQ(ierr);
1227c4762a1bSJed Brown   ierr = MatDestroy(&user->MQ);CHKERRQ(ierr);
1228c4762a1bSJed Brown   if (!user->use_ptap) {
1229c4762a1bSJed Brown     ierr = MatDestroy(&user->Div);CHKERRQ(ierr);
1230c4762a1bSJed Brown     ierr = MatDestroy(&user->Divwork);CHKERRQ(ierr);
1231c4762a1bSJed Brown   } else {
1232c4762a1bSJed Brown     ierr = MatDestroy(&user->Diag);CHKERRQ(ierr);
1233c4762a1bSJed Brown   }
1234c4762a1bSJed Brown   if (user->use_lrc) {
1235c4762a1bSJed Brown     ierr = MatDestroy(&user->Ones);CHKERRQ(ierr);
1236c4762a1bSJed Brown   }
1237c4762a1bSJed Brown 
1238c4762a1bSJed Brown   ierr = MatDestroy(&user->Grad);CHKERRQ(ierr);
1239c4762a1bSJed Brown   ierr = MatDestroy(&user->Av);CHKERRQ(ierr);
1240c4762a1bSJed Brown   ierr = MatDestroy(&user->Avwork);CHKERRQ(ierr);
1241c4762a1bSJed Brown   ierr = MatDestroy(&user->L);CHKERRQ(ierr);
1242c4762a1bSJed Brown   ierr = MatDestroy(&user->Js);CHKERRQ(ierr);
1243c4762a1bSJed Brown   ierr = MatDestroy(&user->Jd);CHKERRQ(ierr);
1244c4762a1bSJed Brown   ierr = MatDestroy(&user->JsBlock);CHKERRQ(ierr);
1245c4762a1bSJed Brown   ierr = MatDestroy(&user->JsInv);CHKERRQ(ierr);
1246c4762a1bSJed Brown 
1247c4762a1bSJed Brown   ierr = VecDestroy(&user->x);CHKERRQ(ierr);
1248c4762a1bSJed Brown   ierr = VecDestroy(&user->u);CHKERRQ(ierr);
1249c4762a1bSJed Brown   ierr = VecDestroy(&user->uwork);CHKERRQ(ierr);
1250c4762a1bSJed Brown   ierr = VecDestroy(&user->utrue);CHKERRQ(ierr);
1251c4762a1bSJed Brown   ierr = VecDestroy(&user->y);CHKERRQ(ierr);
1252c4762a1bSJed Brown   ierr = VecDestroy(&user->ywork);CHKERRQ(ierr);
1253c4762a1bSJed Brown   ierr = VecDestroy(&user->ytrue);CHKERRQ(ierr);
1254c4762a1bSJed Brown   ierr = VecDestroy(&user->c);CHKERRQ(ierr);
1255c4762a1bSJed Brown   ierr = VecDestroy(&user->cwork);CHKERRQ(ierr);
1256c4762a1bSJed Brown   ierr = VecDestroy(&user->ur);CHKERRQ(ierr);
1257c4762a1bSJed Brown   ierr = VecDestroy(&user->q);CHKERRQ(ierr);
1258c4762a1bSJed Brown   ierr = VecDestroy(&user->d);CHKERRQ(ierr);
1259c4762a1bSJed Brown   ierr = VecDestroy(&user->dwork);CHKERRQ(ierr);
1260c4762a1bSJed Brown   ierr = VecDestroy(&user->lwork);CHKERRQ(ierr);
1261c4762a1bSJed Brown   ierr = VecDestroy(&user->S);CHKERRQ(ierr);
1262c4762a1bSJed Brown   ierr = VecDestroy(&user->Swork);CHKERRQ(ierr);
1263c4762a1bSJed Brown   ierr = VecDestroy(&user->Sdiag);CHKERRQ(ierr);
1264c4762a1bSJed Brown   ierr = VecDestroy(&user->Ywork);CHKERRQ(ierr);
1265c4762a1bSJed Brown   ierr = VecDestroy(&user->Twork);CHKERRQ(ierr);
1266c4762a1bSJed Brown   ierr = VecDestroy(&user->Av_u);CHKERRQ(ierr);
1267c4762a1bSJed Brown   ierr = VecDestroy(&user->js_diag);CHKERRQ(ierr);
1268c4762a1bSJed Brown   ierr = ISDestroy(&user->s_is);CHKERRQ(ierr);
1269c4762a1bSJed Brown   ierr = ISDestroy(&user->d_is);CHKERRQ(ierr);
1270c4762a1bSJed Brown   ierr = VecDestroy(&user->suby);CHKERRQ(ierr);
1271c4762a1bSJed Brown   ierr = VecDestroy(&user->subd);CHKERRQ(ierr);
1272c4762a1bSJed Brown   ierr = VecDestroy(&user->subq);CHKERRQ(ierr);
1273c4762a1bSJed Brown   ierr = VecScatterDestroy(&user->state_scatter);CHKERRQ(ierr);
1274c4762a1bSJed Brown   ierr = VecScatterDestroy(&user->design_scatter);CHKERRQ(ierr);
1275c4762a1bSJed Brown   for (i=0;i<user->ns;i++) {
1276c4762a1bSJed Brown     ierr = VecScatterDestroy(&user->yi_scatter[i]);CHKERRQ(ierr);
1277c4762a1bSJed Brown     ierr = VecScatterDestroy(&user->di_scatter[i]);CHKERRQ(ierr);
1278c4762a1bSJed Brown   }
1279c4762a1bSJed Brown   ierr = PetscFree(user->yi_scatter);CHKERRQ(ierr);
1280c4762a1bSJed Brown   ierr = PetscFree(user->di_scatter);CHKERRQ(ierr);
1281c4762a1bSJed Brown   if (user->use_lrc) {
1282c4762a1bSJed Brown     ierr = PetscFree(user->ones);CHKERRQ(ierr);
1283c4762a1bSJed Brown     ierr = MatDestroy(&user->Ones);CHKERRQ(ierr);
1284c4762a1bSJed Brown   }
1285c4762a1bSJed Brown   PetscFunctionReturn(0);
1286c4762a1bSJed Brown }
1287c4762a1bSJed Brown 
1288c4762a1bSJed Brown PetscErrorCode EllipticMonitor(Tao tao, void *ptr)
1289c4762a1bSJed Brown {
1290c4762a1bSJed Brown   PetscErrorCode ierr;
1291c4762a1bSJed Brown   Vec            X;
1292c4762a1bSJed Brown   PetscReal      unorm,ynorm;
1293c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
1294c4762a1bSJed Brown 
1295c4762a1bSJed Brown   PetscFunctionBegin;
1296*a82e8c82SStefano Zampini   ierr = TaoGetSolution(tao,&X);CHKERRQ(ierr);
1297c4762a1bSJed Brown   ierr = Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);CHKERRQ(ierr);
1298c4762a1bSJed Brown   ierr = VecAXPY(user->ywork,-1.0,user->ytrue);CHKERRQ(ierr);
1299c4762a1bSJed Brown   ierr = VecAXPY(user->uwork,-1.0,user->utrue);CHKERRQ(ierr);
1300c4762a1bSJed Brown   ierr = VecNorm(user->uwork,NORM_2,&unorm);CHKERRQ(ierr);
1301c4762a1bSJed Brown   ierr = VecNorm(user->ywork,NORM_2,&ynorm);CHKERRQ(ierr);
1302c4762a1bSJed Brown   ierr = PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);CHKERRQ(ierr);
1303c4762a1bSJed Brown   PetscFunctionReturn(0);
1304c4762a1bSJed Brown }
1305c4762a1bSJed Brown 
1306c4762a1bSJed Brown /*TEST
1307c4762a1bSJed Brown 
1308c4762a1bSJed Brown    build:
1309c4762a1bSJed Brown       requires: !complex
1310c4762a1bSJed Brown 
1311c4762a1bSJed Brown    test:
1312c4762a1bSJed Brown       args: -tao_cmonitor -ns 1 -tao_type lcl -tao_gatol 1.e-3 -tao_max_it 11
1313c4762a1bSJed Brown       requires: !single
1314c4762a1bSJed Brown 
1315c4762a1bSJed Brown    test:
1316c4762a1bSJed Brown       suffix: 2
1317c4762a1bSJed Brown       args: -tao_cmonitor -tao_type lcl -tao_max_it 11 -use_ptap -use_lrc -ns 1 -tao_gatol 1.e-3
1318c4762a1bSJed Brown       requires: !single
1319c4762a1bSJed Brown 
1320c4762a1bSJed Brown TEST*/
1321