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" 42*640f4f53SPatrick 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" 67bbce1389SJed Brown #define KSPSPECEST "specest" 682eac72dbSBarry Smith 698ba1e511SMatthew Knepley /* Logging support */ 70014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 715358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 728ba1e511SMatthew Knepley 73014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 7419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 75014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 8123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 822eac72dbSBarry Smith 83140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 84014dd563SJed Brown PETSC_EXTERN PetscBool KSPRegisterAllCalled; 85607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegisterAll(void); 86bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 87607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void); 8830de9b25SBarry Smith 8919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(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 KSPGetComputeEigenvalues(KSP,PetscBool *); 101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(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*); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 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); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 148ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 151d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *); 152d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1534b0e389bSBarry Smith 154*640f4f53SPatrick Sanan /*E 155*640f4f53SPatrick Sanan 156*640f4f53SPatrick Sanan KSPFCGTruncationStrategy - Define how many stored directions are used to orthogonalize in FCG 157*640f4f53SPatrick Sanan 158*640f4f53SPatrick Sanan KSP_FCG_TRUNCATION_STRATEGY_STANDARD uses all (up to mmax) stored directions 159*640f4f53SPatrick Sanan KSP_FCG_TRUNCATION_STRATEGY_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 160*640f4f53SPatrick Sanan 161*640f4f53SPatrick Sanan .seealso : KSPFCG,KSPFCGSetTruncationStrategy(),KSPFCGGetTruncationStrategy() 162*640f4f53SPatrick Sanan 163*640f4f53SPatrick Sanan E*/ 164*640f4f53SPatrick Sanan typedef enum {KSP_FCG_TRUNCATION_STRATEGY_STANDARD,KSP_FCG_TRUNCATION_STRATEGY_NOTAY} KSPFCGTruncationStrategy; 165*640f4f53SPatrick Sanan 166*640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 167*640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 168*640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 169*640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 170*640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationStrategy(KSP,KSPFCGTruncationStrategy); 171*640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationStrategy(KSP,KSPFCGTruncationStrategy*); 172*640f4f53SPatrick Sanan 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1769f236934SBarry Smith 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1821d73ed98SBarry Smith 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 1851d73ed98SBarry Smith 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 18958ad63f4SBarry Smith 190b49fd9e1SBarry Smith /*E 191b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 192b49fd9e1SBarry Smith 193b49fd9e1SBarry Smith Level: advanced 194b49fd9e1SBarry Smith 195bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 196e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 197b49fd9e1SBarry Smith 198b49fd9e1SBarry Smith E*/ 19978d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2006a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2011f7e983dSSatish Balay /*MC 2021957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2038c5b8ba0SBarry Smith 2048c5b8ba0SBarry Smith Level: advanced 2058c5b8ba0SBarry Smith 2068c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2078c5b8ba0SBarry Smith 208bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 209e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2108c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2118c5b8ba0SBarry Smith M*/ 2128c5b8ba0SBarry Smith 2131f7e983dSSatish Balay /*MC 2148c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2158c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2168c5b8ba0SBarry Smith poor orthogonality. 2178c5b8ba0SBarry Smith 2188c5b8ba0SBarry Smith Level: advanced 2198c5b8ba0SBarry Smith 2208c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2218c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2228c5b8ba0SBarry Smith 223bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 224e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2258c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2268c5b8ba0SBarry Smith M*/ 2278c5b8ba0SBarry Smith 2281f7e983dSSatish Balay /*MC 2298c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2308c5b8ba0SBarry Smith 2318c5b8ba0SBarry Smith Level: advanced 2328c5b8ba0SBarry Smith 2338c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2348c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2358c5b8ba0SBarry Smith 2368c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2378c5b8ba0SBarry Smith 238bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 239e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2408c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2418c5b8ba0SBarry Smith M*/ 2428c5b8ba0SBarry Smith 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 244014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 24508480c60SBarry Smith 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 249c38d4ed2SBarry Smith 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 253121fd945SKris Buschelman 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 257e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 258d9492815SBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2612eac72dbSBarry Smith 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 266390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 267390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 269816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 272e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 273e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 274e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 27684cb2905SBarry Smith 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 279c01c455dSBarry Smith 28023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 28123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 2867287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 2877287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 2881eb62cbbSBarry Smith 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 2931f7f0c4fSBarry Smith 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 29555849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 296ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 29755849f57SBarry Smith 29855849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 2991eb62cbbSBarry Smith 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 302db9b2ab1SHong Zhang 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 304014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 30583ab6a24SBarry Smith 30628ce4d24SBarry Smith /*E 3078a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3088a4b9c5cSBarry Smith test routines. 3098a4b9c5cSBarry Smith 3108a4b9c5cSBarry Smith Level: advanced 3118a4b9c5cSBarry Smith 312a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3131718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 314a3f661c8SBarry Smith 3158a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3168a4b9c5cSBarry Smith 31794b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3181718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3198a4b9c5cSBarry Smith E*/ 3209e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3219e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 322014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 323a21b2a99SBarry Smith 3241f7e983dSSatish Balay /*MC 3259793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3268c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3278c5b8ba0SBarry Smith be based on a norm of a residual etc. 3288c5b8ba0SBarry Smith 3298c5b8ba0SBarry Smith Level: advanced 3308c5b8ba0SBarry Smith 3311957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3328c5b8ba0SBarry Smith 333ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3348c5b8ba0SBarry Smith M*/ 3358c5b8ba0SBarry Smith 3361f7e983dSSatish Balay /*MC 3371957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3388c5b8ba0SBarry Smith convergence test routine. 3398c5b8ba0SBarry Smith 3408c5b8ba0SBarry Smith Level: advanced 3418c5b8ba0SBarry Smith 3429793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3438c5b8ba0SBarry Smith M*/ 3448c5b8ba0SBarry Smith 3451f7e983dSSatish Balay /*MC 346ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3478c5b8ba0SBarry Smith convergence test routine. 3488c5b8ba0SBarry Smith 3498c5b8ba0SBarry Smith Level: advanced 3508c5b8ba0SBarry Smith 3519793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3528c5b8ba0SBarry Smith M*/ 3538c5b8ba0SBarry Smith 3541f7e983dSSatish Balay /*MC 355ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 356*640f4f53SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG 3578c5b8ba0SBarry Smith 3588c5b8ba0SBarry Smith Level: advanced 3598c5b8ba0SBarry Smith 3609793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3618c5b8ba0SBarry Smith M*/ 3628c5b8ba0SBarry Smith 363014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 364014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 365014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 366014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 367014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3688a4b9c5cSBarry Smith 3698a4b9c5cSBarry Smith /*E 3701957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 37128ce4d24SBarry Smith 37228ce4d24SBarry Smith Level: beginner 37328ce4d24SBarry Smith 3744d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 37528ce4d24SBarry Smith 3764d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3774d0a8057SBarry Smith 3784d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3794d0a8057SBarry Smith any of the values here also change them that array of names. 38086c02ca4SBarry Smith 381c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 38228ce4d24SBarry Smith E*/ 383d15094e1SBarry Smith typedef enum {/* converged */ 3849ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 3859ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 386d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 387d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 388b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3898031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3908031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 391329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 392af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 393d15094e1SBarry Smith /* diverged */ 394b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 395d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 396d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 397d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 398b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 399b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 400b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4014d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4026aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 403d15094e1SBarry Smith 404d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 405014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 406d15094e1SBarry Smith 407c838673bSBarry Smith /*MC 408c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 409c838673bSBarry Smith 410c838673bSBarry Smith Level: beginner 411c838673bSBarry Smith 412c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 413c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 414c838673bSBarry Smith 2-norm of the residual for right preconditioning 415c838673bSBarry Smith 416c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 417c838673bSBarry Smith 418c838673bSBarry Smith M*/ 419c838673bSBarry Smith 420c838673bSBarry Smith /*MC 421c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 422c838673bSBarry Smith 423c838673bSBarry Smith Level: beginner 424c838673bSBarry Smith 425c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 426c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 427c838673bSBarry Smith 2-norm of the residual for right preconditioning 428c838673bSBarry Smith 429c838673bSBarry Smith Level: beginner 430c838673bSBarry Smith 431c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 432c838673bSBarry Smith 433c838673bSBarry Smith M*/ 434c838673bSBarry Smith 435c838673bSBarry Smith /*MC 436c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 437c838673bSBarry Smith 438c838673bSBarry Smith Level: beginner 439c838673bSBarry Smith 440c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 441c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 442c838673bSBarry Smith 2-norm of the residual for right preconditioning 443c838673bSBarry Smith 444c838673bSBarry Smith Level: beginner 445c838673bSBarry Smith 446c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 447c838673bSBarry Smith 448c838673bSBarry Smith M*/ 449c838673bSBarry Smith 450c838673bSBarry Smith /*MC 451c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 452c838673bSBarry Smith reached 453c838673bSBarry Smith 454c838673bSBarry Smith Level: beginner 455c838673bSBarry Smith 456c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 457c838673bSBarry Smith 458c838673bSBarry Smith M*/ 459c838673bSBarry Smith 460c838673bSBarry Smith /*MC 4618c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4620059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 4634d0a8057SBarry Smith test routine is set in KSP. 464c838673bSBarry Smith 465c838673bSBarry Smith Level: beginner 466c838673bSBarry Smith 467c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 468c838673bSBarry Smith 469c838673bSBarry Smith M*/ 470c838673bSBarry Smith 471c838673bSBarry Smith /*MC 472c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4733014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4743014e516SBarry Smith preconditioner. 475c838673bSBarry Smith 476c838673bSBarry Smith Level: beginner 477c838673bSBarry Smith 478c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 479c838673bSBarry Smith 480c838673bSBarry Smith M*/ 481c838673bSBarry Smith 482c838673bSBarry Smith /*MC 483c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 484c838673bSBarry Smith method could not continue to enlarge the Krylov space. 485c838673bSBarry Smith 486c838673bSBarry Smith Level: beginner 487c838673bSBarry Smith 488c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 489c838673bSBarry Smith 490c838673bSBarry Smith M*/ 491c838673bSBarry Smith 492c838673bSBarry Smith /*MC 493c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 494c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 495c838673bSBarry Smith 496c838673bSBarry Smith Level: beginner 497c838673bSBarry Smith 498c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 499c838673bSBarry Smith 500c838673bSBarry Smith M*/ 501c838673bSBarry Smith 502c838673bSBarry Smith /*MC 503c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 504c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 505c838673bSBarry Smith be positive definite 506c838673bSBarry Smith 507c838673bSBarry Smith Level: beginner 508c838673bSBarry Smith 5092401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 510c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 511c838673bSBarry Smith 512c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 513c838673bSBarry Smith 514c838673bSBarry Smith M*/ 515c838673bSBarry Smith 516c838673bSBarry Smith /*MC 517c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 518c838673bSBarry Smith while the KSPSolve() is still running. 519c838673bSBarry Smith 520c838673bSBarry Smith Level: beginner 521c838673bSBarry Smith 522c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 523c838673bSBarry Smith 524c838673bSBarry Smith M*/ 525c838673bSBarry Smith 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 5288de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 5308de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 5318de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 5328de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5338de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5340059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 536abef13c0SSatish Balay 5378ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5388ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5398ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5408ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5418ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5428ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5438ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5448ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5458ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 5468ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 5478ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 5488ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 5498ea1b3e6SJed Brown 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 551d4fbbf0eSBarry Smith 55228ce4d24SBarry Smith /*E 55328ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 55428ce4d24SBarry Smith 55528ce4d24SBarry Smith Level: beginner 55628ce4d24SBarry Smith 55728ce4d24SBarry Smith .seealso: KSPCGSetType() 55828ce4d24SBarry Smith E*/ 5599dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5606a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 56128ce4d24SBarry Smith 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5648031f4adStmunson 565014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 568fcae7a14Stmunson 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 570014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 572e559a7feSLois Curfman McInnes 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5788031f4adStmunson 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5801d6018f0SLisandro Dalcin 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5833369ce9aSBarry Smith 5849804daf3SBarry Smith #include <petscdrawtypes.h> 5851d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 5861d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 5871d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 5922f2e5d10SKris Buschelman 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 59503919abeSBarry Smith 596f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 597ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 598f8a50e2bSBarry Smith 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 605f8a50e2bSBarry Smith 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 609f8a50e2bSBarry Smith 610470b340bSDmitry Karpeev /*E 611470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 612470b340bSDmitry Karpeev 613470b340bSDmitry Karpeev Level: intermediate 614470b340bSDmitry Karpeev 615470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 616470b340bSDmitry Karpeev E*/ 617470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 618470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 619470b340bSDmitry Karpeev 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 622d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 623bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 624aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 625bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 626470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 627470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 6285bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 6295a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 630470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 631470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 6323f22127dSBarry Smith 633014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 634014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 635014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 638fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 63923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 640fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 64123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 64223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 643014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 644014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 645fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 646fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6476c699258SBarry Smith 6482eac72dbSBarry Smith #endif 649