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