]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSMultReconstructor.cxx
modifications to be compliant with -Weffc++ option (L. Gaudichet)
[u/mrichter/AliRoot.git] / ITS / AliITSMultReconstructor.cxx
1 //____________________________________________________________________
2 // 
3 // AliITSMultReconstructor - find clusters in the pixels (theta and
4 // phi) and tracklets.
5 // 
6 // These can be used to extract charged particles multiplcicity from the ITS.
7 //
8 // A tracklet consist of two ITS clusters, one in the first pixel
9 // layer and one in the second. The clusters are associates if the 
10 // differencies in Phi (azimuth) and Zeta (longitudinal) are inside 
11 // a fiducial volume. In case of multiple candidates it is selected the
12 // candidate with minimum distance in Phi. 
13 // The parameter AssociationChoice allows to control if two clusters 
14 // in layer 2 can be associated to the same cluster in layer 1 or not.
15 //
16 // -----------------------------------------------------------------
17 // 
18 // NOTE: The cuts on phi and zeta depends on the interacting system (p-p  
19 //  or Pb-Pb). Please, check the file AliITSMultReconstructor.h and be 
20 //  sure that SetPhiWindow and SetZetaWindow are defined accordingly.
21 // 
22 //  
23 //  
24 //
25 //____________________________________________________________________
26
27 #include "AliITSMultReconstructor.h"
28
29 #include "TTree.h"
30 #include "TH1F.h"
31 #include "TH2F.h"
32
33 #include "AliITSRecPoint.h"
34 #include "AliITSgeom.h"
35 #include "AliLog.h"
36
37 //____________________________________________________________________
38 ClassImp(AliITSMultReconstructor)
39
40
41 //____________________________________________________________________
42 AliITSMultReconstructor::AliITSMultReconstructor() {
43   // Method to reconstruct the charged particles multiplicity with the 
44   // SPD (tracklets).
45
46   fGeometry =0;
47
48   SetHistOn();
49   SetPhiWindow();
50   SetZetaWindow();
51   SetOnlyOneTrackletPerC2();
52
53   fClustersLay1       = new Float_t*[300000];
54   fClustersLay2       = new Float_t*[300000];
55   fTracklets          = new Float_t*[300000];
56   fAssociationFlag    = new Bool_t[300000];
57
58   for(Int_t i=0; i<300000; i++) {
59     fClustersLay1[i]       = new Float_t[3];
60     fClustersLay2[i]       = new Float_t[3];
61     fTracklets[i]          = new Float_t[3];
62     fAssociationFlag[i]    = kFALSE;
63   }
64
65   // definition of histograms
66   fhClustersDPhiAcc   = new TH1F("dphiacc",  "dphi",  100,-0.1,0.1);
67   fhClustersDPhiAcc->SetDirectory(0);
68   fhClustersDThetaAcc = new TH1F("dthetaacc","dtheta",100,-0.1,0.1);
69   fhClustersDThetaAcc->SetDirectory(0);
70   fhClustersDZetaAcc = new TH1F("dzetaacc","dzeta",100,-1.,1.);
71   fhClustersDZetaAcc->SetDirectory(0);
72
73   fhDPhiVsDZetaAcc = new TH2F("dphiVsDzetaacc","",100,-1.,1.,100,-0.1,0.1);
74   fhDPhiVsDZetaAcc->SetDirectory(0);
75   fhDPhiVsDThetaAcc = new TH2F("dphiVsDthetaAcc","",100,-0.1,0.1,100,-0.1,0.1);
76   fhDPhiVsDThetaAcc->SetDirectory(0);
77
78   fhClustersDPhiAll   = new TH1F("dphiall",  "dphi",  100,-0.5,0.5);
79   fhClustersDPhiAll->SetDirectory(0);
80   fhClustersDThetaAll = new TH1F("dthetaall","dtheta",100,-0.5,0.5);
81   fhClustersDThetaAll->SetDirectory(0);
82   fhClustersDZetaAll = new TH1F("dzetaall","dzeta",100,-5.,5.);
83   fhClustersDZetaAll->SetDirectory(0);
84
85   fhDPhiVsDZetaAll = new TH2F("dphiVsDzetaall","",100,-5.,5.,100,-0.5,0.5);
86   fhDPhiVsDZetaAll->SetDirectory(0);
87   fhDPhiVsDThetaAll = new TH2F("dphiVsDthetaAll","",100,-0.5,0.5,100,-0.5,0.5);
88   fhDPhiVsDThetaAll->SetDirectory(0);
89
90   fhetaTracklets  = new TH1F("etaTracklets",  "eta",  100,-2.,2.);
91   fhphiTracklets  = new TH1F("phiTracklets",  "phi",  100,-3.14159,3.14159);
92   fhetaClustersLay1  = new TH1F("etaClustersLay1",  "etaCl1",  100,-2.,2.);
93   fhphiClustersLay1  = new TH1F("phiClustersLay1", "phiCl1", 100,-3.141,3.141);
94
95 }
96
97 //______________________________________________________________________
98 AliITSMultReconstructor::AliITSMultReconstructor(const AliITSMultReconstructor &mr) : TObject(mr) {
99   // Copy constructor
100   // Copies are not allowed. The method is protected to avoid misuse.
101   Error("AliITSMultReconstructor","Copy constructor not allowed\n");
102 }
103
104 //______________________________________________________________________
105 AliITSMultReconstructor& AliITSMultReconstructor::operator=(const AliITSMultReconstructor& /* mr */){
106   // Assignment operator
107   // Assignment is not allowed. The method is protected to avoid misuse.
108   Error("= operator","Assignment operator not allowed\n");
109   return *this;
110 }
111
112 //______________________________________________________________________
113 AliITSMultReconstructor::~AliITSMultReconstructor(){
114   // Destructor
115   if(fhClustersDPhiAcc)delete fhClustersDPhiAcc;
116   if(fhClustersDThetaAcc)delete fhClustersDThetaAcc;
117   if(fhClustersDZetaAcc)delete fhClustersDZetaAcc;
118   if(fhClustersDPhiAll)delete fhClustersDPhiAll;
119   if(fhClustersDThetaAll)delete fhClustersDThetaAll;
120   if(fhClustersDZetaAll)delete fhClustersDZetaAll;
121   if(fhDPhiVsDThetaAll)delete fhDPhiVsDThetaAll;
122   if(fhDPhiVsDThetaAcc)delete fhDPhiVsDThetaAcc;
123   if(fhDPhiVsDZetaAll)delete fhDPhiVsDZetaAll;
124   if(fhDPhiVsDZetaAcc)delete fhDPhiVsDZetaAcc;
125   if(fhetaTracklets)delete fhetaTracklets;
126   if(fhphiTracklets)delete fhphiTracklets;
127   if(fhetaClustersLay1)delete fhetaClustersLay1;
128   if(fhphiClustersLay1)delete fhphiClustersLay1;
129   if(fClustersLay1){
130     for(Int_t i=0; i<300000; i++) {
131       delete [] fClustersLay1[i];
132       delete [] fClustersLay2[i];
133       delete [] fTracklets[i];
134     }
135     delete fClustersLay1;
136     delete  fClustersLay2;
137     delete fTracklets;
138   }
139   if(fAssociationFlag)delete fAssociationFlag;
140 }
141
142 //____________________________________________________________________
143 void
144 AliITSMultReconstructor::Reconstruct(TTree* clusterTree, Float_t* vtx, Float_t* /* vtxRes*/) {
145   //
146   // - calls LoadClusterArray that finds the position of the clusters
147   //   (in global coord) 
148   // - convert the cluster coordinates to theta, phi (seen from the
149   //   interaction vertex). The third coordinate is used for ....
150   // - makes an array of tracklets 
151   //   
152   // After this method has been called, the clusters of the two layers
153   // and the tracklets can be retrieved by calling the Get'er methods.
154
155   // reset counters
156   fNClustersLay1 = 0;
157   fNClustersLay2 = 0;
158   fNTracklets = 0; 
159
160   // loading the clusters 
161   LoadClusterArrays(clusterTree);
162
163   // find the tracklets
164   AliDebug(1,"Looking for tracklets... ");  
165
166   //###########################################################
167   // Loop on layer 1 : finding theta, phi and z 
168   for (Int_t iC1=0; iC1<fNClustersLay1; iC1++) {    
169     Float_t x = fClustersLay1[iC1][0] - vtx[0];
170     Float_t y = fClustersLay1[iC1][1] - vtx[1];
171     Float_t z = fClustersLay1[iC1][2] - vtx[2];
172
173     Float_t r    = TMath::Sqrt(TMath::Power(x,2) +
174                                TMath::Power(y,2) +
175                                TMath::Power(z,2));
176     
177     fClustersLay1[iC1][0] = TMath::ACos(z/r);  // Store Theta
178     fClustersLay1[iC1][1] = TMath::ATan2(x,y);  // Store Phi
179     fClustersLay1[iC1][2] = z/r;               // Store scaled z 
180     if (fHistOn) {
181       Float_t eta=fClustersLay1[iC1][0];
182       eta= TMath::Tan(eta/2.);
183       eta=-TMath::Log(eta);
184       fhetaClustersLay1->Fill(eta);    
185       fhphiClustersLay1->Fill(fClustersLay1[iC1][1]);    
186     }      
187 }
188   
189   // Loop on layer 2 : finding theta, phi and r   
190   for (Int_t iC2=0; iC2<fNClustersLay2; iC2++) {    
191     Float_t x = fClustersLay2[iC2][0] - vtx[0];
192     Float_t y = fClustersLay2[iC2][1] - vtx[1];
193     Float_t z = fClustersLay2[iC2][2] - vtx[2];
194    
195     Float_t r    = TMath::Sqrt(TMath::Power(x,2) +
196                                TMath::Power(y,2) +
197                                TMath::Power(z,2));
198     
199     fClustersLay2[iC2][0] = TMath::ACos(z/r);  // Store Theta
200     fClustersLay2[iC2][1] = TMath::ATan2(x,y);  // Store Phi
201     fClustersLay2[iC2][2] = z;                 // Store z
202
203  // this only needs to be initialized for the fNClustersLay2 first associations
204     fAssociationFlag[iC2] = kFALSE;
205   }  
206   
207   //###########################################################
208   // Loop on layer 1 
209   for (Int_t iC1=0; iC1<fNClustersLay1; iC1++) {    
210
211     // reset of variables for multiple candidates
212     Int_t  iC2WithBestDist = 0;     // reset 
213     Float_t distmin        = 100.;  // just to put a huge number! 
214     Float_t dPhimin        = 0.;  // Used for histograms only! 
215     Float_t dThetamin      = 0.;  // Used for histograms only! 
216     Float_t dZetamin       = 0.;  // Used for histograms only! 
217     
218     // Loop on layer 2 
219     for (Int_t iC2=0; iC2<fNClustersLay2; iC2++) {      
220       
221       // The following excludes double associations
222       if (!fAssociationFlag[iC2]) {
223         
224         // find the difference in angles
225         Float_t dTheta = fClustersLay2[iC2][0] - fClustersLay1[iC1][0];
226         Float_t dPhi   = fClustersLay2[iC2][1] - fClustersLay1[iC1][1];
227         
228         // find the difference in z (between linear projection from layer 1
229         // and the actual point: Dzeta= z1/r1*r2 -z2)   
230         Float_t r2   = fClustersLay2[iC2][2]/TMath::Cos(fClustersLay2[iC2][0]);
231         Float_t dZeta  = fClustersLay1[iC1][2]*r2 - fClustersLay2[iC2][2]; 
232
233         if (fHistOn) {
234           fhClustersDPhiAll->Fill(dPhi);    
235           fhClustersDThetaAll->Fill(dTheta);    
236           fhClustersDZetaAll->Fill(dZeta);    
237           fhDPhiVsDThetaAll->Fill(dTheta, dPhi);
238           fhDPhiVsDZetaAll->Fill(dZeta, dPhi);
239         }
240         // make "elliptical" cut in Phi and Zeta! 
241         Float_t d = TMath::Sqrt(TMath::Power(dPhi/fPhiWindow,2) + TMath::Power(dZeta/fZetaWindow,2));
242
243         if (d>1) continue;      
244         
245         //look for the minimum distance: the minimum is in iC2WithBestDist
246         if (TMath::Sqrt(dZeta*dZeta+(r2*dPhi*r2*dPhi)) < distmin ) {
247           distmin=TMath::Sqrt(dZeta*dZeta + (r2*dPhi*r2*dPhi));
248           dPhimin = dPhi;
249           dThetamin = dTheta;
250           dZetamin = dZeta; 
251           iC2WithBestDist = iC2;
252         }
253       } 
254     } // end of loop over clusters in layer 2 
255     
256     if (distmin<100) { // This means that a cluster in layer 2 was found that mathes with iC1
257
258       if (fHistOn) {
259         fhClustersDPhiAcc->Fill(dPhimin);    
260         fhClustersDThetaAcc->Fill(dThetamin);    
261         fhClustersDZetaAcc->Fill(dZetamin);    
262         fhDPhiVsDThetaAcc->Fill(dThetamin, dPhimin);
263         fhDPhiVsDZetaAcc->Fill(dZetamin, dPhimin);
264       }
265       
266       if (fOnlyOneTrackletPerC2) fAssociationFlag[iC2WithBestDist] = kTRUE; // flag the association
267       
268       // store the tracklet
269       
270       // use the theta from the clusters in the first layer 
271       fTracklets[fNTracklets][0] = fClustersLay1[iC1][0];
272       // use the phi from the clusters in the first layer 
273       fTracklets[fNTracklets][1] = fClustersLay1[iC1][1];
274       // Store the difference between phi1 and phi2
275       fTracklets[fNTracklets][2] = fClustersLay1[iC1][1] - fClustersLay2[iC2WithBestDist][1];       
276   
277       if (fHistOn) {
278         Float_t eta=fTracklets[fNTracklets][0];
279         eta= TMath::Tan(eta/2.);
280         eta=-TMath::Log(eta);
281         fhetaTracklets->Fill(eta);    
282         fhphiTracklets->Fill(fTracklets[fNTracklets][1]);    
283       }
284       
285       AliDebug(1,Form(" Adding tracklet candidate %d ", fNTracklets));
286       AliDebug(1,Form(" Cl. %d of Layer 1 and %d of Layer 2", iC1, 
287                       iC2WithBestDist));
288       fNTracklets++;
289     }
290
291     // Delete the following else if you do not want to save Clusters! 
292
293     else { // This means that the cluster has not been associated 
294
295       // store the cluster
296        
297       fTracklets[fNTracklets][0] = fClustersLay1[iC1][0];
298       fTracklets[fNTracklets][1] = fClustersLay1[iC1][1];
299       // Store a flag. This will indicate that the "tracklet" 
300       // was indeed a single cluster! 
301       fTracklets[fNTracklets][2] = -999999.;       
302       AliDebug(1,Form(" Adding a single cluster %d (cluster %d  of layer 1)", 
303                       fNTracklets, iC1));
304       fNTracklets++;
305     }
306
307   } // end of loop over clusters in layer 1
308   
309   AliDebug(1,Form("%d tracklets found", fNTracklets));
310 }
311
312 //____________________________________________________________________
313 void
314 AliITSMultReconstructor::LoadClusterArrays(TTree* itsClusterTree) {
315   // This method
316   // - gets the clusters from the cluster tree 
317   // - convert them into global coordinates 
318   // - store them in the internal arrays
319   
320   AliDebug(1,"Loading clusters ...");
321   
322   fNClustersLay1 = 0;
323   fNClustersLay2 = 0;
324   
325   TClonesArray* itsClusters = new TClonesArray("AliITSRecPoint");
326   TBranch* itsClusterBranch=itsClusterTree->GetBranch("ITSRecPoints");
327
328   itsClusterBranch->SetAddress(&itsClusters);
329
330   Int_t nItsSubs = (Int_t)itsClusterTree->GetEntries();  
331  
332   // loop over the its subdetectors
333   for (Int_t iIts=0; iIts < nItsSubs; iIts++) {
334     
335     if (!itsClusterTree->GetEvent(iIts)) 
336       continue;
337     
338     Int_t nClusters = itsClusters->GetEntriesFast();
339     
340     // stuff needed to get the global coordinates
341     Double_t rot[9];   fGeometry->GetRotMatrix(iIts,rot);
342     Int_t lay,lad,det; fGeometry->GetModuleId(iIts,lay,lad,det);
343     Float_t tx,ty,tz;  fGeometry->GetTrans(lay,lad,det,tx,ty,tz);
344     
345     // Below:
346     // "alpha" is the angle from the global X-axis to the
347     //         local GEANT X'-axis  ( rot[0]=cos(alpha) and rot[1]=sin(alpha) )
348     // "phi" is the angle from the global X-axis to the
349     //       local cluster X"-axis
350     
351     Double_t alpha   = TMath::ATan2(rot[1],rot[0])+TMath::Pi();
352     Double_t itsPhi = TMath::Pi()/2+alpha;
353     
354     if (lay==1) itsPhi+=TMath::Pi();
355     Double_t cp=TMath::Cos(itsPhi), sp=TMath::Sin(itsPhi);
356     Double_t r=tx*cp+ty*sp;
357     
358     // loop over clusters
359     while(nClusters--) {
360       AliITSRecPoint* cluster = (AliITSRecPoint*)itsClusters->UncheckedAt(nClusters);   
361       
362       if (cluster->GetLayer()>1) 
363         continue;            
364       
365       Float_t x = r*cp - cluster->GetY()*sp;
366       Float_t y = r*sp + cluster->GetY()*cp;
367       Float_t z = cluster->GetZ();      
368       
369       if (cluster->GetLayer()==0) {
370         fClustersLay1[fNClustersLay1][0] = x;
371         fClustersLay1[fNClustersLay1][1] = y;
372         fClustersLay1[fNClustersLay1][2] = z;
373         fNClustersLay1++;
374       }
375       if (cluster->GetLayer()==1) {     
376         fClustersLay2[fNClustersLay2][0] = x;
377         fClustersLay2[fNClustersLay2][1] = y;
378         fClustersLay2[fNClustersLay2][2] = z;
379         fNClustersLay2++;
380       }
381       
382     }// end of cluster loop
383   } // end of its "subdetector" loop  
384   
385   AliDebug(1,Form("(clusters in layer 1 : %d,  layer 2: %d)",fNClustersLay1,fNClustersLay2));
386 }
387 //____________________________________________________________________
388 void
389 AliITSMultReconstructor::SaveHists() {
390   // This method save the histograms on the output file
391   // (only if fHistOn is TRUE). 
392   
393   if (!fHistOn)
394     return;
395
396   fhClustersDPhiAll->Write();
397   fhClustersDThetaAll->Write();
398   fhClustersDZetaAll->Write();
399   fhDPhiVsDThetaAll->Write();
400   fhDPhiVsDZetaAll->Write();
401
402   fhClustersDPhiAcc->Write();
403   fhClustersDThetaAcc->Write();
404   fhClustersDZetaAcc->Write();
405   fhDPhiVsDThetaAcc->Write();
406   fhDPhiVsDZetaAcc->Write();
407
408   fhetaTracklets->Write();
409   fhphiTracklets->Write();
410   fhetaClustersLay1->Write();
411   fhphiClustersLay1->Write();
412 }