1 /* $Id: petscksp.h,v 1.104 2001/07/09 22:52:08 buschelm Exp buschelm $ */ 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 #define KSP_COOKIE PETSC_COOKIE+8 10 11 /*S 12 KSP - Abstract PETSc object that manages all Krylov methods 13 14 Level: beginner 15 16 Concepts: Krylov methods 17 18 .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, SLES 19 S*/ 20 typedef struct _p_KSP* KSP; 21 22 /*E 23 KSPType - String with the name of a PETSc Krylov method or the creation function 24 with an optional dynamic library name, for example 25 http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 26 27 Level: beginner 28 29 .seealso: KSPSetType(), KSP 30 E*/ 31 #define KSPRICHARDSON "richardson" 32 #define KSPCHEBYCHEV "chebychev" 33 #define KSPCG "cg" 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 typedef char * KSPType; 48 49 EXTERN int KSPCreate(MPI_Comm,KSP *); 50 EXTERN int KSPSetType(KSP,KSPType); 51 EXTERN int KSPSetUp(KSP); 52 EXTERN int KSPSolve(KSP,int *); 53 EXTERN int KSPSolveTranspose(KSP,int *); 54 EXTERN int KSPDestroy(KSP); 55 56 extern PetscFList KSPList; 57 EXTERN int KSPRegisterAll(char *); 58 EXTERN int KSPRegisterDestroy(void); 59 60 EXTERN int KSPRegister(char*,char*,char*,int(*)(KSP)); 61 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 62 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 63 #else 64 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 65 #endif 66 67 EXTERN int KSPGetType(KSP,KSPType *); 68 EXTERN int KSPSetPreconditionerSide(KSP,PCSide); 69 EXTERN int KSPGetPreconditionerSide(KSP,PCSide*); 70 EXTERN int KSPGetTolerances(KSP,double*,double*,double*,int*); 71 EXTERN int KSPSetTolerances(KSP,double,double,double,int); 72 EXTERN int KSPSetComputeResidual(KSP,PetscTruth); 73 EXTERN int KSPSetUsePreconditionedResidual(KSP,PetscTruth); 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,double*); 83 EXTERN int KSPGetIterationNumber(KSP,int*); 84 85 EXTERN int KSPSetPC(KSP,PC); 86 EXTERN int KSPGetPC(KSP,PC*); 87 88 EXTERN int KSPSetAvoidNorms(KSP,PetscTruth); 89 90 EXTERN int KSPSetMonitor(KSP,int (*)(KSP,int,double,void*),void *,int (*)(void*)); 91 EXTERN int KSPClearMonitor(KSP); 92 EXTERN int KSPGetMonitorContext(KSP,void **); 93 EXTERN int KSPGetResidualHistory(KSP,double **,int *); 94 EXTERN int KSPSetResidualHistory(KSP,double *,int,PetscTruth); 95 96 97 EXTERN int KSPBuildSolution(KSP,Vec,Vec *); 98 EXTERN int KSPBuildResidual(KSP,Vec,Vec,Vec *); 99 100 EXTERN int KSPRichardsonSetScale(KSP,double); 101 EXTERN int KSPChebychevSetEigenvalues(KSP,double,double); 102 EXTERN int KSPComputeExtremeSingularValues(KSP,double*,double*); 103 EXTERN int KSPComputeEigenvalues(KSP,int,double*,double*,int *); 104 EXTERN int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*); 105 106 EXTERN int KSPGMRESSetRestart(KSP,int); 107 EXTERN int KSPGMRESSetHapTol(KSP,double); 108 EXTERN int KSPGMRESSetPreAllocateVectors(KSP); 109 EXTERN int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int)); 110 EXTERN int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int); 111 EXTERN int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 112 EXTERN int KSPGMRESIROrthogonalization(KSP,int); 113 114 EXTERN int KSPFGMRESModifyPCNoChange(KSP,int,int,double,void*); 115 EXTERN int KSPFGMRESModifyPCSLES(KSP,int,int,double,void*); 116 EXTERN int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,double,void*),void*,int(*)(void*)); 117 118 EXTERN int KSPQCGSetTrustRegionRadius(KSP,double); 119 EXTERN int KSPQCGGetQuadratic(KSP,double*); 120 EXTERN int KSPQCGGetTrialStepNorm(KSP,double*); 121 122 EXTERN int KSPSetFromOptions(KSP); 123 EXTERN int KSPSetTypeFromOptions(KSP); 124 EXTERN int KSPAddOptionsChecker(int (*)(KSP)); 125 126 EXTERN int KSPSingularValueMonitor(KSP,int,double,void *); 127 EXTERN int KSPDefaultMonitor(KSP,int,double,void *); 128 EXTERN int KSPTrueMonitor(KSP,int,double,void *); 129 EXTERN int KSPDefaultSMonitor(KSP,int,double,void *); 130 EXTERN int KSPVecViewMonitor(KSP,int,double,void *); 131 EXTERN int KSPGMRESKrylovMonitor(KSP,int,double,void *); 132 133 EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec); 134 EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*); 135 EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 136 137 EXTERN int KSPSetOptionsPrefix(KSP,char*); 138 EXTERN int KSPAppendOptionsPrefix(KSP,char*); 139 EXTERN int KSPGetOptionsPrefix(KSP,char**); 140 141 EXTERN int KSPView(KSP,PetscViewer); 142 143 /*E 144 KSPConvergedReason - reason a Krylov method was said to 145 have converged or diverged 146 147 Level: beginner 148 149 Notes: this must match finclude/petscksp.h 150 151 .seealso: SLESSolve(), KSPSolve(), KSPGetConvergedReason() 152 E*/ 153 typedef enum {/* converged */ 154 KSP_CONVERGED_RTOL = 2, 155 KSP_CONVERGED_ATOL = 3, 156 KSP_CONVERGED_ITS = 4, 157 KSP_CONVERGED_QCG_NEG_CURVE = 5, 158 KSP_CONVERGED_QCG_CONSTRAINED = 6, 159 KSP_CONVERGED_STEP_LENGTH = 7, 160 /* diverged */ 161 KSP_DIVERGED_ITS = -3, 162 KSP_DIVERGED_DTOL = -4, 163 KSP_DIVERGED_BREAKDOWN = -5, 164 KSP_DIVERGED_BREAKDOWN_BICG = -6, 165 KSP_DIVERGED_NONSYMMETRIC = -7, 166 KSP_DIVERGED_INDEFINITE_PC = -8, 167 168 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 169 170 EXTERN int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double,KSPConvergedReason*,void*),void *); 171 EXTERN int KSPGetConvergenceContext(KSP,void **); 172 EXTERN int KSPDefaultConverged(KSP,int,double,KSPConvergedReason*,void *); 173 EXTERN int KSPSkipConverged(KSP,int,double,KSPConvergedReason*,void *); 174 EXTERN int KSPGetConvergedReason(KSP,KSPConvergedReason *); 175 176 EXTERN int KSPComputeExplicitOperator(KSP,Mat *); 177 178 /*E 179 KSPCGType - Determines what type of CG to use 180 181 Level: beginner 182 183 .seealso: KSPCGSetType() 184 E*/ 185 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType; 186 187 EXTERN int KSPCGSetType(KSP,KSPCGType); 188 189 EXTERN int PCPreSolve(PC,KSP); 190 EXTERN int PCPostSolve(PC,KSP); 191 192 EXTERN int KSPLGMonitorCreate(char*,char*,int,int,int,int,PetscDrawLG*); 193 EXTERN int KSPLGMonitor(KSP,int,double,void*); 194 EXTERN int KSPLGMonitorDestroy(PetscDrawLG); 195 EXTERN int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,PetscDrawLG*); 196 EXTERN int KSPLGTrueMonitor(KSP,int,double,void*); 197 EXTERN int KSPLGTrueMonitorDestroy(PetscDrawLG); 198 199 #endif 200 201 202