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