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