Belle II Software development
CDCDedxElectronCollectorModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <reconstruction/modules/CDCDedxElectronCollector/CDCDedxElectronCollectorModule.h>
10
11
12#include <TTree.h>
13#include <TH1D.h>
14#include <TH1I.h>
15#include <TMath.h>
16
17using namespace Belle2;
18//-----------------------------------------------------------------
19// Register the Module
20//-----------------------------------------------------------------
21REG_MODULE(CDCDedxElectronCollector);
22
23//-----------------------------------------------------------------
24// Implementation
25//-----------------------------------------------------------------
27{
28 // Set module properties
29 setDescription("A collector module for CDC dE/dx electron calibrations");
30
31 // Parameter definitions
32 addParam("cleanupCuts", m_cuts, "boolean to apply cleanup cuts", true);
33 addParam("maxHits", m_maxHits, "maximum number of hits per track ", int(100));
34 addParam("setEoP", m_setEoP, "Set E over p Cut values. ", double(0.25));
35 addParam("isCosth", m_isCosth, "true for adding costh tree branch. ", false);
36 addParam("isMom", m_isMom, "true for adding momentum tree branch. ", false);
37 addParam("isCharge", m_isCharge, "true for charge dedx tree branch. ", false);
38 addParam("isRun", m_isRun, "true for adding run number tree branch. ", false);
39 addParam("isWire", m_isWire, "true for adding wires tree branch. ", false);
40 addParam("isLayer", m_isLayer, "true for adding layers tree branch. ", false);
41 addParam("isDoca", m_isDoca, "true for adding doca tree branch. ", false);
42 addParam("isEnta", m_isEnta, "true for adding enta tree branch. ", false);
43 addParam("isInjTime", m_isInjTime, "true for adding time var tree branch. ", false);
44 addParam("isDocaRS", m_isDocaRS, "true for adding doca tree branch. ", false);
45 addParam("isEntaRS", m_isEntaRS, "true for adding enta tree branch. ", false);
46 addParam("isDedxhit", m_isDedxhit, "true for adding dedxhit tree branch. ", false);
47 addParam("isADCcorr", m_isADCcorr, "true for adding adc tree branch. ", false);
48 addParam("isBhabha", m_isBhabha, "true for bhabha events", true);
49 addParam("isRadee", m_isRadee, "true for radee events", false);
50 addParam("isTrgSel", m_isTrgSel, "true to enable trigger sel inside module", false);
51}
52
53//-----------------------------------------------------------------
54// Create ROOT objects
55//-----------------------------------------------------------------
56
58{
59 m_trgResult.isOptional();
60 m_dedxTracks.isRequired();
61 m_tracks.isRequired();
62 m_trackFitResults.isRequired();
63
64 // Data object creation
65 auto means = new TH1D("means", "CDC dE/dx truncated means", 250, 0, 2.5);
66 auto ttree = new TTree("tree", "Tree with dE/dx information");
67
68 auto hestats = new TH1I("hestats", "Event Stats", 6, -0.5, 5.5);
69 hestats->SetFillColor(kYellow);
70 hestats->GetXaxis()->SetBinLabel(1, "all");
71 hestats->GetXaxis()->SetBinLabel(2, "notrig");
72 hestats->GetXaxis()->SetBinLabel(3, "noskim");
73 hestats->GetXaxis()->SetBinLabel(4, "wrongskim");
74 hestats->GetXaxis()->SetBinLabel(5, "unclean");
75 hestats->GetXaxis()->SetBinLabel(6, "selected");
76
77 auto htstats = new TH1I("htstats", "track Stats", 7, -0.5, 6.5);
78 htstats->SetFillColor(kYellow);
79 htstats->GetXaxis()->SetBinLabel(1, "alltrk");
80 htstats->GetXaxis()->SetBinLabel(2, "vtx");
81 htstats->GetXaxis()->SetBinLabel(3, "inCDC");
82 htstats->GetXaxis()->SetBinLabel(4, "whits");
83 htstats->GetXaxis()->SetBinLabel(5, "weop");
84 htstats->GetXaxis()->SetBinLabel(6, "radee");
85 htstats->GetXaxis()->SetBinLabel(7, "selected");
86
87 if (m_isInjTime) {
88 ttree->Branch<double>("injtime", &m_injTime);
89 ttree->Branch<double>("injring", &m_injRing);
90 ttree->Branch<double>("costh", &m_costh);
91 ttree->Branch<int>("nhits", &m_nhits);
92 ttree->Branch<double>("p", &m_p);
93 }
94
95 ttree->Branch<double>("dedx", &m_dedx);
96 if (m_isCosth)ttree->Branch<double>("costh", &m_costh);
97 if (m_isMom)ttree->Branch<double>("p", &m_p);
98 if (m_isCharge)ttree->Branch<int>("charge", &m_charge);
99 if (m_isRun)ttree->Branch<int>("run", &m_run);
100
101 if (m_isWire)ttree->Branch("wire", &m_wire);
102 if (m_isLayer)ttree->Branch("layer", &m_layer);
103 if (m_isDoca)ttree->Branch("doca", &m_doca);
104 if (m_isEnta)ttree->Branch("enta", &m_enta);
105 if (m_isDocaRS)ttree->Branch("docaRS", &m_docaRS);
106 if (m_isEntaRS)ttree->Branch("entaRS", &m_entaRS);
107 if (m_isDedxhit)ttree->Branch("dedxhit", &m_dedxhit);
108 if (m_isADCcorr)ttree->Branch("adccorr", &m_adccorr);
109
110 // Collector object registration
111 registerObject<TH1D>("means", means);
112 registerObject<TTree>("tree", ttree);
113 registerObject<TH1I>("hestats", hestats);
114 registerObject<TH1I>("htstats", htstats);
115}
116
117//-----------------------------------------------------------------
118// Fill ROOT objects
119//-----------------------------------------------------------------
121{
122
123 auto hestats = getObjectPtr<TH1I>("hestats");
124 hestats->Fill(0);
125
126 if (m_isTrgSel) {
127 if (!m_trgResult.isValid()) {
128 B2WARNING("SoftwareTriggerResult required to select bhabha/radee event is not found");
129 hestats->Fill(1);
130 return;
131 }
132
133 //release05: bhabha_all is grand skim = bhabha+bhabhaecl+radee
134 const std::map<std::string, int>& fresults = m_trgResult->getResults();
135 if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end() and
136 fresults.find("software_trigger_cut&skim&accept_radee") == fresults.end()) {
137 B2WARNING("Can't find required bhabha/radee trigger identifiers");
138 hestats->Fill(2);
139 return;
140 }
141
142 const bool eBhabha = (m_trgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
144
145 const bool eRadBhabha = (m_trgResult->getResult("software_trigger_cut&skim&accept_radee") ==
147
148 if (!m_isBhabha && !m_isRadee) {
149 B2WARNING("requested not-supported event type: going back");
150 hestats->Fill(3);
151 return;
152 } else if (m_isBhabha && !m_isRadee && !eBhabha) {
153 B2WARNING("requested bhabha only but event not found: going back");
154 hestats->Fill(3);
155 return;
156 } else if (m_isRadee && !m_isBhabha && !eRadBhabha) {
157 B2WARNING("requested radee only but event not found: going back");
158 hestats->Fill(3);
159 return;
160 }
161 } else {
162 hestats->GetXaxis()->SetBinLabel(2, "inact1");
163 hestats->GetXaxis()->SetBinLabel(3, "inact2");
164 hestats->GetXaxis()->SetBinLabel(4, "inact3");
165 }
166
167 StoreObjPtr<EventMetaData> eventMetaDataPtr;
168 int run = eventMetaDataPtr->getRun();
169 if (m_isRun)m_run = run;
170 int nTracks = m_dedxTracks.getEntries();
171 if (nTracks >= 4) {
172 B2WARNING("too many tracks: unclean bhabha or radee event: " << nTracks);
173 hestats->Fill(4);
174 return;
175 }
176
177 hestats->Fill(5);
178
179 //Collector object access
180 auto tree = getObjectPtr<TTree>("tree");
181 auto htstats = getObjectPtr<TH1I>("htstats");
182 auto hmeans = getObjectPtr<TH1D>("means");
183
184 for (int idedx = 0; idedx < nTracks; idedx++) {
185
186 CDCDedxTrack* dedxTrack = m_dedxTracks[idedx];
187 if (!dedxTrack) {
188 B2WARNING("No dedx track: Going back: " << idedx);
189 continue;
190 }
191
192 const Track* track = dedxTrack->getRelatedFrom<Track>();
193 if (!track) {
194 B2WARNING("No track: Going back: " << idedx);
195 continue;
196 }
197
198 const TrackFitResult* fitResult = track->getTrackFitResultWithClosestMass(Const::pion);
199 if (!fitResult) {
200 B2WARNING("No related fit for this track...");
201 continue;
202 }
203
204 m_dedx = dedxTrack->getDedxNoSat();
205 m_p = dedxTrack->getMomentum();
206 m_costh = dedxTrack->getCosTheta();
207 m_charge = fitResult->getChargeSign();
208 m_injTime = dedxTrack->getInjectionTime();
209 m_injRing = dedxTrack->getInjectionRing();
210 m_nhits = dedxTrack->size();
211
212 htstats->Fill(0);
213
214 if (m_cuts) {
215 // apply cleanup cuts
216 if (fabs(fitResult->getD0()) >= 1.0)continue;
217 if (fabs(fitResult->getZ0()) >= 1.0) continue;
218 htstats->Fill(1);
219
220 //if outside CDC
221 if (m_costh < TMath::Cos(150.0 * TMath::DegToRad()))continue; //-0.866
222 if (m_costh > TMath::Cos(17.0 * TMath::DegToRad())) continue; //0.95
223 htstats->Fill(2);
224
225 if (m_nhits > m_maxHits) continue;
226
227 //making some cuts based on acceptance
228 if (m_costh > -0.55 && m_costh < 0.820) {
229 if (dedxTrack->getNLayerHits() < 25)continue; //all CDC layer available here
230 } else {
231 if (m_costh <= -0.62 || m_costh >= 0.880) {
232 if (dedxTrack->getNLayerHits() < 10)continue; //less layer available here
233 if (m_costh > 0 && dedxTrack->getNLayerHits() < 13)continue;
234 } else {
235 if (dedxTrack->getNLayerHits() < 18)continue;
236 }
237 }
238 htstats->Fill(3);
239
240 const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
241 if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
242 double TrkEoverP = (eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons)) / (fitResult->getMomentum().R());
243 if (abs(TrkEoverP - 1.0) > m_setEoP)continue;
244 }
245 htstats->Fill(4);
246 }
247
248 //if dealing with radee here (do a safe side cleanup)
249 if (m_isRadee) {
250 if (nTracks != 2)continue; //exactly 2 tracks
251 bool goodradee = false;
252 //checking if dedx of other track is restricted
253 //will not do too much as radee is clean enough
254 for (int jdedx = 0; jdedx < nTracks; jdedx++) {
255 CDCDedxTrack* dedxOtherTrack = m_dedxTracks[abs(jdedx - 1)];
256 if (!dedxOtherTrack)continue;
257 if (abs(dedxOtherTrack->getDedxNoSat() - 1.0) > 0.25)continue; //loose for uncalibrated
258 goodradee = true;
259 break;
260 }
261 if (!goodradee)continue;
262 htstats->Fill(5);
263 }
264
265
266 // Make sure to remove all the data in vectors from the previous track
267 if (m_isWire)m_wire.clear();
268 if (m_isLayer)m_layer.clear();
269 if (m_isDoca)m_doca.clear();
270 if (m_isEnta)m_enta.clear();
271 if (m_isDocaRS)m_docaRS.clear();
272 if (m_isEntaRS)m_entaRS.clear();
273 if (m_isDedxhit)m_dedxhit.clear();
274 if (m_isADCcorr)m_adccorr.clear();
275
276 // Simple numbers don't need to be cleared
277 // make sure to use the truncated mean without the hadron saturation correction
278
279 for (int i = 0; i < m_nhits; ++i) {
280 // if (m_DBWireGains->getWireGain(dedxTrack->getWire(i)) == 0)continue;
281 if (m_isWire)m_wire.push_back(dedxTrack->getWire(i));
282 if (m_isLayer)m_layer.push_back(dedxTrack->getHitLayer(i));
283 if (m_isDoca)m_doca.push_back(dedxTrack->getDoca(i));
284 if (m_isEnta)m_enta.push_back(dedxTrack->getEnta(i));
285 if (m_isDocaRS)m_docaRS.push_back(dedxTrack->getDocaRS(i) / dedxTrack->getCellHalfWidth(i));
286 if (m_isEntaRS)m_entaRS.push_back(dedxTrack->getEntaRS(i));
287 if (m_isDedxhit)m_dedxhit.push_back(dedxTrack->getDedx(i));
288 if (m_isADCcorr)m_adccorr.push_back(dedxTrack->getADCCount(i));
289 }
290
291 // Track and/or hit information filled as per config
292 htstats->Fill(6);
293 hmeans->Fill(m_dedx);
294 tree->Fill();
295 }
296}
bool m_isEntaRS
flag to write rescaled enta in tree
bool m_isLayer
flag to write layer number in tree
bool m_isADCcorr
flag to write adc corrected in tree
bool m_isWire
flag to write wire number in tree
std::vector< double > m_adccorr
adc corrected for the hit
std::vector< int > m_wire
hit level information
std::vector< int > m_layer
continuous layer number for the hit
bool m_isRun
flag to write run number in tree
StoreArray< TrackFitResult > m_trackFitResults
Required array for TrackFitResults.
bool m_isTrgSel
flag to enable trigger skim selected in the module (off default)
StoreObjPtr< SoftwareTriggerResult > m_trgResult
required input
std::vector< double > m_dedxhit
dE/dx for the hit
std::vector< double > m_doca
distance of closest approach for the hit
int m_nhits
number of dE/dx hits on the track
virtual void collect() override
Fill ROOT objects.
StoreArray< Track > m_tracks
Required array for Tracks.
bool m_isDedxhit
flag to write dedxhits in tree
virtual void prepare() override
Create and book ROOT objects.
CDCDedxElectronCollectorModule()
Constructor: Sets the description, the properties and the parameters of the module.
std::vector< double > m_entaRS
rescaled entrance angle for the hit
bool m_isInjTime
flag to enable trigger skim (off default)
std::vector< double > m_enta
entrance angle for the hit
StoreArray< CDCDedxTrack > m_dedxTracks
Required array for CDCDedxTracks.
bool m_isDocaRS
flag to write rescaled doca in tree
std::vector< double > m_docaRS
rescaled distance of closest approach for the hit
Debug output for CDCDedxPID module.
Definition: CDCDedxTrack.h:25
int getHitLayer(int i) const
Return the (global) layer number for a hit.
Definition: CDCDedxTrack.h:209
int getADCCount(int i) const
Return the adcCount for this hit.
Definition: CDCDedxTrack.h:218
double getDoca(int i) const
Return the distance of closest approach to the sense wire for this hit.
Definition: CDCDedxTrack.h:227
double getDedx() const
Get dE/dx truncated mean for this track.
Definition: CDCDedxTrack.h:103
double getEntaRS(int i) const
Return rescaled enta value for cell height=width assumption.
Definition: CDCDedxTrack.h:236
double getCellHalfWidth(int i) const
Return the half-width of the CDC cell.
Definition: CDCDedxTrack.h:248
int getNLayerHits() const
Return the number of layer hits for this track.
Definition: CDCDedxTrack.h:174
double getCosTheta() const
Return cos(theta) for this track.
Definition: CDCDedxTrack.h:121
double getInjectionRing() const
Return cos(theta) for this track.
Definition: CDCDedxTrack.h:130
int getWire(int i) const
Return the sensor ID for this hit: wire number for CDC (0-14336)
Definition: CDCDedxTrack.h:206
double getDocaRS(int i) const
Return rescaled doca value for cell height=width assumption.
Definition: CDCDedxTrack.h:233
double getEnta(int i) const
Return the entrance angle in the CDC cell for this hit.
Definition: CDCDedxTrack.h:230
double getInjectionTime() const
Return cos(theta) for this track.
Definition: CDCDedxTrack.h:127
double getDedxNoSat() const
Get dE/dx truncated mean without the saturation correction for this track.
Definition: CDCDedxTrack.h:106
int size() const
Return the number of hits for this track.
Definition: CDCDedxTrack.h:200
double getMomentum() const
Return the track momentum valid in the CDC.
Definition: CDCDedxTrack.h:118
Calibration collector module base class.
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
ECL cluster data.
Definition: ECLCluster.h:27
bool hasHypothesis(EHypothesisBit bitmask) const
Return if specific hypothesis bit is set.
Definition: ECLCluster.h:351
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Definition: ECLCluster.cc:23
@ c_nPhotons
CR is split into n photons (N1)
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Class that bundles various TrackFitResults.
Definition: Track.h:25
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
@ c_accept
Accept this event.
Abstract base class for different kinds of events.