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