]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/AliAnaGammaPhos.cxx
Forgot this one
[u/mrichter/AliRoot.git] / PWG4 / AliAnaGammaPhos.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
16 /* $Id: */
17
18 //_________________________________________________________________________
19 // A basic analysis task to analyse photon detected by PHOS
20 //
21 //*-- Yves Schutz 
22 //////////////////////////////////////////////////////////////////////////////
23
24 #include <TCanvas.h>
25 #include <TChain.h>
26 #include <TFile.h> 
27 #include <TH1.h>
28 #include <TH1F.h>
29 #include <TH1I.h>
30 #include <TLegend.h> 
31 #include <TNtuple.h>
32 #include <TROOT.h> 
33 #include <TVector3.h> 
34
35 #include "AliAnaGammaPhos.h" 
36 #include "AliAnalysisManager.h"
37 #include "AliESDEvent.h" 
38 #include "AliESDCaloCluster.h" 
39 #include "AliAODEvent.h"
40 #include "AliAODHandler.h"
41 #include "AliAODPhoton.h"
42 #include "AliLog.h"
43
44 //______________________________________________________________________________
45 AliAnaGammaPhos::AliAnaGammaPhos() : 
46   fChain(0x0),
47   fDebug(0),
48   fESD(0x0),
49   fAOD(0x0),
50   fAODPhotons(0x0), 
51   fPhotonsInPhos(0),
52   fTreeA(0x0),
53   fPhotonId(1.0),
54   fOutputList(0x0), 
55   fhPHOS(0),
56   fhPHOSEnergy(0),
57   fhPHOSDigits(0),
58   fhPHOSRecParticles(0),
59   fhPHOSPhotons(0),
60   fhPHOSInvariantMass(0),
61   fhPHOSDigitsEvent(0)
62 {
63   //Default constructor
64 }
65 //______________________________________________________________________________
66 AliAnaGammaPhos::AliAnaGammaPhos(const char *name) : 
67   AliAnalysisTask(name,""),  
68   fChain(0x0),
69   fDebug(0),
70   fESD(0x0),
71   fAOD(0x0),
72   fAODPhotons(0x0), 
73   fPhotonsInPhos(0),
74   fTreeA(0x0),
75   fPhotonId(1.0),
76   fOutputList(0x0), 
77   fhPHOS(0),
78   fhPHOSEnergy(0),
79   fhPHOSDigits(0),
80   fhPHOSRecParticles(0),
81   fhPHOSPhotons(0),
82   fhPHOSInvariantMass(0),
83   fhPHOSDigitsEvent(0)
84 {
85   // Constructor.
86   // Input slot #0 
87   DefineInput(0, TChain::Class());
88   // Output slots 
89   DefineOutput(0,  TTree::Class()) ; 
90   DefineOutput(1,  TList::Class()) ; 
91 }
92
93 //______________________________________________________________________________
94 AliAnaGammaPhos::~AliAnaGammaPhos()
95 {
96   // dtor
97   //  fOutputList->Clear() ; 
98   //delete fOutputList ;
99 }
100
101
102 //______________________________________________________________________________
103 void AliAnaGammaPhos::ConnectInputData(const Option_t*)
104 {
105   // Initialisation of branch container and histograms 
106     
107   AliInfo(Form("*** Initialization of %s", GetName())) ; 
108   
109   // Get input data
110   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
111   if (!fChain) {
112     AliError(Form("Input 0 for %s not found\n", GetName()));
113     return ;
114   }
115   
116   fESD = new AliESDEvent();
117   fESD->ReadFromTree(fChain);
118
119 }
120
121 //________________________________________________________________________
122 void AliAnaGammaPhos::CreateOutputObjects()
123 {  
124   // Create the outputs containers
125  
126   OpenFile(0) ;
127   AliAODHandler* handler = (AliAODHandler*) ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());  
128   fTreeA = handler->GetTree() ; 
129   fAOD   = handler->GetAOD();
130   fAODPhotons = fAOD->GetClusters() ; 
131   
132
133   OpenFile(1) ; 
134
135   fhPHOSPos            = new TNtuple("PHOSPos"         , "Position in PHOS"  , "x:y:z");
136   fhPHOS               = new TNtuple("PHOS"            , "PHOS"  , "event:digits:clusters:photons");
137   fhPHOSEnergy         = new TH1D("PHOSEnergy"         , "PHOSEnergy"        , 100, 0., 100. ) ;
138   fhPHOSDigits         = new TH1I("PHOSDigitsCluster"  , "PHOSDigits"        , 20 , 0 , 20  ) ;
139   fhPHOSRecParticles   = new TH1D("PHOSRecParticles"   , "PHOSRecParticles" , 20 , 0., 20. ) ;
140   fhPHOSPhotons        = new TH1I("PHOSPhotons"        , "PHOSPhotons"       , 20 , 0 , 20  ) ;
141   fhPHOSInvariantMass  = new TH1D("PHOSInvariantMass"  , "PHOSInvariantMass" , 400, 0., 400.) ;
142   fhPHOSDigitsEvent    = new TH1I("PHOSDigitsEvent"    , "PHOSDigitsEvent"   , 30 , 0 , 30  ) ;
143   
144   // create output container
145   
146   fOutputList = new TList() ; 
147   fOutputList->SetName(GetName()) ; 
148
149   fOutputList->AddAt(fhPHOSPos,             0) ; 
150   fOutputList->AddAt(fhPHOS,                1) ; 
151   fOutputList->AddAt(fhPHOSEnergy,          2) ; 
152   fOutputList->AddAt(fhPHOSDigits,          3) ; 
153   fOutputList->AddAt(fhPHOSRecParticles,    4) ; 
154   fOutputList->AddAt(fhPHOSPhotons,         5) ; 
155   fOutputList->AddAt(fhPHOSInvariantMass,   6) ; 
156   fOutputList->AddAt(fhPHOSDigitsEvent,     7) ; 
157 }
158
159 //______________________________________________________________________________
160 void AliAnaGammaPhos::Exec(Option_t *) 
161 {
162   // Processing of one event
163     
164   Long64_t entry = fChain->GetReadEntry() ;
165   
166   if (!fESD) {
167     AliError("fESD is not connected to the input!") ; 
168     return ; 
169   }
170   
171   if ( !((entry-1)%100) ) 
172     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
173   
174   //************************  PHOS *************************************
175       
176   Int_t       firstPhosCluster       = fESD->GetFirstPHOSCluster() ;
177   const Int_t numberOfPhosClusters   = fESD->GetNumberOfPHOSClusters() ;
178   
179   TVector3 ** phosVector       = new TVector3*[numberOfPhosClusters] ;
180   Float_t  * phosPhotonsEnergy = new Float_t[numberOfPhosClusters] ;
181   Int_t      phosCluster ; 
182   Int_t      numberOfDigitsInPhos   = 0 ;
183   
184   fPhotonsInPhos  = 0 ;
185   // loop over the PHOS Cluster
186   for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + numberOfPhosClusters ; phosCluster++) {
187     AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
188     if (caloCluster) {
189       Float_t pos[3] ;
190       caloCluster->GetPosition( pos ) ;
191       fhPHOSEnergy->Fill( caloCluster->E() ) ;
192       fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
193       fhPHOSDigits->Fill(entry, caloCluster->GetNumberOfDigits() ) ;
194       numberOfDigitsInPhos += caloCluster->GetNumberOfDigits() ;
195       Float_t * pid = caloCluster->GetPid() ;
196       if(pid[AliPID::kPhoton] > GetPhotonId() ) {
197         phosVector[fPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
198         phosPhotonsEnergy[fPhotonsInPhos]=caloCluster->E() ;
199         //new ((*fAODPhotons)[fPhotonsInPhos++;]) AliAODPhoton ( );
200       }
201     }
202   } //PHOS clusters
203     
204   fhPHOSRecParticles->Fill(numberOfPhosClusters);
205   fhPHOSPhotons->Fill(fPhotonsInPhos);
206   fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
207   fhPHOS->Fill(entry, numberOfDigitsInPhos, numberOfPhosClusters, fPhotonsInPhos) ; 
208
209   // invariant Mass
210   if (fPhotonsInPhos > 1 ) {
211     Int_t phosPhoton1, phosPhoton2 ; 
212     for(phosPhoton1 = 0 ; phosPhoton1 < fPhotonsInPhos ; phosPhoton1++) {
213       for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < fPhotonsInPhos ; phosPhoton2++) {      
214         Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
215                                         ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) ) 
216                                         );
217         fhPHOSInvariantMass->Fill(tempMass*1000.);
218       }
219     }    
220   }
221     
222   PostData(0, fTreeA) ; 
223   PostData(1, fOutputList);
224
225   delete [] phosVector ; 
226   delete [] phosPhotonsEnergy ; 
227   
228 }
229
230
231 //______________________________________________________________________________
232 void AliAnaGammaPhos::Init()
233 {
234   // Intialisation of parameters
235   AliInfo("Doing initialisation") ; 
236   SetPhotonId(0.9) ; 
237 }
238
239 //______________________________________________________________________________
240 void AliAnaGammaPhos::Terminate(Option_t *)
241 {
242   // Processing when the event loop is ended
243
244 //   fOutputList = (TList*)GetOutputData(0);  
245 //   fhPHOSPos            = (TNtuple*)fOutputList->At(0);
246 //   fhPHOS               = (TNtuple*)fOutputList->At(1);
247 //   fhPHOSEnergy         = (TH1D*)fOutputList->At(2);
248 //   fhPHOSDigits         = (TH1I*)fOutputList->At(3);
249 //   fhPHOSRecParticles   = (TH1D*)fOutputList->At(4);
250 //   fhPHOSPhotons        = (TH1I*)fOutputList->At(5);
251 //   fhPHOSInvariantMass  = (TH1D*)fOutputList->At(6);
252 //   fhPHOSDigitsEvent    = (TH1I*)fOutputList->At(7);
253
254   
255   Bool_t problem = kFALSE ; 
256   AliInfo(Form(" *** %s Report:", GetName())) ; 
257   printf("        PHOSEnergy Mean         : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(),         fhPHOSEnergy->GetRMS()         ) ;
258   printf("        PHOSDigits Mean         : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(),         fhPHOSDigits->GetRMS()         ) ;
259   printf("        PHOSRecParticles Mean   : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(),   fhPHOSRecParticles->GetRMS()   ) ;
260   printf("        PHOSPhotons Mean        : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(),        fhPHOSPhotons->GetRMS()        ) ;
261   printf("        PHOSInvariantMass Mean  : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(),  fhPHOSInvariantMass->GetRMS()  ) ;
262   printf("        PHOSDigitsEvent Mean    : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(),    fhPHOSDigitsEvent->GetRMS()    ) ;
263
264   TCanvas  * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
265   cPHOS->Divide(3, 2);
266
267   cPHOS->cd(1) ; 
268   if ( fhPHOSEnergy->GetMaximum() > 0. ) 
269     gPad->SetLogy();
270   fhPHOSEnergy->SetAxisRange(0, 25.);
271   fhPHOSEnergy->SetLineColor(2);
272   fhPHOSEnergy->Draw();
273
274   cPHOS->cd(2) ; 
275   fhPHOSDigits->SetAxisRange(0,25.);
276   fhPHOSDigits->SetLineColor(2);
277   fhPHOSDigits->Draw();
278
279   cPHOS->cd(3) ; 
280   if ( fhPHOSRecParticles->GetMaximum() > 0. ) 
281     gPad->SetLogy();
282   fhPHOSRecParticles->SetAxisRange(0, 25.);
283   fhPHOSRecParticles->SetLineColor(2);
284   fhPHOSRecParticles->Draw();
285
286   cPHOS->cd(4) ; 
287   if ( fhPHOSPhotons->GetMaximum() > 0. ) 
288     gPad->SetLogy();
289   fhPHOSPhotons->SetAxisRange(0,25.);
290   fhPHOSPhotons->SetLineColor(2);
291   fhPHOSPhotons->Draw();
292
293   cPHOS->cd(5) ; 
294   fhPHOSInvariantMass->SetLineColor(2);
295   fhPHOSInvariantMass->Draw();
296  
297   cPHOS->cd(6) ; 
298   if ( fhPHOSDigitsEvent->GetMaximum() > 0. ) 
299     gPad->SetLogy();
300   fhPHOSDigitsEvent->SetAxisRange(0,40.);
301   fhPHOSDigitsEvent->SetLineColor(2);
302   fhPHOSDigitsEvent->Draw();
303  
304   cPHOS->Print("PHOS.eps");
305  
306   char line[1024] ; 
307   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
308   gROOT->ProcessLine(line);
309   sprintf(line, ".!rm -fR *.eps"); 
310   gROOT->ProcessLine(line);
311  
312   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
313
314   char * report ; 
315   if(problem)
316     report="Problems found, please check!!!";  
317   else 
318     report="OK";
319
320   AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report)) ; 
321 }