1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d. Modified from ex13.c for illustrating how to solve DAEs. \n"; 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown u_t = uxx + uyy 5*c4762a1bSJed Brown 0 < x < 1, 0 < y < 1; 6*c4762a1bSJed Brown At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125 7*c4762a1bSJed Brown u(x,y) = 0.0 if r >= .125 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown Boundary conditions: 11*c4762a1bSJed Brown Drichlet BC: 12*c4762a1bSJed Brown At x=0, x=1, y=0, y=1: u = 0.0 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown Neumann BC: 15*c4762a1bSJed Brown At x=0, x=1: du(x,y,t)/dx = 0 16*c4762a1bSJed Brown At y=0, y=1: du(x,y,t)/dy = 0 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor 19*c4762a1bSJed Brown ./ex15 -da_grid_x 40 -da_grid_y 40 -draw_pause .1 -boundary 1 -ts_monitor_draw_solution 20*c4762a1bSJed Brown ./ex15 -da_grid_x 40 -da_grid_y 40 -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9 21*c4762a1bSJed Brown 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 /* 29*c4762a1bSJed Brown User-defined data structures and routines 30*c4762a1bSJed Brown */ 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown /* AppCtx: used by FormIFunction() and FormIJacobian() */ 33*c4762a1bSJed Brown typedef struct { 34*c4762a1bSJed Brown DM da; 35*c4762a1bSJed Brown PetscInt nstencilpts; /* number of stencil points: 5 or 9 */ 36*c4762a1bSJed Brown PetscReal c; 37*c4762a1bSJed Brown PetscInt boundary; /* Type of boundary condition */ 38*c4762a1bSJed Brown PetscBool viewJacobian; 39*c4762a1bSJed Brown } AppCtx; 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown extern PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 42*c4762a1bSJed Brown extern PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 43*c4762a1bSJed Brown extern PetscErrorCode FormInitialSolution(Vec,void*); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown int main(int argc,char **argv) 46*c4762a1bSJed Brown { 47*c4762a1bSJed Brown TS ts; /* nonlinear solver */ 48*c4762a1bSJed Brown Vec u,r; /* solution, residual vectors */ 49*c4762a1bSJed Brown Mat J,Jmf = NULL; /* Jacobian matrices */ 50*c4762a1bSJed Brown PetscErrorCode ierr; 51*c4762a1bSJed Brown DM da; 52*c4762a1bSJed Brown PetscReal dt; 53*c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 54*c4762a1bSJed Brown SNES snes; 55*c4762a1bSJed Brown PetscInt Jtype; /* Jacobian type 56*c4762a1bSJed Brown 0: user provide Jacobian; 57*c4762a1bSJed Brown 1: slow finite difference; 58*c4762a1bSJed Brown 2: fd with coloring; */ 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 61*c4762a1bSJed Brown /* Initialize user application context */ 62*c4762a1bSJed Brown user.da = NULL; 63*c4762a1bSJed Brown user.nstencilpts = 5; 64*c4762a1bSJed Brown user.c = -30.0; 65*c4762a1bSJed Brown user.boundary = 0; /* 0: Drichlet BC; 1: Neumann BC */ 66*c4762a1bSJed Brown user.viewJacobian = PETSC_FALSE; 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nstencilpts",&user.nstencilpts,NULL);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian);CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 73*c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 74*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 75*c4762a1bSJed Brown if (user.nstencilpts == 5) { 76*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr); 77*c4762a1bSJed Brown } else if (user.nstencilpts == 9) { 78*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr); 79*c4762a1bSJed Brown } else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"nstencilpts %d is not supported",user.nstencilpts); 80*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 82*c4762a1bSJed Brown user.da = da; 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 85*c4762a1bSJed Brown Extract global vectors from DMDA; 86*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 87*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&u);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = VecDuplicate(u,&r);CHKERRQ(ierr); 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91*c4762a1bSJed Brown Create timestepping solver context 92*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 93*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = TSSetIFunction(ts,r,FormIFunction,&user);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 102*c4762a1bSJed Brown Set initial conditions 103*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 104*c4762a1bSJed Brown ierr = FormInitialSolution(u,&user);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = TSSetSolution(ts,u);CHKERRQ(ierr); 106*c4762a1bSJed Brown dt = .01; 107*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 108*c4762a1bSJed Brown 109*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 110*c4762a1bSJed Brown Set Jacobian evaluation routine 111*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 112*c4762a1bSJed Brown ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr); 114*c4762a1bSJed Brown Jtype = 0; 115*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL, "-Jtype",&Jtype,NULL);CHKERRQ(ierr); 116*c4762a1bSJed Brown if (Jtype == 0) { /* use user provided Jacobian evaluation routine */ 117*c4762a1bSJed Brown if (user.nstencilpts != 5) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"user Jacobian routine FormIJacobian() does not support nstencilpts=%D",user.nstencilpts); 118*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr); 119*c4762a1bSJed Brown } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */ 120*c4762a1bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = MatCreateSNESMF(snes,&Jmf);CHKERRQ(ierr); 122*c4762a1bSJed Brown if (Jtype == 1) { /* slow finite difference J; */ 123*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefault,NULL);CHKERRQ(ierr); 124*c4762a1bSJed Brown } else if (Jtype == 2) { /* Use coloring to compute finite difference J efficiently */ 125*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefaultColor,0);CHKERRQ(ierr); 126*c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Jtype is not supported"); 127*c4762a1bSJed Brown } 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130*c4762a1bSJed Brown Sets various TS parameters from user options 131*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 132*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135*c4762a1bSJed Brown Solve nonlinear system 136*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137*c4762a1bSJed Brown ierr = TSSolve(ts,u);CHKERRQ(ierr); 138*c4762a1bSJed Brown 139*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140*c4762a1bSJed Brown Free work space. 141*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = MatDestroy(&Jmf);CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown ierr = PetscFinalize(); 150*c4762a1bSJed Brown return ierr; 151*c4762a1bSJed Brown } 152*c4762a1bSJed Brown 153*c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 154*c4762a1bSJed Brown /* 155*c4762a1bSJed Brown FormIFunction = Udot - RHSFunction 156*c4762a1bSJed Brown */ 157*c4762a1bSJed Brown PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 158*c4762a1bSJed Brown { 159*c4762a1bSJed Brown PetscErrorCode ierr; 160*c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 161*c4762a1bSJed Brown DM da = (DM)user->da; 162*c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 163*c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 164*c4762a1bSJed Brown PetscScalar u,uxx,uyy,**uarray,**f,**udot; 165*c4762a1bSJed Brown Vec localU; 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown PetscFunctionBeginUser; 168*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 170*c4762a1bSJed Brown 171*c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx); 172*c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy); 173*c4762a1bSJed Brown if (user->nstencilpts == 9 && hx != hy) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"hx must equal hy when nstencilpts = 9 for this example"); 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown /* 176*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 177*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 178*c4762a1bSJed Brown By placing code between these two statements, computations can be 179*c4762a1bSJed Brown done while messages are in transition. 180*c4762a1bSJed Brown */ 181*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 182*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 183*c4762a1bSJed Brown 184*c4762a1bSJed Brown /* Get pointers to vector data */ 185*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&uarray);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,Udot,&udot);CHKERRQ(ierr); 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /* Get local grid boundaries */ 190*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 193*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 194*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 195*c4762a1bSJed Brown /* Boundary conditions */ 196*c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx-1 || j == My-1) { 197*c4762a1bSJed Brown if (user->boundary == 0) { /* Drichlet BC */ 198*c4762a1bSJed Brown f[j][i] = uarray[j][i]; /* F = U */ 199*c4762a1bSJed Brown } else { /* Neumann BC */ 200*c4762a1bSJed Brown if (i == 0 && j == 0) { /* SW corner */ 201*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j+1][i+1]; 202*c4762a1bSJed Brown } else if (i == Mx-1 && j == 0) { /* SE corner */ 203*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j+1][i-1]; 204*c4762a1bSJed Brown } else if (i == 0 && j == My-1) { /* NW corner */ 205*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j-1][i+1]; 206*c4762a1bSJed Brown } else if (i == Mx-1 && j == My-1) { /* NE corner */ 207*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j-1][i-1]; 208*c4762a1bSJed Brown } else if (i == 0) { /* Left */ 209*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i+1]; 210*c4762a1bSJed Brown } else if (i == Mx-1) { /* Right */ 211*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i-1]; 212*c4762a1bSJed Brown } else if (j == 0) { /* Bottom */ 213*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j+1][i]; 214*c4762a1bSJed Brown } else if (j == My-1) { /* Top */ 215*c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j-1][i]; 216*c4762a1bSJed Brown } 217*c4762a1bSJed Brown } 218*c4762a1bSJed Brown } else { /* Interior */ 219*c4762a1bSJed Brown u = uarray[j][i]; 220*c4762a1bSJed Brown /* 5-point stencil */ 221*c4762a1bSJed Brown uxx = (-2.0*u + uarray[j][i-1] + uarray[j][i+1]); 222*c4762a1bSJed Brown uyy = (-2.0*u + uarray[j-1][i] + uarray[j+1][i]); 223*c4762a1bSJed Brown if (user->nstencilpts == 9) { 224*c4762a1bSJed Brown /* 9-point stencil: assume hx=hy */ 225*c4762a1bSJed Brown uxx = 2.0*uxx/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0; 226*c4762a1bSJed Brown uyy = 2.0*uyy/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0; 227*c4762a1bSJed Brown } 228*c4762a1bSJed Brown f[j][i] = udot[j][i] - (uxx*sx + uyy*sy); 229*c4762a1bSJed Brown } 230*c4762a1bSJed Brown } 231*c4762a1bSJed Brown } 232*c4762a1bSJed Brown 233*c4762a1bSJed Brown /* Restore vectors */ 234*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&uarray);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,Udot,&udot);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = PetscLogFlops(11.0*ym*xm);CHKERRQ(ierr); 239*c4762a1bSJed Brown PetscFunctionReturn(0); 240*c4762a1bSJed Brown } 241*c4762a1bSJed Brown 242*c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 243*c4762a1bSJed Brown /* 244*c4762a1bSJed Brown FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot 245*c4762a1bSJed Brown This routine is not used with option '-use_coloring' 246*c4762a1bSJed Brown */ 247*c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx) 248*c4762a1bSJed Brown { 249*c4762a1bSJed Brown PetscErrorCode ierr; 250*c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym,nc; 251*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ctx; 252*c4762a1bSJed Brown DM da = (DM)user->da; 253*c4762a1bSJed Brown MatStencil col[5],row; 254*c4762a1bSJed Brown PetscScalar vals[5],hx,hy,sx,sy; 255*c4762a1bSJed Brown 256*c4762a1bSJed Brown PetscFunctionBeginUser; 257*c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 258*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx); 261*c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy); 262*c4762a1bSJed Brown 263*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 264*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 265*c4762a1bSJed Brown nc = 0; 266*c4762a1bSJed Brown row.j = j; row.i = i; 267*c4762a1bSJed Brown if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) { 268*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0; 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown } else if (user->boundary > 0 && i == 0) { /* Left Neumann */ 271*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0; 272*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i+1; vals[nc++] = -1.0; 273*c4762a1bSJed Brown } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */ 274*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0; 275*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i-1; vals[nc++] = -1.0; 276*c4762a1bSJed Brown } else if (user->boundary > 0 && j == 0) { /* Bottom Neumann */ 277*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0; 278*c4762a1bSJed Brown col[nc].j = j+1; col[nc].i = i; vals[nc++] = -1.0; 279*c4762a1bSJed Brown } else if (user->boundary > 0 && j == My-1) { /* Top Neumann */ 280*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0; 281*c4762a1bSJed Brown col[nc].j = j-1; col[nc].i = i; vals[nc++] = -1.0; 282*c4762a1bSJed Brown } else { /* Interior */ 283*c4762a1bSJed Brown col[nc].j = j-1; col[nc].i = i; vals[nc++] = -sy; 284*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i-1; vals[nc++] = -sx; 285*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i; vals[nc++] = 2.0*(sx + sy) + a; 286*c4762a1bSJed Brown col[nc].j = j; col[nc].i = i+1; vals[nc++] = -sx; 287*c4762a1bSJed Brown col[nc].j = j+1; col[nc].i = i; vals[nc++] = -sy; 288*c4762a1bSJed Brown } 289*c4762a1bSJed Brown ierr = MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES);CHKERRQ(ierr); 290*c4762a1bSJed Brown } 291*c4762a1bSJed Brown } 292*c4762a1bSJed Brown ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 293*c4762a1bSJed Brown ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 294*c4762a1bSJed Brown if (J != Jpre) { 295*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 296*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 297*c4762a1bSJed Brown } 298*c4762a1bSJed Brown 299*c4762a1bSJed Brown if (user->viewJacobian) { 300*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)Jpre),"Jpre:\n");CHKERRQ(ierr); 301*c4762a1bSJed Brown ierr = MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 302*c4762a1bSJed Brown } 303*c4762a1bSJed Brown PetscFunctionReturn(0); 304*c4762a1bSJed Brown } 305*c4762a1bSJed Brown 306*c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 307*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(Vec U,void *ptr) 308*c4762a1bSJed Brown { 309*c4762a1bSJed Brown AppCtx *user=(AppCtx*)ptr; 310*c4762a1bSJed Brown DM da =user->da; 311*c4762a1bSJed Brown PetscReal c =user->c; 312*c4762a1bSJed Brown PetscErrorCode ierr; 313*c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 314*c4762a1bSJed Brown PetscScalar **u; 315*c4762a1bSJed Brown PetscReal hx,hy,x,y,r; 316*c4762a1bSJed Brown 317*c4762a1bSJed Brown PetscFunctionBeginUser; 318*c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 319*c4762a1bSJed Brown 320*c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 321*c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 322*c4762a1bSJed Brown 323*c4762a1bSJed Brown /* Get pointers to vector data */ 324*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 325*c4762a1bSJed Brown 326*c4762a1bSJed Brown /* Get local grid boundaries */ 327*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 328*c4762a1bSJed Brown 329*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 330*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 331*c4762a1bSJed Brown y = j*hy; 332*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 333*c4762a1bSJed Brown x = i*hx; 334*c4762a1bSJed Brown r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)); 335*c4762a1bSJed Brown if (r < .125) u[j][i] = PetscExpReal(c*r*r*r); 336*c4762a1bSJed Brown else u[j][i] = 0.0; 337*c4762a1bSJed Brown } 338*c4762a1bSJed Brown } 339*c4762a1bSJed Brown 340*c4762a1bSJed Brown /* Restore vectors */ 341*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 342*c4762a1bSJed Brown PetscFunctionReturn(0); 343*c4762a1bSJed Brown } 344*c4762a1bSJed Brown 345*c4762a1bSJed Brown /*TEST 346*c4762a1bSJed Brown 347*c4762a1bSJed Brown test: 348*c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor 349*c4762a1bSJed Brown 350*c4762a1bSJed Brown test: 351*c4762a1bSJed Brown suffix: 2 352*c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown test: 355*c4762a1bSJed Brown suffix: 3 356*c4762a1bSJed Brown requires: !single 357*c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor 358*c4762a1bSJed Brown 359*c4762a1bSJed Brown test: 360*c4762a1bSJed Brown suffix: 4 361*c4762a1bSJed Brown requires: !single 362*c4762a1bSJed Brown nsize: 2 363*c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor 364*c4762a1bSJed Brown 365*c4762a1bSJed Brown test: 366*c4762a1bSJed Brown suffix: 5 367*c4762a1bSJed Brown nsize: 1 368*c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor 369*c4762a1bSJed Brown 370*c4762a1bSJed Brown 371*c4762a1bSJed Brown TEST*/ 372*c4762a1bSJed Brown 373