1c4762a1bSJed Brown static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel.\n\ 2c4762a1bSJed Brown This example employs a user-defined monitoring routine and optionally a user-defined\n\ 3c4762a1bSJed Brown routine to check candidate iterates produced by line search routines.\n\ 4c4762a1bSJed Brown The command line options include:\n\ 5c4762a1bSJed Brown -pre_check_iterates : activate checking of iterates\n\ 6c4762a1bSJed Brown -post_check_iterates : activate checking of iterates\n\ 7c4762a1bSJed Brown -check_tol <tol>: set tolerance for iterate checking\n\ 8c4762a1bSJed Brown -user_precond : activate a (trivial) user-defined preconditioner\n\n"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Include "petscdm.h" so that we can use data management objects (DMs) 12c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 13c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 14c4762a1bSJed Brown file automatically includes: 15c4762a1bSJed Brown petscsys.h - base PETSc routines 16c4762a1bSJed Brown petscvec.h - vectors 17c4762a1bSJed Brown petscmat.h - matrices 18c4762a1bSJed Brown petscis.h - index sets 19c4762a1bSJed Brown petscksp.h - Krylov subspace methods 20c4762a1bSJed Brown petscviewer.h - viewers 21c4762a1bSJed Brown petscpc.h - preconditioners 22c4762a1bSJed Brown petscksp.h - linear solvers 23c4762a1bSJed Brown */ 24c4762a1bSJed Brown 25c4762a1bSJed Brown #include <petscdm.h> 26c4762a1bSJed Brown #include <petscdmda.h> 27c4762a1bSJed Brown #include <petscsnes.h> 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* 30c4762a1bSJed Brown User-defined routines. 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 33c4762a1bSJed Brown PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 34c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec); 35c4762a1bSJed Brown PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *); 36c4762a1bSJed Brown PetscErrorCode PreCheck(SNESLineSearch, Vec, Vec, PetscBool *, void *); 37c4762a1bSJed Brown PetscErrorCode PostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 38c4762a1bSJed Brown PetscErrorCode PostSetSubKSP(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 39c4762a1bSJed Brown PetscErrorCode MatrixFreePreconditioner(PC, Vec, Vec); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown User-defined application context 43c4762a1bSJed Brown */ 44*11486bccSBarry Smith typedef struct 45*11486bccSBarry Smith { 46c4762a1bSJed Brown DM da; /* distributed array */ 47c4762a1bSJed Brown Vec F; /* right-hand-side of PDE */ 48c4762a1bSJed Brown PetscMPIInt rank; /* rank of processor */ 49c4762a1bSJed Brown PetscMPIInt size; /* size of communicator */ 50c4762a1bSJed Brown PetscReal h; /* mesh spacing */ 51c4762a1bSJed Brown PetscBool sjerr; /* if or not to test jacobian domain error */ 52c4762a1bSJed Brown } ApplicationCtx; 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* 55c4762a1bSJed Brown User-defined context for monitoring 56c4762a1bSJed Brown */ 57*11486bccSBarry Smith typedef struct 58*11486bccSBarry Smith { 59c4762a1bSJed Brown PetscViewer viewer; 60c4762a1bSJed Brown } MonitorCtx; 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* 63c4762a1bSJed Brown User-defined context for checking candidate iterates that are 64c4762a1bSJed Brown determined by line search methods 65c4762a1bSJed Brown */ 66*11486bccSBarry Smith typedef struct 67*11486bccSBarry Smith { 68c4762a1bSJed Brown Vec last_step; /* previous iterate */ 69c4762a1bSJed Brown PetscReal tolerance; /* tolerance for changes between successive iterates */ 70c4762a1bSJed Brown ApplicationCtx *user; 71c4762a1bSJed Brown } StepCheckCtx; 72c4762a1bSJed Brown 73*11486bccSBarry Smith typedef struct 74*11486bccSBarry Smith { 75c4762a1bSJed Brown PetscInt its0; /* num of prevous outer KSP iterations */ 76c4762a1bSJed Brown } SetSubKSPCtx; 77c4762a1bSJed Brown 78c4762a1bSJed Brown int main(int argc, char **argv) 79c4762a1bSJed Brown { 80c4762a1bSJed Brown SNES snes; /* SNES context */ 81c4762a1bSJed Brown SNESLineSearch linesearch; /* SNESLineSearch context */ 82c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 83c4762a1bSJed Brown ApplicationCtx ctx; /* user-defined context */ 84c4762a1bSJed Brown Vec x, r, U, F; /* vectors */ 85c4762a1bSJed Brown MonitorCtx monP; /* monitoring context */ 86c4762a1bSJed Brown StepCheckCtx checkP; /* step-checking context */ 87c4762a1bSJed Brown SetSubKSPCtx checkP1; 88c4762a1bSJed Brown PetscBool pre_check, post_check, post_setsubksp; /* flag indicating whether we're checking candidate iterates */ 89c4762a1bSJed Brown PetscScalar xp, *FF, *UU, none = -1.0; 90c4762a1bSJed Brown PetscInt its, N = 5, i, maxit, maxf, xs, xm; 91c4762a1bSJed Brown PetscReal abstol, rtol, stol, norm; 92c0558f20SBarry Smith PetscBool flg, viewinitial = PETSC_FALSE; 93c4762a1bSJed Brown 94327415f7SBarry Smith PetscFunctionBeginUser; 959566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &ctx.rank)); 979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &ctx.size)); 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL)); 99c4762a1bSJed Brown ctx.h = 1.0 / (N - 1); 100c4762a1bSJed Brown ctx.sjerr = PETSC_FALSE; 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_jacobian_domain_error", &ctx.sjerr, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_initial", &viewinitial, NULL)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105c4762a1bSJed Brown Create nonlinear solver context 106c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 111c4762a1bSJed Brown Create vector data structures; set function evaluation routine 112c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 116c4762a1bSJed Brown */ 1179566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, N, 1, 1, NULL, &ctx.da)); 1189566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(ctx.da)); 1199566063dSJacob Faibussowitsch PetscCall(DMSetUp(ctx.da)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown Extract global and local vectors from DMDA; then duplicate for remaining 123c4762a1bSJed Brown vectors that are the same types 124c4762a1bSJed Brown */ 1259566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(ctx.da, &x)); 1269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 127*11486bccSBarry Smith PetscCall(VecDuplicate(x, &F)); 128*11486bccSBarry Smith ctx.F = F; 1299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &U)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* 132c4762a1bSJed Brown Set function evaluation routine and vector. Whenever the nonlinear 133c4762a1bSJed Brown solver needs to compute the nonlinear function, it will call this 134c4762a1bSJed Brown routine. 135c4762a1bSJed Brown - Note that the final routine argument is the user-defined 136c4762a1bSJed Brown context that provides application-specific data for the 137c4762a1bSJed Brown function evaluation routine. 138c4762a1bSJed Brown */ 1399566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, &ctx)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 143c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144c4762a1bSJed Brown 1459566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 1469566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, N, N)); 1479566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 1489566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(J, 3, NULL, 3, NULL)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* 152c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 153c4762a1bSJed Brown routine. Whenever the nonlinear solver needs to compute the 154c4762a1bSJed Brown Jacobian matrix, it will call this routine. 155c4762a1bSJed Brown - Note that the final routine argument is the user-defined 156c4762a1bSJed Brown context that provides application-specific data for the 157c4762a1bSJed Brown Jacobian evaluation routine. 158c4762a1bSJed Brown */ 1599566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &ctx)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* 162c4762a1bSJed Brown Optionally allow user-provided preconditioner 163c4762a1bSJed Brown */ 1649566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-user_precond", &flg)); 165c4762a1bSJed Brown if (flg) { 166c4762a1bSJed Brown KSP ksp; 167c4762a1bSJed Brown PC pc; 1689566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 1699566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1709566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCSHELL)); 1719566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(pc, MatrixFreePreconditioner)); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175c4762a1bSJed Brown Customize nonlinear solver; set runtime options 176c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown Set an optional user-defined monitoring routine 180c4762a1bSJed Brown */ 1819566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, 0, 0, 0, 400, 400, &monP.viewer)); 1829566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(snes, Monitor, &monP, 0)); 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* 185c4762a1bSJed Brown Set names for some vectors to facilitate monitoring (optional) 186c4762a1bSJed Brown */ 1879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution")); 1889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution")); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* 191c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 192c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 193c4762a1bSJed Brown */ 1949566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* 197c4762a1bSJed Brown Set an optional user-defined routine to check the validity of candidate 198c4762a1bSJed Brown iterates that are determined by line search methods 199c4762a1bSJed Brown */ 2009566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 2019566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-post_check_iterates", &post_check)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown if (post_check) { 2049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating post step checking routine\n")); 2059566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPostCheck(linesearch, PostCheck, &checkP)); 2069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &(checkP.last_step))); 207c4762a1bSJed Brown 208c4762a1bSJed Brown checkP.tolerance = 1.0; 209c4762a1bSJed Brown checkP.user = &ctx; 210c4762a1bSJed Brown 2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-check_tol", &checkP.tolerance, NULL)); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 2149566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-post_setsubksp", &post_setsubksp)); 215c4762a1bSJed Brown if (post_setsubksp) { 2169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating post setsubksp\n")); 2179566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPostCheck(linesearch, PostSetSubKSP, &checkP1)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-pre_check_iterates", &pre_check)); 221c4762a1bSJed Brown if (pre_check) { 2229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating pre step checking routine\n")); 2239566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(linesearch, PreCheck, &checkP)); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* 227c4762a1bSJed Brown Print parameters used for convergence testing (optional) ... just 228c4762a1bSJed Brown to demonstrate this routine; this information is also printed with 229c4762a1bSJed Brown the option -snes_view 230c4762a1bSJed Brown */ 2319566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf)); 23263a3b9bcSJacob 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)); 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 235c4762a1bSJed Brown Initialize application: 236c4762a1bSJed Brown Store right-hand-side of PDE and exact solution 237c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 238c4762a1bSJed Brown 239c4762a1bSJed Brown /* 240c4762a1bSJed Brown Get local grid boundaries (for 1-dimensional DMDA): 241c4762a1bSJed Brown xs, xm - starting grid index, width of local grid (no ghost points) 242c4762a1bSJed Brown */ 2439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(ctx.da, &xs, NULL, NULL, &xm, NULL, NULL)); 244c4762a1bSJed Brown 245c4762a1bSJed Brown /* 246c4762a1bSJed Brown Get pointers to vector data 247c4762a1bSJed Brown */ 2489566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(ctx.da, F, &FF)); 2499566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(ctx.da, U, &UU)); 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* 252c4762a1bSJed Brown Compute local vector entries 253c4762a1bSJed Brown */ 254c4762a1bSJed Brown xp = ctx.h * xs; 255c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 256c4762a1bSJed Brown FF[i] = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */ 257c4762a1bSJed Brown UU[i] = xp * xp * xp; 258c4762a1bSJed Brown xp += ctx.h; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* 262c4762a1bSJed Brown Restore vectors 263c4762a1bSJed Brown */ 2649566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(ctx.da, F, &FF)); 2659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(ctx.da, U, &UU)); 266c0558f20SBarry Smith if (viewinitial) { 2679566063dSJacob Faibussowitsch PetscCall(VecView(U, 0)); 2689566063dSJacob Faibussowitsch PetscCall(VecView(F, 0)); 269c0558f20SBarry Smith } 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 272c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 273c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* 276c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 277c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 278c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 279c4762a1bSJed Brown this vector to zero by calling VecSet(). 280c4762a1bSJed Brown */ 2819566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(x)); 2829566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 2839566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 28463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 287c4762a1bSJed Brown Check solution and clean up 288c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* 291c4762a1bSJed Brown Check the error 292c4762a1bSJed Brown */ 2939566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, none, U)); 2949566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 29563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g Iterations %" PetscInt_FMT "\n", (double)norm, its)); 296c4762a1bSJed Brown if (ctx.sjerr) { 297c4762a1bSJed Brown SNESType snestype; 2989566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes, &snestype)); 2999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES Type %s\n", snestype)); 300c4762a1bSJed Brown } 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* 303c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 304c4762a1bSJed Brown are no longer needed. 305c4762a1bSJed Brown */ 3069566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&monP.viewer)); 3079566063dSJacob Faibussowitsch if (post_check) PetscCall(VecDestroy(&checkP.last_step)); 3089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 3129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 3139566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3149566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx.da)); 3159566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 316b122ec5aSJacob Faibussowitsch return 0; 317c4762a1bSJed Brown } 318c4762a1bSJed Brown 319c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 320c4762a1bSJed Brown /* 321c4762a1bSJed Brown FormInitialGuess - Computes initial guess. 322c4762a1bSJed Brown 323c4762a1bSJed Brown Input/Output Parameter: 324c4762a1bSJed Brown . x - the solution vector 325c4762a1bSJed Brown */ 326c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x) 327c4762a1bSJed Brown { 328c4762a1bSJed Brown PetscScalar pfive = .50; 329c4762a1bSJed Brown 330c4762a1bSJed Brown PetscFunctionBeginUser; 3319566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 332c4762a1bSJed Brown PetscFunctionReturn(0); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 335c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 336c4762a1bSJed Brown /* 337c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 338c4762a1bSJed Brown 339c4762a1bSJed Brown Input Parameters: 340c4762a1bSJed Brown . snes - the SNES context 341c4762a1bSJed Brown . x - input vector 342c4762a1bSJed Brown . ctx - optional user-defined context, as set by SNESSetFunction() 343c4762a1bSJed Brown 344c4762a1bSJed Brown Output Parameter: 345c4762a1bSJed Brown . f - function vector 346c4762a1bSJed Brown 347c4762a1bSJed Brown Note: 348c4762a1bSJed Brown The user-defined context can contain any application-specific 349c4762a1bSJed Brown data needed for the function evaluation. 350c4762a1bSJed Brown */ 351c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *ctx) 352c4762a1bSJed Brown { 353c4762a1bSJed Brown ApplicationCtx *user = (ApplicationCtx *)ctx; 354c4762a1bSJed Brown DM da = user->da; 35506cf5b03SBarry Smith PetscScalar *ff, d; 35606cf5b03SBarry Smith const PetscScalar *xx, *FF; 357c4762a1bSJed Brown PetscInt i, M, xs, xm; 358c4762a1bSJed Brown Vec xlocal; 359c4762a1bSJed Brown 360c4762a1bSJed Brown PetscFunctionBeginUser; 3619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &xlocal)); 362c4762a1bSJed Brown /* 363c4762a1bSJed Brown Scatter ghost points to local vector, using the 2-step process 364c4762a1bSJed Brown DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 365c4762a1bSJed Brown By placing code between these two statements, computations can 366c4762a1bSJed Brown be done while messages are in transition. 367c4762a1bSJed Brown */ 3689566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, x, INSERT_VALUES, xlocal)); 3699566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, x, INSERT_VALUES, xlocal)); 370c4762a1bSJed Brown 371c4762a1bSJed Brown /* 372c4762a1bSJed Brown Get pointers to vector data. 373c4762a1bSJed Brown - The vector xlocal includes ghost point; the vectors x and f do 374c4762a1bSJed Brown NOT include ghost points. 375c4762a1bSJed Brown - Using DMDAVecGetArray() allows accessing the values using global ordering 376c4762a1bSJed Brown */ 37795b2e421SBarry Smith PetscCall(DMDAVecGetArrayRead(da, xlocal, (void *)&xx)); 3789566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, f, &ff)); 37995b2e421SBarry Smith PetscCall(DMDAVecGetArrayRead(da, user->F, (void *)&FF)); 380c4762a1bSJed Brown 381c4762a1bSJed Brown /* 382c4762a1bSJed Brown Get local grid boundaries (for 1-dimensional DMDA): 383c4762a1bSJed Brown xs, xm - starting grid index, width of local grid (no ghost points) 384c4762a1bSJed Brown */ 3859566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 3869566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 387c4762a1bSJed Brown 388c4762a1bSJed Brown /* 389c4762a1bSJed Brown Set function values for boundary points; define local interior grid point range: 390c4762a1bSJed Brown xsi - starting interior grid index 391c4762a1bSJed Brown xei - ending interior grid index 392c4762a1bSJed Brown */ 393c4762a1bSJed Brown if (xs == 0) { /* left boundary */ 394c4762a1bSJed Brown ff[0] = xx[0]; 395*11486bccSBarry Smith xs++; 396*11486bccSBarry Smith xm--; 397c4762a1bSJed Brown } 398c4762a1bSJed Brown if (xs + xm == M) { /* right boundary */ 399c4762a1bSJed Brown ff[xs + xm - 1] = xx[xs + xm - 1] - 1.0; 400c4762a1bSJed Brown xm--; 401c4762a1bSJed Brown } 402c4762a1bSJed Brown 403c4762a1bSJed Brown /* 404c4762a1bSJed Brown Compute function over locally owned part of the grid (interior points only) 405c4762a1bSJed Brown */ 406c4762a1bSJed Brown d = 1.0 / (user->h * user->h); 407c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) ff[i] = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - FF[i]; 408c4762a1bSJed Brown 409c4762a1bSJed Brown /* 410c4762a1bSJed Brown Restore vectors 411c4762a1bSJed Brown */ 41295b2e421SBarry Smith PetscCall(DMDAVecRestoreArrayRead(da, xlocal, (void *)&xx)); 4139566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, f, &ff)); 41495b2e421SBarry Smith PetscCall(DMDAVecRestoreArrayRead(da, user->F, (void *)&FF)); 4159566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &xlocal)); 416c4762a1bSJed Brown PetscFunctionReturn(0); 417c4762a1bSJed Brown } 418c4762a1bSJed Brown 419c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 420c4762a1bSJed Brown /* 421c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 422c4762a1bSJed Brown 423c4762a1bSJed Brown Input Parameters: 424c4762a1bSJed Brown . snes - the SNES context 425c4762a1bSJed Brown . x - input vector 426c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 427c4762a1bSJed Brown 428c4762a1bSJed Brown Output Parameters: 429c4762a1bSJed Brown . jac - Jacobian matrix 430c4762a1bSJed Brown . B - optionally different preconditioning matrix 431c4762a1bSJed Brown . flag - flag indicating matrix structure 432c4762a1bSJed Brown */ 433c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *ctx) 434c4762a1bSJed Brown { 435c4762a1bSJed Brown ApplicationCtx *user = (ApplicationCtx *)ctx; 436c4762a1bSJed Brown PetscScalar *xx, d, A[3]; 437c4762a1bSJed Brown PetscInt i, j[3], M, xs, xm; 438c4762a1bSJed Brown DM da = user->da; 439c4762a1bSJed Brown 440c4762a1bSJed Brown PetscFunctionBeginUser; 441c4762a1bSJed Brown if (user->sjerr) { 4429566063dSJacob Faibussowitsch PetscCall(SNESSetJacobianDomainError(snes)); 443c4762a1bSJed Brown PetscFunctionReturn(0); 444c4762a1bSJed Brown } 445c4762a1bSJed Brown /* 446c4762a1bSJed Brown Get pointer to vector data 447c4762a1bSJed Brown */ 4489566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, x, &xx)); 4499566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 450c4762a1bSJed Brown 451c4762a1bSJed Brown /* 452c4762a1bSJed Brown Get range of locally owned matrix 453c4762a1bSJed Brown */ 4549566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 455c4762a1bSJed Brown 456c4762a1bSJed Brown /* 457c4762a1bSJed Brown Determine starting and ending local indices for interior grid points. 458c4762a1bSJed Brown Set Jacobian entries for boundary points. 459c4762a1bSJed Brown */ 460c4762a1bSJed Brown 461c4762a1bSJed Brown if (xs == 0) { /* left boundary */ 462*11486bccSBarry Smith i = 0; 463*11486bccSBarry Smith A[0] = 1.0; 464c4762a1bSJed Brown 4659566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 466*11486bccSBarry Smith xs++; 467*11486bccSBarry Smith xm--; 468c4762a1bSJed Brown } 469c4762a1bSJed Brown if (xs + xm == M) { /* right boundary */ 470c4762a1bSJed Brown i = M - 1; 471c4762a1bSJed Brown A[0] = 1.0; 4729566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 473c4762a1bSJed Brown xm--; 474c4762a1bSJed Brown } 475c4762a1bSJed Brown 476c4762a1bSJed Brown /* 477c4762a1bSJed Brown Interior grid points 478c4762a1bSJed Brown - Note that in this case we set all elements for a particular 479c4762a1bSJed Brown row at once. 480c4762a1bSJed Brown */ 481c4762a1bSJed Brown d = 1.0 / (user->h * user->h); 482c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 483*11486bccSBarry Smith j[0] = i - 1; 484*11486bccSBarry Smith j[1] = i; 485*11486bccSBarry Smith j[2] = i + 1; 486*11486bccSBarry Smith A[0] = A[2] = d; 487*11486bccSBarry Smith A[1] = -2.0 * d + 2.0 * xx[i]; 4889566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 3, j, A, INSERT_VALUES)); 489c4762a1bSJed Brown } 490c4762a1bSJed Brown 491c4762a1bSJed Brown /* 492c4762a1bSJed Brown Assemble matrix, using the 2-step process: 493c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 494c4762a1bSJed Brown By placing code between these two statements, computations can be 495c4762a1bSJed Brown done while messages are in transition. 496c4762a1bSJed Brown 497c4762a1bSJed Brown Also, restore vector. 498c4762a1bSJed Brown */ 499c4762a1bSJed Brown 5009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 5019566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, x, &xx)); 5029566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 503c4762a1bSJed Brown 504c4762a1bSJed Brown PetscFunctionReturn(0); 505c4762a1bSJed Brown } 506c4762a1bSJed Brown 507c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 508c4762a1bSJed Brown /* 509c4762a1bSJed Brown Monitor - Optional user-defined monitoring routine that views the 510c4762a1bSJed Brown current iterate with an x-window plot. Set by SNESMonitorSet(). 511c4762a1bSJed Brown 512c4762a1bSJed Brown Input Parameters: 513c4762a1bSJed Brown snes - the SNES context 514c4762a1bSJed Brown its - iteration number 515c4762a1bSJed Brown norm - 2-norm function value (may be estimated) 516c4762a1bSJed Brown ctx - optional user-defined context for private data for the 517c4762a1bSJed Brown monitor routine, as set by SNESMonitorSet() 518c4762a1bSJed Brown 519c4762a1bSJed Brown Note: 520c4762a1bSJed Brown See the manpage for PetscViewerDrawOpen() for useful runtime options, 521c4762a1bSJed Brown such as -nox to deactivate all x-window output. 522c4762a1bSJed Brown */ 523c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, void *ctx) 524c4762a1bSJed Brown { 525c4762a1bSJed Brown MonitorCtx *monP = (MonitorCtx *)ctx; 526c4762a1bSJed Brown Vec x; 527c4762a1bSJed Brown 528c4762a1bSJed Brown PetscFunctionBeginUser; 52963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iter = %" PetscInt_FMT ",SNES Function norm %g\n", its, (double)fnorm)); 5309566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &x)); 5319566063dSJacob Faibussowitsch PetscCall(VecView(x, monP->viewer)); 532c4762a1bSJed Brown PetscFunctionReturn(0); 533c4762a1bSJed Brown } 534c4762a1bSJed Brown 535c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 536c4762a1bSJed Brown /* 537c4762a1bSJed Brown PreCheck - Optional user-defined routine that checks the validity of 538c4762a1bSJed Brown candidate steps of a line search method. Set by SNESLineSearchSetPreCheck(). 539c4762a1bSJed Brown 540c4762a1bSJed Brown Input Parameters: 541c4762a1bSJed Brown snes - the SNES context 542c4762a1bSJed Brown xcurrent - current solution 543c4762a1bSJed Brown y - search direction and length 544c4762a1bSJed Brown 545c4762a1bSJed Brown Output Parameters: 546c4762a1bSJed Brown y - proposed step (search direction and length) (possibly changed) 547c4762a1bSJed Brown changed_y - tells if the step has changed or not 548c4762a1bSJed Brown */ 549c4762a1bSJed Brown PetscErrorCode PreCheck(SNESLineSearch linesearch, Vec xcurrent, Vec y, PetscBool *changed_y, void *ctx) 550c4762a1bSJed Brown { 551c4762a1bSJed Brown PetscFunctionBeginUser; 552c4762a1bSJed Brown *changed_y = PETSC_FALSE; 553c4762a1bSJed Brown PetscFunctionReturn(0); 554c4762a1bSJed Brown } 555c4762a1bSJed Brown 556c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 557c4762a1bSJed Brown /* 558c4762a1bSJed Brown PostCheck - Optional user-defined routine that checks the validity of 559c4762a1bSJed Brown candidate steps of a line search method. Set by SNESLineSearchSetPostCheck(). 560c4762a1bSJed Brown 561c4762a1bSJed Brown Input Parameters: 562c4762a1bSJed Brown snes - the SNES context 563c4762a1bSJed Brown ctx - optional user-defined context for private data for the 564c4762a1bSJed Brown monitor routine, as set by SNESLineSearchSetPostCheck() 565c4762a1bSJed Brown xcurrent - current solution 566c4762a1bSJed Brown y - search direction and length 567c4762a1bSJed Brown x - the new candidate iterate 568c4762a1bSJed Brown 569c4762a1bSJed Brown Output Parameters: 570c4762a1bSJed Brown y - proposed step (search direction and length) (possibly changed) 571c4762a1bSJed Brown x - current iterate (possibly modified) 572c4762a1bSJed Brown 573c4762a1bSJed Brown */ 574c4762a1bSJed Brown PetscErrorCode PostCheck(SNESLineSearch linesearch, Vec xcurrent, Vec y, Vec x, PetscBool *changed_y, PetscBool *changed_x, void *ctx) 575c4762a1bSJed Brown { 576c4762a1bSJed Brown PetscInt i, iter, xs, xm; 577c4762a1bSJed Brown StepCheckCtx *check; 578c4762a1bSJed Brown ApplicationCtx *user; 579c4762a1bSJed Brown PetscScalar *xa, *xa_last, tmp; 580c4762a1bSJed Brown PetscReal rdiff; 581c4762a1bSJed Brown DM da; 582c4762a1bSJed Brown SNES snes; 583c4762a1bSJed Brown 584c4762a1bSJed Brown PetscFunctionBeginUser; 585c4762a1bSJed Brown *changed_x = PETSC_FALSE; 586c4762a1bSJed Brown *changed_y = PETSC_FALSE; 587c4762a1bSJed Brown 5889566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 589c4762a1bSJed Brown check = (StepCheckCtx *)ctx; 590c4762a1bSJed Brown user = check->user; 5919566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 592c4762a1bSJed Brown 593c4762a1bSJed Brown /* iteration 1 indicates we are working on the second iteration */ 594c4762a1bSJed Brown if (iter > 0) { 595c4762a1bSJed Brown da = user->da; 59663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Checking candidate step at iteration %" PetscInt_FMT " with tolerance %g\n", iter, (double)check->tolerance)); 597c4762a1bSJed Brown 598c4762a1bSJed Brown /* Access local array data */ 5999566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, check->last_step, &xa_last)); 6009566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, x, &xa)); 6019566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 602c4762a1bSJed Brown 603c4762a1bSJed Brown /* 604c4762a1bSJed Brown If we fail the user-defined check for validity of the candidate iterate, 605c4762a1bSJed Brown then modify the iterate as we like. (Note that the particular modification 606c4762a1bSJed Brown below is intended simply to demonstrate how to manipulate this data, not 607c4762a1bSJed Brown as a meaningful or appropriate choice.) 608c4762a1bSJed Brown */ 609c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 610c4762a1bSJed Brown if (!PetscAbsScalar(xa[i])) rdiff = 2 * check->tolerance; 611c4762a1bSJed Brown else rdiff = PetscAbsScalar((xa[i] - xa_last[i]) / xa[i]); 612c4762a1bSJed Brown if (rdiff > check->tolerance) { 613c4762a1bSJed Brown tmp = xa[i]; 614c4762a1bSJed Brown xa[i] = .5 * (xa[i] + xa_last[i]); 615c4762a1bSJed Brown *changed_x = PETSC_TRUE; 61663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Altering entry %" PetscInt_FMT ": x=%g, x_last=%g, diff=%g, x_new=%g\n", i, (double)PetscAbsScalar(tmp), (double)PetscAbsScalar(xa_last[i]), (double)rdiff, (double)PetscAbsScalar(xa[i]))); 617c4762a1bSJed Brown } 618c4762a1bSJed Brown } 6199566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, check->last_step, &xa_last)); 6209566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, x, &xa)); 621c4762a1bSJed Brown } 6229566063dSJacob Faibussowitsch PetscCall(VecCopy(x, check->last_step)); 623c4762a1bSJed Brown PetscFunctionReturn(0); 624c4762a1bSJed Brown } 625c4762a1bSJed Brown 626c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 627c4762a1bSJed Brown /* 628c4762a1bSJed Brown PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used 629c4762a1bSJed Brown e.g, 630c4762a1bSJed Brown mpiexec -n 8 ./ex3 -nox -n 10000 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 4 -sub_ksp_type gmres -sub_ksp_max_it 3 -post_setsubksp -sub_ksp_rtol 1.e-16 631c4762a1bSJed Brown Set by SNESLineSearchSetPostCheck(). 632c4762a1bSJed Brown 633c4762a1bSJed Brown Input Parameters: 634c4762a1bSJed Brown linesearch - the LineSearch context 635c4762a1bSJed Brown xcurrent - current solution 636c4762a1bSJed Brown y - search direction and length 637c4762a1bSJed Brown x - the new candidate iterate 638c4762a1bSJed Brown 639c4762a1bSJed Brown Output Parameters: 640c4762a1bSJed Brown y - proposed step (search direction and length) (possibly changed) 641c4762a1bSJed Brown x - current iterate (possibly modified) 642c4762a1bSJed Brown 643c4762a1bSJed Brown */ 644c4762a1bSJed Brown PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch, Vec xcurrent, Vec y, Vec x, PetscBool *changed_y, PetscBool *changed_x, void *ctx) 645c4762a1bSJed Brown { 646c4762a1bSJed Brown SetSubKSPCtx *check; 647c4762a1bSJed Brown PetscInt iter, its, sub_its, maxit; 648c4762a1bSJed Brown KSP ksp, sub_ksp, *sub_ksps; 649c4762a1bSJed Brown PC pc; 650c4762a1bSJed Brown PetscReal ksp_ratio; 651c4762a1bSJed Brown SNES snes; 652c4762a1bSJed Brown 653c4762a1bSJed Brown PetscFunctionBeginUser; 6549566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 655c4762a1bSJed Brown check = (SetSubKSPCtx *)ctx; 6569566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 6579566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 6589566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 6599566063dSJacob Faibussowitsch PetscCall(PCBJacobiGetSubKSP(pc, NULL, NULL, &sub_ksps)); 660c4762a1bSJed Brown sub_ksp = sub_ksps[0]; 6619566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp, &its)); /* outer KSP iteration number */ 6629566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(sub_ksp, &sub_its)); /* inner KSP iteration number */ 663c4762a1bSJed Brown 664c4762a1bSJed Brown if (iter) { 66563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ...PostCheck snes iteration %" PetscInt_FMT ", ksp_it %" PetscInt_FMT " %" PetscInt_FMT ", subksp_it %" PetscInt_FMT "\n", iter, check->its0, its, sub_its)); 666c4762a1bSJed Brown ksp_ratio = ((PetscReal)(its)) / check->its0; 667c4762a1bSJed Brown maxit = (PetscInt)(ksp_ratio * sub_its + 0.5); 668c4762a1bSJed Brown if (maxit < 2) maxit = 2; 6699566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(sub_ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, maxit)); 67063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ...ksp_ratio %g, new maxit %" PetscInt_FMT "\n\n", (double)ksp_ratio, maxit)); 671c4762a1bSJed Brown } 672c4762a1bSJed Brown check->its0 = its; /* save current outer KSP iteration number */ 673c4762a1bSJed Brown PetscFunctionReturn(0); 674c4762a1bSJed Brown } 675c4762a1bSJed Brown 676c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 677c4762a1bSJed Brown /* 678c4762a1bSJed Brown MatrixFreePreconditioner - This routine demonstrates the use of a 679c4762a1bSJed Brown user-provided preconditioner. This code implements just the null 680c4762a1bSJed Brown preconditioner, which of course is not recommended for general use. 681c4762a1bSJed Brown 682c4762a1bSJed Brown Input Parameters: 683c4762a1bSJed Brown + pc - preconditioner 684c4762a1bSJed Brown - x - input vector 685c4762a1bSJed Brown 686c4762a1bSJed Brown Output Parameter: 687c4762a1bSJed Brown . y - preconditioned vector 688c4762a1bSJed Brown */ 689c4762a1bSJed Brown PetscErrorCode MatrixFreePreconditioner(PC pc, Vec x, Vec y) 690c4762a1bSJed Brown { 6919566063dSJacob Faibussowitsch PetscCall(VecCopy(x, y)); 692c4762a1bSJed Brown return 0; 693c4762a1bSJed Brown } 694c4762a1bSJed Brown 695c4762a1bSJed Brown /*TEST 696c4762a1bSJed Brown 697c4762a1bSJed Brown test: 698c4762a1bSJed Brown args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 699c4762a1bSJed Brown 700c4762a1bSJed Brown test: 701c4762a1bSJed Brown suffix: 2 702c4762a1bSJed Brown nsize: 3 703c4762a1bSJed Brown args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 704c4762a1bSJed Brown 705c4762a1bSJed Brown test: 706c4762a1bSJed Brown suffix: 3 707c4762a1bSJed Brown nsize: 2 708c4762a1bSJed Brown args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 709c4762a1bSJed Brown 710c4762a1bSJed Brown test: 711c4762a1bSJed Brown suffix: 4 712c4762a1bSJed Brown args: -nox -pre_check_iterates -post_check_iterates 713c4762a1bSJed Brown 714c4762a1bSJed Brown test: 715c4762a1bSJed Brown suffix: 5 716c4762a1bSJed Brown requires: double !complex !single 717c4762a1bSJed Brown nsize: 2 718c4762a1bSJed Brown args: -nox -snes_test_jacobian -snes_test_jacobian_view 719c4762a1bSJed Brown 720c4762a1bSJed Brown test: 721c4762a1bSJed Brown suffix: 6 722c4762a1bSJed Brown requires: double !complex !single 723c4762a1bSJed Brown nsize: 4 724c4762a1bSJed Brown args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1 725c4762a1bSJed Brown 726c4762a1bSJed Brown test: 727c4762a1bSJed Brown suffix: 7 728c4762a1bSJed Brown requires: double !complex !single 729c4762a1bSJed Brown nsize: 4 730c4762a1bSJed Brown args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1 731c4762a1bSJed Brown 732c4762a1bSJed Brown test: 733c4762a1bSJed Brown suffix: 8 734c4762a1bSJed Brown requires: double !complex !single 735c4762a1bSJed Brown nsize: 4 736c4762a1bSJed Brown args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1 737c4762a1bSJed Brown 738c4762a1bSJed Brown test: 739c4762a1bSJed Brown suffix: 9 740c4762a1bSJed Brown requires: double !complex !single 741c4762a1bSJed Brown nsize: 4 742c4762a1bSJed Brown args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1 743c4762a1bSJed Brown 744c4762a1bSJed Brown test: 745c4762a1bSJed Brown suffix: 10 746c4762a1bSJed Brown requires: double !complex !single 747c4762a1bSJed Brown nsize: 4 748c4762a1bSJed Brown args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1 749c4762a1bSJed Brown 750c4762a1bSJed Brown test: 751c4762a1bSJed Brown suffix: 11 752c4762a1bSJed Brown requires: double !complex !single 753c4762a1bSJed Brown nsize: 4 75457715debSLisandro Dalcin args: -test_jacobian_domain_error -snes_converged_reason -snes_type ms -snes_ms_type m62 -snes_ms_damping 0.9 -snes_check_jacobian_domain_error 1 755c4762a1bSJed Brown 756c0558f20SBarry Smith test: 757c0558f20SBarry Smith suffix: 12 758c0558f20SBarry Smith args: -view_initial 759c0558f20SBarry Smith filter: grep -v "type:" 760c0558f20SBarry Smith 76141ba4c6cSHeeho Park test: 76241ba4c6cSHeeho Park suffix: 13 76341ba4c6cSHeeho Park requires: double !complex !single 76441ba4c6cSHeeho Park nsize: 4 76541ba4c6cSHeeho Park args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1 76641ba4c6cSHeeho Park 767c4762a1bSJed Brown TEST*/ 768