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", 6, -0.5, 5.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, "selected");
91
92 if (m_isInjTime) {
93 ttree->Branch<double>("injtime", &m_injTime);
94 ttree->Branch<double>("injring", &m_injRing);
95 ttree->Branch<double>("costh", &m_costh);
96 ttree->Branch<int>("nhits", &m_nhits);
97 ttree->Branch<double>("p", &m_p);
98 }
99
100 ttree->Branch<double>("dedx", &m_dedx);
101 if (m_isCosth)ttree->Branch<double>("costh", &m_costh);
102 if (m_isMom)ttree->Branch<double>("p", &m_p);
103 if (m_isPt)ttree->Branch<double>("pt", &m_pt);
104 if (m_isCharge)ttree->Branch<int>("charge", &m_charge);
105 if (m_isRun)ttree->Branch<int>("run", &m_run);
106
107 if (m_isWire)ttree->Branch("wire", &m_wire);
108 if (m_isLayer)ttree->Branch("layer", &m_layer);
109 if (m_isDoca)ttree->Branch("doca", &m_doca);
110 if (m_isEnta)ttree->Branch("enta", &m_enta);
111 if (m_isDocaRS)ttree->Branch("docaRS", &m_docaRS);
112 if (m_isEntaRS)ttree->Branch("entaRS", &m_entaRS);
113 if (m_isDedxhit)ttree->Branch("dedxhit", &m_dedxhit);
114 if (m_isADCcorr)ttree->Branch("adccorr", &m_adccorr);
115
116 if (m_islLayer)ttree->Branch("lLayer", &m_lLayer);
117 if (m_islDedx)ttree->Branch("ldedx", &m_ldedx);
118
119 // Collector object registration
120 registerObject<TH1D>("means", means);
121 registerObject<TTree>("tree", ttree);
122 registerObject<TH1I>("hestats", hestats);
123 registerObject<TH1I>("htstats", htstats);
124}
125
126//-----------------------------------------------------------------
127// Fill ROOT objects
128//-----------------------------------------------------------------
130{
131
132 auto hestats = getObjectPtr<TH1I>("hestats");
133 hestats->Fill(0);
134
135 if (m_isTrgSel) {
136 if (!m_trgResult.isValid()) {
137 B2WARNING("SoftwareTriggerResult required to select bhabha/radee event is not found");
138 hestats->Fill(1);
139 return;
140 }
141
142 //release05: bhabha_all is grand skim = bhabha+bhabhaecl+radee
143 const std::map<std::string, int>& fresults = m_trgResult->getResults();
144 if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end() and
145 fresults.find("software_trigger_cut&skim&accept_bhabha_cdc") == fresults.end()) {
146 B2WARNING("Can't find required bhabha/radee trigger identifiers");
147 hestats->Fill(2);
148 return;
149 }
150
151 const bool eBhabha = (m_trgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
153
154 const bool eRadBhabha = (m_trgResult->getResult("software_trigger_cut&skim&accept_bhabha_cdc") ==
156
157 if (!m_isBhabha && !m_isRadee) {
158 B2WARNING("requested not-supported event type: going back");
159 hestats->Fill(3);
160 return;
161 } else if (m_isBhabha && !m_isRadee && !eBhabha) {
162 B2WARNING("requested bhabha only but event not found: going back");
163 hestats->Fill(3);
164 return;
165 } else if (m_isRadee && !m_isBhabha && !eRadBhabha) {
166 B2WARNING("requested radee only but event not found: going back");
167 hestats->Fill(3);
168 return;
169 }
170 } else {
171 hestats->GetXaxis()->SetBinLabel(2, "inact1");
172 hestats->GetXaxis()->SetBinLabel(3, "inact2");
173 hestats->GetXaxis()->SetBinLabel(4, "inact3");
174 }
175
176 StoreObjPtr<EventMetaData> eventMetaDataPtr;
177 int run = eventMetaDataPtr->getRun();
178 if (m_isRun)m_run = run;
179 int nTracks = m_dedxTracks.getEntries();
180 if (nTracks >= 4) {
181 B2WARNING("too many tracks: unclean bhabha or radee event: " << nTracks);
182 hestats->Fill(4);
183 return;
184 }
185
186 hestats->Fill(5);
187
188 //Collector object access
189 auto tree = getObjectPtr<TTree>("tree");
190 auto htstats = getObjectPtr<TH1I>("htstats");
191 auto hmeans = getObjectPtr<TH1D>("means");
192
193 for (int idedx = 0; idedx < nTracks; idedx++) {
194
195 CDCDedxTrack* dedxTrack = m_dedxTracks[idedx];
196 if (!dedxTrack) {
197 B2WARNING("No dedx track: Going back: " << idedx);
198 continue;
199 }
200
201 const Track* track = dedxTrack->getRelatedFrom<Track>();
202 if (!track) {
203 B2WARNING("No track: Going back: " << idedx);
204 continue;
205 }
206
207 const TrackFitResult* fitResult = track->getTrackFitResultWithClosestMass(Const::pion);
208 if (!fitResult) {
209 B2WARNING("No related fit for this track...");
210 continue;
211 }
212
213 m_dedx = dedxTrack->getDedxNoSat();
214 m_p = dedxTrack->getMomentum();
215 m_costh = dedxTrack->getCosTheta();
216 m_charge = fitResult->getChargeSign();
217 m_pt = fitResult->getTransverseMomentum();
218 m_injTime = dedxTrack->getInjectionTime();
219 m_injRing = dedxTrack->getInjectionRing();
220 m_nhits = dedxTrack->size();
221 m_nlhits = dedxTrack->getNLayerHits();
222
223 htstats->Fill(0);
224
225 if (m_cuts) {
226 // apply cleanup cuts
227 if (std::abs(fitResult->getD0()) >= 1.0)continue;
228 if (std::abs(fitResult->getZ0()) >= 1.0) continue;
229 htstats->Fill(1);
230
231 //if outside CDC
232 if (m_costh < TMath::Cos(150.0 * TMath::DegToRad()))continue; //-0.866
233 if (m_costh > TMath::Cos(17.0 * TMath::DegToRad())) continue; //0.95
234 htstats->Fill(2);
235
236 if (m_nhits > m_maxHits) continue;
237
238 //making some cuts based on acceptance
239 if (m_costh > -0.55 && m_costh < 0.820) {
240 if (dedxTrack->getNLayerHits() < 25)continue; //all CDC layer available here
241 } else {
243 if (dedxTrack->getNLayerHits() < 10)continue; //less layer available here
244 if (m_costh > 0 && dedxTrack->getNLayerHits() < 13)continue;
245 } else {
246 if (dedxTrack->getNLayerHits() < 18)continue;
247 }
248 }
249 htstats->Fill(3);
250
251 const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
252 if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
253 double TrkEoverP = (eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons)) / (fitResult->getMomentum().R());
254 if (std::abs(TrkEoverP - 1.0) > m_setEoP)continue;
255 }
256 htstats->Fill(4);
257 }
258
259 // Make sure to remove all the data in vectors from the previous track
260 if (m_isWire)m_wire.clear();
261 if (m_isLayer)m_layer.clear();
262 if (m_isDoca)m_doca.clear();
263 if (m_isEnta)m_enta.clear();
264 if (m_isDocaRS)m_docaRS.clear();
265 if (m_isEntaRS)m_entaRS.clear();
266 if (m_isDedxhit)m_dedxhit.clear();
267 if (m_isADCcorr)m_adccorr.clear();
268
269 // Simple numbers don't need to be cleared
270 // make sure to use the truncated mean without the hadron saturation correction
271
272 for (int i = 0; i < m_nhits; ++i) {
273 // if (m_DBWireGains->getWireGain(dedxTrack->getWire(i)) == 0)continue;
274 if (m_isWire)m_wire.push_back(dedxTrack->getWire(i));
275 if (m_isLayer)m_layer.push_back(dedxTrack->getHitLayer(i));
276 if (m_isDoca)m_doca.push_back(dedxTrack->getDoca(i));
277 if (m_isEnta)m_enta.push_back(dedxTrack->getEnta(i));
278 if (m_isDocaRS)m_docaRS.push_back(dedxTrack->getDocaRS(i) / dedxTrack->getCellHalfWidth(i));
279 if (m_isEntaRS)m_entaRS.push_back(dedxTrack->getEntaRS(i));
280 if (m_isDedxhit)m_dedxhit.push_back(dedxTrack->getDedx(i));
281 if (m_isADCcorr)m_adccorr.push_back(dedxTrack->getADCCount(i));
282 }
283
284 // Make sure to remove all the data in vectors from the previous track
285 if (m_islLayer)m_lLayer.clear();
286 if (m_islDedx)m_ldedx.clear();
287
288 for (int i = 0; i < m_nlhits; ++i) {
289 if (m_islLayer)m_lLayer.push_back(dedxTrack->getLayer(i));
290 if (m_islDedx)m_ldedx.push_back(dedxTrack->getLayerDedx(i));
291 }
292 // Track and/or hit information filled as per config
293 htstats->Fill(5);
294 hmeans->Fill(m_dedx);
295 tree->Fill();
296 }
297}
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.