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 62c8e378dSBarry Smith #include <petscpc.h> 72eac72dbSBarry Smith 8014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitializePackage(const char[]); 91dbb0a54SBarry Smith 1028ce4d24SBarry Smith /*S 1128ce4d24SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods 1228ce4d24SBarry Smith 1328ce4d24SBarry Smith Level: beginner 1428ce4d24SBarry Smith 1528ce4d24SBarry Smith Concepts: Krylov methods 1628ce4d24SBarry Smith 1794b7f48cSBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP 1828ce4d24SBarry Smith S*/ 19e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 202eac72dbSBarry Smith 2176bdecfbSBarry Smith /*J 2228ce4d24SBarry Smith KSPType - String with the name of a PETSc Krylov method or the creation function 2328ce4d24SBarry Smith with an optional dynamic library name, for example 2428ce4d24SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 2528ce4d24SBarry Smith 2628ce4d24SBarry Smith Level: beginner 2728ce4d24SBarry Smith 2828ce4d24SBarry Smith .seealso: KSPSetType(), KSP 2976bdecfbSBarry Smith J*/ 3019fd82e9SBarry Smith typedef const char* KSPType; 3182bf6240SBarry Smith #define KSPRICHARDSON "richardson" 326c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3382bf6240SBarry Smith #define KSPCG "cg" 342c8fc646SJed Brown #define KSPGROPPCG "groppcg" 352c8fc646SJed Brown #define KSPPIPECG "pipecg" 36df2a969aSvictorle #define KSPCGNE "cgne" 37fcae7a14Stmunson #define KSPNASH "nash" 3880e17431SMatthew Knepley #define KSPSTCG "stcg" 398031f4adStmunson #define KSPGLTR "gltr" 4082bf6240SBarry Smith #define KSPGMRES "gmres" 419dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 429dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 434f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 4461b468f9SJed Brown #define KSPPGMRES "pgmres" 4582bf6240SBarry Smith #define KSPTCQMR "tcqmr" 4682bf6240SBarry Smith #define KSPBCGS "bcgs" 47e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 4818ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 49c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 50f0037002Svictorle #define KSPBCGSL "bcgsl" 5182bf6240SBarry Smith #define KSPCGS "cgs" 5282bf6240SBarry Smith #define KSPTFQMR "tfqmr" 5382bf6240SBarry Smith #define KSPCR "cr" 542c8fc646SJed Brown #define KSPPIPECR "pipecr" 5582bf6240SBarry Smith #define KSPLSQR "lsqr" 5682bf6240SBarry Smith #define KSPPREONLY "preonly" 5782bf6240SBarry Smith #define KSPQCG "qcg" 58c9cf9b11SBarry Smith #define KSPBICG "bicg" 59b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6001c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 6121356919SSatish Balay #define KSPLCD "lcd" 621d6018f0SLisandro Dalcin #define KSPPYTHON "python" 6358ad63f4SBarry Smith #define KSPGCR "gcr" 64bbce1389SJed Brown #define KSPSPECEST "specest" 652eac72dbSBarry Smith 668ba1e511SMatthew Knepley /* Logging support */ 67014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 685358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 698ba1e511SMatthew Knepley 70014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 7119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 72014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 73014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 74014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 75014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 76014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 782eac72dbSBarry Smith 79014dd563SJed Brown PETSC_EXTERN PetscFList KSPList; 80014dd563SJed Brown PETSC_EXTERN PetscBool KSPRegisterAllCalled; 81014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterAll(const char[]); 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegisterDestroy(void); 83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 84f550243cSJed Brown PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(const char[]); 8530de9b25SBarry Smith 8630de9b25SBarry Smith /*MC 8730de9b25SBarry Smith KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 8830de9b25SBarry Smith 8930de9b25SBarry Smith Synopsis: 90*f2ba6396SBarry Smith #include "petscksp.h" 911890ba74SBarry Smith PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP)) 9230de9b25SBarry Smith 9330de9b25SBarry Smith Not Collective 9430de9b25SBarry Smith 9530de9b25SBarry Smith Input Parameters: 9630de9b25SBarry Smith + name_solver - name of a new user-defined solver 9730de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 9830de9b25SBarry Smith . name_create - name of routine to create method context 9930de9b25SBarry Smith - routine_create - routine to create method context 10030de9b25SBarry Smith 10130de9b25SBarry Smith Notes: 10230de9b25SBarry Smith KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 10330de9b25SBarry Smith 10430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 10530de9b25SBarry Smith is ignored. 10630de9b25SBarry Smith 10730de9b25SBarry Smith Sample usage: 10830de9b25SBarry Smith .vb 10930de9b25SBarry Smith KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 11030de9b25SBarry Smith "MySolverCreate",MySolverCreate); 11130de9b25SBarry Smith .ve 11230de9b25SBarry Smith 11330de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 11430de9b25SBarry Smith $ KSPSetType(ksp,"my_solver") 11530de9b25SBarry Smith or at runtime via the option 11630de9b25SBarry Smith $ -ksp_type my_solver 11730de9b25SBarry Smith 11830de9b25SBarry Smith Level: advanced 11930de9b25SBarry Smith 120ab901514SBarry Smith Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 12130de9b25SBarry Smith and others of the form ${any_environmental_variable} occuring in pathname will be 12230de9b25SBarry Smith replaced with appropriate values. 12330de9b25SBarry Smith If your function is not being put into a shared library then use KSPRegister() instead 12430de9b25SBarry Smith 12530de9b25SBarry Smith .keywords: KSP, register 12630de9b25SBarry Smith 12730de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy() 12830de9b25SBarry Smith 12930de9b25SBarry Smith M*/ 130aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 131f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 1326df38c32SLois Curfman McInnes #else 133f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 1346df38c32SLois Curfman McInnes #endif 13582bf6240SBarry Smith 13619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1582eac72dbSBarry Smith 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 161aabeff55SBarry Smith 162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 163014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 164014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 166014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 167014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1684b0e389bSBarry Smith 1690e33f6ddSBarry Smith /* not sure where to put this */ 170014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 171014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 172014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 1752eac72dbSBarry Smith 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1770a71e23eSMatthew Knepley 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1802eac72dbSBarry Smith 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 185ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1904b0e389bSBarry Smith 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1949f236934SBarry Smith 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2001d73ed98SBarry Smith 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2031d73ed98SBarry Smith 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 20758ad63f4SBarry Smith 208b49fd9e1SBarry Smith /*E 209b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 210b49fd9e1SBarry Smith 211b49fd9e1SBarry Smith Level: advanced 212b49fd9e1SBarry Smith 213bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 214e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 215b49fd9e1SBarry Smith 216b49fd9e1SBarry Smith E*/ 21778d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2186a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2191f7e983dSSatish Balay /*MC 2208c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2218c5b8ba0SBarry Smith 2228c5b8ba0SBarry Smith Level: advanced 2238c5b8ba0SBarry Smith 2248c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2258c5b8ba0SBarry Smith 226bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 227e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2288c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2298c5b8ba0SBarry Smith M*/ 2308c5b8ba0SBarry Smith 2311f7e983dSSatish Balay /*MC 2328c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2338c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2348c5b8ba0SBarry Smith poor orthogonality. 2358c5b8ba0SBarry Smith 2368c5b8ba0SBarry Smith Level: advanced 2378c5b8ba0SBarry Smith 2388c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2398c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2408c5b8ba0SBarry Smith 241bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 242e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2438c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2448c5b8ba0SBarry Smith M*/ 2458c5b8ba0SBarry Smith 2461f7e983dSSatish Balay /*MC 2478c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2488c5b8ba0SBarry Smith 2498c5b8ba0SBarry Smith Level: advanced 2508c5b8ba0SBarry Smith 2518c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2528c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2538c5b8ba0SBarry Smith 2548c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2558c5b8ba0SBarry Smith 256bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 257e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2588c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2598c5b8ba0SBarry Smith M*/ 2608c5b8ba0SBarry Smith 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 26308480c60SBarry Smith 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 267c38d4ed2SBarry Smith 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 271121fd945SKris Buschelman 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 275e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 276d9492815SBarry Smith 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2792eac72dbSBarry Smith 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 284390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 285390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 287816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 29484cb2905SBarry Smith 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 299c01c455dSBarry Smith 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3067287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3077287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3081eb62cbbSBarry Smith 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3131f7f0c4fSBarry Smith 314014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 31555849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 31655849f57SBarry Smith 31755849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3181eb62cbbSBarry Smith 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 320014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 321db9b2ab1SHong Zhang 322014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 323014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 32483ab6a24SBarry Smith 32528ce4d24SBarry Smith /*E 3268a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3278a4b9c5cSBarry Smith test routines. 3288a4b9c5cSBarry Smith 3298a4b9c5cSBarry Smith Level: advanced 3308a4b9c5cSBarry Smith 331a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3321718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 333a3f661c8SBarry Smith 3348a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3358a4b9c5cSBarry Smith 33694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3371718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3388a4b9c5cSBarry Smith E*/ 3399e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3409e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 341014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 342a21b2a99SBarry Smith 3431f7e983dSSatish Balay /*MC 3449793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3458c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3468c5b8ba0SBarry Smith be based on a norm of a residual etc. 3478c5b8ba0SBarry Smith 3488c5b8ba0SBarry Smith Level: advanced 3498c5b8ba0SBarry Smith 350085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3518c5b8ba0SBarry Smith 352ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3538c5b8ba0SBarry Smith M*/ 3548c5b8ba0SBarry Smith 3551f7e983dSSatish Balay /*MC 356ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3578c5b8ba0SBarry Smith convergence test routine. 3588c5b8ba0SBarry Smith 3598c5b8ba0SBarry Smith Level: advanced 3608c5b8ba0SBarry Smith 3619793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3628c5b8ba0SBarry Smith M*/ 3638c5b8ba0SBarry Smith 3641f7e983dSSatish Balay /*MC 365ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3668c5b8ba0SBarry Smith convergence test routine. 3678c5b8ba0SBarry Smith 3688c5b8ba0SBarry Smith Level: advanced 3698c5b8ba0SBarry Smith 3709793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3718c5b8ba0SBarry Smith M*/ 3728c5b8ba0SBarry Smith 3731f7e983dSSatish Balay /*MC 374ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 375a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3768c5b8ba0SBarry Smith 3778c5b8ba0SBarry Smith Level: advanced 3788c5b8ba0SBarry Smith 3799793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3808c5b8ba0SBarry Smith M*/ 3818c5b8ba0SBarry Smith 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 3878a4b9c5cSBarry Smith 3888a4b9c5cSBarry Smith /*E 38928ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 39028ce4d24SBarry Smith have converged or diverged 39128ce4d24SBarry Smith 39228ce4d24SBarry Smith Level: beginner 39328ce4d24SBarry Smith 3944d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 39528ce4d24SBarry Smith 3964d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3974d0a8057SBarry Smith 3984d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3994d0a8057SBarry Smith any of the values here also change them that array of names. 40086c02ca4SBarry Smith 401c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 40228ce4d24SBarry Smith E*/ 403d15094e1SBarry Smith typedef enum {/* converged */ 4049ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4059ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 406d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 407d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 408b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4098031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4108031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 411329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 412af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 413d15094e1SBarry Smith /* diverged */ 414b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 415d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 416d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 417d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 418b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 419b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 420b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4216aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 4226aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 423d15094e1SBarry Smith 424d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 425014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 426d15094e1SBarry Smith 427c838673bSBarry Smith /*MC 428c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 429c838673bSBarry Smith 430c838673bSBarry Smith Level: beginner 431c838673bSBarry Smith 432c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 433c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 434c838673bSBarry Smith 2-norm of the residual for right preconditioning 435c838673bSBarry Smith 436c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 437c838673bSBarry Smith 438c838673bSBarry Smith M*/ 439c838673bSBarry Smith 440c838673bSBarry Smith /*MC 441c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 442c838673bSBarry Smith 443c838673bSBarry Smith Level: beginner 444c838673bSBarry Smith 445c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 446c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 447c838673bSBarry Smith 2-norm of the residual for right preconditioning 448c838673bSBarry Smith 449c838673bSBarry Smith Level: beginner 450c838673bSBarry Smith 451c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 452c838673bSBarry Smith 453c838673bSBarry Smith M*/ 454c838673bSBarry Smith 455c838673bSBarry Smith /*MC 456c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 457c838673bSBarry Smith 458c838673bSBarry Smith Level: beginner 459c838673bSBarry Smith 460c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 461c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 462c838673bSBarry Smith 2-norm of the residual for right preconditioning 463c838673bSBarry Smith 464c838673bSBarry Smith Level: beginner 465c838673bSBarry Smith 466c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 467c838673bSBarry Smith 468c838673bSBarry Smith M*/ 469c838673bSBarry Smith 470c838673bSBarry Smith /*MC 471c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 472c838673bSBarry Smith reached 473c838673bSBarry Smith 474c838673bSBarry Smith Level: beginner 475c838673bSBarry Smith 476c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 477c838673bSBarry Smith 478c838673bSBarry Smith M*/ 479c838673bSBarry Smith 480c838673bSBarry Smith /*MC 4818c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4828c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4834d0a8057SBarry Smith test routine is set in KSP. 484c838673bSBarry Smith 485c838673bSBarry Smith 486c838673bSBarry Smith Level: beginner 487c838673bSBarry Smith 488c838673bSBarry Smith 489c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 490c838673bSBarry Smith 491c838673bSBarry Smith M*/ 492c838673bSBarry Smith 493c838673bSBarry Smith /*MC 494c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4953014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4963014e516SBarry Smith preconditioner. 497c838673bSBarry Smith 498c838673bSBarry Smith Level: beginner 499c838673bSBarry Smith 500c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 501c838673bSBarry Smith 502c838673bSBarry Smith M*/ 503c838673bSBarry Smith 504c838673bSBarry Smith /*MC 505c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 506c838673bSBarry Smith method could not continue to enlarge the Krylov space. 507c838673bSBarry Smith 508c838673bSBarry Smith 509c838673bSBarry Smith Level: beginner 510c838673bSBarry Smith 511c838673bSBarry Smith 512c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 513c838673bSBarry Smith 514c838673bSBarry Smith M*/ 515c838673bSBarry Smith 516c838673bSBarry Smith /*MC 517c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 518c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 519c838673bSBarry Smith 520c838673bSBarry Smith Level: beginner 521c838673bSBarry Smith 522c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 523c838673bSBarry Smith 524c838673bSBarry Smith M*/ 525c838673bSBarry Smith 526c838673bSBarry Smith /*MC 527c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 528c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 529c838673bSBarry Smith be positive definite 530c838673bSBarry Smith 531c838673bSBarry Smith Level: beginner 532c838673bSBarry Smith 5332401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 534c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 535c838673bSBarry Smith 536c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 537c838673bSBarry Smith 538c838673bSBarry Smith M*/ 539c838673bSBarry Smith 540c838673bSBarry Smith /*MC 541c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 542c838673bSBarry Smith while the KSPSolve() is still running. 543c838673bSBarry Smith 544c838673bSBarry Smith Level: beginner 545c838673bSBarry Smith 546c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 547c838673bSBarry Smith 548c838673bSBarry Smith M*/ 549c838673bSBarry Smith 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 560abef13c0SSatish Balay 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 562d4fbbf0eSBarry Smith 56328ce4d24SBarry Smith /*E 56428ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 56528ce4d24SBarry Smith 56628ce4d24SBarry Smith Level: beginner 56728ce4d24SBarry Smith 56828ce4d24SBarry Smith .seealso: KSPCGSetType() 56928ce4d24SBarry Smith E*/ 5709dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5716a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 57228ce4d24SBarry Smith 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5758031f4adStmunson 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 579fcae7a14Stmunson 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 583e559a7feSLois Curfman McInnes 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5898031f4adStmunson 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5911d6018f0SLisandro Dalcin 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5943369ce9aSBarry Smith 5951d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 5961d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 5971d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6022f2e5d10SKris Buschelman 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 60503919abeSBarry Smith 606f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 607ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 608f8a50e2bSBarry Smith 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 615f8a50e2bSBarry Smith 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 619f8a50e2bSBarry Smith 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 622d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 6231324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat); 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 625014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 626014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 6273f22127dSBarry Smith 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat); 629fcfd50ebSBarry Smith 630014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 631014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 632014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 633014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 634014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 635fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 637fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 640014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 641014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 642fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 643fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6446c699258SBarry Smith 6452eac72dbSBarry Smith #endif 646