Belle II Software development
HitLevelInfoWriter.h
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#pragma once
10
11#include <mdst/dataobjects/Track.h>
12#include <mdst/dataobjects/TrackFitResult.h>
13#include <mdst/dataobjects/ECLCluster.h>
14#include <mdst/dataobjects/KLMCluster.h>
15#include <mdst/dataobjects/EventLevelTriggerTimeInfo.h>
16
17#include <framework/datastore/StoreArray.h>
18#include <framework/datastore/StoreObjPtr.h>
19#include <framework/database/DBObjPtr.h>
20#include <framework/core/Module.h>
21#include <framework/gearbox/Const.h>
22
23#include <cdc/dataobjects/CDCDedxTrack.h>
24#include <cdc/dbobjects/CDCDedxScaleFactor.h>
25#include <cdc/dbobjects/CDCDedxWireGain.h>
26#include <cdc/dbobjects/CDCDedxRunGain.h>
27#include <cdc/dbobjects/CDCDedxCosineCor.h>
28#include <cdc/dbobjects/CDCDedx2DCell.h>
29#include <cdc/dbobjects/CDCDedx1DCell.h>
30#include <cdc/dbobjects/CDCDedxInjectionTime.h>
31#include <cdc/dbobjects/CDCDedxADCNonLinearity.h>
32#include <cdc/dbobjects/CDCDedxCosineEdge.h>
33#include <cdc/dbobjects/CDCDedxHadronCor.h>
34
35#include <vector>
36
37#include <TFile.h>
38#include <TTree.h>
39
40class TH2F;
41
42namespace Belle2 {
47
48 class CDCDedxTrack;
49
53
54 public:
55
58
61
63 virtual void initialize() override;
64
67 virtual void event() override;
68
70 virtual void terminate() override;
71
73 void bookOutput(std::string filename);
74
78 void recalculateDedx(CDCDedxTrack* dedxTrack, std::map<int, std::vector<double>>& l_var,
79 double (&cdcChi)[Const::ChargedStable::c_SetSize]);
80
84 double GetCorrection(int& adc, int layer, int wireID, double doca, double enta, double costheta, double ring, double time) const;
85
89 void HadronCorrection(double costheta, double& dedx) const;
90
94 double D2I(const double cosTheta, const double D) const;
95
99 double I2D(const double cosTheta, const double I) const;
100
108 void calculateMeans(double* mean, double* truncatedMean, double* truncatedMeanErr, const std::vector<double>& dedx) const;
109
113 void saveChiValue(double(&chi)[Const::ChargedStable::c_SetSize], CDCDedxTrack* dedxTrack, double dedx) const;
114
115 private:
116
118 std::vector<std::string> m_strParticleList;
119 std::vector<std::string> m_filename;
120 std::vector<TFile* > m_file;
121 std::vector<TTree* > m_tree;
122
129
131 void fillTrack(const TrackFitResult* fitResult);
132
134 void fillDedx(CDCDedxTrack* dedxTrack);
135
137 void clearEntries();
138
139 // event level information (from emd)
140 int m_expID{ -1};
141 int m_runID{ -1};
142 int m_eventID{ -1};
143
144 double m_injring{ -1.0};
145 double m_injtime{ -1.0};
146
147 // track level information (from tfr)
148 double m_d0{0.};
149 double m_z0{0.};
150 double m_dz{ -1.};
151 double m_dr{ -1.};
152 double m_dphi{ -1.};
153 double m_vx0{0.};
154 double m_vy0{0.};
155 double m_vz0{0.};
156 double m_tanlambda{ -1.};
157 double m_phi0{ -1.};
158 double m_chi2{ -1.};
159
160 double m_nCDChits{ -1.};
161 int m_inCDC{ -1};
162 int m_trackID{ -1};
163 double m_length{ -1.};
164 int m_charge{0};
165 double m_cosTheta{ -2.};
166 double m_pCDC{ -1.};
167 double m_p{ -1.};
168 double m_pt{ -1.};
169 double m_phi{ -1.};
170 double m_ioasym{ -1.};
171 double m_theta{ -2.};
172
173 // track level Mc
174 double m_PDG{ -1.};
175 // double m_motherPDG; /**< MC PID of mother particle */
176 // double m_pTrue; /**< MC true momentum */
177 // double m_trackDist; /**< the total distance traveled by the track */
178
179 // track level dE/dx measurements
180 double m_mean{ -1.};
181 double m_trunc{ -1.};
182 double m_truncNoSat{ -1.};
183 double m_error{ -1.};
184
185 // other dec specific information
186 double m_eop{ -1.};
187 double m_e{ -1.};
188 double m_e1_9{ -1.};
189 double m_e9_21{ -1.};
190 double m_klmLayers{ -1.};
191 double m_eclsnHits{ -1.};
192
193 // calibration constants
194 double m_scale{ -1.};
195 double m_cosCor{ -1.};
196 double m_cosEdgeCor{ -1.};
197 double m_runGain{ -1.};
198 double m_timeGain{ -1.};
199 double m_timeReso{ -1.};
200
201 // hadron cal and PID related variables
202 double m_chieOld{ -1.};
203 double m_chimuOld{ -1.};
204 double m_chipiOld{ -1.};
205 double m_chikOld{ -1.};
206 double m_chipOld{ -1.};
207 double m_chidOld{ -1.};
208
209 double m_chie{ -1.};
210 double m_chimu{ -1.};
211 double m_chipi{ -1.};
212 double m_chik{ -1.};
213 double m_chip{ -1.};
214 double m_chid{ -1.};
215
216 double m_prese{ -1.};
217 double m_presmu{ -1.};
218 double m_prespi{ -1.};
219 double m_presk{ -1.};
220 double m_presp{ -1.};
221 double m_presd{ -1.};
222
223 double m_pmeane{ -1.};
224 double m_pmeanmu{ -1.};
225 double m_pmeanpi{ -1.};
226 double m_pmeank{ -1.};
227 double m_pmeanp{ -1.};
228 double m_pmeand{ -1.};
229
230 static const int kMaxHits = 200;
231 // layer level information
232 int l_nhits{ -1};
233 int l_nhitsused{ -1};
236 int l_layer[kMaxHits] = {};
237 double l_path[kMaxHits] = {};
238 double l_dedx[kMaxHits] = {};
239
240 // hit level information (references on nhits)
241 int h_nhits{ -1};
242 int h_lwire[kMaxHits] = {};
243 int h_wire[kMaxHits] = {};
244 int h_layer[kMaxHits] = {};
245
246 double h_path[kMaxHits] = {};
247 double h_dedx[kMaxHits] = {};
248 double h_adcraw[kMaxHits] = {};
249 double h_adccorr[kMaxHits] = {};
250 double h_doca[kMaxHits] = {};
251 double h_ndoca[kMaxHits] = {};
252 double h_ndocaRS[kMaxHits] = {};
253 double h_enta[kMaxHits] = {};
254 double h_entaRS[kMaxHits] = {};
255 double h_driftT[kMaxHits] = {};
256 double h_driftD[kMaxHits] = {};
257 double h_facnladc[kMaxHits] = {};
258 double h_wireGain[kMaxHits] = {};
259 double h_twodCor[kMaxHits] = {};
260 double h_onedCor[kMaxHits] = {};
261
262 //Tracking variables for extra hits
267
268 // parameters: calibration constants
279
280 std::vector<double> m_hadronpars;
281
283 //Flag to enable and disable set of variables
288 };
289
290} // Belle2 namespace
Debug output for CDCDedxPID module.
static const unsigned int c_SetSize
Number of elements (for use in array bounds etc.)
Definition Const.h:615
Class for accessing objects in the database.
Definition DBObjPtr.h:21
double h_ndoca[kMaxHits]
normalized distance of closest approach
double h_driftD[kMaxHits]
drift distance
double h_onedCor[kMaxHits]
calibration 1D cleanup correction
int l_nhits
the total number of layer hits for this Track
double m_cosEdgeCor
calibration cosine edge correction
double m_presk
pred reso for kaon hypothesis
double m_phi0
Angle of the transverse momentum in the r-phi plane.
double m_d0
Signed distance to the POCA in the r-phi plane.
double m_presp
pred reso for proton hypothesis
double D2I(const double cosTheta, const double D) const
hadron saturation parameterization part 2
int l_wirelongesthit[kMaxHits]
the wire number of longest hit in this layer
double m_ioasym
asymmetry in increasing vs decreasing layer numbers per track
StoreObjPtr< EventLevelTriggerTimeInfo > m_TTDInfo
Store Object Ptr: EventLevelTriggerTimeInfo.
std::vector< double > m_hadronpars
hadron saturation parameters
double m_p
momentum from tracking
double m_e9_21
ratio of energies of the central 3x3 crystal vs 5x5 crystals
int m_trackID
ID number of the Track.
DBObjPtr< CDCDedxRunGain > m_DBRunGain
Run gain DB object.
double m_vx0
X coordinate of track POCA to origin.
virtual ~HitLevelInfoWriterModule()
Destructor.
StoreArray< KLMCluster > m_klmClusters
Required array of input KLMClusters.
double m_pmeand
pred mean for deuteron hypothesis
virtual void initialize() override
Initialize the module.
DBObjPtr< CDCDedxHadronCor > m_DBHadronCor
hadron saturation parameters
double h_enta[kMaxHits]
entrance angle
double m_chipOld
chi value for proton hypothesis
std::string m_strOutputBaseName
Base name for the output ROOT files.
virtual void event() override
This method is called for each event.
double h_WeightKaonHypo[kMaxHits]
weight for kaon hypothesis from KalmanFitterInfo
DBObjPtr< CDCDedxCosineEdge > m_DBCosEdgeCor
cosine edge calibration
HitLevelInfoWriterModule()
Default constructor.
double m_pmeanp
pred mean for proton hypothesis
double I2D(const double cosTheta, const double I) const
hadron saturation parameterization part 1
double h_WeightProtonHypo[kMaxHits]
weight for proton hypothesis from KalmanFitterInfo
double m_chid
modified chi value for deuteron hypothesis
double m_chikOld
chi value for kaon hypothesis
double m_dz
vertex or POCA in case of tracks z in respect to IPs
bool m_isHitLevel
Flag to switch on/off hit level info.
DBObjPtr< CDCDedxADCNonLinearity > m_DBNonlADC
hadron saturation non linearity
double m_pmeank
pred mean for kaon hypothesis
double h_facnladc[kMaxHits]
calibration hit gain
DBObjPtr< CDCDedx1DCell > m_DB1DCell
1D correction DB object
void saveChiValue(double(&chi)[Const::ChargedStable::c_SetSize], CDCDedxTrack *dedxTrack, double dedx) const
for all particles, save chi values into 'chi' chi array of chi values to be modified
StoreArray< TrackFitResult > m_trackFitResults
Required array of input TrackFitResults.
int h_nhits
the number of good hits for this Track
std::vector< std::string > m_filename
full names of the output ROOT files
virtual void terminate() override
End of the event processing.
double h_adccorr[kMaxHits]
charge per hit corr by nonlinear ADC
std::vector< TTree * > m_tree
output ROOT trees
double m_presmu
pred reso for muon hypothesis
double m_cosCor
calibration cosine correction
double m_tanlambda
Slope of the track in the r-z plane.
int h_layer[kMaxHits]
layer number
double h_adcraw[kMaxHits]
charge per hit
int m_expID
experiment in which this Track was found
void fillTrack(const TrackFitResult *fitResult)
Fill the TTree with the information from the track fit.
double h_path[kMaxHits]
path length in cell
double m_chimu
modified chi value for muon hypothesis
double m_vy0
Y coordinate of track POCA to origin.
double m_chik
modified chi value for kaon hypothesis
void recalculateDedx(CDCDedxTrack *dedxTrack, std::map< int, std::vector< double > > &l_var, double(&cdcChi)[Const::ChargedStable::c_SetSize])
Function to recalculate the dedx with latest constants.
double m_chi2
chi^2 from track fit
void HadronCorrection(double costheta, double &dedx) const
Function to apply the hadron correction.
void calculateMeans(double *mean, double *truncatedMean, double *truncatedMeanErr, const std::vector< double > &dedx) const
Save arithmetic and truncated mean for the 'dedx' values.
double h_wireGain[kMaxHits]
calibration hit gain
double m_theta
cos(theta) for the track
double h_twodCor[kMaxHits]
calibration 2D correction
double m_chie
modified chi value for electron hypothesis
DBObjPtr< CDCDedx2DCell > m_DB2DCell
2D correction DB object
double m_runGain
calibration run gain
int m_eventID
event in which this Track was found
double m_timeGain
calibration injection time gain
void bookOutput(std::string filename)
Create the output TFiles and TTrees.
double m_chidOld
chi value for deuteron hypothesis
bool m_isExtraVar
Flag to switch on/off extra level info and some available w/ release/5 only.
void clearEntries()
Clear the arrays before filling an event.
std::vector< TFile * > m_file
output ROOT files
double m_dr
track d0 relative to IP
double m_prese
pred reso for electron hypothesis
DBObjPtr< CDCDedxCosineCor > m_DBCosineCor
Electron saturation correction DB object.
int h_foundByTrackFinder[kMaxHits]
the 'found by track finder' flag for the given hit
double m_timeReso
calibration injection time reso
StoreArray< Track > m_tracks
Required array of input Tracks.
double l_dedx[kMaxHits]
dE/dx for this layer
double l_path[kMaxHits]
distance travelled in this layer
int l_layer[kMaxHits]
layer number
double m_chipiOld
chi value for pion hypothesis
double m_chip
modified chi value for proton hypothesis
bool m_isCorrection
Flag to switch on/off corrections.
double m_pmeane
pred mean for electron hypothesis
double m_cosTheta
cos(theta) for the track
int l_nhitsused
the total number of layer hits used for this Track
StoreArray< ECLCluster > m_eclClusters
Required array of input ECLClusters.
double h_dedx[kMaxHits]
charge per path length
double m_chieOld
chi value for electron hypothesis
double h_driftT[kMaxHits]
drift time
double m_dphi
POCA in degrees in respect to IP.
DBObjPtr< CDCDedxWireGain > m_DBWireGains
Wire gain DB object.
double m_eop
energy over momentum in the calorimeter
bool m_isDeadwire
write only active wires
double m_presd
pred reso for deuteron hypothesis
double m_scale
calibration scale factor
double m_klmLayers
number of klm layers with hits
int l_nhitscombined[kMaxHits]
the number of hits combined this layer
double m_error
standard deviation of the truncated mean
double m_pCDC
momentum valid in CDC
double m_truncNoSat
dE/dx averaged, truncated mean, with corrections (not hadron)
double m_length
total path length of the Track
double m_injring
HER injection status.
double m_nCDChits
Number of CDC hits associated to the track.
static const int kMaxHits
default hit level index
double h_ndocaRS[kMaxHits]
normalized +RS distance of closest approach
double m_vz0
Z coordinate of track POCA to origin.
double m_trunc
dE/dx averaged, truncated mean, with corrections
double GetCorrection(int &adc, int layer, int wireID, double doca, double enta, double costheta, double ring, double time) const
Function to get the correction factor.
int m_inCDC
frack is CDC acceptance or not
double h_doca[kMaxHits]
distance of closest approach
double m_chimuOld
chi value for muon hypothesis
double m_z0
z coordinate of the POCA
int m_charge
the charge for this Track
double h_WeightPionHypo[kMaxHits]
weight for pion hypothesis from KalmanFitterInfo
int m_runID
run in which this Track was found
double m_e
energy in the calorimeter
StoreArray< CDCDedxTrack > m_dedxTracks
Required array of CDCDedxTracks.
double m_injtime
time since last injection in micro seconds
double m_pmeanmu
pred mean for muon hypothesis
double h_entaRS[kMaxHits]
normalized + RS distance of entrance angle
bool m_isRelative
Flag to switch on/off relative constants.
double m_prespi
pred reso for pion hypothesis
DBObjPtr< CDCDedxInjectionTime > m_DBInjectTime
time gain/reso DB object
int h_wire[kMaxHits]
sense wire ID
DBObjPtr< CDCDedxScaleFactor > m_DBScaleFactor
Scale factor to make electrons ~1.
int h_lwire[kMaxHits]
sense wire within layer
double m_chipi
modified chi value for pion hypothesis
std::vector< std::string > m_strParticleList
Vector of ParticleLists to write out.
double m_pt
transverse momentum from tracking
double m_pmeanpi
pred mean for pion hypothesis
void fillDedx(CDCDedxTrack *dedxTrack)
Fill the TTree with the information from a CDCDedxTrack object.
double m_e1_9
ratio of energies of the central 1 crystal vs 3x3 crystals
Module()
Constructor.
Definition Module.cc:30
Accessor to arrays stored in the data store.
Definition StoreArray.h:113
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.
Abstract base class for different kinds of events.