1a80e93e8SEmil static const char help[] = "Time-dependent Brusselator reaction-diffusion PDE in 1d formulated as a PDAE. Demonstrates solving PDEs with algebraic constraints (PDAE).\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown u_t - alpha u_xx = A + u^2 v - (B+1) u 4c4762a1bSJed Brown v_t - alpha v_xx = B u - u^2 v 5c4762a1bSJed Brown 0 < x < 1; 6c4762a1bSJed Brown A = 1, B = 3, alpha = 1/50 7c4762a1bSJed Brown 8c4762a1bSJed Brown Initial conditions: 9c4762a1bSJed Brown u(x,0) = 1 + sin(2 pi x) 10c4762a1bSJed Brown v(x,0) = 3 11c4762a1bSJed Brown 12c4762a1bSJed Brown Boundary conditions: 13c4762a1bSJed Brown u(0,t) = u(1,t) = 1 14c4762a1bSJed Brown v(0,t) = v(1,t) = 3 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown 17c4762a1bSJed Brown #include <petscdm.h> 18c4762a1bSJed Brown #include <petscdmda.h> 19c4762a1bSJed Brown #include <petscts.h> 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscScalar u,v; 23c4762a1bSJed Brown } Field; 24c4762a1bSJed Brown 25c4762a1bSJed Brown typedef struct _User *User; 26c4762a1bSJed Brown struct _User { 27c4762a1bSJed Brown PetscReal A,B; /* Reaction coefficients */ 28c4762a1bSJed Brown PetscReal alpha; /* Diffusion coefficient */ 29c4762a1bSJed Brown PetscReal uleft,uright; /* Dirichlet boundary conditions */ 30c4762a1bSJed Brown PetscReal vleft,vright; /* Dirichlet boundary conditions */ 31c4762a1bSJed Brown }; 32c4762a1bSJed Brown 33c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 34c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 35c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 36c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 37c4762a1bSJed Brown 38c4762a1bSJed Brown int main(int argc,char **argv) 39c4762a1bSJed Brown { 40c4762a1bSJed Brown TS ts; /* nonlinear solver */ 41c4762a1bSJed Brown Vec X; /* solution, residual vectors */ 42c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 43c4762a1bSJed Brown PetscInt steps,mx; 44c4762a1bSJed Brown PetscErrorCode ierr; 45c4762a1bSJed Brown DM da; 46c4762a1bSJed Brown PetscReal ftime,hx,dt; 47c4762a1bSJed Brown struct _User user; /* user-defined work context */ 48c4762a1bSJed Brown TSConvergedReason reason; 49c4762a1bSJed Brown 50*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 51c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 52c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 53c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 59c4762a1bSJed Brown Extract global vectors from DMDA; 60c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 615f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&X)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Initialize user application context */ 641e1ea65dSPierre Jolivet ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");CHKERRQ(ierr); 65c4762a1bSJed Brown { 66c4762a1bSJed Brown user.A = 1; 67c4762a1bSJed Brown user.B = 3; 68c4762a1bSJed Brown user.alpha = 0.02; 69c4762a1bSJed Brown user.uleft = 1; 70c4762a1bSJed Brown user.uright = 1; 71c4762a1bSJed Brown user.vleft = 3; 72c4762a1bSJed Brown user.vright = 3; 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Create timestepping solver context 85c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 865f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSARKIMEX)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,FormIFunction,&user)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATAIJ)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&J)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,J,J,FormIJacobian,&user)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown ftime = 10.0; 975f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,ftime)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101c4762a1bSJed Brown Set initial conditions 102c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(FormInitialSolution(ts,X,&user)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,X)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(X,&mx)); 106c4762a1bSJed Brown hx = 1.0/(PetscReal)(mx/2-1); 107c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */ 1085f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 111c4762a1bSJed Brown Set runtime options 112c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1135f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Solve nonlinear system 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1185f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,X)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetConvergedReason(ts,&reason)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps)); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125c4762a1bSJed Brown Free work space. 126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 131*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 132*b122ec5aSJacob Faibussowitsch return 0; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 136c4762a1bSJed Brown { 137c4762a1bSJed Brown User user = (User)ptr; 138c4762a1bSJed Brown DM da; 139c4762a1bSJed Brown DMDALocalInfo info; 140c4762a1bSJed Brown PetscInt i; 141c4762a1bSJed Brown Field *x,*xdot,*f; 142c4762a1bSJed Brown PetscReal hx; 143c4762a1bSJed Brown Vec Xloc; 144c4762a1bSJed Brown 145c4762a1bSJed Brown PetscFunctionBeginUser; 1465f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 148c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* 151c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 152c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 153c4762a1bSJed Brown By placing code between these two statements, computations can be 154c4762a1bSJed Brown done while messages are in transition. 155c4762a1bSJed Brown */ 1565f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&Xloc)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* Get pointers to vector data */ 1615f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xloc,&x)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xdot,&xdot)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 166c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 167c4762a1bSJed Brown if (i == 0) { 168c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->uleft); 169c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->vleft); 170c4762a1bSJed Brown } else if (i == info.mx-1) { 171c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->uright); 172c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->vright); 173c4762a1bSJed Brown } else { 174c4762a1bSJed Brown f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx; 175c4762a1bSJed Brown f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx; 176c4762a1bSJed Brown } 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* Restore vectors */ 1805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xloc,&x)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xdot,&xdot)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&Xloc)); 184c4762a1bSJed Brown PetscFunctionReturn(0); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 188c4762a1bSJed Brown { 189c4762a1bSJed Brown User user = (User)ptr; 190c4762a1bSJed Brown DM da; 191c4762a1bSJed Brown DMDALocalInfo info; 192c4762a1bSJed Brown PetscInt i; 193c4762a1bSJed Brown PetscReal hx; 194c4762a1bSJed Brown Field *x,*f; 195c4762a1bSJed Brown 196c4762a1bSJed Brown PetscFunctionBeginUser; 1975f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 199c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Get pointers to vector data */ 2025f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X,&x)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 206c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 207c4762a1bSJed Brown PetscScalar u = x[i].u,v = x[i].v; 208c4762a1bSJed Brown f[i].u = hx*(user->A + u*u*v - (user->B+1)*u); 209c4762a1bSJed Brown f[i].v = hx*(user->B*u - u*u*v); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Restore vectors */ 2135f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X,&x)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 215c4762a1bSJed Brown PetscFunctionReturn(0); 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 219c4762a1bSJed Brown /* 220c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 221c4762a1bSJed Brown */ 222c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 223c4762a1bSJed Brown { 224c4762a1bSJed Brown User user = (User)ptr; 225c4762a1bSJed Brown DMDALocalInfo info; 226c4762a1bSJed Brown PetscInt i; 227c4762a1bSJed Brown PetscReal hx; 228c4762a1bSJed Brown DM da; 229c4762a1bSJed Brown Field *x,*xdot; 230c4762a1bSJed Brown 231c4762a1bSJed Brown PetscFunctionBeginUser; 2325f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 234c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* Get pointers to vector data */ 2375f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X,&x)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xdot,&xdot)); 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 241c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 242c4762a1bSJed Brown if (i == 0 || i == info.mx-1) { 243c4762a1bSJed Brown const PetscInt row = i,col = i; 244c4762a1bSJed Brown const PetscScalar vals[2][2] = {{hx,0},{0,hx}}; 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES)); 246c4762a1bSJed Brown } else { 247c4762a1bSJed Brown const PetscInt row = i,col[] = {i-1,i,i+1}; 248c4762a1bSJed Brown const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx; 249c4762a1bSJed Brown const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}}, 250c4762a1bSJed Brown {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}}; 2515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES)); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown } 254c4762a1bSJed Brown 255c4762a1bSJed Brown /* Restore vectors */ 2565f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X,&x)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xdot,&xdot)); 258c4762a1bSJed Brown 2595f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 261c4762a1bSJed Brown if (J != Jpre) { 2625f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown PetscFunctionReturn(0); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 269c4762a1bSJed Brown { 270c4762a1bSJed Brown User user = (User)ctx; 271c4762a1bSJed Brown DM da; 272c4762a1bSJed Brown PetscInt i; 273c4762a1bSJed Brown DMDALocalInfo info; 274c4762a1bSJed Brown Field *x; 275c4762a1bSJed Brown PetscReal hx; 276c4762a1bSJed Brown 277c4762a1bSJed Brown PetscFunctionBeginUser; 2785f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 280c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 281c4762a1bSJed Brown 282c4762a1bSJed Brown /* Get pointers to vector data */ 2835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,X,&x)); 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 286c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 287c4762a1bSJed Brown PetscReal xi = i*hx; 288c4762a1bSJed Brown x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi); 289c4762a1bSJed Brown x[i].v = user->vleft*(1.-xi) + user->vright*xi; 290c4762a1bSJed Brown } 2915f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,X,&x)); 292c4762a1bSJed Brown PetscFunctionReturn(0); 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown /*TEST 296c4762a1bSJed Brown 297c4762a1bSJed Brown test: 298c4762a1bSJed Brown args: -nox -da_grid_x 20 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none 299c4762a1bSJed Brown 300c4762a1bSJed Brown TEST*/ 301