Belle II Software  release-08-01-10
DQMEventProcessorBase.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 <tracking/dqmUtils/DQMEventProcessorBase.h>
10 #include <tracking/dqmUtils/DQMHistoModuleBase.h>
11 
12 #include <framework/datastore/StoreArray.h>
13 #include <vxd/geometry/GeoTools.h>
14 #include <vxd/geometry/SensorInfoBase.h>
15 
16 using namespace Belle2;
17 
19 {
21  if (!recoTracks.isOptional()) {
22  B2DEBUG(22, "Missing recoTracks array in event() for " + m_histoModule->getName() + " module.");
23  return;
24  }
26  if (!tracks.isOptional()) {
27  B2DEBUG(22, "Missing recoTracks array in event() for " + m_histoModule->getName() + " module.");
28  return;
29  }
30 
31  try {
32  m_iTrack = 0;
33  m_iTrackVXD = 0;
34  m_iTrackCDC = 0;
35  m_iTrackVXDCDC = 0;
36 
37  for (const Track& track : tracks) {
38  ProcessTrack(track);
39  }
40 
42  } catch (const std::exception& e) {
43  B2WARNING("Exception " + std::string(e.what()) + " in " + m_histoModule->getName() + " module!");
44  } catch (...) {
45  B2WARNING("Unhandled exception in " + m_histoModule->getName() + " module!");
46  }
47 }
48 
50 {
51  auto recoTracksVector = track.getRelationsTo<RecoTrack>(m_recoTracksStoreArrayName);
52  if (!recoTracksVector.size())
53  return;
54 
55  m_recoTrack = recoTracksVector[0];
56 
57  int nPXDClusters = 0;
58  if (!m_runningOnHLT) {
59  RelationVector<PXDCluster> pxdClusters = DataStore::getRelationsWithObj<PXDCluster>(m_recoTrack);
60  nPXDClusters = (int)pxdClusters.size();
61  }
62  RelationVector<SVDCluster> svdClusters = DataStore::getRelationsWithObj<SVDCluster>(m_recoTrack);
63  int nSVDClusters = (int)svdClusters.size();
64  RelationVector<CDCHit> cdcHits = DataStore::getRelationsWithObj<CDCHit>(m_recoTrack);
65  int nCDCHits = (int)cdcHits.size();
66 
67  // This method allways returns TrackFitResult so there's no need to check if it's not nullptr
68  auto trackFitResult = track.getTrackFitResultWithClosestMass(Const::pion);
69 
70  TString message = ConstructMessage(trackFitResult, nPXDClusters, nSVDClusters, nCDCHits);
71  B2DEBUG(20, message.Data());
72 
73  FillTrackFitResult(trackFitResult);
74  m_histoModule->FillHitNumbers(nPXDClusters, nSVDClusters, nCDCHits);
75 
78 
79  m_iTrack++;
80  if (((nPXDClusters > 0) || (nSVDClusters > 0)) && (nCDCHits == 0))
81  m_iTrackVXD++;
82  if (((nPXDClusters == 0) && (nSVDClusters == 0)) && (nCDCHits > 0))
83  m_iTrackCDC++;
84  if (((nPXDClusters > 0) || (nSVDClusters > 0)) && (nCDCHits > 0))
86 }
87 
88 TString DQMEventProcessorBase::ConstructMessage(const TrackFitResult* trackFitResult, int nPXDClusters, int nSVDClusters,
89  int nCDCHits)
90 {
91  return Form("%s: track %3i, Mom: %f, %f, %f, Pt: %f, Mag: %f, Hits: PXD %i SVD %i CDC %i Suma %i\n",
92  m_histoModule->getName().c_str(),
93  m_iTrack,
94  (float)trackFitResult->getMomentum().x(),
95  (float)trackFitResult->getMomentum().y(),
96  (float)trackFitResult->getMomentum().z(),
97  (float)trackFitResult->getMomentum().Rho(),
98  (float)trackFitResult->getMomentum().R(),
99  nPXDClusters, nSVDClusters, nCDCHits, nPXDClusters + nSVDClusters + nCDCHits
100  );
101 }
102 
104 {
105  if (m_runningOnHLT) {
106  m_histoModule->FillMomentumAngles(trackFitResult);
107  m_histoModule->FillMomentumCoordinates(trackFitResult);
108  }
110 }
111 
113 {
114  // function wasFitSuccessful already checked if TrackFitStatus is not nullptr so it's not necessary to do so
116 
117  m_isNotFirstHit = false;
118 
119  for (auto recoHitInfo : m_recoTrack->getRecoHitInformations(true)) {
120  ProcessRecoHit(recoHitInfo);
121  }
122 }
123 
125 {
126  if (!recoHitInfo) {
127  B2DEBUG(20, "Missing genfit::pxd recoHitInfo in event() for " + m_histoModule->getName() + " module.");
128  return;
129  }
130 
131  if (!recoHitInfo->useInFit())
132  return;
133 
134  bool isPXD = recoHitInfo->getTrackingDetector() == RecoHitInformation::c_PXD;
135  bool isSVD = recoHitInfo->getTrackingDetector() == RecoHitInformation::c_SVD;
136  if (!isPXD && !isSVD)
137  return;
138 
139  auto fitterInfo = m_recoTrack->getCreatedTrackPoint(recoHitInfo)->getFitterInfo();
140  if (!fitterInfo)
141  return;
142 
143  m_rawSensorResidual = new TVectorT<double>(fitterInfo->getResidual(0, false).getState());
144 
145  if (isPXD && ! m_runningOnHLT) {
146  ProcessPXDRecoHit(recoHitInfo);
147  } else if (isSVD) {
148  ProcessSVDRecoHit(recoHitInfo);
149  }
150 
151  delete m_rawSensorResidual;
152 }
153 
155 {
156  m_position.SetX(recoHitInfo->getRelatedTo<PXDCluster>()->getU());
157  m_position.SetY(recoHitInfo->getRelatedTo<PXDCluster>()->getV());
158  m_residual_um.SetX(m_rawSensorResidual->GetMatrixArray()[0] / Unit::um);
159  m_residual_um.SetY(m_rawSensorResidual->GetMatrixArray()[1] / Unit::um);
160 
161  m_sensorID = recoHitInfo->getRelatedTo<PXDCluster>()->getSensorID();
163 
165 
168 
169  if (m_produce1Dres)
171  if (m_produce2Dres)
173 
175 }
176 
178 {
179  if (recoHitInfo->getRelatedTo<SVDCluster>()->isUCluster()) {
180  m_position.SetX(recoHitInfo->getRelatedTo<SVDCluster>()->getPosition());
181  m_residual_um.SetX(m_rawSensorResidual->GetMatrixArray()[0] / Unit::um);
182  } else {
183  m_position.SetY(recoHitInfo->getRelatedTo<SVDCluster>()->getPosition());
184  m_residual_um.SetY(m_rawSensorResidual->GetMatrixArray()[0] / Unit::um);
185  }
186 
187  m_sensorID = recoHitInfo->getRelatedTo<SVDCluster>()->getSensorID();
188  if (m_sensorIDPrev == m_sensorID) {
190 
191  if (! m_runningOnHLT) {
193 
196 
197  if (m_produce1Dres)
199  if (m_produce2Dres)
201  }
202 
204  }
205 
207 }
208 
210 {
211  auto sensorInfo = &VXD::GeoCache::get(m_sensorID);
212  m_globalResidual_um = sensorInfo->vectorToGlobal(m_residual_um, true);
213  ROOT::Math::XYZVector globalPosition = sensorInfo->pointToGlobal(m_position, true);
214 
215  m_phi_deg = globalPosition.Phi() / Unit::deg;
216  m_theta_deg = globalPosition.Theta() / Unit::deg;
217 
219 
220  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
221  m_layerIndex = gTools->getLayerIndex(m_layerNumber);
222  m_correlationIndex = m_layerIndex - gTools->getFirstLayer();
223  m_sensorIndex = gTools->getSensorIndex(m_sensorID);
224 }
225 
227 {
228  if (m_isNotFirstHit && ((m_layerNumber - m_layerNumberPrev) == 1)) {
230  } else {
231  m_isNotFirstHit = true;
232  }
233 
235 }
236 
238 {
242 }
243 
244 bool DQMEventProcessorBase::IsNotYang(int ladderNumber, int layerNumber)
245 {
246  switch (layerNumber) {
247  case 1:
248  return ladderNumber < 5 || ladderNumber > 8;
249  case 2:
250  return ladderNumber < 7 || ladderNumber > 12;
251  default:
252  return true;
253  }
254 }
255 
256 bool DQMEventProcessorBase::IsNotMat(int ladderNumber, int layerNumber)
257 {
258  switch (layerNumber) {
259  case 3:
260  return ladderNumber < 3 || ladderNumber > 5;
261  case 4:
262  return ladderNumber < 4 || ladderNumber > 8;
263  case 5:
264  return ladderNumber < 5 || ladderNumber > 10;
265  case 6:
266  return ladderNumber < 6 || ladderNumber > 13;
267  default:
268  return true;
269  }
270 }
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
virtual void ProcessSVDRecoHit(RecoHitInformation *recoHitInfo)
Compute position in a SVD way which means we need two consecutive hits to be from the same sensor to ...
ROOT::Math::XYZVector m_position
local coordinates of the hit position (u, v, w)
bool m_produce1Dres
if true, produce 1D Track residuals plots for each VXD sensor
static bool IsNotMat(int ladderNumber, int layerNumber)
Returns true if sensor with given ladderNumber and layerNumber isn't in the Mat half-shell,...
virtual void ProcessTrack(const Track &track)
Find RecoTrack for given track.
float m_phiPrev_deg
global phi in degrees of the previous hit
int m_layerIndex
index of the layer of the hit
std::string m_tracksStoreArrayName
StoreArray name where Tracks are written.
virtual void FillCommonHistograms()
Fill histograms which are common for PXD and SVD hit.
int m_layerNumberPrev
number of the layer of the previous hit
virtual void ProcessPXDRecoHit(RecoHitInformation *recoHitInfo)
Compute position in a PXD way.
int m_layerNumber
number of the layer of the hit
virtual void FillTrackFitResult(const TrackFitResult *trackFitResult)
Fill histograms with values derived from TrackFitResult.
int m_iTrackVXD
index of track where are VXD hits and aren't CDC hits (with valid TrackFitResult and related RecoTrac...
int m_iTrackCDC
index of track where are CDC hits and aren't VXD hits (with valid TrackFitResult and related RecoTrac...
bool m_produce2Dres
if true, produce 2D Track residuals plots for each VXD sensor
ROOT::Math::XYZVector m_residual_um
unbiased residual for the hit in micrometers in local coordinates (u, v, w)
bool m_isNotFirstHit
Determines if the hit is not the first hit in the current track.
int m_correlationIndex
index of the layer of the previous hit
int m_sensorIndex
index of the sensor of the hit
virtual void Run()
Call this to start processing the event data and filling histograms.
virtual void ComputeCommonVariables()
Compute variables which are common for PXD and SVD hit.
virtual void ProcessSuccessfulFit()
Continue track processing by calling ProcessRecoHit function on each RecoHitInformation in given Reco...
virtual void ProcessRecoHit(RecoHitInformation *recoHitInfo)
Compute unbiased residual and the calls ProcesPXDRecoHit or ProcessSVDRecoHit.
float m_theta_deg
global theta in degrees of the hit
TVectorT< double > * m_rawSensorResidual
unbiased residual for the hit obtained from the sensor so its length is different for PXD and SVD sen...
float m_phi_deg
global phi in degrees of the hit
virtual void SetCommonPrevVariables()
Set the value of -Prev values which are common for PXD and SVD hit.
bool m_runningOnHLT
true if the DQM is run on HLT
static bool IsNotYang(int ladderNumber, int layerNumber)
Returns true if sensor with given ladderNumber and layerNumber isn't in the Yang half-shell,...
std::string m_recoTracksStoreArrayName
StoreArray name where RecoTracks are written.
float m_thetaPrev_deg
global theta in degrees of the previous hit
VxdID m_sensorID
ID of the current sensor.
virtual TString ConstructMessage(const TrackFitResult *trackFitResult, int nPXDClusters, int nSVDClusters, int nCDCHits)
Make debug message with information about RecoTrack.
int m_iTrackVXDCDC
index of track where are both VXD hits and CDC hits (with valid TrackFitResult and related RecoTrack)
RecoTrack * m_recoTrack
RecoTrack related to currently processed Track.
DQMHistoModuleBase * m_histoModule
DQM histogram module on which the Fill- functions are called to fill histograms.
int m_iTrack
index of track (with valid TrackFitResult and related RecoTrack)
VxdID m_sensorIDPrev
ID of the prewious sensor.
ROOT::Math::XYZVector m_globalResidual_um
unbiased residual for the hit in micrometers in global coordinates (x, y, z)
virtual void FillHitNumbers(int nPXD, int nSVD, int nCDC)
Fill histograms with numbers of hits.
virtual void FillHalfShellsSVD(const B2Vector3D &globalResidual_um, bool isNotMat)
Fill histograms with unbiased residuals for half-shells for SVD sensors.
virtual void FillUB2DResidualsSensor(const B2Vector3D &residual_um, int sensorIndex)
Fill 2D histograms with unbiased residuals for individual sensors.
virtual void FillUBResidualsPXD(const B2Vector3D &residual_um)
Fill histograms with unbiased residuals in PXD sensors.
virtual void FillTrackFitStatus(const genfit::FitStatus *tfs)
Fill histograms which require FitStatus.
virtual void FillTRClusterHitmap(float phi_deg, float theta_deg, int layerIndex)
Fill cluster hitmap in IP angle range.
virtual void FillHelixParametersAndCorrelations(const TrackFitResult *tfr)
Fill histograms with helix parameters and their correlations.
virtual void FillMomentumCoordinates(const TrackFitResult *tfr)
Fill histograms with track momentum Pt.
virtual void FillHalfShellsPXD(const B2Vector3D &globalResidual_um, bool isNotYang)
Fill histograms with unbiased residuals for half-shells for PXD sensors.
virtual void FillTrackIndexes(int iTrack, int iTrackVXD, int iTrackCDC, int iTrackVXDCDC)
Fill histograms with track indexes.
virtual void FillTRClusterCorrelations(float phi_deg, float phiPrev_deg, float theta_deg, float thetaPrev_deg, int correlationIndex)
Fill histograms with correlations between neighbor layers.
virtual void FillMomentumAngles(const TrackFitResult *tfr)
Fill histograms with track momentum Pt.
virtual void FillUBResidualsSVD(const B2Vector3D &residual_um)
Fill histograms with unbiased residuals in SVD sensors.
virtual void FillUB1DResidualsSensor(const B2Vector3D &residual_um, int sensorIndex)
Fill 1D histograms with unbiased residuals for individual sensors.
const std::string & getName() const
Returns the name of the module.
Definition: Module.h:187
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
float getV() const
Get v coordinate of hit position.
Definition: PXDCluster.h:136
float getU() const
Get u coordinate of hit position.
Definition: PXDCluster.h:131
This class stores additional information to every CDC/SVD/PXD hit stored in a RecoTrack.
RecoHitDetector getTrackingDetector() const
Get the detector this hit comes from.
bool useInFit() const
Get the flag, whether this his should be used in a fit or not.
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
const genfit::FitStatus * getTrackFitStatus(const genfit::AbsTrackRep *representation=nullptr) const
Return the track fit status for the given representation or for the cardinal one. You are not allowed...
Definition: RecoTrack.h:621
bool wasFitSuccessful(const genfit::AbsTrackRep *representation=nullptr) const
Returns true if the last fit with the given representation was successful.
Definition: RecoTrack.cc:336
const genfit::TrackPoint * getCreatedTrackPoint(const RecoHitInformation *recoHitInformation) const
Get a pointer to the TrackPoint that was created from this hit.
Definition: RecoTrack.cc:230
std::vector< RecoHitInformation * > getRecoHitInformations(bool getSorted=false) const
Return a list of all RecoHitInformations associated with the RecoTrack.
Definition: RecoTrack.cc:557
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
bool isUCluster() const
Get the direction of strips.
Definition: SVDCluster.h:110
float getPosition(double v=0) const
Get the coordinate of reconstructed hit.
Definition: SVDCluster.h:117
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Values of the result of a track fit with a given particle hypothesis.
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
static const double deg
degree to radians
Definition: Unit.h:109
static const double um
[micrometers]
Definition: Unit.h:71
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition: GeoCache.h:147
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
Definition: TrackPoint.cc:170
Abstract base class for different kinds of events.