]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/AliAnaGammaIsolCut.cxx
Various fixes in order to compile the DA source code
[u/mrichter/AliRoot.git] / PWG4 / AliAnaGammaIsolCut.cxx
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 /* $Id$ */
16
17 /* History of cvs commits:
18  *
19  * $Log$
20  * Revision 1.3  2007/03/08 10:24:32  schutz
21  * Coding convention
22  *
23  * Revision 1.2  2007/02/09 18:40:40  schutz
24  * B\bNew version from Gustavo
25  *
26  * Revision 1.1  2007/01/25 17:24:20  schutz
27  * new class
28  *
29  * Revision 1.1  2007/01/23 17:17:29  schutz
30  * New Gamma package
31  *
32  *
33  */
34
35 //_________________________________________________________________________
36 // Class for the analysis of gamma correlations (gamma-jet, 
37 // gamma-hadron and isolation cut.
38 // This class makes isolation cut analysis for 2 IC methods 
39 //(cone pt sum and particle pt threshold), for different cone sizes 
40 //and pt thresholds
41 //
42 //  Class created from old AliPHOSGammaJet 
43 //  (see AliRoot versions previous Release 4-09)
44 //
45 //*-- Author: Gustavo Conesa (LNF-INFN) 
46 //////////////////////////////////////////////////////////////////////////////
47
48
49 // --- ROOT system ---
50
51 #include <TFile.h>
52 #include <TParticle.h>
53 #include <TH2.h>
54
55 #include "AliAnaGammaIsolCut.h" 
56 #include "AliESD.h"
57 #include "AliESDtrack.h"
58 #include "AliESDCaloCluster.h"
59 #include "Riostream.h"
60 #include "AliLog.h"
61
62 ClassImp(AliAnaGammaIsolCut)
63
64 //____________________________________________________________________________
65   AliAnaGammaIsolCut::AliAnaGammaIsolCut(const char *name) : 
66     AliAnaGammaDirect(name),
67     fOutputContainer(new TObjArray(100)),  
68     fNCones(0),fNPtThres(0),   fConeSizes(),  fPtThresholds(), fhPtCandidate(0), fhPtThresIsolated(), fhPtSumIsolated()
69 {
70   //Ctor
71
72   //Init parameters
73   InitParameters();
74   
75   // Input slot #0 works with an Ntuple
76   DefineInput(0, TChain::Class());
77   // Output slot #0 writes into a TH1 container
78   DefineOutput(0,  TObjArray::Class()) ; 
79   
80 }
81
82
83 //____________________________________________________________________________
84 AliAnaGammaIsolCut::AliAnaGammaIsolCut(const AliAnaGammaIsolCut & ic) : 
85   AliAnaGammaDirect(ic),
86   fOutputContainer(ic.fOutputContainer), 
87   fNCones(ic.fNCones),fNPtThres(ic.fNPtThres), fConeSizes(),fPtThresholds(), 
88   fhPtCandidate(ic. fhPtCandidate), fhPtThresIsolated(), fhPtSumIsolated()
89 {
90   // cpy ctor
91   SetName (ic.GetName()) ; 
92   SetTitle(ic.GetTitle()) ; 
93
94   for(Int_t i = 0; i < fNCones ; i++){
95     fConeSizes[i] =  ic.fConeSizes[i];
96     fhPtSumIsolated[i] = ic.fhPtSumIsolated[i]; 
97       for(Int_t j = 0; j < fNPtThres ; j++)
98         fhPtThresIsolated[i][j] = ic.fhPtThresIsolated[i][j]; 
99     }
100
101   for(Int_t i = 0; i < fNPtThres ; i++)
102     fPtThresholds[i]=   ic.fPtThresholds[i];
103   
104
105 }
106
107 //_________________________________________________________________________
108 AliAnaGammaIsolCut & AliAnaGammaIsolCut::operator = (const AliAnaGammaIsolCut & source)
109 {
110   //assignment operator
111   if(&source == this) return *this;
112
113   fOutputContainer = source.fOutputContainer ;
114   fNCones = source.fNCones ;
115   fNPtThres = source.fNPtThres ; 
116   fhPtCandidate = source. fhPtCandidate ;
117
118   for(Int_t i = 0; i < fNCones ; i++){
119     fConeSizes[i] =  source.fConeSizes[i];
120     fhPtSumIsolated[i] = source.fhPtSumIsolated[i] ;
121       for(Int_t j = 0; j < fNPtThres ; j++)
122         fhPtThresIsolated[i][j] = source.fhPtThresIsolated[i][j] ;
123   }
124   
125   for(Int_t i = 0; i < fNPtThres ; i++)
126     fPtThresholds[i]=   source.fPtThresholds[i];
127   
128
129 return *this;
130
131 }
132
133 //____________________________________________________________________________
134 AliAnaGammaIsolCut::~AliAnaGammaIsolCut() 
135 {
136   // Remove all pointers
137   fOutputContainer->Clear() ; 
138   delete fOutputContainer ;
139
140   delete fhPtCandidate ;
141   delete [] fhPtThresIsolated ;
142   delete [] fhPtSumIsolated ;
143
144 }
145
146 //______________________________________________________________________________
147 void AliAnaGammaIsolCut::ConnectInputData(const Option_t*)
148 {
149   // Initialisation of branch container and histograms 
150   AliAnaGammaDirect::ConnectInputData("");
151
152 }
153
154 //____________________________________________________
155 void AliAnaGammaIsolCut::CreateOutputObjects()
156 {  
157
158   // Create histograms to be saved in output file and 
159   // store them in fOutputContainer
160
161   fOutputContainer = new TObjArray(100) ; 
162
163   //Isolation cut histograms
164   fhPtCandidate  = new TH1F
165     ("PtCandidate","p_{T} of candidate particles for isolation",240,0,120); 
166   fhPtCandidate->SetXTitle("p_{T} (GeV/c)");
167   fOutputContainer->Add(fhPtCandidate) ;
168   
169   char name[128];
170   char title[128];
171   for(Int_t icone = 0; icone<fNCones; icone++){
172     sprintf(name,"PtSumIsolated_Cone_%d",icone);
173     sprintf(title,"Candidate cone sum p_{T} for cone size %d vs candidate p_{T}",icone);
174     fhPtSumIsolated[icone]  = new TH2F(name, title,240,0,120,120,0,10);
175     fhPtSumIsolated[icone]->SetYTitle("#Sigma p_{T} (GeV/c)");
176     fhPtSumIsolated[icone]->SetXTitle("p_{T} (GeV/c)");
177     fOutputContainer->Add(fhPtSumIsolated[icone]) ; 
178     
179     for(Int_t ipt = 0; ipt<fNPtThres;ipt++){ 
180       sprintf(name,"PtThresIsol_Cone_%d_Pt%d",icone,ipt);
181       sprintf(title,"Isolated candidate p_{T} distribution for cone size %d and p_{T}^{th} %d",icone,ipt);
182       fhPtThresIsolated[icone][ipt]  = new TH1F(name, title,240,0,120);
183       fhPtThresIsolated[icone][ipt]->SetXTitle("p_{T} (GeV/c)");
184       fOutputContainer->Add(fhPtThresIsolated[icone][ipt]) ; 
185     }//icone loop
186   }//ipt loop
187
188 }
189
190 //____________________________________________________________________________
191 void AliAnaGammaIsolCut::Exec(Option_t *) 
192 {
193   
194   // Processing of one event
195   //Get ESDs
196
197   Long64_t entry = GetChain()->GetReadEntry() ;
198   if (!GetESD()) {
199     AliError("fESD is not connected to the input!") ; 
200     return ; 
201   }
202   if (GetPrintInfo()) 
203     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(GetChain()))->GetFile()->GetName(), entry)) ; 
204
205   //CreateTLists with arrays of TParticles. Filled with particles only relevant for the analysis.
206
207   TClonesArray * particleList = new TClonesArray("TParticle",1000); // All particles refitted in CTS and detected in EMCAL (jet)
208   TClonesArray * plCTS         = new TClonesArray("TParticle",1000); // All particles refitted in Central Tracking System (ITS+TPC)
209   TClonesArray * plNe         = new TClonesArray("TParticle",1000);   // All particles measured in Jet Calorimeter (EMCAL)
210   TClonesArray * plPHOS     = new TClonesArray("TParticle",1000);  // All particles measured in PHOS as Gamma calorimeter
211   TClonesArray * plEMCAL   = new TClonesArray("TParticle",1000);  // All particles measured in EMCAL as Gamma calorimeter
212   
213   
214   AliDebug(2, "Fill particle lists");
215   
216   //Fill particle lists 
217   CreateParticleList(particleList, plCTS,plEMCAL,plPHOS); 
218   
219   if(GetCalorimeter() == "PHOS")
220     plNe = plPHOS;
221   if(GetCalorimeter() == "EMCAL")
222     plNe = plEMCAL;
223   
224   //Isolation Cut Analysis for both methods and different pt cuts and cones
225   
226   for(Int_t ipr = 0; ipr < plNe->GetEntries() ; ipr ++ ){
227     TParticle * pCandidate = dynamic_cast<TParticle *>(plNe->At(ipr)) ;
228     
229     if(pCandidate->Pt() > GetMinGammaPt()){
230       
231       Bool_t  icPtThres   = kFALSE;
232       Bool_t  icPtSum     = kFALSE;
233       
234       Float_t ptC             = pCandidate->Pt() ;
235       
236       fhPtCandidate->Fill(ptC);
237       
238       for(Int_t icone = 0; icone<fNCones; icone++){
239         SetConeSize(fConeSizes[icone]) ;
240         Float_t coneptsum = 0 ;
241         for(Int_t ipt = 0; ipt<fNPtThres;ipt++){ 
242           SetPtThreshold(fPtThresholds[ipt]) ;
243           MakeIsolationCut(plCTS,plNe, pCandidate, ipr, icPtThres, icPtSum,coneptsum);
244           AliDebug(4,Form("Candidate pt %f, pt in cone %f, Isolated? ICPt %d, ICSum %d",
245                           pCandidate->Pt(), coneptsum, icPtThres, icPtSum));
246           
247           fhPtThresIsolated[icone][ipt]->Fill(ptC); 
248         }//pt thresh loop
249         fhPtSumIsolated[icone]->Fill(ptC,coneptsum) ;
250       }//cone size loop
251     }//min pt candidate
252   }//candidate loop
253   
254   AliDebug(2, "End of analysis, delete pointers");
255   
256   particleList->Delete() ; 
257   plCTS->Delete() ;
258   plNe->Delete() ;
259   plPHOS->Delete() ;
260   plEMCAL->Delete() ;
261
262   delete plNe ;
263   delete plCTS ;
264   //delete plPHOS ;
265   //delete plEMCAL ;
266   delete particleList ;
267
268   PostData(0, fOutputContainer);
269 }    
270
271   //____________________________________________________________________________
272 void AliAnaGammaIsolCut::InitParameters()
273 {
274   // Initialisation of branch container 
275   AliAnaGammaDirect::InitParameters();
276
277   fNCones           = 4 ; 
278   fNPtThres         = 4 ; 
279   fConeSizes[0] = 0.1; fConeSizes[1] = 0.2; fConeSizes[2] = 0.3; fConeSizes[3] = 0.4;
280   fPtThresholds[0]=1.; fPtThresholds[1]=2.; fPtThresholds[2]=3.; fPtThresholds[3]=4.;
281 }
282
283 void AliAnaGammaIsolCut::Print(const Option_t * opt) const
284 {
285
286   //Print some relevant parameters set for the analysis
287   if(! opt)
288     return;
289
290   Info("Print", "%s %s", GetName(), GetTitle() ) ;
291   printf("N Cone Sizes               =     %d\n", fNCones) ; 
292   printf("N pT thresholds           =     %d\n", fNPtThres) ;
293   printf("Cone Sizes                  =    \n") ;
294   for(Int_t i = 0; i < fNCones; i++)
295     printf("   %f;",  fConeSizes[i]) ;
296   printf("    \n") ;
297   for(Int_t i = 0; i < fNPtThres; i++)
298     printf("   %f;",  fPtThresholds[i]) ;
299   printf("    \n") ;
300   
301
302
303 void AliAnaGammaIsolCut::Terminate(Option_t *)
304 {
305   // The Terminate() function is the last function to be called during
306   // a query. It always runs on the client, it can be used to present
307   // the results graphically or save the results to file.
308  
309 }