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*); 136e113a28aSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetErrorIfNotConverged(KSP,PetscTruth); 137e113a28aSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetErrorIfNotConverged(KSP,PetscTruth *); 138a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*); 139dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth); 140a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*); 141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth); 142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *); 143dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *); 144dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*); 145dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*); 146dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace); 147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*); 148906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1492eac72dbSBarry Smith 150dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC); 151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*); 152aabeff55SBarry Smith 153a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 154a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP); 155dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **); 156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 157dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth); 1584b0e389bSBarry Smith 1590e33f6ddSBarry Smith /* not sure where to put this */ 160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*); 161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 1642eac72dbSBarry Smith 1650a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *); 1660a71e23eSMatthew Knepley 167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *); 168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *); 1692eac72dbSBarry Smith 170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal); 171b61692a6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetSelfScale(KSP,PetscTruth); 172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 174dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1764b0e389bSBarry Smith 177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt); 178e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetRestart(KSP, PetscInt*); 179dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal); 1809f236934SBarry Smith 181dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP); 182dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 183bef9c725SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 184dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 185dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1861d73ed98SBarry Smith 187dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt); 188dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP); 1891d73ed98SBarry Smith 19058ad63f4SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetRestart(KSP,PetscInt); 191e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRGetRestart(KSP,PetscInt*); 1927ae9e5e1SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 19358ad63f4SBarry Smith 194b49fd9e1SBarry Smith /*E 195b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 196b49fd9e1SBarry Smith 197b49fd9e1SBarry Smith Level: advanced 198b49fd9e1SBarry Smith 199bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 200e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 201b49fd9e1SBarry Smith 202b49fd9e1SBarry Smith E*/ 20378d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2049dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[]; 2051f7e983dSSatish Balay /*MC 2068c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2078c5b8ba0SBarry Smith 2088c5b8ba0SBarry Smith Level: advanced 2098c5b8ba0SBarry Smith 2108c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2118c5b8ba0SBarry Smith 212bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 213e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2148c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2158c5b8ba0SBarry Smith M*/ 2168c5b8ba0SBarry Smith 2171f7e983dSSatish Balay /*MC 2188c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2198c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2208c5b8ba0SBarry Smith poor orthogonality. 2218c5b8ba0SBarry Smith 2228c5b8ba0SBarry Smith Level: advanced 2238c5b8ba0SBarry Smith 2248c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2258c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2268c5b8ba0SBarry Smith 227bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 228e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2298c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2308c5b8ba0SBarry Smith M*/ 2318c5b8ba0SBarry Smith 2321f7e983dSSatish Balay /*MC 2338c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2348c5b8ba0SBarry Smith 2358c5b8ba0SBarry Smith Level: advanced 2368c5b8ba0SBarry Smith 2378c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2388c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2398c5b8ba0SBarry Smith 2408c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2418c5b8ba0SBarry Smith 242bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 243e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2448c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2458c5b8ba0SBarry Smith M*/ 2468c5b8ba0SBarry Smith 247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 248e928de20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 24908480c60SBarry Smith 250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 252dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 253c38d4ed2SBarry Smith 254dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal); 255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*); 256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*); 257121fd945SKris Buschelman 258d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal); 259d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth); 260*3cfecb52SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,PetscInt); 261d9492815SBarry Smith 262dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP); 263dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2642eac72dbSBarry Smith 265a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 266a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 26790984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultLSQR(KSP,PetscInt,PetscReal,void *); 2681a2f9f99SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 2697517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 270a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 271a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 272a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 27384cb2905SBarry Smith 274dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec); 275dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*); 276dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 277dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 278c01c455dSBarry Smith 279dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure); 280dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 281906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*); 282dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]); 283dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]); 2842dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]); 2851eb62cbbSBarry Smith 286dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth); 287dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*); 288dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth); 289dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*); 2901f7f0c4fSBarry Smith 291dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer); 2921eb62cbbSBarry Smith 293db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec); 294db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*); 295db9b2ab1SHong Zhang 29628ce4d24SBarry Smith /*E 2978a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 2988a4b9c5cSBarry Smith test routines. 2998a4b9c5cSBarry Smith 3008a4b9c5cSBarry Smith Level: advanced 3018a4b9c5cSBarry Smith 302a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3031718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 304a3f661c8SBarry Smith 3058a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3068a4b9c5cSBarry Smith 30794b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3081718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3098a4b9c5cSBarry Smith E*/ 310ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3119dcbbd2bSBarry Smith extern const char *KSPNormTypes[]; 3121f7e983dSSatish Balay /*MC 313ce9499c7SBarry Smith KSP_NORM_NO - Do not compute a norm during the Krylov process. This will 3148c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3158c5b8ba0SBarry Smith be based on a norm of a residual etc. 3168c5b8ba0SBarry Smith 3178c5b8ba0SBarry Smith Level: advanced 3188c5b8ba0SBarry Smith 319085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3208c5b8ba0SBarry Smith 321ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3228c5b8ba0SBarry Smith M*/ 3238c5b8ba0SBarry Smith 3241f7e983dSSatish Balay /*MC 325ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3268c5b8ba0SBarry Smith convergence test routine. 3278c5b8ba0SBarry Smith 3288c5b8ba0SBarry Smith Level: advanced 3298c5b8ba0SBarry Smith 330ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3318c5b8ba0SBarry Smith M*/ 3328c5b8ba0SBarry Smith 3331f7e983dSSatish Balay /*MC 334ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3358c5b8ba0SBarry Smith convergence test routine. 3368c5b8ba0SBarry Smith 3378c5b8ba0SBarry Smith Level: advanced 3388c5b8ba0SBarry Smith 339ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3408c5b8ba0SBarry Smith M*/ 3418c5b8ba0SBarry Smith 3421f7e983dSSatish Balay /*MC 343ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 344a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3458c5b8ba0SBarry Smith 3468c5b8ba0SBarry Smith Level: advanced 3478c5b8ba0SBarry Smith 348ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3498c5b8ba0SBarry Smith M*/ 3508c5b8ba0SBarry Smith 351dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType); 352863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*); 3533c5daeb9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetCheckNormIteration(KSP,PetscInt); 3541e7308c7SRichard Tran Mills EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetLagNorm(KSP,PetscTruth); 3558a4b9c5cSBarry Smith 3568a4b9c5cSBarry Smith /*E 35728ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 35828ce4d24SBarry Smith have converged or diverged 35928ce4d24SBarry Smith 36028ce4d24SBarry Smith Level: beginner 36128ce4d24SBarry Smith 3624d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 36328ce4d24SBarry Smith 3644d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3654d0a8057SBarry Smith 3664d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3674d0a8057SBarry Smith any of the values here also change them that array of names. 36886c02ca4SBarry Smith 369c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 37028ce4d24SBarry Smith E*/ 371d15094e1SBarry Smith typedef enum {/* converged */ 372d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 373d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 374b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 3758031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 3768031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 377329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 378af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 379d15094e1SBarry Smith /* diverged */ 380b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 381d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 382d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 383d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 384b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 385b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 386b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 3876aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 3886aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 389d15094e1SBarry Smith 390d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 3919dcbbd2bSBarry Smith extern const char **KSPConvergedReasons; 392d15094e1SBarry Smith 393c838673bSBarry Smith /*MC 394c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 395c838673bSBarry Smith 396c838673bSBarry Smith Level: beginner 397c838673bSBarry Smith 398c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 399c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 400c838673bSBarry Smith 2-norm of the residual for right preconditioning 401c838673bSBarry Smith 402c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 403c838673bSBarry Smith 404c838673bSBarry Smith M*/ 405c838673bSBarry Smith 406c838673bSBarry Smith /*MC 407c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 408c838673bSBarry Smith 409c838673bSBarry Smith Level: beginner 410c838673bSBarry Smith 411c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 412c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 413c838673bSBarry Smith 2-norm of the residual for right preconditioning 414c838673bSBarry Smith 415c838673bSBarry Smith Level: beginner 416c838673bSBarry Smith 417c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 418c838673bSBarry Smith 419c838673bSBarry Smith M*/ 420c838673bSBarry Smith 421c838673bSBarry Smith /*MC 422c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 423c838673bSBarry Smith 424c838673bSBarry Smith Level: beginner 425c838673bSBarry Smith 426c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 427c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 428c838673bSBarry Smith 2-norm of the residual for right preconditioning 429c838673bSBarry Smith 430c838673bSBarry Smith Level: beginner 431c838673bSBarry Smith 432c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 433c838673bSBarry Smith 434c838673bSBarry Smith M*/ 435c838673bSBarry Smith 436c838673bSBarry Smith /*MC 437c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 438c838673bSBarry Smith reached 439c838673bSBarry Smith 440c838673bSBarry Smith Level: beginner 441c838673bSBarry Smith 442c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 443c838673bSBarry Smith 444c838673bSBarry Smith M*/ 445c838673bSBarry Smith 446c838673bSBarry Smith /*MC 4478c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4488c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4494d0a8057SBarry Smith test routine is set in KSP. 450c838673bSBarry Smith 451c838673bSBarry Smith 452c838673bSBarry Smith Level: beginner 453c838673bSBarry Smith 454c838673bSBarry Smith 455c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 456c838673bSBarry Smith 457c838673bSBarry Smith M*/ 458c838673bSBarry Smith 459c838673bSBarry Smith /*MC 460c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 461c838673bSBarry Smith method could not continue to enlarge the Krylov space. 462c838673bSBarry Smith 463c838673bSBarry Smith Level: beginner 464c838673bSBarry Smith 465c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 466c838673bSBarry Smith 467c838673bSBarry Smith M*/ 468c838673bSBarry Smith 469c838673bSBarry Smith /*MC 470c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 471c838673bSBarry Smith method could not continue to enlarge the Krylov space. 472c838673bSBarry Smith 473c838673bSBarry Smith 474c838673bSBarry Smith Level: beginner 475c838673bSBarry Smith 476c838673bSBarry Smith 477c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 478c838673bSBarry Smith 479c838673bSBarry Smith M*/ 480c838673bSBarry Smith 481c838673bSBarry Smith /*MC 482c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 483c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 484c838673bSBarry Smith 485c838673bSBarry Smith Level: beginner 486c838673bSBarry Smith 487c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 488c838673bSBarry Smith 489c838673bSBarry Smith M*/ 490c838673bSBarry Smith 491c838673bSBarry Smith /*MC 492c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 493c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 494c838673bSBarry Smith be positive definite 495c838673bSBarry Smith 496c838673bSBarry Smith Level: beginner 497c838673bSBarry Smith 4982401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 499c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 500c838673bSBarry Smith 501c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 502c838673bSBarry Smith 503c838673bSBarry Smith M*/ 504c838673bSBarry Smith 505c838673bSBarry Smith /*MC 506c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 507c838673bSBarry Smith while the KSPSolve() is still running. 508c838673bSBarry Smith 509c838673bSBarry Smith Level: beginner 510c838673bSBarry Smith 511c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 512c838673bSBarry Smith 513c838673bSBarry Smith M*/ 514c838673bSBarry Smith 515971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 516dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **); 517dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 51890984118SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 519971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedDestroy(void *); 520971273eeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedCreate(void **); 52194b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP); 5225f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP); 523dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 524dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *); 525abef13c0SSatish Balay 526dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *); 527d4fbbf0eSBarry Smith 52828ce4d24SBarry Smith /*E 52928ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 53028ce4d24SBarry Smith 53128ce4d24SBarry Smith Level: beginner 53228ce4d24SBarry Smith 53328ce4d24SBarry Smith .seealso: KSPCGSetType() 53428ce4d24SBarry Smith E*/ 5359dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5369dcbbd2bSBarry Smith extern const char *KSPCGTypes[]; 53728ce4d24SBarry Smith 538dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType); 539608df7e2SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGUseSingleReduction(KSP,PetscTruth); 5408031f4adStmunson 541fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHSetRadius(KSP,PetscReal); 542fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetNormD(KSP,PetscReal *); 543fcae7a14Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPNASHGetObjFcn(KSP,PetscReal *); 544fcae7a14Stmunson 5450d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal); 54611f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *); 547d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *); 548e559a7feSLois Curfman McInnes 5498031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal); 5508031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *); 551d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *); 552d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *); 5533c851360Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetLambda(KSP,PetscReal *); 5548031f4adStmunson 5551d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPPythonSetType(KSP,const char[]); 5561d6018f0SLisandro Dalcin 557dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP); 558dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP); 5593369ce9aSBarry Smith 560a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 561a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*); 562a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG); 5637517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 5647517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 5657517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG); 566b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 567b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 568b271bb04SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGRangeDestroy(PetscDrawLG); 5692f2e5d10SKris Buschelman 5706891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 5716891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 57203919abeSBarry Smith 573f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 574436851aaSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscTruth monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 575f8a50e2bSBarry Smith 576f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 577f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessDestroy(KSPFischerGuess); 578f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessReset(KSPFischerGuess); 579f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessUpdate(KSPFischerGuess,Vec); 580f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 581436851aaSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFischerGuessSetFromOptions(KSPFischerGuess); 582f8a50e2bSBarry Smith 583f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 584f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFischerGuess(KSP,KSPFischerGuess); 585f8a50e2bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetFischerGuess(KSP,KSPFischerGuess*); 586f8a50e2bSBarry Smith 587084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 5883f22127dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatSchurComplementGetKSP(Mat,KSP*); 589084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 590084e4875SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 591dfe085dbSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 5923f22127dSBarry Smith 5936c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDM(KSP,DM); 594f22f69f0SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDMActive(KSP,PetscTruth); 5956c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDM(KSP,DM*); 5966c699258SBarry Smith 597e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 5982eac72dbSBarry Smith #endif 599