AAA: Acoustics of a Fish-Tank

 July 29, 2018    Antonin Novak
research    research_acoustics    research_bioacoustics    

A. Novak, M. Bruneau & P. Lotton (2018), “Small-Sized Rectangular Liquid-Filled Acoustical Tank Excitation: A Modal Approach Including Leakage Through the Walls”, Acta Acustica united with Acustica. Vol. 104(4), pp. 586-596.


The acoustics of a water-filled rectangular fish-tank is very different from other 3D spaces known in classical acoustics and room acoustics. Since the impedance of the water (similar to the impedance of the glass walls surrounding the tank) is much higher than the impedance of the air outside, the coefficient of reflection is much closer to $R \simeq -1$ (close to Dirichlet boundary condition). Consequently, the low frequencies are completely absorbed (reflected from the walls with opposite phase). In this paper, a precise but straightforward analytical model is derived from scratch. The goal of the model is to predict the Green's function between any two points in the tank within a small amount of time (< 1s) so that the model could be used effectively in applications calculating the Green's function several times, such as localization algorithms. The following figure shows the measured Frequency Response Function (blue curve) between source (coordinates [27, 8, 8] cm) and receiver (coordinates [23.5, 15 14] cm) in a tank with dimensions 58.6 cm x 28.9 cm x 18.4 cm. The green dashed-dotted line shows the basic Dirichlet model ($R=-1$), the red dashed line shows the Green's function of the model. Computation of the Green's function (all the depicted frequencies) took 0.4 seconds.


In the following figure, several Frequency Response Functions (a) are depicted for different positions of the receiver (b). It is clear from the results (a) and (c), that while the low-frequency content is very quickly attenuated, the value of the pressure field at the resonant frequency remains high and almost independent of the receiver position.



In the following audio examples we show how the sound is propagated in a water tank.

The original audio file:


Position 10 cm from the source:


Position 20 cm from the source:


Position 30 cm from the source:


Position 40 cm from the source:



The original audio file is a sample of the song "The Letter" originally written by Joe Cocker and recorded by the French group Les Commis du Mans.


The following Matlab code calculates the Green Function (Frequency Response Function) inside a fish-tank with given dimensions L = [Lx Ly Lz] between the source position r0 = [x0 y0 z0] and receiver position r = [x y z]. The function small_tank_FRF_function() is shown below.


clc; close all; clearvars;

%% define parameters
L = [0.4906  0.2386  0.141]; % tank dimensions [Lx Ly Lz] [m] 
h = 0.0027;                  % wall thickness
r0 = [0.33  0.13  0.11];     % source position [x0 y0 z0] [m] 
r =  [0.20  0.15  0.11];     % hydrophone position [x y z] [m] 
number_of_modes = 15;        % number of modes

%% define frequency range
f1 = 1e2; f2 = 10e3;
N_of_frequencies = 5000;
f = logspace(log10(f1),log10(f2),N_of_frequencies);

%% calculate the FRF
FRF = small_tank_FRF_function(L,r0,r,h,number_of_modes,2*pi*f);

%% show results
figure; plot(f/1000,20*log10(abs(FRF)),'linewidth',2)
xlabel('Frequency [kHz]'); ylabel('Magnitude [dB]');
xlim([0 10]); ylim([-70 30]);



small_tank_FRF_function.m


function FRF = small_tank_FRF_function(L,r0,r,h,number_of_modes,omega)
%% model of a small fish-tank
%
% Calculates the Green Function (Frequency Response Function) inside 
% a fish-tank with given dimensions L = [Lx Ly Lz] between the source 
% position r0 = [x0 y0 z0] and receiver position r = [x y z]. 
%
% FRF = small_tank_FRF_function(L,r0,r,h,number_of_modes,omega)
%
% inputs:    L               ... fish-tank dimensions [Lx, Ly, Lz]
%            r0              ... source position [x0 y0 z0]
%            r               ... receiver position [x y z]
%            h               ... wall thickness
%            number_of_modes ... number of modes
%            omega           ... angular frequency (2*pi*f)
%
% output:    FRF             ... Frequency Repsonse Function
%
% [Novak et al. (2018)] 
% A. Novak, M. Bruneau & P. Lotton (2018), "Small-Sized Rectangular 
%      Liquid-Filled Acoustical Tank Excitation: A Modal Approach Including 
%      Leakage Through the Walls", Acta Acustica united with Acustica. 
%      Vol. 104(4), pp. 586-596.
%
% Antonin Novak (antonin.novak@univ-lemans.fr), 25/10/2018
% Laboratoire d'Acoustique de l'Université du Mans (LAUM, UMR CNRS 6613), 
% 72085 Le Mans, France
%
% https://ant-novak.com
%
%%

% wall properties
c_w = 5600; rho_w = 2500; Zw = c_w * rho_w;
% liquid properties
c_l = 1488; rho_l = 998; Zl = c_l * rho_l;
% air properties
c_a = 340;  rho_a = 1.255; Za = c_a * rho_a;

% wave number
kl = omega./c_l; % liquid

%% 

psum = 0;

