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" 7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN 82eac72dbSBarry Smith 9dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitializePackage(const char[]); 101dbb0a54SBarry Smith 1128ce4d24SBarry Smith /*S 1228ce4d24SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods 1328ce4d24SBarry Smith 1428ce4d24SBarry Smith Level: beginner 1528ce4d24SBarry Smith 1628ce4d24SBarry Smith Concepts: Krylov methods 1728ce4d24SBarry Smith 1894b7f48cSBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP 1928ce4d24SBarry Smith S*/ 20e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 212eac72dbSBarry Smith 2228ce4d24SBarry Smith /*E 2328ce4d24SBarry Smith KSPType - String with the name of a PETSc Krylov method or the creation function 2428ce4d24SBarry Smith with an optional dynamic library name, for example 2528ce4d24SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 2628ce4d24SBarry Smith 2728ce4d24SBarry Smith Level: beginner 2828ce4d24SBarry Smith 2928ce4d24SBarry Smith .seealso: KSPSetType(), KSP 3028ce4d24SBarry Smith E*/ 31a313700dSBarry Smith #define KSPType char* 3282bf6240SBarry Smith #define KSPRICHARDSON "richardson" 3382bf6240SBarry Smith #define KSPCHEBYCHEV "chebychev" 3482bf6240SBarry Smith #define KSPCG "cg" 35df2a969aSvictorle #define KSPCGNE "cgne" 36fcae7a14Stmunson #define KSPNASH "nash" 3780e17431SMatthew Knepley #define KSPSTCG "stcg" 388031f4adStmunson #define KSPGLTR "gltr" 3982bf6240SBarry Smith #define KSPGMRES "gmres" 409dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 419dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 4282bf6240SBarry Smith #define KSPTCQMR "tcqmr" 4382bf6240SBarry Smith #define KSPBCGS "bcgs" 44e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 45f0037002Svictorle #define KSPBCGSL "bcgsl" 4682bf6240SBarry Smith #define KSPCGS "cgs" 4782bf6240SBarry Smith #define KSPTFQMR "tfqmr" 4882bf6240SBarry Smith #define KSPCR "cr" 4982bf6240SBarry Smith #define KSPLSQR "lsqr" 5082bf6240SBarry Smith #define KSPPREONLY "preonly" 5182bf6240SBarry Smith #define KSPQCG "qcg" 52c9cf9b11SBarry Smith #define KSPBICG "bicg" 53b4ac9ba4SBarry Smith #define KSPMINRES "minres" 5401c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 5521356919SSatish Balay #define KSPLCD "lcd" 561d6018f0SLisandro Dalcin #define KSPPYTHON "python" 579f6b3dcaSBarry Smith #define KSPBROYDEN "broyden" 582eac72dbSBarry Smith 598ba1e511SMatthew Knepley /* Logging support */ 60dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE; 618ba1e511SMatthew Knepley 62dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *); 63a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,const KSPType); 64dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP); 65dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP); 66dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec); 67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec); 68dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP); 692eac72dbSBarry Smith 70b0a32e0cSBarry Smith extern PetscFList KSPList; 71b022a5c1SBarry Smith extern PetscTruth KSPRegisterAllCalled; 72dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]); 73dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void); 74dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 7530de9b25SBarry Smith 7630de9b25SBarry Smith /*MC 7730de9b25SBarry Smith KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 7830de9b25SBarry Smith 7930de9b25SBarry Smith Synopsis: 80d360dc6fSBarry Smith PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP)) 8130de9b25SBarry Smith 8230de9b25SBarry Smith Not Collective 8330de9b25SBarry Smith 8430de9b25SBarry Smith Input Parameters: 8530de9b25SBarry Smith + name_solver - name of a new user-defined solver 8630de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 8730de9b25SBarry Smith . name_create - name of routine to create method context 8830de9b25SBarry Smith - routine_create - routine to create method context 8930de9b25SBarry Smith 9030de9b25SBarry Smith Notes: 9130de9b25SBarry Smith KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 9230de9b25SBarry Smith 9330de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 9430de9b25SBarry Smith is ignored. 9530de9b25SBarry Smith 9630de9b25SBarry Smith Sample usage: 9730de9b25SBarry Smith .vb 9830de9b25SBarry Smith KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 9930de9b25SBarry Smith "MySolverCreate",MySolverCreate); 10030de9b25SBarry Smith .ve 10130de9b25SBarry Smith 10230de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 10330de9b25SBarry Smith $ KSPSetType(ksp,"my_solver") 10430de9b25SBarry Smith or at runtime via the option 10530de9b25SBarry Smith $ -ksp_type my_solver 10630de9b25SBarry Smith 10730de9b25SBarry Smith Level: advanced 10830de9b25SBarry Smith 109ab901514SBarry Smith Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 11030de9b25SBarry Smith and others of the form ${any_environmental_variable} occuring in pathname will be 11130de9b25SBarry Smith replaced with appropriate values. 11230de9b25SBarry Smith If your function is not being put into a shared library then use KSPRegister() instead 11330de9b25SBarry Smith 11430de9b25SBarry Smith .keywords: KSP, register 11530de9b25SBarry Smith 11630de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy() 11730de9b25SBarry Smith 11830de9b25SBarry Smith M*/ 119aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 120f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 1216df38c32SLois Curfman McInnes #else 122f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 1236df38c32SLois Curfman McInnes #endif 12482bf6240SBarry Smith 125a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,const KSPType *); 126dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide); 127dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*); 128dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 129dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 130dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth); 131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *); 132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth); 133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*); 134a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*); 135dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth); 136a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*); 137dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth); 138dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *); 139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *); 140dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*); 141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*); 142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace); 143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*); 144906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1452eac72dbSBarry Smith 146dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC); 147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*); 148aabeff55SBarry Smith 149a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 150a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP); 151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **); 152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth); 1544b0e389bSBarry Smith 1550e33f6ddSBarry Smith /* not sure where to put this */ 156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*); 157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 1602eac72dbSBarry Smith 1610a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *); 1620a71e23eSMatthew Knepley 163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *); 164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *); 1652eac72dbSBarry Smith 166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal); 167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1714b0e389bSBarry Smith 172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt); 173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal); 1749f236934SBarry Smith 175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP); 176dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 178dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1791d73ed98SBarry Smith 180dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt); 181dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP); 1821d73ed98SBarry Smith 183b49fd9e1SBarry Smith /*E 184b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 185b49fd9e1SBarry Smith 186b49fd9e1SBarry Smith Level: advanced 187b49fd9e1SBarry Smith 188b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 1898c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 190b49fd9e1SBarry Smith 191b49fd9e1SBarry Smith E*/ 19278d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 1939dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[]; 1941f7e983dSSatish Balay /*MC 1958c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 1968c5b8ba0SBarry Smith 1978c5b8ba0SBarry Smith Level: advanced 1988c5b8ba0SBarry Smith 1998c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2008c5b8ba0SBarry Smith 2018c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 2028c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2038c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2048c5b8ba0SBarry Smith M*/ 2058c5b8ba0SBarry Smith 2061f7e983dSSatish Balay /*MC 2078c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2088c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2098c5b8ba0SBarry Smith poor orthogonality. 2108c5b8ba0SBarry Smith 2118c5b8ba0SBarry Smith Level: advanced 2128c5b8ba0SBarry Smith 2138c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2148c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2158c5b8ba0SBarry Smith 2168c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 2178c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2188c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2198c5b8ba0SBarry Smith M*/ 2208c5b8ba0SBarry Smith 2211f7e983dSSatish Balay /*MC 2228c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2238c5b8ba0SBarry Smith 2248c5b8ba0SBarry Smith Level: advanced 2258c5b8ba0SBarry Smith 2268c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2278c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2288c5b8ba0SBarry Smith 2298c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2308c5b8ba0SBarry Smith 2318c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 2328c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2338c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2348c5b8ba0SBarry Smith M*/ 2358c5b8ba0SBarry Smith 236dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 23708480c60SBarry Smith 238dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 239dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 240dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 241c38d4ed2SBarry Smith 242dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal); 243dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*); 244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*); 245121fd945SKris Buschelman 246d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal); 247d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth); 248d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int); 249d9492815SBarry Smith 250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP); 251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2522eac72dbSBarry Smith 253a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 254a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 255*90984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *); 2561a2f9f99SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 2577517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 258a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 259a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 260a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 26184cb2905SBarry Smith 262dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec); 263dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*); 264dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 265dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 266c01c455dSBarry Smith 267dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure); 268dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 269906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*); 270dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]); 271dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]); 2722dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]); 2731eb62cbbSBarry Smith 274dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth); 275dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*); 276dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth); 277dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*); 2781f7f0c4fSBarry Smith 279dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer); 2801eb62cbbSBarry Smith 281db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec); 282db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*); 283db9b2ab1SHong Zhang 28428ce4d24SBarry Smith /*E 2858a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 2868a4b9c5cSBarry Smith test routines. 2878a4b9c5cSBarry Smith 2888a4b9c5cSBarry Smith Level: advanced 2898a4b9c5cSBarry Smith 290a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 291a3f661c8SBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 292a3f661c8SBarry Smith 2938a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 2948a4b9c5cSBarry Smith 29594b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 296a3f661c8SBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 2978a4b9c5cSBarry Smith E*/ 298ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 2999dcbbd2bSBarry Smith extern const char *KSPNormTypes[]; 3001f7e983dSSatish Balay /*MC 301ce9499c7SBarry Smith KSP_NORM_NO - Do not compute a norm during the Krylov process. This will 3028c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3038c5b8ba0SBarry Smith be based on a norm of a residual etc. 3048c5b8ba0SBarry Smith 3058c5b8ba0SBarry Smith Level: advanced 3068c5b8ba0SBarry Smith 307085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3088c5b8ba0SBarry Smith 309ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3108c5b8ba0SBarry Smith M*/ 3118c5b8ba0SBarry Smith 3121f7e983dSSatish Balay /*MC 313ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3148c5b8ba0SBarry Smith convergence test routine. 3158c5b8ba0SBarry Smith 3168c5b8ba0SBarry Smith Level: advanced 3178c5b8ba0SBarry Smith 318ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3198c5b8ba0SBarry Smith M*/ 3208c5b8ba0SBarry Smith 3211f7e983dSSatish Balay /*MC 322ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3238c5b8ba0SBarry Smith convergence test routine. 3248c5b8ba0SBarry Smith 3258c5b8ba0SBarry Smith Level: advanced 3268c5b8ba0SBarry Smith 327ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3288c5b8ba0SBarry Smith M*/ 3298c5b8ba0SBarry Smith 3301f7e983dSSatish Balay /*MC 331ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 332a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3338c5b8ba0SBarry Smith 3348c5b8ba0SBarry Smith Level: advanced 3358c5b8ba0SBarry Smith 336ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3378c5b8ba0SBarry Smith M*/ 3388c5b8ba0SBarry Smith 339dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType); 340863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*); 3413c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt); 3421e7308c7SRichard Tran Mills EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetLagNorm(KSP,PetscTruth); 3438a4b9c5cSBarry Smith 3448a4b9c5cSBarry Smith /*E 34528ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 34628ce4d24SBarry Smith have converged or diverged 34728ce4d24SBarry Smith 34828ce4d24SBarry Smith Level: beginner 34928ce4d24SBarry Smith 3504d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 35128ce4d24SBarry Smith 3524d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3534d0a8057SBarry Smith 3544d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3554d0a8057SBarry Smith any of the values here also change them that array of names. 35686c02ca4SBarry Smith 357c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 35828ce4d24SBarry Smith E*/ 359d15094e1SBarry Smith typedef enum {/* converged */ 360d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 361d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 362b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3638031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3648031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 365329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 366af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 367d15094e1SBarry Smith /* diverged */ 368b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 369d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 370d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 371d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 372b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 373b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 374b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 3756aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 3766aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 377d15094e1SBarry Smith 378d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 3799dcbbd2bSBarry Smith extern const char **KSPConvergedReasons; 380d15094e1SBarry Smith 381c838673bSBarry Smith /*MC 382c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 383c838673bSBarry Smith 384c838673bSBarry Smith Level: beginner 385c838673bSBarry Smith 386c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 387c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 388c838673bSBarry Smith 2-norm of the residual for right preconditioning 389c838673bSBarry Smith 390c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 391c838673bSBarry Smith 392c838673bSBarry Smith M*/ 393c838673bSBarry Smith 394c838673bSBarry Smith /*MC 395c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 396c838673bSBarry Smith 397c838673bSBarry Smith Level: beginner 398c838673bSBarry Smith 399c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 400c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 401c838673bSBarry Smith 2-norm of the residual for right preconditioning 402c838673bSBarry Smith 403c838673bSBarry Smith Level: beginner 404c838673bSBarry Smith 405c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 406c838673bSBarry Smith 407c838673bSBarry Smith M*/ 408c838673bSBarry Smith 409c838673bSBarry Smith /*MC 410c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 411c838673bSBarry Smith 412c838673bSBarry Smith Level: beginner 413c838673bSBarry Smith 414c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 415c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 416c838673bSBarry Smith 2-norm of the residual for right preconditioning 417c838673bSBarry Smith 418c838673bSBarry Smith Level: beginner 419c838673bSBarry Smith 420c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 421c838673bSBarry Smith 422c838673bSBarry Smith M*/ 423c838673bSBarry Smith 424c838673bSBarry Smith /*MC 425c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 426c838673bSBarry Smith reached 427c838673bSBarry Smith 428c838673bSBarry Smith Level: beginner 429c838673bSBarry Smith 430c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 431c838673bSBarry Smith 432c838673bSBarry Smith M*/ 433c838673bSBarry Smith 434c838673bSBarry Smith /*MC 4358c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4368c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4374d0a8057SBarry Smith test routine is set in KSP. 438c838673bSBarry Smith 439c838673bSBarry Smith 440c838673bSBarry Smith Level: beginner 441c838673bSBarry Smith 442c838673bSBarry Smith 443c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 444c838673bSBarry Smith 445c838673bSBarry Smith M*/ 446c838673bSBarry Smith 447c838673bSBarry Smith /*MC 448c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 449c838673bSBarry Smith method could not continue to enlarge the Krylov space. 450c838673bSBarry Smith 451c838673bSBarry Smith Level: beginner 452c838673bSBarry Smith 453c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 454c838673bSBarry Smith 455c838673bSBarry Smith M*/ 456c838673bSBarry Smith 457c838673bSBarry Smith /*MC 458c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 459c838673bSBarry Smith method could not continue to enlarge the Krylov space. 460c838673bSBarry Smith 461c838673bSBarry Smith 462c838673bSBarry Smith Level: beginner 463c838673bSBarry Smith 464c838673bSBarry Smith 465c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 466c838673bSBarry Smith 467c838673bSBarry Smith M*/ 468c838673bSBarry Smith 469c838673bSBarry Smith /*MC 470c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 471c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 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 480c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 481c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 482c838673bSBarry Smith be positive definite 483c838673bSBarry Smith 484c838673bSBarry Smith Level: beginner 485c838673bSBarry Smith 4862401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 487c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 488c838673bSBarry Smith 489c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 490c838673bSBarry Smith 491c838673bSBarry Smith M*/ 492c838673bSBarry Smith 493c838673bSBarry Smith /*MC 494c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 495c838673bSBarry Smith while the KSPSolve() is still running. 496c838673bSBarry Smith 497c838673bSBarry Smith Level: beginner 498c838673bSBarry Smith 499c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 500c838673bSBarry Smith 501c838673bSBarry Smith M*/ 502c838673bSBarry Smith 503971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 504dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **); 505dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 506*90984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 507971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *); 508971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **); 50994b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP); 5105f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP); 511dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 512dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *); 513abef13c0SSatish Balay 514dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *); 515d4fbbf0eSBarry Smith 51628ce4d24SBarry Smith /*E 51728ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 51828ce4d24SBarry Smith 51928ce4d24SBarry Smith Level: beginner 52028ce4d24SBarry Smith 52128ce4d24SBarry Smith .seealso: KSPCGSetType() 52228ce4d24SBarry Smith E*/ 5239dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5249dcbbd2bSBarry Smith extern const char *KSPCGTypes[]; 52528ce4d24SBarry Smith 526dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType); 527608df7e2SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGUseSingleReduction(KSP,PetscTruth); 5288031f4adStmunson 529fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal); 530fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *); 531fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *); 532fcae7a14Stmunson 5330d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal); 53411f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *); 535d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *); 536e559a7feSLois Curfman McInnes 5378031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal); 5388031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *); 539d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *); 540d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *); 5413c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *); 5428031f4adStmunson 5431d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPPythonSetType(KSP,const char[]); 5441d6018f0SLisandro Dalcin 545dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP); 546dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP); 5473369ce9aSBarry Smith 548a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 549a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 550a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG); 5517517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 5527517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 5537517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG); 554b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 555b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 556b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeDestroy(PetscDrawLG); 5572f2e5d10SKris Buschelman 5586891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 5596891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 56003919abeSBarry Smith 561e884886fSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatMFFDKSPMonitor(KSP,PetscInt,PetscReal,void *); 562e884886fSBarry Smith 563f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 564436851aaSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscTruth monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 565f8a50e2bSBarry Smith 566f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 567f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess); 568f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess); 569f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec); 570f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 571436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess); 572f8a50e2bSBarry Smith 573f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 574f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess); 575f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*); 576f8a50e2bSBarry Smith 577084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 5783f22127dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetKSP(Mat,KSP*); 579084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 580084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 581dfe085dbSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 5823f22127dSBarry Smith 583e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 5842eac72dbSBarry Smith #endif 585