xref: /petsc/include/petscksp.h (revision 329f5518e9d4bb7ce96c0c5576cc53785c973973)
1*329f5518SBarry Smith /* $Id: ksp.h,v 1.87 1999/12/08 22:17:32 balay Exp bsmith $ */
2f26ada1bSBarry Smith /*
3f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
4f26ada1bSBarry Smith */
588d459dfSBarry Smith #ifndef __KSP_H
688d459dfSBarry Smith #define __KSP_H
7aabeff55SBarry Smith #include "pc.h"
82eac72dbSBarry Smith 
99e25ed09SBarry Smith #define KSP_COOKIE  PETSC_COOKIE+8
10f0479e8cSBarry Smith 
11e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
122eac72dbSBarry Smith 
1382bf6240SBarry Smith #define KSPRICHARDSON "richardson"
1482bf6240SBarry Smith #define KSPCHEBYCHEV  "chebychev"
1582bf6240SBarry Smith #define KSPCG         "cg"
1682bf6240SBarry Smith #define KSPGMRES      "gmres"
1782bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
1882bf6240SBarry Smith #define KSPBCGS       "bcgs"
1982bf6240SBarry Smith #define KSPCGS        "cgs"
2082bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
2182bf6240SBarry Smith #define KSPCR         "cr"
2282bf6240SBarry Smith #define KSPLSQR       "lsqr"
2382bf6240SBarry Smith #define KSPPREONLY    "preonly"
2482bf6240SBarry Smith #define KSPQCG        "qcg"
25c9cf9b11SBarry Smith #define KSPBICG       "bicg"
26c38d4ed2SBarry Smith #define KSPFGMRES     "fgmres"
2782bf6240SBarry Smith typedef char * KSPType;
282eac72dbSBarry Smith 
296b5873e3SBarry Smith extern int KSPCreate(MPI_Comm,KSP *);
304b0e389bSBarry Smith extern int KSPSetType(KSP,KSPType);
318ed539a5SBarry Smith extern int KSPSetUp(KSP);
328ed539a5SBarry Smith extern int KSPSolve(KSP,int *);
337c922b88SBarry Smith extern int KSPSolveTranspose(KSP,int *);
348ed539a5SBarry Smith extern int KSPDestroy(KSP);
352eac72dbSBarry Smith 
36488ecbafSBarry Smith extern FList KSPList;
3782bf6240SBarry Smith extern int KSPRegisterAll(char *);
38cf256101SBarry Smith extern int KSPRegisterDestroy(void);
392eac72dbSBarry Smith 
40f1af5d2fSBarry Smith extern int KSPRegister(char*,char*,char*,int(*)(KSP));
41aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
42f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
436df38c32SLois Curfman McInnes #else
44f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
456df38c32SLois Curfman McInnes #endif
4682bf6240SBarry Smith 
4782bf6240SBarry Smith extern int KSPGetType(KSP,KSPType *);
48aad2872bSLois Curfman McInnes extern int KSPSetPreconditionerSide(KSP,PCSide);
49aad2872bSLois Curfman McInnes extern int KSPGetPreconditionerSide(KSP,PCSide*);
50ae622e99SLois Curfman McInnes extern int KSPGetTolerances(KSP,double*,double*,double*,int*);
51f9a5357eSLois Curfman McInnes extern int KSPSetTolerances(KSP,double,double,double,int);
52d4fbbf0eSBarry Smith extern int KSPSetComputeResidual(KSP,PetscTruth);
538ed539a5SBarry Smith extern int KSPSetUsePreconditionedResidual(KSP);
54c39b5c85SLois Curfman McInnes extern int KSPSetInitialGuessNonzero(KSP);
55af85ccbcSBarry Smith extern int KSPGetInitialGuessNonzero(KSP,PetscTruth *);
5690f02eecSBarry Smith extern int KSPSetComputeEigenvalues(KSP);
57d4fbbf0eSBarry Smith extern int KSPSetComputeSingularValues(KSP);
588ed539a5SBarry Smith extern int KSPSetRhs(KSP,Vec);
598ed539a5SBarry Smith extern int KSPGetRhs(KSP,Vec *);
608ed539a5SBarry Smith extern int KSPSetSolution(KSP,Vec);
618ed539a5SBarry Smith extern int KSPGetSolution(KSP,Vec *);
6225fce39dSBarry Smith extern int KSPGetResidualNorm(KSP,double*);
63597479d4SBarry Smith extern int KSPGetIterationNumber(KSP,int*);
642eac72dbSBarry Smith 
6577c4ece6SBarry Smith extern int KSPSetPC(KSP,PC);
6677c4ece6SBarry Smith extern int KSPGetPC(KSP,PC*);
67aabeff55SBarry Smith 
68fcb7de8cSBarry Smith extern int KSPSetAvoidNorms(KSP);
69fcb7de8cSBarry Smith 
70b8a78c4aSBarry Smith extern int KSPSetMonitor(KSP,int (*)(KSP,int,double,void*),void *,int (*)(void*));
715cd90555SBarry Smith extern int KSPClearMonitor(KSP);
728ed539a5SBarry Smith extern int KSPGetMonitorContext(KSP,void **);
730462333dSBarry Smith extern int KSPGetResidualHistory(KSP,double **,int *);
742c9e95e8SBarry Smith extern int KSPSetResidualHistory(KSP,double *,int,PetscTruth);
754b0e389bSBarry Smith 
762eac72dbSBarry Smith 
778ed539a5SBarry Smith extern int KSPBuildSolution(KSP,Vec,Vec *);
788ed539a5SBarry Smith extern int KSPBuildResidual(KSP,Vec,Vec,Vec *);
792eac72dbSBarry Smith 
808ed539a5SBarry Smith extern int KSPRichardsonSetScale(KSP,double);
818ed539a5SBarry Smith extern int KSPChebychevSetEigenvalues(KSP,double,double);
82d4fbbf0eSBarry Smith extern int KSPComputeExtremeSingularValues(KSP,double*,double*);
8388d459dfSBarry Smith extern int KSPComputeEigenvalues(KSP,int,double*,double*,int *);
84d4fbbf0eSBarry Smith extern int KSPComputeEigenvaluesExplicitly(KSP,int,double*,double*);
854b0e389bSBarry Smith 
868ed539a5SBarry Smith extern int KSPGMRESSetRestart(KSP,int);
87b4fd4287SBarry Smith extern int KSPGMRESSetPreAllocateVectors(KSP);
8877c4ece6SBarry Smith extern int KSPGMRESSetOrthogonalization(KSP,int (*)(KSP,int));
89bcd2baecSBarry Smith extern int KSPGMRESUnmodifiedGramSchmidtOrthogonalization(KSP,int);
90bcd2baecSBarry Smith extern int KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int);
91bcd2baecSBarry Smith extern int KSPGMRESIROrthogonalization(KSP,int);
9208480c60SBarry Smith 
93*329f5518SBarry Smith extern int KSPFGMRESModifyPCNoChange(KSP,int,int,double,void*);
94*329f5518SBarry Smith extern int KSPFGMRESModifyPCSLES(KSP,int,int,double,void*);
95*329f5518SBarry Smith extern int KSPFGMRESSetModifyPC(KSP,int (*)(KSP,int,int,double,void*),void*,int(*)(void*));
96c38d4ed2SBarry Smith 
978ed539a5SBarry Smith extern int KSPSetFromOptions(KSP);
9815091d37SBarry Smith extern int KSPSetTypeFromOptions(KSP);
99639f9d9dSBarry Smith extern int KSPAddOptionsChecker(int (*)(KSP));
1002eac72dbSBarry Smith 
1018f4c8dbaSBarry Smith extern int KSPSingularValueMonitor(KSP,int,double,void *);
1028ed539a5SBarry Smith extern int KSPDefaultMonitor(KSP,int,double,void *);
103af6b99e9SBarry Smith extern int KSPTrueMonitor(KSP,int,double,void *);
1048c48e012SBarry Smith extern int KSPDefaultSMonitor(KSP,int,double,void *);
105a302ce19SBarry Smith extern int KSPVecViewMonitor(KSP,int,double,void *);
10621988a88SBarry Smith extern int KSPGMRESKrylovMonitor(KSP,int,double,void *);
10784cb2905SBarry Smith 
1082eac72dbSBarry Smith 
109c01c455dSBarry Smith extern int KSPResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec);
11077c4ece6SBarry Smith extern int KSPUnwindPreconditioner(KSP,Vec,Vec);
111c01c455dSBarry Smith extern int KSPDefaultBuildSolution(KSP,Vec,Vec*);
112c01c455dSBarry Smith extern int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
113c01c455dSBarry Smith 
1148ed539a5SBarry Smith extern int KSPPrintHelp(KSP);
1158ed539a5SBarry Smith 
1161eb62cbbSBarry Smith extern int KSPSetOptionsPrefix(KSP,char*);
1170985bb82SSatish Balay extern int KSPAppendOptionsPrefix(KSP,char*);
1180985bb82SSatish Balay extern int KSPGetOptionsPrefix(KSP,char**);
1191eb62cbbSBarry Smith 
1208ed539a5SBarry Smith extern int KSPView(KSP,Viewer);
1211eb62cbbSBarry Smith 
122d15094e1SBarry Smith typedef enum {/* converged */
123d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
124d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
125b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
126*329f5518SBarry Smith               KSP_CONVERGED_QCG_NEGATIVE_CURVE =  5,
127*329f5518SBarry Smith               KSP_CONVERGED_QCG_CONSTRAINED    =  6,
128*329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
129d15094e1SBarry Smith               /* diverged */
130d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
131d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
132d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
133d15094e1SBarry Smith 
134d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
135d15094e1SBarry Smith 
136abef13c0SSatish Balay extern int KSPSetConvergenceTest(KSP,int (*)(KSP,int,double,KSPConvergedReason*,void*),void *);
137abef13c0SSatish Balay extern int KSPGetConvergenceContext(KSP,void **);
138abef13c0SSatish Balay extern int KSPDefaultConverged(KSP,int,double,KSPConvergedReason*,void *);
139abef13c0SSatish Balay extern int KSPSkipConverged(KSP,int,double,KSPConvergedReason*,void *);
140*329f5518SBarry Smith extern int KSPGetConvergedReason(KSP,KSPConvergedReason *);
141abef13c0SSatish Balay 
142d4fbbf0eSBarry Smith extern int KSPComputeExplicitOperator(KSP,Mat *);
143d4fbbf0eSBarry Smith 
1446d4a8577SBarry Smith typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType;
1456d4a8577SBarry Smith extern int KSPCGSetType(KSP,KSPCGType);
146e559a7feSLois Curfman McInnes 
1473369ce9aSBarry Smith extern int PCPreSolve(PC,KSP);
1483369ce9aSBarry Smith extern int PCPostSolve(PC,KSP);
1493369ce9aSBarry Smith 
150d7e8b826SBarry Smith extern int KSPLGMonitorCreate(char*,char*,int,int,int,int,DrawLG*);
1511eb62cbbSBarry Smith extern int KSPLGMonitor(KSP,int,double,void*);
152d7e8b826SBarry Smith extern int KSPLGMonitorDestroy(DrawLG);
1533b2fbd54SBarry Smith extern int KSPLGTrueMonitorCreate(MPI_Comm,char*,char*,int,int,int,int,DrawLG*);
154af6b99e9SBarry Smith extern int KSPLGTrueMonitor(KSP,int,double,void*);
155af6b99e9SBarry Smith extern int KSPLGTrueMonitorDestroy(DrawLG);
1561eb62cbbSBarry Smith 
1572eac72dbSBarry Smith #endif
1582eac72dbSBarry Smith 
1592eac72dbSBarry Smith 
160