1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] ="Solves a time-dependent nonlinear PDE with lower and upper bounds on the interior grid points. Uses implicit\n\ 3c4762a1bSJed Brown timestepping. Runtime options include:\n\ 4c4762a1bSJed Brown -M <xg>, where <xg> = number of grid points\n\ 5c4762a1bSJed Brown -debug : Activate debugging printouts\n\ 6c4762a1bSJed Brown -nox : Deactivate x-window graphics\n\ 7c4762a1bSJed Brown -ul : lower bound\n\ 8c4762a1bSJed Brown -uh : upper bound\n\n"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 12c4762a1bSJed Brown Concepts: TS^Variational inequality nonlinear solver 13c4762a1bSJed Brown Processors: n 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* ------------------------------------------------------------------------ 17c4762a1bSJed Brown 18c4762a1bSJed Brown This is a variation of ex2.c to solve the PDE 19c4762a1bSJed Brown 20c4762a1bSJed Brown u * u_xx 21c4762a1bSJed Brown u_t = --------- 22c4762a1bSJed Brown 2*(t+1)^2 23c4762a1bSJed Brown 24c4762a1bSJed Brown with box constraints on the interior grid points 25c4762a1bSJed Brown ul <= u(t,x) <= uh with x != 0,1 26c4762a1bSJed Brown on the domain 0 <= x <= 1, with boundary conditions 27c4762a1bSJed Brown u(t,0) = t + 1, u(t,1) = 2*t + 2, 28c4762a1bSJed Brown and initial condition 29c4762a1bSJed Brown u(0,x) = 1 + x*x. 30c4762a1bSJed Brown 31c4762a1bSJed Brown The exact solution is: 32c4762a1bSJed Brown u(t,x) = (1 + x*x) * (1 + t) 33c4762a1bSJed Brown 34c4762a1bSJed Brown We use by default the backward Euler method. 35c4762a1bSJed Brown 36c4762a1bSJed Brown ------------------------------------------------------------------------- */ 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* 39c4762a1bSJed Brown Include "petscts.h" to use the PETSc timestepping routines. Note that 40c4762a1bSJed Brown this file automatically includes "petscsys.h" and other lower-level 41c4762a1bSJed Brown PETSc include files. 42c4762a1bSJed Brown 43c4762a1bSJed Brown Include the "petscdmda.h" to allow us to use the distributed array data 44c4762a1bSJed Brown structures to manage the parallel grid. 45c4762a1bSJed Brown */ 46c4762a1bSJed Brown #include <petscts.h> 47c4762a1bSJed Brown #include <petscdm.h> 48c4762a1bSJed Brown #include <petscdmda.h> 49c4762a1bSJed Brown #include <petscdraw.h> 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown User-defined application context - contains data needed by the 53c4762a1bSJed Brown application-provided callback routines. 54c4762a1bSJed Brown */ 55c4762a1bSJed Brown typedef struct { 56c4762a1bSJed Brown MPI_Comm comm; /* communicator */ 57c4762a1bSJed Brown DM da; /* distributed array data structure */ 58c4762a1bSJed Brown Vec localwork; /* local ghosted work vector */ 59c4762a1bSJed Brown Vec u_local; /* local ghosted approximate solution vector */ 60c4762a1bSJed Brown Vec solution; /* global exact solution vector */ 61c4762a1bSJed Brown PetscInt m; /* total number of grid points */ 62c4762a1bSJed Brown PetscReal h; /* mesh width: h = 1/(m-1) */ 63c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 64c4762a1bSJed Brown } AppCtx; 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown User-defined routines, provided below. 68c4762a1bSJed Brown */ 69c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*); 70c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 71c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 72c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 73c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*); 74c4762a1bSJed Brown extern PetscErrorCode SetBounds(Vec,Vec,PetscScalar,PetscScalar,AppCtx*); 75c4762a1bSJed Brown 76c4762a1bSJed Brown int main(int argc,char **argv) 77c4762a1bSJed Brown { 78c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 79c4762a1bSJed Brown TS ts; /* timestepping context */ 80c4762a1bSJed Brown Mat A; /* Jacobian matrix data structure */ 81c4762a1bSJed Brown Vec u; /* approximate solution vector */ 82c4762a1bSJed Brown Vec r; /* residual vector */ 83c4762a1bSJed Brown PetscInt time_steps_max = 1000; /* default max timesteps */ 84c4762a1bSJed Brown PetscErrorCode ierr; 85c4762a1bSJed Brown PetscReal dt; 86c4762a1bSJed Brown PetscReal time_total_max = 100.0; /* default max total time */ 87c4762a1bSJed Brown Vec xl,xu; /* Lower and upper bounds on variables */ 88c4762a1bSJed Brown PetscScalar ul=0.0,uh = 3.0; 89c4762a1bSJed Brown PetscBool mymonitor; 90c4762a1bSJed Brown PetscReal bounds[] = {1.0, 3.3}; 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 93c4762a1bSJed Brown Initialize program and set problem parameters 94c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 95c4762a1bSJed Brown 96c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,bounds)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown appctx.comm = PETSC_COMM_WORLD; 100c4762a1bSJed Brown appctx.m = 60; 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&appctx.m,NULL)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ul",&ul,NULL)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-uh",&uh,NULL)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor)); 106c4762a1bSJed Brown appctx.h = 1.0/(appctx.m-1.0); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Create vector data structures 110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* 113c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 114c4762a1bSJed Brown and to set up the ghost point communication pattern. There are M 115c4762a1bSJed Brown total grid values spread equally among all the processors. 116c4762a1bSJed Brown */ 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(appctx.da)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(appctx.da)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 123c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 124c4762a1bSJed Brown have the same types. 125c4762a1bSJed Brown */ 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(appctx.da,&u)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(appctx.da,&appctx.u_local)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* 130c4762a1bSJed Brown Create local work vector for use in evaluating right-hand-side function; 131c4762a1bSJed Brown create global work vector for storing exact solution. 132c4762a1bSJed Brown */ 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(appctx.u_local,&appctx.localwork)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.solution)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* Create residual vector */ 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&r)); 138c4762a1bSJed Brown /* Create lower and upper bound vectors */ 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&xl)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&xu)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetBounds(xl,xu,ul,uh,&appctx)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144c4762a1bSJed Brown Create timestepping solver context; set callback routine for 145c4762a1bSJed Brown right-hand-side function evaluation. 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147c4762a1bSJed Brown 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,r,RHSFunction,&appctx)); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153c4762a1bSJed Brown Set optional user-defined monitoring routine 154c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 155c4762a1bSJed Brown 156c4762a1bSJed Brown if (mymonitor) { 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,Monitor,&appctx,NULL)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 161c4762a1bSJed Brown For nonlinear problems, the user can provide a Jacobian evaluation 162c4762a1bSJed Brown routine (or use a finite differencing approximation). 163c4762a1bSJed Brown 164c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine. 165c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 166c4762a1bSJed Brown 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m)); 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 174c4762a1bSJed Brown Set solution vector and initial timestep 175c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 176c4762a1bSJed Brown 177c4762a1bSJed Brown dt = appctx.h/2.0; 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Customize timestepping solver: 182c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 183c4762a1bSJed Brown - Set timestepping duration info 184c4762a1bSJed Brown Then set runtime options, which can override these defaults. 185c4762a1bSJed Brown For example, 186c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 187c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 189c4762a1bSJed Brown 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSBEULER)); 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,time_steps_max)); 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,time_total_max)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 194c4762a1bSJed Brown /* Set lower and upper bound on the solution vector for each time step */ 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSVISetVariableBounds(ts,xl,xu)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199c4762a1bSJed Brown Solve the problem 200c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* 203c4762a1bSJed Brown Evaluate initial conditions 204c4762a1bSJed Brown */ 205*5f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(u,&appctx)); 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* 208c4762a1bSJed Brown Run the timestepping solver 209c4762a1bSJed Brown */ 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,u)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 213c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 214c4762a1bSJed Brown are no longer needed. 215c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 216c4762a1bSJed Brown 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xl)); 219*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xu)); 220*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 222*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 223*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&appctx.da)); 224*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.localwork)); 225*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.solution)); 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.u_local)); 227c4762a1bSJed Brown 228c4762a1bSJed Brown /* 229c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 230c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 231c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 232c4762a1bSJed Brown options are chosen (e.g., -log_view). 233c4762a1bSJed Brown */ 234c4762a1bSJed Brown ierr = PetscFinalize(); 235c4762a1bSJed Brown return ierr; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 240c4762a1bSJed Brown 241c4762a1bSJed Brown Input Parameters: 242c4762a1bSJed Brown u - uninitialized solution vector (global) 243c4762a1bSJed Brown appctx - user-defined application context 244c4762a1bSJed Brown 245c4762a1bSJed Brown Output Parameter: 246c4762a1bSJed Brown u - vector with solution at initial time (global) 247c4762a1bSJed Brown */ 248c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 249c4762a1bSJed Brown { 250c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h,x; 251c4762a1bSJed Brown PetscInt i,mybase,myend; 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown Determine starting point of each processor's range of 255c4762a1bSJed Brown grid values. 256c4762a1bSJed Brown */ 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(u,&mybase,&myend)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* 260c4762a1bSJed Brown Get a pointer to vector data. 261c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 262c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 263c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 264c4762a1bSJed Brown the array. 265c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 266c4762a1bSJed Brown C version. See the users manual for details. 267c4762a1bSJed Brown */ 268*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u,&u_localptr)); 269c4762a1bSJed Brown 270c4762a1bSJed Brown /* 271c4762a1bSJed Brown We initialize the solution array by simply writing the solution 272c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 273c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 274c4762a1bSJed Brown */ 275c4762a1bSJed Brown for (i=mybase; i<myend; i++) { 276c4762a1bSJed Brown x = h*(PetscReal)i; /* current location in global grid */ 277c4762a1bSJed Brown u_localptr[i-mybase] = 1.0 + x*x; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* 281c4762a1bSJed Brown Restore vector 282c4762a1bSJed Brown */ 283*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u,&u_localptr)); 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* 286c4762a1bSJed Brown Print debugging information if desired 287c4762a1bSJed Brown */ 288c4762a1bSJed Brown if (appctx->debug) { 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"initial guess vector\n")); 290*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown return 0; 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 297c4762a1bSJed Brown /* 298c4762a1bSJed Brown SetBounds - Sets the lower and uper bounds on the interior points 299c4762a1bSJed Brown 300c4762a1bSJed Brown Input parameters: 301c4762a1bSJed Brown xl - vector of lower bounds 302c4762a1bSJed Brown xu - vector of upper bounds 303c4762a1bSJed Brown ul - constant lower bound for all variables 304c4762a1bSJed Brown uh - constant upper bound for all variables 305c4762a1bSJed Brown appctx - Application context 306c4762a1bSJed Brown */ 307c4762a1bSJed Brown PetscErrorCode SetBounds(Vec xl, Vec xu, PetscScalar ul, PetscScalar uh,AppCtx *appctx) 308c4762a1bSJed Brown { 309c4762a1bSJed Brown PetscScalar *l,*u; 310c4762a1bSJed Brown PetscMPIInt rank,size; 311c4762a1bSJed Brown PetscInt localsize; 312c4762a1bSJed Brown 313c4762a1bSJed Brown PetscFunctionBeginUser; 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xl,ul)); 315*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xu,uh)); 316*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(xl,&localsize)); 317*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xl,&l)); 318*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xu,&u)); 319c4762a1bSJed Brown 320*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(appctx->comm,&rank)); 321*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(appctx->comm,&size)); 322dd400576SPatrick Sanan if (rank == 0) { 323c4762a1bSJed Brown l[0] = -PETSC_INFINITY; 324c4762a1bSJed Brown u[0] = PETSC_INFINITY; 325c4762a1bSJed Brown } 326c4762a1bSJed Brown if (rank == size-1) { 327c4762a1bSJed Brown l[localsize-1] = -PETSC_INFINITY; 328c4762a1bSJed Brown u[localsize-1] = PETSC_INFINITY; 329c4762a1bSJed Brown } 330*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xl,&l)); 331*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xu,&u)); 332c4762a1bSJed Brown PetscFunctionReturn(0); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 335c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 336c4762a1bSJed Brown /* 337c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 338c4762a1bSJed Brown 339c4762a1bSJed Brown Input Parameters: 340c4762a1bSJed Brown t - current time 341c4762a1bSJed Brown solution - vector in which exact solution will be computed 342c4762a1bSJed Brown appctx - user-defined application context 343c4762a1bSJed Brown 344c4762a1bSJed Brown Output Parameter: 345c4762a1bSJed Brown solution - vector with the newly computed exact solution 346c4762a1bSJed Brown */ 347c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 348c4762a1bSJed Brown { 349c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,x; 350c4762a1bSJed Brown PetscInt i,mybase,myend; 351c4762a1bSJed Brown 352c4762a1bSJed Brown /* 353c4762a1bSJed Brown Determine starting and ending points of each processor's 354c4762a1bSJed Brown range of grid values 355c4762a1bSJed Brown */ 356*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(solution,&mybase,&myend)); 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* 359c4762a1bSJed Brown Get a pointer to vector data. 360c4762a1bSJed Brown */ 361*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(solution,&s_localptr)); 362c4762a1bSJed Brown 363c4762a1bSJed Brown /* 364c4762a1bSJed Brown Simply write the solution directly into the array locations. 365c4762a1bSJed Brown Alternatively, we could use VecSetValues() or VecSetValuesLocal(). 366c4762a1bSJed Brown */ 367c4762a1bSJed Brown for (i=mybase; i<myend; i++) { 368c4762a1bSJed Brown x = h*(PetscReal)i; 369c4762a1bSJed Brown s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown 372c4762a1bSJed Brown /* 373c4762a1bSJed Brown Restore vector 374c4762a1bSJed Brown */ 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(solution,&s_localptr)); 376c4762a1bSJed Brown return 0; 377c4762a1bSJed Brown } 378c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 379c4762a1bSJed Brown /* 380c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 381c4762a1bSJed Brown each timestep. This example plots the solution and computes the 382c4762a1bSJed Brown error in two different norms. 383c4762a1bSJed Brown 384c4762a1bSJed Brown Input Parameters: 385c4762a1bSJed Brown ts - the timestep context 386c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 387c4762a1bSJed Brown initial condition) 388c4762a1bSJed Brown time - the current time 389c4762a1bSJed Brown u - the solution at this timestep 390c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 391c4762a1bSJed Brown In this case we use the application context which contains 392c4762a1bSJed Brown information about the problem size, workspace and the exact 393c4762a1bSJed Brown solution. 394c4762a1bSJed Brown */ 395c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 396c4762a1bSJed Brown { 397c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 398c4762a1bSJed Brown PetscReal en2,en2s,enmax; 399c4762a1bSJed Brown PetscDraw draw; 400c4762a1bSJed Brown 401c4762a1bSJed Brown /* 402c4762a1bSJed Brown We use the default X windows viewer 403c4762a1bSJed Brown PETSC_VIEWER_DRAW_(appctx->comm) 404c4762a1bSJed Brown that is associated with the current communicator. This saves 405c4762a1bSJed Brown the effort of calling PetscViewerDrawOpen() to create the window. 406c4762a1bSJed Brown Note that if we wished to plot several items in separate windows we 407c4762a1bSJed Brown would create each viewer with PetscViewerDrawOpen() and store them in 408c4762a1bSJed Brown the application context, appctx. 409c4762a1bSJed Brown 410c4762a1bSJed Brown PetscReal buffering makes graphics look better. 411c4762a1bSJed Brown */ 412*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm),0,&draw)); 413*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 414*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_DRAW_(appctx->comm))); 415c4762a1bSJed Brown 416c4762a1bSJed Brown /* 417c4762a1bSJed Brown Compute the exact solution at this timestep 418c4762a1bSJed Brown */ 419*5f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(time,appctx->solution,appctx)); 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* 422c4762a1bSJed Brown Print debugging information if desired 423c4762a1bSJed Brown */ 424c4762a1bSJed Brown if (appctx->debug) { 425*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"Computed solution vector\n")); 426*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 427*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"Exact solution vector\n")); 428*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD)); 429c4762a1bSJed Brown } 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* 432c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 433c4762a1bSJed Brown */ 434*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(appctx->solution,-1.0,u)); 435*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_2,&en2)); 436c4762a1bSJed Brown en2s = PetscSqrtReal(appctx->h)*en2; /* scale the 2-norm by the grid spacing */ 437*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_MAX,&enmax)); 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* 440c4762a1bSJed Brown PetscPrintf() causes only the first processor in this 441c4762a1bSJed Brown communicator to print the timestep information. 442c4762a1bSJed Brown */ 443*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"Timestep %D: time = %g,2-norm error = %g, max norm error = %g\n",step,(double)time,(double)en2s,(double)enmax)); 444c4762a1bSJed Brown 445c4762a1bSJed Brown /* 446c4762a1bSJed Brown Print debugging information if desired 447c4762a1bSJed Brown */ 448c4762a1bSJed Brown /* if (appctx->debug) { 449*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"Error vector\n")); 450*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD)); 451c4762a1bSJed Brown } */ 452c4762a1bSJed Brown return 0; 453c4762a1bSJed Brown } 454c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 455c4762a1bSJed Brown /* 456c4762a1bSJed Brown RHSFunction - User-provided routine that evalues the right-hand-side 457c4762a1bSJed Brown function of the ODE. This routine is set in the main program by 458c4762a1bSJed Brown calling TSSetRHSFunction(). We compute: 459c4762a1bSJed Brown global_out = F(global_in) 460c4762a1bSJed Brown 461c4762a1bSJed Brown Input Parameters: 462c4762a1bSJed Brown ts - timesteping context 463c4762a1bSJed Brown t - current time 464c4762a1bSJed Brown global_in - vector containing the current iterate 465c4762a1bSJed Brown ctx - (optional) user-provided context for function evaluation. 466c4762a1bSJed Brown In this case we use the appctx defined above. 467c4762a1bSJed Brown 468c4762a1bSJed Brown Output Parameter: 469c4762a1bSJed Brown global_out - vector containing the newly evaluated function 470c4762a1bSJed Brown */ 471c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx) 472c4762a1bSJed Brown { 473c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 474c4762a1bSJed Brown DM da = appctx->da; /* distributed array */ 475c4762a1bSJed Brown Vec local_in = appctx->u_local; /* local ghosted input vector */ 476c4762a1bSJed Brown Vec localwork = appctx->localwork; /* local ghosted work vector */ 477c4762a1bSJed Brown PetscInt i,localsize; 478c4762a1bSJed Brown PetscMPIInt rank,size; 479c4762a1bSJed Brown PetscScalar *copyptr,sc; 480c4762a1bSJed Brown const PetscScalar *localptr; 481c4762a1bSJed Brown 482c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 483c4762a1bSJed Brown Get ready for local function computations 484c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 485c4762a1bSJed Brown /* 486c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 487c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 488c4762a1bSJed Brown By placing code between these two statements, computations can be 489c4762a1bSJed Brown done while messages are in transition. 490c4762a1bSJed Brown */ 491*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in)); 492*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in)); 493c4762a1bSJed Brown 494c4762a1bSJed Brown /* 495c4762a1bSJed Brown Access directly the values in our local INPUT work array 496c4762a1bSJed Brown */ 497*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(local_in,&localptr)); 498c4762a1bSJed Brown 499c4762a1bSJed Brown /* 500c4762a1bSJed Brown Access directly the values in our local OUTPUT work array 501c4762a1bSJed Brown */ 502*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(localwork,©ptr)); 503c4762a1bSJed Brown 504c4762a1bSJed Brown sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t)); 505c4762a1bSJed Brown 506c4762a1bSJed Brown /* 507c4762a1bSJed Brown Evaluate our function on the nodes owned by this processor 508c4762a1bSJed Brown */ 509*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(local_in,&localsize)); 510c4762a1bSJed Brown 511c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 512c4762a1bSJed Brown Compute entries for the locally owned part 513c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 514c4762a1bSJed Brown 515c4762a1bSJed Brown /* 516c4762a1bSJed Brown Handle boundary conditions: This is done by using the boundary condition 517c4762a1bSJed Brown u(t,boundary) = g(t,boundary) 518c4762a1bSJed Brown for some function g. Now take the derivative with respect to t to obtain 519c4762a1bSJed Brown u_{t}(t,boundary) = g_{t}(t,boundary) 520c4762a1bSJed Brown 521c4762a1bSJed Brown In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1 522c4762a1bSJed Brown and u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2 523c4762a1bSJed Brown */ 524*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(appctx->comm,&rank)); 525*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(appctx->comm,&size)); 526dd400576SPatrick Sanan if (rank == 0) copyptr[0] = 1.0; 527c4762a1bSJed Brown if (rank == size-1) copyptr[localsize-1] = (t < .5) ? 2.0 : 0.0; 528c4762a1bSJed Brown 529c4762a1bSJed Brown /* 530c4762a1bSJed Brown Handle the interior nodes where the PDE is replace by finite 531c4762a1bSJed Brown difference operators. 532c4762a1bSJed Brown */ 533c4762a1bSJed Brown for (i=1; i<localsize-1; i++) copyptr[i] = localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]); 534c4762a1bSJed Brown 535c4762a1bSJed Brown /* 536c4762a1bSJed Brown Restore vectors 537c4762a1bSJed Brown */ 538*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(local_in,&localptr)); 539*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(localwork,©ptr)); 540c4762a1bSJed Brown 541c4762a1bSJed Brown /* 542c4762a1bSJed Brown Insert values from the local OUTPUT vector into the global 543c4762a1bSJed Brown output vector 544c4762a1bSJed Brown */ 545*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out)); 546*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out)); 547c4762a1bSJed Brown 548c4762a1bSJed Brown /* Print debugging information if desired */ 549c4762a1bSJed Brown /* if (appctx->debug) { 550*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"RHS function vector\n")); 551*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(global_out,PETSC_VIEWER_STDOUT_WORLD)); 552c4762a1bSJed Brown } */ 553c4762a1bSJed Brown 554c4762a1bSJed Brown return 0; 555c4762a1bSJed Brown } 556c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 557c4762a1bSJed Brown /* 558c4762a1bSJed Brown RHSJacobian - User-provided routine to compute the Jacobian of 559c4762a1bSJed Brown the nonlinear right-hand-side function of the ODE. 560c4762a1bSJed Brown 561c4762a1bSJed Brown Input Parameters: 562c4762a1bSJed Brown ts - the TS context 563c4762a1bSJed Brown t - current time 564c4762a1bSJed Brown global_in - global input vector 565c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 566c4762a1bSJed Brown 567c4762a1bSJed Brown Output Parameters: 568c4762a1bSJed Brown AA - Jacobian matrix 569c4762a1bSJed Brown BB - optionally different preconditioning matrix 570c4762a1bSJed Brown str - flag indicating matrix structure 571c4762a1bSJed Brown 572c4762a1bSJed Brown Notes: 573c4762a1bSJed Brown RHSJacobian computes entries for the locally owned part of the Jacobian. 574c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 575c4762a1bSJed Brown contiguous chunks of rows across the processors. 576c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 577c4762a1bSJed Brown locally (but any non-local elements will be sent to the 578c4762a1bSJed Brown appropriate processor during matrix assembly). 579c4762a1bSJed Brown - Always specify global row and columns of matrix entries when 580c4762a1bSJed Brown using MatSetValues(). 581c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 582c4762a1bSJed Brown - Note that MatSetValues() uses 0-based row and column numbers 583c4762a1bSJed Brown in Fortran as well as in C. 584c4762a1bSJed Brown */ 585c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat B,void *ctx) 586c4762a1bSJed Brown { 587c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 588c4762a1bSJed Brown Vec local_in = appctx->u_local; /* local ghosted input vector */ 589c4762a1bSJed Brown DM da = appctx->da; /* distributed array */ 590c4762a1bSJed Brown PetscScalar v[3],sc; 591c4762a1bSJed Brown const PetscScalar *localptr; 592c4762a1bSJed Brown PetscInt i,mstart,mend,mstarts,mends,idx[3],is; 593c4762a1bSJed Brown 594c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 595c4762a1bSJed Brown Get ready for local Jacobian computations 596c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 597c4762a1bSJed Brown /* 598c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 599c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 600c4762a1bSJed Brown By placing code between these two statements, computations can be 601c4762a1bSJed Brown done while messages are in transition. 602c4762a1bSJed Brown */ 603*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in)); 604*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in)); 605c4762a1bSJed Brown 606c4762a1bSJed Brown /* 607c4762a1bSJed Brown Get pointer to vector data 608c4762a1bSJed Brown */ 609*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(local_in,&localptr)); 610c4762a1bSJed Brown 611c4762a1bSJed Brown /* 612c4762a1bSJed Brown Get starting and ending locally owned rows of the matrix 613c4762a1bSJed Brown */ 614*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&mstarts,&mends)); 615c4762a1bSJed Brown mstart = mstarts; mend = mends; 616c4762a1bSJed Brown 617c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 618c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 619c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 620c4762a1bSJed Brown contiguous chunks of rows across the processors. 621c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 622c4762a1bSJed Brown locally (but any non-local elements will be sent to the 623c4762a1bSJed Brown appropriate processor during matrix assembly). 624c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 625c4762a1bSJed Brown - We can set matrix entries either using either 626c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(). 627c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 628c4762a1bSJed Brown 629c4762a1bSJed Brown /* 630c4762a1bSJed Brown Set matrix rows corresponding to boundary data 631c4762a1bSJed Brown */ 632c4762a1bSJed Brown if (mstart == 0) { 633c4762a1bSJed Brown v[0] = 0.0; 634*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&mstart,1,&mstart,v,INSERT_VALUES)); 635c4762a1bSJed Brown mstart++; 636c4762a1bSJed Brown } 637c4762a1bSJed Brown if (mend == appctx->m) { 638c4762a1bSJed Brown mend--; 639c4762a1bSJed Brown v[0] = 0.0; 640*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&mend,1,&mend,v,INSERT_VALUES)); 641c4762a1bSJed Brown } 642c4762a1bSJed Brown 643c4762a1bSJed Brown /* 644c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 645c4762a1bSJed Brown matrix one row at a time. 646c4762a1bSJed Brown */ 647c4762a1bSJed Brown sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t)); 648c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 649c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 650c4762a1bSJed Brown is = i - mstart + 1; 651c4762a1bSJed Brown v[0] = sc*localptr[is]; 652c4762a1bSJed Brown v[1] = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]); 653c4762a1bSJed Brown v[2] = sc*localptr[is]; 654*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,3,idx,v,INSERT_VALUES)); 655c4762a1bSJed Brown } 656c4762a1bSJed Brown 657c4762a1bSJed Brown /* 658c4762a1bSJed Brown Restore vector 659c4762a1bSJed Brown */ 660*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(local_in,&localptr)); 661c4762a1bSJed Brown 662c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 663c4762a1bSJed Brown Complete the matrix assembly process and set some options 664c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 665c4762a1bSJed Brown /* 666c4762a1bSJed Brown Assemble matrix, using the 2-step process: 667c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 668c4762a1bSJed Brown Computations can be done while messages are in transition 669c4762a1bSJed Brown by placing code between these two statements. 670c4762a1bSJed Brown */ 671*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 672*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 673c4762a1bSJed Brown 674c4762a1bSJed Brown /* 675c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 676c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 677c4762a1bSJed Brown */ 678*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 679c4762a1bSJed Brown 680c4762a1bSJed Brown return 0; 681c4762a1bSJed Brown } 682c4762a1bSJed Brown 683c4762a1bSJed Brown /*TEST 684c4762a1bSJed Brown 685c4762a1bSJed Brown test: 686c4762a1bSJed Brown args: -snes_type vinewtonrsls -ts_type glee -mymonitor -ts_max_steps 10 -nox 687c4762a1bSJed Brown requires: !single 688c4762a1bSJed Brown 689c4762a1bSJed Brown TEST*/ 690