xref: /petsc/src/tao/unconstrained/tutorials/burgers_spectral.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] ="Solves a simple data assimilation problem with one dimensional Burger's equation using TSAdjoint\n\n";
3 
4 /*
5 
6     Not yet tested in parallel
7 
8 */
9 /*
10    Concepts: TS^time-dependent nonlinear problems
11    Concepts: TS^Burger's equation
12    Concepts: adjoints
13    Processors: n
14 */
15 
16 /* ------------------------------------------------------------------------
17 
18    This program uses the one-dimensional Burger's equation
19        u_t = mu*u_xx - u u_x,
20    on the domain 0 <= x <= 1, with periodic boundary conditions
21 
22    to demonstrate solving a data assimilation problem of finding the initial conditions
23    to produce a given solution at a fixed time.
24 
25    The operators are discretized with the spectral element method
26 
27    See the paper PDE-CONSTRAINED OPTIMIZATION WITH SPECTRAL ELEMENTS USING PETSC AND TAO
28    by OANA MARIN, EMIL CONSTANTINESCU, AND BARRY SMITH for details on the exact solution
29    used
30 
31   ------------------------------------------------------------------------- */
32 
33 #include <petsctao.h>
34 #include <petscts.h>
35 #include <petscdt.h>
36 #include <petscdraw.h>
37 #include <petscdmda.h>
38 
39 /*
40    User-defined application context - contains data needed by the
41    application-provided call-back routines.
42 */
43 
44 typedef struct {
45   PetscInt  n;                /* number of nodes */
46   PetscReal *nodes;           /* GLL nodes */
47   PetscReal *weights;         /* GLL weights */
48 } PetscGLL;
49 
50 typedef struct {
51   PetscInt    N;              /* grid points per elements*/
52   PetscInt    E;              /* number of elements */
53   PetscReal   tol_L2,tol_max; /* error norms */
54   PetscInt    steps;          /* number of timesteps */
55   PetscReal   Tend;           /* endtime */
56   PetscReal   mu;             /* viscosity */
57   PetscReal   L;              /* total length of domain */
58   PetscReal   Le;
59   PetscReal   Tadj;
60 } PetscParam;
61 
62 typedef struct {
63   Vec         obj;               /* desired end state */
64   Vec         grid;              /* total grid */
65   Vec         grad;
66   Vec         ic;
67   Vec         curr_sol;
68   Vec         true_solution;     /* actual initial conditions for the final solution */
69 } PetscData;
70 
71 typedef struct {
72   Vec         grid;              /* total grid */
73   Vec         mass;              /* mass matrix for total integration */
74   Mat         stiff;             /* stifness matrix */
75   Mat         keptstiff;
76   Mat         grad;
77   PetscGLL    gll;
78 } PetscSEMOperators;
79 
80 typedef struct {
81   DM                da;                /* distributed array data structure */
82   PetscSEMOperators SEMop;
83   PetscParam        param;
84   PetscData         dat;
85   TS                ts;
86   PetscReal         initial_dt;
87 } AppCtx;
88 
89 /*
90    User-defined routines
91 */
92 extern PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
93 extern PetscErrorCode RHSMatrixLaplaciangllDM(TS,PetscReal,Vec,Mat,Mat,void*);
94 extern PetscErrorCode RHSMatrixAdvectiongllDM(TS,PetscReal,Vec,Mat,Mat,void*);
95 extern PetscErrorCode InitialConditions(Vec,AppCtx*);
96 extern PetscErrorCode TrueSolution(Vec,AppCtx*);
97 extern PetscErrorCode ComputeObjective(PetscReal,Vec,AppCtx*);
98 extern PetscErrorCode MonitorError(Tao,void*);
99 extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
100 extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
101 
102 int main(int argc,char **argv)
103 {
104   AppCtx         appctx;                 /* user-defined application context */
105   Tao            tao;
106   Vec            u;                      /* approximate solution vector */
107   PetscErrorCode ierr;
108   PetscInt       i, xs, xm, ind, j, lenglob;
109   PetscReal      x, *wrk_ptr1, *wrk_ptr2;
110   MatNullSpace   nsp;
111   PetscMPIInt    size;
112 
113    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114      Initialize program and set problem parameters
115      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116   PetscFunctionBegin;
117 
118   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
119 
120   /*initialize parameters */
121   appctx.param.N    = 10;  /* order of the spectral element */
122   appctx.param.E    = 10;  /* number of elements */
123   appctx.param.L    = 4.0;  /* length of the domain */
124   appctx.param.mu   = 0.01; /* diffusion coefficient */
125   appctx.initial_dt = 5e-3;
126   appctx.param.steps = PETSC_MAX_INT;
127   appctx.param.Tend  = 4;
128 
129   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&appctx.param.N,NULL));
130   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-E",&appctx.param.E,NULL));
131   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-Tend",&appctx.param.Tend,NULL));
132   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&appctx.param.mu,NULL));
133   appctx.param.Le = appctx.param.L/appctx.param.E;
134 
135   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
136   PetscCheck((appctx.param.E % size) == 0,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of elements must be divisible by number of processes");
137 
138   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139      Create GLL data structures
140      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141   CHKERRQ(PetscMalloc2(appctx.param.N,&appctx.SEMop.gll.nodes,appctx.param.N,&appctx.SEMop.gll.weights));
142   CHKERRQ(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights));
143   appctx.SEMop.gll.n = appctx.param.N;
144   lenglob  = appctx.param.E*(appctx.param.N-1);
145 
146   /*
147      Create distributed array (DMDA) to manage parallel grid and vectors
148      and to set up the ghost point communication pattern.  There are E*(Nl-1)+1
149      total grid values spread equally among all the processors, except first and last
150   */
151 
152   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,lenglob,1,1,NULL,&appctx.da));
153   CHKERRQ(DMSetFromOptions(appctx.da));
154   CHKERRQ(DMSetUp(appctx.da));
155 
156   /*
157      Extract global and local vectors from DMDA; we use these to store the
158      approximate solution.  Then duplicate these for remaining vectors that
159      have the same types.
160   */
161 
162   CHKERRQ(DMCreateGlobalVector(appctx.da,&u));
163   CHKERRQ(VecDuplicate(u,&appctx.dat.ic));
164   CHKERRQ(VecDuplicate(u,&appctx.dat.true_solution));
165   CHKERRQ(VecDuplicate(u,&appctx.dat.obj));
166   CHKERRQ(VecDuplicate(u,&appctx.SEMop.grid));
167   CHKERRQ(VecDuplicate(u,&appctx.SEMop.mass));
168   CHKERRQ(VecDuplicate(u,&appctx.dat.curr_sol));
169 
170   CHKERRQ(DMDAGetCorners(appctx.da,&xs,NULL,NULL,&xm,NULL,NULL));
171   CHKERRQ(DMDAVecGetArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1));
172   CHKERRQ(DMDAVecGetArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2));
173 
174   /* Compute function over the locally owned part of the grid */
175 
176     xs=xs/(appctx.param.N-1);
177     xm=xm/(appctx.param.N-1);
178 
179   /*
180      Build total grid and mass over entire mesh (multi-elemental)
181   */
182 
183   for (i=xs; i<xs+xm; i++) {
184     for (j=0; j<appctx.param.N-1; j++) {
185       x = (appctx.param.Le/2.0)*(appctx.SEMop.gll.nodes[j]+1.0)+appctx.param.Le*i;
186       ind=i*(appctx.param.N-1)+j;
187       wrk_ptr1[ind]=x;
188       wrk_ptr2[ind]=.5*appctx.param.Le*appctx.SEMop.gll.weights[j];
189       if (j==0) wrk_ptr2[ind]+=.5*appctx.param.Le*appctx.SEMop.gll.weights[j];
190     }
191   }
192   CHKERRQ(DMDAVecRestoreArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1));
193   CHKERRQ(DMDAVecRestoreArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2));
194 
195   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196    Create matrix data structure; set matrix evaluation routine.
197    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198   CHKERRQ(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE));
199   CHKERRQ(DMCreateMatrix(appctx.da,&appctx.SEMop.stiff));
200   CHKERRQ(DMCreateMatrix(appctx.da,&appctx.SEMop.grad));
201   /*
202    For linear problems with a time-dependent f(u,t) in the equation
203    u_t = f(u,t), the user provides the discretized right-hand-side
204    as a time-dependent matrix.
205    */
206   CHKERRQ(RHSMatrixLaplaciangllDM(appctx.ts,0.0,u,appctx.SEMop.stiff,appctx.SEMop.stiff,&appctx));
207   CHKERRQ(RHSMatrixAdvectiongllDM(appctx.ts,0.0,u,appctx.SEMop.grad,appctx.SEMop.grad,&appctx));
208    /*
209        For linear problems with a time-dependent f(u,t) in the equation
210        u_t = f(u,t), the user provides the discretized right-hand-side
211        as a time-dependent matrix.
212     */
213 
214   CHKERRQ(MatDuplicate(appctx.SEMop.stiff,MAT_COPY_VALUES,&appctx.SEMop.keptstiff));
215 
216   /* attach the null space to the matrix, this probably is not needed but does no harm */
217   CHKERRQ(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp));
218   CHKERRQ(MatSetNullSpace(appctx.SEMop.stiff,nsp));
219   CHKERRQ(MatSetNullSpace(appctx.SEMop.keptstiff,nsp));
220   CHKERRQ(MatNullSpaceTest(nsp,appctx.SEMop.stiff,NULL));
221   CHKERRQ(MatNullSpaceDestroy(&nsp));
222   /* attach the null space to the matrix, this probably is not needed but does no harm */
223   CHKERRQ(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp));
224   CHKERRQ(MatSetNullSpace(appctx.SEMop.grad,nsp));
225   CHKERRQ(MatNullSpaceTest(nsp,appctx.SEMop.grad,NULL));
226   CHKERRQ(MatNullSpaceDestroy(&nsp));
227 
228   /* Create the TS solver that solves the ODE and its adjoint; set its options */
229   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&appctx.ts));
230   CHKERRQ(TSSetProblemType(appctx.ts,TS_NONLINEAR));
231   CHKERRQ(TSSetType(appctx.ts,TSRK));
232   CHKERRQ(TSSetDM(appctx.ts,appctx.da));
233   CHKERRQ(TSSetTime(appctx.ts,0.0));
234   CHKERRQ(TSSetTimeStep(appctx.ts,appctx.initial_dt));
235   CHKERRQ(TSSetMaxSteps(appctx.ts,appctx.param.steps));
236   CHKERRQ(TSSetMaxTime(appctx.ts,appctx.param.Tend));
237   CHKERRQ(TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP));
238   CHKERRQ(TSSetTolerances(appctx.ts,1e-7,NULL,1e-7,NULL));
239   CHKERRQ(TSSetFromOptions(appctx.ts));
240   /* Need to save initial timestep user may have set with -ts_dt so it can be reset for each new TSSolve() */
241   CHKERRQ(TSGetTimeStep(appctx.ts,&appctx.initial_dt));
242   CHKERRQ(TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx));
243   CHKERRQ(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx));
244 
245   /* Set Objective and Initial conditions for the problem and compute Objective function (evolution of true_solution to final time */
246   CHKERRQ(InitialConditions(appctx.dat.ic,&appctx));
247   CHKERRQ(TrueSolution(appctx.dat.true_solution,&appctx));
248   CHKERRQ(ComputeObjective(appctx.param.Tend,appctx.dat.obj,&appctx));
249 
250   CHKERRQ(TSSetSaveTrajectory(appctx.ts));
251   CHKERRQ(TSSetFromOptions(appctx.ts));
252 
253   /* Create TAO solver and set desired solution method  */
254   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
255   CHKERRQ(TaoSetMonitor(tao,MonitorError,&appctx,NULL));
256   CHKERRQ(TaoSetType(tao,TAOBQNLS));
257   CHKERRQ(TaoSetSolution(tao,appctx.dat.ic));
258   /* Set routine for function and gradient evaluation  */
259   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&appctx));
260   /* Check for any TAO command line options  */
261   CHKERRQ(TaoSetTolerances(tao,1e-8,PETSC_DEFAULT,PETSC_DEFAULT));
262   CHKERRQ(TaoSetFromOptions(tao));
263   CHKERRQ(TaoSolve(tao));
264 
265   CHKERRQ(TaoDestroy(&tao));
266   CHKERRQ(MatDestroy(&appctx.SEMop.stiff));
267   CHKERRQ(MatDestroy(&appctx.SEMop.keptstiff));
268   CHKERRQ(MatDestroy(&appctx.SEMop.grad));
269   CHKERRQ(VecDestroy(&u));
270   CHKERRQ(VecDestroy(&appctx.dat.ic));
271   CHKERRQ(VecDestroy(&appctx.dat.true_solution));
272   CHKERRQ(VecDestroy(&appctx.dat.obj));
273   CHKERRQ(VecDestroy(&appctx.SEMop.grid));
274   CHKERRQ(VecDestroy(&appctx.SEMop.mass));
275   CHKERRQ(VecDestroy(&appctx.dat.curr_sol));
276   CHKERRQ(PetscFree2(appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights));
277   CHKERRQ(DMDestroy(&appctx.da));
278   CHKERRQ(TSDestroy(&appctx.ts));
279 
280   /*
281      Always call PetscFinalize() before exiting a program.  This routine
282        - finalizes the PETSc libraries as well as MPI
283        - provides summary and diagnostic information if certain runtime
284          options are chosen (e.g., -log_summary).
285   */
286   ierr = PetscFinalize();
287   return ierr;
288 }
289 
290 /* --------------------------------------------------------------------- */
291 /*
292    InitialConditions - Computes the initial conditions for the Tao optimization solve (these are also initial conditions for the first TSSolve()
293 
294                        The routine TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function
295 
296    Input Parameter:
297    u - uninitialized solution vector (global)
298    appctx - user-defined application context
299 
300    Output Parameter:
301    u - vector with solution at initial time (global)
302 */
303 PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
304 {
305   PetscScalar       *s;
306   const PetscScalar *xg;
307   PetscInt          i,xs,xn;
308 
309   PetscFunctionBegin;
310   CHKERRQ(DMDAVecGetArray(appctx->da,u,&s));
311   CHKERRQ(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
312   CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL));
313   for (i=xs; i<xs+xn; i++) {
314     s[i]=2.0*appctx->param.mu*PETSC_PI*PetscSinScalar(PETSC_PI*xg[i])/(2.0+PetscCosScalar(PETSC_PI*xg[i]))+0.25*PetscExpReal(-4.0*PetscPowReal(xg[i]-2.0,2.0));
315   }
316   CHKERRQ(DMDAVecRestoreArray(appctx->da,u,&s));
317   CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
318   PetscFunctionReturn(0);
319 }
320 
321 /*
322    TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function.
323 
324              InitialConditions() computes the initial conditions for the beginning of the Tao iterations
325 
326    Input Parameter:
327    u - uninitialized solution vector (global)
328    appctx - user-defined application context
329 
330    Output Parameter:
331    u - vector with solution at initial time (global)
332 */
333 PetscErrorCode TrueSolution(Vec u,AppCtx *appctx)
334 {
335   PetscScalar       *s;
336   const PetscScalar *xg;
337   PetscInt          i,xs,xn;
338 
339   PetscFunctionBegin;
340   CHKERRQ(DMDAVecGetArray(appctx->da,u,&s));
341   CHKERRQ(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
342   CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL));
343   for (i=xs; i<xs+xn; i++) {
344     s[i]=2.0*appctx->param.mu*PETSC_PI*PetscSinScalar(PETSC_PI*xg[i])/(2.0+PetscCosScalar(PETSC_PI*xg[i]));
345   }
346   CHKERRQ(DMDAVecRestoreArray(appctx->da,u,&s));
347   CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
348   PetscFunctionReturn(0);
349 }
350 /* --------------------------------------------------------------------- */
351 /*
352    Sets the desired profile for the final end time
353 
354    Input Parameters:
355    t - final time
356    obj - vector storing the desired profile
357    appctx - user-defined application context
358 
359 */
360 PetscErrorCode ComputeObjective(PetscReal t,Vec obj,AppCtx *appctx)
361 {
362   PetscScalar       *s;
363   const PetscScalar *xg;
364   PetscInt          i, xs,xn;
365 
366   PetscFunctionBegin;
367   CHKERRQ(DMDAVecGetArray(appctx->da,obj,&s));
368   CHKERRQ(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
369   CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL));
370   for (i=xs; i<xs+xn; i++) {
371     s[i]=2.0*appctx->param.mu*PETSC_PI*PetscSinScalar(PETSC_PI*xg[i])*PetscExpScalar(-PETSC_PI*PETSC_PI*t*appctx->param.mu)\
372               /(2.0+PetscExpScalar(-PETSC_PI*PETSC_PI*t*appctx->param.mu)*PetscCosScalar(PETSC_PI*xg[i]));
373   }
374   CHKERRQ(DMDAVecRestoreArray(appctx->da,obj,&s));
375   CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
376   PetscFunctionReturn(0);
377 }
378 
379 PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
380 {
381   AppCtx          *appctx = (AppCtx*)ctx;
382 
383   PetscFunctionBegin;
384   CHKERRQ(MatMult(appctx->SEMop.grad,globalin,globalout)); /* grad u */
385   CHKERRQ(VecPointwiseMult(globalout,globalin,globalout)); /* u grad u */
386   CHKERRQ(VecScale(globalout, -1.0));
387   CHKERRQ(MatMultAdd(appctx->SEMop.keptstiff,globalin,globalout,globalout));
388   PetscFunctionReturn(0);
389 }
390 
391 /*
392 
393       K is the discretiziation of the Laplacian
394       G is the discretization of the gradient
395 
396       Computes Jacobian of      K u + diag(u) G u   which is given by
397               K   + diag(u)G + diag(Gu)
398 */
399 PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A, Mat B,void *ctx)
400 {
401   AppCtx         *appctx = (AppCtx*)ctx;
402   Vec            Gglobalin;
403 
404   PetscFunctionBegin;
405   /*    A = diag(u) G */
406 
407   CHKERRQ(MatCopy(appctx->SEMop.grad,A,SAME_NONZERO_PATTERN));
408   CHKERRQ(MatDiagonalScale(A,globalin,NULL));
409 
410   /*    A  = A + diag(Gu) */
411   CHKERRQ(VecDuplicate(globalin,&Gglobalin));
412   CHKERRQ(MatMult(appctx->SEMop.grad,globalin,Gglobalin));
413   CHKERRQ(MatDiagonalSet(A,Gglobalin,ADD_VALUES));
414   CHKERRQ(VecDestroy(&Gglobalin));
415 
416   /*   A  = K - A    */
417   CHKERRQ(MatScale(A,-1.0));
418   CHKERRQ(MatAXPY(A,1.0,appctx->SEMop.keptstiff,SAME_NONZERO_PATTERN));
419   PetscFunctionReturn(0);
420 }
421 
422 /* --------------------------------------------------------------------- */
423 
424 /*
425    RHSMatrixLaplacian - User-provided routine to compute the right-hand-side
426    matrix for the heat equation.
427 
428    Input Parameters:
429    ts - the TS context
430    t - current time  (ignored)
431    X - current solution (ignored)
432    dummy - optional user-defined context, as set by TSetRHSJacobian()
433 
434    Output Parameters:
435    AA - Jacobian matrix
436    BB - optionally different matrix from which the preconditioner is built
437    str - flag indicating matrix structure
438 
439 */
440 PetscErrorCode RHSMatrixLaplaciangllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx)
441 {
442   PetscReal      **temp;
443   PetscReal      vv;
444   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
445   PetscInt       i,xs,xn,l,j;
446   PetscInt       *rowsDM;
447 
448   PetscFunctionBegin;
449   /*
450    Creates the element stiffness matrix for the given gll
451    */
452   CHKERRQ(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
453   /* workaround for clang analyzer warning: Division by zero */
454   PetscCheck(appctx->param.N > 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Spectral element order should be > 1");
455 
456   /* scale by the size of the element */
457   for (i=0; i<appctx->param.N; i++) {
458     vv=-appctx->param.mu*2.0/appctx->param.Le;
459     for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv;
460   }
461 
462   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
463   CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL));
464 
465   xs   = xs/(appctx->param.N-1);
466   xn   = xn/(appctx->param.N-1);
467 
468   CHKERRQ(PetscMalloc1(appctx->param.N,&rowsDM));
469   /*
470    loop over local elements
471    */
472   for (j=xs; j<xs+xn; j++) {
473     for (l=0; l<appctx->param.N; l++) {
474       rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l;
475     }
476     CHKERRQ(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES));
477   }
478   CHKERRQ(PetscFree(rowsDM));
479   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
480   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
481   CHKERRQ(VecReciprocal(appctx->SEMop.mass));
482   CHKERRQ(MatDiagonalScale(A,appctx->SEMop.mass,0));
483   CHKERRQ(VecReciprocal(appctx->SEMop.mass));
484 
485   CHKERRQ(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
486   PetscFunctionReturn(0);
487 }
488 
489 /*
490    RHSMatrixAdvection - User-provided routine to compute the right-hand-side
491    matrix for the Advection equation.
492 
493    Input Parameters:
494    ts - the TS context
495    t - current time
496    global_in - global input vector
497    dummy - optional user-defined context, as set by TSetRHSJacobian()
498 
499    Output Parameters:
500    AA - Jacobian matrix
501    BB - optionally different preconditioning matrix
502    str - flag indicating matrix structure
503 
504 */
505 PetscErrorCode RHSMatrixAdvectiongllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx)
506 {
507   PetscReal      **temp;
508   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
509   PetscInt       xs,xn,l,j;
510   PetscInt       *rowsDM;
511 
512   PetscFunctionBegin;
513   /*
514    Creates the advection matrix for the given gll
515    */
516   CHKERRQ(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
517   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
518 
519   CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL));
520 
521   xs   = xs/(appctx->param.N-1);
522   xn   = xn/(appctx->param.N-1);
523 
524   CHKERRQ(PetscMalloc1(appctx->param.N,&rowsDM));
525   for (j=xs; j<xs+xn; j++) {
526     for (l=0; l<appctx->param.N; l++) {
527       rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l;
528     }
529     CHKERRQ(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES));
530   }
531   CHKERRQ(PetscFree(rowsDM));
532   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
533   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
534 
535   CHKERRQ(VecReciprocal(appctx->SEMop.mass));
536   CHKERRQ(MatDiagonalScale(A,appctx->SEMop.mass,0));
537   CHKERRQ(VecReciprocal(appctx->SEMop.mass));
538   CHKERRQ(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
539   PetscFunctionReturn(0);
540 }
541 /* ------------------------------------------------------------------ */
542 /*
543    FormFunctionGradient - Evaluates the function and corresponding gradient.
544 
545    Input Parameters:
546    tao - the Tao context
547    IC   - the input vector
548    ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradient()
549 
550    Output Parameters:
551    f   - the newly evaluated function
552    G   - the newly evaluated gradient
553 
554    Notes:
555 
556           The forward equation is
557               M u_t = F(U)
558           which is converted to
559                 u_t = M^{-1} F(u)
560           in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is
561                  M^{-1} J
562           where J is the Jacobian of F. Now the adjoint equation is
563                 M v_t = J^T v
564           but TSAdjoint does not solve this since it can only solve the transposed system for the
565           Jacobian the user provided. Hence TSAdjoint solves
566                  w_t = J^T M^{-1} w  (where w = M v)
567           since there is no way to indicate the mass matrix as a separate entity to TS. Thus one
568           must be careful in initializing the "adjoint equation" and using the result. This is
569           why
570               G = -2 M(u(T) - u_d)
571           below (instead of -2(u(T) - u_d) and why the result is
572               G = G/appctx->SEMop.mass (that is G = M^{-1}w)
573           below (instead of just the result of the "adjoint solve").
574 
575 */
576 PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
577 {
578   AppCtx             *appctx = (AppCtx*)ctx;     /* user-defined application context */
579   Vec                temp;
580   PetscInt           its;
581   PetscReal          ff, gnorm, cnorm, xdiff,errex;
582   TaoConvergedReason reason;
583 
584   PetscFunctionBegin;
585   CHKERRQ(TSSetTime(appctx->ts,0.0));
586   CHKERRQ(TSSetStepNumber(appctx->ts,0));
587   CHKERRQ(TSSetTimeStep(appctx->ts,appctx->initial_dt));
588   CHKERRQ(VecCopy(IC,appctx->dat.curr_sol));
589 
590   CHKERRQ(TSSolve(appctx->ts,appctx->dat.curr_sol));
591 
592   CHKERRQ(VecWAXPY(G,-1.0,appctx->dat.curr_sol,appctx->dat.obj));
593 
594   /*
595      Compute the L2-norm of the objective function, cost function is f
596   */
597   CHKERRQ(VecDuplicate(G,&temp));
598   CHKERRQ(VecPointwiseMult(temp,G,G));
599   CHKERRQ(VecDot(temp,appctx->SEMop.mass,f));
600 
601   /* local error evaluation   */
602   CHKERRQ(VecWAXPY(temp,-1.0,appctx->dat.ic,appctx->dat.true_solution));
603   CHKERRQ(VecPointwiseMult(temp,temp,temp));
604   /* for error evaluation */
605   CHKERRQ(VecDot(temp,appctx->SEMop.mass,&errex));
606   CHKERRQ(VecDestroy(&temp));
607   errex  = PetscSqrtReal(errex);
608 
609   /*
610      Compute initial conditions for the adjoint integration. See Notes above
611   */
612 
613   CHKERRQ(VecScale(G, -2.0));
614   CHKERRQ(VecPointwiseMult(G,G,appctx->SEMop.mass));
615   CHKERRQ(TSSetCostGradients(appctx->ts,1,&G,NULL));
616   CHKERRQ(TSAdjointSolve(appctx->ts));
617   CHKERRQ(VecPointwiseDivide(G,G,appctx->SEMop.mass));
618 
619   CHKERRQ(TaoGetSolutionStatus(tao, &its, &ff, &gnorm, &cnorm, &xdiff, &reason));
620   PetscFunctionReturn(0);
621 }
622 
623 PetscErrorCode MonitorError(Tao tao,void *ctx)
624 {
625   AppCtx         *appctx = (AppCtx*)ctx;
626   Vec            temp;
627   PetscReal      nrm;
628 
629   PetscFunctionBegin;
630   CHKERRQ(VecDuplicate(appctx->dat.ic,&temp));
631   CHKERRQ(VecWAXPY(temp,-1.0,appctx->dat.ic,appctx->dat.true_solution));
632   CHKERRQ(VecPointwiseMult(temp,temp,temp));
633   CHKERRQ(VecDot(temp,appctx->SEMop.mass,&nrm));
634   CHKERRQ(VecDestroy(&temp));
635   nrm  = PetscSqrtReal(nrm);
636   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error for initial conditions %g\n",(double)nrm));
637   PetscFunctionReturn(0);
638 }
639 
640 /*TEST
641 
642     build:
643       requires: !complex
644 
645     test:
646       args: -tao_max_it 5 -tao_gatol 1.e-4
647       requires: !single
648 
649     test:
650       suffix: 2
651       nsize: 2
652       args: -tao_max_it 5 -tao_gatol 1.e-4
653       requires: !single
654 
655 TEST*/
656