_ode.py 49 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409
  1. # Authors: Pearu Peterson, Pauli Virtanen, John Travers
  2. """
  3. First-order ODE integrators.
  4. User-friendly interface to various numerical integrators for solving a
  5. system of first order ODEs with prescribed initial conditions::
  6. d y(t)[i]
  7. --------- = f(t,y(t))[i],
  8. d t
  9. y(t=0)[i] = y0[i],
  10. where::
  11. i = 0, ..., len(y0) - 1
  12. class ode
  13. ---------
  14. A generic interface class to numeric integrators. It has the following
  15. methods::
  16. integrator = ode(f, jac=None)
  17. integrator = integrator.set_integrator(name, **params)
  18. integrator = integrator.set_initial_value(y0, t0=0.0)
  19. integrator = integrator.set_f_params(*args)
  20. integrator = integrator.set_jac_params(*args)
  21. y1 = integrator.integrate(t1, step=False, relax=False)
  22. flag = integrator.successful()
  23. class complex_ode
  24. -----------------
  25. This class has the same generic interface as ode, except it can handle complex
  26. f, y and Jacobians by transparently translating them into the equivalent
  27. real-valued system. It supports the real-valued solvers (i.e., not zvode) and is
  28. an alternative to ode with the zvode solver, sometimes performing better.
  29. """
  30. # XXX: Integrators must have:
  31. # ===========================
  32. # cvode - C version of vode and vodpk with many improvements.
  33. # Get it from http://www.netlib.org/ode/cvode.tar.gz.
  34. # To wrap cvode to Python, one must write the extension module by
  35. # hand. Its interface is too much 'advanced C' that using f2py
  36. # would be too complicated (or impossible).
  37. #
  38. # How to define a new integrator:
  39. # ===============================
  40. #
  41. # class myodeint(IntegratorBase):
  42. #
  43. # runner = <odeint function> or None
  44. #
  45. # def __init__(self,...): # required
  46. # <initialize>
  47. #
  48. # def reset(self,n,has_jac): # optional
  49. # # n - the size of the problem (number of equations)
  50. # # has_jac - whether user has supplied its own routine for Jacobian
  51. # <allocate memory,initialize further>
  52. #
  53. # def run(self,f,jac,y0,t0,t1,f_params,jac_params): # required
  54. # # this method is called to integrate from t=t0 to t=t1
  55. # # with initial condition y0. f and jac are user-supplied functions
  56. # # that define the problem. f_params,jac_params are additional
  57. # # arguments
  58. # # to these functions.
  59. # <calculate y1>
  60. # if <calculation was unsuccessful>:
  61. # self.success = 0
  62. # return t1,y1
  63. #
  64. # # In addition, one can define step() and run_relax() methods (they
  65. # # take the same arguments as run()) if the integrator can support
  66. # # these features (see IntegratorBase doc strings).
  67. #
  68. # if myodeint.runner:
  69. # IntegratorBase.integrator_classes.append(myodeint)
  70. __all__ = ['ode', 'complex_ode']
  71. import re
  72. import types
  73. import warnings
  74. import numpy as np
  75. from numpy import asarray, array, zeros, isscalar, real, imag
  76. from . import _vode
  77. from . import _dop
  78. from ._odepack import lsoda as lsoda_step
  79. # ------------------------------------------------------------------------------
  80. # User interface
  81. # ------------------------------------------------------------------------------
  82. class ode:
  83. """
  84. A generic interface class to numeric integrators.
  85. Solve an equation system :math:`y'(t) = f(t,y)` with (optional) ``jac = df/dy``.
  86. *Note*: The first two arguments of ``f(t, y, ...)`` are in the
  87. opposite order of the arguments in the system definition function used
  88. by `scipy.integrate.odeint`.
  89. Parameters
  90. ----------
  91. f : callable ``f(t, y, *f_args)``
  92. Right-hand side of the differential equation. t is a scalar,
  93. ``y.shape == (n,)``.
  94. ``f_args`` is set by calling ``set_f_params(*args)``.
  95. `f` should return a scalar, array or list (not a tuple).
  96. jac : callable ``jac(t, y, *jac_args)``, optional
  97. Jacobian of the right-hand side, ``jac[i,j] = d f[i] / d y[j]``.
  98. ``jac_args`` is set by calling ``set_jac_params(*args)``.
  99. Attributes
  100. ----------
  101. t : float
  102. Current time.
  103. y : ndarray
  104. Current variable values.
  105. See also
  106. --------
  107. odeint : an integrator with a simpler interface based on lsoda from ODEPACK
  108. quad : for finding the area under a curve
  109. Notes
  110. -----
  111. Available integrators are listed below. They can be selected using
  112. the `set_integrator` method.
  113. "vode"
  114. Real-valued Variable-coefficient Ordinary Differential Equation
  115. solver, with fixed-leading-coefficient implementation. It provides
  116. implicit Adams method (for non-stiff problems) and a method based on
  117. backward differentiation formulas (BDF) (for stiff problems).
  118. Source: http://www.netlib.org/ode/vode.f
  119. This integrator accepts the following parameters in `set_integrator`
  120. method of the `ode` class:
  121. - atol : float or sequence
  122. absolute tolerance for solution
  123. - rtol : float or sequence
  124. relative tolerance for solution
  125. - lband : None or int
  126. - uband : None or int
  127. Jacobian band width, jac[i,j] != 0 for i-lband <= j <= i+uband.
  128. Setting these requires your jac routine to return the jacobian
  129. in packed format, jac_packed[i-j+uband, j] = jac[i,j]. The
  130. dimension of the matrix must be (lband+uband+1, len(y)).
  131. - method: 'adams' or 'bdf'
  132. Which solver to use, Adams (non-stiff) or BDF (stiff)
  133. - with_jacobian : bool
  134. This option is only considered when the user has not supplied a
  135. Jacobian function and has not indicated (by setting either band)
  136. that the Jacobian is banded. In this case, `with_jacobian` specifies
  137. whether the iteration method of the ODE solver's correction step is
  138. chord iteration with an internally generated full Jacobian or
  139. functional iteration with no Jacobian.
  140. - nsteps : int
  141. Maximum number of (internally defined) steps allowed during one
  142. call to the solver.
  143. - first_step : float
  144. - min_step : float
  145. - max_step : float
  146. Limits for the step sizes used by the integrator.
  147. - order : int
  148. Maximum order used by the integrator,
  149. order <= 12 for Adams, <= 5 for BDF.
  150. "zvode"
  151. Complex-valued Variable-coefficient Ordinary Differential Equation
  152. solver, with fixed-leading-coefficient implementation. It provides
  153. implicit Adams method (for non-stiff problems) and a method based on
  154. backward differentiation formulas (BDF) (for stiff problems).
  155. Source: http://www.netlib.org/ode/zvode.f
  156. This integrator accepts the same parameters in `set_integrator`
  157. as the "vode" solver.
  158. .. note::
  159. When using ZVODE for a stiff system, it should only be used for
  160. the case in which the function f is analytic, that is, when each f(i)
  161. is an analytic function of each y(j). Analyticity means that the
  162. partial derivative df(i)/dy(j) is a unique complex number, and this
  163. fact is critical in the way ZVODE solves the dense or banded linear
  164. systems that arise in the stiff case. For a complex stiff ODE system
  165. in which f is not analytic, ZVODE is likely to have convergence
  166. failures, and for this problem one should instead use DVODE on the
  167. equivalent real system (in the real and imaginary parts of y).
  168. "lsoda"
  169. Real-valued Variable-coefficient Ordinary Differential Equation
  170. solver, with fixed-leading-coefficient implementation. It provides
  171. automatic method switching between implicit Adams method (for non-stiff
  172. problems) and a method based on backward differentiation formulas (BDF)
  173. (for stiff problems).
  174. This integrator uses the C translation of the original Fortran 77 ODEPACK
  175. library, which can be found at http://www.netlib.org/odepack
  176. This integrator accepts the following parameters in `set_integrator`
  177. method of the `ode` class:
  178. - atol : float or sequence
  179. absolute tolerance for solution
  180. - rtol : float or sequence
  181. relative tolerance for solution
  182. - lband : None or int
  183. - uband : None or int
  184. Jacobian band width, jac[i,j] != 0 for i-lband <= j <= i+uband.
  185. Setting these requires your jac routine to return the jacobian
  186. in packed format, jac_packed[i-j+uband, j] = jac[i,j].
  187. - with_jacobian : bool
  188. *Not used.*
  189. - nsteps : int
  190. Maximum number of (internally defined) steps allowed during one
  191. call to the solver.
  192. - first_step : float
  193. - min_step : float
  194. - max_step : float
  195. Limits for the step sizes used by the integrator.
  196. - max_order_ns : int
  197. Maximum order used in the nonstiff case (default 12).
  198. - max_order_s : int
  199. Maximum order used in the stiff case (default 5).
  200. - max_hnil : int
  201. Maximum number of messages reporting too small step size (t + h = t)
  202. (default 0)
  203. - ixpr : int
  204. Whether to generate extra printing at method switches (default False).
  205. "dopri5"
  206. This is an explicit runge-kutta method of order (4)5 due to Dormand &
  207. Prince (with stepsize control and dense output).
  208. Authors:
  209. E. Hairer and G. Wanner
  210. Universite de Geneve, Dept. de Mathematiques
  211. CH-1211 Geneve 24, Switzerland
  212. e-mail: ernst.hairer@math.unige.ch, gerhard.wanner@math.unige.ch
  213. This code is described in [HNW93]_.
  214. This integrator accepts the following parameters in set_integrator()
  215. method of the ode class:
  216. - atol : float or sequence
  217. absolute tolerance for solution
  218. - rtol : float or sequence
  219. relative tolerance for solution
  220. - nsteps : int
  221. Maximum number of (internally defined) steps allowed during one
  222. call to the solver.
  223. - first_step : float
  224. - max_step : float
  225. - safety : float
  226. Safety factor on new step selection (default 0.9)
  227. - ifactor : float
  228. - dfactor : float
  229. Maximum factor to increase/decrease step size by in one step
  230. - beta : float
  231. Beta parameter for stabilised step size control.
  232. - verbosity : int
  233. Switch for printing messages (< 0 for no messages).
  234. "dop853"
  235. This is an explicit runge-kutta method of order 8(5,3) due to Dormand
  236. & Prince (with stepsize control and dense output).
  237. Options and references the same as "dopri5".
  238. Examples
  239. --------
  240. A problem to integrate and the corresponding jacobian:
  241. >>> from scipy.integrate import ode
  242. >>>
  243. >>> y0, t0 = [1.0j, 2.0], 0
  244. >>>
  245. >>> def f(t, y, arg1):
  246. ... return [1j*arg1*y[0] + y[1], -arg1*y[1]**2]
  247. >>> def jac(t, y, arg1):
  248. ... return [[1j*arg1, 1], [0, -arg1*2*y[1]]]
  249. The integration:
  250. >>> r = ode(f, jac).set_integrator('zvode', method='bdf')
  251. >>> r.set_initial_value(y0, t0).set_f_params(2.0).set_jac_params(2.0)
  252. >>> t1 = 10
  253. >>> dt = 1
  254. >>> while r.successful() and r.t < t1:
  255. ... print(r.t+dt, r.integrate(r.t+dt))
  256. 1 [-0.71038232+0.23749653j 0.40000271+0.j ]
  257. 2.0 [0.19098503-0.52359246j 0.22222356+0.j ]
  258. 3.0 [0.47153208+0.52701229j 0.15384681+0.j ]
  259. 4.0 [-0.61905937+0.30726255j 0.11764744+0.j ]
  260. 5.0 [0.02340997-0.61418799j 0.09523835+0.j ]
  261. 6.0 [0.58643071+0.339819j 0.08000018+0.j ]
  262. 7.0 [-0.52070105+0.44525141j 0.06896565+0.j ]
  263. 8.0 [-0.15986733-0.61234476j 0.06060616+0.j ]
  264. 9.0 [0.64850462+0.15048982j 0.05405414+0.j ]
  265. 10.0 [-0.38404699+0.56382299j 0.04878055+0.j ]
  266. References
  267. ----------
  268. .. [HNW93] E. Hairer, S.P. Norsett and G. Wanner, Solving Ordinary
  269. Differential Equations i. Nonstiff Problems. 2nd edition.
  270. Springer Series in Computational Mathematics,
  271. Springer-Verlag (1993)
  272. """
  273. # generic type compatibility with scipy-stubs
  274. __class_getitem__ = classmethod(types.GenericAlias)
  275. def __init__(self, f, jac=None):
  276. self.stiff = 0
  277. self.f = f
  278. self.jac = jac
  279. self.f_params = ()
  280. self.jac_params = ()
  281. self._y = []
  282. @property
  283. def y(self):
  284. return self._y
  285. def set_initial_value(self, y, t=0.0):
  286. """Set initial conditions y(t) = y."""
  287. if isscalar(y):
  288. y = [y]
  289. n_prev = len(self._y)
  290. if not n_prev:
  291. self.set_integrator('') # find first available integrator
  292. # NOTE: The C code modifies y in place, hence the copy.
  293. self._y = asarray(y, self._integrator.scalar).copy()
  294. self.t = t
  295. self._integrator.reset(len(self._y), self.jac is not None)
  296. return self
  297. def set_integrator(self, name, **integrator_params):
  298. """
  299. Set integrator by name.
  300. Parameters
  301. ----------
  302. name : str
  303. Name of the integrator.
  304. **integrator_params
  305. Additional parameters for the integrator.
  306. """
  307. integrator = find_integrator(name)
  308. if integrator is None:
  309. # FIXME: this really should be raise an exception. Will that break
  310. # any code?
  311. message = f'No integrator name match with {name!r} or is not available.'
  312. warnings.warn(message, stacklevel=2)
  313. else:
  314. self._integrator = integrator(**integrator_params)
  315. if not len(self._y):
  316. self.t = 0.0
  317. self._y = array([0.0], self._integrator.scalar)
  318. self._integrator.reset(len(self._y), self.jac is not None)
  319. return self
  320. def integrate(self, t, step=False, relax=False):
  321. """Find y=y(t), set y as an initial condition, and return y.
  322. Parameters
  323. ----------
  324. t : float
  325. The endpoint of the integration step.
  326. step : bool
  327. If True, and if the integrator supports the step method,
  328. then perform a single integration step and return.
  329. This parameter is provided in order to expose internals of
  330. the implementation, and should not be changed from its default
  331. value in most cases.
  332. relax : bool
  333. If True and if the integrator supports the run_relax method,
  334. then integrate until t_1 >= t and return. ``relax`` is not
  335. referenced if ``step=True``.
  336. This parameter is provided in order to expose internals of
  337. the implementation, and should not be changed from its default
  338. value in most cases.
  339. Returns
  340. -------
  341. y : float
  342. The integrated value at t
  343. """
  344. if step and self._integrator.supports_step:
  345. mth = self._integrator.step
  346. elif relax and self._integrator.supports_run_relax:
  347. mth = self._integrator.run_relax
  348. else:
  349. mth = self._integrator.run
  350. try:
  351. self._y, self.t = mth(self.f, self.jac or (lambda: None),
  352. self._y, self.t, t,
  353. self.f_params, self.jac_params)
  354. except SystemError as e:
  355. # f2py issue with tuple returns, see ticket 1187.
  356. raise ValueError(
  357. 'Function to integrate must not return a tuple.'
  358. ) from e
  359. return self._y
  360. def successful(self):
  361. """Check if integration was successful."""
  362. try:
  363. self._integrator
  364. except AttributeError:
  365. self.set_integrator('')
  366. return self._integrator.success == 1
  367. def get_return_code(self):
  368. """Extracts the return code for the integration to enable better control
  369. if the integration fails.
  370. In general, a return code > 0 implies success, while a return code < 0
  371. implies failure.
  372. Notes
  373. -----
  374. This section describes possible return codes and their meaning, for available
  375. integrators that can be selected by `set_integrator` method.
  376. "vode"
  377. =========== =======
  378. Return Code Message
  379. =========== =======
  380. 2 Integration successful.
  381. -1 Excess work done on this call. (Perhaps wrong MF.)
  382. -2 Excess accuracy requested. (Tolerances too small.)
  383. -3 Illegal input detected. (See printed message.)
  384. -4 Repeated error test failures. (Check all input.)
  385. -5 Repeated convergence failures. (Perhaps bad Jacobian
  386. supplied or wrong choice of MF or tolerances.)
  387. -6 Error weight became zero during problem. (Solution
  388. component i vanished, and ATOL or ATOL(i) = 0.)
  389. =========== =======
  390. "zvode"
  391. =========== =======
  392. Return Code Message
  393. =========== =======
  394. 2 Integration successful.
  395. -1 Excess work done on this call. (Perhaps wrong MF.)
  396. -2 Excess accuracy requested. (Tolerances too small.)
  397. -3 Illegal input detected. (See printed message.)
  398. -4 Repeated error test failures. (Check all input.)
  399. -5 Repeated convergence failures. (Perhaps bad Jacobian
  400. supplied or wrong choice of MF or tolerances.)
  401. -6 Error weight became zero during problem. (Solution
  402. component i vanished, and ATOL or ATOL(i) = 0.)
  403. =========== =======
  404. "dopri5"
  405. =========== =======
  406. Return Code Message
  407. =========== =======
  408. 1 Integration successful.
  409. 2 Integration successful (interrupted by solout).
  410. -1 Input is not consistent.
  411. -2 Larger nsteps is needed.
  412. -3 Step size becomes too small.
  413. -4 Problem is probably stiff (interrupted).
  414. =========== =======
  415. "dop853"
  416. =========== =======
  417. Return Code Message
  418. =========== =======
  419. 1 Integration successful.
  420. 2 Integration successful (interrupted by solout).
  421. -1 Input is not consistent.
  422. -2 Larger nsteps is needed.
  423. -3 Step size becomes too small.
  424. -4 Problem is probably stiff (interrupted).
  425. =========== =======
  426. "lsoda"
  427. =========== =======
  428. Return Code Message
  429. =========== =======
  430. 2 Integration successful.
  431. -1 Excess work done on this call (perhaps wrong Dfun type).
  432. -2 Excess accuracy requested (tolerances too small).
  433. -3 Illegal input detected (internal error).
  434. -4 Repeated error test failures (internal error).
  435. -5 Repeated convergence failures (perhaps bad Jacobian or tolerances).
  436. -6 Error weight became zero during problem.
  437. -7 Internal workspace insufficient to finish (internal error).
  438. =========== =======
  439. """
  440. try:
  441. self._integrator
  442. except AttributeError:
  443. self.set_integrator('')
  444. return self._integrator.istate
  445. def set_f_params(self, *args):
  446. """Set extra parameters for user-supplied function f."""
  447. self.f_params = args
  448. return self
  449. def set_jac_params(self, *args):
  450. """Set extra parameters for user-supplied function jac."""
  451. self.jac_params = args
  452. return self
  453. def set_solout(self, solout):
  454. """
  455. Set callable to be called at every successful integration step.
  456. Parameters
  457. ----------
  458. solout : callable
  459. ``solout(t, y)`` is called at each internal integrator step,
  460. t is a scalar providing the current independent position
  461. y is the current solution ``y.shape == (n,)``
  462. solout should return -1 to stop integration
  463. otherwise it should return None or 0
  464. """
  465. if self._integrator.supports_solout:
  466. self._integrator.set_solout(solout)
  467. if self._y is not None:
  468. self._integrator.reset(len(self._y), self.jac is not None)
  469. else:
  470. raise ValueError("selected integrator does not support solout,"
  471. " choose another one")
  472. def _transform_banded_jac(bjac):
  473. """
  474. Convert a real matrix of the form (for example)
  475. [0 0 A B] [0 0 0 B]
  476. [0 0 C D] [0 0 A D]
  477. [E F G H] to [0 F C H]
  478. [I J K L] [E J G L]
  479. [I 0 K 0]
  480. That is, every other column is shifted up one.
  481. """
  482. # Shift every other column.
  483. newjac = zeros((bjac.shape[0] + 1, bjac.shape[1]))
  484. newjac[1:, ::2] = bjac[:, ::2]
  485. newjac[:-1, 1::2] = bjac[:, 1::2]
  486. return newjac
  487. class complex_ode(ode):
  488. """
  489. A wrapper of ode for complex systems.
  490. This functions similarly as `ode`, but re-maps a complex-valued
  491. equation system to a real-valued one before using the integrators.
  492. Parameters
  493. ----------
  494. f : callable ``f(t, y, *f_args)``
  495. Rhs of the equation. t is a scalar, ``y.shape == (n,)``.
  496. ``f_args`` is set by calling ``set_f_params(*args)``.
  497. jac : callable ``jac(t, y, *jac_args)``
  498. Jacobian of the rhs, ``jac[i,j] = d f[i] / d y[j]``.
  499. ``jac_args`` is set by calling ``set_f_params(*args)``.
  500. Attributes
  501. ----------
  502. t : float
  503. Current time.
  504. y : ndarray
  505. Current variable values.
  506. Examples
  507. --------
  508. For usage examples, see `ode`.
  509. """
  510. def __init__(self, f, jac=None):
  511. self.cf = f
  512. self.cjac = jac
  513. self.tmp_derivative = None # Work array for derivatives in _wrap
  514. if jac is None:
  515. ode.__init__(self, self._wrap, None)
  516. else:
  517. ode.__init__(self, self._wrap, self._wrap_jac)
  518. def _wrap(self, t, y, *f_args):
  519. f = self.cf(*((t, y[::2] + 1j * y[1::2]) + f_args))
  520. # self.tmp_derivative is a real-valued array containing the interleaved
  521. # real and imaginary parts of f (the derivative).
  522. # IMPORTANT: Must NOT use self.tmp here, as it may alias with y!
  523. self.tmp_derivative[::2] = real(f)
  524. self.tmp_derivative[1::2] = imag(f)
  525. return self.tmp_derivative
  526. def _wrap_jac(self, t, y, *jac_args):
  527. # jac is the complex Jacobian computed by the user-defined function.
  528. jac = self.cjac(*((t, y[::2] + 1j * y[1::2]) + jac_args))
  529. # jac_tmp is the real version of the complex Jacobian. Each complex
  530. # entry in jac, say 2+3j, becomes a 2x2 block of the form
  531. # [2 -3]
  532. # [3 2]
  533. jac_tmp = zeros((2 * jac.shape[0], 2 * jac.shape[1]))
  534. jac_tmp[1::2, 1::2] = jac_tmp[::2, ::2] = real(jac)
  535. jac_tmp[1::2, ::2] = imag(jac)
  536. jac_tmp[::2, 1::2] = -jac_tmp[1::2, ::2]
  537. ml = getattr(self._integrator, 'ml', None)
  538. mu = getattr(self._integrator, 'mu', None)
  539. if ml is not None or mu is not None:
  540. # Jacobian is banded. The user's Jacobian function has computed
  541. # the complex Jacobian in packed format. The corresponding
  542. # real-valued version has every other column shifted up.
  543. jac_tmp = _transform_banded_jac(jac_tmp)
  544. return jac_tmp
  545. @property
  546. def y(self):
  547. return self._y[::2] + 1j * self._y[1::2]
  548. def set_integrator(self, name, **integrator_params):
  549. """
  550. Set integrator by name.
  551. Parameters
  552. ----------
  553. name : str
  554. Name of the integrator
  555. **integrator_params
  556. Additional parameters for the integrator.
  557. """
  558. if name == 'zvode':
  559. raise ValueError("zvode must be used with ode, not complex_ode")
  560. lband = integrator_params.get('lband')
  561. uband = integrator_params.get('uband')
  562. if lband is not None or uband is not None:
  563. # The Jacobian is banded. Override the user-supplied bandwidths
  564. # (which are for the complex Jacobian) with the bandwidths of
  565. # the corresponding real-valued Jacobian wrapper of the complex
  566. # Jacobian.
  567. integrator_params['lband'] = 2 * (lband or 0) + 1
  568. integrator_params['uband'] = 2 * (uband or 0) + 1
  569. return ode.set_integrator(self, name, **integrator_params)
  570. def set_initial_value(self, y, t=0.0):
  571. """Set initial conditions y(t) = y."""
  572. y = asarray(y)
  573. self.tmp = zeros(y.size * 2, 'float')
  574. self.tmp[::2] = real(y)
  575. self.tmp[1::2] = imag(y)
  576. # Create separate work array for derivatives to avoid aliasing issues
  577. self.tmp_derivative = zeros(y.size * 2, 'float')
  578. return ode.set_initial_value(self, self.tmp, t)
  579. def integrate(self, t, step=False, relax=False):
  580. """Find y=y(t), set y as an initial condition, and return y.
  581. Parameters
  582. ----------
  583. t : float
  584. The endpoint of the integration step.
  585. step : bool
  586. If True, and if the integrator supports the step method,
  587. then perform a single integration step and return.
  588. This parameter is provided in order to expose internals of
  589. the implementation, and should not be changed from its default
  590. value in most cases.
  591. relax : bool
  592. If True and if the integrator supports the run_relax method,
  593. then integrate until t_1 >= t and return. ``relax`` is not
  594. referenced if ``step=True``.
  595. This parameter is provided in order to expose internals of
  596. the implementation, and should not be changed from its default
  597. value in most cases.
  598. Returns
  599. -------
  600. y : float
  601. The integrated value at t
  602. """
  603. y = ode.integrate(self, t, step, relax)
  604. return y[::2] + 1j * y[1::2]
  605. def set_solout(self, solout):
  606. """
  607. Set callable to be called at every successful integration step.
  608. Parameters
  609. ----------
  610. solout : callable
  611. ``solout(t, y)`` is called at each internal integrator step,
  612. t is a scalar providing the current independent position
  613. y is the current solution ``y.shape == (n,)``
  614. solout should return -1 to stop integration
  615. otherwise it should return None or 0
  616. """
  617. if self._integrator.supports_solout:
  618. self._integrator.set_solout(solout, complex=True)
  619. else:
  620. raise TypeError("selected integrator does not support solouta, "
  621. "choose another one")
  622. # ------------------------------------------------------------------------------
  623. # ODE integrators
  624. # ------------------------------------------------------------------------------
  625. def find_integrator(name):
  626. for cl in IntegratorBase.integrator_classes:
  627. if re.match(name, cl.__name__, re.I):
  628. return cl
  629. return None
  630. class IntegratorConcurrencyError(RuntimeError):
  631. """
  632. Failure due to concurrent usage of an integrator that can be used
  633. only for a single problem at a time.
  634. """
  635. def __init__(self, name):
  636. msg = (f"Integrator `{name}` can be used to solve only a single problem "
  637. "at a time. If you want to integrate multiple problems, "
  638. "consider using a different integrator (see `ode.set_integrator`)")
  639. RuntimeError.__init__(self, msg)
  640. class IntegratorBase:
  641. runner = None # runner is None => integrator is not available
  642. success = None # success==1 if integrator was called successfully
  643. istate = None # istate > 0 means success, istate < 0 means failure
  644. supports_run_relax = None
  645. supports_step = None
  646. supports_solout = False
  647. integrator_classes = []
  648. scalar = float
  649. # generic type compatibility with scipy-stubs
  650. __class_getitem__ = classmethod(types.GenericAlias)
  651. def acquire_new_handle(self):
  652. # Some of the integrators have internal state (ancient
  653. # Fortran...), and so only one instance can use them at a time.
  654. # We keep track of this, and fail when concurrent usage is tried.
  655. self.__class__.active_global_handle += 1
  656. self.handle = self.__class__.active_global_handle
  657. def check_handle(self):
  658. if self.handle is not self.__class__.active_global_handle:
  659. raise IntegratorConcurrencyError(self.__class__.__name__)
  660. def reset(self, n, has_jac):
  661. """Prepare integrator for call: allocate memory, set flags, etc.
  662. n - number of equations.
  663. has_jac - if user has supplied function for evaluating Jacobian.
  664. """
  665. def run(self, f, jac, y0, t0, t1, f_params, jac_params):
  666. """Integrate from t=t0 to t=t1 using y0 as an initial condition.
  667. Return 2-tuple (y1,t1) where y1 is the result and t=t1
  668. defines the stoppage coordinate of the result.
  669. """
  670. raise NotImplementedError('all integrators must define '
  671. 'run(f, jac, t0, t1, y0, f_params, jac_params)')
  672. def step(self, f, jac, y0, t0, t1, f_params, jac_params):
  673. """Make one integration step and return (y1,t1)."""
  674. raise NotImplementedError(f'{self.__class__.__name__} '
  675. 'does not support step() method')
  676. def run_relax(self, f, jac, y0, t0, t1, f_params, jac_params):
  677. """Integrate from t=t0 to t>=t1 and return (y1,t)."""
  678. raise NotImplementedError(f'{self.__class__.__name__} '
  679. 'does not support run_relax() method')
  680. # XXX: __str__ method for getting visual state of the integrator
  681. class vode(IntegratorBase):
  682. runner = getattr(_vode, 'dvode', None)
  683. messages = {-1: 'Excess work done on this call. (Perhaps wrong MF.)',
  684. -2: 'Excess accuracy requested. (Tolerances too small.)',
  685. -3: 'Illegal input detected. (See printed message.)',
  686. -4: 'Repeated error test failures. (Check all input.)',
  687. -5: 'Repeated convergence failures. (Perhaps bad'
  688. ' Jacobian supplied or wrong choice of MF or tolerances.)',
  689. -6: 'Error weight became zero during problem. (Solution'
  690. ' component i vanished, and ATOL or ATOL(i) = 0.)'
  691. }
  692. supports_run_relax = 1
  693. supports_step = 1
  694. def __init__(self,
  695. method='adams',
  696. with_jacobian=False,
  697. rtol=1e-6, atol=1e-12,
  698. lband=None, uband=None,
  699. order=12,
  700. nsteps=500,
  701. max_step=0.0, # corresponds to infinite
  702. min_step=0.0,
  703. first_step=0.0, # determined by solver
  704. ):
  705. if re.match(method, r'adams', re.I):
  706. self.meth = 1
  707. elif re.match(method, r'bdf', re.I):
  708. self.meth = 2
  709. else:
  710. raise ValueError(f'Unknown integration method {method}')
  711. self.with_jacobian = with_jacobian
  712. self.rtol = rtol
  713. self.atol = atol
  714. self.mu = uband
  715. self.ml = lband
  716. self.order = order
  717. self.nsteps = nsteps
  718. self.max_step = max_step
  719. self.min_step = min_step
  720. self.first_step = first_step
  721. self.success = 1
  722. # State persistence arrays for VODE internal state
  723. # These retain the solver state between calls
  724. self.state_doubles = zeros(51, dtype=np.float64) # VODE_STATE_DOUBLE_SIZE
  725. self.state_ints = zeros(41, dtype=np.int32) # VODE_STATE_INT_SIZE
  726. def _determine_mf_and_set_bands(self, has_jac):
  727. """
  728. Determine the `MF` parameter (Method Flag) for the Fortran subroutine `dvode`.
  729. In the Fortran code, the legal values of `MF` are:
  730. 10, 11, 12, 13, 14, 15, 20, 21, 22, 23, 24, 25,
  731. -11, -12, -14, -15, -21, -22, -24, -25
  732. but this Python wrapper does not use negative values.
  733. Returns
  734. mf = 10*self.meth + miter
  735. self.meth is the linear multistep method:
  736. self.meth == 1: method="adams"
  737. self.meth == 2: method="bdf"
  738. miter is the correction iteration method:
  739. miter == 0: Functional iteration; no Jacobian involved.
  740. miter == 1: Chord iteration with user-supplied full Jacobian.
  741. miter == 2: Chord iteration with internally computed full Jacobian.
  742. miter == 3: Chord iteration with internally computed diagonal Jacobian.
  743. miter == 4: Chord iteration with user-supplied banded Jacobian.
  744. miter == 5: Chord iteration with internally computed banded Jacobian.
  745. Side effects: If either self.mu or self.ml is not None and the other is None,
  746. then the one that is None is set to 0.
  747. """
  748. jac_is_banded = self.mu is not None or self.ml is not None
  749. if jac_is_banded:
  750. if self.mu is None:
  751. self.mu = 0
  752. if self.ml is None:
  753. self.ml = 0
  754. # has_jac is True if the user provided a Jacobian function.
  755. if has_jac:
  756. if jac_is_banded:
  757. miter = 4
  758. else:
  759. miter = 1
  760. else:
  761. if jac_is_banded:
  762. if self.ml == self.mu == 0:
  763. miter = 3 # Chord iteration with internal diagonal Jacobian.
  764. else:
  765. miter = 5 # Chord iteration with internal banded Jacobian.
  766. else:
  767. # self.with_jacobian is set by the user in
  768. # the call to ode.set_integrator.
  769. if self.with_jacobian:
  770. miter = 2 # Chord iteration with internal full Jacobian.
  771. else:
  772. miter = 0 # Functional iteration; no Jacobian involved.
  773. mf = 10 * self.meth + miter
  774. return mf
  775. def reset(self, n, has_jac):
  776. mf = self._determine_mf_and_set_bands(has_jac)
  777. if mf == 10:
  778. lrw = 20 + 16 * n
  779. elif mf in [11, 12]:
  780. lrw = 22 + 16 * n + 2 * n * n
  781. elif mf == 13:
  782. lrw = 22 + 17 * n
  783. elif mf in [14, 15]:
  784. lrw = 22 + 18 * n + (3 * self.ml + 2 * self.mu) * n
  785. elif mf == 20:
  786. lrw = 20 + 9 * n
  787. elif mf in [21, 22]:
  788. lrw = 22 + 9 * n + 2 * n * n
  789. elif mf == 23:
  790. lrw = 22 + 10 * n
  791. elif mf in [24, 25]:
  792. lrw = 22 + 11 * n + (3 * self.ml + 2 * self.mu) * n
  793. else:
  794. raise ValueError(f'Unexpected mf={mf}')
  795. if mf % 10 in [0, 3]:
  796. liw = 30
  797. else:
  798. liw = 30 + n
  799. rwork = zeros((lrw,), float)
  800. rwork[4] = self.first_step
  801. rwork[5] = self.max_step
  802. rwork[6] = self.min_step
  803. self.rwork = rwork
  804. iwork = zeros((liw,), dtype=np.int32)
  805. if self.ml is not None:
  806. iwork[0] = self.ml
  807. if self.mu is not None:
  808. iwork[1] = self.mu
  809. iwork[4] = self.order
  810. iwork[5] = self.nsteps
  811. iwork[6] = 2 # mxhnil
  812. self.iwork = iwork
  813. self.call_args = [self.rtol, self.atol, 1, 1,
  814. self.rwork, self.iwork, mf]
  815. self.success = 1
  816. self.initialized = False
  817. # Zero state arrays on reset to avoid contamination from previous problems.
  818. # State persistence works within a single integration (istate=2), but between
  819. # different problems (different n etc.), state needs to be cleared.
  820. self.state_doubles.fill(0.0)
  821. self.state_ints.fill(0)
  822. def run(self, f, jac, y0, t0, t1, f_params, jac_params):
  823. # Note: For banded Jacobians, the user provides the compressed format
  824. # (ml + mu + 1, n), and the C code handles padding to the expanded
  825. # format (ml + 2*mu + 1, n) internally. No Python wrapper needed.
  826. # VODE C wrapper signature:
  827. # dvode(f, jac, y0, t0, t1, rtol, atol, itask, istate, rwork, iwork, mf,
  828. # f_params, jac_params, state_doubles, state_ints)
  829. args = ((f, jac, y0, t0, t1) + tuple(self.call_args) +
  830. (f_params, jac_params, self.state_doubles, self.state_ints))
  831. y1, t, istate = self.runner(*args)
  832. self.istate = istate
  833. if istate < 0:
  834. unexpected_istate_msg = f'Unexpected istate={istate:d}'
  835. warnings.warn(f'{self.__class__.__name__:s}: '
  836. f'{self.messages.get(istate, unexpected_istate_msg):s}',
  837. stacklevel=2)
  838. self.success = 0
  839. else:
  840. self.call_args[3] = 2 # upgrade istate from 1 to 2
  841. self.istate = 2
  842. return y1, t
  843. def step(self, *args):
  844. itask = self.call_args[2]
  845. self.call_args[2] = 2
  846. r = self.run(*args)
  847. self.call_args[2] = itask
  848. return r
  849. def run_relax(self, *args):
  850. itask = self.call_args[2]
  851. self.call_args[2] = 3
  852. r = self.run(*args)
  853. self.call_args[2] = itask
  854. return r
  855. if vode.runner is not None:
  856. IntegratorBase.integrator_classes.append(vode)
  857. class zvode(vode):
  858. runner = getattr(_vode, 'zvode', None)
  859. supports_run_relax = 1
  860. supports_step = 1
  861. scalar = complex
  862. __class_getitem__ = None
  863. def __init__(self, *args, **kwargs):
  864. super().__init__(*args, **kwargs)
  865. # Override state array sizes for ZVODE (53 doubles vs 51 for VODE)
  866. self.state_doubles = zeros(53, dtype=np.float64) # ZVODE_STATE_DOUBLE_SIZE
  867. self.state_ints = zeros(41, dtype=np.int32) # ZVODE_STATE_INT_SIZE
  868. def reset(self, n, has_jac):
  869. mf = self._determine_mf_and_set_bands(has_jac)
  870. if mf in (10,):
  871. lzw = 15 * n
  872. elif mf in (11, 12):
  873. lzw = 15 * n + 2 * n ** 2
  874. elif mf in (-11, -12):
  875. lzw = 15 * n + n ** 2
  876. elif mf in (13,):
  877. lzw = 16 * n
  878. elif mf in (14, 15):
  879. lzw = 17 * n + (3 * self.ml + 2 * self.mu) * n
  880. elif mf in (-14, -15):
  881. lzw = 16 * n + (2 * self.ml + self.mu) * n
  882. elif mf in (20,):
  883. lzw = 8 * n
  884. elif mf in (21, 22):
  885. lzw = 8 * n + 2 * n ** 2
  886. elif mf in (-21, -22):
  887. lzw = 8 * n + n ** 2
  888. elif mf in (23,):
  889. lzw = 9 * n
  890. elif mf in (24, 25):
  891. lzw = 10 * n + (3 * self.ml + 2 * self.mu) * n
  892. elif mf in (-24, -25):
  893. lzw = 9 * n + (2 * self.ml + self.mu) * n
  894. lrw = 20 + n
  895. if mf % 10 in (0, 3):
  896. liw = 30
  897. else:
  898. liw = 30 + n
  899. zwork = zeros((lzw,), complex)
  900. self.zwork = zwork
  901. rwork = zeros((lrw,), float)
  902. rwork[4] = self.first_step
  903. rwork[5] = self.max_step
  904. rwork[6] = self.min_step
  905. self.rwork = rwork
  906. iwork = zeros((liw,), np.int32)
  907. if self.ml is not None:
  908. iwork[0] = self.ml
  909. if self.mu is not None:
  910. iwork[1] = self.mu
  911. iwork[4] = self.order
  912. iwork[5] = self.nsteps
  913. iwork[6] = 2 # mxhnil
  914. self.iwork = iwork
  915. self.call_args = [self.rtol, self.atol, 1, 1,
  916. self.zwork, self.rwork, self.iwork, mf]
  917. self.success = 1
  918. self.initialized = False
  919. # Zero state arrays on reset to avoid contamination from previous problems.
  920. # State persistence works within a single integration (istate=2), but between
  921. # different problems (different n etc.), state needs to be cleared.
  922. self.state_doubles.fill(0.0)
  923. self.state_ints.fill(0)
  924. if zvode.runner is not None:
  925. IntegratorBase.integrator_classes.append(zvode)
  926. class dopri5(IntegratorBase):
  927. runner = getattr(_dop, 'dopri5', None)
  928. name = 'dopri5'
  929. supports_solout = True
  930. messages = {1: 'computation successful',
  931. 2: 'computation successful (interrupted by solout)',
  932. -1: 'input is not consistent',
  933. -2: 'larger nsteps is needed',
  934. -3: 'step size becomes too small',
  935. -4: 'problem is probably stiff (interrupted)',
  936. }
  937. __class_getitem__ = None
  938. def __init__(self,
  939. rtol=1e-6, atol=1e-12,
  940. nsteps=500,
  941. max_step=0.0,
  942. first_step=0.0, # determined by solver
  943. safety=0.9,
  944. ifactor=10.0,
  945. dfactor=0.2,
  946. beta=0.0,
  947. method=None,
  948. verbosity=-1, # no messages if negative
  949. ):
  950. self.rtol = rtol
  951. self.atol = atol
  952. self.nsteps = nsteps
  953. self.max_step = max_step
  954. self.first_step = first_step
  955. self.safety = safety
  956. self.ifactor = ifactor
  957. self.dfactor = dfactor
  958. self.beta = beta
  959. self.verbosity = verbosity
  960. self.success = 1
  961. self.set_solout(None)
  962. def set_solout(self, solout, complex=False):
  963. self.solout = solout
  964. self.solout_cmplx = complex
  965. if solout is None:
  966. self.iout = 0
  967. else:
  968. self.iout = 1
  969. def reset(self, n, has_jac):
  970. work = zeros((8 * n + 21,), float)
  971. work[1] = self.safety
  972. work[2] = self.dfactor
  973. work[3] = self.ifactor
  974. work[4] = self.beta
  975. work[5] = self.max_step
  976. work[6] = self.first_step
  977. self.work = work
  978. self.iwork = zeros((21,), dtype=np.int32)
  979. self.call_args = [self.rtol, self.atol, self._solout,
  980. self.iout, self.work, self.iwork,
  981. self.nsteps, self.verbosity]
  982. self.success = 1
  983. def run(self, f, jac, y0, t0, t1, f_params, jac_params):
  984. x, y, istate = self.runner(*((f, t0, y0, t1) +
  985. tuple(self.call_args) + (f_params,)))
  986. self.istate = istate
  987. if istate < 0:
  988. unexpected_istate_msg = f'Unexpected istate={istate:d}'
  989. warnings.warn(f'{self.__class__.__name__:s}: '
  990. f'{self.messages.get(istate, unexpected_istate_msg):s}',
  991. stacklevel=2)
  992. self.success = 0
  993. return y, x
  994. def _solout(self, x, y):
  995. if self.solout is not None:
  996. if self.solout_cmplx:
  997. y = y[::2] + 1j * y[1::2]
  998. return self.solout(x, y)
  999. else:
  1000. return 1
  1001. if dopri5.runner is not None:
  1002. IntegratorBase.integrator_classes.append(dopri5)
  1003. class dop853(dopri5):
  1004. runner = getattr(_dop, 'dopri853', None)
  1005. name = 'dop853'
  1006. def __init__(self,
  1007. rtol=1e-6, atol=1e-12,
  1008. nsteps=500,
  1009. max_step=0.0,
  1010. first_step=0.0, # determined by solver
  1011. safety=0.9,
  1012. ifactor=6.0,
  1013. dfactor=0.3,
  1014. beta=0.0,
  1015. method=None,
  1016. verbosity=-1, # no messages if negative
  1017. ):
  1018. super().__init__(rtol, atol, nsteps, max_step, first_step, safety,
  1019. ifactor, dfactor, beta, method, verbosity)
  1020. def reset(self, n, has_jac):
  1021. work = zeros((11 * n + 21,), float)
  1022. work[1] = self.safety
  1023. work[2] = self.dfactor
  1024. work[3] = self.ifactor
  1025. work[4] = self.beta
  1026. work[5] = self.max_step
  1027. work[6] = self.first_step
  1028. self.work = work
  1029. self.iwork = zeros((21,), dtype=np.int32)
  1030. self.call_args = [self.rtol, self.atol, self._solout,
  1031. self.iout, self.work, self.iwork,
  1032. self.nsteps, self.verbosity]
  1033. self.success = 1
  1034. if dop853.runner is not None:
  1035. IntegratorBase.integrator_classes.append(dop853)
  1036. class lsoda(IntegratorBase):
  1037. runner = lsoda_step # Use low-level lsoda wrapper
  1038. messages = {
  1039. 2: "Integration successful.",
  1040. -1: "Excess work done on this call (perhaps wrong Dfun type).",
  1041. -2: "Excess accuracy requested (tolerances too small).",
  1042. -3: "Illegal input detected (internal error).",
  1043. -4: "Repeated error test failures (internal error).",
  1044. -5: "Repeated convergence failures (perhaps bad Jacobian or tolerances).",
  1045. -6: "Error weight became zero during problem.",
  1046. -7: "Internal workspace insufficient to finish (internal error)."
  1047. }
  1048. __class_getitem__ = None
  1049. def __init__(self,
  1050. with_jacobian=False,
  1051. rtol=1e-6, atol=1e-12,
  1052. lband=None, uband=None,
  1053. nsteps=500,
  1054. max_step=0.0, # corresponds to infinite
  1055. min_step=0.0,
  1056. first_step=0.0, # determined by solver
  1057. ixpr=0,
  1058. max_hnil=0,
  1059. max_order_ns=12,
  1060. max_order_s=5,
  1061. method=None
  1062. ):
  1063. self.with_jacobian = with_jacobian
  1064. self.rtol = rtol
  1065. self.atol = atol
  1066. self.mu = uband
  1067. self.ml = lband
  1068. self.max_order_ns = max_order_ns
  1069. self.max_order_s = max_order_s
  1070. self.nsteps = nsteps
  1071. self.max_step = max_step
  1072. self.min_step = min_step
  1073. self.first_step = first_step
  1074. self.ixpr = ixpr
  1075. self.max_hnil = max_hnil
  1076. self.success = 1
  1077. self.initialized = False
  1078. # State persistence arrays for LSODA internal state
  1079. # These retain the solver state between calls
  1080. self.state_doubles = zeros(240, dtype=np.float64) # LSODA_STATE_DOUBLE_SIZE
  1081. self.state_ints = zeros(48, dtype=np.int32) # LSODA_STATE_INT_SIZE
  1082. def reset(self, n, has_jac):
  1083. # Zero state arrays on reset to avoid contamination from previous steps.
  1084. # State persistence works within a single integration (istate=2), but between
  1085. # different problems (different n etc.), state needs to be cleared.
  1086. self.state_doubles.fill(0.0)
  1087. self.state_ints.fill(0)
  1088. # Calculate parameters for lsoda subroutine.
  1089. # jt values: 1=user full, 2=FD full, 4=user banded, 5=FD banded (3=invalid)
  1090. if has_jac:
  1091. if self.mu is None and self.ml is None:
  1092. jt = 1 # User-supplied full Jacobian
  1093. else:
  1094. if self.mu is None:
  1095. self.mu = 0
  1096. if self.ml is None:
  1097. self.ml = 0
  1098. jt = 4 # User-supplied banded Jacobian
  1099. else:
  1100. if self.mu is None and self.ml is None:
  1101. jt = 2 # Internally generated full Jacobian (finite differences)
  1102. else:
  1103. if self.mu is None:
  1104. self.mu = 0
  1105. if self.ml is None:
  1106. self.ml = 0
  1107. jt = 5 # Internally generated banded Jacobian (finite differences)
  1108. # Calculate work array sizes
  1109. lrn = 20 + (self.max_order_ns + 4) * n
  1110. if jt in [1, 2]:
  1111. lrs = 22 + (self.max_order_s + 4) * n + n * n
  1112. elif jt in [4, 5]:
  1113. lrs = 22 + (self.max_order_s + 5 + 2 * self.ml + self.mu) * n
  1114. else:
  1115. raise ValueError(f'Unexpected jt={jt}')
  1116. lrw = max(lrn, lrs)
  1117. liw = 20 + n
  1118. # Create and initialize work arrays
  1119. rwork = zeros((lrw,), float)
  1120. rwork[4] = self.first_step
  1121. rwork[5] = self.max_step
  1122. rwork[6] = self.min_step
  1123. self.rwork = rwork
  1124. iwork = zeros((liw,), dtype=np.int32)
  1125. if self.ml is not None:
  1126. iwork[0] = self.ml
  1127. if self.mu is not None:
  1128. iwork[1] = self.mu
  1129. iwork[4] = self.ixpr
  1130. iwork[5] = self.nsteps
  1131. iwork[6] = self.max_hnil
  1132. iwork[7] = self.max_order_ns
  1133. iwork[8] = self.max_order_s
  1134. self.iwork = iwork
  1135. self.call_args = [self.rtol, self.atol, 1, 1,
  1136. self.rwork, self.iwork, jt]
  1137. self.success = 1
  1138. def run(self, f, jac, y0, t0, t1, f_params, jac_params):
  1139. # Prepare arguments for low-level lsoda wrapper
  1140. rtol = self.call_args[0]
  1141. atol = self.call_args[1]
  1142. itask = self.call_args[2]
  1143. istate = self.call_args[3]
  1144. rwork = self.call_args[4]
  1145. iwork = self.call_args[5]
  1146. jt = self.call_args[6]
  1147. # Note: For banded Jacobians, the user provides the compressed format
  1148. # (ml + mu + 1, n), and the C code handles padding to the expanded
  1149. # format (2*ml + mu + 1, n) internally. No Python wrapper needed.
  1150. # Signature:
  1151. #
  1152. # lsoda(fun, y0, t, tout, rtol, atol, itask, istate, rwork, iwork,
  1153. # jac, jt, f_params, tfirst, jac_params, state_doubles, state_ints)
  1154. #
  1155. # "state_doubles" and "state_ints" are arrays for passing internal
  1156. # state between Python-C code.
  1157. y1, t, istate = self.runner(
  1158. f, y0, t0, t1, rtol, atol, itask, istate, rwork, iwork,
  1159. jac, jt, f_params, 1, jac_params, # tfirst=1 for (t, y) signature
  1160. self.state_doubles, self.state_ints
  1161. )
  1162. self.istate = istate
  1163. if istate < 0:
  1164. unexpected_istate_msg = f'Unexpected istate={istate:d}'
  1165. warnings.warn(f'{self.__class__.__name__:s}: '
  1166. f'{self.messages.get(istate, unexpected_istate_msg):s}',
  1167. stacklevel=2)
  1168. self.success = 0
  1169. else:
  1170. self.call_args[3] = 2 # upgrade istate from 1 to 2
  1171. self.istate = 2
  1172. self.success = 1
  1173. return y1, t
  1174. def step(self, *args):
  1175. itask = self.call_args[2]
  1176. self.call_args[2] = 2
  1177. r = self.run(*args)
  1178. self.call_args[2] = itask
  1179. return r
  1180. def run_relax(self, *args):
  1181. itask = self.call_args[2]
  1182. self.call_args[2] = 3
  1183. r = self.run(*args)
  1184. self.call_args[2] = itask
  1185. return r
  1186. if lsoda.runner:
  1187. IntegratorBase.integrator_classes.append(lsoda)