]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliEMCALQATask.cxx
New method to clone current raw-data event and create a single-event raw-reader....
[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 // An analysis task to check the EMCAL photon data in simulated data
21 // An analysis task to check the EMCAL photon data in simulated data
22 // An analysis task to check the EMCAL photon data in simulated data
23 // An analysis task to check the EMCAL photon data in simulated data
24 //
25 //*-- Yves Schutz 
26 //////////////////////////////////////////////////////////////////////////////
27
28 #include <TCanvas.h>
29 #include <TChain.h>
30 #include <TFile.h> 
31 #include <TH1.h>
32 #include <TNtuple.h>
33 #include <TROOT.h>
34 #include <TVector3.h> 
35 #include <TString.h> 
36
37 #include "AliEMCALQATask.h" 
38 #include "AliESD.h" 
39 #include "AliLog.h"
40
41 //______________________________________________________________________________
42 AliEMCALQATask::AliEMCALQATask(const char *name) : 
43   AliAnalysisTask(name,""),  
44   fChain(0),
45   fESD(0), 
46   fOutputContainer(0), 
47   fhEMCALPos(0),
48   fhEMCAL(0),
49   fhEMCALEnergy(0),
50   fhEMCALDigits(0),
51   fhEMCALRecParticles(0),
52   fhEMCALPhotons(0),
53   fhEMCALInvariantMass(0),
54   fhEMCALDigitsEvent(0)
55 {
56   // Constructor.
57   // Input slot #0 works with an Ntuple
58   DefineInput(0, TChain::Class());
59   // Output slot #0 writes into a TH1 container
60   DefineOutput(0,  TObjArray::Class()) ; 
61 }
62
63 //____________________________________________________________________________
64 AliEMCALQATask::AliEMCALQATask(const AliEMCALQATask& ta) :
65   AliAnalysisTask(ta.GetName(),""),  
66   fChain(ta.fChain),
67   fESD(ta.fESD), 
68   fOutputContainer(ta.fOutputContainer), 
69   fhEMCALPos(ta.fhEMCALPos),
70   fhEMCAL(ta.fhEMCAL),
71   fhEMCALEnergy(ta.fhEMCALEnergy),
72   fhEMCALDigits(ta.fhEMCALDigits),
73   fhEMCALRecParticles(ta.fhEMCALRecParticles),
74   fhEMCALPhotons(ta.fhEMCALPhotons),
75   fhEMCALInvariantMass(ta.fhEMCALInvariantMass),
76   fhEMCALDigitsEvent(ta.fhEMCALDigitsEvent)
77
78   // cpy ctor
79 }
80
81 //_____________________________________________________________________________
82 AliEMCALQATask& AliEMCALQATask::operator = (const AliEMCALQATask& ap)
83 {
84 // assignment operator
85
86   this->~AliEMCALQATask();
87   new(this) AliEMCALQATask(ap);
88   return *this;
89 }
90
91 //______________________________________________________________________________
92 AliEMCALQATask::~AliEMCALQATask()
93 {
94   // dtor
95   fOutputContainer->Clear() ; 
96   delete fOutputContainer ; 
97
98   delete fhEMCALPos ;
99   delete fhEMCAL ;
100   delete fhEMCALEnergy ;
101   delete fhEMCALDigits ;
102   delete fhEMCALRecParticles ;
103   delete fhEMCALPhotons ;
104   delete fhEMCALInvariantMass ;
105   delete fhEMCALDigitsEvent ;
106 }
107
108 //______________________________________________________________________________
109 void AliEMCALQATask::ConnectInputData(const Option_t*)
110 {
111   // Initialisation of branch container and histograms 
112     
113   AliInfo(Form("*** Initialization of %s", GetName())) ; 
114   
115   // Get input data
116   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
117   if (!fChain) {
118     AliError(Form("Input 0 for %s not found\n", GetName()));
119     return ;
120   }
121   
122   // One should first check if the branch address was taken by some other task
123   char ** address = (char **)GetBranchAddress(0, "ESD");
124   if (address) {
125     fESD = (AliESD*)(*address);
126   } else {
127     fESD = new AliESD();
128     SetBranchAddress(0, "ESD", &fESD);
129   }
130 }
131   
132 //______________________________________________________________________________
133 void AliEMCALQATask::CreateOutputObjects()  
134 {
135 // create histograms  
136
137   OpenFile(0) ; 
138
139   fhEMCALPos           = new TNtuple("EMCALPos"        , "Position in EMCAL" , "x:y:z");
140   fhEMCAL              = new TNtuple("EMCAL"           , "EMCAL" , "event:digits:clusters:photons");
141   fhEMCALEnergy        = new TH1D("EMCALEnergy"        , "EMCALEnergy"       , 1000, 0., 10. ) ;
142   fhEMCALDigits        = new TH1I("EMCALDigitsCluster" , "EMCALDigits"       , 20 , 0 , 20  ) ;
143   fhEMCALRecParticles  = new TH1D("EMCALRecParticles"  , "EMCALRecParticles", 20 , 0., 20. ) ;
144   fhEMCALPhotons       = new TH1I("EMCALPhotons"       , "EMCALPhotons"      , 20 , 0 , 20  ) ;
145   fhEMCALInvariantMass = new TH1D("EMCALInvariantMass" , "EMCALInvariantMass", 400, 0., 400.) ;
146   fhEMCALDigitsEvent   = new TH1I("EMCALDigitsEvent"   , "EMCALDigitsEvent"  , 30 , 0 , 30  ) ;
147   
148   // create output container
149   
150   fOutputContainer = new TObjArray(8) ; 
151   fOutputContainer->SetName(GetName()) ; 
152
153   fOutputContainer->AddAt(fhEMCALPos,            0) ; 
154   fOutputContainer->AddAt(fhEMCAL,               1) ; 
155   fOutputContainer->AddAt(fhEMCALEnergy,         2) ; 
156   fOutputContainer->AddAt(fhEMCALDigits,         3) ; 
157   fOutputContainer->AddAt(fhEMCALRecParticles,   4) ; 
158   fOutputContainer->AddAt(fhEMCALPhotons,        5) ; 
159   fOutputContainer->AddAt(fhEMCALInvariantMass,  6) ; 
160   fOutputContainer->AddAt(fhEMCALDigitsEvent,    7) ; 
161 }
162
163 //______________________________________________________________________________
164 void AliEMCALQATask::Exec(Option_t *) 
165 {
166   // Processing of one event
167  
168   Long64_t entry = fChain->GetReadEntry() ;
169  
170   if (!fESD) {
171     AliError("fESD is not connected to the input!") ; 
172     return ; 
173   }
174   
175   if ( !((entry-1)%100) ) 
176     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
177   
178   //************************  EMCAL *************************************
179   Int_t       firstEmcalCluster       = fESD->GetFirstEMCALCluster() ;
180   const Int_t kNumberOfEmcalClusters   = fESD->GetNumberOfEMCALClusters() ;
181
182   TVector3 ** emcalVector        = new TVector3*[kNumberOfEmcalClusters] ;
183   Float_t   * emcalPhotonsEnergy = new Float_t[kNumberOfEmcalClusters] ;
184   Int_t      emcalCluster ; 
185   Int_t      numberOfEmcalClustersv1 = 0 ; 
186   Int_t      numberOfPhotonsInEmcal  = 0 ;
187   Int_t      numberOfDigitsInEmcal   = 0 ;  
188
189
190   // loop over all the EMCAL Cluster
191   for(emcalCluster = firstEmcalCluster ; emcalCluster < firstEmcalCluster + kNumberOfEmcalClusters ; emcalCluster++) {
192     AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(emcalCluster) ;
193     if (caloCluster) {
194       Float_t pos[3] ;
195       if(caloCluster->GetClusterType() == AliESDCaloCluster::kEMCALClusterv1) {  
196         caloCluster->GetPosition(pos) ;
197         fhEMCALPos->Fill(pos[0],pos[1],pos[2]) ;
198         fhEMCALEnergy->Fill(caloCluster->E()) ;
199         fhEMCALDigits->Fill(entry, caloCluster->GetNumberOfDigits()) ;
200         numberOfEmcalClustersv1++ ;
201         numberOfDigitsInEmcal += caloCluster->GetNumberOfDigits() ;    
202         // Float_t * pid = clus->GetPid() ;
203         // if(pid[AliPID::kPhoton]>0.9){
204         emcalVector[numberOfPhotonsInEmcal] = new TVector3(pos[0],pos[1],pos[2]) ;
205         emcalPhotonsEnergy[numberOfPhotonsInEmcal] = caloCluster->E() ;
206         numberOfPhotonsInEmcal++ ; 
207       }
208     }
209   } // EMCAL clusters loop
210
211   fhEMCALRecParticles->Fill(numberOfEmcalClustersv1);
212   fhEMCALPhotons->Fill(numberOfPhotonsInEmcal);
213   fhEMCALDigitsEvent->Fill(numberOfDigitsInEmcal);
214   fhEMCAL->Fill(entry, numberOfDigitsInEmcal, numberOfEmcalClustersv1, numberOfPhotonsInEmcal) ; 
215   
216   // invariant Mass
217   if (numberOfPhotonsInEmcal > 1 ) {
218     Int_t emcalPhoton1, emcalPhoton2 ; 
219     for(emcalPhoton1 = 0 ; emcalPhoton1 < numberOfPhotonsInEmcal ; emcalPhoton1++) {
220       for(emcalPhoton2 = emcalPhoton1 + 1 ; emcalPhoton2 < numberOfPhotonsInEmcal ; emcalPhoton2++) {
221         Float_t tempMass = TMath::Sqrt( 2 * emcalPhotonsEnergy[emcalPhoton1] * emcalPhotonsEnergy[emcalPhoton2] * 
222                                         ( 1 - TMath::Cos( emcalVector[emcalPhoton1]->Angle(*emcalVector[emcalPhoton2])) 
223                                           )
224                                         );
225         fhEMCALInvariantMass->Fill(tempMass*1000.);
226       }
227     }    
228   }
229   
230   PostData(0, fOutputContainer);
231   
232   delete [] emcalVector ; 
233   delete [] emcalPhotonsEnergy ;
234 }
235
236 //______________________________________________________________________________
237 void AliEMCALQATask::Terminate(Option_t *)
238 {
239   // Processing when the event loop is ended
240   fOutputContainer = (TObjArray*)GetOutputData(0);
241   fhEMCALEnergy = (TH1D*)fOutputContainer->At(2);
242   fhEMCALDigits = (TH1I*)fOutputContainer->At(3);
243   fhEMCALRecParticles = (TH1D*)fOutputContainer->At(4);
244   fhEMCALPhotons = (TH1I*)fOutputContainer->At(5);
245   fhEMCALInvariantMass = (TH1D*)fOutputContainer->At(6);
246   fhEMCALDigitsEvent = (TH1I*)fOutputContainer->At(7);
247  
248   Bool_t problem = kFALSE ; 
249   AliInfo(Form(" *** %s Report:", GetName())) ; 
250   printf("        EMCALEnergy Mean        : %5.3f , RMS : %5.3f \n", fhEMCALEnergy->GetMean(),        fhEMCALEnergy->GetRMS()        ) ;
251   printf("        EMCALDigits Mean        : %5.3f , RMS : %5.3f \n", fhEMCALDigits->GetMean(),        fhEMCALDigits->GetRMS()        ) ;
252   printf("        EMCALRecParticles Mean  : %5.3f , RMS : %5.3f \n", fhEMCALRecParticles->GetMean(),  fhEMCALRecParticles->GetRMS()  ) ;
253   printf("        EMCALPhotons Mean       : %5.3f , RMS : %5.3f \n", fhEMCALPhotons->GetMean(),       fhEMCALPhotons->GetRMS()       ) ;
254   printf("        EMCALInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhEMCALInvariantMass->GetMean(), fhEMCALInvariantMass->GetRMS() ) ;
255   printf("        EMCALDigitsEvent Mean   : %5.3f , RMS : %5.3f \n", fhEMCALDigitsEvent->GetMean(),   fhEMCALDigitsEvent->GetRMS()   ) ;
256
257   TCanvas  * cEMCAL = new TCanvas("EMCAL", "EMCAL ESD Test", 400, 10, 600, 700);
258   cEMCAL->Divide(3, 2) ; 
259
260   cEMCAL->cd(1) ; 
261   if ( fhEMCALEnergy->GetMaximum() > 0. ) 
262     gPad->SetLogy();
263   fhEMCALEnergy->SetAxisRange(0, 25.);
264   fhEMCALEnergy->SetXTitle("Energy (GeV)");
265   fhEMCALEnergy->Draw();
266   
267   cEMCAL->cd(2) ; 
268   if ( fhEMCALDigits->GetMaximum() > 0. ) 
269     gPad->SetLogy();
270   fhEMCALDigits->SetAxisRange(0, 25.);
271   fhEMCALDigits->SetXTitle("DigitsPerCluster");
272   fhEMCALDigits->Draw();
273  
274   cEMCAL->cd(3) ; 
275   if ( fhEMCALRecParticles->GetMaximum() > 0. ) 
276      gPad->SetLogy();
277   fhEMCALRecParticles->SetAxisRange(0, 25.);
278   fhEMCALRecParticles->SetXTitle("RecParticles");
279   fhEMCALRecParticles->Draw();
280  
281   cEMCAL->cd(4) ; 
282   if ( fhEMCALPhotons->GetMaximum() > 0. ) 
283     gPad->SetLogy();
284   fhEMCALPhotons->SetAxisRange(0, 25.);
285   fhEMCALPhotons->SetXTitle("Photons");
286   fhEMCALPhotons->Draw();
287  
288   cEMCAL->cd(5) ; 
289   fhEMCALInvariantMass->SetXTitle("InvariantMass (MeV/c²)");
290   fhEMCALInvariantMass->Draw();
291  
292   cEMCAL->cd(6) ; 
293   if ( fhEMCALDigitsEvent->GetMaximum() > 0. ) 
294     gPad->SetLogy();
295   fhEMCALDigitsEvent->SetAxisRange(0, 40.);
296   fhEMCALDigitsEvent->SetXTitle("DigitsPerEvent");
297   fhEMCALDigitsEvent->Draw();
298  
299   cEMCAL->Print("EMCAL.eps");
300  
301   char line[1024] ; 
302   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
303   gROOT->ProcessLine(line);
304   sprintf(line, ".!rm -fR *.eps"); 
305   gROOT->ProcessLine(line);
306  
307   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
308
309   TString report ; 
310   if(problem)
311     report="Problems found, please check!!!";  
312   else 
313     report="OK";
314
315   AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ; 
316 }