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" 5858ad63f4SBarry Smith #define KSPGCR "gcr" 59*7945016fSBarry Smith #define KSPNGMRES "ngmres" 602eac72dbSBarry Smith 618ba1e511SMatthew Knepley /* Logging support */ 620700a824SBarry Smith extern PetscClassId PETSCKSP_DLLEXPORT KSP_CLASSID; 638ba1e511SMatthew Knepley 64dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *); 65a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,const KSPType); 66dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP); 67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP); 68dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec); 69dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec); 70dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP); 712eac72dbSBarry Smith 72b0a32e0cSBarry Smith extern PetscFList KSPList; 73b022a5c1SBarry Smith extern PetscTruth KSPRegisterAllCalled; 74dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]); 75dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void); 76dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 7730de9b25SBarry Smith 7830de9b25SBarry Smith /*MC 7930de9b25SBarry Smith KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 8030de9b25SBarry Smith 8130de9b25SBarry Smith Synopsis: 821890ba74SBarry Smith PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP)) 8330de9b25SBarry Smith 8430de9b25SBarry Smith Not Collective 8530de9b25SBarry Smith 8630de9b25SBarry Smith Input Parameters: 8730de9b25SBarry Smith + name_solver - name of a new user-defined solver 8830de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 8930de9b25SBarry Smith . name_create - name of routine to create method context 9030de9b25SBarry Smith - routine_create - routine to create method context 9130de9b25SBarry Smith 9230de9b25SBarry Smith Notes: 9330de9b25SBarry Smith KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 9430de9b25SBarry Smith 9530de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 9630de9b25SBarry Smith is ignored. 9730de9b25SBarry Smith 9830de9b25SBarry Smith Sample usage: 9930de9b25SBarry Smith .vb 10030de9b25SBarry Smith KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 10130de9b25SBarry Smith "MySolverCreate",MySolverCreate); 10230de9b25SBarry Smith .ve 10330de9b25SBarry Smith 10430de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 10530de9b25SBarry Smith $ KSPSetType(ksp,"my_solver") 10630de9b25SBarry Smith or at runtime via the option 10730de9b25SBarry Smith $ -ksp_type my_solver 10830de9b25SBarry Smith 10930de9b25SBarry Smith Level: advanced 11030de9b25SBarry Smith 111ab901514SBarry Smith Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 11230de9b25SBarry Smith and others of the form ${any_environmental_variable} occuring in pathname will be 11330de9b25SBarry Smith replaced with appropriate values. 11430de9b25SBarry Smith If your function is not being put into a shared library then use KSPRegister() instead 11530de9b25SBarry Smith 11630de9b25SBarry Smith .keywords: KSP, register 11730de9b25SBarry Smith 11830de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy() 11930de9b25SBarry Smith 12030de9b25SBarry Smith M*/ 121aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 122f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 1236df38c32SLois Curfman McInnes #else 124f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 1256df38c32SLois Curfman McInnes #endif 12682bf6240SBarry Smith 127a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,const KSPType *); 1281718446fSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPCSide(KSP,PCSide); 129b037da10SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPCSide(KSP,PCSide*); 130dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth); 133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *); 134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth); 135dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*); 136a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*); 137dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth); 138a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*); 139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth); 140dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *); 141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *); 142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*); 143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*); 144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace); 145dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*); 146906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1472eac72dbSBarry Smith 148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC); 149dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*); 150aabeff55SBarry Smith 151a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 152a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP); 153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **); 154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 155dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth); 1564b0e389bSBarry Smith 1570e33f6ddSBarry Smith /* not sure where to put this */ 158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*); 159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 1622eac72dbSBarry Smith 1630a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *); 1640a71e23eSMatthew Knepley 165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *); 166dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *); 1672eac72dbSBarry Smith 168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal); 169b61692a6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetSelfScale(KSP,PetscTruth); 170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1744b0e389bSBarry Smith 175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt); 176dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal); 1779f236934SBarry Smith 178dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP); 179dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 180dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 181dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1821d73ed98SBarry Smith 183dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt); 184dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP); 1851d73ed98SBarry Smith 18658ad63f4SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetRestart(KSP,PetscInt); 1877ae9e5e1SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 18858ad63f4SBarry Smith 189b49fd9e1SBarry Smith /*E 190b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 191b49fd9e1SBarry Smith 192b49fd9e1SBarry Smith Level: advanced 193b49fd9e1SBarry Smith 194b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 1958c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 196b49fd9e1SBarry Smith 197b49fd9e1SBarry Smith E*/ 19878d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 1999dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[]; 2001f7e983dSSatish Balay /*MC 2018c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2028c5b8ba0SBarry Smith 2038c5b8ba0SBarry Smith Level: advanced 2048c5b8ba0SBarry Smith 2058c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2068c5b8ba0SBarry Smith 2078c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 2088c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2098c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2108c5b8ba0SBarry Smith M*/ 2118c5b8ba0SBarry Smith 2121f7e983dSSatish Balay /*MC 2138c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2148c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2158c5b8ba0SBarry Smith poor orthogonality. 2168c5b8ba0SBarry Smith 2178c5b8ba0SBarry Smith Level: advanced 2188c5b8ba0SBarry Smith 2198c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2208c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2218c5b8ba0SBarry Smith 2228c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 2238c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2248c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2258c5b8ba0SBarry Smith M*/ 2268c5b8ba0SBarry Smith 2271f7e983dSSatish Balay /*MC 2288c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2298c5b8ba0SBarry Smith 2308c5b8ba0SBarry Smith Level: advanced 2318c5b8ba0SBarry Smith 2328c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2338c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2348c5b8ba0SBarry Smith 2358c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2368c5b8ba0SBarry Smith 2378c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 2388c5b8ba0SBarry Smith KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2398c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2408c5b8ba0SBarry Smith M*/ 2418c5b8ba0SBarry Smith 242dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 24308480c60SBarry Smith 244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 245dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 246dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 247c38d4ed2SBarry Smith 248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal); 249dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*); 250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*); 251121fd945SKris Buschelman 252d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal); 253d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth); 254d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int); 255d9492815SBarry Smith 256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP); 257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2582eac72dbSBarry Smith 259a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 260a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 26190984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *); 2621a2f9f99SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 2637517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 264a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 265a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 266a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 26784cb2905SBarry Smith 268dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec); 269dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*); 270dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 271dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 272c01c455dSBarry Smith 273dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure); 274dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 275906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*); 276dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]); 277dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]); 2782dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]); 2791eb62cbbSBarry Smith 280dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth); 281dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*); 282dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth); 283dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*); 2841f7f0c4fSBarry Smith 285dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer); 2861eb62cbbSBarry Smith 287db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec); 288db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*); 289db9b2ab1SHong Zhang 29028ce4d24SBarry Smith /*E 2918a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 2928a4b9c5cSBarry Smith test routines. 2938a4b9c5cSBarry Smith 2948a4b9c5cSBarry Smith Level: advanced 2958a4b9c5cSBarry Smith 296a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 2971718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 298a3f661c8SBarry Smith 2998a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3008a4b9c5cSBarry Smith 30194b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3021718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3038a4b9c5cSBarry Smith E*/ 304ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3059dcbbd2bSBarry Smith extern const char *KSPNormTypes[]; 3061f7e983dSSatish Balay /*MC 307ce9499c7SBarry Smith KSP_NORM_NO - Do not compute a norm during the Krylov process. This will 3088c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3098c5b8ba0SBarry Smith be based on a norm of a residual etc. 3108c5b8ba0SBarry Smith 3118c5b8ba0SBarry Smith Level: advanced 3128c5b8ba0SBarry Smith 313085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3148c5b8ba0SBarry Smith 315ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3168c5b8ba0SBarry Smith M*/ 3178c5b8ba0SBarry Smith 3181f7e983dSSatish Balay /*MC 319ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3208c5b8ba0SBarry Smith convergence test routine. 3218c5b8ba0SBarry Smith 3228c5b8ba0SBarry Smith Level: advanced 3238c5b8ba0SBarry Smith 324ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3258c5b8ba0SBarry Smith M*/ 3268c5b8ba0SBarry Smith 3271f7e983dSSatish Balay /*MC 328ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3298c5b8ba0SBarry Smith convergence test routine. 3308c5b8ba0SBarry Smith 3318c5b8ba0SBarry Smith Level: advanced 3328c5b8ba0SBarry Smith 333ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3348c5b8ba0SBarry Smith M*/ 3358c5b8ba0SBarry Smith 3361f7e983dSSatish Balay /*MC 337ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 338a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3398c5b8ba0SBarry Smith 3408c5b8ba0SBarry Smith Level: advanced 3418c5b8ba0SBarry Smith 342ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3438c5b8ba0SBarry Smith M*/ 3448c5b8ba0SBarry Smith 345dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType); 346863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*); 3473c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt); 3481e7308c7SRichard Tran Mills EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetLagNorm(KSP,PetscTruth); 3498a4b9c5cSBarry Smith 3508a4b9c5cSBarry Smith /*E 35128ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 35228ce4d24SBarry Smith have converged or diverged 35328ce4d24SBarry Smith 35428ce4d24SBarry Smith Level: beginner 35528ce4d24SBarry Smith 3564d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 35728ce4d24SBarry Smith 3584d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3594d0a8057SBarry Smith 3604d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3614d0a8057SBarry Smith any of the values here also change them that array of names. 36286c02ca4SBarry Smith 363c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 36428ce4d24SBarry Smith E*/ 365d15094e1SBarry Smith typedef enum {/* converged */ 366d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 367d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 368b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3698031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3708031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 371329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 372af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 373d15094e1SBarry Smith /* diverged */ 374b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 375d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 376d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 377d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 378b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 379b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 380b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 3816aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 3826aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 383d15094e1SBarry Smith 384d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 3859dcbbd2bSBarry Smith extern const char **KSPConvergedReasons; 386d15094e1SBarry Smith 387c838673bSBarry Smith /*MC 388c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 389c838673bSBarry Smith 390c838673bSBarry Smith Level: beginner 391c838673bSBarry Smith 392c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 393c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 394c838673bSBarry Smith 2-norm of the residual for right preconditioning 395c838673bSBarry Smith 396c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 397c838673bSBarry Smith 398c838673bSBarry Smith M*/ 399c838673bSBarry Smith 400c838673bSBarry Smith /*MC 401c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 402c838673bSBarry Smith 403c838673bSBarry Smith Level: beginner 404c838673bSBarry Smith 405c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 406c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 407c838673bSBarry Smith 2-norm of the residual for right preconditioning 408c838673bSBarry Smith 409c838673bSBarry Smith Level: beginner 410c838673bSBarry Smith 411c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 412c838673bSBarry Smith 413c838673bSBarry Smith M*/ 414c838673bSBarry Smith 415c838673bSBarry Smith /*MC 416c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 417c838673bSBarry Smith 418c838673bSBarry Smith Level: beginner 419c838673bSBarry Smith 420c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 421c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 422c838673bSBarry Smith 2-norm of the residual for right preconditioning 423c838673bSBarry Smith 424c838673bSBarry Smith Level: beginner 425c838673bSBarry Smith 426c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 427c838673bSBarry Smith 428c838673bSBarry Smith M*/ 429c838673bSBarry Smith 430c838673bSBarry Smith /*MC 431c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 432c838673bSBarry Smith reached 433c838673bSBarry Smith 434c838673bSBarry Smith Level: beginner 435c838673bSBarry Smith 436c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 437c838673bSBarry Smith 438c838673bSBarry Smith M*/ 439c838673bSBarry Smith 440c838673bSBarry Smith /*MC 4418c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4428c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4434d0a8057SBarry Smith test routine is set in KSP. 444c838673bSBarry Smith 445c838673bSBarry Smith 446c838673bSBarry Smith Level: beginner 447c838673bSBarry Smith 448c838673bSBarry Smith 449c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 450c838673bSBarry Smith 451c838673bSBarry Smith M*/ 452c838673bSBarry Smith 453c838673bSBarry Smith /*MC 454c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 455c838673bSBarry Smith method could not continue to enlarge the Krylov space. 456c838673bSBarry Smith 457c838673bSBarry Smith Level: beginner 458c838673bSBarry Smith 459c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 460c838673bSBarry Smith 461c838673bSBarry Smith M*/ 462c838673bSBarry Smith 463c838673bSBarry Smith /*MC 464c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 465c838673bSBarry Smith method could not continue to enlarge the Krylov space. 466c838673bSBarry Smith 467c838673bSBarry Smith 468c838673bSBarry Smith Level: beginner 469c838673bSBarry Smith 470c838673bSBarry Smith 471c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 472c838673bSBarry Smith 473c838673bSBarry Smith M*/ 474c838673bSBarry Smith 475c838673bSBarry Smith /*MC 476c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 477c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 478c838673bSBarry Smith 479c838673bSBarry Smith Level: beginner 480c838673bSBarry Smith 481c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 482c838673bSBarry Smith 483c838673bSBarry Smith M*/ 484c838673bSBarry Smith 485c838673bSBarry Smith /*MC 486c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 487c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 488c838673bSBarry Smith be positive definite 489c838673bSBarry Smith 490c838673bSBarry Smith Level: beginner 491c838673bSBarry Smith 4922401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 493c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 494c838673bSBarry Smith 495c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 496c838673bSBarry Smith 497c838673bSBarry Smith M*/ 498c838673bSBarry Smith 499c838673bSBarry Smith /*MC 500c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 501c838673bSBarry Smith while the KSPSolve() is still running. 502c838673bSBarry Smith 503c838673bSBarry Smith Level: beginner 504c838673bSBarry Smith 505c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 506c838673bSBarry Smith 507c838673bSBarry Smith M*/ 508c838673bSBarry Smith 509971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 510dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **); 511dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 51290984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 513971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *); 514971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **); 51594b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP); 5165f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP); 517dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 518dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *); 519abef13c0SSatish Balay 520dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *); 521d4fbbf0eSBarry Smith 52228ce4d24SBarry Smith /*E 52328ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 52428ce4d24SBarry Smith 52528ce4d24SBarry Smith Level: beginner 52628ce4d24SBarry Smith 52728ce4d24SBarry Smith .seealso: KSPCGSetType() 52828ce4d24SBarry Smith E*/ 5299dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5309dcbbd2bSBarry Smith extern const char *KSPCGTypes[]; 53128ce4d24SBarry Smith 532dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType); 533608df7e2SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGUseSingleReduction(KSP,PetscTruth); 5348031f4adStmunson 535fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal); 536fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *); 537fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *); 538fcae7a14Stmunson 5390d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal); 54011f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *); 541d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *); 542e559a7feSLois Curfman McInnes 5438031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal); 5448031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *); 545d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *); 546d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *); 5473c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *); 5488031f4adStmunson 5491d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPPythonSetType(KSP,const char[]); 5501d6018f0SLisandro Dalcin 551dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP); 552dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP); 5533369ce9aSBarry Smith 554a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 555a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 556a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG); 5577517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 5587517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 5597517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG); 560b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 561b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 562b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeDestroy(PetscDrawLG); 5632f2e5d10SKris Buschelman 5646891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 5656891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 56603919abeSBarry Smith 567f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 568436851aaSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscTruth monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 569f8a50e2bSBarry Smith 570f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 571f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess); 572f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess); 573f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec); 574f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 575436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess); 576f8a50e2bSBarry Smith 577f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 578f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess); 579f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*); 580f8a50e2bSBarry Smith 581084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 5823f22127dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetKSP(Mat,KSP*); 583084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 584084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 585dfe085dbSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 5863f22127dSBarry Smith 5876c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDM(KSP,DM); 588f22f69f0SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDMActive(KSP,PetscTruth); 5896c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDM(KSP,DM*); 5906c699258SBarry Smith 591e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 5922eac72dbSBarry Smith #endif 593