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