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" 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); 75*c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 8223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 832eac72dbSBarry Smith 84140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 85014dd563SJed Brown PETSC_EXTERN PetscBool KSPRegisterAllCalled; 86607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegisterAll(void); 87bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 88607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void); 8930de9b25SBarry Smith 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 93*c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 101*c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 103*c60c7ad4SBarry 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*); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 1102a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1112a7a6963SBarry Smith PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1122eac72dbSBarry Smith 113d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 114d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 115d4211eb9SBarry Smith 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 118aabeff55SBarry Smith 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1254b0e389bSBarry Smith 126fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 127fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 128fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 129fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 130fa0ddf94SBarry Smith 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 136b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 137b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 138b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 139b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1410a71e23eSMatthew Knepley 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1442eac72dbSBarry Smith 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 149ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(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[]); 1544b0e389bSBarry Smith 155640f4f53SPatrick Sanan /*E 156640f4f53SPatrick Sanan 157bd242c34SBarry Smith KSPFCGTruncationType - Define how stored directions are used to orthogonalize in FCG 158640f4f53SPatrick Sanan 159bd242c34SBarry Smith KSP_FCG_TRUNCATION_TYPE_STANDARD uses all (up to mmax) stored directions 160bd242c34SBarry Smith KSP_FCG_TRUNCATION_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 161640f4f53SPatrick Sanan 162bd242c34SBarry Smith .seealso : KSPFCG,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 163640f4f53SPatrick Sanan 164640f4f53SPatrick Sanan E*/ 165bd242c34SBarry Smith typedef enum {KSP_FCG_TRUNCATION_TYPE_STANDARD,KSP_FCG_TRUNCATION_TYPE_NOTAY} KSPFCGTruncationType; 166bd242c34SBarry Smith PETSC_EXTERN const char *const KSPFCGTruncationTypes[]; 167640f4f53SPatrick Sanan 168640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 169640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 170640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 171640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 172bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCGTruncationType); 173bd242c34SBarry Smith PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCGTruncationType*); 174640f4f53SPatrick Sanan 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1789f236934SBarry Smith 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1841d73ed98SBarry Smith 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 1871d73ed98SBarry Smith 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 19158ad63f4SBarry Smith 192b49fd9e1SBarry Smith /*E 193b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 194b49fd9e1SBarry Smith 195b49fd9e1SBarry Smith Level: advanced 196b49fd9e1SBarry Smith 197bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 198e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 199b49fd9e1SBarry Smith 200b49fd9e1SBarry Smith E*/ 20178d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2026a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2031f7e983dSSatish Balay /*MC 2041957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2058c5b8ba0SBarry Smith 2068c5b8ba0SBarry Smith Level: advanced 2078c5b8ba0SBarry Smith 2088c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2098c5b8ba0SBarry Smith 210bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 211e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2128c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2138c5b8ba0SBarry Smith M*/ 2148c5b8ba0SBarry Smith 2151f7e983dSSatish Balay /*MC 2168c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2178c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2188c5b8ba0SBarry Smith poor orthogonality. 2198c5b8ba0SBarry Smith 2208c5b8ba0SBarry Smith Level: advanced 2218c5b8ba0SBarry Smith 2228c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2238c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2248c5b8ba0SBarry Smith 225bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 226e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2278c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2288c5b8ba0SBarry Smith M*/ 2298c5b8ba0SBarry Smith 2301f7e983dSSatish Balay /*MC 2318c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2328c5b8ba0SBarry Smith 2338c5b8ba0SBarry Smith Level: advanced 2348c5b8ba0SBarry Smith 2358c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2368c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2378c5b8ba0SBarry Smith 2388c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2398c5b8ba0SBarry Smith 240bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 241e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2428c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2438c5b8ba0SBarry Smith M*/ 2448c5b8ba0SBarry Smith 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 24708480c60SBarry Smith 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 251c38d4ed2SBarry Smith 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 255121fd945SKris Buschelman 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 259e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 260d9492815SBarry Smith 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2632eac72dbSBarry Smith 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 268390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 269390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 271816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 274e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 275e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 276e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 27884cb2905SBarry Smith 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 281c01c455dSBarry Smith 28223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 28323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 2887287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 2897287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 2901eb62cbbSBarry Smith 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 2951f7f0c4fSBarry Smith 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 29755849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 298ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 29955849f57SBarry Smith 30055849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3011eb62cbbSBarry Smith 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 304db9b2ab1SHong Zhang 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 30783ab6a24SBarry Smith 30828ce4d24SBarry Smith /*E 3098a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3108a4b9c5cSBarry Smith test routines. 3118a4b9c5cSBarry Smith 3128a4b9c5cSBarry Smith Level: advanced 3138a4b9c5cSBarry Smith 314a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3151718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 316a3f661c8SBarry Smith 3178a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3188a4b9c5cSBarry Smith 31994b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3201718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3218a4b9c5cSBarry Smith E*/ 3229e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3239e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 324014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 325a21b2a99SBarry Smith 3261f7e983dSSatish Balay /*MC 3279793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3288c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3298c5b8ba0SBarry Smith be based on a norm of a residual etc. 3308c5b8ba0SBarry Smith 3318c5b8ba0SBarry Smith Level: advanced 3328c5b8ba0SBarry Smith 3331957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3348c5b8ba0SBarry Smith 335ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3368c5b8ba0SBarry Smith M*/ 3378c5b8ba0SBarry Smith 3381f7e983dSSatish Balay /*MC 3391957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3408c5b8ba0SBarry Smith convergence test routine. 3418c5b8ba0SBarry Smith 3428c5b8ba0SBarry Smith Level: advanced 3438c5b8ba0SBarry Smith 3449793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3458c5b8ba0SBarry Smith M*/ 3468c5b8ba0SBarry Smith 3471f7e983dSSatish Balay /*MC 348ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3498c5b8ba0SBarry Smith convergence test routine. 3508c5b8ba0SBarry Smith 3518c5b8ba0SBarry Smith Level: advanced 3528c5b8ba0SBarry Smith 3539793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3548c5b8ba0SBarry Smith M*/ 3558c5b8ba0SBarry Smith 3561f7e983dSSatish Balay /*MC 357ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 358640f4f53SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG 3598c5b8ba0SBarry Smith 3608c5b8ba0SBarry Smith Level: advanced 3618c5b8ba0SBarry Smith 3629793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3638c5b8ba0SBarry Smith M*/ 3648c5b8ba0SBarry Smith 365014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 366014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 367014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 368014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3708a4b9c5cSBarry Smith 3718a4b9c5cSBarry Smith /*E 3721957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 37328ce4d24SBarry Smith 37428ce4d24SBarry Smith Level: beginner 37528ce4d24SBarry Smith 3764d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 37728ce4d24SBarry Smith 3784d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3794d0a8057SBarry Smith 3804d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3814d0a8057SBarry Smith any of the values here also change them that array of names. 38286c02ca4SBarry Smith 383c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 38428ce4d24SBarry Smith E*/ 385d15094e1SBarry Smith typedef enum {/* converged */ 3869ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 3879ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 388d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 389d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 390b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3918031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3928031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 393329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 394af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 395d15094e1SBarry Smith /* diverged */ 396b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 397d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 398d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 399d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 400b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 401b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 402b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4034d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4046aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 405d15094e1SBarry Smith 406d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 407014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 408d15094e1SBarry Smith 409c838673bSBarry Smith /*MC 410c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 411c838673bSBarry Smith 412c838673bSBarry Smith Level: beginner 413c838673bSBarry Smith 414c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 415c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 416c838673bSBarry Smith 2-norm of the residual for right preconditioning 417c838673bSBarry Smith 418c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 419c838673bSBarry Smith 420c838673bSBarry Smith M*/ 421c838673bSBarry Smith 422c838673bSBarry Smith /*MC 423c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 424c838673bSBarry Smith 425c838673bSBarry Smith Level: beginner 426c838673bSBarry Smith 427c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 428c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 429c838673bSBarry Smith 2-norm of the residual for right preconditioning 430c838673bSBarry Smith 431c838673bSBarry Smith Level: beginner 432c838673bSBarry Smith 433c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 434c838673bSBarry Smith 435c838673bSBarry Smith M*/ 436c838673bSBarry Smith 437c838673bSBarry Smith /*MC 438c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 439c838673bSBarry Smith 440c838673bSBarry Smith Level: beginner 441c838673bSBarry Smith 442c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 443c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 444c838673bSBarry Smith 2-norm of the residual for right preconditioning 445c838673bSBarry Smith 446c838673bSBarry Smith Level: beginner 447c838673bSBarry Smith 448c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 449c838673bSBarry Smith 450c838673bSBarry Smith M*/ 451c838673bSBarry Smith 452c838673bSBarry Smith /*MC 453c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 454c838673bSBarry Smith reached 455c838673bSBarry Smith 456c838673bSBarry Smith Level: beginner 457c838673bSBarry Smith 458c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 459c838673bSBarry Smith 460c838673bSBarry Smith M*/ 461c838673bSBarry Smith 462c838673bSBarry Smith /*MC 4638c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4640059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 4654d0a8057SBarry Smith test routine is set in KSP. 466c838673bSBarry Smith 467c838673bSBarry Smith Level: beginner 468c838673bSBarry Smith 469c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 470c838673bSBarry Smith 471c838673bSBarry Smith M*/ 472c838673bSBarry Smith 473c838673bSBarry Smith /*MC 474c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4753014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4763014e516SBarry Smith preconditioner. 477c838673bSBarry Smith 478c838673bSBarry Smith Level: beginner 479c838673bSBarry Smith 480c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 481c838673bSBarry Smith 482c838673bSBarry Smith M*/ 483c838673bSBarry Smith 484c838673bSBarry Smith /*MC 485c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 486c838673bSBarry Smith method could not continue to enlarge the Krylov space. 487c838673bSBarry Smith 488c838673bSBarry Smith Level: beginner 489c838673bSBarry Smith 490c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 491c838673bSBarry Smith 492c838673bSBarry Smith M*/ 493c838673bSBarry Smith 494c838673bSBarry Smith /*MC 495c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 496c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 497c838673bSBarry Smith 498c838673bSBarry Smith Level: beginner 499c838673bSBarry Smith 500c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 501c838673bSBarry Smith 502c838673bSBarry Smith M*/ 503c838673bSBarry Smith 504c838673bSBarry Smith /*MC 505c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 506c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 507c838673bSBarry Smith be positive definite 508c838673bSBarry Smith 509c838673bSBarry Smith Level: beginner 510c838673bSBarry Smith 5112401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 512c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 513c838673bSBarry Smith 514c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 515c838673bSBarry Smith 516c838673bSBarry Smith M*/ 517c838673bSBarry Smith 518c838673bSBarry Smith /*MC 519c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 520c838673bSBarry Smith while the KSPSolve() is still running. 521c838673bSBarry Smith 522c838673bSBarry Smith Level: beginner 523c838673bSBarry Smith 524c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 525c838673bSBarry Smith 526c838673bSBarry Smith M*/ 527c838673bSBarry Smith 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 5308de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 5328de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 5338de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 5348de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5358de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5360059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 538abef13c0SSatish Balay 5398ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5408ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5418ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5428ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5438ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5448ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5458ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5468ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5478ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 5488ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 5498ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 5508ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 5518ea1b3e6SJed Brown 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 553d4fbbf0eSBarry Smith 55428ce4d24SBarry Smith /*E 55528ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 55628ce4d24SBarry Smith 55728ce4d24SBarry Smith Level: beginner 55828ce4d24SBarry Smith 55928ce4d24SBarry Smith .seealso: KSPCGSetType() 56028ce4d24SBarry Smith E*/ 5619dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5626a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 56328ce4d24SBarry Smith 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 565014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5668031f4adStmunson 567014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 570fcae7a14Stmunson 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 574e559a7feSLois Curfman McInnes 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5808031f4adStmunson 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5821d6018f0SLisandro Dalcin 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5853369ce9aSBarry Smith 5869804daf3SBarry Smith #include <petscdrawtypes.h> 5871d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 5881d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 5891d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 5942f2e5d10SKris Buschelman 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 59703919abeSBarry Smith 598f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 599ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 600f8a50e2bSBarry Smith 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 605014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 607f8a50e2bSBarry Smith 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 611f8a50e2bSBarry Smith 612470b340bSDmitry Karpeev /*E 613470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 614470b340bSDmitry Karpeev 615470b340bSDmitry Karpeev Level: intermediate 616470b340bSDmitry Karpeev 617470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 618470b340bSDmitry Karpeev E*/ 619470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 620470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 621470b340bSDmitry Karpeev 622014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 624d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 625bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 626aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 627bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 628470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 629470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 6305bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 6315a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 632470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 633470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 6343f22127dSBarry Smith 635014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 640fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 64123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 642fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 64323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 64423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 645014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 646014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 647fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 648fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6496c699258SBarry Smith 65002b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 65102b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexProjectField(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec); 65202b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexProjectFunction(DM, void (**)(const PetscReal [], PetscScalar *, void *), void **, InsertMode, Vec); 6532eac72dbSBarry Smith #endif 654