Previous: , Up: Numerical Integration   [Contents][Index]


23.3 Functions of Multiple Variables

Octave does not have built-in functions for computing the integral of functions of multiple variables directly. It is possible, however, to compute the integral of a function of multiple variables using the existing functions for one-dimensional integrals.

To illustrate how the integration can be performed, we will integrate the function

f(x, y) = sin(pi*x*y)*sqrt(x*y)

for x and y between 0 and 1.

The first approach creates a function that integrates f with respect to x, and then integrates that function with respect to y. Because quad is written in Fortran it cannot be called recursively. This means that quad cannot integrate a function that calls quad, and hence cannot be used to perform the double integration. Any of the other integrators, however, can be used which is what the following code demonstrates.

function q = g(y)
  q = ones (size (y));
  for i = 1:length (y)
    f = @(x) sin (pi*x.*y(i)) .* sqrt (x.*y(i));
    q(i) = quadgk (f, 0, 1);
  endfor
endfunction

I = quadgk ("g", 0, 1)
      ⇒ 0.30022

The above process can be simplified with the dblquad and triplequad functions for integrals over two and three variables. For example:

I = dblquad (@(x, y) sin (pi*x.*y) .* sqrt (x.*y), 0, 1, 0, 1)
      ⇒ 0.30022
: dblquad (f, xa, xb, ya, yb)
: dblquad (f, xa, xb, ya, yb, tol)
: dblquad (f, xa, xb, ya, yb, tol, quadf)
: dblquad (f, xa, xb, ya, yb, tol, quadf, …)

Numerically evaluate the double integral of f.

f is a function handle, inline function, or string containing the name of the function to evaluate. The function f must have the form z = f(x,y) where x is a vector and y is a scalar. It should return a vector of the same length and orientation as x.

xa, ya and xb, yb are the lower and upper limits of integration for x and y respectively. The underlying integrator determines whether infinite bounds are accepted.

The optional argument tol defines the absolute tolerance used to integrate each sub-integral. The default value is 1e^{-6}.

The optional argument quadf specifies which underlying integrator function to use. Any choice but quad is available and the default is quadcc.

Additional arguments, are passed directly to f. To use the default value for tol or quadf one may pass ':' or an empty matrix ([]).

See also: triplequad, quad, quadv, quadl, quadgk, quadcc, trapz.

: triplequad (f, xa, xb, ya, yb, za, zb)
: triplequad (f, xa, xb, ya, yb, za, zb, tol)
: triplequad (f, xa, xb, ya, yb, za, zb, tol, quadf)
: triplequad (f, xa, xb, ya, yb, za, zb, tol, quadf, …)

Numerically evaluate the triple integral of f.

f is a function handle, inline function, or string containing the name of the function to evaluate. The function f must have the form w = f(x,y,z) where either x or y is a vector and the remaining inputs are scalars. It should return a vector of the same length and orientation as x or y.

xa, ya, za and xb, yb, zb are the lower and upper limits of integration for x, y, and z respectively. The underlying integrator determines whether infinite bounds are accepted.

The optional argument tol defines the absolute tolerance used to integrate each sub-integral. The default value is 1e-6.

The optional argument quadf specifies which underlying integrator function to use. Any choice but quad is available and the default is quadcc.

Additional arguments, are passed directly to f. To use the default value for tol or quadf one may pass ':' or an empty matrix ([]).

See also: dblquad, quad, quadv, quadl, quadgk, quadcc, trapz.

The above mentioned approach works, but is fairly slow, and that problem increases exponentially with the dimensionality of the integral. Another possible solution is to use Orthogonal Collocation as described in the previous section (see Orthogonal Collocation). The integral of a function f(x,y) for x and y between 0 and 1 can be approximated using n points by the sum over i=1:n and j=1:n of q(i)*q(j)*f(r(i),r(j)), where q and r is as returned by colloc (n). The generalization to more than two variables is straight forward. The following code computes the studied integral using n=8 points.

f = @(x,y) sin (pi*x*y') .* sqrt (x*y');
n = 8;
[t, ~, ~, q] = colloc (n);
I = q'*f(t,t)*q;
      ⇒ 0.30022

It should be noted that the number of points determines the quality of the approximation. If the integration needs to be performed between a and b, instead of 0 and 1, then a change of variables is needed.


Previous: , Up: Numerical Integration   [Contents][Index]