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 PetscReal dt; 85c4762a1bSJed Brown PetscReal time_total_max = 100.0; /* default max total time */ 86c4762a1bSJed Brown Vec xl,xu; /* Lower and upper bounds on variables */ 87c4762a1bSJed Brown PetscScalar ul=0.0,uh = 3.0; 88c4762a1bSJed Brown PetscBool mymonitor; 89c4762a1bSJed Brown PetscReal bounds[] = {1.0, 3.3}; 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown Initialize program and set problem parameters 93c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 94c4762a1bSJed Brown 95*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,bounds)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown appctx.comm = PETSC_COMM_WORLD; 99c4762a1bSJed Brown appctx.m = 60; 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&appctx.m,NULL)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ul",&ul,NULL)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-uh",&uh,NULL)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor)); 105c4762a1bSJed Brown appctx.h = 1.0/(appctx.m-1.0); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108c4762a1bSJed Brown Create vector data structures 109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* 112c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 113c4762a1bSJed Brown and to set up the ghost point communication pattern. There are M 114c4762a1bSJed Brown total grid values spread equally among all the processors. 115c4762a1bSJed Brown */ 1165f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(appctx.da)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(appctx.da)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown Extract global and local vectors from DMDA; we use these to store the 122c4762a1bSJed Brown approximate solution. Then duplicate these for remaining vectors that 123c4762a1bSJed Brown have the same types. 124c4762a1bSJed Brown */ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(appctx.da,&u)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(appctx.da,&appctx.u_local)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* 129c4762a1bSJed Brown Create local work vector for use in evaluating right-hand-side function; 130c4762a1bSJed Brown create global work vector for storing exact solution. 131c4762a1bSJed Brown */ 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(appctx.u_local,&appctx.localwork)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&appctx.solution)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Create residual vector */ 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&r)); 137c4762a1bSJed Brown /* Create lower and upper bound vectors */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&xl)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&xu)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(SetBounds(xl,xu,ul,uh,&appctx)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 143c4762a1bSJed Brown Create timestepping solver context; set callback routine for 144c4762a1bSJed Brown right-hand-side function evaluation. 145c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 146c4762a1bSJed Brown 1475f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,r,RHSFunction,&appctx)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152c4762a1bSJed Brown Set optional user-defined monitoring routine 153c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154c4762a1bSJed Brown 155c4762a1bSJed Brown if (mymonitor) { 1565f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,Monitor,&appctx,NULL)); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown For nonlinear problems, the user can provide a Jacobian evaluation 161c4762a1bSJed Brown routine (or use a finite differencing approximation). 162c4762a1bSJed Brown 163c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine. 164c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 165c4762a1bSJed Brown 1665f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx)); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 173c4762a1bSJed Brown Set solution vector and initial timestep 174c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 175c4762a1bSJed Brown 176c4762a1bSJed Brown dt = appctx.h/2.0; 1775f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 180c4762a1bSJed Brown Customize timestepping solver: 181c4762a1bSJed Brown - Set the solution method to be the Backward Euler method. 182c4762a1bSJed Brown - Set timestepping duration info 183c4762a1bSJed Brown Then set runtime options, which can override these defaults. 184c4762a1bSJed Brown For example, 185c4762a1bSJed Brown -ts_max_steps <maxsteps> -ts_max_time <maxtime> 186c4762a1bSJed Brown to override the defaults set by TSSetMaxSteps()/TSSetMaxTime(). 187c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 188c4762a1bSJed Brown 1895f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSBEULER)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,time_steps_max)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,time_total_max)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 193c4762a1bSJed Brown /* Set lower and upper bound on the solution vector for each time step */ 1945f80ce2aSJacob Faibussowitsch CHKERRQ(TSVISetVariableBounds(ts,xl,xu)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 198c4762a1bSJed Brown Solve the problem 199c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* 202c4762a1bSJed Brown Evaluate initial conditions 203c4762a1bSJed Brown */ 2045f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(u,&appctx)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown Run the timestepping solver 208c4762a1bSJed Brown */ 2095f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,u)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 212c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 213c4762a1bSJed Brown are no longer needed. 214c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 215c4762a1bSJed Brown 2165f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xl)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xu)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&appctx.da)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.localwork)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.solution)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&appctx.u_local)); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* 228c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 229c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 230c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 231c4762a1bSJed Brown options are chosen (e.g., -log_view). 232c4762a1bSJed Brown */ 233*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 234*b122ec5aSJacob Faibussowitsch return 0; 235c4762a1bSJed Brown } 236c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown InitialConditions - Computes the solution at the initial time. 239c4762a1bSJed Brown 240c4762a1bSJed Brown Input Parameters: 241c4762a1bSJed Brown u - uninitialized solution vector (global) 242c4762a1bSJed Brown appctx - user-defined application context 243c4762a1bSJed Brown 244c4762a1bSJed Brown Output Parameter: 245c4762a1bSJed Brown u - vector with solution at initial time (global) 246c4762a1bSJed Brown */ 247c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx) 248c4762a1bSJed Brown { 249c4762a1bSJed Brown PetscScalar *u_localptr,h = appctx->h,x; 250c4762a1bSJed Brown PetscInt i,mybase,myend; 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* 253c4762a1bSJed Brown Determine starting point of each processor's range of 254c4762a1bSJed Brown grid values. 255c4762a1bSJed Brown */ 2565f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(u,&mybase,&myend)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* 259c4762a1bSJed Brown Get a pointer to vector data. 260c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 261c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 262c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 263c4762a1bSJed Brown the array. 264c4762a1bSJed Brown - Note that the Fortran interface to VecGetArray() differs from the 265c4762a1bSJed Brown C version. See the users manual for details. 266c4762a1bSJed Brown */ 2675f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u,&u_localptr)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* 270c4762a1bSJed Brown We initialize the solution array by simply writing the solution 271c4762a1bSJed Brown directly into the array locations. Alternatively, we could use 272c4762a1bSJed Brown VecSetValues() or VecSetValuesLocal(). 273c4762a1bSJed Brown */ 274c4762a1bSJed Brown for (i=mybase; i<myend; i++) { 275c4762a1bSJed Brown x = h*(PetscReal)i; /* current location in global grid */ 276c4762a1bSJed Brown u_localptr[i-mybase] = 1.0 + x*x; 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown Restore vector 281c4762a1bSJed Brown */ 2825f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u,&u_localptr)); 283c4762a1bSJed Brown 284c4762a1bSJed Brown /* 285c4762a1bSJed Brown Print debugging information if desired 286c4762a1bSJed Brown */ 287c4762a1bSJed Brown if (appctx->debug) { 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"initial guess vector\n")); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown return 0; 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 296c4762a1bSJed Brown /* 297c4762a1bSJed Brown SetBounds - Sets the lower and uper bounds on the interior points 298c4762a1bSJed Brown 299c4762a1bSJed Brown Input parameters: 300c4762a1bSJed Brown xl - vector of lower bounds 301c4762a1bSJed Brown xu - vector of upper bounds 302c4762a1bSJed Brown ul - constant lower bound for all variables 303c4762a1bSJed Brown uh - constant upper bound for all variables 304c4762a1bSJed Brown appctx - Application context 305c4762a1bSJed Brown */ 306c4762a1bSJed Brown PetscErrorCode SetBounds(Vec xl, Vec xu, PetscScalar ul, PetscScalar uh,AppCtx *appctx) 307c4762a1bSJed Brown { 308c4762a1bSJed Brown PetscScalar *l,*u; 309c4762a1bSJed Brown PetscMPIInt rank,size; 310c4762a1bSJed Brown PetscInt localsize; 311c4762a1bSJed Brown 312c4762a1bSJed Brown PetscFunctionBeginUser; 3135f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xl,ul)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xu,uh)); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(xl,&localsize)); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xl,&l)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xu,&u)); 318c4762a1bSJed Brown 3195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(appctx->comm,&rank)); 3205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(appctx->comm,&size)); 321dd400576SPatrick Sanan if (rank == 0) { 322c4762a1bSJed Brown l[0] = -PETSC_INFINITY; 323c4762a1bSJed Brown u[0] = PETSC_INFINITY; 324c4762a1bSJed Brown } 325c4762a1bSJed Brown if (rank == size-1) { 326c4762a1bSJed Brown l[localsize-1] = -PETSC_INFINITY; 327c4762a1bSJed Brown u[localsize-1] = PETSC_INFINITY; 328c4762a1bSJed Brown } 3295f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xl,&l)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xu,&u)); 331c4762a1bSJed Brown PetscFunctionReturn(0); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown 334c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 335c4762a1bSJed Brown /* 336c4762a1bSJed Brown ExactSolution - Computes the exact solution at a given time. 337c4762a1bSJed Brown 338c4762a1bSJed Brown Input Parameters: 339c4762a1bSJed Brown t - current time 340c4762a1bSJed Brown solution - vector in which exact solution will be computed 341c4762a1bSJed Brown appctx - user-defined application context 342c4762a1bSJed Brown 343c4762a1bSJed Brown Output Parameter: 344c4762a1bSJed Brown solution - vector with the newly computed exact solution 345c4762a1bSJed Brown */ 346c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx) 347c4762a1bSJed Brown { 348c4762a1bSJed Brown PetscScalar *s_localptr,h = appctx->h,x; 349c4762a1bSJed Brown PetscInt i,mybase,myend; 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* 352c4762a1bSJed Brown Determine starting and ending points of each processor's 353c4762a1bSJed Brown range of grid values 354c4762a1bSJed Brown */ 3555f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(solution,&mybase,&myend)); 356c4762a1bSJed Brown 357c4762a1bSJed Brown /* 358c4762a1bSJed Brown Get a pointer to vector data. 359c4762a1bSJed Brown */ 3605f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(solution,&s_localptr)); 361c4762a1bSJed Brown 362c4762a1bSJed Brown /* 363c4762a1bSJed Brown Simply write the solution directly into the array locations. 364c4762a1bSJed Brown Alternatively, we could use VecSetValues() or VecSetValuesLocal(). 365c4762a1bSJed Brown */ 366c4762a1bSJed Brown for (i=mybase; i<myend; i++) { 367c4762a1bSJed Brown x = h*(PetscReal)i; 368c4762a1bSJed Brown s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown 371c4762a1bSJed Brown /* 372c4762a1bSJed Brown Restore vector 373c4762a1bSJed Brown */ 3745f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(solution,&s_localptr)); 375c4762a1bSJed Brown return 0; 376c4762a1bSJed Brown } 377c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 378c4762a1bSJed Brown /* 379c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 380c4762a1bSJed Brown each timestep. This example plots the solution and computes the 381c4762a1bSJed Brown error in two different norms. 382c4762a1bSJed Brown 383c4762a1bSJed Brown Input Parameters: 384c4762a1bSJed Brown ts - the timestep context 385c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 386c4762a1bSJed Brown initial condition) 387c4762a1bSJed Brown time - the current time 388c4762a1bSJed Brown u - the solution at this timestep 389c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 390c4762a1bSJed Brown In this case we use the application context which contains 391c4762a1bSJed Brown information about the problem size, workspace and the exact 392c4762a1bSJed Brown solution. 393c4762a1bSJed Brown */ 394c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 395c4762a1bSJed Brown { 396c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 397c4762a1bSJed Brown PetscReal en2,en2s,enmax; 398c4762a1bSJed Brown PetscDraw draw; 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* 401c4762a1bSJed Brown We use the default X windows viewer 402c4762a1bSJed Brown PETSC_VIEWER_DRAW_(appctx->comm) 403c4762a1bSJed Brown that is associated with the current communicator. This saves 404c4762a1bSJed Brown the effort of calling PetscViewerDrawOpen() to create the window. 405c4762a1bSJed Brown Note that if we wished to plot several items in separate windows we 406c4762a1bSJed Brown would create each viewer with PetscViewerDrawOpen() and store them in 407c4762a1bSJed Brown the application context, appctx. 408c4762a1bSJed Brown 409c4762a1bSJed Brown PetscReal buffering makes graphics look better. 410c4762a1bSJed Brown */ 4115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(PETSC_VIEWER_DRAW_(appctx->comm),0,&draw)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_DRAW_(appctx->comm))); 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* 416c4762a1bSJed Brown Compute the exact solution at this timestep 417c4762a1bSJed Brown */ 4185f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(time,appctx->solution,appctx)); 419c4762a1bSJed Brown 420c4762a1bSJed Brown /* 421c4762a1bSJed Brown Print debugging information if desired 422c4762a1bSJed Brown */ 423c4762a1bSJed Brown if (appctx->debug) { 4245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"Computed solution vector\n")); 4255f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD)); 4265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"Exact solution vector\n")); 4275f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD)); 428c4762a1bSJed Brown } 429c4762a1bSJed Brown 430c4762a1bSJed Brown /* 431c4762a1bSJed Brown Compute the 2-norm and max-norm of the error 432c4762a1bSJed Brown */ 4335f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(appctx->solution,-1.0,u)); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_2,&en2)); 435c4762a1bSJed Brown en2s = PetscSqrtReal(appctx->h)*en2; /* scale the 2-norm by the grid spacing */ 4365f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(appctx->solution,NORM_MAX,&enmax)); 437c4762a1bSJed Brown 438c4762a1bSJed Brown /* 439c4762a1bSJed Brown PetscPrintf() causes only the first processor in this 440c4762a1bSJed Brown communicator to print the timestep information. 441c4762a1bSJed Brown */ 4425f80ce2aSJacob 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)); 443c4762a1bSJed Brown 444c4762a1bSJed Brown /* 445c4762a1bSJed Brown Print debugging information if desired 446c4762a1bSJed Brown */ 447c4762a1bSJed Brown /* if (appctx->debug) { 4485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"Error vector\n")); 4495f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD)); 450c4762a1bSJed Brown } */ 451c4762a1bSJed Brown return 0; 452c4762a1bSJed Brown } 453c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 454c4762a1bSJed Brown /* 455c4762a1bSJed Brown RHSFunction - User-provided routine that evalues the right-hand-side 456c4762a1bSJed Brown function of the ODE. This routine is set in the main program by 457c4762a1bSJed Brown calling TSSetRHSFunction(). We compute: 458c4762a1bSJed Brown global_out = F(global_in) 459c4762a1bSJed Brown 460c4762a1bSJed Brown Input Parameters: 461c4762a1bSJed Brown ts - timesteping context 462c4762a1bSJed Brown t - current time 463c4762a1bSJed Brown global_in - vector containing the current iterate 464c4762a1bSJed Brown ctx - (optional) user-provided context for function evaluation. 465c4762a1bSJed Brown In this case we use the appctx defined above. 466c4762a1bSJed Brown 467c4762a1bSJed Brown Output Parameter: 468c4762a1bSJed Brown global_out - vector containing the newly evaluated function 469c4762a1bSJed Brown */ 470c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx) 471c4762a1bSJed Brown { 472c4762a1bSJed Brown AppCtx *appctx = (AppCtx*) ctx; /* user-defined application context */ 473c4762a1bSJed Brown DM da = appctx->da; /* distributed array */ 474c4762a1bSJed Brown Vec local_in = appctx->u_local; /* local ghosted input vector */ 475c4762a1bSJed Brown Vec localwork = appctx->localwork; /* local ghosted work vector */ 476c4762a1bSJed Brown PetscInt i,localsize; 477c4762a1bSJed Brown PetscMPIInt rank,size; 478c4762a1bSJed Brown PetscScalar *copyptr,sc; 479c4762a1bSJed Brown const PetscScalar *localptr; 480c4762a1bSJed Brown 481c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 482c4762a1bSJed Brown Get ready for local function computations 483c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 484c4762a1bSJed Brown /* 485c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 486c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 487c4762a1bSJed Brown By placing code between these two statements, computations can be 488c4762a1bSJed Brown done while messages are in transition. 489c4762a1bSJed Brown */ 4905f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in)); 4915f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in)); 492c4762a1bSJed Brown 493c4762a1bSJed Brown /* 494c4762a1bSJed Brown Access directly the values in our local INPUT work array 495c4762a1bSJed Brown */ 4965f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(local_in,&localptr)); 497c4762a1bSJed Brown 498c4762a1bSJed Brown /* 499c4762a1bSJed Brown Access directly the values in our local OUTPUT work array 500c4762a1bSJed Brown */ 5015f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(localwork,©ptr)); 502c4762a1bSJed Brown 503c4762a1bSJed Brown sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t)); 504c4762a1bSJed Brown 505c4762a1bSJed Brown /* 506c4762a1bSJed Brown Evaluate our function on the nodes owned by this processor 507c4762a1bSJed Brown */ 5085f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(local_in,&localsize)); 509c4762a1bSJed Brown 510c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 511c4762a1bSJed Brown Compute entries for the locally owned part 512c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 513c4762a1bSJed Brown 514c4762a1bSJed Brown /* 515c4762a1bSJed Brown Handle boundary conditions: This is done by using the boundary condition 516c4762a1bSJed Brown u(t,boundary) = g(t,boundary) 517c4762a1bSJed Brown for some function g. Now take the derivative with respect to t to obtain 518c4762a1bSJed Brown u_{t}(t,boundary) = g_{t}(t,boundary) 519c4762a1bSJed Brown 520c4762a1bSJed Brown In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1 521c4762a1bSJed Brown and u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2 522c4762a1bSJed Brown */ 5235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(appctx->comm,&rank)); 5245f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(appctx->comm,&size)); 525dd400576SPatrick Sanan if (rank == 0) copyptr[0] = 1.0; 526c4762a1bSJed Brown if (rank == size-1) copyptr[localsize-1] = (t < .5) ? 2.0 : 0.0; 527c4762a1bSJed Brown 528c4762a1bSJed Brown /* 529c4762a1bSJed Brown Handle the interior nodes where the PDE is replace by finite 530c4762a1bSJed Brown difference operators. 531c4762a1bSJed Brown */ 532c4762a1bSJed Brown for (i=1; i<localsize-1; i++) copyptr[i] = localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]); 533c4762a1bSJed Brown 534c4762a1bSJed Brown /* 535c4762a1bSJed Brown Restore vectors 536c4762a1bSJed Brown */ 5375f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(local_in,&localptr)); 5385f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(localwork,©ptr)); 539c4762a1bSJed Brown 540c4762a1bSJed Brown /* 541c4762a1bSJed Brown Insert values from the local OUTPUT vector into the global 542c4762a1bSJed Brown output vector 543c4762a1bSJed Brown */ 5445f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out)); 5455f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out)); 546c4762a1bSJed Brown 547c4762a1bSJed Brown /* Print debugging information if desired */ 548c4762a1bSJed Brown /* if (appctx->debug) { 5495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(appctx->comm,"RHS function vector\n")); 5505f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(global_out,PETSC_VIEWER_STDOUT_WORLD)); 551c4762a1bSJed Brown } */ 552c4762a1bSJed Brown 553c4762a1bSJed Brown return 0; 554c4762a1bSJed Brown } 555c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 556c4762a1bSJed Brown /* 557c4762a1bSJed Brown RHSJacobian - User-provided routine to compute the Jacobian of 558c4762a1bSJed Brown the nonlinear right-hand-side function of the ODE. 559c4762a1bSJed Brown 560c4762a1bSJed Brown Input Parameters: 561c4762a1bSJed Brown ts - the TS context 562c4762a1bSJed Brown t - current time 563c4762a1bSJed Brown global_in - global input vector 564c4762a1bSJed Brown dummy - optional user-defined context, as set by TSetRHSJacobian() 565c4762a1bSJed Brown 566c4762a1bSJed Brown Output Parameters: 567c4762a1bSJed Brown AA - Jacobian matrix 568c4762a1bSJed Brown BB - optionally different preconditioning matrix 569c4762a1bSJed Brown str - flag indicating matrix structure 570c4762a1bSJed Brown 571c4762a1bSJed Brown Notes: 572c4762a1bSJed Brown RHSJacobian computes entries for the locally owned part of the Jacobian. 573c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 574c4762a1bSJed Brown contiguous chunks of rows across the processors. 575c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 576c4762a1bSJed Brown locally (but any non-local elements will be sent to the 577c4762a1bSJed Brown appropriate processor during matrix assembly). 578c4762a1bSJed Brown - Always specify global row and columns of matrix entries when 579c4762a1bSJed Brown using MatSetValues(). 580c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 581c4762a1bSJed Brown - Note that MatSetValues() uses 0-based row and column numbers 582c4762a1bSJed Brown in Fortran as well as in C. 583c4762a1bSJed Brown */ 584c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat B,void *ctx) 585c4762a1bSJed Brown { 586c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 587c4762a1bSJed Brown Vec local_in = appctx->u_local; /* local ghosted input vector */ 588c4762a1bSJed Brown DM da = appctx->da; /* distributed array */ 589c4762a1bSJed Brown PetscScalar v[3],sc; 590c4762a1bSJed Brown const PetscScalar *localptr; 591c4762a1bSJed Brown PetscInt i,mstart,mend,mstarts,mends,idx[3],is; 592c4762a1bSJed Brown 593c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 594c4762a1bSJed Brown Get ready for local Jacobian computations 595c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 596c4762a1bSJed Brown /* 597c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 598c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 599c4762a1bSJed Brown By placing code between these two statements, computations can be 600c4762a1bSJed Brown done while messages are in transition. 601c4762a1bSJed Brown */ 6025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in)); 6035f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in)); 604c4762a1bSJed Brown 605c4762a1bSJed Brown /* 606c4762a1bSJed Brown Get pointer to vector data 607c4762a1bSJed Brown */ 6085f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(local_in,&localptr)); 609c4762a1bSJed Brown 610c4762a1bSJed Brown /* 611c4762a1bSJed Brown Get starting and ending locally owned rows of the matrix 612c4762a1bSJed Brown */ 6135f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&mstarts,&mends)); 614c4762a1bSJed Brown mstart = mstarts; mend = mends; 615c4762a1bSJed Brown 616c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 617c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 618c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 619c4762a1bSJed Brown contiguous chunks of rows across the processors. 620c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 621c4762a1bSJed Brown locally (but any non-local elements will be sent to the 622c4762a1bSJed Brown appropriate processor during matrix assembly). 623c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 624c4762a1bSJed Brown - We can set matrix entries either using either 625c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(). 626c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 627c4762a1bSJed Brown 628c4762a1bSJed Brown /* 629c4762a1bSJed Brown Set matrix rows corresponding to boundary data 630c4762a1bSJed Brown */ 631c4762a1bSJed Brown if (mstart == 0) { 632c4762a1bSJed Brown v[0] = 0.0; 6335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&mstart,1,&mstart,v,INSERT_VALUES)); 634c4762a1bSJed Brown mstart++; 635c4762a1bSJed Brown } 636c4762a1bSJed Brown if (mend == appctx->m) { 637c4762a1bSJed Brown mend--; 638c4762a1bSJed Brown v[0] = 0.0; 6395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&mend,1,&mend,v,INSERT_VALUES)); 640c4762a1bSJed Brown } 641c4762a1bSJed Brown 642c4762a1bSJed Brown /* 643c4762a1bSJed Brown Set matrix rows corresponding to interior data. We construct the 644c4762a1bSJed Brown matrix one row at a time. 645c4762a1bSJed Brown */ 646c4762a1bSJed Brown sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t)); 647c4762a1bSJed Brown for (i=mstart; i<mend; i++) { 648c4762a1bSJed Brown idx[0] = i-1; idx[1] = i; idx[2] = i+1; 649c4762a1bSJed Brown is = i - mstart + 1; 650c4762a1bSJed Brown v[0] = sc*localptr[is]; 651c4762a1bSJed Brown v[1] = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]); 652c4762a1bSJed Brown v[2] = sc*localptr[is]; 6535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,3,idx,v,INSERT_VALUES)); 654c4762a1bSJed Brown } 655c4762a1bSJed Brown 656c4762a1bSJed Brown /* 657c4762a1bSJed Brown Restore vector 658c4762a1bSJed Brown */ 6595f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(local_in,&localptr)); 660c4762a1bSJed Brown 661c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 662c4762a1bSJed Brown Complete the matrix assembly process and set some options 663c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 664c4762a1bSJed Brown /* 665c4762a1bSJed Brown Assemble matrix, using the 2-step process: 666c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd() 667c4762a1bSJed Brown Computations can be done while messages are in transition 668c4762a1bSJed Brown by placing code between these two statements. 669c4762a1bSJed Brown */ 6705f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 6715f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 672c4762a1bSJed Brown 673c4762a1bSJed Brown /* 674c4762a1bSJed Brown Set and option to indicate that we will never add a new nonzero location 675c4762a1bSJed Brown to the matrix. If we do, it will generate an error. 676c4762a1bSJed Brown */ 6775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 678c4762a1bSJed Brown 679c4762a1bSJed Brown return 0; 680c4762a1bSJed Brown } 681c4762a1bSJed Brown 682c4762a1bSJed Brown /*TEST 683c4762a1bSJed Brown 684c4762a1bSJed Brown test: 685c4762a1bSJed Brown args: -snes_type vinewtonrsls -ts_type glee -mymonitor -ts_max_steps 10 -nox 686c4762a1bSJed Brown requires: !single 687c4762a1bSJed Brown 688c4762a1bSJed Brown TEST*/ 689