Belle II Software release-09-00-03
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 ROOT::Math::XYZVector pos(0, 0, 0);
121 ROOT::Math::XYZVector bfield = BFieldManager::getFieldInTesla(pos);
122 if (bfield.Z() > 0.5) {
123 m_bField = true;
124 B2INFO("CDCCalibrationCollector: Magnetic field is ON");
125 } else {
126 m_bField = false;
127 B2INFO("CDCCalibrationCollector: Magnetic field is OFF");
128 }
129 B2INFO("BField at (0,0,0) = " << bfield.R());
130
131}
132
134{
136
137 /* CDCHit distribution */
138 // make evt t0 incase we dont use evt t0
139 evtT0 = 0;
140
141 //reject events don't have eventT0
143 // event with is fail to extract t0 will be exclude from analysis
144 if (m_eventTimeStoreObject.isValid() && m_eventTimeStoreObject->hasEventT0()) {
145 evtT0 = m_eventTimeStoreObject->getEventT0();
146 getObjectPtr<TH1F>("hEventT0")->Fill(evtT0);
147 } else {
148 return;
149 }
150 }
151 // Collects the WireID and Layer of every hit in this event
152 // Used in wire efficiency building
153 std::vector<unsigned short> wiresInCDCTrack;
154
155 for (CDCTrack& cdcTrack : *m_CDCTracks) {
156 for (CDCRecoHit3D& cdcHit : cdcTrack) {
157 unsigned short eWireID = cdcHit.getWire().getEWire();
158 wiresInCDCTrack.push_back(eWireID);
159 }
160 }
161 // WireID collection finished
162
163 const int nTr = m_Tracks.getEntries();
164 const int nHits = m_CDCHits.getEntries();
165 const int nWires = 14336;
166 float oc = static_cast<float>(nHits) / static_cast<float>(nWires);
167 getObjectPtr<TH1F>("hNTracks")->Fill(nTr);
168 getObjectPtr<TH1F>("hOccupancy")->Fill(oc);
169
170 for (int i = 0; i < nTr; ++i) {
171 const Belle2::Track* b2track = m_Tracks[i];
173 if (!fitresult) {
174 B2WARNING("No track fit result found.");
175 continue;
176 }
177
179 if (!recoTrack) {
180 B2WARNING("Can not access RecoTrack of this Belle2::Track");
181 continue;
182 }
183 const genfit::FitStatus* fs = recoTrack->getTrackFitStatus();
184 if (!fs) continue;
185 ndf = fs->getNdf();
186 if (!m_bField) {
187 ndf += 1;
188 }
189
190 getObjectPtr<TH1F>("hPval")->Fill(Pval);
191 getObjectPtr<TH1F>("hNDF")->Fill(ndf);
192 B2DEBUG(99, "ndf = " << ndf);
193 B2DEBUG(99, "Pval = " << Pval);
194
195 if (ndf < m_minimumNDF) continue;
196 double Chi2 = fs->getChi2();
197 Pval = std::max(0., ROOT::Math::chisquared_cdf_c(Chi2, ndf));
198 //store track parameters
199
200 d0 = fitresult->getD0();
201 z0 = fitresult->getZ0();
202 tanL = fitresult->getTanLambda();
203 omega = fitresult->getOmega();
204 phi0 = fitresult->getPhi0() * 180 / M_PI;
205
206 // Rejection of suspicious cosmic tracks.
207 // phi0 of cosmic track must be negative in our definition!
208 if (m_isCosmic == true && phi0 > 0.0) continue;
209
210 //cut at Pt
211 if (fitresult->getTransverseMomentum() < m_minimumPt) continue;
212 if (!m_effStudy) { // all harvest to fill the tree if collecting calibration info
213 try {
214 harvest(recoTrack);
215 } catch (const genfit::Exception& e) {
216 B2ERROR("Exception when harvest information from recotrack: " << e.what());
217 }
218 }
219 if (m_effStudy) { // call buildEfficiencies for efficiency study
220 // Request tracks coming from IP
221 if (fitresult->getD0() > 2 || fitresult->getZ0() > 5) continue;
222 const Helix helixFit = fitresult->getHelix();
223 buildEfficiencies(wiresInCDCTrack, helixFit);
224 }
225 }
226
227}
228
230{
231}
232
234{
235 B2DEBUG(99, "start collect hit");
236 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
238
239 for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
240 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
241 if (!tp) continue;
242 lay = hit->getICLayer();
243 IWire = hit->getIWire();
244 adc = hit->getADCCount();
245 unsigned short tdc = hit->getTDCCount();
246 WireID wireid(lay, IWire);
247 const genfit::KalmanFitterInfo* kfi = tp->getKalmanFitterInfo();
248 if (!kfi) {B2DEBUG(199, "No Fitter Info: Layer " << hit->getICLayer()); continue;}
249 for (unsigned int iMeas = 0; iMeas < kfi->getNumMeasurements(); ++iMeas) {
250 if ((kfi->getWeights().at(iMeas)) > 0.5) {
251 // int boardID = cdcgeo.getBoardID(WireID(lay,IWire));
252 const genfit::MeasuredStateOnPlane& mop = kfi->getFittedState();
253 const TVector3 pocaOnWire = mop.getPlane()->getO();//Local wire position
254 const TVector3 pocaOnTrack = mop.getPlane()->getU();//residual direction
255 const TVector3 pocaMom = mop.getMom();
256 alpha = cdcgeo.getAlpha(pocaOnWire, pocaMom) ;
257 theta = cdcgeo.getTheta(pocaMom);
258 x_mea = kfi->getMeasurementOnPlane(iMeas)->getState()(0);
259 x_b = kfi->getFittedState(true).getState()(3);// x fit biased
260 x_u = kfi->getFittedState(false).getState()(3);//x fit unbiased
261 weight = kfi->getWeights().at(iMeas);
262
263 int lr;
264 if (x_u > 0) lr = 1;
265 else lr = 0;
266
267 //Convert to outgoing
268 if (fabs(alpha) > M_PI / 2) {
269 x_b *= -1;
270 x_u *= -1;
271 }
272 x_mea = copysign(x_mea, x_b);
273 lr = cdcgeo.getOutgoingLR(lr, alpha);
275 alpha = cdcgeo.getOutgoingAlpha(alpha);
276
277 B2DEBUG(99, "x_unbiased " << x_u << " |left_right " << lr);
278 if (m_calExpectedDriftTime) { t_fit = cdcgeo.getDriftTime(abs(x_u), lay, lr, alpha, theta);}
279 alpha *= 180 / M_PI;
280 theta *= 180 / M_PI;
281 //estimate drift time
282 t = tdcTrans->getDriftTime(tdc, wireid, mop.getTime(), pocaOnWire.Z(), adc);
283 getObjectPtr<TTree>("tree")->Fill();
284 } //NDF
285 // }//end of if isU
286 }//end of for
287 }//end of for tp
288}//end of func
289
290const CDCWire& CDCCalibrationCollectorModule::getIntersectingWire(const ROOT::Math::XYZVector& xyz, const CDCWireLayer& layer,
291 const Helix& helixFit) const
292{
293 Vector3D crosspoint;
294 if (layer.isAxial())
295 crosspoint = Vector3D(xyz);
296 else {
297 const CDCWire& oneWire = layer.getWire(1);
298 double newR = oneWire.getWirePos2DAtZ(xyz.Z()).norm();
299 double arcLength = helixFit.getArcLength2DAtCylindricalR(newR);
300 ROOT::Math::XYZVector xyzOnWire = B2Vector3D(helixFit.getPositionAtArcLength2D(arcLength));
301 crosspoint = Vector3D(xyzOnWire);
302 }
303
304 const CDCWire& wire = layer.getClosestWire(crosspoint);
305
306 return wire;
307}
308
309void CDCCalibrationCollectorModule::buildEfficiencies(std::vector<unsigned short> wireHits, const Helix helixFit)
310{
312 for (const CDCWireLayer& wireLayer : wireTopology.getWireLayers()) {
313 const double radiusofLayer = wireLayer.getRefCylindricalR();
314 //simple extrapolation of fit
315 const double arcLength = helixFit.getArcLength2DAtCylindricalR(radiusofLayer);
316 const ROOT::Math::XYZVector xyz = B2Vector3D(helixFit.getPositionAtArcLength2D(arcLength));
317 if (!xyz.X()) continue;
318 const CDCWire& wireIntersected = getIntersectingWire(xyz, wireLayer, helixFit);
319 unsigned short crossedWire = wireIntersected.getEWire();
320 unsigned short crossedCWire = wireIntersected.getNeighborCW()->getEWire();
321 unsigned short crossedCCWire = wireIntersected.getNeighborCCW()->getEWire();
322
323 if (find(wireHits.begin(), wireHits.end(), crossedWire) != wireHits.end()
324 || find(wireHits.begin(), wireHits.end(), crossedCWire) != wireHits.end()
325 || find(wireHits.begin(), wireHits.end(), crossedCCWire) != wireHits.end())
326 isFound = true;
327 else
328 isFound = false;
329
330 wireID = wireIntersected.getIWire();
331 layerID = wireIntersected.getICLayer();
332 z = xyz.Z();
333 getObjectPtr<TTree>("efftree")->Fill();
334 }
335}
336
337
338
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:61
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 representating a sense wire layer in the central drift chamber.
Definition: CDCWireLayer.h:42
Class representating 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 continious 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).
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.