Skip to content

Commit

Permalink
Create new class to contain DCT methods
Browse files Browse the repository at this point in the history
Issue #500 Create a new class  to contain DCT methods. Shift across methods from MFCCStuff.cs. Make this a static class but it could/should eventually be used to contain the parameters required to perform the DCT.

A bunch of other methods where affected by this shift of DCT methods.
  • Loading branch information
towsey committed Jun 30, 2021
1 parent 6b90e4f commit 4fc5bdc
Show file tree
Hide file tree
Showing 9 changed files with 87 additions and 49 deletions.
2 changes: 1 addition & 1 deletion src/AnalysisPrograms/Recognizers/Frogs/LitoriaBicolor.cs
Original file line number Diff line number Diff line change
Expand Up @@ -193,7 +193,7 @@ public static Tuple<BaseSonogram, double[,], double[], List<AcousticEvent>, Imag
int dctLength = (int)Math.Round(framesPerSecond * dctDuration);

// set up the cosine coefficients
double[,] cosines = MFCCStuff.Cosines(dctLength, dctLength);
double[,] cosines = DctMethods.Cosines(dctLength, dctLength);

int upperBandMinBin = (int)Math.Round(lbConfig.UpperBandMinHz / freqBinWidth) + 1;
int upperBandMaxBin = (int)Math.Round(lbConfig.UpperBandMaxHz / freqBinWidth) + 1;
Expand Down
2 changes: 1 addition & 1 deletion src/AnalysisPrograms/Recognizers/Frogs/LitoriaRothii.cs
Original file line number Diff line number Diff line change
Expand Up @@ -136,7 +136,7 @@ public override RecognizerResults Recognize(AudioRecording recording, Config con
int dctLength = (int)Math.Round(framesPerSecond * dctDuration);

// set up the cosine coefficients
double[,] cosines = MFCCStuff.Cosines(dctLength, dctLength);
double[,] cosines = DctMethods.Cosines(dctLength, dctLength);

BaseSonogram sonogram = new SpectrogramStandard(sonoConfig, recording.WavReader);
int rowCount = sonogram.Data.GetLength(0);
Expand Down
74 changes: 74 additions & 0 deletions src/AudioAnalysisTools/DSP/DctMethods.cs
Original file line number Diff line number Diff line change
@@ -0,0 +1,74 @@
// <copyright file="DctMethods.cs" company="QutEcoacoustics">
// All code in this file and all associated files are the copyright and property of the QUT Ecoacoustics Research Group (formerly MQUTeR, and formerly QUT Bioacoustics Research Group).
// </copyright>

namespace AudioAnalysisTools.DSP
{
using System;
using System.Collections.Generic;
using System.Linq;
using System.Text;
using System.Threading.Tasks;
using TowseyLibrary;

public static class DctMethods
{
/// <summary>
/// Returns a matrix of cosine basis functions.
/// These are prepared prior to performing a DCT, Discrete Cosine Transform.
/// The rows k = 0 to coeffCount are the basis functions.
/// The columns, m = 0 to M where M = signalLength or the length of the required DCT.
/// The value of m/M ranges from 0 to 1.0.
/// The value of Pi*m/M ranges from 0 to Pi radians.
/// The value of k*Pi*m/M ranges from 0 to k*Pi radians. WHen k=2, 2Pi radians corresponds to one rotation.
/// </summary>
/// <param name="signalLength">The length of the signal to be processed. e.g. the frequency bin count or filter bank count or ...</param>
/// <param name="coeffCount">The number of basis funcitons = the rquired number of DCT coefficients.</param>
public static double[,] Cosines(int signalLength, int coeffCount)
{
double[,] cosines = new double[coeffCount + 1, signalLength]; //get an extra coefficient because do not want DC coeff at [0].
for (int k = 0; k < coeffCount + 1; k++)
{
double kPiOnM = k * Math.PI / signalLength;

// for each spectral bin
for (int m = 0; m < signalLength; m++)
{
cosines[k, m] = Math.Cos(kPiOnM * (m + 0.5)); //can also be Cos(kPiOnM * (m - 0.5)
}
}

return cosines;
}

//following two lines write matrix of cos values for checking.
//string fPath = @"C:\SensorNetworks\Sonograms\cosines.txt";
//FileTools.WriteMatrix2File_Formatted(cosines, fPath, "F3");

//following two lines write bmp image of cos values for checking.
//string fPath = @"C:\SensorNetworks\Output\cosines.bmp";
//ImageTools.DrawMatrix(cosines, fPath);

public static double[] DoDct(double[] vector, double[,] cosines, int lowerDctBound)
{
var dctArray = DataTools.SubtractMean(vector);
int dctLength = dctArray.Length;
double[] dctCoeff = MFCCStuff.CalculateCeptrum(dctArray, cosines);

// convert to absolute values because not interested in negative values due to phase.
for (int i = 0; i < dctLength; i++)
{
dctCoeff[i] = Math.Abs(dctCoeff[i]);
}

// remove lower coefficients from consideration because they dominate
for (int i = 0; i < lowerDctBound; i++)
{
dctCoeff[i] = 0.0;
}

dctCoeff = DataTools.normalise2UnitLength(dctCoeff);
return dctCoeff;
}
}
}
44 changes: 4 additions & 40 deletions src/AudioAnalysisTools/DSP/MFCCStuff.cs
Original file line number Diff line number Diff line change
Expand Up @@ -477,21 +477,13 @@ public static double InverseHerzTranform(double m, double c, double div)
int binCount = spectra.GetLength(1); //number of filters in filter bank

//set up the cosine coefficients. Need one extra to compensate for DC coeff.
double[,] cosines = Cosines(binCount, coeffCount + 1);

//following two lines write matrix of cos values for checking.
//string fPath = @"path\cosines.txt";
//FileTools.WriteMatrix2File_Formatted(cosines, fPath, "F3");

//following two lines write bmp image of cos values for checking.
//string fPath = @"path\cosines.bmp";
//ImageTools.DrawMatrix(cosines, fPath);
double[,] cosines = DctMethods.Cosines(binCount, coeffCount + 1);

double[,] op = new double[frameCount, coeffCount];
for (int i = 0; i < frameCount; i++)
{
double[] spectrum = DataTools.GetRow(spectra, i); //transfer matrix row=i to vector
double[] cepstrum = DCT(spectrum, cosines);
double[] cepstrum = CalculateCeptrum(spectrum, cosines);

for (int j = 0; j < coeffCount; j++)
{
Expand All @@ -512,7 +504,7 @@ public static double InverseHerzTranform(double m, double c, double div)
for (int i = 0; i < frameCount; i++)
{
double[] spectrum = DataTools.GetRow(spectra, i); //transfer matrix row=i to vector
double[] cepstrum = DCT(spectrum, cosines);
double[] cepstrum = CalculateCeptrum(spectrum, cosines);

for (int j = 0; j < coeffCount; j++)
{
Expand All @@ -523,35 +515,7 @@ public static double InverseHerzTranform(double m, double c, double div)
return op;
}

/// <summary>
/// Returns a matrix of cosine basis functions.
/// These are prepared prior to performing a DCT, Discrete Cosine Transform.
/// The rows k = 0 to coeffCount are the basis functions.
/// The columns, m = 0 to M where M = signalLength or the length of the required DCT.
/// The value of m/M ranges from 0 to 1.0.
/// The value of Pi*m/M ranges from 0 to Pi radians.
/// The value of k*Pi*m/M ranges from 0 to k*Pi radians. WHen k=2, 2Pi radians corresponds to one rotation.
/// </summary>
/// <param name="signalLength">The length of the signal to be processed. e.g. the frequency bin count or filter bank count or ...</param>
/// <param name="coeffCount">The number of basis funcitons = the rquired number of DCT coefficients.</param>
public static double[,] Cosines(int signalLength, int coeffCount)
{
double[,] cosines = new double[coeffCount + 1, signalLength]; //get an extra coefficient because do not want DC coeff at [0].
for (int k = 0; k < coeffCount + 1; k++)
{
double kPiOnM = k * Math.PI / signalLength;

// for each spectral bin
for (int m = 0; m < signalLength; m++)
{
cosines[k, m] = Math.Cos(kPiOnM * (m + 0.5)); //can also be Cos(kPiOnM * (m - 0.5)
}
}

return cosines;
}

public static double[] DCT(double[] spectrum, double[,] cosines)
public static double[] CalculateCeptrum(double[] spectrum, double[,] cosines)
{
int length = spectrum.Length;
int coeffCount = cosines.GetLength(0);
Expand Down
2 changes: 1 addition & 1 deletion src/AudioAnalysisTools/FindMatchingEvents.cs
Original file line number Diff line number Diff line change
Expand Up @@ -596,7 +596,7 @@ public static Tuple<double[]> Execute_MFCC_XCOR(double[,] target, double dynamic
//set up the matrix of cosine coefficients
int coeffCount = 12; //only use first 12 coefficients.
int binCount = target.GetLength(1); //number of filters in filter bank
double[,] cosines = MFCCStuff.Cosines(binCount, coeffCount + 1); //set up the cosine coefficients
double[,] cosines = DctMethods.Cosines(binCount, coeffCount + 1); //set up the cosine coefficients

//adjust target's dynamic range to that set by user
target = SNR.SetDynamicRange(target, 3.0, dynamicRange); //set event's dynamic range
Expand Down
4 changes: 2 additions & 2 deletions src/AudioAnalysisTools/Harmonics/HarmonicAnalysis.cs
Original file line number Diff line number Diff line change
Expand Up @@ -272,7 +272,7 @@ public static Tuple<double[], double[,]> DetectHarmonicsUsingFormantGap(double[,
int cols = matrix.GetLength(1);
double[,] hits = new double[rows, cols];

double[,] cosines = MFCCStuff.Cosines(dctLength, dctLength); //set up the cosine coefficients
double[,] cosines = DctMethods.Cosines(dctLength, dctLength); //set up the cosine coefficients

for (int r = 0; r < rows - dctLength; r++)
{
Expand All @@ -293,7 +293,7 @@ public static Tuple<double[], double[,]> DetectHarmonicsUsingFormantGap(double[,

// DataTools.writeBarGraph(array);

double[] dct = MFCCStuff.DCT(array, cosines);
double[] dct = MFCCStuff.CalculateCeptrum(array, cosines);
for (int i = 0; i < dctLength; i++)
{
dct[i] = Math.Abs(dct[i]); //convert to absolute values
Expand Down
4 changes: 2 additions & 2 deletions src/AudioAnalysisTools/Harmonics/HarmonicParameters.cs
Original file line number Diff line number Diff line change
Expand Up @@ -201,7 +201,7 @@ public static Tuple<double[], double[], int[]> DetectHarmonicsInSpectrogramData(
var binCount = m.GetLength(1);

//set up the cosine coefficients
double[,] cosines = MFCCStuff.Cosines(binCount, binCount);
double[,] cosines = DctMethods.Cosines(binCount, binCount);

// set up time-frame arrays to store decibels, formant intensity and max index.
var dBArray = new double[rowCount];
Expand Down Expand Up @@ -265,7 +265,7 @@ public static Tuple<double[], double[], int[]> DetectHarmonicsInSpectrogramData(
// set the first four values in the returned DCT coefficients to 0.
// We require a minimum of three formants, that is, two harmonic intervals.
int lowerDctBound = 4;
var dctCoefficients = Oscillations2012.DoDct(normXr, cosines, lowerDctBound);
var dctCoefficients = DctMethods.DoDct(normXr, cosines, lowerDctBound);
int indexOfMaxValue = DataTools.GetMaxIndex(dctCoefficients);
intensity[t] = dctCoefficients[indexOfMaxValue];
maxIndexArray[t] = indexOfMaxValue;
Expand Down
2 changes: 1 addition & 1 deletion src/AudioAnalysisTools/Ocillations/Oscillations2014.cs
Original file line number Diff line number Diff line change
Expand Up @@ -916,7 +916,7 @@ public static double[] GetSpectralIndex_Osc(AudioRecording recordingSegment, int
public static void GetOscillationUsingDct(double[] array, double framesPerSecond, double[,] cosines, out double oscilFreq, out double period, out double intenisty)
{
var modifiedArray = DataTools.SubtractMean(array);
var dctCoeff = MFCCStuff.DCT(modifiedArray, cosines);
var dctCoeff = MFCCStuff.CalculateCeptrum(modifiedArray, cosines);

// convert to absolute values because not interested in negative values due to phase.
for (int i = 0; i < dctCoeff.Length; i++)
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -136,7 +136,7 @@ public void TestHarmonicsAlgorithmOn1000HertzHarmonic()
public void TestCosinesMatrixForDct()
{
// get an 8 x 8 matrix.
double[,] cosineBasisFunctions = MFCCStuff.Cosines(8, 8);
double[,] cosineBasisFunctions = DctMethods.Cosines(8, 8);

//following line writes matrix of cos values for checking.
var outputDir = new FileInfo(Path.Join(this.TestOutputDirectory.FullName, "Cosines.csv"));
Expand Down

0 comments on commit 4fc5bdc

Please sign in to comment.