Analyzer.h
1 #ifndef UCORRELATOR_ANALYZER_H
2 #define UCORRELATOR_ANALYZER_H
3 
4 #include "Correlator.h"
5 #include "WaveformCombiner.h"
6 #include "UCorrelatorGUI.h"
7 #include <vector>
8 #include "TH2D.h"
9 #include "simpleStructs.h"
10 #include "AnitaEventSummary.h"
11 #include "ResponseManager.h"
12 #include "UsefulAnitaEvent.h"
13 
14 class FilteredAnitaEvent;
15 class UsefulAdu5Pat;
16 class AnalysisWaveform;
17 class TPad;
18 
19 namespace FFTtools
20 {
21  class DigitalFilter;
22 }
23 
24 namespace UCorrelator
25 {
26 
27  class AnalysisConfig;
28 
34  class Analyzer
35  {
36 
37  public:
40  Analyzer(const AnalysisConfig * cfg = 0, bool interactive = false);
41 
43  virtual ~Analyzer();
44 
45 
46 
48  void analyze(const FilteredAnitaEvent * event, AnitaEventSummary *summary, const TruthAnitaEvent * truth = 0);
49 
52  Correlator * getCorrelator() { return &corr; }
53 
56  WaveformCombiner * getWaveformCombiner() { return &wfcomb; }
57 
60  WaveformCombiner * getXPolWaveformCombiner() { return &wfcomb_xpol; }
61 
62 
64  const gui::Map * getCorrelationMap(AnitaPol::AnitaPol_t pol) const { return correlation_maps[pol] ; }
65 
67  const gui::Map * getZoomedCorrelationMap(AnitaPol::AnitaPol_t pol, int i) const { return zoomed_correlation_maps[pol][i]; }
68 
70  const AnalysisWaveform * getCoherent(AnitaPol::AnitaPol_t pol, int i, bool filtered = false) const { return coherent[pol][filtered ? 1 : 0][i]; }
71 
73  const AnalysisWaveform * getDeconvolved(AnitaPol::AnitaPol_t pol, int i, bool filtered = false) const { return deconvolved[pol][filtered ? 1: 0][i]; }
74 
76  const AnalysisWaveform * getCoherentXpol(AnitaPol::AnitaPol_t pol, int i, bool filtered = false) const { return coherent_xpol[pol][filtered ? 1: 0][i]; }
77 
79  const AnalysisWaveform * getDeconvolvedXpol(AnitaPol::AnitaPol_t pol, int i, bool filtered = false) const { return deconvolved_xpol[pol][filtered ? 1: 0][i]; }
80 
82  const TGraphAligned * getCoherentPower(AnitaPol::AnitaPol_t pol, int i, bool filtered = false) const { return coherent_power[pol][filtered ? 1: 0][i]; }
83 
85  const TGraphAligned * getDeconvolvedPower(AnitaPol::AnitaPol_t pol, int i, bool filtered = false) const { return deconvolved_power[pol][filtered ? 1: 0][i]; }
86 
87  AnitaResponse::ResponseManager * getResponseManager() { return &responses; }
88 
89  /* Return the summary. This only works in interactive mode. */
90  const AnitaEventSummary* getSummary() const { return &last; }
91 
92 
95  void drawSummary(TPad *chpol = 0, TPad * cvpol = 0, int draw_filtered = 0) const;
96 
97  double getRoughPhi(AnitaPol::AnitaPol_t pol, int i) const { return rough_peaks[pol][i].first; }
98  double getRoughTheta(AnitaPol::AnitaPol_t pol, int i) const { return -rough_peaks[pol][i].second; }
99 
100  void clearInteractiveMemory(double frac = 0.5) const;
101 
103  void setDisallowedAntennas(ULong64_t hpol=0, ULong64_t vpol=0) {disallowedAnts[0] = hpol; disallowedAnts[1] = vpol; }
104 
106  void setExcludePhiRange(double phiMin, double phiMax) {phiRange[0] = phiMin; phiRange[1] = phiMax; exclude = true; }
108  void setExcludeThetaRange(double thetaMin, double thetaMax) {thetaRange[0] = thetaMin; thetaRange[1] = thetaMax; exclude = true; }
110  void setExcludeThetaPhiRange(double thetaMin, double thetaMax, double phiMin, double phiMax) {thetaRange[0] = thetaMin; thetaRange[1] = thetaMax; phiRange[0] = phiMin; phiRange[1] = phiMax; exclude = true; }
111 
113  void setPhiRange(double phiMin, double phiMax) {phiRange[0] = phiMin; phiRange[1] = phiMax; exclude = false; }
115  void setThetaRange(double thetaMin, double thetaMax) {thetaRange[0] = thetaMin; thetaRange[1] = thetaMax; exclude = false; }
117  void setThetaPhiRange(double thetaMin, double thetaMax, double phiMin, double phiMax) {thetaRange[0] = thetaMin; thetaRange[1] = thetaMax; phiRange[0] = phiMin; phiRange[1] = phiMax; exclude = false; }
119  void setTrackSource(double setLon, double setLat, double setAlt, double setdTheta = 2.5, double setdPhi = 5., bool blockOut = false) {sourceLon = setLon; sourceLat = setLat; sourceAlt = setAlt; dTheta = setdTheta; dPhi = setdPhi; exclude = blockOut; trackSource = true; }
121  void setTrackSun(double setdTheta = 2.5, double setdPhi = 5., bool blockOut = false) {dTheta = setdTheta; dPhi = setdPhi; exclude = blockOut; trackSun = true; }
123  void setTrackWAIS(double setdTheta = 2.5, double setdPhi = 5., bool blockOut = false) {sourceLon = AnitaLocations::getWaisLongitude(); sourceLat = AnitaLocations::getWaisLatitude(); sourceAlt = AnitaLocations::getWaisAltitude(); dTheta = setdTheta; dPhi = setdPhi; exclude = blockOut; trackSource = true; }
125  void setTrackLDB(double setdTheta = 2.5, double setdPhi = 5., bool blockOut = false) {sourceLon = AnitaLocations::LONGITUDE_LDB; sourceLat = AnitaLocations::LATITUDE_LDB; sourceAlt = AnitaLocations::ALTITUDE_LDB; dTheta = setdTheta; dPhi = setdPhi; exclude = blockOut; trackSource = true; }
127  void setExtraFilters(FilterStrategy* extra);
130 
131  private:
132 
133  void fillWaveformInfo(const AnalysisWaveform * wf, const AnalysisWaveform * xpol_wf, const TGraph * power, AnitaEventSummary::WaveformInfo * info, AnitaPol::AnitaPol_t pol, double rms, double vpp = 0);
134  void fillPointingInfo(double rough_phi, double rough_theta, AnitaEventSummary::PointingHypothesis * point,
135  UsefulAdu5Pat * pat, double hwAngle, UShort_t triggered_sectors, UShort_t masked_sectors, UShort_t triggered_sectors_xpol, UShort_t masked_sectors_xpol);
136  void fillFlags(const FilteredAnitaEvent * fae, AnitaEventSummary* summary, UsefulAdu5Pat * pat);
137  void fillChannelInfo(const FilteredAnitaEvent* event, AnitaEventSummary* summary);
138 
139  gui::Map* correlation_maps[2];
140  std::vector<gui::Map*> zoomed_correlation_maps[2];
141  std::vector<AnalysisWaveform *> coherent[2][2];
142  std::vector<AnalysisWaveform *> deconvolved[2][2];
143  std::vector<TGraphAligned *> coherent_power[2][2];
144  std::vector<TGraphAligned *> deconvolved_power[2][2];
145  std::vector<AnalysisWaveform *> coherent_xpol[2][2];
146  std::vector<AnalysisWaveform *> deconvolved_xpol[2][2];
147  std::vector<TGraphAligned *> coherent_power_xpol[2][2];
148  std::vector<TGraphAligned *> deconvolved_power_xpol[2][2];
149  std::vector<std::pair<double,double> > rough_peaks[2];
150  AnitaEventSummary last; //used in interactive mode by drawSummary
151  TGraph* avg_spectra[2];
152 
153  const AnalysisConfig * cfg;
154  Correlator corr;
156  WaveformCombiner wfcomb;
157  WaveformCombiner wfcomb_xpol;
158  WaveformCombiner wfcomb_filtered;
159  WaveformCombiner wfcomb_xpol_filtered;
160  TH2D *zoomed;
161  double maprms;
162  FFTtools::DigitalFilter * power_filter;
163  bool interactive;
164  bool interactive_deconvolved;
165  bool interactive_xpol_deconvolved;
166  ULong64_t disallowedAnts[2];
167  double phiRange[2];
168  double thetaRange[2];
169  bool exclude;
170  bool trackSource;
171  double sourceLon;
172  double sourceLat;
173  double sourceAlt;
174  bool trackSun;
175  double dPhi;
176  double dTheta;
177 
178  mutable std::vector<TObject*> delete_list;
179  };
180 
181 }
182 
183 
184 #endif
void setTrackLDB(double setdTheta=2.5, double setdPhi=5., bool blockOut=false)
Definition: Analyzer.h:125
void analyze(const FilteredAnitaEvent *event, AnitaEventSummary *summary, const TruthAnitaEvent *truth=0)
Definition: Analyzer.cc:177
const AnalysisWaveform * getDeconvolvedXpol(AnitaPol::AnitaPol_t pol, int i, bool filtered=false) const
Definition: Analyzer.h:79
void setPhiRange(double phiMin, double phiMax)
Definition: Analyzer.h:113
Correlator * getCorrelator()
Definition: Analyzer.h:52
const gui::Map * getCorrelationMap(AnitaPol::AnitaPol_t pol) const
Definition: Analyzer.h:64
Analyzer(const AnalysisConfig *cfg=0, bool interactive=false)
Definition: Analyzer.cc:62
void setDisallowedAntennas(ULong64_t hpol=0, ULong64_t vpol=0)
Definition: Analyzer.h:103
void setExtraFiltersDeconvolved(FilterStrategy *extra)
Definition: Analyzer.cc:1257
const gui::Map * getZoomedCorrelationMap(AnitaPol::AnitaPol_t pol, int i) const
Definition: Analyzer.h:67
void setThetaRange(double thetaMin, double thetaMax)
Definition: Analyzer.h:115
A filter strategy defines the sets of filters that are used and provides some introspection abilities...
void setTrackSource(double setLon, double setLat, double setAlt, double setdTheta=2.5, double setdPhi=5., bool blockOut=false)
Definition: Analyzer.h:119
WaveformCombiner * getWaveformCombiner()
Definition: Analyzer.h:56
void setExcludePhiRange(double phiMin, double phiMax)
Definition: Analyzer.h:106
void setExcludeThetaPhiRange(double thetaMin, double thetaMax, double phiMin, double phiMax)
Definition: Analyzer.h:110
void setTrackWAIS(double setdTheta=2.5, double setdPhi=5., bool blockOut=false)
Definition: Analyzer.h:123
const Double_t LATITUDE_LDB
Latitude at LDB.
WaveformCombiner * getXPolWaveformCombiner()
Definition: Analyzer.h:60
TruthAnitaEvent – The Truth ANITA Event.
As UsefulAnitaEvent is to RawAnitaEvent, UsefulAdu5Pat is to Adu5Pat. Well not quite as useful but yo...
Definition: UsefulAdu5Pat.h:39
const Double_t LONGITUDE_LDB
Longitude at LDB.
Stores information about a coherently summed waveform (filtered/unfiltered/deconvolved) The coherent ...
const AnalysisWaveform * getDeconvolved(AnitaPol::AnitaPol_t pol, int i, bool filtered=false) const
Definition: Analyzer.h:73
void setExtraFilters(FilterStrategy *extra)
Definition: Analyzer.cc:1251
void setThetaPhiRange(double thetaMin, double thetaMax, double phiMin, double phiMax)
Definition: Analyzer.h:117
const AnalysisWaveform * getCoherent(AnitaPol::AnitaPol_t pol, int i, bool filtered=false) const
Definition: Analyzer.h:70
Common analysis format between UCorrelator and Acclaim.
const TGraphAligned * getCoherentPower(AnitaPol::AnitaPol_t pol, int i, bool filtered=false) const
Definition: Analyzer.h:82
void setTrackSun(double setdTheta=2.5, double setdPhi=5., bool blockOut=false)
Definition: Analyzer.h:121
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.
void setExcludeThetaRange(double thetaMin, double thetaMax)
Definition: Analyzer.h:108
const AnalysisWaveform * getCoherentXpol(AnitaPol::AnitaPol_t pol, int i, bool filtered=false) const
Definition: Analyzer.h:76
This class is intended to be the main storage vessel for ANITA waveforms. It is similar in principle ...
void drawSummary(TPad *chpol=0, TPad *cvpol=0, int draw_filtered=0) const
Definition: Analyzer.cc:944
const TGraphAligned * getDeconvolvedPower(AnitaPol::AnitaPol_t pol, int i, bool filtered=false) const
Definition: Analyzer.h:85
const Double_t ALTITUDE_LDB
Altitude at LDB.