Sophie

Sophie

distrib > Fedora > 18 > x86_64 > by-pkgid > ed9405268ea90acd52c7b758a571e546 > files > 7

cminpack-devel-1.3.1-1.fc18.i686.rpm

<HTML><HEAD><TITLE>Manpage of HYBRD_</TITLE>
</HEAD><BODY>
<H1>HYBRD_</H1>
Section: C Library Functions (3)<BR>Updated: March 8, 2002<BR><A HREF="#index">Index</A>
<A HREF="index.html">Return to Main Contents</A><HR>


<A NAME="lbAB">&nbsp;</A>
<H2>NAME</H2>

hybrd_, hybrd1_ - find a zero of a system of nonlinear function
<A NAME="lbAC">&nbsp;</A>
<H2>SYNOPSIS</H2>

<B>#include &lt;<A HREF="file:/usr/include/minpack.h">minpack.h</A>&gt;</B>

<P>


<DL COMPACT>
<DT>
<B>void hybrd1_ ( </B>

<B>void (*</B><I>fcn</I><B>)(</B>

<B>int *</B><I>n</I><B>,</B>

<B>double *</B><I>x</I><B>,</B>

<DD>
<B>double *</B><I>fvec</I><B>,</B>

<B>int *</B><I>iflag</I><B>),</B>

<DL COMPACT><DT><DD>
<B>int *</B><I>n</I><B>,</B>

<B>double *</B><I>x</I><B>,</B>

<B>double *</B><I>fvec</I><B>,</B>

<DD>
<B>double *</B><I>tol</I><B>,</B>

<B>int *</B><I>info</I><B>,</B>

<B>double *</B><I>wa</I><B>,</B>

<DD>
<B>int *</B><I>lwa</I><B>);</B>

</DL>

<DT>
<B>void hybrd_</B>

<B>( void (*</B><I>fcn</I><B>)(</B>

<B>int * </B><I>n</I><B>,</B>

<B>double *</B><I>x</I><B>,</B>

<DD>
<B>double *</B><I>fvec</I><B>,</B>

<B>int *</B><I>iflag</I><B>),</B>

<DL COMPACT><DT><DD>
<B>int *</B><I>n</I><B>,</B>

<B>double *</B><I>x</I><B>,</B>

<B>double *</B><I>fvec</I><B>,</B>

<DD>
<B>double *</B><I>xtol</I><B>,</B>

<B>int *</B><I>maxfev</I><B>,</B>

<B>int *</B><I>ml</I><B>,</B>

<B>int *</B><I>mu</I><B>,</B>

<DD>
<B>double *</B><I>epsfcn</I><B>,</B>

<B>double *</B><I>diag</I><B>,</B>

<B>int *</B><I>mode</I><B>,</B>

<B>double *</B><I>factor</I><B>,</B>

<B>int *</B><I>nprint</I><B>,</B>

<B>int *</B><I>info</I><B>,</B>

<DD>
<B>int *</B><I>nfev</I><B>,</B>

<B>double *</B><I>fjac</I><B>,</B>

<B>int *</B><I>ldfjac</I><B>,</B>

<DD>
<B>double *</B><I>r</I><B>,</B>

<B>int *</B><I>lr</I><B>,</B>

<B>double *</B><I>qtf</I><B>,</B>

<DD>
<B>double *</B><I>wa1</I><B>,</B>

<B>double *</B><I>wa2</I><B>,</B>

<B>double *</B><I>wa3</I><B>,</B>

<B>double *</B><I>wa4</I><B>);</B>

</DL>



<DD>
</DL>
<A NAME="lbAD">&nbsp;</A>
<H2>DESCRIPTION</H2>

<DD>The purpose of <B>hybrd_</B> is to find a zero of a system of
<I>n</I> nonlinear functions in <I>n</I> 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.
<P>

<B>hybrd1_</B> serves the same function but has a simplified calling
sequence.
<BR>

<A NAME="lbAE">&nbsp;</A>
<H3>Language notes</H3>

<B>hybrd_</B> and <B>hybrd1_</B> are written in FORTRAN. If calling from
C, keep these points in mind:
<DL COMPACT>
<DT>Name mangling.<DD>
With <B>g77</B> version 2.95 or 3.0, all the function names end in an
underscore.  This may change with future versions of <B>g77</B>.
<DT>Compile with <B>g77</B>.<DD>
Even if your program is all C code, you should link with <B>g77</B>
so it will pull in the FORTRAN libraries automatically.  It's easiest
just to use <B>g77</B> to do all the compiling.  (It handles C just fine.)
<DT>Call by reference.<DD>
All function parameters must be pointers.
<DT>Column-major arrays.<DD>
Suppose a function returns an array with 5 rows and 3 columns in an
array <I>z</I> and in the call you have declared a leading dimension of
7.  The FORTRAN and equivalent C references are:
<P>
<PRE>
        z(1,1)          z[0]
        z(2,1)          z[1]
        z(5,1)          z[4]
        z(1,2)          z[7]
        z(1,3)          z[14]
        z(i,j)          z[(i-1) + (j-1)*7]
