GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
ExampleEventHelper.cxx
Go to the documentation of this file.
2
4{
5 // try and get the cycle length if we have a PPG provided
6 // only necessary for the first worker, this is shared with all other workers
7 if(Ppg() != nullptr) {
8 // the ODB cycle length is in microseconds!
10 if(slot == 0) {
11 std::stringstream str;
12 str << "Got ODB cycle length " << fCycleLength << " us = " << fCycleLength / 1e6 << " s" << std::endl;
13 std::cout << str.str();
14 }
15 } else if(slot == 0) {
16 std::stringstream str;
17 str << DRED << "No ppg provided, can't fill cycle spectra!" << RESET_COLOR << std::endl;
18 std::cout << str.str();
19 }
20
21 // some variables to easily change range and binning for multiple histograms at once
22 int energyBins = 10000;
23 double lowEnergy = 0.;
24 double highEnergy = 2000.;
25
26 fH2[slot]["zdsMultGriffinMult"] = new TH2I("zdsMultGriffinMult", "ZDS multiplicity vs. GRIFFIN multiplicity (unsuppressed)", 65, -0.5, 64.5, 10, -0.5, 9.5);
27
28 // unsuppressed spectra
29 fH1[slot]["griffinE"] = new TH1F("griffinE", Form("Unsuppressed griffin energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
30 // suppressed spectra
31 fH1[slot]["griffinESupp"] = new TH1F("griffinESupp", Form("Suppressed griffin energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
32 // unsuppressed addback spectra
33 fH1[slot]["griffinEAddback"] = new TH1F("griffinEAddback", Form("Unsuppressed griffin addback energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
34 // suppressed addback spectra
35 fH1[slot]["griffinESuppAddback"] = new TH1F("griffinESuppAddback", Form("Suppressed griffin addback energy;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
36 fH1[slot]["griffinESuppAddbackBeta"] = new TH1F("griffinESuppAddbackBeta", Form("Suppressed griffin addback energy w/ #beta-tag;energy [keV];counts/%.1f keV", (highEnergy - lowEnergy) / energyBins), energyBins, lowEnergy, highEnergy);
37 fH2[slot]["griffinESuppAddbackMatrixBeta"] = new TH2F("griffinESuppAddbackMatrixBeta", "Suppressed griffin addback energy matrix w/ #beta-tag;energy [keV];energy [keV]", energyBins / 5, lowEnergy, highEnergy, energyBins / 5, lowEnergy, highEnergy);
38 fH2[slot]["griffinESuppAddbackMatrixBetaBg"] = new TH2F("griffinESuppAddbackMatrixBetaBg", "Suppressed griffin addback energy matrix w/ #beta-tag (time random BG);energy [keV];energy [keV]", energyBins / 5, lowEnergy, highEnergy, energyBins / 5, lowEnergy, highEnergy);
39
40 // timing spectra
41 fH2[slot]["griffinZdsTS"] = new TH2F("griffinZdsTS", "GRIFFIN crystal vs. GRIFFIN-ZDS timestamp difference (suppressed addback);#DeltaTS_{GRIFFIN-ZDS}", 200, -1000., 1000., 64, 0.5, 64.5);
42 fH2[slot]["griffinZdsTime"] = new TH2F("griffinZdsTime", "GRIFFIN crystal vs. GRIFFIN-ZDS timing (suppressed addback);#Deltat_{GRIFFIN-ZDS}", 2000, -1000., 1000., 64, 0.5, 64.5);
43 fH2[slot]["griffinGriffinTS"] = new TH2F("griffinGriffinTS", "GRIFFIN crystal vs. GRIFFIN-GRIFFIN timestamp difference (suppressed addback);#DeltaTS_{GRIFFIN-GRIFFIN}", 2000, -1000., 1000., 64, 0.5, 64.5);
44 fH2[slot]["griffinGriffinTime"] = new TH2F("griffinGriffinTime", "GRIFFIN crystal vs. GRIFFIN-GRIFFIN timing (suppressed addback);#Deltat_{GRIFFIN-GRIFFIN}", 2000, -1000., 1000., 64, 0.5, 64.5);
45
46 // cycle spectra
47 if(fCycleLength > 0) {
48 fH2[slot]["griffinCycle"] = new TH2F("griffinCycle", "GRIFFIN suppressed addback energy w/ #beta-tag vs. time in cycle;time in cycle [s];energy [keV]", 100 * fCycleLength / 1e6, 0., fCycleLength / 1e6, energyBins / 5, lowEnergy, highEnergy);
49 fH1[slot]["zdsCycle"] = new TH1F("zdsCycle", "ZDS hits in cycle;time in cycle [s]", 100 * fCycleLength / 1e6, 0., fCycleLength / 1e6);
50 }
51}
52
53// Coincidences Gates
55{ // Griffin-Griffin
56 return -250. <= h2->GetTime() - h1->GetTime() && h2->GetTime() - h1->GetTime() <= 250.;
57}
59{
60 return (-500. <= h1->GetTime() - h2->GetTime() && h1->GetTime() - h2->GetTime() <= -250.) || (250. <= h1->GetTime() - h2->GetTime() && h1->GetTime() - h2->GetTime() <= 500.);
61}
63{ // Griffin-Zds
64 return -200. <= h1->GetTime() - h2->GetTime() && h1->GetTime() - h2->GetTime() <= 20.;
65}
67{
68 return (-310. <= h1->GetTime() - h2->GetTime() && h1->GetTime() - h2->GetTime() <= -200.) || (20. <= h1->GetTime() - h2->GetTime() && h1->GetTime() - h2->GetTime() <= 130.);
69}
70
71// TODO: Change the function arguments to match the detectors you want to use and the declaration in the header file!
72void ExampleEventHelper::Exec(unsigned int slot, TGriffin& grif, TGriffinBgo& grifBgo, TZeroDegree& zds) // NOLINT
73{
74 // we use .at() here instead of [] so that we get meaningful error message if a histogram we try to fill wasn't created
75 // e.g. because of a typo
76
77 // multiplicities
78 fH2[slot].at("zdsMultGriffinMult")->Fill(grif.GetMultiplicity(), zds.GetMultiplicity());
79
80 // loop over unsuppressed griffin hits
81 for(int g = 0; g < grif.GetMultiplicity(); ++g) {
82 auto grif1 = grif.GetGriffinHit(g);
83 fH1[slot].at("griffinE")->Fill(grif1->GetEnergy());
84 }
85
86 // loop over suppressed griffin hits
87 for(int g = 0; g < grif.GetSuppressedMultiplicity(&grifBgo); ++g) {
88 auto grif1 = grif.GetSuppressedHit(g);
89 fH1[slot].at("griffinESupp")->Fill(grif1->GetEnergy());
90 }
91
92 // loop over unsuppressed griffin addback hits
93 for(int g = 0; g < grif.GetAddbackMultiplicity(); ++g) {
94 auto grif1 = grif.GetAddbackHit(g);
95 fH1[slot].at("griffinEAddback")->Fill(grif1->GetEnergy());
96 }
97
98 // loop over suppressed griffin addback hits
99 for(int g = 0; g < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++g) {
100 auto grif1 = grif.GetSuppressedAddbackHit(g);
101 fH1[slot].at("griffinESuppAddback")->Fill(grif1->GetEnergy());
102 // we use a flag to check if there is any beta in coincidence with this gamma-ray and only fill afterwards
103 // otherwise we could double-count gamma-rays that are in coincidence with more than one beta
104 // e.g. because of double triggering of SCEPTAR/ZDS where we get a second trigger from the tail of the real signal
105 bool promptBeta = false;
106 for(int z = 0; z < zds.GetMultiplicity(); ++z) {
107 auto zds1 = zds.GetZeroDegreeHit(z);
108 fH2[slot].at("griffinZdsTS")->Fill(grif1->GetTimeStampNs() - zds1->GetTimeStampNs(), grif1->GetArrayNumber());
109 fH2[slot].at("griffinZdsTime")->Fill(grif1->GetTime() - zds1->GetTime(), grif1->GetArrayNumber());
110 if(PromptCoincidence(grif1, zds1)) promptBeta = true;
111 }
112 if(promptBeta) {
113 fH1[slot].at("griffinESuppAddbackBeta")->Fill(grif1->GetEnergy());
114 if(fCycleLength > 0.) {
115 fH2[slot].at("griffinCycle")->Fill(std::fmod(grif1->GetTime(), fCycleLength), grif1->GetEnergy());
116 }
117 for(int g2 = g + 1; g2 < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++g2) {
118 auto grif2 = grif.GetSuppressedAddbackHit(g2);
119 fH2[slot].at("griffinGriffinTS")->Fill(grif1->GetTimeStampNs() - grif2->GetTimeStampNs(), grif1->GetArrayNumber());
120 fH2[slot].at("griffinGriffinTime")->Fill(grif1->GetTime() - grif2->GetTime(), grif1->GetArrayNumber());
121 if(PromptCoincidence(grif1, grif2)) {
122 // fill twice to get a symmetric matrix
123 fH2[slot].at("griffinESuppAddbackMatrixBeta")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
124 fH2[slot].at("griffinESuppAddbackMatrixBeta")->Fill(grif2->GetEnergy(), grif1->GetEnergy());
125 } else if(TimeRandom(grif1, grif2)) {
126 // fill twice to get a symmetric matrix
127 fH2[slot].at("griffinESuppAddbackMatrixBetaBg")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
128 fH2[slot].at("griffinESuppAddbackMatrixBetaBg")->Fill(grif2->GetEnergy(), grif1->GetEnergy());
129 }
130 }
131 }
132 }
133 if(fCycleLength > 0.) {
134 for(int z = 0; z < zds.GetMultiplicity(); ++z) {
135 auto zds1 = zds.GetZeroDegreeHit(z);
136 fH1[slot].at("zdsCycle")->Fill(std::fmod(zds1->GetTime(), fCycleLength));
137 }
138 }
139}
140
141void ExampleEventHelper::EndOfSort(std::shared_ptr<std::map<std::string, TList>>& list)
142{
143 auto coincident = static_cast<TH2*>(list->at("").FindObject(fH2[0].at("griffinESuppAddbackMatrixBeta")));
144 if(coincident == nullptr) {
145 std::cout << "Failed to find griffinESuppAddbackMatrixBeta histogram in list:" << std::endl;
146 list->at("").Print();
147 return;
148 }
149 auto timeRandom = static_cast<TH2*>(list->at("").FindObject(fH2[0].at("griffinESuppAddbackMatrixBetaBg")));
150 if(timeRandom == nullptr) {
151 std::cout << "Failed to find griffinESuppAddbackMatrixBetaBg histogram in list:" << std::endl;
152 list->at("").Print();
153 return;
154 }
155
156 auto corrected = static_cast<TH2*>(coincident->Clone("griffinESuppAddbackMatrixBetaCorr"));
157 // coinc = -250 - 250 = 500 wide, bg = -500 - -250 plus 250 - 500 = 500 wide
158 corrected->Add(timeRandom, -1.);
159 list->at("").Add(corrected);
160}
bool PromptCoincidence(TGriffinHit *one, TGriffinHit *two)
bool TimeRandom(TGriffinHit *h1, TGriffinHit *h2)
bool TimeRandom(TGriffinHit *h1, TGriffinHit *h2)
bool PromptCoincidence(TGriffinHit *h1, TGriffinHit *h2)
#define DRED
Definition Globals.h:18
#define RESET_COLOR
Definition Globals.h:5
void EndOfSort(std::shared_ptr< std::map< std::string, TList > > &list) override
This method gets called at the end of Finalize()
void Exec(unsigned int slot, TGriffin &grif, TGriffinBgo &grifBgo, TZeroDegree &zds)
void CreateHistograms(unsigned int slot) override
Virtual helper function that the user uses to create their histograms.
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:73
std::vector< TGRSIMap< std::string, TH2 * > > fH2
Definition TGRSIHelper.h:49
std::vector< TGRSIMap< std::string, TH1 * > > fH1
Definition TGRSIHelper.h:48
TPPG * Ppg()
Definition TGRSIHelper.h:43
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
Long64_t OdbCycleLength() const
Definition TPPG.h:171
TZeroDegreeHit * GetZeroDegreeHit(const int &i) const
Definition TZeroDegree.h:36