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 */ 59c4762a1bSJed Brown Mat stiff; /* stifness 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 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; 95c4762a1bSJed Brown 96327415f7SBarry Smith PetscFunctionBeginUser; 979566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /*initialize parameters */ 100c4762a1bSJed Brown appctx.param.N = 10; /* order of the spectral element */ 101c4762a1bSJed Brown appctx.param.E = 10; /* number of elements */ 102c4762a1bSJed Brown appctx.param.L = 4.0; /* length of the domain */ 103c4762a1bSJed Brown appctx.param.mu = 0.01; /* diffusion coefficient */ 104c4762a1bSJed Brown appctx.initial_dt = 5e-3; 105c4762a1bSJed Brown appctx.param.steps = PETSC_MAX_INT; 106c4762a1bSJed Brown appctx.param.Tend = 4; 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &appctx.param.N, NULL)); 1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-E", &appctx.param.E, NULL)); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-Tend", &appctx.param.Tend, NULL)); 1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &appctx.param.mu, NULL)); 112c4762a1bSJed Brown appctx.param.Le = appctx.param.L / appctx.param.E; 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1153c633725SBarry Smith PetscCheck((appctx.param.E % size) == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Number of elements must be divisible by number of processes"); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118c4762a1bSJed Brown Create GLL data structures 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1209566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(appctx.param.N, &appctx.SEMop.gll.nodes, appctx.param.N, &appctx.SEMop.gll.weights)); 1219566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights)); 122c4762a1bSJed Brown appctx.SEMop.gll.n = appctx.param.N; 123c4762a1bSJed Brown lenglob = appctx.param.E * (appctx.param.N - 1); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 127c4762a1bSJed Brown and to set up the ghost point communication pattern. There are E*(Nl-1)+1 128c4762a1bSJed Brown total grid values spread equally among all the processors, except first and last 129c4762a1bSJed Brown */ 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, lenglob, 1, 1, NULL, &appctx.da)); 1329566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(appctx.da)); 1339566063dSJacob Faibussowitsch PetscCall(DMSetUp(appctx.da)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* 136c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 137c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 138c4762a1bSJed Brown have the same types. 139c4762a1bSJed Brown */ 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(appctx.da, &appctx.dat.curr_sol)); 1429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.dat.curr_sol, &appctx.SEMop.grid)); 1439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.dat.curr_sol, &appctx.SEMop.mass)); 144c4762a1bSJed Brown 1459566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx.da, &xs, NULL, NULL, &xm, NULL, NULL)); 1469566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1)); 1479566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 150c4762a1bSJed Brown 151c4762a1bSJed Brown xs = xs / (appctx.param.N - 1); 152c4762a1bSJed Brown xm = xm / (appctx.param.N - 1); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* 155c4762a1bSJed Brown Build total grid and mass over entire mesh (multi-elemental) 156c4762a1bSJed Brown */ 157c4762a1bSJed Brown 158c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 159c4762a1bSJed Brown for (j = 0; j < appctx.param.N - 1; j++) { 160c4762a1bSJed Brown x = (appctx.param.Le / 2.0) * (appctx.SEMop.gll.nodes[j] + 1.0) + appctx.param.Le * i; 161c4762a1bSJed Brown ind = i * (appctx.param.N - 1) + j; 162c4762a1bSJed Brown wrk_ptr1[ind] = x; 163c4762a1bSJed Brown wrk_ptr2[ind] = .5 * appctx.param.Le * appctx.SEMop.gll.weights[j]; 164c4762a1bSJed Brown if (j == 0) wrk_ptr2[ind] += .5 * appctx.param.Le * appctx.SEMop.gll.weights[j]; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 1679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.grid, &wrk_ptr1)); 1689566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx.da, appctx.SEMop.mass, &wrk_ptr2)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 171c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1739566063dSJacob Faibussowitsch PetscCall(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE)); 1749566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.stiff)); 1759566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(appctx.da, &appctx.SEMop.grad)); 176c4762a1bSJed Brown /* 177c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 178c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 179c4762a1bSJed Brown as a time-dependent matrix. 180c4762a1bSJed Brown */ 1819566063dSJacob Faibussowitsch PetscCall(RHSMatrixLaplaciangllDM(appctx.ts, 0.0, appctx.dat.curr_sol, appctx.SEMop.stiff, appctx.SEMop.stiff, &appctx)); 1829566063dSJacob Faibussowitsch PetscCall(RHSMatrixAdvectiongllDM(appctx.ts, 0.0, appctx.dat.curr_sol, appctx.SEMop.grad, appctx.SEMop.grad, &appctx)); 183c4762a1bSJed Brown /* 184c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 185c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 186c4762a1bSJed Brown as a time-dependent matrix. 187c4762a1bSJed Brown */ 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(MatDuplicate(appctx.SEMop.stiff, MAT_COPY_VALUES, &appctx.SEMop.keptstiff)); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 1929566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp)); 1939566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.stiff, nsp)); 1949566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.keptstiff, nsp)); 1959566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.stiff, NULL)); 1969566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp)); 197c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 1989566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nsp)); 1999566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(appctx.SEMop.grad, nsp)); 2009566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nsp, appctx.SEMop.grad, NULL)); 2019566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* Create the TS solver that solves the ODE and its adjoint; set its options */ 2049566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts)); 2059566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(appctx.ts, TS_NONLINEAR)); 2069566063dSJacob Faibussowitsch PetscCall(TSSetType(appctx.ts, TSRK)); 2079566063dSJacob Faibussowitsch PetscCall(TSSetDM(appctx.ts, appctx.da)); 2089566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts, 0.0)); 2099566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx.ts, appctx.initial_dt)); 2109566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(appctx.ts, appctx.param.steps)); 2119566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(appctx.ts, appctx.param.Tend)); 2129566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP)); 2139566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(appctx.ts, 1e-7, NULL, 1e-7, NULL)); 2149566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(appctx.ts)); 2159566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts)); 2169566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(appctx.ts, NULL, RHSFunction, &appctx)); 2179566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts, appctx.SEMop.stiff, appctx.SEMop.stiff, RHSJacobian, &appctx)); 218c4762a1bSJed Brown 219c4762a1bSJed Brown /* Set Initial conditions for the problem */ 2209566063dSJacob Faibussowitsch PetscCall(TrueSolution(appctx.ts, 0, appctx.dat.curr_sol, &appctx)); 221c4762a1bSJed Brown 2229566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(appctx.ts, (PetscErrorCode(*)(TS, PetscReal, Vec, void *))TrueSolution, &appctx)); 2239566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx.ts, 0.0)); 2249566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(appctx.ts, 0)); 225c4762a1bSJed Brown 2269566063dSJacob Faibussowitsch PetscCall(TSSolve(appctx.ts, appctx.dat.curr_sol)); 227c4762a1bSJed Brown 2289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.stiff)); 2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.keptstiff)); 2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.SEMop.grad)); 2319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.grid)); 2329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.SEMop.mass)); 2339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.dat.curr_sol)); 2349566063dSJacob Faibussowitsch PetscCall(PetscFree2(appctx.SEMop.gll.nodes, appctx.SEMop.gll.weights)); 2359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&appctx.da)); 2369566063dSJacob Faibussowitsch PetscCall(TSDestroy(&appctx.ts)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 240c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 241c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 242d75802c7SJacob Faibussowitsch options are chosen (e.g., -log_view). 243c4762a1bSJed Brown */ 2449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 245b122ec5aSJacob Faibussowitsch return 0; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* 249c4762a1bSJed Brown TrueSolution() computes the true solution for the PDE 250c4762a1bSJed Brown 251c4762a1bSJed Brown Input Parameter: 252c4762a1bSJed Brown u - uninitialized solution vector (global) 253c4762a1bSJed Brown appctx - user-defined application context 254c4762a1bSJed Brown 255c4762a1bSJed Brown Output Parameter: 256c4762a1bSJed Brown u - vector with solution at initial time (global) 257c4762a1bSJed Brown */ 258d71ae5a4SJacob Faibussowitsch PetscErrorCode TrueSolution(TS ts, PetscReal t, Vec u, AppCtx *appctx) 259d71ae5a4SJacob Faibussowitsch { 260c4762a1bSJed Brown PetscScalar *s; 261c4762a1bSJed Brown const PetscScalar *xg; 262c4762a1bSJed Brown PetscInt i, xs, xn; 263c4762a1bSJed Brown 2643ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2659566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, u, &s)); 2669566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 2679566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 268c4762a1bSJed Brown for (i = xs; i < xs + xn; i++) { 269c4762a1bSJed 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)); 270c4762a1bSJed Brown } 2719566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, u, &s)); 2729566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, appctx->SEMop.grid, (void *)&xg)); 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec globalin, Vec globalout, void *ctx) 277d71ae5a4SJacob Faibussowitsch { 278c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 279c4762a1bSJed Brown 2807510d9b0SBarry Smith PetscFunctionBeginUser; 2819566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.grad, globalin, globalout)); /* grad u */ 2829566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(globalout, globalin, globalout)); /* u grad u */ 2839566063dSJacob Faibussowitsch PetscCall(VecScale(globalout, -1.0)); 2849566063dSJacob Faibussowitsch PetscCall(MatMultAdd(appctx->SEMop.keptstiff, globalin, globalout, globalout)); 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* 289c4762a1bSJed Brown 290c4762a1bSJed Brown K is the discretiziation of the Laplacian 291c4762a1bSJed Brown G is the discretization of the gradient 292c4762a1bSJed Brown 293c4762a1bSJed Brown Computes Jacobian of K u + diag(u) G u which is given by 294c4762a1bSJed Brown K + diag(u)G + diag(Gu) 295c4762a1bSJed Brown */ 296d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec globalin, Mat A, Mat B, void *ctx) 297d71ae5a4SJacob Faibussowitsch { 298c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 299c4762a1bSJed Brown Vec Gglobalin; 300c4762a1bSJed Brown 3017510d9b0SBarry Smith PetscFunctionBeginUser; 302c4762a1bSJed Brown /* A = diag(u) G */ 303c4762a1bSJed Brown 3049566063dSJacob Faibussowitsch PetscCall(MatCopy(appctx->SEMop.grad, A, SAME_NONZERO_PATTERN)); 3059566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, globalin, NULL)); 306c4762a1bSJed Brown 307c4762a1bSJed Brown /* A = A + diag(Gu) */ 3089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(globalin, &Gglobalin)); 3099566063dSJacob Faibussowitsch PetscCall(MatMult(appctx->SEMop.grad, globalin, Gglobalin)); 3109566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, Gglobalin, ADD_VALUES)); 3119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Gglobalin)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* A = K - A */ 3149566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); 3159566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, 0.0, appctx->SEMop.keptstiff, SAME_NONZERO_PATTERN)); 3163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 317c4762a1bSJed Brown } 318c4762a1bSJed Brown 319*46233b44SBarry Smith #include <petscblaslapack.h> 320c4762a1bSJed Brown /* 321c4762a1bSJed Brown Matrix free operation of 1d Laplacian and Grad for GLL spectral elements 322c4762a1bSJed Brown */ 323d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Laplacian(Mat A, Vec x, Vec y) 324d71ae5a4SJacob Faibussowitsch { 325c4762a1bSJed Brown AppCtx *appctx; 326c4762a1bSJed Brown PetscReal **temp, vv; 327c4762a1bSJed Brown PetscInt i, j, xs, xn; 328c4762a1bSJed Brown Vec xlocal, ylocal; 329c4762a1bSJed Brown const PetscScalar *xl; 330c4762a1bSJed Brown PetscScalar *yl; 331c4762a1bSJed Brown PetscBLASInt _One = 1, n; 332c4762a1bSJed Brown PetscScalar _DOne = 1; 333c4762a1bSJed Brown 3343ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 3359566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &appctx)); 3369566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &xlocal)); 3379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(appctx->da, x, INSERT_VALUES, xlocal)); 3389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(appctx->da, x, INSERT_VALUES, xlocal)); 3399566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &ylocal)); 3409566063dSJacob Faibussowitsch PetscCall(VecSet(ylocal, 0.0)); 3419566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 342c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) { 343c4762a1bSJed Brown vv = -appctx->param.mu * 2.0 / appctx->param.Le; 344c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv; 345c4762a1bSJed Brown } 3469566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, xlocal, (void *)&xl)); 3479566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, ylocal, &yl)); 3489566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 3499566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(appctx->param.N, &n)); 35048a46eb9SPierre 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)); 3519566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, xlocal, (void *)&xl)); 3529566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, ylocal, &yl)); 3539566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 3549566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 3559566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(appctx->da, ylocal, ADD_VALUES, y)); 3569566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(appctx->da, ylocal, ADD_VALUES, y)); 3579566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &xlocal)); 3589566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &ylocal)); 3599566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(y, y, appctx->SEMop.mass)); 3603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_Advection(Mat A, Vec x, Vec y) 364d71ae5a4SJacob Faibussowitsch { 365c4762a1bSJed Brown AppCtx *appctx; 366c4762a1bSJed Brown PetscReal **temp; 367c4762a1bSJed Brown PetscInt j, xs, xn; 368c4762a1bSJed Brown Vec xlocal, ylocal; 369c4762a1bSJed Brown const PetscScalar *xl; 370c4762a1bSJed Brown PetscScalar *yl; 371c4762a1bSJed Brown PetscBLASInt _One = 1, n; 372c4762a1bSJed Brown PetscScalar _DOne = 1; 373c4762a1bSJed Brown 3743ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 3759566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &appctx)); 3769566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &xlocal)); 3779566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(appctx->da, x, INSERT_VALUES, xlocal)); 3789566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(appctx->da, x, INSERT_VALUES, xlocal)); 3799566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(appctx->da, &ylocal)); 3809566063dSJacob Faibussowitsch PetscCall(VecSet(ylocal, 0.0)); 3819566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 3829566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(appctx->da, xlocal, (void *)&xl)); 3839566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(appctx->da, ylocal, &yl)); 3849566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 3859566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(appctx->param.N, &n)); 38648a46eb9SPierre 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)); 3879566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(appctx->da, xlocal, (void *)&xl)); 3889566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(appctx->da, ylocal, &yl)); 3899566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 3909566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 3919566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(appctx->da, ylocal, ADD_VALUES, y)); 3929566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(appctx->da, ylocal, ADD_VALUES, y)); 3939566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &xlocal)); 3949566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(appctx->da, &ylocal)); 3959566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(y, y, appctx->SEMop.mass)); 3969566063dSJacob Faibussowitsch PetscCall(VecScale(y, -1.0)); 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* 401c4762a1bSJed Brown RHSMatrixLaplacian - User-provided routine to compute the right-hand-side 402c4762a1bSJed Brown matrix for the Laplacian operator 403c4762a1bSJed Brown 404c4762a1bSJed Brown Input Parameters: 405c4762a1bSJed Brown ts - the TS context 406c4762a1bSJed Brown t - current time (ignored) 407c4762a1bSJed Brown X - current solution (ignored) 408c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 409c4762a1bSJed Brown 410c4762a1bSJed Brown Output Parameters: 411c4762a1bSJed Brown AA - Jacobian matrix 412c4762a1bSJed Brown BB - optionally different matrix from which the preconditioner is built 413c4762a1bSJed Brown str - flag indicating matrix structure 414c4762a1bSJed Brown 415c4762a1bSJed Brown */ 416d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixLaplaciangllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx) 417d71ae5a4SJacob Faibussowitsch { 418c4762a1bSJed Brown PetscReal **temp; 419c4762a1bSJed Brown PetscReal vv; 420c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 421c4762a1bSJed Brown PetscInt i, xs, xn, l, j; 422c4762a1bSJed Brown PetscInt *rowsDM; 423c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 424c4762a1bSJed Brown 4253ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-gll_mf", &flg, NULL)); 427c4762a1bSJed Brown 428c4762a1bSJed Brown if (!flg) { 429c4762a1bSJed Brown /* 430c4762a1bSJed Brown Creates the element stiffness matrix for the given gll 431c4762a1bSJed Brown */ 4329566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 433a5b23f4aSJose E. Roman /* workaround for clang analyzer warning: Division by zero */ 4343c633725SBarry Smith PetscCheck(appctx->param.N > 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spectral element order should be > 1"); 435c4762a1bSJed Brown 436c4762a1bSJed Brown /* scale by the size of the element */ 437c4762a1bSJed Brown for (i = 0; i < appctx->param.N; i++) { 438c4762a1bSJed Brown vv = -appctx->param.mu * 2.0 / appctx->param.Le; 439c4762a1bSJed Brown for (j = 0; j < appctx->param.N; j++) temp[i][j] = temp[i][j] * vv; 440c4762a1bSJed Brown } 441c4762a1bSJed Brown 4429566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 4439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 444c4762a1bSJed Brown 445c4762a1bSJed Brown xs = xs / (appctx->param.N - 1); 446c4762a1bSJed Brown xn = xn / (appctx->param.N - 1); 447c4762a1bSJed Brown 4489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM)); 449c4762a1bSJed Brown /* 450c4762a1bSJed Brown loop over local elements 451c4762a1bSJed Brown */ 452c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) { 453ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l; 4549566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES)); 455c4762a1bSJed Brown } 4569566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM)); 4579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 4589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 4599566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 4609566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0)); 4619566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 462c4762a1bSJed Brown 4639566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 464c4762a1bSJed Brown } else { 4659566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL)); 4669566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 4679566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, appctx)); 4689566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Laplacian)); 469c4762a1bSJed Brown } 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 471c4762a1bSJed Brown } 472c4762a1bSJed Brown 473c4762a1bSJed Brown /* 474c4762a1bSJed Brown RHSMatrixAdvection - User-provided routine to compute the right-hand-side 475c4762a1bSJed Brown matrix for the Advection (gradient) operator. 476c4762a1bSJed Brown 477c4762a1bSJed Brown Input Parameters: 478c4762a1bSJed Brown ts - the TS context 479c4762a1bSJed Brown t - current time 480c4762a1bSJed Brown global_in - global input vector 481c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 482c4762a1bSJed Brown 483c4762a1bSJed Brown Output Parameters: 484c4762a1bSJed Brown AA - Jacobian matrix 485c4762a1bSJed Brown BB - optionally different preconditioning matrix 486c4762a1bSJed Brown str - flag indicating matrix structure 487c4762a1bSJed Brown 488c4762a1bSJed Brown */ 489d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSMatrixAdvectiongllDM(TS ts, PetscReal t, Vec X, Mat A, Mat BB, void *ctx) 490d71ae5a4SJacob Faibussowitsch { 491c4762a1bSJed Brown PetscReal **temp; 492c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 493c4762a1bSJed Brown PetscInt xs, xn, l, j; 494c4762a1bSJed Brown PetscInt *rowsDM; 495c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 496c4762a1bSJed Brown 4973ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-gll_mf", &flg, NULL)); 499c4762a1bSJed Brown 500c4762a1bSJed Brown if (!flg) { 501c4762a1bSJed Brown /* 502c4762a1bSJed Brown Creates the advection matrix for the given gll 503c4762a1bSJed Brown */ 5049566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 5059566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 5069566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(appctx->da, &xs, NULL, NULL, &xn, NULL, NULL)); 507c4762a1bSJed Brown xs = xs / (appctx->param.N - 1); 508c4762a1bSJed Brown xn = xn / (appctx->param.N - 1); 509c4762a1bSJed Brown 5109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(appctx->param.N, &rowsDM)); 511c4762a1bSJed Brown for (j = xs; j < xs + xn; j++) { 512ad540459SPierre Jolivet for (l = 0; l < appctx->param.N; l++) rowsDM[l] = 1 + (j - xs) * (appctx->param.N - 1) + l; 5139566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, appctx->param.N, rowsDM, appctx->param.N, rowsDM, &temp[0][0], ADD_VALUES)); 514c4762a1bSJed Brown } 5159566063dSJacob Faibussowitsch PetscCall(PetscFree(rowsDM)); 5169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 518c4762a1bSJed Brown 5199566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 5209566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, appctx->SEMop.mass, 0)); 5219566063dSJacob Faibussowitsch PetscCall(VecReciprocal(appctx->SEMop.mass)); 522c4762a1bSJed Brown 5239566063dSJacob Faibussowitsch PetscCall(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n, appctx->SEMop.gll.nodes, appctx->SEMop.gll.weights, &temp)); 524c4762a1bSJed Brown } else { 5259566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSHELL)); 5269566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 5279566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, appctx)); 5289566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Advection)); 529c4762a1bSJed Brown } 5303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 531c4762a1bSJed Brown } 532c4762a1bSJed Brown 533c4762a1bSJed Brown /*TEST 534c4762a1bSJed Brown 535c4762a1bSJed Brown build: 536c4762a1bSJed Brown requires: !complex 537c4762a1bSJed Brown 538c4762a1bSJed Brown test: 539c4762a1bSJed Brown suffix: 1 540c4762a1bSJed Brown requires: !single 541c4762a1bSJed Brown 542c4762a1bSJed Brown test: 543c4762a1bSJed Brown suffix: 2 544c4762a1bSJed Brown nsize: 5 545c4762a1bSJed Brown requires: !single 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: 3 549c4762a1bSJed Brown requires: !single 550c4762a1bSJed Brown args: -ts_view -ts_type beuler -gll_mf -pc_type none -ts_max_steps 5 -ts_monitor_error 551c4762a1bSJed Brown 552c4762a1bSJed Brown test: 553c4762a1bSJed Brown suffix: 4 554c4762a1bSJed Brown requires: !single 555c4762a1bSJed Brown args: -ts_view -ts_type beuler -pc_type none -ts_max_steps 5 -ts_monitor_error 556c4762a1bSJed Brown 557c4762a1bSJed Brown TEST*/ 558