changeset 4115:fc2048d4cd21

[project @ 2002-10-22 21:28:42 by jwe]
author jwe
date Tue, 22 Oct 2002 21:28:42 +0000
parents a32457362437
children 644c20b7b9e1
files doc/interpreter/diffeq.txi src/ChangeLog src/DLD-FUNCTIONS/daspk.cc src/DLD-FUNCTIONS/dasrt.cc src/DLD-FUNCTIONS/dassl.cc src/DLD-FUNCTIONS/lsode.cc src/DLD-FUNCTIONS/odessa.cc src/help.cc src/version.h
diffstat 9 files changed, 474 insertions(+), 100 deletions(-) [+]
line wrap: on
line diff
--- a/doc/interpreter/diffeq.txi	Thu Oct 17 21:39:48 2002 +0000
+++ b/doc/interpreter/diffeq.txi	Tue Oct 22 21:28:42 2002 +0000
@@ -42,6 +42,8 @@
 
 @DOCSTRING(lsode)
 
+@DOCSTRING(lsode_options)
+
 Here is an example of solving a set of three differential equations using
 @code{lsode}.  Given the function
 
@@ -85,7 +87,13 @@
 @end group
 @end example
 
-@DOCSTRING(lsode_options)
+Octave also includes @sc{Odessa}, a modification of @sc{Lsode} that
+allows for the computation of parametric sensitivity information
+simultaneously with the solution of the set of ODEs.
+
+@DOCSTRING(odessa)
+
+@DOCSTRING(odessa_options)
 
 See Alan C. Hindmarsh, @cite{ODEPACK, A Systematized Collection of ODE
 Solvers}, in Scientific Computing, R. S. Stepleman, editor, (1983) for
@@ -94,7 +102,7 @@
 @node Differential-Algebraic Equations,  , Ordinary Differential Equations, Differential Equations
 @section Differential-Algebraic Equations
 
-The function @code{dassl} can be used to solve DAEs of the form
+The function @code{daspk} can be used to solve DAEs of the form
 @iftex
 @tex
 $$
@@ -110,11 +118,19 @@
 @end ifinfo
 
 @noindent
-using Petzold's DAE solver @sc{Dassl}.
+using Petzold's DAE solver @sc{Daspk}.
+
+@DOCSTRING(daspk)
+
+@DOCSTRING(daspk_options)
 
-@DOCSTRING(dassl)
+Octave also includes @sc{Dassl}, an earlier version of @var{Daspk},
+and @var{dasrt}, which can be used to solve DAEs with constraints
+(stopping conditions).
 
-@DOCSTRING(dassl_options)
+@DOCSTRING(dasrt)
+
+@DOCSTRING(dasrt_options)
 
 See K. E. Brenan, et al., @cite{Numerical Solution of Initial-Value
 Problems in Differential-Algebraic Equations}, North-Holland (1989) for
--- a/src/ChangeLog	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/ChangeLog	Tue Oct 22 21:28:42 2002 +0000
@@ -1,3 +1,9 @@
+2002-10-22  John W. Eaton  <jwe@bevo.che.wisc.edu>
+
+	* version.h (OCTAVE_VERSION): Now 2.1.37.
+	(OCTAVE_CONTRIB_STATEMENT): New macro.
+	(OCTAVE_STARTUP_MESSAGE): Use it.
+
 2002-10-17  John W. Eaton  <jwe@bevo.che.wisc.edu>
 
 	* data.cc (fill_matrix): If nargin is zero, use val, not 0.0.
