minpack_module Module

Minpack routines for solving a set of nonlinear equations. The two main routines here are hybrj (user-provided Jacobian) and hybrd (estimates the Jacobian using finite differences).

License

*** Original Minpack License ***

 Minpack Copyright Notice (1999) University of Chicago.  All rights reserved

 Redistribution and use in source and binary forms, with or
 without modification, are permitted provided that the
 following conditions are met:

 1. Redistributions of source code must retain the above
 copyright notice, this list of conditions and the following
 disclaimer.

 2. Redistributions in binary form must reproduce the above
 copyright notice, this list of conditions and the following
 disclaimer in the documentation and/or other materials
 provided with the distribution.

 3. The end-user documentation included with the
 redistribution, if any, must include the following
 acknowledgment:

    "This product includes software developed by the
    University of Chicago, as Operator of Argonne National
    Laboratory.

 Alternately, this acknowledgment may appear in the software
 itself, if and wherever such third-party acknowledgments
 normally appear.

 4. WARRANTY DISCLAIMER. THE SOFTWARE IS SUPPLIED "AS IS"
 WITHOUT WARRANTY OF ANY KIND. THE COPYRIGHT HOLDER, THE
 UNITED STATES, THE UNITED STATES DEPARTMENT OF ENERGY, AND
 THEIR EMPLOYEES: (1) DISCLAIM ANY WARRANTIES, EXPRESS OR
 IMPLIED, INCLUDING BUT NOT LIMITED TO ANY IMPLIED WARRANTIES
 OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE, TITLE
 OR NON-INFRINGEMENT, (2) DO NOT ASSUME ANY LEGAL LIABILITY
 OR RESPONSIBILITY FOR THE ACCURACY, COMPLETENESS, OR
 USEFULNESS OF THE SOFTWARE, (3) DO NOT REPRESENT THAT USE OF
 THE SOFTWARE WOULD NOT INFRINGE PRIVATELY OWNED RIGHTS, (4)
 DO NOT WARRANT THAT THE SOFTWARE WILL FUNCTION
 UNINTERRUPTED, THAT IT IS ERROR-FREE OR THAT ANY ERRORS WILL
 BE CORRECTED.

 5. LIMITATION OF LIABILITY. IN NO EVENT WILL THE COPYRIGHT
 HOLDER, THE UNITED STATES, THE UNITED STATES DEPARTMENT OF
 ENERGY, OR THEIR EMPLOYEES: BE LIABLE FOR ANY INDIRECT,
 INCIDENTAL, CONSEQUENTIAL, SPECIAL OR PUNITIVE DAMAGES OF
 ANY KIND OR NATURE, INCLUDING BUT NOT LIMITED TO LOSS OF
 PROFITS OR LOSS OF DATA, FOR ANY REASON WHATSOEVER, WHETHER
 SUCH LIABILITY IS ASSERTED ON THE BASIS OF CONTRACT, TORT
 (INCLUDING NEGLIGENCE OR STRICT LIABILITY), OR OTHERWISE,
 EVEN IF ANY OF SAID PARTIES HAS BEEN WARNED OF THE
 POSSIBILITY OF SUCH LOSS OR DAMAGES.

*** Modifications ***

Modifications for the Fortran Astrodynamics Toolkit are covered under the following license.

History

  • Argonne National Laboratory. minpack project. march 1980. burton s. garbow, kenneth e. hillstrom, jorge j. more, john l. nazareth
  • Jacob Williams, Jan 2016, extensive refactoring into modern Fortran.

Uses

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

Used by

  • module~~minpack_module~~UsedByGraph module~minpack_module minpack_module module~fortran_astrodynamics_toolkit fortran_astrodynamics_toolkit module~fortran_astrodynamics_toolkit->module~minpack_module proc~halo_to_rv_diffcorr halo_orbit_module::halo_to_rv_diffcorr proc~halo_to_rv_diffcorr->module~minpack_module

Abstract Interfaces

