mirror of
				https://github.com/simon987/wavelib.git
				synced 2025-11-04 10:26:52 +00:00 
			
		
		
		
	CCommit : CWT/ICWT added
This commit is contained in:
		
							parent
							
								
									dc38e9f39b
								
							
						
					
					
						commit
						69bf4dc0b0
					
				@ -1,7 +1,7 @@
 | 
			
		||||
wavelib
 | 
			
		||||
=======
 | 
			
		||||
 | 
			
		||||
C Implementation of Wavelet Transform (DWT,SWT and MODWT) and Packet Transform ( Full Tree Decomposition and Best Basis DPWT).
 | 
			
		||||
C Implementation of Wavelet Transform (DWT,SWT and MODWT) and Packet Transform ( Full Tree Decomposition and Best Basis DWPT).
 | 
			
		||||
 | 
			
		||||
Discrete Wavelet Transform Methods Implemented
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
@ -14,6 +14,16 @@ extern "C" {
 | 
			
		||||
#define fft_type double
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
#ifndef cplx_type
 | 
			
		||||
#define cplx_type double
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
typedef struct cplx_t {
 | 
			
		||||
	cplx_type re;
 | 
			
		||||
	cplx_type im;
 | 
			
		||||
} cplx_data;
 | 
			
		||||
 | 
			
		||||
typedef struct wave_set* wave_object;
 | 
			
		||||
 | 
			
		||||
wave_object wave_init(char* wname);
 | 
			
		||||
@ -153,6 +163,34 @@ struct wpt_set{
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
typedef struct cwt_set* cwt_object;
 | 
			
		||||
 | 
			
		||||
cwt_object cwt_init(char* wave, double param, int siglength,double dt, int J);
 | 
			
		||||
 | 
			
		||||
struct cwt_set{
 | 
			
		||||
	char wave[10];// Wavelet - morl/morlet,paul,dog/dgauss
 | 
			
		||||
	int siglength;// Length of Input Data
 | 
			
		||||
	int J;// Total Number of Scales
 | 
			
		||||
	double s0;// Smallest scale. It depends on the sampling rate. s0 <= 2 * dt for most wavelets
 | 
			
		||||
	double dt;// Sampling Rate
 | 
			
		||||
	double dj;// Separation between scales. eg., scale = s0 * 2 ^ ( [0:N-1] *dj ) or scale = s0 *[0:N-1] * dj
 | 
			
		||||
	char type[10];// Scale Type - Power or Linear
 | 
			
		||||
	int pow;// Base of Power in case type = pow. Typical value is pow = 2
 | 
			
		||||
	int sflag;
 | 
			
		||||
	int pflag;
 | 
			
		||||
	int npad;
 | 
			
		||||
	int mother;
 | 
			
		||||
	double m;// Wavelet parameter param
 | 
			
		||||
	double smean;// Input Signal mean
 | 
			
		||||
 | 
			
		||||
	cplx_data *output;
 | 
			
		||||
	double *scale;
 | 
			
		||||
	double *period;
 | 
			
		||||
	double *coi;
 | 
			
		||||
	double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
void dwt(wt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void idwt(wt_object wt, double *dwtop);
 | 
			
		||||
@ -189,6 +227,18 @@ int getDWPTNodelength(wpt_object wt, int X);
 | 
			
		||||
 | 
			
		||||
void getDWPTCoeffs(wpt_object wt, int X, int Y, double *coeffs, int N);
 | 
			
		||||
 | 
			
		||||
void setCWTScales(cwt_object wt, double s0, double dj, char *type, int power);
 | 
			
		||||
 | 
			
		||||
void setCWTScaleVector(cwt_object wt, double *scale, int J, double s0, double dj);
 | 
			
		||||
 | 
			
		||||
void setCWTPadding(cwt_object wt, int pad);
 | 
			
		||||
 | 
			
		||||
void cwt(cwt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void icwt(cwt_object wt, double *cwtop);
 | 
			
		||||
 | 
			
		||||
int getCWTScaleLength(int N);
 | 
			
		||||
 | 
			
		||||
void wave_summary(wave_object obj);
 | 
			
		||||
 | 
			
		||||
void wt_summary(wt_object wt);
 | 
			
		||||
@ -197,6 +247,8 @@ void wtree_summary(wtree_object wt);
 | 
			
		||||
 | 
			
		||||
void wpt_summary(wpt_object wt);
 | 
			
		||||
 | 
			
		||||
void cwt_summary(cwt_object wt);;
 | 
			
		||||
 | 
			
		||||
void wave_free(wave_object object);
 | 
			
		||||
 | 
			
		||||
void wt_free(wt_object object);
 | 
			
		||||
@ -205,6 +257,8 @@ void wtree_free(wtree_object object);
 | 
			
		||||
 | 
			
		||||
void wpt_free(wpt_object object);
 | 
			
		||||
 | 
			
		||||
void cwt_free(cwt_object object);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
@ -3,17 +3,23 @@
 | 
			
		||||
include_directories(${CMAKE_CURRENT_SOURCE_DIR})
 | 
			
		||||
 | 
			
		||||
set(SOURCE_FILES    conv.c
 | 
			
		||||
					cwt.c
 | 
			
		||||
					cwtmath.c
 | 
			
		||||
					hsfft.c
 | 
			
		||||
					real.c
 | 
			
		||||
					wavefilt.c
 | 
			
		||||
					wavefunc.c
 | 
			
		||||
					wavelib.c
 | 
			
		||||
					wtmath.c
 | 
			
		||||
                    )
 | 
			
		||||
 | 
			
		||||
set(HEADER_FILES    conv.h
 | 
			
		||||
					cwt.h
 | 
			
		||||
					cwtmath.h
 | 
			
		||||
					hsfft.h
 | 
			
		||||
					real.h
 | 
			
		||||
					wavefilt.h
 | 
			
		||||
					wavefunc.h
 | 
			
		||||
					wavelib.h
 | 
			
		||||
					wtmath.h
 | 
			
		||||
                    )
 | 
			
		||||
 | 
			
		||||
							
								
								
									
										420
									
								
								src/cwt.c
									
									
									
									
									
										Executable file
									
								
							
							
						
						
									
										420
									
								
								src/cwt.c
									
									
									
									
									
										Executable file
									
								
							@ -0,0 +1,420 @@
 | 
			
		||||
/*
 | 
			
		||||
  Copyright (c) 2015, Rafat Hussain
 | 
			
		||||
*/
 | 
			
		||||
/*
 | 
			
		||||
This code is a C translation ( with some modifications) of Wavelet Software provided by
 | 
			
		||||
C. Torrence and G. Compo, and is available at URL: http://atoc.colorado.edu/research/wavelets/''.
 | 
			
		||||
*/
 | 
			
		||||
 | 
			
		||||
#include "cwt.h"
 | 
			
		||||
 | 
			
		||||
static int factorial2(int N) {
 | 
			
		||||
	int factorial,i;
 | 
			
		||||
 | 
			
		||||
	factorial = 1;
 | 
			
		||||
 | 
			
		||||
	for (i = 1; i <= N;++i) {
 | 
			
		||||
		factorial *= i;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	return factorial;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
static double factorial3(int N) {
 | 
			
		||||
	int i;
 | 
			
		||||
	double factorial;
 | 
			
		||||
 | 
			
		||||
	factorial = 1;
 | 
			
		||||
 | 
			
		||||
	for (i = 1; i <= N; ++i) {
 | 
			
		||||
		factorial *= i;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	return factorial;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
double factorial(int N) {
 | 
			
		||||
	if (N > 40) {
 | 
			
		||||
		printf("This program is only valid for N <= 40 \n");
 | 
			
		||||
		return -1.0;
 | 
			
		||||
	}
 | 
			
		||||
	double fact[41] = { 1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800, 479001600, 6227020800, 87178291200, 1307674368000,
 | 
			
		||||
		20922789888000, 355687428096000, 6402373705728000, 121645100408832000, 2432902008176640000, 51090942171709440000.0, 1124000727777607680000.0,
 | 
			
		||||
		25852016738884976640000.0, 620448401733239439360000.0, 15511210043330985984000000.0, 403291461126605635584000000.0, 10888869450418352160768000000.0,
 | 
			
		||||
		304888344611713860501504000000.0, 8841761993739701954543616000000.0, 265252859812191058636308480000000.0, 8222838654177922817725562880000000.0,
 | 
			
		||||
		263130836933693530167218012160000000.0, 8683317618811886495518194401280000000.0, 295232799039604140847618609643520000000.0, 10333147966386144929666651337523200000000.0,
 | 
			
		||||
		371993326789901217467999448150835200000000.0, 13763753091226345046315979581580902400000000.0, 523022617466601111760007224100074291200000000.0,
 | 
			
		||||
		20397882081197443358640281739902897356800000000.0, 815915283247897734345611269596115894272000000000.0 };
 | 
			
		||||
 | 
			
		||||
	return fact[N];
 | 
			
		||||
 | 
			
		||||
}
 | 
			
		||||
static void wave_function(int nk, double dt,int mother, double param,double scale1, double *kwave, double pi,double *period1,
 | 
			
		||||
	double *coi1, fft_data *daughter) {
 | 
			
		||||
 | 
			
		||||
	double norm, expnt, fourier_factor;
 | 
			
		||||
	int k, m;
 | 
			
		||||
	double temp;
 | 
			
		||||
	int sign,re;
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
	if (mother == 0) {
 | 
			
		||||
		//MORLET
 | 
			
		||||
		if (param < 0.0) {
 | 
			
		||||
			param = 6.0;
 | 
			
		||||
		}
 | 
			
		||||
		norm = sqrt(2.0*pi*scale1 / dt)*pow(pi,-0.25);
 | 
			
		||||
 | 
			
		||||
		for (k = 1; k <= nk / 2 + 1; ++k) {
 | 
			
		||||
			temp = (scale1*kwave[k-1] - param);
 | 
			
		||||
			expnt = -0.5 * temp * temp;
 | 
			
		||||
			daughter[k - 1].re = norm * exp(expnt);
 | 
			
		||||
			daughter[k - 1].im = 0.0;
 | 
			
		||||
		}
 | 
			
		||||
		for (k = nk / 2 + 2; k <= nk; ++k) {
 | 
			
		||||
			daughter[k - 1].re = daughter[k - 1].im = 0.0;
 | 
			
		||||
		}
 | 
			
		||||
		fourier_factor = (4.0*pi) / (param + sqrt(2.0 + param*param));
 | 
			
		||||
		*period1 = scale1*fourier_factor;
 | 
			
		||||
		*coi1 = fourier_factor / sqrt(2.0);
 | 
			
		||||
	}
 | 
			
		||||
	else if (mother == 1) {
 | 
			
		||||
		// PAUL
 | 
			
		||||
		if (param < 0.0) {
 | 
			
		||||
			param = 4.0;
 | 
			
		||||
		}
 | 
			
		||||
		m = (int)param;
 | 
			
		||||
		norm = sqrt(2.0*pi*scale1 / dt)*(pow(2.0,(double)m) / sqrt((double)(m*factorial(2 * m - 1))));
 | 
			
		||||
		for (k = 1; k <= nk / 2 + 1; ++k) {
 | 
			
		||||
			temp = scale1 * kwave[k - 1];
 | 
			
		||||
			expnt = - temp;
 | 
			
		||||
			daughter[k - 1].re = norm * pow(temp,(double)m) * exp(expnt);
 | 
			
		||||
			daughter[k - 1].im = 0.0;
 | 
			
		||||
		}
 | 
			
		||||
		for (k = nk / 2 + 2; k <= nk; ++k) {
 | 
			
		||||
			daughter[k - 1].re = daughter[k - 1].im = 0.0;
 | 
			
		||||
		}
 | 
			
		||||
		fourier_factor = (4.0*pi) / (2.0 * m + 1.0);
 | 
			
		||||
		*period1 = scale1*fourier_factor;
 | 
			
		||||
		*coi1 = fourier_factor * sqrt(2.0);
 | 
			
		||||
	}
 | 
			
		||||
	else if (mother == 2) {
 | 
			
		||||
		if (param < 0.0) {
 | 
			
		||||
			param = 2.0;
 | 
			
		||||
		}
 | 
			
		||||
		m = (int)param;
 | 
			
		||||
 | 
			
		||||
		if (m % 2 == 0) {
 | 
			
		||||
			re = 1;
 | 
			
		||||
		}
 | 
			
		||||
		else {
 | 
			
		||||
			re = 0;
 | 
			
		||||
		}
 | 
			
		||||
 | 
			
		||||
		if (m % 4 == 0 || m % 4 == 1) {
 | 
			
		||||
			sign = -1;
 | 
			
		||||
		}
 | 
			
		||||
		else {
 | 
			
		||||
			sign = 1;
 | 
			
		||||
		}
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
		norm = sqrt(2.0*pi*scale1 / dt)*sqrt(1.0 / gamma(m + 0.50));
 | 
			
		||||
		norm *= sign;
 | 
			
		||||
 | 
			
		||||
		if (re == 1) {
 | 
			
		||||
			for (k = 1; k <= nk; ++k) {
 | 
			
		||||
				temp = scale1 * kwave[k - 1];
 | 
			
		||||
				daughter[k - 1].re = norm*pow(temp,(double)m)*exp(-0.50*pow(temp,2.0));
 | 
			
		||||
				daughter[k - 1].im = 0.0;
 | 
			
		||||
			}
 | 
			
		||||
		}
 | 
			
		||||
		else if (re == 0) {
 | 
			
		||||
			for (k = 1; k <= nk; ++k) {
 | 
			
		||||
				temp = scale1 * kwave[k - 1];
 | 
			
		||||
				daughter[k - 1].re = 0.0;
 | 
			
		||||
				daughter[k - 1].im = norm*pow(temp, (double)m)*exp(-0.50*pow(temp, 2.0));
 | 
			
		||||
			}
 | 
			
		||||
		}
 | 
			
		||||
		fourier_factor = (2.0*pi) * sqrt(2.0 / (2.0 * m + 1.0));
 | 
			
		||||
		*period1 = scale1*fourier_factor;
 | 
			
		||||
		*coi1 = fourier_factor / sqrt(2.0);
 | 
			
		||||
	}
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void cwavelet(double *y, int N, double dt, int mother, double param, double s0, double dj, int jtot, int npad,
 | 
			
		||||
	double *wave, double *scale, double *period, double *coi) {
 | 
			
		||||
 | 
			
		||||
	int i, j, k, iter;
 | 
			
		||||
	double ymean, freq1, pi, period1, coi1;
 | 
			
		||||
	double tmp1, tmp2;
 | 
			
		||||
	double scale1;
 | 
			
		||||
	double *kwave;
 | 
			
		||||
	fft_object obj, iobj;
 | 
			
		||||
	fft_data *ypad, *yfft,*daughter;
 | 
			
		||||
 | 
			
		||||
	pi = 4.0 * atan(1.0);
 | 
			
		||||
 | 
			
		||||
	if (npad < N) {
 | 
			
		||||
		printf("npad must be >= N \n");
 | 
			
		||||
		exit(-1);
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	obj = fft_init(npad, 1);
 | 
			
		||||
	iobj = fft_init(npad, -1);
 | 
			
		||||
 | 
			
		||||
	ypad = (fft_data*)malloc(sizeof(fft_data)* npad);
 | 
			
		||||
	yfft = (fft_data*)malloc(sizeof(fft_data)* npad);
 | 
			
		||||
	daughter = (fft_data*)malloc(sizeof(fft_data)* npad);
 | 
			
		||||
	kwave = (double*)malloc(sizeof(double)* npad);
 | 
			
		||||
 | 
			
		||||
	ymean = 0.0;
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		ymean += y[i];
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	ymean /= N;
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		ypad[i].re = y[i] - ymean;
 | 
			
		||||
		ypad[i].im = 0.0;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = N; i < npad; ++i) {
 | 
			
		||||
		ypad[i].re = ypad[i].im = 0.0;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
	// Find FFT of the input y (ypad)
 | 
			
		||||
 | 
			
		||||
	fft_exec(obj, ypad, yfft);
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < npad; ++i) {
 | 
			
		||||
		yfft[i].re /= (double) npad;
 | 
			
		||||
		yfft[i].im /= (double) npad;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
	//Construct the wavenumber array
 | 
			
		||||
 | 
			
		||||
	freq1 = 2.0*pi / ((double)npad*dt);
 | 
			
		||||
	kwave[0] = 0.0;
 | 
			
		||||
 | 
			
		||||
	for (i = 1; i < npad / 2 + 1; ++i) {
 | 
			
		||||
		kwave[i] = i * freq1;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = npad / 2 + 1; i < npad; ++i) {
 | 
			
		||||
		kwave[i] = -kwave[npad - i ];
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	
 | 
			
		||||
	// Main loop
 | 
			
		||||
 | 
			
		||||
	for (j = 1; j <= jtot; ++j) {
 | 
			
		||||
		scale1 = scale[j - 1];// = s0*pow(2.0, (double)(j - 1)*dj);
 | 
			
		||||
		wave_function(npad, dt, mother, param, scale1, kwave, pi,&period1,&coi1, daughter);
 | 
			
		||||
		period[j - 1] = period1;
 | 
			
		||||
		for (k = 0; k < npad; ++k) {
 | 
			
		||||
			tmp1 = daughter[k].re * yfft[k].re - daughter[k].im * yfft[k].im;
 | 
			
		||||
			tmp2 = daughter[k].re * yfft[k].im + daughter[k].im * yfft[k].re;
 | 
			
		||||
			daughter[k].re = tmp1;
 | 
			
		||||
			daughter[k].im = tmp2;
 | 
			
		||||
		}
 | 
			
		||||
		fft_exec(iobj, daughter, ypad);
 | 
			
		||||
		iter = 2 * (j - 1) * N;
 | 
			
		||||
		for (i = 0; i < N; ++i) {
 | 
			
		||||
			wave[iter + 2 * i] = ypad[i].re;
 | 
			
		||||
			wave[iter + 2 * i + 1] = ypad[i].im;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
 | 
			
		||||
	for (i = 1; i <= (N + 1) / 2; ++i) {
 | 
			
		||||
		coi[i - 1] = coi1 * dt * ((double)i - 1.0);
 | 
			
		||||
		coi[N - i] = coi[i - 1];
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
 | 
			
		||||
	
 | 
			
		||||
	free(kwave);
 | 
			
		||||
	free(ypad);
 | 
			
		||||
	free(yfft);
 | 
			
		||||
	free(daughter);
 | 
			
		||||
 | 
			
		||||
	free_fft(obj);
 | 
			
		||||
	free_fft(iobj);
 | 
			
		||||
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void psi0(int mother, double param,double *val,int *real) {
 | 
			
		||||
	double pi,coeff;
 | 
			
		||||
	int m,sign;
 | 
			
		||||
 | 
			
		||||
	m = (int)param;
 | 
			
		||||
	pi = 4.0 * atan(1.0);
 | 
			
		||||
 | 
			
		||||
	if (mother == 0) {
 | 
			
		||||
		// Morlet
 | 
			
		||||
		*val = 1.0 / sqrt(sqrt(pi));
 | 
			
		||||
		*real = 1;
 | 
			
		||||
	}
 | 
			
		||||
	else if (mother == 1) {
 | 
			
		||||
		//Paul
 | 
			
		||||
		if (m % 2 == 0) {
 | 
			
		||||
			*real = 1;
 | 
			
		||||
		}
 | 
			
		||||
		else {
 | 
			
		||||
			*real = 0;
 | 
			
		||||
		}
 | 
			
		||||
 | 
			
		||||
		if (m % 4 == 0 || m % 4 == 1) {
 | 
			
		||||
			sign = 1;
 | 
			
		||||
		}
 | 
			
		||||
		else {
 | 
			
		||||
			sign = -1;
 | 
			
		||||
		}
 | 
			
		||||
		*val = sign * pow(2.0, (double)m) * factorial(m) / (sqrt(pi * factorial(2 * m)));
 | 
			
		||||
 | 
			
		||||
	}
 | 
			
		||||
	else if (mother == 2) {
 | 
			
		||||
		// D.O.G
 | 
			
		||||
		*real = 1;
 | 
			
		||||
 | 
			
		||||
		if (m % 2 == 0) {
 | 
			
		||||
			if (m % 4 == 0) {
 | 
			
		||||
				sign = -1;
 | 
			
		||||
			}
 | 
			
		||||
			else {
 | 
			
		||||
				sign = 1;
 | 
			
		||||
			}
 | 
			
		||||
			coeff = sign * pow(2.0, (double)m / 2) / gamma(0.5);
 | 
			
		||||
			*val = coeff * gamma(((double)m + 1.0) / 2.0) / sqrt(gamma(m + 0.50));
 | 
			
		||||
		}
 | 
			
		||||
		else {
 | 
			
		||||
			*val = 0;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
static int maxabs(double *array,int N) {
 | 
			
		||||
	double maxval,temp;
 | 
			
		||||
	int i,index;
 | 
			
		||||
	maxval = 0.0;
 | 
			
		||||
	index = -1;
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		temp = fabs(array[i]);
 | 
			
		||||
		if (temp >= maxval) {
 | 
			
		||||
			maxval = temp;
 | 
			
		||||
			index = i;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	return index;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
double cdelta(int mother, double param, double psi0 ) {
 | 
			
		||||
	int N,i,j,iter;
 | 
			
		||||
	double *delta, *scale,*period,*wave,*coi,*mval;
 | 
			
		||||
	double den,cdel;
 | 
			
		||||
	double subscale,dt,dj,s0;
 | 
			
		||||
	int jtot;
 | 
			
		||||
	int maxarr;
 | 
			
		||||
 | 
			
		||||
	subscale = 8.0;
 | 
			
		||||
	dt = 0.25;
 | 
			
		||||
	if (mother == 0) {
 | 
			
		||||
		N = 16;
 | 
			
		||||
		s0 = dt/4;
 | 
			
		||||
	}
 | 
			
		||||
	else if (mother == 1) {
 | 
			
		||||
		N = 16;
 | 
			
		||||
		s0 = dt / 4.0;
 | 
			
		||||
	}
 | 
			
		||||
	else if (mother == 2)
 | 
			
		||||
	{
 | 
			
		||||
		s0 = dt/8.0;
 | 
			
		||||
		N = 256;
 | 
			
		||||
		if (param == 2.0) {
 | 
			
		||||
			subscale = 16.0;
 | 
			
		||||
			s0 = dt / 16.0;
 | 
			
		||||
			N = 2048;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	dj = 1.0 / subscale;
 | 
			
		||||
	jtot = 16 * (int) subscale;
 | 
			
		||||
 | 
			
		||||
	delta = (double*)malloc(sizeof(double)* N);
 | 
			
		||||
	wave = (double*)malloc(sizeof(double)* 2 * N * jtot);
 | 
			
		||||
	coi = (double*)malloc(sizeof(double)* N);
 | 
			
		||||
	scale = (double*)malloc(sizeof(double)* jtot);
 | 
			
		||||
	period = (double*)malloc(sizeof(double)* jtot);
 | 
			
		||||
	mval = (double*)malloc(sizeof(double)* N);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
	delta[0] = 1;
 | 
			
		||||
 | 
			
		||||
	for (i = 1; i < N; ++i) {
 | 
			
		||||
		delta[i] = 0;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < jtot; ++i) {
 | 
			
		||||
		scale[i] = s0*pow(2.0, (double)(i)*dj);
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	cwavelet(delta, N, dt, mother, param, s0, dj, jtot, N, wave, scale, period, coi);
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		mval[i] = 0;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (j = 0; j < jtot; ++j) {
 | 
			
		||||
		iter = 2 * j * N;
 | 
			
		||||
		den = sqrt(scale[j]);
 | 
			
		||||
		for (i = 0; i < N; ++i) {
 | 
			
		||||
			mval[i] += wave[iter + 2 * i]/den;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
	maxarr = maxabs(mval, N);
 | 
			
		||||
 | 
			
		||||
	cdel = sqrt(dt) * dj * mval[maxarr] / psi0;
 | 
			
		||||
 | 
			
		||||
	free(delta);
 | 
			
		||||
	free(wave);
 | 
			
		||||
 | 
			
		||||
	free(scale);
 | 
			
		||||
	free(period);
 | 
			
		||||
	free(coi);
 | 
			
		||||
	free(mval);
 | 
			
		||||
 | 
			
		||||
	return cdel;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void icwavelet(double *wave, int N, double *scale,int jtot,double dt,double dj,double cdelta,double psi0,double *oup) {
 | 
			
		||||
	int i, j,iter;
 | 
			
		||||
	double den, coeff;
 | 
			
		||||
 | 
			
		||||
	coeff = sqrt(dt) * dj / (cdelta *psi0);
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		oup[i] = 0.0;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (j = 0; j < jtot; ++j) {
 | 
			
		||||
		iter = 2 * j * N;
 | 
			
		||||
		den = sqrt(scale[j]);
 | 
			
		||||
		for (i = 0; i < N; ++i) {
 | 
			
		||||
			oup[i] += wave[iter + 2 * i] / den;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		oup[i] *= coeff;
 | 
			
		||||
	}
 | 
			
		||||
}
 | 
			
		||||
							
								
								
									
										29
									
								
								src/cwt.h
									
									
									
									
									
										Executable file
									
								
							
							
						
						
									
										29
									
								
								src/cwt.h
									
									
									
									
									
										Executable file
									
								
							@ -0,0 +1,29 @@
 | 
			
		||||
#ifndef CWT_H_
 | 
			
		||||
#define CWT_H_
 | 
			
		||||
 | 
			
		||||
#include "wavefunc.h"
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
extern "C" {
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
void cwavelet(double *y, int N, double dt, int mother, double param, double s0, double dj, int jtot, int npad,
 | 
			
		||||
		double *wave, double *scale, double *period, double *coi);
 | 
			
		||||
 | 
			
		||||
void psi0(int mother, double param, double *val, int *real);
 | 
			
		||||
 | 
			
		||||
double factorial(int N);
 | 
			
		||||
 | 
			
		||||
double cdelta(int mother, double param, double psi0);
 | 
			
		||||
 | 
			
		||||
void icwavelet(double *wave, int N, double *scale, int jtot, double dt, double dj, double cdelta, double psi0, double *oup);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
}
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#endif /* WAVELIB_H_ */
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
							
								
								
									
										296
									
								
								src/cwtmath.c
									
									
									
									
									
										Executable file
									
								
							
							
						
						
									
										296
									
								
								src/cwtmath.c
									
									
									
									
									
										Executable file
									
								
							@ -0,0 +1,296 @@
 | 
			
		||||
#include "cwtmath.h"
 | 
			
		||||
 | 
			
		||||
static void nsfft_fd(fft_object obj, fft_data *inp, fft_data *oup,double lb,double ub,double *w) {
 | 
			
		||||
	int M,N,i,j,L;
 | 
			
		||||
	double delta,den,theta,tempr,tempi,plb;
 | 
			
		||||
	double *temp1,*temp2;
 | 
			
		||||
 | 
			
		||||
	N = obj->N;
 | 
			
		||||
	L = N/2;
 | 
			
		||||
	//w = (double*)malloc(sizeof(double)*N);
 | 
			
		||||
	
 | 
			
		||||
	M = divideby(N, 2);
 | 
			
		||||
	
 | 
			
		||||
	if (M == 0) {
 | 
			
		||||
		printf("The Non-Standard FFT Length must be a power of 2");
 | 
			
		||||
		exit(1);
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	temp1 = (double*)malloc(sizeof(double)*L);
 | 
			
		||||
	temp2 = (double*)malloc(sizeof(double)*L);
 | 
			
		||||
	
 | 
			
		||||
	delta = (ub - lb)/ N;
 | 
			
		||||
	j = -N;
 | 
			
		||||
	den = 2 * (ub-lb);
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		w[i] = (double)j/den;
 | 
			
		||||
		j += 2;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	fft_exec(obj,inp,oup);
 | 
			
		||||
	
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < L; ++i) {
 | 
			
		||||
		temp1[i] = oup[i].re;
 | 
			
		||||
		temp2[i] = oup[i].im;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N - L; ++i) {
 | 
			
		||||
		oup[i].re = oup[i + L].re;
 | 
			
		||||
		oup[i].im = oup[i + L].im;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < L; ++i) {
 | 
			
		||||
		oup[N - L + i].re = temp1[i];
 | 
			
		||||
		oup[N - L + i].im = temp2[i];
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	plb = PI2 * lb;
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		tempr = oup[i].re;
 | 
			
		||||
		tempi = oup[i].im;
 | 
			
		||||
		theta = w[i] * plb;
 | 
			
		||||
		
 | 
			
		||||
		oup[i].re = delta * (tempr*cos(theta) + tempi*sin(theta));
 | 
			
		||||
		oup[i].im = delta * (tempi*cos(theta) - tempr*sin(theta));
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	
 | 
			
		||||
	//free(w);
 | 
			
		||||
	free(temp1);
 | 
			
		||||
	free(temp2);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
static void nsfft_bk(fft_object obj, fft_data *inp, fft_data *oup,double lb,double ub,double *t) {
 | 
			
		||||
	int M,N,i,j,L;
 | 
			
		||||
	double *w;
 | 
			
		||||
	double delta,den,plb,theta;
 | 
			
		||||
	double *temp1,*temp2;
 | 
			
		||||
	fft_data *inpt;
 | 
			
		||||
 | 
			
		||||
	N = obj->N;
 | 
			
		||||
	L = N/2;
 | 
			
		||||
	
 | 
			
		||||
	M = divideby(N, 2);
 | 
			
		||||
	
 | 
			
		||||
	if (M == 0) {
 | 
			
		||||
		printf("The Non-Standard FFT Length must be a power of 2");
 | 
			
		||||
		exit(1);
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	temp1 = (double*)malloc(sizeof(double)*L);
 | 
			
		||||
	temp2 = (double*)malloc(sizeof(double)*L);
 | 
			
		||||
	w = (double*)malloc(sizeof(double)*N);
 | 
			
		||||
	inpt = (fft_data*) malloc (sizeof(fft_data) * N);
 | 
			
		||||
	
 | 
			
		||||
	delta = (ub - lb)/ N;
 | 
			
		||||
	j = -N;
 | 
			
		||||
	den = 2 * (ub-lb);
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		w[i] = (double)j/den;
 | 
			
		||||
		j += 2;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	plb = PI2 * lb;
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		theta = w[i] * plb;
 | 
			
		||||
		
 | 
			
		||||
		inpt[i].re = (inp[i].re*cos(theta) - inp[i].im*sin(theta))/delta;
 | 
			
		||||
		inpt[i].im = (inp[i].im*cos(theta) + inp[i].re*sin(theta))/delta;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	for (i = 0; i < L; ++i) {
 | 
			
		||||
		temp1[i] = inpt[i].re;
 | 
			
		||||
		temp2[i] = inpt[i].im;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N - L; ++i) {
 | 
			
		||||
		inpt[i].re = inpt[i + L].re;
 | 
			
		||||
		inpt[i].im = inpt[i + L].im;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < L; ++i) {
 | 
			
		||||
		inpt[N - L + i].re = temp1[i];
 | 
			
		||||
		inpt[N - L + i].im = temp2[i];
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	fft_exec(obj,inpt,oup);
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		t[i] = lb + i*delta;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	free(w);
 | 
			
		||||
	free(temp1);
 | 
			
		||||
	free(temp2);
 | 
			
		||||
	free(inpt);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void nsfft_exec(fft_object obj, fft_data *inp, fft_data *oup,double lb,double ub,double *w) {
 | 
			
		||||
	if (obj->sgn == 1) {
 | 
			
		||||
		nsfft_fd(obj,inp,oup,lb,ub,w);
 | 
			
		||||
	} else if (obj->sgn == -1) {
 | 
			
		||||
		nsfft_bk(obj,inp,oup,lb,ub,w);
 | 
			
		||||
	}
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
static double fix(double x) {
 | 
			
		||||
	// Rounds to the integer nearest to zero 
 | 
			
		||||
	if (x >= 0.) {
 | 
			
		||||
		return floor(x);
 | 
			
		||||
	} else {
 | 
			
		||||
		return ceil(x);
 | 
			
		||||
	}
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
int nint(double N) {
 | 
			
		||||
	int i;
 | 
			
		||||
 | 
			
		||||
	i = (int)(N + 0.49999);
 | 
			
		||||
 | 
			
		||||
	return i;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
double gamma(double x) {
 | 
			
		||||
	/*
 | 
			
		||||
	 * This C program code is based on  W J Cody's fortran code.
 | 
			
		||||
	 * http://www.netlib.org/specfun/gamma
 | 
			
		||||
	 * 
 | 
			
		||||
	 * References:
 | 
			
		||||
   "An Overview of Software Development for Special Functions",
 | 
			
		||||
	W. J. Cody, Lecture Notes in Mathematics, 506,
 | 
			
		||||
	Numerical Analysis Dundee, 1975, G. A. Watson (ed.),
 | 
			
		||||
	Springer Verlag, Berlin, 1976.
 | 
			
		||||
 | 
			
		||||
   Computer Approximations, Hart, Et. Al., Wiley and sons, New York, 1968.
 | 
			
		||||
   */ 
 | 
			
		||||
   
 | 
			
		||||
	// numerator and denominator coefficients for 1 <= x <= 2
 | 
			
		||||
	
 | 
			
		||||
	double y,oup,fact,sum,y2,yi,z,nsum,dsum;
 | 
			
		||||
	int swi,n,i;
 | 
			
		||||
	
 | 
			
		||||
	double spi = 0.9189385332046727417803297;
 | 
			
		||||
    double pi  = 3.1415926535897932384626434;
 | 
			
		||||
	double xmax = 171.624e+0;
 | 
			
		||||
	double xinf = 1.79e308;
 | 
			
		||||
	double eps = 2.22e-16;
 | 
			
		||||
	double xninf = 1.79e-308;
 | 
			
		||||
	
 | 
			
		||||
	double num[8] = { -1.71618513886549492533811e+0,
 | 
			
		||||
                        2.47656508055759199108314e+1,
 | 
			
		||||
                       -3.79804256470945635097577e+2,
 | 
			
		||||
                        6.29331155312818442661052e+2,
 | 
			
		||||
                        8.66966202790413211295064e+2,
 | 
			
		||||
                       -3.14512729688483675254357e+4,
 | 
			
		||||
                       -3.61444134186911729807069e+4,
 | 
			
		||||
                        6.64561438202405440627855e+4 };
 | 
			
		||||
                       
 | 
			
		||||
    double den[8] = { -3.08402300119738975254353e+1,
 | 
			
		||||
                        3.15350626979604161529144e+2,
 | 
			
		||||
                       -1.01515636749021914166146e+3,
 | 
			
		||||
                       -3.10777167157231109440444e+3,
 | 
			
		||||
                        2.25381184209801510330112e+4,
 | 
			
		||||
                        4.75584627752788110767815e+3,
 | 
			
		||||
                       -1.34659959864969306392456e+5,
 | 
			
		||||
                       -1.15132259675553483497211e+5 }; 
 | 
			
		||||
                       
 | 
			
		||||
    // Coefficients for Hart's Minimax approximation x >= 12       
 | 
			
		||||
    
 | 
			
		||||
    
 | 
			
		||||
	double c[7] = { -1.910444077728e-03,
 | 
			
		||||
                        8.4171387781295e-04,
 | 
			
		||||
                       -5.952379913043012e-04,
 | 
			
		||||
                        7.93650793500350248e-04,
 | 
			
		||||
                       -2.777777777777681622553e-03,
 | 
			
		||||
                        8.333333333333333331554247e-02,
 | 
			
		||||
                        5.7083835261e-03 };            
 | 
			
		||||
    
 | 
			
		||||
    y = x;             
 | 
			
		||||
    swi = 0; 
 | 
			
		||||
    fact = 1.0;
 | 
			
		||||
    n = 0;      
 | 
			
		||||
    
 | 
			
		||||
    
 | 
			
		||||
    if ( y < 0.) {
 | 
			
		||||
		// Negative x
 | 
			
		||||
		y = -x;
 | 
			
		||||
		yi = fix(y);
 | 
			
		||||
		oup = y - yi;
 | 
			
		||||
		
 | 
			
		||||
		if (oup != 0.0) {
 | 
			
		||||
			if (yi != fix(yi * .5) * 2.) {
 | 
			
		||||
				swi = 1;
 | 
			
		||||
			}
 | 
			
		||||
			fact = -pi / sin(pi * oup);
 | 
			
		||||
			y += 1.;
 | 
			
		||||
		} else {
 | 
			
		||||
			return xinf;
 | 
			
		||||
		}
 | 
			
		||||
	}      
 | 
			
		||||
	
 | 
			
		||||
	if (y < eps) {
 | 
			
		||||
		if (y >= xninf) {
 | 
			
		||||
			oup = 1.0/y;
 | 
			
		||||
		} else {
 | 
			
		||||
			return xinf;
 | 
			
		||||
		}
 | 
			
		||||
		
 | 
			
		||||
	} else if (y < 12.) {
 | 
			
		||||
		yi = y;
 | 
			
		||||
		if ( y < 1.) {
 | 
			
		||||
			z = y;
 | 
			
		||||
			y += 1.;
 | 
			
		||||
		} else {
 | 
			
		||||
			n = ( int ) y - 1;
 | 
			
		||||
            y -= ( double ) n;
 | 
			
		||||
            z = y - 1.0;
 | 
			
		||||
		}
 | 
			
		||||
		nsum = 0.;
 | 
			
		||||
		dsum = 1.;
 | 
			
		||||
		for (i = 0; i < 8; ++i) {
 | 
			
		||||
			nsum = (nsum + num[i]) * z;
 | 
			
		||||
			dsum = dsum * z + den[i];
 | 
			
		||||
		}
 | 
			
		||||
		oup = nsum / dsum + 1.;
 | 
			
		||||
		
 | 
			
		||||
		if (yi < y) {
 | 
			
		||||
	   
 | 
			
		||||
			oup /= yi;
 | 
			
		||||
		} else if (yi > y) {
 | 
			
		||||
	    
 | 
			
		||||
			for (i = 0; i < n; ++i) {
 | 
			
		||||
				oup *= y;
 | 
			
		||||
				y += 1.;
 | 
			
		||||
			}
 | 
			
		||||
		
 | 
			
		||||
		}      
 | 
			
		||||
	
 | 
			
		||||
	} else {
 | 
			
		||||
		if (y <= xmax) {
 | 
			
		||||
			y2 = y * y;
 | 
			
		||||
			sum = c[6];
 | 
			
		||||
			for (i = 0; i < 6; ++i) {
 | 
			
		||||
				sum = sum / y2 + c[i];
 | 
			
		||||
			}
 | 
			
		||||
			sum = sum / y - y + spi;
 | 
			
		||||
			sum += (y - .5) * log(y);
 | 
			
		||||
			oup = exp(sum);
 | 
			
		||||
		} else {
 | 
			
		||||
			return(xinf);
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	if (swi) {
 | 
			
		||||
		oup = -oup;
 | 
			
		||||
	}
 | 
			
		||||
    if (fact != 1.) {
 | 
			
		||||
		oup = fact / oup; 
 | 
			
		||||
	}       
 | 
			
		||||
	
 | 
			
		||||
	return oup;
 | 
			
		||||
}
 | 
			
		||||
							
								
								
									
										22
									
								
								src/cwtmath.h
									
									
									
									
									
										Executable file
									
								
							
							
						
						
									
										22
									
								
								src/cwtmath.h
									
									
									
									
									
										Executable file
									
								
							@ -0,0 +1,22 @@
 | 
			
		||||
#ifndef CWTMATH_H_
 | 
			
		||||
#define CWTMATH_H_
 | 
			
		||||
 | 
			
		||||
#include "wtmath.h"
 | 
			
		||||
#include "hsfft.h"
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
extern "C" {
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
void nsfft_exec(fft_object obj, fft_data *inp, fft_data *oup,double lb,double ub,double *w);// lb -lower bound, ub - upper bound, w - time or frequency grid (Size N)
 | 
			
		||||
 | 
			
		||||
double gamma(double x);
 | 
			
		||||
 | 
			
		||||
int nint(double N);
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
}
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#endif /* WAVELIB_H_ */
 | 
			
		||||
							
								
								
									
										210
									
								
								src/wavefunc.c
									
									
									
									
									
										Executable file
									
								
							
							
						
						
									
										210
									
								
								src/wavefunc.c
									
									
									
									
									
										Executable file
									
								
							@ -0,0 +1,210 @@
 | 
			
		||||
#include "wavefunc.h"
 | 
			
		||||
 | 
			
		||||
void meyer(int N,double lb,double ub,double *phi,double *psi,double *tgrid) {
 | 
			
		||||
	int M,i;
 | 
			
		||||
	double *w;
 | 
			
		||||
	double delta,j;
 | 
			
		||||
	double theta,x,x2,x3,x4,v,cs,sn;
 | 
			
		||||
	double wf;
 | 
			
		||||
	fft_data *phiw,*psiw,*oup;
 | 
			
		||||
	fft_object obj;
 | 
			
		||||
	
 | 
			
		||||
	M = divideby(N, 2);
 | 
			
		||||
	
 | 
			
		||||
	if (M == 0) {
 | 
			
		||||
		printf("Size of Wavelet must be a power of 2");
 | 
			
		||||
		exit(1);
 | 
			
		||||
	}
 | 
			
		||||
	if (lb >= ub) {
 | 
			
		||||
		printf("upper bound must be greater than lower bound");
 | 
			
		||||
		exit(1);
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	obj = fft_init(N,-1);
 | 
			
		||||
	w = (double*)malloc(sizeof(double)*N);
 | 
			
		||||
	phiw = (fft_data*) malloc (sizeof(fft_data) * N);
 | 
			
		||||
	psiw = (fft_data*) malloc (sizeof(fft_data) * N);
 | 
			
		||||
	oup = (fft_data*) malloc (sizeof(fft_data) * N);
 | 
			
		||||
	
 | 
			
		||||
	delta = 2 * (ub-lb) / PI2;
 | 
			
		||||
	
 | 
			
		||||
	j = (double) N;
 | 
			
		||||
	j *= -1.0;
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		w[i] = j / delta;
 | 
			
		||||
		j += 2.0; 
 | 
			
		||||
		psiw[i].re = psiw[i].im = 0.0;
 | 
			
		||||
		phiw[i].re = phiw[i].im = 0.0;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		wf = fabs(w[i]);
 | 
			
		||||
		if (wf <= PI2/3.0) {
 | 
			
		||||
			phiw[i].re = 1.0;
 | 
			
		||||
		}
 | 
			
		||||
		if (wf > PI2/3.0 && wf <= 2 * PI2 / 3.0) {
 | 
			
		||||
			x = (3 * wf / PI2) - 1.0;
 | 
			
		||||
			x2 = x*x;
 | 
			
		||||
			x3 = x2 * x;
 | 
			
		||||
			x4 = x3 *x;
 | 
			
		||||
			v = x4 *(35 - 84*x + 70*x2 - 20*x3);
 | 
			
		||||
			theta = v * PI2 / 4.0;
 | 
			
		||||
			cs = cos(theta);
 | 
			
		||||
			sn = sin(theta);
 | 
			
		||||
			
 | 
			
		||||
			phiw[i].re = cs;
 | 
			
		||||
			psiw[i].re = cos(w[i]/2.0) * sn;
 | 
			
		||||
			psiw[i].im = sin(w[i]/2.0) * sn;
 | 
			
		||||
		}
 | 
			
		||||
		if (wf > 2.0 * PI2/3.0 && wf <= 4 * PI2 / 3.0) {
 | 
			
		||||
			x = (1.5 * wf / PI2) - 1.0;
 | 
			
		||||
			x2 = x*x;
 | 
			
		||||
			x3 = x2 * x;
 | 
			
		||||
			x4 = x3 *x;
 | 
			
		||||
			v = x4 *(35 - 84*x + 70*x2 - 20*x3);
 | 
			
		||||
			theta = v * PI2 / 4.0;
 | 
			
		||||
			cs = cos(theta);
 | 
			
		||||
			
 | 
			
		||||
			psiw[i].re = cos(w[i]/2.0) * cs;
 | 
			
		||||
			psiw[i].im = sin(w[i]/2.0) * cs;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
	nsfft_exec(obj,phiw,oup,lb,ub,tgrid);
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		phi[i] = oup[i].re/N;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	nsfft_exec(obj,psiw,oup,lb,ub,tgrid);
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
		psi[i] = oup[i].re/N;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
	free(oup);
 | 
			
		||||
	free(phiw);
 | 
			
		||||
	free(psiw);
 | 
			
		||||
	free(w);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void gauss(int N,int p,double lb,double ub,double *psi,double *t) {
 | 
			
		||||
  double delta,num,den,t2,t4;
 | 
			
		||||
	int i;
 | 
			
		||||
	
 | 
			
		||||
	if (lb >= ub) {
 | 
			
		||||
		printf("upper bound must be greater than lower bound");
 | 
			
		||||
		exit(1);
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	t[0] = lb;
 | 
			
		||||
	t[N-1] = ub;
 | 
			
		||||
	delta = (ub - lb) / (N-1);
 | 
			
		||||
	for(i = 1; i < N-1;++i) {
 | 
			
		||||
		t[i] = lb + delta * i;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	den = sqrt(gamma(p+0.5));
 | 
			
		||||
	
 | 
			
		||||
	if ((p+1)%2 == 0) {
 | 
			
		||||
		num = 1.0;
 | 
			
		||||
	} else {
 | 
			
		||||
		num = -1.0;
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	num /= den;
 | 
			
		||||
	
 | 
			
		||||
	//printf("\n%g\n",num);
 | 
			
		||||
	
 | 
			
		||||
	if (p == 1) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			psi[i] = -t[i] * exp(- t[i] * t[i]/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 2) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			psi[i] = (-1.0 + t2) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 3) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			psi[i] = t[i] * (3.0 - t2) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 4) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			psi[i] = (t2 * t2 - 6.0 * t2 + 3.0) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 5) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			psi[i] = t[i] * (-t2 * t2 + 10.0 * t2 - 15.0) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 6) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			psi[i] = (t2 * t2 * t2 - 15.0 * t2 * t2 + 45.0 * t2 - 15.0) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 7) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			psi[i] = t[i] * (-t2 * t2 * t2 + 21.0 * t2 * t2 - 105.0 * t2 + 105.0) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 8) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			t4 = t2 * t2;
 | 
			
		||||
			psi[i] = (t4 * t4 - 28.0 * t4 * t2 + 210.0 * t4 - 420.0 * t2 + 105.0) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 9) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			t4 = t2 * t2;
 | 
			
		||||
			psi[i] = t[i] * (- t4 * t4 + 36.0 * t4 * t2 - 378.0 * t4 + 1260.0 * t2 - 945.0) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else if (p == 10) {
 | 
			
		||||
		for(i = 0; i < N;++i) {
 | 
			
		||||
			t2 = t[i] * t[i];
 | 
			
		||||
			t4 = t2 * t2;
 | 
			
		||||
			psi[i] = (t4 * t4 * t2 - 45.0 * t4 * t4 + 630.0 * t4 * t2 - 3150.0 * t4 + 4725.0 * t2 - 945.0) * exp(- t2/2.0) * num; 
 | 
			
		||||
		}
 | 
			
		||||
	} else {
 | 
			
		||||
	  printf("\n The Gaussian Derivative Wavelet is only available for Derivatives 1 to 10");
 | 
			
		||||
	  exit(1);
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	
 | 
			
		||||
	
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void mexhat(int N,double lb,double ub,double *psi,double *t) {
 | 
			
		||||
  gauss(N,2,lb,ub,psi,t);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void morlet(int N,double lb,double ub,double *psi,double *t) {
 | 
			
		||||
       int i;
 | 
			
		||||
       double delta;
 | 
			
		||||
 | 
			
		||||
  	if (lb >= ub) {
 | 
			
		||||
		printf("upper bound must be greater than lower bound");
 | 
			
		||||
		exit(1);
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	t[0] = lb;
 | 
			
		||||
	t[N-1] = ub;
 | 
			
		||||
	delta = (ub - lb) / (N-1);
 | 
			
		||||
	for(i = 1; i < N-1;++i) {
 | 
			
		||||
		t[i] = lb + delta * i;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for(i = 0; i < N;++i) {
 | 
			
		||||
	  psi[i] = exp(- t[i] * t[i] / 2.0) * cos(5 * t[i]);
 | 
			
		||||
	}
 | 
			
		||||
}
 | 
			
		||||
							
								
								
									
										24
									
								
								src/wavefunc.h
									
									
									
									
									
										Executable file
									
								
							
							
						
						
									
										24
									
								
								src/wavefunc.h
									
									
									
									
									
										Executable file
									
								
							@ -0,0 +1,24 @@
 | 
			
		||||
#ifndef WAVEFUNC_H_
 | 
			
		||||
#define WAVEFUNC_H_
 | 
			
		||||
 | 
			
		||||
#include "cwtmath.h"
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
extern "C" {
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
void meyer(int N,double lb,double ub,double *phi,double *psi,double *tgrid);
 | 
			
		||||
 | 
			
		||||
void gauss(int N,int p,double lb,double ub,double *psi,double *t);
 | 
			
		||||
 | 
			
		||||
void mexhat(int N,double lb,double ub,double *psi,double *t);
 | 
			
		||||
 | 
			
		||||
void morlet(int N,double lb,double ub,double *psi,double *t);  
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
}
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#endif /* WAVEFUNC_H_ */
 | 
			
		||||
							
								
								
									
										4913
									
								
								src/wavelib.c
									
									
									
									
									
								
							
							
						
						
									
										4913
									
								
								src/wavelib.c
									
									
									
									
									
								
							
										
											
												File diff suppressed because it is too large
												Load Diff
											
										
									
								
							
							
								
								
									
										402
									
								
								src/wavelib.h
									
									
									
									
									
								
							
							
						
						
									
										402
									
								
								src/wavelib.h
									
									
									
									
									
								
							@ -1,176 +1,230 @@
 | 
			
		||||
/*
 | 
			
		||||
Copyright (c) 2014, Rafat Hussain
 | 
			
		||||
*/
 | 
			
		||||
#ifndef WAVELIB_H_
 | 
			
		||||
#define WAVELIB_H_
 | 
			
		||||
 | 
			
		||||
#include "wtmath.h"
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
extern "C" {
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
#if defined(_MSC_VER)
 | 
			
		||||
#pragma warning(disable : 4200)
 | 
			
		||||
#pragma warning(disable : 4996)
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
typedef struct wave_set* wave_object;
 | 
			
		||||
 | 
			
		||||
wave_object wave_init(char* wname);
 | 
			
		||||
 | 
			
		||||
struct wave_set{
 | 
			
		||||
	char wname[50];
 | 
			
		||||
	int filtlength;// When all filters are of the same length. [Matlab uses zero-padding to make all filters of the same length]
 | 
			
		||||
	int lpd_len;// Default filtlength = lpd_len = lpr_len = hpd_len = hpr_len
 | 
			
		||||
	int hpd_len;
 | 
			
		||||
	int lpr_len;
 | 
			
		||||
	int hpr_len;
 | 
			
		||||
	double *lpd;
 | 
			
		||||
	double *hpd;
 | 
			
		||||
	double *lpr;
 | 
			
		||||
	double *hpr;
 | 
			
		||||
    double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
typedef struct wt_set* wt_object;
 | 
			
		||||
 | 
			
		||||
wt_object wt_init(wave_object wave,char* method, int siglength, int J);
 | 
			
		||||
 | 
			
		||||
struct wt_set{
 | 
			
		||||
	wave_object wave;
 | 
			
		||||
	conv_object cobj;
 | 
			
		||||
	char method[10];
 | 
			
		||||
	int siglength;// Length of the original signal.
 | 
			
		||||
	int outlength;// Length of the output DWT vector
 | 
			
		||||
	int lenlength;// Length of the Output Dimension Vector "length"
 | 
			
		||||
	int J; // Number of decomposition Levels
 | 
			
		||||
	int MaxIter;// Maximum Iterations J <= MaxIter
 | 
			
		||||
	int even;// even = 1 if signal is of even length. even = 0 otherwise
 | 
			
		||||
	char ext[10];// Type of Extension used - "per" or "sym"
 | 
			
		||||
	char cmethod[10]; // Convolution Method - "direct" or "FFT"
 | 
			
		||||
 | 
			
		||||
	int N; //
 | 
			
		||||
	int cfftset;
 | 
			
		||||
	int zpad;
 | 
			
		||||
	int length[102];
 | 
			
		||||
	double *output;
 | 
			
		||||
    double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
typedef struct wtree_set* wtree_object;
 | 
			
		||||
 | 
			
		||||
wtree_object wtree_init(wave_object wave, int siglength, int J);
 | 
			
		||||
 | 
			
		||||
struct wtree_set{
 | 
			
		||||
	wave_object wave;
 | 
			
		||||
	conv_object cobj;
 | 
			
		||||
	char method[10];
 | 
			
		||||
	int siglength;// Length of the original signal.
 | 
			
		||||
	int outlength;// Length of the output DWT vector
 | 
			
		||||
	int lenlength;// Length of the Output Dimension Vector "length"
 | 
			
		||||
	int J; // Number of decomposition Levels
 | 
			
		||||
	int MaxIter;// Maximum Iterations J <= MaxIter
 | 
			
		||||
	int even;// even = 1 if signal is of even length. even = 0 otherwise
 | 
			
		||||
	char ext[10];// Type of Extension used - "per" or "sym"
 | 
			
		||||
 | 
			
		||||
	int N; //
 | 
			
		||||
	int nodes;
 | 
			
		||||
	int cfftset;
 | 
			
		||||
	int zpad;
 | 
			
		||||
	int length[102];
 | 
			
		||||
	double *output;
 | 
			
		||||
	int *nodelength;
 | 
			
		||||
	int *coeflength;
 | 
			
		||||
    double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
typedef struct wpt_set* wpt_object;
 | 
			
		||||
 | 
			
		||||
wpt_object wpt_init(wave_object wave, int siglength, int J);
 | 
			
		||||
 | 
			
		||||
struct wpt_set{
 | 
			
		||||
	wave_object wave;
 | 
			
		||||
	conv_object cobj;
 | 
			
		||||
	int siglength;// Length of the original signal.
 | 
			
		||||
	int outlength;// Length of the output DWT vector
 | 
			
		||||
	int lenlength;// Length of the Output Dimension Vector "length"
 | 
			
		||||
	int J; // Number of decomposition Levels
 | 
			
		||||
	int MaxIter;// Maximum Iterations J <= MaxIter
 | 
			
		||||
	int even;// even = 1 if signal is of even length. even = 0 otherwise
 | 
			
		||||
	char ext[10];// Type of Extension used - "per" or "sym"
 | 
			
		||||
	char entropy[20];
 | 
			
		||||
	double eparam;
 | 
			
		||||
 | 
			
		||||
	int N; //
 | 
			
		||||
	int nodes;
 | 
			
		||||
	int length[102];
 | 
			
		||||
	double *output;
 | 
			
		||||
	double *costvalues;
 | 
			
		||||
	double *basisvector;
 | 
			
		||||
	int *nodeindex;
 | 
			
		||||
	int *numnodeslevel;
 | 
			
		||||
	int *coeflength;
 | 
			
		||||
    double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
void dwt(wt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void idwt(wt_object wt, double *dwtop);
 | 
			
		||||
 | 
			
		||||
void wtree(wtree_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void dwpt(wpt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void idwpt(wpt_object wt, double *dwtop);
 | 
			
		||||
 | 
			
		||||
void swt(wt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void iswt(wt_object wt, double *swtop);
 | 
			
		||||
 | 
			
		||||
void modwt(wt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void imodwt(wt_object wt, double *dwtop);
 | 
			
		||||
 | 
			
		||||
void setDWTExtension(wt_object wt, char *extension);
 | 
			
		||||
 | 
			
		||||
void setWTREEExtension(wtree_object wt, char *extension);
 | 
			
		||||
 | 
			
		||||
void setDWPTExtension(wpt_object wt, char *extension);
 | 
			
		||||
 | 
			
		||||
void setDWPTEntropy(wpt_object wt, char *entropy, double eparam);
 | 
			
		||||
 | 
			
		||||
void setWTConv(wt_object wt, char *cmethod);
 | 
			
		||||
 | 
			
		||||
int getWTREENodelength(wtree_object wt, int X);
 | 
			
		||||
 | 
			
		||||
void getWTREECoeffs(wtree_object wt, int X, int Y, double *coeffs, int N);
 | 
			
		||||
 | 
			
		||||
int getDWPTNodelength(wpt_object wt, int X);
 | 
			
		||||
 | 
			
		||||
void getDWPTCoeffs(wpt_object wt, int X, int Y, double *coeffs, int N);
 | 
			
		||||
 | 
			
		||||
void wave_summary(wave_object obj);
 | 
			
		||||
 | 
			
		||||
void wt_summary(wt_object wt);
 | 
			
		||||
 | 
			
		||||
void wtree_summary(wtree_object wt);
 | 
			
		||||
 | 
			
		||||
void wpt_summary(wpt_object wt);
 | 
			
		||||
 | 
			
		||||
void wave_free(wave_object object);
 | 
			
		||||
 | 
			
		||||
void wt_free(wt_object object);
 | 
			
		||||
 | 
			
		||||
void wtree_free(wtree_object object);
 | 
			
		||||
 | 
			
		||||
void wpt_free(wpt_object object);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
}
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#endif /* WAVELIB_H_ */
 | 
			
		||||
*/
 | 
			
		||||
#ifndef WAVELIB_H_
 | 
			
		||||
#define WAVELIB_H_
 | 
			
		||||
 | 
			
		||||
#include "wtmath.h"
 | 
			
		||||
#include "cwt.h"
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
extern "C" {
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
#if defined(_MSC_VER)
 | 
			
		||||
#pragma warning(disable : 4200)
 | 
			
		||||
#pragma warning(disable : 4996)
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
#ifndef cplx_type
 | 
			
		||||
#define cplx_type double
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
typedef struct cplx_t {
 | 
			
		||||
	cplx_type re;
 | 
			
		||||
	cplx_type im;
 | 
			
		||||
} cplx_data;
 | 
			
		||||
 | 
			
		||||
typedef struct wave_set* wave_object;
 | 
			
		||||
 | 
			
		||||
wave_object wave_init(char* wname);
 | 
			
		||||
 | 
			
		||||
struct wave_set{
 | 
			
		||||
	char wname[50];
 | 
			
		||||
	int filtlength;// When all filters are of the same length. [Matlab uses zero-padding to make all filters of the same length]
 | 
			
		||||
	int lpd_len;// Default filtlength = lpd_len = lpr_len = hpd_len = hpr_len
 | 
			
		||||
	int hpd_len;
 | 
			
		||||
	int lpr_len;
 | 
			
		||||
	int hpr_len;
 | 
			
		||||
	double *lpd;
 | 
			
		||||
	double *hpd;
 | 
			
		||||
	double *lpr;
 | 
			
		||||
	double *hpr;
 | 
			
		||||
	double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
typedef struct wt_set* wt_object;
 | 
			
		||||
 | 
			
		||||
wt_object wt_init(wave_object wave,char* method, int siglength, int J);
 | 
			
		||||
 | 
			
		||||
struct wt_set{
 | 
			
		||||
	wave_object wave;
 | 
			
		||||
	conv_object cobj;
 | 
			
		||||
	char method[10];
 | 
			
		||||
	int siglength;// Length of the original signal.
 | 
			
		||||
	int outlength;// Length of the output DWT vector
 | 
			
		||||
	int lenlength;// Length of the Output Dimension Vector "length"
 | 
			
		||||
	int J; // Number of decomposition Levels
 | 
			
		||||
	int MaxIter;// Maximum Iterations J <= MaxIter
 | 
			
		||||
	int even;// even = 1 if signal is of even length. even = 0 otherwise
 | 
			
		||||
	char ext[10];// Type of Extension used - "per" or "sym"
 | 
			
		||||
	char cmethod[10]; // Convolution Method - "direct" or "FFT"
 | 
			
		||||
 | 
			
		||||
	int N; //
 | 
			
		||||
	int cfftset;
 | 
			
		||||
	int zpad;
 | 
			
		||||
	int length[102];
 | 
			
		||||
	double *output;
 | 
			
		||||
	double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
typedef struct wtree_set* wtree_object;
 | 
			
		||||
 | 
			
		||||
wtree_object wtree_init(wave_object wave, int siglength, int J);
 | 
			
		||||
 | 
			
		||||
struct wtree_set{
 | 
			
		||||
	wave_object wave;
 | 
			
		||||
	conv_object cobj;
 | 
			
		||||
	char method[10];
 | 
			
		||||
	int siglength;// Length of the original signal.
 | 
			
		||||
	int outlength;// Length of the output DWT vector
 | 
			
		||||
	int lenlength;// Length of the Output Dimension Vector "length"
 | 
			
		||||
	int J; // Number of decomposition Levels
 | 
			
		||||
	int MaxIter;// Maximum Iterations J <= MaxIter
 | 
			
		||||
	int even;// even = 1 if signal is of even length. even = 0 otherwise
 | 
			
		||||
	char ext[10];// Type of Extension used - "per" or "sym"
 | 
			
		||||
 | 
			
		||||
	int N; //
 | 
			
		||||
	int nodes;
 | 
			
		||||
	int cfftset;
 | 
			
		||||
	int zpad;
 | 
			
		||||
	int length[102];
 | 
			
		||||
	double *output;
 | 
			
		||||
	int *nodelength;
 | 
			
		||||
	int *coeflength;
 | 
			
		||||
	double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
typedef struct wpt_set* wpt_object;
 | 
			
		||||
 | 
			
		||||
wpt_object wpt_init(wave_object wave, int siglength, int J);
 | 
			
		||||
 | 
			
		||||
struct wpt_set{
 | 
			
		||||
	wave_object wave;
 | 
			
		||||
	conv_object cobj;
 | 
			
		||||
	int siglength;// Length of the original signal.
 | 
			
		||||
	int outlength;// Length of the output DWT vector
 | 
			
		||||
	int lenlength;// Length of the Output Dimension Vector "length"
 | 
			
		||||
	int J; // Number of decomposition Levels
 | 
			
		||||
	int MaxIter;// Maximum Iterations J <= MaxIter
 | 
			
		||||
	int even;// even = 1 if signal is of even length. even = 0 otherwise
 | 
			
		||||
	char ext[10];// Type of Extension used - "per" or "sym"
 | 
			
		||||
	char entropy[20];
 | 
			
		||||
	double eparam;
 | 
			
		||||
 | 
			
		||||
	int N; //
 | 
			
		||||
	int nodes;
 | 
			
		||||
	int length[102];
 | 
			
		||||
	double *output;
 | 
			
		||||
	double *costvalues;
 | 
			
		||||
	double *basisvector;
 | 
			
		||||
	int *nodeindex;
 | 
			
		||||
	int *numnodeslevel;
 | 
			
		||||
	int *coeflength;
 | 
			
		||||
	double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
typedef struct cwt_set* cwt_object;
 | 
			
		||||
 | 
			
		||||
cwt_object cwt_init(char* wave, double param, int siglength,double dt, int J);
 | 
			
		||||
 | 
			
		||||
struct cwt_set{
 | 
			
		||||
	char wave[10];// Wavelet - morl/morlet,paul,dog/dgauss
 | 
			
		||||
	int siglength;// Length of Input Data
 | 
			
		||||
	int J;// Total Number of Scales
 | 
			
		||||
	double s0;// Smallest scale. It depends on the sampling rate. s0 <= 2 * dt for most wavelets
 | 
			
		||||
	double dt;// Sampling Rate
 | 
			
		||||
	double dj;// Separation between scales. eg., scale = s0 * 2 ^ ( [0:N-1] *dj ) or scale = s0 *[0:N-1] * dj
 | 
			
		||||
	char type[10];// Scale Type - Power or Linear
 | 
			
		||||
	int pow;// Base of Power in case type = pow. Typical value is pow = 2
 | 
			
		||||
	int sflag;
 | 
			
		||||
	int pflag;
 | 
			
		||||
	int npad;
 | 
			
		||||
	int mother;
 | 
			
		||||
	double m;// Wavelet parameter param
 | 
			
		||||
	double smean;// Input Signal mean
 | 
			
		||||
 | 
			
		||||
	cplx_data *output;
 | 
			
		||||
	double *scale;
 | 
			
		||||
	double *period;
 | 
			
		||||
	double *coi;
 | 
			
		||||
	double params[0];
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
void dwt(wt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void idwt(wt_object wt, double *dwtop);
 | 
			
		||||
 | 
			
		||||
void wtree(wtree_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void dwpt(wpt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void idwpt(wpt_object wt, double *dwtop);
 | 
			
		||||
 | 
			
		||||
void swt(wt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void iswt(wt_object wt, double *swtop);
 | 
			
		||||
 | 
			
		||||
void modwt(wt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void imodwt(wt_object wt, double *dwtop);
 | 
			
		||||
 | 
			
		||||
void setDWTExtension(wt_object wt, char *extension);
 | 
			
		||||
 | 
			
		||||
void setWTREEExtension(wtree_object wt, char *extension);
 | 
			
		||||
 | 
			
		||||
void setDWPTExtension(wpt_object wt, char *extension);
 | 
			
		||||
 | 
			
		||||
void setDWPTEntropy(wpt_object wt, char *entropy, double eparam);
 | 
			
		||||
 | 
			
		||||
void setWTConv(wt_object wt, char *cmethod);
 | 
			
		||||
 | 
			
		||||
int getWTREENodelength(wtree_object wt, int X);
 | 
			
		||||
 | 
			
		||||
void getWTREECoeffs(wtree_object wt, int X, int Y, double *coeffs, int N);
 | 
			
		||||
 | 
			
		||||
int getDWPTNodelength(wpt_object wt, int X);
 | 
			
		||||
 | 
			
		||||
void getDWPTCoeffs(wpt_object wt, int X, int Y, double *coeffs, int N);
 | 
			
		||||
 | 
			
		||||
void setCWTScales(cwt_object wt, double s0, double dj, char *type, int power);
 | 
			
		||||
 | 
			
		||||
void setCWTScaleVector(cwt_object wt, double *scale, int J, double s0, double dj);
 | 
			
		||||
 | 
			
		||||
void setCWTPadding(cwt_object wt, int pad);
 | 
			
		||||
 | 
			
		||||
void cwt(cwt_object wt, double *inp);
 | 
			
		||||
 | 
			
		||||
void icwt(cwt_object wt, double *cwtop);
 | 
			
		||||
 | 
			
		||||
int getCWTScaleLength(int N);
 | 
			
		||||
 | 
			
		||||
void wave_summary(wave_object obj);
 | 
			
		||||
 | 
			
		||||
void wt_summary(wt_object wt);
 | 
			
		||||
 | 
			
		||||
void wtree_summary(wtree_object wt);
 | 
			
		||||
 | 
			
		||||
void wpt_summary(wpt_object wt);
 | 
			
		||||
 | 
			
		||||
void cwt_summary(cwt_object wt);
 | 
			
		||||
 | 
			
		||||
void wave_free(wave_object object);
 | 
			
		||||
 | 
			
		||||
void wt_free(wt_object object);
 | 
			
		||||
 | 
			
		||||
void wtree_free(wtree_object object);
 | 
			
		||||
 | 
			
		||||
void wpt_free(wpt_object object);
 | 
			
		||||
 | 
			
		||||
void cwt_free(cwt_object object);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#ifdef __cplusplus
 | 
			
		||||
}
 | 
			
		||||
#endif
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
#endif /* WAVELIB_H_ */
 | 
			
		||||
 | 
			
		||||
							
								
								
									
										109
									
								
								test/cwttest.c
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										109
									
								
								test/cwttest.c
									
									
									
									
									
										Normal file
									
								
							@ -0,0 +1,109 @@
 | 
			
		||||
#include <stdio.h>
 | 
			
		||||
#include <stdlib.h>
 | 
			
		||||
#include <string.h>
 | 
			
		||||
#include <math.h>
 | 
			
		||||
#include "../header/wavelib.h"
 | 
			
		||||
 | 
			
		||||
int main() {
 | 
			
		||||
	int i, N, J,subscale,a0,iter,nd,k;
 | 
			
		||||
	double *inp,*oup;
 | 
			
		||||
	double dt, dj,s0, param,mn;
 | 
			
		||||
	double td,tn,den, num, recon_mean, recon_var;
 | 
			
		||||
	cwt_object wt;
 | 
			
		||||
 | 
			
		||||
	FILE *ifp;
 | 
			
		||||
	double temp[1200];
 | 
			
		||||
 | 
			
		||||
	char *wave = "morlet";// Set Morlet wavelet. Other options "paul" and "dog"
 | 
			
		||||
	char *type = "pow";
 | 
			
		||||
 | 
			
		||||
	N = 504;
 | 
			
		||||
	param = 6.0;
 | 
			
		||||
	subscale = 4;
 | 
			
		||||
	dt = 0.25;
 | 
			
		||||
	s0 = dt;
 | 
			
		||||
	dj = 1.0 / (double)subscale;
 | 
			
		||||
	J = 11 * subscale; // Total Number of scales
 | 
			
		||||
	a0 = 2;//power
 | 
			
		||||
 | 
			
		||||
	ifp = fopen("sst_nino3.dat", "r");
 | 
			
		||||
	i = 0;
 | 
			
		||||
	if (!ifp) {
 | 
			
		||||
		printf("Cannot Open File");
 | 
			
		||||
		exit(100);
 | 
			
		||||
	}
 | 
			
		||||
	while (!feof(ifp)) {
 | 
			
		||||
		fscanf(ifp, "%lf \n", &temp[i]);
 | 
			
		||||
		i++;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	fclose(ifp);
 | 
			
		||||
 | 
			
		||||
	wt = cwt_init(wave, param, N,dt, J);
 | 
			
		||||
 | 
			
		||||
	inp = (double*)malloc(sizeof(double)* N);
 | 
			
		||||
	oup = (double*)malloc(sizeof(double)* N);
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		inp[i] = temp[i] ;
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	setCWTScales(wt, s0, dj, type, a0);
 | 
			
		||||
 | 
			
		||||
	cwt(wt, inp);
 | 
			
		||||
 | 
			
		||||
	printf("\n MEAN %g \n", wt->smean);
 | 
			
		||||
	
 | 
			
		||||
	mn = 0.0;
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		mn += sqrt(wt->output[i].re * wt->output[i].re + wt->output[i].im * wt->output[i].im);
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	cwt_summary(wt);
 | 
			
		||||
 | 
			
		||||
	printf("\n abs mean %g \n", mn / N);
 | 
			
		||||
	
 | 
			
		||||
	printf("\n\n");
 | 
			
		||||
	printf("Let CWT w = w(j, n/2 - 1) where n = %d\n\n", N);
 | 
			
		||||
	nd = N/2 - 1;
 | 
			
		||||
	
 | 
			
		||||
	printf("%-15s%-15s%-15s%-15s \n","j","Scale","Period","ABS(w)^2");
 | 
			
		||||
	for(k = 0; k < wt->J;++k) {
 | 
			
		||||
		iter = nd + k * N;
 | 
			
		||||
		printf("%-15d%-15lf%-15lf%-15lf \n",k,wt->scale[k],wt->period[k],
 | 
			
		||||
		wt->output[iter].re * wt->output[iter].re + wt->output[iter].im * wt->output[iter].im);
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	icwt(wt, oup);
 | 
			
		||||
 | 
			
		||||
	num = den = recon_var = recon_mean = 0.0;
 | 
			
		||||
	printf("\n\n");
 | 
			
		||||
	printf("Signal Reconstruction\n");
 | 
			
		||||
	printf("%-15s%-15s%-15s \n","i","Input(i)","Output(i)");
 | 
			
		||||
 | 
			
		||||
	for (i = N - 10; i < N; ++i) {
 | 
			
		||||
		printf("%-15d%-15lf%-15lf \n", i,inp[i] , oup[i]);
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		//printf("%g %g \n", oup[i] ,inp[i] - wt->smean);
 | 
			
		||||
		td = inp[i] ;
 | 
			
		||||
		tn = oup[i] - td;
 | 
			
		||||
		num += (tn * tn);
 | 
			
		||||
		den += (td * td);
 | 
			
		||||
		recon_mean += oup[i];
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	recon_var = sqrt(num / N);
 | 
			
		||||
	recon_mean /= N;
 | 
			
		||||
 | 
			
		||||
	printf("\nRMS Error %g \n", sqrt(num) / sqrt(den));
 | 
			
		||||
	printf("\nVariance %g \n", recon_var);
 | 
			
		||||
	printf("\nMean %g \n", recon_mean);
 | 
			
		||||
 | 
			
		||||
	free(inp);
 | 
			
		||||
	free(oup);
 | 
			
		||||
	cwt_free(wt);
 | 
			
		||||
	return 0;
 | 
			
		||||
}
 | 
			
		||||
							
								
								
									
										504
									
								
								test/sst_nino3.dat
									
									
									
									
									
										Executable file
									
								
							
							
						
						
									
										504
									
								
								test/sst_nino3.dat
									
									
									
									
									
										Executable file
									
								
							@ -0,0 +1,504 @@
 | 
			
		||||
-0.15
 | 
			
		||||
-0.30
 | 
			
		||||
-0.14
 | 
			
		||||
-0.41
 | 
			
		||||
-0.46
 | 
			
		||||
-0.66
 | 
			
		||||
-0.50
 | 
			
		||||
-0.80
 | 
			
		||||
-0.95
 | 
			
		||||
-0.72
 | 
			
		||||
-0.31
 | 
			
		||||
-0.71
 | 
			
		||||
-1.04
 | 
			
		||||
-0.77
 | 
			
		||||
-0.86
 | 
			
		||||
-0.84
 | 
			
		||||
-0.41
 | 
			
		||||
-0.49
 | 
			
		||||
-0.48
 | 
			
		||||
-0.72
 | 
			
		||||
-1.21
 | 
			
		||||
-0.80
 | 
			
		||||
 0.16
 | 
			
		||||
 0.46
 | 
			
		||||
 0.40
 | 
			
		||||
 1.00
 | 
			
		||||
 2.17
 | 
			
		||||
 2.50
 | 
			
		||||
 2.34
 | 
			
		||||
 0.80
 | 
			
		||||
 0.14
 | 
			
		||||
-0.06
 | 
			
		||||
-0.34
 | 
			
		||||
-0.71
 | 
			
		||||
-0.34
 | 
			
		||||
-0.73
 | 
			
		||||
-0.48
 | 
			
		||||
-0.11
 | 
			
		||||
 0.22
 | 
			
		||||
 0.51
 | 
			
		||||
 0.51
 | 
			
		||||
 0.25
 | 
			
		||||
-0.10
 | 
			
		||||
-0.33
 | 
			
		||||
-0.42
 | 
			
		||||
-0.23
 | 
			
		||||
-0.53
 | 
			
		||||
-0.44
 | 
			
		||||
-0.30
 | 
			
		||||
 0.15
 | 
			
		||||
 0.09
 | 
			
		||||
 0.19
 | 
			
		||||
-0.06
 | 
			
		||||
 0.25
 | 
			
		||||
 0.30
 | 
			
		||||
 0.81
 | 
			
		||||
 0.26
 | 
			
		||||
 0.10
 | 
			
		||||
 0.34
 | 
			
		||||
 1.01
 | 
			
		||||
-0.31
 | 
			
		||||
-0.90
 | 
			
		||||
-0.73
 | 
			
		||||
-0.92
 | 
			
		||||
-0.73
 | 
			
		||||
-0.31
 | 
			
		||||
-0.03
 | 
			
		||||
 0.12
 | 
			
		||||
 0.37
 | 
			
		||||
 0.82
 | 
			
		||||
 1.22
 | 
			
		||||
 1.83
 | 
			
		||||
 1.60
 | 
			
		||||
 0.34
 | 
			
		||||
-0.72
 | 
			
		||||
-0.87
 | 
			
		||||
-0.85
 | 
			
		||||
-0.40
 | 
			
		||||
-0.39
 | 
			
		||||
-0.65
 | 
			
		||||
 0.07
 | 
			
		||||
 0.67
 | 
			
		||||
 0.39
 | 
			
		||||
 0.03
 | 
			
		||||
-0.17
 | 
			
		||||
-0.76
 | 
			
		||||
-0.87
 | 
			
		||||
-1.36
 | 
			
		||||
-1.10
 | 
			
		||||
-0.99
 | 
			
		||||
-0.78
 | 
			
		||||
-0.93
 | 
			
		||||
-0.87
 | 
			
		||||
-0.44
 | 
			
		||||
-0.34
 | 
			
		||||
-0.50
 | 
			
		||||
-0.39
 | 
			
		||||
-0.04
 | 
			
		||||
 0.42
 | 
			
		||||
 0.62
 | 
			
		||||
 0.17
 | 
			
		||||
 0.23
 | 
			
		||||
 1.03
 | 
			
		||||
 1.54
 | 
			
		||||
 1.09
 | 
			
		||||
 0.01
 | 
			
		||||
 0.12
 | 
			
		||||
-0.27
 | 
			
		||||
-0.47
 | 
			
		||||
-0.41
 | 
			
		||||
-0.37
 | 
			
		||||
-0.36
 | 
			
		||||
-0.39
 | 
			
		||||
 0.43
 | 
			
		||||
 1.05
 | 
			
		||||
 1.58
 | 
			
		||||
 1.25
 | 
			
		||||
 0.86
 | 
			
		||||
 0.60
 | 
			
		||||
 0.21
 | 
			
		||||
 0.19
 | 
			
		||||
-0.23
 | 
			
		||||
-0.29
 | 
			
		||||
 0.18
 | 
			
		||||
 0.12
 | 
			
		||||
 0.71
 | 
			
		||||
 1.42
 | 
			
		||||
 1.59
 | 
			
		||||
 0.93
 | 
			
		||||
-0.25
 | 
			
		||||
-0.66
 | 
			
		||||
-0.95
 | 
			
		||||
-0.47
 | 
			
		||||
 0.06
 | 
			
		||||
 0.70
 | 
			
		||||
 0.81
 | 
			
		||||
 0.78
 | 
			
		||||
 1.43
 | 
			
		||||
 1.22
 | 
			
		||||
 1.05
 | 
			
		||||
 0.44
 | 
			
		||||
-0.35
 | 
			
		||||
-0.67
 | 
			
		||||
-0.84
 | 
			
		||||
-0.66
 | 
			
		||||
-0.45
 | 
			
		||||
-0.12
 | 
			
		||||
-0.20
 | 
			
		||||
-0.16
 | 
			
		||||
-0.47
 | 
			
		||||
-0.52
 | 
			
		||||
-0.79
 | 
			
		||||
-0.80
 | 
			
		||||
-0.62
 | 
			
		||||
-0.86
 | 
			
		||||
-1.29
 | 
			
		||||
-1.04
 | 
			
		||||
-1.05
 | 
			
		||||
-0.75
 | 
			
		||||
-0.81
 | 
			
		||||
-0.90
 | 
			
		||||
-0.25
 | 
			
		||||
 0.62
 | 
			
		||||
 1.22
 | 
			
		||||
 0.96
 | 
			
		||||
 0.21
 | 
			
		||||
-0.11
 | 
			
		||||
-0.25
 | 
			
		||||
-0.24
 | 
			
		||||
-0.43
 | 
			
		||||
 0.23
 | 
			
		||||
 0.67
 | 
			
		||||
 0.78
 | 
			
		||||
 0.41
 | 
			
		||||
 0.98
 | 
			
		||||
 1.28
 | 
			
		||||
 1.45
 | 
			
		||||
 1.02
 | 
			
		||||
 0.03
 | 
			
		||||
-0.59
 | 
			
		||||
-1.34
 | 
			
		||||
-0.99
 | 
			
		||||
-1.49
 | 
			
		||||
-1.74
 | 
			
		||||
-1.33
 | 
			
		||||
-0.55
 | 
			
		||||
-0.51
 | 
			
		||||
-0.36
 | 
			
		||||
-0.99
 | 
			
		||||
 0.32
 | 
			
		||||
 1.04
 | 
			
		||||
 1.41
 | 
			
		||||
 0.99
 | 
			
		||||
 0.66
 | 
			
		||||
 0.50
 | 
			
		||||
 0.22
 | 
			
		||||
 0.71
 | 
			
		||||
-0.16
 | 
			
		||||
 0.38
 | 
			
		||||
 0.00
 | 
			
		||||
-1.11
 | 
			
		||||
-1.04
 | 
			
		||||
 0.05
 | 
			
		||||
-0.64
 | 
			
		||||
-0.34
 | 
			
		||||
-0.50
 | 
			
		||||
-1.85
 | 
			
		||||
-0.94
 | 
			
		||||
-0.78
 | 
			
		||||
 0.29
 | 
			
		||||
 0.27
 | 
			
		||||
 0.69
 | 
			
		||||
-0.06
 | 
			
		||||
-0.83
 | 
			
		||||
-0.80
 | 
			
		||||
-1.02
 | 
			
		||||
-0.96
 | 
			
		||||
-0.09
 | 
			
		||||
 0.62
 | 
			
		||||
 0.87
 | 
			
		||||
 1.03
 | 
			
		||||
 0.70
 | 
			
		||||
-0.10
 | 
			
		||||
-0.31
 | 
			
		||||
 0.04
 | 
			
		||||
-0.46
 | 
			
		||||
 0.04
 | 
			
		||||
 0.24
 | 
			
		||||
-0.08
 | 
			
		||||
-0.28
 | 
			
		||||
 0.06
 | 
			
		||||
 0.05
 | 
			
		||||
-0.31
 | 
			
		||||
 0.11
 | 
			
		||||
 0.27
 | 
			
		||||
 0.26
 | 
			
		||||
 0.04
 | 
			
		||||
 0.12
 | 
			
		||||
 1.11
 | 
			
		||||
 1.53
 | 
			
		||||
 1.23
 | 
			
		||||
 0.17
 | 
			
		||||
-0.18
 | 
			
		||||
-0.56
 | 
			
		||||
 0.05
 | 
			
		||||
 0.41
 | 
			
		||||
 0.22
 | 
			
		||||
 0.04
 | 
			
		||||
-0.19
 | 
			
		||||
-0.46
 | 
			
		||||
-0.65
 | 
			
		||||
-1.06
 | 
			
		||||
-0.54
 | 
			
		||||
 0.14
 | 
			
		||||
 0.25
 | 
			
		||||
-0.21
 | 
			
		||||
-0.73
 | 
			
		||||
-0.43
 | 
			
		||||
 0.48
 | 
			
		||||
 0.26
 | 
			
		||||
 0.05
 | 
			
		||||
 0.11
 | 
			
		||||
-0.27
 | 
			
		||||
-0.08
 | 
			
		||||
-0.10
 | 
			
		||||
 0.29
 | 
			
		||||
-0.15
 | 
			
		||||
-0.28
 | 
			
		||||
-0.55
 | 
			
		||||
-0.44
 | 
			
		||||
-1.40
 | 
			
		||||
-0.55
 | 
			
		||||
-0.69
 | 
			
		||||
 0.58
 | 
			
		||||
 0.37
 | 
			
		||||
 0.42
 | 
			
		||||
 1.83
 | 
			
		||||
 1.23
 | 
			
		||||
 0.65
 | 
			
		||||
 0.41
 | 
			
		||||
 1.03
 | 
			
		||||
 0.64
 | 
			
		||||
-0.07
 | 
			
		||||
 0.98
 | 
			
		||||
 0.36
 | 
			
		||||
-0.30
 | 
			
		||||
-1.33
 | 
			
		||||
-1.39
 | 
			
		||||
-0.94
 | 
			
		||||
 0.34
 | 
			
		||||
-0.00
 | 
			
		||||
-0.15
 | 
			
		||||
 0.06
 | 
			
		||||
 0.39
 | 
			
		||||
 0.36
 | 
			
		||||
-0.49
 | 
			
		||||
-0.53
 | 
			
		||||
 0.35
 | 
			
		||||
 0.07
 | 
			
		||||
-0.24
 | 
			
		||||
 0.20
 | 
			
		||||
-0.22
 | 
			
		||||
-0.68
 | 
			
		||||
-0.44
 | 
			
		||||
 0.02
 | 
			
		||||
-0.22
 | 
			
		||||
-0.30
 | 
			
		||||
-0.59
 | 
			
		||||
 0.10
 | 
			
		||||
-0.02
 | 
			
		||||
-0.27
 | 
			
		||||
-0.60
 | 
			
		||||
-0.48
 | 
			
		||||
-0.37
 | 
			
		||||
-0.53
 | 
			
		||||
-1.35
 | 
			
		||||
-1.22
 | 
			
		||||
-0.99
 | 
			
		||||
-0.34
 | 
			
		||||
-0.79
 | 
			
		||||
-0.24
 | 
			
		||||
 0.02
 | 
			
		||||
 0.69
 | 
			
		||||
 0.78
 | 
			
		||||
 0.17
 | 
			
		||||
-0.17
 | 
			
		||||
-0.29
 | 
			
		||||
-0.27
 | 
			
		||||
 0.31
 | 
			
		||||
 0.44
 | 
			
		||||
 0.38
 | 
			
		||||
 0.24
 | 
			
		||||
-0.13
 | 
			
		||||
-0.89
 | 
			
		||||
-0.76
 | 
			
		||||
-0.71
 | 
			
		||||
-0.37
 | 
			
		||||
-0.59
 | 
			
		||||
-0.63
 | 
			
		||||
-1.47
 | 
			
		||||
-0.40
 | 
			
		||||
-0.18
 | 
			
		||||
-0.37
 | 
			
		||||
-0.43
 | 
			
		||||
-0.06
 | 
			
		||||
 0.61
 | 
			
		||||
 1.33
 | 
			
		||||
 1.19
 | 
			
		||||
 1.13
 | 
			
		||||
 0.31
 | 
			
		||||
 0.14
 | 
			
		||||
 0.03
 | 
			
		||||
 0.21
 | 
			
		||||
 0.15
 | 
			
		||||
-0.22
 | 
			
		||||
-0.02
 | 
			
		||||
 0.03
 | 
			
		||||
-0.17
 | 
			
		||||
 0.12
 | 
			
		||||
-0.35
 | 
			
		||||
-0.06
 | 
			
		||||
 0.38
 | 
			
		||||
-0.45
 | 
			
		||||
-0.32
 | 
			
		||||
-0.33
 | 
			
		||||
-0.49
 | 
			
		||||
-0.14
 | 
			
		||||
-0.56
 | 
			
		||||
-0.18
 | 
			
		||||
 0.46
 | 
			
		||||
 1.09
 | 
			
		||||
 1.04
 | 
			
		||||
 0.23
 | 
			
		||||
-0.99
 | 
			
		||||
-0.59
 | 
			
		||||
-0.92
 | 
			
		||||
-0.28
 | 
			
		||||
 0.52
 | 
			
		||||
 1.31
 | 
			
		||||
 1.45
 | 
			
		||||
 0.61
 | 
			
		||||
-0.11
 | 
			
		||||
-0.18
 | 
			
		||||
-0.39
 | 
			
		||||
-0.39
 | 
			
		||||
-0.36
 | 
			
		||||
-0.50
 | 
			
		||||
-0.81
 | 
			
		||||
-1.10
 | 
			
		||||
-0.29
 | 
			
		||||
 0.57
 | 
			
		||||
 0.68
 | 
			
		||||
 0.78
 | 
			
		||||
 0.78
 | 
			
		||||
 0.63
 | 
			
		||||
 0.98
 | 
			
		||||
 0.49
 | 
			
		||||
-0.42
 | 
			
		||||
-1.34
 | 
			
		||||
-1.20
 | 
			
		||||
-1.18
 | 
			
		||||
-0.65
 | 
			
		||||
-0.42
 | 
			
		||||
-0.97
 | 
			
		||||
-0.28
 | 
			
		||||
 0.77
 | 
			
		||||
 1.77
 | 
			
		||||
 2.22
 | 
			
		||||
 1.05
 | 
			
		||||
-0.67
 | 
			
		||||
-0.99
 | 
			
		||||
-1.52
 | 
			
		||||
-1.17
 | 
			
		||||
-0.22
 | 
			
		||||
-0.04
 | 
			
		||||
-0.45
 | 
			
		||||
-0.46
 | 
			
		||||
-0.75
 | 
			
		||||
-0.70
 | 
			
		||||
-1.38
 | 
			
		||||
-1.15
 | 
			
		||||
-0.01
 | 
			
		||||
 0.97
 | 
			
		||||
 1.10
 | 
			
		||||
 0.68
 | 
			
		||||
-0.02
 | 
			
		||||
-0.04
 | 
			
		||||
 0.47
 | 
			
		||||
 0.30
 | 
			
		||||
-0.55
 | 
			
		||||
-0.51
 | 
			
		||||
-0.09
 | 
			
		||||
-0.01
 | 
			
		||||
 0.34
 | 
			
		||||
 0.61
 | 
			
		||||
 0.58
 | 
			
		||||
 0.33
 | 
			
		||||
 0.38
 | 
			
		||||
 0.10
 | 
			
		||||
 0.18
 | 
			
		||||
-0.30
 | 
			
		||||
-0.06
 | 
			
		||||
-0.28
 | 
			
		||||
 0.12
 | 
			
		||||
 0.58
 | 
			
		||||
 0.89
 | 
			
		||||
 0.93
 | 
			
		||||
 2.39
 | 
			
		||||
 2.44
 | 
			
		||||
 1.92
 | 
			
		||||
 0.64
 | 
			
		||||
-0.24
 | 
			
		||||
 0.27
 | 
			
		||||
-0.13
 | 
			
		||||
-0.16
 | 
			
		||||
-0.54
 | 
			
		||||
-0.13
 | 
			
		||||
-0.37
 | 
			
		||||
-0.78
 | 
			
		||||
-0.22
 | 
			
		||||
 0.03
 | 
			
		||||
 0.25
 | 
			
		||||
 0.31
 | 
			
		||||
 1.03
 | 
			
		||||
 1.10
 | 
			
		||||
 1.05
 | 
			
		||||
 1.11
 | 
			
		||||
 1.28
 | 
			
		||||
 0.57
 | 
			
		||||
-0.55
 | 
			
		||||
-1.16
 | 
			
		||||
-0.99
 | 
			
		||||
-0.38
 | 
			
		||||
 0.01
 | 
			
		||||
-0.29
 | 
			
		||||
 0.09
 | 
			
		||||
 0.46
 | 
			
		||||
 0.57
 | 
			
		||||
 0.24
 | 
			
		||||
 0.39
 | 
			
		||||
 0.49
 | 
			
		||||
 0.86
 | 
			
		||||
 0.51
 | 
			
		||||
 0.95
 | 
			
		||||
 1.25
 | 
			
		||||
 1.33
 | 
			
		||||
-0.00
 | 
			
		||||
 0.34
 | 
			
		||||
 0.66
 | 
			
		||||
 1.11
 | 
			
		||||
 0.34
 | 
			
		||||
 0.48
 | 
			
		||||
 0.56
 | 
			
		||||
 0.39
 | 
			
		||||
-0.17
 | 
			
		||||
 1.04
 | 
			
		||||
 0.77
 | 
			
		||||
 0.12
 | 
			
		||||
-0.35
 | 
			
		||||
-0.22
 | 
			
		||||
 0.08
 | 
			
		||||
-0.08
 | 
			
		||||
-0.18
 | 
			
		||||
-0.06
 | 
			
		||||
@ -98,6 +98,17 @@ double RMS_Error(double *data, double *rec, int N) {
 | 
			
		||||
    return sqrt(sum/((double)N-1));
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
double REL_Error(double *data, double *rec, int N) {
 | 
			
		||||
    int i;
 | 
			
		||||
    double sum1 = 0;
 | 
			
		||||
    double sum2 = 0;
 | 
			
		||||
    for (i = 0; i < N; ++i) {
 | 
			
		||||
        sum1 += (data[i] - rec[i])*(data[i] - rec[i]);
 | 
			
		||||
        sum2 += data[i] * data[i];
 | 
			
		||||
    }
 | 
			
		||||
    return sqrt(sum1)/sqrt(sum2);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void ReconstructionTest()
 | 
			
		||||
{
 | 
			
		||||
 | 
			
		||||
@ -334,6 +345,77 @@ void DWPTReconstructionTest()
 | 
			
		||||
    free(inp);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
void CWTReconstructionTest() {
 | 
			
		||||
	int i, j,N, J,subscale,a0,iter;
 | 
			
		||||
	double *inp,*oup;
 | 
			
		||||
	double dt, dj,s0, pi,t;
 | 
			
		||||
	double val, epsilon;
 | 
			
		||||
	int it1,it2;
 | 
			
		||||
	cwt_object wt;
 | 
			
		||||
	
 | 
			
		||||
 | 
			
		||||
	char *wave[3];
 | 
			
		||||
	wave[0] = (char*) "morl";
 | 
			
		||||
	wave[1] =(char*) "paul";
 | 
			
		||||
	wave[2] = (char*) "dog";
 | 
			
		||||
	double param[30] = {4.5,5,5.5,6,6.5,8,10,13,17,20,
 | 
			
		||||
		4,5,7,8,10,12,13,14,17,20,2,4,6,8,10,12,14,16,18,20};
 | 
			
		||||
	char *type = (char*) "pow";
 | 
			
		||||
	
 | 
			
		||||
	epsilon = 0.01;
 | 
			
		||||
	N = 2048;
 | 
			
		||||
	dt = 0.000125;
 | 
			
		||||
	subscale = 20;
 | 
			
		||||
	dj = 1.0 / (double) subscale;
 | 
			
		||||
	s0 = dt/32;	
 | 
			
		||||
	J = 32 * subscale;
 | 
			
		||||
	a0 = 2;//power
 | 
			
		||||
	
 | 
			
		||||
 | 
			
		||||
	inp = (double*)malloc(sizeof(double)* N);
 | 
			
		||||
	oup = (double*)malloc(sizeof(double)* N);
 | 
			
		||||
	
 | 
			
		||||
	pi = 4.0 * atan(1.0);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
	for (i = 0; i < N; ++i) {
 | 
			
		||||
		t = dt * i;
 | 
			
		||||
		inp[i] = sin(2 * pi * 500 * t) + sin(2 * pi * 1000 * t) + 0.1 * sin(2 * pi * 8 * t);
 | 
			
		||||
		if (i == 1200 || i ==1232) {
 | 
			
		||||
			inp[i] += 5.0;
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
	for(it1 = 0; it1 < 3;++it1) {
 | 
			
		||||
		for(it2 = 0; it2 < 10;++it2) {
 | 
			
		||||
	
 | 
			
		||||
			wt = cwt_init(wave[it1], param[it1*10+it2], N,dt, J);
 | 
			
		||||
 | 
			
		||||
			setCWTScales(wt, s0, dj, type, a0);
 | 
			
		||||
 | 
			
		||||
			cwt(wt, inp);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
			icwt(wt, oup);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
			//printf("\nWavelet : %s Parameter %g Error %g \n", wave[it1],param[it1*10+it2],REL_Error(inp,oup, wt->siglength));
 | 
			
		||||
			if (REL_Error(inp,oup, wt->siglength) > epsilon) {
 | 
			
		||||
				printf("\n ERROR : DWPT Reconstruction Unit Test Failed. Exiting. \n");
 | 
			
		||||
				exit(-1);
 | 
			
		||||
			}
 | 
			
		||||
	
 | 
			
		||||
			cwt_free(wt);
 | 
			
		||||
		}
 | 
			
		||||
	}
 | 
			
		||||
	
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
	free(inp);
 | 
			
		||||
	free(oup);
 | 
			
		||||
	
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
void DBCoefTests()
 | 
			
		||||
{
 | 
			
		||||
@ -569,6 +651,9 @@ int main() {
 | 
			
		||||
	printf("Running DWPT ReconstructionTests ... ");
 | 
			
		||||
	DWPTReconstructionTest();
 | 
			
		||||
	printf("DONE \n");
 | 
			
		||||
	printf("Running CWT ReconstructionTests ... ");
 | 
			
		||||
	CWTReconstructionTest();
 | 
			
		||||
	printf("DONE \n");
 | 
			
		||||
	printf("\n\nUnit Tests Successful\n\n");
 | 
			
		||||
	return 0;
 | 
			
		||||
	
 | 
			
		||||
 | 
			
		||||
							
								
								
									
										
											BIN
										
									
								
								wavelib-doc.pdf
									
									
									
									
									
								
							
							
						
						
									
										
											BIN
										
									
								
								wavelib-doc.pdf
									
									
									
									
									
								
							
										
											Binary file not shown.
										
									
								
							
		Loading…
	
	
			
			x
			
			
		
	
		Reference in New Issue
	
	Block a user