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 79140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList 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: 90f2ba6396SBarry 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*[]); 175*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 176*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 177*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 178*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1800a71e23eSMatthew Knepley 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1832eac72dbSBarry Smith 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 188ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1934b0e389bSBarry Smith 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1979f236934SBarry Smith 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2031d73ed98SBarry Smith 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2061d73ed98SBarry Smith 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 21058ad63f4SBarry Smith 211b49fd9e1SBarry Smith /*E 212b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 213b49fd9e1SBarry Smith 214b49fd9e1SBarry Smith Level: advanced 215b49fd9e1SBarry Smith 216bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 217e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 218b49fd9e1SBarry Smith 219b49fd9e1SBarry Smith E*/ 22078d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2216a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2221f7e983dSSatish Balay /*MC 2238c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2248c5b8ba0SBarry Smith 2258c5b8ba0SBarry Smith Level: advanced 2268c5b8ba0SBarry Smith 2278c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2288c5b8ba0SBarry Smith 229bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 230e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2318c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2328c5b8ba0SBarry Smith M*/ 2338c5b8ba0SBarry Smith 2341f7e983dSSatish Balay /*MC 2358c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2368c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2378c5b8ba0SBarry Smith poor orthogonality. 2388c5b8ba0SBarry Smith 2398c5b8ba0SBarry Smith Level: advanced 2408c5b8ba0SBarry Smith 2418c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2428c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2438c5b8ba0SBarry Smith 244bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 245e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2468c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2478c5b8ba0SBarry Smith M*/ 2488c5b8ba0SBarry Smith 2491f7e983dSSatish Balay /*MC 2508c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2518c5b8ba0SBarry Smith 2528c5b8ba0SBarry Smith Level: advanced 2538c5b8ba0SBarry Smith 2548c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2558c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2568c5b8ba0SBarry Smith 2578c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2588c5b8ba0SBarry Smith 259bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 260e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2618c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2628c5b8ba0SBarry Smith M*/ 2638c5b8ba0SBarry Smith 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 26608480c60SBarry Smith 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 270c38d4ed2SBarry Smith 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 274121fd945SKris Buschelman 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 278e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 279d9492815SBarry Smith 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2822eac72dbSBarry Smith 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 287390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 288390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 290816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 29784cb2905SBarry Smith 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 302c01c455dSBarry Smith 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 304014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3097287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3107287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3111eb62cbbSBarry Smith 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 313014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 314014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 315014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3161f7f0c4fSBarry Smith 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 31855849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 31955849f57SBarry Smith 32055849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3211eb62cbbSBarry Smith 322014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 323014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 324db9b2ab1SHong Zhang 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 326014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 32783ab6a24SBarry Smith 32828ce4d24SBarry Smith /*E 3298a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3308a4b9c5cSBarry Smith test routines. 3318a4b9c5cSBarry Smith 3328a4b9c5cSBarry Smith Level: advanced 3338a4b9c5cSBarry Smith 334a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3351718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 336a3f661c8SBarry Smith 3378a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3388a4b9c5cSBarry Smith 33994b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3401718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3418a4b9c5cSBarry Smith E*/ 3429e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3439e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 344014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 345a21b2a99SBarry Smith 3461f7e983dSSatish Balay /*MC 3479793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3488c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3498c5b8ba0SBarry Smith be based on a norm of a residual etc. 3508c5b8ba0SBarry Smith 3518c5b8ba0SBarry Smith Level: advanced 3528c5b8ba0SBarry Smith 353085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3548c5b8ba0SBarry Smith 355ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3568c5b8ba0SBarry Smith M*/ 3578c5b8ba0SBarry Smith 3581f7e983dSSatish Balay /*MC 359ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3608c5b8ba0SBarry Smith convergence test routine. 3618c5b8ba0SBarry Smith 3628c5b8ba0SBarry Smith Level: advanced 3638c5b8ba0SBarry Smith 3649793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3658c5b8ba0SBarry Smith M*/ 3668c5b8ba0SBarry Smith 3671f7e983dSSatish Balay /*MC 368ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3698c5b8ba0SBarry Smith convergence test routine. 3708c5b8ba0SBarry Smith 3718c5b8ba0SBarry Smith Level: advanced 3728c5b8ba0SBarry Smith 3739793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3748c5b8ba0SBarry Smith M*/ 3758c5b8ba0SBarry Smith 3761f7e983dSSatish Balay /*MC 377ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 378a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3798c5b8ba0SBarry Smith 3808c5b8ba0SBarry Smith Level: advanced 3818c5b8ba0SBarry Smith 3829793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3838c5b8ba0SBarry Smith M*/ 3848c5b8ba0SBarry Smith 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 389014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 3908a4b9c5cSBarry Smith 3918a4b9c5cSBarry Smith /*E 39228ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 39328ce4d24SBarry Smith have converged or diverged 39428ce4d24SBarry Smith 39528ce4d24SBarry Smith Level: beginner 39628ce4d24SBarry Smith 3974d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 39828ce4d24SBarry Smith 3994d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 4004d0a8057SBarry Smith 4014d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4024d0a8057SBarry Smith any of the values here also change them that array of names. 40386c02ca4SBarry Smith 404c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 40528ce4d24SBarry Smith E*/ 406d15094e1SBarry Smith typedef enum {/* converged */ 4079ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4089ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 409d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 410d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 411b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4128031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4138031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 414329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 415af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 416d15094e1SBarry Smith /* diverged */ 417b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 418d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 419d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 420d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 421b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 422b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 423b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4246aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 4256aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 426d15094e1SBarry Smith 427d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 428014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 429d15094e1SBarry Smith 430c838673bSBarry Smith /*MC 431c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 432c838673bSBarry Smith 433c838673bSBarry Smith Level: beginner 434c838673bSBarry Smith 435c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 436c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 437c838673bSBarry Smith 2-norm of the residual for right preconditioning 438c838673bSBarry Smith 439c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 440c838673bSBarry Smith 441c838673bSBarry Smith M*/ 442c838673bSBarry Smith 443c838673bSBarry Smith /*MC 444c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 445c838673bSBarry Smith 446c838673bSBarry Smith Level: beginner 447c838673bSBarry Smith 448c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 449c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 450c838673bSBarry Smith 2-norm of the residual for right preconditioning 451c838673bSBarry Smith 452c838673bSBarry Smith Level: beginner 453c838673bSBarry Smith 454c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 455c838673bSBarry Smith 456c838673bSBarry Smith M*/ 457c838673bSBarry Smith 458c838673bSBarry Smith /*MC 459c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 460c838673bSBarry Smith 461c838673bSBarry Smith Level: beginner 462c838673bSBarry Smith 463c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 464c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 465c838673bSBarry Smith 2-norm of the residual for right preconditioning 466c838673bSBarry Smith 467c838673bSBarry Smith Level: beginner 468c838673bSBarry Smith 469c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 470c838673bSBarry Smith 471c838673bSBarry Smith M*/ 472c838673bSBarry Smith 473c838673bSBarry Smith /*MC 474c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 475c838673bSBarry Smith reached 476c838673bSBarry Smith 477c838673bSBarry Smith Level: beginner 478c838673bSBarry Smith 479c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 480c838673bSBarry Smith 481c838673bSBarry Smith M*/ 482c838673bSBarry Smith 483c838673bSBarry Smith /*MC 4848c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4858c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4864d0a8057SBarry Smith test routine is set in KSP. 487c838673bSBarry Smith 488c838673bSBarry Smith 489c838673bSBarry Smith Level: beginner 490c838673bSBarry Smith 491c838673bSBarry Smith 492c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 493c838673bSBarry Smith 494c838673bSBarry Smith M*/ 495c838673bSBarry Smith 496c838673bSBarry Smith /*MC 497c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4983014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4993014e516SBarry Smith preconditioner. 500c838673bSBarry Smith 501c838673bSBarry Smith Level: beginner 502c838673bSBarry Smith 503c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 504c838673bSBarry Smith 505c838673bSBarry Smith M*/ 506c838673bSBarry Smith 507c838673bSBarry Smith /*MC 508c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 509c838673bSBarry Smith method could not continue to enlarge the Krylov space. 510c838673bSBarry Smith 511c838673bSBarry Smith 512c838673bSBarry Smith Level: beginner 513c838673bSBarry Smith 514c838673bSBarry Smith 515c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 516c838673bSBarry Smith 517c838673bSBarry Smith M*/ 518c838673bSBarry Smith 519c838673bSBarry Smith /*MC 520c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 521c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 522c838673bSBarry Smith 523c838673bSBarry Smith Level: beginner 524c838673bSBarry Smith 525c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 526c838673bSBarry Smith 527c838673bSBarry Smith M*/ 528c838673bSBarry Smith 529c838673bSBarry Smith /*MC 530c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 531c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 532c838673bSBarry Smith be positive definite 533c838673bSBarry Smith 534c838673bSBarry Smith Level: beginner 535c838673bSBarry Smith 5362401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 537c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 538c838673bSBarry Smith 539c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 540c838673bSBarry Smith 541c838673bSBarry Smith M*/ 542c838673bSBarry Smith 543c838673bSBarry Smith /*MC 544c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 545c838673bSBarry Smith while the KSPSolve() is still running. 546c838673bSBarry Smith 547c838673bSBarry Smith Level: beginner 548c838673bSBarry Smith 549c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 550c838673bSBarry Smith 551c838673bSBarry Smith M*/ 552c838673bSBarry Smith 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 563abef13c0SSatish Balay 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 565d4fbbf0eSBarry Smith 56628ce4d24SBarry Smith /*E 56728ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 56828ce4d24SBarry Smith 56928ce4d24SBarry Smith Level: beginner 57028ce4d24SBarry Smith 57128ce4d24SBarry Smith .seealso: KSPCGSetType() 57228ce4d24SBarry Smith E*/ 5739dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5746a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 57528ce4d24SBarry Smith 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5788031f4adStmunson 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 582fcae7a14Stmunson 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 586e559a7feSLois Curfman McInnes 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5928031f4adStmunson 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5941d6018f0SLisandro Dalcin 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5973369ce9aSBarry Smith 5981d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 5991d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 6001d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6052f2e5d10SKris Buschelman 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 60803919abeSBarry Smith 609f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 610ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 611f8a50e2bSBarry Smith 612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 618f8a50e2bSBarry Smith 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 622f8a50e2bSBarry Smith 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 625d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 6261324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat); 627014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 629014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 6303f22127dSBarry Smith 631014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat); 632fcfd50ebSBarry Smith 633014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 634014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 635014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 638fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 640fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 641014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 642014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 643014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 644014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 645fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 646fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6476c699258SBarry Smith 6482eac72dbSBarry Smith #endif 649