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