]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSMultReconstructor.h
updates concerning noise reduction
[u/mrichter/AliRoot.git] / ITS / AliITSMultReconstructor.h
1 #ifndef ALIITSMULTRECONSTRUCTOR_H
2 #define ALIITSMULTRECONSTRUCTOR_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 //_________________________________________________________________________
7 // 
8 //        Implementation of the ITS-SPD trackleter class
9 //
10 // It retrieves clusters in the pixels (theta and phi) and finds tracklets.
11 // These can be used to extract charged particle multiplicity from the ITS.
12 //
13 // A tracklet consists of two ITS clusters, one in the first pixel layer and 
14 // one in the second. The clusters are associated if the differences in 
15 // Phi (azimuth) and Theta (polar angle) are within fiducial windows.
16 // In case of multiple candidates the candidate with minimum
17 // distance is selected. 
18 //_________________________________________________________________________
19 #include "AliTrackleter.h"
20 #include "AliITSsegmentationSPD.h"
21 #include "TMath.h"
22
23 class TBits;
24 class TTree;
25 class TH1F;
26 class TH2F; 
27 class AliITSDetTypeRec;
28 class AliITSgeom;
29 class AliESDEvent;
30 class AliESDtrack;
31 class AliVertex;
32 class AliESDVertex;
33 class AliMultiplicity;
34 class AliRefArray;
35 class AliITSRecPoint;
36
37 class AliITSMultReconstructor : public AliTrackleter
38 {
39 public:
40   //
41   enum {kClTh,kClPh,kClZ,kClMC0,kClMC1,kClMC2,kClNPar};
42   enum {kTrTheta,kTrPhi,kTrDPhi,kTrDTheta,kTrLab1,kTrLab2,kClID1,kClID2,kTrNPar};
43   enum {kSCTh,kSCPh,kSCLab,kSCID,kSCNPar};   
44   enum {kITSTPC,kITSSAP,kITSTPCBit=BIT(kITSTPC),kITSSAPBit=BIT(kITSSAP)}; // RS
45   AliITSMultReconstructor();
46   virtual ~AliITSMultReconstructor();
47
48   void Reconstruct(AliESDEvent* esd, TTree* treeRP);
49   void Reconstruct(TTree* tree, Float_t* vtx, Float_t* vtxRes=0);   // old reconstructor invocation
50   void ReconstructMix(TTree* clusterTree, TTree* clusterTreeMix, const Float_t* vtx, Float_t* vtrRes=0);
51   void FindTracklets(const Float_t* vtx); 
52   void LoadClusterFiredChips(TTree* tree);
53   void FlagClustersInOverlapRegions(Int_t ic1,Int_t ic2);
54   void FlagTrackClusters(Int_t id);
55   void FlagIfSecondary(AliESDtrack* track, const AliVertex* vtx);
56   void FlagV0s(const AliESDVertex *vtx);
57   void ProcessESDTracks();
58   Bool_t  CanBeElectron(const AliESDtrack* trc) const;
59   
60   virtual void CreateMultiplicityObject();
61   //
62   // Following members are set via AliITSRecoParam
63   void SetPhiWindow(Float_t w=0.08)    {fDPhiWindow=w;   fDPhiWindow2 = w*w;}
64   void SetThetaWindow(Float_t w=0.025) {fDThetaWindow=w; fDThetaWindow2=w*w;}
65   void SetPhiShift(Float_t w=0.0045) {fPhiShift=w;}
66   void SetRemoveClustersFromOverlaps(Bool_t b = kFALSE) {fRemoveClustersFromOverlaps = b;}
67   void SetPhiOverlapCut(Float_t w=0.005) {fPhiOverlapCut=w;}
68   void SetZetaOverlapCut(Float_t w=0.05) {fZetaOverlapCut=w;}
69   void SetPhiRotationAngle(Float_t w=0.0) {fPhiRotationAngle=w;}
70
71   Int_t GetNClustersLayer1() const {return fNClustersLay[0];}
72   Int_t GetNClustersLayer2() const {return fNClustersLay[1];}
73   Int_t GetNClustersLayer(Int_t i) const {return fNClustersLay[i];}
74   Int_t GetNTracklets() const {return fNTracklets;}
75   Int_t GetNSingleClusters() const {return fNSingleCluster;}
76   Int_t GetNSingleClustersLr(int lr) const {return lr==0 ? fNSingleCluster-fNSingleClusterSPD2:(GetStoreSPD2SingleCl() ? fNSingleClusterSPD2 : -1) ;}
77   Short_t GetNFiredChips(Int_t layer) const {return fNFiredChips[layer];}
78
79   Float_t* GetClusterLayer1(Int_t n) const {return &fClustersLay[0][n*kClNPar];}
80   Float_t* GetClusterLayer2(Int_t n) const {return &fClustersLay[1][n*kClNPar];}
81   Float_t* GetClusterOfLayer(Int_t lr,Int_t n) const {return &fClustersLay[lr][n*kClNPar];}
82   Int_t    GetClusterCopyIndex(Int_t lr,Int_t n) const {return fClusterCopyIndex[lr] ? fClusterCopyIndex[lr][n] : -1;}
83   AliITSRecPoint* GetRecPoint(Int_t lr, Int_t n) const;
84   
85   Float_t* GetTracklet(Int_t n) const {return fTracklets[n];}
86   Float_t* GetCluster(Int_t n) const {return fSClusters[n];}
87
88   void     SetScaleDThetaBySin2T(Bool_t v=kTRUE)         {fScaleDTBySin2T = v;}
89   Bool_t   GetScaleDThetaBySin2T()               const   {return fScaleDTBySin2T;}
90   //
91   void     SetNStdDev(Float_t f=1.)                      {fNStdDev = f<0.01 ? 0.01 : f; fNStdDevSq=TMath::Sqrt(fNStdDev);}
92   Float_t  GetNStdDev()                          const   {return fNStdDev;}
93   //
94   void SetHistOn(Bool_t b=kFALSE) {fHistOn=b;}
95   void SaveHists();
96   //
97   void   SetBuildRefs(Bool_t v=kTRUE)                    {fBuildRefs = v;}
98   Bool_t GetBuildRefs()                          const   {return fBuildRefs;}
99   //
100   void   SetStoreSPD2SingleCl(Bool_t v)                  {fStoreSPD2SingleCl = v;}
101   Bool_t GetStoreSPD2SingleCl()                  const   {return fStoreSPD2SingleCl;}
102   //
103   AliITSDetTypeRec *GetDetTypeRec() const {return fDetTypeRec;}
104   void SetDetTypeRec(AliITSDetTypeRec *ptr){fDetTypeRec = ptr;}
105   //
106   void    SetCutPxDrSPDin(Float_t v=0.1)             { fCutPxDrSPDin = v;}
107   void    SetCutPxDrSPDout(Float_t v=0.15)           { fCutPxDrSPDout = v;}
108   void    SetCutPxDz(Float_t v=0.2)                  { fCutPxDz = v;}
109   void    SetCutDCArz(Float_t v=0.5)                 { fCutDCArz = v;}
110   void    SetCutMinElectronProbTPC(Float_t v=0.5)    { fCutMinElectronProbTPC = v;}
111   void    SetCutMinElectronProbESD(Float_t v=0.1)    { fCutMinElectronProbESD = v;}
112   void    SetCutMinP(Float_t v=0.05)                 { fCutMinP = v;}
113   void    SetCutMinRGamma(Float_t v=2.)              { fCutMinRGamma = v;}
114   void    SetCutMinRK0(Float_t v=1.)                 { fCutMinRK0 = v;}
115   void    SetCutMinPointAngle(Float_t v=0.98)        { fCutMinPointAngle = v;}
116   void    SetCutMaxDCADauther(Float_t v=0.5)         { fCutMaxDCADauther = v;}
117   void    SetCutMassGamma(Float_t v=0.03)            { fCutMassGamma = v;}
118   void    SetCutMassGammaNSigma(Float_t v=5.)        { fCutMassGammaNSigma = v;}
119   void    SetCutMassK0(Float_t v=0.03)               { fCutMassK0 = v;}
120   void    SetCutMassK0NSigma(Float_t v=5.)           { fCutMassK0NSigma = v;}
121   void    SetCutChi2cGamma(Float_t v=2.)             { fCutChi2cGamma = v;}
122   void    SetCutChi2cK0(Float_t v=2.)                { fCutChi2cK0 = v;}
123   void    SetCutGammaSFromDecay(Float_t v=-10.)      { fCutGammaSFromDecay = v;}
124   void    SetCutK0SFromDecay(Float_t v=-10.)         { fCutK0SFromDecay = v;}
125   void    SetCutMaxDCA(Float_t v=1.)                 { fCutMaxDCA = v;}
126   //
127   Float_t GetCutPxDrSPDin()                    const {return fCutPxDrSPDin;}
128   Float_t GetCutPxDrSPDout()                   const {return fCutPxDrSPDout;}
129   Float_t GetCutPxDz()                         const {return fCutPxDz;}
130   Float_t GetCutDCArz()                        const {return fCutDCArz;}
131   Float_t GetCutMinElectronProbTPC()           const {return fCutMinElectronProbTPC;}
132   Float_t GetCutMinElectronProbESD()           const {return fCutMinElectronProbESD;}
133   Float_t GetCutMinP()                         const {return fCutMinP;}
134   Float_t GetCutMinRGamma()                    const {return fCutMinRGamma;}
135   Float_t GetCutMinRK0()                       const {return fCutMinRK0;}
136   Float_t GetCutMinPointAngle()                const {return fCutMinPointAngle;}
137   Float_t GetCutMaxDCADauther()                const {return fCutMaxDCADauther;}
138   Float_t GetCutMassGamma()                    const {return fCutMassGamma;}
139   Float_t GetCutMassGammaNSigma()              const {return fCutMassGammaNSigma;}
140   Float_t GetCutMassK0()                       const {return fCutMassK0;}
141   Float_t GetCutMassK0NSigma()                 const {return fCutMassK0NSigma;}
142   Float_t GetCutChi2cGamma()                   const {return fCutChi2cGamma;}
143   Float_t GetCutChi2cK0()                      const {return fCutChi2cK0;}
144   Float_t GetCutGammaSFromDecay()              const {return fCutGammaSFromDecay;}
145   Float_t GetCutK0SFromDecay()                 const {return fCutK0SFromDecay;}
146   Float_t GetCutMaxDCA()                       const {return fCutMaxDCA;}
147   //
148   void  InitAux();
149   void  ClusterPos2Angles(const Float_t *vtx);
150   void  ClusterPos2Angles(Float_t *clPar, const Float_t *vtx) const;
151   Int_t AssociateClusterOfL1(Int_t iC1);
152   Int_t StoreTrackletForL2Cluster(Int_t iC2);
153   void  StoreL1Singles();
154   TClonesArray* GetClustersOfLayer(Int_t il)   const {return fClArr[il];}
155   void  LoadClusters()                               {LoadClusterArrays(fTreeRP);}
156   void  SetTreeRP(TTree* rp)                         {fTreeRP    = rp;}
157   void  SetTreeRPMix(TTree* rp=0)                    {fTreeRPMix = rp;}
158   Bool_t AreClustersLoaded()                   const {return fClustersLoaded;}
159   Bool_t GetCreateClustersCopy()               const {return fCreateClustersCopy;}
160   Bool_t IsRecoDone()                          const {return fRecoDone;}
161   void  SetCreateClustersCopy(Bool_t v=kTRUE)        {fCreateClustersCopy=v;}
162   //
163   //  Float_t* GetClustersArray(Int_t lr)          const {return (Float_t*) (lr==0) ? fClustersLay[0]:fClustersLay[1];}
164   Float_t* GetClustersArray(Int_t lr)          const {if(lr==0){return fClustersLay[0];} 
165                                                       else {return fClustersLay[1];}}
166   Int_t*   GetPartnersOfL2()                   const {return (Int_t*)fPartners;}
167   Float_t* GetMinDistsOfL2()                   const {return (Float_t*)fMinDists;}
168   Double_t GetDPhiShift()                      const {return fDPhiShift;}
169   Double_t GetDPhiWindow2()                    const {return fDPhiWindow2;}
170   Double_t GetDThetaWindow2()                  const {return fDThetaWindow2;}
171   Double_t CalcDist(Double_t dphi, Double_t dtheta, Double_t theta) const; 
172   //
173  protected:
174   void   SetClustersLoaded(Bool_t v=kTRUE)           {fClustersLoaded = v;}
175   AliITSMultReconstructor(const AliITSMultReconstructor& mr);
176   AliITSMultReconstructor& operator=(const AliITSMultReconstructor& mr);
177   void              CalcThetaPhi(float dx,float dy,float dz,float &theta,float &phi) const;
178   AliITSDetTypeRec* fDetTypeRec;            //! pointer to DetTypeRec
179   AliESDEvent*      fESDEvent;              //! pointer to ESD event
180   TTree*            fTreeRP;                //! ITS recpoints
181   TTree*            fTreeRPMix;             //! ITS recpoints for mixing
182   AliRefArray*      fUsedClusLay[2][2];     //! RS: clusters usage in ESD tracks
183   //
184   Float_t*      fClustersLay[2];            //! clusters in the SPD layers of ITS 
185   Int_t*        fDetectorIndexClustersLay[2];  //! module index for clusters in ITS layers
186   Int_t*        fClusterCopyIndex[2];         //! when clusters copy is requested, store here the reference on the index
187   Bool_t*       fOverlapFlagClustersLay[2];  //! flag for clusters in the overlap regions in ITS layers
188
189   Float_t**     fTracklets;            //! tracklets 
190   Float_t**     fSClusters;            //! single clusters (unassociated)
191   
192   Int_t         fNClustersLay[2];      // Number of clusters on each layer
193   Int_t         fNTracklets;           // Number of tracklets
194   Int_t         fNSingleCluster;       // Number of unassociated clusters
195   Int_t         fNSingleClusterSPD2;   // Number of unassociated clusters on 2nd lr
196   Short_t       fNFiredChips[2];       // Number of fired chips in the two SPD layers
197   //
198   // Following members are set via AliITSRecoParam
199   //
200   Float_t       fDPhiWindow;                   // Search window in phi
201   Float_t       fDThetaWindow;                 // Search window in theta
202   Float_t       fPhiShift;                     // Phi shift reference value (at 0.5 T) 
203   Bool_t        fRemoveClustersFromOverlaps;   // Option to skip clusters in the overlaps
204   Float_t       fPhiOverlapCut;                // Fiducial window in phi for overlap cut
205   Float_t       fZetaOverlapCut;               // Fiducial window in eta for overlap cut
206   Float_t       fPhiRotationAngle;             // Angle to rotate the inner layer cluster for combinatorial reco only 
207   //
208   Bool_t        fScaleDTBySin2T;               // use in distance definition
209   Float_t       fNStdDev;                      // number of standard deviations to keep
210   Float_t       fNStdDevSq;                    // sqrt of number of standard deviations to keep
211   //
212   // cuts for secondaries identification
213   Float_t       fCutPxDrSPDin;                 // max P*DR for primaries involving at least 1 SPD
214   Float_t       fCutPxDrSPDout;                // max P*DR for primaries not involving any SPD
215   Float_t       fCutPxDz;                      // max P*DZ for primaries
216   Float_t       fCutDCArz;                     // max DR or DZ for primares
217   //
218   // cuts for flagging tracks in V0s
219   Float_t       fCutMinElectronProbTPC;     // min probability for e+/e- PID involving TPC
220   Float_t       fCutMinElectronProbESD;     // min probability for e+/e- PID not involving TPC
221   //
222   Float_t       fCutMinP;                   // min P of V0
223   Float_t       fCutMinRGamma;              // min transv. distance from ESDVertex to V0 for gammas
224   Float_t       fCutMinRK0;                 // min transv. distance from ESDVertex to V0 for K0s
225   Float_t       fCutMinPointAngle;          // min pointing angle cosine
226   Float_t       fCutMaxDCADauther;          // max DCA of daughters at V0
227   Float_t       fCutMassGamma;              // max gamma mass
228   Float_t       fCutMassGammaNSigma;        // max standard deviations from 0 for gamma
229   Float_t       fCutMassK0;                 // max K0 mass difference from PGD value
230   Float_t       fCutMassK0NSigma;           // max standard deviations for K0 mass from PDG value
231   Float_t       fCutChi2cGamma;             // max constrained chi2 cut for gammas
232   Float_t       fCutChi2cK0;                // max constrained chi2 cut for K0s
233   Float_t       fCutGammaSFromDecay;        // min path*P for gammas
234   Float_t       fCutK0SFromDecay;           // min path*P for K0s
235   Float_t       fCutMaxDCA;                 // max DCA for V0 at ESD vertex  
236
237   Bool_t        fHistOn;               // Option to define and fill the histograms 
238
239   TH1F*         fhClustersDPhiAcc;     // Phi2 - Phi1 for tracklets 
240   TH1F*         fhClustersDThetaAcc;   // Theta2 - Theta1 for tracklets 
241   TH1F*         fhClustersDPhiAll;     // Phi2 - Phi1 all the combinations 
242   TH1F*         fhClustersDThetaAll;   // Theta2 - Theta1 all the combinations
243  
244   TH2F*         fhDPhiVsDThetaAll;     // 2D plot for all the combinations  
245   TH2F*         fhDPhiVsDThetaAcc;     // same plot for tracklets 
246
247   TH1F*         fhetaTracklets;        // Pseudorapidity distr. for tracklets 
248   TH1F*         fhphiTracklets;        // Azimuthal (Phi) distr. for tracklets  
249   TH1F*         fhetaClustersLay1;     // Pseudorapidity distr. for Clusters L. 1
250   TH1F*         fhphiClustersLay1;     // Azimuthal (Phi) distr. for Clusters L. 1 
251
252   // temporary stuff for single event trackleting
253   Double_t      fDPhiShift;            // shift in dphi due to the curvature
254   Double_t      fDPhiWindow2;          // phi window^2
255   Double_t      fDThetaWindow2;        // theta window^2
256   Int_t*        fPartners;             //! L2 partners of L1
257   Int_t*        fAssociatedLay1;       //! association flag
258   Float_t*      fMinDists;             //! smallest distances for L2->L1
259   AliRefArray*  fBlackList;            //! blacklisted cluster references
260   Bool_t        fStoreRefs[2][2];      //! which cluster to track refs to store
261   //
262   // this is for the analysis mode only
263   TClonesArray  *fClArr[2];            //! original clusters
264   Bool_t        fCreateClustersCopy;   //  read and clone clusters directly from the tree
265   Bool_t        fClustersLoaded;       // flag of clusters loaded
266   Bool_t        fRecoDone;             // flag that reconstruction is done
267   Bool_t        fBuildRefs;            // build cluster to tracks references
268   Bool_t        fStoreSPD2SingleCl;    // do we store SPD2 singles
269   //
270   AliITSsegmentationSPD fSPDSeg;       // SPD segmentation model
271   //
272   void LoadClusterArrays(TTree* tree, TTree* treeMix=0);
273   void LoadClusterArrays(TTree* tree,int il);
274
275   ClassDef(AliITSMultReconstructor,11)
276 };
277
278 //____________________________________________________________________
279 inline void AliITSMultReconstructor::ClusterPos2Angles(Float_t *clPar, const Float_t *vtx) const
280 {
281   // convert cluster coordinates to angles wrt vertex
282   Float_t x = clPar[kClTh] - vtx[0];
283   Float_t y = clPar[kClPh] - vtx[1];
284   Float_t z = clPar[kClZ]  - vtx[2];
285   Float_t r    = TMath::Sqrt(x*x + y*y + z*z);
286   clPar[kClTh] = TMath::ACos(z/r);                   // Store Theta
287   clPar[kClPh] = TMath::Pi() + TMath::ATan2(-y,-x);  // Store Phi 
288   //
289 }
290
291 //____________________________________________________________________
292 inline Double_t AliITSMultReconstructor::CalcDist(Double_t dphi, Double_t dtheta, Double_t theta) const
293 {
294   // calculate eliptical distance. theta is the angle of cl1, dtheta = tht(cl1)-tht(cl2)
295   dphi = TMath::Abs(dphi) - fDPhiShift;
296   if (fScaleDTBySin2T) {
297     double sinTI = TMath::Sin(theta-dtheta/2);
298     sinTI *= sinTI;
299     dtheta /= sinTI>1.e-6 ? sinTI : 1.e-6;
300   }
301   return dphi*dphi/fDPhiWindow2 + dtheta*dtheta/fDThetaWindow2;
302 }
303
304 //____________________________________________________________________
305 inline void AliITSMultReconstructor::CalcThetaPhi(float x, float y,float z,float &theta,float &phi) const
306 {
307   // get theta and phi in tracklet convention
308   theta = TMath::ACos(z/TMath::Sqrt(x*x + y*y + z*z));
309   phi   = TMath::Pi() + TMath::ATan2(-y,-x);
310 }
311
312
313 #endif