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