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 58d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 59d71ae5a4SJacob Faibussowitsch { 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 73327415f7SBarry Smith PetscFunctionBeginUser; 749566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 77c4762a1bSJed Brown Initialize problem parameters 78c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 79c4762a1bSJed Brown user.param = 6.0; 809566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL)); 81e00437b9SBarry Smith PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range"); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Create nonlinear solver context 85c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 869566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 89c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 90c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 919566063dSJacob 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)); 929566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da)); 939566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da)); 94c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 96c4762a1bSJed Brown vectors that are the same types 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 989566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.da, &x)); 999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 102c4762a1bSJed Brown Set function evaluation routine and vector 103c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1049566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 108c4762a1bSJed Brown 109c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 110c4762a1bSJed Brown routine. User can override with: 111c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 112c4762a1bSJed Brown (unless user explicitly sets preconditioner) 113c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 114c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 115c4762a1bSJed Brown products within Newton-Krylov method 116c4762a1bSJed Brown -fdcoloring : using finite differences with coloring to compute the Jacobian 117c4762a1bSJed Brown 118c4762a1bSJed Brown Note one can use -matfd_coloring wp or ds the only reason for the -fdcoloring_ds option 119c4762a1bSJed Brown below is to test the call to MatFDColoringSetType(). 120c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL)); 1229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-fdcoloring", &coloring, NULL)); 1239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-fdcoloring_ds", &coloring_ds, NULL)); 1249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-fdcoloring_local", &local_coloring, NULL)); 125c4762a1bSJed Brown if (!matrix_free) { 1269566063dSJacob Faibussowitsch PetscCall(DMSetMatType(user.da, MATAIJ)); 1279566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.da, &J)); 128c4762a1bSJed Brown if (coloring) { 129c4762a1bSJed Brown ISColoring iscoloring; 130c4762a1bSJed Brown if (!local_coloring) { 1319566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(user.da, IS_COLORING_GLOBAL, &iscoloring)); 1329566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring)); 1339566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))FormFunction, &user)); 134c4762a1bSJed Brown } else { 1359566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(user.da, IS_COLORING_LOCAL, &iscoloring)); 1369566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring)); 1379566063dSJacob Faibussowitsch PetscCall(MatFDColoringUseDM(J, matfdcoloring)); 1389566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))FormFunctionLocal, &user)); 139c4762a1bSJed Brown } 1401baa6e33SBarry Smith if (coloring_ds) PetscCall(MatFDColoringSetType(matfdcoloring, MATMFFD_DS)); 1419566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 1429566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring)); 1439566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring)); 1449566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 145c4762a1bSJed Brown } else { 1469566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &user)); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151c4762a1bSJed Brown Customize nonlinear solver; set runtime options 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1539566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, user.da)); 1549566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 157c4762a1bSJed Brown Evaluate initial guess 158c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 159c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 160c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 161c4762a1bSJed Brown this vector to zero by calling VecSet(). 162c4762a1bSJed Brown */ 1639566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user, x)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 166c4762a1bSJed Brown Solve nonlinear system 167c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1689566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 1699566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172c4762a1bSJed Brown Explicitly check norm of the residual of the solution 173c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1749566063dSJacob Faibussowitsch PetscCall(FormFunction(snes, x, r, (void *)&user)); 1759566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &fnorm)); 17663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT " fnorm %g\n", its, (double)fnorm)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 179c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 180c4762a1bSJed Brown are no longer needed. 181c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 182c4762a1bSJed Brown 1839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1869566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da)); 1889566063dSJacob Faibussowitsch PetscCall(MatFDColoringDestroy(&matfdcoloring)); 1899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 190b122ec5aSJacob Faibussowitsch return 0; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 193c4762a1bSJed Brown /* 194c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 195c4762a1bSJed Brown 196c4762a1bSJed Brown Input Parameters: 197c4762a1bSJed Brown user - user-defined application context 198c4762a1bSJed Brown X - vector 199c4762a1bSJed Brown 200c4762a1bSJed Brown Output Parameter: 201c4762a1bSJed Brown X - vector 202c4762a1bSJed Brown */ 203d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) 204d71ae5a4SJacob Faibussowitsch { 205c4762a1bSJed Brown PetscInt i, j, k, Mx, My, Mz, xs, ys, zs, xm, ym, zm; 206c4762a1bSJed Brown PetscReal lambda, temp1, hx, hy, hz, tempk, tempj; 207c4762a1bSJed Brown PetscScalar ***x; 208c4762a1bSJed Brown 209c4762a1bSJed Brown PetscFunctionBeginUser; 2109566063dSJacob 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)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown lambda = user->param; 213c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1); 214c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1); 215c4762a1bSJed Brown hz = 1.0 / (PetscReal)(Mz - 1); 216c4762a1bSJed Brown temp1 = lambda / (lambda + 1.0); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* 219c4762a1bSJed Brown Get a pointer to vector data. 220c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 221c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 222c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 223c4762a1bSJed Brown the array. 224c4762a1bSJed Brown */ 2259566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, X, &x)); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* 228c4762a1bSJed Brown Get local grid boundaries (for 3-dimensional DMDA): 229c4762a1bSJed Brown xs, ys, zs - starting grid indices (no ghost points) 230c4762a1bSJed Brown xm, ym, zm - widths of local grid (no ghost points) 231c4762a1bSJed Brown 232c4762a1bSJed Brown */ 2339566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->da, &xs, &ys, &zs, &xm, &ym, &zm)); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* 236c4762a1bSJed Brown Compute initial guess over the locally owned part of the grid 237c4762a1bSJed Brown */ 238c4762a1bSJed Brown for (k = zs; k < zs + zm; k++) { 239c4762a1bSJed Brown tempk = (PetscReal)(PetscMin(k, Mz - k - 1)) * hz; 240c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 241c4762a1bSJed Brown tempj = PetscMin((PetscReal)(PetscMin(j, My - j - 1)) * hy, tempk); 242c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 243c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx - 1 || j == My - 1 || k == Mz - 1) { 244c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */ 245c4762a1bSJed Brown x[k][j][i] = 0.0; 246c4762a1bSJed Brown } else { 247c4762a1bSJed Brown x[k][j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, tempj)); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown } 250c4762a1bSJed Brown } 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown Restore vector 255c4762a1bSJed Brown */ 2569566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, X, &x)); 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 260c4762a1bSJed Brown /* 261c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear function, F(x) on a ghosted patch 262c4762a1bSJed Brown 263c4762a1bSJed Brown Input Parameters: 264c4762a1bSJed Brown . snes - the SNES context 265c4762a1bSJed Brown . localX - input vector, this contains the ghosted region needed 266c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 267c4762a1bSJed Brown 268c4762a1bSJed Brown Output Parameter: 269c4762a1bSJed Brown . F - function vector, this does not contain a ghosted region 270c4762a1bSJed Brown */ 271d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(SNES snes, Vec localX, Vec F, void *ptr) 272d71ae5a4SJacob Faibussowitsch { 273c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 274c4762a1bSJed Brown PetscInt i, j, k, Mx, My, Mz, xs, ys, zs, xm, ym, zm; 275c4762a1bSJed Brown PetscReal two = 2.0, lambda, hx, hy, hz, hxhzdhy, hyhzdhx, hxhydhz, sc; 276c4762a1bSJed Brown PetscScalar u_north, u_south, u_east, u_west, u_up, u_down, u, u_xx, u_yy, u_zz, ***x, ***f; 277c4762a1bSJed Brown DM da; 278c4762a1bSJed Brown 279c4762a1bSJed Brown PetscFunctionBeginUser; 2809566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 2819566063dSJacob 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)); 282c4762a1bSJed Brown 283c4762a1bSJed Brown lambda = user->param; 284c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1); 285c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1); 286c4762a1bSJed Brown hz = 1.0 / (PetscReal)(Mz - 1); 287c4762a1bSJed Brown sc = hx * hy * hz * lambda; 288c4762a1bSJed Brown hxhzdhy = hx * hz / hy; 289c4762a1bSJed Brown hyhzdhx = hy * hz / hx; 290c4762a1bSJed Brown hxhydhz = hx * hy / hz; 291c4762a1bSJed Brown 292c4762a1bSJed Brown /* 293c4762a1bSJed Brown Get pointers to vector data 294c4762a1bSJed Brown */ 2959566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 2969566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 297c4762a1bSJed Brown 298c4762a1bSJed Brown /* 299c4762a1bSJed Brown Get local grid boundaries 300c4762a1bSJed Brown */ 3019566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 302c4762a1bSJed Brown 303c4762a1bSJed Brown /* 304c4762a1bSJed Brown Compute function over the locally owned part of the grid 305c4762a1bSJed Brown */ 306c4762a1bSJed Brown for (k = zs; k < zs + zm; k++) { 307c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 308c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 309c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx - 1 || j == My - 1 || k == Mz - 1) { 310c4762a1bSJed Brown f[k][j][i] = x[k][j][i]; 311c4762a1bSJed Brown } else { 312c4762a1bSJed Brown u = x[k][j][i]; 313c4762a1bSJed Brown u_east = x[k][j][i + 1]; 314c4762a1bSJed Brown u_west = x[k][j][i - 1]; 315c4762a1bSJed Brown u_north = x[k][j + 1][i]; 316c4762a1bSJed Brown u_south = x[k][j - 1][i]; 317c4762a1bSJed Brown u_up = x[k + 1][j][i]; 318c4762a1bSJed Brown u_down = x[k - 1][j][i]; 319c4762a1bSJed Brown u_xx = (-u_east + two * u - u_west) * hyhzdhx; 320c4762a1bSJed Brown u_yy = (-u_north + two * u - u_south) * hxhzdhy; 321c4762a1bSJed Brown u_zz = (-u_up + two * u - u_down) * hxhydhz; 322c4762a1bSJed Brown f[k][j][i] = u_xx + u_yy + u_zz - sc * PetscExpScalar(u); 323c4762a1bSJed Brown } 324c4762a1bSJed Brown } 325c4762a1bSJed Brown } 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* 329c4762a1bSJed Brown Restore vectors 330c4762a1bSJed Brown */ 3319566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 3329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 3339566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * ym * xm)); 3343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 335c4762a1bSJed Brown } 336c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 337c4762a1bSJed Brown /* 338c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x) on the entire domain 339c4762a1bSJed Brown 340c4762a1bSJed Brown Input Parameters: 341c4762a1bSJed Brown . snes - the SNES context 342c4762a1bSJed Brown . X - input vector 343c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 344c4762a1bSJed Brown 345c4762a1bSJed Brown Output Parameter: 346c4762a1bSJed Brown . F - function vector 347c4762a1bSJed Brown */ 348d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) 349d71ae5a4SJacob Faibussowitsch { 350c4762a1bSJed Brown Vec localX; 351c4762a1bSJed Brown DM da; 352c4762a1bSJed Brown 353c4762a1bSJed Brown PetscFunctionBeginUser; 3549566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 3559566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 356c4762a1bSJed Brown 357c4762a1bSJed Brown /* 358c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 359c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 360c4762a1bSJed Brown By placing code between these two statements, computations can be 361c4762a1bSJed Brown done while messages are in transition. 362c4762a1bSJed Brown */ 3639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 3649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 365c4762a1bSJed Brown 3669566063dSJacob Faibussowitsch PetscCall(FormFunctionLocal(snes, localX, F, ptr)); 3679566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 371c4762a1bSJed Brown /* 372c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 373c4762a1bSJed Brown 374c4762a1bSJed Brown Input Parameters: 375c4762a1bSJed Brown . snes - the SNES context 376c4762a1bSJed Brown . x - input vector 377c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian() 378c4762a1bSJed Brown 379c4762a1bSJed Brown Output Parameters: 380c4762a1bSJed Brown . A - Jacobian matrix 381c4762a1bSJed Brown . B - optionally different preconditioning matrix 382c4762a1bSJed Brown 383c4762a1bSJed Brown */ 384d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) 385d71ae5a4SJacob Faibussowitsch { 386c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; /* user-defined application context */ 387c4762a1bSJed Brown Vec localX; 388c4762a1bSJed Brown PetscInt i, j, k, Mx, My, Mz; 389c4762a1bSJed Brown MatStencil col[7], row; 390c4762a1bSJed Brown PetscInt xs, ys, zs, xm, ym, zm; 391c4762a1bSJed Brown PetscScalar lambda, v[7], hx, hy, hz, hxhzdhy, hyhzdhx, hxhydhz, sc, ***x; 392c4762a1bSJed Brown DM da; 393c4762a1bSJed Brown 394c4762a1bSJed Brown PetscFunctionBeginUser; 3959566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 3969566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 3979566063dSJacob 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)); 398c4762a1bSJed Brown 399c4762a1bSJed Brown lambda = user->param; 400c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1); 401c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1); 402c4762a1bSJed Brown hz = 1.0 / (PetscReal)(Mz - 1); 403c4762a1bSJed Brown sc = hx * hy * hz * lambda; 404c4762a1bSJed Brown hxhzdhy = hx * hz / hy; 405c4762a1bSJed Brown hyhzdhx = hy * hz / hx; 406c4762a1bSJed Brown hxhydhz = hx * hy / hz; 407c4762a1bSJed Brown 408c4762a1bSJed Brown /* 409c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 410c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 411c4762a1bSJed Brown By placing code between these two statements, computations can be 412c4762a1bSJed Brown done while messages are in transition. 413c4762a1bSJed Brown */ 4149566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 4159566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 416c4762a1bSJed Brown 417c4762a1bSJed Brown /* 418c4762a1bSJed Brown Get pointer to vector data 419c4762a1bSJed Brown */ 4209566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* 423c4762a1bSJed Brown Get local grid boundaries 424c4762a1bSJed Brown */ 4259566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 426c4762a1bSJed Brown 427c4762a1bSJed Brown /* 428c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 429c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 430c4762a1bSJed Brown contiguous chunks of rows across the processors. 431c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 432c4762a1bSJed Brown locally (but any non-local elements will be sent to the 433c4762a1bSJed Brown appropriate processor during matrix assembly). 434c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 435c4762a1bSJed Brown - We can set matrix entries either using either 436c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above. 437c4762a1bSJed Brown */ 438c4762a1bSJed Brown for (k = zs; k < zs + zm; k++) { 439c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 440c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 4419371c9d4SSatish Balay row.k = k; 4429371c9d4SSatish Balay row.j = j; 4439371c9d4SSatish Balay row.i = i; 444c4762a1bSJed Brown /* boundary points */ 445c4762a1bSJed Brown if (i == 0 || j == 0 || k == 0 || i == Mx - 1 || j == My - 1 || k == Mz - 1) { 446c4762a1bSJed Brown v[0] = 1.0; 4479566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES)); 448c4762a1bSJed Brown } else { 449c4762a1bSJed Brown /* interior grid points */ 4509371c9d4SSatish Balay v[0] = -hxhydhz; 4519371c9d4SSatish Balay col[0].k = k - 1; 4529371c9d4SSatish Balay col[0].j = j; 4539371c9d4SSatish Balay col[0].i = i; 4549371c9d4SSatish Balay v[1] = -hxhzdhy; 4559371c9d4SSatish Balay col[1].k = k; 4569371c9d4SSatish Balay col[1].j = j - 1; 4579371c9d4SSatish Balay col[1].i = i; 4589371c9d4SSatish Balay v[2] = -hyhzdhx; 4599371c9d4SSatish Balay col[2].k = k; 4609371c9d4SSatish Balay col[2].j = j; 4619371c9d4SSatish Balay col[2].i = i - 1; 4629371c9d4SSatish Balay v[3] = 2.0 * (hyhzdhx + hxhzdhy + hxhydhz) - sc * PetscExpScalar(x[k][j][i]); 4639371c9d4SSatish Balay col[3].k = row.k; 4649371c9d4SSatish Balay col[3].j = row.j; 4659371c9d4SSatish Balay col[3].i = row.i; 4669371c9d4SSatish Balay v[4] = -hyhzdhx; 4679371c9d4SSatish Balay col[4].k = k; 4689371c9d4SSatish Balay col[4].j = j; 4699371c9d4SSatish Balay col[4].i = i + 1; 4709371c9d4SSatish Balay v[5] = -hxhzdhy; 4719371c9d4SSatish Balay col[5].k = k; 4729371c9d4SSatish Balay col[5].j = j + 1; 4739371c9d4SSatish Balay col[5].i = i; 4749371c9d4SSatish Balay v[6] = -hxhydhz; 4759371c9d4SSatish Balay col[6].k = k + 1; 4769371c9d4SSatish Balay col[6].j = j; 4779371c9d4SSatish Balay col[6].i = i; 4789566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 7, col, v, INSERT_VALUES)); 479c4762a1bSJed Brown } 480c4762a1bSJed Brown } 481c4762a1bSJed Brown } 482c4762a1bSJed Brown } 4839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 4849566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 485c4762a1bSJed Brown 486c4762a1bSJed Brown /* 487c4762a1bSJed Brown Assemble matrix, using the 2-step process: 488c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 489c4762a1bSJed Brown */ 4909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 4919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 492c4762a1bSJed Brown 493c4762a1bSJed Brown /* 494c4762a1bSJed Brown Normally since the matrix has already been assembled above; this 495*01c1178eSBarry Smith would do nothing. But in the matrix-free mode -snes_mf_operator 496c4762a1bSJed Brown this tells the "matrix-free" matrix that a new linear system solve 497c4762a1bSJed Brown is about to be done. 498c4762a1bSJed Brown */ 499c4762a1bSJed Brown 5009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 5019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 502c4762a1bSJed Brown 503c4762a1bSJed Brown /* 504c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 505c4762a1bSJed Brown matrix. If we do, it will generate an error. 506c4762a1bSJed Brown */ 5079566063dSJacob Faibussowitsch PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 509c4762a1bSJed Brown } 510c4762a1bSJed Brown 511c4762a1bSJed Brown /*TEST 512c4762a1bSJed Brown 513c4762a1bSJed Brown test: 514c4762a1bSJed Brown nsize: 4 515c4762a1bSJed Brown args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 516c4762a1bSJed Brown 517c4762a1bSJed Brown test: 518c4762a1bSJed Brown suffix: 2 519c4762a1bSJed Brown nsize: 4 520c4762a1bSJed Brown args: -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 521c4762a1bSJed Brown 522c4762a1bSJed Brown test: 523c4762a1bSJed Brown suffix: 3 524c4762a1bSJed Brown nsize: 4 525c4762a1bSJed Brown args: -fdcoloring -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 526c4762a1bSJed Brown 527c4762a1bSJed Brown test: 528c4762a1bSJed Brown suffix: 3_ds 529c4762a1bSJed Brown nsize: 4 530c4762a1bSJed Brown args: -fdcoloring -fdcoloring_ds -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 531c4762a1bSJed Brown 532c4762a1bSJed Brown test: 533c4762a1bSJed Brown suffix: 4 534c4762a1bSJed Brown nsize: 4 535c4762a1bSJed Brown args: -fdcoloring_local -fdcoloring -ksp_monitor_short -da_refine 1 536c4762a1bSJed Brown requires: !single 537c4762a1bSJed Brown 53841ba4c6cSHeeho Park test: 53941ba4c6cSHeeho Park suffix: 5 54041ba4c6cSHeeho Park nsize: 4 54141ba4c6cSHeeho Park args: -fdcoloring_local -fdcoloring -ksp_monitor_short -da_refine 1 -snes_type newtontrdc 54241ba4c6cSHeeho Park requires: !single 54341ba4c6cSHeeho Park 5444a221d59SStefano Zampini test: 5454a221d59SStefano Zampini suffix: 6 5464a221d59SStefano Zampini nsize: 4 5474a221d59SStefano Zampini args: -fdcoloring_local -fdcoloring -da_refine 1 -snes_type newtontr -snes_tr_fallback_type dogleg 5484a221d59SStefano Zampini requires: !single 5494a221d59SStefano Zampini 550c4762a1bSJed Brown TEST*/ 551