]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAlignmentTracks.cxx
Adding new option to switch off the transformation if the alignment was done (M.Ivanov)
[u/mrichter/AliRoot.git] / STEER / AliAlignmentTracks.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
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 //-----------------------------------------------------------------
24
25 #include <TChain.h>
26 #include <TFile.h>
27
28 #include "AliAlignmentTracks.h"
29 #include "AliTrackPointArray.h"
30 #include "AliAlignObjAngles.h"
31 #include "AliTrackFitterRieman.h"
32 #include "AliTrackResidualsChi2.h"
33 #include "AliESD.h"
34 #include "AliLog.h"
35
36 ClassImp(AliAlignmentTracks)
37
38 //______________________________________________________________________________
39 AliAlignmentTracks::AliAlignmentTracks():
40   fESDChain(0),
41   fPointsFilename("AliTrackPoints.root"),
42   fPointsFile(0),
43   fPointsTree(0),
44   fLastIndex(0),
45   fArrayIndex(0),
46   fIsIndexBuilt(kFALSE),
47   fMisalignObjs(0),
48   fTrackFitter(0),
49   fMinimizer(0),
50   fDoUpdate(kTRUE)
51 {
52   // Default constructor
53   InitIndex();
54   InitAlignObjs();
55 }
56
57 //______________________________________________________________________________
58 AliAlignmentTracks::AliAlignmentTracks(TChain *esdchain):
59   fESDChain(esdchain),
60   fPointsFilename("AliTrackPoints.root"),
61   fPointsFile(0),
62   fPointsTree(0),
63   fLastIndex(0),
64   fArrayIndex(0),
65   fIsIndexBuilt(kFALSE),
66   fMisalignObjs(0),
67   fTrackFitter(0),
68   fMinimizer(0),
69   fDoUpdate(kTRUE)
70 {
71   // Constructor in the case
72   // the user provides an already
73   // built TChain with ESD trees
74   InitIndex();
75   InitAlignObjs();
76 }
77
78
79 //______________________________________________________________________________
80 AliAlignmentTracks::AliAlignmentTracks(const char *esdfilename, const char *esdtreename):
81   fPointsFilename("AliTrackPoints.root"),
82   fPointsFile(0),
83   fPointsTree(0),
84   fLastIndex(0),
85   fArrayIndex(0),
86   fIsIndexBuilt(kFALSE),
87   fMisalignObjs(0),
88   fTrackFitter(0),
89   fMinimizer(0),
90   fDoUpdate(kTRUE)
91 {
92   // Constructor in the case
93   // the user provides a single ESD file
94   // or a directory containing ESD files
95   fESDChain = new TChain(esdtreename);
96   fESDChain->Add(esdfilename);
97
98   InitIndex();
99   InitAlignObjs();
100 }
101
102 //______________________________________________________________________________
103 AliAlignmentTracks::AliAlignmentTracks(const AliAlignmentTracks &alignment):
104   TObject(alignment)
105 {
106   // Copy constructor
107   // not implemented
108   AliWarning("Copy constructor not implemented!");
109 }
110
111 //______________________________________________________________________________
112 AliAlignmentTracks& AliAlignmentTracks::operator= (const AliAlignmentTracks& alignment)
113 {
114   // Asignment operator
115   // not implemented
116   if(this==&alignment) return *this;
117
118   AliWarning("Asignment operator not implemented!");
119
120   ((TObject *)this)->operator=(alignment);
121
122   return *this;
123 }
124
125 //______________________________________________________________________________
126 AliAlignmentTracks::~AliAlignmentTracks()
127 {
128   // Destructor
129   if (fESDChain) delete fESDChain;
130
131   DeleteIndex();
132   DeleteAlignObjs();
133
134   delete fTrackFitter;
135   delete fMinimizer;
136
137   if (fPointsFile) fPointsFile->Close();
138 }
139
140 //______________________________________________________________________________
141 void AliAlignmentTracks::AddESD(TChain *esdchain)
142 {
143   // Add a chain with ESD files
144   if (fESDChain)
145     fESDChain->Add(esdchain);
146   else
147     fESDChain = esdchain;
148 }
149
150 //______________________________________________________________________________
151 void AliAlignmentTracks::AddESD(const char *esdfilename, const char *esdtreename)
152 {
153   // Add a single file or
154   // a directory to the chain
155   // with the ESD files
156   if (fESDChain)
157     fESDChain->AddFile(esdfilename,TChain::kBigNumber,esdtreename);
158   else {
159     fESDChain = new TChain(esdtreename);
160     fESDChain->Add(esdfilename);
161   }
162 }
163
164 //______________________________________________________________________________
165 void AliAlignmentTracks::ProcessESD()
166 {
167   // Analyzes and filters ESD tracks
168   // Stores the selected track space points
169   // into the output file
170
171   if (!fESDChain) return;
172
173   AliESD *esd = 0;
174   fESDChain->SetBranchAddress("ESD",&esd);
175
176   // Open the output file
177   if (fPointsFilename.Data() == "") {
178     AliWarning("Incorrect output filename!");
179     return;
180   }
181
182   TFile *pointsFile = TFile::Open(fPointsFilename,"RECREATE");
183   if (!pointsFile || !pointsFile->IsOpen()) {
184     AliWarning(Form("Can't open %s !",fPointsFilename.Data()));
185     return;
186   }
187
188   TTree *pointsTree = new TTree("spTree", "Tree with track space point arrays");
189   const AliTrackPointArray *array = 0;
190   pointsTree->Branch("SP","AliTrackPointArray", &array);
191   Int_t ievent = 0;
192   while (fESDChain->GetEntry(ievent++)) {
193     if (!esd) break;
194     Int_t ntracks = esd->GetNumberOfTracks();
195     for (Int_t itrack=0; itrack < ntracks; itrack++) {
196       AliESDtrack * track = esd->GetTrack(itrack);
197       if (!track) continue;
198  
199      //  UInt_t status = AliESDtrack::kITSpid; 
200 //       status|=AliESDtrack::kTPCpid; 
201 //       status|=AliESDtrack::kTRDpid; 
202 //       if ((track->GetStatus() & status) != status) continue;
203
204       if (track->GetP() < 0.5) continue;
205
206       array = track->GetTrackPointArray();
207       pointsTree->Fill();
208     }
209   }
210
211   if (!pointsTree->Write()) {
212     AliWarning("Can't write the tree with track point arrays!");
213     return;
214   }
215
216   pointsFile->Close();
217 }
218
219 //______________________________________________________________________________
220 void AliAlignmentTracks::ProcessESD(TSelector *selector)
221 {
222   AliWarning(Form("ESD processing based on selector is not yet implemented (%p) !",selector));
223 }
224
225 //______________________________________________________________________________
226 void AliAlignmentTracks::BuildIndex()
227 {
228   // Build index of points tree entries
229   // Used for access based on the volume IDs
230   if (fIsIndexBuilt) return;
231
232   fIsIndexBuilt = kTRUE;
233
234   // Dummy object is created in order
235   // to initialize the volume paths
236   AliAlignObjAngles alobj;
237
238   TFile *fPointsFile = TFile::Open(fPointsFilename);
239   if (!fPointsFile || !fPointsFile->IsOpen()) {
240     AliWarning(Form("Can't open %s !",fPointsFilename.Data()));
241     return;
242   }
243   
244   //  AliTrackPointArray* array = new AliTrackPointArray;
245   AliTrackPointArray* array = 0;
246   fPointsTree = (TTree*) fPointsFile->Get("spTree");
247   if (!fPointsTree) {
248     AliWarning("No pointsTree found!");
249     return;
250   }
251   fPointsTree->SetBranchAddress("SP", &array);
252
253   Int_t nArrays = fPointsTree->GetEntries();
254   for (Int_t iArray = 0; iArray < nArrays; iArray++)
255     {
256       fPointsTree->GetEvent(iArray);
257       if (!array) continue;
258       for (Int_t ipoint = 0; ipoint < array->GetNPoints(); ipoint++) {
259         UShort_t volId = array->GetVolumeID()[ipoint];
260         // check if the volId is valid
261         if (!AliAlignObj::GetVolPath(volId)) {
262           AliError(Form("The volume id %d has no default volume path !",
263                         volId));
264           continue;
265         }
266         Int_t modId;
267         Int_t layerId = AliAlignObj::VolUIDToLayer(volId,modId)
268                       - AliAlignObj::kFirstLayer;
269         if (!fArrayIndex[layerId][modId]) {
270           //first entry for this volume
271           fArrayIndex[layerId][modId] = new TArrayI(1000);
272         }
273         else {
274           Int_t size = fArrayIndex[layerId][modId]->GetSize();
275           // If needed allocate new size
276           if (fLastIndex[layerId][modId] >= size)
277             fArrayIndex[layerId][modId]->Set(size + 1000);
278         }
279
280         // Check if the index is already filled
281         Bool_t fillIndex = kTRUE;
282         if (fLastIndex[layerId][modId] != 0) {
283           if ((*fArrayIndex[layerId][modId])[fLastIndex[layerId][modId]-1] == iArray)
284             fillIndex = kFALSE;
285         }
286         // Fill the index array and store last filled index
287         if (fillIndex) {
288           (*fArrayIndex[layerId][modId])[fLastIndex[layerId][modId]] = iArray;
289           fLastIndex[layerId][modId]++;
290         }
291       }
292     }
293
294 }
295
296 //______________________________________________________________________________
297 void AliAlignmentTracks::InitIndex()
298 {
299   // Initialize the index arrays
300   Int_t nLayers = AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer;
301   fLastIndex = new Int_t*[nLayers];
302   fArrayIndex = new TArrayI**[nLayers];
303   for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
304     fLastIndex[iLayer] = new Int_t[AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer)];
305     fArrayIndex[iLayer] = new TArrayI*[AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer)];
306     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer); iModule++) {
307       fLastIndex[iLayer][iModule] = 0;
308       fArrayIndex[iLayer][iModule] = 0;
309     }
310   }
311 }
312
313 //______________________________________________________________________________
314 void AliAlignmentTracks::ResetIndex()
315 {
316   // Reset the value of the last filled index
317   // Do not realocate memory
318
319   fIsIndexBuilt = kFALSE;
320   
321   for (Int_t iLayer = 0; iLayer < AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer; iLayer++) {
322     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer); iModule++) {
323       fLastIndex[iLayer][iModule] = 0;
324     }
325   }
326 }
327
328 //______________________________________________________________________________
329 void AliAlignmentTracks::DeleteIndex()
330 {
331   // Delete the index arrays
332   // Called by the destructor
333   for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
334     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer); iModule++) {
335       if (fArrayIndex[iLayer][iModule]) {
336         delete fArrayIndex[iLayer][iModule];
337         fArrayIndex[iLayer][iModule] = 0;
338       }
339     }
340     delete [] fLastIndex[iLayer];
341     delete [] fArrayIndex[iLayer];
342   }
343   delete [] fLastIndex;
344   delete [] fArrayIndex;
345 }
346
347 //______________________________________________________________________________
348 Bool_t AliAlignmentTracks::ReadAlignObjs(const char *alignObjFileName, const char* arrayName)
349 {
350   // Read alignment object from a file
351   // To be replaced by a call to CDB
352   AliWarning(Form("Method not yet implemented (%s in %s) !",arrayName,alignObjFileName));
353
354   return kFALSE;
355 }
356
357 //______________________________________________________________________________
358 void AliAlignmentTracks::InitAlignObjs()
359 {
360   // Initialize the alignment objects array
361   Int_t nLayers = AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer;
362   fAlignObjs = new AliAlignObj**[nLayers];
363   for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
364     fAlignObjs[iLayer] = new AliAlignObj*[AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer)];
365     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer); iModule++) {
366       UShort_t volid = AliAlignObj::LayerToVolUID(iLayer+ AliAlignObj::kFirstLayer,iModule);
367       fAlignObjs[iLayer][iModule] = new AliAlignObjAngles("",volid,0,0,0,0,0,0);
368     }
369   }
370 }
371
372 //______________________________________________________________________________
373 void AliAlignmentTracks::ResetAlignObjs()
374 {
375   // Reset the alignment objects array
376   for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
377     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer); iModule++)
378       fAlignObjs[iLayer][iModule]->SetPars(0,0,0,0,0,0);
379   }
380 }
381
382 //______________________________________________________________________________
383 void AliAlignmentTracks::DeleteAlignObjs()
384 {
385   // Delete the alignment objects array
386   for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
387     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer); iModule++)
388       if (fAlignObjs[iLayer][iModule])
389         delete fAlignObjs[iLayer][iModule];
390     delete [] fAlignObjs[iLayer];
391   }
392   delete [] fAlignObjs;
393   fAlignObjs = 0;
394 }
395
396 void AliAlignmentTracks::AlignDetector(AliAlignObj::ELayerID firstLayer,
397                                        AliAlignObj::ELayerID lastLayer,
398                                        AliAlignObj::ELayerID layerRangeMin,
399                                        AliAlignObj::ELayerID layerRangeMax,
400                                        Int_t iterations)
401 {
402   // Align detector volumes within
403   // a given layer range
404   // (could be whole detector).
405   // Tracks are fitted only within
406   // the range defined by the user.
407   Int_t nModules = 0;
408   for (Int_t iLayer = firstLayer; iLayer < lastLayer; iLayer++)
409     nModules += AliAlignObj::LayerSize(iLayer);
410   TArrayI volIds(nModules);
411
412   Int_t modnum = 0;
413   for (Int_t iLayer = firstLayer; iLayer < lastLayer; iLayer++) {
414     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer); iModule++) {
415       UShort_t volId = AliAlignObj::LayerToVolUID(iLayer,iModule);
416       volIds.AddAt(volId,modnum);
417       modnum++;
418     }
419   }
420
421   while (iterations > 0) {
422     AlignVolumes(&volIds,0x0,layerRangeMin,layerRangeMax);
423     iterations--;
424   }
425 }
426
427 //______________________________________________________________________________
428 void AliAlignmentTracks::AlignLayer(AliAlignObj::ELayerID layer,
429                                     AliAlignObj::ELayerID layerRangeMin,
430                                     AliAlignObj::ELayerID layerRangeMax,
431                                     Int_t iterations)
432 {
433   // Align detector volumes within
434   // a given layer.
435   // Tracks are fitted only within
436   // the range defined by the user.
437   Int_t nModules = AliAlignObj::LayerSize(layer);
438   TArrayI volIds(nModules);
439   for (Int_t iModule = 0; iModule < nModules; iModule++) {
440     UShort_t volId = AliAlignObj::LayerToVolUID(layer,iModule);
441     volIds.AddAt(volId,iModule);
442   }
443
444   while (iterations > 0) {
445     AlignVolumes(&volIds,0x0,layerRangeMin,layerRangeMax);
446     iterations--;
447   }
448 }
449
450 //______________________________________________________________________________
451 void AliAlignmentTracks::AlignVolume(UShort_t volId, UShort_t volIdFit,
452                                      Int_t iterations)
453 {
454   // Align single detector volume to
455   // another volume.
456   // Tracks are fitted only within
457   // the second volume.
458   TArrayI volIds(1);
459   volIds.AddAt(volId,0);
460   TArrayI volIdsFit(1);
461   volIdsFit.AddAt(volIdFit,0);
462
463   while (iterations > 0) {
464     AlignVolumes(&volIds,&volIdsFit);
465     iterations--;
466   }
467 }
468
469 //______________________________________________________________________________
470 void AliAlignmentTracks::AlignVolumes(const TArrayI *volids, const TArrayI *volidsfit,
471                                      AliAlignObj::ELayerID layerRangeMin,
472                                      AliAlignObj::ELayerID layerRangeMax,
473                                      Int_t iterations)
474 {
475   // Align a set of detector volumes.
476   // Tracks are fitted only within
477   // the range defined by the user
478   // (by layerRangeMin and layerRangeMax)
479   // or within the set of volidsfit
480   // Repeat the procedure 'iterations' times
481
482   Int_t nVolIds = volids->GetSize();
483   if (nVolIds == 0) {
484     AliError("Volume IDs array is empty!");
485     return;
486   }
487
488   // Load only the tracks with at least one
489   // space point in the set of volume (volids)
490   BuildIndex();
491   AliTrackPointArray **points;
492   // Start the iterations
493   while (iterations > 0) {
494     Int_t nArrays = LoadPoints(volids, points);
495     if (nArrays == 0) return;
496
497     AliTrackResiduals *minimizer = CreateMinimizer();
498     minimizer->SetNTracks(nArrays);
499     minimizer->InitAlignObj();
500     AliTrackFitter *fitter = CreateFitter();
501     for (Int_t iArray = 0; iArray < nArrays; iArray++) {
502       if (!points[iArray]) continue;
503       fitter->SetTrackPointArray(points[iArray], kFALSE);
504       if (fitter->Fit(volids,volidsfit,layerRangeMin,layerRangeMax) == kFALSE) continue;
505       AliTrackPointArray *pVolId,*pTrack;
506       fitter->GetTrackResiduals(pVolId,pTrack);
507       minimizer->AddTrackPointArrays(pVolId,pTrack);
508     }
509     minimizer->Minimize();
510
511     // Update the alignment object(s)
512     if (fDoUpdate) for (Int_t iVolId = 0; iVolId < nVolIds; iVolId++) {
513       UShort_t volid = (*volids)[iVolId];
514       Int_t iModule;
515       AliAlignObj::ELayerID iLayer = AliAlignObj::VolUIDToLayer(volid,iModule);
516       AliAlignObj *alignObj = fAlignObjs[iLayer-AliAlignObj::kFirstLayer][iModule];      
517       *alignObj *= *minimizer->GetAlignObj();
518       alignObj->Print("");
519     }
520
521     UnloadPoints(nArrays, points);
522     
523     iterations--;
524   }
525 }
526   
527 //______________________________________________________________________________
528 Int_t AliAlignmentTracks::LoadPoints(const TArrayI *volids, AliTrackPointArray** &points)
529 {
530   // Load track point arrays with at least
531   // one space point in a given set of detector
532   // volumes (array volids).
533   // Use the already created tree index for
534   // fast access.
535
536   if (!fPointsTree) {
537     AliError("Tree with the space point arrays not initialized!");
538     points = 0;
539     return 0;
540   }
541
542   Int_t nVolIds = volids->GetSize();
543   if (nVolIds == 0) {
544     AliError("Volume IDs array is empty!");
545     points = 0;
546     return 0;
547   }
548
549   Int_t nArrays = 0;
550   for (Int_t iVolId = 0; iVolId < nVolIds; iVolId++) {
551     UShort_t volid = (*volids)[iVolId];
552     Int_t iModule;
553     AliAlignObj::ELayerID iLayer = AliAlignObj::VolUIDToLayer(volid,iModule);
554
555     // In case of empty index
556     if (fLastIndex[iLayer-AliAlignObj::kFirstLayer][iModule] == 0) {
557       AliWarning(Form("There are no space-points belonging to the volume which is to be aligned (Volume ID =%d)!",volid));
558       continue;
559     }
560     nArrays += fLastIndex[iLayer-AliAlignObj::kFirstLayer][iModule];
561   }
562
563   if (nArrays == 0) {
564     AliError("There are no space-points belonging to all of the volumes which are to be aligned!");
565     points = 0;
566     return 0;
567   }
568
569   AliTrackPointArray* array = 0;
570   fPointsTree->SetBranchAddress("SP", &array);
571
572   // Allocate the pointer to the space-point arrays
573   points = new AliTrackPointArray*[nArrays];
574   for (Int_t i = 0; i < nArrays; i++) points[i] = 0x0;
575
576   // Init the array used to flag already loaded tree entries
577   Bool_t *indexUsed = new Bool_t[fPointsTree->GetEntries()];
578   for (Int_t i = 0; i < fPointsTree->GetEntries(); i++)
579     indexUsed[i] = kFALSE;
580
581   // Start the loop over the volume ids
582   Int_t iArray = 0;
583   for (Int_t iVolId = 0; iVolId < nVolIds; iVolId++) {
584     UShort_t volid = (*volids)[iVolId];
585     Int_t iModule;
586     AliAlignObj::ELayerID iLayer = AliAlignObj::VolUIDToLayer(volid,iModule);
587
588     Int_t nArraysId = fLastIndex[iLayer-AliAlignObj::kFirstLayer][iModule];
589     TArrayI *index = fArrayIndex[iLayer-AliAlignObj::kFirstLayer][iModule];
590     AliTrackPoint p;
591
592     for (Int_t iArrayId = 0; iArrayId < nArraysId; iArrayId++) {
593
594       // Get tree entry
595       Int_t entry = (*index)[iArrayId];
596       if (indexUsed[entry] == kTRUE) continue;
597       fPointsTree->GetEvent(entry);
598       if (!array) {
599         AliWarning("Wrong space point array index!");
600         continue;
601       }
602       indexUsed[entry] = kTRUE;
603
604       // Get the space-point array
605       Int_t nPoints = array->GetNPoints();
606       points[iArray] = new AliTrackPointArray(nPoints);
607       for (Int_t iPoint = 0; iPoint < nPoints; iPoint++) {
608         array->GetPoint(p,iPoint);
609         Int_t modnum;
610         AliAlignObj::ELayerID layer = AliAlignObj::VolUIDToLayer(p.GetVolumeID(),modnum);
611         // check if the layer id is valid
612         if ((layer < AliAlignObj::kFirstLayer) ||
613             (layer >= AliAlignObj::kLastLayer)) {
614           AliError(Form("Layer index is invalid: %d (%d -> %d) !",
615                         layer,AliAlignObj::kFirstLayer,AliAlignObj::kLastLayer-1));
616           continue;
617         }
618         if ((modnum >= AliAlignObj::LayerSize(layer)) ||
619             (modnum < 0)) {
620           AliError(Form("Module number inside layer %d is invalid: %d (0 -> %d)",
621                         layer,modnum,AliAlignObj::LayerSize(layer)));
622           continue;
623         }
624
625         // Misalignment is introduced here
626         // Switch it off in case of real
627         // alignment job!
628         if (fMisalignObjs) {
629           AliAlignObj *misalignObj = fMisalignObjs[layer-AliAlignObj::kFirstLayer][modnum];
630           if (misalignObj)
631             misalignObj->Transform(p);
632         }
633         // End of misalignment
634
635         AliAlignObj *alignObj = fAlignObjs[layer-AliAlignObj::kFirstLayer][modnum];
636         alignObj->Transform(p);
637         points[iArray]->AddPoint(iPoint,&p);
638       }
639       iArray++;
640     }
641   }
642
643
644   delete [] indexUsed;
645
646   return nArrays;
647 }
648
649 //______________________________________________________________________________
650 void AliAlignmentTracks::UnloadPoints(Int_t n, AliTrackPointArray **points)
651 {
652   // Unload track point arrays for a given
653   // detector volume
654   for (Int_t iArray = 0; iArray < n; iArray++)
655     delete points[iArray];
656   delete [] points;
657 }
658
659 //______________________________________________________________________________
660 AliTrackFitter *AliAlignmentTracks::CreateFitter()
661 {
662   // Check if the user has already supplied
663   // a track fitter object.
664   // If not, create a default one.
665   if (!fTrackFitter)
666     fTrackFitter = new AliTrackFitterRieman;
667
668   return fTrackFitter;
669 }
670
671 //______________________________________________________________________________
672 AliTrackResiduals *AliAlignmentTracks::CreateMinimizer()
673 {
674   // Check if the user has already supplied
675   // a track residuals minimizer object.
676   // If not, create a default one.
677   if (!fMinimizer)
678     fMinimizer = new AliTrackResidualsChi2;
679
680   return fMinimizer;
681 }
682
683 //______________________________________________________________________________
684 Bool_t AliAlignmentTracks::Misalign(const char *misalignObjFileName, const char* arrayName)
685 {
686   // The method reads from a file a set of AliAlignObj which are
687   // then used to apply misalignments directly on the track
688   // space-points. The method is supposed to be used only for
689   // fast development and debugging of the alignment algorithms.
690   // Be careful not to use it in the case of 'real' alignment
691   // scenario since it will bias the results.
692
693   // Initialize the misalignment objects array
694   Int_t nLayers = AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer;
695   fMisalignObjs = new AliAlignObj**[nLayers];
696   for (Int_t iLayer = 0; iLayer < (AliAlignObj::kLastLayer - AliAlignObj::kFirstLayer); iLayer++) {
697     fMisalignObjs[iLayer] = new AliAlignObj*[AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer)];
698     for (Int_t iModule = 0; iModule < AliAlignObj::LayerSize(iLayer + AliAlignObj::kFirstLayer); iModule++)
699       fMisalignObjs[iLayer][iModule] = 0x0;
700   }
701
702   // Open the misliagnment file and load the array with
703   // misalignment objects
704   TFile* inFile = TFile::Open(misalignObjFileName,"READ");
705   if (!inFile || !inFile->IsOpen()) {
706     AliError(Form("Could not open misalignment file %s !",misalignObjFileName));
707     return kFALSE;
708   }
709
710   TClonesArray* array = ((TClonesArray*) inFile->Get(arrayName));
711   if (!array) {
712     AliError(Form("Could not find misalignment array %s in the file %s !",arrayName,misalignObjFileName));
713     inFile->Close();
714     return kFALSE;
715   }
716   inFile->Close();
717
718   // Store the misalignment objects for further usage  
719   Int_t nObjs = array->GetEntriesFast();
720   AliAlignObj::ELayerID layerId; // volume layer
721   Int_t modId; // volume ID inside the layer
722   for(Int_t i=0; i<nObjs; i++)
723     {
724       AliAlignObj* alObj = (AliAlignObj*)array->UncheckedAt(i);
725       alObj->GetVolUID(layerId,modId);
726       fMisalignObjs[layerId-AliAlignObj::kFirstLayer][modId] = alObj;
727     }
728   return kTRUE;
729 }