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" 597945016fSBarry 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); 176*e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetRestart(KSP, PetscInt*); 177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal); 1789f236934SBarry Smith 179dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP); 180dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 181dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 182dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1831d73ed98SBarry Smith 184dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt); 185dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP); 1861d73ed98SBarry Smith 18758ad63f4SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetRestart(KSP,PetscInt); 188*e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRGetRestart(KSP,PetscInt*); 1897ae9e5e1SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 19058ad63f4SBarry Smith 191b49fd9e1SBarry Smith /*E 192b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 193b49fd9e1SBarry Smith 194b49fd9e1SBarry Smith Level: advanced 195b49fd9e1SBarry Smith 196b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 197*e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 198b49fd9e1SBarry Smith 199b49fd9e1SBarry Smith E*/ 20078d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2019dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[]; 2021f7e983dSSatish Balay /*MC 2038c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2048c5b8ba0SBarry Smith 2058c5b8ba0SBarry Smith Level: advanced 2068c5b8ba0SBarry Smith 2078c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2088c5b8ba0SBarry Smith 2098c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 210*e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2118c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2128c5b8ba0SBarry Smith M*/ 2138c5b8ba0SBarry Smith 2141f7e983dSSatish Balay /*MC 2158c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2168c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2178c5b8ba0SBarry Smith poor orthogonality. 2188c5b8ba0SBarry Smith 2198c5b8ba0SBarry Smith Level: advanced 2208c5b8ba0SBarry Smith 2218c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2228c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2238c5b8ba0SBarry Smith 2248c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 225*e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2268c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2278c5b8ba0SBarry Smith M*/ 2288c5b8ba0SBarry Smith 2291f7e983dSSatish Balay /*MC 2308c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2318c5b8ba0SBarry Smith 2328c5b8ba0SBarry Smith Level: advanced 2338c5b8ba0SBarry Smith 2348c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2358c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2368c5b8ba0SBarry Smith 2378c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2388c5b8ba0SBarry Smith 2398c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 240*e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2418c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2428c5b8ba0SBarry Smith M*/ 2438c5b8ba0SBarry Smith 244dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 245*e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 24608480c60SBarry Smith 247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 249dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 250c38d4ed2SBarry Smith 251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal); 252dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*); 253dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*); 254121fd945SKris Buschelman 255d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal); 256d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth); 257d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int); 258d9492815SBarry Smith 259dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP); 260dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2612eac72dbSBarry Smith 262a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 263a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 26490984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *); 2651a2f9f99SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 2667517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 267a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 268a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 269a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 27084cb2905SBarry Smith 271dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec); 272dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*); 273dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 274dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 275c01c455dSBarry Smith 276dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure); 277dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 278906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*); 279dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]); 280dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]); 2812dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]); 2821eb62cbbSBarry Smith 283dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth); 284dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*); 285dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth); 286dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*); 2871f7f0c4fSBarry Smith 288dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer); 2891eb62cbbSBarry Smith 290db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec); 291db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*); 292db9b2ab1SHong Zhang 29328ce4d24SBarry Smith /*E 2948a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 2958a4b9c5cSBarry Smith test routines. 2968a4b9c5cSBarry Smith 2978a4b9c5cSBarry Smith Level: advanced 2988a4b9c5cSBarry Smith 299a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3001718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 301a3f661c8SBarry Smith 3028a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3038a4b9c5cSBarry Smith 30494b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3051718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3068a4b9c5cSBarry Smith E*/ 307ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3089dcbbd2bSBarry Smith extern const char *KSPNormTypes[]; 3091f7e983dSSatish Balay /*MC 310ce9499c7SBarry Smith KSP_NORM_NO - Do not compute a norm during the Krylov process. This will 3118c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3128c5b8ba0SBarry Smith be based on a norm of a residual etc. 3138c5b8ba0SBarry Smith 3148c5b8ba0SBarry Smith Level: advanced 3158c5b8ba0SBarry Smith 316085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3178c5b8ba0SBarry Smith 318ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3198c5b8ba0SBarry Smith M*/ 3208c5b8ba0SBarry Smith 3211f7e983dSSatish Balay /*MC 322ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual 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_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3288c5b8ba0SBarry Smith M*/ 3298c5b8ba0SBarry Smith 3301f7e983dSSatish Balay /*MC 331ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3328c5b8ba0SBarry Smith convergence test routine. 3338c5b8ba0SBarry Smith 3348c5b8ba0SBarry Smith Level: advanced 3358c5b8ba0SBarry Smith 336ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3378c5b8ba0SBarry Smith M*/ 3388c5b8ba0SBarry Smith 3391f7e983dSSatish Balay /*MC 340ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 341a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3428c5b8ba0SBarry Smith 3438c5b8ba0SBarry Smith Level: advanced 3448c5b8ba0SBarry Smith 345ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3468c5b8ba0SBarry Smith M*/ 3478c5b8ba0SBarry Smith 348dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType); 349863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*); 3503c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt); 3511e7308c7SRichard Tran Mills EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetLagNorm(KSP,PetscTruth); 3528a4b9c5cSBarry Smith 3538a4b9c5cSBarry Smith /*E 35428ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 35528ce4d24SBarry Smith have converged or diverged 35628ce4d24SBarry Smith 35728ce4d24SBarry Smith Level: beginner 35828ce4d24SBarry Smith 3594d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 36028ce4d24SBarry Smith 3614d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3624d0a8057SBarry Smith 3634d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3644d0a8057SBarry Smith any of the values here also change them that array of names. 36586c02ca4SBarry Smith 366c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 36728ce4d24SBarry Smith E*/ 368d15094e1SBarry Smith typedef enum {/* converged */ 369d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 370d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 371b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3728031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3738031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 374329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 375af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 376d15094e1SBarry Smith /* diverged */ 377b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 378d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 379d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 380d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 381b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 382b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 383b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 3846aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 3856aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 386d15094e1SBarry Smith 387d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 3889dcbbd2bSBarry Smith extern const char **KSPConvergedReasons; 389d15094e1SBarry Smith 390c838673bSBarry Smith /*MC 391c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 392c838673bSBarry Smith 393c838673bSBarry Smith Level: beginner 394c838673bSBarry Smith 395c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 396c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 397c838673bSBarry Smith 2-norm of the residual for right preconditioning 398c838673bSBarry Smith 399c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 400c838673bSBarry Smith 401c838673bSBarry Smith M*/ 402c838673bSBarry Smith 403c838673bSBarry Smith /*MC 404c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 405c838673bSBarry Smith 406c838673bSBarry Smith Level: beginner 407c838673bSBarry Smith 408c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 409c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 410c838673bSBarry Smith 2-norm of the residual for right preconditioning 411c838673bSBarry Smith 412c838673bSBarry Smith Level: beginner 413c838673bSBarry Smith 414c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 415c838673bSBarry Smith 416c838673bSBarry Smith M*/ 417c838673bSBarry Smith 418c838673bSBarry Smith /*MC 419c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 420c838673bSBarry Smith 421c838673bSBarry Smith Level: beginner 422c838673bSBarry Smith 423c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 424c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 425c838673bSBarry Smith 2-norm of the residual for right preconditioning 426c838673bSBarry Smith 427c838673bSBarry Smith Level: beginner 428c838673bSBarry Smith 429c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 430c838673bSBarry Smith 431c838673bSBarry Smith M*/ 432c838673bSBarry Smith 433c838673bSBarry Smith /*MC 434c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 435c838673bSBarry Smith reached 436c838673bSBarry Smith 437c838673bSBarry Smith Level: beginner 438c838673bSBarry Smith 439c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 440c838673bSBarry Smith 441c838673bSBarry Smith M*/ 442c838673bSBarry Smith 443c838673bSBarry Smith /*MC 4448c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4458c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4464d0a8057SBarry Smith test routine is set in KSP. 447c838673bSBarry Smith 448c838673bSBarry Smith 449c838673bSBarry Smith Level: beginner 450c838673bSBarry Smith 451c838673bSBarry Smith 452c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 453c838673bSBarry Smith 454c838673bSBarry Smith M*/ 455c838673bSBarry Smith 456c838673bSBarry Smith /*MC 457c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 458c838673bSBarry Smith method could not continue to enlarge the Krylov space. 459c838673bSBarry Smith 460c838673bSBarry Smith Level: beginner 461c838673bSBarry Smith 462c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 463c838673bSBarry Smith 464c838673bSBarry Smith M*/ 465c838673bSBarry Smith 466c838673bSBarry Smith /*MC 467c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 468c838673bSBarry Smith method could not continue to enlarge the Krylov space. 469c838673bSBarry Smith 470c838673bSBarry Smith 471c838673bSBarry Smith Level: beginner 472c838673bSBarry Smith 473c838673bSBarry Smith 474c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 475c838673bSBarry Smith 476c838673bSBarry Smith M*/ 477c838673bSBarry Smith 478c838673bSBarry Smith /*MC 479c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 480c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 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_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 490c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 491c838673bSBarry Smith be positive definite 492c838673bSBarry Smith 493c838673bSBarry Smith Level: beginner 494c838673bSBarry Smith 4952401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 496c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 497c838673bSBarry Smith 498c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 499c838673bSBarry Smith 500c838673bSBarry Smith M*/ 501c838673bSBarry Smith 502c838673bSBarry Smith /*MC 503c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 504c838673bSBarry Smith while the KSPSolve() is still running. 505c838673bSBarry Smith 506c838673bSBarry Smith Level: beginner 507c838673bSBarry Smith 508c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 509c838673bSBarry Smith 510c838673bSBarry Smith M*/ 511c838673bSBarry Smith 512971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 513dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **); 514dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 51590984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 516971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *); 517971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **); 51894b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP); 5195f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP); 520dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 521dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *); 522abef13c0SSatish Balay 523dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *); 524d4fbbf0eSBarry Smith 52528ce4d24SBarry Smith /*E 52628ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 52728ce4d24SBarry Smith 52828ce4d24SBarry Smith Level: beginner 52928ce4d24SBarry Smith 53028ce4d24SBarry Smith .seealso: KSPCGSetType() 53128ce4d24SBarry Smith E*/ 5329dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5339dcbbd2bSBarry Smith extern const char *KSPCGTypes[]; 53428ce4d24SBarry Smith 535dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType); 536608df7e2SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGUseSingleReduction(KSP,PetscTruth); 5378031f4adStmunson 538fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal); 539fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *); 540fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *); 541fcae7a14Stmunson 5420d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal); 54311f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *); 544d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *); 545e559a7feSLois Curfman McInnes 5468031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal); 5478031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *); 548d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *); 549d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *); 5503c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *); 5518031f4adStmunson 5521d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPPythonSetType(KSP,const char[]); 5531d6018f0SLisandro Dalcin 554dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP); 555dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP); 5563369ce9aSBarry Smith 557a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 558a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 559a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG); 5607517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 5617517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 5627517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG); 563b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 564b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 565b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeDestroy(PetscDrawLG); 5662f2e5d10SKris Buschelman 5676891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 5686891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 56903919abeSBarry Smith 570f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 571436851aaSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscTruth monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 572f8a50e2bSBarry Smith 573f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 574f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess); 575f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess); 576f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec); 577f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 578436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess); 579f8a50e2bSBarry Smith 580f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 581f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess); 582f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*); 583f8a50e2bSBarry Smith 584084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 5853f22127dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetKSP(Mat,KSP*); 586084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 587084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 588dfe085dbSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 5893f22127dSBarry Smith 5906c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDM(KSP,DM); 591f22f69f0SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDMActive(KSP,PetscTruth); 5926c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDM(KSP,DM*); 5936c699258SBarry Smith 594e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 5952eac72dbSBarry Smith #endif 596