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 PETSC_EXTERN_CXX_BEGIN 9 10 /*S 11 KSP - Abstract PETSc object that manages all Krylov methods 12 13 Level: beginner 14 15 Concepts: Krylov methods 16 17 .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, SLES 18 S*/ 19 typedef struct _p_KSP* KSP; 20 21 /*E 22 KSPType - String with the name of a PETSc Krylov method or the creation function 23 with an optional dynamic library name, for example 24 http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 25 26 Level: beginner 27 28 .seealso: KSPSetType(), KSP 29 E*/ 30 #define KSPRICHARDSON "richardson" 31 #define KSPCHEBYCHEV "chebychev" 32 #define KSPCG "cg" 33 #define KSPCGNE "cgne" 34 #define KSPGMRES "gmres" 35 #define KSPTCQMR "tcqmr" 36 #define KSPBCGS "bcgs" 37 #define KSPCGS "cgs" 38 #define KSPTFQMR "tfqmr" 39 #define KSPCR "cr" 40 #define KSPLSQR "lsqr" 41 #define KSPPREONLY "preonly" 42 #define KSPQCG "qcg" 43 #define KSPBICG "bicg" 44 #define KSPFGMRES "fgmres" 45 #define KSPMINRES "minres" 46 #define KSPSYMMLQ "symmlq" 47 #define KSPLGMRES "lgmres" 48 typedef char * KSPType; 49 50 /* Logging support */ 51 extern int KSP_COOKIE; 52 extern int KSP_GMRESOrthogonalization; 53 54 EXTERN int KSPCreate(MPI_Comm,KSP *); 55 EXTERN int KSPSetType(KSP,KSPType); 56 EXTERN int KSPSetUp(KSP); 57 EXTERN int KSPSolve(KSP,int *); 58 EXTERN int KSPSolveTranspose(KSP,int *); 59 EXTERN int KSPDestroy(KSP); 60 61 extern PetscFList KSPList; 62 EXTERN int KSPRegisterAll(const char[]); 63 EXTERN int KSPRegisterDestroy(void); 64 65 EXTERN int KSPRegister(const char[],const char[],const char[],int(*)(KSP)); 66 67 /*MC 68 KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 69 70 Synopsis: 71 int KSPRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(KSP)) 72 73 Not Collective 74 75 Input Parameters: 76 + name_solver - name of a new user-defined solver 77 . path - path (either absolute or relative) the library containing this solver 78 . name_create - name of routine to create method context 79 - routine_create - routine to create method context 80 81 Notes: 82 KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 83 84 If dynamic libraries are used, then the fourth input argument (routine_create) 85 is ignored. 86 87 Sample usage: 88 .vb 89 KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 90 "MySolverCreate",MySolverCreate); 91 .ve 92 93 Then, your solver can be chosen with the procedural interface via 94 $ KSPSetType(ksp,"my_solver") 95 or at runtime via the option 96 $ -ksp_type my_solver 97 98 Level: advanced 99 100 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, 101 and others of the form ${any_environmental_variable} occuring in pathname will be 102 replaced with appropriate values. 103 If your function is not being put into a shared library then use KSPRegister() instead 104 105 .keywords: KSP, register 106 107 .seealso: KSPRegisterAll(), KSPRegisterDestroy() 108 109 M*/ 110 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 111 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 112 #else 113 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 114 #endif 115 116 EXTERN int KSPGetType(KSP,KSPType *); 117 EXTERN int KSPSetPreconditionerSide(KSP,PCSide); 118 EXTERN int KSPGetPreconditionerSide(KSP,PCSide*); 119 EXTERN int KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*); 120 EXTERN int KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int); 121 EXTERN int KSPSetInitialGuessNonzero(KSP,PetscTruth); 122 EXTERN int KSPGetInitialGuessNonzero(KSP,PetscTruth *); 123 EXTERN int KSPSetInitialGuessKnoll(KSP,PetscTruth); 124 EXTERN int KSPGetInitialGuessKnoll(KSP,PetscTruth*); 125 EXTERN int KSPSetComputeEigenvalues(KSP,PetscTruth); 126 EXTERN int KSPSetComputeSingularValues(KSP,PetscTruth); 127 EXTERN int KSPSetRhs(KSP,Vec); 128 EXTERN int KSPGetRhs(KSP,Vec *); 129 EXTERN int KSPSetSolution(KSP,Vec); 130 EXTERN int KSPGetSolution(KSP,Vec *); 131 EXTERN int KSPGetResidualNorm(KSP,PetscReal*); 132 EXTERN int KSPGetIterationNumber(KSP,int*); 133 134 EXTERN int KSPSetPC(KSP,PC); 135 EXTERN int KSPGetPC(KSP,PC*); 136 137 EXTERN int KSPSetMonitor(KSP,int (*)(KSP,int,PetscReal,void*),void *,int (*)(void*)); 138 EXTERN int KSPClearMonitor(KSP); 139 EXTERN int KSPGetMonitorContext(KSP,void **); 140 EXTERN int KSPGetResidualHistory(KSP,PetscReal*[],int *); 141 EXTERN int KSPSetResidualHistory(KSP,PetscReal[],int,PetscTruth); 142 143 144 EXTERN int KSPBuildSolution(KSP,Vec,Vec *); 145 EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *); 146 147 EXTERN int KSPRichardsonSetScale(KSP,PetscReal); 148 EXTERN int KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 149 EXTERN int KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 150 EXTERN int KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *); 151 EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*); 152 153 #define KSPGMRESSetRestart(ksp,r) PetscTryMethod((ksp),KSPGMRESSetRestart_C,(KSP,int),((ksp),(r))) 154 #define KSPGMRESSetHapTol(ksp,tol) PetscTryMethod((ksp),KSPGMRESSetHapTol_C,(KSP,PetscReal),((ksp),(tol))) 155 156 EXTERN int KSPGMRESSetPreAllocateVectors(KSP); 157 EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int)); 158 EXTERN int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int); 159 EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 160 EXTERN int KSPGMRESIROrthogonalization(KSP,int); 161 162 EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*); 163 EXTERN int KSPFGMRESModifyPCSLES(KSP,int,int,PetscReal,void*); 164 EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,PetscReal,void*),void*,int(*)(void*)); 165 166 EXTERN int KSPQCGSetTrustRegionRadius(KSP,PetscReal); 167 EXTERN int KSPQCGGetQuadratic(KSP,PetscReal*); 168 EXTERN int KSPQCGGetTrialStepNorm(KSP,PetscReal*); 169 170 EXTERN int KSPSetFromOptions(KSP); 171 EXTERN int KSPAddOptionsChecker(int (*)(KSP)); 172 173 EXTERN int KSPSingularValueMonitor(KSP,int,PetscReal,void *); 174 EXTERN int KSPDefaultMonitor(KSP,int,PetscReal,void *); 175 EXTERN int KSPTrueMonitor(KSP,int,PetscReal,void *); 176 EXTERN int KSPDefaultSMonitor(KSP,int,PetscReal,void *); 177 EXTERN int KSPVecViewMonitor(KSP,int,PetscReal,void *); 178 EXTERN int KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *); 179 180 EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec); 181 EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*); 182 EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 183 184 EXTERN int KSPSetOptionsPrefix(KSP,const char[]); 185 EXTERN int KSPAppendOptionsPrefix(KSP,const char[]); 186 EXTERN int KSPGetOptionsPrefix(KSP,char*[]); 187 188 EXTERN int KSPView(KSP,PetscViewer); 189 190 /*E 191 KSPNormType - Norm that is passed in the Krylov convergence 192 test routines. 193 194 Level: advanced 195 196 Notes: this must match finclude/petscksp.h 197 198 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 199 KSPSetConvergenceTest() 200 E*/ 201 typedef enum {KSP_NO_NORM = 0, 202 KSP_PRECONDITIONED_NORM = 1, 203 KSP_UNPRECONDITIONED_NORM = 2, 204 KSP_NATURAL_NORM = 3} KSPNormType; 205 EXTERN int KSPSetNormType(KSP,KSPNormType); 206 207 /*E 208 KSPConvergedReason - reason a Krylov method was said to 209 have converged or diverged 210 211 Level: beginner 212 213 Notes: this must match finclude/petscksp.h 214 215 Developer note: The string versions of these are in 216 src/sles/ksp/interface/itfunc.c called convergedreasons. 217 If these enums are changed you much change those. 218 219 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason() 220 E*/ 221 typedef enum {/* converged */ 222 KSP_CONVERGED_RTOL = 2, 223 KSP_CONVERGED_ATOL = 3, 224 KSP_CONVERGED_ITS = 4, 225 KSP_CONVERGED_QCG_NEG_CURVE = 5, 226 KSP_CONVERGED_QCG_CONSTRAINED = 6, 227 KSP_CONVERGED_STEP_LENGTH = 7, 228 /* diverged */ 229 KSP_DIVERGED_ITS = -3, 230 KSP_DIVERGED_DTOL = -4, 231 KSP_DIVERGED_BREAKDOWN = -5, 232 KSP_DIVERGED_BREAKDOWN_BICG = -6, 233 KSP_DIVERGED_NONSYMMETRIC = -7, 234 KSP_DIVERGED_INDEFINITE_PC = -8, 235 236 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 237 238 EXTERN int KSPSetConvergenceTest(KSP,int (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *); 239 EXTERN int KSPGetConvergenceContext(KSP,void **); 240 EXTERN int KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *); 241 EXTERN int KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *); 242 EXTERN int KSPGetConvergedReason(KSP,KSPConvergedReason *); 243 244 EXTERN int KSPComputeExplicitOperator(KSP,Mat *); 245 246 /*E 247 KSPCGType - Determines what type of CG to use 248 249 Level: beginner 250 251 .seealso: KSPCGSetType() 252 E*/ 253 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType; 254 255 EXTERN int KSPCGSetType(KSP,KSPCGType); 256 257 EXTERN int PCPreSolve(PC,KSP); 258 EXTERN int PCPostSolve(PC,KSP); 259 260 EXTERN int KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 261 EXTERN int KSPLGMonitor(KSP,int,PetscReal,void*); 262 EXTERN int KSPLGMonitorDestroy(PetscDrawLG); 263 EXTERN int KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 264 EXTERN int KSPLGTrueMonitor(KSP,int,PetscReal,void*); 265 EXTERN int KSPLGTrueMonitorDestroy(PetscDrawLG); 266 267 PETSC_EXTERN_CXX_END 268 #endif 269