1*c4762a1bSJed Brown /* 2*c4762a1bSJed Brown Formatted test for TS routines. 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown Solves U_t=F(t,u) 5*c4762a1bSJed Brown Where: 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown [2*u1+u2 8*c4762a1bSJed Brown F(t,u)= [u1+2*u2+u3 9*c4762a1bSJed Brown [ u2+2*u3 10*c4762a1bSJed Brown We can compare the solutions from euler, beuler and SUNDIALS to 11*c4762a1bSJed Brown see what is the difference. 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown */ 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown static char help[] = "Solves a linear ODE. \n\n"; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown #include <petscts.h> 18*c4762a1bSJed Brown #include <petscpc.h> 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 21*c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 22*c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 23*c4762a1bSJed Brown extern PetscErrorCode Initial(Vec,void*); 24*c4762a1bSJed Brown extern PetscErrorCode MyMatMult(Mat,Vec,Vec); 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown extern PetscReal solx(PetscReal); 27*c4762a1bSJed Brown extern PetscReal soly(PetscReal); 28*c4762a1bSJed Brown extern PetscReal solz(PetscReal); 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown int main(int argc,char **argv) 31*c4762a1bSJed Brown { 32*c4762a1bSJed Brown PetscErrorCode ierr; 33*c4762a1bSJed Brown PetscInt time_steps = 100,steps; 34*c4762a1bSJed Brown PetscMPIInt size; 35*c4762a1bSJed Brown Vec global; 36*c4762a1bSJed Brown PetscReal dt,ftime; 37*c4762a1bSJed Brown TS ts; 38*c4762a1bSJed Brown Mat A = 0,S; 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 41*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);CHKERRQ(ierr); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* set initial conditions */ 46*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&global);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = VecSetSizes(global,PETSC_DECIDE,3);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = VecSetFromOptions(global);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = Initial(global,NULL);CHKERRQ(ierr); 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown /* make timestep context */ 52*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,NULL,NULL);CHKERRQ(ierr); 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown dt = 0.001; 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown /* 59*c4762a1bSJed Brown The user provides the RHS and Jacobian 60*c4762a1bSJed Brown */ 61*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,NULL);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = RHSJacobian(ts,0.0,global,A,A,NULL);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL);CHKERRQ(ierr); 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = MatCreateShell(PETSC_COMM_WORLD,3,3,3,3,NULL,&S);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MyMatMult);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,S,A,RHSJacobian,NULL);CHKERRQ(ierr); 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,time_steps);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,1);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = TSSetSolution(ts,global);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown ierr = TSSolve(ts,global);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* free the memories */ 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = VecDestroy(&global);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatDestroy(&S);CHKERRQ(ierr); 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown ierr = PetscFinalize(); 94*c4762a1bSJed Brown return ierr; 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat S,Vec x,Vec y) 98*c4762a1bSJed Brown { 99*c4762a1bSJed Brown PetscErrorCode ierr; 100*c4762a1bSJed Brown const PetscScalar *inptr; 101*c4762a1bSJed Brown PetscScalar *outptr; 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown PetscFunctionBegin; 104*c4762a1bSJed Brown ierr = VecGetArrayRead(x,&inptr);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = VecGetArray(y,&outptr);CHKERRQ(ierr); 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown outptr[0] = 2.0*inptr[0]+inptr[1]; 108*c4762a1bSJed Brown outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; 109*c4762a1bSJed Brown outptr[2] = inptr[1]+2.0*inptr[2]; 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&inptr);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = VecRestoreArray(y,&outptr);CHKERRQ(ierr); 113*c4762a1bSJed Brown PetscFunctionReturn(0); 114*c4762a1bSJed Brown } 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown /* -------------------------------------------------------------------*/ 118*c4762a1bSJed Brown /* this test problem has initial values (1,1,1). */ 119*c4762a1bSJed Brown PetscErrorCode Initial(Vec global,void *ctx) 120*c4762a1bSJed Brown { 121*c4762a1bSJed Brown PetscScalar *localptr; 122*c4762a1bSJed Brown PetscInt i,mybase,myend,locsize; 123*c4762a1bSJed Brown PetscErrorCode ierr; 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown /* determine starting point of each processor */ 126*c4762a1bSJed Brown ierr = VecGetOwnershipRange(global,&mybase,&myend);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = VecGetLocalSize(global,&locsize);CHKERRQ(ierr); 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown /* Initialize the array */ 130*c4762a1bSJed Brown ierr = VecGetArray(global,&localptr);CHKERRQ(ierr); 131*c4762a1bSJed Brown for (i=0; i<locsize; i++) localptr[i] = 1.0; 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown if (mybase == 0) localptr[0]=1.0; 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown ierr = VecRestoreArray(global,&localptr);CHKERRQ(ierr); 136*c4762a1bSJed Brown return 0; 137*c4762a1bSJed Brown } 138*c4762a1bSJed Brown 139*c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx) 140*c4762a1bSJed Brown { 141*c4762a1bSJed Brown VecScatter scatter; 142*c4762a1bSJed Brown IS from,to; 143*c4762a1bSJed Brown PetscInt i,n,*idx; 144*c4762a1bSJed Brown Vec tmp_vec; 145*c4762a1bSJed Brown PetscErrorCode ierr; 146*c4762a1bSJed Brown const PetscScalar *tmp; 147*c4762a1bSJed Brown 148*c4762a1bSJed Brown /* Get the size of the vector */ 149*c4762a1bSJed Brown ierr = VecGetSize(global,&n);CHKERRQ(ierr); 150*c4762a1bSJed Brown 151*c4762a1bSJed Brown /* Set the index sets */ 152*c4762a1bSJed Brown ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr); 153*c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 154*c4762a1bSJed Brown 155*c4762a1bSJed Brown /* Create local sequential vectors */ 156*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);CHKERRQ(ierr); 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown /* Create scatter context */ 159*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = VecScatterCreate(global,from,tmp_vec,to,&scatter);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 164*c4762a1bSJed Brown 165*c4762a1bSJed Brown ierr = VecGetArrayRead(tmp_vec,&tmp);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e %14.6e %14.6e \n", 167*c4762a1bSJed Brown (double)time,(double)PetscRealPart(tmp[0]),(double)PetscRealPart(tmp[1]),(double)PetscRealPart(tmp[2]));CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e %14.6e %14.6e \n", 169*c4762a1bSJed Brown (double)time,(double)PetscRealPart(tmp[0]-solx(time)),(double)PetscRealPart(tmp[1]-soly(time)),(double)PetscRealPart(tmp[2]-solz(time)));CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = VecRestoreArrayRead(tmp_vec,&tmp);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = ISDestroy(&from);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = ISDestroy(&to);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = VecDestroy(&tmp_vec);CHKERRQ(ierr); 176*c4762a1bSJed Brown return 0; 177*c4762a1bSJed Brown } 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 180*c4762a1bSJed Brown { 181*c4762a1bSJed Brown PetscScalar *outptr; 182*c4762a1bSJed Brown const PetscScalar *inptr; 183*c4762a1bSJed Brown PetscInt i,n,*idx; 184*c4762a1bSJed Brown PetscErrorCode ierr; 185*c4762a1bSJed Brown IS from,to; 186*c4762a1bSJed Brown VecScatter scatter; 187*c4762a1bSJed Brown Vec tmp_in,tmp_out; 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /* Get the length of parallel vector */ 190*c4762a1bSJed Brown ierr = VecGetSize(globalin,&n);CHKERRQ(ierr); 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown /* Set the index sets */ 193*c4762a1bSJed Brown ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr); 194*c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 195*c4762a1bSJed Brown 196*c4762a1bSJed Brown /* Create local sequential vectors */ 197*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = VecDuplicate(tmp_in,&tmp_out);CHKERRQ(ierr); 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown /* Create scatter context */ 201*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);CHKERRQ(ierr); 203*c4762a1bSJed Brown ierr = VecScatterCreate(globalin,from,tmp_in,to,&scatter);CHKERRQ(ierr); 204*c4762a1bSJed Brown ierr = VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 205*c4762a1bSJed Brown ierr = VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 206*c4762a1bSJed Brown ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr); 207*c4762a1bSJed Brown 208*c4762a1bSJed Brown /*Extract income array */ 209*c4762a1bSJed Brown ierr = VecGetArrayRead(tmp_in,&inptr);CHKERRQ(ierr); 210*c4762a1bSJed Brown 211*c4762a1bSJed Brown /* Extract outcome array*/ 212*c4762a1bSJed Brown ierr = VecGetArray(tmp_out,&outptr);CHKERRQ(ierr); 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown outptr[0] = 2.0*inptr[0]+inptr[1]; 215*c4762a1bSJed Brown outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; 216*c4762a1bSJed Brown outptr[2] = inptr[1]+2.0*inptr[2]; 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown ierr = VecRestoreArrayRead(tmp_in,&inptr);CHKERRQ(ierr); 219*c4762a1bSJed Brown ierr = VecRestoreArray(tmp_out,&outptr);CHKERRQ(ierr); 220*c4762a1bSJed Brown 221*c4762a1bSJed Brown ierr = VecScatterCreate(tmp_out,from,globalout,to,&scatter);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 223*c4762a1bSJed Brown ierr = VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 224*c4762a1bSJed Brown 225*c4762a1bSJed Brown /* Destroy idx aand scatter */ 226*c4762a1bSJed Brown ierr = ISDestroy(&from);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = ISDestroy(&to);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = VecDestroy(&tmp_in);CHKERRQ(ierr); 230*c4762a1bSJed Brown ierr = VecDestroy(&tmp_out);CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 232*c4762a1bSJed Brown return 0; 233*c4762a1bSJed Brown } 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ctx) 236*c4762a1bSJed Brown { 237*c4762a1bSJed Brown PetscScalar v[3]; 238*c4762a1bSJed Brown const PetscScalar *tmp; 239*c4762a1bSJed Brown PetscInt idx[3],i; 240*c4762a1bSJed Brown PetscErrorCode ierr; 241*c4762a1bSJed Brown 242*c4762a1bSJed Brown idx[0]=0; idx[1]=1; idx[2]=2; 243*c4762a1bSJed Brown ierr = VecGetArrayRead(x,&tmp);CHKERRQ(ierr); 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown i = 0; 246*c4762a1bSJed Brown v[0] = 2.0; v[1] = 1.0; v[2] = 0.0; 247*c4762a1bSJed Brown ierr = MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 248*c4762a1bSJed Brown 249*c4762a1bSJed Brown i = 1; 250*c4762a1bSJed Brown v[0] = 1.0; v[1] = 2.0; v[2] = 1.0; 251*c4762a1bSJed Brown ierr = MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 252*c4762a1bSJed Brown 253*c4762a1bSJed Brown i = 2; 254*c4762a1bSJed Brown v[0] = 0.0; v[1] = 1.0; v[2] = 2.0; 255*c4762a1bSJed Brown ierr = MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr); 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown ierr = MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 258*c4762a1bSJed Brown ierr = MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown if (A != BB) { 261*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 263*c4762a1bSJed Brown } 264*c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&tmp);CHKERRQ(ierr); 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown return 0; 267*c4762a1bSJed Brown } 268*c4762a1bSJed Brown 269*c4762a1bSJed Brown /* 270*c4762a1bSJed Brown The exact solutions 271*c4762a1bSJed Brown */ 272*c4762a1bSJed Brown PetscReal solx(PetscReal t) 273*c4762a1bSJed Brown { 274*c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 275*c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 276*c4762a1bSJed Brown } 277*c4762a1bSJed Brown 278*c4762a1bSJed Brown PetscReal soly(PetscReal t) 279*c4762a1bSJed Brown { 280*c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0) + 281*c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0); 282*c4762a1bSJed Brown } 283*c4762a1bSJed Brown 284*c4762a1bSJed Brown PetscReal solz(PetscReal t) 285*c4762a1bSJed Brown { 286*c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 287*c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 288*c4762a1bSJed Brown } 289*c4762a1bSJed Brown 290*c4762a1bSJed Brown 291*c4762a1bSJed Brown /*TEST 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown test: 294*c4762a1bSJed Brown suffix: euler 295*c4762a1bSJed Brown args: -ts_type euler 296*c4762a1bSJed Brown requires: !single 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown test: 299*c4762a1bSJed Brown suffix: beuler 300*c4762a1bSJed Brown args: -ts_type beuler 301*c4762a1bSJed Brown requires: !single 302*c4762a1bSJed Brown 303*c4762a1bSJed Brown TEST*/ 304*c4762a1bSJed Brown 305