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