1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves the nonlinear system, the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\ 3c4762a1bSJed Brown This example also illustrates the use of matrix coloring. Runtime options include:\n\ 4c4762a1bSJed Brown -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\ 5c4762a1bSJed Brown problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\ 6c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the x-direction\n\ 7c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the y-direction\n\n"; 8c4762a1bSJed Brown 9f6dfbefdSBarry Smith /* 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. 19c4762a1bSJed Brown 20c4762a1bSJed Brown A finite difference approximation with the usual 5-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 The parallel version of this code is snes/tutorials/ex5.c 25c4762a1bSJed Brown 26f6dfbefdSBarry Smith */ 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that 30c4762a1bSJed Brown this file automatically includes: 31c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 32c4762a1bSJed Brown petscmat.h - matrices 33c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 34c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 35c4762a1bSJed Brown petscksp.h - linear solvers 36c4762a1bSJed Brown */ 37c4762a1bSJed Brown 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 PetscInt mx; /* Discretization in x-direction */ 48c4762a1bSJed Brown PetscInt my; /* Discretization in y-direction */ 49c4762a1bSJed Brown } AppCtx; 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown User-defined routines 53c4762a1bSJed Brown */ 54c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 55c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 56c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx *, Vec); 57c4762a1bSJed Brown extern PetscErrorCode ConvergenceTest(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 58c4762a1bSJed Brown extern PetscErrorCode ConvergenceDestroy(void *); 59c4762a1bSJed Brown extern PetscErrorCode postcheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 60c4762a1bSJed Brown 61d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 62d71ae5a4SJacob Faibussowitsch { 63c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 64c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 65c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 66c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 67c4762a1bSJed Brown PetscInt i, its, N, hist_its[50]; 68c4762a1bSJed Brown PetscMPIInt size; 69c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81, bratu_lambda_min = 0., history[50]; 70c4762a1bSJed Brown MatFDColoring fdcoloring; 7178e7fe0eSHong Zhang PetscBool matrix_free = PETSC_FALSE, flg, fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE, pc = PETSC_FALSE, prunejacobian = PETSC_FALSE; 72c4762a1bSJed Brown KSP ksp; 73c4762a1bSJed Brown PetscInt *testarray; 74c4762a1bSJed Brown 75327415f7SBarry Smith PetscFunctionBeginUser; 769566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 78be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* 81c4762a1bSJed Brown Initialize problem parameters 82c4762a1bSJed Brown */ 839371c9d4SSatish Balay user.mx = 4; 849371c9d4SSatish Balay user.my = 4; 859371c9d4SSatish Balay user.param = 6.0; 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL)); 879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL)); 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL)); 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc", &pc, NULL)); 90e00437b9SBarry Smith PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range"); 91c4762a1bSJed Brown N = user.mx * user.my; 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL)); 9378e7fe0eSHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-prune_jacobian", &prunejacobian, NULL)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Create nonlinear solver context 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown if (pc) { 1029566063dSJacob Faibussowitsch PetscCall(SNESSetType(snes, SNESNEWTONTR)); 1039566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck, NULL)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107c4762a1bSJed Brown Create vector data structures; set function evaluation routine 108c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 1119566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, N)); 1129566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 1139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* 116c4762a1bSJed Brown Set function evaluation routine and vector. Whenever the nonlinear 117c4762a1bSJed Brown solver needs to evaluate the nonlinear function, it will call this 118c4762a1bSJed Brown routine. 119c4762a1bSJed Brown - Note that the final routine argument is the user-defined 120c4762a1bSJed Brown context that provides application-specific data for the 121c4762a1bSJed Brown function evaluation routine. 122c4762a1bSJed Brown */ 1239566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 126c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 127c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* 130c4762a1bSJed Brown Create matrix. Here we only approximately preallocate storage space 131c4762a1bSJed Brown for the Jacobian. See the users manual for a discussion of better 132c4762a1bSJed Brown techniques for preallocating matrix memory. 133c4762a1bSJed Brown */ 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL)); 135c4762a1bSJed Brown if (!matrix_free) { 136c4762a1bSJed Brown PetscBool matrix_free_operator = PETSC_FALSE; 1379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL)); 138c4762a1bSJed Brown if (matrix_free_operator) matrix_free = PETSC_FALSE; 139c4762a1bSJed Brown } 14048a46eb9SPierre Jolivet if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* 143c4762a1bSJed Brown This option will cause the Jacobian to be computed via finite differences 144c4762a1bSJed Brown efficiently using a coloring of the columns of the matrix. 145c4762a1bSJed Brown */ 1469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL)); 147e00437b9SBarry Smith PetscCheck(!matrix_free || !fd_coloring, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Use only one of -snes_mf, -snes_fd_coloring options!\nYou can do -snes_mf_operator -snes_fd_coloring"); 148c4762a1bSJed Brown 149c4762a1bSJed Brown if (fd_coloring) { 150c4762a1bSJed Brown ISColoring iscoloring; 151c4762a1bSJed Brown MatColoring mc; 15278e7fe0eSHong Zhang if (prunejacobian) { 15378e7fe0eSHong Zhang /* Initialize x with random nonzero values so that the nonzeros in the Jacobian 15478e7fe0eSHong Zhang can better reflect the sparsity structure of the Jacobian. */ 15578e7fe0eSHong Zhang PetscRandom rctx; 15678e7fe0eSHong Zhang PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 15778e7fe0eSHong Zhang PetscCall(PetscRandomSetInterval(rctx, 1.0, 2.0)); 15878e7fe0eSHong Zhang PetscCall(VecSetRandom(x, rctx)); 15978e7fe0eSHong Zhang PetscCall(PetscRandomDestroy(&rctx)); 16078e7fe0eSHong Zhang } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* 163c4762a1bSJed Brown This initializes the nonzero structure of the Jacobian. This is artificial 164c4762a1bSJed Brown because clearly if we had a routine to compute the Jacobian we won't need 165c4762a1bSJed Brown to use finite differences. 166c4762a1bSJed Brown */ 1679566063dSJacob Faibussowitsch PetscCall(FormJacobian(snes, x, J, J, &user)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* 170c4762a1bSJed Brown Color the matrix, i.e. determine groups of columns that share no common 171a5b23f4aSJose E. Roman rows. These columns in the Jacobian can all be computed simultaneously. 172c4762a1bSJed Brown */ 1739566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(J, &mc)); 1749566063dSJacob Faibussowitsch PetscCall(MatColoringSetType(mc, MATCOLORINGSL)); 1759566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 1769566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc, &iscoloring)); 1779566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown Create the data structure that SNESComputeJacobianDefaultColor() uses 180c4762a1bSJed Brown to compute the actual Jacobians via finite differences. 181c4762a1bSJed Brown */ 1829566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring)); 1839566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))FormFunction, &user)); 1849566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1859566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring)); 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown Tell SNES to use the routine SNESComputeJacobianDefaultColor() 188c4762a1bSJed Brown to compute Jacobians. 189c4762a1bSJed Brown */ 1909566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring)); 1919566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 19278e7fe0eSHong Zhang if (prunejacobian) PetscCall(SNESPruneJacobianColor(snes, J, J)); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown /* 195c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 196c4762a1bSJed Brown routine. Whenever the nonlinear solver needs to compute the 197c4762a1bSJed Brown Jacobian matrix, it will call this routine. 198c4762a1bSJed Brown - Note that the final routine argument is the user-defined 199c4762a1bSJed Brown context that provides application-specific data for the 200c4762a1bSJed Brown Jacobian evaluation routine. 201c4762a1bSJed Brown - The user can override with: 202c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian 203c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 204c4762a1bSJed Brown (unless user explicitly sets preconditioner) 205c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 206c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 207c4762a1bSJed Brown products within Newton-Krylov method 208c4762a1bSJed Brown */ 209c4762a1bSJed Brown else if (!matrix_free) { 2109566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user)); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 214c4762a1bSJed Brown Customize nonlinear solver; set runtime options 215c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* 218c4762a1bSJed Brown Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 219c4762a1bSJed Brown */ 2209566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown Set array that saves the function norms. This array is intended 224c4762a1bSJed Brown when the user wants to save the convergence history for later use 225c4762a1bSJed Brown rather than just to view the function norms via -snes_monitor. 226c4762a1bSJed Brown */ 2279566063dSJacob Faibussowitsch PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* 230c4762a1bSJed Brown Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the 231c4762a1bSJed Brown user provided test before the specialized test. The convergence context is just an array to 232c4762a1bSJed Brown test that it gets properly freed at the end 233c4762a1bSJed Brown */ 234c4762a1bSJed Brown if (use_convergence_test) { 2359566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &testarray)); 2379566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy)); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 241c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 242c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 243c4762a1bSJed Brown /* 244c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 245c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 246c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 247c4762a1bSJed Brown this vector to zero by calling VecSet(). 248c4762a1bSJed Brown */ 2499566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user, x)); 2509566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 2519566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 25263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* 255c4762a1bSJed Brown Print the convergence history. This is intended just to demonstrate 256c4762a1bSJed Brown use of the data attained via SNESSetConvergenceHistory(). 257c4762a1bSJed Brown */ 2589566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg)); 259c4762a1bSJed Brown if (flg) { 26048a46eb9SPierre Jolivet for (i = 0; i < its + 1; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n", i, hist_its[i], (double)history[i])); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 265c4762a1bSJed Brown are no longer needed. 266c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 267c4762a1bSJed Brown 26848a46eb9SPierre Jolivet if (!matrix_free) PetscCall(MatDestroy(&J)); 26948a46eb9SPierre Jolivet if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring)); 2709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2729566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2739566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 274b122ec5aSJacob Faibussowitsch return 0; 275c4762a1bSJed Brown } 276f6dfbefdSBarry Smith 277c4762a1bSJed Brown /* 278c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 279c4762a1bSJed Brown 280c4762a1bSJed Brown Input Parameters: 281c4762a1bSJed Brown user - user-defined application context 282c4762a1bSJed Brown X - vector 283c4762a1bSJed Brown 284c4762a1bSJed Brown Output Parameter: 285c4762a1bSJed Brown X - vector 286c4762a1bSJed Brown */ 287d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) 288d71ae5a4SJacob Faibussowitsch { 289c4762a1bSJed Brown PetscInt i, j, row, mx, my; 290c4762a1bSJed Brown PetscReal lambda, temp1, temp, hx, hy; 291c4762a1bSJed Brown PetscScalar *x; 292c4762a1bSJed Brown 2933ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 294c4762a1bSJed Brown mx = user->mx; 295c4762a1bSJed Brown my = user->my; 296c4762a1bSJed Brown lambda = user->param; 297c4762a1bSJed Brown 298c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx - 1); 299c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my - 1); 300c4762a1bSJed Brown 301c4762a1bSJed Brown /* 302c4762a1bSJed Brown Get a pointer to vector data. 303c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 304c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 305c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 306c4762a1bSJed Brown the array. 307c4762a1bSJed Brown */ 3089566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 309c4762a1bSJed Brown temp1 = lambda / (lambda + 1.0); 310c4762a1bSJed Brown for (j = 0; j < my; j++) { 311c4762a1bSJed Brown temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy; 312c4762a1bSJed Brown for (i = 0; i < mx; i++) { 313c4762a1bSJed Brown row = i + j * mx; 314c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 315c4762a1bSJed Brown x[row] = 0.0; 316c4762a1bSJed Brown continue; 317c4762a1bSJed Brown } 318c4762a1bSJed Brown x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp)); 319c4762a1bSJed Brown } 320c4762a1bSJed Brown } 321c4762a1bSJed Brown 322c4762a1bSJed Brown /* 323c4762a1bSJed Brown Restore vector 324c4762a1bSJed Brown */ 3259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 3263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 327c4762a1bSJed Brown } 328f6dfbefdSBarry Smith 329c4762a1bSJed Brown /* 330c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 331c4762a1bSJed Brown 332c4762a1bSJed Brown Input Parameters: 333c4762a1bSJed Brown . snes - the SNES context 334c4762a1bSJed Brown . X - input vector 335c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 336c4762a1bSJed Brown 337c4762a1bSJed Brown Output Parameter: 338c4762a1bSJed Brown . F - function vector 339c4762a1bSJed Brown */ 340d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) 341d71ae5a4SJacob Faibussowitsch { 342c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 343c4762a1bSJed Brown PetscInt i, j, row, mx, my; 344c4762a1bSJed Brown PetscReal two = 2.0, one = 1.0, lambda, hx, hy, hxdhy, hydhx; 345c4762a1bSJed Brown PetscScalar ut, ub, ul, ur, u, uxx, uyy, sc, *f; 346c4762a1bSJed Brown const PetscScalar *x; 347c4762a1bSJed Brown 3483ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 349c4762a1bSJed Brown mx = user->mx; 350c4762a1bSJed Brown my = user->my; 351c4762a1bSJed Brown lambda = user->param; 352c4762a1bSJed Brown hx = one / (PetscReal)(mx - 1); 353c4762a1bSJed Brown hy = one / (PetscReal)(my - 1); 354c4762a1bSJed Brown sc = hx * hy; 355c4762a1bSJed Brown hxdhy = hx / hy; 356c4762a1bSJed Brown hydhx = hy / hx; 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* 359c4762a1bSJed Brown Get pointers to vector data 360c4762a1bSJed Brown */ 3619566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 3629566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* 365c4762a1bSJed Brown Compute function 366c4762a1bSJed Brown */ 367c4762a1bSJed Brown for (j = 0; j < my; j++) { 368c4762a1bSJed Brown for (i = 0; i < mx; i++) { 369c4762a1bSJed Brown row = i + j * mx; 370c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 371c4762a1bSJed Brown f[row] = x[row]; 372c4762a1bSJed Brown continue; 373c4762a1bSJed Brown } 374c4762a1bSJed Brown u = x[row]; 375c4762a1bSJed Brown ub = x[row - mx]; 376c4762a1bSJed Brown ul = x[row - 1]; 377c4762a1bSJed Brown ut = x[row + mx]; 378c4762a1bSJed Brown ur = x[row + 1]; 379c4762a1bSJed Brown uxx = (-ur + two * u - ul) * hydhx; 380c4762a1bSJed Brown uyy = (-ut + two * u - ub) * hxdhy; 381c4762a1bSJed Brown f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u); 382c4762a1bSJed Brown } 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown /* 386c4762a1bSJed Brown Restore vectors 387c4762a1bSJed Brown */ 3889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 3899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 391c4762a1bSJed Brown } 392f6dfbefdSBarry Smith 393c4762a1bSJed Brown /* 394c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 395c4762a1bSJed Brown 396c4762a1bSJed Brown Input Parameters: 397c4762a1bSJed Brown . snes - the SNES context 398c4762a1bSJed Brown . x - input vector 399c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian() 400c4762a1bSJed Brown 401c4762a1bSJed Brown Output Parameters: 402c4762a1bSJed Brown . A - Jacobian matrix 403c4762a1bSJed Brown . B - optionally different preconditioning matrix 404c4762a1bSJed Brown . flag - flag indicating matrix structure 405c4762a1bSJed Brown */ 406d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) 407d71ae5a4SJacob Faibussowitsch { 408da81f932SPierre Jolivet AppCtx *user = (AppCtx *)ptr; /* user-defined application context */ 409c4762a1bSJed Brown PetscInt i, j, row, mx, my, col[5]; 410c4762a1bSJed Brown PetscScalar two = 2.0, one = 1.0, lambda, v[5], sc; 411c4762a1bSJed Brown const PetscScalar *x; 412c4762a1bSJed Brown PetscReal hx, hy, hxdhy, hydhx; 413c4762a1bSJed Brown 4143ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 415c4762a1bSJed Brown mx = user->mx; 416c4762a1bSJed Brown my = user->my; 417c4762a1bSJed Brown lambda = user->param; 418c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx - 1); 419c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my - 1); 420c4762a1bSJed Brown sc = hx * hy; 421c4762a1bSJed Brown hxdhy = hx / hy; 422c4762a1bSJed Brown hydhx = hy / hx; 423c4762a1bSJed Brown 424c4762a1bSJed Brown /* 425c4762a1bSJed Brown Get pointer to vector data 426c4762a1bSJed Brown */ 4279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 428c4762a1bSJed Brown 429c4762a1bSJed Brown /* 430c4762a1bSJed Brown Compute entries of the Jacobian 431c4762a1bSJed Brown */ 432c4762a1bSJed Brown for (j = 0; j < my; j++) { 433c4762a1bSJed Brown for (i = 0; i < mx; i++) { 434c4762a1bSJed Brown row = i + j * mx; 435c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 4369566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES)); 437c4762a1bSJed Brown continue; 438c4762a1bSJed Brown } 4399371c9d4SSatish Balay v[0] = -hxdhy; 4409371c9d4SSatish Balay col[0] = row - mx; 4419371c9d4SSatish Balay v[1] = -hydhx; 4429371c9d4SSatish Balay col[1] = row - 1; 4439371c9d4SSatish Balay v[2] = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]); 4449371c9d4SSatish Balay col[2] = row; 4459371c9d4SSatish Balay v[3] = -hydhx; 4469371c9d4SSatish Balay col[3] = row + 1; 4479371c9d4SSatish Balay v[4] = -hxdhy; 4489371c9d4SSatish Balay col[4] = row + mx; 4499566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES)); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown /* 454c4762a1bSJed Brown Restore vector 455c4762a1bSJed Brown */ 4569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 457c4762a1bSJed Brown 458c4762a1bSJed Brown /* 459c4762a1bSJed Brown Assemble matrix 460c4762a1bSJed Brown */ 4619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 4629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 463c4762a1bSJed Brown 464c4762a1bSJed Brown if (jac != J) { 4659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 467c4762a1bSJed Brown } 468c4762a1bSJed Brown 4693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 470c4762a1bSJed Brown } 471c4762a1bSJed Brown 472d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, void *ctx) 473d71ae5a4SJacob Faibussowitsch { 4743ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 475c4762a1bSJed Brown *reason = KSP_CONVERGED_ITERATING; 476c4762a1bSJed Brown if (it > 1) { 477c4762a1bSJed Brown *reason = KSP_CONVERGED_ITS; 4789566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n")); 479c4762a1bSJed Brown } 4803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 481c4762a1bSJed Brown } 482c4762a1bSJed Brown 483d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvergenceDestroy(void *ctx) 484d71ae5a4SJacob Faibussowitsch { 4853ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4869566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n")); 4879566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 4883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 489c4762a1bSJed Brown } 490c4762a1bSJed Brown 491d71ae5a4SJacob Faibussowitsch PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx) 492d71ae5a4SJacob Faibussowitsch { 493c4762a1bSJed Brown PetscReal norm; 494c4762a1bSJed Brown Vec tmp; 495c4762a1bSJed Brown 4963ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &tmp)); 4989566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tmp, -1.0, x, w)); 4999566063dSJacob Faibussowitsch PetscCall(VecNorm(tmp, NORM_2, &norm)); 5009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp)); 5019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm)); 5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 503c4762a1bSJed Brown } 504c4762a1bSJed Brown 505c4762a1bSJed Brown /*TEST 506c4762a1bSJed Brown 507c4762a1bSJed Brown build: 508c4762a1bSJed Brown requires: !single 509c4762a1bSJed Brown 510c4762a1bSJed Brown test: 511c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always 512c4762a1bSJed Brown 513c4762a1bSJed Brown test: 514c4762a1bSJed Brown suffix: 2 515c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always 516c4762a1bSJed Brown 517c4762a1bSJed Brown test: 518c4762a1bSJed Brown suffix: 2a 519c4762a1bSJed Brown filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault" 520c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info 521dfd57a17SPierre Jolivet requires: defined(PETSC_USE_INFO) 522c4762a1bSJed Brown 523c4762a1bSJed Brown test: 524c4762a1bSJed Brown suffix: 2b 525c4762a1bSJed Brown filter: grep -i "User provided convergence test" > /dev/null && echo "Found User provided convergence test" 526c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info 527dfd57a17SPierre Jolivet requires: defined(PETSC_USE_INFO) 528c4762a1bSJed Brown 529c4762a1bSJed Brown test: 530c4762a1bSJed Brown suffix: 3 531c4762a1bSJed Brown args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always 532c4762a1bSJed Brown 533c4762a1bSJed Brown test: 534c4762a1bSJed Brown suffix: 4 535c4762a1bSJed Brown args: -pc -par 6.807 -snes_monitor -snes_converged_reason 536c4762a1bSJed Brown 53778e7fe0eSHong Zhang test: 53878e7fe0eSHong Zhang suffix: 5 53978e7fe0eSHong Zhang args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always -prune_jacobian 54078e7fe0eSHong Zhang output_file: output/ex1_3.out 541*ccb5f961SBarry Smith 542*ccb5f961SBarry Smith test: 543*ccb5f961SBarry Smith suffix: 6 544*ccb5f961SBarry Smith args: -snes_monitor draw:image:testfile -viewer_view 545*ccb5f961SBarry Smith 546c4762a1bSJed Brown TEST*/ 547