BmnRoot
Loading...
Searching...
No Matches
BmnSRCTriggersCheck.cxx
Go to the documentation of this file.
2#include "BmnTrigWaveDigit.h"
4 fIsExp = isExp;
5
6 fT0Branch = "T0";
7 fVetoBranch = "VETO";
8 //fBC2Branch = "BC2";
9 fBC1Branch = "TQDC_BC1";
10 fBC2Branch = "TQDC_BC2";
11 fBC3Branch = "TQDC_BC3";
12 fX1LBranch = "TQDC_X1_Left";
13 fX2LBranch = "TQDC_X2_Left";
14 fX1RBranch = "TQDC_X1_Right";
15 fX2RBranch = "TQDC_X2_Right";
16 fY1LBranch = "TQDC_Y1_Left";
17 fY2LBranch = "TQDC_Y2_Left";
18 fY1RBranch = "TQDC_Y1_Right";
19 fY2RBranch = "TQDC_Y2_Right";
20 fBDBranch = "BD";
21 fBmnEventHeaderBranchName = "EventHeader";
22 fBmnEvQualityBranchName = "srcTriggers";
23}
24
26 FairRootManager* ioman = FairRootManager::Instance();
27
28 fBmnEventHeader = (TClonesArray*) ioman->GetObject(fBmnEventHeaderBranchName);
29 fT0Array = (TClonesArray*) ioman->GetObject(fT0Branch.Data());
30 fVetoArray = (TClonesArray*) ioman->GetObject(fVetoBranch.Data());
31 fBC1Array = (TClonesArray*) ioman->GetObject(fBC2Branch.Data());
32 fBC2Array = (TClonesArray*) ioman->GetObject(fBC2Branch.Data());
33 fBC3Array = (TClonesArray*) ioman->GetObject(fBC3Branch.Data());
34 fX1RArray = (TClonesArray*) ioman->GetObject(fX1RBranch.Data());
35 fX1LArray = (TClonesArray*) ioman->GetObject(fX1LBranch.Data());
36 fX2RArray = (TClonesArray*) ioman->GetObject(fX2RBranch.Data());
37 fX2LArray = (TClonesArray*) ioman->GetObject(fX2LBranch.Data());
38 fY1RArray = (TClonesArray*) ioman->GetObject(fY1RBranch.Data());
39 fY1LArray = (TClonesArray*) ioman->GetObject(fY1LBranch.Data());
40 fY2RArray = (TClonesArray*) ioman->GetObject(fY2RBranch.Data());
41 fY2LArray = (TClonesArray*) ioman->GetObject(fY2LBranch.Data());
42 fBDArray = (TClonesArray*) ioman->GetObject(fBDBranch.Data());
43
44 fBmnEvQuality = new TClonesArray(fBmnEvQualityBranchName);
45 ioman->Register(fBmnEvQualityBranchName, "QUALITY", fBmnEvQuality, true);
46}
47
48void BmnSRCTriggersCheck::Exec(Option_t* opt) {
49 if (!fIsExp)
50 return;
51
52 fBmnEvQuality->Delete();
53
54 BmnEventQuality* evQual = new ((*fBmnEvQuality)[fBmnEvQuality->GetEntriesFast()]) BmnEventQuality("GOOD");
55
56 const Int_t kEnergies = 4;
57 const Int_t kTargets = 5;
58 Double_t energies[kEnergies] = {3.5, 4., 4.5, 5.14};
59 TString targets[kTargets] = {"C", "Al", "Cu", "Pb", "C2H4"};
60
61 BmnEventHeader* evHeader = (BmnEventHeader*) fBmnEventHeader->UncheckedAt(0);
62 if (!evHeader)
63 return;
64
65 UniRun* runInfo = UniRun::GetRun(6, evHeader->GetRunId());
66
67 BmnTriggerType trigType = evHeader->GetTrig();
68 if (trigType == kBMNMINBIAS) {
69 // Setup BD-threshold
70 // To be read subsequently from the UniDB directly when processing RUN7, FIXME
71 TString target = *runInfo->GetTargetParticle();
72 Double_t energy = *runInfo->GetEnergy();
73 Int_t bdThresh;
74 if (Abs(energy - energies[0]) < FLT_EPSILON) // 3.5 GeV/n, all targets
75 bdThresh = 2;
76 else if (Abs(energy - energies[1]) < FLT_EPSILON && (target == targets[0] || target == targets[4])) // 4 GeV/n, C, C2H4
77 bdThresh = 2;
78 else if (Abs(energy - energies[1]) < FLT_EPSILON && target != targets[0] && target != targets[4]) // 4 GeV/n, Al, Cu, Pb
79 bdThresh = 3;
80 else if (Abs(energy - energies[2]) < FLT_EPSILON && (target == targets[0] || target == targets[4])) // 4.5 GeV/n, C, C2H4
81 bdThresh = 2;
82 else if (Abs(energy - energies[2]) < FLT_EPSILON && target != targets[0] && target != targets[4]) // 4 GeV/n, Al, Cu, Pb
83 bdThresh = 3;
84 else if (Abs(energy - energies[3]) < FLT_EPSILON) // 5.14 GeV/n, all targets
85 bdThresh = 2;
86 else
87 bdThresh = 0;
88
89 if (fT0Array->GetEntriesFast() != 1 || fBC2Array->GetEntriesFast() != 1 || fVetoArray->GetEntriesFast() != 0 || fBDArray->GetEntriesFast() < bdThresh)
90 evQual->SetIsGoodEvent("BAD");
91 }
92
93 else if (trigType == kBMNBEAM) {
94 if (fVetoArray->GetEntriesFast() > 1 || fBDArray->GetEntriesFast() > 0)
95 evQual->SetIsGoodEvent("BAD");
96 }
97
98 else
99 evQual->SetIsGoodEvent("BAD");
100
101 vector<int> eloss;
102 for (Int_t iDig = 0; iDig < fBC1Array->GetEntriesFast(); ++iDig) {
103 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fBC1Array->At(iDig);
104 eloss.push_back(dig->GetIntegral());
105 }
106 evQual->SetElossBC1(eloss);
107
108 eloss.clear();
109 for (Int_t iDig = 0; iDig < fBC2Array->GetEntriesFast(); ++iDig) {
110 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fBC2Array->At(iDig);
111 eloss.push_back(dig->GetIntegral());
112 }
113 evQual->SetElossBC2(eloss);
114
115 eloss.clear();
116 for (Int_t iDig = 0; iDig < fBC3Array->GetEntriesFast(); ++iDig) {
117 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fBC3Array->At(iDig);
118 eloss.push_back(dig->GetIntegral());
119 }
120 evQual->SetElossBC3(eloss);
121
122 eloss.clear();
123 for (Int_t iDig = 0; iDig < fX1RArray->GetEntriesFast(); ++iDig) {
124 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fX1RArray->At(iDig);
125 eloss.push_back(dig->GetIntegral());
126 }
127 evQual->SetElossX1R(eloss);
128
129 eloss.clear();
130 for (Int_t iDig = 0; iDig < fX1LArray->GetEntriesFast(); ++iDig) {
131 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fX1LArray->At(iDig);
132 eloss.push_back(dig->GetIntegral());
133 }
134 evQual->SetElossX1L(eloss);
135
136 eloss.clear();
137 for (Int_t iDig = 0; iDig < fY1RArray->GetEntriesFast(); ++iDig) {
138 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fY1RArray->At(iDig);
139 eloss.push_back(dig->GetIntegral());
140 }
141 evQual->SetElossY1R(eloss);
142
143 eloss.clear();
144 for (Int_t iDig = 0; iDig < fY1LArray->GetEntriesFast(); ++iDig) {
145 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fY1LArray->At(iDig);
146 eloss.push_back(dig->GetIntegral());
147 }
148 evQual->SetElossY1L(eloss);
149
150 eloss.clear();
151 for (Int_t iDig = 0; iDig < fX2RArray->GetEntriesFast(); ++iDig) {
152 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fX2RArray->At(iDig);
153 eloss.push_back(dig->GetIntegral());
154 }
155 evQual->SetElossX2R(eloss);
156
157 eloss.clear();
158 for (Int_t iDig = 0; iDig < fX2LArray->GetEntriesFast(); ++iDig) {
159 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fX2LArray->At(iDig);
160 eloss.push_back(dig->GetIntegral());
161 }
162 evQual->SetElossX2L(eloss);
163
164 eloss.clear();
165 for (Int_t iDig = 0; iDig < fY2RArray->GetEntriesFast(); ++iDig) {
166 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fY2RArray->At(iDig);
167 eloss.push_back(dig->GetIntegral());
168 }
169 evQual->SetElossY2R(eloss);
170
171 eloss.clear();
172 for (Int_t iDig = 0; iDig < fY2LArray->GetEntriesFast(); ++iDig) {
173 BmnTrigWaveDigit *dig = (BmnTrigWaveDigit*) fY2LArray->At(iDig);
174 eloss.push_back(dig->GetIntegral());
175 }
176 evQual->SetElossY2L(eloss);
177}
178
BmnTriggerType
Definition BmnEnums.h:61
@ kBMNMINBIAS
Definition BmnEnums.h:63
@ kBMNBEAM
Definition BmnEnums.h:62
void SetIsGoodEvent(TString str)
virtual void Exec(Option_t *opt)
virtual InitStatus Init()
int GetIntegral(UInt_t start=0, UInt_t stop=1e9) const
double * GetEnergy()
get energy of the current run
Definition UniRun.h:117
static UniRun * GetRun(int period_number, int run_number)
get run from the database
Definition UniRun.cxx:177
TString GetTargetParticle()
get target particle of the current run
Definition UniRun.h:115