Belle II Software prerelease-10-00-00a
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#include <cmath>
29
30//using namespace std;
31using namespace Belle2;
32using namespace CDC;
33using namespace genfit;
34using namespace TrackFindingCDC;
35
36
37REG_MODULE(CDCCalibrationCollector);
38
39
41{
42 setDescription("Collector module for cdc calibration");
43 setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
44 addParam("recoTracksColName", m_recoTrackArrayName, "Name of collection hold genfit::Track", std::string(""));
45 addParam("bField", m_bField, "If true -> #Params ==5 else #params ==4 for calculate P-Val", false);
46 addParam("calExpectedDriftTime", m_calExpectedDriftTime, "if true module will calculate expected drift time, it take a time",
47 true);
48 addParam("storeTrackParams", m_storeTrackParams, "Store Track Parameter or not, it will be multicount for each hit", false);
49 addParam("eventT0Extraction", m_eventT0Extraction, "use event t0 extract t0 or not", true);
50 addParam("minimumPt", m_minimumPt, "Tracks with transverse momentum smaller than this value will not used", 0.15);
51 addParam("minimumNDF", m_minimumNDF, "Discard tracks whose degree-of-freedom below this value", 5.);
52 addParam("isCosmic", m_isCosmic, "True when we process cosmic events, else False (collision)", m_isCosmic);
53 addParam("effStudy", m_effStudy, "When true module collects info only necessary for wire eff study", false);
54}
55
59
61{
62 m_Tracks.isRequired(m_trackArrayName);
65 m_CDCHits.isRequired(m_cdcHitArrayName);
68 //Store names to speed up creation later
69 m_relRecoTrackTrackName = relRecoTrackTrack.getName();
70
71 if (!m_effStudy) { // by default collects calibration data
72 auto m_tree = new TTree(m_treeName.c_str(), "tree for cdc calibration");
73 m_tree->Branch<Float_t>("x_mea", &x_mea);
74 m_tree->Branch<Float_t>("x_u", &x_u);
75 m_tree->Branch<Float_t>("x_b", &x_b);
76 m_tree->Branch<Float_t>("alpha", &alpha);
77 m_tree->Branch<Float_t>("theta", &theta);
78 m_tree->Branch<Float_t>("t", &t);
79 m_tree->Branch<UShort_t>("adc", &adc);
80 // m_tree->Branch<int>("boardID", &boardID);
81 m_tree->Branch<UChar_t>("lay", &lay);
82 m_tree->Branch<Float_t>("weight", &weight);
83 m_tree->Branch<UShort_t>("IWire", &IWire);
84 m_tree->Branch<Float_t>("Pval", &Pval);
85 m_tree->Branch<Float_t>("ndf", &ndf);
87 m_tree->Branch<Float_t>("d0", &d0);
88 m_tree->Branch<Float_t>("z0", &z0);
89 m_tree->Branch<Float_t>("phi0", &phi0);
90 m_tree->Branch<Float_t>("tanL", &tanL);
91 m_tree->Branch<Float_t>("omega", &omega);
92 }
93
94 if (m_calExpectedDriftTime) { // expected drift time, calculated form xfit
95 m_tree->Branch<Float_t>("t_fit", &t_fit);
96 }
97
98 registerObject<TTree>("tree", m_tree);
99 }
100 if (m_effStudy) { //if m_effStudy is changed to true prepares to only run wire efficiency study
101 auto m_efftree = new TTree(m_effTreeName.c_str(), "tree for wire efficiency");
102 m_efftree->Branch<unsigned short>("layerID", &layerID);
103 m_efftree->Branch<unsigned short>("wireID", &wireID);
104 m_efftree->Branch<float>("z", &z);
105 m_efftree->Branch<bool>("isFound", &isFound);
106
107 registerObject<TTree>("efftree", m_efftree);
108 }
109
110 auto m_hNDF = new TH1F("hNDF", "NDF of fitted track;NDF;Tracks", 71, -1, 70);
111 auto m_hPval = new TH1F("hPval", "p-values of tracks;pVal;Tracks", 1000, 0, 1);
112 auto m_hEventT0 = new TH1F("hEventT0", "Event T0", 1000, -100, 100);
113 auto m_hNTracks = new TH1F("hNTracks", "Number of tracks", 50, 0, 10);
114 auto m_hOccupancy = new TH1F("hOccupancy", "occupancy", 100, 0, 1.0);
115
116 registerObject<TH1F>("hNDF", m_hNDF);
117 registerObject<TH1F>("hPval", m_hPval);
118 registerObject<TH1F>("hEventT0", m_hEventT0);
119 registerObject<TH1F>("hNTracks", m_hNTracks);
120 registerObject<TH1F>("hOccupancy", m_hOccupancy);
121}
122
124{
126
127 /* CDCHit distribution */
128 // make evt t0 incase we dont use evt t0
129 evtT0 = 0;
130
131 //reject events don't have eventT0
133 // event with is fail to extract t0 will be exclude from analysis
134 if (m_eventTimeStoreObject.isValid() && m_eventTimeStoreObject->hasEventT0()) {
135 evtT0 = m_eventTimeStoreObject->getEventT0();
136 getObjectPtr<TH1F>("hEventT0")->Fill(evtT0);
137 } else {
138 return;
139 }
140 }
141 // Collects the WireID and Layer of every hit in this event
142 // Used in wire efficiency building
143 std::vector<unsigned short> wiresInCDCTrack;
144
145 for (CDCTrack& cdcTrack : *m_CDCTracks) {
146 for (CDCRecoHit3D& cdcHit : cdcTrack) {
147 unsigned short eWireID = cdcHit.getWire().getEWire();
148 wiresInCDCTrack.push_back(eWireID);
149 }
150 }
151 // WireID collection finished
152
153 const int nTr = m_Tracks.getEntries();
154 // Skip events which have number of charged tracks <= 1.
155 int nCTracks = 0;
156 for (int i = 0; i < nTr; ++i) {
157 const Belle2::Track* b2track = m_Tracks[i];
159 if (!fitresult) continue;
160
161 short charge = fitresult->getChargeSign();
162 if (std::abs(charge) > 0) {
163 nCTracks++;
164 }
165 }
166
167 if (nCTracks <= 1) {
168 return ;
169 } else {
170 getObjectPtr<TH1F>("hNTracks")->Fill(nCTracks);
171 }
172
173 const int nHits = m_CDCHits.getEntries();
174 const int nWires = 14336;
175 float oc = static_cast<float>(nHits) / static_cast<float>(nWires);
176 getObjectPtr<TH1F>("hOccupancy")->Fill(oc);
177
178 for (int i = 0; i < nTr; ++i) {
179 const Belle2::Track* b2track = m_Tracks[i];
181 if (!fitresult) {
182 B2WARNING("No track fit result found.");
183 continue;
184 }
185
187 if (!recoTrack) {
188 B2WARNING("Can not access RecoTrack of this Belle2::Track");
189 continue;
190 }
191 const genfit::FitStatus* fs = recoTrack->getTrackFitStatus();
192 if (!fs) continue;
193 ndf = fs->getNdf();
194 if (!m_bField) {
195 ndf += 1;
196 }
197
198 getObjectPtr<TH1F>("hPval")->Fill(Pval);
199 getObjectPtr<TH1F>("hNDF")->Fill(ndf);
200 B2DEBUG(99, "ndf = " << ndf);
201 B2DEBUG(99, "Pval = " << Pval);
202
203 if (ndf < m_minimumNDF) continue;
204 double Chi2 = fs->getChi2();
205 Pval = std::max(0., ROOT::Math::chisquared_cdf_c(Chi2, ndf));
206 //store track parameters
207
208 d0 = fitresult->getD0();
209 z0 = fitresult->getZ0();
210 tanL = fitresult->getTanLambda();
211 omega = fitresult->getOmega();
212 phi0 = fitresult->getPhi0() * 180 / M_PI;
213
214 // Rejection of suspicious cosmic tracks.
215 // phi0 of cosmic track must be negative in our definition!
216 if (m_isCosmic == true && phi0 > 0.0) continue;
217
218 //cut at Pt
219 if (fitresult->getTransverseMomentum() < m_minimumPt) continue;
220 if (!m_effStudy) { // all harvest to fill the tree if collecting calibration info
221 try {
222 harvest(recoTrack);
223 } catch (const genfit::Exception& e) {
224 B2ERROR("Exception when harvest information from recotrack: " << e.what());
225 }
226 }
227 if (m_effStudy) { // call buildEfficiencies for efficiency study
228 // Request tracks coming from IP
229 if (fitresult->getD0() > 2 || fitresult->getZ0() > 5) continue;
230 const Helix helixFit = fitresult->getHelix();
231 buildEfficiencies(wiresInCDCTrack, helixFit);
232 }
233 }
234
235}
236
240
242{
243 B2DEBUG(99, "start collect hit");
244 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
246
247 for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
248 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
249 if (!tp) continue;
250 lay = hit->getICLayer();
251 IWire = hit->getIWire();
252 adc = hit->getADCCount();
253 unsigned short tdc = hit->getTDCCount();
254 WireID wireid(lay, IWire);
255 const genfit::KalmanFitterInfo* kfi = tp->getKalmanFitterInfo();
256 if (!kfi) {B2DEBUG(199, "No Fitter Info: Layer " << hit->getICLayer()); continue;}
257 for (unsigned int iMeas = 0; iMeas < kfi->getNumMeasurements(); ++iMeas) {
258 if ((kfi->getWeights().at(iMeas)) > 0.5) {
259 // int boardID = cdcgeo.getBoardID(WireID(lay,IWire));
260 const genfit::MeasuredStateOnPlane& mop = kfi->getFittedState();
261 const TVector3 pocaOnWire = mop.getPlane()->getO();//Local wire position
262 const TVector3 pocaOnTrack = mop.getPlane()->getU();//residual direction
263 const TVector3 pocaMom = mop.getMom();
264 alpha = cdcgeo.getAlpha(pocaOnWire, pocaMom) ;
265 theta = cdcgeo.getTheta(pocaMom);
266 x_mea = kfi->getMeasurementOnPlane(iMeas)->getState()(0);
267 x_b = kfi->getFittedState(true).getState()(3);// x fit biased
268 x_u = kfi->getFittedState(false).getState()(3);//x fit unbiased
269 weight = kfi->getWeights().at(iMeas);
270
271 int lr;
272 if (x_u > 0) lr = 1;
273 else lr = 0;
274
275 //Convert to outgoing
276 if (std::abs(alpha) > M_PI / 2) {
277 x_b *= -1;
278 x_u *= -1;
279 }
280 x_mea = copysign(x_mea, x_b);
281 lr = cdcgeo.getOutgoingLR(lr, alpha);
283 alpha = cdcgeo.getOutgoingAlpha(alpha);
284
285 B2DEBUG(99, "x_unbiased " << x_u << " |left_right " << lr);
286 if (m_calExpectedDriftTime) { t_fit = cdcgeo.getDriftTime(std::abs(x_u), lay, lr, alpha, theta);}
287 alpha *= 180 / M_PI;
288 theta *= 180 / M_PI;
289 //estimate drift time
290 t = tdcTrans->getDriftTime(tdc, wireid, mop.getTime(), pocaOnWire.Z(), adc);
291 getObjectPtr<TTree>("tree")->Fill();
292 } //NDF
293 // }//end of if isU
294 }//end of for
295 }//end of for tp
296}//end of func
297
298const CDCWire& CDCCalibrationCollectorModule::getIntersectingWire(const ROOT::Math::XYZVector& xyz, const CDCWireLayer& layer,
299 const Helix& helixFit) const
300{
301 Vector3D crosspoint;
302 if (layer.isAxial())
303 crosspoint = Vector3D(xyz);
304 else {
305 const CDCWire& oneWire = layer.getWire(1);
306 double newR = oneWire.getWirePos2DAtZ(xyz.Z()).norm();
307 double arcLength = helixFit.getArcLength2DAtCylindricalR(newR);
308 ROOT::Math::XYZVector xyzOnWire = B2Vector3D(helixFit.getPositionAtArcLength2D(arcLength));
309 crosspoint = Vector3D(xyzOnWire);
310 }
311
312 const CDCWire& wire = layer.getClosestWire(crosspoint);
313
314 return wire;
315}
316
317void CDCCalibrationCollectorModule::buildEfficiencies(std::vector<unsigned short> wireHits, const Helix helixFit)
318{
320 for (const CDCWireLayer& wireLayer : wireTopology.getWireLayers()) {
321 const double radiusofLayer = wireLayer.getRefCylindricalR();
322 //simple extrapolation of fit
323 const double arcLength = helixFit.getArcLength2DAtCylindricalR(radiusofLayer);
324 const ROOT::Math::XYZVector xyz = B2Vector3D(helixFit.getPositionAtArcLength2D(arcLength));
325 if (!xyz.X()) continue;
326 const CDCWire& wireIntersected = getIntersectingWire(xyz, wireLayer, helixFit);
327 unsigned short crossedWire = wireIntersected.getEWire();
328 unsigned short crossedCWire = wireIntersected.getNeighborCW()->getEWire();
329 unsigned short crossedCCWire = wireIntersected.getNeighborCCW()->getEWire();
330
331 if (find(wireHits.begin(), wireHits.end(), crossedWire) != wireHits.end()
332 || find(wireHits.begin(), wireHits.end(), crossedCWire) != wireHits.end()
333 || find(wireHits.begin(), wireHits.end(), crossedCCWire) != wireHits.end())
334 isFound = true;
335 else
336 isFound = false;
337
338 wireID = wireIntersected.getIWire();
339 layerID = wireIntersected.getICLayer();
340 z = xyz.Z();
341 getObjectPtr<TTree>("efftree")->Fill();
342 }
343}
344
345
346
StoreObjPtr< EventT0 > m_eventTimeStoreObject
Event t0 object.
std::string m_recoTrackArrayName
Belle2::RecoTrack StoreArray name.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 for 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.
void registerObject(std::string name, T *obj)
Register object with a name, takes ownership, do not access the pointer beyond prepare()
CalibrationCollectorModule()
Constructor. Sets the default prefix for calibration dataobjects.
T * getObjectPtr(std::string name)
Calls the CalibObjManager to get the requested stored collector data.
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
CDCHit UsedCDCHit
Define, use of CDC hits as CDC hits (for symmetry).
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.
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.
Class representing a three dimensional reconstructed hit.
Class representing a sequence of three dimensional reconstructed hits.
Definition CDCTrack.h:41
Class representing a sense wire layer in the central drift chamber.
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 the 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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
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.