AnitaDataset.h
1 #ifndef ANITA_DATASET_H
2 #define ANITA_DATASET_H
3 
28 #include <vector>
29 #include "AnitaConventions.h"
30 #include "TString.h"
31 #include "TRandom3.h"
32 
33 class RawAnitaHeader;
34 class TTree;
35 class Adu5Pat;
36 class UsefulAnitaEvent;
37 class RawAnitaEvent;
39 class PrettyAnitaHk;
40 class TFile;
41 class TurfRate;
42 class TCut;
43 class TEventList;
44 class SurfHk;
45 class TruthAnitaEvent;
46 
48 {
49  public:
50 
51 
52  // All the blinding options.
53  // The plan is to have each one represented by a bit so multiple strategies can be set at the same time.
54  // Any output files should probably save the blinding strategy to the tree.
55  enum BlindingStrategy
56  {
57  kNoBlinding = 0x00,
58  kInsertedVPolEvents = 0x01,
59  kInsertedHPolEvents = 0x02,
60  kRandomizePolarity = 0x04,
61  kDefault = kNoBlinding
62  };
63 
64 
65  // various data directories that might be searched by AnitaDataset
66  // Note that ANITA_ROOT_DATA is reverted to searched if none of the others work
67  enum DataDirectory
68  {
69  ANITA_ROOT_DATA = -1,
70  ANITA_MC_DATA= 0,
71  ANITA1_ROOT_DATA=1,
72  ANITA2_ROOT_DATA=2,
73  ANITA3_ROOT_DATA=3,
74  ANITA4_ROOT_DATA=4
75  };
76 
83  static TString getDescription(BlindingStrategy strat);
84 
85 
91  BlindingStrategy setStrategy(BlindingStrategy strat);
92 
93 
98  BlindingStrategy getStrategy();
99 
105  inline TString describeCurrentStrategy(){
106  return getDescription(getStrategy());
107  }
108 
109 
119  static int getRunContainingEventNumber(UInt_t eventNumber);
120 
121 
123  static const char * getDataDir(DataDirectory dir = ANITA_ROOT_DATA);
124 
126  static int getRunAtTime(double t);
127 
138  static bool isKnownHiCalEvent(UInt_t eventNumber, Int_t anita=AnitaVersion::get());
139 
140 
153  AnitaDataset (int run, bool decimated = false, WaveCalType::WaveCalType_t cal = WaveCalType::kDefault,
154  DataDirectory dir = ANITA_ROOT_DATA , BlindingStrategy strat = AnitaDataset::kDefault);
155 
157  void setCalType(WaveCalType::WaveCalType_t cal) { fCalType = cal; }
158 
160  virtual ~AnitaDataset();
161 
162 
168  bool loadRun(int run, bool decimated = false, DataDirectory dir = ANITA_ROOT_DATA );
169 
173  int getEvent(int eventNumber, bool quiet = false);
174 
176  int getEntry(int entryNumber);
177 
178  // returns the currently selected entry
179  int current() { return fDecimated ? fDecimatedEntry : fWantedEntry; }
180 
182  int next() { return getEntry(fDecimated ? fDecimatedEntry +1 : fWantedEntry+1); }
183 
185  int previous() { return getEntry(fDecimated ? fDecimatedEntry - 1 : fWantedEntry-1); }
186 
188  int first() { return getEntry(0); }
189 
191  int last() { return getEntry(N()-1); }
192 
194  int N() const;
195 
196 
198  int firstEvent();
199 
201  int nextEvent();
202 
204  int previousEvent();
205 
207  int lastEvent();
208 
210  int nthEvent(int i);
211 
213  int nextMinBiasEvent();
214 
216  int previousMinBiasEvent();
217 
221  int setCut(const TCut & cut);
222 
224  int NInCut() const;
225 
227  int firstInCut();
228 
230  int lastInCut();
231 
233  int nextInCut();
234 
236  int previousInCut();
237 
239  int nthInCut(int i);
240 
243  int setPlaylist(const char* playlist);
244 
246  int NInPlaylist() const;
247 
249  int firstInPlaylist();
250 
252  int lastInPlaylist();
253 
255  int nextInPlaylist();
256 
258  int previousInPlaylist();
259 
261  int nthInPlaylist(int i);
262 
267  virtual UsefulAnitaEvent * useful(bool force_reload = false);
268 
269 
271  RawAnitaEvent * raw(bool force_reload = false);
272 
275  CalibratedAnitaEvent * calibrated(bool force_reload = false);
276 
277 
281  Adu5Pat * gps(bool force_reload = false);
282 
284  PrettyAnitaHk * hk(bool force_reload = false);
285 
291  virtual RawAnitaHeader * header(bool force_reload = false);
292 
293 
295  TurfRate * turf(bool force_reload = false);
296 
298  SurfHk * surf(bool force_reload = false);
299 
300 
302  TruthAnitaEvent * truth(bool force_reload = true);
303 
305  std::vector<Double_t>* RcoInfo(bool force_reload = false);
307  std::vector<Double_t>* ClockPhiInfo(bool force_reload = false);
309  std::vector<Double_t>* TempFactorInfo(bool force_reload = false);
311  Int_t ClockProblemInfo(bool force_reload = false);
313  Int_t ClockSpikeInfo(bool force_reload = false);
315  Int_t RFSpikeInfo(bool force_reload = false);
316 
319 
321  int getCurrRun() { return currRun; };
322 
323  /* Wraps the random number generator for polarity inversion so it is derministic regardless of event processing order */
324  bool maybeInvertPolarity(UInt_t eventNumber);
325 
326  /* Where was HiCal at a particular time?*/
327  static void hiCal(UInt_t unixTime, Double_t& longitude, Double_t& latitude, Double_t& altitude);
328 
329  /* Where was hical? Uses the current header realTime*/
330  void hiCal(Double_t& longitude, Double_t& latitude, Double_t& altitude);
331 
332  protected:
333  void unloadRun();
334  TTree * fHeadTree;
335  TTree * fDecimatedHeadTree; //only used when using decimated
336  Long64_t * fIndices;
337  Long64_t fIndex;
338  RawAnitaHeader * fHeader;
339  TTree *fEventTree;
340  TTree* fCalibInfoTree;
341  CalibratedAnitaEvent * fCalEvent;
342  RawAnitaEvent * fRawEvent;
343  UsefulAnitaEvent * fUseful;
344  std::vector<Double_t>* fRcoInfo;
345  std::vector<Double_t>* fClockPhiInfo;
346  std::vector<Double_t>* fTempFactorInfo;
347  Int_t fClockProblemInfo;
348  Int_t fClockSpikeInfo;
349  Int_t fRFSpikeInfo;
350  Bool_t fUsefulDirty;
351  Bool_t fCalDirty; // used only with raw data
352  Bool_t fGpsDirty; // used only with gpsFile data
353  TTree* fGpsTree;
354  Adu5Pat * fGps;
355  TTree* fHkTree;
356  PrettyAnitaHk * fHk;
357 
358  TTree * fTurfTree;
359  TurfRate * fTurf;
360  Bool_t fTurfDirty;
361 
362  TTree * fSurfTree;
363  SurfHk * fSurf;
364  Bool_t fSurfDirty;
365 
366  TTree * fTruthTree;
367  TruthAnitaEvent * fTruth;
368 
369 
370  /* place to store the current run but that can't be confusingly changed */
371  int currRun;
372 
373  Long64_t fWantedEntry;
374  Long64_t fDecimatedEntry;
375  Bool_t fHaveGpsEvent;
376  Bool_t fHaveCalibFile;
377  Bool_t fHaveCalibInfo;
379  std::vector<TFile *> filesToClose;
380  bool fDecimated;
381  TEventList * fCutList;
382  int fCutIndex;
383 
384  static void loadRunToEv(int anita); // read runToEvA*.txt
385  static void loadHiCalGps();
386  int loadPlaylist(const char* playlist);
387  int fPlaylistIndex;
388  std::vector<std::vector<long> > fPlaylist;
389  int getPlaylistRun() { return fPlaylist[fPlaylistIndex][0]; }
390  Long64_t getPlaylistEvent() { return fPlaylist[fPlaylistIndex][1]; }
391 
392 
393  /* Blinding stuff */
394  void zeroBlindPointers();
395  void loadBlindTrees();
396 
397  BlindingStrategy theStrat;
398  TRandom3 fRandy;
399 
401  Int_t needToOverwriteEvent(AnitaPol::AnitaPol_t pol, UInt_t eventNumber);
402  void overwriteHeader(RawAnitaHeader* header, AnitaPol::AnitaPol_t pol, Int_t fakeTreeEntry);
403  void overwriteEvent(UsefulAnitaEvent* useful, AnitaPol::AnitaPol_t pol, Int_t fakeTreeEntry);
404 
405  // fake things
406  TFile* fBlindFile;
409 
412 
413  // Here we have a set of three vectors that should all be the same length as the elements of each match up.
414  // They are filled in loadBlindTrees
415  std::vector<UInt_t> eventsToOverwrite;
416  std::vector<AnitaPol::AnitaPol_t> polarityOfEventToInsert;
417  std::vector<Int_t> fakeTreeEntries;
418 
419  DataDirectory datadir;
420 
421 
422 };
423 
424 // define the bitwise or operator | to combine blinding strategies
425 inline AnitaDataset::BlindingStrategy operator|(AnitaDataset::BlindingStrategy strat1, AnitaDataset::BlindingStrategy strat2){
426  return static_cast<AnitaDataset::BlindingStrategy>(static_cast<UInt_t>(strat1) | static_cast<UInt_t>(strat2));
427 }
428 // define the bitwise and opterator & to decode the blinding strategies
429 inline AnitaDataset::BlindingStrategy operator&(AnitaDataset::BlindingStrategy strat1, AnitaDataset::BlindingStrategy strat2){
430  return static_cast<AnitaDataset::BlindingStrategy>(static_cast<UInt_t>(strat1) & static_cast<UInt_t>(strat2));
431 }
432 
433 
434 
435 #endif
bool loadedBlindTrees
!< for deciding whether to do polarity flipping (eventNumber is used as seed)
Definition: AnitaDataset.h:400
int nthInCut(int i)
static int getRunAtTime(double t)
int nextMinBiasEvent()
std::vector< Double_t > * RcoInfo(bool force_reload=false)
int setCut(const TCut &cut)
int nthInPlaylist(int i)
TTree * fBlindEventTree[AnitaPol::kNotAPol]
!< Pointer to file containing tree of UsefulAnitaEvents to insert
Definition: AnitaDataset.h:407
Int_t ClockSpikeInfo(bool force_reload=false)
int previousInPlaylist()
AnitaDataset(int run, bool decimated=false, WaveCalType::WaveCalType_t cal=WaveCalType::kDefault, DataDirectory dir=ANITA_ROOT_DATA, BlindingStrategy strat=AnitaDataset::kDefault)
CalibratedAnitaEvent * calibrated(bool force_reload=false)
static void hiCal(UInt_t unixTime, Double_t &longitude, Double_t &latitude, Double_t &altitude)
virtual ~AnitaDataset()
std::vector< UInt_t > eventsToOverwrite
!< Pointer to fake header
Definition: AnitaDataset.h:415
Adu5Pat – The ADU5 Position and Attitude Data.
Definition: Adu5Pat.h:26
enum WaveCalType::EWaveCalType WaveCalType_t
The calibration enumeration type.
PrettyAnitaHk – The prettified ANITA Hk.
Definition: PrettyAnitaHk.h:22
void setCalType(WaveCalType::WaveCalType_t cal)
Definition: AnitaDataset.h:157
int getEntry(int entryNumber)
int previousMinBiasEvent()
TruthAnitaEvent * truth(bool force_reload=true)
Int_t needToOverwriteEvent(AnitaPol::AnitaPol_t pol, UInt_t eventNumber)
!< Have we loaded the tree of events to insert?
bool loadRun(int run, bool decimated=false, DataDirectory dir=ANITA_ROOT_DATA)
static bool isKnownHiCalEvent(UInt_t eventNumber, Int_t anita=AnitaVersion::get())
int getEvent(int eventNumber, bool quiet=false)
Int_t RFSpikeInfo(bool force_reload=false)
RawAnitaHeader – The Raw ANITA Event Header.
SurfHk – The raw SURF scaler+threshold data.
Definition: SurfHk.h:24
USeful in for loops.
int N() const
What you should call for analysis work.
virtual UsefulAnitaEvent * useful(bool force_reload=false)
SurfHk * surf(bool force_reload=false)
UsefulAnitaEvent – The Calibrated Useful Anita Event object.
std::vector< Double_t > * TempFactorInfo(bool force_reload=false)
int NInCut() const
int nthEvent(int i)
RawAnitaEvent * raw(bool force_reload=false)
CalibratedAnitaEvent – The Calibrated Calibrated Anita Event object.
TTree * fBlindHeadTree[AnitaPol::kNotAPol]
!< Tree of UsefulAnitaEvents to insert
Definition: AnitaDataset.h:408
TruthAnitaEvent – The Truth ANITA Event.
UsefulAnitaEvent * fBlindEvent[AnitaPol::kNotAPol]
!< Tree of headers to insert
Definition: AnitaDataset.h:410
static TString getDescription(BlindingStrategy strat)
int NInPlaylist() const
BlindingStrategy getStrategy()
Int_t ClockProblemInfo(bool force_reload=false)
static const char * getDataDir(DataDirectory dir=ANITA_ROOT_DATA)
TurfRate – The Turf Rate data.
Definition: TurfRate.h:21
Adu5Pat * gps(bool force_reload=false)
RawAnitaHeader * fBlindHeader[AnitaPol::kNotAPol]
!< Pointer to fake UsefulAnitaEvent
Definition: AnitaDataset.h:411
TRandom3 fRandy
!< Currently selected blinding strategy
Definition: AnitaDataset.h:398
PrettyAnitaHk * hk(bool force_reload=false)
int getCurrRun()
Definition: AnitaDataset.h:321
virtual RawAnitaHeader * header(bool force_reload=false)
int loadPlaylist(const char *playlist)
Where was HiCal?
std::vector< Double_t > * ClockPhiInfo(bool force_reload=false)
TString describeCurrentStrategy()
Definition: AnitaDataset.h:105
BlindingStrategy setStrategy(BlindingStrategy strat)
enum AnitaPol::EAnitaPol AnitaPol_t
Polarisation enumeration.
static int getRunContainingEventNumber(UInt_t eventNumber)
Get the run that contains the eventNumber.
int setPlaylist(const char *playlist)
RawAnitaEvent – The Raw ANITA Event Data.
Definition: RawAnitaEvent.h:22
TurfRate * turf(bool force_reload=false)