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