annotate scripts/general/quadv.m @ 20158:7503499a252b stable

doc: Update docstrings to have one sentence summary as first line. Update scripts in audio, elfun, general, geometry, and image directories. * scripts/audio/@audioplayer/__get_properties__.m, scripts/audio/@audioplayer/audioplayer.m, scripts/audio/@audioplayer/get.m, scripts/audio/@audioplayer/isplaying.m, scripts/audio/@audioplayer/play.m, scripts/audio/@audioplayer/playblocking.m, scripts/audio/@audioplayer/set.m, scripts/audio/@audioplayer/subsasgn.m, scripts/audio/@audioplayer/subsref.m, scripts/audio/@audiorecorder/audiorecorder.m, scripts/audio/@audiorecorder/get.m, scripts/audio/@audiorecorder/getaudiodata.m, scripts/audio/@audiorecorder/getplayer.m, scripts/audio/@audiorecorder/isrecording.m, scripts/audio/@audiorecorder/play.m, scripts/audio/@audiorecorder/record.m, scripts/audio/@audiorecorder/recordblocking.m, scripts/audio/@audiorecorder/set.m, scripts/audio/@audiorecorder/stop.m, scripts/audio/@audiorecorder/subsasgn.m, scripts/audio/@audiorecorder/subsref.m, scripts/audio/lin2mu.m, scripts/audio/mu2lin.m, scripts/audio/record.m, scripts/audio/sound.m, scripts/audio/soundsc.m, scripts/audio/wavread.m, scripts/audio/wavwrite.m, scripts/elfun/cosd.m, scripts/elfun/sind.m, scripts/elfun/tand.m, scripts/general/accumarray.m, scripts/general/accumdim.m, scripts/general/bitcmp.m, scripts/general/bitget.m, scripts/general/bitset.m, scripts/general/blkdiag.m, scripts/general/cart2pol.m, scripts/general/cart2sph.m, scripts/general/cell2mat.m, scripts/general/celldisp.m, scripts/general/chop.m, scripts/general/circshift.m, scripts/general/common_size.m, scripts/general/cplxpair.m, scripts/general/cumtrapz.m, scripts/general/dblquad.m, scripts/general/deal.m, scripts/general/del2.m, scripts/general/display.m, scripts/general/divergence.m, scripts/general/fieldnames.m, scripts/general/flip.m, scripts/general/flipdim.m, scripts/general/fliplr.m, scripts/general/flipud.m, scripts/general/gradient.m, scripts/general/interp3.m, scripts/general/interpft.m, scripts/general/interpn.m, scripts/general/loadobj.m, scripts/general/logspace.m, scripts/general/methods.m, scripts/general/nargchk.m, scripts/general/narginchk.m, scripts/general/nargoutchk.m, scripts/general/nextpow2.m, scripts/general/nthargout.m, scripts/general/num2str.m, scripts/general/pol2cart.m, scripts/general/polyarea.m, scripts/general/postpad.m, scripts/general/prepad.m, scripts/general/profile.m, scripts/general/quadgk.m, scripts/general/quadl.m, scripts/general/quadv.m, scripts/general/randi.m, scripts/general/rat.m, scripts/general/repmat.m, scripts/general/rot90.m, scripts/general/rotdim.m, scripts/general/saveobj.m, scripts/general/shift.m, scripts/general/shiftdim.m, scripts/general/sortrows.m, scripts/general/sph2cart.m, scripts/general/structfun.m, scripts/general/subsindex.m, scripts/general/trapz.m, scripts/general/triplequad.m, scripts/geometry/delaunayn.m, scripts/geometry/dsearch.m, scripts/geometry/dsearchn.m, scripts/geometry/griddata.m, scripts/geometry/griddata3.m, scripts/geometry/griddatan.m, scripts/geometry/inpolygon.m, scripts/geometry/rectint.m, scripts/geometry/tsearchn.m, scripts/geometry/voronoi.m, scripts/geometry/voronoin.m, scripts/help/__unimplemented__.m, scripts/help/doc.m, scripts/help/doc_cache_create.m, scripts/help/get_first_help_sentence.m, scripts/help/help.m, scripts/help/lookfor.m, scripts/help/print_usage.m, scripts/help/type.m, scripts/help/which.m, scripts/image/autumn.m, scripts/image/bone.m, scripts/image/brighten.m, scripts/image/cmpermute.m, scripts/image/colorcube.m, scripts/image/contrast.m, scripts/image/cool.m, scripts/image/copper.m, scripts/image/cubehelix.m, scripts/image/flag.m, scripts/image/gmap40.m, scripts/image/gray.m, scripts/image/gray2ind.m, scripts/image/hot.m, scripts/image/hsv.m, scripts/image/image.m, scripts/image/imagesc.m, scripts/image/imfinfo.m, scripts/image/imformats.m, scripts/image/imread.m, scripts/image/imshow.m, scripts/image/imwrite.m, scripts/image/iscolormap.m, scripts/image/jet.m, scripts/image/lines.m, scripts/image/ntsc2rgb.m, scripts/image/ocean.m, scripts/image/pink.m, scripts/image/prism.m, scripts/image/rainbow.m, scripts/image/rgb2ntsc.m, scripts/image/spinmap.m, scripts/image/spring.m, scripts/image/summer.m, scripts/image/white.m, scripts/image/winter.m: Update docstrings to have one sentence summary as first line. Re-structure to have line lengths <= 80 chars.
author Rik <rik@octave.org>
date Sun, 03 May 2015 09:36:20 -0700
parents 9fc020886ae9
children
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
19697
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 17744
diff changeset
1 ## Copyright (C) 2008-2015 David Bateman
14138
72c96de7a403 maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents: 13851
diff changeset
2 ## Copyright (C) 2012 Alexander Klein
13850
64ab148fcbb9 quadv.m: Fixes for convergence issues with vector-valued functions (patch #7627)
Alexander Klein <alexander.klein@math.uni-giessen.de>
parents: 12642
diff changeset
3 ##
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
4 ## This file is part of Octave.
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
5 ##
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
6 ## Octave is free software; you can redistribute it and/or modify it
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
7 ## under the terms of the GNU General Public License as published by
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
8 ## the Free Software Foundation; either version 3 of the License, or (at
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
9 ## your option) any later version.
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
10 ##
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
11 ## Octave is distributed in the hope that it will be useful, but
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
12 ## WITHOUT ANY WARRANTY; without even the implied warranty of
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
13 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
14 ## General Public License for more details.
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
15 ##
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
16 ## You should have received a copy of the GNU General Public License
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
17 ## along with Octave; see the file COPYING. If not, see
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
18 ## <http://www.gnu.org/licenses/>.
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
19
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
20 ## -*- texinfo -*-
10793
be55736a0783 Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents: 10549
diff changeset
21 ## @deftypefn {Function File} {@var{q} =} quadv (@var{f}, @var{a}, @var{b})
10430
f1567b3e1108 scripts/general/quadv.m: Replace 'quadl' with 'quadv' in help text
Soren Hauberg <hauberg@gmail.com>
parents: 9209
diff changeset
22 ## @deftypefnx {Function File} {@var{q} =} quadv (@var{f}, @var{a}, @var{b}, @var{tol})
f1567b3e1108 scripts/general/quadv.m: Replace 'quadl' with 'quadv' in help text
Soren Hauberg <hauberg@gmail.com>
parents: 9209
diff changeset
23 ## @deftypefnx {Function File} {@var{q} =} quadv (@var{f}, @var{a}, @var{b}, @var{tol}, @var{trace})
f1567b3e1108 scripts/general/quadv.m: Replace 'quadl' with 'quadv' in help text
Soren Hauberg <hauberg@gmail.com>
parents: 9209
diff changeset
24 ## @deftypefnx {Function File} {@var{q} =} quadv (@var{f}, @var{a}, @var{b}, @var{tol}, @var{trace}, @var{p1}, @var{p2}, @dots{})
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
25 ## @deftypefnx {Function File} {[@var{q}, @var{nfun}] =} quadv (@dots{})
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
26 ##
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11563
diff changeset
27 ## Numerically evaluate the integral of @var{f} from @var{a} to @var{b}
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
28 ## using an adaptive Simpson's rule.
20158
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
29 ##
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
30 ## @var{f} is a function handle, inline function, or string containing the name
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
31 ## of the function to evaluate. @code{quadv} is a vectorized version of
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
32 ## @code{quad} and the function defined by @var{f} must accept a scalar or
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
33 ## vector as input and return a scalar, vector, or array as output.
11078
2aec7e3b8553 Fix help string in general/quadv.m
Carlo de Falco <kingcrimson@tiscali.it>
parents: 10793
diff changeset
34 ##
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
35 ## @var{a} and @var{b} are the lower and upper limits of integration. Both
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
36 ## limits must be finite.
12642
f96b9b9f141b doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents: 12612
diff changeset
37 ##
20158
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
38 ## The optional argument @var{tol} defines the tolerance used to stop the
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
39 ## adaptation procedure. The default value is 1e-6.
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
40 ##
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
41 ## The algorithm used by @code{quadv} involves recursively subdividing the
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
42 ## integration interval and applying Simpson's rule on each subinterval.
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
43 ## If @var{trace} is true then after computing each of these partial
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
44 ## integrals display: (1) the total number of function evaluations,
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
45 ## (2) the left end of the subinterval, (3) the length of the subinterval,
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
46 ## (4) the approximation of the integral over the subinterval.
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
47 ##
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
48 ## Additional arguments @var{p1}, etc., are passed directly to the function
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
49 ## @var{f}. To use default values for @var{tol} and @var{trace}, one may pass
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
50 ## empty matrices ([]).
11078
2aec7e3b8553 Fix help string in general/quadv.m
Carlo de Falco <kingcrimson@tiscali.it>
parents: 10793
diff changeset
51 ##
20158
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
52 ## The result of the integration is returned in @var{q}
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
53 ##
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
54 ## @var{nfun} indicates the number of function evaluations that were made.
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
55 ##
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
56 ## Note: @code{quadv} is written in Octave's scripting language and can be
12642
f96b9b9f141b doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents: 12612
diff changeset
57 ## used recursively in @code{dblquad} and @code{triplequad}, unlike the
20158
7503499a252b doc: Update docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19833
diff changeset
58 ## @code{quad} function.
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
59 ## @seealso{quad, quadl, quadgk, quadcc, trapz, dblquad, triplequad}
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
60 ## @end deftypefn
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
61
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
62 function [q, nfun] = quadv (f, a, b, tol, trace, varargin)
13851
372869056e0f quadv.m: style fixes
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 13850
diff changeset
63 ## TODO: Make norm for convergence testing configurable
372869056e0f quadv.m: style fixes
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 13850
diff changeset
64
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
65 if (nargin < 3)
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
66 print_usage ();
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
67 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
68 if (nargin < 4)
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11563
diff changeset
69 tol = [];
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
70 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
71 if (nargin < 5)
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11563
diff changeset
72 trace = [];
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
73 endif
7795
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
74 if (isa (a, "single") || isa (b, "single"))
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
75 myeps = eps ("single");
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
76 else
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
77 myeps = eps;
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
78 endif
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
79 if (isempty (tol))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11563
diff changeset
80 tol = 1e-6;
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
81 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
82 if (isempty (trace))
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11563
diff changeset
83 trace = 0;
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
84 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
85
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
86 ## Split the interval into 3 abscissa, and apply a 3 point Simpson's rule
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
87 c = (a + b) / 2;
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
88 fa = feval (f, a, varargin{:});
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
89 fc = feval (f, c, varargin{:});
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
90 fb = feval (f, b, varargin{:});
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
91 nfun = 3;
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
92
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
93 ## If have edge singularities, move edge point by eps*(b-a) as
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
94 ## discussed in Shampine paper used to implement quadgk
13851
372869056e0f quadv.m: style fixes
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 13850
diff changeset
95 if (any (isinf (fa(:))))
7795
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
96 fa = feval (f, a + myeps * (b-a), varargin{:});
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
97 endif
13851
372869056e0f quadv.m: style fixes
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents: 13850
diff changeset
98 if (any (isinf (fb(:))))
7795
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
99 fb = feval (f, b - myeps * (b-a), varargin{:});
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
100 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
101
11078
2aec7e3b8553 Fix help string in general/quadv.m
Carlo de Falco <kingcrimson@tiscali.it>
parents: 10793
diff changeset
102 h = (b - a);
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
103 q = (b - a) / 6 * (fa + 4 * fc + fb);
11587
c792872f8942 all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents: 11563
diff changeset
104
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
105 [q, nfun, hmin] = simpsonstp (f, a, b, c, fa, fb, fc, q, nfun, abs (h),
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10430
diff changeset
106 tol, trace, varargin{:});
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
107
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
108 if (nfun > 10000)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
109 warning ("maximum iteration count reached");
16969
58188d5a2587 Use isfinite() to replace isinf() | isnan() combination for 30% speed-up.
Rik <rik@octave.org>
parents: 14868
diff changeset
110 elseif (any (! isfinite (q(:))))
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
111 warning ("infinite or NaN function evaluations were returned");
7795
df9519e9990c Handle single precision eps values
David Bateman <dbateman@free.fr>
parents: 7771
diff changeset
112 elseif (hmin < (b - a) * myeps)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
113 warning ("minimum step size reached -- possibly singular integral");
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
114 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
115 endfunction
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
116
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
117 function [q, nfun, hmin] = simpsonstp (f, a, b, c, fa, fb, fc, q0,
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
118 nfun, hmin, tol, trace, varargin)
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
119 if (nfun > 10000)
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
120 q = q0;
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
121 else
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
122 d = (a + c) / 2;
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
123 e = (c + b) / 2;
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
124 fd = feval (f, d, varargin{:});
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
125 fe = feval (f, e, varargin{:});
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
126 nfun += 2;
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
127 q1 = (c - a) / 6 * (fa + 4 * fd + fc);
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
128 q2 = (b - c) / 6 * (fc + 4 * fe + fb);
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
129 q = q1 + q2;
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
130
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
131 if (abs(a - c) < hmin)
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
132 hmin = abs (a - c);
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
133 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
134
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
135 if (trace)
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
136 disp ([nfun, a, b-a, q]);
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
137 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
138
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
139 ## Force at least one adpative step.
13850
64ab148fcbb9 quadv.m: Fixes for convergence issues with vector-valued functions (patch #7627)
Alexander Klein <alexander.klein@math.uni-giessen.de>
parents: 12642
diff changeset
140 ## Not vectorizing q-q0 in the norm provides a more rigid criterion for
64ab148fcbb9 quadv.m: Fixes for convergence issues with vector-valued functions (patch #7627)
Alexander Klein <alexander.klein@math.uni-giessen.de>
parents: 12642
diff changeset
141 ## matrix-valued functions.
64ab148fcbb9 quadv.m: Fixes for convergence issues with vector-valued functions (patch #7627)
Alexander Klein <alexander.klein@math.uni-giessen.de>
parents: 12642
diff changeset
142 if (nfun == 5 || norm (q - q0, Inf) > tol)
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
143 [q1, nfun, hmin] = simpsonstp (f, a, c, d, fa, fc, fd, q1, nfun, hmin,
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10430
diff changeset
144 tol, trace, varargin{:});
12612
16cca721117b doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents: 11587
diff changeset
145 [q2, nfun, hmin] = simpsonstp (f, c, b, e, fc, fb, fe, q2, nfun, hmin,
10549
95c3e38098bf Untabify .m scripts
Rik <code@nomad.inbox5.com>
parents: 10430
diff changeset
146 tol, trace, varargin{:});
8507
cadc73247d65 style fixes
John W. Eaton <jwe@octave.org>
parents: 7795
diff changeset
147 q = q1 + q2;
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
148 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
149 endif
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
150 endfunction
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
151
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
152
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
153 %!assert (quadv (@sin, 0, 2 * pi), 0, 1e-5)
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
154 %!assert (quadv (@sin, 0, pi), 2, 1e-5)
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
155
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
156 ## Handles weak singularities at the edge
14868
5d3a684236b0 maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents: 14363
diff changeset
157 %!assert (quadv (@(x) 1 ./ sqrt (x), 0, 1), 2, 1e-5)
7771
680631e787aa Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff changeset
158
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
159 ## Handles vector-valued functions
11078
2aec7e3b8553 Fix help string in general/quadv.m
Carlo de Falco <kingcrimson@tiscali.it>
parents: 10793
diff changeset
160 %!assert (quadv (@(x) [(sin (x)), (sin (2 * x))], 0, pi), [2, 0], 1e-5)
2aec7e3b8553 Fix help string in general/quadv.m
Carlo de Falco <kingcrimson@tiscali.it>
parents: 10793
diff changeset
161
19833
9fc020886ae9 maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 19697
diff changeset
162 ## Handles matrix-valued functions
13850
64ab148fcbb9 quadv.m: Fixes for convergence issues with vector-valued functions (patch #7627)
Alexander Klein <alexander.klein@math.uni-giessen.de>
parents: 12642
diff changeset
163 %!assert (quadv (@(x) [ x, x, x; x, 1./sqrt(x), x; x, x, x ], 0, 1 ), [0.5, 0.5, 0.5; 0.5, 2, 0.5; 0.5, 0.5, 0.5], 1e-5)
14363
f3d52523cde1 Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
164