--- a/src/DLD-FUNCTIONS/daspk.cc	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/DLD-FUNCTIONS/daspk.cc	Tue Oct 22 21:28:42 2002 +0000
@@ -140,13 +140,33 @@
 
 DEFUN_DLD (daspk, args, nargout,
   "-*- texinfo -*-\n\
-@deftypefn {Loadable Function} {[@var{x}, @var{xdot}] =} daspk (@var{fcn}, @var{x0}, @var{xdot0}, @var{t}, @var{t_crit})\n\
-Return a matrix of states and their first derivatives with respect to\n\
-@var{t}.  Each row in the result matrices correspond to one of the\n\
+@deftypefn {Loadable Function} {[@var{x}, @var{xdot}, @var{istate}, @var{msg}] =} daspk (@var{fcn}, @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})\n\
+Solve the set of differential-algebraic equations\n\
+@tex\n\
+$$ 0 = f (\\dot{x}, x, t) $$\n\
+with\n\
+$$ x(t_0) = x_0, \\dot{x}(t_0) = \\dot{x}_0 $$\n\
+@end tex\n\
+@ifinfo\n\
+\n\
+@example\n\
+0 = f (xdot, x, t)\n\
+@end example\n\
+\n\
+with\n\
+\n\
+@example\n\
+x(t_0) = x_0, xdot(t_0) = xdot_0\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+The solution is returned in the matrices @var{x} and @var{xdot},\n\
+with each row in the result matrices corresponding to one of the\n\
 elements in the vector @var{t}.  The first element of @var{t}\n\
-corresponds to the initial state @var{x0} and derivative @var{xdot0}, so\n\
-that the first row of the output @var{x} is @var{x0} and the first row\n\
-of the output @var{xdot} is @var{xdot0}.\n\
+should be @math{t_0} and correspond to the initial state of the\n\
+system @var{x_0} and its derivative @var{xdot_0}, so that the first\n\
+row of the output @var{x} is @var{x_0} and the first row\n\
+of the output @var{xdot} is @var{xdot_0}.\n\
 \n\
 The first argument, @var{fcn}, is a string that names the function to\n\
 call to compute the vector of residuals for the set of equations.\n\
@@ -157,27 +177,61 @@
 @end example\n\
 \n\
 @noindent\n\
-where @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\
+in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\
 scalar.\n\
 \n\
+If @var{fcn} is a two-element string array, the first element names\n\
+the function @math{f} described above, and the second element names\n\
+a function to compute the modified Jacobian
+\n\
+@tex\n\
+$$\n\
+J = {\\partial f \\over \\partial x}\n\
+  + c {\\partial f \\over \\partial \\dot{x}}\n\
+$$\n\
+@end tex\n\
+@ifinfo\n\
+      df       df\n\
+jac = -- + c ------\n\
+      dx     d xdot\n\
+@example\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+\n\
+The modified Jacobian function must have the form\n\
+\n\
+@example\n\
+\n\
+@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})\n\
+\n\
+@end example\n\
+\n\
 The second and third arguments to @code{daspk} specify the initial\n\
 condition of the states and their derivatives, and the fourth argument\n\
-specifies a vector of output times at which the solution is desired, \n\
+specifies a vector of output times at which the solution is desired,\n\
 including the time corresponding to the initial condition.\n\
 \n\
 The set of initial states and derivatives are not strictly required to\n\
-be consistent.  In practice, however, @sc{Dassl} is not very good at\n\
-determining a consistent set for you, so it is best if you ensure that\n\
-the initial values result in the function evaluating to zero.\n\
+be consistent.  If they are not consistent, you must use the\n\
+@code{daspk_options} function to provide additional information so\n\
+that @code{daspk} can compute a consistent starting point.\n\
 \n\
 The fifth argument is optional, and may be used to specify a set of\n\
 times that the DAE solver should not integrate past.  It is useful for\n\
 avoiding difficulties with singularities and points where there is a\n\
 discontinuity in the derivative.\n\
 \n\
+After a successful computation, the value of @var{istate} will be\n\
+greater than zero (consistent with the Fortran version of @sc{Daspk}).\n\
+\n\
+If the computation is not successful, the value of @var{istate} will be\n\
+less than zero and @var{msg} will contain additional information.\n\
+\n\
 You can use the function @code{daspk_options} to set optional\n\
 parameters for @code{daspk}.\n\
