]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/AliCaloTrackESDReader.cxx
Configuration example for correlations of gamma and jets reconstructed around leading
[u/mrichter/AliRoot.git] / PWG4 / AliCaloTrackESDReader.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16 /* $Id:  $ */
17
18 /* History of cvs commits:
19  *
20  * $Log$
21  *
22  *
23  */
24
25 //_________________________________________________________________________
26 // Class for reading data (ESDs) in order to do prompt gamma 
27 // or other particle identification and correlations
28 //
29 //
30 //*-- Author: Gustavo Conesa (LNF-INFN) 
31 //////////////////////////////////////////////////////////////////////////////
32
33
34 // --- ROOT system ---
35
36 //---- ANALYSIS system ----
37 #include "AliCaloTrackESDReader.h" 
38 #include "AliESDEvent.h"
39 #include "AliESDVertex.h"
40 #include "AliESDCaloCluster.h"
41 #include "AliAODCaloCluster.h"
42 #include "AliAODTrack.h"
43 #include "AliAODEvent.h"
44 #include "AliMCEvent.h"
45 #include "AliLog.h"
46 #include "Riostream.h"
47
48 ClassImp(AliCaloTrackESDReader)
49
50 //____________________________________________________________________________
51 AliCaloTrackESDReader::AliCaloTrackESDReader() : 
52   AliCaloTrackReader()
53 {
54   //Default Ctor
55   
56   //Initialize parameters
57   fDataType=kESD;
58   
59 }
60
61 //____________________________________________________________________________
62 AliCaloTrackESDReader::AliCaloTrackESDReader(const AliCaloTrackESDReader & g) :   
63   AliCaloTrackReader(g)
64 {
65   // cpy ctor
66 }
67
68 //_________________________________________________________________________
69 AliCaloTrackESDReader & AliCaloTrackESDReader::operator = (const AliCaloTrackESDReader & source)
70 {
71   // assignment operator
72
73   if(&source == this) return *this;
74
75   return *this;
76
77 }
78
79 //____________________________________________________________________________
80 void AliCaloTrackESDReader::FillInputCTS() {
81   //Return array with CTS tracks
82
83   fAODCTS = new TClonesArray("AliAODTrack",0);
84
85   Int_t nTracks   = fESD->GetNumberOfTracks() ;
86   Int_t naod = 0;
87   Double_t pos[3];
88   Double_t p[3];
89   Double_t covTr[21];
90   Double_t pid[10];
91   for (Int_t itrack =  0; itrack <  nTracks; itrack++) {////////////// track loop
92     AliESDtrack * track = fESD->GetTrack(itrack) ; // retrieve track from esd
93
94     //We want tracks fitted in the detectors:
95     ULong_t status=AliESDtrack::kTPCrefit;
96     status|=AliESDtrack::kITSrefit;
97     
98     if ( (track->GetStatus() & status) == status) {//Check if the bits we want are set
99     
100       track->GetPxPyPz(p) ;
101       TLorentzVector momentum(p[0],p[1],p[2],0);
102    
103       if(fCTSPtMin < momentum.Pt() &&fFidutialCut->IsInFidutialCut(momentum,"CTS")){
104
105         if(fDebug > 3 && momentum.Pt() > 0.2)printf("FillInputCTS():: Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
106                                                    momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
107
108         track->GetXYZ(pos);
109         track->GetCovarianceXYZPxPyPz(covTr);
110         track->GetESDpid(pid);
111         
112         Float_t impactXY, impactZ;
113         
114         track->GetImpactParameters(impactXY,impactZ);
115         
116         if (impactXY<3) {
117           // track inside the beam pipe
118           
119           AliAODTrack *aodTrack = new((*fAODCTS)[naod++]) 
120             AliAODTrack(track->GetID(), track->GetLabel(), p, kTRUE, pos, kFALSE,covTr, (Short_t)track->GetSign(), track->GetITSClusterMap(), 
121                                               pid,
122                                               0x0,//primary,
123                                               kTRUE, // check if this is right
124                                               kTRUE, // check if this is right
125                                               AliAODTrack::kPrimary, 
126                                               0);
127           
128           aodTrack->SetFlags(track->GetStatus());
129           aodTrack->ConvertAliPIDtoAODPID();
130         }
131         else continue;   // outside the beam pipe: orphan track 
132       }//Pt and Fidutial cut passed. 
133     }// track status
134   }// track loop
135   if(fDebug > 1) printf("FillInputCTS():: aod entries %d\n", fAODCTS->GetEntries());
136 }
137   
138 //____________________________________________________________________________
139 void AliCaloTrackESDReader::FillInputEMCAL() {
140   //Return array with EMCAL clusters in aod format
141
142   fAODEMCAL = new TClonesArray("AliAODCaloCluster",0);
143
144   TRefArray * caloClusters = new TRefArray();
145   fESD->GetEMCALClusters(caloClusters);
146
147   //Get vertex for momentum calculation  
148   Double_t v[3] ; //vertex ;
149   GetVertex(v);
150
151   //Loop to select clusters in fidutial cut and fill container with aodClusters
152   Int_t naod = 0;
153   Float_t pos[3] ;
154 //   Double_t * pid = new Double_t[AliPID::kSPECIESN];
155
156   for (Int_t iclus =  0; iclus <  caloClusters->GetEntries(); iclus++) {
157     AliESDCaloCluster * clus = (AliESDCaloCluster *) caloClusters->At(iclus) ;
158     TLorentzVector momentum ;
159     clus->GetMomentum(momentum, v);      
160     if(fDebug > 3 && momentum.E() > 0.1)printf("FillInputEMCAL():: all clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
161                                                momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta()); 
162     if(fEMCALPtMin < momentum.Pt() && fFidutialCut->IsInFidutialCut(momentum,"EMCAL")){
163       
164       if(fDebug > 2 && momentum.E() > 0.1)printf("FillInputEMCAL():: Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
165                                                  momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
166 //       pid=clus->GetPid();      
167       clus->GetPosition(pos) ;
168 //       printf("Reader PID ESD: ph %0.2f, pi0 %0.2f, el %0.2f, conv el %0.2f,pi %0.2f, k %0.2f, p %0.2f, k0 %0.2f, n %0.2f, mu %0.2f \n",
169 //           pid[AliPID::kPhoton],pid[AliPID::kPi0],
170 //           pid[AliPID::kElectron],pid[AliPID::kEleCon],pid[AliPID::kPion],
171 //           pid[AliPID::kKaon],pid[AliPID::kProton], pid[AliPID::kKaon0],
172 //           pid[AliPID::kNeutron], pid[AliPID::kMuon]);
173       Int_t id = clus->GetID();
174       Int_t nLabel = clus->GetNLabels();
175       Int_t *labels=0x0;
176       if(clus->GetLabels()) labels =  (clus->GetLabels())->GetArray();
177
178       Float_t energy = clus->E();
179       Char_t ttype= AliAODCluster::kEMCALClusterv1;
180       AliAODCaloCluster *caloCluster = new((*fAODEMCAL)[naod++]) 
181         AliAODCaloCluster(id,nLabel,labels,energy, pos, NULL,ttype,0);
182
183       caloCluster->SetPIDFromESD(clus->GetPid());
184       caloCluster->SetCaloCluster(clus->GetDistanceToBadChannel(), clus->GetClusterDisp(), 
185                                   clus->GetM20(), clus->GetM02(),   
186                                   clus->GetEmcCpvDistance(),  clus->GetNExMax(), clus->GetTOF()) ;
187       caloCluster->SetNCells(clus->GetNCells());
188       caloCluster->SetCellsAbsId(clus->GetCellsAbsId());
189       caloCluster->SetCellsAmplitudeFraction(clus->GetCellsAmplitudeFraction());
190
191       TArrayI* matchedT =       clus->GetTracksMatched();
192       if (matchedT && clus->GetTrackMatched() > 0) {    
193         for (Int_t im = 0; im < matchedT->GetSize(); im++) {
194           caloCluster->AddTrackMatched((fESD->GetTrack(im)));
195         }
196       }
197       
198     }//Pt and Fidutial cut passed.
199   }//esd cluster loop
200   
201   if(fDebug > 1) printf("FillInputEMCAL():: aod entries %d\n", fAODEMCAL->GetEntries());
202
203 }
204
205 //____________________________________________________________________________
206 void AliCaloTrackESDReader::FillInputPHOS() {
207   //Return array with PHOS clusters in aod format
208   fAODPHOS = new TClonesArray("AliAODCaloCluster",0);
209
210   TRefArray * caloClusters = new TRefArray();
211   fESD->GetPHOSClusters(caloClusters);
212
213   //Get vertex for momentum calculation  
214   Double_t v[3] ; //vertex ;
215   GetVertex(v);
216
217   //Loop to select clusters in fidutial cut and fill container with aodClusters
218   Int_t naod = 0;
219   Float_t pos[3] ;
220   Double_t * pid = new Double_t[AliPID::kSPECIESN];
221  
222   for (Int_t iclus =  0; iclus <  caloClusters->GetEntries(); iclus++) {
223     AliESDCaloCluster * clus = (AliESDCaloCluster *) caloClusters->At(iclus) ;
224     TLorentzVector momentum ;
225     clus->GetMomentum(momentum, v);      
226     if(fDebug > 3 && momentum.E() > 0.1)printf("FillInputPHOS():: all clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
227                                                 momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
228     if(fPHOSPtMin < momentum.Pt() && fFidutialCut->IsInFidutialCut(momentum,"PHOS")){
229       
230       if(fDebug > 2 && momentum.E() > 0.1)printf("FillInputPHOS():: Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
231                                                   momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
232  
233       pid=clus->GetPid();      
234       // printf("Reader PID ESD: ph %0.2f, pi0 %0.2f, el %0.2f, conv el %0.2f,pi %0.2f, k %0.2f, p %0.2f, k0 %0.2f, n %0.2f, mu %0.2f \n",
235       //         pid[AliPID::kPhoton],pid[AliPID::kPi0],
236       //         pid[AliPID::kElectron],pid[AliPID::kEleCon],pid[AliPID::kPion],
237       //         pid[AliPID::kKaon],pid[AliPID::kProton], pid[AliPID::kKaon0],
238       //         pid[AliPID::kNeutron], pid[AliPID::kMuon]);
239
240       clus->GetPosition(pos) ;
241       Int_t id = clus->GetID();
242       Int_t nLabel = clus->GetNLabels();
243       Int_t *labels=0x0;
244       if(clus->GetLabels()) labels =  (clus->GetLabels())->GetArray();
245       Float_t energy = clus->E();
246  
247       //Phos cluster type
248       Char_t ttype= AliAODCluster::kPHOSNeutral;
249       Float_t wNeutral = pid[AliPID::kNeutron]+ pid[AliPID::kKaon0]+pid[AliPID::kPhoton]+pid[AliPID::kPi0];
250       Float_t wCharged = pid[AliPID::kMuon] + pid[AliPID::kElectron] + pid[AliPID::kEleCon]+ 
251         pid[AliPID::kProton]+pid[AliPID::kKaon]+pid[AliPID::kPion];
252       if( wCharged > wNeutral)  ttype= AliAODCluster::kPHOSCharged;
253       
254       AliAODCaloCluster *caloCluster = new((*fAODPHOS)[naod++]) 
255         AliAODCaloCluster(id,nLabel,labels,energy, pos, NULL, ttype, 0);
256       caloCluster->SetPIDFromESD(clus->GetPid());   
257       caloCluster->SetCaloCluster(clus->GetDistanceToBadChannel(), clus->GetClusterDisp(), 
258                                   clus->GetM20(), clus->GetM02(),  
259                                   clus->GetEmcCpvDistance(),  clus->GetNExMax()) ;
260       caloCluster->SetNCells(clus->GetNCells());
261       caloCluster->SetCellsAbsId(clus->GetCellsAbsId());
262       caloCluster->SetCellsAmplitudeFraction(clus->GetCellsAmplitudeFraction());
263       TArrayI* matchedT =       clus->GetTracksMatched();
264       if (matchedT) {   
265         for (Int_t im = 0; im < matchedT->GetSize(); im++) {
266           caloCluster->AddTrackMatched((fESD->GetTrack(im)));
267         }
268       }
269      
270     }//Pt and Fidutial cut passed.
271   }//esd cluster loop
272   if(fDebug > 1) printf("FillInputPHOS():: aod entries %d\n", fAODPHOS->GetEntries());
273
274
275 }
276
277 //____________________________________________________________________________
278 void AliCaloTrackESDReader::FillInputEMCALCells() {
279   //Return array with EMCAL cells in esd format
280
281   fEMCALCells = (TNamed*) fESD->GetEMCALCells(); 
282
283 }
284
285 //____________________________________________________________________________
286 void AliCaloTrackESDReader::FillInputPHOSCells() {
287   //Return array with PHOS cells in esd format
288
289   fPHOSCells = (TNamed*) fESD->GetPHOSCells(); 
290
291 }
292
293 //____________________________________________________________________________
294 void AliCaloTrackESDReader::GetVertex(Double_t  v[3]) {
295   //Return vertex position
296
297   fESD->GetVertex()->GetXYZ(v) ;
298  
299 }
300
301
302 //____________________________________________________________________________
303 void AliCaloTrackESDReader::SetInputEvent(TObject* esd, TObject* aod, TObject* mc) {
304   // Connect the data pointers
305   SetESD((AliESDEvent*) esd);
306   SetAOD ((AliAODEvent*) aod);
307   SetMC((AliMCEvent*) mc);
308 }