1*c4762a1bSJed Brown static const char help[] = "Call PetscInitialize multiple times.\n"; 2*c4762a1bSJed Brown /* 3*c4762a1bSJed Brown This example is based on the Brusselator tutorial of the same name, but tests multiple calls to PetscInitialize(). 4*c4762a1bSJed Brown This is a bad "convergence study" because it only compares min and max values of the solution rather than comparing 5*c4762a1bSJed Brown norms of the errors. For convergence studies, we recommend invoking PetscInitialize() only once and comparing norms 6*c4762a1bSJed Brown of errors (perhaps estimated using an accurate reference solution). 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods and multiple solves. 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown u_t - alpha u_xx = A + u^2 v - (B+1) u 11*c4762a1bSJed Brown v_t - alpha v_xx = B u - u^2 v 12*c4762a1bSJed Brown 0 < x < 1; 13*c4762a1bSJed Brown A = 1, B = 3, alpha = 1/10 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown Initial conditions: 16*c4762a1bSJed Brown u(x,0) = 1 + sin(2 pi x) 17*c4762a1bSJed Brown v(x,0) = 3 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown Boundary conditions: 20*c4762a1bSJed Brown u(0,t) = u(1,t) = 1 21*c4762a1bSJed Brown v(0,t) = v(1,t) = 3 22*c4762a1bSJed Brown */ 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown #include <petscdm.h> 25*c4762a1bSJed Brown #include <petscdmda.h> 26*c4762a1bSJed Brown #include <petscts.h> 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown typedef struct { 29*c4762a1bSJed Brown PetscScalar u,v; 30*c4762a1bSJed Brown } Field; 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown typedef struct _User *User; 33*c4762a1bSJed Brown struct _User { 34*c4762a1bSJed Brown PetscReal A,B; /* Reaction coefficients */ 35*c4762a1bSJed Brown PetscReal alpha; /* Diffusion coefficient */ 36*c4762a1bSJed Brown PetscReal uleft,uright; /* Dirichlet boundary conditions */ 37*c4762a1bSJed Brown PetscReal vleft,vright; /* Dirichlet boundary conditions */ 38*c4762a1bSJed Brown }; 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 41*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 42*c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 43*c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 44*c4762a1bSJed Brown static int Brusselator(int,char**,PetscInt); 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown int main(int argc,char **argv) 47*c4762a1bSJed Brown { 48*c4762a1bSJed Brown PetscInt cycle; 49*c4762a1bSJed Brown PetscErrorCode ierr; 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown ierr = MPI_Init(&argc,&argv);if (ierr) return ierr; 52*c4762a1bSJed Brown for (cycle=0; cycle<4; cycle++) { 53*c4762a1bSJed Brown ierr = Brusselator(argc,argv,cycle); 54*c4762a1bSJed Brown if (ierr) return 1; 55*c4762a1bSJed Brown } 56*c4762a1bSJed Brown ierr = MPI_Finalize(); 57*c4762a1bSJed Brown return ierr; 58*c4762a1bSJed Brown } 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown PetscErrorCode Brusselator(int argc,char **argv,PetscInt cycle) 61*c4762a1bSJed Brown { 62*c4762a1bSJed Brown TS ts; /* nonlinear solver */ 63*c4762a1bSJed Brown Vec X; /* solution, residual vectors */ 64*c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 65*c4762a1bSJed Brown PetscInt steps,mx; 66*c4762a1bSJed Brown PetscErrorCode ierr; 67*c4762a1bSJed Brown DM da; 68*c4762a1bSJed Brown PetscReal ftime,hx,dt,xmax,xmin; 69*c4762a1bSJed Brown struct _User user; /* user-defined work context */ 70*c4762a1bSJed Brown TSConvergedReason reason; 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75*c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 76*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 77*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 82*c4762a1bSJed Brown Extract global vectors from DMDA; 83*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 84*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&X);CHKERRQ(ierr); 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* Initialize user application context */ 87*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options",""); 88*c4762a1bSJed Brown { 89*c4762a1bSJed Brown user.A = 1; 90*c4762a1bSJed Brown user.B = 3; 91*c4762a1bSJed Brown user.alpha = 0.1; 92*c4762a1bSJed Brown user.uleft = 1; 93*c4762a1bSJed Brown user.uright = 1; 94*c4762a1bSJed Brown user.vleft = 3; 95*c4762a1bSJed Brown user.vright = 3; 96*c4762a1bSJed Brown ierr = PetscOptionsReal("-A","Reaction rate","",user.A,&user.A,NULL);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = PetscOptionsReal("-B","Reaction rate","",user.B,&user.B,NULL);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = PetscOptionsReal("-alpha","Diffusion coefficient","",user.alpha,&user.alpha,NULL);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = PetscOptionsReal("-uleft","Dirichlet boundary condition","",user.uleft,&user.uleft,NULL);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = PetscOptionsReal("-uright","Dirichlet boundary condition","",user.uright,&user.uright,NULL);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = PetscOptionsReal("-vleft","Dirichlet boundary condition","",user.vleft,&user.vleft,NULL);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = PetscOptionsReal("-vright","Dirichlet boundary condition","",user.vright,&user.vright,NULL);CHKERRQ(ierr); 103*c4762a1bSJed Brown } 104*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107*c4762a1bSJed Brown Create timestepping solver context 108*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr); 117*c4762a1bSJed Brown 118*c4762a1bSJed Brown ftime = 1.0; 119*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123*c4762a1bSJed Brown Set initial conditions 124*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125*c4762a1bSJed Brown ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = TSSetSolution(ts,X);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = VecGetSize(X,&mx);CHKERRQ(ierr); 128*c4762a1bSJed Brown hx = 1.0/(PetscReal)(mx/2-1); 129*c4762a1bSJed Brown dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */ 130*c4762a1bSJed Brown dt *= PetscPowRealInt(0.2,cycle); /* Shrink the time step in convergence study. */ 131*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = TSSetTolerances(ts,1e-3*PetscPowRealInt(0.5,cycle),NULL,1e-3*PetscPowRealInt(0.5,cycle),NULL);CHKERRQ(ierr); 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135*c4762a1bSJed Brown Set runtime options 136*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 138*c4762a1bSJed Brown 139*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140*c4762a1bSJed Brown Solve nonlinear system 141*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142*c4762a1bSJed Brown ierr = TSSolve(ts,X);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = VecMin(X,NULL,&xmin);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = VecMax(X,NULL,&xmax);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151*c4762a1bSJed Brown Free work space. 152*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = PetscFinalize(); 158*c4762a1bSJed Brown return ierr; 159*c4762a1bSJed Brown } 160*c4762a1bSJed Brown 161*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 162*c4762a1bSJed Brown { 163*c4762a1bSJed Brown User user = (User)ptr; 164*c4762a1bSJed Brown DM da; 165*c4762a1bSJed Brown DMDALocalInfo info; 166*c4762a1bSJed Brown PetscInt i; 167*c4762a1bSJed Brown Field *x,*xdot,*f; 168*c4762a1bSJed Brown PetscReal hx; 169*c4762a1bSJed Brown Vec Xloc; 170*c4762a1bSJed Brown PetscErrorCode ierr; 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown PetscFunctionBeginUser; 173*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 175*c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown /* 178*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 179*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 180*c4762a1bSJed Brown By placing code between these two statements, computations can be 181*c4762a1bSJed Brown done while messages are in transition. 182*c4762a1bSJed Brown */ 183*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&Xloc);CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 186*c4762a1bSJed Brown 187*c4762a1bSJed Brown /* Get pointers to vector data */ 188*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Xloc,&x);CHKERRQ(ierr); 189*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Xdot,&xdot);CHKERRQ(ierr); 190*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 193*c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 194*c4762a1bSJed Brown if (i == 0) { 195*c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->uleft); 196*c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->vleft); 197*c4762a1bSJed Brown } else if (i == info.mx-1) { 198*c4762a1bSJed Brown f[i].u = hx * (x[i].u - user->uright); 199*c4762a1bSJed Brown f[i].v = hx * (x[i].v - user->vright); 200*c4762a1bSJed Brown } else { 201*c4762a1bSJed Brown f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx; 202*c4762a1bSJed Brown f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx; 203*c4762a1bSJed Brown } 204*c4762a1bSJed Brown } 205*c4762a1bSJed Brown 206*c4762a1bSJed Brown /* Restore vectors */ 207*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Xloc,&x);CHKERRQ(ierr); 208*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Xdot,&xdot);CHKERRQ(ierr); 209*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 210*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&Xloc);CHKERRQ(ierr); 211*c4762a1bSJed Brown PetscFunctionReturn(0); 212*c4762a1bSJed Brown } 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 215*c4762a1bSJed Brown { 216*c4762a1bSJed Brown User user = (User)ptr; 217*c4762a1bSJed Brown DM da; 218*c4762a1bSJed Brown DMDALocalInfo info; 219*c4762a1bSJed Brown PetscInt i; 220*c4762a1bSJed Brown PetscReal hx; 221*c4762a1bSJed Brown Field *x,*f; 222*c4762a1bSJed Brown PetscErrorCode ierr; 223*c4762a1bSJed Brown 224*c4762a1bSJed Brown PetscFunctionBeginUser; 225*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 227*c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 228*c4762a1bSJed Brown 229*c4762a1bSJed Brown /* Get pointers to vector data */ 230*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,X,&x);CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 232*c4762a1bSJed Brown 233*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 234*c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 235*c4762a1bSJed Brown PetscScalar u = x[i].u,v = x[i].v; 236*c4762a1bSJed Brown f[i].u = hx*(user->A + u*u*v - (user->B+1)*u); 237*c4762a1bSJed Brown f[i].v = hx*(user->B*u - u*u*v); 238*c4762a1bSJed Brown } 239*c4762a1bSJed Brown 240*c4762a1bSJed Brown /* Restore vectors */ 241*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,X,&x);CHKERRQ(ierr); 242*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 243*c4762a1bSJed Brown PetscFunctionReturn(0); 244*c4762a1bSJed Brown } 245*c4762a1bSJed Brown 246*c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 247*c4762a1bSJed Brown /* 248*c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 249*c4762a1bSJed Brown */ 250*c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr) 251*c4762a1bSJed Brown { 252*c4762a1bSJed Brown User user = (User)ptr; 253*c4762a1bSJed Brown PetscErrorCode ierr; 254*c4762a1bSJed Brown DMDALocalInfo info; 255*c4762a1bSJed Brown PetscInt i; 256*c4762a1bSJed Brown PetscReal hx; 257*c4762a1bSJed Brown DM da; 258*c4762a1bSJed Brown Field *x,*xdot; 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown PetscFunctionBeginUser; 261*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 263*c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 264*c4762a1bSJed Brown 265*c4762a1bSJed Brown /* Get pointers to vector data */ 266*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,X,&x);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Xdot,&xdot);CHKERRQ(ierr); 268*c4762a1bSJed Brown 269*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 270*c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 271*c4762a1bSJed Brown if (i == 0 || i == info.mx-1) { 272*c4762a1bSJed Brown const PetscInt row = i,col = i; 273*c4762a1bSJed Brown const PetscScalar vals[2][2] = {{hx,0},{0,hx}}; 274*c4762a1bSJed Brown ierr = MatSetValuesBlocked(Jpre,1,&row,1,&col,&vals[0][0],INSERT_VALUES);CHKERRQ(ierr); 275*c4762a1bSJed Brown } else { 276*c4762a1bSJed Brown const PetscInt row = i,col[] = {i-1,i,i+1}; 277*c4762a1bSJed Brown const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx; 278*c4762a1bSJed Brown const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}}, 279*c4762a1bSJed Brown {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}}; 280*c4762a1bSJed Brown ierr = MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);CHKERRQ(ierr); 281*c4762a1bSJed Brown } 282*c4762a1bSJed Brown } 283*c4762a1bSJed Brown 284*c4762a1bSJed Brown /* Restore vectors */ 285*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,X,&x);CHKERRQ(ierr); 286*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Xdot,&xdot);CHKERRQ(ierr); 287*c4762a1bSJed Brown 288*c4762a1bSJed Brown ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 289*c4762a1bSJed Brown ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 290*c4762a1bSJed Brown if (J != Jpre) { 291*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 292*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 293*c4762a1bSJed Brown } 294*c4762a1bSJed Brown PetscFunctionReturn(0); 295*c4762a1bSJed Brown } 296*c4762a1bSJed Brown 297*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 298*c4762a1bSJed Brown { 299*c4762a1bSJed Brown User user = (User)ctx; 300*c4762a1bSJed Brown DM da; 301*c4762a1bSJed Brown PetscInt i; 302*c4762a1bSJed Brown DMDALocalInfo info; 303*c4762a1bSJed Brown Field *x; 304*c4762a1bSJed Brown PetscReal hx; 305*c4762a1bSJed Brown PetscErrorCode ierr; 306*c4762a1bSJed Brown 307*c4762a1bSJed Brown PetscFunctionBeginUser; 308*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 309*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 310*c4762a1bSJed Brown hx = 1.0/(PetscReal)(info.mx-1); 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown /* Get pointers to vector data */ 313*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr); 314*c4762a1bSJed Brown 315*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 316*c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 317*c4762a1bSJed Brown PetscReal xi = i*hx; 318*c4762a1bSJed Brown x[i].u = user->uleft*(1.-xi) + user->uright*xi + PetscSinReal(2.*PETSC_PI*xi); 319*c4762a1bSJed Brown x[i].v = user->vleft*(1.-xi) + user->vright*xi; 320*c4762a1bSJed Brown } 321*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr); 322*c4762a1bSJed Brown PetscFunctionReturn(0); 323*c4762a1bSJed Brown } 324*c4762a1bSJed Brown 325*c4762a1bSJed Brown /*TEST 326*c4762a1bSJed Brown 327*c4762a1bSJed Brown build: 328*c4762a1bSJed Brown requires: c99 329*c4762a1bSJed Brown 330*c4762a1bSJed Brown test: 331*c4762a1bSJed Brown args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3 332*c4762a1bSJed Brown 333*c4762a1bSJed Brown test: 334*c4762a1bSJed Brown suffix: 2 335*c4762a1bSJed Brown args: -ts_exact_final_time INTERPOLATE -snes_rtol 1.e-3 336*c4762a1bSJed Brown 337*c4762a1bSJed Brown TEST*/ 338*c4762a1bSJed Brown 339