Belle II Software  release-06-02-00
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 #include <framework/dataobjects/EventMetaData.h>
12 #include <mdst/dataobjects/ECLCluster.h>
13 
14 #include <TTree.h>
15 #include <TH1D.h>
16 #include <TH1I.h>
17 #include <TMath.h>
18 
19 using namespace Belle2;
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_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  addParam("cleanupCuts", m_cuts, "Boolean to apply cleanup cuts", true);
35  addParam("maxNumHits", m_maxNumHits,
36  "Maximum number of hits per track. If there is more than this the track will not be collected. ", int(100));
37  addParam("fSetEoverP", fSetEoverP, "Set E over p Cut values. ", double(0.25));
38  addParam("IsCosth", IsCosth, "true for adding costh tree branch. ", false);
39  addParam("IsMom", IsMom, "true for adding momentum tree branch. ", false);
40  addParam("IsCharge", IsCharge, "true for charge dedx tree branch. ", false);
41  addParam("IsRun", IsRun, "true for adding run number tree branch. ", false);
42  addParam("IsWire", IsWire, "true for adding wires tree branch. ", false);
43  addParam("IsLayer", IsLayer, "true for adding layers tree branch. ", false);
44  addParam("IsDoca", IsDoca, "true for adding doca tree branch. ", false);
45  addParam("IsEnta", IsEnta, "true for adding enta tree branch. ", false);
46  addParam("IsDocaRS", IsDocaRS, "true for adding doca tree branch. ", false);
47  addParam("IsEntaRS", IsEntaRS, "true for adding enta tree branch. ", false);
48  addParam("Isdedxhit", Isdedxhit, "true for adding dedxhit tree branch. ", false);
49  addParam("IsBhabhaEvt", IsBhabhaEvt, "true for bhabha events", true);
50  addParam("IsRadBhabhaEvt", IsRadBhabhaEvt, "true for radee events", false);
51  addParam("enableTrgSel", enableTrgSel, "true to enable trigger sel inside module", false);
52 }
53 
54 //-----------------------------------------------------------------
55 // Create ROOT objects
56 //-----------------------------------------------------------------
57 
59 {
60  m_TrgResult.isOptional();
61  m_dedxTracks.isRequired();
62  m_tracks.isRequired();
63  m_trackFitResults.isRequired();
64 
65  // Data object creation
66  auto means = new TH1D("means", "CDC dE/dx truncated means", 250, 0, 2.5);
67  auto ttree = new TTree("tree", "Tree with dE/dx information");
68 
69  auto hestats = new TH1I("hestats", "Event Stats", 6, -0.5, 5.5);
70  hestats->SetFillColor(kYellow);
71  hestats->GetXaxis()->SetBinLabel(1, "all");
72  hestats->GetXaxis()->SetBinLabel(2, "notrig");
73  hestats->GetXaxis()->SetBinLabel(3, "noskim");
74  hestats->GetXaxis()->SetBinLabel(4, "wrongskim");
75  hestats->GetXaxis()->SetBinLabel(5, "unclean");
76  hestats->GetXaxis()->SetBinLabel(6, "selected");
77 
78  auto htstats = new TH1I("htstats", "track Stats", 7, -0.5, 6.5);
79  htstats->SetFillColor(kYellow);
80  htstats->GetXaxis()->SetBinLabel(1, "alltrk");
81  htstats->GetXaxis()->SetBinLabel(2, "vtx");
82  htstats->GetXaxis()->SetBinLabel(3, "inCDC");
83  htstats->GetXaxis()->SetBinLabel(4, "whits");
84  htstats->GetXaxis()->SetBinLabel(5, "weop");
85  htstats->GetXaxis()->SetBinLabel(6, "radee");
86  htstats->GetXaxis()->SetBinLabel(7, "selected");
87 
88  ttree->Branch<double>("dedx", &m_dedx);
89  if (IsCosth)ttree->Branch<double>("costh", &m_costh);
90  if (IsMom)ttree->Branch<double>("p", &m_p);
91  if (IsCharge)ttree->Branch<int>("charge", &m_charge);
92  if (IsRun)ttree->Branch<int>("run", &m_run);
93 
94  if (IsWire)ttree->Branch("wire", &m_wire);
95  if (IsLayer)ttree->Branch("layer", &m_layer);
96  if (IsDoca)ttree->Branch("doca", &m_doca);
97  if (IsEnta)ttree->Branch("enta", &m_enta);
98  if (IsDocaRS)ttree->Branch("docaRS", &m_docaRS);
99  if (IsEntaRS)ttree->Branch("entaRS", &m_entaRS);
100  if (Isdedxhit)ttree->Branch("dedxhit", &m_dedxhit);
101 
102  // Collector object registration
103  registerObject<TH1D>("means", means);
104  registerObject<TTree>("tree", ttree);
105  registerObject<TH1I>("hestats", hestats);
106  registerObject<TH1I>("htstats", htstats);
107 }
108 
109 //-----------------------------------------------------------------
110 // Fill ROOT objects
111 //-----------------------------------------------------------------
113 {
114 
115  auto hestats = getObjectPtr<TH1I>("hestats");
116  hestats->Fill(0);
117 
118  if (enableTrgSel) {
119  if (!m_TrgResult.isValid()) {
120  B2WARNING("SoftwareTriggerResult required to select bhabha/radee event is not found");
121  hestats->Fill(1);
122  return;
123  }
124 
125  //release05: bhabha_all is grand skim = bhabha+bhabhaecl+radee
126  const std::map<std::string, int>& fresults = m_TrgResult->getResults();
127  if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end() and
128  fresults.find("software_trigger_cut&skim&accept_radee") == fresults.end()) {
129  B2WARNING("Can't find required bhabha/radee trigger identifiers");
130  hestats->Fill(2);
131  return;
132  }
133 
134  const bool eBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
136 
137  const bool eRadBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_radee") ==
139 
140  if (!IsBhabhaEvt && !IsRadBhabhaEvt) {
141  B2WARNING("requested not-supported event type: going back");
142  hestats->Fill(3);
143  return;
144  } else if (IsBhabhaEvt && !IsRadBhabhaEvt && !eBhabha) {
145  B2WARNING("requested bhabha only but event not found: going back");
146  hestats->Fill(3);
147  return;
148  } else if (IsRadBhabhaEvt && !IsBhabhaEvt && !eRadBhabha) {
149  B2WARNING("requested radee only but event not found: going back");
150  hestats->Fill(3);
151  return;
152  }
153  } else {
154  hestats->GetXaxis()->SetBinLabel(2, "inact1");
155  hestats->GetXaxis()->SetBinLabel(3, "inact2");
156  hestats->GetXaxis()->SetBinLabel(4, "inact3");
157  }
158 
159  StoreObjPtr<EventMetaData> eventMetaDataPtr;
160  int run = eventMetaDataPtr->getRun();
161  if (IsRun)m_run = run;
162 
163  int nTracks = m_dedxTracks.getEntries();
164  if (nTracks >= 4) {
165  B2WARNING("too many tracks: unclean bhabha or radee event: " << nTracks);
166  hestats->Fill(4);
167  return;
168  }
169 
170  hestats->Fill(5);
171 
172  //Collector object access
173  auto tree = getObjectPtr<TTree>("tree");
174  auto htstats = getObjectPtr<TH1I>("htstats");
175  auto hmeans = getObjectPtr<TH1D>("means");
176 
177  for (int idedx = 0; idedx < nTracks; idedx++) {
178 
179  CDCDedxTrack* dedxTrack = m_dedxTracks[idedx];
180  if (!dedxTrack) {
181  B2WARNING("No dedx track: Going back: " << idedx);
182  continue;
183  }
184 
185  const Track* track = dedxTrack->getRelatedFrom<Track>();
186  if (!track) {
187  B2WARNING("No track: Going back: " << idedx);
188  continue;
189  }
190 
191  const TrackFitResult* fitResult = track->getTrackFitResultWithClosestMass(Const::pion);
192  if (!fitResult) {
193  B2WARNING("No related fit for this track...");
194  continue;
195  }
196 
197  m_dedx = dedxTrack->getDedxNoSat();
198  m_p = dedxTrack->getMomentum();
199  m_costh = dedxTrack->getCosTheta();
200  m_charge = fitResult->getChargeSign();
201  htstats->Fill(0);
202 
203  if (m_cuts) {
204  // apply cleanup cuts
205  if (fabs(fitResult->getD0()) >= 1.0)continue;
206  if (fabs(fitResult->getZ0()) >= 1.0) continue;
207  htstats->Fill(1);
208 
209  //if outside CDC
210  if (m_costh < TMath::Cos(150.0 * TMath::DegToRad()))continue; //-0.866
211  if (m_costh > TMath::Cos(17.0 * TMath::DegToRad())) continue; //0.95
212  htstats->Fill(2);
213 
214  //making some cuts based on acceptance
215  if (m_costh > -0.55 && m_costh < 0.820) {
216  if (dedxTrack->getNLayerHits() < 25)continue; //all CDC layer available here
217  } else {
218  if (m_costh <= -0.62 || m_costh >= 0.880) {
219  if (dedxTrack->getNLayerHits() < 10)continue; //less layer available here
220  if (m_costh > 0 && dedxTrack->getNLayerHits() < 13)continue;
221  } else {
222  if (dedxTrack->getNLayerHits() < 18)continue;
223  }
224  }
225  htstats->Fill(3);
226 
227  const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
228  if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
229  double TrkEoverP = (eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons)) / (fitResult->getMomentum().Mag());
230  if (abs(TrkEoverP - 1.0) > fSetEoverP)continue;
231  }
232  htstats->Fill(4);
233 
234  //if dealing with radee here (do a safe side cleanup)
235  if (IsRadBhabhaEvt) {
236  if (nTracks != 2)continue; //exactly 2 tracks
237  bool goodradee = false;
238  //checking if dedx of other track is restricted
239  //will not do too much as radee is clean enough
240  for (int jdedx = 0; jdedx < nTracks; jdedx++) {
241  CDCDedxTrack* dedxOtherTrack = m_dedxTracks[abs(jdedx - 1)];
242  if (!dedxOtherTrack)continue;
243  if (abs(dedxOtherTrack->getDedxNoSat() - 1.0) > 0.25)continue; //loose for uncalibrated
244  goodradee = true;
245  break;
246  }
247  if (!goodradee)continue;
248  htstats->Fill(5);
249  }
250  }
251 
252 
253  // Make sure to remove all the data in vectors from the previous track
254  if (IsWire)m_wire.clear();
255  if (IsLayer)m_layer.clear();
256  if (IsDoca)m_doca.clear();
257  if (IsEnta)m_enta.clear();
258  if (IsDocaRS)m_docaRS.clear();
259  if (IsEntaRS)m_entaRS.clear();
260  if (Isdedxhit)m_dedxhit.clear();
261 
262  // Simple numbers don't need to be cleared
263  // make sure to use the truncated mean without the hadron saturation correction
264  m_nhits = dedxTrack->size();
265  if (m_nhits > m_maxNumHits) continue;
266 
267  for (int i = 0; i < m_nhits; ++i) {
268  // if (m_DBWireGains->getWireGain(dedxTrack->getWire(i)) == 0)continue;
269  if (IsWire)m_wire.push_back(dedxTrack->getWire(i));
270  if (IsLayer)m_layer.push_back(dedxTrack->getHitLayer(i));
271  if (IsDoca)m_doca.push_back(dedxTrack->getDoca(i));
272  if (IsEnta)m_enta.push_back(dedxTrack->getEnta(i));
273  if (IsDocaRS)m_docaRS.push_back(dedxTrack->getDocaRS(i) / dedxTrack->getCellHalfWidth(i));
274  if (IsEntaRS)m_entaRS.push_back(dedxTrack->getEntaRS(i));
275  if (Isdedxhit)m_dedxhit.push_back(dedxTrack->getDedx(i));
276  }
277 
278  // Track and/or hit information filled as per config
279  htstats->Fill(6);
280  hmeans->Fill(m_dedx);
281  tree->Fill();
282  }
283 }
A collector module for CDC dE/dx electron calibrations.
bool IsWire
flag to write wire number in tree
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.
std::vector< double > m_dedxhit
dE/dx for the hit
std::vector< double > m_doca
distance of closest approach for the hit
bool IsDocaRS
flag to write rescaled doca in tree
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 enableTrgSel
flag to enable trigger skim selected in the module (off deafult)
bool IsLayer
flag to write layer number in tree
bool IsRun
flag to write run number in tree
virtual void prepare() override
Create and book ROOT objects.
bool IsEntaRS
flag to write rescaled enta in tree
std::vector< double > m_entaRS
rescaled entrance angle for the hit
std::vector< double > m_enta
entrance angle for the hit
StoreArray< CDCDedxTrack > m_dedxTracks
Required array for CDCDedxTracks.
std::vector< double > m_docaRS
rescaled distance of closest approach for the hit
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
required input
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:194
double getDoca(int i) const
Return the distance of closest approach to the sense wire for this hit.
Definition: CDCDedxTrack.h:212
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:221
double getCellHalfWidth(int i) const
Return the half-width of the CDC cell.
Definition: CDCDedxTrack.h:233
int getNLayerHits() const
Return the number of layer hits for this track.
Definition: CDCDedxTrack.h:159
double getCosTheta() const
Return cos(theta) for this track.
Definition: CDCDedxTrack.h:121
int getWire(int i) const
Return the sensor ID for this hit: wire number for CDC (0-14336)
Definition: CDCDedxTrack.h:191
double getDocaRS(int i) const
Return rescaled doca value for cell height=width assumption.
Definition: CDCDedxTrack.h:218
double getEnta(int i) const
Return the entrance angle in the CDC cell for this hit.
Definition: CDCDedxTrack.h:215
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:185
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:542
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:19
@ c_nPhotons
CR is split into n photons (N1)
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:95
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
Class that bundles various TrackFitResults.
Definition: Track.h:25
#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.