minpack_module Module

Modernized Minpack

Authors

  • argonne national laboratory. minpack project. march 1980. burton s. garbow, kenneth e. hillstrom, jorge j. more.
  • Jacob Williams, Sept 2021, updated to modern standards.

Uses

  • module~~minpack_module~~UsesGraph module~minpack_module minpack_module iso_fortran_env iso_fortran_env module~minpack_module->iso_fortran_env

Used by

  • module~~minpack_module~~UsedByGraph module~minpack_module minpack_module module~minpack_capi minpack_capi module~minpack_capi->module~minpack_module

Contents


Variables

TypeVisibilityAttributesNameInitial
real(kind=wp), public, parameter, dimension(3):: dpmpar =[epsilon(1.0_wp), tiny(1.0_wp), huge(1.0_wp)]

machine constants

real(kind=wp), private, parameter:: epsmch =dpmpar(1)

the machine precision

real(kind=wp), private, parameter:: one =1.0_wp
real(kind=wp), private, parameter:: zero =0.0_wp

Abstract Interfaces

abstract interface

  • public subroutine func(n, x, fvec, iflag)

    user-supplied subroutine for hybrd, hybrd1, and fdjac1

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: n

    the number of variables.

    real(kind=wp), intent(in) :: x(n)

    independent variable vector

    real(kind=wp), intent(out) :: fvec(n)

    value of function at x

    integer, intent(inout) :: iflag

    set to <0 to terminate execution

abstract interface

  • public subroutine func2(m, n, x, fvec, iflag)

    user-supplied subroutine for fdjac2, lmdif, and lmdif1

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: m

    the number of functions.

    integer, intent(in) :: n

    the number of variables.

    real(kind=wp), intent(in) :: x(n)

    independent variable vector

    real(kind=wp), intent(out) :: fvec(m)

    value of function at x

    integer, intent(inout) :: iflag

    the value of iflag should not be changed unless the user wants to terminate execution of lmdif. in this case set iflag to a negative integer.

abstract interface

  • public subroutine fcn_hybrj(n, x, fvec, fjac, ldfjac, iflag)

    user-supplied subroutine for hybrj and hybrj1

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: n

    the number of variables.

    real(kind=wp), intent(in), dimension(n):: x

    independent variable vector

    real(kind=wp), intent(inout), dimension(n):: fvec

    value of function at x

    real(kind=wp), intent(inout), dimension(ldfjac, n):: fjac

    jacobian matrix at x

    integer, intent(in) :: ldfjac

    leading dimension of the array fjac.

    integer, intent(inout) :: iflag

    if iflag = 1 calculate the functions at x and return this vector in fvec. do not alter fjac. if iflag = 2 calculate the jacobian at x and return this matrix in fjac. do not alter fvec.

    the value of iflag should not be changed by fcn unless the user wants to terminate execution of hybrj. in this case set iflag to a negative integer.

abstract interface

  • public subroutine fcn_lmder(m, n, x, fvec, fjac, ldfjac, iflag)

    user-supplied subroutine for lmder and lmder1

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: m

    the number of functions.

    integer, intent(in) :: n

    the number of variables.

    real(kind=wp), intent(in) :: x(n)

    independent variable vector

    real(kind=wp), intent(inout) :: fvec(m)

    value of function at x

    real(kind=wp), intent(inout) :: fjac(ldfjac,n)

    jacobian matrix at x

    integer, intent(in) :: ldfjac

    leading dimension of the array fjac.

    integer, intent(inout) :: iflag

    if iflag = 1 calculate the functions at x and return this vector in fvec. do not alter fjac. if iflag = 2 calculate the jacobian at x and return this matrix in fjac. do not alter fvec.

    the value of iflag should not be changed by fcn unless the user wants to terminate execution of lmder. in this case set iflag to a negative integer.

abstract interface

  • public subroutine fcn_lmstr(m, n, x, fvec, fjrow, iflag)

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: m

    the number of functions.

    integer, intent(in) :: n

    the number of variables.

    real(kind=wp), intent(in) :: x(n)

    independent variable vector

    real(kind=wp), intent(inout) :: fvec(m)

    value of function at x

    real(kind=wp), intent(inout) :: fjrow(n)

    jacobian row

    integer, intent(inout) :: iflag

    if iflag = 1 calculate the functions at x and return this vector in fvec. if iflag = i calculate the (i-1)-st row of the jacobian at x and return this vector in fjrow.

    the value of iflag should not be changed by fcn unless the user wants to terminate execution of lmstr. in this case set iflag to a negative integer.


Functions

public pure function enorm(n, x)

given an n-vector x, this function calculates the euclidean norm of x.

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

a positive integer input variable.

real(kind=wp), intent(in) :: x(n)

an input array of length n.

Return Value real(kind=wp)


Subroutines

public subroutine chkder(m, n, x, Fvec, Fjac, Ldfjac, Xp, Fvecp, Mode, Err)

this subroutine checks the gradients of m nonlinear functions in n variables, evaluated at a point x, for consistency with the functions themselves.

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables.

real(kind=wp), intent(in) :: x(n)

input array

real(kind=wp), intent(in) :: Fvec(m)

an array of length m. on input when mode = 2, fvec must contain the functions evaluated at x.

real(kind=wp), intent(in) :: Fjac(Ldfjac,n)

an m by n array. on input when mode = 2, the rows of fjac must contain the gradients of the respective functions evaluated at x.

integer, intent(in) :: Ldfjac

a positive integer input parameter not less than m which specifies the leading dimension of the array fjac.