</PRE>

</DL>
<A NAME="lbAF">&nbsp;</A>
<H3>Parameters for both functions</H3>

<I>fcn</I> is the name of the user-supplied subroutine which calculates
the functions. In FORTRAN, <I>fcn</I> must be declared in an external
statement in the user calling program, and should be written as
follows:
<P>
<PRE>
subroutine fcn(n,x,fvec,iflag)
integer n,iflag
double precision x(n),fvec(n)
----------
calculate the functions at x and
return this vector in fvec.
---------
return
end
</PRE>

<P>
<P>
In C, <I>fcn</I> should be written as follows:
<P>
<PRE>
  void fcn(int n, double *x, double *fvec, int *iflag)
  {
    /* calculate the functions at x and
       return this vector in fvec. */
  }
</PRE>

<P>
The value of <I>iflag</I> should not be changed by <I>fcn</I> unless
the user wants to terminate execution of <B>hybrd_</B>.
In this case set <I>iflag</I> to a negative integer.
<P>
<I>n</I> is a positive integer input variable set to the number
of functions and variables.
<P>
<I>x</I> is an array of length <I>n</I>. On input <I>x</I> must contain
an initial estimate of the solution vector. On output <I>x</I>
contains the final estimate of the solution vector.
<P>
<I>fvec</I> is an output array of length <I>n</I> which contains
the functions evaluated at the output <I>x</I>.
<BR>

<A NAME="lbAG">&nbsp;</A>
<H3>Parameters for <B>hybrd1_</B></H3>

<P>
<I>tol</I> is a nonnegative input variable. Termination occurs
when the algorithm estimates that the relative error
between <I>x</I> and the solution is at most <I>tol</I>.
<P>
<I>info</I> is an integer output variable. If the user has
terminated execution, <I>info</I> is set to the (negative)
value of <I>iflag</I>. See description of <I>fcn</I>. Otherwise,
<I>info</I> is set as follows.
<P>
<I>info</I> = 0   improper input parameters.
<P>
<I>info</I> = 1   algorithm estimates that the relative error
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;between&nbsp;<I>x</I>&nbsp;and&nbsp;the&nbsp;solution&nbsp;is&nbsp;at&nbsp;most&nbsp;<I>tol</I>.
<P>
<I>info</I> = 2   number of calls to fcn has reached or exceeded
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;200*(<I>n</I>+1).
<P>
<I>info</I> = 3   <I>tol</I> is too small. No further improvement in
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;the&nbsp;approximate&nbsp;solution&nbsp;<I>x</I>&nbsp;is&nbsp;possible.
<P>
<I>info</I> = 4   iteration is not making good progress.
<P>
<I>wa</I> is a work array of length <I>lwa</I>.
<P>
<I>lwa</I> is a positive integer input variable not less than
(<I>n</I>*(3*<I>n</I>+13))/2.
<BR>

<A NAME="lbAH">&nbsp;</A>
<H3>Parameters for <B>hybrd_</B></H3>

