1 /* $Id: ksp.h,v 1.55 1997/06/05 12:59:23 bsmith Exp bsmith $ */ 2 /* 3 Defines the interface functions for the Krylov subspace accelerators. 4 */ 5 #ifndef __KSP_PACKAGE 6 #define __KSP_PACKAGE 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 typedef enum { KSPRICHARDSON, KSPCHEBYCHEV, KSPCG, KSPGMRES, KSPTCQMR, KSPBCGS, 17 KSPCGS, KSPTFQMR, KSPCR, KSPLSQR, KSPPREONLY, KSPQCG, KSPNEW} KSPType; 18 19 extern int KSPCreate(MPI_Comm,KSP *); 20 extern int KSPSetType(KSP,KSPType); 21 extern int KSPSetUp(KSP); 22 extern int KSPSolve(KSP,int *); 23 extern int KSPDestroy(KSP); 24 25 extern int KSPRegisterAll(); 26 extern int KSPRegisterDestroy(); 27 extern int KSPRegister(KSPType,KSPType*,char *,int (*)(KSP)); 28 extern int KSPRegisterAllCalled; 29 30 extern int KSPGetType(KSP, KSPType *,char **); 31 extern int KSPSetPreconditionerSide(KSP,PCSide); 32 extern int KSPGetPreconditionerSide(KSP,PCSide*); 33 extern int KSPGetTolerances(KSP,double*,double*,double*,int*); 34 extern int KSPSetTolerances(KSP,double,double,double,int); 35 extern int KSPSetComputeResidual(KSP,PetscTruth); 36 extern int KSPSetUsePreconditionedResidual(KSP); 37 extern int KSPSetInitialGuessNonzero(KSP); 38 extern int KSPSetComputeEigenvalues(KSP); 39 extern int KSPSetComputeSingularValues(KSP); 40 extern int KSPSetRhs(KSP,Vec); 41 extern int KSPGetRhs(KSP,Vec *); 42 extern int KSPSetSolution(KSP,Vec); 43 extern int KSPGetSolution(KSP,Vec *); 44 45 extern int KSPSetPC(KSP,PC); 46 extern int KSPGetPC(KSP,PC*); 47 48 extern int KSPSetMonitor(KSP,int (*)(KSP,int,double, void*), void *); 49 extern int KSPGetMonitorContext(KSP,void **); 50 extern int KSPSetResidualHistory(KSP, double *,int); 51 52 extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double, void*), void *); 53 extern int KSPGetConvergenceContext(KSP,void **); 54 55 extern int KSPBuildSolution(KSP, Vec,Vec *); 56 extern int KSPBuildResidual(KSP, Vec, Vec,Vec *); 57 58 extern int KSPRichardsonSetScale(KSP , double); 59 extern int KSPChebychevSetEigenvalues(KSP , double, double); 60 extern int KSPComputeExtremeSingularValues(KSP, double*,double*); 61 extern int KSPComputeEigenvalues(KSP,int,double*,double*); 62 extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*); 63 64 extern int KSPGMRESSetRestart(KSP, int); 65 extern int KSPGMRESSetPreAllocateVectors(KSP); 66 extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int)); 67 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int); 68 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalizationLocal(KSP,int); 69 extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 70 extern int KSPGMRESIROrthogonalization(KSP,int); 71 extern int KSPGMRESDGKSOrthogonalization(KSP,int); 72 73 extern int KSPSetFromOptions(KSP); 74 extern int KSPAddOptionsChecker(int (*)(KSP)); 75 76 extern int KSPSingularValueMonitor(KSP,int,double, void * ); 77 extern int KSPDefaultMonitor(KSP,int,double, void *); 78 extern int KSPTrueMonitor(KSP,int,double, void *); 79 extern int KSPDefaultSMonitor(KSP,int,double, void *); 80 81 extern int KSPDefaultConverged(KSP,int,double, void *); 82 83 extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec); 84 extern int KSPUnwindPreconditioner(KSP,Vec,Vec); 85 extern int KSPDefaultBuildSolution(KSP,Vec,Vec*); 86 extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 87 88 extern int KSPPrintHelp(KSP); 89 90 extern int KSPSetOptionsPrefix(KSP,char*); 91 extern int KSPAppendOptionsPrefix(KSP,char*); 92 extern int KSPGetOptionsPrefix(KSP,char**); 93 94 extern int KSPView(KSP,Viewer); 95 96 extern int KSPComputeExplicitOperator(KSP,Mat *); 97 98 typedef enum {KSP_CG_SYMMETRIC=1, KSP_CG_HERMITIAN=2} KSPCGType; 99 extern int KSPCGSetType(KSP,KSPCGType); 100 101 extern int PCPreSolve(PC,KSP); 102 extern int PCPostSolve(PC,KSP); 103 104 extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*); 105 extern int KSPLGMonitor(KSP,int,double,void*); 106 extern int KSPLGMonitorDestroy(DrawLG); 107 extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*); 108 extern int KSPLGTrueMonitor(KSP,int,double,void*); 109 extern int KSPLGTrueMonitorDestroy(DrawLG); 110 111 #endif 112 113 114