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 169*fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 170*fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 171*fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 172*fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 173*fa0ddf94SBarry Smith 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 179b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 180b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 181b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 182b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1840a71e23eSMatthew Knepley 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1872eac72dbSBarry Smith 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 192ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 1974b0e389bSBarry Smith 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 2019f236934SBarry Smith 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2071d73ed98SBarry Smith 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2101d73ed98SBarry Smith 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 21458ad63f4SBarry Smith 215b49fd9e1SBarry Smith /*E 216b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 217b49fd9e1SBarry Smith 218b49fd9e1SBarry Smith Level: advanced 219b49fd9e1SBarry Smith 220bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 221e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 222b49fd9e1SBarry Smith 223b49fd9e1SBarry Smith E*/ 22478d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2256a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2261f7e983dSSatish Balay /*MC 2278c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 2288c5b8ba0SBarry Smith 2298c5b8ba0SBarry Smith Level: advanced 2308c5b8ba0SBarry Smith 2318c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2328c5b8ba0SBarry Smith 233bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 234e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2358c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2368c5b8ba0SBarry Smith M*/ 2378c5b8ba0SBarry Smith 2381f7e983dSSatish Balay /*MC 2398c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2408c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2418c5b8ba0SBarry Smith poor orthogonality. 2428c5b8ba0SBarry Smith 2438c5b8ba0SBarry Smith Level: advanced 2448c5b8ba0SBarry Smith 2458c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2468c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2478c5b8ba0SBarry Smith 248bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 249e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2508c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2518c5b8ba0SBarry Smith M*/ 2528c5b8ba0SBarry Smith 2531f7e983dSSatish Balay /*MC 2548c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2558c5b8ba0SBarry Smith 2568c5b8ba0SBarry Smith Level: advanced 2578c5b8ba0SBarry Smith 2588c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2598c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2608c5b8ba0SBarry Smith 2618c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2628c5b8ba0SBarry Smith 263bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 264e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2658c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2668c5b8ba0SBarry Smith M*/ 2678c5b8ba0SBarry Smith 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 27008480c60SBarry Smith 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 274c38d4ed2SBarry Smith 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 278121fd945SKris Buschelman 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 282e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 283d9492815SBarry Smith 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2862eac72dbSBarry Smith 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 291390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 292390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 294816e40a9SMark F. Adams PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 30184cb2905SBarry Smith 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 304c01c455dSBarry Smith 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3117287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3127287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3131eb62cbbSBarry Smith 314014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 315014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 316014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3181f7f0c4fSBarry Smith 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 32055849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 32155849f57SBarry Smith 32255849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3231eb62cbbSBarry Smith 324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 326db9b2ab1SHong Zhang 327014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 32983ab6a24SBarry Smith 33028ce4d24SBarry Smith /*E 3318a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3328a4b9c5cSBarry Smith test routines. 3338a4b9c5cSBarry Smith 3348a4b9c5cSBarry Smith Level: advanced 3358a4b9c5cSBarry Smith 336a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3371718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 338a3f661c8SBarry Smith 3398a4b9c5cSBarry Smith Notes: this must match finclude/petscksp.h 3408a4b9c5cSBarry Smith 34194b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3421718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3438a4b9c5cSBarry Smith E*/ 3449e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3459e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 346014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 347a21b2a99SBarry Smith 3481f7e983dSSatish Balay /*MC 3499793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3508c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3518c5b8ba0SBarry Smith be based on a norm of a residual etc. 3528c5b8ba0SBarry Smith 3538c5b8ba0SBarry Smith Level: advanced 3548c5b8ba0SBarry Smith 355085a36d4SBarry Smith Note: Some Krylov methods need to compute a residual norm and then this option is ignored 3568c5b8ba0SBarry Smith 357ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3588c5b8ba0SBarry Smith M*/ 3598c5b8ba0SBarry Smith 3601f7e983dSSatish Balay /*MC 361ce9499c7SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 3628c5b8ba0SBarry Smith convergence test routine. 3638c5b8ba0SBarry Smith 3648c5b8ba0SBarry Smith Level: advanced 3658c5b8ba0SBarry Smith 3669793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3678c5b8ba0SBarry Smith M*/ 3688c5b8ba0SBarry Smith 3691f7e983dSSatish Balay /*MC 370ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3718c5b8ba0SBarry Smith convergence test routine. 3728c5b8ba0SBarry Smith 3738c5b8ba0SBarry Smith Level: advanced 3748c5b8ba0SBarry Smith 3759793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3768c5b8ba0SBarry Smith M*/ 3778c5b8ba0SBarry Smith 3781f7e983dSSatish Balay /*MC 379ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 380a3f661c8SBarry Smith convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 3818c5b8ba0SBarry Smith 3828c5b8ba0SBarry Smith Level: advanced 3838c5b8ba0SBarry Smith 3849793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3858c5b8ba0SBarry Smith M*/ 3868c5b8ba0SBarry Smith 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 389014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 391014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3928a4b9c5cSBarry Smith 393b30237c6SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 394b30237c6SBarry Smith 3958a4b9c5cSBarry Smith /*E 39628ce4d24SBarry Smith KSPConvergedReason - reason a Krylov method was said to 39728ce4d24SBarry Smith have converged or diverged 39828ce4d24SBarry Smith 39928ce4d24SBarry Smith Level: beginner 40028ce4d24SBarry Smith 4014d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 40228ce4d24SBarry Smith 4034d0a8057SBarry Smith Developer notes: this must match finclude/petscksp.h 4044d0a8057SBarry Smith 4054d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4064d0a8057SBarry Smith any of the values here also change them that array of names. 40786c02ca4SBarry Smith 408c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 40928ce4d24SBarry Smith E*/ 410d15094e1SBarry Smith typedef enum {/* converged */ 4119ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4129ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 413d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 414d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 415b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4168031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4178031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 418329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 419af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 420d15094e1SBarry Smith /* diverged */ 421b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 422d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 423d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 424d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 425b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 426b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 427b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4284d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4296aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 430d15094e1SBarry Smith 431d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 432014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 433d15094e1SBarry Smith 434c838673bSBarry Smith /*MC 435c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 436c838673bSBarry Smith 437c838673bSBarry Smith Level: beginner 438c838673bSBarry Smith 439c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 440c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 441c838673bSBarry Smith 2-norm of the residual for right preconditioning 442c838673bSBarry Smith 443c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 444c838673bSBarry Smith 445c838673bSBarry Smith M*/ 446c838673bSBarry Smith 447c838673bSBarry Smith /*MC 448c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 449c838673bSBarry Smith 450c838673bSBarry Smith Level: beginner 451c838673bSBarry Smith 452c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 453c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 454c838673bSBarry Smith 2-norm of the residual for right preconditioning 455c838673bSBarry Smith 456c838673bSBarry Smith Level: beginner 457c838673bSBarry Smith 458c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 459c838673bSBarry Smith 460c838673bSBarry Smith M*/ 461c838673bSBarry Smith 462c838673bSBarry Smith /*MC 463c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 464c838673bSBarry Smith 465c838673bSBarry Smith Level: beginner 466c838673bSBarry Smith 467c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 468c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 469c838673bSBarry Smith 2-norm of the residual for right preconditioning 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 478c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 479c838673bSBarry Smith reached 480c838673bSBarry Smith 481c838673bSBarry Smith Level: beginner 482c838673bSBarry Smith 483c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 484c838673bSBarry Smith 485c838673bSBarry Smith M*/ 486c838673bSBarry Smith 487c838673bSBarry Smith /*MC 4888c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4898c9406f6SLisandro Dalcin the preconditioner is applied. Also used when the KSPSkipConverged() convergence 4904d0a8057SBarry Smith test routine is set in KSP. 491c838673bSBarry Smith 492c838673bSBarry Smith 493c838673bSBarry Smith Level: beginner 494c838673bSBarry Smith 495c838673bSBarry Smith 496c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 497c838673bSBarry Smith 498c838673bSBarry Smith M*/ 499c838673bSBarry Smith 500c838673bSBarry Smith /*MC 501c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 5023014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 5033014e516SBarry Smith preconditioner. 504c838673bSBarry Smith 505c838673bSBarry Smith Level: beginner 506c838673bSBarry Smith 507c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 508c838673bSBarry Smith 509c838673bSBarry Smith M*/ 510c838673bSBarry Smith 511c838673bSBarry Smith /*MC 512c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 513c838673bSBarry Smith method could not continue to enlarge the Krylov space. 514c838673bSBarry Smith 515c838673bSBarry Smith 516c838673bSBarry Smith Level: beginner 517c838673bSBarry Smith 518c838673bSBarry Smith 519c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 520c838673bSBarry Smith 521c838673bSBarry Smith M*/ 522c838673bSBarry Smith 523c838673bSBarry Smith /*MC 524c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 525c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 526c838673bSBarry Smith 527c838673bSBarry Smith Level: beginner 528c838673bSBarry Smith 529c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 530c838673bSBarry Smith 531c838673bSBarry Smith M*/ 532c838673bSBarry Smith 533c838673bSBarry Smith /*MC 534c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 535c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 536c838673bSBarry Smith be positive definite 537c838673bSBarry Smith 538c838673bSBarry Smith Level: beginner 539c838673bSBarry Smith 5402401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 541c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 542c838673bSBarry Smith 543c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 544c838673bSBarry Smith 545c838673bSBarry Smith M*/ 546c838673bSBarry Smith 547c838673bSBarry Smith /*MC 548c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 549c838673bSBarry Smith while the KSPSolve() is still running. 550c838673bSBarry Smith 551c838673bSBarry Smith Level: beginner 552c838673bSBarry Smith 553c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 554c838673bSBarry Smith 555c838673bSBarry Smith M*/ 556c838673bSBarry Smith 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *); 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **); 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 565014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 567abef13c0SSatish Balay 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 569d4fbbf0eSBarry Smith 57028ce4d24SBarry Smith /*E 57128ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 57228ce4d24SBarry Smith 57328ce4d24SBarry Smith Level: beginner 57428ce4d24SBarry Smith 57528ce4d24SBarry Smith .seealso: KSPCGSetType() 57628ce4d24SBarry Smith E*/ 5779dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 5786a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 57928ce4d24SBarry Smith 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 5828031f4adStmunson 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 586fcae7a14Stmunson 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 590e559a7feSLois Curfman McInnes 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 5968031f4adStmunson 597014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 5981d6018f0SLisandro Dalcin 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 6013369ce9aSBarry Smith 6029804daf3SBarry Smith #include <petscdrawtypes.h> 6031d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 6041d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 6051d1e9da1SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 606014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6102f2e5d10SKris Buschelman 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 612014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 61303919abeSBarry Smith 614f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 615ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 616f8a50e2bSBarry Smith 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 622014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 623f8a50e2bSBarry Smith 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 626014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 627f8a50e2bSBarry Smith 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 629014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 630d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 6311324cf18SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat); 632014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 633014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 634014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 6353f22127dSBarry Smith 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 641fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 643fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 644014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 645014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 646014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 647014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 648fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 649fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6506c699258SBarry Smith 6512eac72dbSBarry Smith #endif 652