static char help[] = "Solves u`` + u^{2} = f with Newton-like methods. Using\n\ matrix-free techniques with user-provided explicit preconditioner matrix.\n\n"; #include extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); extern PetscErrorCode FormJacobianNoMatrix(SNES,Vec,Mat,Mat,void*); extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); extern PetscErrorCode FormFunctioni(void *,PetscInt,Vec,PetscScalar *); extern PetscErrorCode OtherFunctionForDifferencing(void*,Vec,Vec); extern PetscErrorCode FormInitialGuess(SNES,Vec); extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*); typedef struct { PetscViewer viewer; } MonitorCtx; typedef struct { PetscBool variant; } AppCtx; int main(int argc,char **argv) { SNES snes; /* SNES context */ SNESType type = SNESNEWTONLS; /* default nonlinear solution method */ Vec x,r,F,U; /* vectors */ Mat J,B; /* Jacobian matrix-free, explicit preconditioner */ AppCtx user; /* user-defined work context */ PetscScalar h,xp = 0.0,v; PetscInt its,n = 5,i; PetscBool puremf = PETSC_FALSE; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); PetscCall(PetscOptionsHasName(NULL,NULL,"-variant",&user.variant)); h = 1.0/(n-1); /* Set up data structures */ PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x)); PetscCall(PetscObjectSetName((PetscObject)x,"Approximate Solution")); PetscCall(VecDuplicate(x,&r)); PetscCall(VecDuplicate(x,&F)); PetscCall(VecDuplicate(x,&U)); PetscCall(PetscObjectSetName((PetscObject)U,"Exact Solution")); /* create explicit matrix preconditioner */ PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&B)); /* Store right-hand-side of PDE and exact solution */ for (i=0; ivariant) { PetscCall(MatMFFDSetBase(jac,x,NULL)); } PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); return 0; } PetscErrorCode FormJacobianNoMatrix(SNES snes,Vec x,Mat jac,Mat B,void *dummy) { AppCtx *user = (AppCtx*) dummy; if (user->variant) { PetscCall(MatMFFDSetBase(jac,x,NULL)); } PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); return 0; } /* -------------------- User-defined monitor ----------------------- */ PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *dummy) { MonitorCtx *monP = (MonitorCtx*) dummy; Vec x; MPI_Comm comm; PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); PetscCall(PetscFPrintf(comm,stdout,"iter = %D, SNES Function norm %g \n",its,(double)fnorm)); PetscCall(SNESGetSolution(snes,&x)); PetscCall(VecView(x,monP->viewer)); return 0; } /*TEST test: args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short test: suffix: 2 args: -variant -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short output_file: output/ex7_1.out # uses AIJ matrix to define diagonal matrix for Jacobian preconditioning test: suffix: 3 args: -variant -pc_type jacobi -snes_view -ksp_monitor # uses MATMFFD matrix to define diagonal matrix for Jacobian preconditioning test: suffix: 4 args: -variant -pc_type jacobi -puremf -snes_view -ksp_monitor TEST*/