xref: /petsc/src/tao/unconstrained/tutorials/spectraladjointassimilation.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Solves a simple data assimilation problem with one dimensional advection diffusion 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 /* ------------------------------------------------------------------------
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    This program uses the one-dimensional advection-diffusion equation),
13c4762a1bSJed Brown        u_t = mu*u_xx - a u_x,
14c4762a1bSJed Brown    on the domain 0 <= x <= 1, with periodic boundary conditions
15c4762a1bSJed Brown 
16c4762a1bSJed Brown    to demonstrate solving a data assimilation problem of finding the initial conditions
17c4762a1bSJed Brown    to produce a given solution at a fixed time.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown    The operators are discretized with the spectral element method
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   ------------------------------------------------------------------------- */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown /*
24c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this file
25c4762a1bSJed Brown    automatically includes:
26c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h  - vectors
27c4762a1bSJed Brown      petscmat.h  - matrices
28c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
29c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h   - preconditioners
30c4762a1bSJed Brown      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
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   a;              /* advection speed */
58c4762a1bSJed Brown   PetscReal   L;              /* total length of domain */
59c4762a1bSJed Brown   PetscReal   Le;
60c4762a1bSJed Brown   PetscReal   Tadj;
61c4762a1bSJed Brown } PetscParam;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown typedef struct {
64c4762a1bSJed Brown   Vec         reference;               /* desired end state */
65c4762a1bSJed Brown   Vec         grid;              /* total grid */
66c4762a1bSJed Brown   Vec         grad;
67c4762a1bSJed Brown   Vec         ic;
68c4762a1bSJed Brown   Vec         curr_sol;
69c4762a1bSJed Brown   Vec         joe;
70c4762a1bSJed Brown   Vec         true_solution;     /* actual initial conditions for the final solution */
71c4762a1bSJed Brown } PetscData;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown typedef struct {
74c4762a1bSJed Brown   Vec         grid;              /* total grid */
75c4762a1bSJed Brown   Vec         mass;              /* mass matrix for total integration */
76c4762a1bSJed Brown   Mat         stiff;             /* stifness matrix */
77c4762a1bSJed Brown   Mat         advec;
78c4762a1bSJed Brown   Mat         keptstiff;
79c4762a1bSJed Brown   PetscGLL    gll;
80c4762a1bSJed Brown } PetscSEMOperators;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown typedef struct {
83c4762a1bSJed Brown   DM                da;                /* distributed array data structure */
84c4762a1bSJed Brown   PetscSEMOperators SEMop;
85c4762a1bSJed Brown   PetscParam        param;
86c4762a1bSJed Brown   PetscData         dat;
87c4762a1bSJed Brown   TS                ts;
88c4762a1bSJed Brown   PetscReal         initial_dt;
89c4762a1bSJed Brown   PetscReal         *solutioncoefficients;
90c4762a1bSJed Brown   PetscInt          ncoeff;
91c4762a1bSJed Brown } AppCtx;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown /*
94c4762a1bSJed Brown    User-defined routines
95c4762a1bSJed Brown */
96c4762a1bSJed Brown extern PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
97c4762a1bSJed Brown extern PetscErrorCode RHSLaplacian(TS,PetscReal,Vec,Mat,Mat,void*);
98c4762a1bSJed Brown extern PetscErrorCode RHSAdvection(TS,PetscReal,Vec,Mat,Mat,void*);
99c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
100c4762a1bSJed Brown extern PetscErrorCode ComputeReference(TS,PetscReal,Vec,AppCtx*);
101c4762a1bSJed Brown extern PetscErrorCode MonitorError(Tao,void*);
102c4762a1bSJed Brown extern PetscErrorCode MonitorDestroy(void**);
103c4762a1bSJed Brown extern PetscErrorCode ComputeSolutionCoefficients(AppCtx*);
104c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
105c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
106c4762a1bSJed Brown 
107c4762a1bSJed Brown int main(int argc,char **argv)
108c4762a1bSJed Brown {
109c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
110c4762a1bSJed Brown   Tao            tao;
111c4762a1bSJed Brown   Vec            u;                      /* approximate solution vector */
112c4762a1bSJed Brown   PetscInt       i, xs, xm, ind, j, lenglob;
113c4762a1bSJed Brown   PetscReal      x, *wrk_ptr1, *wrk_ptr2;
114c4762a1bSJed Brown   MatNullSpace   nsp;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Initialize program and set problem parameters
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119c4762a1bSJed Brown   PetscFunctionBegin;
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /*initialize parameters */
124c4762a1bSJed Brown   appctx.param.N    = 10;  /* order of the spectral element */
125c4762a1bSJed Brown   appctx.param.E    = 8;  /* number of elements */
126c4762a1bSJed Brown   appctx.param.L    = 1.0;  /* length of the domain */
127c4762a1bSJed Brown   appctx.param.mu   = 0.00001; /* diffusion coefficient */
128c4762a1bSJed Brown   appctx.param.a    = 0.0;     /* advection speed */
129c4762a1bSJed Brown   appctx.initial_dt = 1e-4;
130c4762a1bSJed Brown   appctx.param.steps = PETSC_MAX_INT;
131c4762a1bSJed Brown   appctx.param.Tend  = 0.01;
132c4762a1bSJed Brown   appctx.ncoeff      = 2;
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&appctx.param.N,NULL));
1359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-E",&appctx.param.E,NULL));
1369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ncoeff",&appctx.ncoeff,NULL));
1379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-Tend",&appctx.param.Tend,NULL));
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&appctx.param.mu,NULL));
1399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-a",&appctx.param.a,NULL));
140c4762a1bSJed Brown   appctx.param.Le = appctx.param.L/appctx.param.E;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143c4762a1bSJed Brown      Create GLL data structures
144c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(appctx.param.N,&appctx.SEMop.gll.nodes,appctx.param.N,&appctx.SEMop.gll.weights));
1469566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights));
147c4762a1bSJed Brown   appctx.SEMop.gll.n = appctx.param.N;
148c4762a1bSJed Brown   lenglob  = appctx.param.E*(appctx.param.N-1);
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /*
151c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
152c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are E*(Nl-1)+1
153c4762a1bSJed Brown      total grid values spread equally among all the processors, except first and last
154c4762a1bSJed Brown   */
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,lenglob,1,1,NULL,&appctx.da));
1579566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(appctx.da));
1589566063dSJacob Faibussowitsch   PetscCall(DMSetUp(appctx.da));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /*
161c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
162c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
163c4762a1bSJed Brown      have the same types.
164c4762a1bSJed Brown   */
165c4762a1bSJed Brown 
1669566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(appctx.da,&u));
1679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.dat.ic));
1689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.dat.true_solution));
1699566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.dat.reference));
1709566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.SEMop.grid));
1719566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.SEMop.mass));
1729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.dat.curr_sol));
1739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&appctx.dat.joe));
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx.da,&xs,NULL,NULL,&xm,NULL,NULL));
1769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1));
1779566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
180c4762a1bSJed Brown 
181c4762a1bSJed Brown     xs=xs/(appctx.param.N-1);
182c4762a1bSJed Brown     xm=xm/(appctx.param.N-1);
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /*
185c4762a1bSJed Brown      Build total grid and mass over entire mesh (multi-elemental)
186c4762a1bSJed Brown   */
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
189c4762a1bSJed Brown     for (j=0; j<appctx.param.N-1; j++) {
190c4762a1bSJed Brown       x = (appctx.param.Le/2.0)*(appctx.SEMop.gll.nodes[j]+1.0)+appctx.param.Le*i;
191c4762a1bSJed Brown       ind=i*(appctx.param.N-1)+j;
192c4762a1bSJed Brown       wrk_ptr1[ind]=x;
193c4762a1bSJed Brown       wrk_ptr2[ind]=.5*appctx.param.Le*appctx.SEMop.gll.weights[j];
194c4762a1bSJed Brown       if (j==0) wrk_ptr2[ind]+=.5*appctx.param.Le*appctx.SEMop.gll.weights[j];
195c4762a1bSJed Brown     }
196c4762a1bSJed Brown   }
1979566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1));
1989566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown    Create matrix data structure; set matrix evaluation routine.
202c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2039566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE));
2049566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(appctx.da,&appctx.SEMop.stiff));
2059566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(appctx.da,&appctx.SEMop.advec));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /*
208c4762a1bSJed Brown    For linear problems with a time-dependent f(u,t) in the equation
209c4762a1bSJed Brown    u_t = f(u,t), the user provides the discretized right-hand-side
210c4762a1bSJed Brown    as a time-dependent matrix.
211c4762a1bSJed Brown    */
2129566063dSJacob Faibussowitsch   PetscCall(RHSLaplacian(appctx.ts,0.0,u,appctx.SEMop.stiff,appctx.SEMop.stiff,&appctx));
2139566063dSJacob Faibussowitsch   PetscCall(RHSAdvection(appctx.ts,0.0,u,appctx.SEMop.advec,appctx.SEMop.advec,&appctx));
2149566063dSJacob Faibussowitsch   PetscCall(MatAXPY(appctx.SEMop.stiff,-1.0,appctx.SEMop.advec,DIFFERENT_NONZERO_PATTERN));
2159566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(appctx.SEMop.stiff,MAT_COPY_VALUES,&appctx.SEMop.keptstiff));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* attach the null space to the matrix, this probably is not needed but does no harm */
2189566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp));
2199566063dSJacob Faibussowitsch   PetscCall(MatSetNullSpace(appctx.SEMop.stiff,nsp));
2209566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceTest(nsp,appctx.SEMop.stiff,NULL));
2219566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&nsp));
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* Create the TS solver that solves the ODE and its adjoint; set its options */
2249566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&appctx.ts));
2259566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(appctx.ts,(PetscErrorCode (*)(TS,PetscReal,Vec, void *))ComputeReference,&appctx));
2269566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(appctx.ts,TS_LINEAR));
2279566063dSJacob Faibussowitsch   PetscCall(TSSetType(appctx.ts,TSRK));
2289566063dSJacob Faibussowitsch   PetscCall(TSSetDM(appctx.ts,appctx.da));
2299566063dSJacob Faibussowitsch   PetscCall(TSSetTime(appctx.ts,0.0));
2309566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(appctx.ts,appctx.initial_dt));
2319566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(appctx.ts,appctx.param.steps));
2329566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(appctx.ts,appctx.param.Tend));
2339566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP));
2349566063dSJacob Faibussowitsch   PetscCall(TSSetTolerances(appctx.ts,1e-7,NULL,1e-7,NULL));
2359566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(appctx.ts));
236c4762a1bSJed Brown   /* Need to save initial timestep user may have set with -ts_dt so it can be reset for each new TSSolve() */
2379566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(appctx.ts,&appctx.initial_dt));
2389566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(appctx.ts,NULL,TSComputeRHSFunctionLinear,&appctx));
2399566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,TSComputeRHSJacobianConstant,&appctx));
2409566063dSJacob Faibussowitsch   /*  PetscCall(TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx));
2419566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx)); */
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* Set random initial conditions as initial guess, compute analytic reference solution and analytic (true) initial conditions */
2449566063dSJacob Faibussowitsch   PetscCall(ComputeSolutionCoefficients(&appctx));
2459566063dSJacob Faibussowitsch   PetscCall(InitialConditions(appctx.dat.ic,&appctx));
2469566063dSJacob Faibussowitsch   PetscCall(ComputeReference(appctx.ts,appctx.param.Tend,appctx.dat.reference,&appctx));
2479566063dSJacob Faibussowitsch   PetscCall(ComputeReference(appctx.ts,0.0,appctx.dat.true_solution,&appctx));
248c4762a1bSJed Brown 
249f32d6360SSatish Balay   /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */
2509566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(appctx.ts));
2519566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(appctx.ts));
252f32d6360SSatish Balay 
253c4762a1bSJed Brown   /* Create TAO solver and set desired solution method  */
2549566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
2559566063dSJacob Faibussowitsch   PetscCall(TaoSetMonitor(tao,MonitorError,&appctx,MonitorDestroy));
2569566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBQNLS));
2579566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,appctx.dat.ic));
258c4762a1bSJed Brown   /* Set routine for function and gradient evaluation  */
2599566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&appctx));
260c4762a1bSJed Brown   /* Check for any TAO command line options  */
2619566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao,1e-8,PETSC_DEFAULT,PETSC_DEFAULT));
2629566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
2639566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
264c4762a1bSJed Brown 
2659566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
2669566063dSJacob Faibussowitsch   PetscCall(PetscFree(appctx.solutioncoefficients));
2679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.SEMop.advec));
2689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.SEMop.stiff));
2699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.SEMop.keptstiff));
2709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.ic));
2729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.joe));
2739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.true_solution));
2749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.reference));
2759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.SEMop.grid));
2769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.SEMop.mass));
2779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.dat.curr_sol));
2789566063dSJacob Faibussowitsch   PetscCall(PetscFree2(appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights));
2799566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&appctx.da));
2809566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&appctx.ts));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /*
283c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
284c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
285c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
286c4762a1bSJed Brown          options are chosen (e.g., -log_summary).
287c4762a1bSJed Brown   */
2889566063dSJacob Faibussowitsch     PetscCall(PetscFinalize());
289b122ec5aSJacob Faibussowitsch     return 0;
290c4762a1bSJed Brown }
291c4762a1bSJed Brown 
292c4762a1bSJed Brown /*
293c4762a1bSJed Brown     Computes the coefficients for the analytic solution to the PDE
294c4762a1bSJed Brown */
295c4762a1bSJed Brown PetscErrorCode ComputeSolutionCoefficients(AppCtx *appctx)
296c4762a1bSJed Brown {
297c4762a1bSJed Brown   PetscRandom       rand;
298c4762a1bSJed Brown   PetscInt          i;
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   PetscFunctionBegin;
3019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(appctx->ncoeff,&appctx->solutioncoefficients));
3029566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
3039566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand,.9,1.0));
304c4762a1bSJed Brown   for (i=0; i<appctx->ncoeff; i++) {
3059566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand,&appctx->solutioncoefficients[i]));
306c4762a1bSJed Brown   }
3079566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
308c4762a1bSJed Brown   PetscFunctionReturn(0);
309c4762a1bSJed Brown }
310c4762a1bSJed Brown 
311c4762a1bSJed Brown /* --------------------------------------------------------------------- */
312c4762a1bSJed Brown /*
313c4762a1bSJed Brown    InitialConditions - Computes the (random) initial conditions for the Tao optimization solve (these are also initial conditions for the first TSSolve()
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 */
322c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
323c4762a1bSJed Brown {
324c4762a1bSJed Brown   PetscScalar       *s;
325c4762a1bSJed Brown   const PetscScalar *xg;
326c4762a1bSJed Brown   PetscInt          i,j,lenglob;
327c4762a1bSJed Brown   PetscReal         sum,val;
328c4762a1bSJed Brown   PetscRandom       rand;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   PetscFunctionBegin;
3319566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
3329566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand,.9,1.0));
3339566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx->da,u,&s));
3349566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
335c4762a1bSJed Brown   lenglob  = appctx->param.E*(appctx->param.N-1);
336c4762a1bSJed Brown   for (i=0; i<lenglob; i++) {
337c4762a1bSJed Brown     s[i]= 0;
338c4762a1bSJed Brown     for (j=0; j<appctx->ncoeff; j++) {
3399566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand,&val));
340c4762a1bSJed Brown       s[i] += val*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]);
341c4762a1bSJed Brown     }
342c4762a1bSJed Brown   }
3439566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
3449566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx->da,u,&s));
3459566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
346c4762a1bSJed Brown   /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */
3479566063dSJacob Faibussowitsch   PetscCall(VecSum(u,&sum));
3489566063dSJacob Faibussowitsch   PetscCall(VecShift(u,-sum/lenglob));
349c4762a1bSJed Brown   PetscFunctionReturn(0);
350c4762a1bSJed Brown }
351c4762a1bSJed Brown 
352c4762a1bSJed Brown /*
353c4762a1bSJed Brown    TrueSolution() computes the true solution for the Tao optimization solve which means they are the initial conditions for the objective function.
354c4762a1bSJed Brown 
355a5b23f4aSJose E. Roman              InitialConditions() computes the initial conditions for the beginning of the Tao iterations
356c4762a1bSJed Brown 
357c4762a1bSJed Brown    Input Parameter:
358c4762a1bSJed Brown    u - uninitialized solution vector (global)
359c4762a1bSJed Brown    appctx - user-defined application context
360c4762a1bSJed Brown 
361c4762a1bSJed Brown    Output Parameter:
362c4762a1bSJed Brown    u - vector with solution at initial time (global)
363c4762a1bSJed Brown */
364c4762a1bSJed Brown PetscErrorCode TrueSolution(Vec u,AppCtx *appctx)
365c4762a1bSJed Brown {
366c4762a1bSJed Brown   PetscScalar       *s;
367c4762a1bSJed Brown   const PetscScalar *xg;
368c4762a1bSJed Brown   PetscInt          i,j,lenglob;
369c4762a1bSJed Brown   PetscReal         sum;
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   PetscFunctionBegin;
3729566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx->da,u,&s));
3739566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
374c4762a1bSJed Brown   lenglob  = appctx->param.E*(appctx->param.N-1);
375c4762a1bSJed Brown   for (i=0; i<lenglob; i++) {
376c4762a1bSJed Brown     s[i]= 0;
377c4762a1bSJed Brown     for (j=0; j<appctx->ncoeff; j++) {
378c4762a1bSJed Brown       s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*xg[i]);
379c4762a1bSJed Brown     }
380c4762a1bSJed Brown   }
3819566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx->da,u,&s));
3829566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
383c4762a1bSJed Brown   /* make sure initial conditions do not contain the constant functions, since with periodic boundary conditions the constant functions introduce a null space */
3849566063dSJacob Faibussowitsch   PetscCall(VecSum(u,&sum));
3859566063dSJacob Faibussowitsch   PetscCall(VecShift(u,-sum/lenglob));
386c4762a1bSJed Brown   PetscFunctionReturn(0);
387c4762a1bSJed Brown }
388c4762a1bSJed Brown /* --------------------------------------------------------------------- */
389c4762a1bSJed Brown /*
390c4762a1bSJed Brown    Sets the desired profile for the final end time
391c4762a1bSJed Brown 
392c4762a1bSJed Brown    Input Parameters:
393c4762a1bSJed Brown    t - final time
394c4762a1bSJed Brown    obj - vector storing the desired profile
395c4762a1bSJed Brown    appctx - user-defined application context
396c4762a1bSJed Brown 
397c4762a1bSJed Brown */
398c4762a1bSJed Brown PetscErrorCode ComputeReference(TS ts,PetscReal t,Vec obj,AppCtx *appctx)
399c4762a1bSJed Brown {
400c4762a1bSJed Brown   PetscScalar       *s,tc;
401c4762a1bSJed Brown   const PetscScalar *xg;
402c4762a1bSJed Brown   PetscInt          i, j,lenglob;
403c4762a1bSJed Brown 
404c4762a1bSJed Brown   PetscFunctionBegin;
4059566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(appctx->da,obj,&s));
4069566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
407c4762a1bSJed Brown   lenglob  = appctx->param.E*(appctx->param.N-1);
408c4762a1bSJed Brown   for (i=0; i<lenglob; i++) {
409c4762a1bSJed Brown     s[i]= 0;
410c4762a1bSJed Brown     for (j=0; j<appctx->ncoeff; j++) {
411c4762a1bSJed Brown       tc    = -appctx->param.mu*(j+1)*(j+1)*4.0*PETSC_PI*PETSC_PI*t;
412c4762a1bSJed Brown       s[i] += appctx->solutioncoefficients[j]*PetscSinScalar(2*(j+1)*PETSC_PI*(xg[i] + appctx->param.a*t))*PetscExpReal(tc);
413c4762a1bSJed Brown     }
414c4762a1bSJed Brown   }
4159566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(appctx->da,obj,&s));
4169566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg));
417c4762a1bSJed Brown   PetscFunctionReturn(0);
418c4762a1bSJed Brown }
419c4762a1bSJed Brown 
420c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
421c4762a1bSJed Brown {
422c4762a1bSJed Brown   AppCtx          *appctx = (AppCtx*)ctx;
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   PetscFunctionBegin;
4259566063dSJacob Faibussowitsch   PetscCall(MatMult(appctx->SEMop.keptstiff,globalin,globalout));
426c4762a1bSJed Brown   PetscFunctionReturn(0);
427c4762a1bSJed Brown }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A, Mat B,void *ctx)
430c4762a1bSJed Brown {
431c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   PetscFunctionBegin;
4349566063dSJacob Faibussowitsch   PetscCall(MatCopy(appctx->SEMop.keptstiff,A,DIFFERENT_NONZERO_PATTERN));
435c4762a1bSJed Brown   PetscFunctionReturn(0);
436c4762a1bSJed Brown }
437c4762a1bSJed Brown 
438c4762a1bSJed Brown /* --------------------------------------------------------------------- */
439c4762a1bSJed Brown 
440c4762a1bSJed Brown /*
441c4762a1bSJed Brown    RHSLaplacian -   matrix for diffusion
442c4762a1bSJed Brown 
443c4762a1bSJed Brown    Input Parameters:
444c4762a1bSJed Brown    ts - the TS context
445c4762a1bSJed Brown    t - current time  (ignored)
446c4762a1bSJed Brown    X - current solution (ignored)
447c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
448c4762a1bSJed Brown 
449c4762a1bSJed Brown    Output Parameters:
450c4762a1bSJed Brown    AA - Jacobian matrix
451c4762a1bSJed Brown    BB - optionally different matrix from which the preconditioner is built
452c4762a1bSJed Brown    str - flag indicating matrix structure
453c4762a1bSJed Brown 
454c4762a1bSJed Brown    Scales by the inverse of the mass matrix (perhaps that should be pulled out)
455c4762a1bSJed Brown 
456c4762a1bSJed Brown */
457c4762a1bSJed Brown PetscErrorCode RHSLaplacian(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx)
458c4762a1bSJed Brown {
459c4762a1bSJed Brown   PetscReal      **temp;
460c4762a1bSJed Brown   PetscReal      vv;
461c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
462c4762a1bSJed Brown   PetscInt       i,xs,xn,l,j;
463c4762a1bSJed Brown   PetscInt       *rowsDM;
464c4762a1bSJed Brown 
465c4762a1bSJed Brown   PetscFunctionBegin;
466c4762a1bSJed Brown   /*
467c4762a1bSJed Brown    Creates the element stiffness matrix for the given gll
468c4762a1bSJed Brown    */
4699566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
470c4762a1bSJed Brown 
471c4762a1bSJed Brown   /* scale by the size of the element */
472c4762a1bSJed Brown   for (i=0; i<appctx->param.N; i++) {
473c4762a1bSJed Brown     vv=-appctx->param.mu*2.0/appctx->param.Le;
474c4762a1bSJed Brown     for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv;
475c4762a1bSJed Brown   }
476c4762a1bSJed Brown 
4779566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
4789566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL));
479c4762a1bSJed Brown 
4803c859ba3SBarry Smith   PetscCheck(appctx->param.N-1 >= 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2");
481c4762a1bSJed Brown   xs   = xs/(appctx->param.N-1);
482c4762a1bSJed Brown   xn   = xn/(appctx->param.N-1);
483c4762a1bSJed Brown 
4849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(appctx->param.N,&rowsDM));
485c4762a1bSJed Brown   /*
486c4762a1bSJed Brown    loop over local elements
487c4762a1bSJed Brown    */
488c4762a1bSJed Brown   for (j=xs; j<xs+xn; j++) {
489c4762a1bSJed Brown     for (l=0; l<appctx->param.N; l++) {
490c4762a1bSJed Brown       rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l;
491c4762a1bSJed Brown     }
4929566063dSJacob Faibussowitsch     PetscCall(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES));
493c4762a1bSJed Brown   }
4949566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowsDM));
4959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
4969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
4979566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
4989566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,appctx->SEMop.mass,0));
4999566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
500c4762a1bSJed Brown 
5019566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
502c4762a1bSJed Brown   PetscFunctionReturn(0);
503c4762a1bSJed Brown }
504c4762a1bSJed Brown 
505c4762a1bSJed Brown /*
506c4762a1bSJed Brown     Almost identical to Laplacian
507c4762a1bSJed Brown 
508c4762a1bSJed Brown     Note that the element matrix is NOT scaled by the size of element like the Laplacian term.
509c4762a1bSJed Brown  */
510c4762a1bSJed Brown PetscErrorCode RHSAdvection(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx)
511c4762a1bSJed Brown {
512c4762a1bSJed Brown   PetscReal      **temp;
513c4762a1bSJed Brown   PetscReal      vv;
514c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
515c4762a1bSJed Brown   PetscInt       i,xs,xn,l,j;
516c4762a1bSJed Brown   PetscInt       *rowsDM;
517c4762a1bSJed Brown 
518c4762a1bSJed Brown   PetscFunctionBegin;
519c4762a1bSJed Brown   /*
520c4762a1bSJed Brown    Creates the element stiffness matrix for the given gll
521c4762a1bSJed Brown    */
5229566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
523c4762a1bSJed Brown 
524c4762a1bSJed Brown   /* scale by the size of the element */
525c4762a1bSJed Brown   for (i=0; i<appctx->param.N; i++) {
526c4762a1bSJed Brown     vv = -appctx->param.a;
527c4762a1bSJed Brown     for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv;
528c4762a1bSJed Brown   }
529c4762a1bSJed Brown 
5309566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
5319566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL));
532c4762a1bSJed Brown 
5333c859ba3SBarry Smith   PetscCheck(appctx->param.N-1 >= 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Polynomial order must be at least 2");
534c4762a1bSJed Brown   xs   = xs/(appctx->param.N-1);
535c4762a1bSJed Brown   xn   = xn/(appctx->param.N-1);
536c4762a1bSJed Brown 
5379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(appctx->param.N,&rowsDM));
538c4762a1bSJed Brown   /*
539c4762a1bSJed Brown    loop over local elements
540c4762a1bSJed Brown    */
541c4762a1bSJed Brown   for (j=xs; j<xs+xn; j++) {
542c4762a1bSJed Brown     for (l=0; l<appctx->param.N; l++) {
543c4762a1bSJed Brown       rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l;
544c4762a1bSJed Brown     }
5459566063dSJacob Faibussowitsch     PetscCall(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES));
546c4762a1bSJed Brown   }
5479566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowsDM));
5489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
5499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
5509566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
5519566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,appctx->SEMop.mass,0));
5529566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(appctx->SEMop.mass));
553c4762a1bSJed Brown 
5549566063dSJacob Faibussowitsch   PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp));
555c4762a1bSJed Brown   PetscFunctionReturn(0);
556c4762a1bSJed Brown }
557c4762a1bSJed Brown 
558c4762a1bSJed Brown /* ------------------------------------------------------------------ */
559c4762a1bSJed Brown /*
560c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
561c4762a1bSJed Brown 
562c4762a1bSJed Brown    Input Parameters:
563c4762a1bSJed Brown    tao - the Tao context
564c4762a1bSJed Brown    ic   - the input vector
565a82e8c82SStefano Zampini    ctx - optional user-defined context, as set when calling TaoSetObjectiveAndGradient()
566c4762a1bSJed Brown 
567c4762a1bSJed Brown    Output Parameters:
568c4762a1bSJed Brown    f   - the newly evaluated function
569c4762a1bSJed Brown    G   - the newly evaluated gradient
570c4762a1bSJed Brown 
571c4762a1bSJed Brown    Notes:
572c4762a1bSJed Brown 
573c4762a1bSJed Brown           The forward equation is
574c4762a1bSJed Brown               M u_t = F(U)
575c4762a1bSJed Brown           which is converted to
576c4762a1bSJed Brown                 u_t = M^{-1} F(u)
577c4762a1bSJed Brown           in the user code since TS has no direct way of providing a mass matrix. The Jacobian of this is
578c4762a1bSJed Brown                  M^{-1} J
579c4762a1bSJed Brown           where J is the Jacobian of F. Now the adjoint equation is
580c4762a1bSJed Brown                 M v_t = J^T v
581c4762a1bSJed Brown           but TSAdjoint does not solve this since it can only solve the transposed system for the
582c4762a1bSJed Brown           Jacobian the user provided. Hence TSAdjoint solves
583c4762a1bSJed Brown                  w_t = J^T M^{-1} w  (where w = M v)
584a5b23f4aSJose E. Roman           since there is no way to indicate the mass matrix as a separate entity to TS. Thus one
585c4762a1bSJed Brown           must be careful in initializing the "adjoint equation" and using the result. This is
586c4762a1bSJed Brown           why
587c4762a1bSJed Brown               G = -2 M(u(T) - u_d)
588c4762a1bSJed Brown           below (instead of -2(u(T) - u_d)
589c4762a1bSJed Brown 
590c4762a1bSJed Brown */
591c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec ic,PetscReal *f,Vec G,void *ctx)
592c4762a1bSJed Brown {
593c4762a1bSJed Brown   AppCtx           *appctx = (AppCtx*)ctx;     /* user-defined application context */
594c4762a1bSJed Brown   Vec               temp;
595c4762a1bSJed Brown 
596c4762a1bSJed Brown   PetscFunctionBegin;
5979566063dSJacob Faibussowitsch   PetscCall(TSSetTime(appctx->ts,0.0));
5989566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(appctx->ts,0));
5999566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(appctx->ts,appctx->initial_dt));
6009566063dSJacob Faibussowitsch   PetscCall(VecCopy(ic,appctx->dat.curr_sol));
601c4762a1bSJed Brown 
6029566063dSJacob Faibussowitsch   PetscCall(TSSolve(appctx->ts,appctx->dat.curr_sol));
6039566063dSJacob Faibussowitsch   PetscCall(VecCopy(appctx->dat.curr_sol,appctx->dat.joe));
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   /*     Compute the difference between the current ODE solution and target ODE solution */
6069566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(G,-1.0,appctx->dat.curr_sol,appctx->dat.reference));
607c4762a1bSJed Brown 
608c4762a1bSJed Brown   /*     Compute the objective/cost function   */
6099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(G,&temp));
6109566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(temp,G,G));
6119566063dSJacob Faibussowitsch   PetscCall(VecDot(temp,appctx->SEMop.mass,f));
6129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&temp));
613c4762a1bSJed Brown 
614c4762a1bSJed Brown   /*     Compute initial conditions for the adjoint integration. See Notes above  */
6159566063dSJacob Faibussowitsch   PetscCall(VecScale(G, -2.0));
6169566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(G,G,appctx->SEMop.mass));
6179566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(appctx->ts,1,&G,NULL));
618c4762a1bSJed Brown 
6199566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(appctx->ts));
6209566063dSJacob Faibussowitsch   /* PetscCall(VecPointwiseDivide(G,G,appctx->SEMop.mass));*/
621c4762a1bSJed Brown   PetscFunctionReturn(0);
622c4762a1bSJed Brown }
623c4762a1bSJed Brown 
624c4762a1bSJed Brown PetscErrorCode MonitorError(Tao tao,void *ctx)
625c4762a1bSJed Brown {
626c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
627c4762a1bSJed Brown   Vec            temp,grad;
628c4762a1bSJed Brown   PetscReal      nrm;
629c4762a1bSJed Brown   PetscInt       its;
630c4762a1bSJed Brown   PetscReal      fct,gnorm;
631c4762a1bSJed Brown 
632c4762a1bSJed Brown   PetscFunctionBegin;
6339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(appctx->dat.ic,&temp));
6349566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(temp,-1.0,appctx->dat.ic,appctx->dat.true_solution));
6359566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(temp,temp,temp));
6369566063dSJacob Faibussowitsch   PetscCall(VecDot(temp,appctx->SEMop.mass,&nrm));
637c4762a1bSJed Brown   nrm   = PetscSqrtReal(nrm);
6389566063dSJacob Faibussowitsch   PetscCall(TaoGetGradient(tao,&grad,NULL,NULL));
6399566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(temp,temp,temp));
6409566063dSJacob Faibussowitsch   PetscCall(VecDot(temp,appctx->SEMop.mass,&gnorm));
641c4762a1bSJed Brown   gnorm = PetscSqrtReal(gnorm);
6429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&temp));
6439566063dSJacob Faibussowitsch   PetscCall(TaoGetIterationNumber(tao,&its));
6449566063dSJacob Faibussowitsch   PetscCall(TaoGetSolutionStatus(tao,NULL,&fct,NULL,NULL,NULL,NULL));
645c4762a1bSJed Brown   if (!its) {
6469566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%% Iteration Error Objective Gradient-norm\n"));
6479566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"history = [\n"));
648c4762a1bSJed Brown   }
649*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%3" PetscInt_FMT " %g %g %g\n",its,(double)nrm,(double)fct,(double)gnorm));
650c4762a1bSJed Brown   PetscFunctionReturn(0);
651c4762a1bSJed Brown }
652c4762a1bSJed Brown 
653c4762a1bSJed Brown PetscErrorCode MonitorDestroy(void **ctx)
654c4762a1bSJed Brown {
655c4762a1bSJed Brown   PetscFunctionBegin;
6569566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"];\n"));
657c4762a1bSJed Brown   PetscFunctionReturn(0);
658c4762a1bSJed Brown }
659c4762a1bSJed Brown 
660c4762a1bSJed Brown /*TEST
661c4762a1bSJed Brown 
662c4762a1bSJed Brown    build:
663c4762a1bSJed Brown      requires: !complex
664c4762a1bSJed Brown 
665c4762a1bSJed Brown    test:
666c4762a1bSJed Brown      requires: !single
667c4762a1bSJed Brown      args:  -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none
668c4762a1bSJed Brown 
669c4762a1bSJed Brown    test:
670c4762a1bSJed Brown      suffix: cn
671c4762a1bSJed Brown      requires: !single
672c4762a1bSJed Brown      args:  -ts_type cn -ts_dt .003 -pc_type lu -E 10 -N 8 -ncoeff 5 -tao_bqnls_mat_lmvm_scale_type none
673c4762a1bSJed Brown 
674c4762a1bSJed Brown    test:
675c4762a1bSJed Brown      suffix: 2
676c4762a1bSJed Brown      requires: !single
677c4762a1bSJed Brown      args:  -ts_adapt_dt_max 3.e-3 -E 10 -N 8 -ncoeff 5  -a .1 -tao_bqnls_mat_lmvm_scale_type none
678c4762a1bSJed Brown 
679c4762a1bSJed Brown TEST*/
680