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" 34*2c8fc646SJed Brown #define KSPGROPPCG "groppcg" 35*2c8fc646SJed 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" 54*2c8fc646SJed 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: 901890ba74SBarry Smith PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP)) 9130de9b25SBarry Smith 9230de9b25SBarry Smith Not Collective 9330de9b25SBarry Smith 9430de9b25SBarry Smith Input Parameters: 9530de9b25SBarry Smith + name_solver - name of a new user-defined solver 9630de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 9730de9b25SBarry Smith . name_create - name of routine to create method context 9830de9b25SBarry Smith - routine_create - routine to create method context 9930de9b25SBarry Smith 10030de9b25SBarry Smith Notes: 10130de9b25SBarry Smith KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 10230de9b25SBarry Smith 10330de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 10430de9b25SBarry Smith is ignored. 10530de9b25SBarry Smith 10630de9b25SBarry Smith Sample usage: 10730de9b25SBarry Smith .vb 10830de9b25SBarry Smith KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 10930de9b25SBarry Smith "MySolverCreate",MySolverCreate); 11030de9b25SBarry Smith .ve 11130de9b25SBarry Smith 11230de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 11330de9b25SBarry Smith $ KSPSetType(ksp,"my_solver") 11430de9b25SBarry Smith or at runtime via the option 11530de9b25SBarry Smith $ -ksp_type my_solver 11630de9b25SBarry Smith 11730de9b25SBarry Smith Level: advanced 11830de9b25SBarry Smith 119ab901514SBarry Smith Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 12030de9b25SBarry Smith and others of the form ${any_environmental_variable} occuring in pathname will be 12130de9b25SBarry Smith replaced with appropriate values. 12230de9b25SBarry Smith If your function is not being put into a shared library then use KSPRegister() instead 12330de9b25SBarry Smith 12430de9b25SBarry Smith .keywords: KSP, register 12530de9b25SBarry Smith 12630de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy() 12730de9b25SBarry Smith 12830de9b25SBarry Smith M*/ 129aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 130f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 1316df38c32SLois Curfman McInnes #else 132f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 1336df38c32SLois Curfman McInnes #endif 13482bf6240SBarry Smith 13519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1572eac72dbSBarry Smith 158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 160aabeff55SBarry Smith 161014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 163014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 164014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 166014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1674b0e389bSBarry Smith 1680e33f6ddSBarry Smith /* not sure where to put this */ 169014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 170014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 171014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 172014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 1742eac72dbSBarry Smith 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1760a71e23eSMatthew Knepley 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1792eac72dbSBarry Smith 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1884b0e389bSBarry Smith 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1929f236934SBarry Smith 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1981d73ed98SBarry Smith 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2011d73ed98SBarry Smith 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 20558ad63f4SBarry Smith 206b49fd9e1SBarry Smith /*E 207b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 208b49fd9e1SBarry Smith 209b49fd9e1SBarry Smith Level: advanced 210b49fd9e1SBarry Smith 211bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 212e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 213b49fd9e1SBarry Smith 214b49fd9e1SBarry Smith E*/ 21578d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2166a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2171f7e983dSSatish Balay /*MC 2188c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2198c5b8ba0SBarry Smith 2208c5b8ba0SBarry Smith Level: advanced 2218c5b8ba0SBarry Smith 2228c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2238c5b8ba0SBarry Smith 224bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 225e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2268c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2278c5b8ba0SBarry Smith M*/ 2288c5b8ba0SBarry Smith 2291f7e983dSSatish Balay /*MC 2308c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2318c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2328c5b8ba0SBarry Smith poor orthogonality. 2338c5b8ba0SBarry Smith 2348c5b8ba0SBarry Smith Level: advanced 2358c5b8ba0SBarry Smith 2368c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2378c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2388c5b8ba0SBarry Smith 239bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 240e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2418c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2428c5b8ba0SBarry Smith M*/ 2438c5b8ba0SBarry Smith 2441f7e983dSSatish Balay /*MC 2458c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2468c5b8ba0SBarry Smith 2478c5b8ba0SBarry Smith Level: advanced 2488c5b8ba0SBarry Smith 2498c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2508c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2518c5b8ba0SBarry Smith 2528c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2538c5b8ba0SBarry Smith 254bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 255e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2568c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2578c5b8ba0SBarry Smith M*/ 2588c5b8ba0SBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 26108480c60SBarry Smith 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 265c38d4ed2SBarry Smith 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 269121fd945SKris Buschelman 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 273d9492815SBarry Smith 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2762eac72dbSBarry Smith 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 281390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 282390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 284816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 29184cb2905SBarry Smith 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 296c01c455dSBarry Smith 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3037287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3047287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3051eb62cbbSBarry Smith 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3101f7f0c4fSBarry Smith 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 31255849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 31355849f57SBarry Smith 31455849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3151eb62cbbSBarry Smith 316014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 318db9b2ab1SHong Zhang 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 320014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 32183ab6a24SBarry Smith 32228ce4d24SBarry Smith /*E 3238a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3248a4b9c5cSBarry Smith test routines. 3258a4b9c5cSBarry Smith 3268a4b9c5cSBarry Smith Level: advanced 3278a4b9c5cSBarry Smith 328a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3291718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 330a3f661c8SBarry Smith 3318a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3328a4b9c5cSBarry Smith 33394b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3341718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3358a4b9c5cSBarry Smith E*/ 3369e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3379e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 338014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 339a21b2a99SBarry Smith 3401f7e983dSSatish Balay /*MC 3419793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3428c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3438c5b8ba0SBarry Smith be based on a norm of a residual etc. 3448c5b8ba0SBarry Smith 3458c5b8ba0SBarry Smith Level: advanced 3468c5b8ba0SBarry Smith 347085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3488c5b8ba0SBarry Smith 349ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3508c5b8ba0SBarry Smith M*/ 3518c5b8ba0SBarry Smith 3521f7e983dSSatish Balay /*MC 353ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3548c5b8ba0SBarry Smith convergence test routine. 3558c5b8ba0SBarry Smith 3568c5b8ba0SBarry Smith Level: advanced 3578c5b8ba0SBarry Smith 3589793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3598c5b8ba0SBarry Smith M*/ 3608c5b8ba0SBarry Smith 3611f7e983dSSatish Balay /*MC 362ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3638c5b8ba0SBarry Smith convergence test routine. 3648c5b8ba0SBarry Smith 3658c5b8ba0SBarry Smith Level: advanced 3668c5b8ba0SBarry Smith 3679793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3688c5b8ba0SBarry Smith M*/ 3698c5b8ba0SBarry Smith 3701f7e983dSSatish Balay /*MC 371ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 372a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3738c5b8ba0SBarry Smith 3748c5b8ba0SBarry Smith Level: advanced 3758c5b8ba0SBarry Smith 3769793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3778c5b8ba0SBarry Smith M*/ 3788c5b8ba0SBarry Smith 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 3848a4b9c5cSBarry Smith 3858a4b9c5cSBarry Smith /*E 38628ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 38728ce4d24SBarry Smith have converged or diverged 38828ce4d24SBarry Smith 38928ce4d24SBarry Smith Level: beginner 39028ce4d24SBarry Smith 3914d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 39228ce4d24SBarry Smith 3934d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3944d0a8057SBarry Smith 3954d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3964d0a8057SBarry Smith any of the values here also change them that array of names. 39786c02ca4SBarry Smith 398c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 39928ce4d24SBarry Smith E*/ 400d15094e1SBarry Smith typedef enum {/* converged */ 4019ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4029ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 403d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 404d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 405b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4068031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4078031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 408329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 409af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 410d15094e1SBarry Smith /* diverged */ 411b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 412d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 413d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 414d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 415b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 416b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 417b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4186aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 4196aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 420d15094e1SBarry Smith 421d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 422014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 423d15094e1SBarry Smith 424c838673bSBarry Smith /*MC 425c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*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 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 434c838673bSBarry Smith 435c838673bSBarry Smith M*/ 436c838673bSBarry Smith 437c838673bSBarry Smith /*MC 438c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 439c838673bSBarry Smith 440c838673bSBarry Smith Level: beginner 441c838673bSBarry Smith 442c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 443c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 444c838673bSBarry Smith 2-norm of the residual for right preconditioning 445c838673bSBarry Smith 446c838673bSBarry Smith Level: beginner 447c838673bSBarry Smith 448c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 449c838673bSBarry Smith 450c838673bSBarry Smith M*/ 451c838673bSBarry Smith 452c838673bSBarry Smith /*MC 453c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 454c838673bSBarry Smith 455c838673bSBarry Smith Level: beginner 456c838673bSBarry Smith 457c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 458c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 459c838673bSBarry Smith 2-norm of the residual for right preconditioning 460c838673bSBarry Smith 461c838673bSBarry Smith Level: beginner 462c838673bSBarry Smith 463c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 464c838673bSBarry Smith 465c838673bSBarry Smith M*/ 466c838673bSBarry Smith 467c838673bSBarry Smith /*MC 468c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 469c838673bSBarry Smith reached 470c838673bSBarry Smith 471c838673bSBarry Smith Level: beginner 472c838673bSBarry Smith 473c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 474c838673bSBarry Smith 475c838673bSBarry Smith M*/ 476c838673bSBarry Smith 477c838673bSBarry Smith /*MC 4788c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4798c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4804d0a8057SBarry Smith test routine is set in KSP. 481c838673bSBarry Smith 482c838673bSBarry Smith 483c838673bSBarry Smith Level: beginner 484c838673bSBarry Smith 485c838673bSBarry Smith 486c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 487c838673bSBarry Smith 488c838673bSBarry Smith M*/ 489c838673bSBarry Smith 490c838673bSBarry Smith /*MC 491c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4923014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4933014e516SBarry Smith preconditioner. 494c838673bSBarry Smith 495c838673bSBarry Smith Level: beginner 496c838673bSBarry Smith 497c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 498c838673bSBarry Smith 499c838673bSBarry Smith M*/ 500c838673bSBarry Smith 501c838673bSBarry Smith /*MC 502c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 503c838673bSBarry Smith method could not continue to enlarge the Krylov space. 504c838673bSBarry Smith 505c838673bSBarry Smith 506c838673bSBarry Smith Level: beginner 507c838673bSBarry Smith 508c838673bSBarry Smith 509c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 510c838673bSBarry Smith 511c838673bSBarry Smith M*/ 512c838673bSBarry Smith 513c838673bSBarry Smith /*MC 514c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 515c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 516c838673bSBarry Smith 517c838673bSBarry Smith Level: beginner 518c838673bSBarry Smith 519c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 520c838673bSBarry Smith 521c838673bSBarry Smith M*/ 522c838673bSBarry Smith 523c838673bSBarry Smith /*MC 524c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 525c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 526c838673bSBarry Smith be positive definite 527c838673bSBarry Smith 528c838673bSBarry Smith Level: beginner 529c838673bSBarry Smith 5302401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 531c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 532c838673bSBarry Smith 533c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 534c838673bSBarry Smith 535c838673bSBarry Smith M*/ 536c838673bSBarry Smith 537c838673bSBarry Smith /*MC 538c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 539c838673bSBarry Smith while the KSPSolve() is still running. 540c838673bSBarry Smith 541c838673bSBarry Smith Level: beginner 542c838673bSBarry Smith 543c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 544c838673bSBarry Smith 545c838673bSBarry Smith M*/ 546c838673bSBarry Smith 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 557abef13c0SSatish Balay 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 559d4fbbf0eSBarry Smith 56028ce4d24SBarry Smith /*E 56128ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 56228ce4d24SBarry Smith 56328ce4d24SBarry Smith Level: beginner 56428ce4d24SBarry Smith 56528ce4d24SBarry Smith .seealso: KSPCGSetType() 56628ce4d24SBarry Smith E*/ 5679dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5686a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 56928ce4d24SBarry Smith 570014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5728031f4adStmunson 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 576fcae7a14Stmunson 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 580e559a7feSLois Curfman McInnes 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5868031f4adStmunson 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5881d6018f0SLisandro Dalcin 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5913369ce9aSBarry Smith 5921d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 5931d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 5941d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 597014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 5992f2e5d10SKris Buschelman 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 60203919abeSBarry Smith 603f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 604ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 605f8a50e2bSBarry Smith 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 612f8a50e2bSBarry Smith 613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 616f8a50e2bSBarry Smith 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 619d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 6201324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 622014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 6243f22127dSBarry Smith 625014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat); 626fcfd50ebSBarry Smith 627014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 629014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 630014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 631014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 632fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 633014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 634fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 635014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 639fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 640fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6416c699258SBarry Smith 6422eac72dbSBarry Smith #endif 643