1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 3d.\n\ 3*c4762a1bSJed Brown We solve the Bratu (SFI - solid fuel ignition) problem in a 3D rectangular\n\ 4*c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\ 5*c4762a1bSJed Brown The command line options include:\n\ 6*c4762a1bSJed Brown -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\ 7*c4762a1bSJed Brown problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n"; 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown /*T 10*c4762a1bSJed Brown Concepts: SNES^parallel Bratu example 11*c4762a1bSJed Brown Concepts: DMDA^using distributed arrays; 12*c4762a1bSJed Brown Processors: n 13*c4762a1bSJed Brown T*/ 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown Solid Fuel Ignition (SFI) problem. This problem is modeled by 20*c4762a1bSJed Brown the partial differential equation 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown with boundary conditions 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown u = 0 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown A finite difference approximation with the usual 7-point stencil 29*c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a nonlinear 30*c4762a1bSJed Brown system of equations. 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ------------------------------------------------------------------------- */ 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown /* 36*c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 37*c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 38*c4762a1bSJed Brown file automatically includes: 39*c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 40*c4762a1bSJed Brown petscmat.h - matrices 41*c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 42*c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 43*c4762a1bSJed Brown petscksp.h - linear solvers 44*c4762a1bSJed Brown */ 45*c4762a1bSJed Brown #include <petscdm.h> 46*c4762a1bSJed Brown #include <petscdmda.h> 47*c4762a1bSJed Brown #include <petscsnes.h> 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown /* 51*c4762a1bSJed Brown User-defined application context - contains data needed by the 52*c4762a1bSJed Brown application-provided call-back routines, FormJacobian() and 53*c4762a1bSJed Brown FormFunction(). 54*c4762a1bSJed Brown */ 55*c4762a1bSJed Brown typedef struct { 56*c4762a1bSJed Brown PetscReal param; /* test problem parameter */ 57*c4762a1bSJed Brown DM da; /* distributed array data structure */ 58*c4762a1bSJed Brown } AppCtx; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown /* 61*c4762a1bSJed Brown User-defined routines 62*c4762a1bSJed Brown */ 63*c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(SNES,Vec,Vec,void*); 64*c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 65*c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx*,Vec); 66*c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown int main(int argc,char **argv) 69*c4762a1bSJed Brown { 70*c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 71*c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 72*c4762a1bSJed Brown Mat J = NULL; /* Jacobian matrix */ 73*c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 74*c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 75*c4762a1bSJed Brown MatFDColoring matfdcoloring = NULL; 76*c4762a1bSJed Brown PetscBool matrix_free = PETSC_FALSE,coloring = PETSC_FALSE, coloring_ds = PETSC_FALSE,local_coloring = PETSC_FALSE; 77*c4762a1bSJed Brown PetscErrorCode ierr; 78*c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.,fnorm; 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 81*c4762a1bSJed Brown Initialize program 82*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 87*c4762a1bSJed Brown Initialize problem parameters 88*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 89*c4762a1bSJed Brown user.param = 6.0; 90*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);CHKERRQ(ierr); 91*c4762a1bSJed Brown if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda is out of range"); 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 94*c4762a1bSJed Brown Create nonlinear solver context 95*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 96*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99*c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 100*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 101*c4762a1bSJed Brown ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,4,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&user.da);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = DMSetFromOptions(user.da);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = DMSetUp(user.da);CHKERRQ(ierr); 104*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105*c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 106*c4762a1bSJed Brown vectors that are the same types 107*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108*c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.da,&x);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 112*c4762a1bSJed Brown Set function evaluation routine and vector 113*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 114*c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction,(void*)&user);CHKERRQ(ierr); 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117*c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 120*c4762a1bSJed Brown routine. User can override with: 121*c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 122*c4762a1bSJed Brown (unless user explicitly sets preconditioner) 123*c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 124*c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 125*c4762a1bSJed Brown products within Newton-Krylov method 126*c4762a1bSJed Brown -fdcoloring : using finite differences with coloring to compute the Jacobian 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown Note one can use -matfd_coloring wp or ds the only reason for the -fdcoloring_ds option 129*c4762a1bSJed Brown below is to test the call to MatFDColoringSetType(). 130*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-fdcoloring",&coloring,NULL);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-fdcoloring_ds",&coloring_ds,NULL);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-fdcoloring_local",&local_coloring,NULL);CHKERRQ(ierr); 135*c4762a1bSJed Brown if (!matrix_free) { 136*c4762a1bSJed Brown ierr = DMSetMatType(user.da,MATAIJ);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = DMCreateMatrix(user.da,&J);CHKERRQ(ierr); 138*c4762a1bSJed Brown if (coloring) { 139*c4762a1bSJed Brown ISColoring iscoloring; 140*c4762a1bSJed Brown if (!local_coloring) { 141*c4762a1bSJed Brown ierr = DMCreateColoring(user.da,IS_COLORING_GLOBAL,&iscoloring);CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = MatFDColoringCreate(J,iscoloring,&matfdcoloring);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunction,&user);CHKERRQ(ierr); 144*c4762a1bSJed Brown } else { 145*c4762a1bSJed Brown ierr = DMCreateColoring(user.da,IS_COLORING_LOCAL,&iscoloring);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = MatFDColoringCreate(J,iscoloring,&matfdcoloring);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = MatFDColoringUseDM(J,matfdcoloring);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunctionLocal,&user);CHKERRQ(ierr); 149*c4762a1bSJed Brown } 150*c4762a1bSJed Brown if (coloring_ds) { 151*c4762a1bSJed Brown ierr = MatFDColoringSetType(matfdcoloring,MATMFFD_DS);CHKERRQ(ierr); 152*c4762a1bSJed Brown } 153*c4762a1bSJed Brown ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = MatFDColoringSetUp(J,iscoloring,matfdcoloring);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 157*c4762a1bSJed Brown } else { 158*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian,&user);CHKERRQ(ierr); 159*c4762a1bSJed Brown } 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 163*c4762a1bSJed Brown Customize nonlinear solver; set runtime options 164*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 165*c4762a1bSJed Brown ierr = SNESSetDM(snes,user.da);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 167*c4762a1bSJed Brown 168*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 169*c4762a1bSJed Brown Evaluate initial guess 170*c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 171*c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 172*c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 173*c4762a1bSJed Brown this vector to zero by calling VecSet(). 174*c4762a1bSJed Brown */ 175*c4762a1bSJed Brown ierr = FormInitialGuess(&user,x);CHKERRQ(ierr); 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 178*c4762a1bSJed Brown Solve nonlinear system 179*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 180*c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 181*c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 182*c4762a1bSJed Brown 183*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184*c4762a1bSJed Brown Explicitly check norm of the residual of the solution 185*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 186*c4762a1bSJed Brown ierr = FormFunction(snes,x,r,(void*)&user);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = VecNorm(r,NORM_2,&fnorm);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D fnorm %g\n",its,(double)fnorm);CHKERRQ(ierr); 189*c4762a1bSJed Brown 190*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 192*c4762a1bSJed Brown are no longer needed. 193*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 194*c4762a1bSJed Brown 195*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = DMDestroy(&user.da);CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr); 201*c4762a1bSJed Brown ierr = PetscFinalize(); 202*c4762a1bSJed Brown return ierr; 203*c4762a1bSJed Brown } 204*c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 205*c4762a1bSJed Brown /* 206*c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 207*c4762a1bSJed Brown 208*c4762a1bSJed Brown Input Parameters: 209*c4762a1bSJed Brown user - user-defined application context 210*c4762a1bSJed Brown X - vector 211*c4762a1bSJed Brown 212*c4762a1bSJed Brown Output Parameter: 213*c4762a1bSJed Brown X - vector 214*c4762a1bSJed Brown */ 215*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 216*c4762a1bSJed Brown { 217*c4762a1bSJed Brown PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm; 218*c4762a1bSJed Brown PetscErrorCode ierr; 219*c4762a1bSJed Brown PetscReal lambda,temp1,hx,hy,hz,tempk,tempj; 220*c4762a1bSJed Brown PetscScalar ***x; 221*c4762a1bSJed Brown 222*c4762a1bSJed Brown PetscFunctionBeginUser; 223*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 224*c4762a1bSJed Brown 225*c4762a1bSJed Brown lambda = user->param; 226*c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 227*c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 228*c4762a1bSJed Brown hz = 1.0/(PetscReal)(Mz-1); 229*c4762a1bSJed Brown temp1 = lambda/(lambda + 1.0); 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown /* 232*c4762a1bSJed Brown Get a pointer to vector data. 233*c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 234*c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 235*c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 236*c4762a1bSJed Brown the array. 237*c4762a1bSJed Brown */ 238*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da,X,&x);CHKERRQ(ierr); 239*c4762a1bSJed Brown 240*c4762a1bSJed Brown /* 241*c4762a1bSJed Brown Get local grid boundaries (for 3-dimensional DMDA): 242*c4762a1bSJed Brown xs, ys, zs - starting grid indices (no ghost points) 243*c4762a1bSJed Brown xm, ym, zm - widths of local grid (no ghost points) 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown */ 246*c4762a1bSJed Brown ierr = DMDAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 247*c4762a1bSJed Brown 248*c4762a1bSJed Brown /* 249*c4762a1bSJed Brown Compute initial guess over the locally owned part of the grid 250*c4762a1bSJed Brown */ 251*c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 252*c4762a1bSJed Brown tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz; 253*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 254*c4762a1bSJed Brown tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk); 255*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 256*c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) { 257*c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */ 258*c4762a1bSJed Brown x[k][j][i] = 0.0; 259*c4762a1bSJed Brown } else { 260*c4762a1bSJed Brown x[k][j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj)); 261*c4762a1bSJed Brown } 262*c4762a1bSJed Brown } 263*c4762a1bSJed Brown } 264*c4762a1bSJed Brown } 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown /* 267*c4762a1bSJed Brown Restore vector 268*c4762a1bSJed Brown */ 269*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da,X,&x);CHKERRQ(ierr); 270*c4762a1bSJed Brown PetscFunctionReturn(0); 271*c4762a1bSJed Brown } 272*c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 273*c4762a1bSJed Brown /* 274*c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear function, F(x) on a ghosted patch 275*c4762a1bSJed Brown 276*c4762a1bSJed Brown Input Parameters: 277*c4762a1bSJed Brown . snes - the SNES context 278*c4762a1bSJed Brown . localX - input vector, this contains the ghosted region needed 279*c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 280*c4762a1bSJed Brown 281*c4762a1bSJed Brown Output Parameter: 282*c4762a1bSJed Brown . F - function vector, this does not contain a ghosted region 283*c4762a1bSJed Brown */ 284*c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(SNES snes,Vec localX,Vec F,void *ptr) 285*c4762a1bSJed Brown { 286*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 287*c4762a1bSJed Brown PetscErrorCode ierr; 288*c4762a1bSJed Brown PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm; 289*c4762a1bSJed Brown PetscReal two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc; 290*c4762a1bSJed Brown PetscScalar u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f; 291*c4762a1bSJed Brown DM da; 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown PetscFunctionBeginUser; 294*c4762a1bSJed Brown ierr = SNESGetDM(snes,&da);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 296*c4762a1bSJed Brown 297*c4762a1bSJed Brown lambda = user->param; 298*c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 299*c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 300*c4762a1bSJed Brown hz = 1.0/(PetscReal)(Mz-1); 301*c4762a1bSJed Brown sc = hx*hy*hz*lambda; 302*c4762a1bSJed Brown hxhzdhy = hx*hz/hy; 303*c4762a1bSJed Brown hyhzdhx = hy*hz/hx; 304*c4762a1bSJed Brown hxhydhz = hx*hy/hz; 305*c4762a1bSJed Brown 306*c4762a1bSJed Brown /* 307*c4762a1bSJed Brown Get pointers to vector data 308*c4762a1bSJed Brown */ 309*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr); 310*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown /* 313*c4762a1bSJed Brown Get local grid boundaries 314*c4762a1bSJed Brown */ 315*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 316*c4762a1bSJed Brown 317*c4762a1bSJed Brown /* 318*c4762a1bSJed Brown Compute function over the locally owned part of the grid 319*c4762a1bSJed Brown */ 320*c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 321*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 322*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 323*c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) { 324*c4762a1bSJed Brown f[k][j][i] = x[k][j][i]; 325*c4762a1bSJed Brown } else { 326*c4762a1bSJed Brown u = x[k][j][i]; 327*c4762a1bSJed Brown u_east = x[k][j][i+1]; 328*c4762a1bSJed Brown u_west = x[k][j][i-1]; 329*c4762a1bSJed Brown u_north = x[k][j+1][i]; 330*c4762a1bSJed Brown u_south = x[k][j-1][i]; 331*c4762a1bSJed Brown u_up = x[k+1][j][i]; 332*c4762a1bSJed Brown u_down = x[k-1][j][i]; 333*c4762a1bSJed Brown u_xx = (-u_east + two*u - u_west)*hyhzdhx; 334*c4762a1bSJed Brown u_yy = (-u_north + two*u - u_south)*hxhzdhy; 335*c4762a1bSJed Brown u_zz = (-u_up + two*u - u_down)*hxhydhz; 336*c4762a1bSJed Brown f[k][j][i] = u_xx + u_yy + u_zz - sc*PetscExpScalar(u); 337*c4762a1bSJed Brown } 338*c4762a1bSJed Brown } 339*c4762a1bSJed Brown } 340*c4762a1bSJed Brown } 341*c4762a1bSJed Brown 342*c4762a1bSJed Brown /* 343*c4762a1bSJed Brown Restore vectors 344*c4762a1bSJed Brown */ 345*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr); 346*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 347*c4762a1bSJed Brown ierr = PetscLogFlops(11.0*ym*xm);CHKERRQ(ierr); 348*c4762a1bSJed Brown PetscFunctionReturn(0); 349*c4762a1bSJed Brown } 350*c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 351*c4762a1bSJed Brown /* 352*c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x) on the entire domain 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown Input Parameters: 355*c4762a1bSJed Brown . snes - the SNES context 356*c4762a1bSJed Brown . X - input vector 357*c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 358*c4762a1bSJed Brown 359*c4762a1bSJed Brown Output Parameter: 360*c4762a1bSJed Brown . F - function vector 361*c4762a1bSJed Brown */ 362*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr) 363*c4762a1bSJed Brown { 364*c4762a1bSJed Brown PetscErrorCode ierr; 365*c4762a1bSJed Brown Vec localX; 366*c4762a1bSJed Brown DM da; 367*c4762a1bSJed Brown 368*c4762a1bSJed Brown PetscFunctionBeginUser; 369*c4762a1bSJed Brown ierr = SNESGetDM(snes,&da);CHKERRQ(ierr); 370*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr); 371*c4762a1bSJed Brown 372*c4762a1bSJed Brown /* 373*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 374*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 375*c4762a1bSJed Brown By placing code between these two statements, computations can be 376*c4762a1bSJed Brown done while messages are in transition. 377*c4762a1bSJed Brown */ 378*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 379*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 380*c4762a1bSJed Brown 381*c4762a1bSJed Brown ierr = FormFunctionLocal(snes,localX,F,ptr);CHKERRQ(ierr); 382*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 383*c4762a1bSJed Brown PetscFunctionReturn(0); 384*c4762a1bSJed Brown } 385*c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 386*c4762a1bSJed Brown /* 387*c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 388*c4762a1bSJed Brown 389*c4762a1bSJed Brown Input Parameters: 390*c4762a1bSJed Brown . snes - the SNES context 391*c4762a1bSJed Brown . x - input vector 392*c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian() 393*c4762a1bSJed Brown 394*c4762a1bSJed Brown Output Parameters: 395*c4762a1bSJed Brown . A - Jacobian matrix 396*c4762a1bSJed Brown . B - optionally different preconditioning matrix 397*c4762a1bSJed Brown 398*c4762a1bSJed Brown */ 399*c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 400*c4762a1bSJed Brown { 401*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; /* user-defined application context */ 402*c4762a1bSJed Brown Vec localX; 403*c4762a1bSJed Brown PetscErrorCode ierr; 404*c4762a1bSJed Brown PetscInt i,j,k,Mx,My,Mz; 405*c4762a1bSJed Brown MatStencil col[7],row; 406*c4762a1bSJed Brown PetscInt xs,ys,zs,xm,ym,zm; 407*c4762a1bSJed Brown PetscScalar lambda,v[7],hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc,***x; 408*c4762a1bSJed Brown DM da; 409*c4762a1bSJed Brown 410*c4762a1bSJed Brown PetscFunctionBeginUser; 411*c4762a1bSJed Brown ierr = SNESGetDM(snes,&da);CHKERRQ(ierr); 412*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr); 413*c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr); 414*c4762a1bSJed Brown 415*c4762a1bSJed Brown lambda = user->param; 416*c4762a1bSJed Brown hx = 1.0/(PetscReal)(Mx-1); 417*c4762a1bSJed Brown hy = 1.0/(PetscReal)(My-1); 418*c4762a1bSJed Brown hz = 1.0/(PetscReal)(Mz-1); 419*c4762a1bSJed Brown sc = hx*hy*hz*lambda; 420*c4762a1bSJed Brown hxhzdhy = hx*hz/hy; 421*c4762a1bSJed Brown hyhzdhx = hy*hz/hx; 422*c4762a1bSJed Brown hxhydhz = hx*hy/hz; 423*c4762a1bSJed Brown 424*c4762a1bSJed Brown /* 425*c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 426*c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 427*c4762a1bSJed Brown By placing code between these two statements, computations can be 428*c4762a1bSJed Brown done while messages are in transition. 429*c4762a1bSJed Brown */ 430*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 431*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 432*c4762a1bSJed Brown 433*c4762a1bSJed Brown /* 434*c4762a1bSJed Brown Get pointer to vector data 435*c4762a1bSJed Brown */ 436*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr); 437*c4762a1bSJed Brown 438*c4762a1bSJed Brown /* 439*c4762a1bSJed Brown Get local grid boundaries 440*c4762a1bSJed Brown */ 441*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 442*c4762a1bSJed Brown 443*c4762a1bSJed Brown /* 444*c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 445*c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 446*c4762a1bSJed Brown contiguous chunks of rows across the processors. 447*c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 448*c4762a1bSJed Brown locally (but any non-local elements will be sent to the 449*c4762a1bSJed Brown appropriate processor during matrix assembly). 450*c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 451*c4762a1bSJed Brown - We can set matrix entries either using either 452*c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above. 453*c4762a1bSJed Brown */ 454*c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 455*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 456*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 457*c4762a1bSJed Brown row.k = k; row.j = j; row.i = i; 458*c4762a1bSJed Brown /* boundary points */ 459*c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) { 460*c4762a1bSJed Brown v[0] = 1.0; 461*c4762a1bSJed Brown ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr); 462*c4762a1bSJed Brown } else { 463*c4762a1bSJed Brown /* interior grid points */ 464*c4762a1bSJed Brown v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j; col[0].i = i; 465*c4762a1bSJed Brown v[1] = -hxhzdhy; col[1].k=k; col[1].j=j-1;col[1].i = i; 466*c4762a1bSJed Brown v[2] = -hyhzdhx; col[2].k=k; col[2].j=j; col[2].i = i-1; 467*c4762a1bSJed Brown v[3] = 2.0*(hyhzdhx+hxhzdhy+hxhydhz)-sc*PetscExpScalar(x[k][j][i]);col[3].k=row.k;col[3].j=row.j;col[3].i = row.i; 468*c4762a1bSJed Brown v[4] = -hyhzdhx; col[4].k=k; col[4].j=j; col[4].i = i+1; 469*c4762a1bSJed Brown v[5] = -hxhzdhy; col[5].k=k; col[5].j=j+1;col[5].i = i; 470*c4762a1bSJed Brown v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j; col[6].i = i; 471*c4762a1bSJed Brown ierr = MatSetValuesStencil(jac,1,&row,7,col,v,INSERT_VALUES);CHKERRQ(ierr); 472*c4762a1bSJed Brown } 473*c4762a1bSJed Brown } 474*c4762a1bSJed Brown } 475*c4762a1bSJed Brown } 476*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr); 477*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 478*c4762a1bSJed Brown 479*c4762a1bSJed Brown /* 480*c4762a1bSJed Brown Assemble matrix, using the 2-step process: 481*c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 482*c4762a1bSJed Brown */ 483*c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 484*c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 485*c4762a1bSJed Brown 486*c4762a1bSJed Brown /* 487*c4762a1bSJed Brown Normally since the matrix has already been assembled above; this 488*c4762a1bSJed Brown would do nothing. But in the matrix free mode -snes_mf_operator 489*c4762a1bSJed Brown this tells the "matrix-free" matrix that a new linear system solve 490*c4762a1bSJed Brown is about to be done. 491*c4762a1bSJed Brown */ 492*c4762a1bSJed Brown 493*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 494*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 495*c4762a1bSJed Brown 496*c4762a1bSJed Brown /* 497*c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 498*c4762a1bSJed Brown matrix. If we do, it will generate an error. 499*c4762a1bSJed Brown */ 500*c4762a1bSJed Brown ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 501*c4762a1bSJed Brown PetscFunctionReturn(0); 502*c4762a1bSJed Brown } 503*c4762a1bSJed Brown 504*c4762a1bSJed Brown 505*c4762a1bSJed Brown 506*c4762a1bSJed Brown /*TEST 507*c4762a1bSJed Brown 508*c4762a1bSJed Brown test: 509*c4762a1bSJed Brown nsize: 4 510*c4762a1bSJed Brown args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 511*c4762a1bSJed Brown 512*c4762a1bSJed Brown test: 513*c4762a1bSJed Brown suffix: 2 514*c4762a1bSJed Brown nsize: 4 515*c4762a1bSJed Brown args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 516*c4762a1bSJed Brown 517*c4762a1bSJed Brown test: 518*c4762a1bSJed Brown suffix: 3 519*c4762a1bSJed Brown nsize: 4 520*c4762a1bSJed Brown args: -fdcoloring -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 521*c4762a1bSJed Brown 522*c4762a1bSJed Brown test: 523*c4762a1bSJed Brown suffix: 3_ds 524*c4762a1bSJed Brown nsize: 4 525*c4762a1bSJed Brown args: -fdcoloring -fdcoloring_ds -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 526*c4762a1bSJed Brown 527*c4762a1bSJed Brown test: 528*c4762a1bSJed Brown suffix: 4 529*c4762a1bSJed Brown nsize: 4 530*c4762a1bSJed Brown args: -fdcoloring_local -fdcoloring -ksp_monitor_short -da_refine 1 531*c4762a1bSJed Brown requires: !single 532*c4762a1bSJed Brown 533*c4762a1bSJed Brown TEST*/ 534