1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCKSP_H 50a835dfdSSatish Balay #define __PETSCKSP_H 62c8e378dSBarry Smith #include <petscpc.h> 72eac72dbSBarry Smith 8607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 91dbb0a54SBarry Smith 1028ce4d24SBarry Smith /*S 118f6c3df8SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 128f6c3df8SBarry Smith linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 1328ce4d24SBarry Smith 1428ce4d24SBarry Smith Level: beginner 1528ce4d24SBarry Smith 1628ce4d24SBarry Smith Concepts: Krylov methods 1728ce4d24SBarry Smith 188f6c3df8SBarry Smith Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a 198f6c3df8SBarry Smith KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver). 208f6c3df8SBarry Smith 218f6c3df8SBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy() 2228ce4d24SBarry Smith S*/ 23e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 242eac72dbSBarry Smith 2576bdecfbSBarry Smith /*J 268f6c3df8SBarry Smith KSPType - String with the name of a PETSc Krylov method. 2728ce4d24SBarry Smith 2828ce4d24SBarry Smith Level: beginner 2928ce4d24SBarry Smith 308f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 3176bdecfbSBarry Smith J*/ 3219fd82e9SBarry Smith typedef const char* KSPType; 3382bf6240SBarry Smith #define KSPRICHARDSON "richardson" 346c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3582bf6240SBarry Smith #define KSPCG "cg" 362c8fc646SJed Brown #define KSPGROPPCG "groppcg" 372c8fc646SJed Brown #define KSPPIPECG "pipecg" 38df2a969aSvictorle #define KSPCGNE "cgne" 39fcae7a14Stmunson #define KSPNASH "nash" 4080e17431SMatthew Knepley #define KSPSTCG "stcg" 418031f4adStmunson #define KSPGLTR "gltr" 42640f4f53SPatrick Sanan #define KSPFCG "fcg" 43*390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 4482bf6240SBarry Smith #define KSPGMRES "gmres" 459dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 469dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 474f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 4861b468f9SJed Brown #define KSPPGMRES "pgmres" 4982bf6240SBarry Smith #define KSPTCQMR "tcqmr" 5082bf6240SBarry Smith #define KSPBCGS "bcgs" 51e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5218ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 53c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 54f0037002Svictorle #define KSPBCGSL "bcgsl" 5582bf6240SBarry Smith #define KSPCGS "cgs" 5682bf6240SBarry Smith #define KSPTFQMR "tfqmr" 5782bf6240SBarry Smith #define KSPCR "cr" 582c8fc646SJed Brown #define KSPPIPECR "pipecr" 5982bf6240SBarry Smith #define KSPLSQR "lsqr" 6082bf6240SBarry Smith #define KSPPREONLY "preonly" 6182bf6240SBarry Smith #define KSPQCG "qcg" 62c9cf9b11SBarry Smith #define KSPBICG "bicg" 63b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6401c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 6521356919SSatish Balay #define KSPLCD "lcd" 661d6018f0SLisandro Dalcin #define KSPPYTHON "python" 6758ad63f4SBarry Smith #define KSPGCR "gcr" 68e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 69e4d80e07Szianekhodja #define KSPCGLS "cgls" 702eac72dbSBarry Smith 718ba1e511SMatthew Knepley /* Logging support */ 72014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 735358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 748ba1e511SMatthew Knepley 75014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 7619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 77c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 8423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 8525c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 862eac72dbSBarry Smith 87140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 88bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 8930de9b25SBarry Smith 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 93c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 1017d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool ); 102c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 104c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 109734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1102a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1112a7a6963SBarry Smith PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1122eac72dbSBarry Smith 113d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 114d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 115d4211eb9SBarry Smith 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 118aabeff55SBarry Smith 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1254b0e389bSBarry Smith 126fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 127fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 128fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 129fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 130fa0ddf94SBarry Smith 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 136b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 137b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 138b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 139b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1410a71e23eSMatthew Knepley 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1442eac72dbSBarry Smith 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 14858450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 14982fe98a5SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool); 150ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 15158450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 153d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *); 154d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1557d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]); 1564b0e389bSBarry Smith 157640f4f53SPatrick Sanan /*E 158640f4f53SPatrick Sanan 15906137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 160640f4f53SPatrick Sanan 16106137d0aSPatrick Sanan KSP_FCD_TRUNCATION uses all (up to mmax) stored directions 16206137d0aSPatrick Sanan KSP_FCD_TRUNCATION_RESTART uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 163640f4f53SPatrick Sanan 1642b26979fSBarry Smith Level: intermediate 16506137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 166640f4f53SPatrick Sanan 167640f4f53SPatrick Sanan E*/ 16806137d0aSPatrick Sanan typedef enum {KSP_FCD_TRUNCATION,KSP_FCD_TRUNCATION_RESTART} KSPFCDTruncationType; 16906137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 17006137d0aSPatrick Sanan 17106137d0aSPatrick Sanan #define KSPFCDGetNumOldDirections(ctx,i,mi) ({ \ 17206137d0aSPatrick Sanan if(ctx->trunctype == KSP_FCD_TRUNCATION_RESTART){ \ 17306137d0aSPatrick Sanan mi = ((i-1) % ctx->mmax)+1; \ 17406137d0aSPatrick Sanan if (mi==1 && i!=1) ++(ctx->n_search_space_resets); \ 17506137d0aSPatrick Sanan } else if (ctx->trunctype == KSP_FCD_TRUNCATION) \ 17606137d0aSPatrick Sanan mi = ctx->mmax; \ 17706137d0aSPatrick Sanan else { \ 17806137d0aSPatrick Sanan SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy");CHKERRQ(ierr); \ 17906137d0aSPatrick Sanan } \ 18006137d0aSPatrick Sanan }) 181640f4f53SPatrick Sanan 182640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 183640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 184640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 185640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 18606137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 18706137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 188640f4f53SPatrick Sanan 189*390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 190*390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 191*390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 192*390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 193*390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 194*390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 195*390d8e47SPatrick Sanan 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1999f236934SBarry Smith 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2051d73ed98SBarry Smith 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2081d73ed98SBarry Smith 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 21258ad63f4SBarry Smith 213b49fd9e1SBarry Smith /*E 214b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 215b49fd9e1SBarry Smith 216b49fd9e1SBarry Smith Level: advanced 217b49fd9e1SBarry Smith 218bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 219e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 220b49fd9e1SBarry Smith 221b49fd9e1SBarry Smith E*/ 22278d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2236a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2241f7e983dSSatish Balay /*MC 2251957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2268c5b8ba0SBarry Smith 2278c5b8ba0SBarry Smith Level: advanced 2288c5b8ba0SBarry Smith 2298c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2308c5b8ba0SBarry Smith 231bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 232e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2338c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2348c5b8ba0SBarry Smith M*/ 2358c5b8ba0SBarry Smith 2361f7e983dSSatish Balay /*MC 2378c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2388c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2398c5b8ba0SBarry Smith poor orthogonality. 2408c5b8ba0SBarry Smith 2418c5b8ba0SBarry Smith Level: advanced 2428c5b8ba0SBarry Smith 2438c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2448c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2458c5b8ba0SBarry Smith 246bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 247e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2488c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2498c5b8ba0SBarry Smith M*/ 2508c5b8ba0SBarry Smith 2511f7e983dSSatish Balay /*MC 2528c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2538c5b8ba0SBarry Smith 2548c5b8ba0SBarry Smith Level: advanced 2558c5b8ba0SBarry Smith 2568c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2578c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2588c5b8ba0SBarry Smith 2598c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2608c5b8ba0SBarry Smith 261bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 262e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2638c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2648c5b8ba0SBarry Smith M*/ 2658c5b8ba0SBarry Smith 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 26808480c60SBarry Smith 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 272c38d4ed2SBarry Smith 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 276121fd945SKris Buschelman 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 280e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 281d9492815SBarry Smith 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2842eac72dbSBarry Smith 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 289390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 290390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 292816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 295e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 296e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 297e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 29984cb2905SBarry Smith 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 302c01c455dSBarry Smith 30323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 30423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3097287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3107287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3111eb62cbbSBarry Smith 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 313014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 314014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 315014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3161f7f0c4fSBarry Smith 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 31855849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 319685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 3202a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 3212a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 32255849f57SBarry Smith 32355849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3241eb62cbbSBarry Smith 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 326014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 327db9b2ab1SHong Zhang 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 329014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 33068ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 33183ab6a24SBarry Smith 33228ce4d24SBarry Smith /*E 3338a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3348a4b9c5cSBarry Smith test routines. 3358a4b9c5cSBarry Smith 3368a4b9c5cSBarry Smith Level: advanced 3378a4b9c5cSBarry Smith 338a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3391718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 340a3f661c8SBarry Smith 341af0996ceSBarry Smith Notes: this must match petsc/finclude/petscksp.h 3428a4b9c5cSBarry Smith 34394b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3441718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3458a4b9c5cSBarry Smith E*/ 3469e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3479e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 348014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 349a21b2a99SBarry Smith 3501f7e983dSSatish Balay /*MC 3519793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3528c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3538c5b8ba0SBarry Smith be based on a norm of a residual etc. 3548c5b8ba0SBarry Smith 3558c5b8ba0SBarry Smith Level: advanced 3568c5b8ba0SBarry Smith 3571957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3588c5b8ba0SBarry Smith 359ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3608c5b8ba0SBarry Smith M*/ 3618c5b8ba0SBarry Smith 3621f7e983dSSatish Balay /*MC 3631957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3648c5b8ba0SBarry Smith convergence test routine. 3658c5b8ba0SBarry Smith 3668c5b8ba0SBarry Smith Level: advanced 3678c5b8ba0SBarry Smith 3689793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3698c5b8ba0SBarry Smith M*/ 3708c5b8ba0SBarry Smith 3711f7e983dSSatish Balay /*MC 372ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3738c5b8ba0SBarry Smith convergence test routine. 3748c5b8ba0SBarry Smith 3758c5b8ba0SBarry Smith Level: advanced 3768c5b8ba0SBarry Smith 3779793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3788c5b8ba0SBarry Smith M*/ 3798c5b8ba0SBarry Smith 3801f7e983dSSatish Balay /*MC 381ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 382*390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 3838c5b8ba0SBarry Smith 3848c5b8ba0SBarry Smith Level: advanced 3858c5b8ba0SBarry Smith 3869793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3878c5b8ba0SBarry Smith M*/ 3888c5b8ba0SBarry Smith 389014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 391014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 392014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 393014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3948a4b9c5cSBarry Smith 3958a4b9c5cSBarry Smith /*E 3961957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 39728ce4d24SBarry Smith 39828ce4d24SBarry Smith Level: beginner 39928ce4d24SBarry Smith 4004d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 40128ce4d24SBarry Smith 402af0996ceSBarry Smith Developer notes: this must match petsc/finclude/petscksp.h 4034d0a8057SBarry Smith 4044d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4054d0a8057SBarry Smith any of the values here also change them that array of names. 40686c02ca4SBarry Smith 407c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 40828ce4d24SBarry Smith E*/ 409d15094e1SBarry Smith typedef enum {/* converged */ 4109ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4119ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 412d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 413d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 414b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4158031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4168031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 417329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 418af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 419d15094e1SBarry Smith /* diverged */ 420b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 421d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 422d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 423d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 424b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 425b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 426b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4274d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4286aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 429422a814eSBarry Smith KSP_DIVERGED_PCSETUP_FAILED = -11, 430d15094e1SBarry Smith 431d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 432014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 433d15094e1SBarry Smith 434c838673bSBarry Smith /*MC 435c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 436c838673bSBarry Smith 437c838673bSBarry Smith Level: beginner 438c838673bSBarry Smith 439c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 440c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 441c838673bSBarry Smith 2-norm of the residual for right preconditioning 442c838673bSBarry Smith 443c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 444c838673bSBarry Smith 445c838673bSBarry Smith M*/ 446c838673bSBarry Smith 447c838673bSBarry Smith /*MC 448c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 449c838673bSBarry Smith 450c838673bSBarry Smith Level: beginner 451c838673bSBarry Smith 452c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 453c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 454c838673bSBarry Smith 2-norm of the residual for right preconditioning 455c838673bSBarry Smith 456c838673bSBarry Smith Level: beginner 457c838673bSBarry Smith 458c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 459c838673bSBarry Smith 460c838673bSBarry Smith M*/ 461c838673bSBarry Smith 462c838673bSBarry Smith /*MC 463c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 464c838673bSBarry Smith 465c838673bSBarry Smith Level: beginner 466c838673bSBarry Smith 467c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 468c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 469c838673bSBarry Smith 2-norm of the residual for right preconditioning 470c838673bSBarry Smith 471c838673bSBarry Smith Level: beginner 472c838673bSBarry Smith 473c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 474c838673bSBarry Smith 475c838673bSBarry Smith M*/ 476c838673bSBarry Smith 477c838673bSBarry Smith /*MC 478c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 479c838673bSBarry Smith reached 480c838673bSBarry Smith 481c838673bSBarry Smith Level: beginner 482c838673bSBarry Smith 483c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 484c838673bSBarry Smith 485c838673bSBarry Smith M*/ 486c838673bSBarry Smith 487c838673bSBarry Smith /*MC 4888c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4890059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 4904d0a8057SBarry Smith test routine is set in KSP. 491c838673bSBarry Smith 492c838673bSBarry Smith Level: beginner 493c838673bSBarry Smith 494c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 495c838673bSBarry Smith 496c838673bSBarry Smith M*/ 497c838673bSBarry Smith 498c838673bSBarry Smith /*MC 499c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 5003014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 5013014e516SBarry Smith preconditioner. 502c838673bSBarry Smith 503c838673bSBarry Smith Level: beginner 504c838673bSBarry Smith 505c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 506c838673bSBarry Smith 507c838673bSBarry Smith M*/ 508c838673bSBarry Smith 509c838673bSBarry Smith /*MC 510c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 511c838673bSBarry Smith method could not continue to enlarge the Krylov space. 512c838673bSBarry Smith 513c838673bSBarry Smith Level: beginner 514c838673bSBarry Smith 515c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 516c838673bSBarry Smith 517c838673bSBarry Smith M*/ 518c838673bSBarry Smith 519c838673bSBarry Smith /*MC 520c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 521c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 522c838673bSBarry Smith 523c838673bSBarry Smith Level: beginner 524c838673bSBarry Smith 525c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 526c838673bSBarry Smith 527c838673bSBarry Smith M*/ 528c838673bSBarry Smith 529c838673bSBarry Smith /*MC 530c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 531c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 532c838673bSBarry Smith be positive definite 533c838673bSBarry Smith 534c838673bSBarry Smith Level: beginner 535c838673bSBarry Smith 5362401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 537c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 538c838673bSBarry Smith 539c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 540c838673bSBarry Smith 541c838673bSBarry Smith M*/ 542c838673bSBarry Smith 543c838673bSBarry Smith /*MC 544c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 545c838673bSBarry Smith while the KSPSolve() is still running. 546c838673bSBarry Smith 547c838673bSBarry Smith Level: beginner 548c838673bSBarry Smith 549c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 550c838673bSBarry Smith 551c838673bSBarry Smith M*/ 552c838673bSBarry Smith 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 5558de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 5578de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 5588de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 5598de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5608de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5610059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 563abef13c0SSatish Balay 5648ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5658ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5668ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5678ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5688ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5698ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5708ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5718ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5728ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 5738ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 5748ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 5758ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 5768ea1b3e6SJed Brown 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 578d4fbbf0eSBarry Smith 57928ce4d24SBarry Smith /*E 58028ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 58128ce4d24SBarry Smith 58228ce4d24SBarry Smith Level: beginner 58328ce4d24SBarry Smith 58428ce4d24SBarry Smith .seealso: KSPCGSetType() 58528ce4d24SBarry Smith E*/ 5869dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5876a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 58828ce4d24SBarry Smith 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5918031f4adStmunson 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 595fcae7a14Stmunson 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 597014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 599e559a7feSLois Curfman McInnes 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 6058031f4adStmunson 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 6071d6018f0SLisandro Dalcin 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 6103369ce9aSBarry Smith 6119804daf3SBarry Smith #include <petscdrawtypes.h> 612d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 613d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 614d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 615d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6172f2e5d10SKris Buschelman 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 62003919abeSBarry Smith 621f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 622ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 623f8a50e2bSBarry Smith 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 626014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 627014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 629014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 630f8a50e2bSBarry Smith 631014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 632014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 633014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 634f8a50e2bSBarry Smith 635470b340bSDmitry Karpeev /*E 636470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 637470b340bSDmitry Karpeev 638470b340bSDmitry Karpeev Level: intermediate 639470b340bSDmitry Karpeev 640470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 641470b340bSDmitry Karpeev E*/ 642470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 643470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 644470b340bSDmitry Karpeev 645014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 646014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 647d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 648bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 649aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 650bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 651470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 652470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 6535bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 6545a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 655470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 656470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 6573f22127dSBarry Smith 658014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 659014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 660014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 661014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 662014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 663fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 66423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 665fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 66623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 66723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 668014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 669014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 670fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 671fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6726c699258SBarry Smith 67302b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 67402b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexProjectField(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec); 6752eac72dbSBarry Smith #endif 676