GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TimeWalkSelector.C
Go to the documentation of this file.
1#define TimeWalkSelector_cxx
2// The class definition in TimeWalkSelector.h has been generated automatically
3#include "TimeWalkSelector.h"
4
6{
7 // Define Histograms
8 fH2["time_eng_walk"] = new TH2F("time_eng_walk", "time_eng_walk", 120, -20, 100, 4000, 0, 4000);
9 fH2["time_eng"] = new TH2F("time_eng", "time_eng", 400, -200, 200, 4000, 0, 4000);
10 fH2["time_eng_nogate"] = new TH2F("time_eng_nogate", "time_eng_nogate", 400, -200, 200, 4000, 0, 4000);
11
12 fH2["time_engcfd"] = new TH2F("time_engcfd", "time_eng including CFD correction", 4000, -2000, 2000, 4000, 0, 4000);
13 fH2["time_engcfdnopu"] = new TH2F("time_engcfdnopu", "time_eng including CFD correction", 4000, -2000, 2000, 4000, 0, 4000);
14
15 fH2["time_eng_walk_bg"] = new TH2F("time_eng_walk_bg", "time_eng_walk_bg", 120, -20, 100, 4000, 0, 4000);
16 fH2["time_eng_bg"] = new TH2F("time_eng_bg", "time_eng_bg", 120, -20, 100, 4000, 0, 4000);
17
18 fH2["time_eng_walk_bgnoPU"] = new TH2F("time_eng_walk_bgnoPU", "time_eng_walk_bgnoPU", 120, -20, 100, 4000, 0, 4000);
19 fH2["time_eng_bgnoPU"] = new TH2F("time_eng_bgnoPU", "time_eng_bgnoPU", 120, -20, 100, 4000, 0, 4000);
20
21 fH2["time_eng_walk_noPU"] = new TH2F("time_eng_walk_noPU", "time_eng_walk_noPU", 120, -20, 100, 4000, 0, 4000);
22 fH2["time_eng_noPU"] = new TH2F("time_eng_noPU", "time_eng_noPU", 120, -20, 100, 4000, 0, 4000);
23 fH2["eng_time_noPU"] = new TH2F("eng_time_noPU", "eng_time_noPU", 8000, 0, 8000, 120, -20, 100);
24
25 fH2["time_eng_walk_high"] = new TH2F("time_eng_walk_high", "time_eng_walk_high", 120, -20, 100, 4000, 0, 4000);
26 fH2["time_eng_high"] = new TH2F("time_eng_high", "time_eng_high", 120, -20, 100, 4000, 0, 4000);
27
28 fH2["time_eng_walk2"] = new TH2F("eng_time_walk", "eng_time_walk", 4000, 0, 4000, 120, -20, 100);
29 fH2["time_eng2"] = new TH2F("eng_teng", "eng_time", 4000, 0, 4000, 120, -20, 100);
30
31 fH2["time_eng_walk_no0"] = new TH2F("time_eng_walk_no0", "time_eng_walk_no0", 120, -20, 100, 4000, 0, 4000);
32 fH2["time_eng_no0"] = new TH2F("time_eng_no0", "time_eng_no0", 120, -20, 100, 4000, 0, 4000);
33
34 fH2["time_eng_walk_beta"] = new TH2F("time_eng_walk_beta", "time_eng_walk_beta", 120, -20, 100, 4000, 0, 4000);
35 fH2["time_eng_beta"] = new TH2F("time_eng_beta", "time_eng_beta", 120, -20, 100, 4000, 0, 4000);
36
37 fH3["time_eng_chan"] = new TH3F("time_engd_chan", "time_engd_chan", 120, -20, 100, 2000, 0, 2000, 64, 0, 64);
38 fH3["time_eng_chan2"] = new TH3F("time_engd_chan2", "time_engd_chan2", 120, -20, 100, 2000, 0, 2000, 64, 0, 64);
39 fH3["time_eng_det"] = new TH3F("time_engd_det", "time_engd_det", 120, -20, 100, 2000, 0, 2000, 17, 0, 17);
40 fH3["time_eng_det2"] = new TH3F("time_engd_det2", "time_engd_det2", 120, -20, 100, 2000, 0, 2000, 17, 0, 17);
41
42 fH2["kValueChan"] = new TH2F("kValueChan", "kValueChan", 800, 0, 800, 65, 0, 65);
43 fH2["kValueTDiff"] = new TH2F("kValueTDiff", "kValueTDiff", 400, -200, 200, 800, 0, 800);
44 fH2["kValueTDiff_samechan"] = new TH2F("kValueTDiff_samechan", "kValueTDiff_samechan", 400, -200, 200, 800, 0, 800);
45 fH2["kValueTDiff_nogate"] = new TH2F("kValueTDiff_nogate", "kValueTDiff_nogate", 400, -200, 200, 800, 0, 800);
46
47 // Send histograms to Output list to be added and written.
48 for(auto it : fH1) {
49 GetOutputList()->Add(it.second);
50 }
51 for(auto it : fH2) {
52 GetOutputList()->Add(it.second);
53 }
54 for(auto it : fH3) {
55 GetOutputList()->Add(it.second);
56 }
57 for(auto it : fHSparse) {
58 GetOutputList()->Add(it.second);
59 }
60}
61
63{
64 // Check if hits are less then 300 ns apart.
65 return std::fabs(g->GetTime() - s->GetTime()) < 300.;
66}
67
69{
70 // Check if hits are less then 500 ns apart.
71 return std::fabs(g1->GetTime() - g2->GetTime()) < 500.;
72}
73
75{
76 // Loop over all Griffin Hits
77 for(auto i = 0; i < fGrif->GetMultiplicity(); ++i) {
78 auto grif1 = fGrif->GetGriffinHit(i);
79 fH2.at("kValueChan")->Fill(grif1->GetKValue(), grif1->GetArrayNumber());
80 // Loop over all sceptar hits
81 for(auto j = 0; j < fScep->GetMultiplicity(); ++j) {
82 auto scep = fScep->GetSceptarHit(j);
83 fH2.at("time_eng_bg")->Fill(grif1->GetTimeStamp() - scep->GetTimeStamp(), grif1->GetEnergy());
84 fH2.at("time_eng_walk_bg")->Fill((grif1->GetTime() - scep->GetTime()) / 10., grif1->GetEnergy());
85 if(grif1->GetKValue() == 700) {
86 fH2.at("time_eng_bgnoPU")->Fill(grif1->GetTimeStamp() - scep->GetTimeStamp(), grif1->GetEnergy());
87 fH2.at("time_eng_walk_bgnoPU")->Fill((grif1->GetTime() - scep->GetTime()) / 10., grif1->GetEnergy());
88 }
89 }
90 for(auto j = 0; j < fGrif->GetMultiplicity(); ++j) {
91 if(i == j) continue;
92 auto grif2 = fGrif->GetGriffinHit(j);
93 long timediff = grif2->GetTimeStamp() - grif1->GetTimeStamp();
94 double timediff_walk = (grif2->GetTime() - grif1->GetTime()) / 10.;
95 fH2.at("time_eng_nogate")->Fill(timediff, grif2->GetEnergy());
96 fH2.at("kValueTDiff_nogate")->Fill(timediff, grif2->GetKValue());
97 if(grif1->GetEnergy() > 0.) {
98 fH2.at("time_eng")->Fill(timediff, grif2->GetEnergy());
99 fH2.at("time_eng_walk")->Fill(timediff_walk, grif2->GetEnergy());
100 fH3.at("time_eng_chan")->Fill(timediff, grif2->GetEnergy(), grif2->GetArrayNumber()); // was crystal1 + 4*detector2???
101 fH3.at("time_eng_chan2")->Fill(timediff, grif2->GetEnergy(), grif1->GetArrayNumber()); // was crystal2 + 4*detector1???
102 fH3.at("time_eng_det")->Fill(timediff, grif2->GetEnergy(), grif2->GetDetector());
103 fH3.at("time_eng_det2")->Fill(timediff, grif2->GetEnergy(), grif1->GetDetector());
104 fH2.at("time_engcfd")->Fill(timediff_walk * 10., grif2->GetEnergy());
105 if(grif1->GetKValue() == 700 && grif2->GetKValue() == 700) {
106 fH2.at("time_engcfdnopu")->Fill(timediff_walk * 10., grif2->GetEnergy());
107 fH2.at("time_eng_noPU")->Fill(timediff, grif2->GetEnergy());
108 fH2.at("time_eng_walk_noPU")->Fill(timediff_walk, grif2->GetEnergy());
109 fH2.at("eng_time_noPU")->Fill(grif2->GetEnergy(), timediff);
110 }
111 fH2.at("kValueTDiff")->Fill(timediff, grif2->GetKValue());
112 if(grif1->GetAddress() == grif2->GetAddress()) {
113 fH2.at("kValueTDiff")->Fill(timediff, grif2->GetKValue());
114 }
115 if(grif1->GetAddress() != 0 && grif2->GetAddress() != 0) {
116 fH2.at("time_eng_no0")->Fill(timediff, grif2->GetEnergy());
117 fH2.at("time_eng_walk_no0")->Fill(timediff_walk, grif2->GetEnergy());
118 }
119 fH2.at("time_eng2")->Fill(grif2->GetEnergy(), timediff);
120 fH2.at("time_eng_walk2")->Fill(grif2->GetEnergy(), timediff_walk);
121 }
122 if(grif1->GetEnergy() < 500.) {
123 fH2.at("time_eng_high")->Fill(timediff, grif2->GetEnergy());
124 fH2.at("time_eng_walk_high")->Fill(timediff_walk, grif2->GetEnergy());
125 }
126 }
127 }
128 for(auto i = 0; i < fScep->GetMultiplicity(); ++i) {
129 auto scep1 = fScep->GetSceptarHit(i);
130 for(auto j = 0; j < fScep->GetMultiplicity(); ++j) {
131 if(i == j) continue;
132 auto scep2 = fScep->GetSceptarHit(j);
133 fH2.at("time_eng_beta")->Fill(scep2->GetTimeStamp() - scep1->GetTimeStamp(), scep2->GetEnergy());
134 fH2.at("time_eng_walk_beta")->Fill((scep2->GetTime() - scep1->GetTime()) / 10., scep2->GetEnergy());
135 }
136 }
137}
const Double_t s
bool PromptCoincidence(TGriffinHit *g, TSceptarHit *s)
virtual Long64_t GetTimeStamp(Option_t *="") const
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
TGRSIMap< std::string, TH1 * > fH1
TGRSIMap< std::string, TH2 * > fH2
TGRSIMap< std::string, THnSparseF * > fHSparse
TList * GetOutputList() const override
this does the same as TSelector::GetOutputList()
TGRSIMap< std::string, TH3 * > fH3
Short_t GetMultiplicity() const override
Definition TGriffin.h:52
TGriffinHit * GetGriffinHit(const int &i)
!
Definition TGriffin.h:47
TSceptarHit * GetSceptarHit(const int &i) const
Definition TSceptar.h:36