1c4762a1bSJed Brown static char help[] = "Solves one dimensional Burger's equation compares with exact solution\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown Not yet tested in parallel
5c4762a1bSJed Brown */
6c4762a1bSJed Brown
7c4762a1bSJed Brown /* ------------------------------------------------------------------------
8c4762a1bSJed Brown
9c4762a1bSJed Brown This program uses the one-dimensional Burger's equation
10c4762a1bSJed Brown u_t = mu*u_xx - u u_x,
11c4762a1bSJed Brown on the domain 0 <= x <= 1, with periodic boundary conditions
12c4762a1bSJed Brown
13c4762a1bSJed Brown The operators are discretized with the spectral element method
14c4762a1bSJed Brown
15c4762a1bSJed Brown See the paper PDE-CONSTRAINED OPTIMIZATION WITH SPECTRAL ELEMENTS USING PETSC AND TAO
16c4762a1bSJed Brown by OANA MARIN, EMIL CONSTANTINESCU, AND BARRY SMITH for details on the exact solution
17c4762a1bSJed Brown used
18c4762a1bSJed Brown
19c4762a1bSJed Brown See src/tao/unconstrained/tutorials/burgers_spectral.c
20c4762a1bSJed Brown
21c4762a1bSJed Brown ------------------------------------------------------------------------- */
22c4762a1bSJed Brown
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown #include <petscdt.h>
25c4762a1bSJed Brown #include <petscdraw.h>
26c4762a1bSJed Brown #include <petscdmda.h>
27c4762a1bSJed Brown
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown User-defined application context - contains data needed by the
30c4762a1bSJed Brown application-provided call-back routines.
31c4762a1bSJed Brown */
32c4762a1bSJed Brown
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown PetscInt n; /* number of nodes */
35c4762a1bSJed Brown PetscReal *nodes; /* GLL nodes */
36c4762a1bSJed Brown PetscReal *weights; /* GLL weights */
37c4762a1bSJed Brown } PetscGLL;
38c4762a1bSJed Brown
39c4762a1bSJed Brown typedef struct {
40c4762a1bSJed Brown PetscInt N; /* grid points per elements*/
41c4762a1bSJed Brown PetscInt E; /* number of elements */
42c4762a1bSJed Brown PetscReal tol_L2, tol_max; /* error norms */
43c4762a1bSJed Brown PetscInt steps; /* number of timesteps */
44c4762a1bSJed Brown PetscReal Tend; /* endtime */
45c4762a1bSJed Brown PetscReal mu; /* viscosity */
46c4762a1bSJed Brown PetscReal L; /* total length of domain */
47c4762a1bSJed Brown PetscReal Le;
48c4762a1bSJed Brown PetscReal Tadj;
49c4762a1bSJed Brown } PetscParam;
50c4762a1bSJed Brown
51c4762a1bSJed Brown typedef struct {
52c4762a1bSJed Brown Vec grid; /* total grid */
53c4762a1bSJed Brown Vec curr_sol;
54c4762a1bSJed Brown } PetscData;
55c4762a1bSJed Brown
56c4762a1bSJed Brown typedef struct {
57c4762a1bSJed Brown Vec grid; /* total grid */
58c4762a1bSJed Brown Vec mass; /* mass matrix for total integration */
59943f6d28SPierre Jolivet Mat stiff; /* stiffness matrix */
60c4762a1bSJed Brown Mat keptstiff;
61c4762a1bSJed Brown Mat grad;
62c4762a1bSJed Brown PetscGLL gll;
63c4762a1bSJed Brown } PetscSEMOperators;
64c4762a1bSJed Brown
65c4762a1bSJed Brown typedef struct {
66c4762a1bSJed Brown DM da; /* distributed array data structure */
67c4762a1bSJed Brown PetscSEMOperators SEMop;
68c4762a1bSJed Brown PetscParam param;
69c4762a1bSJed Brown PetscData dat;
70c4762a1bSJed Brown TS ts;
71c4762a1bSJed Brown PetscReal initial_dt;
72c4762a1bSJed Brown } AppCtx;
73c4762a1bSJed Brown
74c4762a1bSJed Brown /*
75c4762a1bSJed Brown User-defined routines
76c4762a1bSJed Brown */
77c4762a1bSJed Brown extern PetscErrorCode RHSMatrixLaplaciangllDM(TS, PetscReal, Vec, Mat, Mat, void *);
78c4762a1bSJed Brown extern PetscErrorCode RHSMatrixAdvectiongllDM(TS, PetscReal, Vec, Mat, Mat, void *);
79c4762a1bSJed Brown extern PetscErrorCode TrueSolution(TS, PetscReal, Vec, AppCtx *);
80c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *);
81c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *);
82c4762a1bSJed Brown
main(int argc,char ** argv)83d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
84d71ae5a4SJacob Faibussowitsch {
85c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */
86c4762a1bSJed Brown PetscInt i, xs, xm, ind, j, lenglob;
87c4762a1bSJed Brown PetscReal x, *wrk_ptr1, *wrk_ptr2;
88c4762a1bSJed Brown MatNullSpace nsp;
89c4762a1bSJed Brown PetscMPIInt size;
90c4762a1bSJed Brown
91c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown Initialize program and set problem parameters
93c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
947510d9b0SBarry Smith PetscFunctionBeginUser;
95c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
96c4762a1bSJed Brown
97c4762a1bSJed Brown /*initialize parameters */
98c4762a1bSJed Brown appctx.param.N = 10; /* order of the spectral element */
99c4762a1bSJed Brown appctx.param.E = 10; /* number of elements */
100c4762a1bSJed Brown appctx.param.L = 4.0; /* length of the domain */
101c4762a1bSJed Brown appctx.param.mu = 0.01; /* diffusion coefficient */
102c4762a1bSJed Brown appctx.initial_dt = 5e-3;
1031690c2aeSBarry Smith appctx.param.steps = PETSC_INT_MAX;
104c4762a1bSJed Brown appctx.param.Tend = 4;
105c4762a1bSJed Brown
1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &appctx.param.N, NULL));
1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-E", &appctx.param.E, NULL));
1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-Tend", &appctx.param.Tend, NULL));
1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &appctx.param.mu, NULL));
110c4762a1bSJed Brown appctx.param.Le = appctx.param.L / appctx.param.E;
111c4762a1bSJed Brown
1129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1133c633725SBarry Smith PetscCheck((appctx.param.E % size) == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of elements must be divisible by number of processes");
114c4762a1bSJed Brown
115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown Create GLL data structures
117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(appctx.param.N, &appctx.SEMop.gll.nodes, appctx.param.N, &appctx.SEMop.gll.weights));
1199566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights));
120c4762a1bSJed Brown appctx.SEMop.gll.n = appctx.param.N;
121c4762a1bSJed Brown lenglob = appctx.param.E * (appctx.param.N - 1);
122c4762a1bSJed Brown
123c4762a1bSJed Brown /*
124c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
125c4762a1bSJed Brown and to set up the ghost point communication pattern. There are E*(Nl-1)+1
126c4762a1bSJed Brown total grid values spread equally among all the processors, except first and last
127c4762a1bSJed Brown */
128c4762a1bSJed Brown
1299566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, lenglob, 1, 1, NULL, &appctx.da));
1309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(appctx.da));
1319566063dSJacob Faibussowitsch PetscCall(DMSetUp(appctx.da));
132c4762a1bSJed Brown
133c4762a1bSJed Brown /*
134c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the
135c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that
136c4762a1bSJed Brown have the same types.
137c4762a1bSJed Brown */
138c4762a1bSJed Brown
1399566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(appctx.da, &appctx.dat.curr_sol));
1409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.dat.curr_sol, &appctx.SEMop.grid));
1419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.dat.curr_sol, &appctx.SEMop.mass));
142c4762a1bSJed Brown
1439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx.da, &xs, NULL, NULL, &xm, NULL, NULL));
1449566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1));
1459566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2));
146c4762a1bSJed Brown
147c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
148c4762a1bSJed Brown
149c4762a1bSJed Brown xs = xs / (appctx.param.N - 1);
150c4762a1bSJed Brown xm = xm / (appctx.param.N - 1);
151c4762a1bSJed Brown
152c4762a1bSJed Brown /*
153c4762a1bSJed Brown Build total grid and mass over entire mesh (multi-elemental)
154c4762a1bSJed Brown */
155c4762a1bSJed Brown
156c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
157c4762a1bSJed Brown for (j = 0; j < appctx.param.N - 1; j++) {
158c4762a1bSJed Brown x = (appctx.param.Le / 2.0) * (appctx.SEMop.gll.nodes[j] + 1.0) + appctx.param.Le * i;
159c4762a1bSJed Brown ind = i * (appctx.param.N - 1) + j;
160c4762a1bSJed Brown wrk_ptr1[ind] = x;
161c4762a1bSJed Brown wrk_ptr2[ind] = .5 * appctx.param.Le * appctx.SEMop.gll.weights[j];
162c4762a1bSJed Brown if (j == 0) wrk_ptr2[ind] += .5 * appctx.param.Le * appctx.SEMop.gll.weights[j];
163c4762a1bSJed Brown }
164c4762a1bSJed Brown }
1659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1));
1669566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2));
167c4762a1bSJed Brown
168c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine.
170c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1719566063dSJacob Faibussowitsch PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE));
1729566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.stiff));
1739566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.grad));
174c4762a1bSJed Brown /*
175c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation
176dd8e379bSPierre Jolivet u_t = f(u,t), the user provides the discretized right-hand side
177c4762a1bSJed Brown as a time-dependent matrix.
178c4762a1bSJed Brown */
1799566063dSJacob Faibussowitsch PetscCall(RHSMatrixLaplaciangllDM(appctx.ts, 0.0, appctx.dat.curr_sol, appctx.SEMop.stiff, appctx.SEMop.stiff, &appctx));
1809566063dSJacob Faibussowitsch PetscCall(RHSMatrixAdvectiongllDM(appctx.ts, 0.0, appctx.dat.curr_sol, appctx.SEMop.grad, appctx.SEMop.grad, &appctx));
181c4762a1bSJed Brown /*
182c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation
183dd8e379bSPierre Jolivet u_t = f(u,t), the user provides the discretized right-hand side
184c4762a1bSJed Brown as a time-dependent matrix.
185c4762a1bSJed Brown */
186c4762a1bSJed Brown
1879566063dSJacob Faibussowitsch PetscCall(MatDuplicate(appctx.SEMop.stiff, MAT_COPY_VALUES, &appctx.SEMop.keptstiff));
188c4762a1bSJed Brown
189c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */
1909566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp));
1919566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.stiff, nsp));
1929566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.keptstiff, nsp));
1939566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.stiff, NULL));
1949566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp));
195c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */
1969566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp));
1979566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.grad, nsp));
1989566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.grad, NULL));
1999566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp));
200c4762a1bSJed Brown
201c4762a1bSJed Brown /* Create the TS solver that solves the ODE and its adjoint; set its options */
2029566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts));
2039566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(appctx.ts, TS_NONLINEAR));
2049566063dSJacob Faibussowitsch PetscCall(TSSetType(appctx.ts, TSRK));
2059566063dSJacob Faibussowitsch PetscCall(TSSetDM(appctx.ts, appctx.da));
2069566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts, 0.0));
2079566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx.ts, appctx.initial_dt));
2089566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(appctx.ts, appctx.param.steps));
2099566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(appctx.ts, appctx.param.Tend));
2109566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP));
2119566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(appctx.ts, 1e-7, NULL, 1e-7, NULL));
2129566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(appctx.ts));
2139566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts));
2149566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(appctx.ts, NULL, RHSFunction, &appctx));
2159566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts, appctx.SEMop.stiff, appctx.SEMop.stiff, RHSJacobian, &appctx));
216c4762a1bSJed Brown
217c4762a1bSJed Brown /* Set Initial conditions for the problem */
2189566063dSJacob Faibussowitsch PetscCall(TrueSolution(appctx.ts, 0, appctx.dat.curr_sol, &appctx));
219c4762a1bSJed Brown
2209566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(appctx.ts, (PetscErrorCode (*)(TS, PetscReal, Vec, void *))TrueSolution, &appctx));
2219566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts, 0.0));
2229566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(appctx.ts, 0));
223c4762a1bSJed Brown
2249566063dSJacob Faibussowitsch PetscCall(TSSolve(appctx.ts, appctx.dat.curr_sol));
225c4762a1bSJed Brown
2269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.stiff));
2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.keptstiff));
2289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.grad));
2299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.grid));
2309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.mass));
2319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.curr_sol));
2329566063dSJacob Faibussowitsch PetscCall(PetscFree2(appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights));
2339566063dSJacob Faibussowitsch PetscCall(DMDestroy(&appctx.da));
2349566063dSJacob Faibussowitsch PetscCall(TSDestroy(&appctx.ts));
235c4762a1bSJed Brown
236c4762a1bSJed Brown /*
237c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine
238c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI
239c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime
240d75802c7SJacob Faibussowitsch options are chosen (e.g., -log_view).
241c4762a1bSJed Brown */
2429566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
243b122ec5aSJacob Faibussowitsch return 0;
244c4762a1bSJed Brown }
245c4762a1bSJed Brown
246c4762a1bSJed Brown /*
247c4762a1bSJed Brown TrueSolution() computes the true solution for the PDE
248c4762a1bSJed Brown
249c4762a1bSJed Brown Input Parameter:
250c4762a1bSJed Brown u - uninitialized solution vector (global)
251c4762a1bSJed Brown appctx - user-defined application context
252c4762a1bSJed Brown
253c4762a1bSJed Brown Output Parameter:
254c4762a1bSJed Brown u - vector with solution at initial time (global)
255c4762a1bSJed Brown */
TrueSolution(TS ts,PetscReal t,Vec u,AppCtx * appctx)256d71ae5a4SJacob Faibussowitsch PetscErrorCode TrueSolution(TS ts, PetscReal t, Vec u, AppCtx *appctx)
257d71ae5a4SJacob Faibussowitsch {
258c4762a1bSJed Brown PetscScalar *s;
259c4762a1bSJed Brown const PetscScalar *xg;
260c4762a1bSJed Brown PetscInt i, xs, xn;
261c4762a1bSJed Brown
2623ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
2639566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, u, &s));
2649566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
2659566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
266c4762a1bSJed Brown for (i = xs; i < xs + xn; i++) {
267c4762a1bSJed 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));
268c4762a1bSJed Brown }
2699566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, u, &s));
2709566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg));
2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
272c4762a1bSJed Brown }
273c4762a1bSJed Brown
RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,PetscCtx ctx)274*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, PetscCtx ctx)
275d71ae5a4SJacob Faibussowitsch {
276c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
277c4762a1bSJed Brown
2787510d9b0SBarry Smith PetscFunctionBeginUser;
2799566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.grad, globalin, globalout)); /* grad u */
2809566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(globalout, globalin, globalout)); /* u grad u */
2819566063dSJacob Faibussowitsch PetscCall(VecScale(globalout, -1.0));
2829566063dSJacob Faibussowitsch PetscCall(MatMultAdd(appctx->SEMop.keptstiff, globalin, globalout, globalout));
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
284c4762a1bSJed Brown }
285c4762a1bSJed Brown
286c4762a1bSJed Brown /*
287c4762a1bSJed Brown
288c4762a1bSJed Brown K is the discretiziation of the Laplacian
289c4762a1bSJed Brown G is the discretization of the gradient
290c4762a1bSJed Brown
291c4762a1bSJed Brown Computes Jacobian of K u + diag(u) G u which is given by
292c4762a1bSJed Brown K + diag(u)G + diag(Gu)
293c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A,Mat B,PetscCtx ctx)294*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec globalin, Mat A, Mat B, PetscCtx ctx)
295d71ae5a4SJacob Faibussowitsch {
296c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
297c4762a1bSJed Brown Vec Gglobalin;
298c4762a1bSJed Brown
2997510d9b0SBarry Smith PetscFunctionBeginUser;
300c4762a1bSJed Brown /* A = diag(u) G */
301c4762a1bSJed Brown
3029566063dSJacob Faibussowitsch PetscCall(MatCopy(appctx->SEMop.grad, A, SAME_NONZERO_PATTERN));
3039566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, globalin, NULL));
304c4762a1bSJed Brown
305c4762a1bSJed Brown /* A = A + diag(Gu) */
3069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(globalin, &Gglobalin));
3079566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.grad, globalin, Gglobalin));
3089566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, Gglobalin, ADD_VALUES));
3099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Gglobalin));
310c4762a1bSJed Brown
311c4762a1bSJed Brown /* A = K - A */
3129566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0));
3139566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, 0.0, appctx->SEMop.keptstiff, SAME_NONZERO_PATTERN));
3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
315c4762a1bSJed Brown }
316c4762a1bSJed Brown
31746233b44SBarry Smith #include <petscblaslapack.h>
318c4762a1bSJed Brown /*
319c4762a1bSJed Brown Matrix free operation of 1d Laplacian and Grad for GLL spectral elements
320c4762a1bSJed Brown */
MatMult_Laplacian(Mat A,Vec x,Vec y)321d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Laplacian(Mat A, Vec x, Vec y)
322d71ae5a4SJacob Faibussowitsch {
323c4762a1bSJed Brown AppCtx *appctx;
324c4762a1bSJed Brown PetscReal **temp, vv;
325c4762a1bSJed Brown PetscInt i, j, xs, xn;
326c4762a1bSJed Brown Vec xlocal, ylocal;
327c4762a1bSJed Brown const PetscScalar *xl;
328c4762a1bSJed Brown PetscScalar *yl;
329c4762a1bSJed Brown PetscBLASInt _One = 1, n;
330c4762a1bSJed Brown PetscScalar _DOne = 1;
331c4762a1bSJed Brown
3323ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
3339566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &appctx));
3349566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &xlocal));
3359566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(appctx->da, x, INSERT_VALUES, xlocal));
3369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(appctx->da, x, INSERT_VALUES, xlocal));
3379566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &ylocal));
3389566063dSJacob Faibussowitsch PetscCall(VecSet(ylocal, 0.0));
3399566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
340c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) {
341c4762a1bSJed Brown vv = -appctx->param.mu * 2.0 / appctx->param.Le;
342c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv;
343c4762a1bSJed Brown }
3449566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, xlocal, (void *)&xl));
3459566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, ylocal, &yl));
3469566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
3479566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(appctx->param.N, &n));
34848a46eb9SPierre Jolivet 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));
3499566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, xlocal, (void *)&xl));
3509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, ylocal, &yl));
3519566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
3529566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0));
3539566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(appctx->da, ylocal, ADD_VALUES, y));
3549566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(appctx->da, ylocal, ADD_VALUES, y));
3559566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &xlocal));
3569566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &ylocal));
3579566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(y, y, appctx->SEMop.mass));
3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
359c4762a1bSJed Brown }
360c4762a1bSJed Brown
MatMult_Advection(Mat A,Vec x,Vec y)361d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Advection(Mat A, Vec x, Vec y)
362d71ae5a4SJacob Faibussowitsch {
363c4762a1bSJed Brown AppCtx *appctx;
364c4762a1bSJed Brown PetscReal **temp;
365c4762a1bSJed Brown PetscInt j, xs, xn;
366c4762a1bSJed Brown Vec xlocal, ylocal;
367c4762a1bSJed Brown const PetscScalar *xl;
368c4762a1bSJed Brown PetscScalar *yl;
369c4762a1bSJed Brown PetscBLASInt _One = 1, n;
370c4762a1bSJed Brown PetscScalar _DOne = 1;
371c4762a1bSJed Brown
3723ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
3739566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &appctx));
3749566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &xlocal));
3759566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(appctx->da, x, INSERT_VALUES, xlocal));
3769566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(appctx->da, x, INSERT_VALUES, xlocal));
3779566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &ylocal));
3789566063dSJacob Faibussowitsch PetscCall(VecSet(ylocal, 0.0));
3799566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
3809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, xlocal, (void *)&xl));
3819566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, ylocal, &yl));
3829566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
3839566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(appctx->param.N, &n));
38448a46eb9SPierre Jolivet 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));
3859566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, xlocal, (void *)&xl));
3869566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, ylocal, &yl));
3879566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
3889566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0));
3899566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(appctx->da, ylocal, ADD_VALUES, y));
3909566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(appctx->da, ylocal, ADD_VALUES, y));
3919566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &xlocal));
3929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &ylocal));
3939566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(y, y, appctx->SEMop.mass));
3949566063dSJacob Faibussowitsch PetscCall(VecScale(y, -1.0));
3953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
396c4762a1bSJed Brown }
397c4762a1bSJed Brown
398c4762a1bSJed Brown /*
399c4762a1bSJed Brown RHSMatrixLaplacian - User-provided routine to compute the right-hand-side
400c4762a1bSJed Brown matrix for the Laplacian operator
401c4762a1bSJed Brown
402c4762a1bSJed Brown Input Parameters:
403c4762a1bSJed Brown ts - the TS context
404c4762a1bSJed Brown t - current time (ignored)
405c4762a1bSJed Brown X - current solution (ignored)
406c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian()
407c4762a1bSJed Brown
408c4762a1bSJed Brown Output Parameters:
409c4762a1bSJed Brown AA - Jacobian matrix
410c4762a1bSJed Brown BB - optionally different matrix from which the preconditioner is built
411c4762a1bSJed Brown
412c4762a1bSJed Brown */
RHSMatrixLaplaciangllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,PetscCtx ctx)413*2a8381b2SBarry Smith PetscErrorCode RHSMatrixLaplaciangllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, PetscCtx ctx)
414d71ae5a4SJacob Faibussowitsch {
415c4762a1bSJed Brown PetscReal **temp;
416c4762a1bSJed Brown PetscReal vv;
417c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
418c4762a1bSJed Brown PetscInt i, xs, xn, l, j;
419c4762a1bSJed Brown PetscInt *rowsDM;
420c4762a1bSJed Brown PetscBool flg = PETSC_FALSE;
421c4762a1bSJed Brown
4223ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-gll_mf", &flg, NULL));
424c4762a1bSJed Brown
425c4762a1bSJed Brown if (!flg) {
426c4762a1bSJed Brown /*
427c4762a1bSJed Brown Creates the element stiffness matrix for the given gll
428c4762a1bSJed Brown */
4299566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
430a5b23f4aSJose E. Roman /* workaround for clang analyzer warning: Division by zero */
4313c633725SBarry Smith PetscCheck(appctx->param.N > 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spectral element order should be > 1");
432c4762a1bSJed Brown
433c4762a1bSJed Brown /* scale by the size of the element */
434c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) {
435c4762a1bSJed Brown vv = -appctx->param.mu * 2.0 / appctx->param.Le;
436c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv;
437c4762a1bSJed Brown }
438c4762a1bSJed Brown
4399566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
4409566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
441c4762a1bSJed Brown
442c4762a1bSJed Brown xs = xs / (appctx->param.N - 1);
443c4762a1bSJed Brown xn = xn / (appctx->param.N - 1);
444c4762a1bSJed Brown
4459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM));
446c4762a1bSJed Brown /*
447c4762a1bSJed Brown loop over local elements
448c4762a1bSJed Brown */
449c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) {
450ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l;
4519566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES));
452c4762a1bSJed Brown }
4539566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM));
4549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4569566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
4579566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0));
4589566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
459c4762a1bSJed Brown
4609566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
461c4762a1bSJed Brown } else {
4629566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL));
4639566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
4649566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, appctx));
46557d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Laplacian));
466c4762a1bSJed Brown }
4673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
468c4762a1bSJed Brown }
469c4762a1bSJed Brown
470c4762a1bSJed Brown /*
471c4762a1bSJed Brown RHSMatrixAdvection - User-provided routine to compute the right-hand-side
472c4762a1bSJed Brown matrix for the Advection (gradient) operator.
473c4762a1bSJed Brown
474c4762a1bSJed Brown Input Parameters:
475c4762a1bSJed Brown ts - the TS context
476c4762a1bSJed Brown t - current time
477c4762a1bSJed Brown global_in - global input vector
478c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian()
479c4762a1bSJed Brown
480c4762a1bSJed Brown Output Parameters:
481c4762a1bSJed Brown AA - Jacobian matrix
4827addb90fSBarry Smith BB - optionally different matrix used to construct the preconditioner
483c4762a1bSJed Brown
484c4762a1bSJed Brown */
RHSMatrixAdvectiongllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,PetscCtx ctx)485*2a8381b2SBarry Smith PetscErrorCode RHSMatrixAdvectiongllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, PetscCtx ctx)
486d71ae5a4SJacob Faibussowitsch {
487c4762a1bSJed Brown PetscReal **temp;
488c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
489c4762a1bSJed Brown PetscInt xs, xn, l, j;
490c4762a1bSJed Brown PetscInt *rowsDM;
491c4762a1bSJed Brown PetscBool flg = PETSC_FALSE;
492c4762a1bSJed Brown
4933ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
4949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-gll_mf", &flg, NULL));
495c4762a1bSJed Brown
496c4762a1bSJed Brown if (!flg) {
497c4762a1bSJed Brown /*
498c4762a1bSJed Brown Creates the advection matrix for the given gll
499c4762a1bSJed Brown */
5009566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
5019566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
5029566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL));
503c4762a1bSJed Brown xs = xs / (appctx->param.N - 1);
504c4762a1bSJed Brown xn = xn / (appctx->param.N - 1);
505c4762a1bSJed Brown
5069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM));
507c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) {
508ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l;
5099566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES));
510c4762a1bSJed Brown }
5119566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM));
5129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
514c4762a1bSJed Brown
5159566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
5169566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0));
5179566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass));
518c4762a1bSJed Brown
5199566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp));
520c4762a1bSJed Brown } else {
5219566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL));
5229566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
5239566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, appctx));
52457d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Advection));
525c4762a1bSJed Brown }
5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
527c4762a1bSJed Brown }
528c4762a1bSJed Brown
529c4762a1bSJed Brown /*TEST
530c4762a1bSJed Brown
531c4762a1bSJed Brown build:
532c4762a1bSJed Brown requires: !complex
533c4762a1bSJed Brown
534c4762a1bSJed Brown test:
535c4762a1bSJed Brown suffix: 1
536c4762a1bSJed Brown requires: !single
5373886731fSPierre Jolivet output_file: output/empty.out
538c4762a1bSJed Brown
539c4762a1bSJed Brown test:
540c4762a1bSJed Brown suffix: 2
541c4762a1bSJed Brown nsize: 5
542c4762a1bSJed Brown requires: !single
5433886731fSPierre Jolivet output_file: output/empty.out
544c4762a1bSJed Brown
545c4762a1bSJed Brown test:
546c4762a1bSJed Brown suffix: 3
547c4762a1bSJed Brown requires: !single
548c4762a1bSJed Brown args: -ts_view -ts_type beuler -gll_mf -pc_type none -ts_max_steps 5 -ts_monitor_error
549c4762a1bSJed Brown
550c4762a1bSJed Brown test:
551c4762a1bSJed Brown suffix: 4
552c4762a1bSJed Brown requires: !single
553c4762a1bSJed Brown args: -ts_view -ts_type beuler -pc_type none -ts_max_steps 5 -ts_monitor_error
554c4762a1bSJed Brown
555c4762a1bSJed Brown TEST*/
556