1 /* $Id: ksp.h,v 1.88 2000/01/11 21:04:04 bsmith Exp bsmith $ */ 2 /* 3 Defines the interface functions for the Krylov subspace accelerators. 4 */ 5 #ifndef __KSP_H 6 #define __KSP_H 7 #include "pc.h" 8 9 #define KSP_COOKIE PETSC_COOKIE+8 10 11 typedef struct _p_KSP* KSP; 12 13 #define KSPRICHARDSON "richardson" 14 #define KSPCHEBYCHEV "chebychev" 15 #define KSPCG "cg" 16 #define KSPGMRES "gmres" 17 #define KSPTCQMR "tcqmr" 18 #define KSPBCGS "bcgs" 19 #define KSPCGS "cgs" 20 #define KSPTFQMR "tfqmr" 21 #define KSPCR "cr" 22 #define KSPLSQR "lsqr" 23 #define KSPPREONLY "preonly" 24 #define KSPQCG "qcg" 25 #define KSPBICG "bicg" 26 #define KSPFGMRES "fgmres" 27 #define KSPMINRES "minres" 28 typedef char * KSPType; 29 30 extern int KSPCreate(MPI_Comm,KSP *); 31 extern int KSPSetType(KSP,KSPType); 32 extern int KSPSetUp(KSP); 33 extern int KSPSolve(KSP,int *); 34 extern int KSPSolveTranspose(KSP,int *); 35 extern int KSPDestroy(KSP); 36 37 extern FList KSPList; 38 extern int KSPRegisterAll(char *); 39 extern int KSPRegisterDestroy(void); 40 41 extern int KSPRegister(char*,char*,char*,int(*)(KSP)); 42 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 43 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 44 #else 45 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 46 #endif 47 48 extern int KSPGetType(KSP,KSPType *); 49 extern int KSPSetPreconditionerSide(KSP,PCSide); 50 extern int KSPGetPreconditionerSide(KSP,PCSide*); 51 extern int KSPGetTolerances(KSP,double*,double*,double*,int*); 52 extern int KSPSetTolerances(KSP,double,double,double,int); 53 extern int KSPSetComputeResidual(KSP,PetscTruth); 54 extern int KSPSetUsePreconditionedResidual(KSP); 55 extern int KSPSetInitialGuessNonzero(KSP); 56 extern int KSPGetInitialGuessNonzero(KSP,PetscTruth *); 57 extern int KSPSetComputeEigenvalues(KSP); 58 extern int KSPSetComputeSingularValues(KSP); 59 extern int KSPSetRhs(KSP,Vec); 60 extern int KSPGetRhs(KSP,Vec *); 61 extern int KSPSetSolution(KSP,Vec); 62 extern int KSPGetSolution(KSP,Vec *); 63 extern int KSPGetResidualNorm(KSP,double*); 64 extern int KSPGetIterationNumber(KSP,int*); 65 66 extern int KSPSetPC(KSP,PC); 67 extern int KSPGetPC(KSP,PC*); 68 69 extern int KSPSetAvoidNorms(KSP); 70 71 extern int KSPSetMonitor(KSP,int (*)(KSP,int,double,void*),void *,int (*)(void*)); 72 extern int KSPClearMonitor(KSP); 73 extern int KSPGetMonitorContext(KSP,void **); 74 extern int KSPGetResidualHistory(KSP,double **,int *); 75 extern int KSPSetResidualHistory(KSP,double *,int,PetscTruth); 76 77 78 extern int KSPBuildSolution(KSP,Vec,Vec *); 79 extern int KSPBuildResidual(KSP,Vec,Vec,Vec *); 80 81 extern int KSPRichardsonSetScale(KSP,double); 82 extern int KSPChebychevSetEigenvalues(KSP,double,double); 83 extern int KSPComputeExtremeSingularValues(KSP,double*,double*); 84 extern int KSPComputeEigenvalues(KSP,int,double*,double*,int *); 85 extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*); 86 87 extern int KSPGMRESSetRestart(KSP,int); 88 extern int KSPGMRESSetPreAllocateVectors(KSP); 89 extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int)); 90 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int); 91 extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 92 extern int KSPGMRESIROrthogonalization(KSP,int); 93 94 extern int KSPFGMRESModifyPCNoChange(KSP,int,int,double,void*); 95 extern int KSPFGMRESModifyPCSLES(KSP,int,int,double,void*); 96 extern int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,double,void*),void*,int(*)(void*)); 97 98 extern int KSPSetFromOptions(KSP); 99 extern int KSPSetTypeFromOptions(KSP); 100 extern int KSPAddOptionsChecker(int (*)(KSP)); 101 102 extern int KSPSingularValueMonitor(KSP,int,double,void *); 103 extern int KSPDefaultMonitor(KSP,int,double,void *); 104 extern int KSPTrueMonitor(KSP,int,double,void *); 105 extern int KSPDefaultSMonitor(KSP,int,double,void *); 106 extern int KSPVecViewMonitor(KSP,int,double,void *); 107 extern int KSPGMRESKrylovMonitor(KSP,int,double,void *); 108 109 110 extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec); 111 extern int KSPUnwindPreconditioner(KSP,Vec,Vec); 112 extern int KSPDefaultBuildSolution(KSP,Vec,Vec*); 113 extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 114 115 extern int KSPPrintHelp(KSP); 116 117 extern int KSPSetOptionsPrefix(KSP,char*); 118 extern int KSPAppendOptionsPrefix(KSP,char*); 119 extern int KSPGetOptionsPrefix(KSP,char**); 120 121 extern int KSPView(KSP,Viewer); 122 123 typedef enum {/* converged */ 124 KSP_CONVERGED_RTOL = 2, 125 KSP_CONVERGED_ATOL = 3, 126 KSP_CONVERGED_ITS = 4, 127 KSP_CONVERGED_QCG_NEGATIVE_CURVE = 5, 128 KSP_CONVERGED_QCG_CONSTRAINED = 6, 129 KSP_CONVERGED_STEP_LENGTH = 7, 130 /* diverged */ 131 KSP_DIVERGED_ITS = -3, 132 KSP_DIVERGED_DTOL = -4, 133 KSP_DIVERGED_BREAKDOWN = -5, 134 KSP_DIVERGED_BREAKDOWN_BICG = -6, 135 KSP_DIVERGED_NONSYMMETRIC = -7, 136 KSP_DIVERGED_INDEFINITE_PC = -8, 137 138 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 139 140 extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double,KSPConvergedReason*,void*),void *); 141 extern int KSPGetConvergenceContext(KSP,void **); 142 extern int KSPDefaultConverged(KSP,int,double,KSPConvergedReason*,void *); 143 extern int KSPSkipConverged(KSP,int,double,KSPConvergedReason*,void *); 144 extern int KSPGetConvergedReason(KSP,KSPConvergedReason *); 145 146 extern int KSPComputeExplicitOperator(KSP,Mat *); 147 148 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType; 149 extern int KSPCGSetType(KSP,KSPCGType); 150 151 extern int PCPreSolve(PC,KSP); 152 extern int PCPostSolve(PC,KSP); 153 154 extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*); 155 extern int KSPLGMonitor(KSP,int,double,void*); 156 extern int KSPLGMonitorDestroy(DrawLG); 157 extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*); 158 extern int KSPLGTrueMonitor(KSP,int,double,void*); 159 extern int KSPLGTrueMonitorDestroy(DrawLG); 160 161 #endif 162 163 164