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