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).
*** 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.
function for hybrd.
calculate the functions at x
and
return this vector in fvec
.
Type | Intent | Optional | 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 |
function for hybrj
Type | Intent | Optional | 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 |
Replacement for the original Minpack routine.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | i |
given an n-vector x, this function calculates the euclidean norm of x.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n |
size of |
||
real(kind=wp), | intent(in), | dimension(n) | :: | x |
input array |
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.
Type | Intent | Optional | 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) |
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.
Type | Intent | Optional | 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) |
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.
Type | Intent | Optional | 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 |
||
real(kind=wp), | intent(out) | :: | fvec(n) |
an output array of length |
||
real(kind=wp), | intent(in) | :: | xtol |
a nonnegative input variable. termination
occurs when the relative error between two consecutive
iterates is at most |
||
integer, | intent(in) | :: | maxfev |
a positive integer input variable. termination
occurs when the number of calls to |
||
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
|
||
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
|
||
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 |
||
real(kind=wp), | intent(inout) | :: | diag(n) |
an array of length |
||
integer, | intent(in) | :: | mode |
if |
||
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
|
||
integer, | intent(in) | :: | nprint |
an integer input variable that enables controlled
printing of iterates if it is positive. in this case,
|
||
integer, | intent(out) | :: | info |
an integer output variable. if the user has
terminated execution, |
||
integer, | intent(out) | :: | nfev |
output variable set to the number of calls to |
||
real(kind=wp), | intent(out) | :: | fjac(ldfjac,n) |
array which contains the
orthogonal matrix |
||
integer, | intent(in) | :: | ldfjac |
a positive integer input variable not less than |
||
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 |
||
real(kind=wp), | intent(out) | :: | qtf(n) |
an output array of length |
||
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 |
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.
Type | Intent | Optional | 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 |
|
real(kind=wp), | intent(out), | dimension(n) | :: | fvec |
an output array of length |
|
real(kind=wp), | intent(in) | :: | tol |
a nonnegative input variable. termination occurs
when the algorithm estimates that the relative error
between |
||
integer, | intent(out) | :: | info |
an integer output variable. if the user has
terminated execution, info is set to the (negative)
value of |
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.
Type | Intent | Optional | 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) |
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.
Type | Intent | Optional | 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 |
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.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer | :: | m | ||||
integer | :: | n | ||||
real(kind=wp) | :: | q(ldq,m) | ||||
integer | :: | ldq | ||||
real(kind=wp) | :: | wa(m) |
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
Type | Intent | Optional | 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) |
given an m by n matrix a, this subroutine computes aq where q is the product of 2(n - 1) transformations
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer | :: | m | ||||
integer | :: | n | ||||
real(kind=wp) | :: | a(lda,n) | ||||
integer | :: | lda | ||||
real(kind=wp) | :: | v(n) | ||||
real(kind=wp) | :: | w(n) |
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
Type | Intent | Optional | 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 |