xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/petscksp.pxi (revision 3220ff8572602716d60bdda8b51773ebceb3c8ea)
1cdef extern from * nogil:
2
3    ctypedef const char* PetscKSPType "KSPType"
4    PetscKSPType KSPRICHARDSON
5    PetscKSPType KSPCHEBYSHEV
6    PetscKSPType KSPCG
7    PetscKSPType KSPGROPPCG
8    PetscKSPType KSPPIPECG
9    PetscKSPType KSPPIPECGRR
10    PetscKSPType KSPPIPELCG
11    PetscKSPType KSPPIPEPRCG
12    PetscKSPType KSPPIPECG2
13    PetscKSPType KSPCGNE
14    PetscKSPType KSPNASH
15    PetscKSPType KSPSTCG
16    PetscKSPType KSPGLTR
17    PetscKSPType KSPFCG
18    PetscKSPType KSPPIPEFCG
19    PetscKSPType KSPGMRES
20    PetscKSPType KSPPIPEFGMRES
21    PetscKSPType   KSPFGMRES
22    PetscKSPType   KSPLGMRES
23    PetscKSPType   KSPDGMRES
24    PetscKSPType   KSPPGMRES
25    PetscKSPType KSPTCQMR
26    PetscKSPType KSPBCGS
27    PetscKSPType   KSPIBCGS
28    PetscKSPType   KSPQMRCGS
29    PetscKSPType   KSPFBCGS
30    PetscKSPType   KSPFBCGSR
31    PetscKSPType   KSPBCGSL
32    PetscKSPType   KSPPIPEBCGS
33    PetscKSPType KSPCGS
34    PetscKSPType KSPTFQMR
35    PetscKSPType KSPCR
36    PetscKSPType KSPPIPECR
37    PetscKSPType KSPLSQR
38    PetscKSPType KSPPREONLY
39    PetscKSPType   KSPNONE
40    PetscKSPType KSPQCG
41    PetscKSPType KSPBICG
42    PetscKSPType KSPMINRES
43    PetscKSPType KSPSYMMLQ
44    PetscKSPType KSPLCD
45    PetscKSPType KSPPYTHON
46    PetscKSPType KSPGCR
47    PetscKSPType KSPPIPEGCR
48    PetscKSPType KSPTSIRM
49    PetscKSPType KSPCGLS
50    PetscKSPType KSPFETIDP
51    PetscKSPType KSPHPDDM
52
53    ctypedef enum PetscKSPNormType "KSPNormType":
54        KSP_NORM_DEFAULT
55        KSP_NORM_NONE
56        KSP_NORM_PRECONDITIONED
57        KSP_NORM_UNPRECONDITIONED
58        KSP_NORM_NATURAL
59
60    ctypedef enum PetscKSPConvergedReason "KSPConvergedReason":
61        # iterating
62        KSP_CONVERGED_ITERATING
63        # converged
64        KSP_CONVERGED_RTOL_NORMAL_EQUATIONS
65        KSP_CONVERGED_ATOL_NORMAL_EQUATIONS
66        KSP_CONVERGED_RTOL
67        KSP_CONVERGED_ATOL
68        KSP_CONVERGED_ITS
69        KSP_CONVERGED_NEG_CURVE
70        KSP_CONVERGED_STEP_LENGTH
71        KSP_CONVERGED_HAPPY_BREAKDOWN
72        # diverged
73        KSP_DIVERGED_NULL
74        KSP_DIVERGED_MAX_IT "KSP_DIVERGED_ITS"
75        KSP_DIVERGED_DTOL
76        KSP_DIVERGED_BREAKDOWN
77        KSP_DIVERGED_BREAKDOWN_BICG
78        KSP_DIVERGED_NONSYMMETRIC
79        KSP_DIVERGED_INDEFINITE_PC
80        KSP_DIVERGED_NANORINF
81        KSP_DIVERGED_INDEFINITE_MAT
82        KSP_DIVERGED_PC_FAILED
83
84    ctypedef enum PetscKSPDMActive "KSPDMActive":
85        KSP_DMACTIVE_OPERATOR = 1
86        KSP_DMACTIVE_RHS = 2
87        KSP_DMACTIVE_INITIAL_GUESS = 4
88        KSP_DMACTIVE_ALL = 7
89
90    ctypedef enum PetscKSPHPDDMType "KSPHPDDMType":
91        KSP_HPDDM_TYPE_GMRES
92        KSP_HPDDM_TYPE_BGMRES
93        KSP_HPDDM_TYPE_CG
94        KSP_HPDDM_TYPE_BCG
95        KSP_HPDDM_TYPE_GCRODR
96        KSP_HPDDM_TYPE_BGCRODR
97        KSP_HPDDM_TYPE_BFBCG
98        KSP_HPDDM_TYPE_PREONLY
99
100    ctypedef PetscErrorCode (*PetscKSPCtxDel)(void*)
101
102    ctypedef PetscErrorCode (*PetscKSPConvergedFunction)(PetscKSP,
103                                                         PetscInt,
104                                                         PetscReal,
105                                                         PetscKSPConvergedReason*,
106                                                         void*) except PETSC_ERR_PYTHON
107
108    ctypedef PetscErrorCode (*PetscKSPMonitorFunction)(PetscKSP,
109                                                       PetscInt,
110                                                       PetscReal,
111                                                       void*) except PETSC_ERR_PYTHON
112
113    ctypedef PetscErrorCode (*PetscKSPComputeRHSFunction)(PetscKSP,
114                                                          PetscVec,
115                                                          void*) except PETSC_ERR_PYTHON
116
117    ctypedef PetscErrorCode (*PetscKSPComputeOpsFunction)(PetscKSP,
118                                                          PetscMat,
119                                                          PetscMat,
120                                                          void*) except PETSC_ERR_PYTHON
121
122    ctypedef PetscErrorCode (*PetscKSPPreSolveFunction)(PetscKSP,
123                                                        PetscVec,
124                                                        PetscVec,
125                                                        void*) except PETSC_ERR_PYTHON
126
127    ctypedef PetscErrorCode (*PetscKSPPostSolveFunction)(PetscKSP,
128                                                         PetscVec,
129                                                         PetscVec,
130                                                         void*) except PETSC_ERR_PYTHON
131
132    PetscErrorCode KSPCreate(MPI_Comm, PetscKSP*)
133    PetscErrorCode KSPDestroy(PetscKSP*)
134    PetscErrorCode KSPView(PetscKSP, PetscViewer)
135
136    PetscErrorCode KSPSetType(PetscKSP, PetscKSPType)
137    PetscErrorCode KSPGetType(PetscKSP, PetscKSPType*)
138
139    PetscErrorCode KSPSetOptionsPrefix(PetscKSP, char[])
140    PetscErrorCode KSPAppendOptionsPrefix(PetscKSP, char[])
141    PetscErrorCode KSPGetOptionsPrefix(PetscKSP, char*[])
142    PetscErrorCode KSPSetFromOptions(PetscKSP)
143
144    PetscErrorCode KSPSetTolerances(PetscKSP, PetscReal, PetscReal, PetscReal, PetscInt)
145    PetscErrorCode KSPGetTolerances(PetscKSP, PetscReal*, PetscReal*, PetscReal*, PetscInt*)
146    PetscErrorCode KSPSetNormType(PetscKSP, PetscKSPNormType)
147    PetscErrorCode KSPGetNormType(PetscKSP, PetscKSPNormType*)
148    PetscErrorCode KSPSetPCSide(PetscKSP, PetscPCSide)
149    PetscErrorCode KSPGetPCSide(PetscKSP, PetscPCSide*)
150    PetscErrorCode KSPSetSupportedNorm(PetscKSP, PetscKSPNormType, PetscPCSide, PetscInt)
151
152    PetscErrorCode KSPSetConvergenceTest(PetscKSP, PetscKSPConvergedFunction, void*, PetscKSPCtxDel)
153    PetscErrorCode KSPSetResidualHistory(PetscKSP, PetscReal[], PetscInt, PetscBool)
154    PetscErrorCode KSPGetResidualHistory(PetscKSP, PetscReal*[], PetscInt*)
155    PetscErrorCode KSPLogResidualHistory(PetscKSP, PetscReal)
156    PetscErrorCode KSPConvergedDefaultCreate(void**)
157    PetscErrorCode KSPConvergedDefaultDestroy(void*)
158    PetscErrorCode KSPConvergedDefault(PetscKSP, PetscInt, PetscReal, PetscKSPConvergedReason*, void*) except PETSC_ERR_PYTHON
159    PetscErrorCode KSPLSQRConvergedDefault(PetscKSP, PetscInt, PetscReal, PetscKSPConvergedReason*, void*) except PETSC_ERR_PYTHON
160    PetscErrorCode KSPConvergedSkip(PetscKSP, PetscInt, PetscReal, PetscKSPConvergedReason*, void*) except PETSC_ERR_PYTHON
161
162    PetscErrorCode KSPMonitorSet(PetscKSP, PetscKSPMonitorFunction, void*, PetscKSPCtxDel)
163    PetscErrorCode KSPMonitorCancel(PetscKSP)
164    PetscErrorCode KSPMonitor(PetscKSP, PetscInt, PetscReal)
165
166    PetscErrorCode KSPSetInitialGuessNonzero(PetscKSP, PetscBool)
167    PetscErrorCode KSPGetInitialGuessNonzero(PetscKSP, PetscBool*)
168    PetscErrorCode KSPSetInitialGuessKnoll(PetscKSP, PetscBool)
169    PetscErrorCode KSPGetInitialGuessKnoll(PetscKSP, PetscBool*)
170    PetscErrorCode KSPSetUseFischerGuess(PetscKSP, PetscInt, PetscInt)
171
172    PetscErrorCode KSPGetComputeEigenvalues(PetscKSP, PetscBool*)
173    PetscErrorCode KSPSetComputeEigenvalues(PetscKSP, PetscBool)
174    PetscErrorCode KSPGetComputeSingularValues(PetscKSP, PetscBool*)
175    PetscErrorCode KSPSetComputeSingularValues(PetscKSP, PetscBool)
176
177    PetscErrorCode KSPSetComputeRHS(PetscKSP, PetscKSPComputeRHSFunction, void*)
178    PetscErrorCode KSPSetComputeOperators(PetscKSP, PetscKSPComputeOpsFunction, void*)
179    PetscErrorCode KSPSetOperators(PetscKSP, PetscMat, PetscMat)
180    PetscErrorCode KSPGetOperators(PetscKSP, PetscMat*, PetscMat*)
181    PetscErrorCode KSPGetOperatorsSet(PetscKSP, PetscBool*, PetscBool*)
182
183    PetscErrorCode KSPSetPC(PetscKSP, PetscPC)
184    PetscErrorCode KSPGetPC(PetscKSP, PetscPC*)
185
186    PetscErrorCode KSPGetDM(PetscKSP, PetscDM*)
187    PetscErrorCode KSPSetDM(PetscKSP, PetscDM)
188    PetscErrorCode KSPSetDMActive(PetscKSP, PetscKSPDMActive, PetscBool)
189
190    PetscErrorCode KSPSetUp(PetscKSP)
191    PetscErrorCode KSPReset(PetscKSP)
192    PetscErrorCode KSPSetUpOnBlocks(PetscKSP)
193    PetscErrorCode KSPSetPreSolve(PetscKSP, PetscKSPPreSolveFunction, void*)
194    PetscErrorCode KSPSetPostSolve(PetscKSP, PetscKSPPostSolveFunction, void*)
195    PetscErrorCode KSPSolve(PetscKSP, PetscVec, PetscVec)
196    PetscErrorCode KSPSolveTranspose(PetscKSP, PetscVec, PetscVec)
197    PetscErrorCode KSPMatSolve(PetscKSP, PetscMat, PetscMat)
198    PetscErrorCode KSPMatSolveTranspose(PetscKSP, PetscMat, PetscMat)
199
200    PetscErrorCode KSPGetRhs(PetscKSP, PetscVec*)
201    PetscErrorCode KSPGetSolution(PetscKSP, PetscVec*)
202    PetscErrorCode KSPGetConvergedReason(PetscKSP, PetscKSPConvergedReason*)
203    PetscErrorCode KSPGetIterationNumber(PetscKSP, PetscInt*)
204    PetscErrorCode KSPGetResidualNorm(PetscKSP, PetscReal*)
205    PetscErrorCode KSPSetErrorIfNotConverged(PetscKSP, PetscBool)
206    PetscErrorCode KSPGetErrorIfNotConverged(PetscKSP, PetscBool*)
207
208    PetscErrorCode KSPBuildSolution(PetscKSP, PetscVec, PetscVec*)
209    PetscErrorCode KSPBuildSolutionDefault(PetscKSP, PetscVec, PetscVec*)
210    PetscErrorCode KSPBuildResidual(PetscKSP, PetscVec, PetscVec, PetscVec*)
211    PetscErrorCode KSPBuildResidualDefault(PetscKSP, PetscVec, PetscVec, PetscVec*)
212
213    PetscErrorCode KSPSetDiagonalScale(PetscKSP, PetscBool)
214    PetscErrorCode KSPGetDiagonalScale(PetscKSP, PetscBool*)
215    PetscErrorCode KSPSetDiagonalScaleFix(PetscKSP, PetscBool)
216    PetscErrorCode KSPGetDiagonalScaleFix(PetscKSP, PetscBool*)
217
218    PetscErrorCode KSPComputeExplicitOperator(PetscKSP, PetscMat*)
219    PetscErrorCode KSPComputeEigenvalues(PetscKSP, PetscInt, PetscReal[], PetscReal[], PetscInt*)
220    PetscErrorCode KSPComputeExtremeSingularValues(PetscKSP, PetscReal*, PetscReal*)
221
222    PetscErrorCode KSPCreateVecs(PetscKSP, PetscInt, PetscVec**, PetscInt, PetscVec**)
223
224    PetscErrorCode KSPGMRESSetRestart(PetscKSP, PetscInt)
225
226    PetscErrorCode KSPCGGetObjFcn(PetscKSP, PetscReal*)
227
228    PetscErrorCode KSPPythonSetType(PetscKSP, char[])
229    PetscErrorCode KSPPythonGetType(PetscKSP, char*[])
230
231    PetscErrorCode KSPHPDDMSetType(PetscKSP, PetscKSPHPDDMType)
232    PetscErrorCode KSPHPDDMGetType(PetscKSP, PetscKSPHPDDMType*)
233
234cdef extern from * nogil: # custom.h
235    PetscErrorCode KSPSetIterationNumber(PetscKSP, PetscInt)
236    PetscErrorCode KSPSetResidualNorm(PetscKSP, PetscReal)
237    PetscErrorCode KSPConvergenceTestCall(PetscKSP, PetscInt, PetscReal, PetscKSPConvergedReason*)
238    PetscErrorCode KSPSetConvergedReason(PetscKSP, PetscKSPConvergedReason)
239    PetscErrorCode KSPAddConvergenceTest(PetscKSP, PetscKSPConvergedFunction, PetscBool)
240
241# -----------------------------------------------------------------------------
242
243cdef inline KSP ref_KSP(PetscKSP ksp):
244    cdef KSP ob = <KSP> KSP()
245    ob.ksp = ksp
246    CHKERR(PetscINCREF(ob.obj))
247    return ob
248
249# -----------------------------------------------------------------------------
250
251cdef PetscErrorCode KSP_Converged(
252    PetscKSP  ksp,
253    PetscInt  its,
254    PetscReal rnm,
255    PetscKSPConvergedReason *r,
256    void      *ctx,
257   ) except PETSC_ERR_PYTHON with gil:
258    cdef KSP Ksp = ref_KSP(ksp)
259    (converged, args, kargs) = Ksp.get_attr('__converged__')
260    reason = converged(Ksp, toInt(its), toReal(rnm), *args, **kargs)
261    if   reason is None:  r[0] = KSP_CONVERGED_ITERATING
262    elif reason is False: r[0] = KSP_CONVERGED_ITERATING
263    elif reason is True:  r[0] = KSP_CONVERGED_ITS # XXX ?
264    else:                 r[0] = reason
265    return PETSC_SUCCESS
266
267# -----------------------------------------------------------------------------
268
269cdef PetscErrorCode KSP_Monitor(
270    PetscKSP  ksp,
271    PetscInt  its,
272    PetscReal rnm,
273    void      *ctx,
274   ) except PETSC_ERR_PYTHON with gil:
275    cdef KSP Ksp = ref_KSP(ksp)
276    cdef object monitorlist = Ksp.get_attr('__monitor__')
277    if monitorlist is None: return PETSC_SUCCESS
278    for (monitor, args, kargs) in monitorlist:
279        monitor(Ksp, toInt(its), toReal(rnm), *args, **kargs)
280    return PETSC_SUCCESS
281
282# -----------------------------------------------------------------------------
283
284cdef PetscErrorCode KSP_ComputeRHS(
285    PetscKSP ksp,
286    PetscVec rhs,
287    void     *ctx,
288   ) except PETSC_ERR_PYTHON with gil:
289    cdef KSP Ksp = ref_KSP(ksp)
290    cdef Vec Rhs = ref_Vec(rhs)
291    cdef object context = Ksp.get_attr('__rhs__')
292    if context is None and ctx != NULL: context = <object>ctx
293    assert context is not None and type(context) is tuple # sanity check
294    (computerhs, args, kargs) = context
295    computerhs(Ksp, Rhs, *args, **kargs)
296    return PETSC_SUCCESS
297
298cdef PetscErrorCode KSP_ComputeOps(
299    PetscKSP ksp,
300    PetscMat A,
301    PetscMat B,
302    void     *ctx,
303   ) except PETSC_ERR_PYTHON with gil:
304    cdef KSP Ksp  = ref_KSP(ksp)
305    cdef Mat Amat = ref_Mat(A)
306    cdef Mat Bmat = ref_Mat(B)
307    cdef object context = Ksp.get_attr('__operators__')
308    if context is None and ctx != NULL: context = <object>ctx
309    assert context is not None and type(context) is tuple # sanity check
310    (computeops, args, kargs) = context
311    computeops(Ksp, Amat, Bmat, *args, **kargs)
312    return PETSC_SUCCESS
313
314# -----------------------------------------------------------------------------
315
316cdef PetscErrorCode KSP_PreSolve(
317    PetscKSP ksp,
318    PetscVec rhs,
319    PetscVec x,
320    void* ctx,
321   ) except PETSC_ERR_PYTHON with gil:
322    cdef KSP Ksp = ref_KSP(ksp)
323    cdef Vec Rhs = ref_Vec(rhs)
324    cdef Vec X = ref_Vec(x)
325    cdef object context = Ksp.get_attr('__presolve__')
326    if context is None and ctx != NULL: context = <object>ctx
327    assert context is not None and type(context) is tuple
328    (presolve, args, kargs) = context
329    presolve(Ksp, Rhs, X, *args, **kargs)
330    return PETSC_SUCCESS
331
332cdef PetscErrorCode KSP_PostSolve(
333    PetscKSP ksp,
334    PetscVec rhs,
335    PetscVec x,
336    void* ctx,
337   ) except PETSC_ERR_PYTHON with gil:
338    cdef KSP Ksp = ref_KSP(ksp)
339    cdef Vec Rhs = ref_Vec(rhs)
340    cdef Vec X = ref_Vec(x)
341    cdef object context = Ksp.get_attr('__postsolve__')
342    if context is None and ctx != NULL: context = <object>ctx
343    assert context is not None and type(context) is tuple
344    (postsolve, args, kargs) = context
345    postsolve(Ksp, Rhs, X, *args, **kargs)
346    return PETSC_SUCCESS
347
348# -----------------------------------------------------------------------------
349