1 /* $Id: petscksp.h,v 1.107 2001/08/06 21:16:38 bsmith Exp $ */ 2 /* 3 Defines the interface functions for the Krylov subspace accelerators. 4 */ 5 #ifndef __PETSCKSP_H 6 #define __PETSCKSP_H 7 #include "petscpc.h" 8 9 /*S 10 KSP - Abstract PETSc object that manages all Krylov methods 11 12 Level: beginner 13 14 Concepts: Krylov methods 15 16 .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, SLES 17 S*/ 18 typedef struct _p_KSP* KSP; 19 20 /*E 21 KSPType - String with the name of a PETSc Krylov method or the creation function 22 with an optional dynamic library name, for example 23 http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 24 25 Level: beginner 26 27 .seealso: KSPSetType(), KSP 28 E*/ 29 #define KSPRICHARDSON "richardson" 30 #define KSPCHEBYCHEV "chebychev" 31 #define KSPCG "cg" 32 #define KSPGMRES "gmres" 33 #define KSPTCQMR "tcqmr" 34 #define KSPBCGS "bcgs" 35 #define KSPCGS "cgs" 36 #define KSPTFQMR "tfqmr" 37 #define KSPCR "cr" 38 #define KSPLSQR "lsqr" 39 #define KSPPREONLY "preonly" 40 #define KSPQCG "qcg" 41 #define KSPBICG "bicg" 42 #define KSPFGMRES "fgmres" 43 #define KSPMINRES "minres" 44 #define KSPSYMMLQ "symmlq" 45 typedef char * KSPType; 46 47 /* Logging support */ 48 extern int KSP_COOKIE; 49 extern int KSP_GMRESOrthogonalization; 50 51 EXTERN int KSPCreate(MPI_Comm,KSP *); 52 EXTERN int KSPSetType(KSP,KSPType); 53 EXTERN int KSPSetUp(KSP); 54 EXTERN int KSPSolve(KSP,int *); 55 EXTERN int KSPSolveTranspose(KSP,int *); 56 EXTERN int KSPDestroy(KSP); 57 58 extern PetscFList KSPList; 59 EXTERN int KSPRegisterAll(char *); 60 EXTERN int KSPRegisterDestroy(void); 61 62 EXTERN int KSPRegister(char*,char*,char*,int(*)(KSP)); 63 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 64 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 65 #else 66 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 67 #endif 68 69 EXTERN int KSPGetType(KSP,KSPType *); 70 EXTERN int KSPSetPreconditionerSide(KSP,PCSide); 71 EXTERN int KSPGetPreconditionerSide(KSP,PCSide*); 72 EXTERN int KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*); 73 EXTERN int KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int); 74 EXTERN int KSPSetInitialGuessNonzero(KSP,PetscTruth); 75 EXTERN int KSPGetInitialGuessNonzero(KSP,PetscTruth *); 76 EXTERN int KSPSetComputeEigenvalues(KSP,PetscTruth); 77 EXTERN int KSPSetComputeSingularValues(KSP,PetscTruth); 78 EXTERN int KSPSetRhs(KSP,Vec); 79 EXTERN int KSPGetRhs(KSP,Vec *); 80 EXTERN int KSPSetSolution(KSP,Vec); 81 EXTERN int KSPGetSolution(KSP,Vec *); 82 EXTERN int KSPGetResidualNorm(KSP,PetscReal*); 83 EXTERN int KSPGetIterationNumber(KSP,int*); 84 85 EXTERN int KSPSetPC(KSP,PC); 86 EXTERN int KSPGetPC(KSP,PC*); 87 88 EXTERN int KSPSetMonitor(KSP,int (*)(KSP,int,PetscReal,void*),void *,int (*)(void*)); 89 EXTERN int KSPClearMonitor(KSP); 90 EXTERN int KSPGetMonitorContext(KSP,void **); 91 EXTERN int KSPGetResidualHistory(KSP,PetscReal **,int *); 92 EXTERN int KSPSetResidualHistory(KSP,PetscReal *,int,PetscTruth); 93 94 95 EXTERN int KSPBuildSolution(KSP,Vec,Vec *); 96 EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *); 97 98 EXTERN int KSPRichardsonSetScale(KSP,PetscReal); 99 EXTERN int KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 100 EXTERN int KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 101 EXTERN int KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *); 102 EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*); 103 104 #define KSPGMRESSetRestart(ksp,r) PetscTryMethod((ksp),KSPGMRESSetRestart_C,(KSP,int),((ksp),(r))) 105 #define KSPGMRESSetHapTol(ksp,tol) PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol))) 106 107 EXTERN int KSPGMRESSetPreAllocateVectors(KSP); 108 EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int)); 109 EXTERN int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int); 110 EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 111 EXTERN int KSPGMRESIROrthogonalization(KSP,int); 112 113 EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*); 114 EXTERN int KSPFGMRESModifyPCSLES(KSP,int,int,PetscReal,void*); 115 EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int(*)(void*)); 116 117 EXTERN int KSPQCGSetTrustRegionRadius(KSP,PetscReal); 118 EXTERN int KSPQCGGetQuadratic(KSP,PetscReal*); 119 EXTERN int KSPQCGGetTrialStepNorm(KSP,PetscReal*); 120 121 EXTERN int KSPSetFromOptions(KSP); 122 EXTERN int KSPAddOptionsChecker(int (*)(KSP)); 123 124 EXTERN int KSPSingularValueMonitor(KSP,int,PetscReal,void *); 125 EXTERN int KSPDefaultMonitor(KSP,int,PetscReal,void *); 126 EXTERN int KSPTrueMonitor(KSP,int,PetscReal,void *); 127 EXTERN int KSPDefaultSMonitor(KSP,int,PetscReal,void *); 128 EXTERN int KSPVecViewMonitor(KSP,int,PetscReal,void *); 129 EXTERN int KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *); 130 131 EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec); 132 EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*); 133 EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 134 135 EXTERN int KSPSetOptionsPrefix(KSP,char*); 136 EXTERN int KSPAppendOptionsPrefix(KSP,char*); 137 EXTERN int KSPGetOptionsPrefix(KSP,char**); 138 139 EXTERN int KSPView(KSP,PetscViewer); 140 141 /*E 142 KSPNormType - Norm that is passed in the Krylov convergence 143 test routines. 144 145 Level: advanced 146 147 Notes: this must match finclude/petscksp.h 148 149 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 150 KSPSetConvergenceTest() 151 E*/ 152 typedef enum {KSP_NO_NORM = 0, 153 KSP_PRECONDITIONED_NORM = 1, 154 KSP_UNPRECONDITIONED_NORM = 2, 155 KSP_NATURAL_NORM = 3} KSPNormType; 156 EXTERN int KSPSetNormType(KSP,KSPNormType); 157 158 /*E 159 KSPConvergedReason - reason a Krylov method was said to 160 have converged or diverged 161 162 Level: beginner 163 164 Notes: this must match finclude/petscksp.h 165 166 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason() 167 E*/ 168 typedef enum {/* converged */ 169 KSP_CONVERGED_RTOL = 2, 170 KSP_CONVERGED_ATOL = 3, 171 KSP_CONVERGED_ITS = 4, 172 KSP_CONVERGED_QCG_NEG_CURVE = 5, 173 KSP_CONVERGED_QCG_CONSTRAINED = 6, 174 KSP_CONVERGED_STEP_LENGTH = 7, 175 /* diverged */ 176 KSP_DIVERGED_ITS = -3, 177 KSP_DIVERGED_DTOL = -4, 178 KSP_DIVERGED_BREAKDOWN = -5, 179 KSP_DIVERGED_BREAKDOWN_BICG = -6, 180 KSP_DIVERGED_NONSYMMETRIC = -7, 181 KSP_DIVERGED_INDEFINITE_PC = -8, 182 183 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 184 185 EXTERN int KSPSetConvergenceTest(KSP,int (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *); 186 EXTERN int KSPGetConvergenceContext(KSP,void **); 187 EXTERN int KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *); 188 EXTERN int KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *); 189 EXTERN int KSPGetConvergedReason(KSP,KSPConvergedReason *); 190 191 EXTERN int KSPComputeExplicitOperator(KSP,Mat *); 192 193 /*E 194 KSPCGType - Determines what type of CG to use 195 196 Level: beginner 197 198 .seealso: KSPCGSetType() 199 E*/ 200 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType; 201 202 EXTERN int KSPCGSetType(KSP,KSPCGType); 203 204 EXTERN int PCPreSolve(PC,KSP); 205 EXTERN int PCPostSolve(PC,KSP); 206 207 EXTERN int KSPLGMonitorCreate(char*,char*,int,int,int,int,PetscDrawLG*); 208 EXTERN int KSPLGMonitor(KSP,int,PetscReal,void*); 209 EXTERN int KSPLGMonitorDestroy(PetscDrawLG); 210 EXTERN int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,PetscDrawLG*); 211 EXTERN int KSPLGTrueMonitor(KSP,int,PetscReal,void*); 212 EXTERN int KSPLGTrueMonitorDestroy(PetscDrawLG); 213 214 #endif 215 216 217