real(kind=wp), intent(out) :: Xp(n)

an array of length n. on output when mode = 1, xp is set to a neighboring point of x.

real(kind=wp), intent(in) :: Fvecp(m)

an array of length m. on input when mode = 2, fvecp must contain the functions evaluated at xp.

integer, intent(in) :: Mode

an integer input variable set to 1 on the first call and 2 on the second. other values of mode are equivalent to mode = 1.

the user must call chkder twice, first with mode = 1 and then with mode = 2.

  • mode = 1. on input, x must contain the point of evaluation. on output, xp is set to a neighboring point.

  • mode = 2. on input, fvec must contain the functions and the rows of fjac must contain the gradients of the respective functions each evaluated at x, and fvecp must contain the functions evaluated at xp. on output, err contains measures of correctness of the respective gradients.

real(kind=wp), intent(out) :: Err(m)

an array of length m. on output when mode = 2, err contains measures of correctness of the respective gradients. if there is no severe loss of significance, then if err(i) is 1.0 the i-th gradient is correct, while if err(i) is 0.0 the i-th gradient is incorrect. for values of err between 0.0 and 1.0, the categorization is less certain. in general, a value of err(i) greater than 0.5 indicates that the i-th gradient is probably correct, while a value of err(i) less than 0.5 indicates that the i-th gradient is probably incorrect.

public subroutine dogleg(n, r, Lr, Diag, Qtb, Delta, x, Wa1, Wa2)

given an m by n matrix a, an n by n nonsingular diagonal matrix d, an m-vector b, and a positive number delta, the problem is to determine the convex combination x of the gauss-newton and scaled gradient directions that minimizes (ax - b) in the least squares sense, subject to the restriction that the euclidean norm of dx be at most delta.

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

a positive integer input variable set to the order of r.

real(kind=wp), intent(in) :: r(Lr)

an input array of length lr which must contain the upper triangular matrix r stored by rows.

integer, intent(in) :: Lr

a positive integer input variable not less than (n*(n+1))/2.

real(kind=wp), intent(in) :: Diag(n)

an input array of length n which must contain the diagonal elements of the matrix d.

real(kind=wp), intent(in) :: Qtb(n)

an input array of length n which must contain the first n elements of the vector (q transpose)*b.

real(kind=wp), intent(in) :: Delta

a positive input variable which specifies an upper bound on the euclidean norm of d*x.

real(kind=wp), intent(out) :: x(n)

an output array of length n which contains the desired convex combination of the gauss-newton direction and the scaled gradient direction.

real(kind=wp), intent(inout) :: Wa1(n)

work arrays of length n

real(kind=wp), intent(inout) :: Wa2(n)

work arrays of length n

public subroutine fdjac1(fcn, n, x, Fvec, Fjac, Ldfjac, Iflag, Ml, Mu, Epsfcn, Wa1, Wa2)

this subroutine computes a forward-difference approximation to the n by n jacobian matrix associated with a specified problem of n functions in n variables. if the jacobian has a banded form, then function evaluations are saved by only approximating the nonzero terms.

Arguments

TypeIntentOptionalAttributesName
procedure(func) :: fcn

the user-supplied subroutine which calculates the functions.

integer, intent(in) :: n

a positive integer input variable set to the number of functions and variables.

real(kind=wp), intent(inout) :: x(n)

an input array of length n.

real(kind=wp), intent(in) :: Fvec(n)

an input array of length n which must contain the functions evaluated at x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

an output n by n array which contains the approximation to the jacobian matrix evaluated at x.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than n which specifies the leading dimension of the array fjac.

integer, intent(inout) :: Iflag

an integer variable which can be used to terminate the execution of fdjac1. see description of func.

integer, intent(in) :: Ml

a nonnegative integer input variable which specifies the number of subdiagonals within the band of the jacobian matrix. if the jacobian is not banded, set ml to at least n - 1.

integer, intent(in) :: Mu

a nonnegative integer input variable which specifies the number of superdiagonals within the band of the jacobian matrix. if the jacobian is not banded, set mu to at least n - 1.

real(kind=wp), intent(in) :: Epsfcn

an input variable used in determining a suitable step length for the forward-difference approximation. this approximation assumes that the relative errors in the functions are of the order of epsfcn. if epsfcn is less than the machine precision, it is assumed that the relative errors in the functions are of the order of the machine precision.

