1 /* $Id: ksp.h,v 1.80 1999/04/21 18:19:32 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 "petsc.h" 8 #include "vec.h" 9 #include "mat.h" 10 #include "pc.h" 11 12 #define KSP_COOKIE PETSC_COOKIE+8 13 14 typedef struct _p_KSP* KSP; 15 16 #define KSPRICHARDSON "richardson" 17 #define KSPCHEBYCHEV "chebychev" 18 #define KSPCG "cg" 19 #define KSPGMRES "gmres" 20 #define KSPTCQMR "tcqmr" 21 #define KSPBCGS "bcgs" 22 #define KSPCGS "cgs" 23 #define KSPTFQMR "tfqmr" 24 #define KSPCR "cr" 25 #define KSPLSQR "lsqr" 26 #define KSPPREONLY "preonly" 27 #define KSPQCG "qcg" 28 #define KSPBICG "bicg" 29 typedef char * KSPType; 30 31 extern int KSPCreate(MPI_Comm,KSP *); 32 extern int KSPSetType(KSP,KSPType); 33 extern int KSPSetUp(KSP); 34 extern int KSPSolve(KSP,int *); 35 extern int KSPSolveTrans(KSP,int *); 36 extern int KSPDestroy(KSP); 37 38 extern FList KSPList; 39 extern int KSPRegisterAll(char *); 40 extern int KSPRegisterDestroy(void); 41 42 extern int KSPRegister_Private(char*,char*,char*,int(*)(KSP)); 43 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 44 #define KSPRegister(a,b,c,d) KSPRegister_Private(a,b,c,0) 45 #else 46 #define KSPRegister(a,b,c,d) KSPRegister_Private(a,b,c,d) 47 #endif 48 49 extern int KSPGetType(KSP, KSPType *); 50 extern int KSPSetPreconditionerSide(KSP,PCSide); 51 extern int KSPGetPreconditionerSide(KSP,PCSide*); 52 extern int KSPGetTolerances(KSP,double*,double*,double*,int*); 53 extern int KSPSetTolerances(KSP,double,double,double,int); 54 extern int KSPSetComputeResidual(KSP,PetscTruth); 55 extern int KSPSetUsePreconditionedResidual(KSP); 56 extern int KSPSetInitialGuessNonzero(KSP); 57 extern int KSPGetInitialGuessNonzero(KSP,PetscTruth *); 58 extern int KSPSetComputeEigenvalues(KSP); 59 extern int KSPSetComputeSingularValues(KSP); 60 extern int KSPSetRhs(KSP,Vec); 61 extern int KSPGetRhs(KSP,Vec *); 62 extern int KSPSetSolution(KSP,Vec); 63 extern int KSPGetSolution(KSP,Vec *); 64 extern int KSPGetResidualNorm(KSP,double*); 65 extern int KSPGetIterationNumber(KSP,int*); 66 67 extern int KSPSetPC(KSP,PC); 68 extern int KSPGetPC(KSP,PC*); 69 70 extern int KSPSetAvoidNorms(KSP); 71 72 extern int KSPSetMonitor(KSP,int (*)(KSP,int,double, void*), void *,int (*)(void*)); 73 extern int KSPClearMonitor(KSP); 74 extern int KSPGetMonitorContext(KSP,void **); 75 extern int KSPGetResidualHistory(KSP, double **, int *); 76 extern int KSPSetResidualHistory(KSP, double *,int,PetscTruth); 77 78 extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double, void*), void *); 79 extern int KSPGetConvergenceContext(KSP,void **); 80 81 extern int KSPBuildSolution(KSP, Vec,Vec *); 82 extern int KSPBuildResidual(KSP, Vec, Vec,Vec *); 83 84 extern int KSPRichardsonSetScale(KSP , double); 85 extern int KSPChebychevSetEigenvalues(KSP , double, double); 86 extern int KSPComputeExtremeSingularValues(KSP, double*,double*); 87 extern int KSPComputeEigenvalues(KSP,int,double*,double*,int *); 88 extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*); 89 90 extern int KSPGMRESSetRestart(KSP, int); 91 extern int KSPGMRESSetPreAllocateVectors(KSP); 92 extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int)); 93 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int); 94 extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 95 extern int KSPGMRESIROrthogonalization(KSP,int); 96 extern int KSPGMRESPrestartSet(KSP,int); 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 extern int KSPDefaultConverged(KSP,int,double, void *); 110 extern int KSPSkipConverged(KSP,int,double, void *); 111 112 extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec); 113 extern int KSPUnwindPreconditioner(KSP,Vec,Vec); 114 extern int KSPDefaultBuildSolution(KSP,Vec,Vec*); 115 extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 116 117 extern int KSPPrintHelp(KSP); 118 119 extern int KSPSetOptionsPrefix(KSP,char*); 120 extern int KSPAppendOptionsPrefix(KSP,char*); 121 extern int KSPGetOptionsPrefix(KSP,char**); 122 123 extern int KSPView(KSP,Viewer); 124 125 extern int KSPComputeExplicitOperator(KSP,Mat *); 126 127 typedef enum {KSP_CG_SYMMETRIC=1, KSP_CG_HERMITIAN=2} KSPCGType; 128 extern int KSPCGSetType(KSP,KSPCGType); 129 130 extern int PCPreSolve(PC,KSP); 131 extern int PCPostSolve(PC,KSP); 132 133 extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*); 134 extern int KSPLGMonitor(KSP,int,double,void*); 135 extern int KSPLGMonitorDestroy(DrawLG); 136 extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*); 137 extern int KSPLGTrueMonitor(KSP,int,double,void*); 138 extern int KSPLGTrueMonitorDestroy(DrawLG); 139 140 #endif 141 142 143