1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 //-----------------------------------------------------------------
17 // Implementation of the alignment steering class
18 // It provides an access to the track space points
19 // written along the esd tracks. The class enables
20 // the user to plug any track fitter (deriving from
21 // AliTrackFitter class) and minimization fo the
22 // track residual sums (deriving from the AliTrackResiduals).
23 //-----------------------------------------------------------------
28 #include "AliAlignmentTracks.h"
29 #include "AliTrackPointArray.h"
30 #include "AliAlignObjAngles.h"
31 #include "AliTrackFitterRieman.h"
32 #include "AliTrackResidualsChi2.h"
36 ClassImp(AliAlignmentTracks)
38 //______________________________________________________________________________
39 AliAlignmentTracks::AliAlignmentTracks():
41 fPointsFilename("AliTrackPoints.root"),
46 fIsIndexBuilt(kFALSE),
51 // Default constructor
56 //______________________________________________________________________________
57 AliAlignmentTracks::AliAlignmentTracks(TChain *esdchain):
59 fPointsFilename("AliTrackPoints.root"),
64 fIsIndexBuilt(kFALSE),
69 // Constructor in the case
70 // the user provides an already
71 // built TChain with ESD trees
77 //______________________________________________________________________________
78 AliAlignmentTracks::AliAlignmentTracks(const char *esdfilename, const char *esdtreename):
79 fPointsFilename("AliTrackPoints.root"),
84 fIsIndexBuilt(kFALSE),
89 // Constructor in the case
90 // the user provides a single ESD file
91 // or a directory containing ESD files
92 fESDChain = new TChain(esdtreename);
93 fESDChain->Add(esdfilename);
99 //______________________________________________________________________________
100 AliAlignmentTracks::AliAlignmentTracks(const AliAlignmentTracks &alignment):
105 AliWarning("Copy constructor not implemented!");
108 //______________________________________________________________________________
109 AliAlignmentTracks& AliAlignmentTracks::operator= (const AliAlignmentTracks& alignment)
111 // Asignment operator
113 if(this==&alignment) return *this;
115 AliWarning("Asignment operator not implemented!");
117 ((TObject *)this)->operator=(alignment);
122 //______________________________________________________________________________
123 AliAlignmentTracks::~AliAlignmentTracks()
126 if (fESDChain) delete fESDChain;
134 if (fPointsFile) fPointsFile->Close();
137 //______________________________________________________________________________
138 void AliAlignmentTracks::AddESD(TChain *esdchain)
140 // Add a chain with ESD files
142 fESDChain->Add(esdchain);
144 fESDChain = esdchain;
147 //______________________________________________________________________________
148 void AliAlignmentTracks::AddESD(const char *esdfilename, const char *esdtreename)
150 // Add a single file or
151 // a directory to the chain
152 // with the ESD files
154 fESDChain->AddFile(esdfilename,TChain::kBigNumber,esdtreename);
156 fESDChain = new TChain(esdtreename);
157 fESDChain->Add(esdfilename);
161 //______________________________________________________________________________
162 void AliAlignmentTracks::ProcessESD()
164 // Analyzes and filters ESD tracks
165 // Stores the selected track space points
166 // into the output file
168 if (!fESDChain) return;
171 fESDChain->SetBranchAddress("ESD",&esd);
173 // Open the output file
174 if (fPointsFilename.Data() == "") {
175 AliWarning("Incorrect output filename!");
179 TFile *pointsFile = TFile::Open(fPointsFilename,"RECREATE");
180 if (!pointsFile || !pointsFile->IsOpen()) {
181 AliWarning(Form("Can't open %s !",fPointsFilename.Data()));
185 TTree *pointsTree = new TTree("spTree", "Tree with track space point arrays");
186 AliTrackPointArray *array = 0;
187 pointsTree->Branch("SP","AliTrackPointArray", &array);
190 while (fESDChain->GetEntry(ievent++)) {
192 Int_t ntracks = esd->GetNumberOfTracks();
193 for (Int_t itrack=0; itrack < ntracks; itrack++) {
194 AliESDtrack * track = esd->GetTrack(itrack);
195 if (!track) continue;
197 UInt_t status = AliESDtrack::kITSpid;
198 status|=AliESDtrack::kTPCpid;
199 status|=AliESDtrack::kTRDpid;
200 if ((track->GetStatus() & status) != status) continue;
202 if (track->GetP() < 0.5) continue;
204 array = track->GetTrackPointArray();
209 if (!pointsTree->Write()) {
210 AliWarning("Can't write the tree with track point arrays!");
217 //______________________________________________________________________________
218 void AliAlignmentTracks::ProcessESD(TSelector *selector)
220 AliWarning(Form("ESD processing based on selector is not yet implemented (%p) !",selector));
223 //______________________________________________________________________________
224 void AliAlignmentTracks::BuildIndex()
226 // Build index of points tree entries
227 // Used for access based on the volume IDs
228 if (fIsIndexBuilt) return;
230 fIsIndexBuilt = kTRUE;
232 TFile *fPointsFile = TFile::Open(fPointsFilename);
233 if (!fPointsFile || !fPointsFile->IsOpen()) {
234 AliWarning(Form("Can't open %s !",fPointsFilename.Data()));
238 // AliTrackPointArray* array = new AliTrackPointArray;
239 AliTrackPointArray* array = 0;
240 fPointsTree = (TTree*) fPointsFile->Get("spTree");
242 AliWarning("No pointsTree found!");
245 fPointsTree->SetBranchAddress("SP", &array);
247 Int_t nArrays = fPointsTree->GetEntries();
248 for (Int_t iArray = 0; iArray < nArrays; iArray++)
250 fPointsTree->GetEvent(iArray);
251 if (!array) continue;
252 for (Int_t ipoint = 0; ipoint < array->GetNPoints(); ipoint++) {
253 UShort_t volId = array->GetVolumeID()[ipoint];
255 Int_t layerId = AliAlignObj::VolUIDToLayer(volId,modId)
256 - AliAlignObj::kFirstLayer;
257 if (!fArrayIndex[layerId][modId]) {
258 //first entry for this volume
259 fArrayIndex[layerId][modId] = new TArrayI(1000);
262 Int_t size = fArrayIndex[layerId][modId]->GetSize();
263 // If needed allocate new size
264 if (fLastIndex[layerId][modId] >= size)
265 fArrayIndex[layerId][modId]->Set(size + 1000);
268 // Check if the index is already filled
269 Bool_t fillIndex = kTRUE;
270 if (fLastIndex[layerId][modId] != 0) {
271 if ((*fArrayIndex[layerId][modId])[fLastIndex[layerId][modId]-1] == iArray)
274 // Fill the index array and store last filled index
276 (*fArrayIndex[layerId][modId])[fLastIndex[layerId][modId]] = iArray;
277 fLastIndex[layerId][modId]++;
284 // //______________________________________________________________________________
285 // void AliAlignmentTracks::BuildIndexLayer(AliAlignObj::ELayerID layer)
289 // //______________________________________________________________________________
290 // void AliAlignmentTracks::BuildIndexVolume(UShort_t volid)
294 //______________________________________________________________________________
295 void AliAlignmentTracks::InitIndex()
297 // Initialize the index arrays
298 Int_t nLayers = AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer;
299 fLastIndex = new Int_t*[nLayers];
300 fArrayIndex = new TArrayI**[nLayers];
301 for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
302 fLastIndex[iLayer] = new Int_t[AliAlignObj::LayerSize(iLayer)];
303 fArrayIndex[iLayer] = new TArrayI*[AliAlignObj::LayerSize(iLayer)];
304 for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++) {
305 fLastIndex[iLayer][iModule] = 0;
306 fArrayIndex[iLayer][iModule] = 0;
311 //______________________________________________________________________________
312 void AliAlignmentTracks::ResetIndex()
314 // Reset the value of the last filled index
315 // Do not realocate memory
317 fIsIndexBuilt = kFALSE;
319 for (Int_t iLayer = 0; iLayer < AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer; iLayer++) {
320 for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++) {
321 fLastIndex[iLayer][iModule] = 0;
326 //______________________________________________________________________________
327 void AliAlignmentTracks::DeleteIndex()
329 // Delete the index arrays
330 // Called by the destructor
331 for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
332 for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++) {
333 if (fArrayIndex[iLayer][iModule]) {
334 delete fArrayIndex[iLayer][iModule];
335 fArrayIndex[iLayer][iModule] = 0;
338 delete [] fLastIndex[iLayer];
339 delete [] fArrayIndex[iLayer];
341 delete [] fLastIndex;
342 delete [] fArrayIndex;
345 //______________________________________________________________________________
346 Bool_t AliAlignmentTracks::ReadAlignObjs(const char *alignObjFileName, const char* arrayName)
348 // Read alignment object from a file
349 // To be replaced by a call to CDB
350 AliWarning(Form("Method not yet implemented (%s in %s) !",arrayName,alignObjFileName));
355 //______________________________________________________________________________
356 void AliAlignmentTracks::InitAlignObjs()
358 // Initialize the alignment objects array
359 Int_t nLayers = AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer;
360 fAlignObjs = new AliAlignObj**[nLayers];
361 for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
362 fAlignObjs[iLayer] = new AliAlignObj*[AliAlignObj::LayerSize(iLayer)];
363 for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++) {
364 UShort_t volid = AliAlignObj::LayerToVolUID(iLayer+ AliAlignObj::kFirstLayer,iModule);
365 fAlignObjs[iLayer][iModule] = new AliAlignObjAngles("",volid,0,0,0,0,0,0);
370 //______________________________________________________________________________
371 void AliAlignmentTracks::ResetAlignObjs()
373 // Reset the alignment objects array
374 for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
375 for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++)
376 fAlignObjs[iLayer][iModule]->SetPars(0,0,0,0,0,0);
380 //______________________________________________________________________________
381 void AliAlignmentTracks::DeleteAlignObjs()
383 // Delete the alignment objects array
384 for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
385 for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++)
386 if (fAlignObjs[iLayer][iModule])
387 delete fAlignObjs[iLayer][iModule];
388 delete [] fAlignObjs[iLayer];
390 delete [] fAlignObjs;
394 //______________________________________________________________________________
395 void AliAlignmentTracks::Align(Int_t iterations)
397 // This method is just an example
398 // how one can user AlignLayer and
399 // AlignVolume methods to construct
400 // a custom alignment procedure
401 while (iterations > 0) {
403 AlignLayer(AliAlignObj::kTPC1);
404 AlignLayer(AliAlignObj::kSSD2);
405 AlignLayer(AliAlignObj::kSSD1);
406 AlignLayer(AliAlignObj::kSDD2);
407 AlignLayer(AliAlignObj::kSDD1);
408 AlignLayer(AliAlignObj::kSPD2);
409 AlignLayer(AliAlignObj::kSPD1);
411 AlignLayer(AliAlignObj::kSPD2);
412 AlignLayer(AliAlignObj::kSDD1);
413 AlignLayer(AliAlignObj::kSDD2);
414 AlignLayer(AliAlignObj::kSSD1);
415 AlignLayer(AliAlignObj::kSSD2);
416 AlignLayer(AliAlignObj::kTPC1);
417 AlignLayer(AliAlignObj::kTPC2);
418 AlignLayer(AliAlignObj::kTRD1);
419 AlignLayer(AliAlignObj::kTRD2);
420 AlignLayer(AliAlignObj::kTRD3);
421 AlignLayer(AliAlignObj::kTRD4);
422 AlignLayer(AliAlignObj::kTRD5);
423 AlignLayer(AliAlignObj::kTRD6);
424 AlignLayer(AliAlignObj::kTOF);
426 AlignLayer(AliAlignObj::kTRD6);
427 AlignLayer(AliAlignObj::kTRD5);
428 AlignLayer(AliAlignObj::kTRD4);
429 AlignLayer(AliAlignObj::kTRD3);
430 AlignLayer(AliAlignObj::kTRD2);
431 AlignLayer(AliAlignObj::kTRD1);
432 AlignLayer(AliAlignObj::kTPC2);
436 //______________________________________________________________________________
437 void AliAlignmentTracks::AlignLayer(AliAlignObj::ELayerID layer,
438 AliAlignObj::ELayerID layerRangeMin,
439 AliAlignObj::ELayerID layerRangeMax,
442 // Align detector volumes within
444 // Tracks are fitted only within
445 // the range defined by the user.
446 Int_t nModules = AliAlignObj::LayerSize(layer - AliAlignObj::kFirstLayer);
447 while (iterations > 0) {
448 for (Int_t iModule = 0; iModule < nModules; iModule++) {
449 UShort_t volId = AliAlignObj::LayerToVolUID(layer,iModule);
450 AlignVolume(volId,layerRangeMin,layerRangeMax);
456 //______________________________________________________________________________
457 void AliAlignmentTracks::AlignVolume(UShort_t volid, UShort_t volidfit,
458 AliAlignObj::ELayerID layerRangeMin,
459 AliAlignObj::ELayerID layerRangeMax,
462 // Align a single detector volume.
463 // Tracks are fitted only within
464 // the range defined by the user
465 // (by layerRangeMin and layerRangeMax)
466 // or within the volid2
467 // Repeat the procedure 'iterations' times
469 // First take the alignment object to be updated
471 AliAlignObj::ELayerID iLayer = AliAlignObj::VolUIDToLayer(volid,iModule);
472 AliAlignObj *alignObj = fAlignObjs[iLayer-AliAlignObj::kFirstLayer][iModule];
474 // Then load only the tracks with at least one
475 // space point in the volume (volid)
477 AliTrackPointArray **points;
478 // Start the iterations
479 while (iterations > 0) {
480 Int_t nArrays = LoadPoints(volid, points);
481 if (nArrays == 0) return;
483 AliTrackResiduals *minimizer = CreateMinimizer();
484 minimizer->SetNTracks(nArrays);
485 minimizer->SetAlignObj(alignObj);
486 AliTrackFitter *fitter = CreateFitter();
487 for (Int_t iArray = 0; iArray < nArrays; iArray++) {
488 fitter->SetTrackPointArray(points[iArray], kFALSE);
489 fitter->Fit(volid,volidfit,layerRangeMin,layerRangeMax);
490 AliTrackPointArray *pVolId,*pTrack;
491 fitter->GetTrackResiduals(pVolId,pTrack);
492 minimizer->AddTrackPointArrays(pVolId,pTrack);
494 minimizer->Minimize();
495 *alignObj *= *minimizer->GetAlignObj();
497 UnloadPoints(nArrays, points);
503 //______________________________________________________________________________
504 Int_t AliAlignmentTracks::LoadPoints(UShort_t volid, AliTrackPointArray** &points)
506 // Load track point arrays with at least
507 // one space point in a given detector
509 // Use the already created tree index for
513 AliWarning("Tree with the space point arrays not initialized!");
519 AliAlignObj::ELayerID iLayer = AliAlignObj::VolUIDToLayer(volid,iModule);
520 Int_t nArrays = fLastIndex[iLayer-AliAlignObj::kFirstLayer][iModule];
522 // In case of empty index
528 AliTrackPointArray* array = 0;
529 fPointsTree->SetBranchAddress("SP", &array);
531 points = new AliTrackPointArray*[nArrays];
532 TArrayI *index = fArrayIndex[iLayer-AliAlignObj::kFirstLayer][iModule];
534 for (Int_t iArray = 0; iArray < nArrays; iArray++) {
535 fPointsTree->GetEvent((*index)[iArray]);
537 AliWarning("Wrong space point array index!");
540 Int_t nPoints = array->GetNPoints();
541 points[iArray] = new AliTrackPointArray(nPoints);
542 for (Int_t iPoint = 0; iPoint < nPoints; iPoint++) {
543 array->GetPoint(p,iPoint);
545 AliAlignObj::ELayerID layer = AliAlignObj::VolUIDToLayer(p.GetVolumeID(),modnum);
547 // Misalignment is introduced here
548 // Switch it off in case of real
551 AliAlignObj *misalignObj = fMisalignObjs[layer-AliAlignObj::kFirstLayer][modnum];
553 misalignObj->Transform(p);
555 // End of misalignment
557 AliAlignObj *alignObj = fAlignObjs[layer-AliAlignObj::kFirstLayer][modnum];
558 alignObj->Transform(p);
559 points[iArray]->AddPoint(iPoint,&p);
566 //______________________________________________________________________________
567 void AliAlignmentTracks::UnloadPoints(Int_t n, AliTrackPointArray **points)
569 // Unload track point arrays for a given
571 for (Int_t iArray = 0; iArray < n; iArray++)
572 delete points[iArray];
576 //______________________________________________________________________________
577 AliTrackFitter *AliAlignmentTracks::CreateFitter()
579 // Check if the user has already supplied
580 // a track fitter object.
581 // If not, create a default one.
583 fTrackFitter = new AliTrackFitterRieman;
588 //______________________________________________________________________________
589 AliTrackResiduals *AliAlignmentTracks::CreateMinimizer()
591 // Check if the user has already supplied
592 // a track residuals minimizer object.
593 // If not, create a default one.
595 fMinimizer = new AliTrackResidualsChi2;
600 //______________________________________________________________________________
601 Bool_t AliAlignmentTracks::Misalign(const char *misalignObjFileName, const char* arrayName)
603 // The method reads from a file a set of AliAlignObj which are
604 // then used to apply misalignments directly on the track
605 // space-points. The method is supposed to be used only for
606 // fast development and debugging of the alignment algorithms.
607 // Be careful not to use it in the case of 'real' alignment
608 // scenario since it will bias the results.
610 // Initialize the misalignment objects array
611 Int_t nLayers = AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer;
612 fMisalignObjs = new AliAlignObj**[nLayers];
613 for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
614 fMisalignObjs[iLayer] = new AliAlignObj*[AliAlignObj::LayerSize(iLayer)];
615 for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++)
616 fMisalignObjs[iLayer][iModule] = 0x0;
619 // Open the misliagnment file and load the array with
620 // misalignment objects
621 TFile* inFile = TFile::Open(misalignObjFileName,"READ");
622 if (!inFile || !inFile->IsOpen()) {
623 AliError(Form("Could not open misalignment file %s !",misalignObjFileName));
627 TClonesArray* array = ((TClonesArray*) inFile->Get(arrayName));
629 AliError(Form("Could not find misalignment array %s in the file %s !",arrayName,misalignObjFileName));
635 // Store the misalignment objects for further usage
636 Int_t nObjs = array->GetEntriesFast();
637 AliAlignObj::ELayerID layerId; // volume layer
638 Int_t modId; // volume ID inside the layer
639 for(Int_t i=0; i<nObjs; i++)
641 AliAlignObj* alObj = (AliAlignObj*)array->UncheckedAt(i);
642 alObj->GetVolUID(layerId,modId);
643 fMisalignObjs[layerId-AliAlignObj::kFirstLayer][modId] = alObj;