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