static char help[] = "Newton method to solve u'' + u^{2} = f, sequentially.\n\ This example employs a user-defined monitoring routine.\n\n"; /*T Concepts: SNES^basic uniprocessor example Concepts: SNES^setting a user-defined monitoring routine Processors: 1 T*/ /* Include "petscdraw.h" so that we can use PETSc drawing routines. Include "petscsnes.h" so that we can use SNES solvers. Note that this file automatically includes: petsc.h - base PETSc routines petscvec.h - vectors petscsys.h - system routines petscmat.h - matrices petscis.h - index sets petscksp.h - Krylov subspace methods petscviewer.h - viewers petscpc.h - preconditioners petscksp.h - linear solvers */ #include "petscsnes.h" /* User-defined routines */ extern PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); extern PetscErrorCode FormInitialGuess(Vec); extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void *); /* User-defined context for monitoring */ typedef struct { PetscViewer viewer; } MonitorCtx; #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc,char **argv) { SNES snes; /* SNES context */ Vec x,r,F,U; /* vectors */ Mat J; /* Jacobian matrix */ MonitorCtx monP; /* monitoring context */ PetscErrorCode ierr; PetscInt its,n = 5,i,maxit,maxf; PetscMPIInt size; PetscScalar h,xp,v,none = -1.0; PetscReal abstol,rtol,stol,norm; PetscInitialize(&argc,&argv,(char *)0,help); ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); if (size != 1) SETERRQ(1,"This is a uniprocessor example only!"); ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr); h = 1.0/(n-1); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create nonlinear solver context - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create vector data structures; set function evaluation routine - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Note that we form 1 vector from scratch and then duplicate as needed. */ ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr); ierr = VecSetFromOptions(x);CHKERRQ(ierr); ierr = VecDuplicate(x,&r);CHKERRQ(ierr); ierr = VecDuplicate(x,&F);CHKERRQ(ierr); ierr = VecDuplicate(x,&U);CHKERRQ(ierr); /* Set function evaluation routine and vector */ ierr = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create matrix data structure; set Jacobian evaluation routine - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr); ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); ierr = MatSetFromOptions(J);CHKERRQ(ierr); /* Set Jacobian matrix data structure and default Jacobian evaluation routine. User can override with: -snes_fd : default finite differencing approximation of Jacobian -snes_mf : matrix-free Newton-Krylov method with no preconditioning (unless user explicitly sets preconditioner) -snes_mf_operator : form preconditioning matrix as set by the user, but use matrix-free approx for Jacobian-vector products within Newton-Krylov method */ ierr = SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Customize nonlinear solver; set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Set an optional user-defined monitoring routine */ ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer);CHKERRQ(ierr); ierr = SNESMonitorSet(snes,Monitor,&monP,0);CHKERRQ(ierr); /* Set names for some vectors to facilitate monitoring (optional) */ ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");CHKERRQ(ierr); /* Set SNES/KSP/KSP/PC runtime options, e.g., -snes_view -snes_monitor -ksp_type -pc_type */ ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); /* Print parameters used for convergence testing (optional) ... just to demonstrate this routine; this information is also printed with the option -snes_view */ ierr = SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"atol=%G, rtol=%G, stol=%G, maxit=%D, maxf=%D\n",abstol,rtol,stol,maxit,maxf);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Initialize application: Store right-hand-side of PDE and exact solution - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ xp = 0.0; for (i=0; iviewer);CHKERRQ(ierr); return 0; }