1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 /* History of cvs commits:
20 * Revision 1.2 2007/02/09 18:40:40 schutz
21 * B
\bNew version from Gustavo
23 * Revision 1.1 2007/01/25 17:24:20 schutz
26 * Revision 1.1 2007/01/23 17:17:29 schutz
32 //_________________________________________________________________________
33 // Class for the analysis of gamma correlations (gamma-jet,
34 // gamma-hadron and isolation cut.
35 // This class makes isolation cut analysis for 2 IC methods
36 //(cone pt sum and particle pt threshold), for different cone sizes
39 // Class created from old AliPHOSGammaJet
40 // (see AliRoot versions previous Release 4-09)
42 //*-- Author: Gustavo Conesa (LNF-INFN)
43 //////////////////////////////////////////////////////////////////////////////
46 // --- ROOT system ---
49 #include <TParticle.h>
52 #include "AliAnaGammaIsolCut.h"
54 #include "AliESDtrack.h"
55 #include "AliESDCaloCluster.h"
56 #include "Riostream.h"
59 ClassImp(AliAnaGammaIsolCut)
61 //____________________________________________________________________________
62 AliAnaGammaIsolCut::AliAnaGammaIsolCut(const char *name) :
63 AliAnaGammaDirect(name),
64 fOutputContainer(new TObjArray(100)),
65 fNCones(0),fNPtThres(0), fConeSizes(), fPtThresholds(), fhPtCandidate(0), fhPtThresIsolated(), fhPtSumIsolated()
72 TList * list = gDirectory->GetListOfKeys() ;
76 for (index = 0 ; index < list->GetSize()-1 ; index++) {
77 //-1 to avoid GammaJet Task
78 h = dynamic_cast<TH2F*>(gDirectory->Get(list->At(index)->GetName())) ;
79 fOutputContainer->Add(h) ;
82 // Input slot #0 works with an Ntuple
83 DefineInput(0, TChain::Class());
84 // Output slot #0 writes into a TH1 container
85 DefineOutput(0, TObjArray::Class()) ;
90 //____________________________________________________________________________
91 AliAnaGammaIsolCut::AliAnaGammaIsolCut(const AliAnaGammaIsolCut & ic) :
92 AliAnaGammaDirect(ic),
93 fOutputContainer(ic.fOutputContainer),
94 fNCones(ic.fNCones),fNPtThres(ic.fNPtThres), fConeSizes(),fPtThresholds(),
95 fhPtCandidate(ic. fhPtCandidate), fhPtThresIsolated(), fhPtSumIsolated()
98 SetName (ic.GetName()) ;
99 SetTitle(ic.GetTitle()) ;
101 for(Int_t i = 0; i < fNCones ; i++){
102 fConeSizes[i] = ic.fConeSizes[i];
103 fhPtSumIsolated[i] = ic.fhPtSumIsolated[i];
104 for(Int_t j = 0; j < fNPtThres ; j++)
105 fhPtThresIsolated[i][j] = ic.fhPtThresIsolated[i][j];
108 for(Int_t i = 0; i < fNPtThres ; i++)
109 fPtThresholds[i]= ic.fPtThresholds[i];
114 //_________________________________________________________________________
115 AliAnaGammaIsolCut & AliAnaGammaIsolCut::operator = (const AliAnaGammaIsolCut & source)
117 //assignment operator
118 if(&source == this) return *this;
120 fOutputContainer = source.fOutputContainer ;
121 fNCones = source.fNCones ;
122 fNPtThres = source.fNPtThres ;
123 fhPtCandidate = source. fhPtCandidate ;
125 for(Int_t i = 0; i < fNCones ; i++){
126 fConeSizes[i] = source.fConeSizes[i];
127 fhPtSumIsolated[i] = source.fhPtSumIsolated[i] ;
128 for(Int_t j = 0; j < fNPtThres ; j++)
129 fhPtThresIsolated[i][j] = source.fhPtThresIsolated[i][j] ;
132 for(Int_t i = 0; i < fNPtThres ; i++)
133 fPtThresholds[i]= source.fPtThresholds[i];
140 //____________________________________________________________________________
141 AliAnaGammaIsolCut::~AliAnaGammaIsolCut()
143 // Remove all pointers
144 fOutputContainer->Clear() ;
145 delete fOutputContainer ;
147 delete fhPtCandidate ;
148 delete [] fhPtThresIsolated ;
149 delete [] fhPtSumIsolated ;
153 //______________________________________________________________________________
154 void AliAnaGammaIsolCut::ConnectInputData(const Option_t*)
156 // Initialisation of branch container and histograms
157 AliAnaGammaDirect::ConnectInputData("");
161 //____________________________________________________
162 void AliAnaGammaIsolCut::CreateOutputObjects()
165 // Create histograms to be saved in output file and
166 // store them in fOutputContainer
168 fOutputContainer = new TObjArray(100) ;
170 //Isolation cut histograms
171 fhPtCandidate = new TH1F
172 ("PtCandidate","p_{T} of candidate particles for isolation",240,0,120);
173 fhPtCandidate->SetXTitle("p_{T} (GeV/c)");
174 fOutputContainer->Add(fhPtCandidate) ;
178 for(Int_t icone = 0; icone<fNCones; icone++){
179 sprintf(name,"PtSumIsolated_Cone_%d",icone);
180 sprintf(title,"Candidate cone sum p_{T} for cone size %d vs candidate p_{T}",icone);
181 fhPtSumIsolated[icone] = new TH2F(name, title,240,0,120,120,0,10);
182 fhPtSumIsolated[icone]->SetYTitle("#Sigma p_{T} (GeV/c)");
183 fhPtSumIsolated[icone]->SetXTitle("p_{T} (GeV/c)");
184 fOutputContainer->Add(fhPtSumIsolated[icone]) ;
186 for(Int_t ipt = 0; ipt<fNPtThres;ipt++){
187 sprintf(name,"PtThresIsol_Cone_%d_Pt%d",icone,ipt);
188 sprintf(title,"Isolated candidate p_{T} distribution for cone size %d and p_{T}^{th} %d",icone,ipt);
189 fhPtThresIsolated[icone][ipt] = new TH1F(name, title,240,0,120);
190 fhPtThresIsolated[icone][ipt]->SetXTitle("p_{T} (GeV/c)");
191 fOutputContainer->Add(fhPtThresIsolated[icone][ipt]) ;
197 //____________________________________________________________________________
198 void AliAnaGammaIsolCut::Exec(Option_t *)
201 // Processing of one event
204 Long64_t entry = GetChain()->GetReadEntry() ;
206 AliError("fESD is not connected to the input!") ;
210 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(GetChain()))->GetFile()->GetName(), entry)) ;
212 //CreateTLists with arrays of TParticles. Filled with particles only relevant for the analysis.
214 TClonesArray * particleList = new TClonesArray("TParticle",1000); // All particles refitted in CTS and detected in EMCAL (jet)
215 TClonesArray * plCTS = new TClonesArray("TParticle",1000); // All particles refitted in Central Tracking System (ITS+TPC)
216 TClonesArray * plNe = new TClonesArray("TParticle",1000); // All particles measured in Jet Calorimeter (EMCAL)
217 TClonesArray * plPHOS = new TClonesArray("TParticle",1000); // All particles measured in PHOS as Gamma calorimeter
218 TClonesArray * plEMCAL = new TClonesArray("TParticle",1000); // All particles measured in EMCAL as Gamma calorimeter
221 AliDebug(2, "Fill particle lists");
223 //Fill particle lists
224 CreateParticleList(particleList, plCTS,plEMCAL,plPHOS);
226 if(GetCalorimeter() == "PHOS")
228 if(GetCalorimeter() == "EMCAL")
231 //Isolation Cut Analysis for both methods and different pt cuts and cones
233 for(Int_t ipr = 0; ipr < plNe->GetEntries() ; ipr ++ ){
234 TParticle * pCandidate = dynamic_cast<TParticle *>(plNe->At(ipr)) ;
236 if(pCandidate->Pt() > GetMinGammaPt()){
238 Bool_t icPtThres = kFALSE;
239 Bool_t icPtSum = kFALSE;
241 Float_t ptC = pCandidate->Pt() ;
243 fhPtCandidate->Fill(ptC);
245 for(Int_t icone = 0; icone<fNCones; icone++){
246 SetConeSize(fConeSizes[icone]) ;
247 Float_t coneptsum = 0 ;
248 for(Int_t ipt = 0; ipt<fNPtThres;ipt++){
249 SetPtThreshold(fPtThresholds[ipt]) ;
250 MakeIsolationCut(plCTS,plNe, pCandidate, ipr, icPtThres, icPtSum,coneptsum);
251 AliDebug(4,Form("Candidate pt %f, pt in cone %f, Isolated? ICPt %d, ICSum %d",
252 pCandidate->Pt(), coneptsum, icPtThres, icPtSum));
254 fhPtThresIsolated[icone][ipt]->Fill(ptC);
256 fhPtSumIsolated[icone]->Fill(ptC,coneptsum) ;
261 AliDebug(2, "End of analysis, delete pointers");
263 particleList->Delete() ;
273 delete particleList ;
275 PostData(0, fOutputContainer);
278 //____________________________________________________________________________
279 void AliAnaGammaIsolCut::InitParameters()
281 // Initialisation of branch container
282 AliAnaGammaDirect::InitParameters();
286 fConeSizes[0] = 0.1; fConeSizes[1] = 0.2; fConeSizes[2] = 0.3; fConeSizes[3] = 0.4;
287 fPtThresholds[0]=1.; fPtThresholds[1]=2.; fPtThresholds[2]=3.; fPtThresholds[3]=4.;
290 void AliAnaGammaIsolCut::Print(const Option_t * opt) const
293 //Print some relevant parameters set for the analysis
297 Info("Print", "%s %s", GetName(), GetTitle() ) ;
298 printf("N Cone Sizes = %d\n", fNCones) ;
299 printf("N pT thresholds = %d\n", fNPtThres) ;
300 printf("Cone Sizes = \n") ;
301 for(Int_t i = 0; i < fNCones; i++)
302 printf(" %f;", fConeSizes[i]) ;
304 for(Int_t i = 0; i < fNPtThres; i++)
305 printf(" %f;", fPtThresholds[i]) ;
310 void AliAnaGammaIsolCut::Terminate(Option_t *)
312 // The Terminate() function is the last function to be called during
313 // a query. It always runs on the client, it can be used to present
314 // the results graphically or save the results to file.