real(kind=wp), intent(inout) :: Wa1(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa2(n)

work array of length n. if ml + mu + 1 is at least n, then the jacobian is considered dense, and wa2 is not referenced.

public subroutine fdjac2(fcn, m, n, x, Fvec, Fjac, Ldfjac, Iflag, Epsfcn, Wa)

this subroutine computes a forward-difference approximation to the m by n jacobian matrix associated with a specified problem of m functions in n variables.

Arguments

TypeIntentOptionalAttributesName
procedure(func2) :: fcn

the user-supplied subroutine which calculates the functions.

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

an input array of length n.

real(kind=wp), intent(in) :: Fvec(m)

an input array of length m which must contain the functions evaluated at x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

an output m by n array which contains the approximation to the jacobian matrix evaluated at x.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than m which specifies the leading dimension of the array fjac.

integer, intent(inout) :: Iflag

an integer variable which can be used to terminate the execution of fdjac2. see description of func2.

real(kind=wp), intent(in) :: Epsfcn

an input variable used in determining a suitable step length for the forward-difference approximation. this approximation assumes that the relative errors in the functions are of the order of epsfcn. if epsfcn is less than the machine precision, it is assumed that the relative errors in the functions are of the order of the machine precision.

real(kind=wp), intent(inout) :: Wa(m)

a work array of length m.

public subroutine hybrd(fcn, n, x, fvec, xtol, maxfev, ml, mu, epsfcn, diag, mode, factor, nprint, info, nfev, fjac, ldfjac, r, lr, qtf, wa1, wa2, wa3, wa4)

the purpose of hybrd is to find a zero of a system of n nonlinear functions in n variables by a modification of the powell hybrid method. the user must provide a subroutine which calculates the functions. the jacobian is then calculated by a forward-difference approximation.

Arguments

TypeIntentOptionalAttributesName
procedure(func) :: fcn

user-supplied subroutine which calculates the functions

integer, intent(in) :: n

a positive integer input variable set to the number of functions and variables.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: fvec(n)

an output array of length n which contains the functions evaluated at the output x.

real(kind=wp), intent(in) :: xtol

a nonnegative input variable. termination occurs when the relative error between two consecutive iterates is at most xtol.

integer, intent(in) :: maxfev

a positive integer input variable. termination occurs when the number of calls to fcn is at least maxfev by the end of an iteration.

integer, intent(in) :: ml

a nonnegative integer input variable which specifies the number of subdiagonals within the band of the jacobian matrix. if the jacobian is not banded, set ml to at least n - 1.

integer, intent(in) :: mu

a nonnegative integer input variable which specifies the number of superdiagonals within the band of the jacobian matrix. if the jacobian is not banded, set mu to at leastn - 1.

real(kind=wp), intent(in) :: epsfcn

an input variable used in determining a suitable step length for the forward-difference approximation. this approximation assumes that the relative errors in the functions are of the order of epsfcn. if epsfcn is less than the machine precision, it is assumed that the relative errors in the functions are of the order of the machine precision.

real(kind=wp), intent(inout) :: diag(n)

an array of length n. if mode = 1 (see below), diag is internally set. if mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.

integer, intent(in) :: mode

if mode = 1, the variables will be scaled internally. if mode = 2, the scaling is specified by the input diag. other values of mode are equivalent to mode = 1.

real(kind=wp), intent(in) :: factor

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 (.1,100.). 100. is a generally recommended value.

integer, intent(in) :: nprint

an integer input variable that enables controlled printing of iterates if it is positive. in this case, fcn is called with iflag = 0 at the beginning of the first iteration and every nprint iterations thereafter and immediately prior to return, with x and fvec available for printing. if nprint is not positive, no special calls of fcn with iflag = 0 are made.

integer, intent(out) :: info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 relative error between two consecutive iterates is at most xtol.
  • info = 2 number of calls to fcn has reached or exceeded maxfev.
  • info = 3 xtol is too small. no further improvement in the approximate solution x is possible.
  • info = 4 iteration is not making good progress, as measured by the improvement from the last five jacobian evaluations.
  • info = 5 iteration is not making good progress, as measured by the improvement from the last ten iterations.
integer, intent(out) :: nfev

output variable set to the number of calls to fcn.

real(kind=wp), intent(out) :: fjac(ldfjac,n)

array which contains the orthogonal matrix q produced by the QR factorization of the final approximate jacobian.

integer, intent(in) :: ldfjac

a positive integer input variable not less than n which specifies the leading dimension of the array fjac.

real(kind=wp), intent(out) :: r(lr)

an output array which contains the upper triangular matrix produced by the QR factorization of the final approximate jacobian, stored rowwise.

integer, intent(in) :: lr

a positive integer input variable not less than (n*(n+1))/2.

real(kind=wp), intent(out) :: qtf(n)

an output array of length n which contains the vector (q transpose)*fvec.

real(kind=wp), intent(inout) :: wa1(n)

work array

real(kind=wp), intent(inout) :: wa2(n)

work array

real(kind=wp), intent(inout) :: wa3(n)

work array

real(kind=wp), intent(inout) :: wa4(n)

work array

public subroutine hybrd1(fcn, n, x, fvec, tol, info, Wa, Lwa)

the purpose of hybrd1 is to find a zero of a system of n nonlinear functions in n variables by a modification of the powell hybrid method. this is done by using the more general nonlinear equation solver hybrd. the user must provide a subroutine which calculates the functions. the jacobian is then calculated by a forward-difference approximation.

Arguments

TypeIntentOptionalAttributesName
procedure(func) :: fcn

user-supplied subroutine which calculates the functions

integer, intent(in) :: n

a positive integer input variable set to the number of functions and variables.

real(kind=wp), intent(inout), dimension(n):: x

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.

real(kind=wp), intent(out), dimension(n):: fvec

an output array of length n which contains the functions evaluated at the output x.

real(kind=wp), intent(in) :: tol

a nonnegative input variable. termination occurs when the algorithm estimates that the relative error between x and the solution is at most tol.

integer, intent(out) :: info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 algorithm estimates that the relative error between x and the solution is at most tol.
  • info = 2 number of calls to fcn has reached or exceeded 200*(n+1).
  • info = 3 tol is too small. no further improvement in the approximate solution x is possible.
  • info = 4 iteration is not making good progress.
real(kind=wp), intent(inout) :: Wa(Lwa)

a work array of length lwa.

integer, intent(in) :: Lwa

a positive integer input variable not less than (n(3n+13))/2.

public subroutine hybrj(fcn, n, x, Fvec, Fjac, Ldfjac, Xtol, Maxfev, Diag, Mode, Factor, Nprint, Info, Nfev, Njev, r, Lr, Qtf, Wa1, Wa2, Wa3, Wa4)

the purpose of hybrj is to find a zero of a system of n nonlinear functions in n variables by a modification of the powell hybrid method. the user must provide a subroutine which calculates the functions and the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_hybrj) :: fcn

