xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/petscts.pxi (revision 00946a4d4715aa5f9fdac01930be79cdd37d52b1)
1cdef extern from * nogil:
2
3    ctypedef const char* PetscTSType "TSType"
4    PetscTSType TSEULER
5    PetscTSType TSBEULER
6    PetscTSType TSBASICSYMPLECTIC
7    PetscTSType TSPSEUDO
8    PetscTSType TSCN
9    PetscTSType TSSUNDIALS
10    PetscTSType TSRK
11    PetscTSType TSPYTHON
12    PetscTSType TSTHETA
13    PetscTSType TSALPHA
14    PetscTSType TSALPHA2
15    PetscTSType TSGLLE
16    PetscTSType TSGLEE
17    PetscTSType TSSSP
18    PetscTSType TSARKIMEX
19    PetscTSType TSDIRK
20    PetscTSType TSROSW
21    PetscTSType TSEIMEX
22    PetscTSType TSMIMEX
23    PetscTSType TSBDF
24    PetscTSType TSRADAU5
25    PetscTSType TSMPRK
26    PetscTSType TSDISCGRAD
27
28    ctypedef enum PetscTSProblemType "TSProblemType":
29        TS_LINEAR
30        TS_NONLINEAR
31
32    ctypedef enum PetscTSEquationType "TSEquationType":
33        TS_EQ_UNSPECIFIED
34        TS_EQ_EXPLICIT
35        TS_EQ_ODE_EXPLICIT
36        TS_EQ_DAE_SEMI_EXPLICIT_INDEX1
37        TS_EQ_DAE_SEMI_EXPLICIT_INDEX2
38        TS_EQ_DAE_SEMI_EXPLICIT_INDEX3
39        TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI
40        TS_EQ_IMPLICIT
41        TS_EQ_ODE_IMPLICIT
42        TS_EQ_DAE_IMPLICIT_INDEX1
43        TS_EQ_DAE_IMPLICIT_INDEX2
44        TS_EQ_DAE_IMPLICIT_INDEX3
45        TS_EQ_DAE_IMPLICIT_INDEXHI
46
47    ctypedef enum PetscTSConvergedReason "TSConvergedReason":
48        # iterating
49        TS_CONVERGED_ITERATING
50        # converged
51        TS_CONVERGED_TIME
52        TS_CONVERGED_ITS
53        TS_CONVERGED_USER
54        TS_CONVERGED_EVENT
55        # diverged
56        TS_DIVERGED_NONLINEAR_SOLVE
57        TS_DIVERGED_STEP_REJECTED
58
59    ctypedef enum PetscTSExactFinalTimeOption "TSExactFinalTimeOption":
60        TS_EXACTFINALTIME_UNSPECIFIED
61        TS_EXACTFINALTIME_STEPOVER
62        TS_EXACTFINALTIME_INTERPOLATE
63        TS_EXACTFINALTIME_MATCHSTEP
64
65    ctypedef PetscErrorCode PetscTSCtxDel(void*)
66
67    ctypedef PetscErrorCode (*PetscTSFunctionFunction)(PetscTS,
68                                                       PetscReal,
69                                                       PetscVec,
70                                                       PetscVec,
71                                                       void*) except PETSC_ERR_PYTHON
72
73    ctypedef PetscErrorCode (*PetscTSJacobianFunction)(PetscTS,
74                                                       PetscReal,
75                                                       PetscVec,
76                                                       PetscMat,
77                                                       PetscMat,
78                                                       void*) except PETSC_ERR_PYTHON
79
80    ctypedef PetscErrorCode (*PetscTSIFunctionFunction)(PetscTS,
81                                                        PetscReal,
82                                                        PetscVec,
83                                                        PetscVec,
84                                                        PetscVec,
85                                                        void*) except PETSC_ERR_PYTHON
86    ctypedef PetscErrorCode (*PetscTSIJacobianFunction)(PetscTS,
87                                                        PetscReal,
88                                                        PetscVec,
89                                                        PetscVec,
90                                                        PetscReal,
91                                                        PetscMat,
92                                                        PetscMat,
93                                                        void*) except PETSC_ERR_PYTHON
94
95    ctypedef PetscErrorCode (*PetscTSIJacobianPFunction)(PetscTS,
96                                                         PetscReal,
97                                                         PetscVec,
98                                                         PetscVec,
99                                                         PetscReal,
100                                                         PetscMat,
101                                                         void*) except PETSC_ERR_PYTHON
102
103    ctypedef PetscErrorCode (*PetscTSI2FunctionFunction)(PetscTS,
104                                                         PetscReal,
105                                                         PetscVec,
106                                                         PetscVec,
107                                                         PetscVec,
108                                                         PetscVec,
109                                                         void*) except PETSC_ERR_PYTHON
110    ctypedef PetscErrorCode (*PetscTSI2JacobianFunction)(PetscTS,
111                                                         PetscReal,
112                                                         PetscVec,
113                                                         PetscVec,
114                                                         PetscVec,
115                                                         PetscReal,
116                                                         PetscReal,
117                                                         PetscMat,
118                                                         PetscMat,
119                                                         void*) except PETSC_ERR_PYTHON
120
121    ctypedef PetscErrorCode (*PetscTSMonitorFunction)(PetscTS,
122                                                      PetscInt,
123                                                      PetscReal,
124                                                      PetscVec,
125                                                      void*) except PETSC_ERR_PYTHON
126
127    ctypedef PetscErrorCode (*PetscTSPreStepFunction)  (PetscTS) except PETSC_ERR_PYTHON
128    ctypedef PetscErrorCode (*PetscTSPostStepFunction) (PetscTS) except PETSC_ERR_PYTHON
129
130    PetscErrorCode TSCreate(MPI_Comm comm, PetscTS*)
131    PetscErrorCode TSClone(PetscTS, PetscTS*)
132    PetscErrorCode TSDestroy(PetscTS*)
133    PetscErrorCode TSView(PetscTS, PetscViewer)
134    PetscErrorCode TSLoad(PetscTS, PetscViewer)
135
136    PetscErrorCode TSSetProblemType(PetscTS, PetscTSProblemType)
137    PetscErrorCode TSGetProblemType(PetscTS, PetscTSProblemType*)
138    PetscErrorCode TSSetEquationType(PetscTS, PetscTSEquationType)
139    PetscErrorCode TSGetEquationType(PetscTS, PetscTSEquationType*)
140    PetscErrorCode TSSetType(PetscTS, PetscTSType)
141    PetscErrorCode TSGetType(PetscTS, PetscTSType*)
142
143    PetscErrorCode TSSetOptionsPrefix(PetscTS, char[])
144    PetscErrorCode TSAppendOptionsPrefix(PetscTS, char[])
145    PetscErrorCode TSGetOptionsPrefix(PetscTS, char*[])
146    PetscErrorCode TSSetFromOptions(PetscTS)
147
148    PetscErrorCode TSSetSolution(PetscTS, PetscVec)
149    PetscErrorCode TSGetSolution(PetscTS, PetscVec*)
150    PetscErrorCode TS2SetSolution(PetscTS, PetscVec, PetscVec)
151    PetscErrorCode TS2GetSolution(PetscTS, PetscVec*, PetscVec*)
152
153    PetscErrorCode TSGetRHSFunction(PetscTS, PetscVec*, PetscTSFunctionFunction*, void*)
154    PetscErrorCode TSGetRHSJacobian(PetscTS, PetscMat*, PetscMat*, PetscTSJacobianFunction*, void**)
155    PetscErrorCode TSSetRHSFunction(PetscTS, PetscVec, PetscTSFunctionFunction, void*)
156    PetscErrorCode TSSetRHSJacobian(PetscTS, PetscMat, PetscMat, PetscTSJacobianFunction, void*)
157    PetscErrorCode TSSetIFunction(PetscTS, PetscVec, PetscTSIFunctionFunction, void*)
158    PetscErrorCode TSSetIJacobian(PetscTS, PetscMat, PetscMat, PetscTSIJacobianFunction, void*)
159    PetscErrorCode TSSetIJacobianP(PetscTS, PetscMat, PetscTSIJacobianPFunction, void*)
160    PetscErrorCode TSGetIFunction(PetscTS, PetscVec*, PetscTSIFunctionFunction*, void*)
161    PetscErrorCode TSGetIJacobian(PetscTS, PetscMat*, PetscMat*, PetscTSIJacobianFunction*, void**)
162    PetscErrorCode TSSetI2Function(PetscTS, PetscVec, PetscTSI2FunctionFunction, void*)
163    PetscErrorCode TSSetI2Jacobian(PetscTS, PetscMat, PetscMat, PetscTSI2JacobianFunction, void*)
164    PetscErrorCode TSGetI2Function(PetscTS, PetscVec*, PetscTSI2FunctionFunction*, void**)
165    PetscErrorCode TSGetI2Jacobian(PetscTS, PetscMat*, PetscMat*, PetscTSI2JacobianFunction*, void**)
166    PetscErrorCode TSRHSSplitSetIS(PetscTS, char[], PetscIS)
167    PetscErrorCode TSRHSSplitSetRHSFunction(PetscTS, char[], PetscVec, PetscTSFunctionFunction, void*)
168    PetscErrorCode TSRHSSplitSetIFunction(PetscTS, char[], PetscVec, PetscTSIFunctionFunction, void*)
169    PetscErrorCode TSRHSSplitSetIJacobian(PetscTS, char[], PetscMat, PetscMat, PetscTSIJacobianFunction, void*)
170
171    PetscErrorCode TSGetKSP(PetscTS, PetscKSP*)
172    PetscErrorCode TSGetSNES(PetscTS, PetscSNES*)
173
174    PetscErrorCode TSGetDM(PetscTS, PetscDM*)
175    PetscErrorCode TSSetDM(PetscTS, PetscDM)
176
177    PetscErrorCode TSComputeRHSFunction(PetscTS, PetscReal, PetscVec, PetscVec)
178    PetscErrorCode TSComputeRHSFunctionLinear(PetscTS, PetscReal, PetscVec, PetscVec, void*)
179    PetscErrorCode TSComputeRHSJacobian(PetscTS, PetscReal, PetscVec, PetscMat, PetscMat)
180    PetscErrorCode TSComputeRHSJacobianConstant(PetscTS, PetscReal, PetscVec, PetscMat, PetscMat, void*)
181    PetscErrorCode TSComputeIFunction(PetscTS, PetscReal, PetscVec, PetscVec, PetscVec, PetscBool)
182    PetscErrorCode TSComputeIJacobian(PetscTS, PetscReal, PetscVec, PetscVec, PetscReal, PetscMat, PetscMat, PetscBool)
183    PetscErrorCode TSComputeIJacobianP(PetscTS, PetscReal, PetscVec, PetscVec, PetscReal, PetscMat, PetscBool)
184    PetscErrorCode TSComputeI2Function(PetscTS, PetscReal, PetscVec, PetscVec, PetscVec, PetscVec)
185    PetscErrorCode TSComputeI2Jacobian(PetscTS, PetscReal, PetscVec, PetscVec, PetscVec, PetscReal, PetscReal, PetscMat, PetscMat)
186
187    PetscErrorCode TSSetTime(PetscTS, PetscReal)
188    PetscErrorCode TSGetTime(PetscTS, PetscReal*)
189    PetscErrorCode TSGetPrevTime(PetscTS, PetscReal*)
190    PetscErrorCode TSGetSolveTime(PetscTS, PetscReal*)
191    PetscErrorCode TSSetTimeStep(PetscTS, PetscReal)
192    PetscErrorCode TSGetTimeStep(PetscTS, PetscReal*)
193    PetscErrorCode TSSetStepNumber(PetscTS, PetscInt)
194    PetscErrorCode TSGetStepNumber(PetscTS, PetscInt*)
195    PetscErrorCode TSSetMaxSteps(PetscTS, PetscInt)
196    PetscErrorCode TSGetMaxSteps(PetscTS, PetscInt*)
197    PetscErrorCode TSSetMaxTime(PetscTS, PetscReal)
198    PetscErrorCode TSGetMaxTime(PetscTS, PetscReal*)
199    PetscErrorCode TSSetExactFinalTime(PetscTS, PetscTSExactFinalTimeOption)
200    PetscErrorCode TSSetTimeSpan(PetscTS, PetscInt, PetscReal*)
201    PetscErrorCode TSSetEvaluationTimes(PetscTS, PetscInt, PetscReal*)
202    PetscErrorCode TSGetEvaluationTimes(PetscTS, PetscInt*, const PetscReal**)
203    PetscErrorCode TSGetEvaluationSolutions(PetscTS, PetscInt*, const PetscReal**, PetscVec**)
204
205    PetscErrorCode TSSetConvergedReason(PetscTS, PetscTSConvergedReason)
206    PetscErrorCode TSGetConvergedReason(PetscTS, PetscTSConvergedReason*)
207    PetscErrorCode TSGetSNESIterations(PetscTS, PetscInt*)
208    PetscErrorCode TSGetKSPIterations(PetscTS, PetscInt*)
209    PetscErrorCode TSGetStepRejections(PetscTS, PetscInt*)
210    PetscErrorCode TSSetMaxStepRejections(PetscTS, PetscInt)
211    PetscErrorCode TSGetSNESFailures(PetscTS, PetscInt*)
212    PetscErrorCode TSSetMaxSNESFailures(PetscTS, PetscInt)
213    PetscErrorCode TSSetErrorIfStepFails(PetscTS, PetscBool)
214    PetscErrorCode TSSetTolerances(PetscTS, PetscReal, PetscVec, PetscReal, PetscVec)
215    PetscErrorCode TSGetTolerances(PetscTS, PetscReal*, PetscVec*, PetscReal*, PetscVec*)
216
217    PetscErrorCode TSMonitorSet(PetscTS, PetscTSMonitorFunction, void*, PetscTSCtxDel*)
218    PetscErrorCode TSMonitorCancel(PetscTS)
219    PetscErrorCode TSMonitor(PetscTS, PetscInt, PetscReal, PetscVec)
220
221    ctypedef PetscErrorCode (*PetscTSIndicator)(PetscTS, PetscReal, PetscVec, PetscReal[], void*) except PETSC_ERR_PYTHON
222    ctypedef PetscErrorCode (*PetscTSPostEvent)(PetscTS, PetscInt, PetscInt[], PetscReal, PetscVec, PetscBool, void*) except PETSC_ERR_PYTHON
223
224    PetscErrorCode TSSetEventHandler(PetscTS, PetscInt, PetscInt[], PetscBool[], PetscTSIndicator, PetscTSPostEvent, void*)
225    PetscErrorCode TSSetEventTolerances(PetscTS, PetscReal, PetscReal[])
226    PetscErrorCode TSGetNumEvents(PetscTS, PetscInt*)
227
228    ctypedef PetscErrorCode (*PetscTSAdjointR)(PetscTS, PetscReal, PetscVec, PetscVec, void*) except PETSC_ERR_PYTHON
229    ctypedef PetscErrorCode (*PetscTSAdjointDRDY)(PetscTS, PetscReal, PetscVec, PetscVec[], void*) except PETSC_ERR_PYTHON
230    ctypedef PetscErrorCode (*PetscTSAdjointDRDP)(PetscTS, PetscReal, PetscVec, PetscVec[], void*) except PETSC_ERR_PYTHON
231    ctypedef PetscErrorCode (*PetscTSRHSJacobianP)(PetscTS, PetscReal, PetscVec, PetscMat, void*) except PETSC_ERR_PYTHON
232
233    PetscErrorCode TSSetSaveTrajectory(PetscTS)
234    PetscErrorCode TSRemoveTrajectory(PetscTS)
235    PetscErrorCode TSSetCostGradients(PetscTS, PetscInt, PetscVec*, PetscVec*)
236    PetscErrorCode TSGetCostGradients(PetscTS, PetscInt*, PetscVec**, PetscVec**)
237    PetscErrorCode TSCreateQuadratureTS(PetscTS, PetscBool, PetscTS*)
238    PetscErrorCode TSGetQuadratureTS(PetscTS, PetscBool*, PetscTS*)
239    PetscErrorCode TSGetCostIntegral(PetscTS, PetscVec*)
240
241    PetscErrorCode TSSetRHSJacobianP(PetscTS, PetscMat, PetscTSRHSJacobianP, void*)
242    PetscErrorCode TSComputeRHSJacobianP(PetscTS, PetscReal, PetscVec, PetscMat)
243
244    PetscErrorCode TSAdjointSolve(PetscTS)
245    PetscErrorCode TSAdjointSetSteps(PetscTS, PetscInt)
246    PetscErrorCode TSAdjointStep(PetscTS)
247    PetscErrorCode TSAdjointSetUp(PetscTS)
248    PetscErrorCode TSAdjointReset(PetscTS)
249    PetscErrorCode TSAdjointComputeDRDPFunction(PetscTS, PetscReal, PetscVec, PetscVec*)
250    PetscErrorCode TSAdjointComputeDRDYFunction(PetscTS, PetscReal, PetscVec, PetscVec*)
251    PetscErrorCode TSAdjointCostIntegral(PetscTS)
252
253    PetscErrorCode TSForwardSetSensitivities(PetscTS, PetscInt, PetscVec*, PetscInt, PetscVec*)
254    PetscErrorCode TSForwardGetSensitivities(PetscTS, PetscInt*, PetscVec**, PetscInt*, PetscVec**)
255    PetscErrorCode TSForwardSetIntegralGradients(PetscTS, PetscInt, PetscVec *, PetscVec *)
256    PetscErrorCode TSForwardGetIntegralGradients(PetscTS, PetscInt*, PetscVec **, PetscVec **)
257    PetscErrorCode TSForwardSetRHSJacobianP(PetscTS, PetscVec*, PetscTSCostIntegrandFunction, void*)
258    PetscErrorCode TSForwardComputeRHSJacobianP(PetscTS, PetscReal, PetscVec, PetscVec*)
259    PetscErrorCode TSForwardSetUp(PetscTS)
260    PetscErrorCode TSForwardCostIntegral(PetscTS)
261    PetscErrorCode TSForwardStep(PetscTS)
262
263    PetscErrorCode TSSetPreStep(PetscTS, PetscTSPreStepFunction)
264    PetscErrorCode TSSetPostStep(PetscTS, PetscTSPostStepFunction)
265
266    PetscErrorCode TSSetUp(PetscTS)
267    PetscErrorCode TSReset(PetscTS)
268    PetscErrorCode TSStep(PetscTS)
269    PetscErrorCode TSRestartStep(PetscTS)
270    PetscErrorCode TSRollBack(PetscTS)
271    PetscErrorCode TSSolve(PetscTS, PetscVec)
272    PetscErrorCode TSInterpolate(PetscTS, PetscReal, PetscVec)
273    PetscErrorCode TSPreStage(PetscTS, PetscReal)
274    PetscErrorCode TSPostStage(PetscTS, PetscReal, PetscInt, PetscVec*)
275
276    PetscErrorCode TSThetaSetTheta(PetscTS, PetscReal)
277    PetscErrorCode TSThetaGetTheta(PetscTS, PetscReal*)
278    PetscErrorCode TSThetaSetEndpoint(PetscTS, PetscBool)
279    PetscErrorCode TSThetaGetEndpoint(PetscTS, PetscBool*)
280
281    PetscErrorCode TSAlphaSetRadius(PetscTS, PetscReal)
282    PetscErrorCode TSAlphaSetParams(PetscTS, PetscReal, PetscReal, PetscReal)
283    PetscErrorCode TSAlphaGetParams(PetscTS, PetscReal*, PetscReal*, PetscReal*)
284
285    ctypedef const char* PetscTSRKType "TSRKType"
286    PetscTSRKType TSRK1FE
287    PetscTSRKType TSRK2A
288    PetscTSRKType TSRK2B
289    PetscTSRKType TSRK3
290    PetscTSRKType TSRK3BS
291    PetscTSRKType TSRK4
292    PetscTSRKType TSRK5F
293    PetscTSRKType TSRK5DP
294    PetscTSRKType TSRK5BS
295    PetscTSRKType TSRK6VR
296    PetscTSRKType TSRK7VR
297    PetscTSRKType TSRK8VR
298
299    PetscErrorCode TSRKGetType(PetscTS, PetscTSRKType*)
300    PetscErrorCode TSRKSetType(PetscTS, PetscTSRKType)
301
302    ctypedef const char* PetscTSARKIMEXType "TSARKIMEXType"
303    PetscTSARKIMEXType TSARKIMEX1BEE
304    PetscTSARKIMEXType TSARKIMEXA2
305    PetscTSARKIMEXType TSARKIMEXL2
306    PetscTSARKIMEXType TSARKIMEXARS122
307    PetscTSARKIMEXType TSARKIMEX2C
308    PetscTSARKIMEXType TSARKIMEX2D
309    PetscTSARKIMEXType TSARKIMEX2E
310    PetscTSARKIMEXType TSARKIMEXPRSSP2
311    PetscTSARKIMEXType TSARKIMEX3
312    PetscTSARKIMEXType TSARKIMEXBPR3
313    PetscTSARKIMEXType TSARKIMEXARS443
314    PetscTSARKIMEXType TSARKIMEX4
315    PetscTSARKIMEXType TSARKIMEX5
316
317    PetscErrorCode TSARKIMEXGetType(PetscTS, PetscTSRKType*)
318    PetscErrorCode TSARKIMEXSetType(PetscTS, PetscTSRKType)
319    PetscErrorCode TSARKIMEXSetFullyImplicit(PetscTS, PetscBool)
320    PetscErrorCode TSARKIMEXSetFastSlowSplit(PetscTS, PetscBool)
321
322    ctypedef const char* PetscTSDIRKType "TSDIRKType"
323    PetscTSDIRKType TSDIRKS212
324    PetscTSDIRKType TSDIRKES122SAL
325    PetscTSDIRKType TSDIRKES213SAL
326    PetscTSDIRKType TSDIRKES324SAL
327    PetscTSDIRKType TSDIRKES325SAL
328    PetscTSDIRKType TSDIRK657A
329    PetscTSDIRKType TSDIRKES648SA
330    PetscTSDIRKType TSDIRK658A
331    PetscTSDIRKType TSDIRKS659A
332    PetscTSDIRKType TSDIRK7510SAL
333    PetscTSDIRKType TSDIRKES7510SA
334    PetscTSDIRKType TSDIRK759A
335    PetscTSDIRKType TSDIRKS7511SAL
336    PetscTSDIRKType TSDIRK8614A
337    PetscTSDIRKType TSDIRK8616SAL
338    PetscTSDIRKType TSDIRKES8516SAL
339
340    PetscErrorCode TSDIRKGetType(PetscTS, PetscTSDIRKType*)
341    PetscErrorCode TSDIRKSetType(PetscTS, PetscTSDIRKType)
342
343    PetscErrorCode TSPythonSetType(PetscTS, char[])
344    PetscErrorCode TSPythonGetType(PetscTS, char*[])
345
346cdef extern from * nogil:
347    struct _p_TSAdapt
348    ctypedef _p_TSAdapt *PetscTSAdapt "TSAdapt"
349    PetscErrorCode TSGetAdapt(PetscTS, PetscTSAdapt*)
350    PetscErrorCode TSAdaptGetStepLimits(PetscTSAdapt, PetscReal*, PetscReal*)
351    PetscErrorCode TSAdaptSetStepLimits(PetscTSAdapt, PetscReal, PetscReal)
352    PetscErrorCode TSAdaptCheckStage(PetscTSAdapt, PetscTS, PetscReal, PetscVec, PetscBool*)
353
354cdef extern from * nogil: # custom.h
355    PetscErrorCode TSSetTimeStepNumber(PetscTS, PetscInt)
356
357# -----------------------------------------------------------------------------
358
359cdef inline TS ref_TS(PetscTS ts):
360    cdef TS ob = <TS> TS()
361    ob.ts = ts
362    CHKERR(PetscINCREF(ob.obj))
363    return ob
364
365# -----------------------------------------------------------------------------
366
367cdef PetscErrorCode TS_RHSFunction(
368    PetscTS   ts,
369    PetscReal t,
370    PetscVec  x,
371    PetscVec  f,
372    void      *ctx,
373   ) except PETSC_ERR_PYTHON with gil:
374    cdef TS  Ts   = ref_TS(ts)
375    cdef Vec Xvec = ref_Vec(x)
376    cdef Vec Fvec = ref_Vec(f)
377    cdef object context = Ts.get_attr('__rhsfunction__')
378    if context is None and ctx != NULL: context = <object>ctx
379    assert context is not None and type(context) is tuple # sanity check
380    (function, args, kargs) = context
381    function(Ts, toReal(t), Xvec, Fvec, *args, **kargs)
382    return PETSC_SUCCESS
383
384cdef PetscErrorCode TS_RHSJacobian(
385    PetscTS   ts,
386    PetscReal t,
387    PetscVec  x,
388    PetscMat  J,
389    PetscMat  P,
390    void      *ctx,
391   ) except PETSC_ERR_PYTHON with gil:
392    cdef TS  Ts   = ref_TS(ts)
393    cdef Vec Xvec = ref_Vec(x)
394    cdef Mat Jmat = ref_Mat(J)
395    cdef Mat Pmat = ref_Mat(P)
396    cdef object context = Ts.get_attr('__rhsjacobian__')
397    if context is None and ctx != NULL: context = <object>ctx
398    assert context is not None and type(context) is tuple # sanity check
399    (jacobian, args, kargs) = context
400    jacobian(Ts, toReal(t), Xvec, Jmat, Pmat, *args, **kargs)
401    return PETSC_SUCCESS
402
403cdef PetscErrorCode TS_RHSJacobianP(
404    PetscTS   ts,
405    PetscReal t,
406    PetscVec  x,
407    PetscMat  J,
408    void      *ctx,
409   ) except PETSC_ERR_PYTHON with gil:
410    cdef TS  Ts   = ref_TS(ts)
411    cdef Vec Xvec = ref_Vec(x)
412    cdef Mat Jmat = ref_Mat(J)
413    cdef object context = Ts.get_attr('__rhsjacobianp__')
414    if context is None and ctx != NULL: context = <object>ctx
415    assert context is not None and type(context) is tuple # sanity check
416    (jacobianp, args, kargs) = context
417    jacobianp(Ts, toReal(t), Xvec, Jmat, *args, **kargs)
418    return PETSC_SUCCESS
419
420# -----------------------------------------------------------------------------
421
422cdef PetscErrorCode TS_IFunction(
423    PetscTS   ts,
424    PetscReal t,
425    PetscVec  x,
426    PetscVec  xdot,
427    PetscVec  f,
428    void      *ctx,
429   ) except PETSC_ERR_PYTHON with gil:
430    cdef TS  Ts    = ref_TS(ts)
431    cdef Vec Xvec  = ref_Vec(x)
432    cdef Vec XDvec = ref_Vec(xdot)
433    cdef Vec Fvec  = ref_Vec(f)
434    cdef object context = Ts.get_attr('__ifunction__')
435    if context is None and ctx != NULL: context = <object>ctx
436    assert context is not None and type(context) is tuple # sanity check
437    (function, args, kargs) = context
438    function(Ts, toReal(t), Xvec, XDvec, Fvec, *args, **kargs)
439    return PETSC_SUCCESS
440
441cdef PetscErrorCode TS_IJacobian(
442    PetscTS   ts,
443    PetscReal t,
444    PetscVec  x,
445    PetscVec  xdot,
446    PetscReal a,
447    PetscMat  J,
448    PetscMat  P,
449    void      *ctx,
450   ) except PETSC_ERR_PYTHON with gil:
451    cdef TS   Ts    = ref_TS(ts)
452    cdef Vec  Xvec  = ref_Vec(x)
453    cdef Vec  XDvec = ref_Vec(xdot)
454    cdef Mat  Jmat  = ref_Mat(J)
455    cdef Mat  Pmat  = ref_Mat(P)
456    cdef object context = Ts.get_attr('__ijacobian__')
457    if context is None and ctx != NULL: context = <object>ctx
458    assert context is not None and type(context) is tuple # sanity check
459    (jacobian, args, kargs) = context
460    jacobian(Ts, toReal(t), Xvec, XDvec, toReal(a), Jmat, Pmat, *args, **kargs)
461    return PETSC_SUCCESS
462
463cdef PetscErrorCode TS_IJacobianP(
464    PetscTS   ts,
465    PetscReal t,
466    PetscVec  x,
467    PetscVec  xdot,
468    PetscReal a,
469    PetscMat  J,
470    void      *ctx,
471   ) except PETSC_ERR_PYTHON with gil:
472    cdef TS   Ts    = ref_TS(ts)
473    cdef Vec  Xvec  = ref_Vec(x)
474    cdef Vec  XDvec = ref_Vec(xdot)
475    cdef Mat  Jmat  = ref_Mat(J)
476    cdef object context = Ts.get_attr('__ijacobianp__')
477    if context is None and ctx != NULL: context = <object>ctx
478    assert context is not None and type(context) is tuple # sanity check
479    (jacobian, args, kargs) = context
480    jacobian(Ts, toReal(t), Xvec, XDvec, toReal(a), Jmat, *args, **kargs)
481    return PETSC_SUCCESS
482
483cdef PetscErrorCode TS_I2Function(
484    PetscTS   ts,
485    PetscReal t,
486    PetscVec  x,
487    PetscVec  xdot,
488    PetscVec  xdotdot,
489    PetscVec  f,
490    void      *ctx,
491   ) except PETSC_ERR_PYTHON with gil:
492    cdef TS  Ts    = ref_TS(ts)
493    cdef Vec Xvec  = ref_Vec(x)
494    cdef Vec XDvec = ref_Vec(xdot)
495    cdef Vec XDDvec = ref_Vec(xdotdot)
496    cdef Vec Fvec  = ref_Vec(f)
497    cdef object context = Ts.get_attr('__i2function__')
498    if context is None and ctx != NULL: context = <object>ctx
499    assert context is not None and type(context) is tuple # sanity check
500    (function, args, kargs) = context
501    function(Ts, toReal(t), Xvec, XDvec, XDDvec, Fvec, *args, **kargs)
502    return PETSC_SUCCESS
503
504cdef PetscErrorCode TS_I2Jacobian(
505    PetscTS   ts,
506    PetscReal t,
507    PetscVec  x,
508    PetscVec  xdot,
509    PetscVec  xdotdot,
510    PetscReal v,
511    PetscReal a,
512    PetscMat  J,
513    PetscMat  P,
514    void      *ctx,
515   ) except PETSC_ERR_PYTHON with gil:
516    cdef TS   Ts    = ref_TS(ts)
517    cdef Vec  Xvec  = ref_Vec(x)
518    cdef Vec  XDvec = ref_Vec(xdot)
519    cdef Vec  XDDvec = ref_Vec(xdotdot)
520    cdef Mat  Jmat  = ref_Mat(J)
521    cdef Mat  Pmat  = ref_Mat(P)
522    cdef object context = Ts.get_attr('__i2jacobian__')
523    if context is None and ctx != NULL: context = <object>ctx
524    assert context is not None and type(context) is tuple # sanity check
525    (jacobian, args, kargs) = context
526    jacobian(Ts, toReal(t), Xvec, XDvec, XDDvec, toReal(v), toReal(a), Jmat, Pmat, *args, **kargs)
527    return PETSC_SUCCESS
528
529# -----------------------------------------------------------------------------
530
531cdef PetscErrorCode TS_Monitor(
532    PetscTS   ts,
533    PetscInt  step,
534    PetscReal time,
535    PetscVec  u,
536    void      *ctx,
537   ) except PETSC_ERR_PYTHON with gil:
538    cdef TS  Ts = ref_TS(ts)
539    cdef Vec Vu = ref_Vec(u)
540    cdef object monitorlist = Ts.get_attr('__monitor__')
541    if monitorlist is None: return PETSC_SUCCESS
542    for (monitor, args, kargs) in monitorlist:
543        monitor(Ts, toInt(step), toReal(time), Vu, *args, **kargs)
544    return PETSC_SUCCESS
545
546# -----------------------------------------------------------------------------
547
548cdef PetscErrorCode TS_Indicator(
549    PetscTS   ts,
550    PetscReal time,
551    PetscVec  u,
552    PetscReal fvalue[],
553    void      *ctx,
554   ) except PETSC_ERR_PYTHON with gil:
555    cdef TS  Ts = ref_TS(ts)
556    cdef Vec Vu = ref_Vec(u)
557    cdef object context = Ts.get_attr('__indicator__')
558    if context is None: return PETSC_SUCCESS
559    (indicator, args, kargs) = context
560    cdef PetscInt nevents = 0
561    CHKERR(TSGetNumEvents(ts, &nevents))
562    cdef npy_intp s = <npy_intp> nevents
563    fvalue_array = PyArray_SimpleNewFromData(1, &s, NPY_PETSC_REAL, fvalue)
564    indicator(Ts, toReal(time), Vu, fvalue_array, *args, **kargs)
565    return PETSC_SUCCESS
566
567cdef PetscErrorCode TS_PostEvent(
568    PetscTS   ts,
569    PetscInt  nevents_zero,
570    PetscInt  events_zero[],
571    PetscReal time,
572    PetscVec  u,
573    PetscBool forward,
574    void      *ctx,
575   ) except PETSC_ERR_PYTHON with gil:
576    cdef TS  Ts = ref_TS(ts)
577    cdef Vec Vu = ref_Vec(u)
578    cdef object context = Ts.get_attr('__postevent__')
579    if context is None: return PETSC_SUCCESS
580    (postevent, args, kargs) = context
581    cdef npy_intp s = <npy_intp> nevents_zero
582    events_zero_array = PyArray_SimpleNewFromData(1, &s, NPY_PETSC_INT, events_zero)
583    postevent(Ts, events_zero_array, toReal(time), Vu, toBool(forward), *args, **kargs)
584    return PETSC_SUCCESS
585
586cdef PetscErrorCode TS_PreStep(
587    PetscTS ts,
588   ) except PETSC_ERR_PYTHON with gil:
589    cdef TS Ts = ref_TS(ts)
590    (prestep, args, kargs) = Ts.get_attr('__prestep__')
591    prestep(Ts, *args, **kargs)
592    return PETSC_SUCCESS
593
594cdef PetscErrorCode TS_PostStep(
595    PetscTS ts,
596   ) except PETSC_ERR_PYTHON with gil:
597    cdef TS Ts = ref_TS(ts)
598    (poststep, args, kargs) = Ts.get_attr('__poststep__')
599    poststep(Ts, *args, **kargs)
600    return PETSC_SUCCESS
601
602# -----------------------------------------------------------------------------
603