Belle II Software development
CDCCalibrationCollector.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 "cdc/modules/cdcCalibrationCollector/CDCCalibrationCollector.h"
10#include <cdc/translators/RealisticTDCCountTranslator.h>
11#include <framework/datastore/RelationArray.h>
12
13#include <tracking/trackFindingCDC/eventdata/hits/CDCWireHit.h>
14#include <tracking/trackFindingCDC/topology/CDCWireTopology.h>
15
16#include <genfit/TrackPoint.h>
17#include <genfit/KalmanFitterInfo.h>
18#include <genfit/MeasurementOnPlane.h>
19#include <genfit/MeasuredStateOnPlane.h>
20
21#include <Math/ProbFuncMathCore.h>
22
23#include <cdc/dataobjects/WireID.h>
24#include <cdc/geometry/CDCGeometryPar.h>
25
26#include <TH1F.h>
27
28//using namespace std;
29using namespace Belle2;
30using namespace CDC;
31using namespace genfit;
32using namespace TrackFindingCDC;
33
34
35REG_MODULE(CDCCalibrationCollector);
36
37
39{
40 setDescription("Collector module for cdc calibration");
41 setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
42 addParam("recoTracksColName", m_recoTrackArrayName, "Name of collection hold genfit::Track", std::string(""));
43 addParam("bField", m_bField, "If true -> #Params ==5 else #params ==4 for calculate P-Val", false);
44 addParam("calExpectedDriftTime", m_calExpectedDriftTime, "if true module will calculate expected drift time, it take a time",
45 true);
46 addParam("storeTrackParams", m_storeTrackParams, "Store Track Parameter or not, it will be multicount for each hit", false);
47 addParam("eventT0Extraction", m_eventT0Extraction, "use event t0 extract t0 or not", true);
48 addParam("minimumPt", m_minimumPt, "Tracks with tranverse momentum smaller than this value will not used", 0.15);
49 addParam("minimumNDF", m_minimumNDF, "Discard tracks whose degree-of-freedom below this value", 5.);
50 addParam("isCosmic", m_isCosmic, "True when we process cosmic events, else False (collision)", m_isCosmic);
51 addParam("effStudy", m_effStudy, "When true module collects info only necessary for wire eff study", false);
52}
53
55{
56}
57
59{
60 m_Tracks.isRequired(m_trackArrayName);
63 m_CDCHits.isRequired(m_cdcHitArrayName);
66 //Store names to speed up creation later
67 m_relRecoTrackTrackName = relRecoTrackTrack.getName();
68
69 if (!m_effStudy) { // by default collects calibration data
70 auto m_tree = new TTree(m_treeName.c_str(), "tree for cdc calibration");
71 m_tree->Branch<Float_t>("x_mea", &x_mea);
72 m_tree->Branch<Float_t>("x_u", &x_u);
73 m_tree->Branch<Float_t>("x_b", &x_b);
74 m_tree->Branch<Float_t>("alpha", &alpha);
75 m_tree->Branch<Float_t>("theta", &theta);
76 m_tree->Branch<Float_t>("t", &t);
77 m_tree->Branch<UShort_t>("adc", &adc);
78 // m_tree->Branch<int>("boardID", &boardID);
79 m_tree->Branch<UChar_t>("lay", &lay);
80 m_tree->Branch<Float_t>("weight", &weight);
81 m_tree->Branch<UShort_t>("IWire", &IWire);
82 m_tree->Branch<Float_t>("Pval", &Pval);
83 m_tree->Branch<Float_t>("ndf", &ndf);
85 m_tree->Branch<Float_t>("d0", &d0);
86 m_tree->Branch<Float_t>("z0", &z0);
87 m_tree->Branch<Float_t>("phi0", &phi0);
88 m_tree->Branch<Float_t>("tanL", &tanL);
89 m_tree->Branch<Float_t>("omega", &omega);
90 }
91
92 if (m_calExpectedDriftTime) { // expected drift time, calculated form xfit
93 m_tree->Branch<Float_t>("t_fit", &t_fit);
94 }
95
96 registerObject<TTree>("tree", m_tree);
97 }
98 if (m_effStudy) { //if m_effStudy is changed to true prepares to only run wire efficiency study
99 auto m_efftree = new TTree(m_effTreeName.c_str(), "tree for wire efficiency");
100 m_efftree->Branch<unsigned short>("layerID", &layerID);
101 m_efftree->Branch<unsigned short>("wireID", &wireID);
102 m_efftree->Branch<float>("z", &z);
103 m_efftree->Branch<bool>("isFound", &isFound);
104
105 registerObject<TTree>("efftree", m_efftree);
106 }
107
108 auto m_hNDF = new TH1F("hNDF", "NDF of fitted track;NDF;Tracks", 71, -1, 70);
109 auto m_hPval = new TH1F("hPval", "p-values of tracks;pVal;Tracks", 1000, 0, 1);
110 auto m_hEventT0 = new TH1F("hEventT0", "Event T0", 1000, -100, 100);
111 auto m_hNTracks = new TH1F("hNTracks", "Number of tracks", 50, 0, 10);
112 auto m_hOccupancy = new TH1F("hOccupancy", "occupancy", 100, 0, 1.0);
113
114 registerObject<TH1F>("hNDF", m_hNDF);
115 registerObject<TH1F>("hPval", m_hPval);
116 registerObject<TH1F>("hEventT0", m_hEventT0);
117 registerObject<TH1F>("hNTracks", m_hNTracks);
118 registerObject<TH1F>("hOccupancy", m_hOccupancy);
119}
120
122{
124
125 /* CDCHit distribution */
126 // make evt t0 incase we dont use evt t0
127 evtT0 = 0;
128
129 //reject events don't have eventT0
131 // event with is fail to extract t0 will be exclude from analysis
132 if (m_eventTimeStoreObject.isValid() && m_eventTimeStoreObject->hasEventT0()) {
133 evtT0 = m_eventTimeStoreObject->getEventT0();
134 getObjectPtr<TH1F>("hEventT0")->Fill(evtT0);
135 } else {
136 return;
137 }
138 }
139 // Collects the WireID and Layer of every hit in this event
140 // Used in wire efficiency building
141 std::vector<unsigned short> wiresInCDCTrack;
142
143 for (CDCTrack& cdcTrack : *m_CDCTracks) {
144 for (CDCRecoHit3D& cdcHit : cdcTrack) {
145 unsigned short eWireID = cdcHit.getWire().getEWire();
146 wiresInCDCTrack.push_back(eWireID);
147 }
148 }
149 // WireID collection finished
150
151 const int nTr = m_Tracks.getEntries();
152 // Skip events which have number of charged tracks <= 1.
153 int nCTracks = 0;
154 for (int i = 0; i < nTr; ++i) {
155 const Belle2::Track* b2track = m_Tracks[i];
157 if (!fitresult) continue;
158
159 short charge = fitresult->getChargeSign();
160 if (fabs(charge) > 0) {
161 nCTracks++;
162 }
163 }
164
165 if (nCTracks <= 1) {
166 return ;
167 } else {
168 getObjectPtr<TH1F>("hNTracks")->Fill(nCTracks);
169 }
170
171 const int nHits = m_CDCHits.getEntries();
172 const int nWires = 14336;
173 float oc = static_cast<float>(nHits) / static_cast<float>(nWires);
174 getObjectPtr<TH1F>("hOccupancy")->Fill(oc);
175
176 for (int i = 0; i < nTr; ++i) {
177 const Belle2::Track* b2track = m_Tracks[i];
179 if (!fitresult) {
180 B2WARNING("No track fit result found.");
181 continue;
182 }
183
185 if (!recoTrack) {
186 B2WARNING("Can not access RecoTrack of this Belle2::Track");
187 continue;
188 }
189 const genfit::FitStatus* fs = recoTrack->getTrackFitStatus();
190 if (!fs) continue;
191 ndf = fs->getNdf();
192 if (!m_bField) {
193 ndf += 1;
194 }
195
196 getObjectPtr<TH1F>("hPval")->Fill(Pval);
197 getObjectPtr<TH1F>("hNDF")->Fill(ndf);
198 B2DEBUG(99, "ndf = " << ndf);
199 B2DEBUG(99, "Pval = " << Pval);
200
201 if (ndf < m_minimumNDF) continue;
202 double Chi2 = fs->getChi2();
203 Pval = std::max(0., ROOT::Math::chisquared_cdf_c(Chi2, ndf));
204 //store track parameters
205
206 d0 = fitresult->getD0();
207 z0 = fitresult->getZ0();
208 tanL = fitresult->getTanLambda();
209 omega = fitresult->getOmega();
210 phi0 = fitresult->getPhi0() * 180 / M_PI;
211
212 // Rejection of suspicious cosmic tracks.
213 // phi0 of cosmic track must be negative in our definition!
214 if (m_isCosmic == true && phi0 > 0.0) continue;
215
216 //cut at Pt
217 if (fitresult->getTransverseMomentum() < m_minimumPt) continue;
218 if (!m_effStudy) { // all harvest to fill the tree if collecting calibration info
219 try {
220 harvest(recoTrack);
221 } catch (const genfit::Exception& e) {
222 B2ERROR("Exception when harvest information from recotrack: " << e.what());
223 }
224 }
225 if (m_effStudy) { // call buildEfficiencies for efficiency study
226 // Request tracks coming from IP
227 if (fitresult->getD0() > 2 || fitresult->getZ0() > 5) continue;
228 const Helix helixFit = fitresult->getHelix();
229 buildEfficiencies(wiresInCDCTrack, helixFit);
230 }
231 }
232
233}
234
236{
237}
238
240{
241 B2DEBUG(99, "start collect hit");
242 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
244
245 for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
246 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
247 if (!tp) continue;
248 lay = hit->getICLayer();
249 IWire = hit->getIWire();
250 adc = hit->getADCCount();
251 unsigned short tdc = hit->getTDCCount();
252 WireID wireid(lay, IWire);
253 const genfit::KalmanFitterInfo* kfi = tp->getKalmanFitterInfo();
254 if (!kfi) {B2DEBUG(199, "No Fitter Info: Layer " << hit->getICLayer()); continue;}
255 for (unsigned int iMeas = 0; iMeas < kfi->getNumMeasurements(); ++iMeas) {
256 if ((kfi->getWeights().at(iMeas)) > 0.5) {
257 // int boardID = cdcgeo.getBoardID(WireID(lay,IWire));
258 const genfit::MeasuredStateOnPlane& mop = kfi->getFittedState();
259 const TVector3 pocaOnWire = mop.getPlane()->getO();//Local wire position
260 const TVector3 pocaOnTrack = mop.getPlane()->getU();//residual direction
261 const TVector3 pocaMom = mop.getMom();
262 alpha = cdcgeo.getAlpha(pocaOnWire, pocaMom) ;
263 theta = cdcgeo.getTheta(pocaMom);
264 x_mea = kfi->getMeasurementOnPlane(iMeas)->getState()(0);
265 x_b = kfi->getFittedState(true).getState()(3);// x fit biased
266 x_u = kfi->getFittedState(false).getState()(3);//x fit unbiased
267 weight = kfi->getWeights().at(iMeas);
268
269 int lr;
270 if (x_u > 0) lr = 1;
271 else lr = 0;
272
273 //Convert to outgoing
274 if (fabs(alpha) > M_PI / 2) {
275 x_b *= -1;
276 x_u *= -1;
277 }
278 x_mea = copysign(x_mea, x_b);
279 lr = cdcgeo.getOutgoingLR(lr, alpha);
281 alpha = cdcgeo.getOutgoingAlpha(alpha);
282
283 B2DEBUG(99, "x_unbiased " << x_u << " |left_right " << lr);
284 if (m_calExpectedDriftTime) { t_fit = cdcgeo.getDriftTime(abs(x_u), lay, lr, alpha, theta);}
285 alpha *= 180 / M_PI;
286 theta *= 180 / M_PI;
287 //estimate drift time
288 t = tdcTrans->getDriftTime(tdc, wireid, mop.getTime(), pocaOnWire.Z(), adc);
289 getObjectPtr<TTree>("tree")->Fill();
290 } //NDF
291 // }//end of if isU
292 }//end of for
293 }//end of for tp
294}//end of func
295
296const CDCWire& CDCCalibrationCollectorModule::getIntersectingWire(const ROOT::Math::XYZVector& xyz, const CDCWireLayer& layer,
297 const Helix& helixFit) const
298{
299 Vector3D crosspoint;
300 if (layer.isAxial())
301 crosspoint = Vector3D(xyz);
302 else {
303 const CDCWire& oneWire = layer.getWire(1);
304 double newR = oneWire.getWirePos2DAtZ(xyz.Z()).norm();
305 double arcLength = helixFit.getArcLength2DAtCylindricalR(newR);
306 ROOT::Math::XYZVector xyzOnWire = B2Vector3D(helixFit.getPositionAtArcLength2D(arcLength));
307 crosspoint = Vector3D(xyzOnWire);
308 }
309
310 const CDCWire& wire = layer.getClosestWire(crosspoint);
311
312 return wire;
313}
314
315void CDCCalibrationCollectorModule::buildEfficiencies(std::vector<unsigned short> wireHits, const Helix helixFit)
316{
318 for (const CDCWireLayer& wireLayer : wireTopology.getWireLayers()) {
319 const double radiusofLayer = wireLayer.getRefCylindricalR();
320 //simple extrapolation of fit
321 const double arcLength = helixFit.getArcLength2DAtCylindricalR(radiusofLayer);
322 const ROOT::Math::XYZVector xyz = B2Vector3D(helixFit.getPositionAtArcLength2D(arcLength));
323 if (!xyz.X()) continue;
324 const CDCWire& wireIntersected = getIntersectingWire(xyz, wireLayer, helixFit);
325 unsigned short crossedWire = wireIntersected.getEWire();
326 unsigned short crossedCWire = wireIntersected.getNeighborCW()->getEWire();
327 unsigned short crossedCCWire = wireIntersected.getNeighborCCW()->getEWire();
328
329 if (find(wireHits.begin(), wireHits.end(), crossedWire) != wireHits.end()
330 || find(wireHits.begin(), wireHits.end(), crossedCWire) != wireHits.end()
331 || find(wireHits.begin(), wireHits.end(), crossedCCWire) != wireHits.end())
332 isFound = true;
333 else
334 isFound = false;
335
336 wireID = wireIntersected.getIWire();
337 layerID = wireIntersected.getICLayer();
338 z = xyz.Z();
339 getObjectPtr<TTree>("efftree")->Fill();
340 }
341}
342
343
344
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
StoreObjPtr< EventT0 > m_eventTimeStoreObject
Event t0 object.
std::string m_recoTrackArrayName
Belle2::RecoTrack StoreArray nam.e.
bool m_calExpectedDriftTime
Calculate expected drift time from x_fit or not.
double m_minimumNDF
minimum NDF required for track
std::string m_effTreeName
Name of efficiency tree for the output file.
std::string m_cdcTrackVectorName
Belle2::CDCTrack vectorpointer name.
bool m_storeTrackParams
Store Track parameter or not.
TrackFindingCDC::StoreWrappedObjPtr< std::vector< TrackFindingCDC::CDCTrack > > m_CDCTracks
CDC tracks.
StoreArray< TrackFitResult > m_TrackFitResults
Track fit results.
void harvest(Belle2::RecoTrack *track)
collect hit information of fitted track.
std::string m_cdcHitArrayName
Belle2::CDCHit StoreArray name.
Float_t theta
Entrance Polar angle of hit (degree).
std::string m_relRecoTrackTrackName
Relation between RecoTrack and Belle2:Track.
unsigned short wireID
wireID for hit-level wire monitoring
void collect() override
Event action, collect information for calibration.
bool m_isCosmic
true when we process cosmic events, else false (collision).
Float_t t_fit
Drift time calculated from x_fit.
bool m_eventT0Extraction
use Event T0 extract t0 or not.
void buildEfficiencies(std::vector< unsigned short > wireHits, const Helix helixFit)
fills efficiency objects
bool m_bField
fit incase no magnetic Field of not, if false, NDF=4 in cal P-value
std::string m_trackArrayName
Belle2::Track StoreArray name.
void prepare() override
Initializes the Module.
Float_t x_u
X_fit for unbiased track fit.
std::string m_trackFitResultArrayName
Belle2::TrackFitResult StoreArray name.
unsigned short layerID
layerID for hit-level wire monitoring
Float_t x_mea
measure drift length (signed by left right).
const TrackFindingCDC::CDCWire & getIntersectingWire(const ROOT::Math::XYZVector &xyz, const TrackFindingCDC::CDCWireLayer &layer, const Helix &helixFit) const
extrapolates the helix fit to a given layer and finds the wire which it would be hitting
std::string m_treeName
Name of tree for the output file.
Float_t alpha
Entrance Azimuthal angle of hit (degree).
double m_minimumPt
minimum pt required for track
bool m_effStudy
When true module collects info only necessary for wire eff study.
float z
z of hit fot hit-level wire monitoring
bool isFound
flag for a hit that has been found near a track as expected by extrapolation
The Class for CDC Geometry Parameters.
double getTheta(const B2Vector3D &momentum) const
Returns track incident angle (theta in rad.).
double getAlpha(const B2Vector3D &posOnWire, const B2Vector3D &momentum) const
Returns track incident angle in rphi plane (alpha in rad.).
double getOutgoingAlpha(const double alpha) const
Converts incoming- to outgoing-alpha.
unsigned short getOutgoingLR(const unsigned short lr, const double alpha) const
Converts incoming-lr to outgoing-lr.
double getOutgoingTheta(const double alpha, const double theta) const
Converts incoming- to outgoing-theta.
double getDriftTime(double dist, unsigned short layer, unsigned short lr, double alpha, double theta) const
Return the drift time to the sense wire.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Helix parameter class.
Definition: Helix.h:48
Translator mirroring the realistic Digitization.
double getDriftTime(unsigned short tdcCount, const WireID &wireID, double timeOfFlightEstimator, double z, unsigned short adcCount) override
Get Drift time.
Calibration collector module base class.
static const ChargedStable muon
muon particle
Definition: Const.h:660
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
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
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
const std::string & getName() const
Return name under which the object is saved in the DataStore.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Class representing a three dimensional reconstructed hit.
Definition: CDCRecoHit3D.h:52
Class representing a sequence of three dimensional reconstructed hits.
Definition: CDCTrack.h:41
Class representing a sense wire layer in the central drift chamber.
Definition: CDCWireLayer.h:42
Class representing the sense wire arrangement in the whole of the central drift chamber.
const std::vector< Belle2::TrackFindingCDC::CDCWireLayer > & getWireLayers() const
Getter for the underlying storing layer vector.
static CDCWireTopology & getInstance()
Getter for the singleton instance of the wire topology.
Class representing a sense wire in the central drift chamber.
Definition: CDCWire.h:58
Vector2D getWirePos2DAtZ(const double z) const
Gives the xy projected position of the wire at the given z coordinate.
Definition: CDCWire.h:192
MayBePtr< const CDCWire > getNeighborCCW() const
Gives the closest neighbor in the counterclockwise direction - always exists.
Definition: CDCWire.cc:159
IWire getIWire() const
Getter for the wire id within its layer.
Definition: CDCWire.h:146
unsigned short getEWire() const
Getter for the encoded wire number.
Definition: CDCWire.h:131
ILayer getICLayer() const
Getter for the continuous layer id ranging from 0 - 55.
Definition: CDCWire.h:150
MayBePtr< const CDCWire > getNeighborCW() const
Gives the closest neighbor in the clockwise direction - always exists.
Definition: CDCWire.cc:164
double norm() const
Calculates the length of the vector.
Definition: Vector2D.h:175
Values of the result of a track fit with a given particle hypothesis.
Helix getHelix() const
Conversion to framework Helix (without covariance).
short getChargeSign() const
Return track charge (1 or -1).
double getOmega() const
Getter for omega.
double getD0() const
Getter for d0.
double getTransverseMomentum() const
Getter for the absolute value of the transverse momentum at the perigee.
double getTanLambda() const
Getter for tanLambda.
double getZ0() const
Getter for z0.
double getPhi0() const
Getter for phi0.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
Definition: Track.cc:104
Class to identify a wire inside the CDC.
Definition: WireID.h:34
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
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
HepGeom::Vector3D< double > Vector3D
3D Vector
Definition: Cell.h:34
Abstract base class for different kinds of events.