GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
DirectoryHelper.cxx
Go to the documentation of this file.
1#include "DirectoryHelper.hh"
2
3void DirectoryHelper::CreateHistograms(unsigned int slot)
4{
5 // some variables to easily change range and binning for multiple histograms at once
6 int energyBins = 10000;
7 double lowEnergy = 0.;
8 double highEnergy = 2000.;
9
10 // unsuppressed spectra
11 fH1[slot]["unsuppressed/singles/griffinE"] = new TH1F("griffinE", Form("Unsuppressed griffin energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
12 // suppressed spectra
13 fH1[slot]["suppressed/singles/griffinESupp"] = new TH1F("griffinESupp", Form("Suppressed griffin energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
14 // unsuppressed addback spectra
15 fH1[slot]["unsuppressed/addback/griffinEAddback"] = new TH1F("griffinEAddback", Form("Unsuppressed griffin addback energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
16 // suppressed addback spectra
17 fH1[slot]["suppressed/addback/griffinESuppAddback"] = new TH1F("griffinESuppAddback", Form("Suppressed griffin addback energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
18
19 // initialize the two arrays to keep time
20 fLastTS[slot].resize(64, 0.);
21 fLastSuppressedTS[slot].resize(64, 0.);
22 fLastTime[slot].resize(64, 0.);
23 fLastSuppressedTime[slot].resize(64, 0.);
24 fLastTSNoPileup[slot].resize(64, 0.);
25 fLastSuppressedTSNoPileup[slot].resize(64, 0.);
26 fLastTimeNoPileup[slot].resize(64, 0.);
27 fLastSuppressedTimeNoPileup[slot].resize(64, 0.);
28
29 // timing spectra
30 fH2[slot]["griffinDeadTS"] = new TH2F("griffinDeadTS", "timestamp difference between consecutive hits in a griffin channel", 200, 0., 2000., 64, 0.5, 64.5);
31 fH2[slot]["griffinSuppressedDeadTS"] = new TH2F("griffinSuppressedDeadTS", "timestamp difference between consecutive suppressed hits in a griffin channel", 200, 0., 2000., 64, 0.5, 64.5);
32 fH2[slot]["griffinDeadTime"] = new TH2F("griffinDeadTime", "time difference between consecutive hits in a griffin channel", 2000, 0., 2000., 64, 0.5, 64.5);
33 fH2[slot]["griffinSuppressedDeadTime"] = new TH2F("griffinSuppressedDeadTime", "time difference between consecutive suppressed hits in a griffin channel", 2000, 0., 2000., 64, 0.5, 64.5);
34 fH2[slot]["griffinDeadTSNoPileup"] = new TH2F("griffinDeadTSNoPileup", "timestamp difference between consecutive hits in a griffin channel w/o pileups", 200, 0., 2000., 64, 0.5, 64.5);
35 fH2[slot]["griffinSuppressedDeadTSNoPileup"] = new TH2F("griffinSuppressedDeadTSNoPileup", "timestamp difference between consecutive suppressed hits in a griffin channel w/o pileups", 200, 0., 2000., 64, 0.5, 64.5);
36 fH2[slot]["griffinDeadTimeNoPileup"] = new TH2F("griffinDeadTimeNoPileup", "time difference between consecutive hits in a griffin channel w/o pileups", 2000, 0., 2000., 64, 0.5, 64.5);
37 fH2[slot]["griffinSuppressedDeadTimeNoPileup"] = new TH2F("griffinSuppressedDeadTimeNoPileup", "time difference between consecutive suppressed hits in a griffin channel w/o pileups", 2000, 0., 2000., 64, 0.5, 64.5);
38}
39
40// Coincidences Gates
42{ // Griffin-Griffin
43 return -250. <= h2->GetTime() - h1->GetTime() && h2->GetTime() - h1->GetTime() <= 250.;
44}
46{
47 return (-500. <= h1->GetTime() - h2->GetTime() && h1->GetTime() - h2->GetTime() <= -250.) || (250. <= h1->GetTime() - h2->GetTime() && h1->GetTime() - h2->GetTime() <= 500.);
48}
49
50void DirectoryHelper::Exec(unsigned int slot, TGriffin& grif, TGriffinBgo& grifBgo)
51{
52 // we use .at() here instead of [] so that we get meaningful error message if a histogram we try to fill wasn't created
53 // e.g. because of a typo
54
55 // loop over unsuppressed griffin hits
56 for(int g = 0; g < grif.GetMultiplicity(); ++g) {
57 auto grif1 = grif.GetGriffinHit(g);
58 fH1[slot].at("unsuppressed/singles/griffinE")->Fill(grif1->GetEnergy());
59 if(grif1->GetArrayNumber() <= 64) {
60 if(fLastTS[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinDeadTS")->Fill(grif1->GetTimeStampNs() - fLastTS[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
61 fLastTS[slot][grif1->GetArrayNumber()] = grif1->GetTimeStampNs();
62 if(fLastTime[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinDeadTime")->Fill(grif1->GetTime() - fLastTime[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
63 fLastTime[slot][grif1->GetArrayNumber()] = grif1->GetTime();
64 if(grif1->GetKValue() != 379) {
65 if(fLastTSNoPileup[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinDeadTSNoPileup")->Fill(grif1->GetTimeStampNs() - fLastTSNoPileup[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
66 fLastTSNoPileup[slot][grif1->GetArrayNumber()] = grif1->GetTimeStampNs();
67 if(fLastTimeNoPileup[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinDeadTimeNoPileup")->Fill(grif1->GetTime() - fLastTimeNoPileup[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
68 fLastTimeNoPileup[slot][grif1->GetArrayNumber()] = grif1->GetTime();
69 }
70 }
71 }
72
73 // loop over suppressed griffin hits
74 for(int g = 0; g < grif.GetSuppressedMultiplicity(&grifBgo); ++g) {
75 auto grif1 = grif.GetSuppressedHit(g);
76 fH1[slot].at("suppressed/singles/griffinESupp")->Fill(grif1->GetEnergy());
77 if(grif1->GetArrayNumber() <= 64) {
78 if(fLastSuppressedTS[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinSuppressedDeadTS")->Fill(grif1->GetTimeStampNs() - fLastSuppressedTS[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
79 fLastSuppressedTS[slot][grif1->GetArrayNumber()] = grif1->GetTimeStampNs();
80 if(fLastSuppressedTime[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinSuppressedDeadTime")->Fill(grif1->GetTime() - fLastSuppressedTime[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
81 fLastSuppressedTime[slot][grif1->GetArrayNumber()] = grif1->GetTime();
82 if(grif1->GetKValue() != 379) {
83 if(fLastSuppressedTSNoPileup[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinSuppressedDeadTSNoPileup")->Fill(grif1->GetTimeStampNs() - fLastSuppressedTSNoPileup[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
84 fLastSuppressedTSNoPileup[slot][grif1->GetArrayNumber()] = grif1->GetTimeStampNs();
85 if(fLastSuppressedTimeNoPileup[slot][grif1->GetArrayNumber()] != 0) fH2[slot].at("griffinSuppressedDeadTimeNoPileup")->Fill(grif1->GetTime() - fLastSuppressedTimeNoPileup[slot][grif1->GetArrayNumber()], grif1->GetArrayNumber());
86 fLastSuppressedTimeNoPileup[slot][grif1->GetArrayNumber()] = grif1->GetTime();
87 }
88 }
89 }
90
91 // loop over unsuppressed griffin addback hits
92 for(int g = 0; g < grif.GetAddbackMultiplicity(); ++g) {
93 auto grif1 = grif.GetAddbackHit(g);
94 fH1[slot].at("unsuppressed/addback/griffinEAddback")->Fill(grif1->GetEnergy());
95 }
96
97 // loop over suppressed griffin addback hits
98 for(int g = 0; g < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++g) {
99 auto grif1 = grif.GetSuppressedAddbackHit(g);
100 fH1[slot].at("suppressed/addback/griffinESuppAddback")->Fill(grif1->GetEnergy());
101 }
102}
103
104void DirectoryHelper::EndOfSort(std::shared_ptr<std::map<std::string, TList>>& list)
105{
106}
bool TimeRandom(TGriffinHit *h1, TGriffinHit *h2)
bool PromptCoincidence(TGriffinHit *h1, TGriffinHit *h2)
std::map< std::vector< double > > fLastSuppressedTime
std::map< std::vector< double > > fLastSuppressedTS
void Exec(unsigned int slot, TGriffin &grif, TGriffinBgo &grifBgo)
std::map< std::vector< double > > fLastTime
void EndOfSort(std::shared_ptr< std::map< std::string, TList > > &list) override
This method gets called at the end of Finalize()
std::map< std::vector< double > > fLastTSNoPileup
std::map< std::vector< double > > fLastSuppressedTimeNoPileup
void CreateHistograms(unsigned int slot) override
Virtual helper function that the user uses to create their histograms.
std::map< std::vector< double > > fLastTimeNoPileup
std::map< std::vector< double > > fLastTS
std::map< std::vector< double > > fLastSuppressedTSNoPileup
virtual Double_t GetTime(const ETimeFlag &correct_flag=ETimeFlag::kAll, Option_t *opt="") const
Returns a time value to the nearest nanosecond!
std::vector< TGRSIMap< std::string, TH2 * > > fH2
Definition TGRSIHelper.h:49
std::vector< TGRSIMap< std::string, TH1 * > > fH1
Definition TGRSIHelper.h:48
TGriffinHit * GetAddbackHit(const int &i)
Definition TGriffin.h:87
Short_t GetAddbackMultiplicity()
Definition TGriffin.h:84
TGriffinHit * GetSuppressedAddbackHit(const int &i)
Definition TGriffin.h:122
Short_t GetSuppressedMultiplicity(const TBgo *bgo)
Definition TGriffin.h:111
Short_t GetSuppressedAddbackMultiplicity(const TBgo *bgo)
Definition TGriffin.h:119
Short_t GetMultiplicity() const override
Definition TGriffin.h:52
TGriffinHit * GetSuppressedHit(const int &i)
!
Definition TGriffin.h:108
TGriffinHit * GetGriffinHit(const int &i)
!
Definition TGriffin.h:47