1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCKSP_H 50a835dfdSSatish Balay #define __PETSCKSP_H 62c8e378dSBarry Smith #include <petscpc.h> 72eac72dbSBarry Smith 8607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 91dbb0a54SBarry Smith 1028ce4d24SBarry Smith /*S 118f6c3df8SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 128f6c3df8SBarry Smith linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 1328ce4d24SBarry Smith 1428ce4d24SBarry Smith Level: beginner 1528ce4d24SBarry Smith 1628ce4d24SBarry Smith Concepts: Krylov methods 1728ce4d24SBarry Smith 188f6c3df8SBarry Smith Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a 198f6c3df8SBarry Smith KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver). 208f6c3df8SBarry Smith 218f6c3df8SBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy() 2228ce4d24SBarry Smith S*/ 23e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 242eac72dbSBarry Smith 2576bdecfbSBarry Smith /*J 268f6c3df8SBarry Smith KSPType - String with the name of a PETSc Krylov method. 2728ce4d24SBarry Smith 2828ce4d24SBarry Smith Level: beginner 2928ce4d24SBarry Smith 308f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 3176bdecfbSBarry Smith J*/ 3219fd82e9SBarry Smith typedef const char* KSPType; 3382bf6240SBarry Smith #define KSPRICHARDSON "richardson" 346c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3582bf6240SBarry Smith #define KSPCG "cg" 362c8fc646SJed Brown #define KSPGROPPCG "groppcg" 372c8fc646SJed Brown #define KSPPIPECG "pipecg" 38df2a969aSvictorle #define KSPCGNE "cgne" 39fcae7a14Stmunson #define KSPNASH "nash" 4080e17431SMatthew Knepley #define KSPSTCG "stcg" 418031f4adStmunson #define KSPGLTR "gltr" 42640f4f53SPatrick Sanan #define KSPFCG "fcg" 4382bf6240SBarry Smith #define KSPGMRES "gmres" 449dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 459dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 464f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 4761b468f9SJed Brown #define KSPPGMRES "pgmres" 4882bf6240SBarry Smith #define KSPTCQMR "tcqmr" 4982bf6240SBarry Smith #define KSPBCGS "bcgs" 50e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5118ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 52c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 53f0037002Svictorle #define KSPBCGSL "bcgsl" 5482bf6240SBarry Smith #define KSPCGS "cgs" 5582bf6240SBarry Smith #define KSPTFQMR "tfqmr" 5682bf6240SBarry Smith #define KSPCR "cr" 572c8fc646SJed Brown #define KSPPIPECR "pipecr" 5882bf6240SBarry Smith #define KSPLSQR "lsqr" 5982bf6240SBarry Smith #define KSPPREONLY "preonly" 6082bf6240SBarry Smith #define KSPQCG "qcg" 61c9cf9b11SBarry Smith #define KSPBICG "bicg" 62b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6301c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 6421356919SSatish Balay #define KSPLCD "lcd" 651d6018f0SLisandro Dalcin #define KSPPYTHON "python" 6658ad63f4SBarry Smith #define KSPGCR "gcr" 67e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 68e4d80e07Szianekhodja #define KSPCGLS "cgls" 692eac72dbSBarry Smith 708ba1e511SMatthew Knepley /* Logging support */ 71014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 725358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 738ba1e511SMatthew Knepley 74014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 7519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 76c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 80014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 8323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 8425c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 852eac72dbSBarry Smith 86140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 87bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 8830de9b25SBarry Smith 89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 92c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 1007d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool ); 101c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 103c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 108734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1092a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1102a7a6963SBarry Smith PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1112eac72dbSBarry Smith 112d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 113d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 114d4211eb9SBarry Smith 115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 117aabeff55SBarry Smith 118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1244b0e389bSBarry Smith 125fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 126fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 127fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 128fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 129fa0ddf94SBarry Smith 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 135b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 136b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 137b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 138b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1400a71e23eSMatthew Knepley 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1432eac72dbSBarry Smith 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 14758450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 14882fe98a5SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool); 149ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 15058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 152d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *); 153d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1547d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]); 1554b0e389bSBarry Smith 156640f4f53SPatrick Sanan /*E 157640f4f53SPatrick Sanan 158*06137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 159640f4f53SPatrick Sanan 160*06137d0aSPatrick Sanan KSP_FCD_TRUNCATION uses all (up to mmax) stored directions 161*06137d0aSPatrick Sanan KSP_FCD_TRUNCATION_RESTART uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 162640f4f53SPatrick Sanan 1632b26979fSBarry Smith Level: intermediate 164*06137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 165640f4f53SPatrick Sanan 166640f4f53SPatrick Sanan E*/ 167*06137d0aSPatrick Sanan typedef enum {KSP_FCD_TRUNCATION,KSP_FCD_TRUNCATION_RESTART} KSPFCDTruncationType; 168*06137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 169*06137d0aSPatrick Sanan 170*06137d0aSPatrick Sanan #define KSPFCDGetNumOldDirections(ctx,i,mi) ({ \ 171*06137d0aSPatrick Sanan if(ctx->trunctype == KSP_FCD_TRUNCATION_RESTART){ \ 172*06137d0aSPatrick Sanan mi = ((i-1) % ctx->mmax)+1; \ 173*06137d0aSPatrick Sanan if (mi==1 && i!=1) ++(ctx->n_search_space_resets); \ 174*06137d0aSPatrick Sanan } else if (ctx->trunctype == KSP_FCD_TRUNCATION) \ 175*06137d0aSPatrick Sanan mi = ctx->mmax; \ 176*06137d0aSPatrick Sanan else { \ 177*06137d0aSPatrick Sanan SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy");CHKERRQ(ierr); \ 178*06137d0aSPatrick Sanan } \ 179*06137d0aSPatrick Sanan }) 180640f4f53SPatrick Sanan 181640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 182640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 183640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 184640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 185*06137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 186*06137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 187640f4f53SPatrick Sanan 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1919f236934SBarry Smith 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1971d73ed98SBarry Smith 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2001d73ed98SBarry Smith 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 20458ad63f4SBarry Smith 205b49fd9e1SBarry Smith /*E 206b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 207b49fd9e1SBarry Smith 208b49fd9e1SBarry Smith Level: advanced 209b49fd9e1SBarry Smith 210bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 211e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 212b49fd9e1SBarry Smith 213b49fd9e1SBarry Smith E*/ 21478d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2156a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2161f7e983dSSatish Balay /*MC 2171957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2188c5b8ba0SBarry Smith 2198c5b8ba0SBarry Smith Level: advanced 2208c5b8ba0SBarry Smith 2218c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2228c5b8ba0SBarry Smith 223bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 224e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2258c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2268c5b8ba0SBarry Smith M*/ 2278c5b8ba0SBarry Smith 2281f7e983dSSatish Balay /*MC 2298c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2308c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2318c5b8ba0SBarry Smith poor orthogonality. 2328c5b8ba0SBarry Smith 2338c5b8ba0SBarry Smith Level: advanced 2348c5b8ba0SBarry Smith 2358c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2368c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2378c5b8ba0SBarry Smith 238bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 239e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2408c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2418c5b8ba0SBarry Smith M*/ 2428c5b8ba0SBarry Smith 2431f7e983dSSatish Balay /*MC 2448c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2458c5b8ba0SBarry Smith 2468c5b8ba0SBarry Smith Level: advanced 2478c5b8ba0SBarry Smith 2488c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2498c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2508c5b8ba0SBarry Smith 2518c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2528c5b8ba0SBarry Smith 253bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 254e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2558c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2568c5b8ba0SBarry Smith M*/ 2578c5b8ba0SBarry Smith 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 26008480c60SBarry Smith 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 264c38d4ed2SBarry Smith 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 268121fd945SKris Buschelman 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 272e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 273d9492815SBarry Smith 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2762eac72dbSBarry Smith 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 281390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 282390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 284816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 287e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 288e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 289e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 29184cb2905SBarry Smith 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 294c01c455dSBarry Smith 29523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 29623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3017287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3027287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3031eb62cbbSBarry Smith 304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3081f7f0c4fSBarry Smith 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 31055849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 311685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 3122a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 3132a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 31455849f57SBarry Smith 31555849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3161eb62cbbSBarry Smith 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 319db9b2ab1SHong Zhang 320014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 321014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 32268ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 32383ab6a24SBarry Smith 32428ce4d24SBarry Smith /*E 3258a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3268a4b9c5cSBarry Smith test routines. 3278a4b9c5cSBarry Smith 3288a4b9c5cSBarry Smith Level: advanced 3298a4b9c5cSBarry Smith 330a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3311718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 332a3f661c8SBarry Smith 333af0996ceSBarry Smith Notes: this must match petsc/finclude/petscksp.h 3348a4b9c5cSBarry Smith 33594b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3361718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3378a4b9c5cSBarry Smith E*/ 3389e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3399e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 340014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 341a21b2a99SBarry Smith 3421f7e983dSSatish Balay /*MC 3439793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3448c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3458c5b8ba0SBarry Smith be based on a norm of a residual etc. 3468c5b8ba0SBarry Smith 3478c5b8ba0SBarry Smith Level: advanced 3488c5b8ba0SBarry Smith 3491957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3508c5b8ba0SBarry Smith 351ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3528c5b8ba0SBarry Smith M*/ 3538c5b8ba0SBarry Smith 3541f7e983dSSatish Balay /*MC 3551957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3568c5b8ba0SBarry Smith convergence test routine. 3578c5b8ba0SBarry Smith 3588c5b8ba0SBarry Smith Level: advanced 3598c5b8ba0SBarry Smith 3609793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3618c5b8ba0SBarry Smith M*/ 3628c5b8ba0SBarry Smith 3631f7e983dSSatish Balay /*MC 364ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3658c5b8ba0SBarry Smith convergence test routine. 3668c5b8ba0SBarry Smith 3678c5b8ba0SBarry Smith Level: advanced 3688c5b8ba0SBarry Smith 3699793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3708c5b8ba0SBarry Smith M*/ 3718c5b8ba0SBarry Smith 3721f7e983dSSatish Balay /*MC 373ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 374640f4f53SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG 3758c5b8ba0SBarry Smith 3768c5b8ba0SBarry Smith Level: advanced 3778c5b8ba0SBarry Smith 3789793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3798c5b8ba0SBarry Smith M*/ 3808c5b8ba0SBarry Smith 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3868a4b9c5cSBarry Smith 3878a4b9c5cSBarry Smith /*E 3881957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 38928ce4d24SBarry Smith 39028ce4d24SBarry Smith Level: beginner 39128ce4d24SBarry Smith 3924d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 39328ce4d24SBarry Smith 394af0996ceSBarry Smith Developer notes: this must match petsc/finclude/petscksp.h 3954d0a8057SBarry Smith 3964d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3974d0a8057SBarry Smith any of the values here also change them that array of names. 39886c02ca4SBarry Smith 399c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 40028ce4d24SBarry Smith E*/ 401d15094e1SBarry Smith typedef enum {/* converged */ 4029ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4039ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 404d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 405d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 406b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4078031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4088031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 409329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 410af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 411d15094e1SBarry Smith /* diverged */ 412b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 413d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 414d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 415d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 416b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 417b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 418b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4194d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4206aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 421422a814eSBarry Smith KSP_DIVERGED_PCSETUP_FAILED = -11, 422d15094e1SBarry Smith 423d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 424014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 425d15094e1SBarry Smith 426c838673bSBarry Smith /*MC 427c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 428c838673bSBarry Smith 429c838673bSBarry Smith Level: beginner 430c838673bSBarry Smith 431c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 432c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 433c838673bSBarry Smith 2-norm of the residual for right preconditioning 434c838673bSBarry Smith 435c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 436c838673bSBarry Smith 437c838673bSBarry Smith M*/ 438c838673bSBarry Smith 439c838673bSBarry Smith /*MC 440c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 441c838673bSBarry Smith 442c838673bSBarry Smith Level: beginner 443c838673bSBarry Smith 444c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 445c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 446c838673bSBarry Smith 2-norm of the residual for right preconditioning 447c838673bSBarry Smith 448c838673bSBarry Smith Level: beginner 449c838673bSBarry Smith 450c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 451c838673bSBarry Smith 452c838673bSBarry Smith M*/ 453c838673bSBarry Smith 454c838673bSBarry Smith /*MC 455c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 456c838673bSBarry Smith 457c838673bSBarry Smith Level: beginner 458c838673bSBarry Smith 459c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 460c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 461c838673bSBarry Smith 2-norm of the residual for right preconditioning 462c838673bSBarry Smith 463c838673bSBarry Smith Level: beginner 464c838673bSBarry Smith 465c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 466c838673bSBarry Smith 467c838673bSBarry Smith M*/ 468c838673bSBarry Smith 469c838673bSBarry Smith /*MC 470c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 471c838673bSBarry Smith reached 472c838673bSBarry Smith 473c838673bSBarry Smith Level: beginner 474c838673bSBarry Smith 475c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 476c838673bSBarry Smith 477c838673bSBarry Smith M*/ 478c838673bSBarry Smith 479c838673bSBarry Smith /*MC 4808c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4810059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 4824d0a8057SBarry Smith test routine is set in KSP. 483c838673bSBarry Smith 484c838673bSBarry Smith Level: beginner 485c838673bSBarry Smith 486c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 487c838673bSBarry Smith 488c838673bSBarry Smith M*/ 489c838673bSBarry Smith 490c838673bSBarry Smith /*MC 491c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4923014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4933014e516SBarry Smith preconditioner. 494c838673bSBarry Smith 495c838673bSBarry Smith Level: beginner 496c838673bSBarry Smith 497c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 498c838673bSBarry Smith 499c838673bSBarry Smith M*/ 500c838673bSBarry Smith 501c838673bSBarry Smith /*MC 502c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 503c838673bSBarry Smith method could not continue to enlarge the Krylov space. 504c838673bSBarry Smith 505c838673bSBarry Smith Level: beginner 506c838673bSBarry Smith 507c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 508c838673bSBarry Smith 509c838673bSBarry Smith M*/ 510c838673bSBarry Smith 511c838673bSBarry Smith /*MC 512c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 513c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 514c838673bSBarry Smith 515c838673bSBarry Smith Level: beginner 516c838673bSBarry Smith 517c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 518c838673bSBarry Smith 519c838673bSBarry Smith M*/ 520c838673bSBarry Smith 521c838673bSBarry Smith /*MC 522c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 523c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 524c838673bSBarry Smith be positive definite 525c838673bSBarry Smith 526c838673bSBarry Smith Level: beginner 527c838673bSBarry Smith 5282401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 529c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 530c838673bSBarry Smith 531c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 532c838673bSBarry Smith 533c838673bSBarry Smith M*/ 534c838673bSBarry Smith 535c838673bSBarry Smith /*MC 536c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 537c838673bSBarry Smith while the KSPSolve() is still running. 538c838673bSBarry Smith 539c838673bSBarry Smith Level: beginner 540c838673bSBarry Smith 541c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 542c838673bSBarry Smith 543c838673bSBarry Smith M*/ 544c838673bSBarry Smith 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 5478de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 5498de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 5508de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 5518de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5528de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5530059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 555abef13c0SSatish Balay 5568ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5578ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5588ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5598ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5608ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5618ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5628ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5638ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5648ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 5658ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 5668ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 5678ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 5688ea1b3e6SJed Brown 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 570d4fbbf0eSBarry Smith 57128ce4d24SBarry Smith /*E 57228ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 57328ce4d24SBarry Smith 57428ce4d24SBarry Smith Level: beginner 57528ce4d24SBarry Smith 57628ce4d24SBarry Smith .seealso: KSPCGSetType() 57728ce4d24SBarry Smith E*/ 5789dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5796a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 58028ce4d24SBarry Smith 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5838031f4adStmunson 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 587fcae7a14Stmunson 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 591e559a7feSLois Curfman McInnes 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5978031f4adStmunson 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5991d6018f0SLisandro Dalcin 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 6023369ce9aSBarry Smith 6039804daf3SBarry Smith #include <petscdrawtypes.h> 604d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 605d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 606d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 607d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6092f2e5d10SKris Buschelman 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 61203919abeSBarry Smith 613f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 614ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 615f8a50e2bSBarry Smith 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 622f8a50e2bSBarry Smith 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 626f8a50e2bSBarry Smith 627470b340bSDmitry Karpeev /*E 628470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 629470b340bSDmitry Karpeev 630470b340bSDmitry Karpeev Level: intermediate 631470b340bSDmitry Karpeev 632470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 633470b340bSDmitry Karpeev E*/ 634470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 635470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 636470b340bSDmitry Karpeev 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 639d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 640bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 641aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 642bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 643470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 644470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 6455bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 6465a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 647470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 648470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 6493f22127dSBarry Smith 650014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 651014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 652014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 653014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 654014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 655fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 65623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 657fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 65823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 65923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 660014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 661014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 662fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 663fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6646c699258SBarry Smith 66502b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 66602b02e71SToby 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); 6672eac72dbSBarry Smith #endif 668