abstract interface

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

    function for hybrd. calculate the functions at x and return this vector in fvec.

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: n
    real(kind=wp), intent(in), dimension(n) :: x
    real(kind=wp), intent(out), dimension(n) :: fvec
    integer, intent(inout) :: iflag

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

abstract interface

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

    function for hybrj

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: n
    real(kind=wp), intent(in), dimension(n) :: x
    real(kind=wp), intent(out), dimension(n) :: fvec
    real(kind=wp), intent(out), dimension(ldfjac,n) :: fjac
    integer, intent(in) :: ldfjac
    integer, intent(inout) :: iflag

Functions

public function dpmpar(i)

Replacement for the original Minpack routine.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: i

Return Value real(kind=wp)

public function enorm(n, x)

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

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

size of x

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

input array

Return Value real(kind=wp)


Subroutines

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

Type IntentOptional Attributes Name
integer :: n
real(kind=wp) :: r(lr)
integer :: lr
real(kind=wp) :: diag(n)
real(kind=wp) :: qtb(n)
real(kind=wp) :: delta
real(kind=wp) :: x(n)
real(kind=wp) :: wa1(n)
real(kind=wp) :: wa2(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.

Read more…

Arguments

Type IntentOptional Attributes Name
procedure(fcn_hybrd) :: fcn
integer :: n
real(kind=wp) :: x(n)
real(kind=wp) :: fvec(n)
real(kind=wp) :: fjac(ldfjac,n)
integer :: ldfjac
integer :: iflag
integer :: ml
integer :: mu
real(kind=wp) :: epsfcn
real(kind=wp) :: wa1(n)
real(kind=wp) :: wa2(n)

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.

Read more…

Arguments

Type IntentOptional Attributes Name
procedure(fcn_hybrd) :: 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)

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

Type IntentOptional Attributes Name
procedure(fcn_hybrd) :: 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.

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.

Read more…

Arguments

Type IntentOptional Attributes Name
procedure(fcn_hybrj) :: fcn
integer :: n
real(kind=wp) :: x(n)
real(kind=wp) :: fvec(n)
real(kind=wp) :: fjac(ldfjac,n)
integer :: ldfjac
real(kind=wp) :: xtol
integer :: maxfev
real(kind=wp) :: diag(n)
integer :: mode
real(kind=wp) :: factor
integer :: nprint
integer :: info
integer :: nfev
integer :: njev
real(kind=wp) :: r(lr)
integer :: lr
real(kind=wp) :: qtf(n)
real(kind=wp) :: wa1(n)
real(kind=wp) :: wa2(n)
real(kind=wp) :: wa3(n)
real(kind=wp) :: wa4(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.

Read more…

Arguments

Type IntentOptional Attributes Name
procedure(fcn_hybrj) :: fcn
integer :: n
real(kind=wp) :: x(n)
real(kind=wp) :: fvec(n)
real(kind=wp) :: fjac(ldfjac,n)
integer :: ldfjac
real(kind=wp) :: tol
integer :: info
real(kind=wp) :: wa(lwa)
integer :: lwa

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.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: m
integer :: n
real(kind=wp) :: q(ldq,m)
integer :: ldq
real(kind=wp) :: wa(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

Type IntentOptional Attributes Name
integer :: m
integer :: n
real(kind=wp) :: a(lda,n)
integer :: lda
logical :: pivot
integer :: ipvt(lipvt)
integer :: lipvt
real(kind=wp) :: rdiag(n)
real(kind=wp) :: acnorm(n)
real(kind=wp) :: wa(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

Type IntentOptional Attributes Name
integer :: m
integer :: n
real(kind=wp) :: a(lda,n)
integer :: lda
real(kind=wp) :: v(n)
real(kind=wp) :: w(n)

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

Type IntentOptional Attributes Name
integer :: m
integer :: n
real(kind=wp) :: s(ls)
integer :: ls
real(kind=wp) :: u(m)
real(kind=wp) :: v(n)
real(kind=wp) :: w(m)
logical :: sing