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 /*T 10c4762a1bSJed Brown Concepts: SNES^sequential Bratu example 11c4762a1bSJed Brown Processors: 1 12c4762a1bSJed Brown T*/ 13c4762a1bSJed Brown 14c4762a1bSJed Brown /* ------------------------------------------------------------------------ 15c4762a1bSJed Brown 16c4762a1bSJed Brown Solid Fuel Ignition (SFI) problem. This problem is modeled by 17c4762a1bSJed Brown the partial differential equation 18c4762a1bSJed Brown 19c4762a1bSJed Brown -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 20c4762a1bSJed Brown 21c4762a1bSJed Brown with boundary conditions 22c4762a1bSJed Brown 23c4762a1bSJed Brown u = 0 for x = 0, x = 1, y = 0, y = 1. 24c4762a1bSJed Brown 25c4762a1bSJed Brown A finite difference approximation with the usual 5-point stencil 26c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a nonlinear 27c4762a1bSJed Brown system of equations. 28c4762a1bSJed Brown 29c4762a1bSJed Brown The parallel version of this code is snes/tutorials/ex5.c 30c4762a1bSJed Brown 31c4762a1bSJed Brown ------------------------------------------------------------------------- */ 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* 34c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that 35c4762a1bSJed Brown this file automatically includes: 36c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 37c4762a1bSJed Brown petscmat.h - matrices 38c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 39c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 40c4762a1bSJed Brown petscksp.h - linear solvers 41c4762a1bSJed Brown */ 42c4762a1bSJed Brown 43c4762a1bSJed Brown #include <petscsnes.h> 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* 46c4762a1bSJed Brown User-defined application context - contains data needed by the 47c4762a1bSJed Brown application-provided call-back routines, FormJacobian() and 48c4762a1bSJed Brown FormFunction(). 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown typedef struct { 51c4762a1bSJed Brown PetscReal param; /* test problem parameter */ 52c4762a1bSJed Brown PetscInt mx; /* Discretization in x-direction */ 53c4762a1bSJed Brown PetscInt my; /* Discretization in y-direction */ 54c4762a1bSJed Brown } AppCtx; 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* 57c4762a1bSJed Brown User-defined routines 58c4762a1bSJed Brown */ 59c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 60c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 61c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx*,Vec); 62c4762a1bSJed Brown extern PetscErrorCode ConvergenceTest(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 63c4762a1bSJed Brown extern PetscErrorCode ConvergenceDestroy(void*); 64c4762a1bSJed Brown extern PetscErrorCode postcheck(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 65c4762a1bSJed Brown 66c4762a1bSJed Brown int main(int argc,char **argv) 67c4762a1bSJed Brown { 68c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 69c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 70c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 71c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 72c4762a1bSJed Brown PetscErrorCode ierr; 73c4762a1bSJed Brown PetscInt i,its,N,hist_its[50]; 74c4762a1bSJed Brown PetscMPIInt size; 75c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81,bratu_lambda_min = 0.,history[50]; 76c4762a1bSJed Brown MatFDColoring fdcoloring; 77c4762a1bSJed Brown PetscBool matrix_free = PETSC_FALSE,flg,fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE,pc = PETSC_FALSE; 78c4762a1bSJed Brown KSP ksp; 79c4762a1bSJed Brown PetscInt *testarray; 80c4762a1bSJed Brown 81c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 82ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 83c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* 86c4762a1bSJed Brown Initialize problem parameters 87c4762a1bSJed Brown */ 88c4762a1bSJed Brown user.mx = 4; user.my = 4; user.param = 6.0; 89c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-pc",&pc,NULL);CHKERRQ(ierr); 93c4762a1bSJed Brown if (user.param >= bratu_lambda_max || user.param <= bratu_lambda_min) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda is out of range"); 94c4762a1bSJed Brown N = user.mx*user.my; 95c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-use_convergence_test",&use_convergence_test,NULL);CHKERRQ(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 98c4762a1bSJed Brown Create nonlinear solver context 99c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 100c4762a1bSJed Brown 101c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 102c4762a1bSJed Brown 103c4762a1bSJed Brown if (pc) { 104c4762a1bSJed Brown ierr = SNESSetType(snes,SNESNEWTONTR);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = SNESNewtonTRSetPostCheck(snes, postcheck,NULL);CHKERRQ(ierr); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Create vector data structures; set function evaluation routine 110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 111c4762a1bSJed Brown 112c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = VecSetSizes(x,PETSC_DECIDE,N);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* 118c4762a1bSJed Brown Set function evaluation routine and vector. Whenever the nonlinear 119c4762a1bSJed Brown solver needs to evaluate the nonlinear function, it will call this 120c4762a1bSJed Brown routine. 121c4762a1bSJed Brown - Note that the final routine argument is the user-defined 122c4762a1bSJed Brown context that provides application-specific data for the 123c4762a1bSJed Brown function evaluation routine. 124c4762a1bSJed Brown */ 125c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction,(void*)&user);CHKERRQ(ierr); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* 132c4762a1bSJed Brown Create matrix. Here we only approximately preallocate storage space 133c4762a1bSJed Brown for the Jacobian. See the users manual for a discussion of better 134c4762a1bSJed Brown techniques for preallocating matrix memory. 135c4762a1bSJed Brown */ 136c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL);CHKERRQ(ierr); 137c4762a1bSJed Brown if (!matrix_free) { 138c4762a1bSJed Brown PetscBool matrix_free_operator = PETSC_FALSE; 139c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-snes_mf_operator",&matrix_free_operator,NULL);CHKERRQ(ierr); 140c4762a1bSJed Brown if (matrix_free_operator) matrix_free = PETSC_FALSE; 141c4762a1bSJed Brown } 142c4762a1bSJed Brown if (!matrix_free) { 143c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&J);CHKERRQ(ierr); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* 147c4762a1bSJed Brown This option will cause the Jacobian to be computed via finite differences 148c4762a1bSJed Brown efficiently using a coloring of the columns of the matrix. 149c4762a1bSJed Brown */ 150c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-snes_fd_coloring",&fd_coloring,NULL);CHKERRQ(ierr); 151c4762a1bSJed Brown if (matrix_free && fd_coloring) SETERRQ(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"); 152c4762a1bSJed Brown 153c4762a1bSJed Brown if (fd_coloring) { 154c4762a1bSJed Brown ISColoring iscoloring; 155c4762a1bSJed Brown MatColoring mc; 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* 158c4762a1bSJed Brown This initializes the nonzero structure of the Jacobian. This is artificial 159c4762a1bSJed Brown because clearly if we had a routine to compute the Jacobian we won't need 160c4762a1bSJed Brown to use finite differences. 161c4762a1bSJed Brown */ 162c4762a1bSJed Brown ierr = FormJacobian(snes,x,J,J,&user);CHKERRQ(ierr); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* 165c4762a1bSJed Brown Color the matrix, i.e. determine groups of columns that share no common 166c4762a1bSJed Brown rows. These columns in the Jacobian can all be computed simulataneously. 167c4762a1bSJed Brown */ 168c4762a1bSJed Brown ierr = MatColoringCreate(J,&mc);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = MatColoringSetType(mc,MATCOLORINGSL);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 173c4762a1bSJed Brown /* 174c4762a1bSJed Brown Create the data structure that SNESComputeJacobianDefaultColor() uses 175c4762a1bSJed Brown to compute the actual Jacobians via finite differences. 176c4762a1bSJed Brown */ 177c4762a1bSJed Brown ierr = MatFDColoringCreate(J,iscoloring,&fdcoloring);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))FormFunction,&user);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = MatFDColoringSetUp(J,iscoloring,fdcoloring);CHKERRQ(ierr); 181c4762a1bSJed Brown /* 182c4762a1bSJed Brown Tell SNES to use the routine SNESComputeJacobianDefaultColor() 183c4762a1bSJed Brown to compute Jacobians. 184c4762a1bSJed Brown */ 185c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown /* 189c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 190c4762a1bSJed Brown routine. Whenever the nonlinear solver needs to compute the 191c4762a1bSJed Brown Jacobian matrix, it will call this routine. 192c4762a1bSJed Brown - Note that the final routine argument is the user-defined 193c4762a1bSJed Brown context that provides application-specific data for the 194c4762a1bSJed Brown Jacobian evaluation routine. 195c4762a1bSJed Brown - The user can override with: 196c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian 197c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 198c4762a1bSJed Brown (unless user explicitly sets preconditioner) 199c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 200c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 201c4762a1bSJed Brown products within Newton-Krylov method 202c4762a1bSJed Brown */ 203c4762a1bSJed Brown else if (!matrix_free) { 204c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user);CHKERRQ(ierr); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 208c4762a1bSJed Brown Customize nonlinear solver; set runtime options 209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 213c4762a1bSJed Brown */ 214c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* 217c4762a1bSJed Brown Set array that saves the function norms. This array is intended 218c4762a1bSJed Brown when the user wants to save the convergence history for later use 219c4762a1bSJed Brown rather than just to view the function norms via -snes_monitor. 220c4762a1bSJed Brown */ 221c4762a1bSJed Brown ierr = SNESSetConvergenceHistory(snes,history,hist_its,50,PETSC_TRUE);CHKERRQ(ierr); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* 224c4762a1bSJed Brown Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the 225c4762a1bSJed Brown user provided test before the specialized test. The convergence context is just an array to 226c4762a1bSJed Brown test that it gets properly freed at the end 227c4762a1bSJed Brown */ 228c4762a1bSJed Brown if (use_convergence_test) { 229c4762a1bSJed Brown ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = PetscMalloc1(5,&testarray);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = KSPSetConvergenceTest(ksp,ConvergenceTest,testarray,ConvergenceDestroy);CHKERRQ(ierr); 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 235c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 236c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 239c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 240c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 241c4762a1bSJed Brown this vector to zero by calling VecSet(). 242c4762a1bSJed Brown */ 243c4762a1bSJed Brown ierr = FormInitialGuess(&user,x);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr); 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* 249c4762a1bSJed Brown Print the convergence history. This is intended just to demonstrate 250c4762a1bSJed Brown use of the data attained via SNESSetConvergenceHistory(). 251c4762a1bSJed Brown */ 252c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-print_history",&flg);CHKERRQ(ierr); 253c4762a1bSJed Brown if (flg) { 254c4762a1bSJed Brown for (i=0; i<its+1; i++) { 255c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"iteration %D: Linear iterations %D Function norm = %g\n",i,hist_its[i],(double)history[i]);CHKERRQ(ierr); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown } 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 260c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 261c4762a1bSJed Brown are no longer needed. 262c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 263c4762a1bSJed Brown 264c4762a1bSJed Brown if (!matrix_free) { 265c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown if (fd_coloring) { 268c4762a1bSJed Brown ierr = MatFDColoringDestroy(&fdcoloring);CHKERRQ(ierr); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = PetscFinalize(); 274c4762a1bSJed Brown return ierr; 275c4762a1bSJed Brown } 276c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 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 */ 287c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 288c4762a1bSJed Brown { 289c4762a1bSJed Brown PetscInt i,j,row,mx,my; 290c4762a1bSJed Brown PetscErrorCode ierr; 291c4762a1bSJed Brown PetscReal lambda,temp1,temp,hx,hy; 292c4762a1bSJed Brown PetscScalar *x; 293c4762a1bSJed Brown 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 */ 308c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 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 */ 325c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 326c4762a1bSJed Brown return 0; 327c4762a1bSJed Brown } 328c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 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 */ 340c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr) 341c4762a1bSJed Brown { 342c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 343c4762a1bSJed Brown PetscInt i,j,row,mx,my; 344c4762a1bSJed Brown PetscErrorCode ierr; 345c4762a1bSJed Brown PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx; 346c4762a1bSJed Brown PetscScalar ut,ub,ul,ur,u,uxx,uyy,sc,*f; 347c4762a1bSJed Brown const PetscScalar *x; 348c4762a1bSJed Brown 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 */ 361c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 362c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 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 */ 388c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 389c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 390c4762a1bSJed Brown return 0; 391c4762a1bSJed Brown } 392c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 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 */ 406c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 407c4762a1bSJed Brown { 408c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; /* user-defined applicatin context */ 409c4762a1bSJed Brown PetscInt i,j,row,mx,my,col[5]; 410c4762a1bSJed Brown PetscErrorCode ierr; 411c4762a1bSJed Brown PetscScalar two = 2.0,one = 1.0,lambda,v[5],sc; 412c4762a1bSJed Brown const PetscScalar *x; 413c4762a1bSJed Brown PetscReal hx,hy,hxdhy,hydhx; 414c4762a1bSJed Brown 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 */ 427c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 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) { 436c4762a1bSJed Brown ierr = MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);CHKERRQ(ierr); 437c4762a1bSJed Brown continue; 438c4762a1bSJed Brown } 439c4762a1bSJed Brown v[0] = -hxdhy; col[0] = row - mx; 440c4762a1bSJed Brown v[1] = -hydhx; col[1] = row - 1; 441c4762a1bSJed Brown v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row; 442c4762a1bSJed Brown v[3] = -hydhx; col[3] = row + 1; 443c4762a1bSJed Brown v[4] = -hxdhy; col[4] = row + mx; 444c4762a1bSJed Brown ierr = MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr); 445c4762a1bSJed Brown } 446c4762a1bSJed Brown } 447c4762a1bSJed Brown 448c4762a1bSJed Brown /* 449c4762a1bSJed Brown Restore vector 450c4762a1bSJed Brown */ 451c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 452c4762a1bSJed Brown 453c4762a1bSJed Brown /* 454c4762a1bSJed Brown Assemble matrix 455c4762a1bSJed Brown */ 456c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 457c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 458c4762a1bSJed Brown 459c4762a1bSJed Brown if (jac != J) { 460c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 461c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 462c4762a1bSJed Brown } 463c4762a1bSJed Brown 464c4762a1bSJed Brown return 0; 465c4762a1bSJed Brown } 466c4762a1bSJed Brown 467c4762a1bSJed Brown PetscErrorCode ConvergenceTest(KSP ksp,PetscInt it,PetscReal nrm,KSPConvergedReason *reason,void *ctx) 468c4762a1bSJed Brown { 469c4762a1bSJed Brown PetscErrorCode ierr; 470c4762a1bSJed Brown 471c4762a1bSJed Brown PetscFunctionBegin; 472c4762a1bSJed Brown *reason = KSP_CONVERGED_ITERATING; 473c4762a1bSJed Brown if (it > 1) { 474c4762a1bSJed Brown *reason = KSP_CONVERGED_ITS; 475c4762a1bSJed Brown ierr = PetscInfo(NULL,"User provided convergence test returning after 2 iterations\n");CHKERRQ(ierr); 476c4762a1bSJed Brown } 477c4762a1bSJed Brown PetscFunctionReturn(0); 478c4762a1bSJed Brown } 479c4762a1bSJed Brown 480c4762a1bSJed Brown PetscErrorCode ConvergenceDestroy(void* ctx) 481c4762a1bSJed Brown { 482c4762a1bSJed Brown PetscErrorCode ierr; 483c4762a1bSJed Brown 484c4762a1bSJed Brown PetscFunctionBegin; 485c4762a1bSJed Brown ierr = PetscInfo(NULL,"User provided convergence destroy called\n");CHKERRQ(ierr); 486c4762a1bSJed Brown ierr = PetscFree(ctx);CHKERRQ(ierr); 487c4762a1bSJed Brown PetscFunctionReturn(0); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown 490c4762a1bSJed Brown PetscErrorCode postcheck(SNES snes,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx) 491c4762a1bSJed Brown { 492c4762a1bSJed Brown PetscErrorCode ierr; 493c4762a1bSJed Brown PetscReal norm; 494c4762a1bSJed Brown Vec tmp; 495c4762a1bSJed Brown 496c4762a1bSJed Brown PetscFunctionBegin; 497c4762a1bSJed Brown ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr); 498c4762a1bSJed Brown ierr = VecWAXPY(tmp,-1.0,x,w);CHKERRQ(ierr); 499c4762a1bSJed Brown ierr = VecNorm(tmp,NORM_2,&norm);CHKERRQ(ierr); 500c4762a1bSJed Brown ierr = VecDestroy(&tmp);CHKERRQ(ierr); 501c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of search step %g\n",(double)norm);CHKERRQ(ierr); 502c4762a1bSJed Brown PetscFunctionReturn(0); 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 521*dfd57a17SPierre 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 527*dfd57a17SPierre 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 537c4762a1bSJed Brown TEST*/ 538c4762a1bSJed Brown 539