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