-@end deftypefn")
+@end deftypefn\n\
+@seealso{dassl}")
 {
   octave_value_list retval;
 
--- a/src/DLD-FUNCTIONS/dasrt.cc	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/DLD-FUNCTIONS/dasrt.cc	Tue Oct 22 21:28:42 2002 +0000
@@ -43,7 +43,7 @@
 
 #include "DASRT-opts.cc"
 
-// Global pointers for user defined function required by dassl.
+// Global pointers for user defined function required by dasrt.
 static octave_function *dasrt_f;
 static octave_function *dasrt_j;
 static octave_function *dasrt_cf;
@@ -235,79 +235,134 @@
 
 DEFUN_DLD (dasrt, args, nargout,
   "-*- texinfo -*-\n\
-@deftypefn {Loadable Function} {[@var{x}, @var{xdot}, @var{t}] =} dasrt (@var{fj} [, @var{g}], @var{x_0}, @var{xdot_0}, @var{t_out} [, @var{t_crit}])\n\
-Solve a system of differential/algebraic equations with functional\n\
-stopping criteria.\n\
+@deftypefn {Loadable Function} {[@var{x}, @var{xdot}, @var{t_out}, @var{istat}, @var{msg}] =} dasrt (@var{fcn} [, @var{g}], @var{x_0}, @var{xdot_0}, @var{t} [, @var{t_crit}])\n\
+Solve the set of differential-algebraic equations\n\
+@tex\n\
+$$ 0 = f (\\dot{x}, x, t) $$\n\
+with\n\
+$$ x(t_0) = x_0, \\dot{x}(t_0) = \\dot{x}_0 $$\n\
+@end tex\n\
+@ifinfo\n\
 \n\
-The function to be integrated must be of the form:\n\
 @example\n\
-@var{res} = f (@var{x}, @var{xdot}, @var{t}) = 0\n\
+0 = f (xdot, x, t)\n\
+@end example\n\
+\n\
+with\n\
+\n\
+@example\n\
+x(t_0) = x_0, xdot(t_0) = xdot_0\n\
 @end example\n\
 \n\
-The stopping condition must be of the form:\n\
+@end ifinfo\n\
+with functional stopping criteria (root solving).\n\
+\n\
+The solution is returned in the matrices @var{x} and @var{xdot},\n\
+with each row in the result matrices corresponding to one of the\n\
+elements in the vector @var{t_out}.  The first element of @var{t}\n\
+should be @math{t_0} and correspond to the initial state of the\n\
+system @var{x_0} and its derivative @var{xdot_0}, so that the first\n\
+row of the output @var{x} is @var{x_0} and the first row\n\
+of the output @var{xdot} is @var{xdot_0}.\n\
+\n\
+The vector @var{t} provides an upper limit on the length of the\n\
+integration.  If the stopping condition is met, the vector\n\
+@var{t_out} will be shorter than @var{t}, and the final element of\n\
+@var{t_out} will be the point at which the stopping condition was met,\n\
+and may not correspond to any element of the vector @var{t}.\n\
+\n\
+The first argument, @var{fcn}, is a string that names the function to\n\
+call to compute the vector of residuals for the set of equations.\n\
+It must have the form\n\
 \n\
 @example\n\
-@var{res} = g (@var{x}, @var{t}) = 0\n\
-@end example\n\
-\n\
-The Jacobian (if included) must be of the form:\n\
-\n\
-@example\n\
-@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{cj})\n\
-   =  df/dx + cj*df/dxdot\n\
+@var{res} = f (@var{x}, @var{xdot}, @var{t})\n\
 @end example\n\
 \n\
 @noindent\n\
-The following inputs are entered:\n\
+in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\
+scalar.\n\
+\n\
+If @var{fcn} is a two-element string array, the first element names\n\
+the function @math{f} described above, and the second element names\n\
+a function to compute the modified Jacobian
 \n\
-@table @var\n\
-@item fj\n\
-The string vector containing either @var{f} or both @var{f} and @var{j}.\n\
-\n\
-@item f\n\
-The function to be integrated.\n\
+@tex\n\
+$$\n\
+J = {\\partial f \\over \\partial x}\n\
+  + c {\\partial f \\over \\partial \\dot{x}}\n\
+$$\n\
+@end tex\n\
+@ifinfo\n\
 \n\
-@item g\n\
-The function with the stopping conditions.\n\
+@example\n\
+      df       df\n\
+jac = -- + c ------\n\
+      dx     d xdot\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+\n\
+The modified Jacobian function must have the form\n\
+\n\
+@example\n\
 \n\
-@item j\n\
-The optional Jacobian function.  If not included, it will be approximated\n\
-by finite differences.\n\
+@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})\n\
+\n\
+@end example\n\
 \n\
