CrossCorrelator.h
1 /* -*- C++ -*-.*********************************************************************************************
2  Author: Ben Strutt
3  Email: strutt@physics.ucla.edu
4 
5  Description:
6  A Cross Correlator to interact with the ROOTified ANITA-3 data and do some interferometry.
7 ***********************************************************************************************************/
8 
9 #ifndef CROSSCORRELATOR_H
10 #define CROSSCORRELATOR_H
11 
12 // Ryan things
13 // #include "UsefulAnitaEvent.h"
14 #include "FilteredAnitaEvent.h"
15 #include "AnalysisWaveform.h"
16 #include "AnitaEventCalibrator.h"
17 #include "AnitaGeomTool.h"
18 #include "UsefulAdu5Pat.h"
19 #include "FFTtools.h"
20 
21 #include "AnitaEventSummary.h"
22 
23 // My things
24 #include "RootTools.h"
25 #include "FancyFFTs.h"
26 
27 // ROOT things
28 #include "TGraph.h"
29 #include "TH2D.h"
30 #include "TROOT.h" // for gDirectory pointer?
31 
32 // standard c++ things
33 #include <iostream>
34 
35 
36 
37 // Anita & Geometry definitions
38 // #define PEAK_THETA_DEG_RANGE 180
39 
40 #define SPEED_OF_LIGHT 2.99792458e8
41 #define SPEED_OF_LIGHT_NS 0.299792458
42 
43 
44 // Number of phi-sectors to cross correlate between
45 #define DELTA_PHI_SECT 2
46 
47 // Offline reconstruction definitions
48 #define NUM_COMBOS 336
49 
50 // Typical number of samples in waveform
51 #define UPSAMPLE_FACTOR 6
52 #define NOMINAL_SAMPLING_DELTAT (1./2.6f)
53 #define PAD_FACTOR 2
54 #define GET_NUM_FREQS(n)((n)/2+1)
55 
56 namespace Acclaim
57 {
58 
65 
66  public:
67 
68  //--------------------------------------------------------------------------------------------------------
69  // Public member functions
70  //--------------------------------------------------------------------------------------------------------
71 
73  virtual ~CrossCorrelator();
74 
75  virtual void correlateEvent(const FilteredAnitaEvent* realEvent);
76  virtual void correlateEvent(const FilteredAnitaEvent* realEvent, AnitaPol::AnitaPol_t pol);
77 
78  void doUpsampledCrossCorrelations(AnitaPol::AnitaPol_t pol, Int_t phiSector);
79 
80  TGraph* getCrossCorrelationGraph(AnitaPol::AnitaPol_t pol, Int_t ant1, Int_t ant2) const;
81  TGraph* getUpsampledCrossCorrelationGraph(AnitaPol::AnitaPol_t pol, Int_t ant1, Int_t ant2) const;
82  TGraph* getCrossCorrelationGraphWorker(Int_t numSamps, AnitaPol::AnitaPol_t pol, Int_t ant1, Int_t ant2) const;
83 
84 
85  virtual Double_t getCrossCorrelation(AnitaPol::AnitaPol_t pol, Int_t combo, Double_t deltaT) const;
86  Double_t getTimeOfMaximumUpsampledCrossCorrelation(AnitaPol::AnitaPol_t pol, Int_t ant1, Int_t ant2) const;
87 
88  void initializeVariables();
89  void getNormalizedInterpolatedTGraphs(const FilteredAnitaEvent* realEvent, AnitaPol::AnitaPol_t pol, bool raw = false);
90  void renormalizeFourierDomain(AnitaPol::AnitaPol_t pol, Int_t ant);
91  void doFFTs(AnitaPol::AnitaPol_t pol);
92  void doCrossCorrelations(AnitaPol::AnitaPol_t pol);
93 
94 
95  Bool_t useCombo(Int_t ant1, Int_t ant2, Int_t phiSector, Int_t deltaPhiSect);
96  void fillCombosToUse();
97  void do5PhiSectorCombinatorics();
98 
99  Double_t getInterpolatedUpsampledCorrelationValue(AnitaPol::AnitaPol_t pol, Int_t combo, Double_t deltaT) const;
100 
101 
102 
103 
104  //--------------------------------------------------------------------------------------------------------
105  // Public member variables
106  //--------------------------------------------------------------------------------------------------------
107 
108  std::complex<Double_t> ffts[AnitaPol::kNotAPol][NUM_SEAVEYS][GET_NUM_FREQS(NUM_SAMP*PAD_FACTOR)];
109  Double_t crossCorrelations[AnitaPol::kNotAPol][NUM_COMBOS][NUM_SAMP*PAD_FACTOR];
110  Double_t crossCorrelationsUpsampled[AnitaPol::kNotAPol][NUM_COMBOS][NUM_SAMP*PAD_FACTOR*UPSAMPLE_FACTOR*PAD_FACTOR];
111 
112  Double_t fVolts[AnitaPol::kNotAPol][NUM_SEAVEYS][NUM_SAMP*PAD_FACTOR];
113  Double_t startTimes[AnitaPol::kNotAPol][NUM_SEAVEYS];
114 
115  Double_t interpRMS[AnitaPol::kNotAPol][NUM_SEAVEYS];
116  Double_t interpRMS2[AnitaPol::kNotAPol][NUM_SEAVEYS];
117  Int_t comboIndices[NUM_SEAVEYS][NUM_SEAVEYS];
118 
121  Double_t correlationDeltaT;
122  Int_t numSamples;
124  Int_t numCombos;
125 
126  std::vector<Int_t> comboToAnt1s;
127  std::vector<Int_t> comboToAnt2s;
128  std::vector<Int_t> combosToUseGlobal[NUM_PHI];
129 
130 
134 
135  private:
136 
137  //--------------------------------------------------------------------------------------------------------
138  // Private member variables
139  //--------------------------------------------------------------------------------------------------------
140 
141 
142  };
143 
144 
145 
154 
155  public:
156  TemplateCorrelator(Int_t run, UInt_t eventNumber);
157  virtual ~TemplateCorrelator();
158  void initTemplate(const FilteredAnitaEvent* fEv);
159  void initTemplate(Int_t run, UInt_t eventNumber);
160  virtual void correlateEvent(const FilteredAnitaEvent* fEv);
161  virtual void correlateEvent(const FilteredAnitaEvent* fEv, AnitaPol::AnitaPol_t pol);
162 
163  TGraph* getCrossCorrelationGraph(AnitaPol::AnitaPol_t pol, Int_t ant) const;
164 
165  Double_t getPeakCorrelation(AnitaPol::AnitaPol_t pol, Double_t minOffset=-100, Double_t maxOffset=100, Double_t stepSize=NOMINAL_SAMPLING_DELTAT) const;
166 
167  virtual Double_t getCrossCorrelation(AnitaPol::AnitaPol_t pol, Int_t ant, Double_t deltaT) const;
168 
169  protected:
170  // double templateAllChannelRMS;
171 
172 };
173 
174 }
175 
176 #endif
Int_t comboIndices[NUM_SEAVEYS][NUM_SEAVEYS]
Array mapping ant1+ant2 to combo index.
virtual void correlateEvent(const FilteredAnitaEvent *fEv)
std::vector< Int_t > comboToAnt1s
Vector mapping combo index to ant1.
Double_t correlationDeltaT
nominalSamplingDeltaT/UPSAMPLE_FACTOR, deltaT of for interpolation.
Int_t multiplyTopRingByMinusOne
For showing how I&#39;m an idiot with respect to compiling the ANITA-3 prioritizer.
std::vector< Int_t > comboToAnt2s
Vector mapping combo index to ant2.
Int_t kDeltaPhiSect
Specifies how many phi-sectors around peak use in reconstruction.
std::complex< Double_t > ffts[AnitaPol::kNotAPol][NUM_SEAVEYS][((NUM_SAMP *2)/2+1)]
FFTs of evenly resampled waveforms.
TemplateCorrelator(Int_t run, UInt_t eventNumber)
Double_t interpRMS2[AnitaPol::kNotAPol][NUM_SEAVEYS]
RMS of interpolated TGraphs with extra zero padding.
Int_t numCombos
Number of possible antenna pairs, counted during initialization. Should equal NUM_COMBOS.
Double_t interpRMS[AnitaPol::kNotAPol][NUM_SEAVEYS]
RMS of interpolated TGraphs.
Namespace which wraps everything in the library.
USeful in for loops.
TGraph * getCrossCorrelationGraph(AnitaPol::AnitaPol_t pol, Int_t ant) const
UInt_t eventNumber[AnitaPol::kNotAPol]
For tracking event number.
void initTemplate(const FilteredAnitaEvent *fEv)
Double_t nominalSamplingDeltaT
ANITA-3 => 1./2.6 ns, deltaT for evenly resampling.
std::vector< Int_t > combosToUseGlobal[NUM_PHI]
Depends on L3 trigger for global image.
Int_t numSamples
Number of samples in waveform after padding.
Cross-correlate all channels of an event with a template.
A class to take in FiteredAnitaEvents and cross-correlate nearby channels.
Int_t numSamplesUpsampled
Number of samples in waveform after padding and up sampling.
Double_t crossCorrelationsUpsampled[AnitaPol::kNotAPol][336][NUM_SAMP *2 *6 *2]
Upsampled cross correlations.
Double_t crossCorrelations[AnitaPol::kNotAPol][336][NUM_SAMP *2]
Cross correlations.
This class is intended to store all the necessary data about an ANITA event for filtering and analysi...
enum AnitaPol::EAnitaPol AnitaPol_t
Polarisation enumeration.
Int_t kOnlyThisCombo
For debugging, only fill histograms with one particular antenna pair.
Double_t fVolts[AnitaPol::kNotAPol][NUM_SEAVEYS][NUM_SAMP *2]
Hold the filtered waveforms for padding...