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" 60bbce1389SJed Brown #define KSPSPECEST "specest" 612eac72dbSBarry Smith 628ba1e511SMatthew Knepley /* Logging support */ 630700a824SBarry Smith extern PetscClassId PETSCKSP_DLLEXPORT KSP_CLASSID; 648ba1e511SMatthew Knepley 65dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *); 66a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,const KSPType); 67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP); 68dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP); 69dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec); 70dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec); 71dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP); 722eac72dbSBarry Smith 73b0a32e0cSBarry Smith extern PetscFList KSPList; 74ace3abfcSBarry Smith extern PetscBool KSPRegisterAllCalled; 75dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]); 76dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void); 77dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 7830de9b25SBarry Smith 7930de9b25SBarry Smith /*MC 8030de9b25SBarry Smith KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 8130de9b25SBarry Smith 8230de9b25SBarry Smith Synopsis: 831890ba74SBarry Smith PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP)) 8430de9b25SBarry Smith 8530de9b25SBarry Smith Not Collective 8630de9b25SBarry Smith 8730de9b25SBarry Smith Input Parameters: 8830de9b25SBarry Smith + name_solver - name of a new user-defined solver 8930de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 9030de9b25SBarry Smith . name_create - name of routine to create method context 9130de9b25SBarry Smith - routine_create - routine to create method context 9230de9b25SBarry Smith 9330de9b25SBarry Smith Notes: 9430de9b25SBarry Smith KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 9530de9b25SBarry Smith 9630de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 9730de9b25SBarry Smith is ignored. 9830de9b25SBarry Smith 9930de9b25SBarry Smith Sample usage: 10030de9b25SBarry Smith .vb 10130de9b25SBarry Smith KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 10230de9b25SBarry Smith "MySolverCreate",MySolverCreate); 10330de9b25SBarry Smith .ve 10430de9b25SBarry Smith 10530de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 10630de9b25SBarry Smith $ KSPSetType(ksp,"my_solver") 10730de9b25SBarry Smith or at runtime via the option 10830de9b25SBarry Smith $ -ksp_type my_solver 10930de9b25SBarry Smith 11030de9b25SBarry Smith Level: advanced 11130de9b25SBarry Smith 112ab901514SBarry Smith Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 11330de9b25SBarry Smith and others of the form ${any_environmental_variable} occuring in pathname will be 11430de9b25SBarry Smith replaced with appropriate values. 11530de9b25SBarry Smith If your function is not being put into a shared library then use KSPRegister() instead 11630de9b25SBarry Smith 11730de9b25SBarry Smith .keywords: KSP, register 11830de9b25SBarry Smith 11930de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy() 12030de9b25SBarry Smith 12130de9b25SBarry Smith M*/ 122aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 123f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 1246df38c32SLois Curfman McInnes #else 125f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 1266df38c32SLois Curfman McInnes #endif 12782bf6240SBarry Smith 128a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,const KSPType *); 1291718446fSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPCSide(KSP,PCSide); 130b037da10SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPCSide(KSP,PCSide*); 131dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 133ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscBool ); 134ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscBool *); 135ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscBool ); 136ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscBool *); 137ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetErrorIfNotConverged(KSP,PetscBool ); 138ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetErrorIfNotConverged(KSP,PetscBool *); 139ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscBool *); 140ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscBool ); 141ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscBool *); 142ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscBool ); 143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *); 144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *); 145dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*); 146dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*); 147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace); 148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*); 149906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1502eac72dbSBarry Smith 151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC); 152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*); 153aabeff55SBarry Smith 154a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 155a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP); 156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **); 157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 158ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1594b0e389bSBarry Smith 1600e33f6ddSBarry Smith /* not sure where to put this */ 161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*); 162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 1652eac72dbSBarry Smith 1660a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *); 1670a71e23eSMatthew Knepley 168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *); 169dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *); 1702eac72dbSBarry Smith 171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal); 172ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetSelfScale(KSP,PetscBool ); 173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 174dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 176dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1774b0e389bSBarry Smith 178dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt); 179e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetRestart(KSP, PetscInt*); 180dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal); 1819f236934SBarry Smith 182dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP); 183dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 184bef9c725SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 185dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 186dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1871d73ed98SBarry Smith 188dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt); 189dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP); 1901d73ed98SBarry Smith 19158ad63f4SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetRestart(KSP,PetscInt); 192e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRGetRestart(KSP,PetscInt*); 1937ae9e5e1SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 19458ad63f4SBarry Smith 195b49fd9e1SBarry Smith /*E 196b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 197b49fd9e1SBarry Smith 198b49fd9e1SBarry Smith Level: advanced 199b49fd9e1SBarry Smith 200bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 201e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 202b49fd9e1SBarry Smith 203b49fd9e1SBarry Smith E*/ 20478d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2059dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[]; 2061f7e983dSSatish Balay /*MC 2078c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2088c5b8ba0SBarry Smith 2098c5b8ba0SBarry Smith Level: advanced 2108c5b8ba0SBarry Smith 2118c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2128c5b8ba0SBarry Smith 213bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 214e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2158c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2168c5b8ba0SBarry Smith M*/ 2178c5b8ba0SBarry Smith 2181f7e983dSSatish Balay /*MC 2198c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2208c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2218c5b8ba0SBarry Smith poor orthogonality. 2228c5b8ba0SBarry Smith 2238c5b8ba0SBarry Smith Level: advanced 2248c5b8ba0SBarry Smith 2258c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2268c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2278c5b8ba0SBarry Smith 228bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 229e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2308c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2318c5b8ba0SBarry Smith M*/ 2328c5b8ba0SBarry Smith 2331f7e983dSSatish Balay /*MC 2348c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2358c5b8ba0SBarry Smith 2368c5b8ba0SBarry Smith Level: advanced 2378c5b8ba0SBarry Smith 2388c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2398c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2408c5b8ba0SBarry Smith 2418c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2428c5b8ba0SBarry Smith 243bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 244e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2458c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2468c5b8ba0SBarry Smith M*/ 2478c5b8ba0SBarry Smith 248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 249e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 25008480c60SBarry Smith 251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 252dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 253dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 254c38d4ed2SBarry Smith 255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal); 256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*); 257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*); 258121fd945SKris Buschelman 259d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal); 260ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscBool ); 2613cfecb52SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,PetscInt); 262d9492815SBarry Smith 263dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP); 264dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2652eac72dbSBarry Smith 266a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 267a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 26890984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *); 2691a2f9f99SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 2707517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 271a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 272a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 273a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 27484cb2905SBarry Smith 275dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec); 276dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*); 277dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 278dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 279c01c455dSBarry Smith 280dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure); 281dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 282ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 283dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]); 284dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]); 2852dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]); 2861eb62cbbSBarry Smith 287ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscBool ); 288ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscBool *); 289ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscBool ); 290ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscBool *); 2911f7f0c4fSBarry Smith 292dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer); 2931eb62cbbSBarry Smith 294db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec); 295db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*); 296db9b2ab1SHong Zhang 29728ce4d24SBarry Smith /*E 2988a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 2998a4b9c5cSBarry Smith test routines. 3008a4b9c5cSBarry Smith 3018a4b9c5cSBarry Smith Level: advanced 3028a4b9c5cSBarry Smith 303a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3041718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 305a3f661c8SBarry Smith 3068a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3078a4b9c5cSBarry Smith 30894b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3091718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3108a4b9c5cSBarry Smith E*/ 311ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3129dcbbd2bSBarry Smith extern const char *KSPNormTypes[]; 3131f7e983dSSatish Balay /*MC 314ce9499c7SBarry Smith KSP_NORM_NO - Do not compute a norm during the Krylov process. This will 3158c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3168c5b8ba0SBarry Smith be based on a norm of a residual etc. 3178c5b8ba0SBarry Smith 3188c5b8ba0SBarry Smith Level: advanced 3198c5b8ba0SBarry Smith 320085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3218c5b8ba0SBarry Smith 322ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3238c5b8ba0SBarry Smith M*/ 3248c5b8ba0SBarry Smith 3251f7e983dSSatish Balay /*MC 326ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3278c5b8ba0SBarry Smith convergence test routine. 3288c5b8ba0SBarry Smith 3298c5b8ba0SBarry Smith Level: advanced 3308c5b8ba0SBarry Smith 331ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3328c5b8ba0SBarry Smith M*/ 3338c5b8ba0SBarry Smith 3341f7e983dSSatish Balay /*MC 335ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3368c5b8ba0SBarry Smith convergence test routine. 3378c5b8ba0SBarry Smith 3388c5b8ba0SBarry Smith Level: advanced 3398c5b8ba0SBarry Smith 340ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3418c5b8ba0SBarry Smith M*/ 3428c5b8ba0SBarry Smith 3431f7e983dSSatish Balay /*MC 344ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 345a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3468c5b8ba0SBarry Smith 3478c5b8ba0SBarry Smith Level: advanced 3488c5b8ba0SBarry Smith 349ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3508c5b8ba0SBarry Smith M*/ 3518c5b8ba0SBarry Smith 352dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType); 353863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*); 3543c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt); 355ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetLagNorm(KSP,PetscBool ); 3568a4b9c5cSBarry Smith 3578a4b9c5cSBarry Smith /*E 35828ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 35928ce4d24SBarry Smith have converged or diverged 36028ce4d24SBarry Smith 36128ce4d24SBarry Smith Level: beginner 36228ce4d24SBarry Smith 3634d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 36428ce4d24SBarry Smith 3654d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3664d0a8057SBarry Smith 3674d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3684d0a8057SBarry Smith any of the values here also change them that array of names. 36986c02ca4SBarry Smith 370c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 37128ce4d24SBarry Smith E*/ 372d15094e1SBarry Smith typedef enum {/* converged */ 3739ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 3749ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 375d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 376d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 377b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3788031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3798031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 380329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 381af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 382d15094e1SBarry Smith /* diverged */ 383b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 384d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 385d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 386d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 387b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 388b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 389b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 3906aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 3916aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 392d15094e1SBarry Smith 393d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 3949dcbbd2bSBarry Smith extern const char **KSPConvergedReasons; 395d15094e1SBarry Smith 396c838673bSBarry Smith /*MC 397c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 398c838673bSBarry Smith 399c838673bSBarry Smith Level: beginner 400c838673bSBarry Smith 401c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 402c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 403c838673bSBarry Smith 2-norm of the residual for right preconditioning 404c838673bSBarry Smith 405c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 406c838673bSBarry Smith 407c838673bSBarry Smith M*/ 408c838673bSBarry Smith 409c838673bSBarry Smith /*MC 410c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 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_DTOL - norm(r) >= dtol*norm(b) 426c838673bSBarry Smith 427c838673bSBarry Smith Level: beginner 428c838673bSBarry Smith 429c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 430c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 431c838673bSBarry Smith 2-norm of the residual for right preconditioning 432c838673bSBarry Smith 433c838673bSBarry Smith Level: beginner 434c838673bSBarry Smith 435c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 436c838673bSBarry Smith 437c838673bSBarry Smith M*/ 438c838673bSBarry Smith 439c838673bSBarry Smith /*MC 440c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 441c838673bSBarry Smith reached 442c838673bSBarry Smith 443c838673bSBarry Smith Level: beginner 444c838673bSBarry Smith 445c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 446c838673bSBarry Smith 447c838673bSBarry Smith M*/ 448c838673bSBarry Smith 449c838673bSBarry Smith /*MC 4508c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4518c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4524d0a8057SBarry Smith test routine is set in KSP. 453c838673bSBarry Smith 454c838673bSBarry Smith 455c838673bSBarry Smith Level: beginner 456c838673bSBarry Smith 457c838673bSBarry Smith 458c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 459c838673bSBarry Smith 460c838673bSBarry Smith M*/ 461c838673bSBarry Smith 462c838673bSBarry Smith /*MC 463c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 464*3014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 465*3014e516SBarry Smith preconditioner. 466c838673bSBarry Smith 467c838673bSBarry Smith Level: beginner 468c838673bSBarry Smith 469c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 470c838673bSBarry Smith 471c838673bSBarry Smith M*/ 472c838673bSBarry Smith 473c838673bSBarry Smith /*MC 474c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 475c838673bSBarry Smith method could not continue to enlarge the Krylov space. 476c838673bSBarry Smith 477c838673bSBarry Smith 478c838673bSBarry Smith Level: beginner 479c838673bSBarry Smith 480c838673bSBarry Smith 481c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 482c838673bSBarry Smith 483c838673bSBarry Smith M*/ 484c838673bSBarry Smith 485c838673bSBarry Smith /*MC 486c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 487c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 488c838673bSBarry Smith 489c838673bSBarry Smith Level: beginner 490c838673bSBarry Smith 491c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 492c838673bSBarry Smith 493c838673bSBarry Smith M*/ 494c838673bSBarry Smith 495c838673bSBarry Smith /*MC 496c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 497c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 498c838673bSBarry Smith be positive definite 499c838673bSBarry Smith 500c838673bSBarry Smith Level: beginner 501c838673bSBarry Smith 5022401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 503c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 504c838673bSBarry Smith 505c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 506c838673bSBarry Smith 507c838673bSBarry Smith M*/ 508c838673bSBarry Smith 509c838673bSBarry Smith /*MC 510c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 511c838673bSBarry Smith while the KSPSolve() is still running. 512c838673bSBarry Smith 513c838673bSBarry Smith Level: beginner 514c838673bSBarry Smith 515c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 516c838673bSBarry Smith 517c838673bSBarry Smith M*/ 518c838673bSBarry Smith 519971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 520dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **); 521dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 52290984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 523971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *); 524971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **); 52594b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP); 5265f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP); 527dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 528dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *); 529abef13c0SSatish Balay 530dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *); 531d4fbbf0eSBarry Smith 53228ce4d24SBarry Smith /*E 53328ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 53428ce4d24SBarry Smith 53528ce4d24SBarry Smith Level: beginner 53628ce4d24SBarry Smith 53728ce4d24SBarry Smith .seealso: KSPCGSetType() 53828ce4d24SBarry Smith E*/ 5399dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5409dcbbd2bSBarry Smith extern const char *KSPCGTypes[]; 54128ce4d24SBarry Smith 542dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType); 543ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGUseSingleReduction(KSP,PetscBool ); 5448031f4adStmunson 545fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal); 546fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *); 547fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *); 548fcae7a14Stmunson 5490d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal); 55011f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *); 551d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *); 552e559a7feSLois Curfman McInnes 5538031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal); 5548031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *); 555d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *); 556d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *); 5573c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *); 5588031f4adStmunson 5591d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPPythonSetType(KSP,const char[]); 5601d6018f0SLisandro Dalcin 561dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP); 562dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP); 5633369ce9aSBarry Smith 564a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 565a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 566a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG); 5677517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 5687517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 5697517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG); 570b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 571b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 572b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeDestroy(PetscDrawLG); 5732f2e5d10SKris Buschelman 5746891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 5756891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 57603919abeSBarry Smith 577f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 578ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 579f8a50e2bSBarry Smith 580f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 581f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess); 582f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess); 583f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec); 584f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 585436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess); 586f8a50e2bSBarry Smith 587f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 588f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess); 589f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*); 590f8a50e2bSBarry Smith 591cfb32e20SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 5923f22127dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetKSP(Mat,KSP*); 593cfb32e20SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 594cfb32e20SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 595cfb32e20SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 5963f22127dSBarry Smith 5976c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDM(KSP,DM); 598ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDMActive(KSP,PetscBool ); 5996c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDM(KSP,DM*); 6006c699258SBarry Smith 601e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 6022eac72dbSBarry Smith #endif 603