-@item x_0\n\
-The initial state.\n\
+The optional second argument names a function that defines the\n\
+constraint functions whose roots are desired during the integration.\n\
+This function must have the form\n\
 \n\
-@item xdot_0\n\
-The time derivative of the initial state.\n\
+@example\n\
+@var{g_out} = g (@var{x}, @var{t})\n\
+@end example\n\
 \n\
-@item t_out\n\
-The times at which the solution should be returned.  This vector should\n\
-include the initial time a the first value.\n\
+and return a vector of the constraint function values.\n\
+If the value of any of the constraint functions changes sign, @sc{Dasrt}\n\
+will attempt to stop the integration at the point of the sign change.\n\
+\n\
+If the name of the constraint function is omitted, @code{dasrt} solves\n\
+the saem problem as @code{daspk} or @code{dassl}.
 \n\
-@item t_crit\n\
-The times at which the function is non-smooth or poorly behaved.\n\
+Note that because of numerical errors in the constraint functions\n\
+due to roundoff and integration error, @sc{Dasrt} may return false\n\
+roots, or return the same root at two or more nearly equal values of\n\
+@var{T}.  If such false roots are suspected, the user should consider\n\
+smaller error tolerances or higher precision in the evaluation of the\n\
+constraint functions.\n\
 \n\
-@end table\n\
-\n\
-@noindent\n\
-The following outputs are returned:\n\
+If a root of some constraint function defines the end of the problem,\n\
+the input to @sc{Dasrt} should nevertheless allow integration to a\n\
+point slightly past that root, so that @sc{Dasrt} can locate the root\n\
+by interpolation.\n\
 \n\
-@table @var\n\
-@item x\n\
-The states at each time in @var{t}.\n\
+The third and fourth arguments to @code{dasrt} specify the initial\n\
+condition of the states and their derivatives, and the fourth argument\n\
+specifies a vector of output times at which the solution is desired,\n\
+including the time corresponding to the initial condition.\n\
+\n\
+The set of initial states and derivatives are not strictly required to\n\
+be consistent.  In practice, however, @sc{Dassl} is not very good at\n\
+determining a consistent set for you, so it is best if you ensure that\n\
+the initial values result in the function evaluating to zero.\n\
 \n\
-@item xdot\n\
-The time derivative of the states at each time in @var{t}.\n\
+The sixth argument is optional, and may be used to specify a set of\n\
+times that the DAE solver should not integrate past.  It is useful for\n\
+avoiding difficulties with singularities and points where there is a\n\
+discontinuity in the derivative.\n\
 \n\
-@item t\n\
-All the times in the requested output time vector up to the stopping\n\
-criteria.  The time at which the stopping criteria is achieved is returned\n\
-as the last time in the vector.\n\
-@end table\n\
+After a successful computation, the value of @var{istate} will be\n\
+greater than zero (consistent with the Fortran version of @sc{Dassl}).\n\
+\n\
+If the computation is not successful, the value of @var{istate} will be\n\
+less than zero and @var{msg} will contain additional information.\n\
 \n\
 You can use the function @code{dasrt_options} to set optional\n\
 parameters for @code{dasrt}.\n\
