Project
Loading...
Searching...
No Matches
TracksToRecords.cxx
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
13
14#include <TTreeReader.h>
15#include <TTreeReaderValue.h>
16
17#include "Framework/InputSpec.h"
18#include "Framework/Logger.h"
20#include "MFTBase/Geometry.h"
22
24
25using namespace o2::mft;
26
28
29//__________________________________________________________________________
31 : Aligner(),
32 mRunNumber(0),
33 mBz(0),
34 mNumberTFs(0),
35 mNumberOfClusterChainROFs(0),
36 mNumberOfTrackChainROFs(0),
37 mCounterLocalEquationFailed(0),
38 mCounterSkippedTracks(0),
39 mCounterUsedTracks(0),
40 mGlobalDerivatives(std::vector<double>(mNumberOfGlobalParam)),
41 mLocalDerivatives(std::vector<double>(mNumberOfTrackParam)),
42 mMinNumberClusterCut(6),
43 mWeightRecord(1.),
44 mDictionary(nullptr),
45 mAlignPoint(new AlignPointHelper()),
46 mWithControl(false),
47 mNEntriesAutoSave(10000),
48 mRecordWriter(new o2::fwdalign::MilleRecordWriter()),
49 mWithConstraintsRecWriter(false),
50 mConstraintsRecWriter(nullptr),
51 mMillepede(new o2::fwdalign::MillePede2())
52{
55 }
56 // initialise the content of each array
59 LOGF(debug, "TracksToRecords instantiated");
60}
61
62//__________________________________________________________________________
64{
67 }
68 if (mMillepede) {
69 delete mMillepede;
70 }
71 if (mRecordWriter) {
72 delete mRecordWriter;
73 }
74 if (mAlignPoint) {
75 delete mAlignPoint;
76 }
77 if (mDictionary) {
78 mDictionary = nullptr;
79 }
80 LOGF(debug, "TracksToRecords destroyed");
81}
82
83//__________________________________________________________________________
85{
86 if (mIsInitDone) {
87 return;
88 }
89 if (mDictionary == nullptr) {
90 LOGF(fatal, "TracksToRecords::init() failed because no cluster dictionary is defined");
91 mIsInitDone = false;
92 return;
93 }
94
98
103 }
104
106
110 mResCut,
112
113 LOG(info) << "-------------- TracksToRecords configured with -----------------";
114 LOGF(info, "Chi2CutNStdDev = %d", mChi2CutNStdDev);
115 LOGF(info, "ResidualCutInitial = %.3f", mResCutInitial);
116 LOGF(info, "ResidualCut = %.3f", mResCut);
117 LOGF(info, "MinNumberClusterCut = %d", mMinNumberClusterCut);
118 LOGF(info, "mStartFac = %.3f", mStartFac);
119 LOGF(info,
120 "Allowed variation: dx = %.3f, dy = %.3f, dz = %.3f, dRz = %.4f",
121 mAllowVar[0], mAllowVar[1], mAllowVar[3], mAllowVar[2]);
122 LOG(info) << "-----------------------------------------------------------";
123
124 // set allowed variations for all parameters
125 for (int chipId = 0; chipId < mNumberOfSensors; ++chipId) {
126 for (Int_t iPar = 0; iPar < mNDofPerSensor; ++iPar) {
127 mMillepede->SetParSigma(chipId * mNDofPerSensor + iPar, mAllowVar[iPar]);
128 }
129 }
130
131 // set iterations
132 if (mStartFac > 1) {
134 }
135
136 mIsInitDone = true;
137 LOGF(info, "TracksToRecords init done");
138}
139
140//__________________________________________________________________________
142{
143 mNumberTFs++; // TF Counter
144
145 // get tracks
146 mMFTTracks = ctx.inputs().get<gsl::span<o2::mft::TrackMFT>>("tracks");
147 mMFTTracksROF = ctx.inputs().get<gsl::span<o2::itsmft::ROFRecord>>("tracksrofs");
148 mMFTTrackClusIdx = ctx.inputs().get<gsl::span<int>>("trackClIdx");
149
150 // get clusters
151 mMFTClusters = ctx.inputs().get<gsl::span<o2::itsmft::CompClusterExt>>("compClusters");
152 mMFTClustersROF = ctx.inputs().get<gsl::span<o2::itsmft::ROFRecord>>("clustersrofs");
153 mMFTClusterPatterns = ctx.inputs().get<gsl::span<unsigned char>>("patterns");
157}
158
159//__________________________________________________________________________
161{
162 if (!mIsInitDone) {
163 LOGF(fatal, "TracksToRecords::processRecoTracks() aborted because init was not done !");
164 return;
165 }
167 LOGF(fatal, "TracksToRecords::processRecoTracks() aborted because uninitialised mRecordWriter !");
168 return;
169 }
170
171 LOG(info) << "TracksToRecords::processRecoTracks() - start";
172
173 int nCounterAllTracks = 0;
174
175 for (auto& oneTrack : mMFTTracks) { // track loop
176
177 LOGF(debug, "Processing track # %5d", nCounterAllTracks);
178
179 // Skip the track if not enough clusters
180 auto ncls = oneTrack.getNumberOfPoints();
181 if (ncls < mMinNumberClusterCut) {
182 nCounterAllTracks++;
184 continue;
185 }
186
187 // Skip presumably quite low momentum track
188 if (!oneTrack.isLTF()) {
189 nCounterAllTracks++;
191 continue;
192 }
193
194 auto offset = oneTrack.getExternalClusterIndexOffset();
195
197
198 // Store the initial track parameters
199 auto track = oneTrack;
202
203 bool isTrackUsed = true;
204
205 for (int icls = 0; icls < ncls; ++icls) { // cluster loop
206
208
209 // Store measured positions
210 auto clsEntry = mMFTTrackClusIdx[offset + icls];
211 auto localCluster = mMFTClustersLocal[clsEntry];
212 auto globalCluster = mMFTClustersGlobal[clsEntry];
213 mAlignPoint->setMeasuredPosition(localCluster, globalCluster);
214 if (!mAlignPoint->isClusterOk()) {
215 LOGF(warning, "TracksToRecords::processRecoTracks() - will not use track # %5d with at least a bad cluster", nCounterAllTracks);
217 isTrackUsed = false;
218 break;
219 }
220
221 // Propagate track to the current z plane of this cluster
222 track.propagateParamToZlinear(mAlignPoint->getGlobalMeasuredPosition().Z());
223
224 // Store reco positions
226
227 // compute residuals
230
231 // Compute derivatives
234
235 // Set local equations
236 bool success = true;
237 success &= setLocalEquationX();
238 success &= setLocalEquationY();
239 success &= setLocalEquationZ();
240 isTrackUsed &= success;
241 if (mWithControl && success) {
243 }
244 if (!success) {
245 LOGF(error, "TracksToRecords::processRecoTracks() - track %i h %d d %d l %d s %4d lMpos x %.2e y %.2e z %.2e gMpos x %.2e y %.2e z %.2e gRpos x %.2e y %.2e z %.2e",
250 }
251
252 } // end of loop on clusters
253
254 if (isTrackUsed) {
257 const bool doPrint = false;
258 mRecordWriter->fillRecordTree(doPrint); // save record data
260 }
261 nCounterAllTracks++;
262 } // end of loop on tracks
263 LOG(info) << "TracksToRecords::processRecoTracks() - end";
264}
265
266//__________________________________________________________________________
267void TracksToRecords::processROFs(TChain* mfttrackChain, TChain* mftclusterChain)
268{
269 if (!mIsInitDone) {
270 LOGF(fatal, "TracksToRecords::processROFs() aborted because init was not done !");
271 return;
272 }
273
275 LOGF(fatal, "TracksToRecords::processROFs() aborted because uninitialised mRecordWriter !");
276 return;
277 }
278
279 LOG(info) << "TracksToRecords::processROFs() - start";
280
281 TTreeReader mftTrackChainReader(mfttrackChain);
282 TTreeReader mftClusterChainReader(mftclusterChain);
283 std::vector<unsigned char>::iterator pattIterator;
284
285 TTreeReaderValue<std::vector<o2::mft::TrackMFT>> mftTracks =
286 {mftTrackChainReader, "MFTTrack"};
287 TTreeReaderValue<std::vector<o2::itsmft::ROFRecord>> mftTracksROF =
288 {mftTrackChainReader, "MFTTracksROF"};
289 TTreeReaderValue<std::vector<int>> mftTrackClusIdx =
290 {mftTrackChainReader, "MFTTrackClusIdx"};
291
292 TTreeReaderValue<std::vector<o2::itsmft::CompClusterExt>> mftClusters =
293 {mftClusterChainReader, "MFTClusterComp"};
294 TTreeReaderValue<std::vector<o2::itsmft::ROFRecord>> mftClustersROF =
295 {mftClusterChainReader, "MFTClustersROF"};
296 TTreeReaderValue<std::vector<unsigned char>> mftClusterPatterns =
297 {mftClusterChainReader, "MFTClusterPatt"};
298
299 int nCounterAllTracks = 0;
300
301 while (mftTrackChainReader.Next() && mftClusterChainReader.Next()) {
302
303 mNumberOfTrackChainROFs += (*mftTracksROF).size();
304 mNumberOfClusterChainROFs += (*mftClustersROF).size();
306
307 pattIterator = (*mftClusterPatterns).begin();
309 *mftClusters, pattIterator, mMFTClustersLocal, mMFTClustersGlobal);
310
311 //______________________________________________________
312 for (auto& oneTrack : *mftTracks) { // track loop
313
314 LOGF(debug, "Processing track # %5d", nCounterAllTracks);
315
316 // Skip the track if not enough clusters
317 auto ncls = oneTrack.getNumberOfPoints();
318 if (ncls < mMinNumberClusterCut) {
319 nCounterAllTracks++;
321 continue;
322 }
323
324 // Skip presumably quite low momentum track
325 if (!oneTrack.isLTF()) {
326 nCounterAllTracks++;
328 continue;
329 }
330
331 auto offset = oneTrack.getExternalClusterIndexOffset();
332
334
335 // Store the initial track parameters
338
339 bool isTrackUsed = true;
340
341 for (int icls = 0; icls < ncls; ++icls) { // cluster loop
342
344
345 // Store measured positions
346 auto clsEntry = (*mftTrackClusIdx)[offset + icls];
347 auto localCluster = mMFTClustersLocal[clsEntry];
348 auto globalCluster = mMFTClustersGlobal[clsEntry];
349 mAlignPoint->setMeasuredPosition(localCluster, globalCluster);
350 if (!mAlignPoint->isClusterOk()) {
351 LOGF(warning, "TracksToRecords::processROFs() - will not use track # %5d with at least a bad cluster", nCounterAllTracks);
353 isTrackUsed = false;
354 break;
355 }
356
357 // Propagate track to the current z plane of this cluster
358 oneTrack.propagateParamToZlinear(mAlignPoint->getGlobalMeasuredPosition().Z());
359
360 // Store reco positions
362
363 // compute residuals
366
367 // Compute derivatives
370
371 // Set local equations
372 bool success = true;
373 success &= setLocalEquationX();
374 success &= setLocalEquationY();
375 success &= setLocalEquationZ();
376 isTrackUsed &= success;
377 if (mWithControl && success) {
379 }
380 if (!success) {
381 LOGF(error, "TracksToRecords::processROFs() - track %i h %d d %d l %d s %4d lMpos x %.2e y %.2e z %.2e gMpos x %.2e y %.2e z %.2e gRpos x %.2e y %.2e z %.2e",
386 }
387
388 } // end of loop on clusters
389
390 if (isTrackUsed) {
391 // copy track record
394 const bool doPrint = false;
395 mRecordWriter->fillRecordTree(doPrint); // save record data
397 }
398 nCounterAllTracks++;
399 } // end of loop on tracks
400
401 } // end of loop on TChain reader
402
403 LOG(info) << "TracksToRecords::processROFs() - end";
404}
405
406//__________________________________________________________________________
408{
409 LOGF(info, "TracksToRecords processRecoTracks() summary: ");
411 LOGF(info,
412 "n ROFs = %d, used tracks = %d, skipped tracks = %d, local equations failed = %d",
415 } else {
416 LOGF(info,
417 "n TFs = %d, used tracks = %d, skipped tracks = %d, local equations failed = %d",
420 }
421}
422
423//__________________________________________________________________________
434
435//__________________________________________________________________________
437{
438 if (mRecordWriter) {
439 mRecordWriter->terminate(); // write record tree and close output file
440 }
441 if (mWithControl) {
443 }
444}
445
446//__________________________________________________________________________
457
458//__________________________________________________________________________
468
469//__________________________________________________________________________
471{
472 // index [0 .. 3] for {dX0, dTx, dY0, dTz}
473
474 bool success = false;
477 success = true;
478 } else {
479 LOGF(error,
480 "AlignHelper::setLocalDerivative() - index %d >= %d",
482 }
483 return success;
484}
485
486//__________________________________________________________________________
488{
489 // index [0 .. 3] for {dDeltaX, dDeltaY, dDeltaRz, dDeltaZ}
490
491 bool success = false;
494 success = true;
495 } else {
496 LOGF(error,
497 "AlignHelper::setGlobalDerivative() - index %d >= %d",
499 }
500 return success;
501}
502
503//__________________________________________________________________________
505{
506 bool success = false;
507 std::fill(mLocalDerivatives.begin(), mLocalDerivatives.end(), 0.);
508 success = true;
509 return success;
510}
511
512//__________________________________________________________________________
514{
515 bool success = false;
516 std::fill(mGlobalDerivatives.begin(), mGlobalDerivatives.end(), 0.);
517 success = true;
518 return success;
519}
520
521//__________________________________________________________________________
523{
524
526 LOGF(error,
527 "TracksToRecords::setLocalEquationX() - no align point coordinates set !");
528 return false;
529 }
531 return false;
532 }
534 return false;
535 }
536
537 bool success = true;
538
539 // clean slate for the local equation for this measurement
540
541 success &= resetGlocalDerivative();
542 success &= resetLocalDerivative();
543
544 // local derivatives
545 // index [0 .. 3] for {dX0, dTx, dY0, dTz}
546
551
552 // global derivatives
553 // index [0 .. 3] for {dDeltaX, dDeltaY, dDeltaRz, dDeltaZ}
554
555 Int_t chipId = mAlignPoint->getSensorId();
560
561 if (success) {
562 if (mCounterUsedTracks < 5) {
563 LOGF(debug,
564 "TracksToRecords::setLocalEquationX(): track %i sr %4d local %.3e %.3e %.3e %.3e, global %.3e %.3e %.3e %.3e X %.3e sigma %.3e",
565 mCounterUsedTracks, chipId,
573 }
579 } else {
581 }
582
583 return success;
584}
585
586//__________________________________________________________________________
588{
590 LOGF(error,
591 "TracksToRecords::setLocalEquationY() - no align point coordinates set !");
592 return false;
593 }
595 return false;
596 }
598 return false;
599 }
600
601 bool success = true;
602
603 // clean slate for the local equation for this measurement
604
605 success &= resetGlocalDerivative();
606 success &= resetLocalDerivative();
607
608 // local derivatives
609 // index [0 .. 3] for {dX0, dTx, dY0, dTz}
610
615
616 // global derivatives
617 // index [0 .. 3] for {dDeltaX, dDeltaY, dDeltaRz, dDeltaZ}
618
619 Int_t chipId = mAlignPoint->getSensorId();
624
625 if (success) {
626 if (mCounterUsedTracks < 5) {
627 LOGF(debug,
628 "TracksToRecords::setLocalEquationY(): track %i sr %4d local %.3e %.3e %.3e %.3e, global %.3e %.3e %.3e %.3e Y %.3e sigma %.3e",
629 mCounterUsedTracks, chipId,
637 }
643 } else {
645 }
646
647 return success;
648}
649
650//__________________________________________________________________________
652{
654 LOGF(error,
655 "TracksToRecords::setLocalEquationZ() - no align point coordinates set !");
656 return false;
657 }
659 return false;
660 }
662 return false;
663 }
664
665 bool success = true;
666
667 // clean slate for the local equation for this measurement
668
669 success &= resetGlocalDerivative();
670 success &= resetLocalDerivative();
671
672 // local derivatives
673 // index [0 .. 3] for {dX0, dTx, dY0, dTz}
674
679
680 // global derivatives
681 // index [0 .. 3] for {dDeltaX, dDeltaY, dDeltaRz, dDeltaZ}
682
683 Int_t chipId = mAlignPoint->getSensorId();
688
689 if (success) {
690 if (mCounterUsedTracks < 5) {
691 LOGF(debug,
692 "TracksToRecords::setLocalEquationZ(): track %i sr %4d local %.3e %.3e %.3e %.3e, global %.3e %.3e %.3e %.3e Z %.3e sigma %.3e",
693 mCounterUsedTracks, chipId,
701 }
707 } else {
709 }
710
711 return success;
712}
Class handling both virtual segmentation and real volumes.
Class to store the data of single track processing.
ClassImp(o2::mft::TracksToRecords)
Class responsible to create records from tracks (and attached clusters) to feed alignment.
std::ostringstream debug
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
int InitMille(int nGlo, const int nLoc, const int lNStdDev=-1, const double lResCut=-1., const double lResCutInit=-1., const std::vector< int > &regroup={})
init all
void SetParSigma(int i, double par)
Definition MillePede2.h:278
void SetConstraintsRecWriter(o2::fwdalign::MilleRecordWriter *myP)
Definition MillePede2.h:268
void SetLocalEquation(std::vector< double > &dergb, std::vector< double > &derlc, const double lMeas, const double lSigma)
assing derivs of loc.eq.
int SetIterations(const double lChi2CutFac)
Number of iterations is calculated from lChi2CutFac.
void SetRecordWriter(o2::fwdalign::MilleRecordWriter *myP)
Definition MillePede2.h:267
void setCyclicAutoSave(const long nEntries)
Set the number of entries to be used by TTree::AutoSave()
bool isInitOk() const
check if init went well
void setDataFileName(TString fname)
choose data records filename
void init()
init output file and tree
void setRecordRun(int run)
assign run
void terminate()
write tree and close output file
void changeDataBranchName(const bool isConstraintsRec=true)
choose data records filename
void setRecordWeight(double wgh)
assign weight
void fillRecordTree(const bool doPrint=false)
fill tree
o2::fwdalign::MillePedeRecord * getRecord()
return the record
void setCyclicAutoSave(const long nEntries)
Set the number of entries to be used by TTree::AutoSave()
void init()
init output file and tree
void fill(o2::mft::AlignPointHelper *aPoint, const int iTrack=0, const bool doPrint=false)
fill the tree from an align point
void terminate()
write tree and close output file
Container of a single alignment point and methods to fill it.
o2::math_utils::Point3D< double > getLocalMeasuredPositionSigma() const
GlobalDerivative globalDerivativeX() const
void recordTrackInitialParam(o2::mft::TrackMFT &mftTrack)
store the track parameters at the initial z0 plane
GlobalDerivative globalDerivativeZ() const
o2::math_utils::Point3D< double > getLocalResidual() const
void setClusterDictionary(const o2::itsmft::TopologyDictionary *d)
set cluster pattern dictionary (needed to compute cluster coordinates)
o2::math_utils::Point3D< double > getLocalMeasuredPosition() const
o2::math_utils::Point3D< double > getGlobalRecoPosition() const
LocalDerivative localDerivativeY() const
void setGlobalRecoPosition(o2::mft::TrackMFT &mftTrack)
LocalDerivative localDerivativeZ() const
void convertCompactClusters(gsl::span< const itsmft::CompClusterExt > clusters, gsl::span< const unsigned char >::iterator &pattIt, std::vector< o2::BaseCluster< double > > &outputLocalClusters, std::vector< o2::BaseCluster< double > > &outputGlobalClusters)
convert compact clusters (pixel coordinates in row, col) from workflow to base clusters with 3D posit...
void setMeasuredPosition(const o2::BaseCluster< double > &localCluster, const o2::BaseCluster< double > &globalCluster)
void computeGlobalDerivatives()
method to call the computation of all three components of the global derivative
void resetAlignPoint()
reset all quantities that define an alignment point to their default value
void computeLocalDerivatives()
method to call the computation of all three compnonents of the local derivative
o2::math_utils::Point3D< double > getGlobalMeasuredPosition() const
LocalDerivative localDerivativeX() const
void resetTrackInitialParam()
reset all track parameters to their default value (zero)
GlobalDerivative globalDerivativeY() const
TString mMilleRecordsFileName
output file name when saving the Mille records
Definition Aligner.h:71
double mResCutInitial
Cut on residual on first iteration.
Definition Aligner.h:69
static constexpr int mNumberOfSensors
Total number of sensors (detection elements) in the MFT.
Definition Aligner.h:64
TString mMilleConstraintsRecFileName
output file name when saving the records of the constraints
Definition Aligner.h:72
static constexpr int mNDofPerSensor
translation in global x, y, z, and rotation Rz around global z-axis
Definition Aligner.h:62
double mResCut
Cut on residual for other iterations.
Definition Aligner.h:70
bool mIsInitDone
boolean to follow the initialisation status
Definition Aligner.h:73
double mStartFac
Initial value for chi2 cut, used to reject outliers i.e. bad tracks with sum(chi2) > Chi2DoFLim(fNStd...
Definition Aligner.h:67
static constexpr int mNumberOfGlobalParam
Number of alignment (= global) parameters.
Definition Aligner.h:65
std::array< double, mNDofPerSensor > mAllowVar
"Encouraged" variation for degrees of freedom {dx, dy, dRz, dz}
Definition Aligner.h:66
static constexpr int mNumberOfTrackParam
Number of track (= local) parameters (X0, Tx, Y0, Ty)
Definition Aligner.h:61
int mChi2CutNStdDev
Number of standard deviations for chi2 cut.
Definition Aligner.h:68
bool setGlobalDerivative(Int_t index, Double_t value)
set array of global derivatives
void processRecoTracks()
use valid tracks (and associated clusters) from the workflow to build Mille records
o2::fwdalign::MilleRecordWriter * mRecordWriter
utility that handles the writing of the data records to a ROOT file
void startConstraintsRecWriter()
init the utility needed to write constraints records
o2::mft::AlignPointHelper * mAlignPoint
Alignment point helper.
o2::fwdalign::MilleRecordWriter * mConstraintsRecWriter
utility that handles the writing of the constraints records
gsl::span< const int > mMFTTrackClusIdx
std::vector< double > mGlobalDerivatives
vector of global derivatives {dDeltaX, dDeltaY, dDeltaRz, dDeltaZ}
long mNEntriesAutoSave
number of entries needed to call AutoSave for the output TTrees
void printProcessTrackSummary()
print a summary status of what happened in processRecoTracks() or processROFs()
bool mWithConstraintsRecWriter
boolean to be set to true if one wants to also write constaints records
bool setLocalEquationY()
set the 2nd component of the local equation vector for a given alignment point
int mCounterSkippedTracks
count how many tracks did not met the cut on the min. nb of clusters
int mNumberOfTrackChainROFs
number of ROFs in the track chain
gsl::span< const o2::mft::TrackMFT > mMFTTracks
void processTimeFrame(o2::framework::ProcessingContext &ctx)
access mft tracks and clusters in the timeframe provided by the workflow
void processROFs(TChain *mfttrackChain, TChain *mftclusterChain)
use mft tracks and clusters provided by ROOT files accessed via TChain to build Mille records
bool resetLocalDerivative()
reset the array of the Local derivative
int mMinNumberClusterCut
Minimum number of clusters in the track to be used for alignment.
std::vector< o2::BaseCluster< double > > mMFTClustersGlobal
MFT clusters in global coordinate system.
bool mWithControl
boolean to set the use of the control tree
void endConstraintsRecWriter()
end the utility used to write constraints records
int mCounterUsedTracks
count how many tracks were used to make Mille records
gsl::span< constunsignedchar >::iterator mPattIt
std::vector< o2::BaseCluster< double > > mMFTClustersLocal
MFT clusters in local coordinate system.
bool setLocalEquationZ()
set the last component of the local equation vector for a given alignment point
gsl::span< const unsigned char > mMFTClusterPatterns
int mNumberTFs
number of timeframes processed
gsl::span< const o2::itsmft::ROFRecord > mMFTClustersROF
gsl::span< const o2::itsmft::CompClusterExt > mMFTClusters
int mNumberOfClusterChainROFs
number of ROFs in the cluster chain
double mWeightRecord
the weight given to a single Mille record in Millepede algorithm
const o2::itsmft::TopologyDictionary * mDictionary
cluster patterns dictionary
void init() override
init Millipede and AlignPointHelper
bool resetGlocalDerivative()
reset the array of the Global derivative
void endRecordWriter()
end the utility used to write data records and its control tree
gsl::span< const o2::itsmft::ROFRecord > mMFTTracksROF
o2::fwdalign::MillePede2 * mMillepede
Millepede2 implementation copied from AliROOT.
~TracksToRecords() override
destructor
void startRecordWriter()
init the utility needed to write data records and its control tree
int mCounterLocalEquationFailed
count how many times we failed to set a local equation
o2::mft::AlignPointControl mPointControl
AlignPointControl handles the control tree.
bool setLocalEquationX()
set the first component of the local equation vector for a given alignment point
bool setLocalDerivative(Int_t index, Double_t value)
set array of local derivatives
std::vector< double > mLocalDerivatives
vector of local derivatives {dX0, dTx, dY0, dTz}
GLuint index
Definition glcorearb.h:781
GLsizei const GLfloat * value
Definition glcorearb.h:819
GLintptr offset
Definition glcorearb.h:660
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Defining DataPointCompositeObject explicitly as copiable.
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"