for mx=1:number_of_modes
    display(['mode n.' num2str(mx)]);
    
    for my=1:number_of_modes
    
        for mz=1:number_of_modes
            
            % Dirichlet (0th order) modal wave-number components (in liquid)
            kmx = mx*pi/L(1);
            kmy = my*pi/L(2);
            kmz = mz*pi/L(3);               
            % (0th order) modulus of the  wave-number (in liquid)
            km = sqrt(kmx.^2 + kmy.^2 + kmz.^2);
            
            % (0th order) modulus of modal wave-numbers in air (kam) and wall (kwm)
            kwm = c_l/c_w*km.*(1-1i*0.022); % 0.022 losses in walls
            kam = c_l/c_a*km;

            % modal wave-number components in wall
            % Eq.(A9) of [Novak et al. (2018)] (for wall)                      
            kwmx = -1i*sqrt(-kwm.^2 + ((kmy).^2 + (kmz).^2));
            kwmy = -1i*sqrt(-kwm.^2 + ((kmx).^2 + (kmz).^2));
            kwmz = -1i*sqrt(-kwm.^2 + ((kmx).^2 + (kmy).^2));  

            % modal wave-number components in air
            % modified Eq.(A9) of [Novak et al. (2018)] for air
            kamx = -1i*sqrt(-kam.^2 + ((kmy).^2 + (kmz).^2));
            kamy = -1i*sqrt(-kam.^2 + ((kmx).^2 + (kmz).^2));
            kamz = -1i*sqrt(-kam.^2 + ((kmx).^2 + (kmy).^2));

            % Eq.(9) of [Novak et al. (2018)]           
            Ealx = Za./Zl .* (kmx./km)./(kamx./kam);
            Ealy = Za./Zl .* (kmy./km)./(kamy./kam);
            Ealz = Za./Zl .* (kmz./km)./(kamz./kam);
            
            % Eq.(9) of [Novak et al. (2018)]           
            Ewlx = Zw./Zl .* (kmx./km)./(kwmx./kwm);
            Ewly = Zw./Zl .* (kmy./km)./(kwmy./kwm);
            Ewlz = Zw./Zl .* (kmz./km)./(kwmz./kwm);
            
            % Eq.(9) of [Novak et al. (2018)]           
            Eawx = Za./Zw .* (kwmx./kwm)./(kamx./kam);
            Eawy = Za./Zw .* (kwmy./kwm)./(kamy./kam);
            Eawz = Za./Zw .* (kwmz./kwm)./(kamz./kam);
            
            % modal specific impedance components  
            % Eq.(8) of [Novak et al. (2018)]
            zetax = (Ealx + 1i*Ewlx.*tan(kwmx*h))./(1 + Eawx.*tan(kwmx*h));
            zetay = (Ealy + 1i*Ewly.*tan(kwmy*h))./(1 + Eawy.*tan(kwmy*h));
            zetaz = (Ealz + 1i*Ewlz.*tan(kwmz*h))./(1 + Eawz.*tan(kwmz*h));
            zetaz2 = Ealz;            
          
            % Eq.(16) of [Novak et al. (2018)]
            kmx = (mx*pi +  1i .* (zetax + zetax))./L(1);
            kmy = (my*pi +  1i .* (zetay + zetay))./L(2);
            kmz = (mz*pi +  1i .* (zetaz + zetaz2))./L(3);            
            km2 = kmx.^2 + kmy.^2 + kmz.^2;
               
            % Eq.(18) of [Novak et al. (2018)]
            Psi0 = sin(kmx*r0(1) - 1i.*zetax) .* sin(kmy*r0(2) - 1i.*zetay) .* sin(kmz*r0(3) - 1i.*zetaz);            
            Psir  = sin(kmx*r(1)  - 1i.*zetax) .* sin(kmy*r(2)  - 1i.*zetay) .* sin(kmz*r(3)  - 1i.*zetaz);            
                       
            % Eq.(25) of [Novak et al. (2018)]
            psum = psum + Psir.*Psi0 ./ (km2-kl.^2);
            
        end
    end
end

FRF = sqrt(2^3./prod(L)) * psum;


Small-sized water-filled tanks with thin (compared to the acoustic wavelength) walls that are surrounded on all sides by air, have been used over the past decades for underwater studies. So far, the analytical modelling of the acoustic field inside such an enclosure considers either the eigenvalue problem with pressure release conditions at the walls or imposes empirical impedance boundary conditions, or even proceeds in analogy with the field within a waveguide. At low frequencies, approximations involving an incompressible fluid or the Laplace equation have been used. Those models have limitations that are always caused by simplification of the boundary conditions. This paper deals with both an analytical formulation expressing the acoustic leakage through the walls (lossy and reacting walls) and modal solutions for the acoustic pressure and acoustic velocity fields when a source emits energy and when it is shut off, providing the transient acoustic response of small tanks quantitatively. Several analytical results are compared with experimental observations.

@article{novak2018small,
author={Novak, Antonin and Bruneau, Michel and Lotton, Pierrick},
    title={Small-Sized Rectangular Liquid-Filled Acoustical Tank Excitation: A Modal Approach Including Leakage Through the Walls},
    journal={Acta Acustica united with Acustica},
    volume={104},
    number={4},
    pages={586-596},
    year={2018},
    publisher={S. Hirzel Verlag}
}