1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 426bd1501SBarry Smith #ifndef PETSCKSP_H 526bd1501SBarry Smith #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 1695452b02SPatrick Sanan Notes: 17db3b2ab5SMatthew Knepley When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a 18db3b2ab5SMatthew Knepley KSPType of KSPPREONLY, meaning that only application of the preconditioner is used as the linear solver. 198f6c3df8SBarry Smith 20bca11509SBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES 2128ce4d24SBarry Smith S*/ 22e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 232eac72dbSBarry Smith 2476bdecfbSBarry Smith /*J 258f6c3df8SBarry Smith KSPType - String with the name of a PETSc Krylov method. 2628ce4d24SBarry Smith 2728ce4d24SBarry Smith Level: beginner 2828ce4d24SBarry Smith 298f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 3076bdecfbSBarry Smith J*/ 3119fd82e9SBarry Smith typedef const char* KSPType; 3282bf6240SBarry Smith #define KSPRICHARDSON "richardson" 336c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3482bf6240SBarry Smith #define KSPCG "cg" 352c8fc646SJed Brown #define KSPGROPPCG "groppcg" 362c8fc646SJed Brown #define KSPPIPECG "pipecg" 37901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 38265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 39b21a8899STyler Chen #define KSPPIPEPRCG "pipeprcg" 40325e8391SManasi #define KSPPIPECG2 "pipecg2" 41df2a969aSvictorle #define KSPCGNE "cgne" 4205de396fSBarry Smith #define KSPNASH "nash" 4305de396fSBarry Smith #define KSPSTCG "stcg" 4405de396fSBarry Smith #define KSPGLTR "gltr" 4505de396fSBarry Smith #define KSPCGNASH PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash" 4605de396fSBarry Smith #define KSPCGSTCG PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg" 4705de396fSBarry Smith #define KSPCGGLTR PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr" 48640f4f53SPatrick Sanan #define KSPFCG "fcg" 49390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 5082bf6240SBarry Smith #define KSPGMRES "gmres" 51483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 529dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 539dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 544f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5561b468f9SJed Brown #define KSPPGMRES "pgmres" 5682bf6240SBarry Smith #define KSPTCQMR "tcqmr" 5782bf6240SBarry Smith #define KSPBCGS "bcgs" 58e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5918ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 60c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 61f0037002Svictorle #define KSPBCGSL "bcgsl" 62f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 6382bf6240SBarry Smith #define KSPCGS "cgs" 6482bf6240SBarry Smith #define KSPTFQMR "tfqmr" 6582bf6240SBarry Smith #define KSPCR "cr" 662c8fc646SJed Brown #define KSPPIPECR "pipecr" 6782bf6240SBarry Smith #define KSPLSQR "lsqr" 6882bf6240SBarry Smith #define KSPPREONLY "preonly" 6982bf6240SBarry Smith #define KSPQCG "qcg" 70c9cf9b11SBarry Smith #define KSPBICG "bicg" 71b4ac9ba4SBarry Smith #define KSPMINRES "minres" 7201c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 7321356919SSatish Balay #define KSPLCD "lcd" 741d6018f0SLisandro Dalcin #define KSPPYTHON "python" 7558ad63f4SBarry Smith #define KSPGCR "gcr" 76fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 77e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 78e4d80e07Szianekhodja #define KSPCGLS "cgls" 79329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 8038cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 812eac72dbSBarry Smith 828ba1e511SMatthew Knepley /* Logging support */ 83014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 84ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 855358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 868ba1e511SMatthew Knepley 87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*); 8819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 89c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 9475f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP,PetscBool); 955ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat); 963e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP,PetscInt); 973e00ce70SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp,PetscInt n) {return KSPSetMatSolveBatchSize(ksp,n);} 983e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP,PetscInt*); 993e00ce70SStefano Zampini PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp,PetscInt *n) {return KSPGetMatSolveBatchSize(ksp,n);} 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 101ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 10323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 1047d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*); 10525c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 106c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec); 1072eac72dbSBarry Smith 108140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 109ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 110798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList; 111798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 112798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 113bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 114798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **)); 11530de9b25SBarry Smith 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 119c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*); 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*); 124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool); 1257d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool); 126c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*); 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool); 128c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 133734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1342a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 13525ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1362eac72dbSBarry Smith 137d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 138d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 139d4211eb9SBarry Smith 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 142aabeff55SBarry Smith 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**)); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 1463ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void*); 14794a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool); 14994a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*); 15094a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool); 1514b0e389bSBarry Smith 152fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 153fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*); 154fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 155fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 156fa0ddf94SBarry Smith 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 158cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP); 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 160014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 161014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 162014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 163285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]); 164b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 165b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 166b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 167b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 168014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*); 1694a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*); 170f3b08a26SMatthew G. Knepley /* 171f3b08a26SMatthew G. Knepley PCMGCoarseList contains the list of coarse space constructor currently registered 172f3b08a26SMatthew G. Knepley These are added with PCMGRegisterCoarseSpaceConstructor() 173f3b08a26SMatthew G. Knepley */ 174f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 175f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)); 176f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)); 177f3b08a26SMatthew G. Knepley 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*); 1802eac72dbSBarry Smith 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 18458450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 185b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 18658450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 188d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*); 189d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1907d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]); 1914b0e389bSBarry Smith 192640f4f53SPatrick Sanan /*E 193640f4f53SPatrick Sanan 19406137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 195640f4f53SPatrick Sanan 19693f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 19793f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 198640f4f53SPatrick Sanan 1992b26979fSBarry Smith Level: intermediate 20006137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 201640f4f53SPatrick Sanan 202640f4f53SPatrick Sanan E*/ 20393f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 20406137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 205640f4f53SPatrick Sanan 206640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 207640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 208640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 209640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 21006137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 21106137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 212640f4f53SPatrick Sanan 213390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 214390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 215390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 216390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 217390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 218390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 219390d8e47SPatrick Sanan 220fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 221fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 222fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 223fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 224fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 225fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 226fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 227fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 22883f0b094SBarry Smith 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 232e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal); 2339f236934SBarry Smith 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2391d73ed98SBarry Smith 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2421d73ed98SBarry Smith 243483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 244483d6965SPatrick Sanan 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 24858ad63f4SBarry Smith 24904d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*); 25004d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC); 25104d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*); 252568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat); 253e82af88dSprj- 254e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat); 255e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*); 2565ea7661aSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); } 257d55ab951SPierre Jolivet /*E 258d55ab951SPierre Jolivet KSPHPDDMType - Type of Krylov method used by KSPHPDDM 259d55ab951SPierre Jolivet 260d55ab951SPierre Jolivet Level: intermediate 261d55ab951SPierre Jolivet 262d55ab951SPierre Jolivet Values: 263d55ab951SPierre Jolivet + KSP_HPDDM_TYPE_GMRES (default) 264d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BGMRES 265d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_CG 266d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BCG 267d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_GCRODR 268d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BGCRODR 269d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BFBCG 270d55ab951SPierre Jolivet - KSP_HPDDM_TYPE_PREONLY 271d55ab951SPierre Jolivet 272d55ab951SPierre Jolivet .seealso: KSPHPDDM, KSPHPDDMSetType() 273d55ab951SPierre Jolivet E*/ 274d55ab951SPierre Jolivet typedef enum { 275d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 276d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 277d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 278d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 279d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 280d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 281d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 282d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 283d55ab951SPierre Jolivet } KSPHPDDMType; 284d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 285d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType); 286d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*); 287e82af88dSprj- 288b49fd9e1SBarry Smith /*E 289b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 290b49fd9e1SBarry Smith 291b49fd9e1SBarry Smith Level: advanced 292b49fd9e1SBarry Smith 293bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 294e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 295b49fd9e1SBarry Smith 296b49fd9e1SBarry Smith E*/ 29778d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2986a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2991f7e983dSSatish Balay /*MC 3001957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 3018c5b8ba0SBarry Smith 3028c5b8ba0SBarry Smith Level: advanced 3038c5b8ba0SBarry Smith 3048c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 3058c5b8ba0SBarry Smith 306bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 307e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 3088c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 3098c5b8ba0SBarry Smith M*/ 3108c5b8ba0SBarry Smith 3111f7e983dSSatish Balay /*MC 3128c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 3138c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 3148c5b8ba0SBarry Smith poor orthogonality. 3158c5b8ba0SBarry Smith 3168c5b8ba0SBarry Smith Level: advanced 3178c5b8ba0SBarry Smith 3188c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 3198c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 3208c5b8ba0SBarry Smith 321bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 322e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 3238c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 3248c5b8ba0SBarry Smith M*/ 3258c5b8ba0SBarry Smith 3261f7e983dSSatish Balay /*MC 3278c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 3288c5b8ba0SBarry Smith 3298c5b8ba0SBarry Smith Level: advanced 3308c5b8ba0SBarry Smith 3318c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 3328c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 3338c5b8ba0SBarry Smith 3348c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 3358c5b8ba0SBarry Smith 336bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 337e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 3388c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 3398c5b8ba0SBarry Smith M*/ 3408c5b8ba0SBarry Smith 341014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 342014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 34308480c60SBarry Smith 344014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 346014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 347c38d4ed2SBarry Smith 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 350014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 351121fd945SKris Buschelman 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 353014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool); 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 355e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 356d9492815SBarry Smith 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 35887d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 3602eac72dbSBarry Smith 361798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],void*); 362798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char *[],int,int,int,int,PetscDrawLG *); 363798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 364798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 365798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 366798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 367798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 368798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 369798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 370798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 371798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 372798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 373798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 374798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 375798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 376798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 377798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 378798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 379798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 380798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 381798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 382fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 383798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 3842bd2df7bSSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorResidual(ksp,n,rnorm,vf);} 3852bd2df7bSSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidual(ksp,n,rnorm,vf);} 3862bd2df7bSSatish Balay PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidualMax(ksp,n,rnorm,vf);} 387798534f6SMatthew G. Knepley 388798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*); 389390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 390390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 391e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 392e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 393e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 39484cb2905SBarry Smith 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 396014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 397c01c455dSBarry Smith 39823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 39923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 400014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*); 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 4041eb62cbbSBarry Smith 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool); 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*); 4091f7f0c4fSBarry Smith 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 41155849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 412fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]); 41319a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer); 414c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP,PetscErrorCode (*)(KSP,void*),void *vctx,PetscErrorCode (*)(void**)); 4151b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 416c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 41794a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer); 4181b2b9847SBarry Smith 41919a666eeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);} 4201b2b9847SBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);} 42155849f57SBarry Smith 42255849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 4231eb62cbbSBarry Smith 424dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool); 4250e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool); 426014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 427884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*); 428798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 429798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 430798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 431db9b2ab1SHong Zhang 432014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 433014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 43468ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 43583ab6a24SBarry Smith 43628ce4d24SBarry Smith /*E 4378a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 4388a4b9c5cSBarry Smith test routines. 4398a4b9c5cSBarry Smith 4408a4b9c5cSBarry Smith Level: advanced 4418a4b9c5cSBarry Smith 442a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 4431718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 444a3f661c8SBarry Smith 44595452b02SPatrick Sanan Notes: 44695452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 4478a4b9c5cSBarry Smith 44894b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 4491718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 4508a4b9c5cSBarry Smith E*/ 4519e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 4529e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 453014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 454a21b2a99SBarry Smith 4551f7e983dSSatish Balay /*MC 4569793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 4578c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 4588c5b8ba0SBarry Smith be based on a norm of a residual etc. 4598c5b8ba0SBarry Smith 4608c5b8ba0SBarry Smith Level: advanced 4618c5b8ba0SBarry Smith 4621957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 4638c5b8ba0SBarry Smith 464ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 4658c5b8ba0SBarry Smith M*/ 4668c5b8ba0SBarry Smith 4671f7e983dSSatish Balay /*MC 4681957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 4698c5b8ba0SBarry Smith convergence test routine. 4708c5b8ba0SBarry Smith 4718c5b8ba0SBarry Smith Level: advanced 4728c5b8ba0SBarry Smith 4739793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 4748c5b8ba0SBarry Smith M*/ 4758c5b8ba0SBarry Smith 4761f7e983dSSatish Balay /*MC 477ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 4788c5b8ba0SBarry Smith convergence test routine. 4798c5b8ba0SBarry Smith 4808c5b8ba0SBarry Smith Level: advanced 4818c5b8ba0SBarry Smith 4829793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 4838c5b8ba0SBarry Smith M*/ 4848c5b8ba0SBarry Smith 4851f7e983dSSatish Balay /*MC 486ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 487390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 4888c5b8ba0SBarry Smith 4898c5b8ba0SBarry Smith Level: advanced 4908c5b8ba0SBarry Smith 4919793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 4928c5b8ba0SBarry Smith M*/ 4938c5b8ba0SBarry Smith 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 4998a4b9c5cSBarry Smith 500f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 5018a4b9c5cSBarry Smith /*E 5021957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 50328ce4d24SBarry Smith 50428ce4d24SBarry Smith Level: beginner 50528ce4d24SBarry Smith 50695452b02SPatrick Sanan Notes: 50795452b02SPatrick Sanan See KSPGetConvergedReason() for explanation of each value 50828ce4d24SBarry Smith 50995452b02SPatrick Sanan Developer Notes: 51095452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 5114d0a8057SBarry Smith 5124d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 5134d0a8057SBarry Smith any of the values here also change them that array of names. 51486c02ca4SBarry Smith 5151b2b9847SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances(), KSPConvergedReasonView() 51628ce4d24SBarry Smith E*/ 517d15094e1SBarry Smith typedef enum {/* converged */ 5189ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 5199ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 520d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 521d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 522b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 5238031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 5248031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 525329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 526af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 527d15094e1SBarry Smith /* diverged */ 528b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 529d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 530d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 531d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 532b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 533b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 534b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 5354d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 5366aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 537c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 538aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 539d15094e1SBarry Smith 540d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 541014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 542d15094e1SBarry Smith 543c838673bSBarry Smith /*MC 544f9fed41fSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called 545c838673bSBarry Smith 546c838673bSBarry Smith Level: beginner 547c838673bSBarry Smith 548c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 549c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 550c838673bSBarry Smith 2-norm of the residual for right preconditioning 551c838673bSBarry Smith 552f9fed41fSBarry Smith See also KSP_CONVERGED_ATOL which may apply before this tolerance. 553f9fed41fSBarry Smith 554f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 555c838673bSBarry Smith 556c838673bSBarry Smith M*/ 557c838673bSBarry Smith 558c838673bSBarry Smith /*MC 559c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 560c838673bSBarry Smith 561c838673bSBarry Smith Level: beginner 562c838673bSBarry Smith 563c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 564c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 565c838673bSBarry Smith 2-norm of the residual for right preconditioning 566c838673bSBarry Smith 567f9fed41fSBarry Smith See also KSP_CONVERGED_RTOL which may apply before this tolerance. 568c838673bSBarry Smith 569f9fed41fSBarry Smith .seealso: KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 570c838673bSBarry Smith 571c838673bSBarry Smith M*/ 572c838673bSBarry Smith 573c838673bSBarry Smith /*MC 574c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 575c838673bSBarry Smith 576c838673bSBarry Smith Level: beginner 577c838673bSBarry Smith 578c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 579c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 580c838673bSBarry Smith 2-norm of the residual for right preconditioning 581c838673bSBarry Smith 582c838673bSBarry Smith Level: beginner 583c838673bSBarry Smith 584f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 585c838673bSBarry Smith 586c838673bSBarry Smith M*/ 587c838673bSBarry Smith 588c838673bSBarry Smith /*MC 589c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 590c838673bSBarry Smith reached 591c838673bSBarry Smith 592c838673bSBarry Smith Level: beginner 593c838673bSBarry Smith 594c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 595c838673bSBarry Smith 596c838673bSBarry Smith M*/ 597c838673bSBarry Smith 598c838673bSBarry Smith /*MC 5998c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 6000059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 6014d0a8057SBarry Smith test routine is set in KSP. 602c838673bSBarry Smith 603c838673bSBarry Smith Level: beginner 604c838673bSBarry Smith 605c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 606c838673bSBarry Smith 607c838673bSBarry Smith M*/ 608c838673bSBarry Smith 609c838673bSBarry Smith /*MC 610c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 6111de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 6121de96524SPierre Jolivet preconditioner. In KSPHPDDM, this is also returned when some search directions within a block 6131de96524SPierre Jolivet are colinear. 614c838673bSBarry Smith 615c838673bSBarry Smith Level: beginner 616c838673bSBarry Smith 617c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 618c838673bSBarry Smith 619c838673bSBarry Smith M*/ 620c838673bSBarry Smith 621c838673bSBarry Smith /*MC 622c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 623c838673bSBarry Smith method could not continue to enlarge the Krylov space. 624c838673bSBarry Smith 625c838673bSBarry Smith Level: beginner 626c838673bSBarry Smith 627c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 628c838673bSBarry Smith 629c838673bSBarry Smith M*/ 630c838673bSBarry Smith 631c838673bSBarry Smith /*MC 632c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 633c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 634c838673bSBarry Smith 635c838673bSBarry Smith Level: beginner 636c838673bSBarry Smith 637c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 638c838673bSBarry Smith 639c838673bSBarry Smith M*/ 640c838673bSBarry Smith 641c838673bSBarry Smith /*MC 642c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 643c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 644c838673bSBarry Smith be positive definite 645c838673bSBarry Smith 646c838673bSBarry Smith Level: beginner 647c838673bSBarry Smith 64895452b02SPatrick Sanan Notes: 64995452b02SPatrick Sanan This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 650c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 651c838673bSBarry Smith 652c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 653c838673bSBarry Smith 654c838673bSBarry Smith M*/ 655c838673bSBarry Smith 656c838673bSBarry Smith /*MC 657c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 6589fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 6599fc87aa7SBarry Smith such as PCFIELDSPLIT. 6609fc87aa7SBarry Smith 6619fc87aa7SBarry Smith Level: beginner 6629fc87aa7SBarry Smith 66395452b02SPatrick Sanan Notes: 66495452b02SPatrick Sanan Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 6659fc87aa7SBarry Smith 6669fc87aa7SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 6679fc87aa7SBarry Smith 6689fc87aa7SBarry Smith M*/ 6699fc87aa7SBarry Smith 6709fc87aa7SBarry Smith /*MC 671c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 672c838673bSBarry Smith while the KSPSolve() is still running. 673c838673bSBarry Smith 674c838673bSBarry Smith Level: beginner 675c838673bSBarry Smith 676c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 677c838673bSBarry Smith 678c838673bSBarry Smith M*/ 679c838673bSBarry Smith 680014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*)); 681df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 682df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 6833ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void*); 6848de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 6853eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 6868de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*); 6878de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**); 6888de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 6898de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 69054b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool); 6910059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 692014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*); 693c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP,const char**); 69494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*); 695abef13c0SSatish Balay 69625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 6978ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 69825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 6998ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 70025ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 7018ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 70225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 7038ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 70425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 7058ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 70625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 7078ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 7088ea1b3e6SJed Brown 7090bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*); 71025ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); } 711d4fbbf0eSBarry Smith 71228ce4d24SBarry Smith /*E 71328ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 71428ce4d24SBarry Smith 71528ce4d24SBarry Smith Level: beginner 71628ce4d24SBarry Smith 71728ce4d24SBarry Smith .seealso: KSPCGSetType() 71828ce4d24SBarry Smith E*/ 7199dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 7206a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 72128ce4d24SBarry Smith 722014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 723014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool); 7248031f4adStmunson 725dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal); 726dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*); 727dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*); 728fcae7a14Stmunson 72905de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*); 73005de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*); 73125ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);} 73225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);} 7338031f4adStmunson 734014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 7351d6018f0SLisandro Dalcin 736*f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC,PetscErrorCode (*)(PC,KSP)); 737014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 738014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 7393369ce9aSBarry Smith 7409804daf3SBarry Smith #include <petscdrawtypes.h> 741014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 7422f2e5d10SKris Buschelman 743014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 744014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 74503919abeSBarry Smith 746ba36735cSStefano Zampini /*S 747ba36735cSStefano Zampini KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 748f8a50e2bSBarry Smith 749ba36735cSStefano Zampini Level: beginner 750f8a50e2bSBarry Smith 751ba36735cSStefano Zampini .seealso: KSPCreate(), KSPSetGuessType(), KSPGuessType 752ba36735cSStefano Zampini S*/ 753ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess; 754ba36735cSStefano Zampini /*J 755ba36735cSStefano Zampini KSPGuessType - String with the name of a PETSc initial guess for Krylov methods. 756ba36735cSStefano Zampini 757ba36735cSStefano Zampini Level: beginner 758ba36735cSStefano Zampini 759ba36735cSStefano Zampini .seealso: KSPGuess 760ba36735cSStefano Zampini J*/ 761ba36735cSStefano Zampini typedef const char* KSPGuessType; 762ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 763ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 7641d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess)); 765ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess); 766ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*); 767ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer); 768ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*); 769ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*); 770ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType); 771ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*); 772ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 773ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec); 774ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec); 775ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 776ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt); 777014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 778ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool); 779ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*); 780f8a50e2bSBarry Smith 781470b340bSDmitry Karpeev /*E 782470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 783470b340bSDmitry Karpeev 784470b340bSDmitry Karpeev Level: intermediate 785470b340bSDmitry Karpeev 786470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 787470b340bSDmitry Karpeev E*/ 7880c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType; 789470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 790470b340bSDmitry Karpeev 791864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 792864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 793864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 794864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_USER = 3} MatLMVMSymBroydenScaleType; 795864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 79692f76d53SAlp Dener 797014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 798014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 799d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 800bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 801aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 802bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 803470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 804470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 8055bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 8065a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 807470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*); 808470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 8093f22127dSBarry Smith 81078e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*); 81178e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*); 81278e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*); 813864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 814864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 815864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 816864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 817864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 818cd929ea3SAlp Dener 819cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 820b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*); 821cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 822cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 823e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 8240ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 825cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 826cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 827cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 828cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 829cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 8302d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 831cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 832cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*); 833cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*); 834cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*); 83592f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 836cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*); 837cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*); 838864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 839864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 840cd929ea3SAlp Dener 841014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 842014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool); 843014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 844014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 845014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 846fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*); 84723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 848fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 84923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 85023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 851014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 852014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 853fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 854fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 8556c699258SBarry Smith 85602b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 8578c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, 858e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 859e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 860e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 861191494d9SMatthew G. Knepley PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec); 862557cf195SMatthew G. Knepley 863557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *); 864557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal); 8652eac72dbSBarry Smith #endif 866