]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliEMCALQATask.cxx
Bug fix. The equippment ID was not given to alirawreader
[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::ConnectInputData(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   // One should first check if the branch address was taken by some other task
92   char ** address = (char **)GetBranchAddress(0, "ESD");
93   if (address) {
94     fESD = (AliESD*)(*address);
95   } else {
96     fESD = new AliESD();
97     SetBranchAddress(0, "ESD", &fESD);
98   }
99 }
100   
101 //______________________________________________________________________________
102 void AliEMCALQATask::CreateOutputObjects()  
103 {
104 // create histograms  
105   fhEMCALPos           = new TNtuple("EMCALPos"        , "Position in EMCAL" , "x:y:z");
106   fhEMCAL              = new TNtuple("EMCAL"           , "EMCAL" , "event:digits:clusters:photons");
107   fhEMCALEnergy        = new TH1D("EMCALEnergy"        , "EMCALEnergy"       , 1000, 0., 10. ) ;
108   fhEMCALDigits        = new TH1I("EMCALDigitsCluster" , "EMCALDigits"       , 20 , 0 , 20  ) ;
109   fhEMCALRecParticles  = new TH1D("EMCALRecParticles"  , "EMCALRecParticles", 20 , 0., 20. ) ;
110   fhEMCALPhotons       = new TH1I("EMCALPhotons"       , "EMCALPhotons"      , 20 , 0 , 20  ) ;
111   fhEMCALInvariantMass = new TH1D("EMCALInvariantMass" , "EMCALInvariantMass", 400, 0., 400.) ;
112   fhEMCALDigitsEvent   = new TH1I("EMCALDigitsEvent"   , "EMCALDigitsEvent"  , 30 , 0 , 30  ) ;
113   
114   // create output container
115   
116   fOutputContainer = new TObjArray(8) ; 
117   fOutputContainer->SetName(GetName()) ; 
118
119   fOutputContainer->AddAt(fhEMCALPos,            0) ; 
120   fOutputContainer->AddAt(fhEMCAL,               1) ; 
121   fOutputContainer->AddAt(fhEMCALEnergy,         2) ; 
122   fOutputContainer->AddAt(fhEMCALDigits,         3) ; 
123   fOutputContainer->AddAt(fhEMCALRecParticles,   4) ; 
124   fOutputContainer->AddAt(fhEMCALPhotons,        5) ; 
125   fOutputContainer->AddAt(fhEMCALInvariantMass,  6) ; 
126   fOutputContainer->AddAt(fhEMCALDigitsEvent,    7) ; 
127 }
128
129 //______________________________________________________________________________
130 void AliEMCALQATask::Exec(Option_t *) 
131 {
132   // Processing of one event
133  
134   Long64_t entry = fChain->GetReadEntry() ;
135  
136   if (!fESD) {
137     AliError("fESD is not connected to the input!") ; 
138     return ; 
139   }
140   
141   if ( !((entry-1)%100) ) 
142     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
143   
144   //************************  EMCAL *************************************
145   Int_t       firstEmcalCluster       = fESD->GetFirstEMCALCluster() ;
146   const Int_t mumberOfEmcalClusters   = fESD->GetNumberOfEMCALClusters() ;
147
148   TVector3 ** emcalVector        = new TVector3*[mumberOfEmcalClusters] ;
149   Float_t   * emcalPhotonsEnergy = new Float_t[mumberOfEmcalClusters] ;
150   Int_t      emcalCluster ; 
151   Int_t      numberOfEmcalClustersv1 = 0 ; 
152   Int_t      numberOfPhotonsInEmcal  = 0 ;
153   Int_t      numberOfDigitsInEmcal   = 0 ;  
154
155
156   // loop over all the EMCAL Cluster
157   for(emcalCluster = firstEmcalCluster ; emcalCluster < firstEmcalCluster + mumberOfEmcalClusters ; emcalCluster++) {
158     AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(emcalCluster) ;
159     if (caloCluster) {
160       Float_t pos[3] ;
161       if(caloCluster->GetClusterType() == AliESDCaloCluster::kClusterv1) {  
162         caloCluster->GetGlobalPosition(pos) ;
163         fhEMCALPos->Fill(pos[0],pos[1],pos[2]) ;
164         fhEMCALEnergy->Fill(caloCluster->GetClusterEnergy()) ;
165         fhEMCALDigits->Fill(entry, caloCluster->GetNumberOfDigits()) ;
166         numberOfEmcalClustersv1++ ;
167         numberOfDigitsInEmcal += caloCluster->GetNumberOfDigits() ;    
168         // Float_t * pid = clus->GetPid() ;
169         // if(pid[AliPID::kPhoton]>0.9){
170         emcalVector[numberOfPhotonsInEmcal] = new TVector3(pos[0],pos[1],pos[2]) ;
171         emcalPhotonsEnergy[numberOfPhotonsInEmcal] = caloCluster->GetClusterEnergy() ;
172         numberOfPhotonsInEmcal++ ; 
173       }
174     }
175   } // EMCAL clusters loop
176
177   fhEMCALRecParticles->Fill(numberOfEmcalClustersv1);
178   fhEMCALPhotons->Fill(numberOfPhotonsInEmcal);
179   fhEMCALDigitsEvent->Fill(numberOfDigitsInEmcal);
180   fhEMCAL->Fill(entry, numberOfDigitsInEmcal, numberOfEmcalClustersv1, numberOfPhotonsInEmcal) ; 
181   
182   // invariant Mass
183   if (numberOfPhotonsInEmcal > 1 ) {
184     Int_t emcalPhoton1, emcalPhoton2 ; 
185     for(emcalPhoton1 = 0 ; emcalPhoton1 < numberOfPhotonsInEmcal ; emcalPhoton1++) {
186       for(emcalPhoton2 = emcalPhoton1 + 1 ; emcalPhoton2 < numberOfPhotonsInEmcal ; emcalPhoton2++) {
187         Float_t tempMass = TMath::Sqrt( 2 * emcalPhotonsEnergy[emcalPhoton1] * emcalPhotonsEnergy[emcalPhoton2] * 
188                                         ( 1 - TMath::Cos( emcalVector[emcalPhoton1]->Angle(*emcalVector[emcalPhoton2])) 
189                                           )
190                                         );
191         fhEMCALInvariantMass->Fill(tempMass*1000.);
192       }
193     }    
194   }
195   
196   PostData(0, fOutputContainer);
197   
198   delete [] emcalVector ; 
199   delete [] emcalPhotonsEnergy ;
200 }
201
202 //______________________________________________________________________________
203 void AliEMCALQATask::Terminate(Option_t *)
204 {
205   // Processing when the event loop is ended
206   fOutputContainer = (TObjArray*)GetOutputData(0);
207   fhEMCALEnergy = (TH1D*)fOutputContainer->At(2);
208   fhEMCALDigits = (TH1I*)fOutputContainer->At(3);
209   fhEMCALRecParticles = (TH1D*)fOutputContainer->At(4);
210   fhEMCALPhotons = (TH1I*)fOutputContainer->At(5);
211   fhEMCALInvariantMass = (TH1D*)fOutputContainer->At(6);
212   fhEMCALDigitsEvent = (TH1I*)fOutputContainer->At(7);
213  
214   Bool_t problem = kFALSE ; 
215   AliInfo(Form(" *** %s Report:", GetName())) ; 
216   printf("        EMCALEnergy Mean        : %5.3f , RMS : %5.3f \n", fhEMCALEnergy->GetMean(),        fhEMCALEnergy->GetRMS()        ) ;
217   printf("        EMCALDigits Mean        : %5.3f , RMS : %5.3f \n", fhEMCALDigits->GetMean(),        fhEMCALDigits->GetRMS()        ) ;
218   printf("        EMCALRecParticles Mean  : %5.3f , RMS : %5.3f \n", fhEMCALRecParticles->GetMean(),  fhEMCALRecParticles->GetRMS()  ) ;
219   printf("        EMCALPhotons Mean       : %5.3f , RMS : %5.3f \n", fhEMCALPhotons->GetMean(),       fhEMCALPhotons->GetRMS()       ) ;
220   printf("        EMCALInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhEMCALInvariantMass->GetMean(), fhEMCALInvariantMass->GetRMS() ) ;
221   printf("        EMCALDigitsEvent Mean   : %5.3f , RMS : %5.3f \n", fhEMCALDigitsEvent->GetMean(),   fhEMCALDigitsEvent->GetRMS()   ) ;
222
223   TCanvas  * cEMCAL = new TCanvas("EMCAL", "EMCAL ESD Test", 400, 10, 600, 700);
224   cEMCAL->Divide(3, 2) ; 
225
226   cEMCAL->cd(1) ; 
227   if ( fhEMCALEnergy->GetMaximum() > 0. ) 
228     gPad->SetLogy();
229   fhEMCALEnergy->SetAxisRange(0, 25.);
230   fhEMCALEnergy->SetXTitle("Energy (GeV)");
231   fhEMCALEnergy->Draw();
232   
233   cEMCAL->cd(2) ; 
234   if ( fhEMCALDigits->GetMaximum() > 0. ) 
235     gPad->SetLogy();
236   fhEMCALDigits->SetAxisRange(0, 25.);
237   fhEMCALDigits->SetXTitle("DigitsPerCluster");
238   fhEMCALDigits->Draw();
239  
240   cEMCAL->cd(3) ; 
241   if ( fhEMCALRecParticles->GetMaximum() > 0. ) 
242      gPad->SetLogy();
243   fhEMCALRecParticles->SetAxisRange(0, 25.);
244   fhEMCALRecParticles->SetXTitle("RecParticles");
245   fhEMCALRecParticles->Draw();
246  
247   cEMCAL->cd(4) ; 
248   if ( fhEMCALPhotons->GetMaximum() > 0. ) 
249     gPad->SetLogy();
250   fhEMCALPhotons->SetAxisRange(0, 25.);
251   fhEMCALPhotons->SetXTitle("Photons");
252   fhEMCALPhotons->Draw();
253  
254   cEMCAL->cd(5) ; 
255   fhEMCALInvariantMass->SetXTitle("InvariantMass (MeV/c²)");
256   fhEMCALInvariantMass->Draw();
257  
258   cEMCAL->cd(6) ; 
259   if ( fhEMCALDigitsEvent->GetMaximum() > 0. ) 
260     gPad->SetLogy();
261   fhEMCALDigitsEvent->SetAxisRange(0, 40.);
262   fhEMCALDigitsEvent->SetXTitle("DigitsPerEvent");
263   fhEMCALDigitsEvent->Draw();
264  
265   cEMCAL->Print("EMCAL.eps");
266  
267   char line[1024] ; 
268   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
269   gROOT->ProcessLine(line);
270   sprintf(line, ".!rm -fR *.eps"); 
271   gROOT->ProcessLine(line);
272  
273   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
274
275   char * report ; 
276   if(problem)
277     report="Problems found, please check!!!";  
278   else 
279     report="OK";
280
281   AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report)) ; 
282 }