]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSESDQA.cxx
effc++ warnings
[u/mrichter/AliRoot.git] / PHOS / AliPHOSESDQA.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 PHOS 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 "AliPHOSESDQA.h" 
37 #include "AliESD.h" 
38 #include "AliLog.h"
39
40 //______________________________________________________________________________
41 AliPHOSESDQA::AliPHOSESDQA(const char *name) : 
42   AliAnalysisTask(name,""),  
43   fChain(0),
44   fESD(0), 
45   fhPHOS(0),
46   fhPHOSEnergy(0),
47   fhPHOSDigits(0),
48   fhPHOSRecParticles(0),
49   fhPHOSPhotons(0),
50   fhPHOSInvariantMass(0),
51   fhPHOSDigitsEvent(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 AliPHOSESDQA::~AliPHOSESDQA()
62 {
63   // dtor
64   fOutputContainer->Clear() ; 
65   delete fOutputContainer ;
66
67   delete fhPHOSPos ;
68   delete fhPHOS ;
69   delete fhPHOSEnergy ;
70   delete fhPHOSDigits ;
71   delete fhPHOSRecParticles ;
72   delete fhPHOSPhotons ;
73   delete fhPHOSInvariantMass ;
74   delete fhPHOSDigitsEvent ;
75 }
76
77
78 //______________________________________________________________________________
79 void AliPHOSESDQA::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 AliPHOSESDQA::CreateOutputObjects()
104 {  
105   // create histograms 
106   
107   OpenFile(0) ; 
108
109   fhPHOSPos            = new TNtuple("PHOSPos"         , "Position in PHOS"  , "x:y:z");
110   fhPHOS               = new TNtuple("PHOS"            , "PHOS"  , "event:digits:clusters:photons");
111   fhPHOSEnergy         = new TH1D("PHOSEnergy"         , "PHOSEnergy"        , 1000, 0., 10. ) ;
112   fhPHOSDigits         = new TH1I("PHOSDigitsCluster"  , "PHOSDigits"        , 20 , 0 , 20  ) ;
113   fhPHOSRecParticles   = new TH1D("PHOSRecParticles"   , "PHOSRecParticles" , 20 , 0., 20. ) ;
114   fhPHOSPhotons        = new TH1I("PHOSPhotons"        , "PHOSPhotons"       , 20 , 0 , 20  ) ;
115   fhPHOSInvariantMass  = new TH1D("PHOSInvariantMass"  , "PHOSInvariantMass" , 400, 0., 400.) ;
116   fhPHOSDigitsEvent    = new TH1I("PHOSDigitsEvent"    , "PHOSDigitsEvent"   , 30 , 0 , 30  ) ;
117   
118   // create output container
119   
120   fOutputContainer = new TObjArray(8) ; 
121   fOutputContainer->SetName(GetName()) ; 
122
123   fOutputContainer->AddAt(fhPHOSPos,             0) ; 
124   fOutputContainer->AddAt(fhPHOS,                1) ; 
125   fOutputContainer->AddAt(fhPHOSEnergy,          2) ; 
126   fOutputContainer->AddAt(fhPHOSDigits,          3) ; 
127   fOutputContainer->AddAt(fhPHOSRecParticles,    4) ; 
128   fOutputContainer->AddAt(fhPHOSPhotons,         5) ; 
129   fOutputContainer->AddAt(fhPHOSInvariantMass,   6) ; 
130   fOutputContainer->AddAt(fhPHOSDigitsEvent,     7) ; 
131 }
132
133 //______________________________________________________________________________
134 void AliPHOSESDQA::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   //************************  PHOS *************************************
149       
150   Int_t       firstPhosCluster       = fESD->GetFirstPHOSCluster() ;
151   const Int_t numberOfPhosClusters   = fESD->GetNumberOfPHOSClusters() ;
152   
153   TVector3 ** phosVector       = new TVector3*[numberOfPhosClusters] ;
154   Float_t  * phosPhotonsEnergy = new Float_t[numberOfPhosClusters] ;
155   Int_t      phosCluster ; 
156   Int_t      numberOfPhotonsInPhos  = 0 ;
157   Int_t      numberOfDigitsInPhos   = 0 ;
158   
159   // loop over the PHOS Cluster
160   for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + numberOfPhosClusters ; phosCluster++) {
161     AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
162     if (caloCluster) {
163       Float_t pos[3] ;
164       caloCluster->GetPosition( pos ) ;
165       fhPHOSEnergy->Fill( caloCluster->E() ) ;
166       fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
167       fhPHOSDigits->Fill(entry, caloCluster->GetNumberOfDigits() ) ;
168       numberOfDigitsInPhos += caloCluster->GetNumberOfDigits() ;
169       Float_t * pid = caloCluster->GetPid() ;
170       if(pid[AliPID::kPhoton] > 0.9) {
171         phosVector[numberOfPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
172         phosPhotonsEnergy[numberOfPhotonsInPhos]=caloCluster->E() ;
173         numberOfPhotonsInPhos++;
174       }
175     }
176   } //PHOS clusters
177     
178   fhPHOSRecParticles->Fill(numberOfPhosClusters);
179   fhPHOSPhotons->Fill(numberOfPhotonsInPhos);
180   fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
181   fhPHOS->Fill(entry, numberOfDigitsInPhos, numberOfPhosClusters, numberOfPhotonsInPhos) ; 
182
183   // invariant Mass
184   if (numberOfPhotonsInPhos > 1 ) {
185     Int_t phosPhoton1, phosPhoton2 ; 
186     for(phosPhoton1 = 0 ; phosPhoton1 < numberOfPhotonsInPhos ; phosPhoton1++) {
187       for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < numberOfPhotonsInPhos ; phosPhoton2++) {      
188         Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
189                                         ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) ) 
190                                         );
191         fhPHOSInvariantMass->Fill(tempMass*1000.);
192       }
193     }    
194   }
195   
196   PostData(0, fOutputContainer);
197
198   delete [] phosVector ; 
199   delete [] phosPhotonsEnergy ; 
200   
201 }
202
203 //______________________________________________________________________________
204 void AliPHOSESDQA::Terminate(Option_t *)
205 {
206   // Processing when the event loop is ended
207
208   fOutputContainer = (TObjArray*)GetOutputData(0);  
209   fhPHOSPos            = (TNtuple*)fOutputContainer->At(0);
210   fhPHOS               = (TNtuple*)fOutputContainer->At(1);
211   fhPHOSEnergy         = (TH1D*)fOutputContainer->At(2);
212   fhPHOSDigits         = (TH1I*)fOutputContainer->At(3);
213   fhPHOSRecParticles   = (TH1D*)fOutputContainer->At(4);
214   fhPHOSPhotons        = (TH1I*)fOutputContainer->At(5);
215   fhPHOSInvariantMass  = (TH1D*)fOutputContainer->At(6);
216   fhPHOSDigitsEvent    = (TH1I*)fOutputContainer->At(7);
217
218   Bool_t problem = kFALSE ; 
219   AliInfo(Form(" *** %s Report:", GetName())) ; 
220   printf("        PHOSEnergy Mean         : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(),         fhPHOSEnergy->GetRMS()         ) ;
221   printf("        PHOSDigits Mean         : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(),         fhPHOSDigits->GetRMS()         ) ;
222   printf("        PHOSRecParticles Mean   : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(),   fhPHOSRecParticles->GetRMS()   ) ;
223   printf("        PHOSPhotons Mean        : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(),        fhPHOSPhotons->GetRMS()        ) ;
224   printf("        PHOSInvariantMass Mean  : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(),  fhPHOSInvariantMass->GetRMS()  ) ;
225   printf("        PHOSDigitsEvent Mean    : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(),    fhPHOSDigitsEvent->GetRMS()    ) ;
226
227   TCanvas  * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
228   cPHOS->Divide(3, 2);
229
230   cPHOS->cd(1) ; 
231   if ( fhPHOSEnergy->GetMaximum() > 0. ) 
232     gPad->SetLogy();
233   fhPHOSEnergy->SetAxisRange(0, 25.);
234   fhPHOSEnergy->SetLineColor(2);
235   fhPHOSEnergy->Draw();
236
237   cPHOS->cd(2) ; 
238   fhPHOSDigits->SetAxisRange(0,25.);
239   fhPHOSDigits->SetLineColor(2);
240   fhPHOSDigits->Draw();
241
242   cPHOS->cd(3) ; 
243   if ( fhPHOSRecParticles->GetMaximum() > 0. ) 
244     gPad->SetLogy();
245   fhPHOSRecParticles->SetAxisRange(0, 25.);
246   fhPHOSRecParticles->SetLineColor(2);
247   fhPHOSRecParticles->Draw();
248
249   cPHOS->cd(4) ; 
250   if ( fhPHOSPhotons->GetMaximum() > 0. ) 
251     gPad->SetLogy();
252   fhPHOSPhotons->SetAxisRange(0,25.);
253   fhPHOSPhotons->SetLineColor(2);
254   fhPHOSPhotons->Draw();
255
256   cPHOS->cd(5) ; 
257   fhPHOSInvariantMass->SetLineColor(2);
258   fhPHOSInvariantMass->Draw();
259  
260   cPHOS->cd(6) ; 
261   if ( fhPHOSDigitsEvent->GetMaximum() > 0. ) 
262     gPad->SetLogy();
263   fhPHOSDigitsEvent->SetAxisRange(0,40.);
264   fhPHOSDigitsEvent->SetLineColor(2);
265   fhPHOSDigitsEvent->Draw();
266  
267   cPHOS->Print("PHOS.eps");
268  
269   char line[1024] ; 
270   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
271   gROOT->ProcessLine(line);
272   sprintf(line, ".!rm -fR *.eps"); 
273   gROOT->ProcessLine(line);
274  
275   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
276
277   TString report ; 
278   if(problem)
279     report="Problems found, please check!!!";  
280   else 
281     report="OK";
282
283   AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ; 
284 }