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" 4282bf6240SBarry Smith #define KSPGMRES "gmres" 439dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 449dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 454f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 4661b468f9SJed Brown #define KSPPGMRES "pgmres" 4782bf6240SBarry Smith #define KSPTCQMR "tcqmr" 4882bf6240SBarry Smith #define KSPBCGS "bcgs" 49e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5018ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 51c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 52f0037002Svictorle #define KSPBCGSL "bcgsl" 5382bf6240SBarry Smith #define KSPCGS "cgs" 5482bf6240SBarry Smith #define KSPTFQMR "tfqmr" 5582bf6240SBarry Smith #define KSPCR "cr" 562c8fc646SJed Brown #define KSPPIPECR "pipecr" 5782bf6240SBarry Smith #define KSPLSQR "lsqr" 5882bf6240SBarry Smith #define KSPPREONLY "preonly" 5982bf6240SBarry Smith #define KSPQCG "qcg" 60c9cf9b11SBarry Smith #define KSPBICG "bicg" 61b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6201c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 6321356919SSatish Balay #define KSPLCD "lcd" 641d6018f0SLisandro Dalcin #define KSPPYTHON "python" 6558ad63f4SBarry Smith #define KSPGCR "gcr" 66bbce1389SJed Brown #define KSPSPECEST "specest" 672eac72dbSBarry Smith 688ba1e511SMatthew Knepley /* Logging support */ 69014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 705358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 718ba1e511SMatthew Knepley 72014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 7319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 74014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 75014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 802eac72dbSBarry Smith 81140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 82014dd563SJed Brown PETSC_EXTERN PetscBool KSPRegisterAllCalled; 83607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegisterAll(void); 84bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 85607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void); 8630de9b25SBarry Smith 8719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 88014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1092eac72dbSBarry Smith 110d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 111d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 112d4211eb9SBarry Smith 113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 115aabeff55SBarry Smith 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1224b0e389bSBarry Smith 123fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 124fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 125fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 126fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 127fa0ddf94SBarry Smith 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 133b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 134b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 135b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 136b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1380a71e23eSMatthew Knepley 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1412eac72dbSBarry Smith 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 146ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1514b0e389bSBarry Smith 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1559f236934SBarry Smith 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1611d73ed98SBarry Smith 162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 163014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 1641d73ed98SBarry Smith 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 166014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 167014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 16858ad63f4SBarry Smith 169b49fd9e1SBarry Smith /*E 170b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 171b49fd9e1SBarry Smith 172b49fd9e1SBarry Smith Level: advanced 173b49fd9e1SBarry Smith 174bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 175e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 176b49fd9e1SBarry Smith 177b49fd9e1SBarry Smith E*/ 17878d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 1796a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 1801f7e983dSSatish Balay /*MC 1818c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 1828c5b8ba0SBarry Smith 1838c5b8ba0SBarry Smith Level: advanced 1848c5b8ba0SBarry Smith 1858c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 1868c5b8ba0SBarry Smith 187bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 188e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 1898c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 1908c5b8ba0SBarry Smith M*/ 1918c5b8ba0SBarry Smith 1921f7e983dSSatish Balay /*MC 1938c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 1948c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 1958c5b8ba0SBarry Smith poor orthogonality. 1968c5b8ba0SBarry Smith 1978c5b8ba0SBarry Smith Level: advanced 1988c5b8ba0SBarry Smith 1998c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2008c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2018c5b8ba0SBarry Smith 202bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 203e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2048c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2058c5b8ba0SBarry Smith M*/ 2068c5b8ba0SBarry Smith 2071f7e983dSSatish Balay /*MC 2088c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2098c5b8ba0SBarry Smith 2108c5b8ba0SBarry Smith Level: advanced 2118c5b8ba0SBarry Smith 2128c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2138c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2148c5b8ba0SBarry Smith 2158c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2168c5b8ba0SBarry Smith 217bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 218e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2198c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2208c5b8ba0SBarry Smith M*/ 2218c5b8ba0SBarry Smith 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 22408480c60SBarry Smith 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 228c38d4ed2SBarry Smith 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 232121fd945SKris Buschelman 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 236e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 237d9492815SBarry Smith 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2402eac72dbSBarry Smith 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 244014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 245390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 246390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 248816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 251e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 252e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 253e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 25584cb2905SBarry Smith 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 258c01c455dSBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 2657287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 2667287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 2671eb62cbbSBarry Smith 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 2721f7f0c4fSBarry Smith 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 27455849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 275ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 27655849f57SBarry Smith 27755849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 2781eb62cbbSBarry Smith 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 281db9b2ab1SHong Zhang 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 28483ab6a24SBarry Smith 28528ce4d24SBarry Smith /*E 2868a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 2878a4b9c5cSBarry Smith test routines. 2888a4b9c5cSBarry Smith 2898a4b9c5cSBarry Smith Level: advanced 2908a4b9c5cSBarry Smith 291a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 2921718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 293a3f661c8SBarry Smith 2948a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 2958a4b9c5cSBarry Smith 29694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 2971718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 2988a4b9c5cSBarry Smith E*/ 2999e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3009e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 301014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 302a21b2a99SBarry Smith 3031f7e983dSSatish Balay /*MC 3049793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3058c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3068c5b8ba0SBarry Smith be based on a norm of a residual etc. 3078c5b8ba0SBarry Smith 3088c5b8ba0SBarry Smith Level: advanced 3098c5b8ba0SBarry Smith 310085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3118c5b8ba0SBarry Smith 312ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3138c5b8ba0SBarry Smith M*/ 3148c5b8ba0SBarry Smith 3151f7e983dSSatish Balay /*MC 316ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3178c5b8ba0SBarry Smith convergence test routine. 3188c5b8ba0SBarry Smith 3198c5b8ba0SBarry Smith Level: advanced 3208c5b8ba0SBarry Smith 3219793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3228c5b8ba0SBarry Smith M*/ 3238c5b8ba0SBarry Smith 3241f7e983dSSatish Balay /*MC 325ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3268c5b8ba0SBarry Smith convergence test routine. 3278c5b8ba0SBarry Smith 3288c5b8ba0SBarry Smith Level: advanced 3298c5b8ba0SBarry Smith 3309793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3318c5b8ba0SBarry Smith M*/ 3328c5b8ba0SBarry Smith 3331f7e983dSSatish Balay /*MC 334ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 335a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3368c5b8ba0SBarry Smith 3378c5b8ba0SBarry Smith Level: advanced 3388c5b8ba0SBarry Smith 3399793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3408c5b8ba0SBarry Smith M*/ 3418c5b8ba0SBarry Smith 342014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 343014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 344014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 346014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3478a4b9c5cSBarry Smith 3488a4b9c5cSBarry Smith /*E 34928ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 35028ce4d24SBarry Smith have converged or diverged 35128ce4d24SBarry Smith 35228ce4d24SBarry Smith Level: beginner 35328ce4d24SBarry Smith 3544d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 35528ce4d24SBarry Smith 3564d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3574d0a8057SBarry Smith 3584d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3594d0a8057SBarry Smith any of the values here also change them that array of names. 36086c02ca4SBarry Smith 361c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 36228ce4d24SBarry Smith E*/ 363d15094e1SBarry Smith typedef enum {/* converged */ 3649ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 3659ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 366d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 367d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 368b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3698031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3708031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 371329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 372af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 373d15094e1SBarry Smith /* diverged */ 374b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 375d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 376d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 377d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 378b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 379b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 380b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 3814d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 3826aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 383d15094e1SBarry Smith 384d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 385014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 386d15094e1SBarry Smith 387c838673bSBarry Smith /*MC 388c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 389c838673bSBarry Smith 390c838673bSBarry Smith Level: beginner 391c838673bSBarry Smith 392c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 393c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 394c838673bSBarry Smith 2-norm of the residual for right preconditioning 395c838673bSBarry Smith 396c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 397c838673bSBarry Smith 398c838673bSBarry Smith M*/ 399c838673bSBarry Smith 400c838673bSBarry Smith /*MC 401c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 402c838673bSBarry Smith 403c838673bSBarry Smith Level: beginner 404c838673bSBarry Smith 405c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 406c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 407c838673bSBarry Smith 2-norm of the residual for right preconditioning 408c838673bSBarry Smith 409c838673bSBarry Smith Level: beginner 410c838673bSBarry Smith 411c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 412c838673bSBarry Smith 413c838673bSBarry Smith M*/ 414c838673bSBarry Smith 415c838673bSBarry Smith /*MC 416c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 417c838673bSBarry Smith 418c838673bSBarry Smith Level: beginner 419c838673bSBarry Smith 420c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 421c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 422c838673bSBarry Smith 2-norm of the residual for right preconditioning 423c838673bSBarry Smith 424c838673bSBarry Smith Level: beginner 425c838673bSBarry Smith 426c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 427c838673bSBarry Smith 428c838673bSBarry Smith M*/ 429c838673bSBarry Smith 430c838673bSBarry Smith /*MC 431c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 432c838673bSBarry Smith reached 433c838673bSBarry Smith 434c838673bSBarry Smith Level: beginner 435c838673bSBarry Smith 436c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 437c838673bSBarry Smith 438c838673bSBarry Smith M*/ 439c838673bSBarry Smith 440c838673bSBarry Smith /*MC 4418c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4420059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 4434d0a8057SBarry Smith test routine is set in KSP. 444c838673bSBarry Smith 445c838673bSBarry Smith 446c838673bSBarry Smith Level: beginner 447c838673bSBarry Smith 448c838673bSBarry Smith 449c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 450c838673bSBarry Smith 451c838673bSBarry Smith M*/ 452c838673bSBarry Smith 453c838673bSBarry Smith /*MC 454c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4553014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4563014e516SBarry Smith preconditioner. 457c838673bSBarry Smith 458c838673bSBarry Smith Level: beginner 459c838673bSBarry Smith 460c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 461c838673bSBarry Smith 462c838673bSBarry Smith M*/ 463c838673bSBarry Smith 464c838673bSBarry Smith /*MC 465c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 466c838673bSBarry Smith method could not continue to enlarge the Krylov space. 467c838673bSBarry Smith 468c838673bSBarry Smith 469c838673bSBarry Smith Level: beginner 470c838673bSBarry Smith 471c838673bSBarry Smith 472c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 473c838673bSBarry Smith 474c838673bSBarry Smith M*/ 475c838673bSBarry Smith 476c838673bSBarry Smith /*MC 477c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 478c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 479c838673bSBarry Smith 480c838673bSBarry Smith Level: beginner 481c838673bSBarry Smith 482c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 483c838673bSBarry Smith 484c838673bSBarry Smith M*/ 485c838673bSBarry Smith 486c838673bSBarry Smith /*MC 487c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 488c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 489c838673bSBarry Smith be positive definite 490c838673bSBarry Smith 491c838673bSBarry Smith Level: beginner 492c838673bSBarry Smith 4932401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 494c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 495c838673bSBarry Smith 496c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 497c838673bSBarry Smith 498c838673bSBarry Smith M*/ 499c838673bSBarry Smith 500c838673bSBarry Smith /*MC 501c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 502c838673bSBarry Smith while the KSPSolve() is still running. 503c838673bSBarry Smith 504c838673bSBarry Smith Level: beginner 505c838673bSBarry Smith 506c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 507c838673bSBarry Smith 508c838673bSBarry Smith M*/ 509c838673bSBarry Smith 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 5128de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 5148de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 5158de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 5168de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5178de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5180059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 520abef13c0SSatish Balay 5218ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5228ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5238ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5248ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5258ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5268ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5278ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5288ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5298ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 5308ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 5318ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 5328ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 5338ea1b3e6SJed Brown 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 535d4fbbf0eSBarry Smith 53628ce4d24SBarry Smith /*E 53728ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 53828ce4d24SBarry Smith 53928ce4d24SBarry Smith Level: beginner 54028ce4d24SBarry Smith 54128ce4d24SBarry Smith .seealso: KSPCGSetType() 54228ce4d24SBarry Smith E*/ 5439dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5446a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 54528ce4d24SBarry Smith 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5488031f4adStmunson 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 552fcae7a14Stmunson 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 556e559a7feSLois Curfman McInnes 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5628031f4adStmunson 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5641d6018f0SLisandro Dalcin 565014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5673369ce9aSBarry Smith 5689804daf3SBarry Smith #include <petscdrawtypes.h> 5691d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 5701d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 5711d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 5762f2e5d10SKris Buschelman 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 57903919abeSBarry Smith 580f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 581ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 582f8a50e2bSBarry Smith 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 589f8a50e2bSBarry Smith 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 593f8a50e2bSBarry Smith 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 596d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 5971324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 600*5a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 6023f22127dSBarry Smith 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 608fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 610fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 612014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 613014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 614014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 615fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 616fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6176c699258SBarry Smith 6182eac72dbSBarry Smith #endif 619