ANITA Analysis Framework
AnitaEventSummary.h
1 #ifndef _ANITA_EVENT_SUMMARY_H_
2 #define _ANITA_EVENT_SUMMARY_H_
3 
4 #include "TObject.h"
5 #include "AnitaConventions.h"
6 #include "Adu5Pat.h"
7 #include <iostream>
8 #include "TVector3.h"
9 
10 class UsefulAdu5Pat;
11 class RawAnitaHeader;
12 class TruthAnitaEvent;
13 
32 class AnitaEventSummary : public TObject
33 {
34  public:
35 
36 
37  //------------------------------------------------------------------------------------
38  /*************************************************************************************
39  * Static members (for array sizes inside class)
40  *************************************************************************************/
41  static const Int_t maxDirectionsPerPol = 5;
42  static const Int_t peaksPerSpectrum = 3;
43  static const Int_t numFracPowerWindows = 5;
44  static const Int_t numBlastPowerBands = 3;
45 
46  //------------------------------------------------------------------------------------
47  /*************************************************************************************
48  * Sub-classes
49  *************************************************************************************/
50 
51 
56  class SourceHypothesis; // Forward class declarations for PointingHypothesis utils
57  class PayloadLocation; // Forward class declarations for PointingHypothesis utils
58 
60  {
61  public:
62  PointingHypothesis() : fContainer(NULL) { ; }
63  Double_t phi;
64  Double_t theta;
65  Double_t value;
66  Double_t snr;
67  Double_t mapRMS;
68  Double_t mapHistoryVal;
69  Double_t hwAngle;
70  Double_t hwAngleXPol;
71  Double_t latitude;
72  Double_t longitude;
73  Double_t altitude;
74  Double_t distanceToSource;
75  Double_t sigma_theta;
76  Double_t sigma_phi;
77  Double_t rho;
78  Double_t chisq;
80  Double_t phi_separation;
81  Double_t dphi_rough;
82  Double_t dtheta_rough;
83  Bool_t triggered;
84  Bool_t triggered_xpol;
85  Bool_t masked;
86  Bool_t masked_xpol;
87  Double_t antennaPeakAverage;
88 
89  // Most basic resolution utility functions in payload coordinates relative to ADU5-aft-fore line
90  Double_t dPhi(Double_t phi) const;
91  Double_t dTheta(Double_t theta, bool different_sign_conventions = false) const;
92 
93  // Peak direction relative to north (N->E->S->W)
94  Double_t bearing() const;
95  Double_t dPhiNorth() const;
96 
97  // Find angle from stored source hypotheses
98  Double_t dPhiWais() const;
99  Double_t dThetaWais(bool different_sign_conventions = true) const;
100  Double_t dPhiSun() const;
101  Double_t dThetaSun(bool different_sign_conventions = true) const;
102  Double_t dPhiLDB() const;
103  Double_t dThetaLDB(bool different_sign_conventions = true) const;
104  Double_t dPhiMC() const;
105  Double_t dThetaMC(bool different_sign_conventions = true) const;
106  Double_t dPhiTagged() const;
107  Double_t dThetaTagged(bool different_sign_conventions = true) const;
108 
109  // Are you within this theta/phi of a stored hypothesis?
110  Bool_t closeToMC(Double_t deltaPhiDeg, Double_t deltaThetaDeg) const;
111  Bool_t closeToWais(Double_t deltaPhiDeg, Double_t deltaThetaDeg) const;
112  Bool_t closeToLDB(Double_t deltaPhiDeg, Double_t deltaThetaDeg) const;
113  Bool_t closeToSun(Double_t deltaPhiDeg, Double_t deltaThetaDeg) const;
114  Bool_t closeToTagged(Double_t deltaPhiDeg, Double_t deltaThetaDeg) const;
115 
116  // peak direction relative to trigger information
117  Double_t minAbsHwAngle() const;
118  Bool_t absHwAngleLessThanAbsHwAngleXPol() const;
119 
120  private:
121  //----------------------------------------------------------------------------------------------------
122  // WARNING! This private nonsense is fragile, and a bit hacky.
123  // The //! comment after the fContainer member means it does not persist in ROOT.
124  // Please do not edit that comment.
125  //----------------------------------------------------------------------------------------------------
126  mutable AnitaEventSummary* fContainer;
127  const AnitaEventSummary* getContainer(const char* funcName) const;
128  Double_t dPhiSource(const SourceHypothesis& source) const; // Won't work inside TTree::Draw due to limitations in TTreeFormula, so are private, use e.g. dPhiWais() instead.
129  Double_t dThetaSource(const SourceHypothesis& source, bool different_sign_conventions = true) const; // Won't work inside TTree::Draw due to limitations in TTreeFormula, so are private, use e.g. dThetaWais() instead.
130  void printEvent() const;
131  friend class AnitaEventSummary;
132 
133  ClassDefNV(PointingHypothesis,15);
134  };
135 
143  {
144 
145  public:
146  WaveformInfo() : fContainer(NULL), fLastEventNumberCache(0), nwMeanCache(-1),
147  nwGradCache(-1), nwInterceptCache(-1), nwChisquareCache(-1) {; }
148  Double_t snr; //[0,100,16] /// Signal to Noise of waveform
149  Double_t peakHilbert;//[0,4096,21] /// peak of hilbert envelope
150  Double_t peakVal; //[0,4096,21] /// peak value
151  Double_t xPolPeakVal; //[0,4096,21] /// Peak of xpol trace
152  Double_t xPolPeakHilbert; //[0,4096,21] /// Peak of xpol hilbert Envelope
153 
154  Double_t I,Q,U,V;
155  Double_t max_dI,max_dQ,max_dU,max_dV,polErr;
157 
158  Double_t totalPower;
159  Double_t totalPowerXpol;
160 
161  //spectrum info
162  Double_t bandwidth[peaksPerSpectrum]; //[0,2,16] /// bandwidth of each peak (implementation defined, may not be comparable between analyses)
163  Double_t peakFrequency[peaksPerSpectrum]; //[0,2,16] /// peak frequency of power spectrum
164  Double_t peakPower[peaksPerSpectrum]; //power within +/- bandwidth of each peak
165  Double_t spectrumSlope; // [-100,100,16] /// Slope of line fit to spectrum (in log-space, so this is spectral-index)
166  Double_t spectrumIntercept; // [-200,200,16] /// Intercept of line fit to spectrum (in log-space)
167 
168  //Shape parameters, computed using hilbert envelope
169  // This should probably taken out into its own class
170  Double_t riseTime_10_90; //[0,128,16] /// Rise time of hilbert env from 10% to 90% of peak
171  Double_t riseTime_10_50; //[0,128,16] /// Rise time of hilbert env from 10% to 50% of peak
172  Double_t fallTime_90_10; //[0,128,16]/// Fall time of hilbert env from 90% to 10% of peak
173  Double_t fallTime_50_10; //[0,128,16] /// Fall time of hilbert env from 50% to 10% of peak
174  Double_t width_50_50; //[0,128,16] /// Width from first envelope crossing of 50 percent of peak to last
175  Double_t width_10_10; //[0,128,16]/// Width from first envelope crossing of 10 percent of peak to last
176  Double_t power_10_10;
177  Double_t power_50_50;
178  Double_t peakTime; //[-128,384,18] /// Time that peak hilbert env occurs
179  Double_t peakMoments[5];
180 
181  Double_t impulsivityMeasure; //[-1,1, 16] /// A number that has something to do with how impulsive it is
182  Double_t bandwidthMeasure;
183  Double_t fracPowerWindowBegins[numFracPowerWindows]; //[0,128,16] /// Narrowest width containing {10%, 20%, 30%, 40%, 50%} of the total power
184  Double_t fracPowerWindowEnds[numFracPowerWindows]; //[0,128,16] /// Narrowest width containing {10%, 20%, 30%, 40%, 50%} of the total power
185 
186  Int_t numAntennasInCoherent;
187 
188  Double_t localMaxToMin; //[0,4096,21] /// Largest value of local max to neighbouring local min (see Acclaim::RootTools::getLocalMaxToMin)
189  Double_t localMaxToMinTime; //[0,100,16] /// Time between local maxima and minima +ve means max is before min, -ve means min is before max
190  Double_t globalMaxToMin; //[0,4096,21] /// Difference between maximum and minimum voltage
191  Double_t globalMaxToMinTime; //[0,128,16] /// Time between maximum and minimum volts, +ve means max is before min, -ve means min is before max
192 
193 
194  //some utilities for polarization info
195  Double_t linearPolFrac() const;
196  Double_t linearPolAngle() const;
197  Double_t circPolFrac() const;
198  Double_t totalPolFrac() const;
199 
200  //utilities for instantaneous polarization info
201  Double_t instantaneousLinearPolFrac() const;
202  Double_t instantaneousLinearPolAngle() const;
203  Double_t instantaneousCircPolFrac() const;
204  Double_t instantaneousTotalPolFrac() const;
205 
206  Double_t standardizedPeakMoment(Int_t i) const;
207  inline Double_t skewness(){return standardizedPeakMoment(3);}
208  inline Double_t kurtosis(){return standardizedPeakMoment(4);}
209 
210  Double_t fracPowerWindowMean() const;
211  Double_t fracPowerWindowGradient() const;
212  Double_t fracPowerWindowIntercept() const;
213  Double_t fracPowerWindowChisquare() const;
214 
215  ClassDefNV(WaveformInfo, 18);
216 
217  private:
218  friend class AnitaEventSummary;
219  void cacheQuantitiesDerivedFromNarrowestWidths() const;
220  mutable AnitaEventSummary* fContainer;
221  mutable UInt_t fLastEventNumberCache;
222  mutable Double_t nwMeanCache;
223  mutable Double_t nwGradCache;
224  mutable Double_t nwInterceptCache;
225  mutable Double_t nwChisquareCache;
226  };
227 
233  {
234 
235  public:
237  ChannelInfo() : pol(AnitaPol::kNotAPol), ant(-1) {; }
238 
239  Double_t rms; //[0,1024,20]
240  Double_t avgPower;
241  Double_t snr; //[0,100,16]/// Signal to Noise of waveform
242  Double_t peakHilbert;//[0,4096,21] /// peak of hilbert envelope
243 
244  Double_t getPhi() const;
245  inline Int_t getAnt() const {return ant;} // could add some errors on -1 here...
246  inline AnitaPol::AnitaPol_t getPol() const {return pol;} // could add some errors on 2 here...
247 
248  private:
249  friend class AnitaEventSummary;
250  // WARNING! THE COMMENTS TRAILING THESE AFFECT WHAT HAPPENS IN ROOT!
251  // The //! comment initializer means ROOT does not store these variables when writing to files.
252  // These variables are filled in the AnitaEventSummary constructor, so each ChannelInfo knows
253  // its location in the ChannelInfo channels[2][48] array in AnitaEventSummary
254  AnitaPol::AnitaPol_t pol;
255  Int_t ant;
256 
257  ClassDefNV(ChannelInfo, 5);
258  };
259 
260 
266  {
267  public:
268  EventFlags() {; }
271  {
272  NONE,
273  WAIS, // is actually Hpol wais in both A3 and A4
274  LDB,
275  SIPLE,
276  TD,
277  WAIS_V, // the Vpol wais in A4
278  HICAL // HICAL2 flag
279  };
280 
281  Int_t isGood;
282  Int_t isRF;
283  Int_t isAdu5Trigger;
284  Int_t isG12Trigger;
285  Int_t isSoftwareTrigger;
286  Int_t isMinBiasTrigger;
287  Int_t isPayloadBlast;
288  Int_t nadirFlag;
289  Int_t strongCWFlag;
290  Int_t isHPolTrigger;
291  Int_t isVPolTrigger;
292  Int_t isStepFunction;
293  Int_t hasGlitch;
294 
295  CalPulser pulser;
296  Bool_t isVarner;
297  Bool_t isVarner2;
298 
303  Double_t meanPower[1+AnitaRing::kNotARing];
304  Double_t medianPower[1+AnitaRing::kNotARing];
305  Double_t meanPowerFiltered[1+AnitaRing::kNotARing];
306  Double_t medianPowerFiltered[1+AnitaRing::kNotARing];
307 
308  Double_t maxBottomToTopRatio[AnitaPol::kNotAPol];
309  Int_t maxBottomToTopRatioSector[AnitaPol::kNotAPol];
310  Double_t minBottomToTopRatio[AnitaPol::kNotAPol];
311  Int_t minBottomToTopRatioSector[AnitaPol::kNotAPol];
312 
313  Int_t nSectorsWhereBottomExceedsTop;
314 
316  Double_t blastFraction;
317 
318  Double_t middleOrBottomPower[numBlastPowerBands];
319  Double_t topPower[numBlastPowerBands];
320  Int_t middleOrBottomAnt[numBlastPowerBands];
321  Int_t middleOrBottomPol[numBlastPowerBands];
322 
323  ClassDefNV(EventFlags,15);
324  };
325 
326 
327 
328 
329 
335  {
336  public:
337  SourceHypothesis() { reset(); }
338  Double_t theta;
339  Double_t phi;
340  Double_t distance;
341 
342  Double_t mapValue[NUM_POLS];
343  Double_t mapHistoryVal[NUM_POLS];
344 
345  void reset();
346  virtual ~SourceHypothesis() { ; }
347 
348  ClassDef(SourceHypothesis,4);
349  };
350 
351 
352 
353 
358  class MCTruth : public SourceHypothesis
359  {
360  public:
361  MCTruth() { reset(); }
362  WaveformInfo wf[AnitaPol::kNotAPol];
363  Double_t weight;
364  Double_t energy;
365  Double_t triggerSNR[2];
366  TVector3 nuDirection; //in earth centered fixed coordinates (i.e. Cartesian coordinates of neutrino far away)
367  double nuTheta, nuPhi; //neutrino (not RF!) direction in payload coordinates -> ray at infinity
368  double interactionTheta, interactionPhi; //neutrino interaction angle to ANITA in payload coordinates
369  double interactionLon, interactionLat, interactionAlt;
370  double rfExitTheta, rfExitPhi; //rf exit
371  double rfExitLon, rfExitLat, rfExitAlt;
372  double balloonLon, balloonLat, balloonAlt; //balloon co-ords
373  void reset();
374  virtual ~MCTruth() { ; }
375 
376  ClassDefNV(MCTruth,9);
377  };
378 
379 
380 
381 
387  {
388  public:
389  PayloadLocation() { reset(); }
390  PayloadLocation(const Adu5Pat* pat);
391 
392  Float_t latitude;
393  Float_t longitude;
394  Float_t altitude;
395  Float_t heading;
396 
397  Float_t prevHeading; //useful for determining rotation rate
398 
399  void reset() { latitude = -999; longitude = -999; altitude = -999; heading = -999; prevHeading = -999;};
400  void update(const Adu5Pat* pat);
401 
406  Adu5Pat pat () const {
407  Adu5Pat pat;
408  pat.longitude = longitude;
409  pat.latitude = latitude;
410  pat.altitude = altitude;
411  pat.heading = heading;
412  return pat;
413  }
414 
415  ClassDefNV(PayloadLocation,2);
416  };
417 
418 
419 
420 
421 
422 
423  //------------------------------------------------------------------------------------
424  /*************************************************************************************
425  * Public member variables
426  *************************************************************************************/
427  Int_t run;
428  UInt_t eventNumber;
429  UInt_t realTime;
430  Int_t nPeaks[AnitaPol::kNotAPol];
431  PointingHypothesis peak[AnitaPol::kNotAPol][maxDirectionsPerPol];
432  WaveformInfo coherent[AnitaPol::kNotAPol][maxDirectionsPerPol];
433  WaveformInfo deconvolved[AnitaPol::kNotAPol][maxDirectionsPerPol];
434  WaveformInfo coherent_filtered[AnitaPol::kNotAPol][maxDirectionsPerPol];
435  WaveformInfo deconvolved_filtered[AnitaPol::kNotAPol][maxDirectionsPerPol];
436  ChannelInfo channels[AnitaPol::kNotAPol][NUM_SEAVEYS];
443 
444 
445  //------------------------------------------------------------------------------------
446  /*************************************************************************************
447  * Public member functions
448  *************************************************************************************/
449  // See source file for doxygen comments
451  AnitaEventSummary(const RawAnitaHeader* header);
452  AnitaEventSummary(const RawAnitaHeader* header, UsefulAdu5Pat* pat, const TruthAnitaEvent * truth = 0);
453  void setTriggerInfomation(const RawAnitaHeader* header);
454  void setSourceInformation(UsefulAdu5Pat* pat, const TruthAnitaEvent * truth = 0);
455  void zeroInternals();
456  Bool_t update() const {resetNonPersistent(); return true;}
457 
458  // Utilities to find interesting entries in AnitaEventSummary
459  AnitaPol::AnitaPol_t highestPol() const;
460  Int_t highestPolAsInt() const;
461  Int_t highestPeakInd() const;
462  const PointingHypothesis& highestPeak() const;
463  const WaveformInfo& highestCoherent() const;
464  const WaveformInfo& highestDeconvolved() const;
465  const WaveformInfo& highestCoherentFiltered() const;
466  const WaveformInfo& highestDeconvolvedFiltered() const;
467 
468  static void setThresholdForMostImpulsive(double threshold);
469  AnitaPol::AnitaPol_t mostImpulsivePol(int whichMetric=0) const;
470  Int_t mostImpulsivePolAsInt(int whichMetric=0) const;
471  Int_t mostImpulsiveInd(int whichMetric=0) const;
472  const PointingHypothesis& mostImpulsivePeak(int whichMetric=0) const;
473  const WaveformInfo& mostImpulsiveCoherent(int whichMetric=0) const;
474  const WaveformInfo& mostImpulsiveDeconvolved(int whichMetric=0) const;
475  const WaveformInfo& mostImpulsiveCoherentFiltered(int whichMetric=0) const;
476  const WaveformInfo& mostImpulsiveDeconvolvedFiltered(int whichMetric=0) const;
477 
478  inline Double_t weight(){return mc.weight > 0 ? mc.weight : 1;}
479  AnitaPol::AnitaPol_t trainingPol() const;
480  Int_t trainingPolAsInt() const;
481  Int_t trainingPeakInd() const;
482  const PointingHypothesis& trainingPeak() const;
483  const WaveformInfo& trainingCoherent() const;
484  const WaveformInfo& trainingDeconvolved() const;
485  const WaveformInfo& trainingCoherentFiltered() const;
486  const WaveformInfo& trainingDeconvolvedFiltered() const;
487  Int_t countChannelAboveThreshold(int threshold=100) const;
488 
489  //------------------------------------------------------------------------------------
490  private:
491 
492  /*************************************************************************************
493  * Private member variables/functions
494  *************************************************************************************/
495  // WARNING! THE COMMENTS TRAILING THESE AFFECT WHAT HAPPENS IN ROOT!
496  // The //! comment initializer means ROOT does not store these variables when writing to files.
497  // We want to keep it this way, since they are only used to cache the results of the utility
498  // funtions to find the highest peak and peak nearest the monte carlo truth info.
499  mutable Int_t fHighestPeakIndex;
500  mutable AnitaPol::AnitaPol_t fHighestPol;
501  mutable Int_t fTrainingPeakIndex;
502  mutable AnitaPol::AnitaPol_t fTrainingPol;
503  mutable Int_t fMostImpulsiveIndex;
504  mutable AnitaPol::AnitaPol_t fMostImpulsivePol;
505  mutable UInt_t fLastEventNumber;
506 
507 
508  void findHighestPeak() const;
509  void findTrainingPeak() const;
510  void findMostImpulsive(int whichMetric) const;
511  void resetNonPersistent() const;
512  const SourceHypothesis* sourceFromTag() const;
513 
514  ClassDefNV(AnitaEventSummary, 41);
515 };
516 
517 #endif
WaveformInfo deconvolved[AnitaPol::kNotAPol][maxDirectionsPerPol]
Summaries of the (unfiltered) coherently summed waveforms, array index correponds to entry in peak[][...
Bool_t masked_xpol
was this in a masked phi sector?
Double_t totalPower
The number of points used in the above estimates.
Double_t hwAngleXPol
angle with respect to triggering phi sector
PointingHypothesis peak[AnitaPol::kNotAPol][maxDirectionsPerPol]
Number of peaks stored in this AnitaEventSummary (might be less than maxDirectionsPerPol) ...
Double_t mapHistoryVal[NUM_POLS]
what the instantaneous map value is at this source hypothesis
Double_t dphi_rough
angular separation from higher value peak in same event. 1000 if highest value event (i...
Double_t dThetaTagged(bool different_sign_conventions=true) const
See AnitaEventSummary::sourceFromTag()
Double_t phi_separation
If an event barely missed the ground, it is useful to see the coordinates at which it would hit if th...
Adu5Pat pat() const
Copy the data from the pat into the object.
Bool_t closeToMC(Double_t deltaPhiDeg, Double_t deltaThetaDeg) const
See AnitaEventSummary::sourceFromTag()
Double_t hwAngle
value of average of the peak location over the past 60 min-bias events
MCTruth mc
Contains location of LDB cal pulser in map coordinates at time of event.
virtual ~SourceHypothesis()
sets all the values to nonsense. Sorry, mapHistoryVal means this is in source now ...
Double_t sigma_theta
on continent, or -9999 if doesn&#39;t intersect
Double_t meanPower[1+AnitaRing::kNotARing]
Double_t impulsivityMeasure
moments about Peak (1st - 5th moments)
SourceHypothesis wais
Contains location of sun in map coordinates at time of event.
AnitaPol::AnitaPol_t mostImpulsivePol(int whichMetric=0) const
value between 0 and 1, will find the brightest peak that is within threshold as impulsive as the most...
static const Int_t numFracPowerWindows
The maximum number of frequency peaks per waveform spectrum.
static const Int_t peaksPerSpectrum
The maximum number of hypotheses storable per polarization */.
WaveformInfo coherent_filtered[AnitaPol::kNotAPol][maxDirectionsPerPol]
Summaries of the (unfiltered) de-dispersed coherently summed waveforms, array index correponds to ent...
Float_t latitude
Slightly more useful constructor.
Double_t max_dI
Integral Stokes Parameters (over the entire waveform)
Bool_t triggered
theta - theta rough
ChannelInfo()
Correct indices are set in the AnitaEventSummary constructor.
Double_t totalPowerXpol
Total power in waveform.
EventFlags flags
Summaries of each channel&#39;s waveform.
Bool_t masked
was this in a triggered xpol phi sector?
Double_t mapHistoryVal
rms of interferometric map
Double_t antennaPeakAverage
was this in a masked phi xpol sector?
Double_t power_50_50
Power enclosed within 10_10 width.
Double_t dPhi(Double_t phi) const
the average of channel peaks in this direction
AnitaEventSummary()
Reduced GPS data.
Double_t theta_adjustment_needed
chisq/ndof of peak finding process, if available (otherwise zero)
Double_t chisq
correlation coefficient between theta and phi
AnitaPol::AnitaPol_t trainingPol() const
Return the weight of the event, always returns 1 for data, the weight from MCTruth otherwise...
Double_t localMaxToMin
Number of antennas used to make this.
Int_t NPointsMaxStokes
instantanteous stokes parameters (computed near max_dI).
SourceHypothesis ldb
Contains location of WAIS divide cal pulser in map coordinates at time of event.
PayloadLocation anitaLocation
Contains summary information about MC truth, if real data then this filled with constant, unphysical values.
Double_t longitude
on continent, or -9999 if doesn&#39;t intersect
Stores information about a coherently summed waveform (filtered/unfiltered/deconvolved) The coherent ...
UInt_t realTime
Event number.
Double_t theta
peak phi, degrees
ChannelInfo channels[AnitaPol::kNotAPol][NUM_SEAVEYS]
Summaries of the filtered, de-dispersed, coherently summed waveforms, array index correponds to entry...
Double_t dtheta_rough
phi - phi rough
Double_t value
peak theta, degrees
Common analysis format between UCorrelator and Acclaim.
Double_t peakTime
Power enclosed within 50_50 width.
SourceHypothesis sun
Flags corresponding the event quality, trigger type, calibration pulser timing, etc.
Int_t nPeaks[AnitaPol::kNotAPol]
Time of the event.
WaveformInfo coherent[AnitaPol::kNotAPol][maxDirectionsPerPol]
Summaries of the event peak directions (indices of all WaveformInfo member arrays match peak index) ...
Bool_t triggered_xpol
was this in a triggered phi sector?
Double_t latitude
angle with respect to triggering phi sector in opposite polarisation
Double_t distanceToSource
on continent, or -9999 if doesn&#39;t intersect
WaveformInfo deconvolved_filtered[AnitaPol::kNotAPol][maxDirectionsPerPol]
Summaries of the filtered, coherently summed waveforms, array index correponds to entry in peak[][]...
void reset()
a history of the interferometric map value for the source location
Double_t altitude
on continent, or -9999 if doesn&#39;t intersect