GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
ExampleTreeHelper.cxx
Go to the documentation of this file.
2
4{
5 // define Trees
6 fTree[slot].emplace("coinc", new TTree("coinc", "coincident events"));
7 fTree[slot]["timebg"] = new TTree("timebg", "time random events");
8
9 fSuppressedAddback2[slot] = new double[3];
10 // set branch addresses for output tree (these can be different for different trees)
11 // four coincident gammas a,b,c,d will be saved as a,b,c, a,b,d, a,c,d, and b,c,d
12 fTree[slot]["coinc"]->Branch("energy", &(fSuppressedAddback2[slot]), "energy[3]/D");
13 // fTree[slot]["coinc"]->Branch("timing", &(fBetaGammaTiming[slot]));
14 fTree[slot]["coinc"]->Branch("mult", &(fGriffinMultiplicity[slot]), "mult/I");
15
16 fTree[slot]["timebg"]->Branch("energy", &(fSuppressedAddback[slot]));
17 fTree[slot]["timebg"]->Branch("timing", &(fBetaGammaTiming[slot]));
18
19 // We can also create histograms at the same time, maybe for some diagnostics or simple checks
20 fH1[slot]["asE"] = new TH1D("asE", "suppressed Addback Singles", 12000, 0, 3000);
21 fH2[slot]["sceptarMult"] = new TH2D("sceptarMult", "Sceptar multiplicity vs. Griffin multiplicity", 64, 0, 64, 40, 0, 40);
22 fH2[slot]["zdsMult"] = new TH2D("zdsMult", "ZeroDegree multiplicity vs. Griffin multiplicity", 64, 0, 64, 10, 0, 10);
23}
24
26{
27 // Check if hits are less then 300 ns apart.
28 return std::fabs(g->GetTime() - z->GetTime()) < 300.;
29}
30
32{
33 // Check if hits are less then 300 ns apart.
34 return std::fabs(g->GetTime() - s->GetTime()) < 300.;
35}
36
38{
39 // Check if hits are less then 500 ns apart.
40 return std::fabs(g1->GetTime() - g2->GetTime()) < 500.;
41}
42
43void ExampleTreeHelper::Exec(unsigned int slot, TGriffin& grif, TGriffinBgo& grifBgo, TZeroDegree& zds, TSceptar& scep)
44{
45 // we could check multiplicities here and skip events where we do not have at least
46 // three suppressed addback energies and a beta, but we want to fill some general
47 // histograms without these cuts.
48
49 // clear the vectors and other variables
50 fSuppressedAddback[slot].resize(3, 0.);
51 fBetaGammaTiming[slot].resize(3, -1e6);
52 for(int i = 0; i < 3; ++i) fSuppressedAddback2[slot][i] = 0.;
53
54 fH2[slot].at("sceptarMult")->Fill(grif.GetSuppressedAddbackMultiplicity(&grifBgo), scep.GetMultiplicity());
55 fH2[slot].at("zdsMult")->Fill(grif.GetSuppressedAddbackMultiplicity(&grifBgo), zds.GetMultiplicity());
56
57 // Loop over all suppressed addback Griffin Hits
59 if(fGriffinMultiplicity[slot] < 3) { return; }
60 for(auto i = 0; i < fGriffinMultiplicity[slot]; ++i) {
61 auto grif1 = grif.GetSuppressedAddbackHit(i);
62 fH1[slot].at("asE")->Fill(grif1->GetEnergy());
63 fSuppressedAddback[slot][0] = grif1->GetEnergy();
64 fSuppressedAddback2[slot][0] = grif1->GetEnergy();
65 for(auto j = i + 1; j < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++j) {
66 auto grif2 = grif.GetSuppressedAddbackHit(j);
67 fSuppressedAddback[slot][1] = grif2->GetEnergy();
68 fSuppressedAddback2[slot][1] = grif2->GetEnergy();
69 for(auto k = j + 1; k < grif.GetSuppressedAddbackMultiplicity(&grifBgo); ++k) {
70 auto grif3 = grif.GetSuppressedAddbackHit(k);
71 fSuppressedAddback[slot][2] = grif3->GetEnergy();
72 fSuppressedAddback2[slot][2] = grif3->GetEnergy();
73 // we now have three suppressed addback hits i, j, and k so now we need a coincident beta-tag
74 bool foundBeta = false;
75 for(auto b = 0; b < zds.GetMultiplicity(); ++b) {
76 auto beta = zds.GetZeroDegreeHit(b);
77 if(b == 0) {
78 // use the time of the first beta as reference in case we don't find a coincident beta
79 fBetaGammaTiming[slot][0] = grif1->GetTime() - beta->GetTime();
80 fBetaGammaTiming[slot][1] = grif2->GetTime() - beta->GetTime();
81 fBetaGammaTiming[slot][2] = grif3->GetTime() - beta->GetTime();
82 }
83 if(PromptCoincidence(grif1, beta) && PromptCoincidence(grif2, beta) && PromptCoincidence(grif3, beta)) {
84 foundBeta = true;
85 fBetaGammaTiming[slot][0] = grif1->GetTime() - beta->GetTime();
86 fBetaGammaTiming[slot][1] = grif2->GetTime() - beta->GetTime();
87 fBetaGammaTiming[slot][2] = grif3->GetTime() - beta->GetTime();
88 break;
89 }
90 }
91 // only check sceptar if we haven't found a zds yet
92 for(auto b = 0; !foundBeta && b < scep.GetMultiplicity(); ++b) {
93 auto beta = scep.GetSceptarHit(b);
94 if(b == 0) {
95 // use the time of the first beta as reference in case we don't find a coincident beta
96 fBetaGammaTiming[slot][0] = grif1->GetTime() - beta->GetTime();
97 fBetaGammaTiming[slot][1] = grif2->GetTime() - beta->GetTime();
98 fBetaGammaTiming[slot][2] = grif3->GetTime() - beta->GetTime();
99 }
100 if(PromptCoincidence(grif1, beta) && PromptCoincidence(grif2, beta) && PromptCoincidence(grif3, beta)) {
101 foundBeta = true;
102 fBetaGammaTiming[slot][0] = grif1->GetTime() - beta->GetTime();
103 fBetaGammaTiming[slot][1] = grif2->GetTime() - beta->GetTime();
104 fBetaGammaTiming[slot][2] = grif3->GetTime() - beta->GetTime();
105 break;
106 }
107 }
108
109 if(foundBeta) {
110 fTree[slot].at("coinc")->Fill();
111 } else {
112 fTree[slot].at("timebg")->Fill();
113 }
114 }
115 }
116 }
117}
bool PromptCoincidence(TGriffinHit *one, TGriffinHit *two)
bool PromptCoincidence(TGriffinHit *g, TZeroDegreeHit *z)
const Double_t s
void CreateHistograms(unsigned int slot) override
Virtual helper function that the user uses to create their histograms.
std::map< unsigned int, int > fGriffinMultiplicity
multiplicity of suppressed addback energies
void Exec(unsigned int slot, TGriffin &grif, TGriffinBgo &grifBgo, TZeroDegree &zds, TSceptar &scep)
std::map< unsigned int, std::vector< double > > fBetaGammaTiming
vector of beta-gamma timing
std::map< unsigned int, std::vector< double > > fSuppressedAddback
vector of suppressed addback energies
std::map< unsigned int, double * > fSuppressedAddback2
vector of suppressed addback energies
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, TTree * > > fTree
Definition TGRSIHelper.h:53
std::vector< TGRSIMap< std::string, TH2 * > > fH2
Definition TGRSIHelper.h:49
std::vector< TGRSIMap< std::string, TH1 * > > fH1
Definition TGRSIHelper.h:48
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