Project
Loading...
Searching...
No Matches
StrangenessTracker.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
15
16#ifndef _ALICEO2_STRANGENESS_TRACKER_
17#define _ALICEO2_STRANGENESS_TRACKER_
18
19#include <gsl/gsl>
20
32
38
41
42#ifdef ENABLE_UPGRADES
44#endif
45
46namespace o2
47{
48namespace strangeness_tracking
49{
50
51enum DauType : int {
54 kBach = 2
55};
56
58
59 std::array<unsigned int, 7> arr;
60};
61
63{
64 public:
80 using MCLabSpan = gsl::span<const o2::MCCompLabel>;
82
83 StrangenessTracker() = default;
85
86 bool loadData(const o2::globaltracking::RecoContainer& recoData);
87 bool matchDecayToITStrack(float decayR, StrangeTrack& strangeTrack, ClusAttachments& structClus, const TrackITS& itsTrack, std::vector<o2::track::TrackParCovF>& daughterTracks, int iThread = 0);
88 void prepareITStracks();
89 void process();
90 void processV0(int iv0, const V0& v0, const V0Index& v0Idx, int iThread = 0);
91 void processCascade(int icasc, const Cascade& casc, const CascadeIndex& cascIdx, const V0& cascV0, int iThread = 0);
92 void process3Body(int i3body, const Decay3Body& dec3body, const Decay3BodyIndex& dec3bodyIdx, int iThread = 0);
93 bool updateTrack(const ITSCluster& clus, o2::track::TrackParCov& track);
94
95 std::vector<ClusAttachments>& getClusAttachments(int iThread = 0) { return mClusAttachments[iThread]; };
96 std::vector<StrangeTrack>& getStrangeTrackVec(int iThread = 0) { return mStrangeTrackVec[iThread]; };
97 std::vector<o2::MCCompLabel>& getStrangeTrackLabels(int iThread = 0) { return mStrangeTrackLabels[iThread]; };
98 size_t getNTracks(int ithread = 0) const { return ithread < (int)mStrangeTrackVec.size() ? mStrangeTrackVec[ithread].size() : 0; }
99
100 float getBz() const { return mBz; }
101 void setBz(float d) { mBz = d; }
105 void setMCTruthOn(bool v) { mMCTruthON = v; }
106 bool getMCTruthOn() const { return mMCTruthON; }
107
108#ifdef ENABLE_UPGRADES
109 void setClusterDictionaryIT3(const o2::its3::TopologyDictionary* d)
110 {
111 mIT3Dict = d;
112 }
113#endif
114
115 void clear()
116 {
117 for (int i = 0; i < mNThreads; i++) {
118 mDaughterTracks[i].clear();
119 mClusAttachments[i].clear();
120 mStrangeTrackVec[i].clear();
121 if (mMCTruthON) {
122 mStrangeTrackLabels[i].clear();
123 }
124 }
125 mTracksIdxTable.clear();
126 mSortedITStracks.clear();
127 mSortedITSindexes.clear();
128 mITSvtxBrackets.clear();
129 mInputITSclusters.clear();
130 mInputClusterSizes.clear();
131 }
132
133 void setupThreads(int nThreads = 1)
134 {
135 mNThreads = nThreads;
136 mFitterV0.resize(nThreads);
137 mFitter3Body.resize(nThreads);
138 mFitter4Body.resize(nThreads);
139 mStrangeTrackVec.resize(nThreads);
140 mClusAttachments.resize(nThreads);
141 mStrangeTrackLabels.resize(nThreads);
142 mDaughterTracks.resize(nThreads);
143 }
144
146 {
147 for (auto& fitter : mFitterV0) {
148 fitter.setBz(mBz);
149 fitter.setUseAbsDCA(true);
150 }
151 for (auto& fitter : mFitter3Body) {
152 fitter.setBz(mBz);
153 fitter.setUseAbsDCA(true);
154 }
155 for (auto& fitter : mFitter4Body) {
156 fitter.setBz(mBz);
157 fitter.setUseAbsDCA(true);
158 }
159 }
160
161 double calcV0alpha(const V0& v0)
162 {
163 std::array<float, 3> momT, momP, momN;
164 v0.getProng(0).getPxPyPzGlo(momP);
165 v0.getProng(1).getPxPyPzGlo(momN);
166 v0.getPxPyPzGlo(momT);
167 float qNeg = momN[0] * momT[0] + momN[1] * momT[1] + momN[2] * momT[2];
168 float qPos = momP[0] * momT[0] + momP[1] * momT[1] + momP[2] * momT[2];
169 return (qPos - qNeg) / (qPos + qNeg);
170 };
171
172 double calcMotherMass(const std::array<float, 3>& pDauFirst, const std::array<float, 3>& pDauSecond, PID pidDauFirst, PID pidDauSecond)
173 {
174 double m2DauFirst = PID::getMass2(pidDauFirst);
175 double m2DauSecond = PID::getMass2(pidDauSecond);
176 double p2DauFirst = (pDauFirst[0] * pDauFirst[0]) + (pDauFirst[1] * pDauFirst[1]) + (pDauFirst[2] * pDauFirst[2]);
177 double p2DauSecond = (pDauSecond[0] * pDauSecond[0]) + (pDauSecond[1] * pDauSecond[1]) + (pDauSecond[2] * pDauSecond[2]);
178 float ePos = std::sqrt(p2DauFirst + m2DauFirst), eNeg = std::sqrt(p2DauSecond + m2DauSecond);
179
180 double e2Mother = (ePos + eNeg) * (ePos + eNeg);
181 double pxMother = (pDauFirst[0] + pDauSecond[0]);
182 double pyMother = (pDauFirst[1] + pDauSecond[1]);
183 double pzMother = (pDauFirst[2] + pDauSecond[2]);
184 double p2Mother = (pxMother * pxMother) + (pyMother * pyMother) + (pzMother * pzMother);
185 return std::sqrt(e2Mother - p2Mother);
186 }
187
188 double calcMotherMass3body(const std::array<float, 3>& pDauFirst, const std::array<float, 3>& pDauSecond, const std::array<float, 3>& pDauThird, PID pidDauFirst, PID pidDauSecond, PID pidDauThird)
189 {
190 double m2DauFirst = PID::getMass2(pidDauFirst);
191 double m2DauSecond = PID::getMass2(pidDauSecond);
192 double m2DauThird = PID::getMass2(pidDauThird);
193 double p2DauFirst = (pDauFirst[0] * pDauFirst[0]) + (pDauFirst[1] * pDauFirst[1]) + (pDauFirst[2] * pDauFirst[2]);
194 double p2DauSecond = (pDauSecond[0] * pDauSecond[0]) + (pDauSecond[1] * pDauSecond[1]) + (pDauSecond[2] * pDauSecond[2]);
195 double p2DauThird = (pDauThird[0] * pDauThird[0]) + (pDauThird[1] * pDauThird[1]) + (pDauThird[2] * pDauThird[2]);
196 float eFirst = std::sqrt(p2DauFirst + m2DauFirst), eSecond = std::sqrt(p2DauSecond + m2DauSecond), eThird = std::sqrt(p2DauThird + m2DauThird);
197
198 double e2Mother = (eFirst + eSecond + eThird) * (eFirst + eSecond + eThird);
199 double pxMother = (pDauFirst[0] + pDauSecond[0] + pDauThird[0]);
200 double pyMother = (pDauFirst[1] + pDauSecond[1] + pDauThird[1]);
201 double pzMother = (pDauFirst[2] + pDauSecond[2] + pDauThird[2]);
202 double p2Mother = (pxMother * pxMother) + (pyMother * pyMother) + (pzMother * pzMother);
203 return std::sqrt(e2Mother - p2Mother);
204 }
205
206 bool recreateV0(const o2::track::TrackParCov& posTrack, const o2::track::TrackParCov& negTrack, V0& newV0, int iThread = 0)
207 {
208 int nCand;
209 try {
210 nCand = mFitterV0[iThread].process(posTrack, negTrack);
211 } catch (std::runtime_error& e) {
212 return false;
213 }
214 if (!nCand || !mFitterV0[iThread].propagateTracksToVertex()) {
215 return false;
216 }
217
218 const auto& v0XYZ = mFitterV0[iThread].getPCACandidatePos();
219
220 auto& propPos = mFitterV0[iThread].getTrack(0, 0);
221 auto& propNeg = mFitterV0[iThread].getTrack(1, 0);
222
223 std::array<float, 3> pP, pN;
224 propPos.getPxPyPzGlo(pP);
225 propNeg.getPxPyPzGlo(pN);
226 std::array<float, 3> pV0 = {pP[0] + pN[0], pP[1] + pN[1], pP[2] + pN[2]};
227 newV0 = V0(v0XYZ, pV0, mFitterV0[iThread].calcPCACovMatrixFlat(0), propPos, propNeg, PID::HyperTriton);
228 return true;
229 };
230
231 std::vector<ITSCluster> getTrackClusters(const TrackITS& itsTrack)
232 {
233 std::vector<ITSCluster> outVec;
234 outVec.reserve(7);
235 auto firstClus = itsTrack.getFirstClusterEntry();
236 auto ncl = itsTrack.getNumberOfClusters();
237 for (int icl = 0; icl < ncl; icl++) {
238 outVec.push_back(mInputITSclusters[mInputITSidxs[firstClus + icl]]);
239 }
240 return outVec;
241 };
242
243 std::vector<int> getTrackClusterSizes(const TrackITS& itsTrack)
244 {
245 std::vector<int> outVec;
246 outVec.reserve(7);
247 auto firstClus = itsTrack.getFirstClusterEntry();
248 auto ncl = itsTrack.getNumberOfClusters();
249 for (int icl = 0; icl < ncl; icl++) {
250 outVec.push_back(mInputClusterSizes[mInputITSidxs[firstClus + icl]]);
251 }
252 return outVec;
253 };
254
255 void getClusterSizesITS(std::vector<int>& clusSizeVec, const gsl::span<const o2::itsmft::CompClusterExt> ITSclus, gsl::span<const unsigned char>::iterator& pattIt, const o2::itsmft::TopologyDictionary* mdict)
256 {
257 for (unsigned int iClus{0}; iClus < ITSclus.size(); ++iClus) {
258 auto& clus = ITSclus[iClus];
259 auto pattID = clus.getPatternID();
260 int npix;
262
263 if (pattID == o2::itsmft::CompCluster::InvalidPatternID || mdict->isGroup(pattID)) {
264 patt.acquirePattern(pattIt);
265 npix = patt.getNPixels();
266 } else {
267
268 npix = mdict->getNpixels(pattID);
269 patt = mdict->getPattern(pattID);
270 }
271 clusSizeVec[iClus] = npix;
272 }
273 // LOG(info) << " Patt Npixel: " << pattVec[0].getNPixels();
274 }
275
276#ifdef ENABLE_UPGRADES
277 void getClusterSizesIT3(std::vector<int>& clusSizeVec, const gsl::span<const o2::itsmft::CompClusterExt> ITSclus, gsl::span<const unsigned char>::iterator& pattIt, const o2::its3::TopologyDictionary* mdict)
278 {
279 for (unsigned int iClus{0}; iClus < ITSclus.size(); ++iClus) {
280 auto& clus = ITSclus[iClus];
281 auto pattID = clus.getPatternID();
282 int npix;
284
285 if (pattID == o2::itsmft::CompCluster::InvalidPatternID || mdict->isGroup(pattID)) {
286 patt.acquirePattern(pattIt);
287 npix = patt.getNPixels();
288 } else {
289
290 npix = mdict->getNpixels(pattID);
291 patt = mdict->getPattern(pattID);
292 }
293 clusSizeVec[iClus] = npix;
294 }
295 // LOG(info) << " Patt Npixel: " << pattVec[0].getNPixels();
296 }
297#endif
298
300 {
301 if (v0.rotate(itsTrack.getParamOut().getAlpha()) && v0.propagateTo(itsTrack.getParamOut().getX(), mBz)) {
302 return v0.getPredictedChi2(itsTrack.getParamOut());
303 }
304 return -100;
305 };
306
307 o2::MCCompLabel getStrangeTrackLabel(const TrackITS& itsTrack, const StrangeTrack& strangeTrack, const ClusAttachments& structClus) // ITS label with fake flag recomputed
308 {
309 bool isFake = false;
310 auto itsTrkLab = mITSTrkLabels[strangeTrack.mITSRef];
311 for (unsigned int iLay = 0; iLay < 7; iLay++) {
312 if (itsTrack.hasHitOnLayer(iLay) && itsTrack.isFakeOnLayer(iLay) && structClus.arr[iLay] == 0) {
313 isFake = true;
314 break;
315 }
316 }
317 itsTrkLab.setFakeFlag(isFake);
318 return itsTrkLab;
319 }
320
321 protected:
322 bool mMCTruthON = false;
323 int mNThreads = 1;
324 gsl::span<const TrackITS> mInputITStracks; // input ITS tracks
325 std::vector<VBracket> mITSvtxBrackets; // time brackets for ITS tracks
326 std::vector<int> mTracksIdxTable; // index table for ITS tracks
327 std::vector<int> mInputClusterSizes; // input cluster sizes
328 std::vector<ITSCluster> mInputITSclusters; // input ITS clusters
329 gsl::span<const int> mInputITSidxs; // input ITS track-cluster indexes
330 gsl::span<const V0> mInputV0tracks; // input V0 of decay daughters
331 gsl::span<const V0Index> mInputV0Indices; // input V0 indices of decay daughters
332 gsl::span<const Cascade> mInputCascadeTracks; // input cascade of decay daughters
333 gsl::span<const CascadeIndex> mInputCascadeIndices; // input cascade indices of decay daughters
334 gsl::span<const Decay3Body> mInput3BodyTracks; // input decay3body of decay daughters
335 gsl::span<const Decay3BodyIndex> mInput3BodyIndices; // input decay3body indices of decay daughters
336 const MCLabContCl* mITSClsLabels = nullptr;
338
339 std::vector<o2::its::TrackITS> mSortedITStracks; // sorted ITS tracks
340 std::vector<int> mSortedITSindexes; // indexes of sorted ITS tracks
341 IndexTableUtils mUtils; // structure for computing eta/phi matching selections
342
343 std::vector<std::vector<StrangeTrack>> mStrangeTrackVec; // structure containing updated mother and daughter tracks (per thread)
344 std::vector<std::vector<ClusAttachments>> mClusAttachments; // # of attached tracks, -1 not attached, 0 for the mother, > 0 for the daughters (per thread)
345 std::vector<std::vector<o2::MCCompLabel>> mStrangeTrackLabels; // vector of MC labels for mother track (per thread)
346
348 float mBz = -5; // Magnetic field
350#ifdef ENABLE_UPGRADES
351 const o2::its3::TopologyDictionary* mIT3Dict = nullptr;
352#endif
353
354 std::vector<DCAFitter2> mFitterV0; // optional DCA Fitter for recreating V0 with hypertriton mass hypothesis (per thread)
355 std::vector<DCAFitter3> mFitter3Body; // optional DCA Fitter for final 3 Body refit (per thread)
356 std::vector<DCAFitter4> mFitter4Body; // optional DCA Fitter for final 4 Body refit (per thread)
357
359
360 std::vector<std::vector<o2::track::TrackParCovF>> mDaughterTracks; // vector of daughter tracks (per thread)
361 ClusAttachments mStructClus; // # of attached tracks, 1 for mother, 2 for daughter
362
364};
365
366} // namespace strangeness_tracking
367} // namespace o2
368
369#endif // _ALICEO2_STRANGENESS_TRACKER_
Definition of the ITSMFT compact cluster.
Defintions for N-prongs secondary vertex fit.
Wrapper container for different reconstructed object types.
Definition of the ClusterTopology class.
particle ids, masses, names class definition
Base track model for the Barrel, params only, w/o covariance.
Definition of the BuildTopologyDictionary class for ITS3.
int32_t i
Definition of the GeometryTGeo class.
Definition of the ITS track.
Extention of GlobalTrackID by flags relevant for verter-track association.
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
void setFakeFlag(bool v=true)
TO BE DONE: extend to generic N body vertex.
Definition Decay3Body.h:26
bool isGroup(int n) const
Returns true if the element corresponds to a group of rare topologies.
int getNpixels(int n) const
Returns the number of fired pixels of the n_th element.
const itsmft::ClusterPattern & getPattern(int n) const
Returns the pattern of the topology.
bool hasHitOnLayer(uint32_t i) const
Definition TrackITS.h:104
bool isFakeOnLayer(uint32_t i) const
Definition TrackITS.h:105
int getFirstClusterEntry() const
Definition TrackITS.h:66
void acquirePattern(iterator &pattIt)
int getNPixels() const
Returns the number of fired pixels.
static constexpr unsigned short InvalidPatternID
Definition CompCluster.h:46
const ClusterPattern & getPattern(int n) const
Returns the pattern of the topology.
int getNpixels(int n) const
Returns the number of fired pixels of the n_th element.
bool isGroup(int n) const
Returns true if the element corresponds to a group of rare topologies.
void setConfigParams(const StrangenessTrackingParamConfig *params)
std::vector< std::vector< o2::track::TrackParCovF > > mDaughterTracks
o2::MCCompLabel getStrangeTrackLabel(const TrackITS &itsTrack, const StrangeTrack &strangeTrack, const ClusAttachments &structClus)
const StrangenessTrackingParamConfig * mStrParams
double calcMotherMass(const std::array< float, 3 > &pDauFirst, const std::array< float, 3 > &pDauSecond, PID pidDauFirst, PID pidDauSecond)
gsl::span< const TrackITS > mInputITStracks
number of threads (externally driven)
int mNThreads
flag availability of MC truth
std::vector< std::vector< StrangeTrack > > mStrangeTrackVec
void processV0(int iv0, const V0 &v0, const V0Index &v0Idx, int iThread=0)
std::vector< StrangeTrack > & getStrangeTrackVec(int iThread=0)
void processCascade(int icasc, const Cascade &casc, const CascadeIndex &cascIdx, const V0 &cascV0, int iThread=0)
std::vector< ITSCluster > getTrackClusters(const TrackITS &itsTrack)
o2::base::PropagatorImpl< float >::MatCorrType mCorrType
std::vector< int > getTrackClusterSizes(const TrackITS &itsTrack)
void getClusterSizesITS(std::vector< int > &clusSizeVec, const gsl::span< const o2::itsmft::CompClusterExt > ITSclus, gsl::span< const unsigned char >::iterator &pattIt, const o2::itsmft::TopologyDictionary *mdict)
std::vector< std::vector< o2::MCCompLabel > > mStrangeTrackLabels
bool matchDecayToITStrack(float decayR, StrangeTrack &strangeTrack, ClusAttachments &structClus, const TrackITS &itsTrack, std::vector< o2::track::TrackParCovF > &daughterTracks, int iThread=0)
void setCorrType(const o2::base::PropagatorImpl< float >::MatCorrType &type)
std::vector< o2::its::TrackITS > mSortedITStracks
input ITS Track MC labels
gsl::span< const Decay3Body > mInput3BodyTracks
bool recreateV0(const o2::track::TrackParCov &posTrack, const o2::track::TrackParCov &negTrack, V0 &newV0, int iThread=0)
std::vector< std::vector< ClusAttachments > > mClusAttachments
std::vector< ClusAttachments > & getClusAttachments(int iThread=0)
double calcMotherMass3body(const std::array< float, 3 > &pDauFirst, const std::array< float, 3 > &pDauSecond, const std::array< float, 3 > &pDauThird, PID pidDauFirst, PID pidDauSecond, PID pidDauThird)
bool updateTrack(const ITSCluster &clus, o2::track::TrackParCov &track)
void process3Body(int i3body, const Decay3Body &dec3body, const Decay3BodyIndex &dec3bodyIdx, int iThread=0)
gsl::span< const Decay3BodyIndex > mInput3BodyIndices
float getMatchingChi2(o2::track::TrackParCovF v0, const TrackITS &itsTrack)
void setClusterDictionaryITS(const o2::itsmft::TopologyDictionary *d)
const o2::itsmft::TopologyDictionary * mITSDict
gsl::span< const CascadeIndex > mInputCascadeIndices
bool loadData(const o2::globaltracking::RecoContainer &recoData)
std::vector< o2::MCCompLabel > & getStrangeTrackLabels(int iThread=0)
MCLabSpan mITSTrkLabels
input ITS Cluster MC labels
gsl::span< const o2::MCCompLabel > MCLabSpan
static constexpr ID HyperTriton
Definition PID.h:113
const GLdouble * v
Definition glcorearb.h:832
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLenum const GLfloat * params
Definition glcorearb.h:272
GLfloat v0
Definition glcorearb.h:811
TrackParCovF TrackParCov
Definition Track.h:33
TrackParametrizationWithError< float > TrackParCovF
Definition Track.h:31
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...