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