xref: /petsc/src/tao/unconstrained/tutorials/burgers_spectral.c (revision dd8e379b23d2ef935d8131fb74f7eb73fef09263)
1c4762a1bSJed Brown static char help[] = "Solves a simple data assimilation problem with one dimensional Burger's equation using TSAdjoint\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown 
5c4762a1bSJed Brown     Not yet tested in parallel
6c4762a1bSJed Brown 
7c4762a1bSJed Brown */
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    This program uses the one-dimensional Burger's equation
12c4762a1bSJed Brown        u_t = mu*u_xx - u u_x,
13c4762a1bSJed Brown    on the domain 0 <= x <= 1, with periodic boundary conditions
14c4762a1bSJed Brown 
15c4762a1bSJed Brown    to demonstrate solving a data assimilation problem of finding the initial conditions
16c4762a1bSJed Brown    to produce a given solution at a fixed time.
17c4762a1bSJed Brown 
18c4762a1bSJed Brown    The operators are discretized with the spectral element method
19c4762a1bSJed Brown 
20c4762a1bSJed Brown    See the paper PDE-CONSTRAINED OPTIMIZATION WITH SPECTRAL ELEMENTS USING PETSC AND TAO
21c4762a1bSJed Brown    by OANA MARIN, EMIL CONSTANTINESCU, AND BARRY SMITH for details on the exact solution
22c4762a1bSJed Brown    used
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   ------------------------------------------------------------------------- */
25c4762a1bSJed Brown 
26c4762a1bSJed Brown #include <petsctao.h>
27c4762a1bSJed Brown #include <petscts.h>
28c4762a1bSJed Brown #include <petscdt.h>
29c4762a1bSJed Brown #include <petscdraw.h>
30c4762a1bSJed Brown #include <petscdmda.h>
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown    User-defined application context - contains data needed by the
34c4762a1bSJed Brown    application-provided call-back routines.
35c4762a1bSJed Brown */
36c4762a1bSJed Brown 
37c4762a1bSJed Brown typedef struct {
38c4762a1bSJed Brown   PetscInt   n;       /* number of nodes */
39c4762a1bSJed Brown   PetscReal *nodes;   /* GLL nodes */
40c4762a1bSJed Brown   PetscReal *weights; /* GLL weights */
41c4762a1bSJed Brown } PetscGLL;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown typedef struct {
44c4762a1bSJed Brown   PetscInt  N;               /* grid points per elements*/
45c4762a1bSJed Brown   PetscInt  E;               /* number of elements */
46c4762a1bSJed Brown   PetscReal tol_L2, tol_max; /* error norms */
47c4762a1bSJed Brown   PetscInt  steps;           /* number of timesteps */
48c4762a1bSJed Brown   PetscReal Tend;            /* endtime */
49c4762a1bSJed Brown   PetscReal mu;              /* viscosity */
50c4762a1bSJed Brown   PetscReal L;               /* total length of domain */
51c4762a1bSJed Brown   PetscReal Le;
52c4762a1bSJed Brown   PetscReal Tadj;
53c4762a1bSJed Brown } PetscParam;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown typedef struct {
56c4762a1bSJed Brown   Vec obj;  /* desired end state */
57c4762a1bSJed Brown   Vec grid; /* total grid */
58c4762a1bSJed Brown   Vec grad;
59c4762a1bSJed Brown   Vec ic;
60c4762a1bSJed Brown   Vec curr_sol;
61c4762a1bSJed Brown   Vec true_solution; /* actual initial conditions for the final solution */
62c4762a1bSJed Brown } PetscData;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown typedef struct {
65c4762a1bSJed Brown   Vec      grid;  /* total grid */
66c4762a1bSJed Brown   Vec      mass;  /* mass matrix for total integration */
67c4762a1bSJed Brown   Mat      stiff; /* stifness matrix */
68c4762a1bSJed Brown   Mat      keptstiff;
69c4762a1bSJed Brown   Mat      grad;
70c4762a1bSJed Brown   PetscGLL gll;
71c4762a1bSJed Brown } PetscSEMOperators;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown typedef struct {
74c4762a1bSJed Brown   DM                da; /* distributed array data structure */
75c4762a1bSJed Brown   PetscSEMOperators SEMop;
76c4762a1bSJed Brown   PetscParam        param;
77c4762a1bSJed Brown   PetscData         dat;
78c4762a1bSJed Brown   TS                ts;
79c4762a1bSJed Brown   PetscReal         initial_dt;
80c4762a1bSJed Brown } AppCtx;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown /*
83c4762a1bSJed Brown    User-defined routines
84c4762a1bSJed Brown */
85c4762a1bSJed Brown extern PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
86c4762a1bSJed Brown extern PetscErrorCode RHSMatrixLaplaciangllDM(TS, PetscReal, Vec, Mat, Mat, void *);
87c4762a1bSJed Brown extern PetscErrorCode RHSMatrixAdvectiongllDM(TS, PetscReal, Vec, Mat, Mat, void *);
88c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec, AppCtx *);
89c4762a1bSJed Brown extern PetscErrorCode TrueSolution(Vec, AppCtx *);
90c4762a1bSJed Brown extern PetscErrorCode ComputeObjective(PetscReal, Vec, AppCtx *);
91c4762a1bSJed Brown extern PetscErrorCode MonitorError(Tao, void *);
92c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
93c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
94c4762a1bSJed Brown 
95d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
96d71ae5a4SJacob Faibussowitsch {
97c4762a1bSJed Brown   AppCtx       appctx; /* user-defined application context */
98c4762a1bSJed Brown   Tao          tao;
99c4762a1bSJed Brown   Vec          u; /* approximate solution vector */
100c4762a1bSJed Brown   PetscInt     i, xs, xm, ind, j, lenglob;
101c4762a1bSJed Brown   PetscReal    x, *wrk_ptr1, *wrk_ptr2;
102c4762a1bSJed Brown   MatNullSpace nsp;
103c4762a1bSJed Brown   PetscMPIInt  size;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown      Initialize program and set problem parameters
107c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108327415f7SBarry Smith   PetscFunctionBeginUser;
1099566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /*initialize parameters */
112c4762a1bSJed Brown   appctx.param.N     = 10;   /* order of the spectral element */
113c4762a1bSJed Brown   appctx.param.E     = 10;   /* number of elements */
114c4762a1bSJed Brown   appctx.param.L     = 4.0;  /* length of the domain */
115c4762a1bSJed Brown   appctx.param.mu    = 0.01; /* diffusion coefficient */
116c4762a1bSJed Brown   appctx.initial_dt  = 5e-3;
117c4762a1bSJed Brown   appctx.param.steps = PETSC_MAX_INT;
118c4762a1bSJed Brown   appctx.param.Tend  = 4;
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &appctx.param.N, NULL));
1219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-E", &appctx.param.E, NULL));
1229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-Tend", &appctx.param.Tend, NULL));
1239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &appctx.param.mu, NULL));
124c4762a1bSJed Brown   appctx.param.Le = appctx.param.L / appctx.param.E;
125c4762a1bSJed Brown 
1269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1273c859ba3SBarry Smith   PetscCheck((appctx.param.E % size) == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of elements must be divisible by number of processes");
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown      Create GLL data structures
131c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(appctx.param.N, &appctx.SEMop.gll.nodes, appctx.param.N, &appctx.SEMop.gll.weights));
1339566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights));
134c4762a1bSJed Brown   appctx.SEMop.gll.n = appctx.param.N;
135c4762a1bSJed Brown   lenglob            = appctx.param.E * (appctx.param.N - 1);
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
139c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are E*(Nl-1)+1
140c4762a1bSJed Brown      total grid values spread equally among all the processors, except first and last
141c4762a1bSJed Brown   */
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, lenglob, 1, 1, NULL, &appctx.da));
1449566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(appctx.da));
1459566063dSJacob Faibussowitsch   PetscCall(DMSetUp(appctx.da));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /*
148c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
149c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
150c4762a1bSJed Brown      have the same types.
151c4762a1bSJed Brown   */
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(appctx.da, &u));
1549566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.dat.ic));
1559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.dat.true_solution));
1569566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.dat.obj));
1579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.SEMop.grid));
1589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.SEMop.mass));
1599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &appctx.dat.curr_sol));
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx.da, &xs, NULL, NULL, &xm, NULL, NULL));
1629566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1));
1639566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   xs = xs / (appctx.param.N - 1);
168c4762a1bSJed Brown   xm = xm / (appctx.param.N - 1);
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /*
171c4762a1bSJed Brown      Build total grid and mass over entire mesh (multi-elemental)
172c4762a1bSJed Brown   */
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
175c4762a1bSJed Brown     for (j = 0; j < appctx.param.N - 1; j++) {
176c4762a1bSJed Brown       x             = (appctx.param.Le / 2.0) * (appctx.SEMop.gll.nodes[j] + 1.0) + appctx.param.Le * i;
177c4762a1bSJed Brown       ind           = i * (appctx.param.N - 1) + j;
178c4762a1bSJed Brown       wrk_ptr1[ind] = x;
179c4762a1bSJed Brown       wrk_ptr2[ind] = .5 * appctx.param.Le * appctx.SEMop.gll.weights[j];
180c4762a1bSJed Brown       if (j == 0) wrk_ptr2[ind] += .5 * appctx.param.Le * appctx.SEMop.gll.weights[j];
181c4762a1bSJed Brown     }
182c4762a1bSJed Brown   }
1839566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1));
1849566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown    Create matrix data structure; set matrix evaluation routine.
188c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1899566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE));
1909566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.stiff));
1919566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.grad));
192c4762a1bSJed Brown   /*
193c4762a1bSJed Brown    For linear problems with a time-dependent f(u,t) in the equation
194*dd8e379bSPierre Jolivet    u_t = f(u,t), the user provides the discretized right-hand side
195c4762a1bSJed Brown    as a time-dependent matrix.
196c4762a1bSJed Brown    */
1979566063dSJacob Faibussowitsch   PetscCall(RHSMatrixLaplaciangllDM(appctx.ts, 0.0, u, appctx.SEMop.stiff, appctx.SEMop.stiff, &appctx));
1989566063dSJacob Faibussowitsch   PetscCall(RHSMatrixAdvectiongllDM(appctx.ts, 0.0, u, appctx.SEMop.grad, appctx.SEMop.grad, &appctx));
199c4762a1bSJed Brown   /*
200c4762a1bSJed Brown        For linear problems with a time-dependent f(u,t) in the equation
201*dd8e379bSPierre Jolivet        u_t = f(u,t), the user provides the discretized right-hand side
202c4762a1bSJed Brown        as a time-dependent matrix.
203c4762a1bSJed Brown     */
204c4762a1bSJed Brown 
2059566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(appctx.SEMop.stiff, MAT_COPY_VALUES, &appctx.SEMop.keptstiff));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /* attach the null space to the matrix, this probably is not needed but does no harm */
2089566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp));
2099566063dSJacob Faibussowitsch   PetscCall(MatSetNullSpace(appctx.SEMop.stiff, nsp));
2109566063dSJacob Faibussowitsch   PetscCall(MatSetNullSpace(appctx.SEMop.keptstiff, nsp));
2119566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.stiff, NULL));
2129566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&nsp));
213c4762a1bSJed Brown   /* attach the null space to the matrix, this probably is not needed but does no harm */
2149566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp));
2159566063dSJacob Faibussowitsch   PetscCall(MatSetNullSpace(appctx.SEMop.grad, nsp));
2169566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.grad, NULL));
2179566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&nsp));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* Create the TS solver that solves the ODE and its adjoint; set its options */
2209566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts));
2219566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(appctx.ts, TS_NONLINEAR));
2229566063dSJacob Faibussowitsch   PetscCall(TSSetType(appctx.ts, TSRK));
2239566063dSJacob Faibussowitsch   PetscCall(TSSetDM(appctx.ts, appctx.da));
2249566063dSJacob Faibussowitsch   PetscCall(TSSetTime(appctx.ts, 0.0));
2259566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(appctx.ts, appctx.initial_dt));
2269566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(appctx.ts, appctx.param.steps));
2279566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(appctx.ts, appctx.param.Tend));
2289566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP));
2299566063dSJacob Faibussowitsch   PetscCall(TSSetTolerances(appctx.ts, 1e-7, NULL, 1e-7, NULL));
2309566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(appctx.ts));
231c4762a1bSJed Brown   /* Need to save initial timestep user may have set with -ts_dt so it can be reset for each new TSSolve() */
2329566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(appctx.ts, &appctx.initial_dt));
2339566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(appctx.ts, NULL, RHSFunction, &appctx));
2349566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(appctx.ts, appctx.SEMop.stiff, appctx.SEMop.stiff, RHSJacobian, &appctx));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* Set Objective and Initial conditions for the problem and compute Objective function (evolution of true_solution to final time */
2379566063dSJacob Faibussowitsch   PetscCall(InitialConditions(appctx.dat.ic, &appctx));
2389566063dSJacob Faibussowitsch   PetscCall(TrueSolution(appctx.dat.true_solution, &appctx));
2399566063dSJacob Faibussowitsch   PetscCall(ComputeObjective(appctx.param.Tend, appctx.dat.obj, &appctx));
240c4762a1bSJed Brown 
2419566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(appctx.ts));
2429566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(appctx.ts));
243f32d6360SSatish Balay 
244c4762a1bSJed Brown   /* Create TAO solver and set desired solution method  */
2459566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
24610978b7dSBarry Smith   PetscCall(TaoMonitorSet(tao, MonitorError, &appctx, NULL));
2479566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBQNLS));
2489566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, appctx.dat.ic));
249c4762a1bSJed Brown   /* Set routine for function and gradient evaluation  */
2509566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&appctx));
251c4762a1bSJed Brown   /* Check for any TAO command line options  */
2529566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT));
2539566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
2549566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
255c4762a1bSJed Brown 
2569566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
2579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.SEMop.stiff));
2589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.SEMop.keptstiff));
2599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.SEMop.grad));
2609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.ic));
2629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.true_solution));
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.obj));
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.SEMop.grid));
2659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.SEMop.mass));
2669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.curr_sol));
2679566063dSJacob Faibussowitsch   PetscCall(PetscFree2(appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights));
2689566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&appctx.da));
2699566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&appctx.ts));
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /*
272c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
273c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
274c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
275d75802c7SJacob Faibussowitsch          options are chosen (e.g., -log_view).
276c4762a1bSJed Brown   */
2779566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
278b122ec5aSJacob Faibussowitsch   return 0;
279c4762a1bSJed Brown }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown /* --------------------------------------------------------------------- */
282c4762a1bSJed Brown /*
283c4762a1bSJed Brown    InitialConditions - Computes the initial conditions for the Tao optimization solve (these are also initial conditions for the first TSSolve()
284c4762a1bSJed Brown 
285c4762a1bSJed Brown                        The routine TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    Input Parameter:
288c4762a1bSJed Brown    u - uninitialized solution vector (global)
289c4762a1bSJed Brown    appctx - user-defined application context
290c4762a1bSJed Brown 
291c4762a1bSJed Brown    Output Parameter:
292c4762a1bSJed Brown    u - vector with solution at initial time (global)
293c4762a1bSJed Brown */
294d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(Vec u, AppCtx *appctx)
295d71ae5a4SJacob Faibussowitsch {
296c4762a1bSJed Brown   PetscScalar       *s;
297c4762a1bSJed Brown   const PetscScalar *xg;
298c4762a1bSJed Brown   PetscInt           i, xs, xn;
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   PetscFunctionBegin;
3019566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx->da, u, &s));
3029566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
3039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
304ad540459SPierre Jolivet   for (i = xs; i < xs + xn; i++) 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));
3059566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx->da, u, &s));
3069566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
3073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
308c4762a1bSJed Brown }
309c4762a1bSJed Brown 
310c4762a1bSJed Brown /*
311c4762a1bSJed Brown    TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function.
312c4762a1bSJed Brown 
313a5b23f4aSJose E. Roman              InitialConditions() computes the initial conditions for the beginning of the Tao iterations
314c4762a1bSJed Brown 
315c4762a1bSJed Brown    Input Parameter:
316c4762a1bSJed Brown    u - uninitialized solution vector (global)
317c4762a1bSJed Brown    appctx - user-defined application context
318c4762a1bSJed Brown 
319c4762a1bSJed Brown    Output Parameter:
320c4762a1bSJed Brown    u - vector with solution at initial time (global)
321c4762a1bSJed Brown */
322d71ae5a4SJacob Faibussowitsch PetscErrorCode TrueSolution(Vec u, AppCtx *appctx)
323d71ae5a4SJacob Faibussowitsch {
324c4762a1bSJed Brown   PetscScalar       *s;
325c4762a1bSJed Brown   const PetscScalar *xg;
326c4762a1bSJed Brown   PetscInt           i, xs, xn;
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   PetscFunctionBegin;
3299566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx->da, u, &s));
3309566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
3319566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
332ad540459SPierre Jolivet   for (i = xs; i < xs + xn; i++) s[i] = 2.0 * appctx->param.mu * PETSC_PI * PetscSinScalar(PETSC_PI * xg[i]) / (2.0 + PetscCosScalar(PETSC_PI * xg[i]));
3339566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx->da, u, &s));
3349566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
336c4762a1bSJed Brown }
337c4762a1bSJed Brown /* --------------------------------------------------------------------- */
338c4762a1bSJed Brown /*
339c4762a1bSJed Brown    Sets the desired profile for the final end time
340c4762a1bSJed Brown 
341c4762a1bSJed Brown    Input Parameters:
342c4762a1bSJed Brown    t - final time
343c4762a1bSJed Brown    obj - vector storing the desired profile
344c4762a1bSJed Brown    appctx - user-defined application context
345c4762a1bSJed Brown 
346c4762a1bSJed Brown */
347d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeObjective(PetscReal t, Vec obj, AppCtx *appctx)
348d71ae5a4SJacob Faibussowitsch {
349c4762a1bSJed Brown   PetscScalar       *s;
350c4762a1bSJed Brown   const PetscScalar *xg;
351c4762a1bSJed Brown   PetscInt           i, xs, xn;
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   PetscFunctionBegin;
3549566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx->da, obj, &s));
3559566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
3569566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
357c4762a1bSJed Brown   for (i = xs; i < xs + xn; i++) {
3589371c9d4SSatish Balay     s[i] = 2.0 * appctx->param.mu * PETSC_PI * PetscSinScalar(PETSC_PI * xg[i]) * PetscExpScalar(-PETSC_PI * PETSC_PI * t * appctx->param.mu) / (2.0 + PetscExpScalar(-PETSC_PI * PETSC_PI * t * appctx->param.mu) * PetscCosScalar(PETSC_PI * xg[i]));
359c4762a1bSJed Brown   }
3609566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx->da, obj, &s));
3619566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
363c4762a1bSJed Brown }
364c4762a1bSJed Brown 
365d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx)
366d71ae5a4SJacob Faibussowitsch {
367c4762a1bSJed Brown   AppCtx *appctx = (AppCtx *)ctx;
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   PetscFunctionBegin;
3709566063dSJacob Faibussowitsch   PetscCall(MatMult(appctx->SEMop.grad, globalin, globalout)); /* grad u */
3719566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(globalout, globalin, globalout)); /* u grad u */
3729566063dSJacob Faibussowitsch   PetscCall(VecScale(globalout, -1.0));
3739566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(appctx->SEMop.keptstiff, globalin, globalout, globalout));
3743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown /*
378c4762a1bSJed Brown 
379c4762a1bSJed Brown       K is the discretiziation of the Laplacian
380c4762a1bSJed Brown       G is the discretization of the gradient
381c4762a1bSJed Brown 
382c4762a1bSJed Brown       Computes Jacobian of      K u + diag(u) G u   which is given by
383c4762a1bSJed Brown               K   + diag(u)G + diag(Gu)
384c4762a1bSJed Brown */
385d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec globalin, Mat A, Mat B, void *ctx)
386d71ae5a4SJacob Faibussowitsch {
387c4762a1bSJed Brown   AppCtx *appctx = (AppCtx *)ctx;
388c4762a1bSJed Brown   Vec     Gglobalin;
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   PetscFunctionBegin;
391c4762a1bSJed Brown   /*    A = diag(u) G */
392c4762a1bSJed Brown 
3939566063dSJacob Faibussowitsch   PetscCall(MatCopy(appctx->SEMop.grad, A, SAME_NONZERO_PATTERN));
3949566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A, globalin, NULL));
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   /*    A  = A + diag(Gu) */
3979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(globalin, &Gglobalin));
3989566063dSJacob Faibussowitsch   PetscCall(MatMult(appctx->SEMop.grad, globalin, Gglobalin));
3999566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(A, Gglobalin, ADD_VALUES));
4009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Gglobalin));
401c4762a1bSJed Brown 
402c4762a1bSJed Brown   /*   A  = K - A    */
4039566063dSJacob Faibussowitsch   PetscCall(MatScale(A, -1.0));
4049566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A, 1.0, appctx->SEMop.keptstiff, SAME_NONZERO_PATTERN));
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
406c4762a1bSJed Brown }
407c4762a1bSJed Brown 
408c4762a1bSJed Brown /* --------------------------------------------------------------------- */
409c4762a1bSJed Brown 
410c4762a1bSJed Brown /*
411c4762a1bSJed Brown    RHSMatrixLaplacian - User-provided routine to compute the right-hand-side
412c4762a1bSJed Brown    matrix for the heat equation.
413c4762a1bSJed Brown 
414c4762a1bSJed Brown    Input Parameters:
415c4762a1bSJed Brown    ts - the TS context
416c4762a1bSJed Brown    t - current time  (ignored)
417c4762a1bSJed Brown    X - current solution (ignored)
418c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
419c4762a1bSJed Brown 
420c4762a1bSJed Brown    Output Parameters:
421c4762a1bSJed Brown    AA - Jacobian matrix
422c4762a1bSJed Brown    BB - optionally different matrix from which the preconditioner is built
423c4762a1bSJed Brown    str - flag indicating matrix structure
424c4762a1bSJed Brown 
425c4762a1bSJed Brown */
426d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixLaplaciangllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx)
427d71ae5a4SJacob Faibussowitsch {
428c4762a1bSJed Brown   PetscReal **temp;
429c4762a1bSJed Brown   PetscReal   vv;
430c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
431c4762a1bSJed Brown   PetscInt    i, xs, xn, l, j;
432c4762a1bSJed Brown   PetscInt   *rowsDM;
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   PetscFunctionBegin;
435c4762a1bSJed Brown   /*
436c4762a1bSJed Brown    Creates the element stiffness matrix for the given gll
437c4762a1bSJed Brown    */
4389566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
439a5b23f4aSJose E. Roman   /* workaround for clang analyzer warning: Division by zero */
4403c859ba3SBarry Smith   PetscCheck(appctx->param.N > 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spectral element order should be > 1");
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   /* scale by the size of the element */
443c4762a1bSJed Brown   for (i = 0; i < appctx->param.N; i++) {
444c4762a1bSJed Brown     vv = -appctx->param.mu * 2.0 / appctx->param.Le;
445c4762a1bSJed Brown     for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv;
446c4762a1bSJed Brown   }
447c4762a1bSJed Brown 
4489566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
4499566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   xs = xs / (appctx->param.N - 1);
452c4762a1bSJed Brown   xn = xn / (appctx->param.N - 1);
453c4762a1bSJed Brown 
4549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(appctx->param.N, &rowsDM));
455c4762a1bSJed Brown   /*
456c4762a1bSJed Brown    loop over local elements
457c4762a1bSJed Brown    */
458c4762a1bSJed Brown   for (j = xs; j < xs + xn; j++) {
459ad540459SPierre Jolivet     for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l;
4609566063dSJacob Faibussowitsch     PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES));
461c4762a1bSJed Brown   }
4629566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowsDM));
4639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4659566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
4669566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0));
4679566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
468c4762a1bSJed Brown 
4699566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471c4762a1bSJed Brown }
472c4762a1bSJed Brown 
473c4762a1bSJed Brown /*
474c4762a1bSJed Brown    RHSMatrixAdvection - User-provided routine to compute the right-hand-side
475c4762a1bSJed Brown    matrix for the Advection equation.
476c4762a1bSJed Brown 
477c4762a1bSJed Brown    Input Parameters:
478c4762a1bSJed Brown    ts - the TS context
479c4762a1bSJed Brown    t - current time
480c4762a1bSJed Brown    global_in - global input vector
481c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
482c4762a1bSJed Brown 
483c4762a1bSJed Brown    Output Parameters:
484c4762a1bSJed Brown    AA - Jacobian matrix
485c4762a1bSJed Brown    BB - optionally different preconditioning matrix
486c4762a1bSJed Brown    str - flag indicating matrix structure
487c4762a1bSJed Brown 
488c4762a1bSJed Brown */
489d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixAdvectiongllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx)
490d71ae5a4SJacob Faibussowitsch {
491c4762a1bSJed Brown   PetscReal **temp;
492c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
493c4762a1bSJed Brown   PetscInt    xs, xn, l, j;
494c4762a1bSJed Brown   PetscInt   *rowsDM;
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   PetscFunctionBegin;
497c4762a1bSJed Brown   /*
498c4762a1bSJed Brown    Creates the advection matrix for the given gll
499c4762a1bSJed Brown    */
5009566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
5019566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
502c4762a1bSJed Brown 
5039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
504c4762a1bSJed Brown 
505c4762a1bSJed Brown   xs = xs / (appctx->param.N - 1);
506c4762a1bSJed Brown   xn = xn / (appctx->param.N - 1);
507c4762a1bSJed Brown 
5089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(appctx->param.N, &rowsDM));
509c4762a1bSJed Brown   for (j = xs; j < xs + xn; j++) {
510ad540459SPierre Jolivet     for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l;
5119566063dSJacob Faibussowitsch     PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES));
512c4762a1bSJed Brown   }
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowsDM));
5149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
516c4762a1bSJed Brown 
5179566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
5189566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0));
5199566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
5209566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
5213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
522c4762a1bSJed Brown }
523c4762a1bSJed Brown /* ------------------------------------------------------------------ */
524c4762a1bSJed Brown /*
525c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
526c4762a1bSJed Brown 
527c4762a1bSJed Brown    Input Parameters:
528c4762a1bSJed Brown    tao - the Tao context
529c4762a1bSJed Brown    IC   - the input vector
530a82e8c82SStefano Zampini    ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradient()
531c4762a1bSJed Brown 
532c4762a1bSJed Brown    Output Parameters:
533c4762a1bSJed Brown    f   - the newly evaluated function
534c4762a1bSJed Brown    G   - the newly evaluated gradient
535c4762a1bSJed Brown 
536c4762a1bSJed Brown    Notes:
537c4762a1bSJed Brown 
538c4762a1bSJed Brown           The forward equation is
539c4762a1bSJed Brown               M u_t = F(U)
540c4762a1bSJed Brown           which is converted to
541c4762a1bSJed Brown                 u_t = M^{-1} F(u)
542c4762a1bSJed Brown           in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is
543c4762a1bSJed Brown                  M^{-1} J
544c4762a1bSJed Brown           where J is the Jacobian of F. Now the adjoint equation is
545c4762a1bSJed Brown                 M v_t = J^T v
546c4762a1bSJed Brown           but TSAdjoint does not solve this since it can only solve the transposed system for the
547c4762a1bSJed Brown           Jacobian the user provided. Hence TSAdjoint solves
548c4762a1bSJed Brown                  w_t = J^T M^{-1} w  (where w = M v)
549a5b23f4aSJose E. Roman           since there is no way to indicate the mass matrix as a separate entity to TS. Thus one
550c4762a1bSJed Brown           must be careful in initializing the "adjoint equation" and using the result. This is
551c4762a1bSJed Brown           why
552c4762a1bSJed Brown               G = -2 M(u(T) - u_d)
553c4762a1bSJed Brown           below (instead of -2(u(T) - u_d) and why the result is
554c4762a1bSJed Brown               G = G/appctx->SEMop.mass (that is G = M^{-1}w)
555c4762a1bSJed Brown           below (instead of just the result of the "adjoint solve").
556c4762a1bSJed Brown 
557c4762a1bSJed Brown */
558d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec IC, PetscReal *f, Vec G, void *ctx)
559d71ae5a4SJacob Faibussowitsch {
560c4762a1bSJed Brown   AppCtx            *appctx = (AppCtx *)ctx; /* user-defined application context */
561c4762a1bSJed Brown   Vec                temp;
562c4762a1bSJed Brown   PetscInt           its;
563c4762a1bSJed Brown   PetscReal          ff, gnorm, cnorm, xdiff, errex;
564c4762a1bSJed Brown   TaoConvergedReason reason;
565c4762a1bSJed Brown 
566c4762a1bSJed Brown   PetscFunctionBegin;
5679566063dSJacob Faibussowitsch   PetscCall(TSSetTime(appctx->ts, 0.0));
5689566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(appctx->ts, 0));
5699566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(appctx->ts, appctx->initial_dt));
5709566063dSJacob Faibussowitsch   PetscCall(VecCopy(IC, appctx->dat.curr_sol));
571c4762a1bSJed Brown 
5729566063dSJacob Faibussowitsch   PetscCall(TSSolve(appctx->ts, appctx->dat.curr_sol));
573c4762a1bSJed Brown 
5749566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(G, -1.0, appctx->dat.curr_sol, appctx->dat.obj));
575c4762a1bSJed Brown 
576c4762a1bSJed Brown   /*
577c4762a1bSJed Brown      Compute the L2-norm of the objective function, cost function is f
578c4762a1bSJed Brown   */
5799566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(G, &temp));
5809566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(temp, G, G));
5819566063dSJacob Faibussowitsch   PetscCall(VecDot(temp, appctx->SEMop.mass, f));
582c4762a1bSJed Brown 
583c4762a1bSJed Brown   /* local error evaluation   */
5849566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(temp, -1.0, appctx->dat.ic, appctx->dat.true_solution));
5859566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(temp, temp, temp));
586c4762a1bSJed Brown   /* for error evaluation */
5879566063dSJacob Faibussowitsch   PetscCall(VecDot(temp, appctx->SEMop.mass, &errex));
5889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&temp));
589c4762a1bSJed Brown   errex = PetscSqrtReal(errex);
590c4762a1bSJed Brown 
591c4762a1bSJed Brown   /*
592c4762a1bSJed Brown      Compute initial conditions for the adjoint integration. See Notes above
593c4762a1bSJed Brown   */
594c4762a1bSJed Brown 
5959566063dSJacob Faibussowitsch   PetscCall(VecScale(G, -2.0));
5969566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(G, G, appctx->SEMop.mass));
5979566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(appctx->ts, 1, &G, NULL));
5989566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(appctx->ts));
5999566063dSJacob Faibussowitsch   PetscCall(VecPointwiseDivide(G, G, appctx->SEMop.mass));
600c4762a1bSJed Brown 
6019566063dSJacob Faibussowitsch   PetscCall(TaoGetSolutionStatus(tao, &its, &ff, &gnorm, &cnorm, &xdiff, &reason));
6023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
603c4762a1bSJed Brown }
604c4762a1bSJed Brown 
605d71ae5a4SJacob Faibussowitsch PetscErrorCode MonitorError(Tao tao, void *ctx)
606d71ae5a4SJacob Faibussowitsch {
607c4762a1bSJed Brown   AppCtx   *appctx = (AppCtx *)ctx;
608c4762a1bSJed Brown   Vec       temp;
609c4762a1bSJed Brown   PetscReal nrm;
610c4762a1bSJed Brown 
611c4762a1bSJed Brown   PetscFunctionBegin;
6129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(appctx->dat.ic, &temp));
6139566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(temp, -1.0, appctx->dat.ic, appctx->dat.true_solution));
6149566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(temp, temp, temp));
6159566063dSJacob Faibussowitsch   PetscCall(VecDot(temp, appctx->SEMop.mass, &nrm));
6169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&temp));
617c4762a1bSJed Brown   nrm = PetscSqrtReal(nrm);
6189566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error for initial conditions %g\n", (double)nrm));
6193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
620c4762a1bSJed Brown }
621c4762a1bSJed Brown 
622c4762a1bSJed Brown /*TEST
623c4762a1bSJed Brown 
624c4762a1bSJed Brown     build:
625c4762a1bSJed Brown       requires: !complex
626c4762a1bSJed Brown 
627c4762a1bSJed Brown     test:
628c4762a1bSJed Brown       args: -tao_max_it 5 -tao_gatol 1.e-4
629c4762a1bSJed Brown       requires: !single
630c4762a1bSJed Brown 
631c4762a1bSJed Brown     test:
632c4762a1bSJed Brown       suffix: 2
633c4762a1bSJed Brown       nsize: 2
634c4762a1bSJed Brown       args: -tao_max_it 5 -tao_gatol 1.e-4
635c4762a1bSJed Brown       requires: !single
636c4762a1bSJed Brown 
637c4762a1bSJed Brown TEST*/
638