Belle II Software development
CDCCRTestModule.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 <framework/core/HistoModule.h>
12
13#include <cdc/dataobjects/WireID.h>
14
15#include <string>
16
17#include <mdst/dataobjects/Track.h>
18#include <mdst/dataobjects/TrackFitResult.h>
19#include <tracking/dataobjects/RecoTrack.h>
20#include <framework/dataobjects/EventT0.h>
21#include <framework/datastore/StoreArray.h>
22#include <framework/geometry/B2Vector3.h>
23#include "TH1.h"
24#include "TH2.h"
25#include "TProfile.h"
26#include "TTree.h"
27
28namespace Belle2 {
34 namespace CDC {
35
40
41 public:
42
47
51 virtual ~CDCCRTestModule();
52
56 void initialize() override;
57
61 void beginRun() override;
62
67 void event() override;
72 void endRun() override;
73
77 void terminate() override;
78
83 void defineHisto() override;
84
85 private:
86
90 TH1* getHist(const char* name, const char* title,
91 int nBins, double x0, double x1)
92 {
93 TH1* h = new TH1D(name, title, nBins, x0, x1);
94 m_allHistos.push_back(h);
95 return h;
96 }
97
101 TProfile* getHistProfile(const char* name, const char* title,
102 int nBins, double x0, double x1)
103 {
104 TProfile* h = new TProfile(name, title, nBins, x0, x1);
105 m_allHistos.push_back(h);
106 return h;
107 }
108
112 TH2* getHist(const char* name, const char* title,
113 int nBinsX, double x0, double x1,
114 int nBinsY, double y0, double y1)
115 {
116 TH2* h = new TH2D(name, title, nBinsX, x0, x1, nBinsY, y0, y1);
117 m_allHistos.push_back(h);
118 return h;
119 }
120
124 TH1* getHist(const std::string& name, const std::string& title,
125 int nBins, double x0, double x1)
126 {
127 return getHist(name.c_str(), title.c_str(), nBins, x0, x1);
128 }
129
133 TProfile* getHistProfile(const std::string& name, const std::string& title,
134 int nBins, double x0, double x1)
135 {
136 return getHistProfile(name.c_str(), title.c_str(), nBins, x0, x1);
137 }
138
142 TH2* getHist(const std::string& name, const std::string& title,
143 int nBinsX, double x0, double x1,
144 int nBinsY, double y0, double y1)
145 {
146 return getHist(name.c_str(), title.c_str(), nBinsX, x0, x1, nBinsY, y0, y1);
147 }
148
149
153 const genfit::SharedPlanePtr constructPlane(const genfit::MeasuredStateOnPlane& state, WireID m_wireID);
154
159
163 void plotResults(Belle2::RecoTrack* track);
164
168 void getHitDistInTrackCand(const RecoTrack* track);//Draw hit distribution from track candidate
169
174
178 void HitEfficiency(const Belle2::RecoTrack* track);
179
183 int getICLayer(int slayer, int ilayer)
184 {
185 if (slayer == 0) {return ilayer;}
186 else {return 8 + (slayer - 1) * 6 + ilayer;}
187 }
188
195
198
201
204
207
208 std::string m_trackArrayName;
209 std::string m_cdcHitArrayName ;
210 std::string m_recoTrackArrayName ;
214 std::vector<TH1*> m_allHistos;
215 std::string m_treeName;
217 TTree* m_tree;
221 TH1* m_hNDF;
222 TH1* m_hNHits;
224 TH1* m_hChi2;
225 TH1* m_hPval;
228 TH1* m_hAlpha;
229 TH1* m_hPhi0;
230 TH1* m_hTheta;
234 TH1* m_hResidualU[56];
235 TH1* m_hEvtT0;
239 TH2* m_hDxDt[56];
240 TProfile* m_hHitEff_soft[56];
247 double res_b;
248 double res_u;
249 double res_b_err;
250 double res_u_err;
251 double weight;
252 double absRes_b;
253 double absRes_u;
254 double alpha;
255 double theta;
256 unsigned short adc;
257 short tdc;
258 double t;
259 double t_fit;
260 double dt_flight;
262 double dt_prop;
263 double evtT0;
264 double Pt;
266 double x_mea;
267 double x_u;
268 double x_b;
269 double x_sim;
270 double z;
271 double z_sim;
272 double z_prop;
273 int lay;
274 int IWire;
275 int lr;
278 double Pval;
279 double TrPval;
280 double ndf;
281 double d0;
282 double z0;
283 double phi0;
284 double tanL;
285 double omega;
286 double m_MinimumPt;
289 std::vector<double> m_TriggerPos;
290 std::vector<double> m_TriggerPlaneDirection;
291 std::vector<double> m_TriggerSize;
292 std::vector<int> m_up;
293 std::vector<int> m_low;
302 bool m_noBFit;
303 bool m_ToP;
304 bool m_ToF;
311 };
312 }
314}
CDC Cosmic test calibration module.
StoreObjPtr< EventT0 > m_eventTimeStoreObject
Event timing.
std::string m_recoTrackArrayName
Belle2::RecoTrack StoreArray nam.e.
double m_MinimumPt
Minimum Transverse momentum of tracks.
TProfile * getHistProfile(const char *name, const char *title, int nBins, double x0, double x1)
Create profile plot.
int trighit
Trigger hit information.
bool m_calExpectedDriftTime
Calculate expected drift time from x_fit or not.
double x_sim
Simulation DriftLength .
TTree * m_tree
Output tree recording the information of each hit.
double res_b
Biased residual.
double alpha
Entrance Azimuthal angle of hit (degree).
void getHitDistInTrackCand(const RecoTrack *track)
Make hit distribution from track candidate.
TH1 * m_hNHits
Number of Hits per track.
bool m_plotResidual
Process track to get the hit information of fitted track.
double res_u_err
Unbiased residual error.
TH1 * m_hNHits_trackcand
Number of Hits per trackCand.
B2Vector3D m_trigHitPos
Trigger position.
void initialize() override
Initializes the Module.
int getICLayer(int slayer, int ilayer)
Convert slayer and ilayer to iclayer.
StoreArray< TrackFitResult > m_TrackFitResults
Track fit results.
double absRes_b
absolute value of biased residual.
bool m_StoreTrackParams
Store Track parameter or not.
void event() override
Event action (main routine).
TProfile * getHistProfile(const std::string &name, const std::string &title, int nBins, double x0, double x1)
Create profile plot.
double dt_flight
Time of flight.
bool m_SmallerOutput
make output smaller by ignore some variable.
void HitEfficiency(const Belle2::RecoTrack *track)
Cal Hit eff.
TH1 * m_hHitDistInTrack[56]
Hit Dist.
std::string m_cdcHitArrayName
Belle2::CDCHit StoreArray name.
void endRun() override
End run action.
double res_u
Unbiased residual.
TH1 * m_hPhi0
Phi0 of ttrack, see Helix.
void terminate() override
Termination action.
TH1 * m_hNormalizedResidualU[56]
Residual distribution normalized with tracking error.
TH1 * m_hNTracks
Number of track fitted, Convergence, not conv, not fit.
TH1 * m_hHitDistInCDCHit[56]
Hit Dist.
double z0
Track Parameter, z0.
bool m_hitEfficiency
calculate hit eff or not, Haven't finished.
std::string m_relRecoTrackTrackName
Relation between RecoTrack and Belle2:Track.
std::vector< int > m_low
lower channel list for each board.
void plotResults(Belle2::RecoTrack *track)
Plot track parameters and related variables.
const genfit::SharedPlanePtr constructPlane(const genfit::MeasuredStateOnPlane &state, WireID m_wireID)
Construct a plane for the hit.
TH2 * m_hNDFNormalizedResidualU[56]
Normalized residual vs.
double Pval
P-value of fitted track.
double z_prop
Propagation Length along the sense wire.
double res_b_err
Biased residual error.
double t
Measurement Drift time.
TH2 * m_hNDFPval
Degree-of-freedom vs Probability histo.
TH1 * m_hTheta
Theta of each Hit.
TH1 * m_hNTracksPerEvent
Number of TrackCand per Event.
void beginRun() override
Begin run action.
std::string m_histogramDirectoryName
subdir where to place the histograms.
void getResidualOfUnFittedLayer(Belle2::RecoTrack *track)
Calculate residual for Layers which didn't use int fitting.
double omega
Track Parameter, omega.
double weight
Weight of hit.
double Pt
Transverse momentum.
double x_b
X_fit for biased track fit.
double t_fit
Drift time calculated from x_fit.
bool m_fillExpertHistos
Fill some histogram for monitoring fit quality.
TH1 * m_hAlpha
Alpha of each Hit.
std::string m_trackArrayName
Belle2::Track StoreArray name.
double z_sim
Z of hit on wire (simulation).
const Belle2::TrackFitResult * fitresult
Track fit result.
double TrPval
P-value of fitted track.
bool m_MakeHitDist
Switch to make histograms of hit distribution.
TH2 * m_hNDFChi2
Chi2 vs degree-of-freedom histo.
virtual ~CDCCRTestModule()
Destructor.
bool m_noBFit
fit incase no magnetic Field of not, if true, NDF=4 in cal P-value
bool m_EstimateResultForUnFittedLayer
Calculate residual for layer that we do not use in track fitting.
int boardID
Electrical Board ID.
double z
Z of hit on wire.
std::vector< TH1 * > m_allHistos
A list of 1d histograms.
TProfile * m_hHitEff_soft[56]
Hit efficience of each layer, software.
TH2 * m_h2DHitDistInTrack
2D Hit Dist..(ICLay vs IWire) have weight>0.5 after fit with DAF
double trigHitPos_x
X-position of track at trigger counter.
std::string m_trackFitResultArrayName
Belle2::TrackFitResult StoreArray name.
double dt_prop
Time of propagation.
bool m_ToP
Enable to correct ToP if true.
double x_u
X_fit for unbiased track fit.
double theta
Entrance Polar angle of hit (degree).
std::vector< int > m_up
upper channel list for each board.
TH1 * m_hNTracksPerEventFitted
Number of TrackCand per Event.
bool m_ToF
Enable to correct ToF if true.
TH1 * m_hPval
Fit Probability histo.
B2Vector3D getTriggerHitPosition(Belle2::RecoTrack *track)
extrapolation track to trigger counter plane (y position).
std::string m_treeName
Name of tree for the output file.
double d0
Track Parameter, d0.
StoreArray< RecoTrack > m_RecoTracks
Tracks.
TH1 * m_hHitDistInTrCand[56]
Hit Dist.
TH2 * getHist(const std::string &name, const std::string &title, int nBinsX, double x0, double x1, int nBinsY, double y0, double y1)
Create 2d-histogram.
StoreArray< Track > m_Tracks
Tracks.
TH2 * m_h2DHitDistInTrCand
2D Hit Dist.
std::vector< double > m_TriggerPlaneDirection
Nominal center position of trigger counter.
TH1 * getHist(const std::string &name, const std::string &title, int nBins, double x0, double x1)
Create 1d-histogram.
TH2 * m_hDxDt[56]
Unbiased x_fit vs.
TH2 * m_hNDFResidualU[56]
Residual vs.
double tanL
Track Parameter, tanL.
unsigned short adc
adc value.
bool m_StoreCDCSimHitInfo
Store CDCSimHit Information.
double phi0
Track Parameter, phi0.
double trigHitPos_z
Z-position of track at trigger counter.
std::vector< double > m_TriggerPos
Nominal center position of trigger counter.
bool m_EventT0Extraction
use Event T0 extract t0 or not.
double ndf
degree of freedom.
double absRes_u
absolute value of unbiased residual.
TH1 * m_hNDF
Number of Degree Freedom.
TH1 * m_hResidualU[56]
Residual distribution (in cm)
TH2 * getHist(const char *name, const char *title, int nBinsX, double x0, double x1, int nBinsY, double y0, double y1)
Create 2d-histogram.
TH2 * m_hTriggerHitZX
Trigger hit image.
std::vector< double > m_TriggerSize
Size of trigger counter (Width x length).
TH2 * m_h2DHitDistInCDCHit
2D Hit Dist.
TH1 * getHist(const char *name, const char *title, int nBins, double x0, double x1)
Create 1D histogram.
double dt_flight_sim
Time of flight (Simulation).
void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
StoreArray< CDCHit > m_CDCHits
CDC hits.
double x_mea
measure drift length (signed by left right).
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
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.
Class to identify a wire inside the CDC.
Definition: WireID.h:34
Abstract base class for different kinds of events.