]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliEMCALQATask.cxx
1c94fe2656c60331111dca3fb1b250e16a354001
[u/mrichter/AliRoot.git] / ESDCheck / AliEMCALQATask.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 // An analysis task to check the EMCAL photon data in simulated data
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 "AliEMCALQATask.h" 
36 #include "AliESD.h" 
37 #include "AliLog.h"
38
39 //______________________________________________________________________________
40 AliEMCALQATask::AliEMCALQATask(const char *name) : 
41   AliAnalysisTask(name,""),  
42   fChain(0),
43   fESD(0), 
44   fhEMCALPos(0),
45   fhEMCAL(0),
46   fhEMCALEnergy(0),
47   fhEMCALDigits(0),
48   fhEMCALRecParticles(0),
49   fhEMCALPhotons(0),
50   fhEMCALInvariantMass(0),
51   fhEMCALDigitsEvent(0)
52 {
53   // Constructor.
54   // Input slot #0 works with an Ntuple
55   DefineInput(0, TChain::Class());
56   // Output slot #0 writes into a TH1 container
57   DefineOutput(0,  TObjArray::Class()) ; 
58 }
59
60 //______________________________________________________________________________
61 AliEMCALQATask::~AliEMCALQATask()
62 {
63   // dtor
64   fOutputContainer->Clear() ; 
65   delete fOutputContainer ; 
66
67   delete fhEMCALPos ;
68   delete fhEMCAL ;
69   delete fhEMCALEnergy ;
70   delete fhEMCALDigits ;
71   delete fhEMCALRecParticles ;
72   delete fhEMCALPhotons ;
73   delete fhEMCALInvariantMass ;
74   delete fhEMCALDigitsEvent ;
75 }
76
77 //______________________________________________________________________________
78 void AliEMCALQATask::Init(const Option_t*)
79 {
80   // Initialisation of branch container and histograms 
81     
82   AliInfo(Form("*** Initialization of %s", GetName())) ; 
83   
84   // Get input data
85   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
86   if (!fChain) {
87     AliError(Form("Input 0 for %s not found\n", GetName()));
88     return ;
89   }
90   
91   if (!fESD) {
92     // One should first check if the branch address was taken by some other task
93     char ** address = (char **)GetBranchAddress(0, "ESD") ;
94     if (address) 
95       fESD = (AliESD *)(*address) ; 
96     if (!fESD) 
97       fChain->SetBranchAddress("ESD", &fESD) ;  
98   }
99   // The output objects will be written to 
100   TDirectory * cdir = gDirectory ; 
101   // Open a file for output #0
102   char outputName[1024] ; 
103   sprintf(outputName, "%s.root", GetName() ) ; 
104   OpenFile(0, outputName , "RECREATE") ; 
105   if (cdir) 
106     cdir->cd() ; 
107   
108   // create histograms 
109   
110   fhEMCALPos           = new TNtuple("EMCALPos"        , "Position in EMCAL" , "x:y:z");
111   fhEMCAL              = new TNtuple("EMCAL"           , "EMCAL" , "event:digits:clusters:photons");
112   fhEMCALEnergy        = new TH1D("EMCALEnergy"        , "EMCALEnergy"       , 1000, 0., 10. ) ;
113   fhEMCALDigits        = new TH1I("EMCALDigitsCluster" , "EMCALDigits"       , 20 , 0 , 20  ) ;
114   fhEMCALRecParticles  = new TH1D("EMCALRecParticles"  , "EMCALRecParticles", 20 , 0., 20. ) ;
115   fhEMCALPhotons       = new TH1I("EMCALPhotons"       , "EMCALPhotons"      , 20 , 0 , 20  ) ;
116   fhEMCALInvariantMass = new TH1D("EMCALInvariantMass" , "EMCALInvariantMass", 400, 0., 400.) ;
117   fhEMCALDigitsEvent   = new TH1I("EMCALDigitsEvent"   , "EMCALDigitsEvent"  , 30 , 0 , 30  ) ;
118   
119   // create output container
120   
121   fOutputContainer = new TObjArray(8) ; 
122   fOutputContainer->SetName(GetName()) ; 
123
124   fOutputContainer->AddAt(fhEMCALPos,            0) ; 
125   fOutputContainer->AddAt(fhEMCAL,               1) ; 
126   fOutputContainer->AddAt(fhEMCALEnergy,         2) ; 
127   fOutputContainer->AddAt(fhEMCALDigits,         3) ; 
128   fOutputContainer->AddAt(fhEMCALRecParticles,   4) ; 
129   fOutputContainer->AddAt(fhEMCALPhotons,        5) ; 
130   fOutputContainer->AddAt(fhEMCALInvariantMass,  6) ; 
131   fOutputContainer->AddAt(fhEMCALDigitsEvent,    7) ; 
132 }
133
134 //______________________________________________________________________________
135 void AliEMCALQATask::Exec(Option_t *) 
136 {
137   // Processing of one event
138  
139   Long64_t entry = fChain->GetReadEntry() ;
140  
141   if (!fESD) {
142     AliError("fESD is not connected to the input!") ; 
143     return ; 
144   }
145   
146   if ( !((entry-1)%100) ) 
147     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
148   
149   //************************  EMCAL *************************************
150   Int_t       firstEmcalCluster       = fESD->GetFirstEMCALCluster() ;
151   const Int_t mumberOfEmcalClusters   = fESD->GetNumberOfEMCALClusters() ;
152
153   TVector3 ** emcalVector        = new TVector3*[mumberOfEmcalClusters] ;
154   Float_t   * emcalPhotonsEnergy = new Float_t[mumberOfEmcalClusters] ;
155   Int_t      emcalCluster ; 
156   Int_t      numberOfEmcalClustersv1 = 0 ; 
157   Int_t      numberOfPhotonsInEmcal  = 0 ;
158   Int_t      numberOfDigitsInEmcal   = 0 ;  
159
160
161   // loop over all the EMCAL Cluster
162   for(emcalCluster = firstEmcalCluster ; emcalCluster < firstEmcalCluster + mumberOfEmcalClusters ; emcalCluster++) {
163     AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(emcalCluster) ;
164     if (caloCluster) {
165       Float_t pos[3] ;
166       fhEMCALPos->Fill(pos[0],pos[1],pos[2]) ;
167       if(caloCluster->GetClusterType() == AliESDCaloCluster::kClusterv1) {  
168         caloCluster->GetGlobalPosition(pos) ;
169         fhEMCALEnergy->Fill(caloCluster->GetClusterEnergy()) ;
170         fhEMCALDigits->Fill(entry, caloCluster->GetNumberOfDigits()) ;
171         numberOfEmcalClustersv1++ ;
172         numberOfDigitsInEmcal += caloCluster->GetNumberOfDigits() ;    
173         // Float_t * pid = clus->GetPid() ;
174         // if(pid[AliPID::kPhoton]>0.9){
175         emcalVector[numberOfPhotonsInEmcal] = new TVector3(pos[0],pos[1],pos[2]) ;
176         emcalPhotonsEnergy[numberOfPhotonsInEmcal] = caloCluster->GetClusterEnergy() ;
177         numberOfPhotonsInEmcal++ ; 
178       }
179     }
180   } // EMCAL clusters loop
181
182   fhEMCALRecParticles->Fill(numberOfEmcalClustersv1);
183   fhEMCALPhotons->Fill(numberOfPhotonsInEmcal);
184   fhEMCALDigitsEvent->Fill(numberOfDigitsInEmcal);
185   fhEMCAL->Fill(entry, numberOfDigitsInEmcal, numberOfEmcalClustersv1, numberOfPhotonsInEmcal) ; 
186   
187   // invariant Mass
188   if (numberOfPhotonsInEmcal > 1 ) {
189     Int_t emcalPhoton1, emcalPhoton2 ; 
190     for(emcalPhoton1 = 0 ; emcalPhoton1 < numberOfPhotonsInEmcal ; emcalPhoton1++) {
191       for(emcalPhoton2 = emcalPhoton1 + 1 ; emcalPhoton2 < numberOfPhotonsInEmcal ; emcalPhoton2++) {
192         Float_t tempMass = TMath::Sqrt( 2 * emcalPhotonsEnergy[emcalPhoton1] * emcalPhotonsEnergy[emcalPhoton2] * 
193                                         ( 1 - TMath::Cos( emcalVector[emcalPhoton1]->Angle(*emcalVector[emcalPhoton2])) 
194                                           )
195                                         );
196         fhEMCALInvariantMass->Fill(tempMass*1000.);
197       }
198     }    
199   }
200   
201   PostData(0, fOutputContainer);
202   
203   delete [] emcalVector ; 
204   delete [] emcalPhotonsEnergy ;
205 }
206
207 //______________________________________________________________________________
208 void AliEMCALQATask::Terminate(Option_t *)
209 {
210   // Processing when the event loop is ended
211   
212   printf("EMCALEnergy Mean        : %5.3f , RMS : %5.3f \n", fhEMCALEnergy->GetMean(),        fhEMCALEnergy->GetRMS()        ) ;
213   printf("EMCALDigits Mean        : %5.3f , RMS : %5.3f \n", fhEMCALDigits->GetMean(),        fhEMCALDigits->GetRMS()        ) ;
214   printf("EMCALRecParticles Mean  : %5.3f , RMS : %5.3f \n", fhEMCALRecParticles->GetMean(),  fhEMCALRecParticles->GetRMS()  ) ;
215   printf("EMCALPhotons Mean       : %5.3f , RMS : %5.3f \n", fhEMCALPhotons->GetMean(),       fhEMCALPhotons->GetRMS()       ) ;
216   printf("EMCALInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhEMCALInvariantMass->GetMean(), fhEMCALInvariantMass->GetRMS() ) ;
217   printf("EMCALDigitsEvent Mean   : %5.3f , RMS : %5.3f \n", fhEMCALDigitsEvent->GetMean(),   fhEMCALDigitsEvent->GetRMS()   ) ;
218
219   TCanvas  * cEMCAL = new TCanvas("EMCAL", "EMCAL ESD Test", 400, 10, 600, 700);
220   cEMCAL->Divide(3, 2) ; 
221
222   cEMCAL->cd(1) ; 
223   gPad->SetLogy();
224   fhEMCALEnergy->SetAxisRange(0, 25.);
225   fhEMCALEnergy->SetXTitle("Energy (GeV)");
226   fhEMCALEnergy->Draw();
227   
228   cEMCAL->cd(2) ; 
229   gPad->SetLogy();
230   fhEMCALDigits->SetAxisRange(0, 25.);
231   fhEMCALDigits->SetXTitle("DigitsPerCluster");
232   fhEMCALDigits->Draw();
233  
234   cEMCAL->cd(3) ; 
235   gPad->SetLogy();
236   fhEMCALRecParticles->SetAxisRange(0, 25.);
237   fhEMCALRecParticles->SetXTitle("RecParticles");
238   fhEMCALRecParticles->Draw();
239  
240   cEMCAL->cd(4) ; 
241   gPad->SetLogy();
242   fhEMCALPhotons->SetAxisRange(0, 25.);
243   fhEMCALPhotons->SetXTitle("Photons");
244   fhEMCALPhotons->Draw();
245  
246   cEMCAL->cd(5) ; 
247   fhEMCALInvariantMass->SetXTitle("InvariantMass (MeV/c²)");
248   fhEMCALInvariantMass->Draw();
249  
250   cEMCAL->cd(6) ; 
251   gPad->SetLogy();
252   fhEMCALDigitsEvent->SetAxisRange(0, 40.);
253   fhEMCALDigitsEvent->SetXTitle("DigitsPerEvent");
254   fhEMCALDigitsEvent->Draw();
255  
256   cEMCAL->Print("EMCAL.eps");
257  
258   char line[1024] ; 
259   sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ; 
260   gROOT->ProcessLine(line);
261   sprintf(line, ".!rm -fR *.eps"); 
262   gROOT->ProcessLine(line);
263  
264   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
265 }