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 /* 6c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 7c4762a1bSJed Brown file automatically includes: 8c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 9c4762a1bSJed Brown petscmat.h - matrices 10c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 11c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 12c4762a1bSJed Brown petscksp.h - linear solvers 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscsnes.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* 18c4762a1bSJed Brown User-defined routines 19c4762a1bSJed Brown */ 20c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 21c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 22c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(Vec); 23c4762a1bSJed Brown 24c4762a1bSJed Brown int main(int argc,char **argv) 25c4762a1bSJed Brown { 26c4762a1bSJed Brown SNES snes; /* SNES context */ 27c4762a1bSJed Brown Vec x,r,F,U; /* vectors */ 28c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 29c4762a1bSJed Brown PetscInt its,n = 5,i,maxit,maxf,lens[3] = {1,2,2}; 30c4762a1bSJed Brown PetscMPIInt size; 31c4762a1bSJed Brown PetscScalar h,xp,v,none = -1.0; 32c4762a1bSJed Brown PetscReal abstol,rtol,stol,norm; 33c4762a1bSJed Brown KSP ksp; 34c4762a1bSJed Brown PC pc; 35c4762a1bSJed Brown 36*327415f7SBarry Smith PetscFunctionBeginUser; 379566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 39be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 41c4762a1bSJed Brown h = 1.0/(n-1); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 44c4762a1bSJed Brown Create nonlinear solver context 45c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 46c4762a1bSJed Brown 479566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 489566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 499566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 509566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCVPBJACOBI)); 51c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 52c4762a1bSJed Brown Create vector data structures; set function evaluation routine 53c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 54c4762a1bSJed Brown /* 55c4762a1bSJed Brown Note that we form 1 vector from scratch and then duplicate as needed. 56c4762a1bSJed Brown */ 579566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 589566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,PETSC_DECIDE,n)); 599566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&F)); 629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&U)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown Set function evaluation routine and vector 66c4762a1bSJed Brown */ 679566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)F)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 70c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 71c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&J)); 749566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n)); 759566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 769566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J,3,NULL)); 779566063dSJacob Faibussowitsch PetscCall(MatSetVariableBlockSizes(J,3,lens)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* 80c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 81c4762a1bSJed Brown routine. User can override with: 82c4762a1bSJed Brown -snes_fd : default finite differencing approximation of Jacobian 83c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 84c4762a1bSJed Brown (unless user explicitly sets preconditioner) 85c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 86c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 87c4762a1bSJed Brown products within Newton-Krylov method 88c4762a1bSJed Brown */ 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,NULL)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 93c4762a1bSJed Brown Customize nonlinear solver; set runtime options 94c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* 97c4762a1bSJed Brown Set names for some vectors to facilitate monitoring (optional) 98c4762a1bSJed Brown */ 999566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x,"Approximate Solution")); 1009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U,"Exact Solution")); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 104c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 105c4762a1bSJed Brown */ 1069566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* 109c4762a1bSJed Brown Print parameters used for convergence testing (optional) ... just 110c4762a1bSJed Brown to demonstrate this routine; this information is also printed with 111c4762a1bSJed Brown the option -snes_view 112c4762a1bSJed Brown */ 1139566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf)); 11463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%" PetscInt_FMT ", maxf=%" PetscInt_FMT "\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117c4762a1bSJed Brown Initialize application: 118c4762a1bSJed Brown Store right-hand-side of PDE and exact solution 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120c4762a1bSJed Brown 121c4762a1bSJed Brown xp = 0.0; 122c4762a1bSJed Brown for (i=0; i<n; i++) { 123c4762a1bSJed Brown v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 1249566063dSJacob Faibussowitsch PetscCall(VecSetValues(F,1,&i,&v,INSERT_VALUES)); 125c4762a1bSJed Brown v = xp*xp*xp; 1269566063dSJacob Faibussowitsch PetscCall(VecSetValues(U,1,&i,&v,INSERT_VALUES)); 127c4762a1bSJed Brown xp += h; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 132c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 133c4762a1bSJed Brown /* 134c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 135c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 136c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 137c4762a1bSJed Brown this vector to zero by calling VecSet(). 138c4762a1bSJed Brown */ 1399566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(x)); 1409566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 1419566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 14263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"number of SNES iterations = %" PetscInt_FMT "\n\n",its)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145c4762a1bSJed Brown Check solution and clean up 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* 149c4762a1bSJed Brown Check the error 150c4762a1bSJed Brown */ 1519566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,none,U)); 1529566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&norm)); 15363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %" PetscInt_FMT "\n",(double)norm,its)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* 156c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 157c4762a1bSJed Brown are no longer needed. 158c4762a1bSJed Brown */ 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); PetscCall(VecDestroy(&r)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); PetscCall(VecDestroy(&F)); 1619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); PetscCall(SNESDestroy(&snes)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 163b122ec5aSJacob Faibussowitsch return 0; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 166c4762a1bSJed Brown /* 167c4762a1bSJed Brown FormInitialGuess - Computes initial guess. 168c4762a1bSJed Brown 169c4762a1bSJed Brown Input/Output Parameter: 170c4762a1bSJed Brown . x - the solution vector 171c4762a1bSJed Brown */ 172c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x) 173c4762a1bSJed Brown { 174c4762a1bSJed Brown PetscScalar pfive = .50; 1759566063dSJacob Faibussowitsch PetscCall(VecSet(x,pfive)); 176c4762a1bSJed Brown return 0; 177c4762a1bSJed Brown } 178c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 179c4762a1bSJed Brown /* 180c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 181c4762a1bSJed Brown 182c4762a1bSJed Brown Input Parameters: 183c4762a1bSJed Brown . snes - the SNES context 184c4762a1bSJed Brown . x - input vector 185c4762a1bSJed Brown . ctx - optional user-defined context, as set by SNESSetFunction() 186c4762a1bSJed Brown 187c4762a1bSJed Brown Output Parameter: 188c4762a1bSJed Brown . f - function vector 189c4762a1bSJed Brown 190c4762a1bSJed Brown Note: 191c4762a1bSJed Brown The user-defined context can contain any application-specific data 192c4762a1bSJed Brown needed for the function evaluation (such as various parameters, work 193c4762a1bSJed Brown vectors, and grid information). In this program the context is just 194c4762a1bSJed Brown a vector containing the right-hand-side of the discretized PDE. 195c4762a1bSJed Brown */ 196c4762a1bSJed Brown 197c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx) 198c4762a1bSJed Brown { 199c4762a1bSJed Brown Vec g = (Vec)ctx; 200c4762a1bSJed Brown const PetscScalar *xx,*gg; 201c4762a1bSJed Brown PetscScalar *ff,d; 202c4762a1bSJed Brown PetscInt i,n; 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* 205c4762a1bSJed Brown Get pointers to vector data. 206c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 207c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 208c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 209c4762a1bSJed Brown the array. 210c4762a1bSJed Brown */ 2119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetArray(f,&ff)); 2139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(g,&gg)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* 216c4762a1bSJed Brown Compute function 217c4762a1bSJed Brown */ 2189566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&n)); 219c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 220c4762a1bSJed Brown ff[0] = xx[0]; 221c4762a1bSJed 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]; 222c4762a1bSJed Brown ff[n-1] = xx[n-1] - 1.0; 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* 225c4762a1bSJed Brown Restore vectors 226c4762a1bSJed Brown */ 2279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 2289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f,&ff)); 2299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(g,&gg)); 230c4762a1bSJed Brown return 0; 231c4762a1bSJed Brown } 232c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 233c4762a1bSJed Brown /* 234c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 235c4762a1bSJed Brown 236c4762a1bSJed Brown Input Parameters: 237c4762a1bSJed Brown . snes - the SNES context 238c4762a1bSJed Brown . x - input vector 239c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 240c4762a1bSJed Brown 241c4762a1bSJed Brown Output Parameters: 242c4762a1bSJed Brown . jac - Jacobian matrix 243c4762a1bSJed Brown . B - optionally different preconditioning matrix 244c4762a1bSJed Brown 245c4762a1bSJed Brown */ 246c4762a1bSJed Brown 247c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 248c4762a1bSJed Brown { 249c4762a1bSJed Brown const PetscScalar *xx; 250c4762a1bSJed Brown PetscScalar A[3],d; 251c4762a1bSJed Brown PetscInt i,n,j[3]; 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown Get pointer to vector data 255c4762a1bSJed Brown */ 2569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* 259c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 260c4762a1bSJed Brown - Note that in this case we set all elements for a particular 261c4762a1bSJed Brown row at once. 262c4762a1bSJed Brown */ 2639566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&n)); 264c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* 267c4762a1bSJed Brown Interior grid points 268c4762a1bSJed Brown */ 269c4762a1bSJed Brown for (i=1; i<n-1; i++) { 270c4762a1bSJed Brown j[0] = i - 1; j[1] = i; j[2] = i + 1; 271c4762a1bSJed Brown A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i]; 2729566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES)); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* 276c4762a1bSJed Brown Boundary points 277c4762a1bSJed Brown */ 278c4762a1bSJed Brown i = 0; A[0] = 1.0; 279c4762a1bSJed Brown 2809566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES)); 281c4762a1bSJed Brown 282c4762a1bSJed Brown i = n-1; A[0] = 1.0; 283c4762a1bSJed Brown 2849566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,1,&i,A,INSERT_VALUES)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* 287c4762a1bSJed Brown Restore vector 288c4762a1bSJed Brown */ 2899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* 292c4762a1bSJed Brown Assemble matrix 293c4762a1bSJed Brown */ 2949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 296c4762a1bSJed Brown if (jac != B) { 2979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 299c4762a1bSJed Brown } 300c4762a1bSJed Brown return 0; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303c4762a1bSJed Brown /*TEST 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown args: -snes_monitor_short -snes_view -ksp_monitor 307c4762a1bSJed Brown 308c4762a1bSJed Brown # this is just a test for SNESKSPTRASPOSEONLY and KSPSolveTranspose to behave properly 309c4762a1bSJed Brown # the solution is wrong on purpose 310c4762a1bSJed Brown test: 311c4762a1bSJed Brown requires: !single !complex 312c4762a1bSJed Brown suffix: transpose_only 313c4762a1bSJed 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 314c4762a1bSJed Brown 315c4762a1bSJed Brown TEST*/ 316