This is a compilation of answers to questions frequently asked by levmar users.


  1. Q1 -- What is levmar?
  2. Q2 -- What types of minimization problems can levmar solve?
  3. Q3 -- Where can I find more information on nonlinear regression/data fitting?
  4. Q4 -- Who is using levmar?
  5. Q5 -- What do I need to use levmar?
  6. Q6 -- How can I compile levmar?
  7. Q7 -- Where are the arguments of the various routines explained in detail?
  8. Q8 -- Are matrices stored in row or column major order?
  9. Q9 -- How can I pass my own data for use by the objective function and its Jacobian?
  10. Q10 -- Can I use levmar for problems whose Jacobian is unavailable or difficult/expensive to compute?
  11. Q11 -- How can I obtain an analytic Jacobian for my minimization problem??
  12. Q12 -- How can I verify the correctness of a Jacobian?
  13. Q13 -- How should I choose the starting point for a minimization?
  14. Q14 -- How can I force an early termination of the minimization?
  15. Q15 -- How can I deal with badly scaled problems?
  16. Q16 -- How can I compile just the single or double variants of levmar's functions?
  17. Q17 -- How can I force the linear solvers free their allocated memory before returning?
  18. Q18 -- What is the covar matrix computed by the various routines?
  19. Q19 -- How can I solve a weighted non-linear least squares problem using levmar?
  20. Q20 -- How can I solve a non-linear least squares problem with some fixed parameters?
  21. Q21 -- How can I solve a problem involving box constraints with infinite endpoints?
  22. Q22 -- What if my constrained optimization problem has an infeasible starting point?
  23. Q23 -- Where can I find an example of using levmar for curve fitting?
  24. Q24 -- Where can I find an example of a real-world problem using levmar?
  25. Q25 -- Why am I getting syntax errors in lmaccess.h with Visual Studio?
  26. Q26 -- Why am I getting link errors for the f2c library?
  27. Q27 -- Why am I getting a link error for MAIN__?
  28. Q28 -- Why am I getting a link error for sgemm/dgemm?
  29. Q29 -- How can I determine what kind of BLAS library I have installed?
  30. Q30 -- Can I use levmar from matlab?
  31. Q31 -- How can I compile levmar's MEX-file?
  32. Q32 -- Can levmar benefit from a multi-core architecture?
  33. Q33 -- Can I use levmar from perl?
  34. Q34 -- Can I use levmar from python?
  35. Q35 -- Can I use levmar from haskell?
  36. Q36 -- Can I use levmar from tcl?
  37. Q37 -- Can I include levmar in a commercial product?
  38. Q38 -- What are the licensing terms of LAPACK/BLAS?

  1. Q1 -- What is levmar? -- [top]

    levmar is a C/C++ implementation of the Levenberg-Marquardt nonlinear least squares minimization algorithm that is distributed under the GNU General Public License.

  2. Q2 -- What types of minimization problems can levmar solve? -- [top]

    levmar supports unconstrained nonlinear least squares minimization with analytic or approximate Jacobians. Such problems typically arise from the solution of systems of nonlinear equations or in data fitting applications. levmar has also some support for constrained nonlinear least squares, allowing linear equation, linear inequality and box inequality constraints to be imposed.

  3. Q3 -- Where can I find more information on nonlinear regression/data fitting? -- [top]

    Nonlinear regression, i.e. the problem of fitting a nonlinear model to your data, is far too broad to be covered here. A useful introduction can be found in these notes from TU Denmark. Another free source is http://www.curvefit.com/; yet another option is to consult sections 15.4 and 15.5 of the Numerical Recipes book.

  4. Q4 -- Who is using levmar? -- [top]

    levmar is being used by thousands of individuals to solve nonlinear regression problems in various disciplines. It has also been integrated into open source and commercial products, the list of which currently includes:

    1. Hugin, a panorama photo stitcher
    2. SciGraphica, a package for scientific data analysis and technical graphics
    3. Tarquin, an analysis tool for NMR spectroscopic data
    4. Teem, a group of libraries for representing, processing, and visualizing scientific raster data
    5. CellML, a language for describing mathematical models
    6. Astrometry.net, an astrometric calibraton service
    7. Stimfit, a software for the analysis of electrophysiological data
    8. StochFit, a software for modeling specular x-ray reflectivity or neutron reflectivity data.
    9. Autopano, a program for creating image panoramas

    Additionally, several commercial license agreements have been signed with private companies, allowing them to use levmar in their products.

    levmar has also been employed in a number of research papers, check Google Scholar for an incomplete list.

  5. Q5 -- What do I need to use levmar? -- [top]

    To be able to use all levmar functions, you should have installed LAPACK or an equivalent vendor library. On Intel processors, MKL is highly recommended. Netlib has an f2c'ed free version at http://www.netlib.org/clapack. levmar users have also reported that the precompiled Windows libraries at Netlib are broken and suggest to build them yourself. Alternatively, precompiled Windows LAPACK/BLAS libraries are offered here, here and here. Precompiled BLAS Windows libraries are provided by OpenBLAS.
    It is also noted that a subset of levmar's functionality is available even without LAPACK. The use (or not) of LAPACK is controlled by the HAVE_LAPACK macro definition in levmar.h.

  6. Q6 -- How can I compile levmar? -- [top]

    First, you have to decide if you want LAPACK or not. In case that you plan to use it, you should make sure that it is installed at your site. The second step is to compile levmar itself. The tarfile contains makefiles for Unix/Linux using gcc (Makefile) or Intel's icc (Makefile.icc) and Windows using Visual Studio/Visual Studio Express (Makefile.vc). Please read the comments within those files for more information.
    There is also the option of not using any of the supplied Makefiles but rather generate platform-specific compilation commands using the supplied CMakeLists.txt file with the CMake build system; please refer to CMake's online documentation for more information on how to use it.

  7. Q7 -- Where are the arguments of the various routines explained in detail? -- [top]

    The source code of each function is preceeded by detailed comments explaining the role of each argument. See the *_core.c files.

  8. Q8 -- Are matrices stored in row or column major order? -- [top]

    Matrices are stored using the normal C convention, i.e. in row major order.

  9. Q9 -- How can I pass my own data for use by the objective function and its Jacobian? -- [top]

    Assume that you want to pass two arrays, one of doubles and one of integers. A quick but not very elegant solution is to use global variables. The simplest way to avoid using globals is to declare a structure as

    struct mydata{
       double dar[XXX];
       int iar[YYY];
    };
    
    where XXX and YYY denote the appropriate array sizes. Then, define a structure variable with
    struct mydata data;
    
    and fill it:
    data.dar[0]=7.0;
    data.iar[0]=-17;
    // etc
    
    Following this, call the appropriate LM routine passing it the address of data as the adata argument, e.g.
    ret=dlevmar_der(func, fjac, ..., (void *)&data);
    
    Your func and fjac routines should interpret the supplied data using type casting:
    struct mydata *dptr;
    
    dptr=(struct mydata *)adata; // adata is passed as void *
    
    // supplied data can now be accessed as dptr->dar[0], etc
    
  10. Q10 -- Can I use levmar for problems whose Jacobian is unavailable or difficult/expensive to compute? -- [top]

    Yes. levmar includes routines that approximate the Jacobian using a combination of finite differences and secant updates. These routines have the suffix _dif, as opposed to those using analytic Jacobians that end in _der. Typically, minimizing with approximate Jacobians takes longer to converge compared to when the Jacobian is supplied analytically (see also Q11). However, cases where the Jacobian is very expensive to compute are exceptions to this rule. In such situations, better performance can be achieved by employing the _dif routines, even if the analytic expression for the Jacobian is known.

  11. Q11 -- How can I obtain an analytic Jacobian for my minimization problem? -- [top]

    Employing analytic Jacobians is generally recommended due to faster convergence. When the objective function is too complicated for hand-coding its Jacobian, the latter can be generated through the use of automatic differentiation tools such as those listed here, that work by systematically applying the chain rule to a given code segment.

  12. Q12 -- How can I verify the correctness of a Jacobian? -- [top]

    levmar includes the functions dlevmar_chkjac() and slevmar_chkjac() which allow you to check the Jacobian of a function for consistency with the function itself. More details can be found in the comments preceeding the source code.

  13. Q13 -- How should I choose the starting point for a minimization? -- [top]

    Selecting the starting point for a particular minimization problem is largely problem-specific. In some cases, your background knowledge of the problem may guide you to select reasonable initial values for the problem's parameters. In other cases, you might be able to define a simplified linear version of the non-linear problem under consideration, whose solution will serve as the starting point for the non-linear problem. You might even use random values for your starting point. Beware, however, that the success of this approach is highly dependent on the nature of your problem: If it is highly non-linear or it has several local minima, it is very likely that you will end up with unexpected surprises.

  14. Q14 -- How can I force an early termination of the minimization? -- [top]

    In certain applications, a minimization problem can be considered to be satisfactory solved as soon as the least squares error becomes smaller than some value. Despite that more iterations can reduce this error further, the reduction is not so critical for the problem at hand and might take more computation time to be achieved. In such cases, by setting opts[3] to a suitable value, the user can instruct levmar to return (with code 6) as soon as the least squares error becomes sufficiently small.

  15. Q15 -- How can I deal with badly scaled problems? -- [top]

    Practical applications of minimization often involve variables that vary greatly in magnitude. Such minimization problems are characterized as being badly scaled and, unless adjusted, may lead to poor convergence. levmar does not include an explicit mechanism for dealing with badly scaled problems since such mechanisms are, to a large extent, problem-specific. Nevertheless, the user can account for badly scaled problems by using her knowledge of the problem domain to rescale the independent variables, i.e. change their units so that they all have approximately the same range. See chapter 7 in Dennis-Schnabel for more details.

  16. Q16 -- How can I compile just the single or double variants of levmar's functions? -- [top]

    levmar offers single and double precision variants of all optimization functions it includes. Starting with version 2.3, it supports the option of compiling only the single or only the double precision variants. This results in a library of smaller size. The choice of compiling the single or double precision variants is controlled by defining the macros LM_SNGL_PREC and LM_DBL_PREC, respectively. By default, both macros are defined.

  17. Q17 -- How can I force the linear solvers free their allocated memory before returning? -- [top]

    The default behavior of the linear solver routines is to retain memory between calls, aiming to increase performance by avoiding excessive free's and malloc's. If this is undesirable, it can be turned off by recompiling levmar without the LINSOLVERS_RETAIN_MEMORY macro definition.

  18. Q18 -- What is the covar matrix computed by the various routines? -- [top]

    Upon return, all routines fill covar with the variance-covariance matrix of the best-fit parameters, which estimates the statistical error on the latter: covar's diagonal (off-diagonal) elements are estimates of the variances (covariances) of the best-fit parameters. More information can be found at the background section of this note.

  19. Q19 -- How can I solve a weighted non-linear least squares problem using levmar? -- [top]

    Currently levmar does not directly offer support for weighted non-linear least squares problems. However, it is quite straightforward to program this for a particular problem, as follows: Let us assume that you want to minimize ||e||_C = e^T*C*e with e being your error vector (e=x-f(p)) and C a symmetric positive-definite weight matrix. Then, you can equivalently minimize s^T*s without weights, where s=W*e=W*x - W*f(p) and W is the Cholesky decomposition of C, i.e. C=W^T*W. This is because s^T*s=e^t*W^t*W*e=e^T*C*e=||e||_C. In the case that all non-diagonal elements of C are zero, then W is simply the diagonal matrix resulting after taking the sqrt of C's diagonal elements. See also the discussion here.
    Functions dlevmar_chol()/slevmar_chol() included in levmar can be employed to compute the Cholesky decomposition of a symmetric positive-definite matrix.

  20. Q20 -- How can I solve a non-linear least squares problem with some fixed parameters? -- [top]

    Occasionally, it is required to solve problems where certain parameters are held to fixed values and the objective function is minimized with respect to the remaining ones only. This can be achieved in the following two ways:

    • The fixed parameters are not included in the parameter vector. The values of the fixed parameters can be communicated to the objective function as user-supplied arguments through the adata mechanism (see Q9) or in global variables.
    • Columns of the Jacobian which correspond to the fixed parameters are explicitly set to zero.
    Although it requires slightly more programming, the first option is the cleanest one.

  21. Q21 -- How can I solve a problem involving box constraints with infinite endpoints? -- [top]

    To deal with non-linear least squares problems for which no lower bound inequality constraints apply to certain variables, the corresponding entries in the lb argument of dlevmar_bc_XXX/slevmar_bc_XXX or dlevmar_blec_XXX/slevmar_blec_XXX should be set to -DBL_MAX/-FLT_MAX (i.e., minus infinity), respectively. Similarly, when there are no upper bound constraints, the entries in the ub argument should be set to DBL_MAX/FLT_MAX.

  22. Q22 -- What if my constrained optimization problem has an infeasible starting point? -- [top]

    It is strongly advised that the starting point to a constrained optimization problem is feasible, i.e. satisfies any equality, inequality and box constraints that should apply. In the opposite case, levmar maps an infeasible point to a feasible one but does so in a very crude manner that might move the resulting feasible point far from the local minimum.

  23. Q23 -- Where can I find an example of using levmar for curve fitting? -- [top]

    An example showing how to use levmar for fitting a three-parameter exponential model to a set of data measurements is in expfit.c; refer to the comments annotating the code for more details.

  24. Q24 -- Where can I find an example of a real-world problem using levmar? -- [top]

    An example of a real-word use of levmar is the homest library which is concerned with estimating a homography, that is a plane projective transformation that arises in computer vision.

  25. Q25 -- Why am I getting syntax errors in lmaccess.h with Visual Studio? -- [top]

    This is most likely to occur when compiling an application that makes use of levmar version 2.4 or earlier under Visual Studio. The problem is caused by an erroneous include path option passed to the compiler, which results in the inadvertent inclusion of an irrelevant header file. Specifically, Visual Studio comes with a header file also called lm.h. If a directive such as #include <lm.h>  is included in a source file while the /I switch to the compiler does not specify the directory where levmar's lm.h resides, then Visual Studio's lm.h gets included resulting in compilation errors. The solution is to ensure that an appropriate include path is specified to the compiler with the /I switch. Starting with levmar version 2.5, lm.h was renamed to levmar.h to avoid such conflicts. For backwards compatibility, lm.h has been retained as a symbolic link to levmar.h. However, its use in new code is discouraged.

  26. Q26 -- Why am I getting link errors for the f2c library? -- [top]

    The f2c library is needed for f2c'ed versions of LAPACK/BLAS. On some systems it is equivalent to the F77 and I77 libraries, whereas on others it is not needed at all. If your linker complains about a missing f2c library, edit CMakeLists.txt and try substituting it with F77 and I77 or remove it altogether.

  27. Q27 -- Why am I getting a link error for MAIN__? -- [top]

    When building lmdemo, certain LAPACK/BLAS libraries cause the linker to fail with an error like undefined reference to `MAIN__', or similar. This is due to these libraries superfluously calling a Fortran main() routine for initialization. There are several workarounds to this problem:

    • If using gcc, pass the option -u MAIN__ to the linker.
    • If compiling under Windows, link against these precompiled libraries.
    • Define a dummy MAIN__ function as int MAIN__(int argc, char **argv) { abort(); return 1; }

    See also the discussion here.

  28. Q28 -- Why am I getting a link error for sgemm/dgemm? -- [top]

    This is a BLAS name mangling problem. Make sure that the definition of LM_BLAS_PREFIX in misc.h matches your system:

    • Define it being empty if you have a F77 BLAS (default).
    • Define it as f2c_ if you have an f2c'ed BLAS.
    • Define it as cblas_ if you have C interface to BLAS.

    For the f2c'ed/CBLAS cases, you should also define LM_BLAS_SUFFIX as being empty in misc.h, instead of the default of being defined as an underscore. See also Q29.

  29. Q29 -- How can I determine what kind of BLAS library I have installed? -- [top]

    On unix, try nm libblas.a | grep gemm . Under Windows, use dumpbin /all blas.lib | more. Inspect the output of those commands for prefixes like f2c_ and cblas_ to the gemm routines.

  30. Q30 -- Can I use levmar from matlab? -- [top]

    Yes. Starting with version 2.2, levmar includes a MEX-file that can be used to interface it with matlab. See levmar.c in the matlab directory of levmar's source distribution. Documentation for using levmar with matlab is provided in file levmar.m. Assuming that the MEX-file has been compiled (see also Q31), file lmdemo.m contains examples of its use with matlab. Try it with the command matlab < lmdemo.m. It has been reported that the MEX-file also works with GNU octave, but I have not tested this myself; if you succeed with it, please let me know.

  31. Q31 -- How can I compile levmar's MEX-file? -- [top]

    The answer depends on your operating system/compiler. The Makefile in the matlab directory works with GCC under Linux, whereas Makefile.w32 works with MSVC under Windows. You might also need to run "mex -setup" from the matlab or OS command prompt to configure your system to build MEX-files. In any case, the provided makefiles can be used as a starting point for compiling on your own system. If building under Windows, you might also consider using gnumex.

  32. Q32 -- Can levmar benefit from a multi-core architecture? -- [top]

    Yes, in some cases. Starting with version 2.6, levmar includes support for a parallel Cholesky solver from the PLASMA linear algebra package for homogeneous multi-core architectures. This parallel solver can improve the performance of large optimization problems, i.e. problems whose number of parameters is in the order of several hundreds. Since the gain in performance is platform-specific, the user is advised to do some testing for determining the speed up achieved by the parallel solver, if any. Obviously, for small optimization problems, the overheads of parallel execution will outweigh the benefits. To use the parallel Cholesky solver, the HAVE_PLASMA macro in levmar.h should be defined prior to compilation and the invocations of AX_EQ_B_PLASMA_CHOL() uncommented in the source code.

  33. Q33 -- Can I use levmar from perl? -- [top]

    Yes. John Lapeyre has written a PDL interface to levmar.

  34. Q34 -- Can I use levmar from python? -- [top]

    Yes. Alastair Tse has written a python interface to levmar.

  35. Q35 -- Can I use levmar from haskell? -- [top]

    Yes. Bas and Roel van Dijk have written a haskell interface to levmar.

  36. Q36 -- Can I use levmar from tcl? -- [top]

    Yes. Andrey Nakin has written a tcl interface to levmar.

  37. Q37 -- Can I include levmar in a commercial product? -- [top]

    levmar is distributed under GPL, which allows you to modify the original code into a new program and then sell copies of the modified program commercially, but only under the terms of the GPL. Thus, you must make the source code available to the users of the program as described in the GPL, and they must be allowed to redistribute and modify it as described in the GPL. Clearly, such requirements are in conflict with using levmar as part of a proprietary application. However, it is possible to obtain a paid license for proprietary commercial use under terms different than those of GPL. Please contact the author for further information on such licensing. More answers to questions related to the GPL license can be found here.

  38. Q38 -- What are the licensing terms of LAPACK/BLAS? -- [top]

    Netlib's LAPACK/BLAS are freely distributable, their exact licensing terms are available here: LAPACK | BLAS. If you are using a vendor-specific flavor (see Q5), please consult its documentation for details.