Mercurial > forge
view main/gsl/src/DDDD_to_D.cc.template @ 11807:31d04ab60d20 octave-forge
gsl: remove duplicated content (thanks to Thomas Weber <tweber@debian.org>)
* DDDD_to_D.cc.template: file content was duplicated on r8745. Removing
duplicated text, and emacs text editor specific tags.
author | carandraug |
---|---|
date | Sun, 16 Jun 2013 15:43:53 +0000 |
parents | c762ea517481 |
children |
line wrap: on
line source
DEFUN_DLD(GSL_OCTAVE_NAME, args, nargout, " -*- texinfo -*-\n\ @deftypefn {Loadable Function} {@var{out} =} GSL_OCTAVE_NAME (@var{x0}, @var{x1}, @var{x2}, @var{x3})\n\ @deftypefnx {Loadable Function} {[@var{out}, @var{err}] =} GSL_OCTAVE_NAME (@dots{})\n\ \n\ GSL_FUNC_DOCSTRING \n\ @var{err} contains an estimate of the absolute error in the value @var{out}.a.\n\ \n\ This function is from the GNU Scientific Library,\n\ see @url{http://www.gnu.org/software/gsl/} for documentation.\n\ @end deftypefn\n\ ") // // // Generated R. Rogers 4/21/2008 // Version 1 Expanded to three argument input and added maintainence hints // Version 2 Expanded to four argument input (M. Helm 2011-10-08) // { int i; dim_vector dv; gsl_set_error_handler (octave_gsl_errorhandler); // Check number of arguments here if((args.length() != 4 )|| (nargout > 2)) { print_usage (); return octave_value(); } // Check argument types here if(!args(0).is_real_type() || !args(1).is_real_type() || !args(2).is_real_type() || !args(3).is_real_type()) { error("The arguments must be real."); print_usage (); return octave_value(); } // Nice combinatorial explosion here // Generate internal variables NDArray x0 = args(0).array_value(); NDArray x1 = args(1).array_value(); NDArray x2 = args(2).array_value(); NDArray x3 = args(3).array_value(); // // Case one; all inputs the same length A A A A if((x0.length() == x1.length() ) && (x0.length()==x2.length()) && (x0.length()==x3.length())) { dv = x0.dims(); NDArray out(dv); int len = x0.length(); // One output argument if(nargout < 2) { for(i = 0; i < len; i++) { out.xelem(i) = GSL_FUNC_NAME (x0.xelem(i), x1.xelem(i),x2.xelem(i), x3.xelem(i)); } return octave_value(out); // Two arguments } else { NDArray err(dv); gsl_sf_result result; octave_value_list retval; for(i = 0; i < len; i++) { GSL_FUNC_NAME_e (x0.xelem(i), x1.xelem(i), x2.xelem(i), x3.xelem(i), &result); out.xelem(i) = result.val; err.xelem(i) = result.err; } retval(1) = octave_value(err); retval(0) = octave_value(out); return retval; } // // Now we start on having only one array and three scalars, A S S S } else if(( x0.length() != 1) && (x1.length() == 1) && (x2.length()==1) && (x3.length()==1)) { dv = x0.dims(); NDArray out(dv); int len = x0.length(); // int x1_int = static_cast<int>(x1.xelem(0)); // int x2_int = static_cast<int>(x2.xelem(0)); double x1_real = x1.xelem(0); double x2_real = x2.xelem(0); double x3_real = x3.xelem(0); // One output argument if(nargout < 2) { for(i = 0; i < len; i++) { out.xelem(i) = GSL_FUNC_NAME (x0.xelem(i),x1_real,x2_real, x3_real); } return octave_value(out); // Two output argument } else { NDArray err(dv); gsl_sf_result result; octave_value_list retval; for(i = 0; i < len; i++) { GSL_FUNC_NAME_e (x0.xelem(i),x1_real,x2_real, x3_real, &result); out.xelem(i) = result.val; err.xelem(i) = result.err; } retval(1) = octave_value(err); retval(0) = octave_value(out); return retval; } // S A S S input form } else if((x0.length() == 1)&& ( x1.length() != 1) && (x2.length()==1) && (x3.length()==1)) { dv = x1.dims(); NDArray out(dv); int len = x1.length(); // int x0_int = static_cast<int>(x0.xelem(0)); // int x2_int = static_cast<int>(x2.xelem(0)); double x0_real = x0.xelem(0); double x2_real = x2.xelem(0); double x3_real = x3.xelem(0); // One output argument if(nargout < 2) { for(i = 0; i < len; i++) { out.xelem(i) = GSL_FUNC_NAME (x0_real,x1.xelem(i),x2_real,x3_real); } return octave_value(out); // Two output argument } else { NDArray err(dv); gsl_sf_result result; octave_value_list retval; for(i = 0; i < len; i++) { GSL_FUNC_NAME_e (x0_real,x1.xelem(i),x2_real, x3_real, &result); out.xelem(i) = result.val; err.xelem(i) = result.err; } retval(1) = octave_value(err); retval(0) = octave_value(out); return retval; } // S S A S input form } else if((x0.length() == 1)&& ( x1.length() == 1) && ( x2.length()!=1) && ( x3.length()==1)) { dv = x2.dims(); NDArray out(dv); int len = x2.length(); // int x0_int = static_cast<int>(x0.xelem(0)); // int x1_int = static_cast<int>(x1.xelem(0)); double x0_real = x0.xelem(0); double x1_real = x1.xelem(0); double x3_real = x3.xelem(0); // One output argument if(nargout < 2) { for(i = 0; i < len; i++) { out.xelem(i) = GSL_FUNC_NAME (x0_real,x1_real,x2.xelem(i),x3_real); } return octave_value(out); // Two output argument } else { NDArray err(dv); gsl_sf_result result; octave_value_list retval; for(i = 0; i < len; i++) { GSL_FUNC_NAME_e (x0_real,x1_real,x2.xelem(i),x3_real, &result); out.xelem(i) = result.val; err.xelem(i) = result.err; } retval(1) = octave_value(err); retval(0) = octave_value(out); return retval; } // S S S A input form } else if((x0.length() == 1)&& ( x1.length() == 1) && ( x2.length()==1) && ( x3.length()!=1)) { dv = x3.dims(); NDArray out(dv); int len = x3.length(); // int x0_int = static_cast<int>(x0.xelem(0)); // int x1_int = static_cast<int>(x1.xelem(0)); double x0_real = x0.xelem(0); double x1_real = x1.xelem(0); double x2_real = x2.xelem(0); // One output argument if(nargout < 2) { for(i = 0; i < len; i++) { out.xelem(i) = GSL_FUNC_NAME (x0_real,x1_real,x2_real,x3.xelem(i)); } return octave_value(out); // Two output argument } else { NDArray err(dv); gsl_sf_result result; octave_value_list retval; for(i = 0; i < len; i++) { GSL_FUNC_NAME_e (x0_real,x1_real,x2_real,x3.xelem(i), &result); out.xelem(i) = result.val; err.xelem(i) = result.err; } retval(1) = octave_value(err); retval(0) = octave_value(out); return retval; } } else { error("All arguments must either have the same size, or three of them must be scalar."); print_usage (); } return octave_value(); }