Package omero :: Package util :: Module mpfit
[hide private]
[frames] | no frames]

Source Code for Module omero.util.mpfit

   1  """ 
   2  Perform Levenberg-Marquardt least-squares minimization, based on MINPACK-1. 
   3   
   4                                                                     AUTHORS 
   5    The original version of this software, called LMFIT, was written in FORTRAN 
   6    as part of the MINPACK-1 package by XXX. 
   7   
   8    Craig Markwardt converted the FORTRAN code to IDL.  The information for the 
   9    IDL version is: 
  10           Craig B. Markwardt, NASA/GSFC Code 662, Greenbelt, MD 20770 
  11           craigm@lheamail.gsfc.nasa.gov 
  12           UPDATED VERSIONs can be found on my WEB PAGE: 
  13                  http://cow.physics.wisc.edu/~craigm/idl/idl.html 
  14           
  15    Mark Rivers created this Python version from Craig's IDL version. 
  16          Mark Rivers, University of Chicago 
  17          Building 434A, Argonne National Laboratory 
  18          9700 South Cass Avenue, Argonne, IL 60439 
  19          rivers@cars.uchicago.edu 
  20          Updated versions can be found at http://cars.uchicago.edu/software 
  21    
  22   Sergey Koposov converted the Mark's Python version from Numeric to numpy 
  23          Sergey Koposov, Max Planck Institute for Astronomy 
  24          Heidelberg, Germany, D-69117 
  25          koposov@mpia.de 
  26          Updated versions can be found at http://code.google.com/p/astrolibpy/source/browse/trunk/ 
  27   
  28                                                                   DESCRIPTION 
  29   
  30   MPFIT uses the Levenberg-Marquardt technique to solve the 
  31   least-squares problem.  In its typical use, MPFIT will be used to 
  32   fit a user-supplied function (the "model") to user-supplied data 
  33   points (the "data") by adjusting a set of parameters.  MPFIT is 
  34   based upon MINPACK-1 (LMDIF.F) by More' and collaborators. 
  35   
  36   For example, a researcher may think that a set of observed data 
  37   points is best modelled with a Gaussian curve.  A Gaussian curve is 
  38   parameterized by its mean, standard deviation and normalization. 
  39   MPFIT will, within certain constraints, find the set of parameters 
  40   which best fits the data.  The fit is "best" in the least-squares 
  41   sense; that is, the sum of the weighted squared differences between 
  42   the model and data is minimized. 
  43   
  44   The Levenberg-Marquardt technique is a particular strategy for 
  45   iteratively searching for the best fit.  This particular 
  46   implementation is drawn from MINPACK-1 (see NETLIB), and is much faster 
  47   and more accurate than the version provided in the Scientific Python package 
  48   in Scientific.Functions.LeastSquares. 
  49   This version allows upper and lower bounding constraints to be placed on each 
  50   parameter, or the parameter can be held fixed. 
  51   
  52   The user-supplied Python function should return an array of weighted 
  53   deviations between model and data.  In a typical scientific problem 
  54   the residuals should be weighted so that each deviate has a 
  55   gaussian sigma of 1.0.  If X represents values of the independent 
  56   variable, Y represents a measurement for each value of X, and ERR 
  57   represents the error in the measurements, then the deviates could 
  58   be calculated as follows: 
  59   
  60     DEVIATES = (Y - F(X)) / ERR 
  61   
  62   where F is the analytical function representing the model.  You are 
  63   recommended to use the convenience functions MPFITFUN and 
  64   MPFITEXPR, which are driver functions that calculate the deviates 
  65   for you.  If ERR are the 1-sigma uncertainties in Y, then 
  66   
  67     TOTAL( DEVIATES^2 ) 
  68   
  69   will be the total chi-squared value.  MPFIT will minimize the 
  70   chi-square value.  The values of X, Y and ERR are passed through 
  71   MPFIT to the user-supplied function via the FUNCTKW keyword. 
  72   
  73   Simple constraints can be placed on parameter values by using the 
  74   PARINFO keyword to MPFIT.  See below for a description of this 
  75   keyword. 
  76   
  77   MPFIT does not perform more general optimization tasks.  See TNMIN 
  78   instead.  MPFIT is customized, based on MINPACK-1, to the 
  79   least-squares minimization problem. 
  80   
  81   
  82                                                             USER FUNCTION 
  83   
  84   The user must define a function which returns the appropriate 
  85   values as specified above.  The function should return the weighted 
  86   deviations between the model and the data.  It should also return a status 
  87   flag and an optional partial derivative array.  For applications which 
  88   use finite-difference derivatives -- the default -- the user 
  89   function should be declared in the following way: 
  90   
  91     def myfunct(p, fjac=None, x=None, y=None, err=None) 
  92          # Parameter values are passed in "p" 
  93          # If fjac==None then partial derivatives should not be 
  94          # computed.  It will always be None if MPFIT is called with default 
  95          # flag. 
  96          model = F(x, p) 
  97          # Non-negative status value means MPFIT should continue, negative means 
  98          # stop the calculation. 
  99          status = 0 
 100          return([status, (y-model)/err] 
 101   
 102   See below for applications with analytical derivatives. 
 103   
 104   The keyword parameters X, Y, and ERR in the example above are 
 105   suggestive but not required.  Any parameters can be passed to 
 106   MYFUNCT by using the functkw keyword to MPFIT.  Use MPFITFUN and 
 107   MPFITEXPR if you need ideas on how to do that.  The function *must* 
 108   accept a parameter list, P. 
 109   
 110   In general there are no restrictions on the number of dimensions in 
 111   X, Y or ERR.  However the deviates *must* be returned in a 
 112   one-dimensional Numeric array of type Float. 
 113   
 114   User functions may also indicate a fatal error condition using the 
 115   status return described above. If status is set to a number between 
 116   -15 and -1 then MPFIT will stop the calculation and return to the caller. 
 117   
 118   
 119                                                          ANALYTIC DERIVATIVES 
 120   
 121   In the search for the best-fit solution, MPFIT by default 
 122   calculates derivatives numerically via a finite difference 
 123   approximation.  The user-supplied function need not calculate the 
 124   derivatives explicitly.  However, if you desire to compute them 
 125   analytically, then the AUTODERIVATIVE=0 keyword must be passed to MPFIT. 
 126   As a practical matter, it is often sufficient and even faster to allow 
 127   MPFIT to calculate the derivatives numerically, and so 
 128   AUTODERIVATIVE=0 is not necessary. 
 129   
 130   If AUTODERIVATIVE=0 is used then the user function must check the parameter 
 131   FJAC, and if FJAC!=None then return the partial derivative array in the 
 132   return list. 
 133     def myfunct(p, fjac=None, x=None, y=None, err=None) 
 134          # Parameter values are passed in "p" 
 135          # If FJAC!=None then partial derivatives must be comptuer. 
 136          # FJAC contains an array of len(p), where each entry 
 137          # is 1 if that parameter is free and 0 if it is fixed. 
 138          model = F(x, p) 
 139          Non-negative status value means MPFIT should continue, negative means 
 140          # stop the calculation. 
 141          status = 0 
 142          if (dojac): 
 143             pderiv = zeros([len(x), len(p)], Float) 
 144             for j in range(len(p)): 
 145                   pderiv[:,j] = FGRAD(x, p, j) 
 146          else: 
 147             pderiv = None 
 148          return([status, (y-model)/err, pderiv] 
 149   
 150   where FGRAD(x, p, i) is a user function which must compute the 
 151   derivative of the model with respect to parameter P[i] at X.  When 
 152   finite differencing is used for computing derivatives (ie, when 
 153   AUTODERIVATIVE=1), or when MPFIT needs only the errors but not the 
 154   derivatives the parameter FJAC=None. 
 155   
 156   Derivatives should be returned in the PDERIV array. PDERIV should be an m x 
 157   n array, where m is the number of data points and n is the number 
 158   of parameters.  dp[i,j] is the derivative at the ith point with 
 159   respect to the jth parameter. 
 160   
 161   The derivatives with respect to fixed parameters are ignored; zero 
 162   is an appropriate value to insert for those derivatives.  Upon 
 163   input to the user function, FJAC is set to a vector with the same 
 164   length as P, with a value of 1 for a parameter which is free, and a 
 165   value of zero for a parameter which is fixed (and hence no 
 166   derivative needs to be calculated). 
 167   
 168   If the data is higher than one dimensional, then the *last* 
 169   dimension should be the parameter dimension.  Example: fitting a 
 170   50x50 image, "dp" should be 50x50xNPAR. 
 171   
 172   
 173                     CONSTRAINING PARAMETER VALUES WITH THE PARINFO KEYWORD 
 174   
 175   The behavior of MPFIT can be modified with respect to each 
 176   parameter to be fitted.  A parameter value can be fixed; simple 
 177   boundary constraints can be imposed; limitations on the parameter 
 178   changes can be imposed; properties of the automatic derivative can 
 179   be modified; and parameters can be tied to one another. 
 180   
 181   These properties are governed by the PARINFO structure, which is 
 182   passed as a keyword parameter to MPFIT. 
 183   
 184   PARINFO should be a list of dictionaries, one list entry for each parameter. 
 185   Each parameter is associated with one element of the array, in 
 186   numerical order.  The dictionary can have the following keys 
 187   (none are required, keys are case insensitive): 
 188   
 189          'value' - the starting parameter value (but see the START_PARAMS 
 190                           parameter for more information). 
 191   
 192          'fixed' - a boolean value, whether the parameter is to be held 
 193                           fixed or not.  Fixed parameters are not varied by 
 194                           MPFIT, but are passed on to MYFUNCT for evaluation. 
 195   
 196          'limited' - a two-element boolean array.  If the first/second 
 197                             element is set, then the parameter is bounded on the 
 198                             lower/upper side.  A parameter can be bounded on both 
 199                             sides.  Both LIMITED and LIMITS must be given 
 200                             together. 
 201   
 202          'limits' - a two-element float array.  Gives the 
 203                            parameter limits on the lower and upper sides, 
 204                            respectively.  Zero, one or two of these values can be 
 205                            set, depending on the values of LIMITED.  Both LIMITED 
 206                            and LIMITS must be given together. 
 207   
 208          'parname' - a string, giving the name of the parameter.  The 
 209                             fitting code of MPFIT does not use this tag in any 
 210                             way.  However, the default iterfunct will print the 
 211                             parameter name if available. 
 212   
 213          'step' - the step size to be used in calculating the numerical 
 214                          derivatives.  If set to zero, then the step size is 
 215                          computed automatically.  Ignored when AUTODERIVATIVE=0. 
 216   
 217          'mpside' - the sidedness of the finite difference when computing 
 218                            numerical derivatives.  This field can take four 
 219                            values: 
 220   
 221                                   0 - one-sided derivative computed automatically 
 222                                   1 - one-sided derivative (f(x+h) - f(x)  )/h 
 223                                  -1 - one-sided derivative (f(x)   - f(x-h))/h 
 224                                   2 - two-sided derivative (f(x+h) - f(x-h))/(2*h) 
 225   
 226                           Where H is the STEP parameter described above.  The 
 227                           "automatic" one-sided derivative method will chose a 
 228                           direction for the finite difference which does not 
 229                           violate any constraints.  The other methods do not 
 230                           perform this check.  The two-sided method is in 
 231                           principle more precise, but requires twice as many 
 232                           function evaluations.  Default: 0. 
 233   
 234          'mpmaxstep' - the maximum change to be made in the parameter 
 235                                   value.  During the fitting process, the parameter 
 236                                   will never be changed by more than this value in 
 237                                   one iteration. 
 238   
 239                                   A value of 0 indicates no maximum.  Default: 0. 
 240   
 241          'tied' - a string expression which "ties" the parameter to other 
 242                          free or fixed parameters.  Any expression involving 
 243                          constants and the parameter array P are permitted. 
 244                          Example: if parameter 2 is always to be twice parameter 
 245                          1 then use the following: parinfo(2).tied = '2 * p(1)'. 
 246                          Since they are totally constrained, tied parameters are 
 247                          considered to be fixed; no errors are computed for them. 
 248                          [ NOTE: the PARNAME can't be used in expressions. ] 
 249   
 250          'mpprint' - if set to 1, then the default iterfunct will print the 
 251                             parameter value.  If set to 0, the parameter value 
 252                             will not be printed.  This tag can be used to 
 253                             selectively print only a few parameter values out of 
 254                             many.  Default: 1 (all parameters printed) 
 255   
 256   
 257   Future modifications to the PARINFO structure, if any, will involve 
 258   adding dictionary tags beginning with the two letters "MP". 
 259   Therefore programmers are urged to avoid using tags starting with 
 260   the same letters; otherwise they are free to include their own 
 261   fields within the PARINFO structure, and they will be ignored. 
 262   
 263   PARINFO Example: 
 264   parinfo = [{'value':0., 'fixed':0, 'limited':[0,0], 'limits':[0.,0.]}]*5 
 265   parinfo[0]['fixed'] = 1 
 266   parinfo[4]['limited'][0] = 1 
 267   parinfo[4]['limits'][0]  = 50. 
 268   values = [5.7, 2.2, 500., 1.5, 2000.] 
 269   for i in range(5): parinfo[i]['value']=values[i] 
 270   
 271   A total of 5 parameters, with starting values of 5.7, 
 272   2.2, 500, 1.5, and 2000 are given.  The first parameter 
 273   is fixed at a value of 5.7, and the last parameter is 
 274   constrained to be above 50. 
 275   
 276   
 277                                                                     EXAMPLE 
 278   
 279     import mpfit 
 280     import numpy.oldnumeric as Numeric 
 281     x = arange(100, float) 
 282     p0 = [5.7, 2.2, 500., 1.5, 2000.] 
 283     y = ( p[0] + p[1]*[x] + p[2]*[x**2] + p[3]*sqrt(x) + 
 284                   p[4]*log(x)) 
 285     fa = {'x':x, 'y':y, 'err':err} 
 286     m = mpfit('myfunct', p0, functkw=fa) 
 287     print 'status = ', m.status 
 288     if (m.status <= 0): print 'error message = ', m.errmsg 
 289     print 'parameters = ', m.params 
 290   
 291     Minimizes sum of squares of MYFUNCT.  MYFUNCT is called with the X, 
 292     Y, and ERR keyword parameters that are given by FUNCTKW.  The 
 293     results can be obtained from the returned object m. 
 294   
 295   
 296                                                          THEORY OF OPERATION 
 297   
 298     There are many specific strategies for function minimization.  One 
 299     very popular technique is to use function gradient information to 
 300     realize the local structure of the function.  Near a local minimum 
 301     the function value can be taylor expanded about x0 as follows: 
 302   
 303            f(x) = f(x0) + f'(x0) . (x-x0) + (1/2) (x-x0) . f''(x0) . (x-x0) 
 304                           -----   ---------------   -------------------------------  (1) 
 305           Order  0th               1st                                     2nd 
 306   
 307     Here f'(x) is the gradient vector of f at x, and f''(x) is the 
 308     Hessian matrix of second derivatives of f at x.  The vector x is 
 309     the set of function parameters, not the measured data vector.  One 
 310     can find the minimum of f, f(xm) using Newton's method, and 
 311     arrives at the following linear equation: 
 312   
 313            f''(x0) . (xm-x0) = - f'(x0)                                                  (2) 
 314   
 315     If an inverse can be found for f''(x0) then one can solve for 
 316     (xm-x0), the step vector from the current position x0 to the new 
 317     projected minimum.  Here the problem has been linearized (ie, the 
 318     gradient information is known to first order).  f''(x0) is 
 319     symmetric n x n matrix, and should be positive definite. 
 320   
 321     The Levenberg - Marquardt technique is a variation on this theme. 
 322     It adds an additional diagonal term to the equation which may aid the 
 323     convergence properties: 
 324   
 325            (f''(x0) + nu I) . (xm-x0) = -f'(x0)                            (2a) 
 326   
 327     where I is the identity matrix.  When nu is large, the overall 
 328     matrix is diagonally dominant, and the iterations follow steepest 
 329     descent.  When nu is small, the iterations are quadratically 
 330     convergent. 
 331   
 332     In principle, if f''(x0) and f'(x0) are known then xm-x0 can be 
 333     determined.  However the Hessian matrix is often difficult or 
 334     impossible to compute.  The gradient f'(x0) may be easier to 
 335     compute, if even by finite difference techniques.  So-called 
 336     quasi-Newton techniques attempt to successively estimate f''(x0) 
 337     by building up gradient information as the iterations proceed. 
 338   
 339     In the least squares problem there are further simplifications 
 340     which assist in solving eqn (2).  The function to be minimized is 
 341     a sum of squares: 
 342   
 343             f = Sum(hi^2)                                                                                 (3) 
 344   
 345     where hi is the ith residual out of m residuals as described 
 346     above.  This can be substituted back into eqn (2) after computing 
 347     the derivatives: 
 348   
 349             f'  = 2 Sum(hi  hi') 
 350             f'' = 2 Sum(hi' hj') + 2 Sum(hi hi'')                                (4) 
 351   
 352     If one assumes that the parameters are already close enough to a 
 353     minimum, then one typically finds that the second term in f'' is 
 354     negligible [or, in any case, is too difficult to compute].  Thus, 
 355     equation (2) can be solved, at least approximately, using only 
 356     gradient information. 
 357   
 358     In matrix notation, the combination of eqns (2) and (4) becomes: 
 359   
 360                  hT' . h' . dx = - hT' . h                                                 (5) 
 361   
 362     Where h is the residual vector (length m), hT is its transpose, h' 
 363     is the Jacobian matrix (dimensions n x m), and dx is (xm-x0).  The 
 364     user function supplies the residual vector h, and in some cases h' 
 365     when it is not found by finite differences (see MPFIT_FDJAC2, 
 366     which finds h and hT').  Even if dx is not the best absolute step 
 367     to take, it does provide a good estimate of the best *direction*, 
 368     so often a line minimization will occur along the dx vector 
 369     direction. 
 370   
 371     The method of solution employed by MINPACK is to form the Q . R 
 372     factorization of h', where Q is an orthogonal matrix such that QT . 
 373     Q = I, and R is upper right triangular.  Using h' = Q . R and the 
 374     ortogonality of Q, eqn (5) becomes 
 375   
 376                  (RT . QT) . (Q . R) . dx = - (RT . QT) . h 
 377                                           RT . R . dx = - RT . QT . h             (6) 
 378                                                    R . dx = - QT . h 
 379   
 380     where the last statement follows because R is upper triangular. 
 381     Here, R, QT and h are known so this is a matter of solving for dx. 
 382     The routine MPFIT_QRFAC provides the QR factorization of h, with 
 383     pivoting, and MPFIT_QRSOLV provides the solution for dx. 
 384   
 385   
 386                                                                   REFERENCES 
 387   
 388     MINPACK-1, Jorge More', available from netlib (www.netlib.org). 
 389     "Optimization Software Guide," Jorge More' and Stephen Wright, 
 390           SIAM, *Frontiers in Applied Mathematics*, Number 14. 
 391     More', Jorge J., "The Levenberg-Marquardt Algorithm: 
 392           Implementation and Theory," in *Numerical Analysis*, ed. Watson, 
 393           G. A., Lecture Notes in Mathematics 630, Springer-Verlag, 1977. 
 394   
 395   
 396                                                     MODIFICATION HISTORY 
 397   
 398     Translated from MINPACK-1 in FORTRAN, Apr-Jul 1998, CM 
 399   Copyright (C) 1997-2002, Craig Markwardt 
 400   This software is provided as is without any warranty whatsoever. 
 401   Permission to use, copy, modify, and distribute modified or 
 402   unmodified copies is granted, provided this copyright and disclaimer 
 403   are included unchanged. 
 404   
 405     Translated from MPFIT (Craig Markwardt's IDL package) to Python, 
 406     August, 2002.  Mark Rivers 
 407     Converted from Numeric to numpy (Sergey Koposov, July 2008) 
 408  """ 
 409   
 410  from numpy import * 
 411  import numpy 
 412  import types 
 413   
 414  #        Original FORTRAN documentation 
 415  #        ********** 
 416  # 
 417  #        subroutine lmdif 
 418  # 
 419  #        the purpose of lmdif is to minimize the sum of the squares of 
 420  #        m nonlinear functions in n variables by a modification of 
 421  #        the levenberg-marquardt algorithm. the user must provide a 
 422  #        subroutine which calculates the functions. the jacobian is 
 423  #        then calculated by a forward-difference approximation. 
 424  # 
 425  #        the subroutine statement is 
 426  # 
 427  #          subroutine lmdif(fcn,m,n,x,fvec,ftol,xtol,gtol,maxfev,epsfcn, 
 428  #                                               diag,mode,factor,nprint,info,nfev,fjac, 
 429  #                                               ldfjac,ipvt,qtf,wa1,wa2,wa3,wa4) 
 430  # 
 431  #        where 
 432  # 
 433  #          fcn is the name of the user-supplied subroutine which 
 434  #                calculates the functions. fcn must be declared 
 435  #                in an external statement in the user calling 
 436  #                program, and should be written as follows. 
 437  # 
 438  #                subroutine fcn(m,n,x,fvec,iflag) 
 439  #                integer m,n,iflag 
 440  #                double precision x(n),fvec(m) 
 441  #                ---------- 
 442  #                calculate the functions at x and 
 443  #                return this vector in fvec. 
 444  #                ---------- 
 445  #                return 
 446  #                end 
 447  # 
 448  #                the value of iflag should not be changed by fcn unless 
 449  #                the user wants to terminate execution of lmdif. 
 450  #                in this case set iflag to a negative integer. 
 451  # 
 452  #          m is a positive integer input variable set to the number 
 453  #                of functions. 
 454  # 
 455  #          n is a positive integer input variable set to the number 
 456  #                of variables. n must not exceed m. 
 457  # 
 458  #          x is an array of length n. on input x must contain 
 459  #                an initial estimate of the solution vector. on output x 
 460  #                contains the final estimate of the solution vector. 
 461  # 
 462  #          fvec is an output array of length m which contains 
 463  #                the functions evaluated at the output x. 
 464  # 
 465  #          ftol is a nonnegative input variable. termination 
 466  #                occurs when both the actual and predicted relative 
 467  #                reductions in the sum of squares are at most ftol. 
 468  #                therefore, ftol measures the relative error desired 
 469  #                in the sum of squares. 
 470  # 
 471  #          xtol is a nonnegative input variable. termination 
 472  #                occurs when the relative error between two consecutive 
 473  #                iterates is at most xtol. therefore, xtol measures the 
 474  #                relative error desired in the approximate solution. 
 475  # 
 476  #          gtol is a nonnegative input variable. termination 
 477  #                occurs when the cosine of the angle between fvec and 
 478  #                any column of the jacobian is at most gtol in absolute 
 479  #                value. therefore, gtol measures the orthogonality 
 480  #                desired between the function vector and the columns 
 481  #                of the jacobian. 
 482  # 
 483  #          maxfev is a positive integer input variable. termination 
 484  #                occurs when the number of calls to fcn is at least 
 485  #                maxfev by the end of an iteration. 
 486  # 
 487  #          epsfcn is an input variable used in determining a suitable 
 488  #                step length for the forward-difference approximation. this 
 489  #                approximation assumes that the relative errors in the 
 490  #                functions are of the order of epsfcn. if epsfcn is less 
 491  #                than the machine precision, it is assumed that the relative 
 492  #                errors in the functions are of the order of the machine 
 493  #                precision. 
 494  # 
 495  #          diag is an array of length n. if mode = 1 (see 
 496  #                below), diag is internally set. if mode = 2, diag 
 497  #                must contain positive entries that serve as 
 498  #                multiplicative scale factors for the variables. 
 499  # 
 500  #          mode is an integer input variable. if mode = 1, the 
 501  #                variables will be scaled internally. if mode = 2, 
 502  #                the scaling is specified by the input diag. other 
 503  #                values of mode are equivalent to mode = 1. 
 504  # 
 505  #          factor is a positive input variable used in determining the 
 506  #                initial step bound. this bound is set to the product of 
 507  #                factor and the euclidean norm of diag*x if nonzero, or else 
 508  #                to factor itself. in most cases factor should lie in the 
 509  #                interval (.1,100.). 100. is a generally recommended value. 
 510  # 
 511  #          nprint is an integer input variable that enables controlled 
 512  #                printing of iterates if it is positive. in this case, 
 513  #                fcn is called with iflag = 0 at the beginning of the first 
 514  #                iteration and every nprint iterations thereafter and 
 515  #                immediately prior to return, with x and fvec available 
 516  #                for printing. if nprint is not positive, no special calls 
 517  #                of fcn with iflag = 0 are made. 
 518  # 
 519  #          info is an integer output variable. if the user has 
 520  #                terminated execution, info is set to the (negative) 
 521  #                value of iflag. see description of fcn. otherwise, 
 522  #                info is set as follows. 
 523  # 
 524  #                info = 0  improper input parameters. 
 525  # 
 526  #                info = 1  both actual and predicted relative reductions 
 527  #                                  in the sum of squares are at most ftol. 
 528  # 
 529  #                info = 2  relative error between two consecutive iterates 
 530  #                                  is at most xtol. 
 531  # 
 532  #                info = 3  conditions for info = 1 and info = 2 both hold. 
 533  # 
 534  #                info = 4  the cosine of the angle between fvec and any 
 535  #                                  column of the jacobian is at most gtol in 
 536  #                                  absolute value. 
 537  # 
 538  #                info = 5  number of calls to fcn has reached or 
 539  #                                  exceeded maxfev. 
 540  # 
 541  #                info = 6  ftol is too small. no further reduction in 
 542  #                                  the sum of squares is possible. 
 543  # 
 544  #                info = 7  xtol is too small. no further improvement in 
 545  #                                  the approximate solution x is possible. 
 546  # 
 547  #                info = 8  gtol is too small. fvec is orthogonal to the 
 548  #                                  columns of the jacobian to machine precision. 
 549  # 
 550  #          nfev is an integer output variable set to the number of 
 551  #                calls to fcn. 
 552  # 
 553  #          fjac is an output m by n array. the upper n by n submatrix 
 554  #                of fjac contains an upper triangular matrix r with 
 555  #                diagonal elements of nonincreasing magnitude such that 
 556  # 
 557  #                               t        t                 t 
 558  #                          p *(jac *jac)*p = r *r, 
 559  # 
 560  #                where p is a permutation matrix and jac is the final 
 561  #                calculated jacobian. column j of p is column ipvt(j) 
 562  #                (see below) of the identity matrix. the lower trapezoidal 
 563  #                part of fjac contains information generated during 
 564  #                the computation of r. 
 565  # 
 566  #          ldfjac is a positive integer input variable not less than m 
 567  #                which specifies the leading dimension of the array fjac. 
 568  # 
 569  #          ipvt is an integer output array of length n. ipvt 
 570  #                defines a permutation matrix p such that jac*p = q*r, 
 571  #                where jac is the final calculated jacobian, q is 
 572  #                orthogonal (not stored), and r is upper triangular 
 573  #                with diagonal elements of nonincreasing magnitude. 
 574  #                column j of p is column ipvt(j) of the identity matrix. 
 575  # 
 576  #          qtf is an output array of length n which contains 
 577  #                the first n elements of the vector (q transpose)*fvec. 
 578  # 
 579  #          wa1, wa2, and wa3 are work arrays of length n. 
 580  # 
 581  #          wa4 is a work array of length m. 
 582  # 
 583  #        subprograms called 
 584  # 
 585  #          user-supplied ...... fcn 
 586  # 
 587  #          minpack-supplied ... dpmpar,enorm,fdjac2,,qrfac 
 588  # 
 589  #          fortran-supplied ... dabs,dmax1,dmin1,dsqrt,mod 
 590  # 
 591  #        argonne national laboratory. minpack project. march 1980. 
 592  #        burton s. garbow, kenneth e. hillstrom, jorge j. more 
 593  # 
 594  #        ********** 
 595   
