static char help[] = "u`` + u^{2} = f. Different matrices for the Jacobian and the preconditioner.\n\ Demonstrates the use of matrix-free Newton-Krylov methods in conjunction\n\ with a user-provided preconditioner. Input arguments are:\n\ -snes_mf : Use matrix-free Newton methods\n\ -user_precond : Employ a user-defined preconditioner. Used only with\n\ matrix-free methods in this example.\n\n"; /*T Concepts: SNES^different matrices for the Jacobian and preconditioner; Concepts: SNES^matrix-free methods Concepts: SNES^user-provided preconditioner; Concepts: matrix-free methods Concepts: user-provided preconditioner; Processors: 1 T*/ /* 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 */ PetscErrorCode FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); PetscErrorCode FormFunction(SNES,Vec,Vec,void*); PetscErrorCode MatrixFreePreconditioner(void*,Vec,Vec); int main(int argc,char **argv) { SNES snes; /* SNES context */ KSP ksp; /* KSP context */ PC pc; /* PC context */ Vec x,r,F; /* vectors */ Mat J,JPrec; /* Jacobian,preconditioner matrices */ PetscErrorCode ierr; PetscInt it,n = 5,i; PetscMPIInt size; PetscInt *Shistit = 0,Khistl = 200,Shistl = 10; PetscReal h,xp = 0.0,*Khist = 0,*Shist = 0; PetscScalar v,pfive = .5; PetscTruth flg; 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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = VecCreate(PETSC_COMM_SELF,&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 = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create matrix data structures; set Jacobian evaluation routine - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,PETSC_NULL,&J);CHKERRQ(ierr); ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,1,PETSC_NULL,&JPrec);CHKERRQ(ierr); /* Note that in this case we create separate matrices for the Jacobian and preconditioner matrix. Both of these are computed in the routine FormJacobian() */ ierr = SNESSetJacobian(snes,J,JPrec,FormJacobian,0);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Customize nonlinear solver; set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ /* Set preconditioner for matrix-free method */ ierr = PetscOptionsHasName(PETSC_NULL,"-snes_mf",&flg);CHKERRQ(ierr); if (flg) { ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); ierr = PetscOptionsHasName(PETSC_NULL,"-user_precond",&flg);CHKERRQ(ierr); if (flg) { /* user-defined precond */ ierr = PCSetType(pc,PCSHELL);CHKERRQ(ierr); ierr = PCShellSetApply(pc,MatrixFreePreconditioner);CHKERRQ(ierr); } else {ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);} } ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); /* Save all the linear residuals for all the Newton steps; this enables us to retain complete convergence history for printing after the conclusion of SNESSolve(). Alternatively, one could use the monitoring options -snes_monitor -ksp_monitor to see this information during the solver's execution; however, such output during the run distorts performance evaluation data. So, the following is a good option when monitoring code performance, for example when using -log_summary. */ ierr = PetscOptionsHasName(PETSC_NULL,"-rhistory",&flg);CHKERRQ(ierr); if (flg) { ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); ierr = PetscMalloc(Khistl*sizeof(PetscReal),&Khist);CHKERRQ(ierr); ierr = KSPSetResidualHistory(ksp,Khist,Khistl,PETSC_FALSE);CHKERRQ(ierr); ierr = PetscMalloc(Shistl*sizeof(PetscReal),&Shist);CHKERRQ(ierr); ierr = PetscMalloc(Shistl*sizeof(PetscInt),&Shistit);CHKERRQ(ierr); ierr = SNESSetConvergenceHistory(snes,Shist,Shistit,Shistl,PETSC_FALSE);CHKERRQ(ierr); } /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Initialize application: Store right-hand-side of PDE and exact solution - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ xp = 0.0; for (i=0; i