Class Lmmin

  • Direct Known Subclasses:
    AlgorithmELSUNCOpt2D.FitOAR2DLMModel, AlgorithmELSUNCOpt3D.FitOAR3DLMModel

    public abstract class Lmmin
    extends java.lang.Object
    To run self tests in another module put: if (testMode) { new FitAll(); setCompleted(false); return; } class FitAll extends Lmmin { public FitAll() { super(); } public void fitToFunction(double x[], double fvec[], int info[]) { return; } }
    • Field Summary

      Fields 
      Modifier and Type Field Description
      private int BARD  
      private int BOX_3D  
      private int BROWN_ALMOST_LINEAR  
      private int BROWN_AND_DENNIS  
      private int CHEBYQUAD  
      private double DBL_EPSILON  
      private double DBL_MAX  
      private double DBL_MIN
      2**-1022 = D1MACH(1).
      protected boolean debugMatrix  
      protected boolean debugMessages  
      protected double[] diag
      diag is an array of length n.
      private int DRAPER24D  
      protected double epsfcn
      step used to calculate the jacobian epsfcn is an input variable used in choosing a step length for the forward-difference approximation.
      private int EQUILIBRIUM_COMBUSTION  
      protected double factor
      initial bound to steps in the outer loop factor is a positive input variable used in determining the initial step bound.
      private double[] fjac
      fjac is an OUTPUT m by n array.
      protected double fnorm
      norm of the residue vector fvec
      private int FREUDENSTEIN_AND_ROTH  
      protected double ftol
      ftol is a nonnegative input variable.
      private double[] fvec
      fvec is an output array of length m which contains the functions evaluated at the output x.
      protected double gtol
      gtol is a nonnegative input variable.
      private int HATFLDB  
      private int HATFLDC  
      private int HELICAL_VALLEY  
      private int HOCK1  
      private int HOCK21_MODIFIED  
      private int HOCK25  
      protected int[] info
      status of minimization info is an integer OUTPUT variable that indicates the termination status of lm_lmdif as follows: info
      private int[] ipvt
      ipvt is an integer OUTPUT array of length n.
      private int JENNRICH_AND_SAMPSON  
      private int KOWALIK_AND_OSBORNE  
      private int LEVMAR_ROSENBROCK  
      private int LINEAR_FULL_RANK  
      private int LINEAR_RANK1  
      private int LINEAR_RANK1_WITH_ZERO_COLUMNS_AND_ROWS  
      private double LM_DWARF  
      private java.lang.String[] lm_infmsg
      Set message texts (indexed by info)
      private double LM_MACHEP
      resolution of arithmetic
      private java.lang.String[] lm_shortmsg  
      private double LM_SQRT_DWARF
      Square should not underflow
      private double LM_SQRT_GIANT
      Square should not overflow
      private double LM_USERTOL  
      private int m
      m is a positive integer input variable set to the number of functions.
      protected int maxcall
      maximum number of iterations
      private int maxfev
      maxfev is a positive integer input variable.
      private int MEYER  
      protected int mode
      mode is an integer input variable.
      private int MODIFIED_ROSENBROCK  
      private int n
      n is a positive integer variable set to the number of variables; n must not exceed m.
      protected int nfev
      actual number of iterations
      private int OSBORNE1  
      private int OSBORNE2  
      private double[] par  
      private int POWELL_2_PARAMETER  
      private int POWELL_SINGULAR  
      protected int printflags
      Bits ored to show where to output diagnostics
      private double[] qtf
      qtf is an OUTPUT array of length n which contains the first n elements of the vector (q transpose)*fvec.
      private int ROSENBROCK  
      private int testCase  
      private boolean testMode  
      private double[] wa1
      wa1, wa2, and wa3 are work arrays of length n.
      private double[] wa2  
      private double[] wa3  
      private double[] wa4
      wa4 is a work array of length m, used among others to hold residuals from evaluate.
      private int WATSON  
      private int WOOD  
      protected double[] x
      x is an array of length n.
      private double[] xSeries  
      protected double xtol
      xtol is a nonnegative input variable.
      private double[] ySeries  
    • Constructor Summary

      Constructors 
      Constructor Description
      Lmmin()  
      Lmmin​(int m, int n)  
      Lmmin​(int m, int n, double[] x)  
    • Method Summary

      All Methods Instance Methods Abstract Methods Concrete Methods 
      Modifier and Type Method Description
      void driver()
      driver.
      void dumpResults()  
      private void fitTestModel()  
      abstract void fitToFunction​(double[] x, double[] fvec, int[] info)  
      private void fitToTestFunction​(double[] par, double[] fvec)
      lm_enorm.
      double[] getParameters()  
      int getStatus()  
      (package private) double lm_enorm​(int n, double[] x, int offset)  
      private void lm_lmdif()  
      private void lm_lmpar​(int n, double[] r, int ldr, int[] ipvt, double[] diag, double[] qtb, double delta, double[] par, double[] x, double[] sdiag, double[] aux, double[] xdi)  
      private void lm_printout_std​(double[] par, double[] fvec, int printflags, int iflag, int iter, int nfev)  
      private void lm_qrfac​(int m, int n, double[] a, int pivot, int[] ipvt, double[] rdiag, double[] acnorm, double[] wa)  
      private void lm_qrsolv​(int n, double[] r, int ldr, int[] ipvt, double[] diag, double[] qtb, double[] x, double[] sdiag, double[] wa)  
      private double shiftedChebyshev​(double x, int n)  
      void statusMessage​(int status)  
      • Methods inherited from class java.lang.Object

        clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
    • Field Detail

      • testMode

        private boolean testMode
      • DBL_EPSILON

        private final double DBL_EPSILON
      • DBL_MIN

        private final double DBL_MIN
        2**-1022 = D1MACH(1).
      • LM_MACHEP

        private final double LM_MACHEP
        resolution of arithmetic
      • LM_DWARF

        private final double LM_DWARF
      • LM_SQRT_DWARF

        private final double LM_SQRT_DWARF
        Square should not underflow
      • LM_SQRT_GIANT

        private final double LM_SQRT_GIANT
        Square should not overflow
      • LM_USERTOL

        private final double LM_USERTOL
      • ftol

        protected double ftol
        ftol is a nonnegative input variable. Termination occurs when both the actual and predicted relative reductions in the sum of squares are at most ftol. Therefore, ftol measures the relative error desired in the sum of squares.
      • xtol

        protected double xtol
        xtol is a nonnegative input variable. Termination occurs when the relative error between two consecutive iterates is at most xtol. Therefore, xtol measures the relative error desired in the approximate solution.
      • gtol

        protected double gtol
        gtol is a nonnegative input variable. Termination occurs when the cosine of the angle between fvec and any column of the jacobian is at most gtol in absolute value. Therefore, gtol measures the orthogonality desired between the function vector and the columns of the jacobian.
      • epsfcn

        protected double epsfcn
        step used to calculate the jacobian epsfcn is an input variable used in choosing a step length for the forward-difference approximation. The relative errors in the functions are assumed to be of the order of epsfcn.
      • factor

        protected double factor
        initial bound to steps in the outer loop factor is a positive input variable used in determining the initial step bound. This bound is set to the product of factor and the euclidean norm of diag*x if nonzero, or else to factor itself. In most cases factor should lie in the interval (0.1,100.0). Generally, the value 100.0 is recommended.
      • maxcall

        protected int maxcall
        maximum number of iterations
      • maxfev

        private int maxfev
        maxfev is a positive integer input variable. Termination occurs when the number of calls to lm_fcn is at least maxfev by the end of an iteration.
      • mode

        protected int mode
        mode is an integer input variable. If mode = 1, the variables will be scaled internally. If mode = 2, the scaling is specified by the input diag.
      • printflags

        protected int printflags
        Bits ored to show where to output diagnostics
      • lm_infmsg

        private final java.lang.String[] lm_infmsg
        Set message texts (indexed by info)
      • lm_shortmsg

        private final java.lang.String[] lm_shortmsg
      • debugMessages

        protected boolean debugMessages
      • debugMatrix

        protected boolean debugMatrix
      • m

        private int m
        m is a positive integer input variable set to the number of functions.
      • n

        private int n
        n is a positive integer variable set to the number of variables; n must not exceed m.
      • x

        protected double[] x
        x is an array of length n. On input x must contain an initial estimate of the solution vector. On output x contains the final estimate of the solution vector.
      • fnorm

        protected double fnorm
        norm of the residue vector fvec
      • nfev

        protected int nfev
        actual number of iterations
      • info

        protected int[] info
        status of minimization info is an integer OUTPUT variable that indicates the termination status of lm_lmdif as follows: info
        • fvec

          private double[] fvec
          fvec is an output array of length m which contains the functions evaluated at the output x.
        • diag

          protected double[] diag
          diag is an array of length n. If mode = 1, diag is internally set. If mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.
        • qtf

          private double[] qtf
          qtf is an OUTPUT array of length n which contains the first n elements of the vector (q transpose)*fvec.
        • fjac

          private double[] fjac
          fjac is an OUTPUT m by n array. The upper n by n submatrix of fjac contains an upper triangular matrix r with diagonal elements of nonincreasing magnitude such that pT*(jacT*jac)*p = rT*r (NOTE: T stands for matrix transposition), where p is a permutation matrix and jac is the final calculated jacobian. Column j of p is column ipvt(j) (see below) of the identity matrix. The lower trapezoidal part of fjac contains information generated during the computation of r.
        • wa1

          private double[] wa1
          wa1, wa2, and wa3 are work arrays of length n.
        • wa2

          private double[] wa2
        • wa3

          private double[] wa3
        • wa4

          private double[] wa4
          wa4 is a work array of length m, used among others to hold residuals from evaluate.
        • ipvt

          private int[] ipvt
          ipvt is an integer OUTPUT array of length n. It defines a permutation matrix p such that jac*p = q*r, where jac is the final calculated jacobian, q is orthogonal (not stored), and r is upper triangular with diagonal elements of nonincreasing magnitude. Column j of p is column ipvt(j) of the identity matrix.
        • par

          private double[] par
        • testCase

          private int testCase
        • xSeries

          private double[] xSeries
        • ySeries

          private double[] ySeries
        • LINEAR_RANK1_WITH_ZERO_COLUMNS_AND_ROWS

          private final int LINEAR_RANK1_WITH_ZERO_COLUMNS_AND_ROWS
          See Also:
          Constant Field Values
        • Constructor Detail

          • Lmmin

            public Lmmin()
          • Lmmin

            public Lmmin​(int m,
                         int n)
            Parameters:
            m - the number of functions
            n - the number of variables; n must not exceed m.
          • Lmmin

            public Lmmin​(int m,
                         int n,
                         double[] x)
            Parameters:
            m - the number of functions
            n - the number of variables; n must not exceed m.
            x -
        • Method Detail

          • statusMessage

            public void statusMessage​(int status)
          • getParameters

            public double[] getParameters()
          • fitTestModel

            private void fitTestModel()
          • dumpResults

            public void dumpResults()
          • getStatus

            public int getStatus()
          • driver

            public void driver()
            driver.
          • lm_lmdif

            private void lm_lmdif()
          • lm_qrfac

            private void lm_qrfac​(int m,
                                  int n,
                                  double[] a,
                                  int pivot,
                                  int[] ipvt,
                                  double[] rdiag,
                                  double[] acnorm,
                                  double[] wa)
          • fitToFunction

            public abstract void fitToFunction​(double[] x,
                                               double[] fvec,
                                               int[] info)
          • lm_printout_std

            private void lm_printout_std​(double[] par,
                                         double[] fvec,
                                         int printflags,
                                         int iflag,
                                         int iter,
                                         int nfev)
          • lm_lmpar

            private void lm_lmpar​(int n,
                                  double[] r,
                                  int ldr,
                                  int[] ipvt,
                                  double[] diag,
                                  double[] qtb,
                                  double delta,
                                  double[] par,
                                  double[] x,
                                  double[] sdiag,
                                  double[] aux,
                                  double[] xdi)
          • lm_qrsolv

            private void lm_qrsolv​(int n,
                                   double[] r,
                                   int ldr,
                                   int[] ipvt,
                                   double[] diag,
                                   double[] qtb,
                                   double[] x,
                                   double[] sdiag,
                                   double[] wa)
          • lm_enorm

            double lm_enorm​(int n,
                            double[] x,
                            int offset)
          • fitToTestFunction

            private void fitToTestFunction​(double[] par,
                                           double[] fvec)
            lm_enorm.
          • shiftedChebyshev

            private double shiftedChebyshev​(double x,
                                            int n)