xref: /petsc/src/ts/tutorials/ex50.c (revision a5b23f4acc7afc99d3844ebd5fb65a81c16e8b8c)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Solves one dimensional Burger's equation compares with exact solution\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 nonlinear problems
11c4762a1bSJed Brown    Concepts: TS^Burger's equation
12c4762a1bSJed Brown    Processors: n
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown /* ------------------------------------------------------------------------
16c4762a1bSJed Brown 
17c4762a1bSJed Brown    This program uses the one-dimensional Burger's equation
18c4762a1bSJed Brown        u_t = mu*u_xx - u u_x,
19c4762a1bSJed Brown    on the domain 0 <= x <= 1, with periodic boundary conditions
20c4762a1bSJed Brown 
21c4762a1bSJed Brown    The operators are discretized with the spectral element method
22c4762a1bSJed Brown 
23c4762a1bSJed Brown    See the paper PDE-CONSTRAINED OPTIMIZATION WITH SPECTRAL ELEMENTS USING PETSC AND TAO
24c4762a1bSJed Brown    by OANA MARIN, EMIL CONSTANTINESCU, AND BARRY SMITH for details on the exact solution
25c4762a1bSJed Brown    used
26c4762a1bSJed Brown 
27c4762a1bSJed Brown    See src/tao/unconstrained/tutorials/burgers_spectral.c
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   ------------------------------------------------------------------------- */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown #include <petscts.h>
32c4762a1bSJed Brown #include <petscdt.h>
33c4762a1bSJed Brown #include <petscdraw.h>
34c4762a1bSJed Brown #include <petscdmda.h>
35c4762a1bSJed Brown 
36c4762a1bSJed Brown /*
37c4762a1bSJed Brown    User-defined application context - contains data needed by the
38c4762a1bSJed Brown    application-provided call-back routines.
39c4762a1bSJed Brown */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown typedef struct {
42c4762a1bSJed Brown   PetscInt  n;                /* number of nodes */
43c4762a1bSJed Brown   PetscReal *nodes;           /* GLL nodes */
44c4762a1bSJed Brown   PetscReal *weights;         /* GLL weights */
45c4762a1bSJed Brown } PetscGLL;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown typedef struct {
48c4762a1bSJed Brown   PetscInt    N;             /* grid points per elements*/
49c4762a1bSJed Brown   PetscInt    E;              /* number of elements */
50c4762a1bSJed Brown   PetscReal   tol_L2,tol_max; /* error norms */
51c4762a1bSJed Brown   PetscInt    steps;          /* number of timesteps */
52c4762a1bSJed Brown   PetscReal   Tend;           /* endtime */
53c4762a1bSJed Brown   PetscReal   mu;             /* viscosity */
54c4762a1bSJed Brown   PetscReal   L;              /* total length of domain */
55c4762a1bSJed Brown   PetscReal   Le;
56c4762a1bSJed Brown   PetscReal   Tadj;
57c4762a1bSJed Brown } PetscParam;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown typedef struct {
60c4762a1bSJed Brown   Vec         grid;              /* total grid */
61c4762a1bSJed Brown   Vec         curr_sol;
62c4762a1bSJed Brown } PetscData;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown typedef struct {
65c4762a1bSJed Brown   Vec         grid;              /* total grid */
66c4762a1bSJed Brown   Vec         mass;              /* mass matrix for total integration */
67c4762a1bSJed Brown   Mat         stiff;             /* stifness matrix */
68c4762a1bSJed Brown   Mat         keptstiff;
69c4762a1bSJed Brown   Mat         grad;
70c4762a1bSJed Brown   PetscGLL    gll;
71c4762a1bSJed Brown } PetscSEMOperators;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown typedef struct {
74c4762a1bSJed Brown   DM                da;                /* distributed array data structure */
75c4762a1bSJed Brown   PetscSEMOperators SEMop;
76c4762a1bSJed Brown   PetscParam        param;
77c4762a1bSJed Brown   PetscData         dat;
78c4762a1bSJed Brown   TS                ts;
79c4762a1bSJed Brown   PetscReal         initial_dt;
80c4762a1bSJed Brown } AppCtx;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown /*
83c4762a1bSJed Brown    User-defined routines
84c4762a1bSJed Brown */
85c4762a1bSJed Brown extern PetscErrorCode RHSMatrixLaplaciangllDM(TS,PetscReal,Vec,Mat,Mat,void*);
86c4762a1bSJed Brown extern PetscErrorCode RHSMatrixAdvectiongllDM(TS,PetscReal,Vec,Mat,Mat,void*);
87c4762a1bSJed Brown extern PetscErrorCode TrueSolution(TS,PetscReal,Vec,AppCtx*);
88c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
89c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown int main(int argc,char **argv)
92c4762a1bSJed Brown {
93c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
94c4762a1bSJed Brown   PetscErrorCode ierr;
95c4762a1bSJed Brown   PetscInt       i, xs, xm, ind, j, lenglob;
96c4762a1bSJed Brown   PetscReal      x, *wrk_ptr1, *wrk_ptr2;
97c4762a1bSJed Brown   MatNullSpace   nsp;
98c4762a1bSJed Brown   PetscMPIInt    size;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown      Initialize program and set problem parameters
102c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103c4762a1bSJed Brown   PetscFunctionBegin;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /*initialize parameters */
108c4762a1bSJed Brown   appctx.param.N    = 10;  /* order of the spectral element */
109c4762a1bSJed Brown   appctx.param.E    = 10;  /* number of elements */
110c4762a1bSJed Brown   appctx.param.L    = 4.0;  /* length of the domain */
111c4762a1bSJed Brown   appctx.param.mu   = 0.01; /* diffusion coefficient */
112c4762a1bSJed Brown   appctx.initial_dt = 5e-3;
113c4762a1bSJed Brown   appctx.param.steps = PETSC_MAX_INT;
114c4762a1bSJed Brown   appctx.param.Tend  = 4;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&appctx.param.N,NULL);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-E",&appctx.param.E,NULL);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-Tend",&appctx.param.Tend,NULL);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&appctx.param.mu,NULL);CHKERRQ(ierr);
120c4762a1bSJed Brown   appctx.param.Le = appctx.param.L/appctx.param.E;
121c4762a1bSJed Brown 
122ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
123c4762a1bSJed Brown   if (appctx.param.E % size) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of elements must be divisible by number of processes");
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Create GLL data structures
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128c4762a1bSJed Brown   ierr = PetscMalloc2(appctx.param.N,&appctx.SEMop.gll.nodes,appctx.param.N,&appctx.SEMop.gll.weights);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = PetscDTGaussLobattoLegendreQuadrature(appctx.param.N,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights);CHKERRQ(ierr);
130c4762a1bSJed Brown   appctx.SEMop.gll.n = appctx.param.N;
131c4762a1bSJed Brown   lenglob  = appctx.param.E*(appctx.param.N-1);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /*
134c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
135c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are E*(Nl-1)+1
136c4762a1bSJed Brown      total grid values spread equally among all the processors, except first and last
137c4762a1bSJed Brown   */
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,lenglob,1,1,NULL,&appctx.da);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = DMSetFromOptions(appctx.da);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = DMSetUp(appctx.da);CHKERRQ(ierr);
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /*
144c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
145c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
146c4762a1bSJed Brown      have the same types.
147c4762a1bSJed Brown   */
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   ierr = DMCreateGlobalVector(appctx.da,&appctx.dat.curr_sol);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = VecDuplicate(appctx.dat.curr_sol,&appctx.SEMop.grid);CHKERRQ(ierr);
151c4762a1bSJed Brown   ierr = VecDuplicate(appctx.dat.curr_sol,&appctx.SEMop.mass);CHKERRQ(ierr);
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   ierr = DMDAGetCorners(appctx.da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
154c4762a1bSJed Brown   ierr = DMDAVecGetArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1);CHKERRQ(ierr);
155c4762a1bSJed Brown   ierr = DMDAVecGetArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2);CHKERRQ(ierr);
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
158c4762a1bSJed Brown 
159c4762a1bSJed Brown     xs=xs/(appctx.param.N-1);
160c4762a1bSJed Brown     xm=xm/(appctx.param.N-1);
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /*
163c4762a1bSJed Brown      Build total grid and mass over entire mesh (multi-elemental)
164c4762a1bSJed Brown   */
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
167c4762a1bSJed Brown     for (j=0; j<appctx.param.N-1; j++) {
168c4762a1bSJed Brown       x = (appctx.param.Le/2.0)*(appctx.SEMop.gll.nodes[j]+1.0)+appctx.param.Le*i;
169c4762a1bSJed Brown       ind=i*(appctx.param.N-1)+j;
170c4762a1bSJed Brown       wrk_ptr1[ind]=x;
171c4762a1bSJed Brown       wrk_ptr2[ind]=.5*appctx.param.Le*appctx.SEMop.gll.weights[j];
172c4762a1bSJed Brown       if (j==0) wrk_ptr2[ind]+=.5*appctx.param.Le*appctx.SEMop.gll.weights[j];
173c4762a1bSJed Brown     }
174c4762a1bSJed Brown   }
175c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1);CHKERRQ(ierr);
176c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2);CHKERRQ(ierr);
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179c4762a1bSJed Brown    Create matrix data structure; set matrix evaluation routine.
180c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181c4762a1bSJed Brown   ierr = DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE);CHKERRQ(ierr);
182c4762a1bSJed Brown   ierr = DMCreateMatrix(appctx.da,&appctx.SEMop.stiff);CHKERRQ(ierr);
183c4762a1bSJed Brown   ierr = DMCreateMatrix(appctx.da,&appctx.SEMop.grad);CHKERRQ(ierr);
184c4762a1bSJed Brown   /*
185c4762a1bSJed Brown    For linear problems with a time-dependent f(u,t) in the equation
186c4762a1bSJed Brown    u_t = f(u,t), the user provides the discretized right-hand-side
187c4762a1bSJed Brown    as a time-dependent matrix.
188c4762a1bSJed Brown    */
189c4762a1bSJed Brown   ierr = RHSMatrixLaplaciangllDM(appctx.ts,0.0,appctx.dat.curr_sol,appctx.SEMop.stiff,appctx.SEMop.stiff,&appctx);CHKERRQ(ierr);
190c4762a1bSJed Brown   ierr = RHSMatrixAdvectiongllDM(appctx.ts,0.0,appctx.dat.curr_sol,appctx.SEMop.grad,appctx.SEMop.grad,&appctx);CHKERRQ(ierr);
191c4762a1bSJed Brown    /*
192c4762a1bSJed Brown        For linear problems with a time-dependent f(u,t) in the equation
193c4762a1bSJed Brown        u_t = f(u,t), the user provides the discretized right-hand-side
194c4762a1bSJed Brown        as a time-dependent matrix.
195c4762a1bSJed Brown     */
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   ierr = MatDuplicate(appctx.SEMop.stiff,MAT_COPY_VALUES,&appctx.SEMop.keptstiff);CHKERRQ(ierr);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* attach the null space to the matrix, this probably is not needed but does no harm */
200c4762a1bSJed Brown   ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = MatSetNullSpace(appctx.SEMop.stiff,nsp);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = MatSetNullSpace(appctx.SEMop.keptstiff,nsp);CHKERRQ(ierr);
203c4762a1bSJed Brown   ierr = MatNullSpaceTest(nsp,appctx.SEMop.stiff,NULL);CHKERRQ(ierr);
204c4762a1bSJed Brown   ierr = MatNullSpaceDestroy(&nsp);CHKERRQ(ierr);
205c4762a1bSJed Brown   /* attach the null space to the matrix, this probably is not needed but does no harm */
206c4762a1bSJed Brown   ierr = MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp);CHKERRQ(ierr);
207c4762a1bSJed Brown   ierr = MatSetNullSpace(appctx.SEMop.grad,nsp);CHKERRQ(ierr);
208c4762a1bSJed Brown   ierr = MatNullSpaceTest(nsp,appctx.SEMop.grad,NULL);CHKERRQ(ierr);
209c4762a1bSJed Brown   ierr = MatNullSpaceDestroy(&nsp);CHKERRQ(ierr);
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* Create the TS solver that solves the ODE and its adjoint; set its options */
212c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&appctx.ts);CHKERRQ(ierr);
213c4762a1bSJed Brown   ierr = TSSetProblemType(appctx.ts,TS_NONLINEAR);CHKERRQ(ierr);
214c4762a1bSJed Brown   ierr = TSSetType(appctx.ts,TSRK);CHKERRQ(ierr);
215c4762a1bSJed Brown   ierr = TSSetDM(appctx.ts,appctx.da);CHKERRQ(ierr);
216c4762a1bSJed Brown   ierr = TSSetTime(appctx.ts,0.0);CHKERRQ(ierr);
217c4762a1bSJed Brown   ierr = TSSetTimeStep(appctx.ts,appctx.initial_dt);CHKERRQ(ierr);
218c4762a1bSJed Brown   ierr = TSSetMaxSteps(appctx.ts,appctx.param.steps);CHKERRQ(ierr);
219c4762a1bSJed Brown   ierr = TSSetMaxTime(appctx.ts,appctx.param.Tend);CHKERRQ(ierr);
220c4762a1bSJed Brown   ierr = TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
221c4762a1bSJed Brown   ierr = TSSetTolerances(appctx.ts,1e-7,NULL,1e-7,NULL);CHKERRQ(ierr);
222c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(appctx.ts);CHKERRQ(ierr);
223c4762a1bSJed Brown   ierr = TSSetFromOptions(appctx.ts);CHKERRQ(ierr);
224c4762a1bSJed Brown   ierr = TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr);
225c4762a1bSJed Brown   ierr = TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx);CHKERRQ(ierr);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* Set Initial conditions for the problem  */
228c4762a1bSJed Brown   ierr = TrueSolution(appctx.ts,0,appctx.dat.curr_sol,&appctx);CHKERRQ(ierr);
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   ierr = TSSetSolutionFunction(appctx.ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void *))TrueSolution,&appctx);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = TSSetTime(appctx.ts,0.0);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = TSSetStepNumber(appctx.ts,0);CHKERRQ(ierr);
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   ierr = TSSolve(appctx.ts,appctx.dat.curr_sol);CHKERRQ(ierr);
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   ierr = MatDestroy(&appctx.SEMop.stiff);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = MatDestroy(&appctx.SEMop.keptstiff);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = MatDestroy(&appctx.SEMop.grad);CHKERRQ(ierr);
239c4762a1bSJed Brown   ierr = VecDestroy(&appctx.SEMop.grid);CHKERRQ(ierr);
240c4762a1bSJed Brown   ierr = VecDestroy(&appctx.SEMop.mass);CHKERRQ(ierr);
241c4762a1bSJed Brown   ierr = VecDestroy(&appctx.dat.curr_sol);CHKERRQ(ierr);
242c4762a1bSJed Brown   ierr = PetscFree2(appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights);CHKERRQ(ierr);
243c4762a1bSJed Brown   ierr = DMDestroy(&appctx.da);CHKERRQ(ierr);
244c4762a1bSJed Brown   ierr = TSDestroy(&appctx.ts);CHKERRQ(ierr);
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /*
247c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
248c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
249c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
250c4762a1bSJed Brown          options are chosen (e.g., -log_summary).
251c4762a1bSJed Brown   */
252c4762a1bSJed Brown     ierr = PetscFinalize();
253c4762a1bSJed Brown     return ierr;
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown /*
257c4762a1bSJed Brown    TrueSolution() computes the true solution for the PDE
258c4762a1bSJed Brown 
259c4762a1bSJed Brown    Input Parameter:
260c4762a1bSJed Brown    u - uninitialized solution vector (global)
261c4762a1bSJed Brown    appctx - user-defined application context
262c4762a1bSJed Brown 
263c4762a1bSJed Brown    Output Parameter:
264c4762a1bSJed Brown    u - vector with solution at initial time (global)
265c4762a1bSJed Brown */
266c4762a1bSJed Brown PetscErrorCode TrueSolution(TS ts, PetscReal t, Vec u,AppCtx *appctx)
267c4762a1bSJed Brown {
268c4762a1bSJed Brown   PetscScalar       *s;
269c4762a1bSJed Brown   const PetscScalar *xg;
270c4762a1bSJed Brown   PetscErrorCode    ierr;
271c4762a1bSJed Brown   PetscInt          i,xs,xn;
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   ierr = DMDAVecGetArray(appctx->da,u,&s);CHKERRQ(ierr);
274c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr);
275c4762a1bSJed Brown   ierr = DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr);
276c4762a1bSJed Brown   for (i=xs; i<xs+xn; i++) {
277c4762a1bSJed Brown     s[i]=2.0*appctx->param.mu*PETSC_PI*PetscSinScalar(PETSC_PI*xg[i])*PetscExpReal(-appctx->param.mu*PETSC_PI*PETSC_PI*t)/(2.0+PetscCosScalar(PETSC_PI*xg[i])*PetscExpReal(-appctx->param.mu*PETSC_PI*PETSC_PI*t));
278c4762a1bSJed Brown   }
279c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(appctx->da,u,&s);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg);CHKERRQ(ierr);
281c4762a1bSJed Brown   return 0;
282c4762a1bSJed Brown }
283c4762a1bSJed Brown 
284c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
285c4762a1bSJed Brown {
286c4762a1bSJed Brown   PetscErrorCode ierr;
287c4762a1bSJed Brown   AppCtx          *appctx = (AppCtx*)ctx;
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   PetscFunctionBegin;
290c4762a1bSJed Brown   ierr = MatMult(appctx->SEMop.grad,globalin,globalout);CHKERRQ(ierr); /* grad u */
291c4762a1bSJed Brown   ierr = VecPointwiseMult(globalout,globalin,globalout);CHKERRQ(ierr); /* u grad u */
292c4762a1bSJed Brown   ierr = VecScale(globalout, -1.0);CHKERRQ(ierr);
293c4762a1bSJed Brown   ierr = MatMultAdd(appctx->SEMop.keptstiff,globalin,globalout,globalout);CHKERRQ(ierr);
294c4762a1bSJed Brown   PetscFunctionReturn(0);
295c4762a1bSJed Brown }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown /*
298c4762a1bSJed Brown 
299c4762a1bSJed Brown       K is the discretiziation of the Laplacian
300c4762a1bSJed Brown       G is the discretization of the gradient
301c4762a1bSJed Brown 
302c4762a1bSJed Brown       Computes Jacobian of      K u + diag(u) G u   which is given by
303c4762a1bSJed Brown               K   + diag(u)G + diag(Gu)
304c4762a1bSJed Brown */
305c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A, Mat B,void *ctx)
306c4762a1bSJed Brown {
307c4762a1bSJed Brown   PetscErrorCode ierr;
308c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
309c4762a1bSJed Brown   Vec            Gglobalin;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   PetscFunctionBegin;
312c4762a1bSJed Brown   /*    A = diag(u) G */
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   ierr = MatCopy(appctx->SEMop.grad,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
315c4762a1bSJed Brown   ierr = MatDiagonalScale(A,globalin,NULL);CHKERRQ(ierr);
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   /*    A  = A + diag(Gu) */
318c4762a1bSJed Brown   ierr = VecDuplicate(globalin,&Gglobalin);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = MatMult(appctx->SEMop.grad,globalin,Gglobalin);CHKERRQ(ierr);
320c4762a1bSJed Brown   ierr = MatDiagonalSet(A,Gglobalin,ADD_VALUES);CHKERRQ(ierr);
321c4762a1bSJed Brown   ierr = VecDestroy(&Gglobalin);CHKERRQ(ierr);
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /*   A  = K - A    */
324c4762a1bSJed Brown   ierr = MatScale(A,-1.0);CHKERRQ(ierr);
325c4762a1bSJed Brown   ierr = MatAXPY(A,0.0,appctx->SEMop.keptstiff,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
326c4762a1bSJed Brown   PetscFunctionReturn(0);
327c4762a1bSJed Brown }
328c4762a1bSJed Brown 
329c4762a1bSJed Brown /* --------------------------------------------------------------------- */
330c4762a1bSJed Brown 
331c4762a1bSJed Brown #include "petscblaslapack.h"
332c4762a1bSJed Brown /*
333c4762a1bSJed Brown      Matrix free operation of 1d Laplacian and Grad for GLL spectral elements
334c4762a1bSJed Brown */
335c4762a1bSJed Brown PetscErrorCode MatMult_Laplacian(Mat A,Vec x,Vec y)
336c4762a1bSJed Brown {
337c4762a1bSJed Brown   AppCtx            *appctx;
338c4762a1bSJed Brown   PetscErrorCode    ierr;
339c4762a1bSJed Brown   PetscReal         **temp,vv;
340c4762a1bSJed Brown   PetscInt          i,j,xs,xn;
341c4762a1bSJed Brown   Vec               xlocal,ylocal;
342c4762a1bSJed Brown   const PetscScalar *xl;
343c4762a1bSJed Brown   PetscScalar       *yl;
344c4762a1bSJed Brown   PetscBLASInt      _One = 1,n;
345c4762a1bSJed Brown   PetscScalar       _DOne = 1;
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   ierr = MatShellGetContext(A,&appctx);CHKERRQ(ierr);
348c4762a1bSJed Brown   ierr = DMGetLocalVector(appctx->da,&xlocal);CHKERRQ(ierr);
349c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(appctx->da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
350c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(appctx->da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
351c4762a1bSJed Brown   ierr = DMGetLocalVector(appctx->da,&ylocal);CHKERRQ(ierr);
352c4762a1bSJed Brown   ierr = VecSet(ylocal,0.0);CHKERRQ(ierr);
353c4762a1bSJed Brown   ierr = PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
354c4762a1bSJed Brown   for (i=0; i<appctx->param.N; i++) {
355c4762a1bSJed Brown     vv =-appctx->param.mu*2.0/appctx->param.Le;
356c4762a1bSJed Brown     for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv;
357c4762a1bSJed Brown   }
358c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(appctx->da,xlocal,(void*)&xl);CHKERRQ(ierr);
359c4762a1bSJed Brown   ierr = DMDAVecGetArray(appctx->da,ylocal,&yl);CHKERRQ(ierr);
360c4762a1bSJed Brown   ierr = DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr);
361c4762a1bSJed Brown   ierr = PetscBLASIntCast(appctx->param.N,&n);CHKERRQ(ierr);
362c4762a1bSJed Brown   for (j=xs; j<xs+xn; j += appctx->param.N-1) {
363c4762a1bSJed Brown     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&_DOne,&temp[0][0],&n,&xl[j],&_One,&_DOne,&yl[j],&_One));
364c4762a1bSJed Brown   }
365c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(appctx->da,xlocal,(void*)&xl);CHKERRQ(ierr);
366c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(appctx->da,ylocal,&yl);CHKERRQ(ierr);
367c4762a1bSJed Brown   ierr = PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
368c4762a1bSJed Brown   ierr = VecSet(y,0.0);CHKERRQ(ierr);
369c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(appctx->da,ylocal,ADD_VALUES,y);CHKERRQ(ierr);
370c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(appctx->da,ylocal,ADD_VALUES,y);CHKERRQ(ierr);
371c4762a1bSJed Brown   ierr = DMRestoreLocalVector(appctx->da,&xlocal);CHKERRQ(ierr);
372c4762a1bSJed Brown   ierr = DMRestoreLocalVector(appctx->da,&ylocal);CHKERRQ(ierr);
373c4762a1bSJed Brown   ierr = VecPointwiseDivide(y,y,appctx->SEMop.mass);CHKERRQ(ierr);
374c4762a1bSJed Brown   return 0;
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown PetscErrorCode MatMult_Advection(Mat A,Vec x,Vec y)
378c4762a1bSJed Brown {
379c4762a1bSJed Brown   AppCtx            *appctx;
380c4762a1bSJed Brown   PetscErrorCode    ierr;
381c4762a1bSJed Brown   PetscReal         **temp;
382c4762a1bSJed Brown   PetscInt          j,xs,xn;
383c4762a1bSJed Brown   Vec               xlocal,ylocal;
384c4762a1bSJed Brown   const PetscScalar *xl;
385c4762a1bSJed Brown   PetscScalar       *yl;
386c4762a1bSJed Brown   PetscBLASInt      _One = 1,n;
387c4762a1bSJed Brown   PetscScalar       _DOne = 1;
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   ierr = MatShellGetContext(A,&appctx);CHKERRQ(ierr);
390c4762a1bSJed Brown   ierr = DMGetLocalVector(appctx->da,&xlocal);CHKERRQ(ierr);
391c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(appctx->da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
392c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(appctx->da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
393c4762a1bSJed Brown   ierr = DMGetLocalVector(appctx->da,&ylocal);CHKERRQ(ierr);
394c4762a1bSJed Brown   ierr = VecSet(ylocal,0.0);CHKERRQ(ierr);
395c4762a1bSJed Brown   ierr = PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
396c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(appctx->da,xlocal,(void*)&xl);CHKERRQ(ierr);
397c4762a1bSJed Brown   ierr = DMDAVecGetArray(appctx->da,ylocal,&yl);CHKERRQ(ierr);
398c4762a1bSJed Brown   ierr = DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr);
399c4762a1bSJed Brown   ierr = PetscBLASIntCast(appctx->param.N,&n);CHKERRQ(ierr);
400c4762a1bSJed Brown   for (j=xs; j<xs+xn; j += appctx->param.N-1) {
401c4762a1bSJed Brown     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&_DOne,&temp[0][0],&n,&xl[j],&_One,&_DOne,&yl[j],&_One));
402c4762a1bSJed Brown   }
403c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(appctx->da,xlocal,(void*)&xl);CHKERRQ(ierr);
404c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(appctx->da,ylocal,&yl);CHKERRQ(ierr);
405c4762a1bSJed Brown   ierr = PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
406c4762a1bSJed Brown   ierr = VecSet(y,0.0);CHKERRQ(ierr);
407c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(appctx->da,ylocal,ADD_VALUES,y);CHKERRQ(ierr);
408c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(appctx->da,ylocal,ADD_VALUES,y);CHKERRQ(ierr);
409c4762a1bSJed Brown   ierr = DMRestoreLocalVector(appctx->da,&xlocal);CHKERRQ(ierr);
410c4762a1bSJed Brown   ierr = DMRestoreLocalVector(appctx->da,&ylocal);CHKERRQ(ierr);
411c4762a1bSJed Brown   ierr = VecPointwiseDivide(y,y,appctx->SEMop.mass);CHKERRQ(ierr);
412c4762a1bSJed Brown   ierr = VecScale(y,-1.0);CHKERRQ(ierr);
413c4762a1bSJed Brown   return 0;
414c4762a1bSJed Brown }
415c4762a1bSJed Brown 
416c4762a1bSJed Brown /*
417c4762a1bSJed Brown    RHSMatrixLaplacian - User-provided routine to compute the right-hand-side
418c4762a1bSJed Brown    matrix for the Laplacian operator
419c4762a1bSJed Brown 
420c4762a1bSJed Brown    Input Parameters:
421c4762a1bSJed Brown    ts - the TS context
422c4762a1bSJed Brown    t - current time  (ignored)
423c4762a1bSJed Brown    X - current solution (ignored)
424c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
425c4762a1bSJed Brown 
426c4762a1bSJed Brown    Output Parameters:
427c4762a1bSJed Brown    AA - Jacobian matrix
428c4762a1bSJed Brown    BB - optionally different matrix from which the preconditioner is built
429c4762a1bSJed Brown    str - flag indicating matrix structure
430c4762a1bSJed Brown 
431c4762a1bSJed Brown */
432c4762a1bSJed Brown PetscErrorCode RHSMatrixLaplaciangllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx)
433c4762a1bSJed Brown {
434c4762a1bSJed Brown   PetscReal      **temp;
435c4762a1bSJed Brown   PetscReal      vv;
436c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
437c4762a1bSJed Brown   PetscErrorCode ierr;
438c4762a1bSJed Brown   PetscInt       i,xs,xn,l,j;
439c4762a1bSJed Brown   PetscInt       *rowsDM;
440c4762a1bSJed Brown   PetscBool      flg = PETSC_FALSE;
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-gll_mf",&flg,NULL);CHKERRQ(ierr);
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   if (!flg) {
445c4762a1bSJed Brown     /*
446c4762a1bSJed Brown      Creates the element stiffness matrix for the given gll
447c4762a1bSJed Brown      */
448c4762a1bSJed Brown     ierr = PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
449*a5b23f4aSJose E. Roman     /* workaround for clang analyzer warning: Division by zero */
450c4762a1bSJed Brown     if (appctx->param.N <= 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Spectral element order should be > 1");
451c4762a1bSJed Brown 
452c4762a1bSJed Brown     /* scale by the size of the element */
453c4762a1bSJed Brown     for (i=0; i<appctx->param.N; i++) {
454c4762a1bSJed Brown       vv=-appctx->param.mu*2.0/appctx->param.Le;
455c4762a1bSJed Brown       for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv;
456c4762a1bSJed Brown     }
457c4762a1bSJed Brown 
458c4762a1bSJed Brown     ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
459c4762a1bSJed Brown     ierr = DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr);
460c4762a1bSJed Brown 
461c4762a1bSJed Brown     xs   = xs/(appctx->param.N-1);
462c4762a1bSJed Brown     xn   = xn/(appctx->param.N-1);
463c4762a1bSJed Brown 
464c4762a1bSJed Brown     ierr = PetscMalloc1(appctx->param.N,&rowsDM);CHKERRQ(ierr);
465c4762a1bSJed Brown     /*
466c4762a1bSJed Brown      loop over local elements
467c4762a1bSJed Brown      */
468c4762a1bSJed Brown     for (j=xs; j<xs+xn; j++) {
469c4762a1bSJed Brown       for (l=0; l<appctx->param.N; l++) {
470c4762a1bSJed Brown         rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l;
471c4762a1bSJed Brown       }
472c4762a1bSJed Brown       ierr = MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES);CHKERRQ(ierr);
473c4762a1bSJed Brown     }
474c4762a1bSJed Brown     ierr = PetscFree(rowsDM);CHKERRQ(ierr);
475c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
476c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
477c4762a1bSJed Brown     ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr);
478c4762a1bSJed Brown     ierr = MatDiagonalScale(A,appctx->SEMop.mass,0);CHKERRQ(ierr);
479c4762a1bSJed Brown     ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr);
480c4762a1bSJed Brown 
481c4762a1bSJed Brown     ierr = PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
482c4762a1bSJed Brown   } else {
483c4762a1bSJed Brown     ierr = MatSetType(A,MATSHELL);CHKERRQ(ierr);
484c4762a1bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
485c4762a1bSJed Brown     ierr = MatShellSetContext(A,appctx);CHKERRQ(ierr);
486c4762a1bSJed Brown     ierr = MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_Laplacian);CHKERRQ(ierr);
487c4762a1bSJed Brown   }
488c4762a1bSJed Brown   return 0;
489c4762a1bSJed Brown }
490c4762a1bSJed Brown 
491c4762a1bSJed Brown /*
492c4762a1bSJed Brown    RHSMatrixAdvection - User-provided routine to compute the right-hand-side
493c4762a1bSJed Brown    matrix for the Advection (gradient) operator.
494c4762a1bSJed Brown 
495c4762a1bSJed Brown    Input Parameters:
496c4762a1bSJed Brown    ts - the TS context
497c4762a1bSJed Brown    t - current time
498c4762a1bSJed Brown    global_in - global input vector
499c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
500c4762a1bSJed Brown 
501c4762a1bSJed Brown    Output Parameters:
502c4762a1bSJed Brown    AA - Jacobian matrix
503c4762a1bSJed Brown    BB - optionally different preconditioning matrix
504c4762a1bSJed Brown    str - flag indicating matrix structure
505c4762a1bSJed Brown 
506c4762a1bSJed Brown */
507c4762a1bSJed Brown PetscErrorCode RHSMatrixAdvectiongllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx)
508c4762a1bSJed Brown {
509c4762a1bSJed Brown   PetscReal      **temp;
510c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
511c4762a1bSJed Brown   PetscErrorCode ierr;
512c4762a1bSJed Brown   PetscInt       xs,xn,l,j;
513c4762a1bSJed Brown   PetscInt       *rowsDM;
514c4762a1bSJed Brown   PetscBool      flg = PETSC_FALSE;
515c4762a1bSJed Brown 
516c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-gll_mf",&flg,NULL);CHKERRQ(ierr);
517c4762a1bSJed Brown 
518c4762a1bSJed Brown   if (!flg) {
519c4762a1bSJed Brown     /*
520c4762a1bSJed Brown      Creates the advection matrix for the given gll
521c4762a1bSJed Brown      */
522c4762a1bSJed Brown     ierr = PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
523c4762a1bSJed Brown     ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
524c4762a1bSJed Brown     ierr = DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr);
525c4762a1bSJed Brown     xs   = xs/(appctx->param.N-1);
526c4762a1bSJed Brown     xn   = xn/(appctx->param.N-1);
527c4762a1bSJed Brown 
528c4762a1bSJed Brown     ierr = PetscMalloc1(appctx->param.N,&rowsDM);CHKERRQ(ierr);
529c4762a1bSJed Brown     for (j=xs; j<xs+xn; j++) {
530c4762a1bSJed Brown       for (l=0; l<appctx->param.N; l++) {
531c4762a1bSJed Brown         rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l;
532c4762a1bSJed Brown       }
533c4762a1bSJed Brown       ierr = MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES);CHKERRQ(ierr);
534c4762a1bSJed Brown     }
535c4762a1bSJed Brown     ierr = PetscFree(rowsDM);CHKERRQ(ierr);
536c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
537c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
538c4762a1bSJed Brown 
539c4762a1bSJed Brown     ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr);
540c4762a1bSJed Brown     ierr = MatDiagonalScale(A,appctx->SEMop.mass,0);CHKERRQ(ierr);
541c4762a1bSJed Brown     ierr = VecReciprocal(appctx->SEMop.mass);CHKERRQ(ierr);
542c4762a1bSJed Brown 
543c4762a1bSJed Brown     ierr = PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp);CHKERRQ(ierr);
544c4762a1bSJed Brown   } else {
545c4762a1bSJed Brown     ierr = MatSetType(A,MATSHELL);CHKERRQ(ierr);
546c4762a1bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
547c4762a1bSJed Brown     ierr = MatShellSetContext(A,appctx);CHKERRQ(ierr);
548c4762a1bSJed Brown     ierr = MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_Advection);CHKERRQ(ierr);
549c4762a1bSJed Brown   }
550c4762a1bSJed Brown   return 0;
551c4762a1bSJed Brown }
552c4762a1bSJed Brown 
553c4762a1bSJed Brown /*TEST
554c4762a1bSJed Brown 
555c4762a1bSJed Brown     build:
556c4762a1bSJed Brown       requires: !complex
557c4762a1bSJed Brown 
558c4762a1bSJed Brown     test:
559c4762a1bSJed Brown       suffix: 1
560c4762a1bSJed Brown       requires: !single
561c4762a1bSJed Brown 
562c4762a1bSJed Brown     test:
563c4762a1bSJed Brown       suffix: 2
564c4762a1bSJed Brown       nsize: 5
565c4762a1bSJed Brown       requires: !single
566c4762a1bSJed Brown 
567c4762a1bSJed Brown     test:
568c4762a1bSJed Brown       suffix: 3
569c4762a1bSJed Brown       requires: !single
570c4762a1bSJed Brown       args: -ts_view  -ts_type beuler -gll_mf -pc_type none -ts_max_steps 5 -ts_monitor_error
571c4762a1bSJed Brown 
572c4762a1bSJed Brown     test:
573c4762a1bSJed Brown       suffix: 4
574c4762a1bSJed Brown       requires: !single
575c4762a1bSJed Brown       args: -ts_view  -ts_type beuler  -pc_type none -ts_max_steps 5 -ts_monitor_error
576c4762a1bSJed Brown 
577c4762a1bSJed Brown TEST*/
578