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