InterferometricMap.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 INTERFEROMETRIC_MAP_H
10 #define INTERFEROMETRIC_MAP_H
11 
12 #define NUM_BINS_THETA 70
13 #define NUM_BINS_PHI 9
14 
15 #define MIN_THETA -55
16 #define MAX_THETA 35
17 #define PHI_RANGE 22.5
18 
19 
20 #define DEGREES_IN_CIRCLE 360
21 
22 #define MAX_NUM_PEAKS 5
23 #define PEAK_PHI_DEG_RANGE 10
24 #define PEAK_THETA_DEG_RANGE 10
25 
26 #define NUM_BINS_THETA_ZOOM 40
27 #define NUM_BINS_PHI_ZOOM 40
28 #define ZOOM_BINS_PER_DEGREE_PHI 4
29 #define ZOOM_BINS_PER_DEGREE_THETA 4
30 
31 #define ZOOM_BIN_SIZE_PHI (1./ZOOM_BINS_PER_DEGREE_PHI)
32 #define ZOOM_BIN_SIZE_THETA (1./ZOOM_BINS_PER_DEGREE_THETA)
33 #define THETA_RANGE_ZOOM (NUM_BINS_THETA_ZOOM*ZOOM_BIN_SIZE_THETA)
34 #define PHI_RANGE_ZOOM (NUM_BINS_PHI_ZOOM*ZOOM_BIN_SIZE_PHI)
35 
36 #define NUM_BINS_PHI_ZOOM_TOTAL (DEGREES_IN_CIRCLE*ZOOM_BINS_PER_DEGREE_PHI + NUM_BINS_PHI_ZOOM)
37 #define NUM_BINS_THETA_ZOOM_TOTAL ((MAX_THETA - MIN_THETA)*ZOOM_BINS_PER_DEGREE_THETA + NUM_BINS_THETA_ZOOM)
38 
39 #define NUM_BINS_QUAD_FIT_PHI 5
40 #define NUM_BINS_QUAD_FIT_THETA 5
41 
42 #include "AnitaConventions.h"
43 
44 #include "TH2D.h"
45 #include "TGraph.h"
46 #include <map>
47 #include "TDecompSVD.h"
48 #include "TGraphInteractive.h"
49 
50 class Adu5Pat;
51 class UsefulAdu5Pat;
52 class TProfile2D;
53 class TruthAnitaEvent;
54 
55 namespace Acclaim
56 {
57 
58 class InterferometryCache;
59 class CrossCorrelator;
60 
61 
62 class InterferometricMap : public TH2D, public GuiParent {
63 
64  public:
66  InterferometricMap(Int_t peakInd, Int_t phiSector, Double_t centrePhi, Double_t phiRange, Double_t centreTheta, Double_t thetaRange);
67  virtual ~InterferometricMap();
68 
70  void findPeakValues(Int_t numPeaks, std::vector<Double_t>& peakValues, std::vector<Double_t>& phiDegs, std::vector<Double_t>& thetaDegs) const;
71  void getPeakInfo(Double_t& peak, Double_t& phiDeg, Double_t& thetaDeg) const{
72  peak = fPeakBinValue;
73  phiDeg = fPeakPhi;
74  thetaDeg = fPeakTheta;
75  }
76  void getPeakInfo(Double_t& peak, Double_t& phiDeg, Double_t& thetaDeg, Double_t& reducedChisquare) const{
77  getPeakInfo(peak, phiDeg, thetaDeg);
78  reducedChisquare = fPeakReducedChisquare;
79  }
80 
81  virtual void Draw(Option_t* opt);
82  virtual void ExecuteEvent(int event, int x, int y);
83 
84  void addGpsInfo(const Adu5Pat* pat);
85  void addTruthInfo(const TruthAnitaEvent* truth);
86  // void addGpsInfo(const UsefulAdu5Pat* usefulPat);
87 
88 
99  TPad* makeProjectionCanvas(TPad* pad);
100 
106  TPad* makeProjectionCanvas(){return makeProjectionCanvas(NULL);} //*MENU*
107 
108  inline Int_t GetNbinsPhi() const {return GetNbinsX();}
109  inline Int_t GetNbinsTheta() const {return GetNbinsY();}
110  inline const TAxis* GetPhiAxis() const {return GetXaxis();}
111  inline const TAxis* GetThetaAxis() const {return GetYaxis();}
112 
113  const TGraphInteractive* getPeakPointGraph(); // for plotting
114  const TGraphInteractive* getEdgeBoxGraph(); // for plotting
115  const TGraphInteractive* getSunGraph(); // for plotting
116  const TGraphInteractive* getTruthGraph(); // for plotting
117 
118  bool isAZoomMap() const {return fIsZoomMap;}
119  Int_t getPeakPhiSector() const {return fPeakPhiSector;}
120  AnitaPol::AnitaPol_t getPol() const {return pol;}
121 
122  virtual void Reset(Option_t* = "");
123 
132  void setResolutionEstimateFromWaveformSNR(double snr);
133 
134 
135  protected:
136 
137  UsefulAdu5Pat* fUsefulPat;
138  void deletePat();
139  Double_t truthLat;
140  Double_t truthLon;
141  Double_t truthAlt;
142 
143  void fitPeakWithQuadratic(Int_t peakPhiBin, Int_t peakThetaBin);
144  void setPeakInfoJustFromPeakBin(Int_t peakPhiBin, Int_t peakThetaBin);
145 
146  TGraph& findOrMakeGraph(TString name);
147 
148  // name and title
150  UInt_t eventNumber;
151  void setNameAndTitle(); // once we have the eventNumber/polarization
152  void setDefaultName(); // before then...
153 
154  bool fThetaAxisInSinTheta;
155  void initializeInternals();
156 
157 
158 
159  // doing the zoomed in maps requires knowing a little more information
160  // isZoomMap = false, and all other = -1 if doing a coarse map
161  bool fIsZoomMap;
162  int fPeakPhiSector;
163  int minThetaBin;
164  int minPhiBin;
165  int peakIndex; // for name and title of zoomed map only
166 
167  // Get the appropriate bin edges for the zoom map
168  void getIndicesOfEdgeBins(const std::vector<double>& binEdges, Double_t lowVal, Double_t highVal, Int_t& lowIndex, Int_t& highIndex);
169 
170  Double_t fPeakBinValue;
171  Int_t fPeakBinPhi;
173 
174  Double_t fPeakValue;
175  Double_t fPeakPhi;
176  Double_t fPeakTheta;
178  Double_t fPeakSigmaPhi;
179  Double_t fPeakSigmaTheta;
180  Double_t fPeakCovariance;
181 
182  Double_t fSigmaTheta;
183  Double_t fSigmaPhi;
184 
186 
187  // static members, may end up elsewhere at some point
188  public:
189  static const std::vector<Double_t>& getCoarseBinEdgesTheta();
190  static const std::vector<Double_t>& getFineBinEdgesTheta();
191  static const std::vector<Double_t>& getCoarseBinEdgesPhi();
192  static const std::vector<Double_t>& getFineBinEdgesPhi();
193  static Double_t getBin0PhiDeg();
194  static Double_t getPhiSectorCenterPhiDeg(int phi);
195  static Int_t getPhiSectorFromPhiRough(double phiRough);
196 
197  private:
198  static const TDecompSVD& getSVD();
199  static void getDeltaBinRangeSVD(Int_t& dPhiBinLow, Int_t& dPhiBinHigh, Int_t& dThetaBinLow, Int_t& dThetaBinHigh);
200 
201 };
202 }
203 
204 
205 
206 #endif // INTERFEROMETRIC_MAP_H
Double_t fPeakSigmaTheta
Width of peak in phi.
static const std::vector< Double_t > & getFineBinEdgesTheta()
void addTruthInfo(const TruthAnitaEvent *truth)
Double_t fPeakReducedChisquare
Location in theta of the peak of the quadratic fit to the region around the peak bin.
Adu5Pat – The ADU5 Position and Attitude Data.
Definition: Adu5Pat.h:26
static const std::vector< Double_t > & getFineBinEdgesPhi()
void setPeakInfoJustFromPeakBin(Int_t peakPhiBin, Int_t peakThetaBin)
static const std::vector< Double_t > & getCoarseBinEdgesTheta()
Double_t fPeakSigmaPhi
Residual of the peak fit... I hope.
ClassDef(InterferometricMap, 4)
Estimate of of the pointing resultion from the waveform SNR.
virtual void Draw(Option_t *opt)
void Fill(AnitaPol::AnitaPol_t pol, CrossCorrelator *cc, InterferometryCache *dtCache)
void fitPeakWithQuadratic(Int_t peakPhiBin, Int_t peakThetaBin)
Fit a 2D quadratic function to the bins around the peak of the map. Only implemented for zoomed maps...
Inherit from this to draw interactive TGraphs on top of you.
Double_t fPeakCovariance
Width of peak theta.
Namespace which wraps everything in the library.
Double_t fSigmaPhi
Estimate of the pointing resolution from the waveform SNR.
Double_t fPeakPhi
Value of the peak of quadratic fit to the region around the peak bin.
Double_t fPeakValue
The bin in theta (counting from 0) containing the peak value.
Double_t fSigmaTheta
Width of peak theta.
TruthAnitaEvent – The Truth ANITA Event.
static Double_t getPhiSectorCenterPhiDeg(int phi)
As UsefulAnitaEvent is to RawAnitaEvent, UsefulAdu5Pat is to Adu5Pat. Well not quite as useful but yo...
Definition: UsefulAdu5Pat.h:39
virtual void ExecuteEvent(int event, int x, int y)
A class to take in FiteredAnitaEvents and cross-correlate nearby channels.
Class to cache the deltaTs or parts of their calculation for fast map making.
Int_t fPeakBinTheta
The bin in phi (counting from 0) containing the peak value.
void setResolutionEstimateFromWaveformSNR(double snr)
static const std::vector< Double_t > & getCoarseBinEdgesPhi()
Int_t fPeakBinPhi
The value of the highest bin in the map.
void addGpsInfo(const Adu5Pat *pat)
enum AnitaPol::EAnitaPol AnitaPol_t
Polarisation enumeration.
Double_t fPeakTheta
Location in phi of the peak of the quadratic fit to the region around the peak bin.