GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
ExampleTreeSelector.C
Go to the documentation of this file.
1#define ExampleTreeSelector_cxx
2// The class definition in ExampleTreeSelector.h has been generated automatically
4
6{
7 // define Trees
8 fTree["coinc"] = new TTree("coinc", "coincident events");
9 fTree["timebg"] = new TTree("timebg", "time random events");
10
11 // set branch addresses for output tree (these can be different for different trees)
12 // we save the entry number as well to detect possible double counting
13 // four coincident gammas a,b,c,d will be saved as a,b,c, a,b,d, a,c,d, and b,c,d
14 fTree["coinc"]->Branch("energy", &fSuppressedAddback);
15 fTree["coinc"]->Branch("timing", &fBetaGammaTiming);
16 fTree["coinc"]->Branch("entry", &fEntry, "entry/L");
17
18 fTree["timebg"]->Branch("energy", &fSuppressedAddback);
19 fTree["timebg"]->Branch("timing", &fBetaGammaTiming);
20 fTree["timebg"]->Branch("entry", &fEntry, "entry/L");
21
22 // We can also create histograms at the same time, maybe for some diagnostics or simple checks
23 fH1["asE"] = new TH1D("asE", "suppressed Addback Singles", 12000, 0, 3000);
24 fH2["sceptarMult"] = new TH2D("sceptarMult", "Sceptar multiplicity vs. Griffin multiplicity", 64, 0, 64, 40, 0, 40);
25 fH2["zdsMult"] = new TH2D("zdsMult", "ZeroDegree multiplicity vs. Griffin multiplicity", 64, 0, 64, 10, 0, 10);
26
27 // Send histograms to Output list to be added and written.
28 for(auto it : fH1) {
29 GetOutputList()->Add(it.second);
30 }
31 for(auto it : fH2) {
32 GetOutputList()->Add(it.second);
33 }
34 for(auto it : fTree) {
35 GetOutputList()->Add(it.second);
36 }
37}
38
40{
41 // Check if hits are less then 300 ns apart.
42 return std::fabs(g->GetTime() - z->GetTime()) < 300.;
43}
44
46{
47 // Check if hits are less then 300 ns apart.
48 return std::fabs(g->GetTime() - s->GetTime()) < 300.;
49}
50
52{
53 // Check if hits are less then 500 ns apart.
54 return std::fabs(g1->GetTime() - g2->GetTime()) < 500.;
55}
56
58{
59 // we could check multiplicities here and skip events where we do not have at least
60 // three suppressed addback energies and a beta, but we want to fill some general
61 // histograms without these cuts.
62
63 // clear the vectors and other variables
64 fSuppressedAddback.resize(3, 0.);
65 fBetaGammaTiming.resize(3, -1e6);
66
69
70 // Loop over all suppressed addback Griffin Hits
71 for(auto i = 0; i < fGrif->GetSuppressedAddbackMultiplicity(fGriffinBgo); ++i) {
72 auto grif1 = fGrif->GetSuppressedAddbackHit(i);
73 fH1.at("asE")->Fill(grif1->GetEnergy());
74 fSuppressedAddback[0] = grif1->GetEnergy();
75 for(auto j = i + 1; j < fGrif->GetSuppressedAddbackMultiplicity(fGriffinBgo); ++j) {
76 auto grif2 = fGrif->GetSuppressedAddbackHit(j);
77 fSuppressedAddback[1] = grif2->GetEnergy();
78 for(auto k = j + 1; k < fGrif->GetSuppressedAddbackMultiplicity(fGriffinBgo); ++k) {
79 auto grif3 = fGrif->GetSuppressedAddbackHit(k);
80 fSuppressedAddback[2] = grif3->GetEnergy();
81 // we now have three suppressed addback hits i, j, and k so now we need a coincident beta-tag
82 bool foundBeta = false;
83 for(auto b = 0; b < fZds->GetMultiplicity(); ++b) {
84 auto beta = fZds->GetZeroDegreeHit(b);
85 if(b == 0) {
86 // use the time of the first beta as reference in case we don't find a coincident beta
87 fBetaGammaTiming[0] = grif1->GetTime() - beta->GetTime();
88 fBetaGammaTiming[1] = grif2->GetTime() - beta->GetTime();
89 fBetaGammaTiming[2] = grif3->GetTime() - beta->GetTime();
90 }
91 if(PromptCoincidence(grif1, beta) && PromptCoincidence(grif2, beta) && PromptCoincidence(grif3, beta)) {
92 foundBeta = true;
93 fBetaGammaTiming[0] = grif1->GetTime() - beta->GetTime();
94 fBetaGammaTiming[1] = grif2->GetTime() - beta->GetTime();
95 fBetaGammaTiming[2] = grif3->GetTime() - beta->GetTime();
96 break;
97 }
98 }
99 // only check sceptar if we haven't found a zds yet
100 for(auto b = 0; !foundBeta && b < fScep->GetMultiplicity(); ++b) {
101 auto beta = fScep->GetSceptarHit(b);
102 if(b == 0) {
103 // use the time of the first beta as reference in case we don't find a coincident beta
104 fBetaGammaTiming[0] = grif1->GetTime() - beta->GetTime();
105 fBetaGammaTiming[1] = grif2->GetTime() - beta->GetTime();
106 fBetaGammaTiming[2] = grif3->GetTime() - beta->GetTime();
107 }
108 if(PromptCoincidence(grif1, beta) && PromptCoincidence(grif2, beta) && PromptCoincidence(grif3, beta)) {
109 foundBeta = true;
110 fBetaGammaTiming[0] = grif1->GetTime() - beta->GetTime();
111 fBetaGammaTiming[1] = grif2->GetTime() - beta->GetTime();
112 fBetaGammaTiming[2] = grif3->GetTime() - beta->GetTime();
113 break;
114 }
115 }
116
117 if(foundBeta) {
118 fTree.at("coinc")->Fill();
119 } else {
120 fTree.at("timebg")->Fill();
121 }
122 }
123 }
124 }
125}
bool PromptCoincidence(TGriffinHit *one, TGriffinHit *two)
bool PromptCoincidence(TGriffinHit *g, TZeroDegreeHit *z)
const Double_t s
TSceptar * fScep
pointer for sceptar branch
std::vector< double > fSuppressedAddback
vector of suppressed addback energies
TGriffinBgo * fGriffinBgo
pointer for griffin BGO branch
std::vector< double > fBetaGammaTiming
vector of beta-gamma timing
TZeroDegree * fZds
pointer for ZDS branch
TGriffin * fGrif
pointer for griffin branch
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
TList * GetOutputList() const override
this does the same as TSelector::GetOutputList()
TGRSIMap< std::string, TTree * > fTree
TGriffinHit * GetSuppressedAddbackHit(const int &i)
Definition TGriffin.h:122
Short_t GetSuppressedAddbackMultiplicity(const TBgo *bgo)
Definition TGriffin.h:119
TSceptarHit * GetSceptarHit(const int &i) const
Definition TSceptar.h:36
TZeroDegreeHit * GetZeroDegreeHit(const int &i) const
Definition TZeroDegree.h:36