1c4762a1bSJed Brown static char help[] = "Solves the nonlinear system, the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\ 2c4762a1bSJed Brown This example also illustrates the use of matrix coloring. Runtime options include:\n\ 3c4762a1bSJed Brown -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\ 4c4762a1bSJed Brown problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\ 5c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the x-direction\n\ 6c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the y-direction\n\n"; 7c4762a1bSJed Brown 8f6dfbefdSBarry Smith /* 9c4762a1bSJed Brown 10c4762a1bSJed Brown Solid Fuel Ignition (SFI) problem. This problem is modeled by 11c4762a1bSJed Brown the partial differential equation 12c4762a1bSJed Brown 13c4762a1bSJed Brown -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 14c4762a1bSJed Brown 15c4762a1bSJed Brown with boundary conditions 16c4762a1bSJed Brown 17c4762a1bSJed Brown u = 0 for x = 0, x = 1, y = 0, y = 1. 18c4762a1bSJed Brown 19c4762a1bSJed Brown A finite difference approximation with the usual 5-point stencil 20c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a nonlinear 21c4762a1bSJed Brown system of equations. 22c4762a1bSJed Brown 23c4762a1bSJed Brown The parallel version of this code is snes/tutorials/ex5.c 24c4762a1bSJed Brown 25f6dfbefdSBarry Smith */ 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that 29c4762a1bSJed Brown this 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 37c4762a1bSJed Brown #include <petscsnes.h> 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* 40c4762a1bSJed Brown User-defined application context - contains data needed by the 41c4762a1bSJed Brown application-provided call-back routines, FormJacobian() and 42c4762a1bSJed Brown FormFunction(). 43c4762a1bSJed Brown */ 44c4762a1bSJed Brown typedef struct { 45c4762a1bSJed Brown PetscReal param; /* test problem parameter */ 46c4762a1bSJed Brown PetscInt mx; /* Discretization in x-direction */ 47c4762a1bSJed Brown PetscInt my; /* Discretization in y-direction */ 48c4762a1bSJed Brown } AppCtx; 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* 51c4762a1bSJed Brown User-defined routines 52c4762a1bSJed Brown */ 53c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 54c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 55c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx *, Vec); 56c4762a1bSJed Brown extern PetscErrorCode ConvergenceTest(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 57c4762a1bSJed Brown extern PetscErrorCode ConvergenceDestroy(void *); 58c4762a1bSJed Brown extern PetscErrorCode postcheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 59c4762a1bSJed Brown 60d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 61d71ae5a4SJacob Faibussowitsch { 62c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 63c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 64c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 65c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 66c4762a1bSJed Brown PetscInt i, its, N, hist_its[50]; 67c4762a1bSJed Brown PetscMPIInt size; 68c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81, bratu_lambda_min = 0., history[50]; 69c4762a1bSJed Brown MatFDColoring fdcoloring; 709be84c52SStefano Zampini PetscBool matrix_free = PETSC_FALSE, flg, fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE, pc = PETSC_FALSE, prunejacobian = PETSC_FALSE, null_appctx = PETSC_TRUE; 71c4762a1bSJed Brown KSP ksp; 72c4762a1bSJed Brown PetscInt *testarray; 73c4762a1bSJed Brown 74327415f7SBarry Smith PetscFunctionBeginUser; 759566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 77be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* 80c4762a1bSJed Brown Initialize problem parameters 81c4762a1bSJed Brown */ 829371c9d4SSatish Balay user.mx = 4; 839371c9d4SSatish Balay user.my = 4; 849371c9d4SSatish Balay user.param = 6.0; 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL)); 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL)); 879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL)); 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc", &pc, NULL)); 89e00437b9SBarry Smith PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range"); 90c4762a1bSJed Brown N = user.mx * user.my; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL)); 9278e7fe0eSHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-prune_jacobian", &prunejacobian, NULL)); 939be84c52SStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-null_appctx", &null_appctx, 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 1069be84c52SStefano Zampini /* Test application context handling from Python */ 1079be84c52SStefano Zampini if (!null_appctx) { PetscCall(SNESSetApplicationContext(snes, (void *)&user)); } 1089be84c52SStefano Zampini 109c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 110c4762a1bSJed Brown Create vector data structures; set function evaluation routine 111c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 1149566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, N)); 1159566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 1169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* 119c4762a1bSJed Brown Set function evaluation routine and vector. Whenever the nonlinear 120c4762a1bSJed Brown solver needs to evaluate the nonlinear function, it will call this 121c4762a1bSJed Brown routine. 122c4762a1bSJed Brown - Note that the final routine argument is the user-defined 123c4762a1bSJed Brown context that provides application-specific data for the 124c4762a1bSJed Brown function evaluation routine. 125c4762a1bSJed Brown */ 1269566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* 133c4762a1bSJed Brown Create matrix. Here we only approximately preallocate storage space 134c4762a1bSJed Brown for the Jacobian. See the users manual for a discussion of better 135c4762a1bSJed Brown techniques for preallocating matrix memory. 136c4762a1bSJed Brown */ 1379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL)); 138c4762a1bSJed Brown if (!matrix_free) { 139c4762a1bSJed Brown PetscBool matrix_free_operator = PETSC_FALSE; 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL)); 141c4762a1bSJed Brown if (matrix_free_operator) matrix_free = PETSC_FALSE; 142c4762a1bSJed Brown } 14348a46eb9SPierre Jolivet if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* 146c4762a1bSJed Brown This option will cause the Jacobian to be computed via finite differences 147c4762a1bSJed Brown efficiently using a coloring of the columns of the matrix. 148c4762a1bSJed Brown */ 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL)); 150e00437b9SBarry 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"); 151c4762a1bSJed Brown 152c4762a1bSJed Brown if (fd_coloring) { 153c4762a1bSJed Brown ISColoring iscoloring; 154c4762a1bSJed Brown MatColoring mc; 15578e7fe0eSHong Zhang if (prunejacobian) { 15678e7fe0eSHong Zhang /* Initialize x with random nonzero values so that the nonzeros in the Jacobian 15778e7fe0eSHong Zhang can better reflect the sparsity structure of the Jacobian. */ 15878e7fe0eSHong Zhang PetscRandom rctx; 15978e7fe0eSHong Zhang PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 16078e7fe0eSHong Zhang PetscCall(PetscRandomSetInterval(rctx, 1.0, 2.0)); 16178e7fe0eSHong Zhang PetscCall(VecSetRandom(x, rctx)); 16278e7fe0eSHong Zhang PetscCall(PetscRandomDestroy(&rctx)); 16378e7fe0eSHong Zhang } 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* 166c4762a1bSJed Brown This initializes the nonzero structure of the Jacobian. This is artificial 167c4762a1bSJed Brown because clearly if we had a routine to compute the Jacobian we won't need 168c4762a1bSJed Brown to use finite differences. 169c4762a1bSJed Brown */ 1709566063dSJacob Faibussowitsch PetscCall(FormJacobian(snes, x, J, J, &user)); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown Color the matrix, i.e. determine groups of columns that share no common 174a5b23f4aSJose E. Roman rows. These columns in the Jacobian can all be computed simultaneously. 175c4762a1bSJed Brown */ 1769566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(J, &mc)); 1779566063dSJacob Faibussowitsch PetscCall(MatColoringSetType(mc, MATCOLORINGSL)); 1789566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 1799566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc, &iscoloring)); 1809566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 181c4762a1bSJed Brown /* 182c4762a1bSJed Brown Create the data structure that SNESComputeJacobianDefaultColor() uses 183c4762a1bSJed Brown to compute the actual Jacobians via finite differences. 184c4762a1bSJed Brown */ 1859566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring)); 1869566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))FormFunction, &user)); 1879566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1889566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring)); 189c4762a1bSJed Brown /* 190c4762a1bSJed Brown Tell SNES to use the routine SNESComputeJacobianDefaultColor() 191c4762a1bSJed Brown to compute Jacobians. 192c4762a1bSJed Brown */ 1939566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring)); 1949566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 19578e7fe0eSHong Zhang if (prunejacobian) PetscCall(SNESPruneJacobianColor(snes, J, J)); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown /* 198c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 199c4762a1bSJed Brown routine. Whenever the nonlinear solver needs to compute the 200c4762a1bSJed Brown Jacobian matrix, it will call this routine. 201c4762a1bSJed Brown - Note that the final routine argument is the user-defined 202c4762a1bSJed Brown context that provides application-specific data for the 203c4762a1bSJed Brown Jacobian evaluation routine. 204c4762a1bSJed Brown - The user can override with: 205c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian 206c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 207c4762a1bSJed Brown (unless user explicitly sets preconditioner) 208c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 209c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 210c4762a1bSJed Brown products within Newton-Krylov method 211c4762a1bSJed Brown */ 212c4762a1bSJed Brown else if (!matrix_free) { 2139566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user)); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 217c4762a1bSJed Brown Customize nonlinear solver; set runtime options 218c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* 221c4762a1bSJed Brown Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 222c4762a1bSJed Brown */ 2239566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* 226c4762a1bSJed Brown Set array that saves the function norms. This array is intended 227c4762a1bSJed Brown when the user wants to save the convergence history for later use 228c4762a1bSJed Brown rather than just to view the function norms via -snes_monitor. 229c4762a1bSJed Brown */ 2309566063dSJacob Faibussowitsch PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE)); 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* 233c4762a1bSJed Brown Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the 234c4762a1bSJed Brown user provided test before the specialized test. The convergence context is just an array to 235c4762a1bSJed Brown test that it gets properly freed at the end 236c4762a1bSJed Brown */ 237c4762a1bSJed Brown if (use_convergence_test) { 2389566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &testarray)); 2409566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy)); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 243c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 244c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 245c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 246c4762a1bSJed Brown /* 247c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 248c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 249c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 250c4762a1bSJed Brown this vector to zero by calling VecSet(). 251c4762a1bSJed Brown */ 2529566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user, x)); 2539566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 2549566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 25563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* 258c4762a1bSJed Brown Print the convergence history. This is intended just to demonstrate 259c4762a1bSJed Brown use of the data attained via SNESSetConvergenceHistory(). 260c4762a1bSJed Brown */ 2619566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg)); 262c4762a1bSJed Brown if (flg) { 26348a46eb9SPierre 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])); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 267c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 268c4762a1bSJed Brown are no longer needed. 269c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 270c4762a1bSJed Brown 27148a46eb9SPierre Jolivet if (!matrix_free) PetscCall(MatDestroy(&J)); 27248a46eb9SPierre Jolivet if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2759566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 277b122ec5aSJacob Faibussowitsch return 0; 278c4762a1bSJed Brown } 279f6dfbefdSBarry Smith 280c4762a1bSJed Brown /* 281c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 282c4762a1bSJed Brown 283c4762a1bSJed Brown Input Parameters: 284c4762a1bSJed Brown user - user-defined application context 285c4762a1bSJed Brown X - vector 286c4762a1bSJed Brown 287c4762a1bSJed Brown Output Parameter: 288c4762a1bSJed Brown X - vector 289c4762a1bSJed Brown */ 290d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) 291d71ae5a4SJacob Faibussowitsch { 292c4762a1bSJed Brown PetscInt i, j, row, mx, my; 293c4762a1bSJed Brown PetscReal lambda, temp1, temp, hx, hy; 294c4762a1bSJed Brown PetscScalar *x; 295c4762a1bSJed Brown 2963ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 297c4762a1bSJed Brown mx = user->mx; 298c4762a1bSJed Brown my = user->my; 299c4762a1bSJed Brown lambda = user->param; 300c4762a1bSJed Brown 301c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx - 1); 302c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my - 1); 303c4762a1bSJed Brown 304c4762a1bSJed Brown /* 305c4762a1bSJed Brown Get a pointer to vector data. 306c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 307c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 308c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 309c4762a1bSJed Brown the array. 310c4762a1bSJed Brown */ 3119566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 312c4762a1bSJed Brown temp1 = lambda / (lambda + 1.0); 313c4762a1bSJed Brown for (j = 0; j < my; j++) { 314c4762a1bSJed Brown temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy; 315c4762a1bSJed Brown for (i = 0; i < mx; i++) { 316c4762a1bSJed Brown row = i + j * mx; 317c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 318c4762a1bSJed Brown x[row] = 0.0; 319c4762a1bSJed Brown continue; 320c4762a1bSJed Brown } 321c4762a1bSJed Brown x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp)); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown } 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* 326c4762a1bSJed Brown Restore vector 327c4762a1bSJed Brown */ 3289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 330c4762a1bSJed Brown } 331f6dfbefdSBarry Smith 332c4762a1bSJed Brown /* 333c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 334c4762a1bSJed Brown 335c4762a1bSJed Brown Input Parameters: 336c4762a1bSJed Brown . snes - the SNES context 337c4762a1bSJed Brown . X - input vector 338c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 339c4762a1bSJed Brown 340c4762a1bSJed Brown Output Parameter: 341c4762a1bSJed Brown . F - function vector 342c4762a1bSJed Brown */ 343d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) 344d71ae5a4SJacob Faibussowitsch { 345c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 346c4762a1bSJed Brown PetscInt i, j, row, mx, my; 347c4762a1bSJed Brown PetscReal two = 2.0, one = 1.0, lambda, hx, hy, hxdhy, hydhx; 348c4762a1bSJed Brown PetscScalar ut, ub, ul, ur, u, uxx, uyy, sc, *f; 349c4762a1bSJed Brown const PetscScalar *x; 350c4762a1bSJed Brown 3513ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 352c4762a1bSJed Brown mx = user->mx; 353c4762a1bSJed Brown my = user->my; 354c4762a1bSJed Brown lambda = user->param; 355c4762a1bSJed Brown hx = one / (PetscReal)(mx - 1); 356c4762a1bSJed Brown hy = one / (PetscReal)(my - 1); 357c4762a1bSJed Brown sc = hx * hy; 358c4762a1bSJed Brown hxdhy = hx / hy; 359c4762a1bSJed Brown hydhx = hy / hx; 360c4762a1bSJed Brown 361c4762a1bSJed Brown /* 362c4762a1bSJed Brown Get pointers to vector data 363c4762a1bSJed Brown */ 3649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 3659566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* 368c4762a1bSJed Brown Compute function 369c4762a1bSJed Brown */ 370c4762a1bSJed Brown for (j = 0; j < my; j++) { 371c4762a1bSJed Brown for (i = 0; i < mx; i++) { 372c4762a1bSJed Brown row = i + j * mx; 373c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 374c4762a1bSJed Brown f[row] = x[row]; 375c4762a1bSJed Brown continue; 376c4762a1bSJed Brown } 377c4762a1bSJed Brown u = x[row]; 378c4762a1bSJed Brown ub = x[row - mx]; 379c4762a1bSJed Brown ul = x[row - 1]; 380c4762a1bSJed Brown ut = x[row + mx]; 381c4762a1bSJed Brown ur = x[row + 1]; 382c4762a1bSJed Brown uxx = (-ur + two * u - ul) * hydhx; 383c4762a1bSJed Brown uyy = (-ut + two * u - ub) * hxdhy; 384c4762a1bSJed Brown f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown } 387c4762a1bSJed Brown 388c4762a1bSJed Brown /* 389c4762a1bSJed Brown Restore vectors 390c4762a1bSJed Brown */ 3919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 3929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 394c4762a1bSJed Brown } 395f6dfbefdSBarry Smith 396c4762a1bSJed Brown /* 397c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 398c4762a1bSJed Brown 399c4762a1bSJed Brown Input Parameters: 400c4762a1bSJed Brown . snes - the SNES context 401c4762a1bSJed Brown . x - input vector 402c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian() 403c4762a1bSJed Brown 404c4762a1bSJed Brown Output Parameters: 405c4762a1bSJed Brown . A - Jacobian matrix 406c4762a1bSJed Brown . B - optionally different preconditioning matrix 407c4762a1bSJed Brown . flag - flag indicating matrix structure 408c4762a1bSJed Brown */ 409d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) 410d71ae5a4SJacob Faibussowitsch { 411da81f932SPierre Jolivet AppCtx *user = (AppCtx *)ptr; /* user-defined application context */ 412c4762a1bSJed Brown PetscInt i, j, row, mx, my, col[5]; 413c4762a1bSJed Brown PetscScalar two = 2.0, one = 1.0, lambda, v[5], sc; 414c4762a1bSJed Brown const PetscScalar *x; 415c4762a1bSJed Brown PetscReal hx, hy, hxdhy, hydhx; 416c4762a1bSJed Brown 4173ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 418c4762a1bSJed Brown mx = user->mx; 419c4762a1bSJed Brown my = user->my; 420c4762a1bSJed Brown lambda = user->param; 421c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx - 1); 422c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my - 1); 423c4762a1bSJed Brown sc = hx * hy; 424c4762a1bSJed Brown hxdhy = hx / hy; 425c4762a1bSJed Brown hydhx = hy / hx; 426c4762a1bSJed Brown 427c4762a1bSJed Brown /* 428c4762a1bSJed Brown Get pointer to vector data 429c4762a1bSJed Brown */ 4309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* 433c4762a1bSJed Brown Compute entries of the Jacobian 434c4762a1bSJed Brown */ 435c4762a1bSJed Brown for (j = 0; j < my; j++) { 436c4762a1bSJed Brown for (i = 0; i < mx; i++) { 437c4762a1bSJed Brown row = i + j * mx; 438c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 4399566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES)); 440c4762a1bSJed Brown continue; 441c4762a1bSJed Brown } 4429371c9d4SSatish Balay v[0] = -hxdhy; 4439371c9d4SSatish Balay col[0] = row - mx; 4449371c9d4SSatish Balay v[1] = -hydhx; 4459371c9d4SSatish Balay col[1] = row - 1; 4469371c9d4SSatish Balay v[2] = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]); 4479371c9d4SSatish Balay col[2] = row; 4489371c9d4SSatish Balay v[3] = -hydhx; 4499371c9d4SSatish Balay col[3] = row + 1; 4509371c9d4SSatish Balay v[4] = -hxdhy; 4519371c9d4SSatish Balay col[4] = row + mx; 4529566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES)); 453c4762a1bSJed Brown } 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 456c4762a1bSJed Brown /* 457c4762a1bSJed Brown Restore vector 458c4762a1bSJed Brown */ 4599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 460c4762a1bSJed Brown 461c4762a1bSJed Brown /* 462c4762a1bSJed Brown Assemble matrix 463c4762a1bSJed Brown */ 4649566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 4659566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 466c4762a1bSJed Brown 467c4762a1bSJed Brown if (jac != J) { 4689566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4699566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 470c4762a1bSJed Brown } 471c4762a1bSJed Brown 4723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 473c4762a1bSJed Brown } 474c4762a1bSJed Brown 475d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, void *ctx) 476d71ae5a4SJacob Faibussowitsch { 4773ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 478c4762a1bSJed Brown *reason = KSP_CONVERGED_ITERATING; 479c4762a1bSJed Brown if (it > 1) { 480c4762a1bSJed Brown *reason = KSP_CONVERGED_ITS; 4819566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n")); 482c4762a1bSJed Brown } 4833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 484c4762a1bSJed Brown } 485c4762a1bSJed Brown 486d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvergenceDestroy(void *ctx) 487d71ae5a4SJacob Faibussowitsch { 4883ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4899566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n")); 4909566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 4913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 492c4762a1bSJed Brown } 493c4762a1bSJed Brown 494d71ae5a4SJacob Faibussowitsch PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx) 495d71ae5a4SJacob Faibussowitsch { 496c4762a1bSJed Brown PetscReal norm; 497c4762a1bSJed Brown Vec tmp; 498c4762a1bSJed Brown 4993ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 5009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &tmp)); 5019566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tmp, -1.0, x, w)); 5029566063dSJacob Faibussowitsch PetscCall(VecNorm(tmp, NORM_2, &norm)); 5039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp)); 5049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm)); 5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 506c4762a1bSJed Brown } 507c4762a1bSJed Brown 508c4762a1bSJed Brown /*TEST 509c4762a1bSJed Brown 510c4762a1bSJed Brown build: 511c4762a1bSJed Brown requires: !single 512c4762a1bSJed Brown 513c4762a1bSJed Brown test: 514c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always 515c4762a1bSJed Brown 516c4762a1bSJed Brown test: 517c4762a1bSJed Brown suffix: 2 518c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always 519c4762a1bSJed Brown 520c4762a1bSJed Brown test: 521c4762a1bSJed Brown suffix: 2a 522c4762a1bSJed Brown filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault" 523c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info 524dfd57a17SPierre Jolivet requires: defined(PETSC_USE_INFO) 525c4762a1bSJed Brown 526c4762a1bSJed Brown test: 527c4762a1bSJed Brown suffix: 2b 528c4762a1bSJed Brown filter: grep -i "User provided convergence test" > /dev/null && echo "Found User provided convergence test" 529c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info 530dfd57a17SPierre Jolivet requires: defined(PETSC_USE_INFO) 531c4762a1bSJed Brown 532c4762a1bSJed Brown test: 533*24fb275aSStefano Zampini suffix: 2c 534*24fb275aSStefano Zampini args: -snes_converged_reason -snes_type newtontr -snes_tr_qn {{same different}separate output} -pc_type mat -snes_view -snes_tr_qn_mat_type lmvmdfp -snes_tr_norm_type infinity 535*24fb275aSStefano Zampini 536*24fb275aSStefano Zampini test: 537c4762a1bSJed Brown suffix: 3 538c4762a1bSJed Brown args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always 539c4762a1bSJed Brown 540c4762a1bSJed Brown test: 541c4762a1bSJed Brown suffix: 4 542c4762a1bSJed Brown args: -pc -par 6.807 -snes_monitor -snes_converged_reason 543c4762a1bSJed Brown 54478e7fe0eSHong Zhang test: 54578e7fe0eSHong Zhang suffix: 5 54678e7fe0eSHong 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 54778e7fe0eSHong Zhang output_file: output/ex1_3.out 548ccb5f961SBarry Smith 549ccb5f961SBarry Smith test: 550ccb5f961SBarry Smith suffix: 6 551ccb5f961SBarry Smith args: -snes_monitor draw:image:testfile -viewer_view 552ccb5f961SBarry Smith 5539be84c52SStefano Zampini test: 5549be84c52SStefano Zampini suffix: python 5559be84c52SStefano Zampini requires: petsc4py 5569be84c52SStefano Zampini args: -python -snes_type python -snes_python_type ex1.py:MySNES -snes_view -null_appctx {{0 1}separate output} 5579be84c52SStefano Zampini localrunfiles: ex1.py 5589be84c52SStefano Zampini 559c4762a1bSJed Brown TEST*/ 560