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" 4382bf6240SBarry Smith #define KSPGMRES "gmres" 449dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 459dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 464f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 4761b468f9SJed Brown #define KSPPGMRES "pgmres" 4882bf6240SBarry Smith #define KSPTCQMR "tcqmr" 4982bf6240SBarry Smith #define KSPBCGS "bcgs" 50e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5118ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 52c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 53f0037002Svictorle #define KSPBCGSL "bcgsl" 5482bf6240SBarry Smith #define KSPCGS "cgs" 5582bf6240SBarry Smith #define KSPTFQMR "tfqmr" 5682bf6240SBarry Smith #define KSPCR "cr" 572c8fc646SJed Brown #define KSPPIPECR "pipecr" 5882bf6240SBarry Smith #define KSPLSQR "lsqr" 5982bf6240SBarry Smith #define KSPPREONLY "preonly" 6082bf6240SBarry Smith #define KSPQCG "qcg" 61c9cf9b11SBarry Smith #define KSPBICG "bicg" 62b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6301c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 6421356919SSatish Balay #define KSPLCD "lcd" 651d6018f0SLisandro Dalcin #define KSPPYTHON "python" 6658ad63f4SBarry Smith #define KSPGCR "gcr" 67e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 68e4d80e07Szianekhodja #define KSPCGLS "cgls" 692eac72dbSBarry Smith 708ba1e511SMatthew Knepley /* Logging support */ 71014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 725358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 738ba1e511SMatthew Knepley 74014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 7519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 76c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 8323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 8425c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 852eac72dbSBarry Smith 86140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 87bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 8830de9b25SBarry Smith 89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 92c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 1007d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool ); 101c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 103c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 108734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1092a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1102a7a6963SBarry 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);} 1112eac72dbSBarry Smith 112d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 113d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 114d4211eb9SBarry Smith 115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 117aabeff55SBarry Smith 118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1244b0e389bSBarry Smith 125fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 126fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 127fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 128fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 129fa0ddf94SBarry Smith 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 135b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 136b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 137b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 138b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1400a71e23eSMatthew Knepley 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1432eac72dbSBarry Smith 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 14758450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 14882fe98a5SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool); 149ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 15058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 152d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *); 153d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1547d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]); 1554b0e389bSBarry Smith 156640f4f53SPatrick Sanan /*E 157640f4f53SPatrick Sanan 158bd242c34SBarry Smith KSPFCGTruncationType - Define how stored directions are used to orthogonalize in FCG 159640f4f53SPatrick Sanan 1609662b6c4SBarry Smith KSP_FCG_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 1619662b6c4SBarry Smith KSP_FCG_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 162640f4f53SPatrick Sanan 1632b26979fSBarry Smith Level: intermediate 1642b26979fSBarry Smith 165bd242c34SBarry Smith .seealso : KSPFCG,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 166640f4f53SPatrick Sanan 167640f4f53SPatrick Sanan E*/ 1689662b6c4SBarry Smith typedef enum {KSP_FCG_TRUNC_TYPE_STANDARD,KSP_FCG_TRUNC_TYPE_NOTAY} KSPFCGTruncationType; 169bd242c34SBarry Smith PETSC_EXTERN const char *const KSPFCGTruncationTypes[]; 170640f4f53SPatrick Sanan 171640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 172640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 173640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 174640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 175bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCGTruncationType); 176bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCGTruncationType*); 177640f4f53SPatrick Sanan 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1819f236934SBarry Smith 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1871d73ed98SBarry Smith 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 1901d73ed98SBarry Smith 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 19458ad63f4SBarry Smith 195b49fd9e1SBarry Smith /*E 196b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 197b49fd9e1SBarry Smith 198b49fd9e1SBarry Smith Level: advanced 199b49fd9e1SBarry Smith 200bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 201e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 202b49fd9e1SBarry Smith 203b49fd9e1SBarry Smith E*/ 20478d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2056a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2061f7e983dSSatish Balay /*MC 2071957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2088c5b8ba0SBarry Smith 2098c5b8ba0SBarry Smith Level: advanced 2108c5b8ba0SBarry Smith 2118c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2128c5b8ba0SBarry Smith 213bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 214e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2158c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2168c5b8ba0SBarry Smith M*/ 2178c5b8ba0SBarry Smith 2181f7e983dSSatish Balay /*MC 2198c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2208c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2218c5b8ba0SBarry Smith poor orthogonality. 2228c5b8ba0SBarry Smith 2238c5b8ba0SBarry Smith Level: advanced 2248c5b8ba0SBarry Smith 2258c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2268c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2278c5b8ba0SBarry Smith 228bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 229e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2308c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2318c5b8ba0SBarry Smith M*/ 2328c5b8ba0SBarry Smith 2331f7e983dSSatish Balay /*MC 2348c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2358c5b8ba0SBarry Smith 2368c5b8ba0SBarry Smith Level: advanced 2378c5b8ba0SBarry Smith 2388c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2398c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2408c5b8ba0SBarry Smith 2418c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2428c5b8ba0SBarry Smith 243bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 244e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2458c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2468c5b8ba0SBarry Smith M*/ 2478c5b8ba0SBarry Smith 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 25008480c60SBarry Smith 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 254c38d4ed2SBarry Smith 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 258121fd945SKris Buschelman 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 262e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 263d9492815SBarry Smith 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2662eac72dbSBarry Smith 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 271390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 272390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 274816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 277e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 278e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 279e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 28184cb2905SBarry Smith 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 284c01c455dSBarry Smith 28523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 28623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 2917287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 2927287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 2931eb62cbbSBarry Smith 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 2981f7f0c4fSBarry Smith 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 30055849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 301685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 3022a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 3032a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 30455849f57SBarry Smith 30555849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3061eb62cbbSBarry Smith 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 309db9b2ab1SHong Zhang 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 31268ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 31383ab6a24SBarry Smith 31428ce4d24SBarry Smith /*E 3158a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3168a4b9c5cSBarry Smith test routines. 3178a4b9c5cSBarry Smith 3188a4b9c5cSBarry Smith Level: advanced 3198a4b9c5cSBarry Smith 320a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3211718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 322a3f661c8SBarry Smith 323af0996ceSBarry Smith Notes: this must match petsc/finclude/petscksp.h 3248a4b9c5cSBarry Smith 32594b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3261718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3278a4b9c5cSBarry Smith E*/ 3289e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3299e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 330014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 331a21b2a99SBarry Smith 3321f7e983dSSatish Balay /*MC 3339793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3348c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3358c5b8ba0SBarry Smith be based on a norm of a residual etc. 3368c5b8ba0SBarry Smith 3378c5b8ba0SBarry Smith Level: advanced 3388c5b8ba0SBarry Smith 3391957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3408c5b8ba0SBarry Smith 341ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3428c5b8ba0SBarry Smith M*/ 3438c5b8ba0SBarry Smith 3441f7e983dSSatish Balay /*MC 3451957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3468c5b8ba0SBarry Smith convergence test routine. 3478c5b8ba0SBarry Smith 3488c5b8ba0SBarry Smith Level: advanced 3498c5b8ba0SBarry Smith 3509793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3518c5b8ba0SBarry Smith M*/ 3528c5b8ba0SBarry Smith 3531f7e983dSSatish Balay /*MC 354ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3558c5b8ba0SBarry Smith convergence test routine. 3568c5b8ba0SBarry Smith 3578c5b8ba0SBarry Smith Level: advanced 3588c5b8ba0SBarry Smith 3599793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3608c5b8ba0SBarry Smith M*/ 3618c5b8ba0SBarry Smith 3621f7e983dSSatish Balay /*MC 363ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 364640f4f53SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG 3658c5b8ba0SBarry Smith 3668c5b8ba0SBarry Smith Level: advanced 3678c5b8ba0SBarry Smith 3689793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3698c5b8ba0SBarry Smith M*/ 3708c5b8ba0SBarry Smith 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3768a4b9c5cSBarry Smith 3778a4b9c5cSBarry Smith /*E 3781957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 37928ce4d24SBarry Smith 38028ce4d24SBarry Smith Level: beginner 38128ce4d24SBarry Smith 3824d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 38328ce4d24SBarry Smith 384af0996ceSBarry Smith Developer notes: this must match petsc/finclude/petscksp.h 3854d0a8057SBarry Smith 3864d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3874d0a8057SBarry Smith any of the values here also change them that array of names. 38886c02ca4SBarry Smith 389c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 39028ce4d24SBarry Smith E*/ 391d15094e1SBarry Smith typedef enum {/* converged */ 3929ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 3939ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 394d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 395d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 396b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3978031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3988031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 399329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 400af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 401d15094e1SBarry Smith /* diverged */ 402b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 403d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 404d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 405d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 406b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 407b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 408b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4094d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4106aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 411422a814eSBarry Smith KSP_DIVERGED_PCSETUP_FAILED = -11, 412d15094e1SBarry Smith 413d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 414014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 415d15094e1SBarry Smith 416c838673bSBarry Smith /*MC 417c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 418c838673bSBarry Smith 419c838673bSBarry Smith Level: beginner 420c838673bSBarry Smith 421c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 422c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 423c838673bSBarry Smith 2-norm of the residual for right preconditioning 424c838673bSBarry Smith 425c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 426c838673bSBarry Smith 427c838673bSBarry Smith M*/ 428c838673bSBarry Smith 429c838673bSBarry Smith /*MC 430c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 431c838673bSBarry Smith 432c838673bSBarry Smith Level: beginner 433c838673bSBarry Smith 434c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 435c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 436c838673bSBarry Smith 2-norm of the residual for right preconditioning 437c838673bSBarry Smith 438c838673bSBarry Smith Level: beginner 439c838673bSBarry Smith 440c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 441c838673bSBarry Smith 442c838673bSBarry Smith M*/ 443c838673bSBarry Smith 444c838673bSBarry Smith /*MC 445c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 446c838673bSBarry Smith 447c838673bSBarry Smith Level: beginner 448c838673bSBarry Smith 449c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 450c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 451c838673bSBarry Smith 2-norm of the residual for right preconditioning 452c838673bSBarry Smith 453c838673bSBarry Smith Level: beginner 454c838673bSBarry Smith 455c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 456c838673bSBarry Smith 457c838673bSBarry Smith M*/ 458c838673bSBarry Smith 459c838673bSBarry Smith /*MC 460c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 461c838673bSBarry Smith reached 462c838673bSBarry Smith 463c838673bSBarry Smith Level: beginner 464c838673bSBarry Smith 465c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 466c838673bSBarry Smith 467c838673bSBarry Smith M*/ 468c838673bSBarry Smith 469c838673bSBarry Smith /*MC 4708c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4710059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 4724d0a8057SBarry Smith test routine is set in KSP. 473c838673bSBarry Smith 474c838673bSBarry Smith Level: beginner 475c838673bSBarry Smith 476c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 477c838673bSBarry Smith 478c838673bSBarry Smith M*/ 479c838673bSBarry Smith 480c838673bSBarry Smith /*MC 481c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4823014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4833014e516SBarry Smith preconditioner. 484c838673bSBarry Smith 485c838673bSBarry Smith Level: beginner 486c838673bSBarry Smith 487c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 488c838673bSBarry Smith 489c838673bSBarry Smith M*/ 490c838673bSBarry Smith 491c838673bSBarry Smith /*MC 492c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 493c838673bSBarry Smith method could not continue to enlarge the Krylov space. 494c838673bSBarry Smith 495c838673bSBarry Smith Level: beginner 496c838673bSBarry Smith 497c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 498c838673bSBarry Smith 499c838673bSBarry Smith M*/ 500c838673bSBarry Smith 501c838673bSBarry Smith /*MC 502c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 503c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 504c838673bSBarry Smith 505c838673bSBarry Smith Level: beginner 506c838673bSBarry Smith 507c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 508c838673bSBarry Smith 509c838673bSBarry Smith M*/ 510c838673bSBarry Smith 511c838673bSBarry Smith /*MC 512c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 513c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 514c838673bSBarry Smith be positive definite 515c838673bSBarry Smith 516c838673bSBarry Smith Level: beginner 517c838673bSBarry Smith 5182401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 519c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 520c838673bSBarry Smith 521c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 522c838673bSBarry Smith 523c838673bSBarry Smith M*/ 524c838673bSBarry Smith 525c838673bSBarry Smith /*MC 526c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 527c838673bSBarry Smith while the KSPSolve() is still running. 528c838673bSBarry Smith 529c838673bSBarry Smith Level: beginner 530c838673bSBarry Smith 531c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 532c838673bSBarry Smith 533c838673bSBarry Smith M*/ 534c838673bSBarry Smith 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 5378de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 5398de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 5408de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 5418de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5428de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5430059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 545abef13c0SSatish Balay 5468ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5478ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5488ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5498ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5508ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5518ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5528ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5538ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5548ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 5558ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 5568ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 5578ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 5588ea1b3e6SJed Brown 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 560d4fbbf0eSBarry Smith 56128ce4d24SBarry Smith /*E 56228ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 56328ce4d24SBarry Smith 56428ce4d24SBarry Smith Level: beginner 56528ce4d24SBarry Smith 56628ce4d24SBarry Smith .seealso: KSPCGSetType() 56728ce4d24SBarry Smith E*/ 5689dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5696a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 57028ce4d24SBarry Smith 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5738031f4adStmunson 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 577fcae7a14Stmunson 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 581e559a7feSLois Curfman McInnes 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5878031f4adStmunson 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5891d6018f0SLisandro Dalcin 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5923369ce9aSBarry Smith 5939804daf3SBarry Smith #include <petscdrawtypes.h> 594d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 595d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 596d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 597d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 5992f2e5d10SKris Buschelman 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 60203919abeSBarry Smith 603f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 604ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 605f8a50e2bSBarry Smith 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 612f8a50e2bSBarry Smith 613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 616f8a50e2bSBarry Smith 617470b340bSDmitry Karpeev /*E 618470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 619470b340bSDmitry Karpeev 620470b340bSDmitry Karpeev Level: intermediate 621470b340bSDmitry Karpeev 622470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 623470b340bSDmitry Karpeev E*/ 624470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 625470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 626470b340bSDmitry Karpeev 627014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 629d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 630bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 631aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 632bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 633470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 634470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 6355bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 6365a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 637470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 638470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 6393f22127dSBarry Smith 640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 643014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 644014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 645fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 64623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 647fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 64823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 64923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 650014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 651014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 652fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 653fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6546c699258SBarry Smith 65502b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 656*bfc4295aSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectField(DM, Vec, 657e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 658e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 659e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 660e29c0834SMatthew G. Knepley PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec); 6612eac72dbSBarry Smith #endif 662