1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\ 3c4762a1bSJed Brown This example tests PCVPBJacobiSetBlocks().\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /*T 6c4762a1bSJed Brown Concepts: SNES^basic uniprocessor example 7c4762a1bSJed Brown Processors: 1 8c4762a1bSJed Brown T*/ 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 12c4762a1bSJed Brown file automatically includes: 13c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 14c4762a1bSJed Brown petscmat.h - matrices 15c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 16c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 17c4762a1bSJed Brown petscksp.h - linear solvers 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown 20c4762a1bSJed Brown #include <petscsnes.h> 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* 23c4762a1bSJed Brown User-defined routines 24c4762a1bSJed Brown */ 25c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 26c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 27c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec); 28c4762a1bSJed Brown 29c4762a1bSJed Brown int main(int argc,char **argv) 30c4762a1bSJed Brown { 31c4762a1bSJed Brown SNES snes; /* SNES context */ 32c4762a1bSJed Brown Vec x,r,F,U; /* vectors */ 33c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 34c4762a1bSJed Brown PetscErrorCode ierr; 35c4762a1bSJed Brown PetscInt its,n = 5,i,maxit,maxf,lens[3] = {1,2,2}; 36c4762a1bSJed Brown PetscMPIInt size; 37c4762a1bSJed Brown PetscScalar h,xp,v,none = -1.0; 38c4762a1bSJed Brown PetscReal abstol,rtol,stol,norm; 39c4762a1bSJed Brown KSP ksp; 40c4762a1bSJed Brown PC pc; 41c4762a1bSJed Brown 42c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 43ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 44*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 45c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 46c4762a1bSJed Brown h = 1.0/(n-1); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 49c4762a1bSJed Brown Create nonlinear solver context 50c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 51c4762a1bSJed Brown 52c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PCSetType(pc,PCVPBJACOBI);CHKERRQ(ierr); 56c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 57c4762a1bSJed Brown Create vector data structures; set function evaluation routine 58c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 59c4762a1bSJed Brown /* 60c4762a1bSJed Brown Note that we form 1 vector from scratch and then duplicate as needed. 61c4762a1bSJed Brown */ 62c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = VecDuplicate(x,&F);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = VecDuplicate(x,&U);CHKERRQ(ierr); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* 70c4762a1bSJed Brown Set function evaluation routine and vector 71c4762a1bSJed Brown */ 72c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 76c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 77c4762a1bSJed Brown 78c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = MatSetFromOptions(J);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(J,3,NULL);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = MatSetVariableBlockSizes(J,3,lens);CHKERRQ(ierr); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* 85c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 86c4762a1bSJed Brown routine. User can override with: 87c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian 88c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 89c4762a1bSJed Brown (unless user explicitly sets preconditioner) 90c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 91c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 92c4762a1bSJed Brown products within Newton-Krylov method 93c4762a1bSJed Brown */ 94c4762a1bSJed Brown 95c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian,NULL);CHKERRQ(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 98c4762a1bSJed Brown Customize nonlinear solver; set runtime options 99c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* 102c4762a1bSJed Brown Set names for some vectors to facilitate monitoring (optional) 103c4762a1bSJed Brown */ 104c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");CHKERRQ(ierr); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* 108c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 109c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 110c4762a1bSJed Brown */ 111c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* 114c4762a1bSJed Brown Print parameters used for convergence testing (optional) ... just 115c4762a1bSJed Brown to demonstrate this routine; this information is also printed with 116c4762a1bSJed Brown the option -snes_view 117c4762a1bSJed Brown */ 118c4762a1bSJed Brown ierr = SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf);CHKERRQ(ierr); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122c4762a1bSJed Brown Initialize application: 123c4762a1bSJed Brown Store right-hand-side of PDE and exact solution 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125c4762a1bSJed Brown 126c4762a1bSJed Brown xp = 0.0; 127c4762a1bSJed Brown for (i=0; i<n; i++) { 128c4762a1bSJed Brown v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 129c4762a1bSJed Brown ierr = VecSetValues(F,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 130c4762a1bSJed Brown v = xp*xp*xp; 131c4762a1bSJed Brown ierr = VecSetValues(U,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 132c4762a1bSJed Brown xp += h; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 140c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 141c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 142c4762a1bSJed Brown this vector to zero by calling VecSet(). 143c4762a1bSJed Brown */ 144c4762a1bSJed Brown ierr = FormInitialGuess(x);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 150c4762a1bSJed Brown Check solution and clean up 151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* 154c4762a1bSJed Brown Check the error 155c4762a1bSJed Brown */ 156c4762a1bSJed Brown ierr = VecAXPY(x,none,U);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);CHKERRQ(ierr); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* 161c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 162c4762a1bSJed Brown are no longer needed. 163c4762a1bSJed Brown */ 164c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); ierr = VecDestroy(&F);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = PetscFinalize(); 168c4762a1bSJed Brown return ierr; 169c4762a1bSJed Brown } 170c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 171c4762a1bSJed Brown /* 172c4762a1bSJed Brown FormInitialGuess - Computes initial guess. 173c4762a1bSJed Brown 174c4762a1bSJed Brown Input/Output Parameter: 175c4762a1bSJed Brown . x - the solution vector 176c4762a1bSJed Brown */ 177c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x) 178c4762a1bSJed Brown { 179c4762a1bSJed Brown PetscErrorCode ierr; 180c4762a1bSJed Brown PetscScalar pfive = .50; 181c4762a1bSJed Brown ierr = VecSet(x,pfive);CHKERRQ(ierr); 182c4762a1bSJed Brown return 0; 183c4762a1bSJed Brown } 184c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 187c4762a1bSJed Brown 188c4762a1bSJed Brown Input Parameters: 189c4762a1bSJed Brown . snes - the SNES context 190c4762a1bSJed Brown . x - input vector 191c4762a1bSJed Brown . ctx - optional user-defined context, as set by SNESSetFunction() 192c4762a1bSJed Brown 193c4762a1bSJed Brown Output Parameter: 194c4762a1bSJed Brown . f - function vector 195c4762a1bSJed Brown 196c4762a1bSJed Brown Note: 197c4762a1bSJed Brown The user-defined context can contain any application-specific data 198c4762a1bSJed Brown needed for the function evaluation (such as various parameters, work 199c4762a1bSJed Brown vectors, and grid information). In this program the context is just 200c4762a1bSJed Brown a vector containing the right-hand-side of the discretized PDE. 201c4762a1bSJed Brown */ 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx) 204c4762a1bSJed Brown { 205c4762a1bSJed Brown Vec g = (Vec)ctx; 206c4762a1bSJed Brown const PetscScalar *xx,*gg; 207c4762a1bSJed Brown PetscScalar *ff,d; 208c4762a1bSJed Brown PetscErrorCode ierr; 209c4762a1bSJed Brown PetscInt i,n; 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Get pointers to vector data. 213c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 214c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 215c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 216c4762a1bSJed Brown the array. 217c4762a1bSJed Brown */ 218c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = VecGetArrayRead(g,&gg);CHKERRQ(ierr); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown Compute function 224c4762a1bSJed Brown */ 225c4762a1bSJed Brown ierr = VecGetSize(x,&n);CHKERRQ(ierr); 226c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 227c4762a1bSJed Brown ff[0] = xx[0]; 228c4762a1bSJed Brown for (i=1; i<n-1; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - gg[i]; 229c4762a1bSJed Brown ff[n-1] = xx[n-1] - 1.0; 230c4762a1bSJed Brown 231c4762a1bSJed Brown /* 232c4762a1bSJed Brown Restore vectors 233c4762a1bSJed Brown */ 234c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = VecRestoreArrayRead(g,&gg);CHKERRQ(ierr); 237c4762a1bSJed Brown return 0; 238c4762a1bSJed Brown } 239c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 240c4762a1bSJed Brown /* 241c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 242c4762a1bSJed Brown 243c4762a1bSJed Brown Input Parameters: 244c4762a1bSJed Brown . snes - the SNES context 245c4762a1bSJed Brown . x - input vector 246c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 247c4762a1bSJed Brown 248c4762a1bSJed Brown Output Parameters: 249c4762a1bSJed Brown . jac - Jacobian matrix 250c4762a1bSJed Brown . B - optionally different preconditioning matrix 251c4762a1bSJed Brown 252c4762a1bSJed Brown */ 253c4762a1bSJed Brown 254c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 255c4762a1bSJed Brown { 256c4762a1bSJed Brown const PetscScalar *xx; 257c4762a1bSJed Brown PetscScalar A[3],d; 258c4762a1bSJed Brown PetscErrorCode ierr; 259c4762a1bSJed Brown PetscInt i,n,j[3]; 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* 262c4762a1bSJed Brown Get pointer to vector data 263c4762a1bSJed Brown */ 264c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* 267c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 268c4762a1bSJed Brown - Note that in this case we set all elements for a particular 269c4762a1bSJed Brown row at once. 270c4762a1bSJed Brown */ 271c4762a1bSJed Brown ierr = VecGetSize(x,&n);CHKERRQ(ierr); 272c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* 275c4762a1bSJed Brown Interior grid points 276c4762a1bSJed Brown */ 277c4762a1bSJed Brown for (i=1; i<n-1; i++) { 278c4762a1bSJed Brown j[0] = i - 1; j[1] = i; j[2] = i + 1; 279c4762a1bSJed Brown A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i]; 280c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr); 281c4762a1bSJed Brown } 282c4762a1bSJed Brown 283c4762a1bSJed Brown /* 284c4762a1bSJed Brown Boundary points 285c4762a1bSJed Brown */ 286c4762a1bSJed Brown i = 0; A[0] = 1.0; 287c4762a1bSJed Brown 288c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr); 289c4762a1bSJed Brown 290c4762a1bSJed Brown i = n-1; A[0] = 1.0; 291c4762a1bSJed Brown 292c4762a1bSJed Brown ierr = MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr); 293c4762a1bSJed Brown 294c4762a1bSJed Brown /* 295c4762a1bSJed Brown Restore vector 296c4762a1bSJed Brown */ 297c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* 300c4762a1bSJed Brown Assemble matrix 301c4762a1bSJed Brown */ 302c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 303c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 304c4762a1bSJed Brown if (jac != B) { 305c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown return 0; 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown /*TEST 312c4762a1bSJed Brown 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown args: -snes_monitor_short -snes_view -ksp_monitor 315c4762a1bSJed Brown 316c4762a1bSJed Brown # this is just a test for SNESKSPTRASPOSEONLY and KSPSolveTranspose to behave properly 317c4762a1bSJed Brown # the solution is wrong on purpose 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown requires: !single !complex 320c4762a1bSJed Brown suffix: transpose_only 321c4762a1bSJed Brown args: -snes_monitor_short -snes_view -ksp_monitor -snes_type ksptransposeonly -pc_type ilu -snes_test_jacobian -snes_test_jacobian_view -ksp_view_rhs -ksp_view_solution -ksp_view_mat_explicit -ksp_view_preconditioned_operator_explicit 322c4762a1bSJed Brown 323c4762a1bSJed Brown TEST*/ 324