<P>
<I>xtol</I> is a nonnegative input variable. Termination
occurs when the relative error between two consecutive
iterates is at most <I>xtol</I>.
<P>
<I>maxfev</I> is a positive integer input variable. Termination
occurs when the number of calls to <I>fcn</I> is at least <I>maxfev</I>
by the end of an iteration.
<P>
<I>ml</I> is 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
<I>ml</I> to at least <I>n</I> - 1.
<P>
<I>mu</I> is 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 <I>n</I> - 1.
<P>
<I>epsfcn</I> is 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 <I>epsfcn</I>. If <I>epsfcn</I> is less
than the machine precision, it is assumed that the relative
errors in the functions are of the order of the machine
precision.
<P>
<I>diag</I> is an array of length <I>n</I>. If <I>mode</I> = 1 (see
below), <I>diag</I> is internally set. If <I>mode</I> = 2, <I>diag</I>
must contain positive entries that serve as
multiplicative scale factors for the variables.
<P>
<I>mode</I> is an integer input variable. If <I>mode</I> = 1, the
variables will be scaled internally. If <I>mode</I> = 2,
the scaling is specified by the input <I>diag</I>. Other
values of mode are equivalent to <I>mode</I> = 1.
<P>
<I>factor</I> is a positive input variable used in determining the
initial step bound. This bound is set to the product of
<I>factor</I> and the euclidean norm of diag*x if nonzero, or else
to <I>factor</I> itself. In most cases factor should lie in the
interval (.1,100.). 100. Is a generally recommended value.
<P>
<I>nprint</I> is an integer input variable that enables controlled
printing of iterates if it is positive. In this case,
<I>fcn</I> is called with <I>iflag</I> = 0 at the beginning of the first
iteration and every nprint iterations thereafter and
immediately prior to return, with <I>x</I> and <I>fvec</I> available
for printing. If <I>nprint</I> is not positive, no special calls
of <I>fcn</I> with <I>iflag</I> = 0 are made.
<P>
<I>info</I> is an integer output variable. If the user has
terminated execution, <I>info</I> is set to the (negative)
value of <I>iflag</I>. See description of <I>fcn</I>. Otherwise,
<I>info</I> is set as follows.
<P>
<I>info</I> = 0   improper input parameters.
<P>
<I>info</I> = 1   relative error between two consecutive iterates
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;is&nbsp;at&nbsp;most&nbsp;<I>xtol</I>.
<P>
<I>info</I> = 2   number of calls to <I>fcn</I> has reached or exceeded
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;<I>maxfev</I>.
<P>
<I>info</I> = 3   <I>xtol</I> is too small. No further improvement in
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;the&nbsp;approximate&nbsp;solution&nbsp;<I>x</I>&nbsp;is&nbsp;possible.
<P>
<I>info</I> = 4   iteration is not making good progress, as
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;measured&nbsp;by&nbsp;the&nbsp;improvement&nbsp;from&nbsp;the&nbsp;last
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;five&nbsp;jacobian&nbsp;evaluations.
<P>
<I>info</I> = 5   iteration is not making good progress, as
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;measured&nbsp;by&nbsp;the&nbsp;improvement&nbsp;from&nbsp;the&nbsp;last
<BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;ten&nbsp;iterations.
<P>
<I>nfev</I> is an integer output variable set to the number of
calls to <I>fcn</I>.
<P>
<I>fjac</I> is an output <I>n</I> by <I>n</I> array which contains the
orthogonal matrix <I>q</I> produced by the <I>qr</I> factorization
of the final approximate jacobian.
<P>
<I>ldfjac</I> is a positive integer input variable not less than <I>n</I>
which specifies the leading dimension of the array <I>fjac</I>.
<P>
<I>r</I> is an output array of length <I>lr</I> which contains the
upper triangular matrix produced by the <I>qr</I> factorization
of the final approximate Jacobian, stored rowwise.
<P>
<I>lr</I> is a positive integer input variable not less than
(<I>n</I>*(<I>n</I>+1))/2.
<P>
<I>qtf</I> is an output array of length <I>n</I> which contains
the vector (q transpose)*<I>fvec</I>.
<P>
<I>wa1</I>, <I>wa2</I>, <I>wa3</I>, and <I>wa4</I> are work arrays of length <I>n</I>.
<P>
<A NAME="lbAI">&nbsp;</A>
<H2>SEE ALSO</H2>

<B><A HREF="hybrj_.html">hybrj</A></B>(3),

<B><A HREF="hybrj_.html">hybrj1</A></B>(3).

<BR>

<P>
<A NAME="lbAJ">&nbsp;</A>
<H2>AUTHORS</H2>

Burton S. Garbow, Kenneth E. Hillstrom, Jorge J. More.
<BR>

This manual page was written by Jim Van Zandt &lt;<A HREF="mailto:jrv@debian.org">jrv@debian.org</A>&gt;,
for the Debian GNU/Linux system (but may be used by others).
<P>

<HR>
<A NAME="index">&nbsp;</A><H2>Index</H2>
<DL>
<DT><A HREF="#lbAB">NAME</A><DD>
<DT><A HREF="#lbAC">SYNOPSIS</A><DD>
<DT><A HREF="#lbAD">DESCRIPTION</A><DD>
<DL>
<DT><A HREF="#lbAE">Language notes</A><DD>
<DT><A HREF="#lbAF">Parameters for both functions</A><DD>
<DT><A HREF="#lbAG">Parameters for <B>hybrd1_</B></A><DD>
<DT><A HREF="#lbAH">Parameters for <B>hybrd_</B></A><DD>
</DL>
<DT><A HREF="#lbAI">SEE ALSO</A><DD>
<DT><A HREF="#lbAJ">AUTHORS</A><DD>
</DL>
<HR>
This document was created by
<A HREF="index.html">man2html</A>,
using the manual pages.<BR>
Time: 10:19:50 GMT, April 20, 2007
</BODY>
</HTML>