wavecatcher-analysis
Loading...
Searching...
No Matches
Filters.cc
Go to the documentation of this file.
1#include "Filters.h"
2
12void Filters::Convolute(double*& result, double* first, double* second, int size) {
13
14 // FFT real -> im
15 auto fft_re_im = [&size](double* orig, double*& re, double*& im) {
16 TVirtualFFT* fft = TVirtualFFT::FFT(1, &size, "R2C ES");
17 fft->SetPoints(orig);
18 fft->Transform();
19 fft->GetPointsComplex(re, im);
20 delete fft;
21 };
22
23 double* refirst = new double[size];
24 double* imfirst = new double[size];
25 double* resecond = new double[size];
26 double* imsecond = new double[size];
27 double* reres = new double[size];
28 double* imres = new double[size];
29
30 fft_re_im(first, refirst, imfirst);
31 fft_re_im(second, resecond, imsecond);
32
33 TComplex cofirst;
34 TComplex cosecond;
35 TComplex cores;
36 double i_size = 1. / static_cast<double>(max(1, size));
37
38 for (int i = 0; i < size; i++) {
39 cofirst(refirst[i], imfirst[i]);
40 cosecond(resecond[i], imsecond[i]);
41
42 cores = cofirst * cosecond * i_size;
43
44 reres[i] = cores.Re();
45 imres[i] = cores.Im();
46 }
47
48
49 //cout << "performing IFFT ... ";
50 TVirtualFFT* fft_back = TVirtualFFT::FFT(1, &size, "C2R ES");
51 fft_back->SetPointsComplex(reres, imres);
52 fft_back->Transform();
53 fft_back->GetPoints(result);
54 delete fft_back;
55 delete[] imres; delete[] reres; delete[] refirst; delete[] imfirst; delete[] resecond; delete[] imsecond;
56}
57
66void Filters::Convolute(double*& first, double* second, int size) {
67 Convolute(first, first, second, size);
68}
69
81void Filters::Deconvolute(double*& result, double* first, double* second, int size, double sigma, double x0, double bin_size) {
82
83 // FFT real -> im
84 auto fft_re_im = [&size](double* orig, double*& re, double*& im) {
85 TVirtualFFT* fft = TVirtualFFT::FFT(1, &size, "R2C ES");
86 fft->SetPoints(orig);
87 fft->Transform();
88 fft->GetPointsComplex(re, im);
89 delete fft;
90 };
91
92 double* refirst = new double[size];
93 double* imfirst = new double[size];
94 double* resecond = new double[size];
95 double* imsecond = new double[size];
96 double* regauss = new double[size];
97 double* imgauss = new double[size];
98 double* reres = new double[size];
99 double* imres = new double[size];
100 double* gauss = new double[size];
101
102 fft_re_im(first, refirst, imfirst);
103 fft_re_im(second, resecond, imsecond);
104 if (sigma > 0.) {
105 for (int i = 0; i < size; i++) {
106 gauss[i] = TMath::Gaus(static_cast<double>(i) * bin_size, x0, sigma);
107 }
108 fft_re_im(gauss, regauss, imgauss);
109 }
110
111 TComplex cofirst;
112 TComplex cosecond;
113 TComplex cogauss;
114 TComplex cores;
115 double i_size = 1. / static_cast<double>(max(1, size));
116
117 for (int i = 0; i < size; i++) {
118 cofirst(refirst[i], imfirst[i]);
119 cosecond(resecond[i], imsecond[i]);
120 if (sigma > 0.) cogauss(regauss[i], imgauss[i]);
121 else cogauss(1., 0.);
122
123 cores = cofirst * cogauss / cosecond * i_size;
124
125 reres[i] = cores.Re();
126 imres[i] = cores.Im();
127 }
128
129 TVirtualFFT* fft_back = TVirtualFFT::FFT(1, &size, "C2R ES");
130 fft_back->SetPointsComplex(reres, imres);
131 fft_back->Transform();
132 fft_back->GetPoints(result);
133 delete fft_back;
134
135 delete[] imres; delete[] reres; delete[] refirst; delete[] imfirst;
136 delete[] resecond; delete[] imsecond; delete[] regauss; delete[] imgauss; delete[] gauss;
137}
138
162void Filters::SmoothArray(double*& ar, int nbins, double sigma, string method, double bin_size, double sigma2) {
163
164 if (method == "Box" || method == "Mean" || method == "Average") BoxFilter(ar, nbins, static_cast<int>(sigma));
165 else if (method == "GausFFT") GausFFTFilter(ar, nbins, sigma, bin_size);
166 else if (method == "Gaus") GausFilter(ar, nbins, sigma, bin_size);
167 else if (method == "Bilateral") BilateralFilter(ar, nbins, sigma, sigma2, bin_size);
168 else if (method == "Bilateral2") Bilateral2Filter(ar, nbins, sigma, sigma2, .05, bin_size);
169 else if (method == "Median") MedianFilter(ar, nbins, static_cast<int>(sigma));
170 else {
171 cout << "Smooth method " << method << " not known. Using default Gaus method." << endl;
172 GausFilter(ar, nbins, sigma, bin_size);
173 }
174}
176
177
182void Filters::SmoothArray(double*& ar, int nbins, double sigma, int method, double bin_size, double sigma2) {
183
184 switch (method) {
185 case 0:
186 SmoothArray(ar, nbins, sigma, "Box");
187 break;
188 case 1:
189 SmoothArray(ar, nbins, sigma, "GausFFT", bin_size);
190 break;
191 case 3:
192 SmoothArray(ar, nbins, sigma, "Bilateral", bin_size, sigma2);
193 break;
194 case 4:
195 SmoothArray(ar, nbins, sigma, "Bilateral2", bin_size, sigma2);
196 break;
197 case 5:
198 SmoothArray(ar, nbins, sigma, "Median");
199 break;
200 default:
201 SmoothArray(ar, nbins, sigma, "Gaus", bin_size);
202 break;
203 }
204}
205
210void Filters::BoxFilter(double*& ar, int nbins, int sigma) {
211 // calculate running average from -sigma until +sigma (sigma = number of bins)
212 if (sigma == 0) sigma = 1.;
213
214 double* result = new double[nbins];
215 double w_sum = .0;
216 int w_nbins = 2 * sigma + 1;
217
218 // init
219 for (int i = 0; i <= min(nbins - 1, sigma); i++) w_sum += ar[i];
220 // start
221 for (int i = 0; i < sigma && i < nbins; i++) {
222 result[i] = w_sum / static_cast<double>(i + sigma + 1);
223 if (i + sigma + 1 < nbins) w_sum += ar[i + sigma + 1];
224 }
225 // center (full window)
226 for (int i = sigma; i < nbins - sigma; i++) {
227 if (i - sigma - 1 >= 0) w_sum -= ar[i - sigma - 1];
228 if (i + sigma < nbins) w_sum += ar[i + sigma];
229 result[i] = w_sum / static_cast<double>(w_nbins);
230 }
231 // end
232 for (int i = max(nbins - sigma, 0); i < nbins; i++) {
233 int nmean = nbins - (i - sigma);
234 result[i] = w_sum / static_cast<double>(nmean);
235 if (i - sigma >= 0) w_sum -= ar[i - sigma];
236 }
237 for (int i = 0; i < nbins; i++) ar[i] = result[i];
238 delete[] result;
239}
240
246void Filters::GausFilter(double*& ar, int nbins, double sigma, double bin_size) {
247 // gauss kernel +/-3*sigma
248 if (sigma <= 0 || bin_size <= 0) return;
249
250 double* artmp = new double[nbins];
251 for (int i = 0; i < nbins; i++) artmp[i] = ar[i];
252
253 // calculate running average from -sigma until +sigma (sigma = number of bins)
254 int nbins_6sigma = static_cast<int>(ceil(6. * sigma / bin_size));
255 if (nbins_6sigma % 2 == 0) nbins_6sigma++;
256
257 if (nbins_6sigma > 1) {
258 int half_window = nbins_6sigma * 0.5;
259 double gauss_offset = static_cast<double>(half_window) * bin_size;
260 double* gauss = new double[nbins_6sigma];
261 for (int i = 0; i < nbins_6sigma; i++) {
262 gauss[i] = TMath::Gaus(static_cast<double>(i) * bin_size, gauss_offset, sigma, true);
263 }
264
265 double res = 0, norm = 0;
266 for (int i = 0; i < nbins; i++) {
267 res = 0., norm = 0.;
268 int start = max(0, half_window - i);
269 int end = min(nbins - i + half_window, nbins_6sigma);
270
271 for (int j = start; j < end; j++) {
272 res += gauss[j] * artmp[i + j - half_window];
273 norm += gauss[j];
274 }
275 if (norm != 0.) ar[i] = res / norm;
276 }
277 delete[] gauss;
278 }
279 delete[] artmp;
280}
281
287void Filters::GausFFTFilter(double*& ar, int nbins, double sigma, double bin_size) {
288 // convolution with gauss clipped at +-5 sigma (very inefficient and slow)
289 double* gauss = new double[nbins];
290
291 double sum = 0.;
292 double five_sigma = 5 * sigma;
293 double i_denom1 = -1. / (2. * sigma * sigma);
294 double i_denom2 = 1. / sigma * 2.506628;
295
296 for (int i = 0; i < nbins; i++) {
297 double position = static_cast<double>(i) * bin_size;
298 if (position < five_sigma) gauss[i] = exp(pow((position - five_sigma), 2) * i_denom1) * i_denom2;
299 else gauss[i] = 0.;
300 sum += gauss[i];
301 }
302
303 double i_sum = (sum == 0.) ? 1. : 1. / sum;
304 for (int i = 0; i < nbins; i++) {
305 gauss[i] *= i_sum;
306 }
307
308 Convolute(ar, gauss, nbins);
309 delete[] gauss;
310}
311
321void Filters::BilateralFilter(double*& ar, int nbins, double sigma_x, double sigma_y, double bin_size) {
322 double* artmp = new double[nbins];
323 for (int i = 0; i < nbins; i++) artmp[i] = ar[i];
324
325 int nbins_3sigma = static_cast<int>(ceil(3 * sigma_x / bin_size));
326
327 for (int i = 0; i < nbins; i++) {
328 double weighted_sum = 0.0;
329 double sum_of_weights = 0.0;
330
331 for (int j = max(0, i - nbins_3sigma); j < min(nbins - 1, i + nbins_3sigma); j++) {
332 // Calculate spatial and intensity distances
333 double spatial_dist = abs(i - j) * bin_size;
334 double intensity_dist = abs(artmp[i] - artmp[j]);
335
336 // weights
337 double w_x = exp(-0.5 * (spatial_dist * spatial_dist) / (sigma_x * sigma_x));
338 double w_y = exp(-0.5 * (intensity_dist * intensity_dist) / (sigma_y * sigma_y));
339 double weight = w_x * w_y;
340
341 weighted_sum += artmp[j] * weight;
342 sum_of_weights += weight;
343 }
344
345 if (sum_of_weights != 0.) ar[i] = weighted_sum / sum_of_weights;
346 }
347 delete[] artmp;
348}
349
360void Filters::Bilateral2Filter(double*& ar, int nbins, int sigma_x, double sigma_y, double sigma_slope, double bin_size) {
361 // Asymmetric smoothing using a gauss of the change of the data relative to the central bin in a window for weighting
362 double* artmp = new double[nbins];
363 for (int i = 0; i < nbins; i++) artmp[i] = ar[i];
364
365 double i_denom = -1. / (2. * sigma_y * sigma_y);
366 int x_shift = static_cast<int>(sigma_x / bin_size);
367
368 double weightedSum = 0.;
369 double sumOfWeights = 0.;
370
371 double* slope = new double[nbins - 1];
372 double* slope_weight = new double[nbins - 1];
373 for (int i = 0; i < nbins - 1; i++) {
374 slope[i] = ar[i + 1] - ar[i];
375 slope_weight[i] = 1. / (1. + pow(.5 * slope[i] / sigma_slope, 2.));
376 }
377
378 for (int i = 0; i < nbins; i++) {
379 weightedSum = 0.;
380 sumOfWeights = 0.;
381 int start = max(0, i - x_shift);
382 int end = min(i + x_shift, nbins - 1);
383
384 for (int j = start; j <= end; j++) {
385 double diff = artmp[i] - artmp[j];
386 double weight = exp(pow(diff, 2) * i_denom);
387 if (j < nbins - 1) weight *= slope_weight[j];
388 weightedSum += weight * artmp[j];
389 sumOfWeights += weight;
390 }
391 if (sumOfWeights != 0) ar[i] = weightedSum / sumOfWeights;
392 }
393 delete[] artmp;
394 delete[] slope;
395 delete[] slope_weight;
396}
397
402void Filters::MedianFilter(double*& ar, int nbins, int window_size) {
403 double* artmp = new double[nbins];
404 for (int i = 0; i < nbins; i++) artmp[i] = ar[i];
405
406 int half_window = window_size * 0.5;
407
408 for (int i = 0; i < nbins; i++) {
409 double* window = new double[window_size];
410 int window_index = 0;
411
412 for (int j = i - half_window; j <= i + half_window; j++) {
413 if (j >= 0 && j < nbins) {
414 window[window_index] = artmp[j];
415 window_index++;
416 }
417 }
418
419 // calculate the median
420 sort(window, window + window_index);
421 double median = (window_index % 2 == 0) ?
422 (window[window_index / 2 - 1] + window[window_index / 2]) * 0.5 :
423 window[window_index / 2];
424
425 ar[i] = median;
426
427 delete[] window;
428 }
429 delete[] artmp;
430}
431
445void Filters::ResponseFilter(double*& ar, int nbins, double sigma1, double sigma2, double factor, double bin_size) {
446
447 double* artmp = new double[nbins];
448 for (int i = 0; i < nbins; i++) artmp[i] = ar[i];
449
450 // shifted difference of two gauss functions (~smoothed derivative)
451 int nbins_23sigma = static_cast<int>(ceil((2. * sigma1 + 3. * sigma2) / bin_size));
452 int nbins_3sigma = static_cast<int>(ceil(3. * sigma2 / bin_size));
453 double* sdog = new double[nbins_23sigma];
454
455 double i_denom1 = 1. / (2. * sigma1 * sigma1);
456 double i_denom2 = 1. / (2. * sigma2 * sigma2);
457 double norm = 0.;
458 for (int i = 0; i < nbins_23sigma; i++) {
459 sdog[i] = exp(-1. * pow(static_cast<double>(i) * bin_size - 3. * sigma2, 2.) * i_denom1) -
460 factor * exp(-1. * pow(static_cast<double>(i) * bin_size - 2. * sigma2, 2.) * i_denom2);
461 if (sdog[i] > 0.) norm += sdog[i];
462 }
463 if (norm == 0.) norm = 1.;
464 double i_norm = 1. / norm;
465
466 for (int i = 0; i < nbins; i++) {
467 double res = 0.;
468 for (int j = -1 * nbins_3sigma; j < nbins_23sigma - nbins_3sigma; j++) {
469 if (i + j >= 0 && i + j < nbins) res += sdog[j + nbins_3sigma] * artmp[i + j];
470 else if (i + j < 0) res += sdog[j + nbins_3sigma] * artmp[0];
471 else res += sdog[j + nbins_3sigma] * artmp[nbins - 1];
472 }
473 ar[i] = res * i_norm;
474 }
475 delete[] artmp;
476 delete[] sdog;
477}
478
492void Filters::SecondOrderUnderdampedFilter(double*& ar, int nbins, double period, double damping, double shift, double bin_size, bool do_plot) {
493
494 double* artmp = new double[nbins];
495 for (int i = 0; i < nbins; i++) artmp[i] = ar[i];
496
497 int shift_n = static_cast<int>(ceil(shift / bin_size));
498 int nbins_response = shift_n + static_cast<int>(3. * period / bin_size);
499
500 if (nbins_response > nbins) {
501 cout << "ERROR: nbins_response > nbins. Reduce shift or period." << endl;
502 nbins_response = nbins - 1;
503 }
504 if (damping <= 0) {
505 cout << "ERROR: damping <= 0. Setting to 1 ns" << endl;
506 damping = 1;
507 }
508
509 if (nbins_response % 2 == 0) nbins_response++;
510 double* souf = new double[nbins_response];
511
512 double exp_factor = -1. / damping;
513 double omega = 2. * M_PI / period;
514 double sum_norm = 0.;
515 for (int i = 0; i < nbins_response; i++) {
516 double position = static_cast<double>(i) * bin_size;
517
518 if (i < shift_n) souf[i] = 0.;
519 else {
520 souf[i] = exp(exp_factor * (position - shift)) * sin(omega * (position - shift));
521 sum_norm += souf[i];
522 }
523 }
524 if (sum_norm == 0.) sum_norm = 1.;
525 double i_sum_norm = 1. / sum_norm;
526
527 for (int i = 0; i < nbins; i++) {
528 double res = 0.;
529 for (int j = max(0, nbins_response / 2 - i); j < min(nbins_response, nbins + nbins_response / 2 - i); j++) {
530 res += souf[nbins_response - j - 1] * artmp[i + j - nbins_response / 2];
531 }
532 ar[i] = res * i_sum_norm;
533 }
534
535 // plot the results
536 if (do_plot) {
537 double* x = new double[nbins_response];
538 double* x_full = new double[nbins];
539 for (int i = 0; i < nbins; i++) {
540 x_full[i] = static_cast<double>(i) * bin_size;
541 if (i < nbins_response) x[i] = static_cast<double>(i) * bin_size;
542 }
543
544 TCanvas* c_response = new TCanvas("response", "response", 800, 600);
545 TGraph* g_original = new TGraph(nbins, x_full, artmp);
546 g_original->SetTitle("original");
547 g_original->SetMarkerStyle(34);
548 g_original->Draw("AP");
549
550 TGraph* g_response = new TGraph(nbins_response, x, souf);
551 g_response->SetTitle("response");
552 g_response->SetLineColor(kBlue);
553 g_response->SetLineWidth(5);
554 g_response->SetLineStyle(7);
555 g_response->Draw("L same");
556
557 TGraph* g_result = new TGraph(nbins, x_full, ar);
558 g_result->SetTitle("result");
559 g_result->SetLineColor(kRed);
560 g_result->SetLineWidth(3);
561 g_result->Draw("L same");
562
563 g_response->Scale(TMath::MaxElement(nbins, ar) / TMath::MaxElement(nbins_response, souf));
564 c_response->Update(); c_response->Modified();
565 gPad->BuildLegend(.5, .7, .9, .9);
566 c_response->SaveAs("response.png");
567
568 delete[] x;
569 delete[] x_full;
570 }
571 delete[] artmp;
572 delete[] souf;
573}
574
578void Filters::SecondOrderUnderdampedFilter(TH1F*& his, int nbins, double period, double damping, double shift, double bin_size, bool do_plot) {
579 double* yvals = Helpers::gety(his);
580 SecondOrderUnderdampedFilter(yvals, nbins, period, damping, shift, bin_size, do_plot);
581 for (int i = 1; i <= his->GetNbinsX(); i++) his->SetBinContent(i, yvals[i - 1]);
582 delete[] yvals;
583}
static void Deconvolute(double *&, double *, double *, int, double, double=0., double=0.)
Helper to perform partial deconvolution of two 1D arrays.
Definition Filters.cc:81
static void BilateralFilter(double *&, int, double, double, double)
Bilateral filter (non-linear) with 3 sigma kernel.
Definition Filters.cc:321
static void MedianFilter(double *&ar, int, int)
Median filter.
Definition Filters.cc:402
static void ResponseFilter(double *&, int, double=.4, double=1.2, double=.25, double=.3125)
Custom filter emulating primitive response function.
Definition Filters.cc:445
static void BoxFilter(double *&, int, int)
Simple running average (box) filter.
Definition Filters.cc:210
static void GausFFTFilter(double *&, int, double, double)
Gaussian smoothing with FFT convolution.
Definition Filters.cc:287
static void SecondOrderUnderdampedFilter(double *&, int, double, double, double, double=.3125, bool=false)
Shifted second order underdamped filter.
Definition Filters.cc:492
static void Convolute(double *&, double *, double *, int)
Helper to perform convolution of two 1D arrays.
Definition Filters.cc:12
static void GausFilter(double *&, int, double, double)
Gaussian smoothing with simple kernel.
Definition Filters.cc:246
static void Bilateral2Filter(double *&, int, int, double, double, double)
Nonlinear filter based on bilateral filter.
Definition Filters.cc:360
static void SmoothArray(double *&, int, double=.625, string="Gaus", double=.3125, double=1.5)
Apply smoothing array of double with length nbins.
Definition Filters.cc:162
static double * gety(HistType *his)
Get array of y values for a histogram.
Definition Helpers.cc:138