_odrpack.py 41 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142
  1. """
  2. Python wrappers for Orthogonal Distance Regression (ODRPACK).
  3. Notes
  4. =====
  5. * Array formats -- FORTRAN stores its arrays in memory column first, i.e., an
  6. array element A(i, j, k) will be next to A(i+1, j, k). In C and, consequently,
  7. NumPy, arrays are stored row first: A[i, j, k] is next to A[i, j, k+1]. For
  8. efficiency and convenience, the input and output arrays of the fitting
  9. function (and its Jacobians) are passed to FORTRAN without transposition.
  10. Therefore, where the ODRPACK documentation says that the X array is of shape
  11. (N, M), it will be passed to the Python function as an array of shape (M, N).
  12. If M==1, the 1-D case, then nothing matters; if M>1, then your
  13. Python functions will be dealing with arrays that are indexed in reverse of
  14. the ODRPACK documentation. No real issue, but watch out for your indexing of
  15. the Jacobians: the i,jth elements (@f_i/@x_j) evaluated at the nth
  16. observation will be returned as jacd[j, i, n]. Except for the Jacobians, it
  17. really is easier to deal with x[0] and x[1] than x[:,0] and x[:,1]. Of course,
  18. you can always use the transpose() function from SciPy explicitly.
  19. * Examples -- See the accompanying file test/test.py for examples of how to set
  20. up fits of your own. Some are taken from the User's Guide; some are from
  21. other sources.
  22. * Models -- Some common models are instantiated in the accompanying module
  23. models.py . Contributions are welcome.
  24. Credits
  25. =======
  26. * Thanks to Arnold Moene and Gerard Vermeulen for fixing some killer bugs.
  27. Robert Kern
  28. robert.kern@gmail.com
  29. """
  30. import os
  31. import numpy
  32. from warnings import warn
  33. from scipy.odr import __odrpack
  34. __all__ = ['odr', 'OdrWarning', 'OdrError', 'OdrStop',
  35. 'Data', 'RealData', 'Model', 'Output', 'ODR',
  36. 'odr_error', 'odr_stop']
  37. odr = __odrpack.odr
  38. class OdrWarning(UserWarning):
  39. """
  40. Warning indicating that the data passed into
  41. ODR will cause problems when passed into 'odr'
  42. that the user should be aware of.
  43. """
  44. pass
  45. class OdrError(Exception):
  46. """
  47. Exception indicating an error in fitting.
  48. This is raised by `~scipy.odr.odr` if an error occurs during fitting.
  49. """
  50. pass
  51. class OdrStop(Exception):
  52. """
  53. Exception stopping fitting.
  54. You can raise this exception in your objective function to tell
  55. `~scipy.odr.odr` to stop fitting.
  56. """
  57. pass
  58. # Backwards compatibility
  59. odr_error = OdrError
  60. odr_stop = OdrStop
  61. __odrpack._set_exceptions(OdrError, OdrStop)
  62. def _conv(obj, dtype=None):
  63. """ Convert an object to the preferred form for input to the odr routine.
  64. """
  65. if obj is None:
  66. return obj
  67. else:
  68. if dtype is None:
  69. obj = numpy.asarray(obj)
  70. else:
  71. obj = numpy.asarray(obj, dtype)
  72. if obj.shape == ():
  73. # Scalar.
  74. return obj.dtype.type(obj)
  75. else:
  76. return obj
  77. def _report_error(info):
  78. """ Interprets the return code of the odr routine.
  79. Parameters
  80. ----------
  81. info : int
  82. The return code of the odr routine.
  83. Returns
  84. -------
  85. problems : list(str)
  86. A list of messages about why the odr() routine stopped.
  87. """
  88. stopreason = ('Blank',
  89. 'Sum of squares convergence',
  90. 'Parameter convergence',
  91. 'Both sum of squares and parameter convergence',
  92. 'Iteration limit reached')[info % 5]
  93. if info >= 5:
  94. # questionable results or fatal error
  95. I = (info//10000 % 10,
  96. info//1000 % 10,
  97. info//100 % 10,
  98. info//10 % 10,
  99. info % 10)
  100. problems = []
  101. if I[0] == 0:
  102. if I[1] != 0:
  103. problems.append('Derivatives possibly not correct')
  104. if I[2] != 0:
  105. problems.append('Error occurred in callback')
  106. if I[3] != 0:
  107. problems.append('Problem is not full rank at solution')
  108. problems.append(stopreason)
  109. elif I[0] == 1:
  110. if I[1] != 0:
  111. problems.append('N < 1')
  112. if I[2] != 0:
  113. problems.append('M < 1')
  114. if I[3] != 0:
  115. problems.append('NP < 1 or NP > N')
  116. if I[4] != 0:
  117. problems.append('NQ < 1')
  118. elif I[0] == 2:
  119. if I[1] != 0:
  120. problems.append('LDY and/or LDX incorrect')
  121. if I[2] != 0:
  122. problems.append('LDWE, LD2WE, LDWD, and/or LD2WD incorrect')
  123. if I[3] != 0:
  124. problems.append('LDIFX, LDSTPD, and/or LDSCLD incorrect')
  125. if I[4] != 0:
  126. problems.append('LWORK and/or LIWORK too small')
  127. elif I[0] == 3:
  128. if I[1] != 0:
  129. problems.append('STPB and/or STPD incorrect')
  130. if I[2] != 0:
  131. problems.append('SCLB and/or SCLD incorrect')
  132. if I[3] != 0:
  133. problems.append('WE incorrect')
  134. if I[4] != 0:
  135. problems.append('WD incorrect')
  136. elif I[0] == 4:
  137. problems.append('Error in derivatives')
  138. elif I[0] == 5:
  139. problems.append('Error occurred in callback')
  140. elif I[0] == 6:
  141. problems.append('Numerical error detected')
  142. return problems
  143. else:
  144. return [stopreason]
  145. class Data:
  146. """
  147. The data to fit.
  148. Parameters
  149. ----------
  150. x : array_like
  151. Observed data for the independent variable of the regression
  152. y : array_like, optional
  153. If array-like, observed data for the dependent variable of the
  154. regression. A scalar input implies that the model to be used on
  155. the data is implicit.
  156. we : array_like, optional
  157. If `we` is a scalar, then that value is used for all data points (and
  158. all dimensions of the response variable).
  159. If `we` is a rank-1 array of length q (the dimensionality of the
  160. response variable), then this vector is the diagonal of the covariant
  161. weighting matrix for all data points.
  162. If `we` is a rank-1 array of length n (the number of data points), then
  163. the i'th element is the weight for the i'th response variable
  164. observation (single-dimensional only).
  165. If `we` is a rank-2 array of shape (q, q), then this is the full
  166. covariant weighting matrix broadcast to each observation.
  167. If `we` is a rank-2 array of shape (q, n), then `we[:,i]` is the
  168. diagonal of the covariant weighting matrix for the i'th observation.
  169. If `we` is a rank-3 array of shape (q, q, n), then `we[:,:,i]` is the
  170. full specification of the covariant weighting matrix for each
  171. observation.
  172. If the fit is implicit, then only a positive scalar value is used.
  173. wd : array_like, optional
  174. If `wd` is a scalar, then that value is used for all data points
  175. (and all dimensions of the input variable). If `wd` = 0, then the
  176. covariant weighting matrix for each observation is set to the identity
  177. matrix (so each dimension of each observation has the same weight).
  178. If `wd` is a rank-1 array of length m (the dimensionality of the input
  179. variable), then this vector is the diagonal of the covariant weighting
  180. matrix for all data points.
  181. If `wd` is a rank-1 array of length n (the number of data points), then
  182. the i'th element is the weight for the ith input variable observation
  183. (single-dimensional only).
  184. If `wd` is a rank-2 array of shape (m, m), then this is the full
  185. covariant weighting matrix broadcast to each observation.
  186. If `wd` is a rank-2 array of shape (m, n), then `wd[:,i]` is the
  187. diagonal of the covariant weighting matrix for the ith observation.
  188. If `wd` is a rank-3 array of shape (m, m, n), then `wd[:,:,i]` is the
  189. full specification of the covariant weighting matrix for each
  190. observation.
  191. fix : array_like of ints, optional
  192. The `fix` argument is the same as ifixx in the class ODR. It is an
  193. array of integers with the same shape as data.x that determines which
  194. input observations are treated as fixed. One can use a sequence of
  195. length m (the dimensionality of the input observations) to fix some
  196. dimensions for all observations. A value of 0 fixes the observation,
  197. a value > 0 makes it free.
  198. meta : dict, optional
  199. Free-form dictionary for metadata.
  200. Notes
  201. -----
  202. Each argument is attached to the member of the instance of the same name.
  203. The structures of `x` and `y` are described in the Model class docstring.
  204. If `y` is an integer, then the Data instance can only be used to fit with
  205. implicit models where the dimensionality of the response is equal to the
  206. specified value of `y`.
  207. The `we` argument weights the effect a deviation in the response variable
  208. has on the fit. The `wd` argument weights the effect a deviation in the
  209. input variable has on the fit. To handle multidimensional inputs and
  210. responses easily, the structure of these arguments has the n'th
  211. dimensional axis first. These arguments heavily use the structured
  212. arguments feature of ODRPACK to conveniently and flexibly support all
  213. options. See the ODRPACK User's Guide for a full explanation of how these
  214. weights are used in the algorithm. Basically, a higher value of the weight
  215. for a particular data point makes a deviation at that point more
  216. detrimental to the fit.
  217. """
  218. def __init__(self, x, y=None, we=None, wd=None, fix=None, meta=None):
  219. self.x = _conv(x)
  220. if not isinstance(self.x, numpy.ndarray):
  221. raise ValueError(("Expected an 'ndarray' of data for 'x', "
  222. "but instead got data of type '{name}'").format(
  223. name=type(self.x).__name__))
  224. self.y = _conv(y)
  225. self.we = _conv(we)
  226. self.wd = _conv(wd)
  227. self.fix = _conv(fix)
  228. self.meta = {} if meta is None else meta
  229. def set_meta(self, **kwds):
  230. """ Update the metadata dictionary with the keywords and data provided
  231. by keywords.
  232. Examples
  233. --------
  234. ::
  235. data.set_meta(lab="Ph 7; Lab 26", title="Ag110 + Ag108 Decay")
  236. """
  237. self.meta.update(kwds)
  238. def __getattr__(self, attr):
  239. """ Dispatch attribute access to the metadata dictionary.
  240. """
  241. if attr in self.meta:
  242. return self.meta[attr]
  243. else:
  244. raise AttributeError("'%s' not in metadata" % attr)
  245. class RealData(Data):
  246. """
  247. The data, with weightings as actual standard deviations and/or
  248. covariances.
  249. Parameters
  250. ----------
  251. x : array_like
  252. Observed data for the independent variable of the regression
  253. y : array_like, optional
  254. If array-like, observed data for the dependent variable of the
  255. regression. A scalar input implies that the model to be used on
  256. the data is implicit.
  257. sx : array_like, optional
  258. Standard deviations of `x`.
  259. `sx` are standard deviations of `x` and are converted to weights by
  260. dividing 1.0 by their squares.
  261. sy : array_like, optional
  262. Standard deviations of `y`.
  263. `sy` are standard deviations of `y` and are converted to weights by
  264. dividing 1.0 by their squares.
  265. covx : array_like, optional
  266. Covariance of `x`
  267. `covx` is an array of covariance matrices of `x` and are converted to
  268. weights by performing a matrix inversion on each observation's
  269. covariance matrix.
  270. covy : array_like, optional
  271. Covariance of `y`
  272. `covy` is an array of covariance matrices and are converted to
  273. weights by performing a matrix inversion on each observation's
  274. covariance matrix.
  275. fix : array_like, optional
  276. The argument and member fix is the same as Data.fix and ODR.ifixx:
  277. It is an array of integers with the same shape as `x` that
  278. determines which input observations are treated as fixed. One can
  279. use a sequence of length m (the dimensionality of the input
  280. observations) to fix some dimensions for all observations. A value
  281. of 0 fixes the observation, a value > 0 makes it free.
  282. meta : dict, optional
  283. Free-form dictionary for metadata.
  284. Notes
  285. -----
  286. The weights `wd` and `we` are computed from provided values as follows:
  287. `sx` and `sy` are converted to weights by dividing 1.0 by their squares.
  288. For example, ``wd = 1./numpy.power(`sx`, 2)``.
  289. `covx` and `covy` are arrays of covariance matrices and are converted to
  290. weights by performing a matrix inversion on each observation's covariance
  291. matrix. For example, ``we[i] = numpy.linalg.inv(covy[i])``.
  292. These arguments follow the same structured argument conventions as wd and
  293. we only restricted by their natures: `sx` and `sy` can't be rank-3, but
  294. `covx` and `covy` can be.
  295. Only set *either* `sx` or `covx` (not both). Setting both will raise an
  296. exception. Same with `sy` and `covy`.
  297. """
  298. def __init__(self, x, y=None, sx=None, sy=None, covx=None, covy=None,
  299. fix=None, meta=None):
  300. if (sx is not None) and (covx is not None):
  301. raise ValueError("cannot set both sx and covx")
  302. if (sy is not None) and (covy is not None):
  303. raise ValueError("cannot set both sy and covy")
  304. # Set flags for __getattr__
  305. self._ga_flags = {}
  306. if sx is not None:
  307. self._ga_flags['wd'] = 'sx'
  308. else:
  309. self._ga_flags['wd'] = 'covx'
  310. if sy is not None:
  311. self._ga_flags['we'] = 'sy'
  312. else:
  313. self._ga_flags['we'] = 'covy'
  314. self.x = _conv(x)
  315. if not isinstance(self.x, numpy.ndarray):
  316. raise ValueError(("Expected an 'ndarray' of data for 'x', "
  317. "but instead got data of type '{name}'").format(
  318. name=type(self.x).__name__))
  319. self.y = _conv(y)
  320. self.sx = _conv(sx)
  321. self.sy = _conv(sy)
  322. self.covx = _conv(covx)
  323. self.covy = _conv(covy)
  324. self.fix = _conv(fix)
  325. self.meta = {} if meta is None else meta
  326. def _sd2wt(self, sd):
  327. """ Convert standard deviation to weights.
  328. """
  329. return 1./numpy.power(sd, 2)
  330. def _cov2wt(self, cov):
  331. """ Convert covariance matrix(-ices) to weights.
  332. """
  333. from scipy.linalg import inv
  334. if len(cov.shape) == 2:
  335. return inv(cov)
  336. else:
  337. weights = numpy.zeros(cov.shape, float)
  338. for i in range(cov.shape[-1]): # n
  339. weights[:,:,i] = inv(cov[:,:,i])
  340. return weights
  341. def __getattr__(self, attr):
  342. lookup_tbl = {('wd', 'sx'): (self._sd2wt, self.sx),
  343. ('wd', 'covx'): (self._cov2wt, self.covx),
  344. ('we', 'sy'): (self._sd2wt, self.sy),
  345. ('we', 'covy'): (self._cov2wt, self.covy)}
  346. if attr not in ('wd', 'we'):
  347. if attr in self.meta:
  348. return self.meta[attr]
  349. else:
  350. raise AttributeError("'%s' not in metadata" % attr)
  351. else:
  352. func, arg = lookup_tbl[(attr, self._ga_flags[attr])]
  353. if arg is not None:
  354. return func(*(arg,))
  355. else:
  356. return None
  357. class Model:
  358. """
  359. The Model class stores information about the function you wish to fit.
  360. It stores the function itself, at the least, and optionally stores
  361. functions which compute the Jacobians used during fitting. Also, one
  362. can provide a function that will provide reasonable starting values
  363. for the fit parameters possibly given the set of data.
  364. Parameters
  365. ----------
  366. fcn : function
  367. fcn(beta, x) --> y
  368. fjacb : function
  369. Jacobian of fcn wrt the fit parameters beta.
  370. fjacb(beta, x) --> @f_i(x,B)/@B_j
  371. fjacd : function
  372. Jacobian of fcn wrt the (possibly multidimensional) input
  373. variable.
  374. fjacd(beta, x) --> @f_i(x,B)/@x_j
  375. extra_args : tuple, optional
  376. If specified, `extra_args` should be a tuple of extra
  377. arguments to pass to `fcn`, `fjacb`, and `fjacd`. Each will be called
  378. by `apply(fcn, (beta, x) + extra_args)`
  379. estimate : array_like of rank-1
  380. Provides estimates of the fit parameters from the data
  381. estimate(data) --> estbeta
  382. implicit : boolean
  383. If TRUE, specifies that the model
  384. is implicit; i.e `fcn(beta, x)` ~= 0 and there is no y data to fit
  385. against
  386. meta : dict, optional
  387. freeform dictionary of metadata for the model
  388. Notes
  389. -----
  390. Note that the `fcn`, `fjacb`, and `fjacd` operate on NumPy arrays and
  391. return a NumPy array. The `estimate` object takes an instance of the
  392. Data class.
  393. Here are the rules for the shapes of the argument and return
  394. arrays of the callback functions:
  395. `x`
  396. if the input data is single-dimensional, then `x` is rank-1
  397. array; i.e., ``x = array([1, 2, 3, ...]); x.shape = (n,)``
  398. If the input data is multi-dimensional, then `x` is a rank-2 array;
  399. i.e., ``x = array([[1, 2, ...], [2, 4, ...]]); x.shape = (m, n)``.
  400. In all cases, it has the same shape as the input data array passed to
  401. `~scipy.odr.odr`. `m` is the dimensionality of the input data,
  402. `n` is the number of observations.
  403. `y`
  404. if the response variable is single-dimensional, then `y` is a
  405. rank-1 array, i.e., ``y = array([2, 4, ...]); y.shape = (n,)``.
  406. If the response variable is multi-dimensional, then `y` is a rank-2
  407. array, i.e., ``y = array([[2, 4, ...], [3, 6, ...]]); y.shape =
  408. (q, n)`` where `q` is the dimensionality of the response variable.
  409. `beta`
  410. rank-1 array of length `p` where `p` is the number of parameters;
  411. i.e. ``beta = array([B_1, B_2, ..., B_p])``
  412. `fjacb`
  413. if the response variable is multi-dimensional, then the
  414. return array's shape is `(q, p, n)` such that ``fjacb(x,beta)[l,k,i] =
  415. d f_l(X,B)/d B_k`` evaluated at the ith data point. If `q == 1`, then
  416. the return array is only rank-2 and with shape `(p, n)`.
  417. `fjacd`
  418. as with fjacb, only the return array's shape is `(q, m, n)`
  419. such that ``fjacd(x,beta)[l,j,i] = d f_l(X,B)/d X_j`` at the ith data
  420. point. If `q == 1`, then the return array's shape is `(m, n)`. If
  421. `m == 1`, the shape is (q, n). If `m == q == 1`, the shape is `(n,)`.
  422. """
  423. def __init__(self, fcn, fjacb=None, fjacd=None,
  424. extra_args=None, estimate=None, implicit=0, meta=None):
  425. self.fcn = fcn
  426. self.fjacb = fjacb
  427. self.fjacd = fjacd
  428. if extra_args is not None:
  429. extra_args = tuple(extra_args)
  430. self.extra_args = extra_args
  431. self.estimate = estimate
  432. self.implicit = implicit
  433. self.meta = meta if meta is not None else {}
  434. def set_meta(self, **kwds):
  435. """ Update the metadata dictionary with the keywords and data provided
  436. here.
  437. Examples
  438. --------
  439. set_meta(name="Exponential", equation="y = a exp(b x) + c")
  440. """
  441. self.meta.update(kwds)
  442. def __getattr__(self, attr):
  443. """ Dispatch attribute access to the metadata.
  444. """
  445. if attr in self.meta:
  446. return self.meta[attr]
  447. else:
  448. raise AttributeError("'%s' not in metadata" % attr)
  449. class Output:
  450. """
  451. The Output class stores the output of an ODR run.
  452. Attributes
  453. ----------
  454. beta : ndarray
  455. Estimated parameter values, of shape (q,).
  456. sd_beta : ndarray
  457. Standard deviations of the estimated parameters, of shape (p,).
  458. cov_beta : ndarray
  459. Covariance matrix of the estimated parameters, of shape (p,p).
  460. delta : ndarray, optional
  461. Array of estimated errors in input variables, of same shape as `x`.
  462. eps : ndarray, optional
  463. Array of estimated errors in response variables, of same shape as `y`.
  464. xplus : ndarray, optional
  465. Array of ``x + delta``.
  466. y : ndarray, optional
  467. Array ``y = fcn(x + delta)``.
  468. res_var : float, optional
  469. Residual variance.
  470. sum_square : float, optional
  471. Sum of squares error.
  472. sum_square_delta : float, optional
  473. Sum of squares of delta error.
  474. sum_square_eps : float, optional
  475. Sum of squares of eps error.
  476. inv_condnum : float, optional
  477. Inverse condition number (cf. ODRPACK UG p. 77).
  478. rel_error : float, optional
  479. Relative error in function values computed within fcn.
  480. work : ndarray, optional
  481. Final work array.
  482. work_ind : dict, optional
  483. Indices into work for drawing out values (cf. ODRPACK UG p. 83).
  484. info : int, optional
  485. Reason for returning, as output by ODRPACK (cf. ODRPACK UG p. 38).
  486. stopreason : list of str, optional
  487. `info` interpreted into English.
  488. Notes
  489. -----
  490. Takes one argument for initialization, the return value from the
  491. function `~scipy.odr.odr`. The attributes listed as "optional" above are
  492. only present if `~scipy.odr.odr` was run with ``full_output=1``.
  493. """
  494. def __init__(self, output):
  495. self.beta = output[0]
  496. self.sd_beta = output[1]
  497. self.cov_beta = output[2]
  498. if len(output) == 4:
  499. # full output
  500. self.__dict__.update(output[3])
  501. self.stopreason = _report_error(self.info)
  502. def pprint(self):
  503. """ Pretty-print important results.
  504. """
  505. print('Beta:', self.beta)
  506. print('Beta Std Error:', self.sd_beta)
  507. print('Beta Covariance:', self.cov_beta)
  508. if hasattr(self, 'info'):
  509. print('Residual Variance:',self.res_var)
  510. print('Inverse Condition #:', self.inv_condnum)
  511. print('Reason(s) for Halting:')
  512. for r in self.stopreason:
  513. print(' %s' % r)
  514. class ODR:
  515. """
  516. The ODR class gathers all information and coordinates the running of the
  517. main fitting routine.
  518. Members of instances of the ODR class have the same names as the arguments
  519. to the initialization routine.
  520. Parameters
  521. ----------
  522. data : Data class instance
  523. instance of the Data class
  524. model : Model class instance
  525. instance of the Model class
  526. Other Parameters
  527. ----------------
  528. beta0 : array_like of rank-1
  529. a rank-1 sequence of initial parameter values. Optional if
  530. model provides an "estimate" function to estimate these values.
  531. delta0 : array_like of floats of rank-1, optional
  532. a (double-precision) float array to hold the initial values of
  533. the errors in the input variables. Must be same shape as data.x
  534. ifixb : array_like of ints of rank-1, optional
  535. sequence of integers with the same length as beta0 that determines
  536. which parameters are held fixed. A value of 0 fixes the parameter,
  537. a value > 0 makes the parameter free.
  538. ifixx : array_like of ints with same shape as data.x, optional
  539. an array of integers with the same shape as data.x that determines
  540. which input observations are treated as fixed. One can use a sequence
  541. of length m (the dimensionality of the input observations) to fix some
  542. dimensions for all observations. A value of 0 fixes the observation,
  543. a value > 0 makes it free.
  544. job : int, optional
  545. an integer telling ODRPACK what tasks to perform. See p. 31 of the
  546. ODRPACK User's Guide if you absolutely must set the value here. Use the
  547. method set_job post-initialization for a more readable interface.
  548. iprint : int, optional
  549. an integer telling ODRPACK what to print. See pp. 33-34 of the
  550. ODRPACK User's Guide if you absolutely must set the value here. Use the
  551. method set_iprint post-initialization for a more readable interface.
  552. errfile : str, optional
  553. string with the filename to print ODRPACK errors to. If the file already
  554. exists, an error will be thrown. The `overwrite` argument can be used to
  555. prevent this. *Do Not Open This File Yourself!*
  556. rptfile : str, optional
  557. string with the filename to print ODRPACK summaries to. If the file
  558. already exists, an error will be thrown. The `overwrite` argument can be
  559. used to prevent this. *Do Not Open This File Yourself!*
  560. ndigit : int, optional
  561. integer specifying the number of reliable digits in the computation
  562. of the function.
  563. taufac : float, optional
  564. float specifying the initial trust region. The default value is 1.
  565. The initial trust region is equal to taufac times the length of the
  566. first computed Gauss-Newton step. taufac must be less than 1.
  567. sstol : float, optional
  568. float specifying the tolerance for convergence based on the relative
  569. change in the sum-of-squares. The default value is eps**(1/2) where eps
  570. is the smallest value such that 1 + eps > 1 for double precision
  571. computation on the machine. sstol must be less than 1.
  572. partol : float, optional
  573. float specifying the tolerance for convergence based on the relative
  574. change in the estimated parameters. The default value is eps**(2/3) for
  575. explicit models and ``eps**(1/3)`` for implicit models. partol must be less
  576. than 1.
  577. maxit : int, optional
  578. integer specifying the maximum number of iterations to perform. For
  579. first runs, maxit is the total number of iterations performed and
  580. defaults to 50. For restarts, maxit is the number of additional
  581. iterations to perform and defaults to 10.
  582. stpb : array_like, optional
  583. sequence (``len(stpb) == len(beta0)``) of relative step sizes to compute
  584. finite difference derivatives wrt the parameters.
  585. stpd : optional
  586. array (``stpd.shape == data.x.shape`` or ``stpd.shape == (m,)``) of relative
  587. step sizes to compute finite difference derivatives wrt the input
  588. variable errors. If stpd is a rank-1 array with length m (the
  589. dimensionality of the input variable), then the values are broadcast to
  590. all observations.
  591. sclb : array_like, optional
  592. sequence (``len(stpb) == len(beta0)``) of scaling factors for the
  593. parameters. The purpose of these scaling factors are to scale all of
  594. the parameters to around unity. Normally appropriate scaling factors
  595. are computed if this argument is not specified. Specify them yourself
  596. if the automatic procedure goes awry.
  597. scld : array_like, optional
  598. array (scld.shape == data.x.shape or scld.shape == (m,)) of scaling
  599. factors for the *errors* in the input variables. Again, these factors
  600. are automatically computed if you do not provide them. If scld.shape ==
  601. (m,), then the scaling factors are broadcast to all observations.
  602. work : ndarray, optional
  603. array to hold the double-valued working data for ODRPACK. When
  604. restarting, takes the value of self.output.work.
  605. iwork : ndarray, optional
  606. array to hold the integer-valued working data for ODRPACK. When
  607. restarting, takes the value of self.output.iwork.
  608. overwrite : bool, optional
  609. If it is True, output files defined by `errfile` and `rptfile` are
  610. overwritten. The default is False.
  611. Attributes
  612. ----------
  613. data : Data
  614. The data for this fit
  615. model : Model
  616. The model used in fit
  617. output : Output
  618. An instance if the Output class containing all of the returned
  619. data from an invocation of ODR.run() or ODR.restart()
  620. """
  621. def __init__(self, data, model, beta0=None, delta0=None, ifixb=None,
  622. ifixx=None, job=None, iprint=None, errfile=None, rptfile=None,
  623. ndigit=None, taufac=None, sstol=None, partol=None, maxit=None,
  624. stpb=None, stpd=None, sclb=None, scld=None, work=None, iwork=None,
  625. overwrite=False):
  626. self.data = data
  627. self.model = model
  628. if beta0 is None:
  629. if self.model.estimate is not None:
  630. self.beta0 = _conv(self.model.estimate(self.data))
  631. else:
  632. raise ValueError(
  633. "must specify beta0 or provide an estimater with the model"
  634. )
  635. else:
  636. self.beta0 = _conv(beta0)
  637. if ifixx is None and data.fix is not None:
  638. ifixx = data.fix
  639. if overwrite:
  640. # remove output files for overwriting.
  641. if rptfile is not None and os.path.exists(rptfile):
  642. os.remove(rptfile)
  643. if errfile is not None and os.path.exists(errfile):
  644. os.remove(errfile)
  645. self.delta0 = _conv(delta0)
  646. # These really are 32-bit integers in FORTRAN (gfortran), even on 64-bit
  647. # platforms.
  648. # XXX: some other FORTRAN compilers may not agree.
  649. self.ifixx = _conv(ifixx, dtype=numpy.int32)
  650. self.ifixb = _conv(ifixb, dtype=numpy.int32)
  651. self.job = job
  652. self.iprint = iprint
  653. self.errfile = errfile
  654. self.rptfile = rptfile
  655. self.ndigit = ndigit
  656. self.taufac = taufac
  657. self.sstol = sstol
  658. self.partol = partol
  659. self.maxit = maxit
  660. self.stpb = _conv(stpb)
  661. self.stpd = _conv(stpd)
  662. self.sclb = _conv(sclb)
  663. self.scld = _conv(scld)
  664. self.work = _conv(work)
  665. self.iwork = _conv(iwork)
  666. self.output = None
  667. self._check()
  668. def _check(self):
  669. """ Check the inputs for consistency, but don't bother checking things
  670. that the builtin function odr will check.
  671. """
  672. x_s = list(self.data.x.shape)
  673. if isinstance(self.data.y, numpy.ndarray):
  674. y_s = list(self.data.y.shape)
  675. if self.model.implicit:
  676. raise OdrError("an implicit model cannot use response data")
  677. else:
  678. # implicit model with q == self.data.y
  679. y_s = [self.data.y, x_s[-1]]
  680. if not self.model.implicit:
  681. raise OdrError("an explicit model needs response data")
  682. self.set_job(fit_type=1)
  683. if x_s[-1] != y_s[-1]:
  684. raise OdrError("number of observations do not match")
  685. n = x_s[-1]
  686. if len(x_s) == 2:
  687. m = x_s[0]
  688. else:
  689. m = 1
  690. if len(y_s) == 2:
  691. q = y_s[0]
  692. else:
  693. q = 1
  694. p = len(self.beta0)
  695. # permissible output array shapes
  696. fcn_perms = [(q, n)]
  697. fjacd_perms = [(q, m, n)]
  698. fjacb_perms = [(q, p, n)]
  699. if q == 1:
  700. fcn_perms.append((n,))
  701. fjacd_perms.append((m, n))
  702. fjacb_perms.append((p, n))
  703. if m == 1:
  704. fjacd_perms.append((q, n))
  705. if p == 1:
  706. fjacb_perms.append((q, n))
  707. if m == q == 1:
  708. fjacd_perms.append((n,))
  709. if p == q == 1:
  710. fjacb_perms.append((n,))
  711. # try evaluating the supplied functions to make sure they provide
  712. # sensible outputs
  713. arglist = (self.beta0, self.data.x)
  714. if self.model.extra_args is not None:
  715. arglist = arglist + self.model.extra_args
  716. res = self.model.fcn(*arglist)
  717. if res.shape not in fcn_perms:
  718. print(res.shape)
  719. print(fcn_perms)
  720. raise OdrError("fcn does not output %s-shaped array" % y_s)
  721. if self.model.fjacd is not None:
  722. res = self.model.fjacd(*arglist)
  723. if res.shape not in fjacd_perms:
  724. raise OdrError(
  725. "fjacd does not output %s-shaped array" % repr((q, m, n)))
  726. if self.model.fjacb is not None:
  727. res = self.model.fjacb(*arglist)
  728. if res.shape not in fjacb_perms:
  729. raise OdrError(
  730. "fjacb does not output %s-shaped array" % repr((q, p, n)))
  731. # check shape of delta0
  732. if self.delta0 is not None and self.delta0.shape != self.data.x.shape:
  733. raise OdrError(
  734. "delta0 is not a %s-shaped array" % repr(self.data.x.shape))
  735. if self.data.x.size == 0:
  736. warn(("Empty data detected for ODR instance. "
  737. "Do not expect any fitting to occur"),
  738. OdrWarning)
  739. def _gen_work(self):
  740. """ Generate a suitable work array if one does not already exist.
  741. """
  742. n = self.data.x.shape[-1]
  743. p = self.beta0.shape[0]
  744. if len(self.data.x.shape) == 2:
  745. m = self.data.x.shape[0]
  746. else:
  747. m = 1
  748. if self.model.implicit:
  749. q = self.data.y
  750. elif len(self.data.y.shape) == 2:
  751. q = self.data.y.shape[0]
  752. else:
  753. q = 1
  754. if self.data.we is None:
  755. ldwe = ld2we = 1
  756. elif len(self.data.we.shape) == 3:
  757. ld2we, ldwe = self.data.we.shape[1:]
  758. else:
  759. # Okay, this isn't precisely right, but for this calculation,
  760. # it's fine
  761. ldwe = 1
  762. ld2we = self.data.we.shape[1]
  763. if self.job % 10 < 2:
  764. # ODR not OLS
  765. lwork = (18 + 11*p + p*p + m + m*m + 4*n*q + 6*n*m + 2*n*q*p +
  766. 2*n*q*m + q*q + 5*q + q*(p+m) + ldwe*ld2we*q)
  767. else:
  768. # OLS not ODR
  769. lwork = (18 + 11*p + p*p + m + m*m + 4*n*q + 2*n*m + 2*n*q*p +
  770. 5*q + q*(p+m) + ldwe*ld2we*q)
  771. if isinstance(self.work, numpy.ndarray) and self.work.shape == (lwork,)\
  772. and self.work.dtype.str.endswith('f8'):
  773. # the existing array is fine
  774. return
  775. else:
  776. self.work = numpy.zeros((lwork,), float)
  777. def set_job(self, fit_type=None, deriv=None, var_calc=None,
  778. del_init=None, restart=None):
  779. """
  780. Sets the "job" parameter is a hopefully comprehensible way.
  781. If an argument is not specified, then the value is left as is. The
  782. default value from class initialization is for all of these options set
  783. to 0.
  784. Parameters
  785. ----------
  786. fit_type : {0, 1, 2} int
  787. 0 -> explicit ODR
  788. 1 -> implicit ODR
  789. 2 -> ordinary least-squares
  790. deriv : {0, 1, 2, 3} int
  791. 0 -> forward finite differences
  792. 1 -> central finite differences
  793. 2 -> user-supplied derivatives (Jacobians) with results
  794. checked by ODRPACK
  795. 3 -> user-supplied derivatives, no checking
  796. var_calc : {0, 1, 2} int
  797. 0 -> calculate asymptotic covariance matrix and fit
  798. parameter uncertainties (V_B, s_B) using derivatives
  799. recomputed at the final solution
  800. 1 -> calculate V_B and s_B using derivatives from last iteration
  801. 2 -> do not calculate V_B and s_B
  802. del_init : {0, 1} int
  803. 0 -> initial input variable offsets set to 0
  804. 1 -> initial offsets provided by user in variable "work"
  805. restart : {0, 1} int
  806. 0 -> fit is not a restart
  807. 1 -> fit is a restart
  808. Notes
  809. -----
  810. The permissible values are different from those given on pg. 31 of the
  811. ODRPACK User's Guide only in that one cannot specify numbers greater than
  812. the last value for each variable.
  813. If one does not supply functions to compute the Jacobians, the fitting
  814. procedure will change deriv to 0, finite differences, as a default. To
  815. initialize the input variable offsets by yourself, set del_init to 1 and
  816. put the offsets into the "work" variable correctly.
  817. """
  818. if self.job is None:
  819. job_l = [0, 0, 0, 0, 0]
  820. else:
  821. job_l = [self.job // 10000 % 10,
  822. self.job // 1000 % 10,
  823. self.job // 100 % 10,
  824. self.job // 10 % 10,
  825. self.job % 10]
  826. if fit_type in (0, 1, 2):
  827. job_l[4] = fit_type
  828. if deriv in (0, 1, 2, 3):
  829. job_l[3] = deriv
  830. if var_calc in (0, 1, 2):
  831. job_l[2] = var_calc
  832. if del_init in (0, 1):
  833. job_l[1] = del_init
  834. if restart in (0, 1):
  835. job_l[0] = restart
  836. self.job = (job_l[0]*10000 + job_l[1]*1000 +
  837. job_l[2]*100 + job_l[3]*10 + job_l[4])
  838. def set_iprint(self, init=None, so_init=None,
  839. iter=None, so_iter=None, iter_step=None, final=None, so_final=None):
  840. """ Set the iprint parameter for the printing of computation reports.
  841. If any of the arguments are specified here, then they are set in the
  842. iprint member. If iprint is not set manually or with this method, then
  843. ODRPACK defaults to no printing. If no filename is specified with the
  844. member rptfile, then ODRPACK prints to stdout. One can tell ODRPACK to
  845. print to stdout in addition to the specified filename by setting the
  846. so_* arguments to this function, but one cannot specify to print to
  847. stdout but not a file since one can do that by not specifying a rptfile
  848. filename.
  849. There are three reports: initialization, iteration, and final reports.
  850. They are represented by the arguments init, iter, and final
  851. respectively. The permissible values are 0, 1, and 2 representing "no
  852. report", "short report", and "long report" respectively.
  853. The argument iter_step (0 <= iter_step <= 9) specifies how often to make
  854. the iteration report; the report will be made for every iter_step'th
  855. iteration starting with iteration one. If iter_step == 0, then no
  856. iteration report is made, regardless of the other arguments.
  857. If the rptfile is None, then any so_* arguments supplied will raise an
  858. exception.
  859. """
  860. if self.iprint is None:
  861. self.iprint = 0
  862. ip = [self.iprint // 1000 % 10,
  863. self.iprint // 100 % 10,
  864. self.iprint // 10 % 10,
  865. self.iprint % 10]
  866. # make a list to convert iprint digits to/from argument inputs
  867. # rptfile, stdout
  868. ip2arg = [[0, 0], # none, none
  869. [1, 0], # short, none
  870. [2, 0], # long, none
  871. [1, 1], # short, short
  872. [2, 1], # long, short
  873. [1, 2], # short, long
  874. [2, 2]] # long, long
  875. if (self.rptfile is None and
  876. (so_init is not None or
  877. so_iter is not None or
  878. so_final is not None)):
  879. raise OdrError(
  880. "no rptfile specified, cannot output to stdout twice")
  881. iprint_l = ip2arg[ip[0]] + ip2arg[ip[1]] + ip2arg[ip[3]]
  882. if init is not None:
  883. iprint_l[0] = init
  884. if so_init is not None:
  885. iprint_l[1] = so_init
  886. if iter is not None:
  887. iprint_l[2] = iter
  888. if so_iter is not None:
  889. iprint_l[3] = so_iter
  890. if final is not None:
  891. iprint_l[4] = final
  892. if so_final is not None:
  893. iprint_l[5] = so_final
  894. if iter_step in range(10):
  895. # 0..9
  896. ip[2] = iter_step
  897. ip[0] = ip2arg.index(iprint_l[0:2])
  898. ip[1] = ip2arg.index(iprint_l[2:4])
  899. ip[3] = ip2arg.index(iprint_l[4:6])
  900. self.iprint = ip[0]*1000 + ip[1]*100 + ip[2]*10 + ip[3]
  901. def run(self):
  902. """ Run the fitting routine with all of the information given and with ``full_output=1``.
  903. Returns
  904. -------
  905. output : Output instance
  906. This object is also assigned to the attribute .output .
  907. """
  908. args = (self.model.fcn, self.beta0, self.data.y, self.data.x)
  909. kwds = {'full_output': 1}
  910. kwd_l = ['ifixx', 'ifixb', 'job', 'iprint', 'errfile', 'rptfile',
  911. 'ndigit', 'taufac', 'sstol', 'partol', 'maxit', 'stpb',
  912. 'stpd', 'sclb', 'scld', 'work', 'iwork']
  913. if self.delta0 is not None and (self.job // 10000) % 10 == 0:
  914. # delta0 provided and fit is not a restart
  915. self._gen_work()
  916. d0 = numpy.ravel(self.delta0)
  917. self.work[:len(d0)] = d0
  918. # set the kwds from other objects explicitly
  919. if self.model.fjacb is not None:
  920. kwds['fjacb'] = self.model.fjacb
  921. if self.model.fjacd is not None:
  922. kwds['fjacd'] = self.model.fjacd
  923. if self.data.we is not None:
  924. kwds['we'] = self.data.we
  925. if self.data.wd is not None:
  926. kwds['wd'] = self.data.wd
  927. if self.model.extra_args is not None:
  928. kwds['extra_args'] = self.model.extra_args
  929. # implicitly set kwds from self's members
  930. for attr in kwd_l:
  931. obj = getattr(self, attr)
  932. if obj is not None:
  933. kwds[attr] = obj
  934. self.output = Output(odr(*args, **kwds))
  935. return self.output
  936. def restart(self, iter=None):
  937. """ Restarts the run with iter more iterations.
  938. Parameters
  939. ----------
  940. iter : int, optional
  941. ODRPACK's default for the number of new iterations is 10.
  942. Returns
  943. -------
  944. output : Output instance
  945. This object is also assigned to the attribute .output .
  946. """
  947. if self.output is None:
  948. raise OdrError("cannot restart: run() has not been called before")
  949. self.set_job(restart=1)
  950. self.work = self.output.work
  951. self.iwork = self.output.iwork
  952. self.maxit = iter
  953. return self.run()