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