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