1c4762a1bSJed Brown static const char help[] = "Call PetscInitialize multiple times.\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown This example is based on the Brusselator tutorial of the same name, but tests multiple calls to PetscInitialize(). 4c4762a1bSJed Brown This is a bad "convergence study" because it only compares min and max values of the solution rather than comparing 5c4762a1bSJed Brown norms of the errors. For convergence studies, we recommend invoking PetscInitialize() only once and comparing norms 6c4762a1bSJed Brown of errors (perhaps estimated using an accurate reference solution). 7c4762a1bSJed Brown 8c4762a1bSJed Brown Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods and multiple solves. 9c4762a1bSJed Brown 10c4762a1bSJed Brown u_t - alpha u_xx = A + u^2 v - (B+1) u 11c4762a1bSJed Brown v_t - alpha v_xx = B u - u^2 v 12c4762a1bSJed Brown 0 < x < 1; 13c4762a1bSJed Brown A = 1, B = 3, alpha = 1/10 14c4762a1bSJed Brown 15c4762a1bSJed Brown Initial conditions: 16c4762a1bSJed Brown u(x,0) = 1 + sin(2 pi x) 17c4762a1bSJed Brown v(x,0) = 3 18c4762a1bSJed Brown 19c4762a1bSJed Brown Boundary conditions: 20c4762a1bSJed Brown u(0,t) = u(1,t) = 1 21c4762a1bSJed Brown v(0,t) = v(1,t) = 3 22c4762a1bSJed Brown */ 23c4762a1bSJed Brown 24c4762a1bSJed Brown #include <petscdm.h> 25c4762a1bSJed Brown #include <petscdmda.h> 26c4762a1bSJed Brown #include <petscts.h> 27c4762a1bSJed Brown 28c4762a1bSJed Brown typedef struct { 29c4762a1bSJed Brown PetscScalar u,v; 30c4762a1bSJed Brown } Field; 31c4762a1bSJed Brown 32c4762a1bSJed Brown typedef struct _User *User; 33c4762a1bSJed Brown struct _User { 34c4762a1bSJed Brown PetscReal A,B; /* Reaction coefficients */ 35c4762a1bSJed Brown PetscReal alpha; /* Diffusion coefficient */ 36c4762a1bSJed Brown PetscReal uleft,uright; /* Dirichlet boundary conditions */ 37c4762a1bSJed Brown PetscReal vleft,vright; /* Dirichlet boundary conditions */ 38c4762a1bSJed Brown }; 39c4762a1bSJed Brown 40c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 41c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 42c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 43c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 44c4762a1bSJed Brown static int Brusselator(int,char**,PetscInt); 45c4762a1bSJed Brown 46c4762a1bSJed Brown int main(int argc,char **argv) 47c4762a1bSJed Brown { 48c4762a1bSJed Brown PetscInt cycle; 49c4762a1bSJed Brown PetscErrorCode ierr; 50c4762a1bSJed Brown 51c4762a1bSJed Brown ierr = MPI_Init(&argc,&argv);if (ierr) return ierr; 52c4762a1bSJed Brown for (cycle=0; cycle<4; cycle++) { 53c4762a1bSJed Brown ierr = Brusselator(argc,argv,cycle); 54c4762a1bSJed Brown if (ierr) return 1; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown ierr = MPI_Finalize(); 57c4762a1bSJed Brown return ierr; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown PetscErrorCode Brusselator(int argc,char **argv,PetscInt cycle) 61c4762a1bSJed Brown { 62c4762a1bSJed Brown TS ts; /* nonlinear solver */ 63c4762a1bSJed Brown Vec X; /* solution, residual vectors */ 64c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 65c4762a1bSJed Brown PetscInt steps,mx; 66c4762a1bSJed Brown PetscErrorCode ierr; 67c4762a1bSJed Brown DM da; 68c4762a1bSJed Brown PetscReal ftime,hx,dt,xmax,xmin; 69c4762a1bSJed Brown struct _User user; /* user-defined work context */ 70c4762a1bSJed Brown TSConvergedReason reason; 71c4762a1bSJed Brown 72*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 76c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 82c4762a1bSJed Brown Extract global vectors from DMDA; 83c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 845f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&X)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Initialize user application context */ 871e1ea65dSPierre Jolivet ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");CHKERRQ(ierr); 88c4762a1bSJed Brown { 89c4762a1bSJed Brown user.A = 1; 90c4762a1bSJed Brown user.B = 3; 91c4762a1bSJed Brown user.alpha = 0.1; 92c4762a1bSJed Brown user.uleft = 1; 93c4762a1bSJed Brown user.uright = 1; 94c4762a1bSJed Brown user.vleft = 3; 95c4762a1bSJed Brown user.vright = 3; 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL)); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107c4762a1bSJed Brown Create timestepping solver context 108c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSARKIMEX)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,FormIFunction,&user)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATAIJ)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&J)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,J,J,FormIJacobian,&user)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown ftime = 1.0; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,ftime)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123c4762a1bSJed Brown Set initial conditions 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(FormInitialSolution(ts,X,&user)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,X)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(X,&mx)); 128c4762a1bSJed Brown hx = 1.0/(PetscReal)(mx/2-1); 129c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */ 130c4762a1bSJed Brown dt *= PetscPowRealInt(0.2,cycle); /* Shrink the time step in convergence study. */ 1315f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Set runtime options 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1375f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140c4762a1bSJed Brown Solve nonlinear system 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1425f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,X)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetConvergedReason(ts,&reason)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(VecMin(X,NULL,&xmin)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(VecMax(X,NULL,&xmax)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after % 3D steps. Range [%6.4f,%6.4f]\n",TSConvergedReasons[reason],(double)ftime,steps,(double)xmin,(double)xmax)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151c4762a1bSJed Brown Free work space. 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 157*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 158*b122ec5aSJacob Faibussowitsch return 0; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown User user = (User)ptr; 164c4762a1bSJed Brown DM da; 165c4762a1bSJed Brown DMDALocalInfo info; 166c4762a1bSJed Brown PetscInt i; 167c4762a1bSJed Brown Field *x,*xdot,*f; 168c4762a1bSJed Brown PetscReal hx; 169c4762a1bSJed Brown Vec Xloc; 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscFunctionBeginUser; 1725f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 174c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* 177c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 178c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 179c4762a1bSJed Brown By placing code between these two statements, computations can be 180c4762a1bSJed Brown done while messages are in transition. 181c4762a1bSJed Brown */ 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&Xloc)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* Get pointers to vector data */ 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xloc,&x)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xdot,&xdot)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 192c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 193c4762a1bSJed Brown if (i == 0) { 194c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->uleft); 195c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->vleft); 196c4762a1bSJed Brown } else if (i == info.mx-1) { 197c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->uright); 198c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->vright); 199c4762a1bSJed Brown } else { 200c4762a1bSJed Brown f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx; 201c4762a1bSJed Brown f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx; 202c4762a1bSJed Brown } 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* Restore vectors */ 2065f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xloc,&x)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xdot,&xdot)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&Xloc)); 210c4762a1bSJed Brown PetscFunctionReturn(0); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown 213c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 214c4762a1bSJed Brown { 215c4762a1bSJed Brown User user = (User)ptr; 216c4762a1bSJed Brown DM da; 217c4762a1bSJed Brown DMDALocalInfo info; 218c4762a1bSJed Brown PetscInt i; 219c4762a1bSJed Brown PetscReal hx; 220c4762a1bSJed Brown Field *x,*f; 221c4762a1bSJed Brown 222c4762a1bSJed Brown PetscFunctionBeginUser; 2235f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 225c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* Get pointers to vector data */ 2285f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X,&x)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 230c4762a1bSJed Brown 231c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 232c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 233c4762a1bSJed Brown PetscScalar u = x[i].u,v = x[i].v; 234c4762a1bSJed Brown f[i].u = hx*(user->A + u*u*v - (user->B+1)*u); 235c4762a1bSJed Brown f[i].v = hx*(user->B*u - u*u*v); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* Restore vectors */ 2395f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X,&x)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 241c4762a1bSJed Brown PetscFunctionReturn(0); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 245c4762a1bSJed Brown /* 246c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 247c4762a1bSJed Brown */ 248c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 249c4762a1bSJed Brown { 250c4762a1bSJed Brown User user = (User)ptr; 251c4762a1bSJed Brown DMDALocalInfo info; 252c4762a1bSJed Brown PetscInt i; 253c4762a1bSJed Brown PetscReal hx; 254c4762a1bSJed Brown DM da; 255c4762a1bSJed Brown Field *x,*xdot; 256c4762a1bSJed Brown 257c4762a1bSJed Brown PetscFunctionBeginUser; 2585f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 260c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* Get pointers to vector data */ 2635f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,X,&x)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Xdot,&xdot)); 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 267c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 268c4762a1bSJed Brown if (i == 0 || i == info.mx-1) { 269c4762a1bSJed Brown const PetscInt row = i,col = i; 270c4762a1bSJed Brown const PetscScalar vals[2][2] = {{hx,0},{0,hx}}; 2715f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES)); 272c4762a1bSJed Brown } else { 273c4762a1bSJed Brown const PetscInt row = i,col[] = {i-1,i,i+1}; 274c4762a1bSJed Brown const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx; 275c4762a1bSJed Brown const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}}, 276c4762a1bSJed Brown {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}}; 2775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES)); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* Restore vectors */ 2825f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,X,&x)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Xdot,&xdot)); 284c4762a1bSJed Brown 2855f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 287c4762a1bSJed Brown if (J != Jpre) { 2885f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown PetscFunctionReturn(0); 292c4762a1bSJed Brown } 293c4762a1bSJed Brown 294c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 295c4762a1bSJed Brown { 296c4762a1bSJed Brown User user = (User)ctx; 297c4762a1bSJed Brown DM da; 298c4762a1bSJed Brown PetscInt i; 299c4762a1bSJed Brown DMDALocalInfo info; 300c4762a1bSJed Brown Field *x; 301c4762a1bSJed Brown PetscReal hx; 302c4762a1bSJed Brown 303c4762a1bSJed Brown PetscFunctionBeginUser; 3045f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 306c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* Get pointers to vector data */ 3095f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,X,&x)); 310c4762a1bSJed Brown 311c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 312c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 313c4762a1bSJed Brown PetscReal xi = i*hx; 314c4762a1bSJed Brown x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi); 315c4762a1bSJed Brown x[i].v = user->vleft*(1.-xi) + user->vright*xi; 316c4762a1bSJed Brown } 3175f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,X,&x)); 318c4762a1bSJed Brown PetscFunctionReturn(0); 319c4762a1bSJed Brown } 320c4762a1bSJed Brown 321c4762a1bSJed Brown /*TEST 322c4762a1bSJed Brown 323c4762a1bSJed Brown test: 324c4762a1bSJed Brown args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3 325c4762a1bSJed Brown 326c4762a1bSJed Brown test: 327c4762a1bSJed Brown suffix: 2 328c4762a1bSJed Brown args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3 329c4762a1bSJed Brown 330c4762a1bSJed Brown TEST*/ 331