1 /* $Id: ksp.h,v 1.63 1998/03/06 00:21:12 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 #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 KSPTRLS "trls" 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 KSPDestroy(KSP); 36 37 extern DLList KSPList; 38 extern int KSPRegisterAll(char *); 39 extern int KSPRegisterDestroy(void); 40 41 42 extern int KSPGetType(KSP, KSPType *); 43 extern int KSPSetPreconditionerSide(KSP,PCSide); 44 extern int KSPGetPreconditionerSide(KSP,PCSide*); 45 extern int KSPGetTolerances(KSP,double*,double*,double*,int*); 46 extern int KSPSetTolerances(KSP,double,double,double,int); 47 extern int KSPSetComputeResidual(KSP,PetscTruth); 48 extern int KSPSetUsePreconditionedResidual(KSP); 49 extern int KSPSetInitialGuessNonzero(KSP); 50 extern int KSPSetComputeEigenvalues(KSP); 51 extern int KSPSetComputeSingularValues(KSP); 52 extern int KSPSetRhs(KSP,Vec); 53 extern int KSPGetRhs(KSP,Vec *); 54 extern int KSPSetSolution(KSP,Vec); 55 extern int KSPGetSolution(KSP,Vec *); 56 extern int KSPGetResidualNorm(KSP,double*); 57 58 extern int KSPSetPC(KSP,PC); 59 extern int KSPGetPC(KSP,PC*); 60 61 extern int KSPSetMonitor(KSP,int (*)(KSP,int,double, void*), void *); 62 extern int KSPGetMonitorContext(KSP,void **); 63 extern int KSPSetResidualHistory(KSP, double *,int); 64 65 extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double, void*), void *); 66 extern int KSPGetConvergenceContext(KSP,void **); 67 68 extern int KSPBuildSolution(KSP, Vec,Vec *); 69 extern int KSPBuildResidual(KSP, Vec, Vec,Vec *); 70 71 extern int KSPRichardsonSetScale(KSP , double); 72 extern int KSPChebychevSetEigenvalues(KSP , double, double); 73 extern int KSPComputeExtremeSingularValues(KSP, double*,double*); 74 extern int KSPComputeEigenvalues(KSP,int,double*,double*); 75 extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*); 76 77 extern int KSPGMRESSetRestart(KSP, int); 78 extern int KSPGMRESSetPreAllocateVectors(KSP); 79 extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int)); 80 extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int); 81 extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 82 extern int KSPGMRESIROrthogonalization(KSP,int); 83 extern int KSPGMRESDGKSOrthogonalization(KSP,int); 84 85 extern int KSPSetFromOptions(KSP); 86 extern int KSPAddOptionsChecker(int (*)(KSP)); 87 88 extern int KSPSingularValueMonitor(KSP,int,double, void * ); 89 extern int KSPDefaultMonitor(KSP,int,double, void *); 90 extern int KSPTrueMonitor(KSP,int,double, void *); 91 extern int KSPDefaultSMonitor(KSP,int,double, void *); 92 93 extern int KSPDefaultConverged(KSP,int,double, void *); 94 95 extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec); 96 extern int KSPUnwindPreconditioner(KSP,Vec,Vec); 97 extern int KSPDefaultBuildSolution(KSP,Vec,Vec*); 98 extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 99 100 extern int KSPPrintHelp(KSP); 101 102 extern int KSPSetOptionsPrefix(KSP,char*); 103 extern int KSPAppendOptionsPrefix(KSP,char*); 104 extern int KSPGetOptionsPrefix(KSP,char**); 105 106 extern int KSPView(KSP,Viewer); 107 108 extern int KSPComputeExplicitOperator(KSP,Mat *); 109 110 typedef enum {KSP_CG_SYMMETRIC=1, KSP_CG_HERMITIAN=2} KSPCGType; 111 extern int KSPCGSetType(KSP,KSPCGType); 112 113 extern int PCPreSolve(PC,KSP); 114 extern int PCPostSolve(PC,KSP); 115 116 extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*); 117 extern int KSPLGMonitor(KSP,int,double,void*); 118 extern int KSPLGMonitorDestroy(DrawLG); 119 extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*); 120 extern int KSPLGTrueMonitor(KSP,int,double,void*); 121 extern int KSPLGTrueMonitorDestroy(DrawLG); 122 123 #endif 124 125 126