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 8607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 91dbb0a54SBarry Smith 1028ce4d24SBarry Smith /*S 118f6c3df8SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 128f6c3df8SBarry Smith linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 1328ce4d24SBarry Smith 1428ce4d24SBarry Smith Level: beginner 1528ce4d24SBarry Smith 1628ce4d24SBarry Smith Concepts: Krylov methods 1728ce4d24SBarry Smith 1895452b02SPatrick Sanan Notes: 1995452b02SPatrick Sanan When a direct solver is used but no Krylov solver is used the KSP object is still used by with a 208f6c3df8SBarry Smith KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver). 218f6c3df8SBarry Smith 22bca11509SBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES 2328ce4d24SBarry Smith S*/ 24e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 252eac72dbSBarry Smith 2676bdecfbSBarry Smith /*J 278f6c3df8SBarry Smith KSPType - String with the name of a PETSc Krylov method. 2828ce4d24SBarry Smith 2928ce4d24SBarry Smith Level: beginner 3028ce4d24SBarry Smith 318f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 3276bdecfbSBarry Smith J*/ 3319fd82e9SBarry Smith typedef const char* KSPType; 3482bf6240SBarry Smith #define KSPRICHARDSON "richardson" 356c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3682bf6240SBarry Smith #define KSPCG "cg" 372c8fc646SJed Brown #define KSPGROPPCG "groppcg" 382c8fc646SJed Brown #define KSPPIPECG "pipecg" 39901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 40265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 41df2a969aSvictorle #define KSPCGNE "cgne" 421155b781STodd Munson #define KSPCGNASH "nash" 431155b781STodd Munson #define KSPCGSTCG "stcg" 441155b781STodd Munson #define KSPCGGLTR "gltr" 45640f4f53SPatrick Sanan #define KSPFCG "fcg" 46390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 4782bf6240SBarry Smith #define KSPGMRES "gmres" 48483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 499dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 509dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 514f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5261b468f9SJed Brown #define KSPPGMRES "pgmres" 5382bf6240SBarry Smith #define KSPTCQMR "tcqmr" 5482bf6240SBarry Smith #define KSPBCGS "bcgs" 55e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5618ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 57c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 58f0037002Svictorle #define KSPBCGSL "bcgsl" 59f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 6082bf6240SBarry Smith #define KSPCGS "cgs" 6182bf6240SBarry Smith #define KSPTFQMR "tfqmr" 6282bf6240SBarry Smith #define KSPCR "cr" 632c8fc646SJed Brown #define KSPPIPECR "pipecr" 6482bf6240SBarry Smith #define KSPLSQR "lsqr" 6582bf6240SBarry Smith #define KSPPREONLY "preonly" 6682bf6240SBarry Smith #define KSPQCG "qcg" 67c9cf9b11SBarry Smith #define KSPBICG "bicg" 68b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6901c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 7021356919SSatish Balay #define KSPLCD "lcd" 711d6018f0SLisandro Dalcin #define KSPPYTHON "python" 7258ad63f4SBarry Smith #define KSPGCR "gcr" 73fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 74e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 75e4d80e07Szianekhodja #define KSPCGLS "cgls" 76329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 772eac72dbSBarry Smith 788ba1e511SMatthew Knepley /* Logging support */ 79014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 80ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 815358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 828ba1e511SMatthew Knepley 83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*); 8419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 85c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*); 86014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 88014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 9223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 9325c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 942eac72dbSBarry Smith 95140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 96ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 97bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 9830de9b25SBarry Smith 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 102c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool); 1087d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool); 109c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*); 110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool); 111c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*); 112014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*); 113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*); 114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 116734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1172a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1182a7a6963SBarry Smith PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1192eac72dbSBarry Smith 120d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 121d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 122d4211eb9SBarry Smith 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 125aabeff55SBarry Smith 126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**)); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1324b0e389bSBarry Smith 133fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 134fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*); 135fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 136fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 137fa0ddf94SBarry Smith 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 143b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 144b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 145b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 146b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*); 1480a71e23eSMatthew Knepley 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*); 1512eac72dbSBarry Smith 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 15558450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 156b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 15758450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 159d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*); 160d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1617d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]); 1624b0e389bSBarry Smith 163640f4f53SPatrick Sanan /*E 164640f4f53SPatrick Sanan 16506137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 166640f4f53SPatrick Sanan 16793f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 16893f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 169640f4f53SPatrick Sanan 1702b26979fSBarry Smith Level: intermediate 17106137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 172640f4f53SPatrick Sanan 173640f4f53SPatrick Sanan E*/ 17493f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 17506137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 176640f4f53SPatrick Sanan 177640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 178640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 179640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 180640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 18106137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 18206137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 183640f4f53SPatrick Sanan 184390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 185390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 186390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 187390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 188390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 189390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 190390d8e47SPatrick Sanan 191fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 192fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 193fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 194fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 195fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 196fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 197fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 198fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 19983f0b094SBarry Smith 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 2039f236934SBarry Smith 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2091d73ed98SBarry Smith 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2121d73ed98SBarry Smith 213483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 214483d6965SPatrick Sanan 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 21858ad63f4SBarry Smith 21904d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*); 22004d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC); 22104d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*); 222568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat); 223b49fd9e1SBarry Smith /*E 224b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 225b49fd9e1SBarry Smith 226b49fd9e1SBarry Smith Level: advanced 227b49fd9e1SBarry Smith 228bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 229e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 230b49fd9e1SBarry Smith 231b49fd9e1SBarry Smith E*/ 23278d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2336a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2341f7e983dSSatish Balay /*MC 2351957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2368c5b8ba0SBarry Smith 2378c5b8ba0SBarry Smith Level: advanced 2388c5b8ba0SBarry Smith 2398c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2408c5b8ba0SBarry Smith 241bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 242e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2438c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2448c5b8ba0SBarry Smith M*/ 2458c5b8ba0SBarry Smith 2461f7e983dSSatish Balay /*MC 2478c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2488c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2498c5b8ba0SBarry Smith poor orthogonality. 2508c5b8ba0SBarry Smith 2518c5b8ba0SBarry Smith Level: advanced 2528c5b8ba0SBarry Smith 2538c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2548c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2558c5b8ba0SBarry Smith 256bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 257e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2588c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2598c5b8ba0SBarry Smith M*/ 2608c5b8ba0SBarry Smith 2611f7e983dSSatish Balay /*MC 2628c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2638c5b8ba0SBarry Smith 2648c5b8ba0SBarry Smith Level: advanced 2658c5b8ba0SBarry Smith 2668c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2678c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2688c5b8ba0SBarry Smith 2698c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2708c5b8ba0SBarry Smith 271bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 272e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2738c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2748c5b8ba0SBarry Smith M*/ 2758c5b8ba0SBarry Smith 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 27808480c60SBarry Smith 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 282c38d4ed2SBarry Smith 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 286121fd945SKris Buschelman 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 290e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 291d9492815SBarry Smith 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2942eac72dbSBarry Smith 295fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 296fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 297fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 298fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 299390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 300390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 301fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 302fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 303fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 304fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 305e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 306e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 307e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*); 309fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*)); 31084cb2905SBarry Smith 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 313c01c455dSBarry Smith 31423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 31523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 316014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*); 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3207287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3217287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3221eb62cbbSBarry Smith 323014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*); 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 326014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*); 3271f7f0c4fSBarry Smith 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 32955849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 330685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 3312a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 3322a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 33355849f57SBarry Smith 33455849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3351eb62cbbSBarry Smith 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 337014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 33863aaf0cbSSatish Balay PETSC_EXTERN PetscErrorCode KSPLSQRGetArnorm(KSP,PetscReal*,PetscReal*,PetscReal*); 339db9b2ab1SHong Zhang 340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 341014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 34268ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 34383ab6a24SBarry Smith 34428ce4d24SBarry Smith /*E 3458a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3468a4b9c5cSBarry Smith test routines. 3478a4b9c5cSBarry Smith 3488a4b9c5cSBarry Smith Level: advanced 3498a4b9c5cSBarry Smith 350a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3511718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 352a3f661c8SBarry Smith 35395452b02SPatrick Sanan Notes: 35495452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 3558a4b9c5cSBarry Smith 35694b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3571718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3588a4b9c5cSBarry Smith E*/ 3599e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3609e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 361014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 362a21b2a99SBarry Smith 3631f7e983dSSatish Balay /*MC 3649793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3658c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3668c5b8ba0SBarry Smith be based on a norm of a residual etc. 3678c5b8ba0SBarry Smith 3688c5b8ba0SBarry Smith Level: advanced 3698c5b8ba0SBarry Smith 3701957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3718c5b8ba0SBarry Smith 372ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3738c5b8ba0SBarry Smith M*/ 3748c5b8ba0SBarry Smith 3751f7e983dSSatish Balay /*MC 3761957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3778c5b8ba0SBarry Smith convergence test routine. 3788c5b8ba0SBarry Smith 3798c5b8ba0SBarry Smith Level: advanced 3808c5b8ba0SBarry Smith 3819793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3828c5b8ba0SBarry Smith M*/ 3838c5b8ba0SBarry Smith 3841f7e983dSSatish Balay /*MC 385ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3868c5b8ba0SBarry Smith convergence test routine. 3878c5b8ba0SBarry Smith 3888c5b8ba0SBarry Smith Level: advanced 3898c5b8ba0SBarry Smith 3909793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3918c5b8ba0SBarry Smith M*/ 3928c5b8ba0SBarry Smith 3931f7e983dSSatish Balay /*MC 394ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 395390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 3968c5b8ba0SBarry Smith 3978c5b8ba0SBarry Smith Level: advanced 3988c5b8ba0SBarry Smith 3999793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 4008c5b8ba0SBarry Smith M*/ 4018c5b8ba0SBarry Smith 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 4078a4b9c5cSBarry Smith 4088a4b9c5cSBarry Smith /*E 4091957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 41028ce4d24SBarry Smith 41128ce4d24SBarry Smith Level: beginner 41228ce4d24SBarry Smith 41395452b02SPatrick Sanan Notes: 41495452b02SPatrick Sanan See KSPGetConvergedReason() for explanation of each value 41528ce4d24SBarry Smith 41695452b02SPatrick Sanan Developer Notes: 41795452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 4184d0a8057SBarry Smith 4194d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4204d0a8057SBarry Smith any of the values here also change them that array of names. 42186c02ca4SBarry Smith 422c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 42328ce4d24SBarry Smith E*/ 424d15094e1SBarry Smith typedef enum {/* converged */ 4259ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4269ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 427d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 428d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 429b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4308031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4318031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 432329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 433af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 434d15094e1SBarry Smith /* diverged */ 435b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 436d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 437d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 438d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 439b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 440b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 441b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4424d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4436aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 444422a814eSBarry Smith KSP_DIVERGED_PCSETUP_FAILED = -11, 445d15094e1SBarry Smith 446d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 447014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 448d15094e1SBarry Smith 449c838673bSBarry Smith /*MC 450f9fed41fSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called 451c838673bSBarry Smith 452c838673bSBarry Smith Level: beginner 453c838673bSBarry Smith 454c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 455c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 456c838673bSBarry Smith 2-norm of the residual for right preconditioning 457c838673bSBarry Smith 458f9fed41fSBarry Smith See also KSP_CONVERGED_ATOL which may apply before this tolerance. 459f9fed41fSBarry Smith 460f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 461c838673bSBarry Smith 462c838673bSBarry Smith M*/ 463c838673bSBarry Smith 464c838673bSBarry Smith /*MC 465c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 466c838673bSBarry Smith 467c838673bSBarry Smith Level: beginner 468c838673bSBarry Smith 469c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 470c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 471c838673bSBarry Smith 2-norm of the residual for right preconditioning 472c838673bSBarry Smith 473f9fed41fSBarry Smith See also KSP_CONVERGED_RTOL which may apply before this tolerance. 474c838673bSBarry Smith 475f9fed41fSBarry Smith .seealso: KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 476c838673bSBarry Smith 477c838673bSBarry Smith M*/ 478c838673bSBarry Smith 479c838673bSBarry Smith /*MC 480c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 481c838673bSBarry Smith 482c838673bSBarry Smith Level: beginner 483c838673bSBarry Smith 484c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 485c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 486c838673bSBarry Smith 2-norm of the residual for right preconditioning 487c838673bSBarry Smith 488c838673bSBarry Smith Level: beginner 489c838673bSBarry Smith 490f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 491c838673bSBarry Smith 492c838673bSBarry Smith M*/ 493c838673bSBarry Smith 494c838673bSBarry Smith /*MC 495c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 496c838673bSBarry Smith reached 497c838673bSBarry Smith 498c838673bSBarry Smith Level: beginner 499c838673bSBarry Smith 500c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 501c838673bSBarry Smith 502c838673bSBarry Smith M*/ 503c838673bSBarry Smith 504c838673bSBarry Smith /*MC 5058c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 5060059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 5074d0a8057SBarry Smith test routine is set in KSP. 508c838673bSBarry Smith 509c838673bSBarry Smith Level: beginner 510c838673bSBarry Smith 511c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 512c838673bSBarry Smith 513c838673bSBarry Smith M*/ 514c838673bSBarry Smith 515c838673bSBarry Smith /*MC 516c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 5173014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 5183014e516SBarry Smith preconditioner. 519c838673bSBarry Smith 520c838673bSBarry Smith Level: beginner 521c838673bSBarry Smith 522c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 523c838673bSBarry Smith 524c838673bSBarry Smith M*/ 525c838673bSBarry Smith 526c838673bSBarry Smith /*MC 527c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 528c838673bSBarry Smith method could not continue to enlarge the Krylov space. 529c838673bSBarry Smith 530c838673bSBarry Smith Level: beginner 531c838673bSBarry Smith 532c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 533c838673bSBarry Smith 534c838673bSBarry Smith M*/ 535c838673bSBarry Smith 536c838673bSBarry Smith /*MC 537c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 538c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 539c838673bSBarry Smith 540c838673bSBarry Smith Level: beginner 541c838673bSBarry Smith 542c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 543c838673bSBarry Smith 544c838673bSBarry Smith M*/ 545c838673bSBarry Smith 546c838673bSBarry Smith /*MC 547c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 548c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 549c838673bSBarry Smith be positive definite 550c838673bSBarry Smith 551c838673bSBarry Smith Level: beginner 552c838673bSBarry Smith 55395452b02SPatrick Sanan Notes: 55495452b02SPatrick Sanan This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 555c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 556c838673bSBarry Smith 557c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 558c838673bSBarry Smith 559c838673bSBarry Smith M*/ 560c838673bSBarry Smith 561c838673bSBarry Smith /*MC 5629fc87aa7SBarry Smith KSP_DIVERGED_PCSETUP_FAILED - It was not possible to build the requested preconditioner. This is usually due to a 5639fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 5649fc87aa7SBarry Smith such as PCFIELDSPLIT. 5659fc87aa7SBarry Smith 5669fc87aa7SBarry Smith Level: beginner 5679fc87aa7SBarry Smith 56895452b02SPatrick Sanan Notes: 56995452b02SPatrick Sanan Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 5709fc87aa7SBarry Smith 5719fc87aa7SBarry Smith 5729fc87aa7SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 5739fc87aa7SBarry Smith 5749fc87aa7SBarry Smith M*/ 5759fc87aa7SBarry Smith 5769fc87aa7SBarry Smith /*MC 577c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 578c838673bSBarry Smith while the KSPSolve() is still running. 579c838673bSBarry Smith 580c838673bSBarry Smith Level: beginner 581c838673bSBarry Smith 582c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 583c838673bSBarry Smith 584c838673bSBarry Smith M*/ 585c838673bSBarry Smith 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*)); 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**); 5888de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 5908de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*); 5918de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**); 5928de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5938de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5940059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*); 596abef13c0SSatish Balay 5978ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5988ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5998ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 6008ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 6018ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 6028ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 6038ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 6048ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 6058ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 6068ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 6078ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 6088ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 6098ea1b3e6SJed Brown 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat*); 611d4fbbf0eSBarry Smith 61228ce4d24SBarry Smith /*E 61328ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 61428ce4d24SBarry Smith 61528ce4d24SBarry Smith Level: beginner 61628ce4d24SBarry Smith 61728ce4d24SBarry Smith .seealso: KSPCGSetType() 61828ce4d24SBarry Smith E*/ 6199dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 6206a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 62128ce4d24SBarry Smith 622014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 6248031f4adStmunson 625dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal); 626dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*); 627dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*); 628fcae7a14Stmunson 629d2e0dd5eSTodd Munson PETSC_EXTERN PetscErrorCode KSPCGGLTRGetMinEig(KSP,PetscReal*); 630d2e0dd5eSTodd Munson PETSC_EXTERN PetscErrorCode KSPCGGLTRGetLambda(KSP,PetscReal*); 6318031f4adStmunson 632014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 6331d6018f0SLisandro Dalcin 634014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 635014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 6363369ce9aSBarry Smith 6379804daf3SBarry Smith #include <petscdrawtypes.h> 638d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 639d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 640d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 641d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6432f2e5d10SKris Buschelman 644014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 645014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 64603919abeSBarry Smith 647ba36735cSStefano Zampini /*S 648ba36735cSStefano Zampini KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 649f8a50e2bSBarry Smith 650ba36735cSStefano Zampini Level: beginner 651f8a50e2bSBarry Smith 652ba36735cSStefano Zampini Concepts: Krylov methods 653ba36735cSStefano Zampini 654ba36735cSStefano Zampini .seealso: KSPCreate(), KSPSetGuessType(), KSPGuessType 655ba36735cSStefano Zampini S*/ 656ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess; 657ba36735cSStefano Zampini /*J 658ba36735cSStefano Zampini KSPGuessType - String with the name of a PETSc initial guess for Krylov methods. 659ba36735cSStefano Zampini 660ba36735cSStefano Zampini Level: beginner 661ba36735cSStefano Zampini 662ba36735cSStefano Zampini .seealso: KSPGuess 663ba36735cSStefano Zampini J*/ 664ba36735cSStefano Zampini typedef const char* KSPGuessType; 665ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 666ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 667ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess); 668ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*); 669ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer); 670ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*); 671ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*); 672ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType); 673ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*); 674ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 675ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec); 676ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec); 677ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 678ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt); 679014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 680ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool); 681ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*); 682f8a50e2bSBarry Smith 683470b340bSDmitry Karpeev /*E 684470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 685470b340bSDmitry Karpeev 686470b340bSDmitry Karpeev Level: intermediate 687470b340bSDmitry Karpeev 688470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 689470b340bSDmitry Karpeev E*/ 6900c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType; 691470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 692470b340bSDmitry Karpeev 693014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 694014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 695d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 696bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 697aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 698bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 699470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 700470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 7015bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 7025a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 703470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*); 704470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 7053f22127dSBarry Smith 706*78e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*); 707*78e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*); 708*78e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*); 709*78e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 710*78e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 711*78e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 712*78e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 713cd929ea3SAlp Dener 714cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 715b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*); 716cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 717cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 718cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetJ0(Mat); 719cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 720cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 721cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 722cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 723cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 724cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 725cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*); 726cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*); 727cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*); 728cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*); 729cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*); 730cd929ea3SAlp Dener 731014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 732014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 733014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 734014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 735014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 736fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*); 73723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 738fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 73923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 74023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 741014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 742014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 743fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 744fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 7456c699258SBarry Smith 74602b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 7478c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, 748e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 749e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 750e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 751191494d9SMatthew G. Knepley PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec); 7522eac72dbSBarry Smith #endif 753