1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] ="Solves one dimensional Burger's equation compares with exact solution\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown 6c4762a1bSJed Brown Not yet tested in parallel 7c4762a1bSJed Brown 8c4762a1bSJed Brown */ 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 11c4762a1bSJed Brown Concepts: TS^Burger's equation 12c4762a1bSJed Brown Processors: n 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown /* ------------------------------------------------------------------------ 16c4762a1bSJed Brown 17c4762a1bSJed Brown This program uses the one-dimensional Burger's equation 18c4762a1bSJed Brown u_t = mu*u_xx - u u_x, 19c4762a1bSJed Brown on the domain 0 <= x <= 1, with periodic boundary conditions 20c4762a1bSJed Brown 21c4762a1bSJed Brown The operators are discretized with the spectral element method 22c4762a1bSJed Brown 23c4762a1bSJed Brown See the paper PDE-CONSTRAINED OPTIMIZATION WITH SPECTRAL ELEMENTS USING PETSC AND TAO 24c4762a1bSJed Brown by OANA MARIN, EMIL CONSTANTINESCU, AND BARRY SMITH for details on the exact solution 25c4762a1bSJed Brown used 26c4762a1bSJed Brown 27c4762a1bSJed Brown See src/tao/unconstrained/tutorials/burgers_spectral.c 28c4762a1bSJed Brown 29c4762a1bSJed Brown ------------------------------------------------------------------------- */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown #include <petscts.h> 32c4762a1bSJed Brown #include <petscdt.h> 33c4762a1bSJed Brown #include <petscdraw.h> 34c4762a1bSJed Brown #include <petscdmda.h> 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* 37c4762a1bSJed Brown User-defined application context - contains data needed by the 38c4762a1bSJed Brown application-provided call-back routines. 39c4762a1bSJed Brown */ 40c4762a1bSJed Brown 41c4762a1bSJed Brown typedef struct { 42c4762a1bSJed Brown PetscInt n; /* number of nodes */ 43c4762a1bSJed Brown PetscReal *nodes; /* GLL nodes */ 44c4762a1bSJed Brown PetscReal *weights; /* GLL weights */ 45c4762a1bSJed Brown } PetscGLL; 46c4762a1bSJed Brown 47c4762a1bSJed Brown typedef struct { 48c4762a1bSJed Brown PetscInt N; /* grid points per elements*/ 49c4762a1bSJed Brown PetscInt E; /* number of elements */ 50c4762a1bSJed Brown PetscReal tol_L2,tol_max; /* error norms */ 51c4762a1bSJed Brown PetscInt steps; /* number of timesteps */ 52c4762a1bSJed Brown PetscReal Tend; /* endtime */ 53c4762a1bSJed Brown PetscReal mu; /* viscosity */ 54c4762a1bSJed Brown PetscReal L; /* total length of domain */ 55c4762a1bSJed Brown PetscReal Le; 56c4762a1bSJed Brown PetscReal Tadj; 57c4762a1bSJed Brown } PetscParam; 58c4762a1bSJed Brown 59c4762a1bSJed Brown typedef struct { 60c4762a1bSJed Brown Vec grid; /* total grid */ 61c4762a1bSJed Brown Vec curr_sol; 62c4762a1bSJed Brown } PetscData; 63c4762a1bSJed Brown 64c4762a1bSJed Brown typedef struct { 65c4762a1bSJed Brown Vec grid; /* total grid */ 66c4762a1bSJed Brown Vec mass; /* mass matrix for total integration */ 67c4762a1bSJed Brown Mat stiff; /* stifness matrix */ 68c4762a1bSJed Brown Mat keptstiff; 69c4762a1bSJed Brown Mat grad; 70c4762a1bSJed Brown PetscGLL gll; 71c4762a1bSJed Brown } PetscSEMOperators; 72c4762a1bSJed Brown 73c4762a1bSJed Brown typedef struct { 74c4762a1bSJed Brown DM da; /* distributed array data structure */ 75c4762a1bSJed Brown PetscSEMOperators SEMop; 76c4762a1bSJed Brown PetscParam param; 77c4762a1bSJed Brown PetscData dat; 78c4762a1bSJed Brown TS ts; 79c4762a1bSJed Brown PetscReal initial_dt; 80c4762a1bSJed Brown } AppCtx; 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* 83c4762a1bSJed Brown User-defined routines 84c4762a1bSJed Brown */ 85c4762a1bSJed Brown extern PetscErrorCode RHSMatrixLaplaciangllDM(TS,PetscReal,Vec,Mat,Mat,void*); 86c4762a1bSJed Brown extern PetscErrorCode RHSMatrixAdvectiongllDM(TS,PetscReal,Vec,Mat,Mat,void*); 87c4762a1bSJed Brown extern PetscErrorCode TrueSolution(TS,PetscReal,Vec,AppCtx*); 88c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 89c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 90c4762a1bSJed Brown 91c4762a1bSJed Brown int main(int argc,char **argv) 92c4762a1bSJed Brown { 93c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 94c4762a1bSJed Brown PetscInt i, xs, xm, ind, j, lenglob; 95c4762a1bSJed Brown PetscReal x, *wrk_ptr1, *wrk_ptr2; 96c4762a1bSJed Brown MatNullSpace nsp; 97c4762a1bSJed Brown PetscMPIInt size; 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100c4762a1bSJed Brown Initialize program and set problem parameters 101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 102c4762a1bSJed Brown PetscFunctionBegin; 103c4762a1bSJed Brown 104*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /*initialize parameters */ 107c4762a1bSJed Brown appctx.param.N = 10; /* order of the spectral element */ 108c4762a1bSJed Brown appctx.param.E = 10; /* number of elements */ 109c4762a1bSJed Brown appctx.param.L = 4.0; /* length of the domain */ 110c4762a1bSJed Brown appctx.param.mu = 0.01; /* diffusion coefficient */ 111c4762a1bSJed Brown appctx.initial_dt = 5e-3; 112c4762a1bSJed Brown appctx.param.steps = PETSC_MAX_INT; 113c4762a1bSJed Brown appctx.param.Tend = 4; 114c4762a1bSJed Brown 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&appctx.param.N,NULL)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-E",&appctx.param.E,NULL)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-Tend",&appctx.param.Tend,NULL)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&appctx.param.mu,NULL)); 119c4762a1bSJed Brown appctx.param.Le = appctx.param.L/appctx.param.E; 120c4762a1bSJed Brown 1215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1223c633725SBarry Smith PetscCheck((appctx.param.E % size) == 0,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of elements must be divisible by number of processes"); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125c4762a1bSJed Brown Create GLL data structures 126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(appctx.param.N,&appctx.SEMop.gll.nodes,appctx.param.N,&appctx.SEMop.gll.weights)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTGaussLobattoLegendreQuadrature(appctx.param.N,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights)); 129c4762a1bSJed Brown appctx.SEMop.gll.n = appctx.param.N; 130c4762a1bSJed Brown lenglob = appctx.param.E*(appctx.param.N-1); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* 133c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 134c4762a1bSJed Brown and to set up the ghost point communication pattern. There are E*(Nl-1)+1 135c4762a1bSJed Brown total grid values spread equally among all the processors, except first and last 136c4762a1bSJed Brown */ 137c4762a1bSJed Brown 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,lenglob,1,1,NULL,&appctx.da)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(appctx.da)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(appctx.da)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* 143c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 144c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 145c4762a1bSJed Brown have the same types. 146c4762a1bSJed Brown */ 147c4762a1bSJed Brown 1485f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(appctx.da,&appctx.dat.curr_sol)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(appctx.dat.curr_sol,&appctx.SEMop.grid)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(appctx.dat.curr_sol,&appctx.SEMop.mass)); 151c4762a1bSJed Brown 1525f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx.da,&xs,NULL,NULL,&xm,NULL,NULL)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 157c4762a1bSJed Brown 158c4762a1bSJed Brown xs=xs/(appctx.param.N-1); 159c4762a1bSJed Brown xm=xm/(appctx.param.N-1); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* 162c4762a1bSJed Brown Build total grid and mass over entire mesh (multi-elemental) 163c4762a1bSJed Brown */ 164c4762a1bSJed Brown 165c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 166c4762a1bSJed Brown for (j=0; j<appctx.param.N-1; j++) { 167c4762a1bSJed Brown x = (appctx.param.Le/2.0)*(appctx.SEMop.gll.nodes[j]+1.0)+appctx.param.Le*i; 168c4762a1bSJed Brown ind=i*(appctx.param.N-1)+j; 169c4762a1bSJed Brown wrk_ptr1[ind]=x; 170c4762a1bSJed Brown wrk_ptr2[ind]=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 171c4762a1bSJed Brown if (j==0) wrk_ptr2[ind]+=.5*appctx.param.Le*appctx.SEMop.gll.weights[j]; 172c4762a1bSJed Brown } 173c4762a1bSJed Brown } 1745f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx.da,appctx.SEMop.grid,&wrk_ptr1)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx.da,appctx.SEMop.mass,&wrk_ptr2)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 178c4762a1bSJed Brown Create matrix data structure; set matrix evaluation routine. 179c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1805f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatrixPreallocateOnly(appctx.da, PETSC_TRUE)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(appctx.da,&appctx.SEMop.stiff)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(appctx.da,&appctx.SEMop.grad)); 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 */ 1885f80ce2aSJacob Faibussowitsch CHKERRQ(RHSMatrixLaplaciangllDM(appctx.ts,0.0,appctx.dat.curr_sol,appctx.SEMop.stiff,appctx.SEMop.stiff,&appctx)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(RHSMatrixAdvectiongllDM(appctx.ts,0.0,appctx.dat.curr_sol,appctx.SEMop.grad,appctx.SEMop.grad,&appctx)); 190c4762a1bSJed Brown /* 191c4762a1bSJed Brown For linear problems with a time-dependent f(u,t) in the equation 192c4762a1bSJed Brown u_t = f(u,t), the user provides the discretized right-hand-side 193c4762a1bSJed Brown as a time-dependent matrix. 194c4762a1bSJed Brown */ 195c4762a1bSJed Brown 1965f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(appctx.SEMop.stiff,MAT_COPY_VALUES,&appctx.SEMop.keptstiff)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 1995f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetNullSpace(appctx.SEMop.stiff,nsp)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetNullSpace(appctx.SEMop.keptstiff,nsp)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceTest(nsp,appctx.SEMop.stiff,NULL)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&nsp)); 204c4762a1bSJed Brown /* attach the null space to the matrix, this probably is not needed but does no harm */ 2055f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&nsp)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetNullSpace(appctx.SEMop.grad,nsp)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceTest(nsp,appctx.SEMop.grad,NULL)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&nsp)); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* Create the TS solver that solves the ODE and its adjoint; set its options */ 2115f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&appctx.ts)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(appctx.ts,TS_NONLINEAR)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(appctx.ts,TSRK)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(appctx.ts,appctx.da)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(appctx.ts,0.0)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(appctx.ts,appctx.initial_dt)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(appctx.ts,appctx.param.steps)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(appctx.ts,appctx.param.Tend)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTolerances(appctx.ts,1e-7,NULL,1e-7,NULL)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSaveTrajectory(appctx.ts)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(appctx.ts)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(appctx.ts,appctx.SEMop.stiff,appctx.SEMop.stiff,RHSJacobian,&appctx)); 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* Set Initial conditions for the problem */ 2275f80ce2aSJacob Faibussowitsch CHKERRQ(TrueSolution(appctx.ts,0,appctx.dat.curr_sol,&appctx)); 228c4762a1bSJed Brown 2295f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolutionFunction(appctx.ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void *))TrueSolution,&appctx)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(appctx.ts,0.0)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetStepNumber(appctx.ts,0)); 232c4762a1bSJed Brown 2335f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(appctx.ts,appctx.dat.curr_sol)); 234c4762a1bSJed Brown 2355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.SEMop.stiff)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.SEMop.keptstiff)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.SEMop.grad)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.SEMop.grid)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.SEMop.mass)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.dat.curr_sol)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(appctx.SEMop.gll.nodes,appctx.SEMop.gll.weights)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&appctx.da)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&appctx.ts)); 244c4762a1bSJed Brown 245c4762a1bSJed Brown /* 246c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 247c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 248c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 249c4762a1bSJed Brown options are chosen (e.g., -log_summary). 250c4762a1bSJed Brown */ 251*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 252*b122ec5aSJacob Faibussowitsch return 0; 253c4762a1bSJed Brown } 254c4762a1bSJed Brown 255c4762a1bSJed Brown /* 256c4762a1bSJed Brown TrueSolution() computes the true solution for the PDE 257c4762a1bSJed Brown 258c4762a1bSJed Brown Input Parameter: 259c4762a1bSJed Brown u - uninitialized solution vector (global) 260c4762a1bSJed Brown appctx - user-defined application context 261c4762a1bSJed Brown 262c4762a1bSJed Brown Output Parameter: 263c4762a1bSJed Brown u - vector with solution at initial time (global) 264c4762a1bSJed Brown */ 265c4762a1bSJed Brown PetscErrorCode TrueSolution(TS ts, PetscReal t, Vec u,AppCtx *appctx) 266c4762a1bSJed Brown { 267c4762a1bSJed Brown PetscScalar *s; 268c4762a1bSJed Brown const PetscScalar *xg; 269c4762a1bSJed Brown PetscInt i,xs,xn; 270c4762a1bSJed Brown 2715f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx->da,u,&s)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 274c4762a1bSJed Brown for (i=xs; i<xs+xn; i++) { 275c4762a1bSJed 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)); 276c4762a1bSJed Brown } 2775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx->da,u,&s)); 2785f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,appctx->SEMop.grid,(void*)&xg)); 279c4762a1bSJed Brown return 0; 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 283c4762a1bSJed Brown { 284c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 285c4762a1bSJed Brown 286c4762a1bSJed Brown PetscFunctionBegin; 2875f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(appctx->SEMop.grad,globalin,globalout)); /* grad u */ 2885f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(globalout,globalin,globalout)); /* u grad u */ 2895f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(globalout, -1.0)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(appctx->SEMop.keptstiff,globalin,globalout,globalout)); 291c4762a1bSJed Brown PetscFunctionReturn(0); 292c4762a1bSJed Brown } 293c4762a1bSJed Brown 294c4762a1bSJed Brown /* 295c4762a1bSJed Brown 296c4762a1bSJed Brown K is the discretiziation of the Laplacian 297c4762a1bSJed Brown G is the discretization of the gradient 298c4762a1bSJed Brown 299c4762a1bSJed Brown Computes Jacobian of K u + diag(u) G u which is given by 300c4762a1bSJed Brown K + diag(u)G + diag(Gu) 301c4762a1bSJed Brown */ 302c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec globalin,Mat A, Mat B,void *ctx) 303c4762a1bSJed Brown { 304c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 305c4762a1bSJed Brown Vec Gglobalin; 306c4762a1bSJed Brown 307c4762a1bSJed Brown PetscFunctionBegin; 308c4762a1bSJed Brown /* A = diag(u) G */ 309c4762a1bSJed Brown 3105f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(appctx->SEMop.grad,A,SAME_NONZERO_PATTERN)); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,globalin,NULL)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* A = A + diag(Gu) */ 3145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(globalin,&Gglobalin)); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(appctx->SEMop.grad,globalin,Gglobalin)); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalSet(A,Gglobalin,ADD_VALUES)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Gglobalin)); 318c4762a1bSJed Brown 319c4762a1bSJed Brown /* A = K - A */ 3205f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,-1.0)); 3215f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A,0.0,appctx->SEMop.keptstiff,SAME_NONZERO_PATTERN)); 322c4762a1bSJed Brown PetscFunctionReturn(0); 323c4762a1bSJed Brown } 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 326c4762a1bSJed Brown 327c4762a1bSJed Brown #include "petscblaslapack.h" 328c4762a1bSJed Brown /* 329c4762a1bSJed Brown Matrix free operation of 1d Laplacian and Grad for GLL spectral elements 330c4762a1bSJed Brown */ 331c4762a1bSJed Brown PetscErrorCode MatMult_Laplacian(Mat A,Vec x,Vec y) 332c4762a1bSJed Brown { 333c4762a1bSJed Brown AppCtx *appctx; 334c4762a1bSJed Brown PetscReal **temp,vv; 335c4762a1bSJed Brown PetscInt i,j,xs,xn; 336c4762a1bSJed Brown Vec xlocal,ylocal; 337c4762a1bSJed Brown const PetscScalar *xl; 338c4762a1bSJed Brown PetscScalar *yl; 339c4762a1bSJed Brown PetscBLASInt _One = 1,n; 340c4762a1bSJed Brown PetscScalar _DOne = 1; 341c4762a1bSJed Brown 3425f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&appctx)); 3435f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(appctx->da,&xlocal)); 3445f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(appctx->da,x,INSERT_VALUES,xlocal)); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(appctx->da,x,INSERT_VALUES,xlocal)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(appctx->da,&ylocal)); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(ylocal,0.0)); 3485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 349c4762a1bSJed Brown for (i=0; i<appctx->param.N; i++) { 350c4762a1bSJed Brown vv =-appctx->param.mu*2.0/appctx->param.Le; 351c4762a1bSJed Brown for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 352c4762a1bSJed Brown } 3535f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(appctx->da,xlocal,(void*)&xl)); 3545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx->da,ylocal,&yl)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBLASIntCast(appctx->param.N,&n)); 357c4762a1bSJed Brown for (j=xs; j<xs+xn; j += appctx->param.N-1) { 358c4762a1bSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&_DOne,&temp[0][0],&n,&xl[j],&_One,&_DOne,&yl[j],&_One)); 359c4762a1bSJed Brown } 3605f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,xlocal,(void*)&xl)); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx->da,ylocal,&yl)); 3625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 3635f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(y,0.0)); 3645f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(appctx->da,ylocal,ADD_VALUES,y)); 3655f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(appctx->da,ylocal,ADD_VALUES,y)); 3665f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(appctx->da,&xlocal)); 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(appctx->da,&ylocal)); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseDivide(y,y,appctx->SEMop.mass)); 369c4762a1bSJed Brown return 0; 370c4762a1bSJed Brown } 371c4762a1bSJed Brown 372c4762a1bSJed Brown PetscErrorCode MatMult_Advection(Mat A,Vec x,Vec y) 373c4762a1bSJed Brown { 374c4762a1bSJed Brown AppCtx *appctx; 375c4762a1bSJed Brown PetscReal **temp; 376c4762a1bSJed Brown PetscInt j,xs,xn; 377c4762a1bSJed Brown Vec xlocal,ylocal; 378c4762a1bSJed Brown const PetscScalar *xl; 379c4762a1bSJed Brown PetscScalar *yl; 380c4762a1bSJed Brown PetscBLASInt _One = 1,n; 381c4762a1bSJed Brown PetscScalar _DOne = 1; 382c4762a1bSJed Brown 3835f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&appctx)); 3845f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(appctx->da,&xlocal)); 3855f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(appctx->da,x,INSERT_VALUES,xlocal)); 3865f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(appctx->da,x,INSERT_VALUES,xlocal)); 3875f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(appctx->da,&ylocal)); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(ylocal,0.0)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(appctx->da,xlocal,(void*)&xl)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(appctx->da,ylocal,&yl)); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 3935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBLASIntCast(appctx->param.N,&n)); 394c4762a1bSJed Brown for (j=xs; j<xs+xn; j += appctx->param.N-1) { 395c4762a1bSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&_DOne,&temp[0][0],&n,&xl[j],&_One,&_DOne,&yl[j],&_One)); 396c4762a1bSJed Brown } 3975f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(appctx->da,xlocal,(void*)&xl)); 3985f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(appctx->da,ylocal,&yl)); 3995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 4005f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(y,0.0)); 4015f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(appctx->da,ylocal,ADD_VALUES,y)); 4025f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(appctx->da,ylocal,ADD_VALUES,y)); 4035f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(appctx->da,&xlocal)); 4045f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(appctx->da,&ylocal)); 4055f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseDivide(y,y,appctx->SEMop.mass)); 4065f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,-1.0)); 407c4762a1bSJed Brown return 0; 408c4762a1bSJed Brown } 409c4762a1bSJed Brown 410c4762a1bSJed Brown /* 411c4762a1bSJed Brown RHSMatrixLaplacian - User-provided routine to compute the right-hand-side 412c4762a1bSJed Brown matrix for the Laplacian operator 413c4762a1bSJed Brown 414c4762a1bSJed Brown Input Parameters: 415c4762a1bSJed Brown ts - the TS context 416c4762a1bSJed Brown t - current time (ignored) 417c4762a1bSJed Brown X - current solution (ignored) 418c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 419c4762a1bSJed Brown 420c4762a1bSJed Brown Output Parameters: 421c4762a1bSJed Brown AA - Jacobian matrix 422c4762a1bSJed Brown BB - optionally different matrix from which the preconditioner is built 423c4762a1bSJed Brown str - flag indicating matrix structure 424c4762a1bSJed Brown 425c4762a1bSJed Brown */ 426c4762a1bSJed Brown PetscErrorCode RHSMatrixLaplaciangllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 427c4762a1bSJed Brown { 428c4762a1bSJed Brown PetscReal **temp; 429c4762a1bSJed Brown PetscReal vv; 430c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 431c4762a1bSJed Brown PetscInt i,xs,xn,l,j; 432c4762a1bSJed Brown PetscInt *rowsDM; 433c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 434c4762a1bSJed Brown 4355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-gll_mf",&flg,NULL)); 436c4762a1bSJed Brown 437c4762a1bSJed Brown if (!flg) { 438c4762a1bSJed Brown /* 439c4762a1bSJed Brown Creates the element stiffness matrix for the given gll 440c4762a1bSJed Brown */ 4415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementLaplacianCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 442a5b23f4aSJose E. Roman /* workaround for clang analyzer warning: Division by zero */ 4433c633725SBarry Smith PetscCheck(appctx->param.N > 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Spectral element order should be > 1"); 444c4762a1bSJed Brown 445c4762a1bSJed Brown /* scale by the size of the element */ 446c4762a1bSJed Brown for (i=0; i<appctx->param.N; i++) { 447c4762a1bSJed Brown vv=-appctx->param.mu*2.0/appctx->param.Le; 448c4762a1bSJed Brown for (j=0; j<appctx->param.N; j++) temp[i][j]=temp[i][j]*vv; 449c4762a1bSJed Brown } 450c4762a1bSJed Brown 4515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 4525f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 453c4762a1bSJed Brown 454c4762a1bSJed Brown xs = xs/(appctx->param.N-1); 455c4762a1bSJed Brown xn = xn/(appctx->param.N-1); 456c4762a1bSJed Brown 4575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(appctx->param.N,&rowsDM)); 458c4762a1bSJed Brown /* 459c4762a1bSJed Brown loop over local elements 460c4762a1bSJed Brown */ 461c4762a1bSJed Brown for (j=xs; j<xs+xn; j++) { 462c4762a1bSJed Brown for (l=0; l<appctx->param.N; l++) { 463c4762a1bSJed Brown rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 464c4762a1bSJed Brown } 4655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES)); 466c4762a1bSJed Brown } 4675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rowsDM)); 4685f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 4695f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 4705f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 4715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,appctx->SEMop.mass,0)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 473c4762a1bSJed Brown 4745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementLaplacianDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 475c4762a1bSJed Brown } else { 4765f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSHELL)); 4775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 4785f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(A,appctx)); 4795f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_Laplacian)); 480c4762a1bSJed Brown } 481c4762a1bSJed Brown return 0; 482c4762a1bSJed Brown } 483c4762a1bSJed Brown 484c4762a1bSJed Brown /* 485c4762a1bSJed Brown RHSMatrixAdvection - User-provided routine to compute the right-hand-side 486c4762a1bSJed Brown matrix for the Advection (gradient) operator. 487c4762a1bSJed Brown 488c4762a1bSJed Brown Input Parameters: 489c4762a1bSJed Brown ts - the TS context 490c4762a1bSJed Brown t - current time 491c4762a1bSJed Brown global_in - global input vector 492c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 493c4762a1bSJed Brown 494c4762a1bSJed Brown Output Parameters: 495c4762a1bSJed Brown AA - Jacobian matrix 496c4762a1bSJed Brown BB - optionally different preconditioning matrix 497c4762a1bSJed Brown str - flag indicating matrix structure 498c4762a1bSJed Brown 499c4762a1bSJed Brown */ 500c4762a1bSJed Brown PetscErrorCode RHSMatrixAdvectiongllDM(TS ts,PetscReal t,Vec X,Mat A,Mat BB,void *ctx) 501c4762a1bSJed Brown { 502c4762a1bSJed Brown PetscReal **temp; 503c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 504c4762a1bSJed Brown PetscInt xs,xn,l,j; 505c4762a1bSJed Brown PetscInt *rowsDM; 506c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 507c4762a1bSJed Brown 5085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-gll_mf",&flg,NULL)); 509c4762a1bSJed Brown 510c4762a1bSJed Brown if (!flg) { 511c4762a1bSJed Brown /* 512c4762a1bSJed Brown Creates the advection matrix for the given gll 513c4762a1bSJed Brown */ 5145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementAdvectionCreate(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 5155f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 5165f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(appctx->da,&xs,NULL,NULL,&xn,NULL,NULL)); 517c4762a1bSJed Brown xs = xs/(appctx->param.N-1); 518c4762a1bSJed Brown xn = xn/(appctx->param.N-1); 519c4762a1bSJed Brown 5205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(appctx->param.N,&rowsDM)); 521c4762a1bSJed Brown for (j=xs; j<xs+xn; j++) { 522c4762a1bSJed Brown for (l=0; l<appctx->param.N; l++) { 523c4762a1bSJed Brown rowsDM[l] = 1+(j-xs)*(appctx->param.N-1)+l; 524c4762a1bSJed Brown } 5255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesLocal(A,appctx->param.N,rowsDM,appctx->param.N,rowsDM,&temp[0][0],ADD_VALUES)); 526c4762a1bSJed Brown } 5275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rowsDM)); 5285f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5295f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 530c4762a1bSJed Brown 5315f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 5325f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,appctx->SEMop.mass,0)); 5335f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(appctx->SEMop.mass)); 534c4762a1bSJed Brown 5355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscGaussLobattoLegendreElementAdvectionDestroy(appctx->SEMop.gll.n,appctx->SEMop.gll.nodes,appctx->SEMop.gll.weights,&temp)); 536c4762a1bSJed Brown } else { 5375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSHELL)); 5385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 5395f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(A,appctx)); 5405f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_Advection)); 541c4762a1bSJed Brown } 542c4762a1bSJed Brown return 0; 543c4762a1bSJed Brown } 544c4762a1bSJed Brown 545c4762a1bSJed Brown /*TEST 546c4762a1bSJed Brown 547c4762a1bSJed Brown build: 548c4762a1bSJed Brown requires: !complex 549c4762a1bSJed Brown 550c4762a1bSJed Brown test: 551c4762a1bSJed Brown suffix: 1 552c4762a1bSJed Brown requires: !single 553c4762a1bSJed Brown 554c4762a1bSJed Brown test: 555c4762a1bSJed Brown suffix: 2 556c4762a1bSJed Brown nsize: 5 557c4762a1bSJed Brown requires: !single 558c4762a1bSJed Brown 559c4762a1bSJed Brown test: 560c4762a1bSJed Brown suffix: 3 561c4762a1bSJed Brown requires: !single 562c4762a1bSJed Brown args: -ts_view -ts_type beuler -gll_mf -pc_type none -ts_max_steps 5 -ts_monitor_error 563c4762a1bSJed Brown 564c4762a1bSJed Brown test: 565c4762a1bSJed Brown suffix: 4 566c4762a1bSJed Brown requires: !single 567c4762a1bSJed Brown args: -ts_view -ts_type beuler -pc_type none -ts_max_steps 5 -ts_monitor_error 568c4762a1bSJed Brown 569c4762a1bSJed Brown TEST*/ 570