Mercurial > hg > octave-lyh
annotate scripts/general/dblquad.m @ 17181:3a23cbde59d5
interpft.m: Fix interpolation to preserve spectral symmetry (bug #39566)
* interpft.m: Fix interpolation to preserve spectral symmetry, be compatible
with Matlab. Add test cases.
author | Mike Miller <mtmiller@ieee.org> |
---|---|
date | Sun, 04 Aug 2013 17:27:40 -0400 |
parents | 5d3a684236b0 |
children | bc924baa2c4e |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
12653
diff
changeset
|
1 ## Copyright (C) 2008-2012 David Bateman |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
2 ## |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
3 ## This file is part of Octave. |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
4 ## |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
5 ## 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
|
6 ## 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
|
7 ## 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
|
8 ## your option) any later version. |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
9 ## |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
10 ## 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
|
11 ## 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
|
12 ## 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
|
13 ## General Public License for more details. |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
14 ## |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
15 ## 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
|
16 ## 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
|
17 ## <http://www.gnu.org/licenses/>. |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
18 |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
19 ## -*- texinfo -*- |
12612
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
20 ## @deftypefn {Function File} {} dblquad (@var{f}, @var{xa}, @var{xb}, @var{ya}, @var{yb}) |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
21 ## @deftypefnx {Function File} {} dblquad (@var{f}, @var{xa}, @var{xb}, @var{ya}, @var{yb}, @var{tol}) |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
22 ## @deftypefnx {Function File} {} dblquad (@var{f}, @var{xa}, @var{xb}, @var{ya}, @var{yb}, @var{tol}, @var{quadf}) |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
23 ## @deftypefnx {Function File} {} dblquad (@var{f}, @var{xa}, @var{xb}, @var{ya}, @var{yb}, @var{tol}, @var{quadf}, @dots{}) |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
24 ## Numerically evaluate the double integral of @var{f}. |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
25 ## @var{f} is a function handle, inline function, or string |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
26 ## containing the name of the function to evaluate. The function @var{f} must |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
27 ## have the form @math{z = f(x,y)} where @var{x} is a vector and @var{y} is a |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
28 ## scalar. It should return a vector of the same length and orientation as |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
29 ## @var{x}. |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
30 ## |
12612
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
31 ## @var{xa}, @var{ya} and @var{xb}, @var{yb} are the lower and upper limits of |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
32 ## integration for x and y respectively. The underlying integrator determines |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
33 ## whether infinite bounds are accepted. |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
34 ## |
12642
f96b9b9f141b
doc: Periodic grammarcheck and spellcheck of documentation.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
35 ## The optional argument @var{tol} defines the absolute tolerance used to |
12612
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
36 ## integrate each sub-integral. The default value is @math{1e^{-6}}. |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
37 ## |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
38 ## The optional argument @var{quadf} specifies which underlying integrator |
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
39 ## function to use. Any choice but @code{quad} is available and the default |
12616
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
40 ## is @code{quadcc}. |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
41 ## |
9051
1bf0ce0930be
Grammar check TexInfo in all .m files
Rik <rdrider0-list@yahoo.com>
parents:
8920
diff
changeset
|
42 ## Additional arguments, are passed directly to @var{f}. To use the default |
12616
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
43 ## value for @var{tol} or @var{quadf} one may pass ':' or an empty matrix ([]). |
12575
d0b799dafede
Grammarcheck files for 3.4.1 release.
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
44 ## @seealso{triplequad, quad, quadv, quadl, quadgk, quadcc, trapz} |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
45 ## @end deftypefn |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
46 |
12616
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
47 function q = dblquad (f, xa, xb, ya, yb, tol = 1e-6, quadf = @quadcc, varargin) |
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
48 |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
49 if (nargin < 5) |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
50 print_usage (); |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
51 endif |
12616
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
52 if (isempty (tol)) |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11562
diff
changeset
|
53 tol = 1e-6; |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
54 endif |
12616
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
55 if (isempty (quadf)) |
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
56 quadf = @quadcc; |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
57 endif |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
58 |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
59 inner = @__dblquad_inner__; |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
60 if (ischar (f)) |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
61 f = @(x,y) feval (f, x, y, varargin{:}); |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
62 varargin = {}; |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
63 endif |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
64 |
8507 | 65 q = feval (quadf, @(y) inner (y, f, xa, xb, tol, quadf, |
10549 | 66 varargin{:}), ya, yb, tol); |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
67 endfunction |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
68 |
8507 | 69 function q = __dblquad_inner__ (y, f, xa, xb, tol, quadf, varargin) |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14363
diff
changeset
|
70 q = zeros (size (y)); |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
71 for i = 1 : length (y) |
8507 | 72 q(i) = feval (quadf, @(x) f(x, y(i), varargin{:}), xa, xb, tol); |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
73 endfor |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
74 endfunction |
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
75 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
76 |
12616
eb4afb6a1a51
dblquad.m,triplequad.m: Switch to quadcc as default integrator.
Rik <octave@nomad.inbox5.com>
parents:
12612
diff
changeset
|
77 %% Nasty integrand to show quadcc off |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
78 %!assert (dblquad (@(x,y) 1 ./ (x+y), 0, 1, 0, 1), 2*log (2), 1e-6) |
7771
680631e787aa
Add quadv, quadgk, dblquad and triplequad functions
David Bateman <dbateman@free.fr>
parents:
diff
changeset
|
79 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
80 %!assert (dblquad (@(x,y) exp (-x.^2 - y.^2) , -1, 1, -1, 1, 1e-6, @quadgk), pi * erf (1).^2, 1e-6) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
81 %!assert (dblquad (@(x,y) exp (-x.^2 - y.^2) , -1, 1, -1, 1, 1e-6, @quadl), pi * erf (1).^2, 1e-6) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
82 %!assert (dblquad (@(x,y) exp (-x.^2 - y.^2) , -1, 1, -1, 1, 1e-6, @quadv), pi * erf (1).^2, 1e-6) |
12612
16cca721117b
doc: Update all documentation for chapter on Numerical Integration
Rik <octave@nomad.inbox5.com>
parents:
12575
diff
changeset
|
83 |