596 -class mpfit:
597 - def __init__(self, fcn, xall=None, functkw={}, parinfo=None, 598 ftol=1.e-10, xtol=1.e-10, gtol=1.e-10, 599 damp=0., maxiter=200, factor=100., nprint=1, 600 iterfunct='default', iterkw={}, nocovar=0, 601 fastnorm=0, rescale=0, autoderivative=1, quiet=0, 602 diag=None, epsfcn=None, debug=0):
603 """ 604 Inputs: 605 fcn: 606 The function to be minimized. The function should return the weighted 607 deviations between the model and the data, as described above. 608 609 xall: 610 An array of starting values for each of the parameters of the model. 611 The number of parameters should be fewer than the number of measurements. 612 613 This parameter is optional if the parinfo keyword is used (but see 614 parinfo). The parinfo keyword provides a mechanism to fix or constrain 615 individual parameters. 616 617 Keywords: 618 619 autoderivative: 620 If this is set, derivatives of the function will be computed 621 automatically via a finite differencing procedure. If not set, then 622 fcn must provide the (analytical) derivatives. 623 Default: set (=1) 624 NOTE: to supply your own analytical derivatives, 625 explicitly pass autoderivative=0 626 627 fastnorm: 628 Set this keyword to select a faster algorithm to compute sum-of-square 629 values internally. For systems with large numbers of data points, the 630 standard algorithm can become prohibitively slow because it cannot be 631 vectorized well. By setting this keyword, MPFIT will run faster, but 632 it will be more prone to floating point overflows and underflows. Thus, setting 633 this keyword may sacrifice some stability in the fitting process. 634 Default: clear (=0) 635 636 ftol: 637 A nonnegative input variable. Termination occurs when both the actual 638 and predicted relative reductions in the sum of squares are at most 639 ftol (and status is accordingly set to 1 or 3). Therefore, ftol 640 measures the relative error desired in the sum of squares. 641 Default: 1E-10 642 643 functkw: 644 A dictionary which contains the parameters to be passed to the 645 user-supplied function specified by fcn via the standard Python 646 keyword dictionary mechanism. This is the way you can pass additional 647 data to your user-supplied function without using global variables. 648 649 Consider the following example: 650 if functkw = {'xval':[1.,2.,3.], 'yval':[1.,4.,9.], 651 'errval':[1.,1.,1.] } 652 then the user supplied function should be declared like this: 653 def myfunct(p, fjac=None, xval=None, yval=None, errval=None): 654 655 Default: {} No extra parameters are passed to the user-supplied 656 function. 657 658 gtol: 659 A nonnegative input variable. Termination occurs when the cosine of 660 the angle between fvec and any column of the jacobian is at most gtol 661 in absolute value (and status is accordingly set to 4). Therefore, 662 gtol measures the orthogonality desired between the function vector 663 and the columns of the jacobian. 664 Default: 1e-10 665 666 iterkw: 667 The keyword arguments to be passed to iterfunct via the dictionary 668 keyword mechanism. This should be a dictionary and is similar in 669 operation to FUNCTKW. 670 Default: {} No arguments are passed. 671 672 iterfunct: 673 The name of a function to be called upon each NPRINT iteration of the 674 MPFIT routine. It should be declared in the following way: 675 def iterfunct(myfunct, p, iter, fnorm, functkw=None, 676 parinfo=None, quiet=0, dof=None, [iterkw keywords here]) 677 # perform custom iteration update 678 679 iterfunct must accept all three keyword parameters (FUNCTKW, PARINFO 680 and QUIET). 681 682 myfunct: The user-supplied function to be minimized, 683 p: The current set of model parameters 684 iter: The iteration number 685 functkw: The arguments to be passed to myfunct. 686 fnorm: The chi-squared value. 687 quiet: Set when no textual output should be printed. 688 dof: The number of degrees of freedom, normally the number of points 689 less the number of free parameters. 690 See below for documentation of parinfo. 691 692 In implementation, iterfunct can perform updates to the terminal or 693 graphical user interface, to provide feedback while the fit proceeds. 694 If the fit is to be stopped for any reason, then iterfunct should return a 695 a status value between -15 and -1. Otherwise it should return None 696 (e.g. no return statement) or 0. 697 In principle, iterfunct should probably not modify the parameter values, 698 because it may interfere with the algorithm's stability. In practice it 699 is allowed. 700 701 Default: an internal routine is used to print the parameter values. 702 703 Set iterfunct=None if there is no user-defined routine and you don't 704 want the internal default routine be called. 705 706 maxiter: 707 The maximum number of iterations to perform. If the number is exceeded, 708 then the status value is set to 5 and MPFIT returns. 709 Default: 200 iterations 710 711 nocovar: 712 Set this keyword to prevent the calculation of the covariance matrix 713 before returning (see COVAR) 714 Default: clear (=0) The covariance matrix is returned 715 716 nprint: 717 The frequency with which iterfunct is called. A value of 1 indicates 718 that iterfunct is called with every iteration, while 2 indicates every 719 other iteration, etc. Note that several Levenberg-Marquardt attempts 720 can be made in a single iteration. 721 Default value: 1 722 723 parinfo 724 Provides a mechanism for more sophisticated constraints to be placed on 725 parameter values. When parinfo is not passed, then it is assumed that 726 all parameters are free and unconstrained. Values in parinfo are never 727 modified during a call to MPFIT. 728 729 See description above for the structure of PARINFO. 730 731 Default value: None All parameters are free and unconstrained. 732 733 quiet: 734 Set this keyword when no textual output should be printed by MPFIT 735 736 damp: 737 A scalar number, indicating the cut-off value of residuals where 738 "damping" will occur. Residuals with magnitudes greater than this 739 number will be replaced by their hyperbolic tangent. This partially 740 mitigates the so-called large residual problem inherent in 741 least-squares solvers (as for the test problem CURVI, 742 http://www.maxthis.com/curviex.htm). 743 A value of 0 indicates no damping. 744 Default: 0 745 746 Note: DAMP doesn't work with autoderivative=0 747 748 xtol: 749 A nonnegative input variable. Termination occurs when the relative error 750 between two consecutive iterates is at most xtol (and status is 751 accordingly set to 2 or 3). Therefore, xtol measures the relative error 752 desired in the approximate solution. 753 Default: 1E-10 754 755 Outputs: 756 757 Returns an object of type mpfit. The results are attributes of this class, 758 e.g. mpfit.status, mpfit.errmsg, mpfit.params, npfit.niter, mpfit.covar. 759 760 .status 761 An integer status code is returned. All values greater than zero can 762 represent success (however .status == 5 may indicate failure to 763 converge). It can have one of the following values: 764 765 -16 766 A parameter or function value has become infinite or an undefined 767 number. This is usually a consequence of numerical overflow in the 768 user's model function, which must be avoided. 769 770 -15 to -1 771 These are error codes that either MYFUNCT or iterfunct may return to 772 terminate the fitting process. Values from -15 to -1 are reserved 773 for the user functions and will not clash with MPFIT. 774 775 0 Improper input parameters. 776 777 1 Both actual and predicted relative reductions in the sum of squares 778 are at most ftol. 779 780 2 Relative error between two consecutive iterates is at most xtol 781 782 3 Conditions for status = 1 and status = 2 both hold. 783 784 4 The cosine of the angle between fvec and any column of the jacobian 785 is at most gtol in absolute value. 786 787 5 The maximum number of iterations has been reached. 788 789 6 ftol is too small. No further reduction in the sum of squares is 790 possible. 791 792 7 xtol is too small. No further improvement in the approximate solution 793 x is possible. 794 795 8 gtol is too small. fvec is orthogonal to the columns of the jacobian 796 to machine precision. 797 798 .fnorm 799 The value of the summed squared residuals for the returned parameter 800 values. 801 802 .covar 803 The covariance matrix for the set of parameters returned by MPFIT. 804 The matrix is NxN where N is the number of parameters. The square root 805 of the diagonal elements gives the formal 1-sigma statistical errors on 806 the parameters if errors were treated "properly" in fcn. 807 Parameter errors are also returned in .perror. 808 809 To compute the correlation matrix, pcor, use this example: 810 cov = mpfit.covar 811 pcor = cov * 0. 812 for i in range(n): 813 for j in range(n): 814 pcor[i,j] = cov[i,j]/sqrt(cov[i,i]*cov[j,j]) 815 816 If nocovar is set or MPFIT terminated abnormally, then .covar is set to 817 a scalar with value None. 818 819 .errmsg 820 A string error or warning message is returned. 821 822 .nfev 823 The number of calls to MYFUNCT performed. 824 825 .niter 826 The number of iterations completed. 827 828 .perror 829 The formal 1-sigma errors in each parameter, computed from the 830 covariance matrix. If a parameter is held fixed, or if it touches a 831 boundary, then the error is reported as zero. 832 833 If the fit is unweighted (i.e. no errors were given, or the weights 834 were uniformly set to unity), then .perror will probably not represent 835 the true parameter uncertainties. 836 837 *If* you can assume that the true reduced chi-squared value is unity -- 838 meaning that the fit is implicitly assumed to be of good quality -- 839 then the estimated parameter uncertainties can be computed by scaling 840 .perror by the measured chi-squared value. 841 842 dof = len(x) - len(mpfit.params) # deg of freedom 843 # scaled uncertainties 844 pcerror = mpfit.perror * sqrt(mpfit.fnorm / dof) 845 846 """ 847 self.niter = 0 848 self.params = None 849 self.covar = None 850 self.perror = None 851 self.status = 0 # Invalid input flag set while we check inputs 852 self.debug = debug 853 self.errmsg = '' 854 self.fastnorm = fastnorm 855 self.nfev = 0 856 self.damp = damp 857 self.machar = machar(double=1) 858 machep = self.machar.machep 859 self.dof=0 860 861 if (fcn==None): 862 self.errmsg = "Usage: parms = mpfit('myfunt', ... )" 863 return 864 865 if (iterfunct == 'default'): iterfunct = self.defiter 866 867 ## Parameter damping doesn't work when user is providing their own 868 ## gradients. 869 if (self.damp != 0) and (autoderivative == 0): 870 self.errmsg = 'ERROR: keywords DAMP and AUTODERIVATIVE are mutually exclusive' 871 return 872 873 ## Parameters can either be stored in parinfo, or x. x takes precedence if it exists 874 if (xall == None) and (parinfo == None): 875 self.errmsg = 'ERROR: must pass parameters in P or PARINFO' 876 return 877 878 ## Be sure that PARINFO is of the right type 879 if (parinfo != None): 880 if (type(parinfo) != types.ListType): 881 self.errmsg = 'ERROR: PARINFO must be a list of dictionaries.' 882 return 883 else: 884 if (type(parinfo[0]) != types.DictionaryType): 885 self.errmsg = 'ERROR: PARINFO must be a list of dictionaries.' 886 return 887 if ((xall != None) and (len(xall) != len(parinfo))): 888 self.errmsg = 'ERROR: number of elements in PARINFO and P must agree' 889 return 890 891 ## If the parameters were not specified at the command line, then 892 ## extract them from PARINFO 893 if (xall == None): 894 xall = self.parinfo(parinfo, 'value') 895 if (xall == None): 896 self.errmsg = 'ERROR: either P or PARINFO(*)["value"] must be supplied.' 897 return 898 899 ## Make sure parameters are Numeric arrays of type Float 900 xall = asarray(xall, float) 901 902 npar = len(xall) 903 self.fnorm = -1. 904 fnorm1 = -1. 905 906 ## TIED parameters? 907 ptied = self.parinfo(parinfo, 'tied', default='', n=npar) 908 self.qanytied = 0 909 for i in range(npar): 910 ptied[i] = ptied[i].strip() 911 if (ptied[i] != ''): self.qanytied = 1 912 self.ptied = ptied 913 914 ## FIXED parameters ? 915 pfixed = self.parinfo(parinfo, 'fixed', default=0, n=npar) 916 pfixed = (pfixed == 1) 917 for i in range(npar): 918 pfixed[i] = pfixed[i] or (ptied[i] != '') ## Tied parameters are also effectively fixed 919 920 ## Finite differencing step, absolute and relative, and sidedness of deriv. 921 step = self.parinfo(parinfo, 'step', default=0., n=npar) 922 dstep = self.parinfo(parinfo, 'relstep', default=0., n=npar) 923 dside = self.parinfo(parinfo, 'mpside', default=0, n=npar) 924 925 ## Maximum and minimum steps allowed to be taken in one iteration 926 maxstep = self.parinfo(parinfo, 'mpmaxstep', default=0., n=npar) 927 minstep = self.parinfo(parinfo, 'mpminstep', default=0., n=npar) 928 qmin = minstep * 0 ## Remove minstep for now!! 929 qmax = maxstep != 0 930 wh = (nonzero(((qmin!=0.) & (qmax!=0.)) & (maxstep < minstep)))[0] 931 if (len(wh) > 0): 932 self.errmsg = 'ERROR: MPMINSTEP is greater than MPMAXSTEP' 933 return 934 wh = (nonzero((qmin!=0.) | (qmax!=0.)))[0] 935 qminmax = len(wh > 0) 936 937 ## Finish up the free parameters 938 ifree = (nonzero(pfixed != 1))[0] 939 nfree = len(ifree) 940 if nfree == 0: 941 self.errmsg = 'ERROR: no free parameters' 942 return 943 944 ## Compose only VARYING parameters 945 self.params = xall ## self.params is the set of parameters to be returned 946 x = take(self.params, ifree) ## x is the set of free parameters 947 948 ## LIMITED parameters ? 949 limited = self.parinfo(parinfo, 'limited', default=[0,0], n=npar) 950 limits = self.parinfo(parinfo, 'limits', default=[0.,0.], n=npar) 951 if (limited != None) and (limits != None): 952 ## Error checking on limits in parinfo 953 wh = (nonzero((limited[:,0] & (xall < limits[:,0])) | 954 (limited[:,1] & (xall > limits[:,1]))))[0] 955 if (len(wh) > 0): 956 self.errmsg = 'ERROR: parameters are not within PARINFO limits' 957 return 958 wh = (nonzero((limited[:,0] & limited[:,1]) & 959 (limits[:,0] >= limits[:,1]) & 960 (pfixed == 0)))[0] 961 if (len(wh) > 0): 962 self.errmsg = 'ERROR: PARINFO parameter limits are not consistent' 963 return 964 965 ## Transfer structure values to local variables 966 qulim = take(limited[:,1], ifree) 967 ulim = take(limits [:,1], ifree) 968 qllim = take(limited[:,0], ifree) 969 llim = take(limits [:,0], ifree) 970 971 wh = (nonzero((qulim!=0.) | (qllim!=0.)))[0] 972 if (len(wh) > 0): qanylim = 1 973 else: qanylim = 0 974 else: 975 ## Fill in local variables with dummy values 976 qulim = zeros(nfree) 977 ulim = x * 0. 978 qllim = qulim 979 llim = x * 0. 980 qanylim = 0 981 982 n = len(x) 983 ## Check input parameters for errors 984 if ((n < 0) or (ftol <= 0) or (xtol <= 0) or (gtol <= 0) 985 or (maxiter < 0) or (factor <= 0)): 986 self.errmsg = 'ERROR: input keywords are inconsistent' 987 return 988 989 if (rescale != 0): 990 self.errmsg = 'ERROR: DIAG parameter scales are inconsistent' 991 if (len(diag) < n): return 992 wh = (nonzero(diag <= 0))[0] 993 if (len(wh) > 0): return 994 self.errmsg = '' 995 996 # Make sure x is a Numeric array of type Float 997 x = asarray(x, float) 998 999 [self.status, fvec] = self.call(fcn, self.params, functkw) 1000 if (self.status < 0): 1001 self.errmsg = 'ERROR: first call to "'+str(fcn)+'" failed' 1002 return 1003 1004 m = len(fvec) 1005 if (m < n): 1006 self.errmsg = 'ERROR: number of parameters must not exceed data' 1007 return 1008 self.dof = m-nfree 1009 self.fnorm = self.enorm(fvec) 1010 1011 ## Initialize Levelberg-Marquardt parameter and iteration counter 1012 1013 par = 0. 1014 self.niter = 1 1015 qtf = x * 0. 1016 self.status = 0 1017 1018 ## Beginning of the outer loop 1019 1020 while(1): 1021 1022 ## If requested, call fcn to enable printing of iterates 1023 put(self.params, ifree, x) 1024 if (self.qanytied): self.params = self.tie(self.params, ptied) 1025 1026 if (nprint > 0) and (iterfunct != None): 1027 if (((self.niter-1) % nprint) == 0): 1028 mperr = 0 1029 xnew0 = self.params.copy() 1030 1031 dof = numpy.max([len(fvec) - len(x), 0]) 1032 status = iterfunct(fcn, self.params, self.niter, self.fnorm**2, 1033 functkw=functkw, parinfo=parinfo, quiet=quiet, 1034 dof=dof, **iterkw) 1035 if (status != None): self.status = status 1036 1037 ## Check for user termination 1038 if (self.status < 0): 1039 self.errmsg = 'WARNING: premature termination by ' + str(iterfunct) 1040 return 1041 1042 ## If parameters were changed (grrr..) then re-tie 1043 if (numpy.max(abs(xnew0-self.params)) > 0): 1044 if (self.qanytied): self.params = self.tie(self.params, ptied) 1045 x = take(self.params, ifree) 1046 1047 1048 ## Calculate the jacobian matrix 1049 self.status = 2 1050 catch_msg = 'calling MPFIT_FDJAC2' 1051 fjac = self.fdjac2(fcn, x, fvec, step, qulim, ulim, dside, 1052 epsfcn=epsfcn, 1053 autoderivative=autoderivative, dstep=dstep, 1054 functkw=functkw, ifree=ifree, xall=self.params) 1055 if (fjac == None): 1056 self.errmsg = 'WARNING: premature termination by FDJAC2' 1057 return 1058 1059 ## Determine if any of the parameters are pegged at the limits 1060 if (qanylim): 1061 catch_msg = 'zeroing derivatives of pegged parameters' 1062 whlpeg = (nonzero(qllim & (x == llim)))[0] 1063 nlpeg = len(whlpeg) 1064 whupeg = (nonzero(qulim & (x == ulim)))[0] 1065 nupeg = len(whupeg) 1066 ## See if any "pegged" values should keep their derivatives 1067 if (nlpeg > 0): 1068 ## Total derivative of sum wrt lower pegged parameters 1069 for i in range(nlpeg): 1070 sum0 = sum(fvec * fjac[:,whlpeg[i]]) 1071 if (sum0 > 0): fjac[:,whlpeg[i]] = 0 1072 if (nupeg > 0): 1073 ## Total derivative of sum wrt upper pegged parameters 1074 for i in range(nupeg): 1075 sum0 = sum(fvec * fjac[:,whupeg[i]]) 1076 if (sum0 < 0): fjac[:,whupeg[i]] = 0 1077 1078 ## Compute the QR factorization of the jacobian 1079 [fjac, ipvt, wa1, wa2] = self.qrfac(fjac, pivot=1) 1080 1081 ## On the first iteration if "diag" is unspecified, scale 1082 ## according to the norms of the columns of the initial jacobian 1083 catch_msg = 'rescaling diagonal elements' 1084 if (self.niter == 1): 1085 if ((rescale==0) or (len(diag) < n)): 1086 diag = wa2.copy() 1087 wh = (nonzero(diag == 0))[0] 1088 put(diag, wh, 1.) 1089 1090 ## On the first iteration, calculate the norm of the scaled x 1091 ## and initialize the step bound delta 1092 wa3 = diag * x 1093 xnorm = self.enorm(wa3) 1094 delta = factor*xnorm 1095 if (delta == 0.): delta = factor 1096 1097 ## Form (q transpose)*fvec and store the first n components in qtf 1098 catch_msg = 'forming (q transpose)*fvec' 1099 wa4 = fvec.copy() 1100 for j in range(n): 1101 lj = ipvt[j] 1102 temp3 = fjac[j,lj] 1103 if (temp3 != 0): 1104 fj = fjac[j:,lj] 1105 wj = wa4[j:] 1106 ## *** optimization wa4(j:*) 1107 wa4[j:] = wj - fj * sum(fj*wj) / temp3 1108 fjac[j,lj] = wa1[j] 1109 qtf[j] = wa4[j] 1110 ## From this point on, only the square matrix, consisting of the 1111 ## triangle of R, is needed. 1112 fjac = fjac[0:n, 0:n] 1113 fjac.shape = [n, n] 1114 temp = fjac.copy() 1115 for i in range(n): 1116 temp[:,i] = fjac[:, ipvt[i]] 1117 fjac = temp.copy() 1118 1119 ## Check for overflow. This should be a cheap test here since FJAC 1120 ## has been reduced to a (small) square matrix, and the test is 1121 ## O(N^2). 1122 #wh = where(finite(fjac) EQ 0, ct) 1123 #if ct GT 0 then goto, FAIL_OVERFLOW 1124 1125 ## Compute the norm of the scaled gradient 1126 catch_msg = 'computing the scaled gradient' 1127 gnorm = 0. 1128 if (self.fnorm != 0): 1129 for j in range(n): 1130 l = ipvt[j] 1131 if (wa2[l] != 0): 1132 sum0 = sum(fjac[0:j+1,j]*qtf[0:j+1])/self.fnorm 1133 gnorm = numpy.max([gnorm,abs(sum0/wa2[l])]) 1134 1135 ## Test for convergence of the gradient norm 1136 if (gnorm <= gtol): 1137 self.status = 4 1138 break 1139 if maxiter == 0: 1140 self.status = 5 1141 break 1142 1143 ## Rescale if necessary 1144 if (rescale == 0): 1145 diag = choose(diag>wa2, (wa2, diag)) 1146 1147 ## Beginning of the inner loop 1148 while(1): 1149 1150 ## Determine the levenberg-marquardt parameter 1151 catch_msg = 'calculating LM parameter (MPFIT_)' 1152 [fjac, par, wa1, wa2] = self.lmpar(fjac, ipvt, diag, qtf, 1153 delta, wa1, wa2, par=par) 1154 ## Store the direction p and x+p. Calculate the norm of p 1155 wa1 = -wa1 1156 1157 if (qanylim == 0) and (qminmax == 0): 1158 ## No parameter limits, so just move to new position WA2 1159 alpha = 1. 1160 wa2 = x + wa1 1161 1162 else: 1163 1164 ## Respect the limits. If a step were to go out of bounds, then 1165 ## we should take a step in the same direction but shorter distance. 1166 ## The step should take us right to the limit in that case. 1167 alpha = 1. 1168 1169 if (qanylim): 1170 ## Do not allow any steps out of bounds 1171 catch_msg = 'checking for a step out of bounds' 1172 if (nlpeg > 0): 1173 put(wa1, whlpeg, clip( 1174 take(wa1, whlpeg), 0., numpy.max(wa1))) 1175 if (nupeg > 0): 1176 put(wa1, whupeg, clip( 1177 take(wa1, whupeg), numpy.min(wa1), 0.)) 1178 1179 dwa1 = abs(wa1) > machep 1180 whl = (nonzero(((dwa1!=0.) & qllim) & ((x + wa1) < llim)))[0] 1181 if (len(whl) > 0): 1182 t = ((take(llim, whl) - take(x, whl)) / 1183 take(wa1, whl)) 1184 alpha = numpy.min([alpha, numpy.min(t)]) 1185 whu = (nonzero(((dwa1!=0.) & qulim) & ((x + wa1) > ulim)))[0] 1186 if (len(whu) > 0): 1187 t = ((take(ulim, whu) - take(x, whu)) / 1188 take(wa1, whu)) 1189 alpha = numpy.min([alpha, numpy.min(t)]) 1190 1191 ## Obey any max step values. 1192 if (qminmax): 1193 nwa1 = wa1 * alpha 1194 whmax = (nonzero((qmax != 0.) & (maxstep > 0)))[0] 1195 if (len(whmax) > 0): 1196 mrat = numpy.max(numpy.abs(nwa1[whmax]) / 1197 numpy.abs(maxstep[ifree[whmax]])) 1198 if (mrat > 1): alpha = alpha / mrat 1199 1200 ## Scale the resulting vector 1201 wa1 = wa1 * alpha 1202 wa2 = x + wa1 1203 1204 ## Adjust the final output values. If the step put us exactly 1205 ## on a boundary, make sure it is exact. 1206 sgnu = (ulim >= 0) * 2. - 1. 1207 sgnl = (llim >= 0) * 2. - 1. 1208 ## Handles case of 1209 ## ... nonzero *LIM ... ...zero * LIM 1210 ulim1 = ulim * (1 - sgnu * machep) - (ulim == 0) * machep 1211 llim1 = llim * (1 + sgnl * machep) + (llim == 0) * machep 1212 wh = (nonzero((qulim!=0) & (wa2 >= ulim1)))[0] 1213 if (len(wh) > 0): put(wa2, wh, take(ulim, wh)) 1214 wh = (nonzero((qllim!=0.) & (wa2 <= llim1)))[0] 1215 if (len(wh) > 0): put(wa2, wh, take(llim, wh)) 1216 # endelse 1217 wa3 = diag * wa1 1218 pnorm = self.enorm(wa3) 1219 1220 ## On the first iteration, adjust the initial step bound 1221 if (self.niter == 1): delta = numpy.min([delta,pnorm]) 1222 1223 put(self.params, ifree, wa2) 1224 1225 ## Evaluate the function at x+p and calculate its norm 1226 mperr = 0 1227 catch_msg = 'calling '+str(fcn) 1228 [self.status, wa4] = self.call(fcn, self.params, functkw) 1229 if (self.status < 0): 1230 self.errmsg = 'WARNING: premature termination by "'+fcn+'"' 1231 return 1232 fnorm1 = self.enorm(wa4) 1233 1234 ## Compute the scaled actual reduction 1235 catch_msg = 'computing convergence criteria' 1236 actred = -1. 1237 if ((0.1 * fnorm1) < self.fnorm): actred = - (fnorm1/self.fnorm)**2 + 1. 1238 1239 ## Compute the scaled predicted reduction and the scaled directional 1240 ## derivative 1241 for j in range(n): 1242 wa3[j] = 0 1243 wa3[0:j+1] = wa3[0:j+1] + fjac[0:j+1,j]*wa1[ipvt[j]] 1244 1245 ## Remember, alpha is the fraction of the full LM step actually 1246 ## taken 1247 temp1 = self.enorm(alpha*wa3)/self.fnorm 1248 temp2 = (sqrt(alpha*par)*pnorm)/self.fnorm 1249 prered = temp1*temp1 + (temp2*temp2)/0.5 1250 dirder = -(temp1*temp1 + temp2*temp2) 1251 1252 ## Compute the ratio of the actual to the predicted reduction. 1253 ratio = 0. 1254 if (prered != 0): ratio = actred/prered 1255 1256 ## Update the step bound 1257 if (ratio <= 0.25): 1258 if (actred >= 0): temp = .5 1259 else: temp = .5*dirder/(dirder + .5*actred) 1260 if ((0.1*fnorm1) >= self.fnorm) or (temp < 0.1): temp = 0.1 1261 delta = temp*numpy.min([delta,pnorm/0.1]) 1262 par = par/temp 1263 else: 1264 if (par == 0) or (ratio >= 0.75): 1265 delta = pnorm/.5 1266 par = .5*par 1267 1268 ## Test for successful iteration 1269 if (ratio >= 0.0001): 1270 ## Successful iteration. Update x, fvec, and their norms 1271 x = wa2 1272 wa2 = diag * x 1273 fvec = wa4 1274 xnorm = self.enorm(wa2) 1275 self.fnorm = fnorm1 1276 self.niter = self.niter + 1 1277 1278 ## Tests for convergence 1279 if ((abs(actred) <= ftol) and (prered <= ftol) 1280 and (0.5 * ratio <= 1)): self.status = 1 1281 if delta <= xtol*xnorm: self.status = 2 1282 if ((abs(actred) <= ftol) and (prered <= ftol) 1283 and (0.5 * ratio <= 1) and (self.status == 2)): self.status = 3 1284 if (self.status != 0): break 1285 1286 ## Tests for termination and stringent tolerances 1287 if (self.niter >= maxiter): self.status = 5 1288 if ((abs(actred) <= machep) and (prered <= machep) 1289 and (0.5*ratio <= 1)): self.status = 6 1290 if delta <= machep*xnorm: self.status = 7 1291 if gnorm <= machep: self.status = 8 1292 if (self.status != 0): break 1293 1294 ## End of inner loop. Repeat if iteration unsuccessful 1295 if (ratio >= 0.0001): break 1296 1297 ## Check for over/underflow 1298 if ~numpy.all(numpy.isfinite(wa1) & numpy.isfinite(wa2) & 1299 numpy.isfinite(x)) or ~numpy.isfinite(ratio): 1300 errmsg = ('''ERROR: parameter or function value(s) have become 1301 'infinite; check model function for over- 'and underflow''') 1302 self.status = -16 1303 break 1304 ##wh = where(finite(wa1) EQ 0 OR finite(wa2) EQ 0 OR finite(x) EQ 0, ct) 1305 ##if ct GT 0 OR finite(ratio) EQ 0 then begin 1306 1307 if (self.status != 0): break; 1308 ## End of outer loop. 1309 1310 catch_msg = 'in the termination phase' 1311 ## Termination, either normal or user imposed. 1312 if (len(self.params) == 0): 1313 return 1314 #if xnew==None: # failed patch WR 1315 # return 1316 if (nfree == 0): self.params = xall.copy() 1317 else: put(self.params, ifree, x) 1318 if (nprint > 0) and (self.status > 0): 1319 catch_msg = 'calling ' + str(fcn) 1320 [status, fvec] = self.call(fcn, self.params, functkw) 1321 catch_msg = 'in the termination phase' 1322 self.fnorm = self.enorm(fvec) 1323 1324 if ((self.fnorm != None) and (fnorm1 != None)): 1325 self.fnorm = numpy.max([self.fnorm, fnorm1]) 1326 self.fnorm = self.fnorm**2. 1327 1328 self.covar = None 1329 self.perror = None 1330 ## (very carefully) set the covariance matrix COVAR 1331 if ((self.status > 0) and (nocovar==0) and (n != None) 1332 and (fjac != None) and (ipvt != None)): 1333 sz = shape(fjac) 1334 if ((n > 0) and (sz[0] >= n) and (sz[1] >= n) 1335 and (len(ipvt) >= n)): 1336 1337 catch_msg = 'computing the covariance matrix' 1338 cv = self.calc_covar(fjac[0:n,0:n], ipvt[0:n]) 1339 cv.shape = [n, n] 1340 nn = len(xall) 1341 1342 ## Fill in actual covariance matrix, accounting for fixed 1343 ## parameters. 1344 self.covar = zeros([nn, nn], dtype=float) 1345 for i in range(n): 1346 self.covar[ifree,ifree[i]] = cv[:,i] 1347 1348 ## Compute errors in parameters 1349 catch_msg = 'computing parameter errors' 1350 self.perror = zeros(nn, dtype=float) 1351 d = diagonal(self.covar) 1352 wh = (nonzero(d >= 0))[0] 1353 if len(wh) > 0: 1354 put(self.perror, wh, sqrt(take(d, wh))) 1355 return
1356 1357
1358 - def __str__(self):
1359 return {'params': self.params, 1360 'niter': self.niter, 1361 'params': self.params, 1362 'covar': self.covar, 1363 'perror': self.perror, 1364 'status': self.status, 1365 'debug': self.debug, 1366 'errmsg': self.errmsg, 1367 'fastnorm': self.fastnorm, 1368 'nfev': self.nfev, 1369 'damp': self.damp 1370 #,'machar':self.machar 1371 }.__str__()
1372 1373 ## Default procedure to be called every iteration. It simply prints 1374 ## the parameter values.
1375 - def defiter(self, fcn, x, iter, fnorm=None, functkw=None, 1376 quiet=0, iterstop=None, parinfo=None, 1377 format=None, pformat='%.10g', dof=1):
1378 1379 if (self.debug): print 'Entering defiter...' 1380 if (quiet): return 1381 if (fnorm == None): 1382 [status, fvec] = self.call(fcn, x, functkw) 1383 fnorm = self.enorm(fvec)**2 1384 1385 ## Determine which parameters to print 1386 nprint = len(x) 1387 print "Iter ", ('%6i' % iter)," CHI-SQUARE = ",('%.10g' % fnorm)," DOF = ", ('%i' % dof) 1388 for i in range(nprint): 1389 if (parinfo != None) and (parinfo[i].has_key('parname')): 1390 p = ' ' + parinfo[i]['parname'] + ' = ' 1391 else: 1392 p = ' P' + str(i) + ' = ' 1393 if (parinfo != None) and (parinfo[i].has_key('mpprint')): 1394 iprint = parinfo[i]['mpprint'] 1395 else: 1396 iprint = 1 1397 if (iprint): 1398 print p + (pformat % x[i]) + ' ' 1399 return(0)
1400 1401 ## DO_ITERSTOP: 1402 ## if keyword_set(iterstop) then begin 1403 ## k = get_kbrd(0) 1404 ## if k EQ string(byte(7)) then begin 1405 ## message, 'WARNING: minimization not complete', /info 1406 ## print, 'Do you want to terminate this procedure? (y/n)', $ 1407 ## format='(A,$)' 1408 ## k = '' 1409 ## read, k 1410 ## if strupcase(strmid(k,0,1)) EQ 'Y' then begin 1411 ## message, 'WARNING: Procedure is terminating.', /info 1412 ## mperr = -1 1413 ## endif 1414 ## endif 1415 ## endif 1416 1417 1418 ## Procedure to parse the parameter values in PARINFO, which is a list of dictionaries
1419 - def parinfo(self, parinfo=None, key='a', default=None, n=0):
1420 if (self.debug): print 'Entering parinfo...' 1421 if (n == 0) and (parinfo != None): n = len(parinfo) 1422 if (n == 0): 1423 values = default 1424 return(values) 1425 values = [] 1426 for i in range(n): 1427 if ((parinfo != None) and (parinfo[i].has_key(key))): 1428 values.append(parinfo[i][key]) 1429 else: 1430 values.append(default) 1431 1432 # Convert to numeric arrays if possible 1433 test = default 1434 if (type(default) == types.ListType): test=default[0] 1435 if (isinstance(test, types.IntType)): 1436 values = asarray(values, int) 1437 elif (isinstance(test, types.FloatType)): 1438 values = asarray(values, float) 1439 return(values)
1440 1441 ## Call user function or procedure, with _EXTRA or not, with 1442 ## derivatives or not.
1443 - def call(self, fcn, x, functkw, fjac=None):
1444 if (self.debug): print 'Entering call...' 1445 if (self.qanytied): x = self.tie(x, self.ptied) 1446 self.nfev = self.nfev + 1 1447 if (fjac == None): 1448 [status, f] = fcn(x, fjac=fjac, **functkw) 1449 if (self.damp > 0): 1450 ## Apply the damping if requested. This replaces the residuals 1451 ## with their hyperbolic tangent. Thus residuals larger than 1452 ## DAMP are essentially clipped. 1453 f = tanh(f/self.damp) 1454 return([status, f]) 1455 else: 1456 return(fcn(x, fjac=fjac, **functkw))
1457 1458
1459 - def enorm(self, vec):
1460 1461 if (self.debug): print 'Entering enorm...' 1462 ## NOTE: it turns out that, for systems that have a lot of data 1463 ## points, this routine is a big computing bottleneck. The extended 1464 ## computations that need to be done cannot be effectively 1465 ## vectorized. The introduction of the FASTNORM configuration 1466 ## parameter allows the user to select a faster routine, which is 1467 ## based on TOTAL() alone. 1468 1469 # Very simple-minded sum-of-squares 1470 if (self.fastnorm): 1471 ans = sqrt(sum(vec*vec)) 1472 else: 1473 agiant = self.machar.rgiant / len(vec) 1474 adwarf = self.machar.rdwarf * len(vec) 1475 1476 ## This is hopefully a compromise between speed and robustness. 1477 ## Need to do this because of the possibility of over- or underflow. 1478 mx = numpy.max(vec) 1479 mn = numpy.min(vec) 1480 mx = max([abs(mx), abs(mn)]) 1481 if mx == 0: return(vec[0]*0.) 1482 if mx > agiant or mx < adwarf: 1483 ans = mx * sqrt(sum((vec/mx)*(vec/mx))) 1484 else: 1485 ans = sqrt(sum(vec*vec)) 1486 1487 return(ans)
1488 1489
1490 - def fdjac2(self, fcn, x, fvec, step=None, ulimited=None, ulimit=None, dside=None, 1491 epsfcn=None, autoderivative=1, 1492 functkw=None, xall=None, ifree=None, dstep=None):
1493 1494 if (self.debug): print 'Entering fdjac2...' 1495 machep = self.machar.machep 1496 if epsfcn == None: epsfcn = machep 1497 if xall == None: xall = x 1498 if ifree == None: ifree = arange(len(xall)) 1499 if step == None: step = x * 0. 1500 nall = len(xall) 1501 1502 eps = sqrt(numpy.max([epsfcn, machep])) 1503 m = len(fvec) 1504 n = len(x) 1505 1506 ## Compute analytical derivative if requested 1507 if (autoderivative == 0): 1508 mperr = 0 1509 fjac = zeros(nall, dtype=float) 1510 Put(fjac, ifree, 1.0) ## Specify which parameters need derivatives 1511 [status, fp] = self.call(fcn, xall, functkw, fjac=fjac) 1512 1513 if len(fjac) != m*nall: 1514 print 'ERROR: Derivative matrix was not computed properly.' 1515 return(None) 1516 1517 ## This definition is consistent with CURVEFIT 1518 ## Sign error found (thanks Jesus Fernandez <fernande@irm.chu-caen.fr>) 1519 fjac.shape = [m,nall] 1520 fjac = -fjac 1521 1522 ## Select only the free parameters 1523 if len(ifree) < nall: 1524 fjac = fjac[:,ifree] 1525 fjac.shape = [m, n] 1526 return(fjac) 1527 1528 fjac = zeros([m, n], dtype=float) 1529 1530 h = eps * abs(x) 1531 1532 ## if STEP is given, use that 1533 ## STEP includes the fixed parameters 1534 if step != None: 1535 stepi = take(step, ifree) 1536 wh = (nonzero(stepi > 0))[0] 1537 if (len(wh) > 0): put(h, wh, take(stepi, wh)) 1538 1539 ## if relative step is given, use that 1540 ## DSTEP includes the fixed parameters 1541 if (len(dstep) > 0): 1542 dstepi = take(dstep, ifree) 1543 wh = (nonzero(dstepi > 0))[0] 1544 if len(wh) > 0: put(h, wh, abs(take(dstepi,wh)*take(x,wh))) 1545 1546 ## In case any of the step values are zero 1547 wh = (nonzero(h == 0))[0] 1548 if len(wh) > 0: put(h, wh, eps) 1549 1550 ## Reverse the sign of the step if we are up against the parameter 1551 ## limit, or if the user requested it. 1552 ## DSIDE includes the fixed parameters (ULIMITED/ULIMIT have only 1553 ## varying ones) 1554 mask = dside[ifree] == -1 1555 if len(ulimited) > 0 and len(ulimit) > 0: 1556 mask = logical_or((mask!=0), logical_and((ulimited!=0),(x > ulimit-h))) 1557 wh = (nonzero(mask))[0] 1558 if len(wh) > 0: put(h, wh, -take(h, wh)) 1559 ## Loop through parameters, computing the derivative for each 1560 for j in range(n): 1561 xp = xall.copy() 1562 xp[ifree[j]] = xp[ifree[j]] + h[j] 1563 [status, fp] = self.call(fcn, xp, functkw) 1564 if (status < 0): return(None) 1565 1566 if abs(dside[ifree[j]]) <= 1: 1567 ## COMPUTE THE ONE-SIDED DERIVATIVE 1568 ## Note optimization fjac(0:*,j) 1569 fjac[0:,j] = (fp-fvec)/h[j] 1570 1571 else: 1572 ## COMPUTE THE TWO-SIDED DERIVATIVE 1573 xp[ifree[j]] = xall[ifree[j]] - h[j] 1574 1575 mperr = 0 1576 [status, fm] = self.call(fcn, xp, functkw) 1577 if (status < 0): return(None) 1578 1579 ## Note optimization fjac(0:*,j) 1580 fjac[0:,j] = (fp-fm)/(2*h[j]) 1581 return(fjac)
1582 1583 1584 1585 # Original FORTRAN documentation 1586 # ********** 1587 # 1588 # subroutine qrfac 1589 # 1590 # this subroutine uses householder transformations with column 1591 # pivoting (optional) to compute a qr factorization of the 1592 # m by n matrix a. that is, qrfac determines an orthogonal 1593 # matrix q, a permutation matrix p, and an upper trapezoidal 1594 # matrix r with diagonal elements of nonincreasing magnitude, 1595 # such that a*p = q*r. the householder transformation for 1596 # column k, k = 1,2,...,min(m,n), is of the form 1597 # 1598 # t 1599 # i - (1/u(k))*u*u 1600 # 1601 # where u has zeros in the first k-1 positions. the form of 1602 # this transformation and the method of pivoting first 1603 # appeared in the corresponding linpack subroutine. 1604 # 1605 # the subroutine statement is 1606 # 1607 # subroutine qrfac(m,n,a,lda,pivot,ipvt,lipvt,rdiag,acnorm,wa) 1608 # 1609 # where 1610 # 1611 # m is a positive integer input variable set to the number 1612 # of rows of a. 1613 # 1614 # n is a positive integer input variable set to the number 1615 # of columns of a. 1616 # 1617 # a is an m by n array. on input a contains the matrix for 1618 # which the qr factorization is to be computed. on output 1619 # the strict upper trapezoidal part of a contains the strict 1620 # upper trapezoidal part of r, and the lower trapezoidal 1621 # part of a contains a factored form of q (the non-trivial 1622 # elements of the u vectors described above). 1623 # 1624 # lda is a positive integer input variable not less than m 1625 # which specifies the leading dimension of the array a. 1626 # 1627 # pivot is a logical input variable. if pivot is set true, 1628 # then column pivoting is enforced. if pivot is set false, 1629 # then no column pivoting is done. 1630 # 1631 # ipvt is an integer output array of length lipvt. ipvt 1632 # defines the permutation matrix p such that a*p = q*r. 1633 # column j of p is column ipvt(j) of the identity matrix. 1634 # if pivot is false, ipvt is not referenced. 1635 # 1636 # lipvt is a positive integer input variable. if pivot is false, 1637 # then lipvt may be as small as 1. if pivot is true, then 1638 # lipvt must be at least n. 1639 # 1640 # rdiag is an output array of length n which contains the 1641 # diagonal elements of r. 1642 # 1643 # acnorm is an output array of length n which contains the 1644 # norms of the corresponding columns of the input matrix a. 1645 # if this information is not needed, then acnorm can coincide 1646 # with rdiag. 1647 # 1648 # wa is a work array of length n. if pivot is false, then wa 1649 # can coincide with rdiag. 1650 # 1651 # subprograms called 1652 # 1653 # minpack-supplied ... dpmpar,enorm 1654 # 1655 # fortran-supplied ... dmax1,dsqrt,min0 1656 # 1657 # argonne national laboratory. minpack project. march 1980. 1658 # burton s. garbow, kenneth e. hillstrom, jorge j. more 1659 # 1660 # ********** 1661 # 1662 # PIVOTING / PERMUTING: 1663 # 1664 # Upon return, A(*,*) is in standard parameter order, A(*,IPVT) is in 1665 # permuted order. 1666 # 1667 # RDIAG is in permuted order. 1668 # ACNORM is in standard parameter order. 1669 # 1670 # 1671 # NOTE: in IDL the factors appear slightly differently than described 1672 # above. The matrix A is still m x n where m >= n. 1673 # 1674 # The "upper" triangular matrix R is actually stored in the strict 1675 # lower left triangle of A under the standard notation of IDL. 1676 # 1677 # The reflectors that generate Q are in the upper trapezoid of A upon 1678 # output. 1679 # 1680 # EXAMPLE: decompose the matrix [[9.,2.,6.],[4.,8.,7.]] 1681 # aa = [[9.,2.,6.],[4.,8.,7.]] 1682 # mpfit_qrfac, aa, aapvt, rdiag, aanorm 1683 # IDL> print, aa 1684 # 1.81818* 0.181818* 0.545455* 1685 # -8.54545+ 1.90160* 0.432573* 1686 # IDL> print, rdiag 1687 # -11.0000+ -7.48166+ 1688 # 1689 # The components marked with a * are the components of the 1690 # reflectors, and those marked with a + are components of R. 1691 # 1692 # To reconstruct Q and R we proceed as follows. First R. 1693 # r = fltarr(m, n) 1694 # for i = 0, n-1 do r(0:i,i) = aa(0:i,i) # fill in lower diag 1695 # r(lindgen(n)*(m+1)) = rdiag 1696 # 1697 # Next, Q, which are composed from the reflectors. Each reflector v 1698 # is taken from the upper trapezoid of aa, and converted to a matrix 1699 # via (I - 2 vT . v / (v . vT)). 1700 # 1701 # hh = ident ## identity matrix 1702 # for i = 0, n-1 do begin 1703 # v = aa(*,i) & if i GT 0 then v(0:i-1) = 0 ## extract reflector 1704 # hh = hh ## (ident - 2*(v # v)/total(v * v)) ## generate matrix 1705 # endfor 1706 # 1707 # Test the result: 1708 # IDL> print, hh ## transpose(r) 1709 # 9.00000 4.00000 1710 # 2.00000 8.00000 1711 # 6.00000 7.00000 1712 # 1713 # Note that it is usually never necessary to form the Q matrix 1714 # explicitly, and MPFIT does not. 1715 1716
1717 - def qrfac(self, a, pivot=0):
1718 1719 if (self.debug): print 'Entering qrfac...' 1720 machep = self.machar.machep 1721 sz = shape(a) 1722 m = sz[0] 1723 n = sz[1] 1724 1725 ## Compute the initial column norms and initialize arrays 1726 acnorm = zeros(n, dtype=float) 1727 for j in range(n): 1728 acnorm[j] = self.enorm(a[:,j]) 1729 rdiag = acnorm.copy() 1730 wa = rdiag.copy() 1731 ipvt = arange(n) 1732 1733 ## Reduce a to r with householder transformations 1734 minmn = numpy.min([m,n]) 1735 for j in range(minmn): 1736 if (pivot != 0): 1737 ## Bring the column of largest norm into the pivot position 1738 rmax = numpy.max(rdiag[j:]) 1739 kmax = (nonzero(rdiag[j:] == rmax))[0] 1740 ct = len(kmax) 1741 kmax = kmax + j 1742 if ct > 0: 1743 kmax = kmax[0] 1744 1745 ## Exchange rows via the pivot only. Avoid actually exchanging 1746 ## the rows, in case there is lots of memory transfer. The 1747 ## exchange occurs later, within the body of MPFIT, after the 1748 ## extraneous columns of the matrix have been shed. 1749 if kmax != j: 1750 temp = ipvt[j] ; ipvt[j] = ipvt[kmax] ; ipvt[kmax] = temp 1751 rdiag[kmax] = rdiag[j] 1752 wa[kmax] = wa[j] 1753 1754 ## Compute the householder transformation to reduce the jth 1755 ## column of A to a multiple of the jth unit vector 1756 lj = ipvt[j] 1757 ajj = a[j:,lj] 1758 ajnorm = self.enorm(ajj) 1759 if ajnorm == 0: break 1760 if a[j,lj] < 0: ajnorm = -ajnorm 1761 1762 ajj = ajj / ajnorm 1763 ajj[0] = ajj[0] + 1 1764 ## *** Note optimization a(j:*,j) 1765 a[j:,lj] = ajj 1766 1767 ## Apply the transformation to the remaining columns 1768 ## and update the norms 1769 1770 ## NOTE to SELF: tried to optimize this by removing the loop, 1771 ## but it actually got slower. Reverted to "for" loop to keep 1772 ## it simple. 1773 if (j+1 < n): 1774 for k in range(j+1, n): 1775 lk = ipvt[k] 1776 ajk = a[j:,lk] 1777 ## *** Note optimization a(j:*,lk) 1778 ## (corrected 20 Jul 2000) 1779 if a[j,lj] != 0: 1780 a[j:,lk] = ajk - ajj * sum(ajk*ajj)/a[j,lj] 1781 if ((pivot != 0) and (rdiag[k] != 0)): 1782 temp = a[j,lk]/rdiag[k] 1783 rdiag[k] = rdiag[k] * sqrt(numpy.max([(1.-temp**2), 0.])) 1784 temp = rdiag[k]/wa[k] 1785 if ((0.05*temp*temp) <= machep): 1786 rdiag[k] = self.enorm(a[j+1:,lk]) 1787 wa[k] = rdiag[k] 1788 rdiag[j] = -ajnorm 1789 return([a, ipvt, rdiag, acnorm])
1790 1791 1792 # Original FORTRAN documentation 1793 # ********** 1794 # 1795 # subroutine qrsolv 1796 # 1797 # given an m by n matrix a, an n by n diagonal matrix d, 1798 # and an m-vector b, the problem is to determine an x which 1799 # solves the system 1800 # 1801 # a*x = b , d*x = 0 , 1802 # 1803 # in the least squares sense. 1804 # 1805 # this subroutine completes the solution of the problem 1806 # if it is provided with the necessary information from the 1807 # factorization, with column pivoting, of a. that is, if 1808 # a*p = q*r, where p is a permutation matrix, q has orthogonal 1809 # columns, and r is an upper triangular matrix with diagonal 1810 # elements of nonincreasing magnitude, then qrsolv expects 1811 # the full upper triangle of r, the permutation matrix p, 1812 # and the first n components of (q transpose)*b. the system 1813 # a*x = b, d*x = 0, is then equivalent to 1814 # 1815 # t t 1816 # r*z = q *b , p *d*p*z = 0 , 1817 # 1818 # where x = p*z. if this system does not have full rank, 1819 # then a least squares solution is obtained. on output qrsolv 1820 # also provides an upper triangular matrix s such that 1821 # 1822 # t t t 1823 # p *(a *a + d*d)*p = s *s . 1824 # 1825 # s is computed within qrsolv and may be of separate interest. 1826 # 1827 # the subroutine statement is 1828 # 1829 # subroutine qrsolv(n,r,ldr,ipvt,diag,qtb,x,sdiag,wa) 1830 # 1831 # where 1832 # 1833 # n is a positive integer input variable set to the order of r. 1834 # 1835 # r is an n by n array. on input the full upper triangle 1836 # must contain the full upper triangle of the matrix r. 1837 # on output the full upper triangle is unaltered, and the 1838 # strict lower triangle contains the strict upper triangle 1839 # (transposed) of the upper triangular matrix s. 1840 # 1841 # ldr is a positive integer input variable not less than n 1842 # which specifies the leading dimension of the array r. 1843 # 1844 # ipvt is an integer input array of length n which defines the 1845 # permutation matrix p such that a*p = q*r. column j of p 1846 # is column ipvt(j) of the identity matrix. 1847 # 1848 # diag is an input array of length n which must contain the 1849 # diagonal elements of the matrix d. 1850 # 1851 # qtb is an input array of length n which must contain the first 1852 # n elements of the vector (q transpose)*b. 1853 # 1854 # x is an output array of length n which contains the least 1855 # squares solution of the system a*x = b, d*x = 0. 1856 # 1857 # sdiag is an output array of length n which contains the 1858 # diagonal elements of the upper triangular matrix s. 1859 # 1860 # wa is a work array of length n. 1861 # 1862 # subprograms called 1863 # 1864 # fortran-supplied ... dabs,dsqrt 1865 # 1866 # argonne national laboratory. minpack project. march 1980. 1867 # burton s. garbow, kenneth e. hillstrom, jorge j. more 1868 # 1869
1870 - def qrsolv(self, r, ipvt, diag, qtb, sdiag):
1871 if (self.debug): print 'Entering qrsolv...' 1872 sz = shape(r) 1873 m = sz[0] 1874 n = sz[1] 1875 1876 ## copy r and (q transpose)*b to preserve input and initialize s. 1877 ## in particular, save the diagonal elements of r in x. 1878 1879 for j in range(n): 1880 r[j:n,j] = r[j,j:n] 1881 x = diagonal(r) 1882 wa = qtb.copy() 1883 1884 ## Eliminate the diagonal matrix d using a givens rotation 1885 for j in range(n): 1886 l = ipvt[j] 1887 if (diag[l] == 0): break 1888 sdiag[j:] = 0 1889 sdiag[j] = diag[l] 1890 1891 ## The transformations to eliminate the row of d modify only a 1892 ## single element of (q transpose)*b beyond the first n, which 1893 ## is initially zero. 1894 1895 qtbpj = 0. 1896 for k in range(j,n): 1897 if (sdiag[k] == 0): break 1898 if (abs(r[k,k]) < abs(sdiag[k])): 1899 cotan = r[k,k]/sdiag[k] 1900 sine = 0.5/sqrt(.25 + .25*cotan*cotan) 1901 cosine = sine*cotan 1902 else: 1903 tang = sdiag[k]/r[k,k] 1904 cosine = 0.5/sqrt(.25 + .25*tang*tang) 1905 sine = cosine*tang 1906 1907 ## Compute the modified diagonal element of r and the 1908 ## modified element of ((q transpose)*b,0). 1909 r[k,k] = cosine*r[k,k] + sine*sdiag[k] 1910 temp = cosine*wa[k] + sine*qtbpj 1911 qtbpj = -sine*wa[k] + cosine*qtbpj 1912 wa[k] = temp 1913 1914 ## Accumulate the transformation in the row of s 1915 if (n > k+1): 1916 temp = cosine*r[k+1:n,k] + sine*sdiag[k+1:n] 1917 sdiag[k+1:n] = -sine*r[k+1:n,k] + cosine*sdiag[k+1:n] 1918 r[k+1:n,k] = temp 1919 sdiag[j] = r[j,j] 1920 r[j,j] = x[j] 1921 1922 ## Solve the triangular system for z. If the system is singular 1923 ## then obtain a least squares solution 1924 nsing = n 1925 wh = (nonzero(sdiag == 0))[0] 1926 if (len(wh) > 0): 1927 nsing = wh[0] 1928 wa[nsing:] = 0 1929 1930 if (nsing >= 1): 1931 wa[nsing-1] = wa[nsing-1]/sdiag[nsing-1] ## Degenerate case 1932 ## *** Reverse loop *** 1933 for j in range(nsing-2,-1,-1): 1934 sum0 = sum(r[j+1:nsing,j]*wa[j+1:nsing]) 1935 wa[j] = (wa[j]-sum0)/sdiag[j] 1936 1937 ## Permute the components of z back to components of x 1938 put(x, ipvt, wa) 1939 return(r, x, sdiag)
1940 1941 1942 1943 1944 # Original FORTRAN documentation 1945 # 1946 # subroutine lmpar 1947 # 1948 # given an m by n matrix a, an n by n nonsingular diagonal 1949 # matrix d, an m-vector b, and a positive number delta, 1950 # the problem is to determine a value for the parameter 1951 # par such that if x solves the system 1952 # 1953 # a*x = b , sqrt(par)*d*x = 0 , 1954 # 1955 # in the least squares sense, and dxnorm is the euclidean 1956 # norm of d*x, then either par is zero and 1957 # 1958 # (dxnorm-delta) .le. 0.1*delta , 1959 # 1960 # or par is positive and 1961 # 1962 # abs(dxnorm-delta) .le. 0.1*delta . 1963 # 1964 # this subroutine completes the solution of the problem 1965 # if it is provided with the necessary information from the 1966 # qr factorization, with column pivoting, of a. that is, if 1967 # a*p = q*r, where p is a permutation matrix, q has orthogonal 1968 # columns, and r is an upper triangular matrix with diagonal 1969 # elements of nonincreasing magnitude, then lmpar expects 1970 # the full upper triangle of r, the permutation matrix p, 1971 # and the first n components of (q transpose)*b. on output 1972 # lmpar also provides an upper triangular matrix s such that 1973 # 1974 # t t t 1975 # p *(a *a + par*d*d)*p = s *s . 1976 # 1977 # s is employed within lmpar and may be of separate interest. 1978 # 1979 # only a few iterations are generally needed for convergence 1980 # of the algorithm. if, however, the limit of 10 iterations 1981 # is reached, then the output par will contain the best 1982 # value obtained so far. 1983 # 1984 # the subroutine statement is 1985 # 1986 # subroutine lmpar(n,r,ldr,ipvt,diag,qtb,delta,par,x,sdiag, 1987 # wa1,wa2) 1988 # 1989 # where 1990 # 1991 # n is a positive integer input variable set to the order of r. 1992 # 1993 # r is an n by n array. on input the full upper triangle 1994 # must contain the full upper triangle of the matrix r. 1995 # on output the full upper triangle is unaltered, and the 1996 # strict lower triangle contains the strict upper triangle 1997 # (transposed) of the upper triangular matrix s. 1998 # 1999 # ldr is a positive integer input variable not less than n 2000 # which specifies the leading dimension of the array r. 2001 # 2002 # ipvt is an integer input array of length n which defines the 2003 # permutation matrix p such that a*p = q*r. column j of p 2004 # is column ipvt(j) of the identity matrix. 2005 # 2006 # diag is an input array of length n which must contain the 2007 # diagonal elements of the matrix d. 2008 # 2009 # qtb is an input array of length n which must contain the first 2010 # n elements of the vector (q transpose)*b. 2011 # 2012 # delta is a positive input variable which specifies an upper 2013 # bound on the euclidean norm of d*x. 2014 # 2015 # par is a nonnegative variable. on input par contains an 2016 # initial estimate of the levenberg-marquardt parameter. 2017 # on output par contains the final estimate. 2018 # 2019 # x is an output array of length n which contains the least 2020 # squares solution of the system a*x = b, sqrt(par)*d*x = 0, 2021 # for the output par. 2022 # 2023 # sdiag is an output array of length n which contains the 2024 # diagonal elements of the upper triangular matrix s. 2025 # 2026 # wa1 and wa2 are work arrays of length n. 2027 # 2028 # subprograms called 2029 # 2030 # minpack-supplied ... dpmpar,enorm,qrsolv 2031 # 2032 # fortran-supplied ... dabs,dmax1,dmin1,dsqrt 2033 # 2034 # argonne national laboratory. minpack project. march 1980. 2035 # burton s. garbow, kenneth e. hillstrom, jorge j. more 2036 # 2037
2038 - def lmpar(self, r, ipvt, diag, qtb, delta, x, sdiag, par=None):
2039 2040 if (self.debug): print 'Entering lmpar...' 2041 dwarf = self.machar.minnum 2042 machep = self.machar.machep 2043 sz = shape(r) 2044 m = sz[0] 2045 n = sz[1] 2046 2047 ## Compute and store in x the gauss-newton direction. If the 2048 ## jacobian is rank-deficient, obtain a least-squares solution 2049 nsing = n 2050 wa1 = qtb.copy() 2051 rthresh = numpy.max(numpy.abs(diagonal(r))) * machep 2052 wh = (nonzero(numpy.abs(diagonal(r)) < rthresh))[0] #patched from lthresh which doesn't exist WR 2053 if len(wh) > 0: 2054 nsing = wh[0] 2055 wa1[wh[0]:] = 0 2056 if nsing >= 1: 2057 ## *** Reverse loop *** 2058 for j in range(nsing-1,-1,-1): 2059 wa1[j] = wa1[j]/r[j,j] 2060 if (j-1 >= 0): 2061 wa1[0:j] = wa1[0:j] - r[0:j,j]*wa1[j] 2062 2063 ## Note: ipvt here is a permutation array 2064 put(x, ipvt, wa1) 2065 2066 ## Initialize the iteration counter. Evaluate the function at the 2067 ## origin, and test for acceptance of the gauss-newton direction 2068 iter = 0 2069 wa2 = diag * x 2070 dxnorm = self.enorm(wa2) 2071 fp = dxnorm - delta 2072 if (fp <= 0.1*delta): 2073 return[r, 0., x, sdiag] 2074 2075 ## If the jacobian is not rank deficient, the newton step provides a 2076 ## lower bound, parl, for the zero of the function. Otherwise set 2077 ## this bound to zero. 2078 2079 parl = 0. 2080 if nsing >= n: 2081 wa1 = take(diag, ipvt)*take(wa2, ipvt)/dxnorm 2082 wa1[0] = wa1[0] / r[0,0] ## Degenerate case 2083 for j in range(1,n): ## Note "1" here, not zero 2084 sum0 = sum(r[0:j,j]*wa1[0:j]) 2085 wa1[j] = (wa1[j] - sum0)/r[j,j] 2086 2087 temp = self.enorm(wa1) 2088 parl = ((fp/delta)/temp)/temp 2089 2090 ## Calculate an upper bound, paru, for the zero of the function 2091 for j in range(n): 2092 sum0 = sum(r[0:j+1,j]*qtb[0:j+1]) 2093 wa1[j] = sum0/diag[ipvt[j]] 2094 gnorm = self.enorm(wa1) 2095 paru = gnorm/delta 2096 if paru == 0: paru = dwarf/numpy.min([delta,0.1]) 2097 2098 ## If the input par lies outside of the interval (parl,paru), set 2099 ## par to the closer endpoint 2100 2101 par = numpy.max([par,parl]) 2102 par = numpy.min([par,paru]) 2103 if par == 0: par = gnorm/dxnorm 2104 2105 ## Beginning of an interation 2106 while(1): 2107 iter = iter + 1 2108 2109 ## Evaluate the function at the current value of par 2110 if par == 0: par = numpy.max([dwarf, paru*0.001]) 2111 temp = sqrt(par) 2112 wa1 = temp * diag 2113 [r, x, sdiag] = self.qrsolv(r, ipvt, wa1, qtb, sdiag) 2114 wa2 = diag*x 2115 dxnorm = self.enorm(wa2) 2116 temp = fp 2117 fp = dxnorm - delta 2118 2119 if ((abs(fp) <= 0.1*delta) or 2120 ((parl == 0) and (fp <= temp) and (temp < 0)) or 2121 (iter == 10)): break; 2122 2123 ## Compute the newton correction 2124 wa1 = take(diag, ipvt)*take(wa2, ipvt)/dxnorm 2125 2126 for j in range(n-1): 2127 wa1[j] = wa1[j]/sdiag[j] 2128 wa1[j+1:n] = wa1[j+1:n] - r[j+1:n,j]*wa1[j] 2129 wa1[n-1] = wa1[n-1]/sdiag[n-1] ## Degenerate case 2130 2131 temp = self.enorm(wa1) 2132 parc = ((fp/delta)/temp)/temp 2133 2134 ## Depending on the sign of the function, update parl or paru 2135 if fp > 0: parl = numpy.max([parl,par]) 2136 if fp < 0: paru = numpy.min([paru,par]) 2137 2138 ## Compute an improved estimate for par 2139 par = numpy.max([parl, par+parc]) 2140 2141 ## End of an iteration 2142 2143 ## Termination 2144 return[r, par, x, sdiag]
2145 2146 2147 ## Procedure to tie one parameter to another.
2148 - def tie(self, p, ptied=None):
2149 if (self.debug): print 'Entering tie...' 2150 if (ptied == None): return 2151 for i in range(len(ptied)): 2152 if ptied[i] == '': continue 2153 cmd = 'p[' + str(i) + '] = ' + ptied[i] 2154 exec(cmd) 2155 return(p)
2156 2157 2158 # Original FORTRAN documentation 2159 # ********** 2160 # 2161 # subroutine covar 2162 # 2163 # given an m by n matrix a, the problem is to determine 2164 # the covariance matrix corresponding to a, defined as 2165 # 2166 # t 2167 # inverse(a *a) . 2168 # 2169 # this subroutine completes the solution of the problem 2170 # if it is provided with the necessary information from the 2171 # qr factorization, with column pivoting, of a. that is, if 2172 # a*p = q*r, where p is a permutation matrix, q has orthogonal 2173 # columns, and r is an upper triangular matrix with diagonal 2174 # elements of nonincreasing magnitude, then covar expects 2175 # the full upper triangle of r and the permutation matrix p. 2176 # the covariance matrix is then computed as 2177 # 2178 # t t 2179 # p*inverse(r *r)*p . 2180 # 2181 # if a is nearly rank deficient, it may be desirable to compute 2182 # the covariance matrix corresponding to the linearly independent 2183 # columns of a. to define the numerical rank of a, covar uses 2184 # the tolerance tol. if l is the largest integer such that 2185 # 2186 # abs(r(l,l)) .gt. tol*abs(r(1,1)) , 2187 # 2188 # then covar computes the covariance matrix corresponding to 2189 # the first l columns of r. for k greater than l, column 2190 # and row ipvt(k) of the covariance matrix are set to zero. 2191 # 2192 # the subroutine statement is 2193 # 2194 # subroutine covar(n,r,ldr,ipvt,tol,wa) 2195 # 2196 # where 2197 # 2198 # n is a positive integer input variable set to the order of r. 2199 # 2200 # r is an n by n array. on input the full upper triangle must 2201 # contain the full upper triangle of the matrix r. on output 2202 # r contains the square symmetric covariance matrix. 2203 # 2204 # ldr is a positive integer input variable not less than n 2205 # which specifies the leading dimension of the array r. 2206 # 2207 # ipvt is an integer input array of length n which defines the 2208 # permutation matrix p such that a*p = q*r. column j of p 2209 # is column ipvt(j) of the identity matrix. 2210 # 2211 # tol is a nonnegative input variable used to define the 2212 # numerical rank of a in the manner described above. 2213 # 2214 # wa is a work array of length n. 2215 # 2216 # subprograms called 2217 # 2218 # fortran-supplied ... dabs 2219 # 2220 # argonne national laboratory. minpack project. august 1980. 2221 # burton s. garbow, kenneth e. hillstrom, jorge j. more 2222 # 2223 # ********** 2224
2225 - def calc_covar(self, rr, ipvt=None, tol=1.e-14):
2226 2227 if (self.debug): print 'Entering calc_covar...' 2228 if rank(rr) != 2: 2229 print 'ERROR: r must be a two-dimensional matrix' 2230 return(-1) 2231 s = shape(rr) 2232 n = s[0] 2233 if s[0] != s[1]: 2234 print 'ERROR: r must be a square matrix' 2235 return(-1) 2236 2237 if (ipvt == None): ipvt = arange(n) 2238 r = rr.copy() 2239 r.shape = [n,n] 2240 2241 ## For the inverse of r in the full upper triangle of r 2242 l = -1 2243 tolr = tol * abs(r[0,0]) 2244 for k in range(n): 2245 if (abs(r[k,k]) <= tolr): break 2246 r[k,k] = 1./r[k,k] 2247 for j in range(k): 2248 temp = r[k,k] * r[j,k] 2249 r[j,k] = 0. 2250 r[0:j+1,k] = r[0:j+1,k] - temp*r[0:j+1,j] 2251 l = k 2252 2253 ## Form the full upper triangle of the inverse of (r transpose)*r 2254 ## in the full upper triangle of r 2255 if l >= 0: 2256 for k in range(l+1): 2257 for j in range(k): 2258 temp = r[j,k] 2259 r[0:j+1,j] = r[0:j+1,j] + temp*r[0:j+1,k] 2260 temp = r[k,k] 2261 r[0:k+1,k] = temp * r[0:k+1,k] 2262 2263 ## For the full lower triangle of the covariance matrix 2264 ## in the strict lower triangle or and in wa 2265 wa = repeat([r[0,0]], n) 2266 for j in range(n): 2267 jj = ipvt[j] 2268 sing = j > l 2269 for i in range(j+1): 2270 if sing: r[i,j] = 0. 2271 ii = ipvt[i] 2272 if ii > jj: r[ii,jj] = r[i,j] 2273 if ii < jj: r[jj,ii] = r[i,j] 2274 wa[jj] = r[j,j] 2275 2276 ## Symmetrize the covariance matrix in r 2277 for j in range(n): 2278 r[0:j+1,j] = r[j,0:j+1] 2279 r[j,j] = wa[j] 2280 2281 return(r)
2282
2283 -class machar:
2284 - def __init__(self, double=1):
2285 if (double == 0): 2286 self.machep = 1.19209e-007 2287 self.maxnum = 3.40282e+038 2288 self.minnum = 1.17549e-038 2289 self.maxgam = 171.624376956302725 2290 else: 2291 self.machep = 2.2204460e-016 2292 self.maxnum = 1.7976931e+308 2293 self.minnum = 2.2250739e-308 2294 self.maxgam = 171.624376956302725 2295 2296 self.maxlog = log(self.maxnum) 2297 self.minlog = log(self.minnum) 2298 self.rdwarf = sqrt(self.minnum*1.5) * 10 2299 self.rgiant = sqrt(self.maxnum) * 0.1
2300