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