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