Package rma.lang

Class RmaMath

java.lang.Object
rma.lang.RmaMath

public class RmaMath extends Object
  • Field Summary

    Fields
    Modifier and Type
    Field
    Description
    static double
     
    static int
     
  • Constructor Summary

    Constructors
    Constructor
    Description
     
  • Method Summary

    Modifier and Type
    Method
    Description
    static double
    balancedDistribution(double target, double[] dist, double[] weight, double[] min, double[] max, double tolfac, boolean[] active, boolean[] minlimited, boolean[] maxlimited)
    Develops a balanced distribution constrained that sums to the given target.
    static double
    balancedDistWithLoss(double target, double[] dist, double[] weight, double[] min, double[] max, double tolfac, boolean[] active, boolean[] minlimited, boolean[] maxlimited, double[] lossfac)
    Develops a balanced distribution constrained that sums to the given target.
    static double[]
    bestFitLinear(double[] x, double[] y, int ncoord, int istart)
    computes the best linear fit to a set of x and y coordinate pairs using least squares, (AT)A m = (AT) c return slope and y intercept in result[0] and result[1] If slope is infinte, result[0] is rma.util.UNDEF_DOUBLE and x intercept is returned in result[1]
    static double[]
    bestFitLinear(int[] x, int[] y, int ncoord, int istart)
    develops the best fit line to a set of x,y points.
    static double
    bisectRoot(RmaSVFunction func, double x0, double x1, int maxiter, double tol)
    Bisection method for finding the root of a function.
    static double
    BrentsMethodSolve(RmaSVFunction function, double lowerLimit, double upperLimit, int maxiter, double errorTol)
    Brent's Method for root finding from http://en.wikipedia.org/wiki/Brent%27s_method Also see Numerical Recipes in C Section 9.3 Van Wijngaarden-Dekker-Brent Method
    static double[][]
    cholesky(double[][] A)
     
    static void
    computeSunriseSunset(int dayOfYear, double latitude, double longitude, TimeZone tz, double[] sunriseContainer, double[] sunsetContainer)
     
    static double[]
    constructAffineTransform(double[] x, double[] y, double[] xp, double[] yp)
    constructs the best fit affine transform parameters from a set of native cooresponding trasformed coordinate pairs.
    static int[]
    decimalToFraction(double decimal)
    This method takes in a decimal value (0 <= X <= 1) and returns a fraction that will create that decimal value.
    static boolean
    equals(double d1, double d2, double tolerance)
    compare to doubles to see if they're equal to within tolerance
    static double
    fpRoot(RmaSVFunction func, double x0, double x1, int maxiter, double tol)
    False Position Bisection method for finding the root of a function.
    static double
    FUNCTION TO CALCULATE A GAUSSIAN DEVIATE, I.E.
    static boolean
    invertMatrix(double[] a)
    Inverts a square matrix using luDComp and luBackSub as described in Numerical Recipies in C (pg.
    static double
    linearInterpolate(double x1, double y1, double x2, double y2, double x)
    Interpolates a Y value along a ling for a given X position
    static double
    log10(double x)
    Returns the base 10 logarithm of a number
    static void
    luBackSub(double[] a, int n, int[] indx, double[] b)
    Backsubstitution routine from Numerical Receipies (pg 44).
    static boolean
    luDComp(double[] a, int n, int[] indx, int[] d)
    LU decomposition from Numerical Recipies (pg 43).
    static void
    main(String[] args)
     
    static boolean
    matAddS(double[] A, double S)
    Adds scalar S to matrix A.
    static boolean
    matAddTo(double[] A, double[] B)
    Add values in matrix A to the values in matrix B, with result returned in matrix B.
    static boolean
    matCharPoly(double[] A, double[] coeff, int n)
    Computes the charasteristic polynomial coefficients of a square matrix A.
    static double
    matDet(double[] A, int n)
    Computes the determinant of a square matrix A.
    static boolean
    matMultAB(double[] A, double[] B, double[] C, int nra, int nca, int ncb)
    Matrix Multiply - A * B = C.
    static boolean
    matMultABt(double[] A, double[] B, double[] C, int nra, int nca, int nrb)
    Matrix Multiply - A * B transpose = C.
    static boolean
    matMultAS(double[] m, double s, int nr, int nc)
    Multiplies the first matrix with a scalar value.
    static boolean
    matMultAtB(double[] A, double[] B, double[] C, int nra, int nca, int ncb)
    Matrix Multiply - A transpose * B = C.
    static int
    matNullity(double[] A, int ra, int ca)
    Computes the nullity of a matrix A.
    static int
    matNullityPSD(double[] A, int n)
    Computes the nullity of positive semi-definite matrix A.
    static int
    matRank(double[] A, int ra, int ca)
    Computes the rank of a matrix A.
    static int
    matRankPSD(double[] A, int n)
    Computes the rank of a positive semi-definite matrix A.
    static boolean
    matSubtTo(double[] A, double[] B)
    Subtracts values in matrix A to the values in matrix B, with result returned in matrix B (B - A => B) Matrices must be of the same length (number of rows and columns).
    static double
    matTrace(double[] A, int n)
    Computes the trace of a square matrix A.
    static boolean
    matTranspose(double[] m, int nr, int nc)
    Transposes the specified matrix.
    static int
    mod(int i1, int i2)
     
    static double[]
    normalize(double[] origValues, double sumValue)
    Takes a list of numbers and normalizes them so they sum to sumValue
    static int
    odeint2(double[] ystart, double x1, double x2, double eps, double h1, double hmin, int[] nok, int[] nbad, RmaDerivFunction derivfunc, int maxsteps, double[] xp, double[][] yp, double dxsav, int[] kount)
    ODE solver that utilizes second order runge kutta with adaptive step size to integrate a set of ODE equations over a given interval.
    static double
    redistributeWithConstraint(double total, double[] origDist, double[] dist, double[] min, double[] max, double tolfac, boolean[] active)
    Develops a balanced distribution constrained that sums to the given target.
    static void
    rk2qcStep(double[] y, double[] dydx, double[] x, double htry, double eps, double[] yscl, double[] hdid, double[] hnext, RmaDerivFunction derivfunc, int maxintv, double[] ysav, double[] dysav, double[] ytmp, double[] dytmp)
    Quality Controlled stepper function calling the second order runge kutta ODE routine.
    static double
    roundToSignificantDigits(double value, int significantDigits)
     
    static void
    rungeKutta2(double[] y0, double[] dydx, double x, double h, RmaDerivFunction derivfunc, double[] y1, double[] dydxm)
    Second Order Runge-Kutta integration step for an array of ODE's

    Methods inherited from class java.lang.Object

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

    • gdISET

      public static int gdISET
    • gdGSET

      public static double gdGSET
  • Constructor Details

    • RmaMath

      public RmaMath()
  • Method Details

    • gasdev

      public static double gasdev()
      FUNCTION TO CALCULATE A GAUSSIAN DEVIATE, I.E. NORMALLLY DISTRIBUTED WHITE NOISE WITH N(0,1) Taken from Numerical Recipies
      Returns:
      gaussian deviate
    • odeint2

      public static int odeint2(double[] ystart, double x1, double x2, double eps, double h1, double hmin, int[] nok, int[] nbad, RmaDerivFunction derivfunc, int maxsteps, double[] xp, double[][] yp, double dxsav, int[] kount)
      ODE solver that utilizes second order runge kutta with adaptive step size to integrate a set of ODE equations over a given interval. Based on Chapter 15.2 of Numerical Recipes in C
      Parameters:
      ystart - array of y values at initial condition
      x1 - starting value of indepedent variable
      x2 - ending value of intependent variable
      eps - error tolerance
      h1 - initial step size
      hmin - minimum allowable step size, may be zero
      nok - number of steps that were successfull without subdivision
      nbad - number of unsucessfull steps requiring further subdivision
      derivfunc - function that computes the derivative of each function
      maxsteps - max number of steps to span the interval
      xp - optional array to hold intermediate steps
      yp - optional array to hold intermediate y values at each step y[eqn][step]
      dxsav - interval at which so save intermediate steps
      Returns:
      zero if sucessful, -1=hit min step size limit, -2=hit max number of steps, -3 invalid input arrays
    • rk2qcStep

      public static void rk2qcStep(double[] y, double[] dydx, double[] x, double htry, double eps, double[] yscl, double[] hdid, double[] hnext, RmaDerivFunction derivfunc, int maxintv, double[] ysav, double[] dysav, double[] ytmp, double[] dytmp)
      Quality Controlled stepper function calling the second order runge kutta ODE routine. Based on Chapter 15.2 of Numerical Recipes in C.
      Parameters:
      y - array of initial y values
      dydx - array of initial dydx values
      x - initial x value
      htry - step size to try
      eps - error test scale factor
      yscl - array of y value error range
      hdid - step actually taken
      hnext - recommended next step
      derivfunc - function that evaluates derivatives
      maxintv - max number of interval subdivisions
      ysav - y value working array
      dysav - dy value working array
      ytmp - y value working array
      dytmp - dy value working array
    • rungeKutta2

      public static void rungeKutta2(double[] y0, double[] dydx, double x, double h, RmaDerivFunction derivfunc, double[] y1, double[] dydxm)
      Second Order Runge-Kutta integration step for an array of ODE's
      Parameters:
      y0 - array of values of dependent variables at x
      dydx - array of values of derivatives at x
      x - initial value of independent variable
      h - integration step
      derivfunc - function that returns dydx at a given x
      y1 - array of estimates of function y=f(x) at end of interval h
    • fpRoot

      public static double fpRoot(RmaSVFunction func, double x0, double x1, int maxiter, double tol)
      False Position Bisection method for finding the root of a function. For more info see Numerical Recipies in C, Section 9.2 Secant Method and False Position Method
      Returns:
    • bisectRoot

      public static double bisectRoot(RmaSVFunction func, double x0, double x1, int maxiter, double tol)
      Bisection method for finding the root of a function. For more info see Numerical Recipies in C, Section 9.2 Secant Method and False Position Method
      Returns:
    • BrentsMethodSolve

      public static double BrentsMethodSolve(RmaSVFunction function, double lowerLimit, double upperLimit, int maxiter, double errorTol)
      Brent's Method for root finding from http://en.wikipedia.org/wiki/Brent%27s_method Also see Numerical Recipes in C Section 9.3 Van Wijngaarden-Dekker-Brent Method
      Returns:
      root of given function
    • balancedDistribution

      public static double balancedDistribution(double target, double[] dist, double[] weight, double[] min, double[] max, double tolfac, boolean[] active, boolean[] minlimited, boolean[] maxlimited)
      Develops a balanced distribution constrained that sums to the given target. Optionally supports a weighted balance, constraints by min and max limits, matching within a tolerance, and limiting balance to only active items. Only the target and dist parameters are required. If other parameter arrays are null, they will be instantiated within the method with default values, which will add computation overhead.
      Parameters:
      target - total value to be distributed
      dist - double array that will contained the distributed values, must be a valid array with length > 0
      weight - optional array identifying balance weighting, if null an even balance is assumed
      min - optional array identifying the minimum distributed values, if null no min checking is performed
      max - optional array identifying the maximum distributed values, if null no max checking is performed
      tolfac - conditional tests in the algorithm use a tolerance value of tolfac*target, if tolfac is not a valid value then tolerance is zero
      active - optional boolean array identifying which distribution elements are included in the balance, if null all elements are active
      minlimited - optional array that will identify which distribution elements were minimum limited, may be null
      maxlimited - optional array that will identify which distribution elements were maximum limited, may be null
      Returns:
      total distribution or RMAConst.UNDEF_DOUBLE if there is an error
    • balancedDistWithLoss

      public static double balancedDistWithLoss(double target, double[] dist, double[] weight, double[] min, double[] max, double tolfac, boolean[] active, boolean[] minlimited, boolean[] maxlimited, double[] lossfac)
      Develops a balanced distribution constrained that sums to the given target. Optionally supports a weighted balance, constraints by min and max limits, matching within a tolerance, loss factors for each item, and limiting balance to only active items. Only the target and dist parameters are required. If other parameter arrays are null, they will be instantiated within the method with default values, which will add computation overhead.
      Parameters:
      target - total value to be distributed
      dist - double array that will contained the distributed values, must be a valid array with length > 0
      weight - optional array identifying balance weighting, if null an even balance is assumed
      min - optional array identifying the minimum distributed values, if null no min checking is performed
      max - optional array identifying the maximum distributed values, if null no max checking is performed
      tolfac - conditional tests in the algorithm use a tolerance value of tolfac*target, if tolfac is not a valid value then tolerance is zero
      active - optional boolean array identifying which distribution elements are included in the balance, if null all elements are active
      minlimited - optional array that will identify which distribution elements were minimum limited, may be null
      maxlimited - optional array that will identify which distribution elements were maximum limited, may be null
      lossfac - optional array that identifies loss factor, may be null
      Returns:
      total distribution or RMAConst.UNDEF_DOUBLE if there is an error
    • redistributeWithConstraint

      public static double redistributeWithConstraint(double total, double[] origDist, double[] dist, double[] min, double[] max, double tolfac, boolean[] active)
      Develops a balanced distribution constrained that sums to the given target. constraints by min and max limits, matching within a tolerance, and limiting balance to only active items. Only the target and dist parameters are required. If other parameter arrays are null, they will be instantiated within the method with default values, which will add computation overhead. Final distribution values are required to be greater than or equal to zero. Min and Max constraints are met if possible. This method was developed for ResSim downsream operation rule parallel system balancing.
      Parameters:
      total - value for final distribution
      origDist - distibution double array
      dist - double array that will contained the distributed values, must be a valid array with length > 0
      min - optional array identifying the minimum distributed values, if null no min checking is performed
      max - optional array identifying the maximum distributed values, if null no max checking is performed
      tolfac - conditional tests in the algorithm use a tolerance value of tolfac*target, if tolfac is not a valid value then tolerance is zero
      active - optional boolean array identifying which distribution elements are included in the balance, if null all elements are active
      Returns:
      total distribution or RMAConst.UNDEF_DOUBLE if there is an error
    • normalize

      public static double[] normalize(double[] origValues, double sumValue)
      Takes a list of numbers and normalizes them so they sum to sumValue
      Parameters:
      origValues - the list of numbers that will be noramlized
      sumValue - the value that all of the normalized values will sum to.
      Returns:
      any array of numbers that have been normalized. Returns null if the input is null. Otherwise returns an array of the exact same size as the input array.
    • linearInterpolate

      public static double linearInterpolate(double x1, double y1, double x2, double y2, double x)
      Interpolates a Y value along a ling for a given X position
      Parameters:
      x1 - - the x coordinate from the first point on the line
      y1 - - the y coordinate from the first point on the line
      x2 - - the x coordinate from the second point on the line
      y2 - - the y coordinate from the second point on the line
      x - - the x position of some point on the line
      Returns:
      - the y coordinate at position (x) on the line.
    • mod

      public static int mod(int i1, int i2)
    • bestFitLinear

      public static double[] bestFitLinear(double[] x, double[] y, int ncoord, int istart)
      computes the best linear fit to a set of x and y coordinate pairs using least squares, (AT)A m = (AT) c return slope and y intercept in result[0] and result[1] If slope is infinte, result[0] is rma.util.UNDEF_DOUBLE and x intercept is returned in result[1]
    • bestFitLinear

      public static double[] bestFitLinear(int[] x, int[] y, int ncoord, int istart)
      develops the best fit line to a set of x,y points.
      Returns:
      result contains slope [0] and intercept [1]
    • constructAffineTransform

      public static double[] constructAffineTransform(double[] x, double[] y, double[] xp, double[] yp)
      constructs the best fit affine transform parameters from a set of native cooresponding trasformed coordinate pairs. All input arrays must have the same length, and there must be at least three coordinate points. Solves for parameters using least squares (AT)A m = (AT) c
      Returns:
      a double array containing the affine tranform matrix parameters {m00,m10,m01,m11, m02, m12}, or null if there is an error
    • luDComp

      public static boolean luDComp(double[] a, int n, int[] indx, int[] d)
      LU decomposition from Numerical Recipies (pg 43). a is a 1-d array containing a square matrix in indexed by [row*n+col], n is the number of rows and columns, indx is an array of dimension n which will be used to store row permutations, d is a int array of dimension 1 that contiains +-1 indicating that row interchanges were odd or even. The LU decomposition is returned in a.
    • luBackSub

      public static void luBackSub(double[] a, int n, int[] indx, double[] b)
      Backsubstitution routine from Numerical Receipies (pg 44). a is the square matrix created by ludComp, n is the number of rows and columns, indx is also output from ludComp, b is input as the right-hand side and becomes the solution.
    • invertMatrix

      public static boolean invertMatrix(double[] a)
      Inverts a square matrix using luDComp and luBackSub as described in Numerical Recipies in C (pg. 45). Note that the matrix is inverted in place so the original matrix lost. a is a 1-d array containing a square matrix indexed by [row*n+col], n is the number of rows and columns.
    • matTranspose

      public static boolean matTranspose(double[] m, int nr, int nc)

      Transposes the specified matrix. Returns null if m is null

      Parameters:
      m -
      nr - - Number of Rows of original matrix.
      nc - - Number of Columns of original matrix.
      Returns:
      true if successful
    • matMultAS

      public static boolean matMultAS(double[] m, double s, int nr, int nc)

      Multiplies the first matrix with a scalar value. Returns null if matrix is null.

      Parameters:
      m -
      s - -- The scale value
      nr - -- Number of Rows
      nc - -- Number of Columns.
    • matMultAB

      public static boolean matMultAB(double[] A, double[] B, double[] C, int nra, int nca, int ncb)
      Matrix Multiply - A * B = C. Matricies are stored in singlely dimensioned arrays [irow * numcols + jcol]. The number of columns in A must equal the number of rows in B. Results are written into the C matrix. Matricies A and B are not changed.
      Parameters:
      A - first matrix (nra rows by nca columns)
      B - second matrix (nca rows by ncb columns)
      C - resulting matrix (nra rows by ncb columns)
      nra - number of rows in matrix A
      nca - number of columns in matrix A
      ncb - number of columns in matrix B
      Returns:
      true if successful
    • matMultAtB

      public static boolean matMultAtB(double[] A, double[] B, double[] C, int nra, int nca, int ncb)
      Matrix Multiply - A transpose * B = C. Matricies are stored in singlely dimensioned arrays [irow * numcols + jcol]. The number of rows in A must equal the number of rows in B. Results are written into the C matrix. Matricies A and B are not changed.
      Parameters:
      A - first matrix (nra rows by nca columns)
      B - second matrix (nca rows by ncb columns)
      C - resulting matrix (nra rows by ncb columns)
      nra - number of rows in matrix A
      nca - number of columns in matrix A
      ncb - number of columns in matrix B
      Returns:
      true if successful
    • matMultABt

      public static boolean matMultABt(double[] A, double[] B, double[] C, int nra, int nca, int nrb)
      Matrix Multiply - A * B transpose = C. Matricies are stored in singlely dimensioned arrays [irow * numcols + jcol]. The number of cols in A must equal the number of cols in B. Results are written into the C matrix. Matricies A and B are not changed.
      Parameters:
      A - first matrix (nra rows by nca columns)
      B - second matrix (nrb rows by nca columns)
      C - resulting matrix (nra rows by nrb columns)
      nra - number of rows in matrix A
      nca - number of columns in matrix A
      nrb - number of rows in matrix B
      Returns:
      true if successful
    • matAddTo

      public static boolean matAddTo(double[] A, double[] B)
      Add values in matrix A to the values in matrix B, with result returned in matrix B. Matricies must be of the same length (number of rows and columns).
      Parameters:
      A -
      B -
      Returns:
      true if successful
    • matSubtTo

      public static boolean matSubtTo(double[] A, double[] B)
      Subtracts values in matrix A to the values in matrix B, with result returned in matrix B (B - A => B) Matrices must be of the same length (number of rows and columns).
      Parameters:
      A -
      B -
      Returns:
      true if successful
    • matAddS

      public static boolean matAddS(double[] A, double S)
      Adds scalar S to matrix A.
      Parameters:
      A -
      S -
      Returns:
      true if successful
    • matTrace

      public static double matTrace(double[] A, int n)
      Computes the trace of a square matrix A.
      Parameters:
      A - , matrix of size n x n
      n - , matrix size = number of rows = number of columns
      Returns:
      trace of A
    • matCharPoly

      public static boolean matCharPoly(double[] A, double[] coeff, int n)
      Computes the charasteristic polynomial coefficients of a square matrix A. Based on the Faddeev-LeVerrier algorithm.
      Parameters:
      A - , matrix of size n x n
      coeff - array (a_n, a_n-1, ... , a_0) of length n+1 of charasteristic polynomial coefficients of A, in particular a_n = 1, a_n-1 = - trace(a), a_0 = (-1)^n det(A)
      n - , matrix size = number of rows = number of columns
      Returns:
      true if successful
    • matDet

      public static double matDet(double[] A, int n)
      Computes the determinant of a square matrix A.
      Parameters:
      A - , matrix of size n x n
      n - , matrix size = number of rows = number of columns
      Returns:
      determinant of A
    • matNullity

      public static int matNullity(double[] A, int ra, int ca)
      Computes the nullity of a matrix A.
      Parameters:
      A - , matrix of size ra x ca
      ra - , number of rows of matrix A
      ca - , number of columns of matrix A
      Returns:
      nullity of A
    • matRank

      public static int matRank(double[] A, int ra, int ca)
      Computes the rank of a matrix A.
      Parameters:
      A - , matrix of size ra x ca
      ra - , number of rows of matrix A
      ca - , number of columns of matrix A
      Returns:
      rank of A Note: based on the idea that rank(A) = rank (At A) = rank(A At)
    • matNullityPSD

      public static int matNullityPSD(double[] A, int n)
      Computes the nullity of positive semi-definite matrix A.
      Parameters:
      A - , positive semi-definite matrix of size n x n
      n - , number of rows and columns of matrix A
      Returns:
      nullity of A
    • matRankPSD

      public static int matRankPSD(double[] A, int n)
      Computes the rank of a positive semi-definite matrix A.
      Parameters:
      A - , positive semi-definite matrix of size n x n
      n - , number of rows and columns of matrix A
      Returns:
      rank of A Note: Re-scaling is necessary because char poly terms grow exponentially Re-scaling is based in the identity det(A)<=(tr(A)/n)^n and similarly for other characteristic polynomial coefficients a_i = O(tr(A)/n)^(n-i) (for the positive semi-definite matrix A)
    • roundToSignificantDigits

      public static double roundToSignificantDigits(double value, int significantDigits)
    • log10

      public static double log10(double x)
      Returns the base 10 logarithm of a number
    • computeSunriseSunset

      public static void computeSunriseSunset(int dayOfYear, double latitude, double longitude, TimeZone tz, double[] sunriseContainer, double[] sunsetContainer)
    • decimalToFraction

      public static int[] decimalToFraction(double decimal)
      This method takes in a decimal value (0 <= X <= 1) and returns a fraction that will create that decimal value. The return falue is int[0] numerator int[1] denominator
    • equals

      public static boolean equals(double d1, double d2, double tolerance)
      compare to doubles to see if they're equal to within tolerance
      Returns:
      true if the numbers are equal within the tolerance
    • cholesky

      public static double[][] cholesky(double[][] A)
    • main

      public static void main(String[] args)