1c4762a1bSJed Brown static const char help[] = "Time-dependent advection-reaction PDE in 1d, demonstrates IMEX methods.\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown u_t + a1*u_x = -k1*u + k2*v + s1 4c4762a1bSJed Brown v_t + a2*v_x = k1*u - k2*v + s2 5c4762a1bSJed Brown 0 < x < 1; 6c4762a1bSJed Brown a1 = 1, k1 = 10^6, s1 = 0, 7c4762a1bSJed Brown a2 = 0, k2 = 2*k1, s2 = 1 8c4762a1bSJed Brown 9c4762a1bSJed Brown Initial conditions: 10c4762a1bSJed Brown u(x,0) = 1 + s2*x 11c4762a1bSJed Brown v(x,0) = k0/k1*u(x,0) + s1/k1 12c4762a1bSJed Brown 13c4762a1bSJed Brown Upstream boundary conditions: 14c4762a1bSJed Brown u(0,t) = 1-sin(12*t)^4 15c4762a1bSJed Brown 16c4762a1bSJed Brown Useful command-line parameters: 17c4762a1bSJed Brown -ts_arkimex_fully_implicit - treats advection implicitly, only relevant with -ts_type arkimex (default) 18c4762a1bSJed Brown -ts_type rosw - use Rosenbrock-W method (linearly implicit IMEX, amortizes assembly and preconditioner setup) 19c4762a1bSJed Brown */ 20c4762a1bSJed Brown 21c4762a1bSJed Brown #include <petscdm.h> 22c4762a1bSJed Brown #include <petscdmda.h> 23c4762a1bSJed Brown #include <petscts.h> 24c4762a1bSJed Brown 25c4762a1bSJed Brown typedef PetscScalar Field[2]; 26c4762a1bSJed Brown 27c4762a1bSJed Brown typedef struct _User *User; 28c4762a1bSJed Brown struct _User { 29c4762a1bSJed Brown PetscReal a[2]; /* Advection speeds */ 30c4762a1bSJed Brown PetscReal k[2]; /* Reaction coefficients */ 31c4762a1bSJed Brown PetscReal s[2]; /* Source terms */ 32c4762a1bSJed Brown }; 33c4762a1bSJed Brown 34c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 35c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 36c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 37c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 38c4762a1bSJed Brown 39c4762a1bSJed Brown int main(int argc,char **argv) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown TS ts; /* time integrator */ 42c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 43c4762a1bSJed Brown SNESLineSearch linesearch; /* line search */ 44c4762a1bSJed Brown Vec X; /* solution, residual vectors */ 45c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 46c4762a1bSJed Brown PetscInt steps,mx; 47c4762a1bSJed Brown PetscErrorCode ierr; 48c4762a1bSJed Brown DM da; 49c4762a1bSJed Brown PetscReal ftime,dt; 50c4762a1bSJed Brown struct _User user; /* user-defined work context */ 51c4762a1bSJed Brown TSConvergedReason reason; 52c4762a1bSJed Brown 53c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 57c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 63c4762a1bSJed Brown Extract global vectors from DMDA; 64c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&X)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Initialize user application context */ 681e1ea65dSPierre Jolivet ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");CHKERRQ(ierr); 69c4762a1bSJed Brown { 70*5f80ce2aSJacob Faibussowitsch user.a[0] = 1; CHKERRQ(PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL)); 71*5f80ce2aSJacob Faibussowitsch user.a[1] = 0; CHKERRQ(PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL)); 72*5f80ce2aSJacob Faibussowitsch user.k[0] = 1e6; CHKERRQ(PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL)); 73*5f80ce2aSJacob Faibussowitsch user.k[1] = 2*user.k[0]; CHKERRQ(PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL)); 74*5f80ce2aSJacob Faibussowitsch user.s[0] = 0; CHKERRQ(PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL)); 75*5f80ce2aSJacob Faibussowitsch user.s[1] = 1; CHKERRQ(PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL)); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 80c4762a1bSJed Brown Create timestepping solver context 81c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSARKIMEX)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,FormIFunction,&user)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATAIJ)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&J)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,J,J,FormIJacobian,&user)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since 92c4762a1bSJed Brown * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use 93c4762a1bSJed Brown * SNESSetType(snes,SNESKSPONLY). */ 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSNES(ts,&snes)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetLineSearch(snes,&linesearch)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown ftime = .1; 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,ftime)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103c4762a1bSJed Brown Set initial conditions 104c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(FormInitialSolution(ts,X,&user)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,X)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(X,&mx)); 108c4762a1bSJed Brown dt = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */ 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 112c4762a1bSJed Brown Set runtime options 113c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117c4762a1bSJed Brown Solve nonlinear system 118c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,X)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetConvergedReason(ts,&reason)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 126c4762a1bSJed Brown Free work space. 127c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 132c4762a1bSJed Brown ierr = PetscFinalize(); 133c4762a1bSJed Brown return ierr; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 137c4762a1bSJed Brown { 138c4762a1bSJed Brown User user = (User)ptr; 139c4762a1bSJed Brown DM da; 140c4762a1bSJed Brown DMDALocalInfo info; 141c4762a1bSJed Brown PetscInt i; 142c4762a1bSJed Brown Field *f; 143c4762a1bSJed Brown const Field *x,*xdot; 144c4762a1bSJed Brown 145c4762a1bSJed Brown PetscFunctionBeginUser; 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Get pointers to vector data */ 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X,(void*)&x)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xdot,(void*)&xdot)); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 155c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 156c4762a1bSJed Brown f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0]; 157c4762a1bSJed Brown f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1]; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* Restore vectors */ 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X,(void*)&x)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown User user = (User)ptr; 170c4762a1bSJed Brown DM da; 171c4762a1bSJed Brown Vec Xloc; 172c4762a1bSJed Brown DMDALocalInfo info; 173c4762a1bSJed Brown PetscInt i,j; 174c4762a1bSJed Brown PetscReal hx; 175c4762a1bSJed Brown Field *f; 176c4762a1bSJed Brown const Field *x; 177c4762a1bSJed Brown 178c4762a1bSJed Brown PetscFunctionBeginUser; 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 181c4762a1bSJed Brown hx = 1.0/(PetscReal)info.mx; 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* 184c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 185c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 186c4762a1bSJed Brown By placing code between these two statements, computations can be 187c4762a1bSJed Brown done while messages are in transition. 188c4762a1bSJed Brown */ 189*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&Xloc)); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* Get pointers to vector data */ 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xloc,(void*)&x)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 198c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 199c4762a1bSJed Brown const PetscReal *a = user->a; 200c4762a1bSJed Brown PetscReal u0t[2]; 201c4762a1bSJed Brown u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12*t),4); 202c4762a1bSJed Brown u0t[1] = 0.0; 203c4762a1bSJed Brown for (j=0; j<2; j++) { 204c4762a1bSJed Brown if (i == 0) f[i][j] = a[j]/hx*(1./3*u0t[j] + 0.5*x[i][j] - x[i+1][j] + 1./6*x[i+2][j]); 205c4762a1bSJed Brown else if (i == 1) f[i][j] = a[j]/hx*(-1./12*u0t[j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]); 206c4762a1bSJed Brown else if (i == info.mx-2) f[i][j] = a[j]/hx*(-1./6*x[i-2][j] + x[i-1][j] - 0.5*x[i][j] - 1./3*x[i+1][j]); 207c4762a1bSJed Brown else if (i == info.mx-1) f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]); 208c4762a1bSJed Brown else f[i][j] = a[j]/hx*(-1./12*x[i-2][j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Restore vectors */ 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xloc,(void*)&x)); 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&Xloc)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 220c4762a1bSJed Brown /* 221c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 222c4762a1bSJed Brown */ 223c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 224c4762a1bSJed Brown { 225c4762a1bSJed Brown User user = (User)ptr; 226c4762a1bSJed Brown DMDALocalInfo info; 227c4762a1bSJed Brown PetscInt i; 228c4762a1bSJed Brown DM da; 229c4762a1bSJed Brown const Field *x,*xdot; 230c4762a1bSJed Brown 231c4762a1bSJed Brown PetscFunctionBeginUser; 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 233*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* Get pointers to vector data */ 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X,(void*)&x)); 237*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xdot,(void*)&xdot)); 238c4762a1bSJed Brown 239c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 240c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 241c4762a1bSJed Brown const PetscReal *k = user->k; 242c4762a1bSJed Brown PetscScalar v[2][2]; 243c4762a1bSJed Brown 244c4762a1bSJed Brown v[0][0] = a + k[0]; v[0][1] = -k[1]; 245c4762a1bSJed Brown v[1][0] = -k[0]; v[1][1] = a+k[1]; 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES)); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* Restore vectors */ 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X,(void*)&x)); 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot)); 252c4762a1bSJed Brown 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 254*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 255c4762a1bSJed Brown if (J != Jpre) { 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown PetscFunctionReturn(0); 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 263c4762a1bSJed Brown { 264c4762a1bSJed Brown User user = (User)ctx; 265c4762a1bSJed Brown DM da; 266c4762a1bSJed Brown PetscInt i; 267c4762a1bSJed Brown DMDALocalInfo info; 268c4762a1bSJed Brown Field *x; 269c4762a1bSJed Brown PetscReal hx; 270c4762a1bSJed Brown 271c4762a1bSJed Brown PetscFunctionBeginUser; 272*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 273*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 274c4762a1bSJed Brown hx = 1.0/(PetscReal)info.mx; 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* Get pointers to vector data */ 277*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,X,&x)); 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 280c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 281c4762a1bSJed Brown PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0; 282c4762a1bSJed Brown x[i][0] = 1 + user->s[1]*r; 283c4762a1bSJed Brown x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik; 284c4762a1bSJed Brown } 285*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,X,&x)); 286c4762a1bSJed Brown PetscFunctionReturn(0); 287c4762a1bSJed Brown } 288c4762a1bSJed Brown 289c4762a1bSJed Brown /*TEST 290c4762a1bSJed Brown 291c4762a1bSJed Brown test: 292c4762a1bSJed Brown args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_arkimex_type 4 -ts_adapt_type none -ts_dt .005 -ts_max_time .1 293c4762a1bSJed Brown requires: !single 294c4762a1bSJed Brown 295c4762a1bSJed Brown test: 296c4762a1bSJed Brown suffix: 2 297c4762a1bSJed Brown args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_dt 1e-3 -ts_adapt_type none -ts_dt .005 -ts_max_time .1 298c4762a1bSJed Brown nsize: 2 299c4762a1bSJed Brown 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown suffix: 3 302c4762a1bSJed Brown args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type ra34pw2 -ts_dt 5e-3 -ts_max_time .1 -ts_adapt_type none 303c4762a1bSJed Brown nsize: 2 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown suffix: 4 307c4762a1bSJed Brown args: -ts_type eimex -da_grid_x 200 -ts_eimex_order_adapt true -ts_dt 0.001 -ts_monitor -ts_max_steps 100 308c4762a1bSJed Brown filter: sed "s/ITS/TIME/g" 309c4762a1bSJed Brown nsize: 2 310c4762a1bSJed Brown 311c4762a1bSJed Brown TEST*/ 312