xref: /petsc/src/tao/pde_constrained/tutorials/hyperbolic.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown #include <petsctao.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();
7a82e8c82SStefano Zampini    Routines: TaoSetSolution();
8a82e8c82SStefano Zampini    Routines: TaoSetObjective();
9a82e8c82SStefano 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: 1
18c4762a1bSJed Brown T*/
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct {
21c4762a1bSJed Brown   PetscInt n; /*  Number of variables */
22c4762a1bSJed Brown   PetscInt m; /*  Number of constraints */
23c4762a1bSJed Brown   PetscInt mx; /*  grid points in each direction */
24c4762a1bSJed Brown   PetscInt nt; /*  Number of time steps */
25c4762a1bSJed Brown   PetscInt ndata; /*  Number of data points per sample */
26c4762a1bSJed Brown   IS       s_is;
27c4762a1bSJed Brown   IS       d_is;
28c4762a1bSJed Brown   VecScatter state_scatter;
29c4762a1bSJed Brown   VecScatter design_scatter;
30c4762a1bSJed Brown   VecScatter *uxi_scatter,*uyi_scatter,*ux_scatter,*uy_scatter,*ui_scatter;
31c4762a1bSJed Brown   VecScatter *yi_scatter;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   Mat       Js,Jd,JsBlockPrec,JsInv,JsBlock;
34c4762a1bSJed Brown   PetscBool jformed,c_formed;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   PetscReal alpha; /*  Regularization parameter */
37c4762a1bSJed Brown   PetscReal gamma;
38c4762a1bSJed Brown   PetscReal ht; /*  Time step */
39c4762a1bSJed Brown   PetscReal T; /*  Final time */
40c4762a1bSJed Brown   Mat Q,QT;
41c4762a1bSJed Brown   Mat L,LT;
42c4762a1bSJed Brown   Mat Div,Divwork,Divxy[2];
43c4762a1bSJed Brown   Mat Grad,Gradxy[2];
44c4762a1bSJed Brown   Mat M;
45c4762a1bSJed Brown   Mat *C,*Cwork;
46c4762a1bSJed Brown   /* Mat Hs,Hd,Hsd; */
47c4762a1bSJed Brown   Vec q;
48c4762a1bSJed Brown   Vec ur; /*  reference */
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   Vec d;
51c4762a1bSJed Brown   Vec dwork;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   Vec y; /*  state variables */
54c4762a1bSJed Brown   Vec ywork;
55c4762a1bSJed Brown   Vec ytrue;
56c4762a1bSJed Brown   Vec *yi,*yiwork,*ziwork;
57c4762a1bSJed Brown   Vec *uxi,*uyi,*uxiwork,*uyiwork,*ui,*uiwork;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   Vec u; /*  design variables */
60c4762a1bSJed Brown   Vec uwork,vwork;
61c4762a1bSJed Brown   Vec utrue;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   Vec js_diag;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   Vec c; /*  constraint vector */
66c4762a1bSJed Brown   Vec cwork;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   Vec lwork;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   KSP      solver;
71c4762a1bSJed Brown   PC       prec;
72c4762a1bSJed Brown   PetscInt block_index;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   PetscInt ksp_its;
75c4762a1bSJed Brown   PetscInt ksp_its_initial;
76c4762a1bSJed Brown } AppCtx;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*);
79c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void*);
80c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*);
81c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void*);
82c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat,void*);
83c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void*);
84c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void*);
85c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
86c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
87c4762a1bSJed Brown PetscErrorCode HyperbolicInitialize(AppCtx *user);
88c4762a1bSJed Brown PetscErrorCode HyperbolicDestroy(AppCtx *user);
89c4762a1bSJed Brown PetscErrorCode HyperbolicMonitor(Tao, void*);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat,Vec,Vec);
92c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat,Vec,Vec);
93c4762a1bSJed Brown PetscErrorCode StateMatBlockMultTranspose(Mat,Vec,Vec);
94c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat,Vec,Vec);
95c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat,Vec);
96c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat,MatDuplicateOption,Mat*);
97c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat,Vec,Vec);
98c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat,Vec,Vec);
99c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC,Vec,Vec);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat,Vec,Vec);
102c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);
103c4762a1bSJed Brown 
104c4762a1bSJed Brown PetscErrorCode Scatter_yi(Vec,Vec*,VecScatter*,PetscInt); /*  y to y1,y2,...,y_nt */
105c4762a1bSJed Brown PetscErrorCode Gather_yi(Vec,Vec*,VecScatter*,PetscInt);
106c4762a1bSJed Brown PetscErrorCode Scatter_uxi_uyi(Vec,Vec*,VecScatter*,Vec*,VecScatter*,PetscInt); /*  u to ux_1,uy_1,ux_2,uy_2,...,u */
107c4762a1bSJed Brown PetscErrorCode Gather_uxi_uyi(Vec,Vec*,VecScatter*,Vec*,VecScatter*,PetscInt);
108c4762a1bSJed Brown 
109c4762a1bSJed Brown static  char help[]="";
110c4762a1bSJed Brown 
111c4762a1bSJed Brown int main(int argc, char **argv)
112c4762a1bSJed Brown {
113c4762a1bSJed Brown   PetscErrorCode     ierr;
114c4762a1bSJed Brown   Vec                x,x0;
115c4762a1bSJed Brown   Tao                tao;
116c4762a1bSJed Brown   AppCtx             user;
117c4762a1bSJed Brown   IS                 is_allstate,is_alldesign;
118c4762a1bSJed Brown   PetscInt           lo,hi,hi2,lo2,ksp_old;
119c4762a1bSJed Brown   PetscInt           ntests = 1;
120c4762a1bSJed Brown   PetscInt           i;
121c4762a1bSJed Brown #if defined(PETSC_USE_LOG)
122c4762a1bSJed Brown   PetscLogStage      stages[1];
123c4762a1bSJed Brown #endif
124c4762a1bSJed Brown 
125*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, (char*)0,help));
126c4762a1bSJed Brown   user.mx = 32;
12776280437SVaclav Hapla   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"hyperbolic example",NULL);CHKERRQ(ierr);
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL));
129c4762a1bSJed Brown   user.nt = 16;
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL));
131c4762a1bSJed Brown   user.ndata = 64;
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL));
133c4762a1bSJed Brown   user.alpha = 10.0;
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL));
135c4762a1bSJed Brown   user.T = 1.0/32.0;
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-Tfinal","Final time","",user.T,&user.T,NULL));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL));
13876280437SVaclav Hapla   ierr = PetscOptionsEnd();CHKERRQ(ierr);
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   user.m = user.mx*user.mx*user.nt; /*  number of constraints */
141c4762a1bSJed Brown   user.n = user.mx*user.mx*3*user.nt; /*  number of variables */
142c4762a1bSJed Brown   user.ht = user.T/user.nt; /*  Time step */
143c4762a1bSJed Brown   user.gamma = user.T*user.ht / (user.mx*user.mx);
144c4762a1bSJed Brown 
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user.u));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user.y));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user.c));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user.u,PETSC_DECIDE,user.n-user.m));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user.y,PETSC_DECIDE,user.m));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user.c,PETSC_DECIDE,user.m));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user.u));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user.y));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user.c));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Create scatters for reduced spaces.
156c4762a1bSJed Brown      If the state vector y and design vector u are partitioned as
157c4762a1bSJed Brown      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
158c4762a1bSJed Brown      then the solution vector x is organized as
159c4762a1bSJed Brown      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
160c4762a1bSJed Brown      The index sets user.s_is and user.d_is correspond to the indices of the
161c4762a1bSJed Brown      state and design variables owned by the current processor.
162c4762a1bSJed Brown   */
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
164c4762a1bSJed Brown 
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(user.y,&lo,&hi));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(user.u,&lo2,&hi2));
167c4762a1bSJed Brown 
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user.s_is));
170c4762a1bSJed Brown 
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user.d_is));
173c4762a1bSJed Brown 
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,hi-lo+hi2-lo2,user.n));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
176c4762a1bSJed Brown 
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,user.s_is,user.y,is_allstate,&user.state_scatter));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,user.d_is,user.u,is_alldesign,&user.design_scatter));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is_alldesign));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is_allstate));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOLCL));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* Set up initial vectors and matrices */
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(HyperbolicInitialize(&user));
188c4762a1bSJed Brown 
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather(x,user.y,user.state_scatter,user.u,user.design_scatter));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&x0));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,x0));
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   /* Set solution vector with an initial guess */
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjective(tao, FormFunction, &user));
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetGradient(tao, NULL, FormGradient, &user));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, &user));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetStateDesignIS(tao,user.s_is,user.d_is));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("Trials",&stages[0]));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(stages[0]));
206c4762a1bSJed Brown   user.ksp_its_initial = user.ksp_its;
207c4762a1bSJed Brown   ksp_old = user.ksp_its;
208c4762a1bSJed Brown   for (i=0; i<ntests; i++) {
2095f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSolve(tao));
2105f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its-ksp_old));
2115f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(x0,x));
2125f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetSolution(tao,x));
213c4762a1bSJed Brown   }
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBarrier((PetscObject)x));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: "));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %D trial(s): ",ntests));
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its));
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: "));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D\n",(user.ksp_its-user.ksp_its_initial)/ntests));
222c4762a1bSJed Brown 
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x0));
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(HyperbolicDestroy(&user));
227*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
228*b122ec5aSJacob Faibussowitsch   return 0;
229c4762a1bSJed Brown }
230c4762a1bSJed Brown /* ------------------------------------------------------------------- */
231c4762a1bSJed Brown /*
232c4762a1bSJed Brown    dwork = Qy - d
233c4762a1bSJed Brown    lwork = L*(u-ur).^2
234c4762a1bSJed Brown    f = 1/2 * (dwork.dork + alpha*y.lwork)
235c4762a1bSJed Brown */
236c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
237c4762a1bSJed Brown {
238c4762a1bSJed Brown   PetscReal      d1=0,d2=0;
239c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   PetscFunctionBegin;
2425f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Q,user->y,user->dwork));
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->dwork,-1.0,user->d));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(user->dwork,user->dwork,&d1));
246c4762a1bSJed Brown 
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(user->uwork,-1.0,user->ur,user->u));
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->uwork,user->uwork,user->uwork));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->L,user->uwork,user->lwork));
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(user->y,user->lwork,&d2));
251c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
252c4762a1bSJed Brown   PetscFunctionReturn(0);
253c4762a1bSJed Brown }
254c4762a1bSJed Brown 
255c4762a1bSJed Brown /* ------------------------------------------------------------------- */
256c4762a1bSJed Brown /*
257c4762a1bSJed Brown     state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2
258c4762a1bSJed Brown     design: g_d = alpha*(L'y).*(u-ur)
259c4762a1bSJed Brown */
260c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
261c4762a1bSJed Brown {
262c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   PetscFunctionBegin;
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Q,user->y,user->dwork));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->dwork,-1.0,user->d));
268c4762a1bSJed Brown 
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->QT,user->dwork,user->ywork));
270c4762a1bSJed Brown 
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->LT,user->y,user->uwork));
2725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(user->vwork,-1.0,user->ur,user->u));
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->uwork,user->vwork,user->uwork));
2745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(user->uwork,user->alpha));
275c4762a1bSJed Brown 
2765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->vwork,user->vwork,user->vwork));
2775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->L,user->vwork,user->lwork));
2785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->ywork,0.5*user->alpha,user->lwork));
279c4762a1bSJed Brown 
2805f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
281c4762a1bSJed Brown   PetscFunctionReturn(0);
282c4762a1bSJed Brown }
283c4762a1bSJed Brown 
284c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
285c4762a1bSJed Brown {
286c4762a1bSJed Brown   PetscReal      d1,d2;
287c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   PetscFunctionBegin;
2905f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Q,user->y,user->dwork));
2925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->dwork,-1.0,user->d));
293c4762a1bSJed Brown 
2945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->QT,user->dwork,user->ywork));
295c4762a1bSJed Brown 
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(user->dwork,user->dwork,&d1));
297c4762a1bSJed Brown 
2985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->LT,user->y,user->uwork));
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(user->vwork,-1.0,user->ur,user->u));
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->uwork,user->vwork,user->uwork));
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(user->uwork,user->alpha));
302c4762a1bSJed Brown 
3035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->vwork,user->vwork,user->vwork));
3045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->L,user->vwork,user->lwork));
3055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->ywork,0.5*user->alpha,user->lwork));
306c4762a1bSJed Brown 
3075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(user->y,user->lwork,&d2));
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   *f = 0.5 * (d1 + user->alpha*d2);
3105f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
311c4762a1bSJed Brown   PetscFunctionReturn(0);
312c4762a1bSJed Brown }
313c4762a1bSJed Brown 
314c4762a1bSJed Brown /* ------------------------------------------------------------------- */
315c4762a1bSJed Brown /* A
316c4762a1bSJed Brown MatShell object
317c4762a1bSJed Brown */
318c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
319c4762a1bSJed Brown {
320c4762a1bSJed Brown   PetscInt       i;
321c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   PetscFunctionBegin;
3245f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
3255f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(user->u,user->ui,user->ui_scatter,user->nt));
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
327c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
3285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCopy(user->Divxy[0],user->C[i],SUBSET_NONZERO_PATTERN));
3295f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCopy(user->Divxy[1],user->Cwork[i],SAME_NONZERO_PATTERN));
330c4762a1bSJed Brown 
3315f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDiagonalScale(user->C[i],NULL,user->uxi[i]));
3325f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]));
3335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(user->C[i],1.0,user->Cwork[i],SUBSET_NONZERO_PATTERN));
3345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatScale(user->C[i],user->ht));
3355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShift(user->C[i],1.0));
336c4762a1bSJed Brown   }
337c4762a1bSJed Brown   PetscFunctionReturn(0);
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown /* ------------------------------------------------------------------- */
341c4762a1bSJed Brown /* B */
342c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
343c4762a1bSJed Brown {
344c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   PetscFunctionBegin;
3475f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
348c4762a1bSJed Brown   PetscFunctionReturn(0);
349c4762a1bSJed Brown }
350c4762a1bSJed Brown 
351c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
352c4762a1bSJed Brown {
353c4762a1bSJed Brown   PetscInt       i;
354c4762a1bSJed Brown   AppCtx         *user;
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   PetscFunctionBegin;
3575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
3585f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
359c4762a1bSJed Brown   user->block_index = 0;
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->JsBlock,user->yi[0],user->yiwork[0]));
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
363c4762a1bSJed Brown     user->block_index = i;
3645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->JsBlock,user->yi[i],user->yiwork[i]));
3655f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->M,user->yi[i-1],user->ziwork[i-1]));
3665f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]));
367c4762a1bSJed Brown   }
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
369c4762a1bSJed Brown   PetscFunctionReturn(0);
370c4762a1bSJed Brown }
371c4762a1bSJed Brown 
372c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
373c4762a1bSJed Brown {
374c4762a1bSJed Brown   PetscInt       i;
375c4762a1bSJed Brown   AppCtx         *user;
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   PetscFunctionBegin;
3785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
3795f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
380c4762a1bSJed Brown 
381c4762a1bSJed Brown   for (i=0; i<user->nt-1; i++) {
382c4762a1bSJed Brown     user->block_index = i;
3835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]));
3845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->M,user->yi[i+1],user->ziwork[i+1]));
3855f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(user->yiwork[i],-1.0,user->ziwork[i+1]));
386c4762a1bSJed Brown   }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   i = user->nt-1;
389c4762a1bSJed Brown   user->block_index = i;
3905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(user->JsBlock,user->yi[i],user->yiwork[i]));
3915f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
392c4762a1bSJed Brown   PetscFunctionReturn(0);
393c4762a1bSJed Brown }
394c4762a1bSJed Brown 
395c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
396c4762a1bSJed Brown {
397c4762a1bSJed Brown   PetscInt       i;
398c4762a1bSJed Brown   AppCtx         *user;
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   PetscFunctionBegin;
4015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
402c4762a1bSJed Brown   i = user->block_index;
4035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->uxiwork[i],X,user->uxi[i]));
4045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->uyiwork[i],X,user->uyi[i]));
4055f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Div,user->uiwork[i],Y));
4075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAYPX(Y,user->ht,X));
408c4762a1bSJed Brown   PetscFunctionReturn(0);
409c4762a1bSJed Brown }
410c4762a1bSJed Brown 
411c4762a1bSJed Brown PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y)
412c4762a1bSJed Brown {
413c4762a1bSJed Brown   PetscInt       i;
414c4762a1bSJed Brown   AppCtx         *user;
415c4762a1bSJed Brown 
416c4762a1bSJed Brown   PetscFunctionBegin;
4175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
418c4762a1bSJed Brown   i = user->block_index;
4195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Grad,X,user->uiwork[i]));
4205f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->uxiwork[i],user->uxi[i],user->uxiwork[i]));
4225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(user->uyiwork[i],user->uyi[i],user->uyiwork[i]));
4235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(Y,1.0,user->uxiwork[i],user->uyiwork[i]));
4245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAYPX(Y,user->ht,X));
425c4762a1bSJed Brown   PetscFunctionReturn(0);
426c4762a1bSJed Brown }
427c4762a1bSJed Brown 
428c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
429c4762a1bSJed Brown {
430c4762a1bSJed Brown   PetscInt       i;
431c4762a1bSJed Brown   AppCtx         *user;
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   PetscFunctionBegin;
4345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
4355f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt));
4365f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_uxi_uyi(X,user->uxiwork,user->uxi_scatter,user->uyiwork,user->uyi_scatter,user->nt));
437c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
4385f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]));
4395f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]));
4405f80ce2aSJacob Faibussowitsch     CHKERRQ(Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->Div,user->uiwork[i],user->ziwork[i]));
4425f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(user->ziwork[i],user->ht));
443c4762a1bSJed Brown   }
4445f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(Y,user->ziwork,user->yi_scatter,user->nt));
445c4762a1bSJed Brown   PetscFunctionReturn(0);
446c4762a1bSJed Brown }
447c4762a1bSJed Brown 
448c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
449c4762a1bSJed Brown {
450c4762a1bSJed Brown   PetscInt       i;
451c4762a1bSJed Brown   AppCtx         *user;
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   PetscFunctionBegin;
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
4555f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt));
4565f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(X,user->yiwork,user->yi_scatter,user->nt));
457c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
4585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->Grad,user->yiwork[i],user->uiwork[i]));
4595f80ce2aSJacob Faibussowitsch     CHKERRQ(Scatter(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4605f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPointwiseMult(user->uxiwork[i],user->yi[i],user->uxiwork[i]));
4615f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPointwiseMult(user->uyiwork[i],user->yi[i],user->uyiwork[i]));
4625f80ce2aSJacob Faibussowitsch     CHKERRQ(Gather(user->uiwork[i],user->uxiwork[i],user->ux_scatter[i],user->uyiwork[i],user->uy_scatter[i]));
4635f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(user->uiwork[i],user->ht));
464c4762a1bSJed Brown   }
4655f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(Y,user->uiwork,user->ui_scatter,user->nt));
466c4762a1bSJed Brown   PetscFunctionReturn(0);
467c4762a1bSJed Brown }
468c4762a1bSJed Brown 
469c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
470c4762a1bSJed Brown {
471c4762a1bSJed Brown   PetscInt       i;
472c4762a1bSJed Brown   AppCtx         *user;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   PetscFunctionBegin;
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(PCShellGetContext(PC_shell,&user));
476c4762a1bSJed Brown   i = user->block_index;
477c4762a1bSJed Brown   if (user->c_formed) {
4785f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y));
479c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
480c4762a1bSJed Brown   PetscFunctionReturn(0);
481c4762a1bSJed Brown }
482c4762a1bSJed Brown 
483c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y)
484c4762a1bSJed Brown {
485c4762a1bSJed Brown   PetscInt       i;
486c4762a1bSJed Brown   AppCtx         *user;
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   PetscFunctionBegin;
4895f80ce2aSJacob Faibussowitsch   CHKERRQ(PCShellGetContext(PC_shell,&user));
490c4762a1bSJed Brown 
491c4762a1bSJed Brown   i = user->block_index;
492c4762a1bSJed Brown   if (user->c_formed) {
4935f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSOR(user->C[i],X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y));
494c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not formed");
495c4762a1bSJed Brown   PetscFunctionReturn(0);
496c4762a1bSJed Brown }
497c4762a1bSJed Brown 
498c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
499c4762a1bSJed Brown {
500c4762a1bSJed Brown   AppCtx         *user;
501c4762a1bSJed Brown   PetscInt       its,i;
502c4762a1bSJed Brown 
503c4762a1bSJed Brown   PetscFunctionBegin;
5045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
505c4762a1bSJed Brown 
506c4762a1bSJed Brown   if (Y == user->ytrue) {
507c4762a1bSJed Brown     /* First solve is done using true solution to set up problem */
5085f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(user->solver,1e-4,1e-20,PETSC_DEFAULT,PETSC_DEFAULT));
509c4762a1bSJed Brown   } else {
5105f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
511c4762a1bSJed Brown   }
5125f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
5135f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt));
5145f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
515c4762a1bSJed Brown 
516c4762a1bSJed Brown   user->block_index = 0;
5175f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolve(user->solver,user->yi[0],user->yiwork[0]));
518c4762a1bSJed Brown 
5195f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(user->solver,&its));
520c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
521c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
5225f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->M,user->yiwork[i-1],user->ziwork[i-1]));
5235f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(user->yi[i],1.0,user->ziwork[i-1]));
524c4762a1bSJed Brown     user->block_index = i;
5255f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSolve(user->solver,user->yi[i],user->yiwork[i]));
526c4762a1bSJed Brown 
5275f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetIterationNumber(user->solver,&its));
528c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
529c4762a1bSJed Brown   }
5305f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
531c4762a1bSJed Brown   PetscFunctionReturn(0);
532c4762a1bSJed Brown }
533c4762a1bSJed Brown 
534c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
535c4762a1bSJed Brown {
536c4762a1bSJed Brown   AppCtx         *user;
537c4762a1bSJed Brown   PetscInt       its,i;
538c4762a1bSJed Brown 
539c4762a1bSJed Brown   PetscFunctionBegin;
5405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
541c4762a1bSJed Brown 
5425f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(X,user->yi,user->yi_scatter,user->nt));
5435f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(Y,user->yiwork,user->yi_scatter,user->nt));
5445f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
545c4762a1bSJed Brown 
546c4762a1bSJed Brown   i = user->nt - 1;
547c4762a1bSJed Brown   user->block_index = i;
5485f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]));
549c4762a1bSJed Brown 
5505f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(user->solver,&its));
551c4762a1bSJed Brown   user->ksp_its = user->ksp_its + its;
552c4762a1bSJed Brown 
553c4762a1bSJed Brown   for (i=user->nt-2; i>=0; i--) {
5545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->M,user->yiwork[i+1],user->ziwork[i+1]));
5555f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(user->yi[i],1.0,user->ziwork[i+1]));
556c4762a1bSJed Brown     user->block_index = i;
5575f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSolveTranspose(user->solver,user->yi[i],user->yiwork[i]));
558c4762a1bSJed Brown 
5595f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetIterationNumber(user->solver,&its));
560c4762a1bSJed Brown     user->ksp_its = user->ksp_its + its;
561c4762a1bSJed Brown   }
5625f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(Y,user->yiwork,user->yi_scatter,user->nt));
563c4762a1bSJed Brown   PetscFunctionReturn(0);
564c4762a1bSJed Brown }
565c4762a1bSJed Brown 
566c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
567c4762a1bSJed Brown {
568c4762a1bSJed Brown   AppCtx         *user;
569c4762a1bSJed Brown 
570c4762a1bSJed Brown   PetscFunctionBegin;
5715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
572c4762a1bSJed Brown 
5735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell));
5745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult));
5755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate));
5765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose));
5775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal));
578c4762a1bSJed Brown   PetscFunctionReturn(0);
579c4762a1bSJed Brown }
580c4762a1bSJed Brown 
581c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
582c4762a1bSJed Brown {
583c4762a1bSJed Brown   AppCtx         *user;
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   PetscFunctionBegin;
5865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(J_shell,&user));
5875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(user->js_diag,X));
588c4762a1bSJed Brown   PetscFunctionReturn(0);
589c4762a1bSJed Brown }
590c4762a1bSJed Brown 
591c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
592c4762a1bSJed Brown {
593c4762a1bSJed Brown   /* con = Ay - q, A = [C(u1)  0     0     ...   0;
594c4762a1bSJed Brown                          -M  C(u2)   0     ...   0;
595c4762a1bSJed Brown                           0   -M   C(u3)   ...   0;
596c4762a1bSJed Brown                                       ...         ;
597c4762a1bSJed Brown                           0    ...      -M C(u_nt)]
598c4762a1bSJed Brown      C(u) = eye + ht*Div*[diag(u1); diag(u2)]       */
599c4762a1bSJed Brown   PetscInt       i;
600c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
601c4762a1bSJed Brown 
602c4762a1bSJed Brown   PetscFunctionBegin;
6035f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter));
6045f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_yi(user->y,user->yi,user->yi_scatter,user->nt));
6055f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
606c4762a1bSJed Brown 
607c4762a1bSJed Brown   user->block_index = 0;
6085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->JsBlock,user->yi[0],user->yiwork[0]));
609c4762a1bSJed Brown 
610c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
611c4762a1bSJed Brown     user->block_index = i;
6125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->JsBlock,user->yi[i],user->yiwork[i]));
6135f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(user->M,user->yi[i-1],user->ziwork[i-1]));
6145f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(user->yiwork[i],-1.0,user->ziwork[i-1]));
615c4762a1bSJed Brown   }
616c4762a1bSJed Brown 
6175f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(C,user->yiwork,user->yi_scatter,user->nt));
6185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(C,-1.0,user->q));
619c4762a1bSJed Brown 
620c4762a1bSJed Brown   PetscFunctionReturn(0);
621c4762a1bSJed Brown }
622c4762a1bSJed Brown 
623c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
624c4762a1bSJed Brown {
625c4762a1bSJed Brown   PetscFunctionBegin;
6265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD));
6275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD));
6285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD));
6295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD));
630c4762a1bSJed Brown   PetscFunctionReturn(0);
631c4762a1bSJed Brown }
632c4762a1bSJed Brown 
633c4762a1bSJed Brown PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
634c4762a1bSJed Brown {
635c4762a1bSJed Brown   PetscInt       i;
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   PetscFunctionBegin;
638c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6395f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD));
6405f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(scatx[i],u,uxi[i],INSERT_VALUES,SCATTER_FORWARD));
6415f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD));
6425f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(scaty[i],u,uyi[i],INSERT_VALUES,SCATTER_FORWARD));
643c4762a1bSJed Brown   }
644c4762a1bSJed Brown   PetscFunctionReturn(0);
645c4762a1bSJed Brown }
646c4762a1bSJed Brown 
647c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
648c4762a1bSJed Brown {
649c4762a1bSJed Brown   PetscFunctionBegin;
6505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE));
6515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE));
6525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE));
6535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE));
654c4762a1bSJed Brown   PetscFunctionReturn(0);
655c4762a1bSJed Brown }
656c4762a1bSJed Brown 
657c4762a1bSJed Brown PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt)
658c4762a1bSJed Brown {
659c4762a1bSJed Brown   PetscInt       i;
660c4762a1bSJed Brown 
661c4762a1bSJed Brown   PetscFunctionBegin;
662c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6635f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE));
6645f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(scatx[i],uxi[i],u,INSERT_VALUES,SCATTER_REVERSE));
6655f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE));
6665f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(scaty[i],uyi[i],u,INSERT_VALUES,SCATTER_REVERSE));
667c4762a1bSJed Brown   }
668c4762a1bSJed Brown   PetscFunctionReturn(0);
669c4762a1bSJed Brown }
670c4762a1bSJed Brown 
671c4762a1bSJed Brown PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
672c4762a1bSJed Brown {
673c4762a1bSJed Brown   PetscInt       i;
674c4762a1bSJed Brown 
675c4762a1bSJed Brown   PetscFunctionBegin;
676c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6775f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD));
6785f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD));
679c4762a1bSJed Brown   }
680c4762a1bSJed Brown   PetscFunctionReturn(0);
681c4762a1bSJed Brown }
682c4762a1bSJed Brown 
683c4762a1bSJed Brown PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
684c4762a1bSJed Brown {
685c4762a1bSJed Brown   PetscInt       i;
686c4762a1bSJed Brown 
687c4762a1bSJed Brown   PetscFunctionBegin;
688c4762a1bSJed Brown   for (i=0; i<nt; i++) {
6895f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE));
6905f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE));
691c4762a1bSJed Brown   }
692c4762a1bSJed Brown   PetscFunctionReturn(0);
693c4762a1bSJed Brown }
694c4762a1bSJed Brown 
695c4762a1bSJed Brown PetscErrorCode HyperbolicInitialize(AppCtx *user)
696c4762a1bSJed Brown {
697c4762a1bSJed Brown   PetscInt       n,i,j,linear_index,istart,iend,iblock,lo,hi;
698c4762a1bSJed Brown   Vec            XX,YY,XXwork,YYwork,yi,uxi,ui,bc;
699c4762a1bSJed Brown   PetscReal      h,sum;
700c4762a1bSJed Brown   PetscScalar    hinv,neg_hinv,quarter=0.25,one=1.0,half_hinv,neg_half_hinv;
701c4762a1bSJed Brown   PetscScalar    vx,vy,zero=0.0;
702c4762a1bSJed Brown   IS             is_from_y,is_to_yi,is_from_u,is_to_uxi,is_to_uyi;
703c4762a1bSJed Brown 
704c4762a1bSJed Brown   PetscFunctionBegin;
705c4762a1bSJed Brown   user->jformed = PETSC_FALSE;
706c4762a1bSJed Brown   user->c_formed = PETSC_FALSE;
707c4762a1bSJed Brown 
708c4762a1bSJed Brown   user->ksp_its = 0;
709c4762a1bSJed Brown   user->ksp_its_initial = 0;
710c4762a1bSJed Brown 
711c4762a1bSJed Brown   n = user->mx * user->mx;
712c4762a1bSJed Brown 
713c4762a1bSJed Brown   h = 1.0/user->mx;
714c4762a1bSJed Brown   hinv = user->mx;
715c4762a1bSJed Brown   neg_hinv = -hinv;
716c4762a1bSJed Brown   half_hinv = hinv / 2.0;
717c4762a1bSJed Brown   neg_half_hinv = neg_hinv / 2.0;
718c4762a1bSJed Brown 
719c4762a1bSJed Brown   /* Generate Grad matrix */
7205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->Grad));
7215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,2*n,n));
7225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Grad));
7235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->Grad,3,NULL,3,NULL));
7245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->Grad,3,NULL));
7255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(user->Grad,&istart,&iend));
726c4762a1bSJed Brown 
727c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
728c4762a1bSJed Brown     if (i<n) {
729c4762a1bSJed Brown       iblock = i / user->mx;
730c4762a1bSJed Brown       j = iblock*user->mx + ((i+user->mx-1) % user->mx);
7315f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES));
732c4762a1bSJed Brown       j = iblock*user->mx + ((i+1) % user->mx);
7335f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
734c4762a1bSJed Brown     }
735c4762a1bSJed Brown     if (i>=n) {
736c4762a1bSJed Brown       j = (i - user->mx) % n;
7375f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(user->Grad,1,&i,1,&j,&half_hinv,INSERT_VALUES));
738c4762a1bSJed Brown       j = (j + 2*user->mx) % n;
7395f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(user->Grad,1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
740c4762a1bSJed Brown     }
741c4762a1bSJed Brown   }
742c4762a1bSJed Brown 
7435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY));
7445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY));
745c4762a1bSJed Brown 
7465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->Gradxy[0]));
7475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->Gradxy[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
7485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Gradxy[0]));
7495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->Gradxy[0],3,NULL,3,NULL));
7505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->Gradxy[0],3,NULL));
7515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(user->Gradxy[0],&istart,&iend));
752c4762a1bSJed Brown 
753c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
754c4762a1bSJed Brown     iblock = i / user->mx;
755c4762a1bSJed Brown     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
7565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->Gradxy[0],1,&i,1,&j,&half_hinv,INSERT_VALUES));
757c4762a1bSJed Brown     j = iblock*user->mx + ((i+1) % user->mx);
7585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->Gradxy[0],1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
7595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->Gradxy[0],1,&i,1,&i,&zero,INSERT_VALUES));
760c4762a1bSJed Brown   }
7615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->Gradxy[0],MAT_FINAL_ASSEMBLY));
7625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->Gradxy[0],MAT_FINAL_ASSEMBLY));
763c4762a1bSJed Brown 
7645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->Gradxy[1]));
7655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->Gradxy[1],PETSC_DECIDE,PETSC_DECIDE,n,n));
7665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Gradxy[1]));
7675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->Gradxy[1],3,NULL,3,NULL));
7685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->Gradxy[1],3,NULL));
7695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(user->Gradxy[1],&istart,&iend));
770c4762a1bSJed Brown 
771c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
772c4762a1bSJed Brown     j = (i + n - user->mx) % n;
7735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->Gradxy[1],1,&i,1,&j,&half_hinv,INSERT_VALUES));
774c4762a1bSJed Brown     j = (j + 2*user->mx) % n;
7755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->Gradxy[1],1,&i,1,&j,&neg_half_hinv,INSERT_VALUES));
7765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->Gradxy[1],1,&i,1,&i,&zero,INSERT_VALUES));
777c4762a1bSJed Brown   }
7785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->Gradxy[1],MAT_FINAL_ASSEMBLY));
7795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->Gradxy[1],MAT_FINAL_ASSEMBLY));
780c4762a1bSJed Brown 
781c4762a1bSJed Brown   /* Generate Div matrix */
7825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div));
7835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(user->Gradxy[0],MAT_INITIAL_MATRIX,&user->Divxy[0]));
7845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(user->Gradxy[1],MAT_INITIAL_MATRIX,&user->Divxy[1]));
785c4762a1bSJed Brown 
786c4762a1bSJed Brown   /* Off-diagonal averaging matrix */
7875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->M));
7885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->M,PETSC_DECIDE,PETSC_DECIDE,n,n));
7895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->M));
7905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->M,4,NULL,4,NULL));
7915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->M,4,NULL));
7925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(user->M,&istart,&iend));
793c4762a1bSJed Brown 
794c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
795c4762a1bSJed Brown     /* kron(Id,Av) */
796c4762a1bSJed Brown     iblock = i / user->mx;
797c4762a1bSJed Brown     j = iblock*user->mx + ((i+user->mx-1) % user->mx);
7985f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
799c4762a1bSJed Brown     j = iblock*user->mx + ((i+1) % user->mx);
8005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
801c4762a1bSJed Brown 
802c4762a1bSJed Brown     /* kron(Av,Id) */
803c4762a1bSJed Brown     j = (i + user->mx) % n;
8045f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
805c4762a1bSJed Brown     j = (i + n - user->mx) % n;
8065f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->M,1,&i,1,&j,&quarter,INSERT_VALUES));
807c4762a1bSJed Brown   }
8085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->M,MAT_FINAL_ASSEMBLY));
8095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->M,MAT_FINAL_ASSEMBLY));
810c4762a1bSJed Brown 
811c4762a1bSJed Brown   /* Generate 2D grid */
8125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&XX));
8135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->q));
8145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(XX,PETSC_DECIDE,n));
8155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user->q,PETSC_DECIDE,n*user->nt));
8165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(XX));
8175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->q));
818c4762a1bSJed Brown 
8195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(XX,&YY));
8205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(XX,&XXwork));
8215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(XX,&YYwork));
8225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(XX,&user->d));
8235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(XX,&user->dwork));
824c4762a1bSJed Brown 
8255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(XX,&istart,&iend));
826c4762a1bSJed Brown   for (linear_index=istart; linear_index<iend; linear_index++) {
827c4762a1bSJed Brown     i = linear_index % user->mx;
828c4762a1bSJed Brown     j = (linear_index-i)/user->mx;
829c4762a1bSJed Brown     vx = h*(i+0.5);
830c4762a1bSJed Brown     vy = h*(j+0.5);
8315f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES));
8325f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES));
833c4762a1bSJed Brown   }
834c4762a1bSJed Brown 
8355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(XX));
8365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(XX));
8375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(YY));
8385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(YY));
839c4762a1bSJed Brown 
840c4762a1bSJed Brown   /* Compute final density function yT
841c4762a1bSJed Brown      yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2))
842c4762a1bSJed Brown      yT = yT / (h^2*sum(yT)) */
8435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(XX,XXwork));
8445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(YY,YYwork));
845c4762a1bSJed Brown 
8465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(XXwork,-0.25));
8475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(YYwork,-0.25));
848c4762a1bSJed Brown 
8495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(XXwork,XXwork,XXwork));
8505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(YYwork,YYwork,YYwork));
851c4762a1bSJed Brown 
8525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(XXwork,user->dwork));
8535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->dwork,1.0,YYwork));
8545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(user->dwork,-30.0));
8555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecExp(user->dwork));
8565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(user->dwork,user->d));
857c4762a1bSJed Brown 
8585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(XX,XXwork));
8595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(YY,YYwork));
860c4762a1bSJed Brown 
8615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(XXwork,-0.75));
8625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(YYwork,-0.75));
863c4762a1bSJed Brown 
8645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(XXwork,XXwork,XXwork));
8655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(YYwork,YYwork,YYwork));
866c4762a1bSJed Brown 
8675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(XXwork,user->dwork));
8685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->dwork,1.0,YYwork));
8695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(user->dwork,-30.0));
8705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecExp(user->dwork));
871c4762a1bSJed Brown 
8725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->d,1.0,user->dwork));
8735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(user->d,1.0));
8745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSum(user->d,&sum));
8755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(user->d,1.0/(h*h*sum)));
876c4762a1bSJed Brown 
877c4762a1bSJed Brown   /* Initial conditions of forward problem */
8785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(XX,&bc));
8795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(XX,XXwork));
8805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(YY,YYwork));
881c4762a1bSJed Brown 
8825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(XXwork,-0.5));
8835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(YYwork,-0.5));
884c4762a1bSJed Brown 
8855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(XXwork,XXwork,XXwork));
8865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(YYwork,YYwork,YYwork));
887c4762a1bSJed Brown 
8885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(bc,1.0,XXwork,YYwork));
8895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(bc,-50.0));
8905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecExp(bc));
8915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecShift(bc,1.0));
8925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSum(bc,&sum));
8935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(bc,1.0/(h*h*sum)));
894c4762a1bSJed Brown 
895c4762a1bSJed Brown   /* Create scatter from y to y_1,y_2,...,y_nt */
896c4762a1bSJed Brown   /*  TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */
8975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(user->nt*user->mx*user->mx,&user->yi_scatter));
8985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&yi));
8995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx));
9005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(yi));
9015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(yi,user->nt,&user->yi));
9025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(yi,user->nt,&user->yiwork));
9035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(yi,user->nt,&user->ziwork));
904c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
9055f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetOwnershipRange(user->yi[i],&lo,&hi));
9065f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi));
9075f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+i*user->mx*user->mx,1,&is_from_y));
9085f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]));
9095f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_to_yi));
9105f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_from_y));
911c4762a1bSJed Brown   }
912c4762a1bSJed Brown 
913c4762a1bSJed Brown   /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */
914c4762a1bSJed Brown   /*  TODO: reorder for better parallelism */
9155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(user->nt*user->mx*user->mx,&user->uxi_scatter));
9165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(user->nt*user->mx*user->mx,&user->uyi_scatter));
9175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(user->nt*user->mx*user->mx,&user->ux_scatter));
9185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(user->nt*user->mx*user->mx,&user->uy_scatter));
9195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*user->nt*user->mx*user->mx,&user->ui_scatter));
9205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&uxi));
9215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&ui));
9225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(uxi,PETSC_DECIDE,user->mx*user->mx));
9235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(ui,PETSC_DECIDE,2*user->mx*user->mx));
9245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(uxi));
9255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(ui));
9265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(uxi,user->nt,&user->uxi));
9275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(uxi,user->nt,&user->uyi));
9285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(uxi,user->nt,&user->uxiwork));
9295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(uxi,user->nt,&user->uyiwork));
9305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(ui,user->nt,&user->ui));
9315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(ui,user->nt,&user->uiwork));
932c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
9335f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetOwnershipRange(user->uxi[i],&lo,&hi));
9345f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi));
9355f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u));
9365f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(user->u,is_from_u,user->uxi[i],is_to_uxi,&user->uxi_scatter[i]));
937c4762a1bSJed Brown 
9385f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_to_uxi));
9395f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_from_u));
940c4762a1bSJed Brown 
9415f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetOwnershipRange(user->uyi[i],&lo,&hi));
9425f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi));
9435f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+(2*i+1)*user->mx*user->mx,1,&is_from_u));
9445f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(user->u,is_from_u,user->uyi[i],is_to_uyi,&user->uyi_scatter[i]));
945c4762a1bSJed Brown 
9465f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_to_uyi));
9475f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_from_u));
948c4762a1bSJed Brown 
9495f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetOwnershipRange(user->uxi[i],&lo,&hi));
9505f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi));
9515f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_from_u));
9525f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(user->ui[i],is_from_u,user->uxi[i],is_to_uxi,&user->ux_scatter[i]));
953c4762a1bSJed Brown 
9545f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_to_uxi));
9555f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_from_u));
956c4762a1bSJed Brown 
9575f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetOwnershipRange(user->uyi[i],&lo,&hi));
9585f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uyi));
9595f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+user->mx*user->mx,1,&is_from_u));
9605f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(user->ui[i],is_from_u,user->uyi[i],is_to_uyi,&user->uy_scatter[i]));
961c4762a1bSJed Brown 
9625f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_to_uyi));
9635f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_from_u));
964c4762a1bSJed Brown 
9655f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetOwnershipRange(user->ui[i],&lo,&hi));
9665f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_uxi));
9675f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+2*i*user->mx*user->mx,1,&is_from_u));
9685f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(user->u,is_from_u,user->ui[i],is_to_uxi,&user->ui_scatter[i]));
969c4762a1bSJed Brown 
9705f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_to_uxi));
9715f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is_from_u));
972c4762a1bSJed Brown   }
973c4762a1bSJed Brown 
974c4762a1bSJed Brown   /* RHS of forward problem */
9755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->M,bc,user->yiwork[0]));
976c4762a1bSJed Brown   for (i=1; i<user->nt; i++) {
9775f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(user->yiwork[i],0.0));
978c4762a1bSJed Brown   }
9795f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_yi(user->q,user->yiwork,user->yi_scatter,user->nt));
980c4762a1bSJed Brown 
981c4762a1bSJed Brown   /* Compute true velocity field utrue */
9825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->u,&user->utrue));
983c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
9845f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(YY,user->uxi[i]));
9855f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(user->uxi[i],150.0*i*user->ht));
9865f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(XX,user->uyi[i]));
9875f80ce2aSJacob Faibussowitsch     CHKERRQ(VecShift(user->uyi[i],-10.0));
9885f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(user->uyi[i],15.0*i*user->ht));
989c4762a1bSJed Brown   }
9905f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_uxi_uyi(user->utrue,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
991c4762a1bSJed Brown 
992c4762a1bSJed Brown   /* Initial guess and reference model */
9935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->utrue,&user->ur));
994c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
9955f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(XX,user->uxi[i]));
9965f80ce2aSJacob Faibussowitsch     CHKERRQ(VecShift(user->uxi[i],i*user->ht));
9975f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(YY,user->uyi[i]));
9985f80ce2aSJacob Faibussowitsch     CHKERRQ(VecShift(user->uyi[i],-i*user->ht));
999c4762a1bSJed Brown   }
10005f80ce2aSJacob Faibussowitsch   CHKERRQ(Gather_uxi_uyi(user->ur,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
1001c4762a1bSJed Brown 
1002c4762a1bSJed Brown   /* Generate regularization matrix L */
10035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->LT));
10045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->LT,PETSC_DECIDE,PETSC_DECIDE,2*n*user->nt,n*user->nt));
10055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->LT));
10065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->LT,1,NULL,1,NULL));
10075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->LT,1,NULL));
10085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(user->LT,&istart,&iend));
1009c4762a1bSJed Brown 
1010c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
1011c4762a1bSJed Brown     iblock = (i+n) / (2*n);
1012c4762a1bSJed Brown     j = i - iblock*n;
10135f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->LT,1,&i,1,&j,&user->gamma,INSERT_VALUES));
1014c4762a1bSJed Brown   }
1015c4762a1bSJed Brown 
10165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->LT,MAT_FINAL_ASSEMBLY));
10175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->LT,MAT_FINAL_ASSEMBLY));
1018c4762a1bSJed Brown 
10195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(user->LT,MAT_INITIAL_MATRIX,&user->L));
1020c4762a1bSJed Brown 
1021c4762a1bSJed Brown   /* Build work vectors and matrices */
10225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->lwork));
10235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetType(user->lwork,VECMPI));
10245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user->lwork,PETSC_DECIDE,user->m));
10255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->lwork));
1026c4762a1bSJed Brown 
10275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork));
1028c4762a1bSJed Brown 
10295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->y,&user->ywork));
10305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->u,&user->uwork));
10315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->u,&user->vwork));
10325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->u,&user->js_diag));
10335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user->c,&user->cwork));
1034c4762a1bSJed Brown 
1035c4762a1bSJed Brown   /* Create matrix-free shell user->Js for computing A*x */
10365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->Js));
10375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult));
10385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate));
10395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose));
10405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal));
1041c4762a1bSJed Brown 
1042c4762a1bSJed Brown   /* Diagonal blocks of user->Js */
10435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlock));
10445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult));
10455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMultTranspose));
1046c4762a1bSJed Brown 
1047c4762a1bSJed Brown   /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1048c4762a1bSJed Brown      D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1049c4762a1bSJed Brown      This is an SOR preconditioner for user->JsBlock. */
10505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,n,n,user,&user->JsBlockPrec));
10515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult));
10525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMultTranspose));
1053c4762a1bSJed Brown 
1054c4762a1bSJed Brown   /* Create a matrix-free shell user->Jd for computing B*x */
10555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->n-user->m,user,&user->Jd));
10565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult));
10575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose));
1058c4762a1bSJed Brown 
1059c4762a1bSJed Brown   /* User-defined routines for computing user->Js\x and user->Js^T\x*/
10605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsInv));
10615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult));
10625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult));
1063c4762a1bSJed Brown 
1064c4762a1bSJed Brown   /* Build matrices for SOR preconditioner */
10655f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter_uxi_uyi(user->u,user->uxi,user->uxi_scatter,user->uyi,user->uyi_scatter,user->nt));
10665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(5*n,&user->C));
10675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*n,&user->Cwork));
1068c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
10695f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(user->Divxy[0],MAT_COPY_VALUES,&user->C[i]));
10705f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(user->Divxy[1],MAT_COPY_VALUES,&user->Cwork[i]));
1071c4762a1bSJed Brown 
10725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDiagonalScale(user->C[i],NULL,user->uxi[i]));
10735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDiagonalScale(user->Cwork[i],NULL,user->uyi[i]));
10745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(user->C[i],1.0,user->Cwork[i],DIFFERENT_NONZERO_PATTERN));
10755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatScale(user->C[i],user->ht));
10765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShift(user->C[i],1.0));
1077c4762a1bSJed Brown   }
1078c4762a1bSJed Brown 
1079c4762a1bSJed Brown   /* Solver options and tolerances */
10805f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&user->solver));
10815f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetType(user->solver,KSPGMRES));
10825f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec));
10835f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500));
10845f80ce2aSJacob Faibussowitsch   /* CHKERRQ(KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500)); */
10855f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(user->solver,&user->prec));
10865f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(user->prec,PCSHELL));
1087c4762a1bSJed Brown 
10885f80ce2aSJacob Faibussowitsch   CHKERRQ(PCShellSetApply(user->prec,StateMatBlockPrecMult));
10895f80ce2aSJacob Faibussowitsch   CHKERRQ(PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMultTranspose));
10905f80ce2aSJacob Faibussowitsch   CHKERRQ(PCShellSetContext(user->prec,user));
1091c4762a1bSJed Brown 
1092c4762a1bSJed Brown   /* Compute true state function yt given ut */
10935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->ytrue));
10945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt));
10955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->ytrue));
1096c4762a1bSJed Brown   user->c_formed = PETSC_TRUE;
10975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(user->utrue,user->u)); /*  Set u=utrue temporarily for StateMatInv */
10985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(user->ytrue,0.0)); /*  Initial guess */
10995f80ce2aSJacob Faibussowitsch   CHKERRQ(StateMatInvMult(user->Js,user->q,user->ytrue));
11005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(user->ur,user->u)); /*  Reset u=ur */
1101c4762a1bSJed Brown 
1102c4762a1bSJed Brown   /* Initial guess y0 for state given u0 */
11035f80ce2aSJacob Faibussowitsch   CHKERRQ(StateMatInvMult(user->Js,user->q,user->y));
1104c4762a1bSJed Brown 
1105c4762a1bSJed Brown   /* Data discretization */
11065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->Q));
11075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->Q,PETSC_DECIDE,PETSC_DECIDE,user->mx*user->mx,user->m));
11085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Q));
11095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->Q,0,NULL,1,NULL));
11105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->Q,1,NULL));
1111c4762a1bSJed Brown 
11125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(user->Q,&istart,&iend));
1113c4762a1bSJed Brown 
1114c4762a1bSJed Brown   for (i=istart; i<iend; i++) {
1115c4762a1bSJed Brown     j = i + user->m - user->mx*user->mx;
11165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(user->Q,1,&i,1,&j,&one,INSERT_VALUES));
1117c4762a1bSJed Brown   }
1118c4762a1bSJed Brown 
11195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->Q,MAT_FINAL_ASSEMBLY));
11205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->Q,MAT_FINAL_ASSEMBLY));
1121c4762a1bSJed Brown 
11225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(user->Q,MAT_INITIAL_MATRIX,&user->QT));
1123c4762a1bSJed Brown 
11245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&XX));
11255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&YY));
11265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&XXwork));
11275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&YYwork));
11285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&yi));
11295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&uxi));
11305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ui));
11315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&bc));
1132c4762a1bSJed Brown 
1133c4762a1bSJed Brown   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
11345f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(user->solver));
1135c4762a1bSJed Brown   PetscFunctionReturn(0);
1136c4762a1bSJed Brown }
1137c4762a1bSJed Brown 
1138c4762a1bSJed Brown PetscErrorCode HyperbolicDestroy(AppCtx *user)
1139c4762a1bSJed Brown {
1140c4762a1bSJed Brown   PetscInt       i;
1141c4762a1bSJed Brown 
1142c4762a1bSJed Brown   PetscFunctionBegin;
11435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Q));
11445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->QT));
11455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Div));
11465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Divwork));
11475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Grad));
11485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->L));
11495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->LT));
11505f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&user->solver));
11515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Js));
11525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Jd));
11535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->JsBlockPrec));
11545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->JsInv));
11555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->JsBlock));
11565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Divxy[0]));
11575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Divxy[1]));
11585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Gradxy[0]));
11595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Gradxy[1]));
11605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->M));
1161c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
11625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&user->C[i]));
11635f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&user->Cwork[i]));
1164c4762a1bSJed Brown   }
11655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->C));
11665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->Cwork));
11675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->u));
11685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->uwork));
11695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->vwork));
11705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->utrue));
11715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->y));
11725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->ywork));
11735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->ytrue));
11745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->yi));
11755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->yiwork));
11765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->ziwork));
11775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->uxi));
11785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->uyi));
11795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->uxiwork));
11805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->uyiwork));
11815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->ui));
11825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroyVecs(user->nt,&user->uiwork));
11835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->c));
11845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->cwork));
11855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->ur));
11865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->q));
11875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->d));
11885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->dwork));
11895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->lwork));
11905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->js_diag));
11915f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&user->s_is));
11925f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&user->d_is));
11935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&user->state_scatter));
11945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&user->design_scatter));
1195c4762a1bSJed Brown   for (i=0; i<user->nt; i++) {
11965f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&user->uxi_scatter[i]));
11975f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&user->uyi_scatter[i]));
11985f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&user->ux_scatter[i]));
11995f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&user->uy_scatter[i]));
12005f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&user->ui_scatter[i]));
12015f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&user->yi_scatter[i]));
1202c4762a1bSJed Brown   }
12035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->uxi_scatter));
12045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->uyi_scatter));
12055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->ux_scatter));
12065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->uy_scatter));
12075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->ui_scatter));
12085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->yi_scatter));
1209c4762a1bSJed Brown   PetscFunctionReturn(0);
1210c4762a1bSJed Brown }
1211c4762a1bSJed Brown 
1212c4762a1bSJed Brown PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr)
1213c4762a1bSJed Brown {
1214c4762a1bSJed Brown   Vec            X;
1215c4762a1bSJed Brown   PetscReal      unorm,ynorm;
1216c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
1217c4762a1bSJed Brown 
1218c4762a1bSJed Brown   PetscFunctionBegin;
12195f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetSolution(tao,&X));
12205f80ce2aSJacob Faibussowitsch   CHKERRQ(Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter));
12215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->ywork,-1.0,user->ytrue));
12225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(user->uwork,-1.0,user->utrue));
12235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(user->uwork,NORM_2,&unorm));
12245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(user->ywork,NORM_2,&ynorm));
12255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm));
1226c4762a1bSJed Brown   PetscFunctionReturn(0);
1227c4762a1bSJed Brown }
1228c4762a1bSJed Brown 
1229c4762a1bSJed Brown /*TEST
1230c4762a1bSJed Brown 
1231c4762a1bSJed Brown    build:
1232c4762a1bSJed Brown       requires: !complex
1233c4762a1bSJed Brown 
1234c4762a1bSJed Brown    test:
1235c4762a1bSJed Brown       requires: !single
1236c4762a1bSJed Brown       args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -tao_gatol 1.e-5
1237c4762a1bSJed Brown 
1238c4762a1bSJed Brown    test:
1239c4762a1bSJed Brown       suffix: guess_pod
1240c4762a1bSJed Brown       requires: !single
1241c4762a1bSJed Brown       args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -ksp_guess_type pod -tao_gatol 1.e-5
1242c4762a1bSJed Brown 
1243c4762a1bSJed Brown TEST*/
1244