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 60a835dfdSSatish Balay #include "petscpc.h" 72eac72dbSBarry Smith 8*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitializePackage(const char[]); 91dbb0a54SBarry Smith 1028ce4d24SBarry Smith /*S 1128ce4d24SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods 1228ce4d24SBarry Smith 1328ce4d24SBarry Smith Level: beginner 1428ce4d24SBarry Smith 1528ce4d24SBarry Smith Concepts: Krylov methods 1628ce4d24SBarry Smith 1794b7f48cSBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP 1828ce4d24SBarry Smith S*/ 19e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 202eac72dbSBarry Smith 2176bdecfbSBarry Smith /*J 2228ce4d24SBarry Smith KSPType - String with the name of a PETSc Krylov method or the creation function 2328ce4d24SBarry Smith with an optional dynamic library name, for example 2428ce4d24SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 2528ce4d24SBarry Smith 2628ce4d24SBarry Smith Level: beginner 2728ce4d24SBarry Smith 2828ce4d24SBarry Smith .seealso: KSPSetType(), KSP 2976bdecfbSBarry Smith J*/ 30a313700dSBarry Smith #define KSPType char* 3182bf6240SBarry Smith #define KSPRICHARDSON "richardson" 326c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3382bf6240SBarry Smith #define KSPCG "cg" 34df2a969aSvictorle #define KSPCGNE "cgne" 35fcae7a14Stmunson #define KSPNASH "nash" 3680e17431SMatthew Knepley #define KSPSTCG "stcg" 378031f4adStmunson #define KSPGLTR "gltr" 3882bf6240SBarry Smith #define KSPGMRES "gmres" 399dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 409dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 414f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 4261b468f9SJed Brown #define KSPPGMRES "pgmres" 4382bf6240SBarry Smith #define KSPTCQMR "tcqmr" 4482bf6240SBarry Smith #define KSPBCGS "bcgs" 45e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 4618ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 4718ac38e6SHong Zhang #define KSPIFBCGS "ifbcgs" 48f0037002Svictorle #define KSPBCGSL "bcgsl" 4982bf6240SBarry Smith #define KSPCGS "cgs" 5082bf6240SBarry Smith #define KSPTFQMR "tfqmr" 5182bf6240SBarry Smith #define KSPCR "cr" 5282bf6240SBarry Smith #define KSPLSQR "lsqr" 5382bf6240SBarry Smith #define KSPPREONLY "preonly" 5482bf6240SBarry Smith #define KSPQCG "qcg" 55c9cf9b11SBarry Smith #define KSPBICG "bicg" 56b4ac9ba4SBarry Smith #define KSPMINRES "minres" 5701c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 5821356919SSatish Balay #define KSPLCD "lcd" 591d6018f0SLisandro Dalcin #define KSPPYTHON "python" 6058ad63f4SBarry Smith #define KSPGCR "gcr" 61bbce1389SJed Brown #define KSPSPECEST "specest" 622eac72dbSBarry Smith 638ba1e511SMatthew Knepley /* Logging support */ 64*014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 658ba1e511SMatthew Knepley 66*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 67*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetType(KSP,const KSPType); 68*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 69*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 70*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 71*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 72*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 73*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 742eac72dbSBarry Smith 75*014dd563SJed Brown PETSC_EXTERN PetscFList KSPList; 76*014dd563SJed Brown PETSC_EXTERN PetscBool KSPRegisterAllCalled; 77*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterAll(const char[]); 78*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterDestroy(void); 79*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 8030de9b25SBarry Smith 8130de9b25SBarry Smith /*MC 8230de9b25SBarry Smith KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 8330de9b25SBarry Smith 8430de9b25SBarry Smith Synopsis: 851890ba74SBarry Smith PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP)) 8630de9b25SBarry Smith 8730de9b25SBarry Smith Not Collective 8830de9b25SBarry Smith 8930de9b25SBarry Smith Input Parameters: 9030de9b25SBarry Smith + name_solver - name of a new user-defined solver 9130de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 9230de9b25SBarry Smith . name_create - name of routine to create method context 9330de9b25SBarry Smith - routine_create - routine to create method context 9430de9b25SBarry Smith 9530de9b25SBarry Smith Notes: 9630de9b25SBarry Smith KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 9730de9b25SBarry Smith 9830de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 9930de9b25SBarry Smith is ignored. 10030de9b25SBarry Smith 10130de9b25SBarry Smith Sample usage: 10230de9b25SBarry Smith .vb 10330de9b25SBarry Smith KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 10430de9b25SBarry Smith "MySolverCreate",MySolverCreate); 10530de9b25SBarry Smith .ve 10630de9b25SBarry Smith 10730de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 10830de9b25SBarry Smith $ KSPSetType(ksp,"my_solver") 10930de9b25SBarry Smith or at runtime via the option 11030de9b25SBarry Smith $ -ksp_type my_solver 11130de9b25SBarry Smith 11230de9b25SBarry Smith Level: advanced 11330de9b25SBarry Smith 114ab901514SBarry Smith Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 11530de9b25SBarry Smith and others of the form ${any_environmental_variable} occuring in pathname will be 11630de9b25SBarry Smith replaced with appropriate values. 11730de9b25SBarry Smith If your function is not being put into a shared library then use KSPRegister() instead 11830de9b25SBarry Smith 11930de9b25SBarry Smith .keywords: KSP, register 12030de9b25SBarry Smith 12130de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy() 12230de9b25SBarry Smith 12330de9b25SBarry Smith M*/ 124aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 125f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 1266df38c32SLois Curfman McInnes #else 127f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 1286df38c32SLois Curfman McInnes #endif 12982bf6240SBarry Smith 130*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetType(KSP,const KSPType *); 131*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 132*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 133*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 134*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 135*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 136*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 137*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 138*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 139*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 140*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 141*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 142*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 143*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 144*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 145*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 146*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 147*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 148*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 149*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 150*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 151*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1522eac72dbSBarry Smith 153*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 154*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 155aabeff55SBarry Smith 156*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 157*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 158*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 159*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 160*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 161*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1624b0e389bSBarry Smith 1630e33f6ddSBarry Smith /* not sure where to put this */ 164*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 165*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 166*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 167*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 168*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 1692eac72dbSBarry Smith 170*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1710a71e23eSMatthew Knepley 172*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 173*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1742eac72dbSBarry Smith 175*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 176*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 177*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 178*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 179*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 180*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 181*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 182*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1834b0e389bSBarry Smith 184*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 185*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 186*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1879f236934SBarry Smith 188*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 189*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 190*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 191*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 192*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1931d73ed98SBarry Smith 194*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 195*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 1961d73ed98SBarry Smith 197*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 198*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 199*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 20058ad63f4SBarry Smith 201b49fd9e1SBarry Smith /*E 202b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 203b49fd9e1SBarry Smith 204b49fd9e1SBarry Smith Level: advanced 205b49fd9e1SBarry Smith 206bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 207e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 208b49fd9e1SBarry Smith 209b49fd9e1SBarry Smith E*/ 21078d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 211*014dd563SJed Brown PETSC_EXTERN const char *KSPGMRESCGSRefinementTypes[]; 2121f7e983dSSatish Balay /*MC 2138c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2148c5b8ba0SBarry Smith 2158c5b8ba0SBarry Smith Level: advanced 2168c5b8ba0SBarry Smith 2178c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2188c5b8ba0SBarry Smith 219bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 220e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2218c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2228c5b8ba0SBarry Smith M*/ 2238c5b8ba0SBarry Smith 2241f7e983dSSatish Balay /*MC 2258c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2268c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2278c5b8ba0SBarry Smith poor orthogonality. 2288c5b8ba0SBarry Smith 2298c5b8ba0SBarry Smith Level: advanced 2308c5b8ba0SBarry Smith 2318c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2328c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2338c5b8ba0SBarry Smith 234bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 235e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2368c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2378c5b8ba0SBarry Smith M*/ 2388c5b8ba0SBarry Smith 2391f7e983dSSatish Balay /*MC 2408c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2418c5b8ba0SBarry Smith 2428c5b8ba0SBarry Smith Level: advanced 2438c5b8ba0SBarry Smith 2448c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2458c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2468c5b8ba0SBarry Smith 2478c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2488c5b8ba0SBarry Smith 249bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 250e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2518c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2528c5b8ba0SBarry Smith M*/ 2538c5b8ba0SBarry Smith 254*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 255*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 25608480c60SBarry Smith 257*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 258*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 259*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 260c38d4ed2SBarry Smith 261*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 262*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 263*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 264121fd945SKris Buschelman 265*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 266*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 267*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 268d9492815SBarry Smith 269*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 270*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2712eac72dbSBarry Smith 272*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 273*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 274*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 275*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 276*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 277*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 278*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 279*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*); 280*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**); 281*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**); 282*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 28384cb2905SBarry Smith 284*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 285*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 286*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 287*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 288c01c455dSBarry Smith 289*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 290*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 291*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 292*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 293*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 294*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 2951eb62cbbSBarry Smith 296*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 297*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 298*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 299*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3001f7f0c4fSBarry Smith 301*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 3021eb62cbbSBarry Smith 303*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 304*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 305db9b2ab1SHong Zhang 306*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 307*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 30883ab6a24SBarry Smith 30928ce4d24SBarry Smith /*E 3108a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3118a4b9c5cSBarry Smith test routines. 3128a4b9c5cSBarry Smith 3138a4b9c5cSBarry Smith Level: advanced 3148a4b9c5cSBarry Smith 315a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3161718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 317a3f661c8SBarry Smith 3188a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3198a4b9c5cSBarry Smith 32094b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3211718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3228a4b9c5cSBarry Smith E*/ 3239e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3249e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 325*014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 326a21b2a99SBarry Smith 3271f7e983dSSatish Balay /*MC 3289793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3298c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3308c5b8ba0SBarry Smith be based on a norm of a residual etc. 3318c5b8ba0SBarry Smith 3328c5b8ba0SBarry Smith Level: advanced 3338c5b8ba0SBarry Smith 334085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3358c5b8ba0SBarry Smith 336ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3378c5b8ba0SBarry Smith M*/ 3388c5b8ba0SBarry Smith 3391f7e983dSSatish Balay /*MC 340ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3418c5b8ba0SBarry Smith convergence test routine. 3428c5b8ba0SBarry Smith 3438c5b8ba0SBarry Smith Level: advanced 3448c5b8ba0SBarry Smith 3459793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3468c5b8ba0SBarry Smith M*/ 3478c5b8ba0SBarry Smith 3481f7e983dSSatish Balay /*MC 349ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3508c5b8ba0SBarry Smith convergence test routine. 3518c5b8ba0SBarry Smith 3528c5b8ba0SBarry Smith Level: advanced 3538c5b8ba0SBarry Smith 3549793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3558c5b8ba0SBarry Smith M*/ 3568c5b8ba0SBarry Smith 3571f7e983dSSatish Balay /*MC 358ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 359a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3608c5b8ba0SBarry Smith 3618c5b8ba0SBarry Smith Level: advanced 3628c5b8ba0SBarry Smith 3639793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3648c5b8ba0SBarry Smith M*/ 3658c5b8ba0SBarry Smith 366*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 367*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 368*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 369*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 370*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 3718a4b9c5cSBarry Smith 3728a4b9c5cSBarry Smith /*E 37328ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 37428ce4d24SBarry Smith have converged or diverged 37528ce4d24SBarry Smith 37628ce4d24SBarry Smith Level: beginner 37728ce4d24SBarry Smith 3784d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 37928ce4d24SBarry Smith 3804d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3814d0a8057SBarry Smith 3824d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3834d0a8057SBarry Smith any of the values here also change them that array of names. 38486c02ca4SBarry Smith 385c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 38628ce4d24SBarry Smith E*/ 387d15094e1SBarry Smith typedef enum {/* converged */ 3889ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 3899ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 390d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 391d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 392b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3938031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3948031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 395329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 396af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 397d15094e1SBarry Smith /* diverged */ 398b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 399d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 400d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 401d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 402b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 403b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 404b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4056aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 4066aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 407d15094e1SBarry Smith 408d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 409*014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 410d15094e1SBarry Smith 411c838673bSBarry Smith /*MC 412c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 413c838673bSBarry Smith 414c838673bSBarry Smith Level: beginner 415c838673bSBarry Smith 416c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 417c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 418c838673bSBarry Smith 2-norm of the residual for right preconditioning 419c838673bSBarry Smith 420c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 421c838673bSBarry Smith 422c838673bSBarry Smith M*/ 423c838673bSBarry Smith 424c838673bSBarry Smith /*MC 425c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 426c838673bSBarry Smith 427c838673bSBarry Smith Level: beginner 428c838673bSBarry Smith 429c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 430c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 431c838673bSBarry Smith 2-norm of the residual for right preconditioning 432c838673bSBarry Smith 433c838673bSBarry Smith Level: beginner 434c838673bSBarry Smith 435c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 436c838673bSBarry Smith 437c838673bSBarry Smith M*/ 438c838673bSBarry Smith 439c838673bSBarry Smith /*MC 440c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 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_ITS - Ran out of iterations before any convergence criteria was 456c838673bSBarry Smith reached 457c838673bSBarry Smith 458c838673bSBarry Smith Level: beginner 459c838673bSBarry Smith 460c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 461c838673bSBarry Smith 462c838673bSBarry Smith M*/ 463c838673bSBarry Smith 464c838673bSBarry Smith /*MC 4658c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4668c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4674d0a8057SBarry Smith test routine is set in KSP. 468c838673bSBarry Smith 469c838673bSBarry Smith 470c838673bSBarry Smith Level: beginner 471c838673bSBarry Smith 472c838673bSBarry Smith 473c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 474c838673bSBarry Smith 475c838673bSBarry Smith M*/ 476c838673bSBarry Smith 477c838673bSBarry Smith /*MC 478c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4793014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4803014e516SBarry Smith preconditioner. 481c838673bSBarry Smith 482c838673bSBarry Smith Level: beginner 483c838673bSBarry Smith 484c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 485c838673bSBarry Smith 486c838673bSBarry Smith M*/ 487c838673bSBarry Smith 488c838673bSBarry Smith /*MC 489c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 490c838673bSBarry Smith method could not continue to enlarge the Krylov space. 491c838673bSBarry Smith 492c838673bSBarry Smith 493c838673bSBarry Smith Level: beginner 494c838673bSBarry Smith 495c838673bSBarry Smith 496c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 497c838673bSBarry Smith 498c838673bSBarry Smith M*/ 499c838673bSBarry Smith 500c838673bSBarry Smith /*MC 501c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 502c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 503c838673bSBarry Smith 504c838673bSBarry Smith Level: beginner 505c838673bSBarry Smith 506c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 507c838673bSBarry Smith 508c838673bSBarry Smith M*/ 509c838673bSBarry Smith 510c838673bSBarry Smith /*MC 511c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 512c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 513c838673bSBarry Smith be positive definite 514c838673bSBarry Smith 515c838673bSBarry Smith Level: beginner 516c838673bSBarry Smith 5172401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 518c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 519c838673bSBarry Smith 520c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 521c838673bSBarry Smith 522c838673bSBarry Smith M*/ 523c838673bSBarry Smith 524c838673bSBarry Smith /*MC 525c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 526c838673bSBarry Smith while the KSPSolve() is still running. 527c838673bSBarry Smith 528c838673bSBarry Smith Level: beginner 529c838673bSBarry Smith 530c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 531c838673bSBarry Smith 532c838673bSBarry Smith M*/ 533c838673bSBarry Smith 534*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 535*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 536*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 537*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 538*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *); 539*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **); 540*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 541*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 542*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 543*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 544abef13c0SSatish Balay 545*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 546d4fbbf0eSBarry Smith 54728ce4d24SBarry Smith /*E 54828ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 54928ce4d24SBarry Smith 55028ce4d24SBarry Smith Level: beginner 55128ce4d24SBarry Smith 55228ce4d24SBarry Smith .seealso: KSPCGSetType() 55328ce4d24SBarry Smith E*/ 5549dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 555*014dd563SJed Brown PETSC_EXTERN const char *KSPCGTypes[]; 55628ce4d24SBarry Smith 557*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 558*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5598031f4adStmunson 560*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 561*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 562*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 563fcae7a14Stmunson 564*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 565*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 566*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 567e559a7feSLois Curfman McInnes 568*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 569*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 570*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 571*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 572*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5738031f4adStmunson 574*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5751d6018f0SLisandro Dalcin 576*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 577*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5783369ce9aSBarry Smith 579*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 580*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 581*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGDestroy(PetscDrawLG*); 582*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 583*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 584*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 585*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 586*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 587*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRangeDestroy(PetscDrawLG*); 5882f2e5d10SKris Buschelman 589*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 590*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 59103919abeSBarry Smith 592f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 593ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 594f8a50e2bSBarry Smith 595*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 596*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 597*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 598*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 599*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 600*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 601f8a50e2bSBarry Smith 602*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 603*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 604*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 605f8a50e2bSBarry Smith 606*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 607*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 608*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 609*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 610*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 6113f22127dSBarry Smith 612*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat); 613fcfd50ebSBarry Smith 614*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 615*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 616*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 617*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 618*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 619*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 620*014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 621*014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 622*014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 623*014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 624*014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6256c699258SBarry Smith 6262eac72dbSBarry Smith #endif 627