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 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. 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 26c4762a1bSJed Brown ------------------------------------------------------------------------- */ 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 61c4762a1bSJed Brown int main(int argc,char **argv) 62c4762a1bSJed Brown { 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; 71c4762a1bSJed Brown PetscBool matrix_free = PETSC_FALSE,flg,fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE,pc = PETSC_FALSE; 72c4762a1bSJed Brown KSP ksp; 73c4762a1bSJed Brown PetscInt *testarray; 74c4762a1bSJed Brown 75*327415f7SBarry 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 */ 83c4762a1bSJed Brown user.mx = 4; user.my = 4; user.param = 6.0; 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL)); 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL)); 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL)); 879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-pc",&pc,NULL)); 88e00437b9SBarry Smith PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda is out of range"); 89c4762a1bSJed Brown N = user.mx*user.my; 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_convergence_test",&use_convergence_test,NULL)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 93c4762a1bSJed Brown Create nonlinear solver context 94c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown if (pc) { 999566063dSJacob Faibussowitsch PetscCall(SNESSetType(snes,SNESNEWTONTR)); 1009566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck,NULL)); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 104c4762a1bSJed Brown Create vector data structures; set function evaluation routine 105c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 1089566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,PETSC_DECIDE,N)); 1099566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 1109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* 113c4762a1bSJed Brown Set function evaluation routine and vector. Whenever the nonlinear 114c4762a1bSJed Brown solver needs to evaluate the nonlinear function, it will call this 115c4762a1bSJed Brown routine. 116c4762a1bSJed Brown - Note that the final routine argument is the user-defined 117c4762a1bSJed Brown context that provides application-specific data for the 118c4762a1bSJed Brown function evaluation routine. 119c4762a1bSJed Brown */ 1209566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)&user)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* 127c4762a1bSJed Brown Create matrix. Here we only approximately preallocate storage space 128c4762a1bSJed Brown for the Jacobian. See the users manual for a discussion of better 129c4762a1bSJed Brown techniques for preallocating matrix memory. 130c4762a1bSJed Brown */ 1319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL)); 132c4762a1bSJed Brown if (!matrix_free) { 133c4762a1bSJed Brown PetscBool matrix_free_operator = PETSC_FALSE; 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_mf_operator",&matrix_free_operator,NULL)); 135c4762a1bSJed Brown if (matrix_free_operator) matrix_free = PETSC_FALSE; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown if (!matrix_free) { 1389566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&J)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* 142c4762a1bSJed Brown This option will cause the Jacobian to be computed via finite differences 143c4762a1bSJed Brown efficiently using a coloring of the columns of the matrix. 144c4762a1bSJed Brown */ 1459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_fd_coloring",&fd_coloring,NULL)); 146e00437b9SBarry 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"); 147c4762a1bSJed Brown 148c4762a1bSJed Brown if (fd_coloring) { 149c4762a1bSJed Brown ISColoring iscoloring; 150c4762a1bSJed Brown MatColoring mc; 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* 153c4762a1bSJed Brown This initializes the nonzero structure of the Jacobian. This is artificial 154c4762a1bSJed Brown because clearly if we had a routine to compute the Jacobian we won't need 155c4762a1bSJed Brown to use finite differences. 156c4762a1bSJed Brown */ 1579566063dSJacob Faibussowitsch PetscCall(FormJacobian(snes,x,J,J,&user)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* 160c4762a1bSJed Brown Color the matrix, i.e. determine groups of columns that share no common 161a5b23f4aSJose E. Roman rows. These columns in the Jacobian can all be computed simultaneously. 162c4762a1bSJed Brown */ 1639566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(J,&mc)); 1649566063dSJacob Faibussowitsch PetscCall(MatColoringSetType(mc,MATCOLORINGSL)); 1659566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 1669566063dSJacob Faibussowitsch PetscCall(MatColoringApply(mc,&iscoloring)); 1679566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown Create the data structure that SNESComputeJacobianDefaultColor() uses 170c4762a1bSJed Brown to compute the actual Jacobians via finite differences. 171c4762a1bSJed Brown */ 1729566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J,iscoloring,&fdcoloring)); 1739566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))FormFunction,&user)); 1749566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1759566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(J,iscoloring,fdcoloring)); 176c4762a1bSJed Brown /* 177c4762a1bSJed Brown Tell SNES to use the routine SNESComputeJacobianDefaultColor() 178c4762a1bSJed Brown to compute Jacobians. 179c4762a1bSJed Brown */ 1809566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring)); 1819566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown /* 184c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 185c4762a1bSJed Brown routine. Whenever the nonlinear solver needs to compute the 186c4762a1bSJed Brown Jacobian matrix, it will call this routine. 187c4762a1bSJed Brown - Note that the final routine argument is the user-defined 188c4762a1bSJed Brown context that provides application-specific data for the 189c4762a1bSJed Brown Jacobian evaluation routine. 190c4762a1bSJed Brown - The user can override with: 191c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian 192c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 193c4762a1bSJed Brown (unless user explicitly sets preconditioner) 194c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 195c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 196c4762a1bSJed Brown products within Newton-Krylov method 197c4762a1bSJed Brown */ 198c4762a1bSJed Brown else if (!matrix_free) { 1999566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user)); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 203c4762a1bSJed Brown Customize nonlinear solver; set runtime options 204c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 208c4762a1bSJed Brown */ 2099566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Set array that saves the function norms. This array is intended 213c4762a1bSJed Brown when the user wants to save the convergence history for later use 214c4762a1bSJed Brown rather than just to view the function norms via -snes_monitor. 215c4762a1bSJed Brown */ 2169566063dSJacob Faibussowitsch PetscCall(SNESSetConvergenceHistory(snes,history,hist_its,50,PETSC_TRUE)); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* 219c4762a1bSJed Brown Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the 220c4762a1bSJed Brown user provided test before the specialized test. The convergence context is just an array to 221c4762a1bSJed Brown test that it gets properly freed at the end 222c4762a1bSJed Brown */ 223c4762a1bSJed Brown if (use_convergence_test) { 2249566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 2259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5,&testarray)); 2269566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp,ConvergenceTest,testarray,ConvergenceDestroy)); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 231c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 232c4762a1bSJed Brown /* 233c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 234c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 235c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 236c4762a1bSJed Brown this vector to zero by calling VecSet(). 237c4762a1bSJed Brown */ 2389566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(&user,x)); 2399566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 2409566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 24163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its)); 242c4762a1bSJed Brown 243c4762a1bSJed Brown /* 244c4762a1bSJed Brown Print the convergence history. This is intended just to demonstrate 245c4762a1bSJed Brown use of the data attained via SNESSetConvergenceHistory(). 246c4762a1bSJed Brown */ 2479566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-print_history",&flg)); 248c4762a1bSJed Brown if (flg) { 249c4762a1bSJed Brown for (i=0; i<its+1; i++) { 25063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n",i,hist_its[i],(double)history[i])); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 255c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 256c4762a1bSJed Brown are no longer needed. 257c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 258c4762a1bSJed Brown 259c4762a1bSJed Brown if (!matrix_free) { 2609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown if (fd_coloring) { 2639566063dSJacob Faibussowitsch PetscCall(MatFDColoringDestroy(&fdcoloring)); 264c4762a1bSJed Brown } 2659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2679566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 269b122ec5aSJacob Faibussowitsch return 0; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 274c4762a1bSJed Brown 275c4762a1bSJed Brown Input Parameters: 276c4762a1bSJed Brown user - user-defined application context 277c4762a1bSJed Brown X - vector 278c4762a1bSJed Brown 279c4762a1bSJed Brown Output Parameter: 280c4762a1bSJed Brown X - vector 281c4762a1bSJed Brown */ 282c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 283c4762a1bSJed Brown { 284c4762a1bSJed Brown PetscInt i,j,row,mx,my; 285c4762a1bSJed Brown PetscReal lambda,temp1,temp,hx,hy; 286c4762a1bSJed Brown PetscScalar *x; 287c4762a1bSJed Brown 288c4762a1bSJed Brown mx = user->mx; 289c4762a1bSJed Brown my = user->my; 290c4762a1bSJed Brown lambda = user->param; 291c4762a1bSJed Brown 292c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx-1); 293c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my-1); 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* 296c4762a1bSJed Brown Get a pointer to vector data. 297c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 298c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 299c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 300c4762a1bSJed Brown the array. 301c4762a1bSJed Brown */ 3029566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 303c4762a1bSJed Brown temp1 = lambda/(lambda + 1.0); 304c4762a1bSJed Brown for (j=0; j<my; j++) { 305c4762a1bSJed Brown temp = (PetscReal)(PetscMin(j,my-j-1))*hy; 306c4762a1bSJed Brown for (i=0; i<mx; i++) { 307c4762a1bSJed Brown row = i + j*mx; 308c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 309c4762a1bSJed Brown x[row] = 0.0; 310c4762a1bSJed Brown continue; 311c4762a1bSJed Brown } 312c4762a1bSJed Brown x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp)); 313c4762a1bSJed Brown } 314c4762a1bSJed Brown } 315c4762a1bSJed Brown 316c4762a1bSJed Brown /* 317c4762a1bSJed Brown Restore vector 318c4762a1bSJed Brown */ 3199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x)); 320c4762a1bSJed Brown return 0; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 323c4762a1bSJed Brown /* 324c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 325c4762a1bSJed Brown 326c4762a1bSJed Brown Input Parameters: 327c4762a1bSJed Brown . snes - the SNES context 328c4762a1bSJed Brown . X - input vector 329c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 330c4762a1bSJed Brown 331c4762a1bSJed Brown Output Parameter: 332c4762a1bSJed Brown . F - function vector 333c4762a1bSJed Brown */ 334c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr) 335c4762a1bSJed Brown { 336c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 337c4762a1bSJed Brown PetscInt i,j,row,mx,my; 338c4762a1bSJed Brown PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx; 339c4762a1bSJed Brown PetscScalar ut,ub,ul,ur,u,uxx,uyy,sc,*f; 340c4762a1bSJed Brown const PetscScalar *x; 341c4762a1bSJed Brown 342c4762a1bSJed Brown mx = user->mx; 343c4762a1bSJed Brown my = user->my; 344c4762a1bSJed Brown lambda = user->param; 345c4762a1bSJed Brown hx = one / (PetscReal)(mx-1); 346c4762a1bSJed Brown hy = one / (PetscReal)(my-1); 347c4762a1bSJed Brown sc = hx*hy; 348c4762a1bSJed Brown hxdhy = hx/hy; 349c4762a1bSJed Brown hydhx = hy/hx; 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* 352c4762a1bSJed Brown Get pointers to vector data 353c4762a1bSJed Brown */ 3549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 3559566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 356c4762a1bSJed Brown 357c4762a1bSJed Brown /* 358c4762a1bSJed Brown Compute function 359c4762a1bSJed Brown */ 360c4762a1bSJed Brown for (j=0; j<my; j++) { 361c4762a1bSJed Brown for (i=0; i<mx; i++) { 362c4762a1bSJed Brown row = i + j*mx; 363c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 364c4762a1bSJed Brown f[row] = x[row]; 365c4762a1bSJed Brown continue; 366c4762a1bSJed Brown } 367c4762a1bSJed Brown u = x[row]; 368c4762a1bSJed Brown ub = x[row - mx]; 369c4762a1bSJed Brown ul = x[row - 1]; 370c4762a1bSJed Brown ut = x[row + mx]; 371c4762a1bSJed Brown ur = x[row + 1]; 372c4762a1bSJed Brown uxx = (-ur + two*u - ul)*hydhx; 373c4762a1bSJed Brown uyy = (-ut + two*u - ub)*hxdhy; 374c4762a1bSJed Brown f[row] = uxx + uyy - sc*lambda*PetscExpScalar(u); 375c4762a1bSJed Brown } 376c4762a1bSJed Brown } 377c4762a1bSJed Brown 378c4762a1bSJed Brown /* 379c4762a1bSJed Brown Restore vectors 380c4762a1bSJed Brown */ 3819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 3829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 383c4762a1bSJed Brown return 0; 384c4762a1bSJed Brown } 385c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 386c4762a1bSJed Brown /* 387c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 388c4762a1bSJed Brown 389c4762a1bSJed Brown Input Parameters: 390c4762a1bSJed Brown . snes - the SNES context 391c4762a1bSJed Brown . x - input vector 392c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian() 393c4762a1bSJed Brown 394c4762a1bSJed Brown Output Parameters: 395c4762a1bSJed Brown . A - Jacobian matrix 396c4762a1bSJed Brown . B - optionally different preconditioning matrix 397c4762a1bSJed Brown . flag - flag indicating matrix structure 398c4762a1bSJed Brown */ 399c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 400c4762a1bSJed Brown { 401c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; /* user-defined applicatin context */ 402c4762a1bSJed Brown PetscInt i,j,row,mx,my,col[5]; 403c4762a1bSJed Brown PetscScalar two = 2.0,one = 1.0,lambda,v[5],sc; 404c4762a1bSJed Brown const PetscScalar *x; 405c4762a1bSJed Brown PetscReal hx,hy,hxdhy,hydhx; 406c4762a1bSJed Brown 407c4762a1bSJed Brown mx = user->mx; 408c4762a1bSJed Brown my = user->my; 409c4762a1bSJed Brown lambda = user->param; 410c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx-1); 411c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my-1); 412c4762a1bSJed Brown sc = hx*hy; 413c4762a1bSJed Brown hxdhy = hx/hy; 414c4762a1bSJed Brown hydhx = hy/hx; 415c4762a1bSJed Brown 416c4762a1bSJed Brown /* 417c4762a1bSJed Brown Get pointer to vector data 418c4762a1bSJed Brown */ 4199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* 422c4762a1bSJed Brown Compute entries of the Jacobian 423c4762a1bSJed Brown */ 424c4762a1bSJed Brown for (j=0; j<my; j++) { 425c4762a1bSJed Brown for (i=0; i<mx; i++) { 426c4762a1bSJed Brown row = i + j*mx; 427c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 4289566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES)); 429c4762a1bSJed Brown continue; 430c4762a1bSJed Brown } 431c4762a1bSJed Brown v[0] = -hxdhy; col[0] = row - mx; 432c4762a1bSJed Brown v[1] = -hydhx; col[1] = row - 1; 433c4762a1bSJed Brown v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row; 434c4762a1bSJed Brown v[3] = -hydhx; col[3] = row + 1; 435c4762a1bSJed Brown v[4] = -hxdhy; col[4] = row + mx; 4369566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES)); 437c4762a1bSJed Brown } 438c4762a1bSJed Brown } 439c4762a1bSJed Brown 440c4762a1bSJed Brown /* 441c4762a1bSJed Brown Restore vector 442c4762a1bSJed Brown */ 4439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 444c4762a1bSJed Brown 445c4762a1bSJed Brown /* 446c4762a1bSJed Brown Assemble matrix 447c4762a1bSJed Brown */ 4489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 4499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 450c4762a1bSJed Brown 451c4762a1bSJed Brown if (jac != J) { 4529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 4539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 456c4762a1bSJed Brown return 0; 457c4762a1bSJed Brown } 458c4762a1bSJed Brown 459c4762a1bSJed Brown PetscErrorCode ConvergenceTest(KSP ksp,PetscInt it,PetscReal nrm,KSPConvergedReason *reason,void *ctx) 460c4762a1bSJed Brown { 461c4762a1bSJed Brown PetscFunctionBegin; 462c4762a1bSJed Brown *reason = KSP_CONVERGED_ITERATING; 463c4762a1bSJed Brown if (it > 1) { 464c4762a1bSJed Brown *reason = KSP_CONVERGED_ITS; 4659566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"User provided convergence test returning after 2 iterations\n")); 466c4762a1bSJed Brown } 467c4762a1bSJed Brown PetscFunctionReturn(0); 468c4762a1bSJed Brown } 469c4762a1bSJed Brown 470c4762a1bSJed Brown PetscErrorCode ConvergenceDestroy(void* ctx) 471c4762a1bSJed Brown { 472c4762a1bSJed Brown PetscFunctionBegin; 4739566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"User provided convergence destroy called\n")); 4749566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 475c4762a1bSJed Brown PetscFunctionReturn(0); 476c4762a1bSJed Brown } 477c4762a1bSJed Brown 478c4762a1bSJed Brown PetscErrorCode postcheck(SNES snes,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx) 479c4762a1bSJed Brown { 480c4762a1bSJed Brown PetscReal norm; 481c4762a1bSJed Brown Vec tmp; 482c4762a1bSJed Brown 483c4762a1bSJed Brown PetscFunctionBegin; 4849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&tmp)); 4859566063dSJacob Faibussowitsch PetscCall(VecWAXPY(tmp,-1.0,x,w)); 4869566063dSJacob Faibussowitsch PetscCall(VecNorm(tmp,NORM_2,&norm)); 4879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp)); 4889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of search step %g\n",(double)norm)); 489c4762a1bSJed Brown PetscFunctionReturn(0); 490c4762a1bSJed Brown } 491c4762a1bSJed Brown 492c4762a1bSJed Brown /*TEST 493c4762a1bSJed Brown 494c4762a1bSJed Brown build: 495c4762a1bSJed Brown requires: !single 496c4762a1bSJed Brown 497c4762a1bSJed Brown test: 498c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always 499c4762a1bSJed Brown 500c4762a1bSJed Brown test: 501c4762a1bSJed Brown suffix: 2 502c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always 503c4762a1bSJed Brown 504c4762a1bSJed Brown test: 505c4762a1bSJed Brown suffix: 2a 506c4762a1bSJed Brown filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault" 507c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info 508dfd57a17SPierre Jolivet requires: defined(PETSC_USE_INFO) 509c4762a1bSJed Brown 510c4762a1bSJed Brown test: 511c4762a1bSJed Brown suffix: 2b 512c4762a1bSJed Brown filter: grep -i "User provided convergence test" > /dev/null && echo "Found User provided convergence test" 513c4762a1bSJed Brown args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info 514dfd57a17SPierre Jolivet requires: defined(PETSC_USE_INFO) 515c4762a1bSJed Brown 516c4762a1bSJed Brown test: 517c4762a1bSJed Brown suffix: 3 518c4762a1bSJed Brown args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always 519c4762a1bSJed Brown 520c4762a1bSJed Brown test: 521c4762a1bSJed Brown suffix: 4 522c4762a1bSJed Brown args: -pc -par 6.807 -snes_monitor -snes_converged_reason 523c4762a1bSJed Brown 524c4762a1bSJed Brown TEST*/ 525