xref: /petsc/include/petscksp.h (revision d8e2560e552b4ed8fbe5358ac1ca98b70eb70dd5)
1f26ada1bSBarry Smith /*
2f26ada1bSBarry Smith    Defines the interface functions for the Krylov subspace accelerators.
3f26ada1bSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCKSP_H
50a835dfdSSatish Balay #define __PETSCKSP_H
60a835dfdSSatish Balay #include "petscpc.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitializePackage(const char[]);
101dbb0a54SBarry Smith 
1128ce4d24SBarry Smith /*S
1228ce4d24SBarry Smith      KSP - Abstract PETSc object that manages all Krylov methods
1328ce4d24SBarry Smith 
1428ce4d24SBarry Smith    Level: beginner
1528ce4d24SBarry Smith 
1628ce4d24SBarry Smith   Concepts: Krylov methods
1728ce4d24SBarry Smith 
1894b7f48cSBarry Smith .seealso:  KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP
1928ce4d24SBarry Smith S*/
20e2a1c21fSSatish Balay typedef struct _p_KSP*     KSP;
212eac72dbSBarry Smith 
2228ce4d24SBarry Smith /*E
2328ce4d24SBarry Smith     KSPType - String with the name of a PETSc Krylov method or the creation function
2428ce4d24SBarry Smith        with an optional dynamic library name, for example
2528ce4d24SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mykspcreate()
2628ce4d24SBarry Smith 
2728ce4d24SBarry Smith    Level: beginner
2828ce4d24SBarry Smith 
2928ce4d24SBarry Smith .seealso: KSPSetType(), KSP
3028ce4d24SBarry Smith E*/
31e5a9bf91SBarry Smith #define KSPType const char*
3282bf6240SBarry Smith #define KSPRICHARDSON "richardson"
3382bf6240SBarry Smith #define KSPCHEBYCHEV  "chebychev"
3482bf6240SBarry Smith #define KSPCG         "cg"
35df2a969aSvictorle #define   KSPCGNE       "cgne"
3680e17431SMatthew Knepley #define   KSPSTCG       "stcg"
378031f4adStmunson #define   KSPGLTR       "gltr"
3882bf6240SBarry Smith #define KSPGMRES      "gmres"
399dcbbd2bSBarry Smith #define   KSPFGMRES     "fgmres"
409dcbbd2bSBarry Smith #define   KSPLGMRES     "lgmres"
4182bf6240SBarry Smith #define KSPTCQMR      "tcqmr"
4282bf6240SBarry Smith #define KSPBCGS       "bcgs"
43f0037002Svictorle #define KSPBCGSL      "bcgsl"
4482bf6240SBarry Smith #define KSPCGS        "cgs"
4582bf6240SBarry Smith #define KSPTFQMR      "tfqmr"
4682bf6240SBarry Smith #define KSPCR         "cr"
4782bf6240SBarry Smith #define KSPLSQR       "lsqr"
4882bf6240SBarry Smith #define KSPPREONLY    "preonly"
4982bf6240SBarry Smith #define KSPQCG        "qcg"
50c9cf9b11SBarry Smith #define KSPBICG       "bicg"
51b4ac9ba4SBarry Smith #define KSPMINRES     "minres"
5201c5daebSSatish Balay #define KSPSYMMLQ     "symmlq"
5321356919SSatish Balay #define KSPLCD        "lcd"
542eac72dbSBarry Smith 
558ba1e511SMatthew Knepley /* Logging support */
56dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE;
578ba1e511SMatthew Knepley 
58dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *);
59f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,KSPType);
60dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP);
61dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP);
62dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec);
63dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec);
64dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP);
652eac72dbSBarry Smith 
66b0a32e0cSBarry Smith extern PetscFList KSPList;
67dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]);
68dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void);
69dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP));
7030de9b25SBarry Smith 
7130de9b25SBarry Smith /*MC
7230de9b25SBarry Smith    KSPRegisterDynamic - Adds a method to the Krylov subspace solver package.
7330de9b25SBarry Smith 
7430de9b25SBarry Smith    Synopsis:
75d360dc6fSBarry Smith    PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP))
7630de9b25SBarry Smith 
7730de9b25SBarry Smith    Not Collective
7830de9b25SBarry Smith 
7930de9b25SBarry Smith    Input Parameters:
8030de9b25SBarry Smith +  name_solver - name of a new user-defined solver
8130de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
8230de9b25SBarry Smith .  name_create - name of routine to create method context
8330de9b25SBarry Smith -  routine_create - routine to create method context
8430de9b25SBarry Smith 
8530de9b25SBarry Smith    Notes:
8630de9b25SBarry Smith    KSPRegisterDynamic() may be called multiple times to add several user-defined solvers.
8730de9b25SBarry Smith 
8830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
8930de9b25SBarry Smith    is ignored.
9030de9b25SBarry Smith 
9130de9b25SBarry Smith    Sample usage:
9230de9b25SBarry Smith .vb
9330de9b25SBarry Smith    KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
9430de9b25SBarry Smith                "MySolverCreate",MySolverCreate);
9530de9b25SBarry Smith .ve
9630de9b25SBarry Smith 
9730de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
9830de9b25SBarry Smith $     KSPSetType(ksp,"my_solver")
9930de9b25SBarry Smith    or at runtime via the option
10030de9b25SBarry Smith $     -ksp_type my_solver
10130de9b25SBarry Smith 
10230de9b25SBarry Smith    Level: advanced
10330de9b25SBarry Smith 
104ab901514SBarry Smith    Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
10530de9b25SBarry Smith           and others of the form ${any_environmental_variable} occuring in pathname will be
10630de9b25SBarry Smith           replaced with appropriate values.
10730de9b25SBarry Smith          If your function is not being put into a shared library then use KSPRegister() instead
10830de9b25SBarry Smith 
10930de9b25SBarry Smith .keywords: KSP, register
11030de9b25SBarry Smith 
11130de9b25SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy()
11230de9b25SBarry Smith 
11330de9b25SBarry Smith M*/
114aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
115f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0)
1166df38c32SLois Curfman McInnes #else
117f1af5d2fSBarry Smith #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d)
1186df38c32SLois Curfman McInnes #endif
11982bf6240SBarry Smith 
120dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *);
121dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide);
122dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*);
123dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*);
124dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt);
125dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth);
126dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *);
127dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth);
128dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*);
129a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*);
130dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth);
131a6e7709dSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*);
132dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth);
133dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *);
134dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *);
135dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*);
136dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*);
137dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace);
138dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*);
139906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**);
1402eac72dbSBarry Smith 
141dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC);
142dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*);
143aabeff55SBarry Smith 
144a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*));
145a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorCancel(KSP);
146dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **);
147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *);
148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth);
1494b0e389bSBarry Smith 
1500e33f6ddSBarry Smith /* not sure where to put this */
151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*);
152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]);
154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]);
1552eac72dbSBarry Smith 
1560a71e23eSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC,KSP *);
1570a71e23eSMatthew Knepley 
158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *);
159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *);
1602eac72dbSBarry Smith 
161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal);
162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal);
163dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*);
164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*);
1664b0e389bSBarry Smith 
167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt);
168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal);
1699f236934SBarry Smith 
170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP);
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt));
172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt);
173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt);
1741d73ed98SBarry Smith 
175dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt);
176dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP);
1771d73ed98SBarry Smith 
178b49fd9e1SBarry Smith /*E
179b49fd9e1SBarry Smith     KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed.
180b49fd9e1SBarry Smith 
181b49fd9e1SBarry Smith    Level: advanced
182b49fd9e1SBarry Smith 
183b49fd9e1SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1848c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization()
185b49fd9e1SBarry Smith 
186b49fd9e1SBarry Smith E*/
18778d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType;
1889dcbbd2bSBarry Smith extern const char *KSPGMRESCGSRefinementTypes[];
1891f7e983dSSatish Balay /*MC
1908c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process
1918c5b8ba0SBarry Smith 
1928c5b8ba0SBarry Smith    Level: advanced
1938c5b8ba0SBarry Smith 
1948c5b8ba0SBarry Smith    Note: Possible unstable, but the fastest to compute
1958c5b8ba0SBarry Smith 
1968c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
1978c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
1988c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
1998c5b8ba0SBarry Smith M*/
2008c5b8ba0SBarry Smith 
2011f7e983dSSatish Balay /*MC
2028c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of
2038c5b8ba0SBarry Smith           iterative refinement if an estimate of the orthogonality of the resulting vectors indicates
2048c5b8ba0SBarry Smith           poor orthogonality.
2058c5b8ba0SBarry Smith 
2068c5b8ba0SBarry Smith    Level: advanced
2078c5b8ba0SBarry Smith 
2088c5b8ba0SBarry Smith    Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to
2098c5b8ba0SBarry Smith      estimate the orthogonality but is more stable.
2108c5b8ba0SBarry Smith 
2118c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2128c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS,
2138c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2148c5b8ba0SBarry Smith M*/
2158c5b8ba0SBarry Smith 
2161f7e983dSSatish Balay /*MC
2178c5b8ba0SBarry Smith     KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process.
2188c5b8ba0SBarry Smith 
2198c5b8ba0SBarry Smith    Level: advanced
2208c5b8ba0SBarry Smith 
2218c5b8ba0SBarry Smith    Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice
2228c5b8ba0SBarry Smith      but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED.
2238c5b8ba0SBarry Smith 
2248c5b8ba0SBarry Smith         You should only use this if you absolutely know that the iterative refinement is needed.
2258c5b8ba0SBarry Smith 
2268c5b8ba0SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(),
2278c5b8ba0SBarry Smith           KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS,
2288c5b8ba0SBarry Smith           KSPGMRESModifiedGramSchmidtOrthogonalization()
2298c5b8ba0SBarry Smith M*/
2308c5b8ba0SBarry Smith 
231dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType);
23208480c60SBarry Smith 
233dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*);
234dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*);
235dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*));
236c38d4ed2SBarry Smith 
237dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal);
238dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*);
239dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*);
240121fd945SKris Buschelman 
241d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetXRes(KSP,PetscReal);
242d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetPol(KSP,PetscTruth);
243d9492815SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBCGSLSetEll(KSP,int);
244d9492815SBarry Smith 
245dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP);
246dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP));
2472eac72dbSBarry Smith 
248a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *);
249a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefault(KSP,PetscInt,PetscReal,void *);
2507517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *);
251a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *);
252a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorSolution(KSP,PetscInt,PetscReal,void *);
253a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *);
25484cb2905SBarry Smith 
255dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec);
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*);
257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
258dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
259c01c455dSBarry Smith 
260dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure);
261dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*);
262906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperatorsSet(KSP,PetscTruth*,PetscTruth*);
263dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]);
264dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]);
2652dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]);
2661eb62cbbSBarry Smith 
267dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth);
268dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*);
269dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth);
270dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*);
2711f7f0c4fSBarry Smith 
272dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer);
2731eb62cbbSBarry Smith 
274db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRSetStandardErrorVec(KSP,Vec);
275db9b2ab1SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLSQRGetStandardErrorVec(KSP,Vec*);
276db9b2ab1SHong Zhang 
27728ce4d24SBarry Smith /*E
2788a4b9c5cSBarry Smith     KSPNormType - Norm that is passed in the Krylov convergence
2798a4b9c5cSBarry Smith        test routines.
2808a4b9c5cSBarry Smith 
2818a4b9c5cSBarry Smith    Level: advanced
2828a4b9c5cSBarry Smith 
2838a4b9c5cSBarry Smith    Notes: this must match finclude/petscksp.h
2848a4b9c5cSBarry Smith 
28594b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(),
2868a4b9c5cSBarry Smith           KSPSetConvergenceTest()
2878a4b9c5cSBarry Smith E*/
288ce9499c7SBarry Smith typedef enum {KSP_NORM_NO = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType;
2899dcbbd2bSBarry Smith extern const char *KSPNormTypes[];
2901f7e983dSSatish Balay /*MC
291ce9499c7SBarry Smith     KSP_NORM_NO - Do not compute a norm during the Krylov process. This will
2928c5b8ba0SBarry Smith           possibly save some computation but means the convergence test cannot
2938c5b8ba0SBarry Smith           be based on a norm of a residual etc.
2948c5b8ba0SBarry Smith 
2958c5b8ba0SBarry Smith    Level: advanced
2968c5b8ba0SBarry Smith 
297085a36d4SBarry Smith     Note: Some Krylov methods need to compute a residual norm and then this option is ignored
2988c5b8ba0SBarry Smith 
299ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL
3008c5b8ba0SBarry Smith M*/
3018c5b8ba0SBarry Smith 
3021f7e983dSSatish Balay /*MC
303ce9499c7SBarry Smith     KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the
3048c5b8ba0SBarry Smith        convergence test routine.
3058c5b8ba0SBarry Smith 
3068c5b8ba0SBarry Smith    Level: advanced
3078c5b8ba0SBarry Smith 
308ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3098c5b8ba0SBarry Smith M*/
3108c5b8ba0SBarry Smith 
3111f7e983dSSatish Balay /*MC
312ce9499c7SBarry Smith     KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the
3138c5b8ba0SBarry Smith        convergence test routine.
3148c5b8ba0SBarry Smith 
3158c5b8ba0SBarry Smith    Level: advanced
3168c5b8ba0SBarry Smith 
317ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest()
3188c5b8ba0SBarry Smith M*/
3198c5b8ba0SBarry Smith 
3201f7e983dSSatish Balay /*MC
321ce9499c7SBarry Smith     KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the
3228c5b8ba0SBarry Smith        convergence test routine.
3238c5b8ba0SBarry Smith 
3248c5b8ba0SBarry Smith    Level: advanced
3258c5b8ba0SBarry Smith 
326ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NO, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest()
3278c5b8ba0SBarry Smith M*/
3288c5b8ba0SBarry Smith 
329dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType);
330863a724eSLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNormType(KSP,KSPNormType*);
3318a4b9c5cSBarry Smith 
3328a4b9c5cSBarry Smith /*E
33328ce4d24SBarry Smith     KSPConvergedReason - reason a Krylov method was said to
33428ce4d24SBarry Smith          have converged or diverged
33528ce4d24SBarry Smith 
33628ce4d24SBarry Smith    Level: beginner
33728ce4d24SBarry Smith 
33828ce4d24SBarry Smith    Notes: this must match finclude/petscksp.h
33928ce4d24SBarry Smith 
34086c02ca4SBarry Smith    Developer note: The string versions of these are in
3414b9ad928SBarry Smith      src/ksp/ksp/interface/itfunc.c called convergedreasons.
34250f1acb2SBarry Smith      If these enums are changed you must change those.
34386c02ca4SBarry Smith 
344c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances()
34528ce4d24SBarry Smith E*/
346d15094e1SBarry Smith typedef enum {/* converged */
347d15094e1SBarry Smith               KSP_CONVERGED_RTOL               =  2,
348d15094e1SBarry Smith               KSP_CONVERGED_ATOL               =  3,
349b335793eSSatish Balay               KSP_CONVERGED_ITS                =  4,
3508031f4adStmunson               KSP_CONVERGED_CG_NEG_CURVE       =  5,
3518031f4adStmunson               KSP_CONVERGED_CG_CONSTRAINED     =  6,
352329f5518SBarry Smith               KSP_CONVERGED_STEP_LENGTH        =  7,
353af3d8002SBarry Smith               KSP_CONVERGED_HAPPY_BREAKDOWN    =  8,
354d15094e1SBarry Smith               /* diverged */
355b3cc6726SBarry Smith               KSP_DIVERGED_NULL                = -2,
356d15094e1SBarry Smith               KSP_DIVERGED_ITS                 = -3,
357d15094e1SBarry Smith               KSP_DIVERGED_DTOL                = -4,
358d15094e1SBarry Smith               KSP_DIVERGED_BREAKDOWN           = -5,
359b4ac9ba4SBarry Smith               KSP_DIVERGED_BREAKDOWN_BICG      = -6,
360b4ac9ba4SBarry Smith               KSP_DIVERGED_NONSYMMETRIC        = -7,
361b4ac9ba4SBarry Smith               KSP_DIVERGED_INDEFINITE_PC       = -8,
3626aee1118SBarry Smith               KSP_DIVERGED_NAN                 = -9,
3636aee1118SBarry Smith               KSP_DIVERGED_INDEFINITE_MAT      = -10,
364d15094e1SBarry Smith 
365d15094e1SBarry Smith               KSP_CONVERGED_ITERATING          =  0} KSPConvergedReason;
3669dcbbd2bSBarry Smith extern const char **KSPConvergedReasons;
367d15094e1SBarry Smith 
368c838673bSBarry Smith /*MC
369c838673bSBarry Smith      KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b)
370c838673bSBarry Smith 
371c838673bSBarry Smith    Level: beginner
372c838673bSBarry Smith 
373c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
374c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
375c838673bSBarry Smith        2-norm of the residual for right preconditioning
376c838673bSBarry Smith 
377c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
378c838673bSBarry Smith 
379c838673bSBarry Smith M*/
380c838673bSBarry Smith 
381c838673bSBarry Smith /*MC
382c838673bSBarry Smith      KSP_CONVERGED_ATOL - norm(r) <= atol
383c838673bSBarry Smith 
384c838673bSBarry Smith    Level: beginner
385c838673bSBarry Smith 
386c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
387c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
388c838673bSBarry Smith        2-norm of the residual for right preconditioning
389c838673bSBarry Smith 
390c838673bSBarry Smith    Level: beginner
391c838673bSBarry Smith 
392c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
393c838673bSBarry Smith 
394c838673bSBarry Smith M*/
395c838673bSBarry Smith 
396c838673bSBarry Smith /*MC
397c838673bSBarry Smith      KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b)
398c838673bSBarry Smith 
399c838673bSBarry Smith    Level: beginner
400c838673bSBarry Smith 
401c838673bSBarry Smith    See KSPNormType and KSPSetNormType() for possible norms that may be used. By default
402c838673bSBarry Smith        for left preconditioning it is the 2-norm of the preconditioned residual, and the
403c838673bSBarry Smith        2-norm of the residual for right preconditioning
404c838673bSBarry Smith 
405c838673bSBarry Smith    Level: beginner
406c838673bSBarry Smith 
407c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
408c838673bSBarry Smith 
409c838673bSBarry Smith M*/
410c838673bSBarry Smith 
411c838673bSBarry Smith /*MC
412c838673bSBarry Smith      KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was
413c838673bSBarry Smith       reached
414c838673bSBarry Smith 
415c838673bSBarry Smith    Level: beginner
416c838673bSBarry Smith 
417c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
418c838673bSBarry Smith 
419c838673bSBarry Smith M*/
420c838673bSBarry Smith 
421c838673bSBarry Smith /*MC
422c838673bSBarry Smith      KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the
423c838673bSBarry Smith            preconditioner is applied.
424c838673bSBarry Smith 
425c838673bSBarry Smith 
426c838673bSBarry Smith    Level: beginner
427c838673bSBarry Smith 
428c838673bSBarry Smith 
429c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
430c838673bSBarry Smith 
431c838673bSBarry Smith M*/
432c838673bSBarry Smith 
433c838673bSBarry Smith /*MC
434c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the
435c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
436c838673bSBarry Smith 
437c838673bSBarry Smith    Level: beginner
438c838673bSBarry Smith 
439c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
440c838673bSBarry Smith 
441c838673bSBarry Smith M*/
442c838673bSBarry Smith 
443c838673bSBarry Smith /*MC
444c838673bSBarry Smith      KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the
445c838673bSBarry Smith           method could not continue to enlarge the Krylov space.
446c838673bSBarry Smith 
447c838673bSBarry Smith 
448c838673bSBarry Smith    Level: beginner
449c838673bSBarry Smith 
450c838673bSBarry Smith 
451c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
452c838673bSBarry Smith 
453c838673bSBarry Smith M*/
454c838673bSBarry Smith 
455c838673bSBarry Smith /*MC
456c838673bSBarry Smith      KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not
457c838673bSBarry Smith         symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry
458c838673bSBarry Smith 
459c838673bSBarry Smith    Level: beginner
460c838673bSBarry Smith 
461c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
462c838673bSBarry Smith 
463c838673bSBarry Smith M*/
464c838673bSBarry Smith 
465c838673bSBarry Smith /*MC
466c838673bSBarry Smith      KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both
467c838673bSBarry Smith         positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to
468c838673bSBarry Smith         be positive definite
469c838673bSBarry Smith 
470c838673bSBarry Smith    Level: beginner
471c838673bSBarry Smith 
4722401956bSBarry Smith      Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force
473c838673bSBarry Smith   the PCICC preconditioner to generate a positive definite preconditioner
474c838673bSBarry Smith 
475c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
476c838673bSBarry Smith 
477c838673bSBarry Smith M*/
478c838673bSBarry Smith 
479c838673bSBarry Smith /*MC
480c838673bSBarry Smith      KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason()
481c838673bSBarry Smith         while the KSPSolve() is still running.
482c838673bSBarry Smith 
483c838673bSBarry Smith    Level: beginner
484c838673bSBarry Smith 
485c838673bSBarry Smith .seealso:  KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances()
486c838673bSBarry Smith 
487c838673bSBarry Smith M*/
488c838673bSBarry Smith 
489dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *);
490dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **);
491dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
49294b02367SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUIRNorm(KSP);
4935f4984edSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConvergedSetUMIRNorm(KSP);
494dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *);
495dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *);
496abef13c0SSatish Balay 
497dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *);
498d4fbbf0eSBarry Smith 
49928ce4d24SBarry Smith /*E
50028ce4d24SBarry Smith     KSPCGType - Determines what type of CG to use
50128ce4d24SBarry Smith 
50228ce4d24SBarry Smith    Level: beginner
50328ce4d24SBarry Smith 
50428ce4d24SBarry Smith .seealso: KSPCGSetType()
50528ce4d24SBarry Smith E*/
5069dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType;
5079dcbbd2bSBarry Smith extern const char *KSPCGTypes[];
50828ce4d24SBarry Smith 
509dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType);
5108031f4adStmunson 
5110d2a7ff2SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGSetRadius(KSP,PetscReal);
51211f224b9Stmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetNormD(KSP,PetscReal *);
513*d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSTCGGetObjFcn(KSP,PetscReal *);
514e559a7feSLois Curfman McInnes 
5158031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRSetRadius(KSP,PetscReal);
5168031f4adStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetNormD(KSP,PetscReal *);
517*d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetObjFcn(KSP,PetscReal *);
518*d8e2560eStmunson EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGLTRGetMinEig(KSP,PetscReal *);
5198031f4adStmunson 
520dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP);
521dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP);
5223369ce9aSBarry Smith 
523a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
524a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLG(KSP,PetscInt,PetscReal,void*);
525a6570f20SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGDestroy(PetscDrawLG);
5267517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
5277517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*);
5287517059bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG);
5292f2e5d10SKris Buschelman 
5300338f08bSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
53103919abeSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec));
53203919abeSBarry Smith 
533e884886fSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT MatMFFDKSPMonitor(KSP,PetscInt,PetscReal,void *);
534e884886fSBarry Smith 
535e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
5362eac72dbSBarry Smith #endif
537