]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/ITS/AliAlignmentDataFilterITS.cxx
11a5996774ad0619d38d7674bdaafbe1d15a7467
[u/mrichter/AliRoot.git] / PWG1 / ITS / AliAlignmentDataFilterITS.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2009, 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 //
18 // AliAnalysisTask to extract from ESD tracks the AliTrackPointArrays
19 // with ITS points for selected tracks. This are the input data for alignment
20 //
21 // Author: A.Dainese, andrea.dainese@pd.infn.it
22 /////////////////////////////////////////////////////////////
23
24 #include <TTree.h>
25 #include <TChain.h>
26 #include <TNtuple.h>
27 #include <TList.h>
28 #include <TH1F.h>
29 #include <TH2F.h>
30 #include <TMap.h>
31 #include <TVector3.h>
32 #include <TGeoManager.h>
33
34 #include "AliLog.h"
35 #include "AliGeomManager.h"
36 #include "AliITSReconstructor.h"
37 #include "AliITSgeomTGeo.h"
38 #include "AliTrackPointArray.h"
39 #include "AliESDInputHandler.h"
40 #include "AliESDVertex.h"
41 #include "AliESDtrack.h"
42 #include "AliESDEvent.h"
43 #include "AliESDfriend.h"
44 #include "AliAnalysisTask.h"
45 #include "AliAnalysisManager.h"
46 #include "AliAlignmentDataFilterITS.h"
47
48
49 ClassImp(AliAlignmentDataFilterITS)
50
51
52 //________________________________________________________________________
53 AliAlignmentDataFilterITS::AliAlignmentDataFilterITS(const char *name):
54 AliAnalysisTask(name,"task"),
55 fOnlySPDFO(kFALSE),
56 fGeometryFileName("geometry.root"),
57 fITSRecoParam(0),
58 fESD(0),
59 fESDfriend(0),
60 fListOfHistos(0),
61 fspTree(0),
62 fHistNevents(0),
63 fHistNpoints(0),
64 fHistPt(0),
65 fHistLayer0(0),
66 fHistLayer1(0),
67 fHistLayer2(0),
68 fHistLayer3(0),
69 fHistLayer4(0),
70 fHistLayer5(0),
71 fntExtra(0),
72 fntCosmicMatching(0)
73 {
74   // Constructor
75
76   // Define input and output slots here
77   // Input slot #0 works with a TChain
78   DefineInput(0, TChain::Class());
79   // Output slot #0 writes into a TTree
80   DefineOutput(0,TTree::Class());  //My private output
81   // Output slot #1 writes into a TList
82   DefineOutput(1,TList::Class());  //My private output
83 }
84
85 //________________________________________________________________________
86 AliAlignmentDataFilterITS::~AliAlignmentDataFilterITS()
87 {
88   // Destructor
89   if (fListOfHistos) {
90     delete fListOfHistos;
91     fListOfHistos = 0;
92   }
93   if (fspTree) {
94     delete fspTree;
95     fspTree = 0;
96   }
97   if (fHistNevents) {
98     delete fHistNevents;
99     fHistNevents = 0;
100   }
101   if (fHistNpoints) {
102     delete fHistNpoints;
103     fHistNpoints = 0;
104   }
105   if (fHistPt) {
106     delete fHistPt;
107     fHistPt = 0;
108   }
109   if (fHistLayer0) {
110     delete fHistLayer0;
111     fHistLayer0 = 0;
112   }
113   if (fHistLayer1) {
114     delete fHistLayer1;
115     fHistLayer1 = 0;
116   }
117   if (fHistLayer2) {
118     delete fHistLayer2;
119     fHistLayer2 = 0;
120   }
121   if (fHistLayer3) {
122     delete fHistLayer3;
123     fHistLayer3 = 0;
124   }
125   if (fHistLayer4) {
126     delete fHistLayer4;
127     fHistLayer4 = 0;
128   }
129   if (fHistLayer5) {
130     delete fHistLayer5;
131     fHistLayer5 = 0;
132   }
133   if (fntExtra) {
134     delete fntExtra;
135     fntExtra = 0;
136   }
137   if (fntCosmicMatching) {
138     delete fntCosmicMatching;
139     fntCosmicMatching = 0;
140   }
141 }  
142
143 //________________________________________________________________________
144 void AliAlignmentDataFilterITS::ConnectInputData(Option_t *) 
145 {
146   // Connect ESD
147   // Called once
148
149   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
150   if(!tree) {
151     printf("ERROR: Could not read chain from input slot 0\n");
152   } else {
153     // Get the OCDB path and the list of OCDB objects used for reco 
154     TMap *cdbMap = (TMap*)(tree->GetTree()->GetUserInfo())->FindObject("cdbMap");
155     TList *cdbList = (TList*)(tree->GetTree()->GetUserInfo())->FindObject("cdbList");
156
157     //cdbList->Print();
158     // write the list to the user info of the output tree
159     if(!fspTree) {
160       printf("ERROR: fspTree does not exist\n");
161     } else {
162       TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());         
163       cdbMapCopy->SetOwner(1);   
164       cdbMapCopy->SetName("cdbMap");     
165       TIter iter1(cdbMap->GetTable());   
166          
167       TPair* pair = 0;   
168       while((pair = dynamic_cast<TPair*> (iter1.Next()))){       
169         TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());    
170         TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());  
171         cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));   
172       }  
173
174       TList *cdbListCopy = new TList();  
175       cdbListCopy->SetOwner(1);  
176       cdbListCopy->SetName("cdbList");   
177       
178       TIter iter2(cdbList);      
179       
180       TObjString* cdbEntry=0;
181       while((cdbEntry =(TObjString*)(iter2.Next()))) {
182         cdbListCopy->Add(new TObjString(*cdbEntry));
183       }  
184       cdbListCopy->Print();
185
186
187       fspTree->GetUserInfo()->Add(cdbMapCopy);   
188       fspTree->GetUserInfo()->Add(cdbListCopy);
189     }
190
191     // Disable all branches and enable only the needed ones
192     tree->SetBranchStatus("fTriggerMask", 1);
193     tree->SetBranchStatus("fSPDVertex*", 1);
194     tree->SetBranchStatus("ESDfriend*", 1);
195     tree->SetBranchAddress("ESDfriend.",&fESDfriend);
196     tree->SetBranchStatus("fSPDMult*", 1);
197     
198     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
199     
200     if(!esdH) {
201       printf("ERROR: Could not get ESDInputHandler\n");
202     } else {
203       fESD = esdH->GetEvent();
204     }
205   }
206   
207   return;
208 }
209
210 //________________________________________________________________________
211 void AliAlignmentDataFilterITS::Init()
212 {
213   // Initialization
214
215   return;
216 }
217
218 //________________________________________________________________________
219 void AliAlignmentDataFilterITS::CreateOutputObjects()
220 {
221   // Create the output container
222   //
223   
224   // load the geometry  
225   if(!gGeoManager) {    
226     printf("AliAlignmentDataFilterITS::CreateOutputObjects(): loading geometry from %s\n",fGeometryFileName.Data());
227     AliGeomManager::LoadGeometry(fGeometryFileName.Data());
228     if(!gGeoManager) { 
229       printf("AliAlignmentDataFilterITS::CreateOutputObjects(): no geometry loaded \n");
230       return;
231     }
232   }
233
234   // Several histograms are more conveniently managed in a TList
235   fListOfHistos = new TList();
236   fListOfHistos->SetOwner();
237
238   fHistNevents = new TH1F("fHistNevents", "Number of processed events; N events; bin",5,-0.5,4.5);
239   fHistNevents->Sumw2();
240   fHistNevents->SetMinimum(0);
241   fListOfHistos->Add(fHistNevents);
242
243   fHistNpoints = new TH1F("fHistNpoints", "Number of AliTrackPoints per track; N points; tracks",25,-0.5,24.5);
244   fHistNpoints->Sumw2();
245   fHistNpoints->SetMinimum(0);
246   fListOfHistos->Add(fHistNpoints);
247
248   fHistPt = new TH1F("fHistPt", "p_{t} of tracks; p_{t} [GeV/c]; tracks",100,0,50);
249   fHistPt->Sumw2();
250   fHistPt->SetMinimum(0);
251   fListOfHistos->Add(fHistPt);
252
253
254   Float_t zmax=14.;
255   Int_t nbinsphi=20,nbinsz=4;
256   fHistLayer0 = new TH2F("fHistLayer0","Points in layer inner SPD; global   #phi; global z [cm]",nbinsphi,-3.14,3.14,nbinsz,-zmax,zmax);
257   fListOfHistos->Add(fHistLayer0);
258   zmax=14.;
259   nbinsphi=40;nbinsz=4;
260   fHistLayer1 = new TH2F("fHistLayer1","Points in layer outer SPD; global   #phi; global z [cm]",nbinsphi,-3.14,3.14,nbinsz,-zmax,zmax);
261   fListOfHistos->Add(fHistLayer1);
262   zmax=22.;
263   nbinsphi=14;nbinsz=6;
264   fHistLayer2 = new TH2F("fHistLayer2","Points in layer inner SDD; global   #phi; global z [cm]",nbinsphi,-3.14,3.14,nbinsz,-zmax,zmax);
265   fListOfHistos->Add(fHistLayer2);
266   zmax=29.5;
267   nbinsphi=22;nbinsz=8;
268   fHistLayer3 = new TH2F("fHistLayer3","Points in layer outer SDD; global   #phi; global z [cm]",nbinsphi,-3.14,3.14,nbinsz,-zmax,zmax);
269   fListOfHistos->Add(fHistLayer3);
270   zmax=45.;
271   nbinsphi=34;nbinsz=23;
272   fHistLayer4 = new TH2F("fHistLayer4","Points in layer inner SSD; global   #phi; global z [cm]",nbinsphi,-3.14,3.14,nbinsz,-zmax,zmax);
273   fListOfHistos->Add(fHistLayer4);
274   zmax=51.;
275   nbinsphi=38;nbinsz=26;
276   fHistLayer5 = new TH2F("fHistLayer5","Points in layer outer SSD; global   #phi; global z [cm]",nbinsphi,-3.14,3.14,nbinsz,-zmax,zmax);
277   fListOfHistos->Add(fHistLayer5);
278
279
280   fntExtra = new TNtuple("fntExtra","extra clusters in ITS","ncls:layer:ladder:volid:phi:x:y:z:xloc:zloc:dxy:dz:d0mu:z0mu");
281   fListOfHistos->Add(fntExtra);
282
283   fntCosmicMatching = new TNtuple("fntCosmicMatching","cosmic tracks matching in ITS","ncls1:ncls2:pt1:pt2:sigmad01:sigmad02:sigmaz01:sigmaz02:dxy:dz:phimu:thetamu:d0mu:z0mu");
284   fListOfHistos->Add(fntCosmicMatching);
285
286   fspTree = new TTree("spTree","Tree with ITS track points");
287   const AliTrackPointArray *array = 0;
288   Float_t curv,curverr;
289   fspTree->Branch("SP","AliTrackPointArray",&array);
290   fspTree->Branch("curv",&curv);
291   fspTree->Branch("curverr",&curverr);
292
293   return;
294 }
295
296 //________________________________________________________________________
297 void AliAlignmentDataFilterITS::Exec(Option_t */*option*/)
298 {
299   // Execute analysis for current event:
300   // write ITS AliTrackPoints for selected tracks to fspTree
301   
302   // check the geometry  
303   if(!gGeoManager) { 
304     printf("AliAlignmentDataFilterITS::Exec(): no geometry loaded \n");
305     return;
306   }
307
308   // check if we have AliITSRecoParam
309   if(!GetRecoParam()) {
310     if(!fITSRecoParam) {
311       printf("AliAlignmentDataFilterITS::Exec(): no AliITSRecoParam\n");
312       return;
313     }
314   }
315
316
317   if(!fESD) {
318     printf("AliAlignmentDataFilterITS::Exec(): no ESD \n");
319     return;
320   } 
321   if(!fESDfriend) {
322     printf("AliAlignmentDataFilterITS::Exec(): no ESDfriend \n");
323     return;
324   } 
325   // attach ESDfriend
326   fESD->SetESDfriend(fESDfriend);
327
328   // Post the data for slot 0
329   fHistNevents->Fill(0);
330
331
332   // write field value to spTree UserInfo
333   if(!((fspTree->GetUserInfo())->FindObject("BzkGauss"))) {
334     Double_t bz=fESD->GetMagneticField();
335     TString bzString; bzString+=bz;
336     TObjString *bzObjString = new TObjString(bzString);
337     TList *bzList = new TList();         
338     bzList->SetOwner(1);         
339     bzList->SetName("BzkGauss");         
340     bzList->Add(bzObjString);
341     fspTree->GetUserInfo()->Add(bzList);
342   }
343
344
345   // Process event as Cosmic or Collision
346   //if(esd->GetEventType()== ???? ) {
347   printf("AliAlignmentDataFilterITS::Exec(): MOVE ASAP TO esd->GetEventType() !\n");
348   if(GetRecoParam()->GetAlignFilterCosmics()) {
349     FilterCosmic(fESD);
350   } else {
351     FilterCollision(fESD);
352   }
353
354   PostData(1,fListOfHistos);
355
356   return;
357 }
358
359 //________________________________________________________________________
360 void AliAlignmentDataFilterITS::FilterCosmic(const AliESDEvent *esd)
361 {
362   // Extract ITS AliTrackPoints for Cosmics (check angular matching
363   // of top and bottom track, merge the two tracks, if requested)
364   //
365
366   // Set branch addresses for space points tree
367   AliTrackPointArray *arrayForTree=0;
368   Float_t curv,curverr;
369   fspTree->SetBranchAddress("SP",&arrayForTree);
370   fspTree->SetBranchAddress("curv",&curv);
371   fspTree->SetBranchAddress("curverr",&curverr);
372
373   TString triggeredClass = fESD->GetFiredTriggerClasses(); 
374   if(fOnlySPDFO && !triggeredClass.Contains("C0SCO-ABCE-NOPF-CENT")) return;
375
376
377   Int_t ntracks = esd->GetNumberOfTracks();
378   if(ntracks<2) return;
379
380   if(esd->GetPrimaryVertexSPD()->GetNContributors()<0) return;
381
382   Double_t vtxpos[3]; esd->GetPrimaryVertexSPD()->GetXYZ(vtxpos);
383
384   Int_t *goodtracksArray = new Int_t[ntracks];
385   Float_t *phiArray = new Float_t[ntracks];
386   Float_t *thetaArray = new Float_t[ntracks];
387   Int_t *nclsArray = new Int_t[ntracks];
388   Int_t ngt=0;
389   Int_t itrack=0;
390   for (itrack=0; itrack < ntracks; itrack++) {
391     AliESDtrack *track = esd->GetTrack(itrack);
392     if (!track) continue;
393
394
395     if(track->GetNcls(0)<GetRecoParam()->GetAlignFilterMinITSPoints()) continue;
396
397     if((GetRecoParam()->GetAlignFilterOnlyITSSATracks()) && track->GetNcls(1)>0) continue;
398     if((GetRecoParam()->GetAlignFilterOnlyITSTPCTracks()) && track->GetNcls(1)==0) continue;
399
400     Float_t phi = track->GetAlpha()+TMath::ASin(track->GetSnp());
401     Float_t theta = 0.5*TMath::Pi()-TMath::ATan(track->GetTgl());
402
403     if(track->Pt()<GetRecoParam()->GetAlignFilterMinPt() || 
404        track->Pt()>GetRecoParam()->GetAlignFilterMaxPt()) continue;
405
406     goodtracksArray[ngt] = itrack;
407     phiArray[ngt]        = phi;
408     thetaArray[ngt]      = theta;
409     nclsArray[ngt]       = track->GetNcls(0);
410     ngt++;
411   }
412
413   if(ngt<2) {
414     delete [] goodtracksArray; goodtracksArray=0;
415     delete [] phiArray; phiArray=0;
416     delete [] thetaArray; thetaArray=0;
417     delete [] nclsArray; nclsArray=0;
418     return;
419   }
420
421   // check matching of the two tracks from the muon
422   Float_t min = 10000000.;
423   Int_t maxCls = 0;
424   Int_t good1 = -1, good2 = -1;
425   for(Int_t itr1=0; itr1<ngt-1; itr1++) {
426     for(Int_t itr2=itr1+1; itr2<ngt; itr2++) {
427       Float_t deltatheta = TMath::Abs(TMath::Pi()-thetaArray[itr1]-thetaArray[itr2]);
428       if(deltatheta>GetRecoParam()->GetAlignFilterMaxMatchingAngle()) continue;
429       Float_t deltaphi = TMath::Abs(TMath::Abs(phiArray[itr1]-phiArray[itr2])-TMath::Pi());
430       if(deltaphi>GetRecoParam()->GetAlignFilterMaxMatchingAngle()) continue;
431       if(nclsArray[itr1]+nclsArray[itr2] > maxCls) {
432         maxCls = nclsArray[itr1]+nclsArray[itr2];
433         min = deltatheta+deltaphi;
434         good1 = goodtracksArray[itr1];
435         good2 = goodtracksArray[itr2];
436       } else if(nclsArray[itr1]+nclsArray[itr2] == maxCls) {
437         if(deltatheta+deltaphi < min) {
438           min = deltatheta+deltaphi;
439           good1 = goodtracksArray[itr1];
440           good2 = goodtracksArray[itr2];
441         }
442       }
443     }
444   }
445   
446   delete [] goodtracksArray; goodtracksArray=0;
447   delete [] phiArray; phiArray=0;
448   delete [] thetaArray; thetaArray=0;
449   delete [] nclsArray; nclsArray=0;
450
451   if(good1<0) return;
452   AliDebug(2,"ok track matching");
453   
454   // track1 will be the inward track (top)
455   // track2 the outward (bottom)
456   AliESDtrack *track1=0; 
457   AliESDtrack *track2=0;
458   AliESDtrack *track = esd->GetTrack(good1);
459   if(track->Py()>0) { 
460     track1 = esd->GetTrack(good1);
461     track2 = esd->GetTrack(good2);
462   } else {
463     track1 = esd->GetTrack(good2);
464     track2 = esd->GetTrack(good1);
465   }
466
467   AliTrackPoint point;
468   const AliTrackPointArray *array=0;
469   Int_t ipt,volId,modId,layerId,lay,lad,det;
470   Int_t jpt=0;
471   Bool_t layerOK[6][2]; 
472   Int_t nclsTrk[2]={0,0};
473
474   for(Int_t l1=0;l1<6;l1++) for(Int_t l2=0;l2<2;l2++) layerOK[l1][l2]=kFALSE;
475     
476   for(itrack=0; itrack<2; itrack++) {
477     if(itrack==0) {
478       track = track1;
479     } else {
480       track = track2;
481     }
482     array = track->GetTrackPointArray();
483     if(!array) {
484       AliWarning("No tracks points avaialble");
485       continue;
486     }
487     for(ipt=0; ipt<array->GetNPoints(); ipt++) {
488       array->GetPoint(point,ipt);
489       volId = point.GetVolumeID();
490       layerId = AliGeomManager::VolUIDToLayer(volId,modId);
491       AliDebug(2,Form("%d %d\n",ipt,layerId-1));
492       if(point.IsExtra() && 
493          (GetRecoParam()->GetAlignFilterSkipExtra())) continue;
494       if(layerId>6) continue;
495       if(!GetRecoParam()->GetAlignFilterUseLayer(layerId-1)) continue;
496       // check minAngleWrtITSModulePlanes
497       Double_t p[3]; track->GetDirection(p);
498       TVector3 pvec(p[0],p[1],p[2]);
499       Double_t rot[9]; AliGeomManager::GetOrigRotation(volId,rot);
500       TVector3 normvec(rot[1],rot[4],rot[7]);
501       Double_t angle = pvec.Angle(normvec);
502       if(angle>0.5*TMath::Pi()) angle = TMath::Pi()-angle;
503       angle = 0.5*TMath::Pi()-angle;
504       if(angle<GetRecoParam()->GetAlignFilterMinAngleWrtModulePlanes()) continue;
505       layerOK[layerId-1][itrack]=kTRUE;
506       jpt++;
507       nclsTrk[itrack]++;
508     }
509   }
510   AliDebug(2,Form("nClsTrk1 %d nClsTrk2 %d\n",nclsTrk[0],nclsTrk[1]));
511     
512   // read ITS cluster maps
513   Int_t map1[6],map2[6];
514   for(Int_t ilay=0;ilay<6;ilay++) {
515     map1[ilay]=0; map2[ilay]=0;
516     if(track1->HasPointOnITSLayer(ilay)) map1[ilay]=1;
517     if(track2->HasPointOnITSLayer(ilay)) map2[ilay]=1;
518   }
519   AliDebug(2,Form("ITS map 1: %d %d %d %d %d %d pt %f\n",map1[0],map1[1],map1[2],map1[3],map1[4],map1[5],track1->Pt()));
520   AliDebug(2,Form("ITS map 2: %d %d %d %d %d %d pt %f\n",map2[0],map2[1],map2[2],map2[3],map2[4],map2[5],track2->Pt()));
521   Int_t idx1[12],idx2[12];
522   track1->GetITSclusters(idx1);
523   track2->GetITSclusters(idx2);
524   AliDebug(2,Form("cls idx 1 %d %d %d %d %d %d %d %d %d %d %d %d\n",idx1[0],idx1[1],idx1[2],idx1[3],idx1[4],idx1[5],idx1[6],idx1[7],idx1[8],idx1[9],idx1[10],idx1[11]));
525   AliDebug(2,Form("cls idx 2 %d %d %d %d %d %d %d %d %d %d %d %d\n",idx2[0],idx2[1],idx2[2],idx2[3],idx2[4],idx2[5],idx2[6],idx2[7],idx2[8],idx2[9],idx2[10],idx2[11]));
526   
527
528   if(jpt<GetRecoParam()->GetAlignFilterMinITSPointsMerged()) return;
529   AliDebug(2,Form(" Total points %d, accepted\n",jpt));  
530   fHistNpoints->Fill(jpt);
531   fHistPt->Fill(0.5*(track1->Pt()+track2->Pt()));
532   
533   Float_t d0z0mu[2];
534   track1->GetDZ(0,0,0,esd->GetMagneticField(),d0z0mu);
535   //printf("d0mu %f  z0mu %f\n",d0z0mu[0],d0z0mu[1]);
536
537   Float_t dzOverlap[2];
538   Float_t curvArray[2],curverrArray[2];
539   Double_t globExtra[3],locExtra[3];
540   if(GetRecoParam()->GetAlignFilterCosmicMergeTracks()) 
541     arrayForTree = new AliTrackPointArray(jpt);
542   
543   jpt=0;
544   for(itrack=0; itrack<2; itrack++) {
545     if(itrack==0) {
546       track = track1;
547     } else {
548       track = track2;
549     }
550     curvArray[itrack] = track->GetC(esd->GetMagneticField());
551     curverrArray[itrack] = TMath::Sqrt(track->GetSigma1Pt2())*track->GetC(esd->GetMagneticField())/track->OneOverPt();
552
553     if(!(GetRecoParam()->GetAlignFilterCosmicMergeTracks())) {
554       jpt=0;
555       arrayForTree = new AliTrackPointArray(nclsTrk[itrack]);
556     }
557     array = track->GetTrackPointArray();
558     for(ipt=0; ipt<array->GetNPoints(); ipt++) {
559       array->GetPoint(point,ipt);
560       volId = point.GetVolumeID();
561       layerId = AliGeomManager::VolUIDToLayer(volId,modId);
562       if(layerId>6 || !layerOK[layerId-1][itrack]) continue;
563       arrayForTree->AddPoint(jpt,&point);
564       jpt++;
565       switch(layerId) {
566       case 1:
567         fHistLayer0->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
568         break;
569       case 2:
570         fHistLayer1->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
571         break;
572       case 3:
573         fHistLayer2->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
574         break;
575       case 4:
576         fHistLayer3->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
577         break;
578       case 5:
579         fHistLayer4->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
580         break;
581       case 6:
582         fHistLayer5->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
583         break;
584       }
585       // Post the data for slot 0
586       if(jpt==1) PostData(1,fListOfHistos); // only if this is the first points
587       if(!point.IsExtra() || 
588          !(GetRecoParam()->GetAlignFilterFillQANtuples())) continue;
589       nclsTrk[itrack]--;
590       for(Int_t ll=1;ll<layerId;ll++) modId+=AliITSgeomTGeo::GetNLadders(ll)*AliITSgeomTGeo::GetNDetectors(ll);
591       AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
592       globExtra[0]=point.GetX();
593       globExtra[1]=point.GetY();
594       globExtra[2]=point.GetZ();
595       AliITSgeomTGeo::GlobalToLocal(lay,lad,det,globExtra,locExtra);
596       //printf("%d %d %d %d %d  %f %f %f\n",volId,modId,lay,lad,det,locExtra[0],locExtra[1],locExtra[2]);
597       track->GetDZ(point.GetX(),point.GetY(),point.GetZ(),esd->GetMagneticField(),dzOverlap);
598       AliTrackPoint pointT;
599       Float_t radius,radiusT,phiv,phivT,thetav,thetavT;
600       for(Int_t lll=0;lll<ipt;lll++) {
601         array->GetPoint(pointT,lll);
602         Int_t layerIdT = AliGeomManager::VolUIDToLayer(pointT.GetVolumeID(),modId);
603         if(layerIdT!=layerId) continue;
604         radius=TMath::Sqrt((point.GetX()-vtxpos[0])*(point.GetX()-vtxpos[0])+(point.GetY()-vtxpos[1])*(point.GetY()-vtxpos[1]));
605         radiusT=TMath::Sqrt((pointT.GetX()-vtxpos[0])*(pointT.GetX()-vtxpos[0])+(pointT.GetY()-vtxpos[1])*(pointT.GetY()-vtxpos[1]));
606         phiv=TMath::ATan2(point.GetY()-vtxpos[1],point.GetX()-vtxpos[0]);
607         phivT=TMath::ATan2(pointT.GetY()-vtxpos[1],pointT.GetX()-vtxpos[0]);
608         if(TMath::Abs(point.GetZ()-vtxpos[2])<0.00001 || TMath::Abs(pointT.GetZ()-vtxpos[2])<0.00001) continue;
609         thetav=TMath::ATan(radius/(point.GetZ()-vtxpos[2]));
610         thetavT=TMath::ATan(radiusT/(pointT.GetZ()-vtxpos[2]));
611         dzOverlap[0]=(Float_t)((phivT-phiv)*0.5*(radiusT+radius));
612         if(TMath::Abs(TMath::Tan(0.5*(thetav+thetavT)))<0.00001) continue;
613         dzOverlap[1]=(Float_t)((pointT.GetZ()-point.GetZ())-(radiusT-radius)/TMath::Tan(0.5*(thetav+thetavT)));
614         fntExtra->Fill((Float_t)nclsTrk[itrack],(Float_t)(layerId-1),lad,volId,TMath::ATan2(point.GetY(),point.GetX()),point.GetX(),point.GetY(),point.GetZ(),locExtra[0],locExtra[2],dzOverlap[0],dzOverlap[1],d0z0mu[0],d0z0mu[1]);
615       }
616     }
617
618     if(!(GetRecoParam()->GetAlignFilterCosmicMergeTracks())) {
619       curv = curvArray[itrack];
620       curverr = curverrArray[itrack];
621       fspTree->Fill();
622     }
623   }
624
625   if(GetRecoParam()->GetAlignFilterCosmicMergeTracks()) {
626     curv = 0.5*(curvArray[0]+curvArray[1]);
627     curverr = 0.5*TMath::Sqrt(curverrArray[0]*curverrArray[0]+curverrArray[1]*curverrArray[1]);
628     fspTree->Fill();
629   }
630   PostData(0,fspTree);
631
632   if(!(GetRecoParam()->GetAlignFilterFillQANtuples())) return; 
633   // fill ntuple with track-to-track matching
634   Float_t phimu,thetamu,phiout,thetaout,dphi,dtheta,rotymu,rotyout,droty;    
635   Float_t d0[2],z0[2];
636   Float_t sigmad0[2],sigmaz0[2];
637   phimu = track1->GetAlpha()+TMath::ASin(track1->GetSnp());
638   thetamu = 0.5*TMath::Pi()-TMath::ATan(track1->GetTgl());
639   phiout = track2->GetAlpha()+TMath::ASin(track2->GetSnp());
640   thetaout = 0.5*TMath::Pi()-TMath::ATan(track2->GetTgl());
641   rotymu = TMath::ATan2(track1->Px(),track1->Pz());
642   rotyout = TMath::ATan2(track2->Px(),track2->Pz());
643
644   dphi = phimu - (phiout+TMath::Pi());
645   dtheta = thetamu - (TMath::Pi()-thetaout);
646   if(rotymu>0) {
647     droty = rotymu - (rotyout+TMath::Pi());
648   } else {
649     droty = rotymu - (rotyout-TMath::Pi());
650   }
651
652   Double_t alpha = TMath::ATan2(track1->Py(),track1->Px());
653
654   track1->Propagate(alpha,0.,esd->GetMagneticField());
655   track2->Propagate(alpha,0.,esd->GetMagneticField());
656   d0[0] = track1->GetY();
657   z0[0] = track1->GetZ();
658   d0[1] = track2->GetY();
659   z0[1] = track2->GetZ();
660   Float_t dxy = -(d0[0]-d0[1]);
661   Float_t dz  = z0[0]-z0[1];
662   sigmad0[0] = TMath::Sqrt(track1->GetSigmaY2());
663   sigmaz0[0] = TMath::Sqrt(track1->GetSigmaZ2());
664   sigmad0[1] = TMath::Sqrt(track2->GetSigmaY2());
665   sigmaz0[1] = TMath::Sqrt(track2->GetSigmaZ2());
666   /*  
667   Double_t xyz1atxl0[3],xyz1atxl1[3],xyz2atxl0[3],xyz2atxl1[3];
668   track1->GetXYZAt(0.,esd->GetMagneticField(),xyz1atxl0);
669   track1->GetXYZAt(1.,esd->GetMagneticField(),xyz1atxl1);
670   track2->GetXYZAt(0.,esd->GetMagneticField(),xyz2atxl0);
671   track2->GetXYZAt(1.,esd->GetMagneticField(),xyz2atxl1);
672   Float_t x1aty0 = (xyz1atxl0[0]*xyz1atxl1[1]-xyz1atxl0[1]*xyz1atxl1[0])/(xyz1atxl1[1]-xyz1atxl0[1]);
673   Float_t x2aty0 = (xyz2atxl0[0]*xyz2atxl1[1]-xyz2atxl0[1]*xyz2atxl1[0])/(xyz2atxl1[1]-xyz2atxl0[1]);
674   Float_t dxaty0 = x1aty0-x2aty0;
675   */
676   fntCosmicMatching->Fill((Float_t)nclsTrk[0],(Float_t)nclsTrk[1],track1->Pt(),track2->Pt(),sigmad0[0],sigmad0[1],sigmaz0[0],sigmaz0[1],dxy,dz,phimu,thetamu,TMath::Abs(d0z0mu[0]),d0z0mu[1]);
677   
678   return;
679 }
680
681 //________________________________________________________________________
682 void AliAlignmentDataFilterITS::FilterCollision(const AliESDEvent *esd)
683 {
684   // Extract ITS AliTrackPoints for Cosmics (check angular matching
685   // of top and bottom track, merge the two tracks, if requested)
686   //
687
688   // Set branch addresses for space points tree
689   AliTrackPointArray *arrayForTree=0;
690   Float_t curv,curverr;
691   fspTree->SetBranchAddress("SP",&arrayForTree);
692   fspTree->SetBranchAddress("curv",&curv);
693   fspTree->SetBranchAddress("curverr",&curverr);
694
695   Int_t ntracks = esd->GetNumberOfTracks();
696
697   if(ntracks==0) return;
698
699   if(esd->GetPrimaryVertexTracks()->GetNContributors()<=0) return;
700
701   Double_t vtxpos[3]; esd->GetPrimaryVertexTracks()->GetXYZ(vtxpos);
702
703   Int_t ncls=0;
704   Double_t pt=-10000.;
705   Double_t d0z0[2],covd0z0[3];
706   const AliTrackPointArray *array = 0;
707
708   for (Int_t itrack=0; itrack < ntracks; itrack++) {
709     AliESDtrack * track = esd->GetTrack(itrack);
710     if (!track) continue;
711
712     if(track->GetNcls(0)<GetRecoParam()->GetAlignFilterMinITSPoints()) continue;
713
714     if((GetRecoParam()->GetAlignFilterOnlyITSSATracks()) && track->GetNcls(1)>0) continue;
715     if((GetRecoParam()->GetAlignFilterOnlyITSTPCTracks()) && track->GetNcls(1)==0) continue;
716
717     if(track->Pt()<GetRecoParam()->GetAlignFilterMinPt() || 
718        track->Pt()>GetRecoParam()->GetAlignFilterMaxPt()) continue;
719
720     pt = track->Pt();
721     ncls = track->GetNcls(0);
722     Double_t maxd=10000.;
723     track->PropagateToDCA(esd->GetPrimaryVertex(),esd->GetMagneticField(),maxd,d0z0,covd0z0);
724
725     // read ITS cluster map
726     Int_t map[6];
727     for(Int_t ilay=0;ilay<6;ilay++) {
728       map[ilay]=0;
729       if(track->HasPointOnITSLayer(ilay)) map[ilay]=1;
730     }
731     AliDebug(2,Form("ITS map : %d %d %d %d %d %d pt %f\n",map[0],map[1],map[2],map[3],map[4],map[5],track->Pt()));
732     Int_t idx[12];
733     track->GetITSclusters(idx);
734     AliDebug(2,Form("cls idx %d %d %d %d %d %d %d %d %d %d %d %d\n",idx[0],idx[1],idx[2],idx[3],idx[4],idx[5],idx[6],idx[7],idx[8],idx[9],idx[10],idx[11]));
735     
736
737     AliTrackPoint point;
738     Int_t ipt,volId,modId,layerId,lay,lad,det;
739     Int_t jpt=0;
740     Bool_t layerOK[6]; for(Int_t l1=0;l1<6;l1++) layerOK[l1]=kFALSE;
741     
742     array = track->GetTrackPointArray();
743     if(!array) continue;
744     for(ipt=0; ipt<array->GetNPoints(); ipt++) {
745       array->GetPoint(point,ipt);
746       volId = point.GetVolumeID();
747       layerId = AliGeomManager::VolUIDToLayer(volId,modId);
748       if(layerId<1 || layerId>6) continue;
749       if(point.IsExtra() && 
750          (GetRecoParam()->GetAlignFilterSkipExtra())) continue;
751       layerOK[layerId-1]=kTRUE;
752       jpt++;
753     }
754
755     if(jpt < GetRecoParam()->GetAlignFilterMinITSPoints()) continue;
756
757     fHistNpoints->Fill(jpt);
758     fHistPt->Fill(pt);
759     PostData(1,fListOfHistos);
760
761     Float_t dzOverlap[2];
762     Double_t globExtra[3],locExtra[3];
763     arrayForTree = new AliTrackPointArray(jpt);
764     jpt=0;
765     array = track->GetTrackPointArray();
766     if(!array) continue;
767     for(ipt=0; ipt<array->GetNPoints(); ipt++) {
768       array->GetPoint(point,ipt);
769       volId = point.GetVolumeID();
770       layerId = AliGeomManager::VolUIDToLayer(volId,modId);
771       if(layerId<1 || layerId>6 || !layerOK[layerId-1]) continue;
772       if(!point.IsExtra() || 
773          !(GetRecoParam()->GetAlignFilterFillQANtuples())) continue;
774       ncls--;
775       for(Int_t ll=1;ll<layerId;ll++) modId+=AliITSgeomTGeo::GetNLadders(ll)*AliITSgeomTGeo::GetNDetectors(ll);
776       AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
777       globExtra[0]=point.GetX();
778       globExtra[1]=point.GetY();
779       globExtra[2]=point.GetZ();
780       AliITSgeomTGeo::GlobalToLocal(lay,lad,det,globExtra,locExtra);
781       track->GetDZ(point.GetX(),point.GetY(),point.GetZ(),esd->GetMagneticField(),dzOverlap);
782       AliTrackPoint pointT;
783       Float_t radius,radiusT,phiv,phivT,thetav,thetavT;
784       for(Int_t lll=0;lll<ipt;lll++) {
785         array->GetPoint(pointT,lll);
786         Int_t layerIdT = AliGeomManager::VolUIDToLayer(pointT.GetVolumeID(),modId);
787         if(layerIdT!=layerId) continue;
788         radius=TMath::Sqrt((point.GetX()-vtxpos[0])*(point.GetX()-vtxpos[0])+(point.GetY()-vtxpos[1])*(point.GetY()-vtxpos[1]));
789         radiusT=TMath::Sqrt((pointT.GetX()-vtxpos[0])*(pointT.GetX()-vtxpos[0])+(pointT.GetY()-vtxpos[1])*(pointT.GetY()-vtxpos[1]));
790         phiv=TMath::ATan2(point.GetY()-vtxpos[1],point.GetX()-vtxpos[0]);
791         phivT=TMath::ATan2(pointT.GetY()-vtxpos[1],pointT.GetX()-vtxpos[0]);
792         if(TMath::Abs(point.GetZ()-vtxpos[2])<0.00001 || TMath::Abs(pointT.GetZ()-vtxpos[2])<0.00001) continue;
793         thetav=TMath::ATan(radius/(point.GetZ()-vtxpos[2]));
794         thetavT=TMath::ATan(radiusT/(pointT.GetZ()-vtxpos[2]));
795         dzOverlap[0]=(Float_t)((phivT-phiv)*0.5*(radiusT+radius));
796         if(TMath::Abs(TMath::Tan(0.5*(thetav+thetavT)))<0.00001) continue;
797         dzOverlap[1]=(Float_t)((pointT.GetZ()-point.GetZ())-(radiusT-radius)/TMath::Tan(0.5*(thetav+thetavT)));
798         fntExtra->Fill((Float_t)ncls,(Float_t)(layerId-1),lad,volId,TMath::ATan2(point.GetY(),point.GetX()),point.GetX(),point.GetY(),point.GetZ(),locExtra[0],locExtra[2],dzOverlap[0],dzOverlap[1],d0z0[0],d0z0[1]);
799       }
800       arrayForTree->AddPoint(jpt,&point);
801       switch(layerId) {
802       case 1:
803         fHistLayer0->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
804         break;
805       case 2:
806         fHistLayer1->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
807         break;
808       case 3:
809         fHistLayer2->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
810         break;
811       case 4:
812         fHistLayer3->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
813         break;
814       case 5:
815         fHistLayer4->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
816         break;
817       case 6:
818         fHistLayer5->Fill(TMath::ATan2(point.GetY(),point.GetX()),point.GetZ());
819         break;
820       }
821       jpt++;
822     }
823
824     curv = track->GetC(esd->GetMagneticField());
825     curverr = TMath::Sqrt(track->GetSigma1Pt2())*track->GetC(esd->GetMagneticField())/track->OneOverPt();
826
827     fspTree->Fill();
828  
829   } // end of tracks loop
830
831   PostData(0,fspTree);
832
833   return;
834 }
835
836 //________________________________________________________________________
837 void AliAlignmentDataFilterITS::Terminate(Option_t */*option*/)
838 {
839   // Terminate analysis
840   //
841   AliDebug(2,"AliITSAlignmentDataFiler: Terminate() \n");
842
843   fspTree = dynamic_cast<TTree*> (GetOutputData(0));
844   if (!fspTree) {     
845     printf("ERROR: fspTree not available\n");
846     return;
847   }
848
849   fListOfHistos = dynamic_cast<TList*> (GetOutputData(1));
850   if (!fListOfHistos) {     
851     printf("ERROR: fListOfHistos not available\n");
852     return;
853   }
854
855   fHistNevents = dynamic_cast<TH1F*>(fListOfHistos->FindObject("fHistNevents"));
856   fHistNpoints = dynamic_cast<TH1F*>(fListOfHistos->FindObject("fHistNpoints"));
857   fHistPt = dynamic_cast<TH1F*>(fListOfHistos->FindObject("fHistPt"));
858   fHistLayer0 = dynamic_cast<TH2F*>(fListOfHistos->FindObject("fHistLayer0"));
859   fHistLayer1 = dynamic_cast<TH2F*>(fListOfHistos->FindObject("fHistLayer1"));
860   fHistLayer2 = dynamic_cast<TH2F*>(fListOfHistos->FindObject("fHistLayer2"));
861   fHistLayer3 = dynamic_cast<TH2F*>(fListOfHistos->FindObject("fHistLayer3"));
862   fHistLayer4 = dynamic_cast<TH2F*>(fListOfHistos->FindObject("fHistLayer4"));
863   fHistLayer5 = dynamic_cast<TH2F*>(fListOfHistos->FindObject("fHistLayer5"));
864   fntExtra = dynamic_cast<TNtuple*>(fListOfHistos->FindObject("fntExtra"));
865   fntCosmicMatching = dynamic_cast<TNtuple*>(fListOfHistos->FindObject("fntCosmicMatching"));
866
867
868
869   return;
870 }
871