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: 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); 184*ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1894b0e389bSBarry Smith 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 1939f236934SBarry Smith 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 1991d73ed98SBarry Smith 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2021d73ed98SBarry Smith 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 20658ad63f4SBarry Smith 207b49fd9e1SBarry Smith /*E 208b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 209b49fd9e1SBarry Smith 210b49fd9e1SBarry Smith Level: advanced 211b49fd9e1SBarry Smith 212bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 213e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 214b49fd9e1SBarry Smith 215b49fd9e1SBarry Smith E*/ 21678d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2176a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2181f7e983dSSatish Balay /*MC 2198c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2208c5b8ba0SBarry Smith 2218c5b8ba0SBarry Smith Level: advanced 2228c5b8ba0SBarry Smith 2238c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2248c5b8ba0SBarry Smith 225bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 226e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2278c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2288c5b8ba0SBarry Smith M*/ 2298c5b8ba0SBarry Smith 2301f7e983dSSatish Balay /*MC 2318c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2328c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2338c5b8ba0SBarry Smith poor orthogonality. 2348c5b8ba0SBarry Smith 2358c5b8ba0SBarry Smith Level: advanced 2368c5b8ba0SBarry Smith 2378c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2388c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2398c5b8ba0SBarry Smith 240bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 241e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2428c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2438c5b8ba0SBarry Smith M*/ 2448c5b8ba0SBarry Smith 2451f7e983dSSatish Balay /*MC 2468c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2478c5b8ba0SBarry Smith 2488c5b8ba0SBarry Smith Level: advanced 2498c5b8ba0SBarry Smith 2508c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2518c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2528c5b8ba0SBarry Smith 2538c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2548c5b8ba0SBarry Smith 255bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 256e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2578c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2588c5b8ba0SBarry Smith M*/ 2598c5b8ba0SBarry Smith 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 26208480c60SBarry Smith 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 266c38d4ed2SBarry Smith 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 270121fd945SKris Buschelman 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 274d9492815SBarry Smith 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2772eac72dbSBarry Smith 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 282390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 283390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 285816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 29284cb2905SBarry Smith 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 297c01c455dSBarry Smith 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3047287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3057287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3061eb62cbbSBarry Smith 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3111f7f0c4fSBarry Smith 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 31355849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 31455849f57SBarry Smith 31555849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3161eb62cbbSBarry Smith 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 319db9b2ab1SHong Zhang 320014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 321014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 32283ab6a24SBarry Smith 32328ce4d24SBarry Smith /*E 3248a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3258a4b9c5cSBarry Smith test routines. 3268a4b9c5cSBarry Smith 3278a4b9c5cSBarry Smith Level: advanced 3288a4b9c5cSBarry Smith 329a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3301718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 331a3f661c8SBarry Smith 3328a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3338a4b9c5cSBarry Smith 33494b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3351718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3368a4b9c5cSBarry Smith E*/ 3379e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3389e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 339014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 340a21b2a99SBarry Smith 3411f7e983dSSatish Balay /*MC 3429793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3438c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3448c5b8ba0SBarry Smith be based on a norm of a residual etc. 3458c5b8ba0SBarry Smith 3468c5b8ba0SBarry Smith Level: advanced 3478c5b8ba0SBarry Smith 348085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3498c5b8ba0SBarry Smith 350ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3518c5b8ba0SBarry Smith M*/ 3528c5b8ba0SBarry Smith 3531f7e983dSSatish Balay /*MC 354ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3558c5b8ba0SBarry Smith convergence test routine. 3568c5b8ba0SBarry Smith 3578c5b8ba0SBarry Smith Level: advanced 3588c5b8ba0SBarry Smith 3599793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3608c5b8ba0SBarry Smith M*/ 3618c5b8ba0SBarry Smith 3621f7e983dSSatish Balay /*MC 363ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3648c5b8ba0SBarry Smith convergence test routine. 3658c5b8ba0SBarry Smith 3668c5b8ba0SBarry Smith Level: advanced 3678c5b8ba0SBarry Smith 3689793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3698c5b8ba0SBarry Smith M*/ 3708c5b8ba0SBarry Smith 3711f7e983dSSatish Balay /*MC 372ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 373a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3748c5b8ba0SBarry Smith 3758c5b8ba0SBarry Smith Level: advanced 3768c5b8ba0SBarry Smith 3779793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3788c5b8ba0SBarry Smith M*/ 3798c5b8ba0SBarry Smith 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 3858a4b9c5cSBarry Smith 3868a4b9c5cSBarry Smith /*E 38728ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 38828ce4d24SBarry Smith have converged or diverged 38928ce4d24SBarry Smith 39028ce4d24SBarry Smith Level: beginner 39128ce4d24SBarry Smith 3924d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 39328ce4d24SBarry Smith 3944d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 3954d0a8057SBarry Smith 3964d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 3974d0a8057SBarry Smith any of the values here also change them that array of names. 39886c02ca4SBarry Smith 399c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 40028ce4d24SBarry Smith E*/ 401d15094e1SBarry Smith typedef enum {/* converged */ 4029ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4039ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 404d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 405d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 406b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4078031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4088031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 409329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 410af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 411d15094e1SBarry Smith /* diverged */ 412b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 413d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 414d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 415d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 416b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 417b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 418b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4196aee1118SBarry Smith KSP_DIVERGED_NAN = -9, 4206aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 421d15094e1SBarry Smith 422d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 423014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 424d15094e1SBarry Smith 425c838673bSBarry Smith /*MC 426c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 427c838673bSBarry Smith 428c838673bSBarry Smith Level: beginner 429c838673bSBarry Smith 430c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 431c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 432c838673bSBarry Smith 2-norm of the residual for right preconditioning 433c838673bSBarry Smith 434c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 435c838673bSBarry Smith 436c838673bSBarry Smith M*/ 437c838673bSBarry Smith 438c838673bSBarry Smith /*MC 439c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 440c838673bSBarry Smith 441c838673bSBarry Smith Level: beginner 442c838673bSBarry Smith 443c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 444c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 445c838673bSBarry Smith 2-norm of the residual for right preconditioning 446c838673bSBarry Smith 447c838673bSBarry Smith Level: beginner 448c838673bSBarry Smith 449c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 450c838673bSBarry Smith 451c838673bSBarry Smith M*/ 452c838673bSBarry Smith 453c838673bSBarry Smith /*MC 454c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 455c838673bSBarry Smith 456c838673bSBarry Smith Level: beginner 457c838673bSBarry Smith 458c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 459c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 460c838673bSBarry Smith 2-norm of the residual for right preconditioning 461c838673bSBarry Smith 462c838673bSBarry Smith Level: beginner 463c838673bSBarry Smith 464c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 465c838673bSBarry Smith 466c838673bSBarry Smith M*/ 467c838673bSBarry Smith 468c838673bSBarry Smith /*MC 469c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 470c838673bSBarry Smith reached 471c838673bSBarry Smith 472c838673bSBarry Smith Level: beginner 473c838673bSBarry Smith 474c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 475c838673bSBarry Smith 476c838673bSBarry Smith M*/ 477c838673bSBarry Smith 478c838673bSBarry Smith /*MC 4798c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4808c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4814d0a8057SBarry Smith test routine is set in KSP. 482c838673bSBarry Smith 483c838673bSBarry Smith 484c838673bSBarry Smith Level: beginner 485c838673bSBarry Smith 486c838673bSBarry Smith 487c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 488c838673bSBarry Smith 489c838673bSBarry Smith M*/ 490c838673bSBarry Smith 491c838673bSBarry Smith /*MC 492c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 4933014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 4943014e516SBarry Smith preconditioner. 495c838673bSBarry Smith 496c838673bSBarry Smith Level: beginner 497c838673bSBarry Smith 498c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 499c838673bSBarry Smith 500c838673bSBarry Smith M*/ 501c838673bSBarry Smith 502c838673bSBarry Smith /*MC 503c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 504c838673bSBarry Smith method could not continue to enlarge the Krylov space. 505c838673bSBarry Smith 506c838673bSBarry Smith 507c838673bSBarry Smith Level: beginner 508c838673bSBarry Smith 509c838673bSBarry Smith 510c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 511c838673bSBarry Smith 512c838673bSBarry Smith M*/ 513c838673bSBarry Smith 514c838673bSBarry Smith /*MC 515c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 516c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 517c838673bSBarry Smith 518c838673bSBarry Smith Level: beginner 519c838673bSBarry Smith 520c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 521c838673bSBarry Smith 522c838673bSBarry Smith M*/ 523c838673bSBarry Smith 524c838673bSBarry Smith /*MC 525c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 526c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 527c838673bSBarry Smith be positive definite 528c838673bSBarry Smith 529c838673bSBarry Smith Level: beginner 530c838673bSBarry Smith 5312401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 532c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 533c838673bSBarry Smith 534c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 535c838673bSBarry Smith 536c838673bSBarry Smith M*/ 537c838673bSBarry Smith 538c838673bSBarry Smith /*MC 539c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 540c838673bSBarry Smith while the KSPSolve() is still running. 541c838673bSBarry Smith 542c838673bSBarry Smith Level: beginner 543c838673bSBarry Smith 544c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 545c838673bSBarry Smith 546c838673bSBarry Smith M*/ 547c838673bSBarry Smith 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 558abef13c0SSatish Balay 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 560d4fbbf0eSBarry Smith 56128ce4d24SBarry Smith /*E 56228ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 56328ce4d24SBarry Smith 56428ce4d24SBarry Smith Level: beginner 56528ce4d24SBarry Smith 56628ce4d24SBarry Smith .seealso: KSPCGSetType() 56728ce4d24SBarry Smith E*/ 5689dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5696a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 57028ce4d24SBarry Smith 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5738031f4adStmunson 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 577fcae7a14Stmunson 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 581e559a7feSLois Curfman McInnes 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5878031f4adStmunson 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5891d6018f0SLisandro Dalcin 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 5923369ce9aSBarry Smith 5931d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 5941d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 5951d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 597014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6002f2e5d10SKris Buschelman 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 60303919abeSBarry Smith 604f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 605ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 606f8a50e2bSBarry Smith 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 613f8a50e2bSBarry Smith 614014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 617f8a50e2bSBarry Smith 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 620d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 6211324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat); 622014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 6253f22127dSBarry Smith 626014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat); 627fcfd50ebSBarry Smith 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 629014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 630014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 631014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 632014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 633fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 634014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 635fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 640fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 641fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6426c699258SBarry Smith 6432eac72dbSBarry Smith #endif 644