-@end deftypefn")
+@end deftypefn\n\
+@seealso{daspk, dasrt, lsode, odessa}")
 {
   octave_value_list retval;
 
--- a/src/DLD-FUNCTIONS/dassl.cc	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/DLD-FUNCTIONS/dassl.cc	Tue Oct 22 21:28:42 2002 +0000
@@ -207,13 +207,33 @@
 
 DEFUN_DLD (dassl, args, nargout,
   "-*- texinfo -*-\n\
-@deftypefn {Loadable Function} {[@var{x}, @var{xdot}] =} dassl (@var{fcn}, @var{x0}, @var{xdot0}, @var{t}, @var{t_crit})\n\
-Return a matrix of states and their first derivatives with respect to\n\
-@var{t}.  Each row in the result matrices correspond to one of the\n\
+@deftypefn {Loadable Function} {[@var{x}, @var{xdot}, @var{istate}, @var{msg}] =} dassl (@var{fcn}, @var{x_0}, @var{xdot_0}, @var{t}, @var{t_crit})\n\
+Solve the set of differential-algebraic equations\n\
+@tex\n\
+$$ 0 = f (\\dot{x}, x, t) $$\n\
+with\n\
+$$ x(t_0) = x_0, \\dot{x}(t_0) = \\dot{x}_0 $$\n\
+@end tex\n\
+@ifinfo\n\
+\n\
+@example\n\
+0 = f (xdot, x, t)\n\
+@end example\n\
+\n\
+with\n\
+\n\
+@example\n\
+x(t_0) = x_0, xdot(t_0) = xdot_0\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+The solution is returned in the matrices @var{x} and @var{xdot},\n\
+with each row in the result matrices corresponding to one of the\n\
 elements in the vector @var{t}.  The first element of @var{t}\n\
-corresponds to the initial state @var{x0} and derivative @var{xdot0}, so\n\
-that the first row of the output @var{x} is @var{x0} and the first row\n\
-of the output @var{xdot} is @var{xdot0}.\n\
+should be @math{t_0} and correspond to the initial state of the\n\
+system @var{x_0} and its derivative @var{xdot_0}, so that the first\n\
+row of the output @var{x} is @var{x_0} and the first row\n\
+of the output @var{xdot} is @var{xdot_0}.\n\
 \n\
 The first argument, @var{fcn}, is a string that names the function to\n\
 call to compute the vector of residuals for the set of equations.\n\
@@ -224,12 +244,39 @@
 @end example\n\
 \n\
 @noindent\n\
-where @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\
+in which @var{x}, @var{xdot}, and @var{res} are vectors, and @var{t} is a\n\
 scalar.\n\
 \n\
+If @var{fcn} is a two-element string array, the first element names\n\
+the function @math{f} described above, and the second element names\n\
+a function to compute the modified Jacobian
+\n\
+@tex\n\
+$$\n\
+J = {\\partial f \\over \\partial x}\n\
+  + c {\\partial f \\over \\partial \\dot{x}}\n\
+$$\n\
+@end tex\n\
+@ifinfo\n\
+      df       df\n\
+jac = -- + c ------\n\
+      dx     d xdot\n\
+@example\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+\n\
+The modified Jacobian function must have the form\n\
+\n\
+@example\n\
+\n\
+@var{jac} = j (@var{x}, @var{xdot}, @var{t}, @var{c})\n\
+\n\
+@end example\n\
+\n\
 The second and third arguments to @code{dassl} specify the initial\n\
 condition of the states and their derivatives, and the fourth argument\n\
-specifies a vector of output times at which the solution is desired, \n\
+specifies a vector of output times at which the solution is desired,\n\
 including the time corresponding to the initial condition.\n\
 \n\
 The set of initial states and derivatives are not strictly required to\n\
@@ -242,9 +289,16 @@
 avoiding difficulties with singularities and points where there is a\n\
 discontinuity in the derivative.\n\
 \n\
+After a successful computation, the value of @var{istate} will be\n\
+greater than zero (consistent with the Fortran version of @sc{Dassl}).\n\
+\n\
+If the computation is not successful, the value of @var{istate} will be\n\
+less than zero and @var{msg} will contain additional information.\n\
+\n\
 You can use the function @code{dassl_options} to set optional\n\
 parameters for @code{dassl}.\n\
-@end deftypefn")
+@end deftypefn\n\
+@seealso{daspk, dasrt, lsode, odessa}")
 {
   octave_value_list retval;
 
--- a/src/DLD-FUNCTIONS/lsode.cc	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/DLD-FUNCTIONS/lsode.cc	Tue Oct 22 21:28:42 2002 +0000
@@ -160,32 +160,111 @@
 
 DEFUN_DLD (lsode, args, nargout,
   "-*- texinfo -*-\n\
-@deftypefn {Loadable Function} {[@var{x}, @var{istate}, @var{msg}]} lsode (@var{fcn}, @var{x0}, @var{t}, @var{t_crit})\n\
-Return a matrix of @var{x} as a function of @var{t}, given the initial\n\
-state of the system @var{x0}.  Each row in the result matrix corresponds\n\
-to one of the elements in the vector @var{t}.  The first element of\n\
-@var{t} corresponds to the initial state @var{x0}, so that the first row\n\
-of the output is @var{x0}.\n\
+@deftypefn {Loadable Function} {[@var{x}, @var{istate}, @var{msg}]} lsode (@var{fcn}, @var{x_0}, @var{t}, @var{t_crit})\n\
+Solve the set of differential equations\n\
+@tex\n\
+$$ {dx \\over dt} = f (x, t) $$\n\
+with\n\
+$$ x(t_0) = x_0 $$\n\
+@end tex\n\
+@ifinfo\n\
+\n\
+@example\n\
+dx\n\
+-- = f(x, t)\n\
+dt
+@end example\n\
+\n\
+with\n\
+\n\
+@example\n\
+x(t_0) = x_0\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+The solution is returned in the matrix @var{x}, with each row\n\
+corresponding to an element of the vector @var{t}.  The first element\n\
+of @var{t} should be @math{t_0} and should correspond to the initial\n\
+state of the system @var{x_0}, so that the first row of the output\n\
+is @var{x_0}.\n\
 \n\
 The first argument, @var{fcn}, is a string that names the function to\n\
 call to compute the vector of right hand sides for the set of equations.\n\
-It must have the form\n\
+The function must have the form\n\
 \n\
 @example\n\
 @var{xdot} = f (@var{x}, @var{t})\n\
 @end example\n\
 \n\
 @noindent\n\
-where @var{xdot} and @var{x} are vectors and @var{t} is a scalar.\n\
+in which @var{xdot} and @var{x} are vectors and @var{t} is a scalar.\n\
+\n\
+If @var{fcn} is a two-element string array, the first element names the\n\
+function @math{f} described above, and the second element names a function\n\
+to compute the Jacobian of @math{f}.  The Jacobian function must have the\n\
+form\n\
+\n\
+@example\n\
+@var{jac} = j (@var{x}, @var{t})\n\
+@end example\n\
+\n\
+in which @var{jac} is the matrix of partial derivatives\n\
+@tex\n\
+$$ J = {\\partial f_i \\over \\partial x_j} = \\left[\\matrix{\n\
+{\\partial f_1 \\over \\partial x_1}\n\
+  & {\\partial f_1 \\over \\partial x_2}\n\
+  & \\cdots\n\
+  & {\\partial f_1 \\over \\partial x_N} \\cr\n\
+{\\partial f_2 \\over \\partial x_1}\n\
+  & {\\partial f_2 \\over \\partial x_2}\n\
+  & \\cdots\n\
+  & {\\partial f_2 \\over \\partial x_N} \\cr\n\
+ \\vdots & \\vdots & \\ddots & \\vdots \\cr\n\
+{\\partial f_3 \\over \\partial x_1}\n\
+  & {\\partial f_3 \\over \\partial x_2}\n\
+  & \\cdots\n\
+  & {\\partial f_3 \\over \\partial x_N} \\cr}\\right]$$\n\
+@end tex\n\
+@ifinfo\n\
+\n\
+@example\n\
+             | df_1  df_1       df_1 |\n\
+             | ----  ----  ...  ---- |\n\
+             | dx_1  dx_2       dx_N |\n\
+             |                       |\n\
+             | df_2  df_2       df_2 |\n\
+             | ----  ----  ...  ---- |\n\
+      df_i   | dx_1  dx_2       dx_N |\n\
+jac = ---- = |                       |\n\
+      dx_j   |  .    .     .    .    |\n\
+             |  .    .      .   .    |\n\
+             |  .    .       .  .    |\n\
+             |                       |\n\
+             | df_N  df_N       df_N |\n\
+             | ----  ----  ...  ---- |\n\
+             | dx_1  dx_2       dx_N |\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+\n\
+The second and third arguments specify the intial state of the system,\n\
+@math{x_0}, and the initial value of the independent variable @math{t_0}.\n\
 \n\
 The fourth argument is optional, and may be used to specify a set of\n\
 times that the ODE solver should not integrate past.  It is useful for\n\
 avoiding difficulties with singularities and points where there is a\n\
 discontinuity in the derivative.\n\
 \n\
+After a successful computation, the value of @var{istate} will be 2\n\
+(consistent with the Fortran version of @sc{Lsode}).\n\
+\n\
+If the computation is not successful, @var{istate} will be something\n\
+other than 2 and @var{msg} will contain additional information.\n\
+\n\
 You can use the function @code{lsode_options} to set optional\n\
 parameters for @code{lsode}.\n\
-@end deftypefn")
+@end deftypefn\n\
+@seealso{daspk, dassl, dasrt, odessa}")
 {
   octave_value_list retval;
 
--- a/src/DLD-FUNCTIONS/odessa.cc	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/DLD-FUNCTIONS/odessa.cc	Tue Oct 22 21:28:42 2002 +0000
@@ -253,19 +253,126 @@
     } \
   while (0)
 
-// --------------------------------
-// Everthing is so great above here
-// --------------------------------
-
 DEFUN_DLD (odessa, args, nargout,
-  "odessa (\"f\", x_0, theta, sx_0, t_out, t_crit)\n\
+  "-*- texinfo -*-\n\
+@deftypefn {Loadable Function} {[@var{x}, @var{sx}, @var{istate}, @var{msg}]} odessa (@var{fcn}, @var{x_0}, @var{p}, @var{sx_0}, @var{t}, @var{t_crit})\n\
+Solve the set of differential equations\n\
+@tex\n\
+$$ {dx \\over dt} = f (x, t; p) $$\n\
+with\n\
+$$ x(t_0) = x_0 $$\n\
+@end tex\n\
+@ifinfo\n\
+\n\
+@example\n\
+dx\n\
+-- = f(x, t; p)\n\
+dt\n\
+@end example\n\
+\n\
+with\n\
+\n\
+@example\n\
+x(t_0) = x_0\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+and simultaneously compute the first-order sensitivity coefficients\n\
+given by\n\
+\n\
+@example\n\
+s'(t) = j(t)*s(t) + df/dp\n\
+@end example\n\
+\n\
+in which\n\
+\n\
+@example\n\
+s(t)  = dx(t)/dp        (sensitivity functions)\n\
+s'(t) = d(dx(t)/dp)/dt\n\
+j(t)  = df(x,t;p)/dx(t) (Jacobian matrix)\n\
+df/dp = df(x,t;p)/dp    (inhomogeneity matrix)\n\
+@end example\n\
+\n\
+The solution is returned in the matrix @var{x}, with each row\n\
+corresponding to an element of the vector @var{t}.  The first element\n\
+of @var{t} should be @math{t_0} and should correspond to the initial\n\
+state of the system @var{x_0}, so that the first row of the output\n\
+is @var{x_0}.\n\
+\n\
+The sensitivities are returned in a list of matrices, @var{sx},\n\
+with each element of the list corresponding to an element of the\n\
+vector @var{t}.\n\
+\n\
+The first argument, @var{fcn}, is a string that names the function to\n\
+call to compute the vector of right hand sides for the set of equations.\n\
+The function must have the form\n\
+\n\
+@example\n\
+@var{xdot} = f (@var{x}, @var{t})\n\
+@end example\n\
 \n\
-The string \"f\" may be substituted for the vector of strings\n\
+@noindent\n\
+in which @var{xdot} and @var{x} are vectors and @var{t} is a scalar.\n\
+\n\
+The @var{fcn} argument may also be an array of strings\n\
+\n\
+@example\n\
+[\"f\"; \"j\"; \"b\"]\n\
+@end example\n\
+\n\
+in which the first element names the function @math{f} described\n\
+above, the second element names a function to compute the Jacobian\n\
+of @math{f}, and the third element names a function to compute the\n
+inhomogeneity matrix.\n\
+\n\
+The Jacobian function must have the form\n\
+\n\
+@example\n\
+@var{jac} = j (@var{x}, @var{t})\n\
+@end example\n\
+\n\
+in which @var{jac} is the matrix of partial derivatives\n\
+@tex\n\
+$$ J = {\\partial f_i \\over \\partial x_j} $$\n\
+@end tex\n\
+@ifinfo\n\
 \n\
-               [\"f\"; \"j\"; \"b\"] \n\
+@example\n\
+      df_i\n\
+jac = ----\n\
+      dx_j\n\
+@end example\n\
+\n\
+@end ifinfo\n\
+\n\
+The function @var{b} must have the form\n\
+\n\
+@example\n\
+@var{dfdp} = b (@var{t}, @var{x}, @var{p})\n\
+@end example\n\
+\n\
+The second argument, @var{x_0}, specifies the intial state of the system.\n\
+\n\
+The third argument, @var{p}, specifies the set of parameters.\n\
 \n\
-You can use the function @code{odessa_options} to set optional\n\
-parameters for @code{odessa}.")
+The fourth argument, @var{sx_0} specifies the initial values of the\n\
+sensitivities.\n\
+\n\
+The sixth argument is optional, and may be used to specify a set of\n\
+times that the ODE solver should not integrate past.  It is useful for\n\
+avoiding difficulties with singularities and points where there is a\n\
+discontinuity in the derivative.\n\
+\n\
+After a successful computation, the value of @var{istate} will be 2\n\
+(consistent with the Fortran version of @sc{Odessa}).\n\
+\n\
+If the computation is not successful, @var{istate} will be something\n\
+other than 2 and @var{msg} will contain additional information.\n\
+\n\
+You can use the function @code{lsode_options} to set optional\n\
+parameters for @code{lsode}.\n\
+@end deftypefn\n\
+@seealso{daspk, dassl, dasrt, lsode}")
 {
   octave_value_list retval;
 
--- a/src/help.cc	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/help.cc	Tue Oct 22 21:28:42 2002 +0000
@@ -613,7 +613,7 @@
 
       OSSTREAM buf;
 
-      buf << "sed -e 's/^[#%]* *//' -e 's/^ *@/@/' | "
+      buf << "sed -e 's/^[#%][#%]* *//' -e 's/^ *@/@/' | "
 	  << Vmakeinfo_prog
 	  << " -D \"VERSION " << OCTAVE_VERSION << "\""
 	  << " -D \"OCTAVEHOME " << OCTAVE_PREFIX << "\""
--- a/src/version.h	Thu Oct 17 21:39:48 2002 +0000
+++ b/src/version.h	Tue Oct 22 21:28:42 2002 +0000
@@ -23,7 +23,7 @@
 #if !defined (octave_version_h)
 #define octave_version_h 1
 
-#define OCTAVE_VERSION "2.1.36"
+#define OCTAVE_VERSION "2.1.37"
 
 #define OCTAVE_COPYRIGHT \
   "Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002 John W. Eaton."
@@ -38,6 +38,10 @@
   "There is ABSOLUTELY NO WARRANTY; not even for MERCHANTIBILITY or\n\
 FITNESS FOR A PARTICULAR PURPOSE."
 
+#define OCTAVE_CONTRIB_STATEMENT \
+  "Please contribute if you find this software useful.  For more\n\
+information, please visit http://www.octave.org/help-wanted.html"
+
 #define OCTAVE_BUGS_STATEMENT \
   "Report bugs to <bug-octave@bevo.che.wisc.edu>."
 
@@ -53,15 +57,14 @@
 #define X_OCTAVE_NAME_VERSION_COPYRIGHT_COPYING_WARRANTY_AND_BUGS(ARG) \
   OCTAVE_NAME_VERSION_COPYRIGHT_COPYING_AND_WARRANTY \
   ARG \
-  "\n\n" \
   OCTAVE_BUGS_STATEMENT
 
 #define OCTAVE_NAME_VERSION_COPYRIGHT_COPYING_WARRANTY_AND_BUGS \
-  X_OCTAVE_NAME_VERSION_COPYRIGHT_COPYING_WARRANTY_AND_BUGS ("")
+  X_OCTAVE_NAME_VERSION_COPYRIGHT_COPYING_WARRANTY_AND_BUGS ("\n\n")
 
 #define OCTAVE_STARTUP_MESSAGE \
   X_OCTAVE_NAME_VERSION_COPYRIGHT_COPYING_WARRANTY_AND_BUGS \
-    ("  For details, type `warranty'.")
+    ("  For details, type `warranty'.\n\n" OCTAVE_CONTRIB_STATEMENT "\n\n")
 
 #endif