the user-supplied subroutine which calculates the functions and the jacobian

integer, intent(in) :: n

a positive integer input variable set to the number of functions and variables.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(n)

an output array of length n which contains the functions evaluated at the output x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

an output n by n array which contains the orthogonal matrix q produced by the qr factorization of the final approximate jacobian.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than n which specifies the leading dimension of the array fjac.

real(kind=wp), intent(in) :: Xtol

a nonnegative input variable. termination occurs when the relative error between two consecutive iterates is at most xtol.

integer, intent(in) :: Maxfev

a positive integer input variable. termination occurs when the number of calls to fcn with iflag = 1 has reached maxfev.

real(kind=wp), intent(inout) :: Diag(n)

an array of length n. if mode = 1 (see below), diag is internally set. if mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.

integer, intent(in) :: Mode

an integer input variable. if mode = 1, the variables will be scaled internally. if mode = 2, the scaling is specified by the input diag. other values of mode are equivalent to mode = 1.

real(kind=wp), intent(in) :: Factor

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 (.1,100.). 100. is a generally recommended value.

integer, intent(in) :: Nprint

an integer input variable that enables controlled printing of iterates if it is positive. in this case, fcn is called with iflag = 0 at the beginning of the first iteration and every nprint iterations thereafter and immediately prior to return, with x and fvec available for printing. fvec and fjac should not be altered. if nprint is not positive, no special calls of fcn with iflag = 0 are made.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 relative error between two consecutive iterates is at most xtol.
  • info = 2 number of calls to fcn with iflag = 1 has reached maxfev.
  • info = 3 xtol is too small. no further improvement in the approximate solution x is possible.
  • info = 4 iteration is not making good progress, as measured by the improvement from the last five jacobian evaluations.
  • info = 5 iteration is not making good progress, as measured by the improvement from the last ten iterations.
integer, intent(out) :: Nfev

an integer output variable set to the number of calls to fcn with iflag = 1.

integer, intent(out) :: Njev

an integer output variable set to the number of calls to fcn with iflag = 2.

real(kind=wp), intent(out) :: r(Lr)

an output array of length lr which contains the upper triangular matrix produced by the qr factorization of the final approximate jacobian, stored rowwise.

integer, intent(in) :: Lr

a positive integer input variable not less than (n*(n+1))/2.

real(kind=wp), intent(out) :: Qtf(n)

an output array of length n which contains the vector (q transpose)*fvec.

