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