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