GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
EfficiencyHelper.cxx
Go to the documentation of this file.
1#include "EfficiencyHelper.hh"
2
3void EfficiencyHelper::CreateHistograms(unsigned int slot)
4{
5 // three histograms for each type: 1D singles, coincident and time random 2D matrix (unsuppressed on y-axis for suppressed data)
6 // all 2D histograms are at 180 degree
7 // 180 degree means for addback that the detectors are opposite, so the crystal angles are 158 degrees and larger!
8 // (for 145 mm distance)
9
10 // unsuppressed spectra
11 fH1[slot]["griffinE"] = new TH1F("griffinE", Form("Unsuppressed griffin energy;energy [keV];counts/%.1f keV", (fHighEnergy - fLowEnergy) / fEnergyBins1D), fEnergyBins1D, fLowEnergy, fHighEnergy);
12 fH2[slot]["griffinGriffinE"] = new TH2F("griffinGriffinE", "Unsuppressed griffin-griffin energy @ 180^{o};energy [keV];energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
13 fH2[slot]["griffinGriffinESum"] = new TH2F("griffinGriffinESum", "Unsuppressed griffin energy vs griffin sum energy @ 180^{o};sum energy [keV];energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
14 // suppressed spectra
15 fH1[slot]["griffinESupp"] = new TH1F("griffinESupp", Form("Suppressed griffin energy;energy [keV];counts/%.1f keV", (fHighEnergy - fLowEnergy) / fEnergyBins1D), fEnergyBins1D, fLowEnergy, fHighEnergy);
16 fH2[slot]["griffinGriffinESuppSum"] = new TH2F("griffinGriffinESuppSum", "Suppressed griffin energy vs griffin sum energy @ 180^{o};sum energy [keV];energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
17 fH2[slot]["griffinGriffinEMixed"] = new TH2F("griffinGriffinEMixed", "Unsuppressed/suppressed griffin-griffin energy @ 180^{o};suppressed energy [keV];unsuppressed energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
18 // unsuppressed addback spectra
19 fH1[slot]["griffinEAddback"] = new TH1F("griffinEAddback", Form("Unsuppressed griffin addback energy;energy [keV];counts/%.1f keV", (fHighEnergy - fLowEnergy) / fEnergyBins1D), fEnergyBins1D, fLowEnergy, fHighEnergy);
20 fH2[slot]["griffinGriffinEAddback"] = new TH2F("griffinGriffinEAddback", "Unsuppressed griffin-griffin addback @ 180^{o};energy [keV];energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
21 fH2[slot]["griffinGriffinEAddbackSum"] = new TH2F("griffinGriffinEAddbackSum", "Unsuppressed griffin addback vs griffin sum energy @ 180^{o};sum energy [keV];energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
22 // suppressed addback spectra
23 fH1[slot]["griffinESuppAddback"] = new TH1F("griffinESuppAddback", Form("Suppressed griffin addback energy;energy [keV];counts/%.1f keV", (fHighEnergy - fLowEnergy) / fEnergyBins1D), fEnergyBins1D, fLowEnergy, fHighEnergy);
24 fH2[slot]["griffinGriffinESuppAddbackSum"] = new TH2F("griffinGriffinESuppAddbackSum", "Suppressed griffin addback vs sum energy @ 180^{o};sum energy [keV];energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
25 fH2[slot]["griffinGriffinEMixedAddback"] = new TH2F("griffinGriffinEMixedAddback", "Unsuppressed/suppressed griffin-griffin addback @ 180^{o};suppressed energy [keV];unsuppressed energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
26 // suppressed single crystal addback spectra
27 fH1[slot]["griffinESingleCrystal"] = new TH1F("griffinESingleCrystal", Form("Suppressed griffin addback energy single crystal method;energy [keV];counts/%.1f keV", (fHighEnergy - fLowEnergy) / fEnergyBins1D), fEnergyBins1D, fLowEnergy, fHighEnergy);
28 fH2[slot]["griffinGriffinESingleCrystalSum"] = new TH2F("griffinGriffinESingleCrystalSum", "Suppressed griffin addback vs sum energy @ 180^{o}, single crystal method;sum energy [keV];energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
29 fH2[slot]["griffinGriffinEMixedSingleCrystal"] = new TH2F("griffinGriffinEMixedSingleCrystal", "Unsuppressed/suppressed griffin-griffin addback @ 180^{o};suppressed energy [keV];unsuppressed energy [keV]", fEnergyBins2D, fLowEnergy, fHighEnergy, fEnergyBins2D, fLowEnergy, fHighEnergy);
30}
31
32void EfficiencyHelper::Exec(unsigned int slot, TGriffin& grif, TGriffinBgo& grifBgo)
33{
34 // we use .at() here instead of [] so that we get meaningful error message if a histogram we try to fill wasn't created
35 // e.g. because of a typo
36
37 // loop over unsuppressed griffin hits
38 for(int g = 0; g < grif.GetMultiplicity(); ++g) {
39 auto grif1 = grif.GetGriffinHit(g);
40 if(Reject(grif1) || !GoodCfd(grif1)) continue;
41 fH1[slot].at("griffinE")->Fill(grif1->GetEnergy());
42 for(int g2 = g + 1; g2 < grif.GetMultiplicity(); ++g2) {
43 auto grif2 = grif.GetGriffinHit(g2);
44 if(Reject(grif2) || !GoodCfd(grif2)) continue;
45 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 179.) {
46 if(PromptCoincidence(grif1, grif2)) {
47 fH2[slot].at("griffinGriffinE")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
48 fH2[slot].at("griffinGriffinE")->Fill(grif2->GetEnergy(), grif1->GetEnergy());
49 fH2[slot].at("griffinGriffinESum")->Fill(grif1->GetEnergy() + grif2->GetEnergy(), grif1->GetEnergy());
50 }
51 }
52 }
53 }
54
55 // loop over suppressed griffin hits
56 for(int g = 0; g < grif.GetSuppressedMultiplicity(&grifBgo); ++g) {
57 auto grif1 = grif.GetSuppressedHit(g);
58 if(Reject(grif1) || !GoodCfd(grif1)) continue;
59 fH1[slot].at("griffinESupp")->Fill(grif1->GetEnergy());
60 for(int g2 = g + 1; g2 < grif.GetSuppressedMultiplicity(&grifBgo); ++g2) {
61 auto grif2 = grif.GetSuppressedHit(g2);
62 if(Reject(grif2) || !GoodCfd(grif2)) continue;
63 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 179.) {
64 if(PromptCoincidence(grif1, grif2)) {
65 fH2[slot].at("griffinGriffinESuppSum")->Fill(grif1->GetEnergy() + grif2->GetEnergy(), grif1->GetEnergy());
66 }
67 }
68 }
69 for(int g2 = 0; g2 < grif.GetMultiplicity(); ++g2) {
70 if(g == g2) continue;
71 auto grif2 = grif.GetGriffinHit(g2);
72 if(Reject(grif2) || !GoodCfd(grif2)) continue;
73 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 179.) {
74 if(PromptCoincidence(grif1, grif2)) {
75 fH2[slot].at("griffinGriffinEMixed")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
76 }
77 }
78 }
79 }
80
81 // loop over unsuppressed griffin addback hits
82 for(int g = 0; g < grif.GetAddbackMultiplicity(); ++g) {
83 auto grif1 = grif.GetAddbackHit(g);
84 if(Reject(grif1) || !GoodCfd(grif1)) continue;
85 fH1[slot].at("griffinEAddback")->Fill(grif1->GetEnergy());
86 for(int g2 = g + 1; g2 < grif.GetAddbackMultiplicity(); ++g2) {
87 auto grif2 = grif.GetAddbackHit(g2);
88 if(Reject(grif2) || !GoodCfd(grif2)) continue;
89 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 157.) {
90 if(PromptCoincidence(grif1, grif2)) {
91 fH2[slot].at("griffinGriffinEAddback")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
92 fH2[slot].at("griffinGriffinEAddback")->Fill(grif2->GetEnergy(), grif1->GetEnergy());
93 fH2[slot].at("griffinGriffinEAddbackSum")->Fill(grif1->GetEnergy() + grif2->GetEnergy(), grif1->GetEnergy());
94 }
95 }
96 }
97 }
98
99 // loop over suppressed griffin addback hits
100 for(int g = 0; g < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++g) {
101 auto grif1 = grif.GetSuppressedAddbackHit(g);
102 if(Reject(grif1) || !GoodCfd(grif1)) continue;
103 fH1[slot].at("griffinESuppAddback")->Fill(grif1->GetEnergy());
104 for(int g2 = g + 1; g2 < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++g2) {
105 auto grif2 = grif.GetSuppressedAddbackHit(g2);
106 if(Reject(grif2) || !GoodCfd(grif2)) continue;
107 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 157.) {
108 if(PromptCoincidence(grif1, grif2)) {
109 fH2[slot].at("griffinGriffinESuppAddbackSum")->Fill(grif1->GetEnergy() + grif2->GetEnergy(), grif1->GetEnergy());
110 }
111 }
112 }
113 for(int g2 = 0; g2 < grif.GetAddbackMultiplicity(); ++g2) {
114 if(g == g2) continue;
115 auto grif2 = grif.GetAddbackHit(g2);
116 if(Reject(grif2) || !GoodCfd(grif2)) continue;
117 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 157.) {
118 if(PromptCoincidence(grif1, grif2)) {
119 fH2[slot].at("griffinGriffinEMixedAddback")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
120 }
121 }
122 }
123 }
124
125 // loop over suppressed griffin addback hits, single crystal method
126 for(int g = 0; g < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++g) {
127 auto grif1 = grif.GetSuppressedAddbackHit(g);
128 if(Reject(grif1) || !GoodCfd(grif1) || grif.GetNSuppressedAddbackFrags(g) > 1) continue;
129 fH1[slot].at("griffinESingleCrystal")->Fill(grif1->GetEnergy());
130 for(int g2 = g + 1; g2 < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++g2) {
131 auto grif2 = grif.GetSuppressedAddbackHit(g2);
132 if(Reject(grif2) || !GoodCfd(grif2) || grif.GetNSuppressedAddbackFrags(g2) > 1) continue;
133 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 179.) {
134 if(PromptCoincidence(grif1, grif2)) {
135 fH2[slot].at("griffinGriffinESingleCrystalSum")->Fill(grif1->GetEnergy() + grif2->GetEnergy(), grif1->GetEnergy());
136 }
137 }
138 }
139 for(int g2 = 0; g2 < grif.GetMultiplicity(); ++g2) {
140 if(g == g2) continue;
141 auto grif2 = grif.GetGriffinHit(g2);
142 if(Reject(grif2) || !GoodCfd(grif2)) continue;
143 if(grif1->GetPosition().Angle(grif2->GetPosition()) / TMath::Pi() * 180. > 179.) {
144 if(PromptCoincidence(grif1, grif2)) {
145 fH2[slot].at("griffinGriffinEMixedSingleCrystal")->Fill(grif1->GetEnergy(), grif2->GetEnergy());
146 }
147 }
148 }
149 }
150}
bool GoodCfd(TDetectorHit *h1)
void CreateHistograms(unsigned int slot) override
Virtual helper function that the user uses to create their histograms.
bool PromptCoincidence(TGriffinHit *h1, TGriffinHit *h2)
void Exec(unsigned int slot, TGriffin &grif, TGriffinBgo &grifBgo)
bool Reject(TGriffinHit *hit)
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
TGriffinHit * GetAddbackHit(const int &i)
Definition TGriffin.cxx:289
Short_t GetAddbackMultiplicity()
Definition TGriffin.cxx:265
UShort_t GetNSuppressedAddbackFrags(const size_t &idx)
Definition TGriffin.cxx:594
TGriffinHit * GetSuppressedAddbackHit(const int &i)
Definition TGriffin.cxx:539
Short_t GetSuppressedMultiplicity(const TBgo *bgo)
Definition TGriffin.cxx:502
Short_t GetSuppressedAddbackMultiplicity(const TBgo *bgo)
Definition TGriffin.cxx:555
TGriffinHit * GetSuppressedHit(const int &i)
!
Definition TGriffin.cxx:486
TGriffinHit * GetGriffinHit(const int &i)
!
Definition TGriffin.cxx:244