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