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 /* ------------------------------------------------------------------------ 11c4762a1bSJed Brown 12c4762a1bSJed Brown This program uses the one-dimensional Burger's equation 13c4762a1bSJed Brown u_t = mu*u_xx - u u_x, 14c4762a1bSJed Brown on the domain 0 <= x <= 1, with periodic boundary conditions 15c4762a1bSJed Brown 16c4762a1bSJed Brown The operators are discretized with the spectral element method 17c4762a1bSJed Brown 18c4762a1bSJed Brown See the paper PDE-CONSTRAINED OPTIMIZATION WITH SPECTRAL ELEMENTS USING PETSC AND TAO 19c4762a1bSJed Brown by OANA MARIN, EMIL CONSTANTINESCU, AND BARRY SMITH for details on the exact solution 20c4762a1bSJed Brown used 21c4762a1bSJed Brown 22c4762a1bSJed Brown See src/tao/unconstrained/tutorials/burgers_spectral.c 23c4762a1bSJed Brown 24c4762a1bSJed Brown ------------------------------------------------------------------------- */ 25c4762a1bSJed Brown 26c4762a1bSJed Brown #include <petscts.h> 27c4762a1bSJed Brown #include <petscdt.h> 28c4762a1bSJed Brown #include <petscdraw.h> 29c4762a1bSJed Brown #include <petscdmda.h> 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* 32c4762a1bSJed Brown User-defined application context - contains data needed by the 33c4762a1bSJed Brown application-provided call-back routines. 34c4762a1bSJed Brown */ 35c4762a1bSJed Brown 36c4762a1bSJed Brown typedef struct { 37c4762a1bSJed Brown PetscInt n; /* number of nodes */ 38c4762a1bSJed Brown PetscReal *nodes; /* GLL nodes */ 39c4762a1bSJed Brown PetscReal *weights; /* GLL weights */ 40c4762a1bSJed Brown } PetscGLL; 41c4762a1bSJed Brown 42c4762a1bSJed Brown typedef struct { 43c4762a1bSJed Brown PetscInt N; /* grid points per elements*/ 44c4762a1bSJed Brown PetscInt E; /* number of elements */ 45c4762a1bSJed Brown PetscReal tol_L2, tol_max; /* error norms */ 46c4762a1bSJed Brown PetscInt steps; /* number of timesteps */ 47c4762a1bSJed Brown PetscReal Tend; /* endtime */ 48c4762a1bSJed Brown PetscReal mu; /* viscosity */ 49c4762a1bSJed Brown PetscReal L; /* total length of domain */ 50c4762a1bSJed Brown PetscReal Le; 51c4762a1bSJed Brown PetscReal Tadj; 52c4762a1bSJed Brown } PetscParam; 53c4762a1bSJed Brown 54c4762a1bSJed Brown typedef struct { 55c4762a1bSJed Brown Vec grid; /* total grid */ 56c4762a1bSJed Brown Vec curr_sol; 57c4762a1bSJed Brown } PetscData; 58c4762a1bSJed Brown 59c4762a1bSJed Brown typedef struct { 60c4762a1bSJed Brown Vec grid; /* total grid */ 61c4762a1bSJed Brown Vec mass; /* mass matrix for total integration */ 62c4762a1bSJed Brown Mat stiff; /* stifness matrix */ 63c4762a1bSJed Brown Mat keptstiff; 64c4762a1bSJed Brown Mat grad; 65c4762a1bSJed Brown PetscGLL gll; 66c4762a1bSJed Brown } PetscSEMOperators; 67c4762a1bSJed Brown 68c4762a1bSJed Brown typedef struct { 69c4762a1bSJed Brown DM da; /* distributed array data structure */ 70c4762a1bSJed Brown PetscSEMOperators SEMop; 71c4762a1bSJed Brown PetscParam param; 72c4762a1bSJed Brown PetscData dat; 73c4762a1bSJed Brown TS ts; 74c4762a1bSJed Brown PetscReal initial_dt; 75c4762a1bSJed Brown } AppCtx; 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* 78c4762a1bSJed Brown User-defined routines 79c4762a1bSJed Brown */ 80c4762a1bSJed Brown extern PetscErrorCode RHSMatrixLaplaciangllDM(TS, PetscReal, Vec, Mat, Mat, void *); 81c4762a1bSJed Brown extern PetscErrorCode RHSMatrixAdvectiongllDM(TS, PetscReal, Vec, Mat, Mat, void *); 82c4762a1bSJed Brown extern PetscErrorCode TrueSolution(TS, PetscReal, Vec, AppCtx *); 83c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *); 84c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *); 85c4762a1bSJed Brown 86d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 87d71ae5a4SJacob Faibussowitsch { 88c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 89c4762a1bSJed Brown PetscInt i, xs, xm, ind, j, lenglob; 90c4762a1bSJed Brown PetscReal x, *wrk_ptr1, *wrk_ptr2; 91c4762a1bSJed Brown MatNullSpace nsp; 92c4762a1bSJed Brown PetscMPIInt size; 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95c4762a1bSJed Brown Initialize program and set problem parameters 96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 977510d9b0SBarry Smith PetscFunctionBeginUser; 98c4762a1bSJed Brown 99327415f7SBarry Smith PetscFunctionBeginUser; 1009566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /*initialize parameters */ 103c4762a1bSJed Brown appctx.param.N = 10; /* order of the spectral element */ 104c4762a1bSJed Brown appctx.param.E = 10; /* number of elements */ 105c4762a1bSJed Brown appctx.param.L = 4.0; /* length of the domain */ 106c4762a1bSJed Brown appctx.param.mu = 0.01; /* diffusion coefficient */ 107c4762a1bSJed Brown appctx.initial_dt = 5e-3; 108c4762a1bSJed Brown appctx.param.steps = PETSC_MAX_INT; 109c4762a1bSJed Brown appctx.param.Tend = 4; 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &appctx.param.N, NULL)); 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-E", &appctx.param.E, NULL)); 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-Tend", &appctx.param.Tend, NULL)); 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &appctx.param.mu, NULL)); 115c4762a1bSJed Brown appctx.param.Le = appctx.param.L / appctx.param.E; 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1183c633725SBarry Smith PetscCheck((appctx.param.E % size) == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of elements must be divisible by number of processes"); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121c4762a1bSJed Brown Create GLL data structures 122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1239566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(appctx.param.N, &appctx.SEMop.gll.nodes, appctx.param.N, &appctx.SEMop.gll.weights)); 1249566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights)); 125c4762a1bSJed Brown appctx.SEMop.gll.n = appctx.param.N; 126c4762a1bSJed Brown lenglob = appctx.param.E * (appctx.param.N - 1); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* 129c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 130c4762a1bSJed Brown and to set up the ghost point communication pattern. There are E*(Nl-1)+1 131c4762a1bSJed Brown total grid values spread equally among all the processors, except first and last 132c4762a1bSJed Brown */ 133c4762a1bSJed Brown 1349566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, lenglob, 1, 1, NULL, &appctx.da)); 1359566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(appctx.da)); 1369566063dSJacob Faibussowitsch PetscCall(DMSetUp(appctx.da)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 140c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 141c4762a1bSJed Brown have the same types. 142c4762a1bSJed Brown */ 143c4762a1bSJed Brown 1449566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(appctx.da, &appctx.dat.curr_sol)); 1459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.dat.curr_sol, &appctx.SEMop.grid)); 1469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.dat.curr_sol, &appctx.SEMop.mass)); 147c4762a1bSJed Brown 1489566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx.da, &xs, NULL, NULL, &xm, NULL, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1)); 1509566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2)); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 153c4762a1bSJed Brown 154c4762a1bSJed Brown xs = xs / (appctx.param.N - 1); 155c4762a1bSJed Brown xm = xm / (appctx.param.N - 1); 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* 158c4762a1bSJed Brown Build total grid and mass over entire mesh (multi-elemental) 159c4762a1bSJed Brown */ 160c4762a1bSJed Brown 161c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 162c4762a1bSJed Brown for (j = 0; j < appctx.param.N - 1; j++) { 163c4762a1bSJed Brown x = (appctx.param.Le / 2.0) * (appctx.SEMop.gll.nodes[j] + 1.0) + appctx.param.Le * i; 164c4762a1bSJed Brown ind = i * (appctx.param.N - 1) + j; 165c4762a1bSJed Brown wrk_ptr1[ind] = x; 166c4762a1bSJed Brown wrk_ptr2[ind] = .5 * appctx.param.Le * appctx.SEMop.gll.weights[j]; 167c4762a1bSJed Brown if (j == 0) wrk_ptr2[ind] += .5 * appctx.param.Le * appctx.SEMop.gll.weights[j]; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown } 1709566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1)); 1719566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 174c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 175c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1769566063dSJacob Faibussowitsch PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE)); 1779566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.stiff)); 1789566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.grad)); 179c4762a1bSJed Brown /* 180c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 181c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 182c4762a1bSJed Brown as a time-dependent matrix. 183c4762a1bSJed Brown */ 1849566063dSJacob Faibussowitsch PetscCall(RHSMatrixLaplaciangllDM(appctx.ts, 0.0, appctx.dat.curr_sol, appctx.SEMop.stiff, appctx.SEMop.stiff, &appctx)); 1859566063dSJacob Faibussowitsch PetscCall(RHSMatrixAdvectiongllDM(appctx.ts, 0.0, appctx.dat.curr_sol, appctx.SEMop.grad, appctx.SEMop.grad, &appctx)); 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 188c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 189c4762a1bSJed Brown as a time-dependent matrix. 190c4762a1bSJed Brown */ 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(MatDuplicate(appctx.SEMop.stiff, MAT_COPY_VALUES, &appctx.SEMop.keptstiff)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 1959566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp)); 1969566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.stiff, nsp)); 1979566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.keptstiff, nsp)); 1989566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.stiff, NULL)); 1999566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp)); 200c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 2019566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp)); 2029566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.grad, nsp)); 2039566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.grad, NULL)); 2049566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* Create the TS solver that solves the ODE and its adjoint; set its options */ 2079566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts)); 2089566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(appctx.ts, TS_NONLINEAR)); 2099566063dSJacob Faibussowitsch PetscCall(TSSetType(appctx.ts, TSRK)); 2109566063dSJacob Faibussowitsch PetscCall(TSSetDM(appctx.ts, appctx.da)); 2119566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts, 0.0)); 2129566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx.ts, appctx.initial_dt)); 2139566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(appctx.ts, appctx.param.steps)); 2149566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(appctx.ts, appctx.param.Tend)); 2159566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP)); 2169566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(appctx.ts, 1e-7, NULL, 1e-7, NULL)); 2179566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(appctx.ts)); 2189566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts)); 2199566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(appctx.ts, NULL, RHSFunction, &appctx)); 2209566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts, appctx.SEMop.stiff, appctx.SEMop.stiff, RHSJacobian, &appctx)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* Set Initial conditions for the problem */ 2239566063dSJacob Faibussowitsch PetscCall(TrueSolution(appctx.ts, 0, appctx.dat.curr_sol, &appctx)); 224c4762a1bSJed Brown 2259566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(appctx.ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))TrueSolution, &appctx)); 2269566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts, 0.0)); 2279566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(appctx.ts, 0)); 228c4762a1bSJed Brown 2299566063dSJacob Faibussowitsch PetscCall(TSSolve(appctx.ts, appctx.dat.curr_sol)); 230c4762a1bSJed Brown 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.stiff)); 2329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.keptstiff)); 2339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.grad)); 2349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.grid)); 2359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.mass)); 2369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.curr_sol)); 2379566063dSJacob Faibussowitsch PetscCall(PetscFree2(appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights)); 2389566063dSJacob Faibussowitsch PetscCall(DMDestroy(&appctx.da)); 2399566063dSJacob Faibussowitsch PetscCall(TSDestroy(&appctx.ts)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* 242c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 243c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 244c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 245c4762a1bSJed Brown options are chosen (e.g., -log_summary). 246c4762a1bSJed Brown */ 2479566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 248b122ec5aSJacob Faibussowitsch return 0; 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* 252c4762a1bSJed Brown TrueSolution() computes the true solution for the PDE 253c4762a1bSJed Brown 254c4762a1bSJed Brown Input Parameter: 255c4762a1bSJed Brown u - uninitialized solution vector (global) 256c4762a1bSJed Brown appctx - user-defined application context 257c4762a1bSJed Brown 258c4762a1bSJed Brown Output Parameter: 259c4762a1bSJed Brown u - vector with solution at initial time (global) 260c4762a1bSJed Brown */ 261d71ae5a4SJacob Faibussowitsch PetscErrorCode TrueSolution(TS ts, PetscReal t, Vec u, AppCtx *appctx) 262d71ae5a4SJacob Faibussowitsch { 263c4762a1bSJed Brown PetscScalar *s; 264c4762a1bSJed Brown const PetscScalar *xg; 265c4762a1bSJed Brown PetscInt i, xs, xn; 266c4762a1bSJed Brown 267*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, u, &s)); 2699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 2709566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 271c4762a1bSJed Brown for (i = xs; i < xs + xn; i++) { 272c4762a1bSJed 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)); 273c4762a1bSJed Brown } 2749566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, u, &s)); 2759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 276*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) 280d71ae5a4SJacob Faibussowitsch { 281c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 282c4762a1bSJed Brown 2837510d9b0SBarry Smith PetscFunctionBeginUser; 2849566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.grad, globalin, globalout)); /* grad u */ 2859566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(globalout, globalin, globalout)); /* u grad u */ 2869566063dSJacob Faibussowitsch PetscCall(VecScale(globalout, -1.0)); 2879566063dSJacob Faibussowitsch PetscCall(MatMultAdd(appctx->SEMop.keptstiff, globalin, globalout, globalout)); 288*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* 292c4762a1bSJed Brown 293c4762a1bSJed Brown K is the discretiziation of the Laplacian 294c4762a1bSJed Brown G is the discretization of the gradient 295c4762a1bSJed Brown 296c4762a1bSJed Brown Computes Jacobian of K u + diag(u) G u which is given by 297c4762a1bSJed Brown K + diag(u)G + diag(Gu) 298c4762a1bSJed Brown */ 299d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec globalin, Mat A, Mat B, void *ctx) 300d71ae5a4SJacob Faibussowitsch { 301c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 302c4762a1bSJed Brown Vec Gglobalin; 303c4762a1bSJed Brown 3047510d9b0SBarry Smith PetscFunctionBeginUser; 305c4762a1bSJed Brown /* A = diag(u) G */ 306c4762a1bSJed Brown 3079566063dSJacob Faibussowitsch PetscCall(MatCopy(appctx->SEMop.grad, A, SAME_NONZERO_PATTERN)); 3089566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, globalin, NULL)); 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* A = A + diag(Gu) */ 3119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(globalin, &Gglobalin)); 3129566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.grad, globalin, Gglobalin)); 3139566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, Gglobalin, ADD_VALUES)); 3149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Gglobalin)); 315c4762a1bSJed Brown 316c4762a1bSJed Brown /* A = K - A */ 3179566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); 3189566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, 0.0, appctx->SEMop.keptstiff, SAME_NONZERO_PATTERN)); 319*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown 322c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 323c4762a1bSJed Brown 324c4762a1bSJed Brown #include "petscblaslapack.h" 325c4762a1bSJed Brown /* 326c4762a1bSJed Brown Matrix free operation of 1d Laplacian and Grad for GLL spectral elements 327c4762a1bSJed Brown */ 328d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Laplacian(Mat A, Vec x, Vec y) 329d71ae5a4SJacob Faibussowitsch { 330c4762a1bSJed Brown AppCtx *appctx; 331c4762a1bSJed Brown PetscReal **temp, vv; 332c4762a1bSJed Brown PetscInt i, j, xs, xn; 333c4762a1bSJed Brown Vec xlocal, ylocal; 334c4762a1bSJed Brown const PetscScalar *xl; 335c4762a1bSJed Brown PetscScalar *yl; 336c4762a1bSJed Brown PetscBLASInt _One = 1, n; 337c4762a1bSJed Brown PetscScalar _DOne = 1; 338c4762a1bSJed Brown 339*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 3409566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &appctx)); 3419566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &xlocal)); 3429566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(appctx->da, x, INSERT_VALUES, xlocal)); 3439566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(appctx->da, x, INSERT_VALUES, xlocal)); 3449566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &ylocal)); 3459566063dSJacob Faibussowitsch PetscCall(VecSet(ylocal, 0.0)); 3469566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 347c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) { 348c4762a1bSJed Brown vv = -appctx->param.mu * 2.0 / appctx->param.Le; 349c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv; 350c4762a1bSJed Brown } 3519566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, xlocal, (void *)&xl)); 3529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, ylocal, &yl)); 3539566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 3549566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(appctx->param.N, &n)); 35548a46eb9SPierre 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)); 3569566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, xlocal, (void *)&xl)); 3579566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, ylocal, &yl)); 3589566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 3599566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 3609566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(appctx->da, ylocal, ADD_VALUES, y)); 3619566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(appctx->da, ylocal, ADD_VALUES, y)); 3629566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &xlocal)); 3639566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &ylocal)); 3649566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(y, y, appctx->SEMop.mass)); 365*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 366c4762a1bSJed Brown } 367c4762a1bSJed Brown 368d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Advection(Mat A, Vec x, Vec y) 369d71ae5a4SJacob Faibussowitsch { 370c4762a1bSJed Brown AppCtx *appctx; 371c4762a1bSJed Brown PetscReal **temp; 372c4762a1bSJed Brown PetscInt j, xs, xn; 373c4762a1bSJed Brown Vec xlocal, ylocal; 374c4762a1bSJed Brown const PetscScalar *xl; 375c4762a1bSJed Brown PetscScalar *yl; 376c4762a1bSJed Brown PetscBLASInt _One = 1, n; 377c4762a1bSJed Brown PetscScalar _DOne = 1; 378c4762a1bSJed Brown 379*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 3809566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &appctx)); 3819566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &xlocal)); 3829566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(appctx->da, x, INSERT_VALUES, xlocal)); 3839566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(appctx->da, x, INSERT_VALUES, xlocal)); 3849566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &ylocal)); 3859566063dSJacob Faibussowitsch PetscCall(VecSet(ylocal, 0.0)); 3869566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 3879566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, xlocal, (void *)&xl)); 3889566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, ylocal, &yl)); 3899566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 3909566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(appctx->param.N, &n)); 39148a46eb9SPierre 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)); 3929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, xlocal, (void *)&xl)); 3939566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, ylocal, &yl)); 3949566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 3959566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 3969566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(appctx->da, ylocal, ADD_VALUES, y)); 3979566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(appctx->da, ylocal, ADD_VALUES, y)); 3989566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &xlocal)); 3999566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &ylocal)); 4009566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(y, y, appctx->SEMop.mass)); 4019566063dSJacob Faibussowitsch PetscCall(VecScale(y, -1.0)); 402*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown 405c4762a1bSJed Brown /* 406c4762a1bSJed Brown RHSMatrixLaplacian - User-provided routine to compute the right-hand-side 407c4762a1bSJed Brown matrix for the Laplacian operator 408c4762a1bSJed Brown 409c4762a1bSJed Brown Input Parameters: 410c4762a1bSJed Brown ts - the TS context 411c4762a1bSJed Brown t - current time (ignored) 412c4762a1bSJed Brown X - current solution (ignored) 413c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 414c4762a1bSJed Brown 415c4762a1bSJed Brown Output Parameters: 416c4762a1bSJed Brown AA - Jacobian matrix 417c4762a1bSJed Brown BB - optionally different matrix from which the preconditioner is built 418c4762a1bSJed Brown str - flag indicating matrix structure 419c4762a1bSJed Brown 420c4762a1bSJed Brown */ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixLaplaciangllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx) 422d71ae5a4SJacob Faibussowitsch { 423c4762a1bSJed Brown PetscReal **temp; 424c4762a1bSJed Brown PetscReal vv; 425c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 426c4762a1bSJed Brown PetscInt i, xs, xn, l, j; 427c4762a1bSJed Brown PetscInt *rowsDM; 428c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 429c4762a1bSJed Brown 430*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-gll_mf", &flg, NULL)); 432c4762a1bSJed Brown 433c4762a1bSJed Brown if (!flg) { 434c4762a1bSJed Brown /* 435c4762a1bSJed Brown Creates the element stiffness matrix for the given gll 436c4762a1bSJed Brown */ 4379566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 438a5b23f4aSJose E. Roman /* workaround for clang analyzer warning: Division by zero */ 4393c633725SBarry Smith PetscCheck(appctx->param.N > 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spectral element order should be > 1"); 440c4762a1bSJed Brown 441c4762a1bSJed Brown /* scale by the size of the element */ 442c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) { 443c4762a1bSJed Brown vv = -appctx->param.mu * 2.0 / appctx->param.Le; 444c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv; 445c4762a1bSJed Brown } 446c4762a1bSJed Brown 4479566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 4489566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 449c4762a1bSJed Brown 450c4762a1bSJed Brown xs = xs / (appctx->param.N - 1); 451c4762a1bSJed Brown xn = xn / (appctx->param.N - 1); 452c4762a1bSJed Brown 4539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM)); 454c4762a1bSJed Brown /* 455c4762a1bSJed Brown loop over local elements 456c4762a1bSJed Brown */ 457c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) { 458ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l; 4599566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES)); 460c4762a1bSJed Brown } 4619566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM)); 4629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 4639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 4649566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 4659566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0)); 4669566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 467c4762a1bSJed Brown 4689566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 469c4762a1bSJed Brown } else { 4709566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL)); 4719566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 4729566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, appctx)); 4739566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Laplacian)); 474c4762a1bSJed Brown } 475*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 476c4762a1bSJed Brown } 477c4762a1bSJed Brown 478c4762a1bSJed Brown /* 479c4762a1bSJed Brown RHSMatrixAdvection - User-provided routine to compute the right-hand-side 480c4762a1bSJed Brown matrix for the Advection (gradient) operator. 481c4762a1bSJed Brown 482c4762a1bSJed Brown Input Parameters: 483c4762a1bSJed Brown ts - the TS context 484c4762a1bSJed Brown t - current time 485c4762a1bSJed Brown global_in - global input vector 486c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 487c4762a1bSJed Brown 488c4762a1bSJed Brown Output Parameters: 489c4762a1bSJed Brown AA - Jacobian matrix 490c4762a1bSJed Brown BB - optionally different preconditioning matrix 491c4762a1bSJed Brown str - flag indicating matrix structure 492c4762a1bSJed Brown 493c4762a1bSJed Brown */ 494d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixAdvectiongllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx) 495d71ae5a4SJacob Faibussowitsch { 496c4762a1bSJed Brown PetscReal **temp; 497c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 498c4762a1bSJed Brown PetscInt xs, xn, l, j; 499c4762a1bSJed Brown PetscInt *rowsDM; 500c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 501c4762a1bSJed Brown 502*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 5039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-gll_mf", &flg, NULL)); 504c4762a1bSJed Brown 505c4762a1bSJed Brown if (!flg) { 506c4762a1bSJed Brown /* 507c4762a1bSJed Brown Creates the advection matrix for the given gll 508c4762a1bSJed Brown */ 5099566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 5109566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 5119566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 512c4762a1bSJed Brown xs = xs / (appctx->param.N - 1); 513c4762a1bSJed Brown xn = xn / (appctx->param.N - 1); 514c4762a1bSJed Brown 5159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM)); 516c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) { 517ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l; 5189566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES)); 519c4762a1bSJed Brown } 5209566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM)); 5219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 523c4762a1bSJed Brown 5249566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 5259566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0)); 5269566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 527c4762a1bSJed Brown 5289566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 529c4762a1bSJed Brown } else { 5309566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL)); 5319566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 5329566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, appctx)); 5339566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Advection)); 534c4762a1bSJed Brown } 535*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown 538c4762a1bSJed Brown /*TEST 539c4762a1bSJed Brown 540c4762a1bSJed Brown build: 541c4762a1bSJed Brown requires: !complex 542c4762a1bSJed Brown 543c4762a1bSJed Brown test: 544c4762a1bSJed Brown suffix: 1 545c4762a1bSJed Brown requires: !single 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: 2 549c4762a1bSJed Brown nsize: 5 550c4762a1bSJed Brown requires: !single 551c4762a1bSJed Brown 552c4762a1bSJed Brown test: 553c4762a1bSJed Brown suffix: 3 554c4762a1bSJed Brown requires: !single 555c4762a1bSJed Brown args: -ts_view -ts_type beuler -gll_mf -pc_type none -ts_max_steps 5 -ts_monitor_error 556c4762a1bSJed Brown 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 4 559c4762a1bSJed Brown requires: !single 560c4762a1bSJed Brown args: -ts_view -ts_type beuler -pc_type none -ts_max_steps 5 -ts_monitor_error 561c4762a1bSJed Brown 562c4762a1bSJed Brown TEST*/ 563