real(kind=wp), intent(inout) :: Wa1(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa2(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa3(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa4(n)

work array of length n.

public subroutine hybrj1(fcn, n, x, Fvec, Fjac, Ldfjac, Tol, Info, Wa, Lwa)

the purpose of hybrj1 is to find a zero of a system of n nonlinear functions in n variables by a modification of the powell hybrid method. this is done by using the more general nonlinear equation solver hybrj. the user must provide a subroutine which calculates the functions and the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_hybrj) :: fcn

the user-supplied subroutine which calculates the functions and the jacobian

integer, intent(in) :: n

a positive integer input variable set to the number of functions and variables.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(n)

an output array of length n which contains the functions evaluated at the output x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

an output n by n array which contains the orthogonal matrix q produced by the qr factorization of the final approximate jacobian.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than n which specifies the leading dimension of the array fjac.

real(kind=wp), intent(in) :: Tol

a nonnegative input variable. termination occurs when the algorithm estimates that the relative error between x and the solution is at most tol.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 algorithm estimates that the relative error between x and the solution is at most tol.
  • info = 2 number of calls to fcn with iflag = 1 has reached 100*(n+1).
  • info = 3 tol is too small. no further improvement in the approximate solution x is possible.
  • info = 4 iteration is not making good progress.
real(kind=wp), intent(inout) :: Wa(Lwa)

a work array of length lwa.

integer, intent(in) :: Lwa

a positive integer input variable not less than (n*(n+13))/2.

public subroutine lmder(fcn, m, n, x, Fvec, Fjac, Ldfjac, Ftol, Xtol, Gtol, Maxfev, Diag, Mode, Factor, Nprint, Info, Nfev, Njev, Ipvt, Qtf, Wa1, Wa2, Wa3, Wa4)

the purpose of lmder is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm. the user must provide a subroutine which calculates the functions and the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_lmder) :: fcn

the user-supplied subroutine which calculates the functions and the jacobian

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(m)

an output array of length m which contains the functions evaluated at the output x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

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

        t     t           t
       p *(jac *jac)*p = r *r,

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.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than m which specifies the leading dimension of the array fjac.

real(kind=wp), intent(in) :: Ftol

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.

real(kind=wp), intent(in) :: Xtol

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.

real(kind=wp), intent(in) :: Gtol

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.

integer, intent(in) :: Maxfev

a positive integer input variable. termination occurs when the number of calls to fcn with iflag = 1 has reached maxfev.

real(kind=wp), intent(inout) :: Diag(n)

an array of length n. if mode = 1 (see below), diag is internally set. if mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.

integer, intent(in) :: Mode

an integer input variable. if mode = 1, the variables will be scaled internally. if mode = 2, the scaling is specified by the input diag. other values of mode are equivalent to mode = 1.

real(kind=wp), intent(in) :: Factor

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 (.1,100.).100. is a generally recommended value.

integer, intent(in) :: Nprint

an integer input variable that enables controlled printing of iterates if it is positive. in this case, fcn is called with iflag = 0 at the beginning of the first iteration and every nprint iterations thereafter and immediately prior to return, with x, fvec, and fjac available for printing. fvec and fjac should not be altered. if nprint is not positive, no special calls of fcn with iflag = 0 are made.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 both actual and predicted relative reductions in the sum of squares are at most ftol.
  • info = 2 relative error between two consecutive iterates is at most xtol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 the cosine of the angle between fvec and any column of the jacobian is at most gtol in absolute value.
  • info = 5 number of calls to fcn with iflag = 1 has reached maxfev.
  • info = 6 ftol is too small. no further reduction in the sum of squares is possible.
  • info = 7 xtol is too small. no further improvement in the approximate solution x is possible.
  • info = 8 gtol is too small. fvec is orthogonal to the columns of the jacobian to machine precision.
integer, intent(out) :: Nfev

an integer output variable set to the number of calls to fcn with iflag = 1.

integer, intent(out) :: Njev

an integer output variable set to the number of calls to fcn with iflag = 2.

integer, intent(out) :: Ipvt(n)

an integer output array of length n. ipvt defines a permutation matrix p such that jacp = qr, 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.

real(kind=wp), intent(out) :: Qtf(n)

an output array of length n which contains the first n elements of the vector (q transpose)*fvec.

real(kind=wp), intent(inout) :: Wa1(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa2(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa3(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa4(m)

work array of length m.

public subroutine lmder1(fcn, m, n, x, Fvec, Fjac, Ldfjac, Tol, Info, Ipvt, Wa, Lwa)

the purpose of lmder1 is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm. this is done by using the more general least-squares solver lmder. the user must provide a subroutine which calculates the functions and the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_lmder) :: fcn

user-supplied subroutine which calculates the functions and the jacobian.

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(m)

an output array of length m which contains the functions evaluated at the output x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

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

        t     t           t
       p *(jac *jac)*p = r *r,

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.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than m which specifies the leading dimension of the array fjac.

real(kind=wp), intent(in) :: Tol

a nonnegative input variable. termination occurs when the algorithm estimates either that the relative error in the sum of squares is at most tol or that the relative error between x and the solution is at most tol.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows.

  • info = 0 improper input parameters.
  • info = 1 algorithm estimates that the relative error in the sum of squares is at most tol.
  • info = 2 algorithm estimates that the relative error between x and the solution is at most tol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 fvec is orthogonal to the columns of the jacobian to machine precision.
  • info = 5 number of calls to fcn with iflag = 1 has reached 100*(n+1).
  • info = 6 tol is too small. no further reduction in the sum of squares is possible.
  • info = 7 tol is too small. no further improvement in the approximate solution x is possible.
integer, intent(out) :: Ipvt(n)

an integer output array of length n. ipvt defines a permutation matrix p such that jacp = qr, 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.

real(kind=wp), intent(inout) :: Wa(Lwa)

a work array of length lwa.

integer, intent(in) :: Lwa

a positive integer input variable not less than 5*n+m.

public subroutine lmdif(fcn, m, n, x, Fvec, Ftol, Xtol, Gtol, Maxfev, Epsfcn, Diag, Mode, Factor, Nprint, Info, Nfev, Fjac, Ldfjac, Ipvt, Qtf, Wa1, Wa2, Wa3, Wa4)

the purpose of lmdif is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm. the user must provide a subroutine which calculates the functions. the jacobian is then calculated by a forward-difference approximation.

Arguments

TypeIntentOptionalAttributesName
procedure(func2) :: fcn

the user-supplied subroutine which calculates the functions.

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(m)

an output array of length m which contains the functions evaluated at the output x.

real(kind=wp), intent(in) :: Ftol

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.

real(kind=wp), intent(in) :: Xtol

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.

real(kind=wp), intent(in) :: Gtol

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.

integer, intent(in) :: Maxfev

a positive integer input variable. termination occurs when the number of calls to fcn is at least maxfev by the end of an iteration.

real(kind=wp), intent(in) :: Epsfcn

an input variable used in determining a suitable step length for the forward-difference approximation. this approximation assumes that the relative errors in the functions are of the order of epsfcn. if epsfcn is less than the machine precision, it is assumed that the relative errors in the functions are of the order of the machine precision.

real(kind=wp), intent(inout) :: Diag(n)

an array of length n. if mode = 1 (see below), diag is internally set. if mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.

integer, intent(in) :: Mode

an integer input variable. if mode = 1, the variables will be scaled internally. if mode = 2, the scaling is specified by the input diag. other values of mode are equivalent to mode = 1.

real(kind=wp), intent(in) :: Factor

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 (.1,100.). 100. is a generally recommended value.

integer, intent(in) :: Nprint

an integer input variable that enables controlled printing of iterates if it is positive. in this case, fcn is called with iflag = 0 at the beginning of the first iteration and every nprint iterations thereafter and immediately prior to return, with x and fvec available for printing. if nprint is not positive, no special calls of fcn with iflag = 0 are made.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 both actual and predicted relative reductions in the sum of squares are at most ftol.
  • info = 2 relative error between two consecutive iterates is at most xtol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 the cosine of the angle between fvec and any column of the jacobian is at most gtol in absolute value.
  • info = 5 number of calls to fcn has reached or exceeded maxfev.
  • info = 6 ftol is too small. no further reduction in the sum of squares is possible.
  • info = 7 xtol is too small. no further improvement in the approximate solution x is possible.
  • info = 8 gtol is too small. fvec is orthogonal to the columns of the jacobian to machine precision.
integer, intent(out) :: Nfev

an integer output variable set to the number of calls to fcn.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

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

        t     t           t
       p *(jac *jac)*p = r *r,

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.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than m which specifies the leading dimension of the array fjac.

integer, intent(out) :: Ipvt(n)

an integer output array of length n. ipvt defines a permutation matrix p such that jacp = qr, 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.

real(kind=wp), intent(out) :: Qtf(n)

an output array of length n which contains the first n elements of the vector (q transpose)*fvec.

real(kind=wp), intent(inout) :: Wa1(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa2(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa3(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa4(m)

work array of length n.

public subroutine lmdif1(fcn, m, n, x, Fvec, Tol, Info, Iwa, Wa, Lwa)

the purpose of lmdif1 is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm. this is done by using the more general least-squares solver lmdif. the user must provide a subroutine which calculates the functions. the jacobian is then calculated by a forward-difference approximation.

Arguments

TypeIntentOptionalAttributesName
procedure(func2) :: fcn

the user-supplied subroutine which calculates the functions.

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(m)

an output array of length m which contains the functions evaluated at the output x.

real(kind=wp), intent(in) :: Tol

a nonnegative input variable. termination occurs when the algorithm estimates either that the relative error in the sum of squares is at most tol or that the relative error between x and the solution is at most tol.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 algorithm estimates that the relative error in the sum of squares is at most tol.
  • info = 2 algorithm estimates that the relative error between x and the solution is at most tol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 fvec is orthogonal to the columns of the jacobian to machine precision.
  • info = 5 number of calls to fcn has reached or exceeded 200*(n+1).
  • info = 6 tol is too small. no further reduction in the sum of squares is possible.
  • info = 7 tol is too small. no further improvement in the approximate solution x is possible.
integer, intent(inout) :: Iwa(n)

an integer work array of length n.

real(kind=wp), intent(inout) :: Wa(Lwa)

a work array of length lwa.

integer, intent(in) :: Lwa

a positive integer input variable not less than mn+5n+m.

public subroutine lmpar(n, r, Ldr, Ipvt, Diag, Qtb, Delta, Par, x, Sdiag, Wa1, Wa2)

given an m by n matrix a, an n by n nonsingular diagonal matrix d, an m-vector b, and a positive number delta, the problem is to determine a value for the parameter par such that if x solves the system

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

a positive integer input variable set to the order of r.

real(kind=wp), intent(inout) :: r(Ldr,n)

an n by n array. on input the full upper triangle must contain the full upper triangle of the matrix r. on output the full upper triangle is unaltered, and the strict lower triangle contains the strict upper triangle (transposed) of the upper triangular matrix s.

integer, intent(in) :: Ldr

a positive integer input variable not less than n which specifies the leading dimension of the array r.

integer, intent(in) :: Ipvt(n)

an integer input array of length n which defines the permutation matrix p such that ap = qr. column j of p is column ipvt(j) of the identity matrix.

real(kind=wp), intent(in) :: Diag(n)

an input array of length n which must contain the diagonal elements of the matrix d.

real(kind=wp), intent(in) :: Qtb(n)

an input array of length n which must contain the first n elements of the vector (q transpose)*b.

real(kind=wp) :: Delta

a positive input variable which specifies an upper bound on the euclidean norm of d*x.

real(kind=wp), intent(inout) :: Par

a nonnegative variable. on input par contains an initial estimate of the levenberg-marquardt parameter. on output par contains the final estimate.

real(kind=wp), intent(out) :: x(n)

an output array of length n which contains the least squares solution of the system ax = b, sqrt(par)d*x = 0, for the output par.

real(kind=wp), intent(out) :: Sdiag(n)

an output array of length n which contains the diagonal elements of the upper triangular matrix s.

real(kind=wp), intent(inout) :: Wa1(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa2(n)

work array of length n.

public subroutine lmstr(fcn, m, n, x, Fvec, Fjac, Ldfjac, Ftol, Xtol, Gtol, Maxfev, Diag, Mode, Factor, Nprint, Info, Nfev, Njev, Ipvt, Qtf, Wa1, Wa2, Wa3, Wa4)

the purpose of lmstr is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm which uses minimal storage. the user must provide a subroutine which calculates the functions and the rows of the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_lmstr) :: fcn

user-supplied subroutine which calculates the functions and the rows of the jacobian.

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(m)

an output array of length m which contains the functions evaluated at the output x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

an output n by n array. the upper triangle of fjac contains an upper triangular matrix r such that

        t     t           t
       p *(jac *jac)*p = r *r,

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 triangular part of fjac contains information generated during the computation of r.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than n which specifies the leading dimension of the array fjac.

real(kind=wp), intent(in) :: Ftol

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.

real(kind=wp), intent(in) :: Xtol

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.

real(kind=wp), intent(in) :: Gtol

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.

integer, intent(in) :: Maxfev

a positive integer input variable. termination occurs when the number of calls to fcn with iflag = 1 has reached maxfev.

real(kind=wp), intent(inout) :: Diag(n)

an array of length n. if mode = 1 (see below), diag is internally set. if mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.

integer, intent(in) :: Mode

an integer input variable. if mode = 1, the variables will be scaled internally. if mode = 2, the scaling is specified by the input diag. other values of mode are equivalent to mode = 1.

real(kind=wp), intent(in) :: Factor

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 (.1,100.). 100. is a generally recommended value.

integer, intent(in) :: Nprint

an integer input variable that enables controlled printing of iterates if it is positive. in this case, fcn is called with iflag = 0 at the beginning of the first iteration and every nprint iterations thereafter and immediately prior to return, with x and fvec available for printing. if nprint is not positive, no special calls of fcn with iflag = 0 are made.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 both actual and predicted relative reductions in the sum of squares are at most ftol.
  • info = 2 relative error between two consecutive iterates is at most xtol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 the cosine of the angle between fvec and any column of the jacobian is at most gtol in absolute value.
  • info = 5 number of calls to fcn with iflag = 1 has reached maxfev.
  • info = 6 ftol is too small. no further reduction in the sum of squares is possible.
  • info = 7 xtol is too small. no further improvement in the approximate solution x is possible.
  • info = 8 gtol is too small. fvec is orthogonal to the columns of the jacobian to machine precision.
integer, intent(out) :: Nfev

an integer output variable set to the number of calls to fcn with iflag = 1.

integer, intent(out) :: Njev

an integer output variable set to the number of calls to fcn with iflag = 2.

integer, intent(out) :: Ipvt(n)

an integer output array of length n. ipvt defines a permutation matrix p such that jacp = qr, where jac is the final calculated jacobian, q is orthogonal (not stored), and r is upper triangular. column j of p is column ipvt(j) of the identity matrix.

real(kind=wp), intent(out) :: Qtf(n)

an output array of length n which contains the first n elements of the vector (q transpose)*fvec.

real(kind=wp), intent(inout) :: Wa1(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa2(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa3(n)

work array of length n.

real(kind=wp), intent(inout) :: Wa4(m)

work array of length m.

public subroutine lmstr1(fcn, m, n, x, Fvec, Fjac, Ldfjac, Tol, Info, Ipvt, Wa, Lwa)

the purpose of lmstr1 is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm which uses minimal storage. this is done by using the more general least-squares solver lmstr. the user must provide a subroutine which calculates the functions and the rows of the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_lmstr) :: fcn

user-supplied subroutine which calculates the functions and the rows of the jacobian.

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

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.

real(kind=wp), intent(out) :: Fvec(m)

an output array of length m which contains the functions evaluated at the output x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

an output n by n array. the upper triangle of fjac contains an upper triangular matrix r such that

        t     t           t
       p *(jac *jac)*p = r *r,

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 triangular part of fjac contains information generated during the computation of r.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than n which specifies the leading dimension of the array fjac.

real(kind=wp), intent(in) :: Tol

a nonnegative input variable. termination occurs when the algorithm estimates either that the relative error in the sum of squares is at most tol or that the relative error between x and the solution is at most tol.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows:

  • info = 0 improper input parameters.
  • info = 1 algorithm estimates that the relative error in the sum of squares is at most tol.
  • info = 2 algorithm estimates that the relative error between x and the solution is at most tol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 fvec is orthogonal to the columns of the jacobian to machine precision.
  • info = 5 number of calls to fcn with iflag = 1 has reached 100*(n+1).
  • info = 6 tol is too small. no further reduction in the sum of squares is possible.
  • info = 7 tol is too small. no further improvement in the approximate solution x is possible.
integer, intent(out) :: Ipvt(n)

an integer output array of length n. ipvt defines a permutation matrix p such that jacp = qr, where jac is the final calculated jacobian, q is orthogonal (not stored), and r is upper triangular. column j of p is column ipvt(j) of the identity matrix.

real(kind=wp), intent(inout) :: Wa(Lwa)

a work array of length lwa.

integer, intent(in) :: Lwa

a positive integer input variable not less than 5*n+m.

public subroutine qform(m, n, q, Ldq, Wa)

this subroutine proceeds from the computed qr factorization of an m by n matrix a to accumulate the m by m orthogonal matrix q from its factored form.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: m

a positive integer input variable set to the number of rows of a and the order of q.

integer, intent(in) :: n

a positive integer input variable set to the number of columns of a.

real(kind=wp), intent(inout) :: q(Ldq,m)

an m by m array. on input the full lower trapezoid in the first min(m,n) columns of q contains the factored form. on output q has been accumulated into a square matrix.

integer, intent(in) :: Ldq

a positive integer input variable not less than m which specifies the leading dimension of the array q.

real(kind=wp), intent(inout) :: Wa(m)

a work array of length m.

public subroutine qrfac(m, n, a, Lda, Pivot, Ipvt, Lipvt, Rdiag, Acnorm, Wa)

this subroutine uses householder transformations with column pivoting (optional) to compute a qr factorization of the m by n matrix a. that is, qrfac determines an orthogonal matrix q, a permutation matrix p, and an upper trapezoidal matrix r with diagonal elements of nonincreasing magnitude, such that ap = qr. the householder transformation for column k, k = 1,2,...,min(m,n), is of the form

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: m

a positive integer input variable set to the number of rows of a.

integer, intent(in) :: n

a positive integer input variable set to the number of columns of a.

real(kind=wp), intent(inout) :: a(Lda,n)

an m by n array. on input a contains the matrix for which the qr factorization is to be computed. on output the strict upper trapezoidal part of a contains the strict upper trapezoidal part of r, and the lower trapezoidal part of a contains a factored form of q (the non-trivial elements of the u vectors described above).

integer, intent(in) :: Lda

a positive integer input variable not less than m which specifies the leading dimension of the array a.

logical, intent(in) :: Pivot

a logical input variable. if pivot is set true, then column pivoting is enforced. if pivot is set false, then no column pivoting is done.

integer, intent(out) :: Ipvt(Lipvt)

an integer output array of length lipvt. ipvt defines the permutation matrix p such that ap = qr. column j of p is column ipvt(j) of the identity matrix. if pivot is false, ipvt is not referenced.

integer, intent(in) :: Lipvt

a positive integer input variable. if pivot is false, then lipvt may be as small as 1. if pivot is true, then lipvt must be at least n.

real(kind=wp), intent(out) :: Rdiag(n)

an output array of length n which contains the diagonal elements of r.

real(kind=wp), intent(out) :: Acnorm(n)

an output array of length n which contains the norms of the corresponding columns of the input matrix a. if this information is not needed, then acnorm can coincide with rdiag.

real(kind=wp), intent(inout) :: Wa(n)

a work array of length n. if pivot is false, then wa can coincide with rdiag.

public subroutine qrsolv(n, r, Ldr, Ipvt, Diag, Qtb, x, Sdiag, Wa)

given an m by n matrix a, an n by n diagonal matrix d, and an m-vector b, the problem is to determine an x which solves the system

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

a positive integer input variable set to the order of r.

real(kind=wp), intent(inout) :: r(Ldr,n)

an n by n array. on input the full upper triangle must contain the full upper triangle of the matrix r. on output the full upper triangle is unaltered, and the strict lower triangle contains the strict upper triangle (transposed) of the upper triangular matrix s.

integer, intent(in) :: Ldr

a positive integer input variable not less than n which specifies the leading dimension of the array r.

integer, intent(in) :: Ipvt(n)

an integer input array of length n which defines the permutation matrix p such that ap = qr. column j of p is column ipvt(j) of the identity matrix.

real(kind=wp), intent(in) :: Diag(n)

an input array of length n which must contain the diagonal elements of the matrix d.

real(kind=wp), intent(in) :: Qtb(n)

an input array of length n which must contain the first n elements of the vector (q transpose)*b.

real(kind=wp), intent(out) :: x(n)

an output array of length n which contains the least squares solution of the system ax = b, dx = 0.

real(kind=wp), intent(out) :: Sdiag(n)

an output array of length n which contains the diagonal elements of the upper triangular matrix s.

real(kind=wp), intent(inout) :: Wa(n)

a work array of length n.

public subroutine r1mpyq(m, n, a, Lda, v, w)

given an m by n matrix a, this subroutine computes aq where q is the product of 2(n - 1) transformations

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: m

a positive integer input variable set to the number of rows of a.

integer, intent(in) :: n

a positive integer input variable set to the number of columns of a.

real(kind=wp), intent(inout) :: a(Lda,n)

an m by n array. on input a must contain the matrix to be postmultiplied by the orthogonal matrix q described above. on output a*q has replaced a.

integer, intent(in) :: Lda

a positive integer input variable not less than m which specifies the leading dimension of the array a.

real(kind=wp), intent(in) :: v(n)

an input array of length n. v(i) must contain the information necessary to recover the givens rotation gv(i) described above.

real(kind=wp), intent(in) :: w(n)

an input array of length n. w(i) must contain the information necessary to recover the givens rotation gw(i) described above.

public subroutine r1updt(m, n, s, Ls, u, v, w, Sing)

given an m by n lower trapezoidal matrix s, an m-vector u, and an n-vector v, the problem is to determine an orthogonal matrix q such that

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: m

a positive integer input variable set to the number of rows of s.

integer, intent(in) :: n

a positive integer input variable set to the number of columns of s. n must not exceed m.

real(kind=wp), intent(inout) :: s(Ls)

an array of length ls. on input s must contain the lower trapezoidal matrix s stored by columns. on output s contains the lower trapezoidal matrix produced as described above.

integer, intent(in) :: Ls

a positive integer input variable not less than (n(2m-n+1))/2.

real(kind=wp), intent(in) :: u(m)

an input array of length m which must contain the vector u.

real(kind=wp), intent(inout) :: v(n)

an array of length n. on input v must contain the vector v. on output v(i) contains the information necessary to recover the givens rotation gv(i) described above.

real(kind=wp), intent(out) :: w(m)

an output array of length m. w(i) contains information necessary to recover the givens rotation gw(i) described above.

logical, intent(out) :: Sing

a logical output variable. sing is set true if any of the diagonal elements of the output s are zero. otherwise sing is set false.

public subroutine rwupdt(n, r, Ldr, w, b, Alpha, Cos, Sin)

given an n by n upper triangular matrix r, this subroutine computes the qr decomposition of the matrix formed when a row is added to r. if the row is specified by the vector w, then rwupdt determines an orthogonal matrix q such that when the n+1 by n matrix composed of r augmented by w is premultiplied by (q transpose), the resulting matrix is upper trapezoidal. the matrix (q transpose) is the product of n transformations

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: n

a positive integer input variable set to the order of r.

real(kind=wp), intent(inout) :: r(Ldr,n)

an n by n array. on input the upper triangular part of r must contain the matrix to be updated. on output r contains the updated triangular matrix.

integer, intent(in) :: Ldr

a positive integer input variable not less than n which specifies the leading dimension of the array r.

real(kind=wp), intent(in) :: w(n)

an input array of length n which must contain the row vector to be added to r.

real(kind=wp), intent(inout) :: b(n)

an array of length n. on input b must contain the first n elements of the vector c. on output b contains the first n elements of the vector (q transpose)*c.

real(kind=wp), intent(inout) :: Alpha

a variable. on input alpha must contain the (n+1)-st element of the vector c. on output alpha contains the (n+1)-st element of the vector (q transpose)*c.

real(kind=wp), intent(out) :: Cos(n)

an output array of length n which contains the cosines of the transforming givens rotations.

real(kind=wp), intent(out) :: Sin(n)

an output array of length n which contains the sines of the transforming givens rotations.