]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliPHOSQATask.cxx
add aliroot macros to look at data from strip modules and from LED reference system
[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 // An analysis task to check the PHOS photon data in simulated data
21 // An analysis task to check the PHOS photon data in simulated data
22 // An analysis task to check the PHOS photon data in simulated data
23 // An analysis task to check the PHOS photon data in simulated data
24 // An analysis task to check the PHOS photon data in simulated data
25 //
26 //*-- Yves Schutz 
27 //////////////////////////////////////////////////////////////////////////////
28
29 #include <TCanvas.h>
30 #include <TChain.h>
31 #include <TFile.h> 
32 #include <TH1.h>
33 #include <TNtuple.h>
34 #include <TROOT.h> 
35 #include <TVector3.h> 
36 #include <TString.h> 
37
38 #include "AliPHOSQATask.h" 
39 #include "AliESD.h" 
40 #include "AliLog.h"
41
42 //______________________________________________________________________________
43 AliPHOSQATask::AliPHOSQATask(const char *name) : 
44   AliAnalysisTask(name,""),  
45   fChain(0),
46   fESD(0), 
47   fOutputContainer(0),
48   fhPHOSPos(0),
49   fhPHOS(0),
50   fhPHOSEnergy(0),
51   fhPHOSDigits(0),
52   fhPHOSRecParticles(0),
53   fhPHOSPhotons(0),
54   fhPHOSInvariantMass(0),
55   fhPHOSDigitsEvent(0)
56 {
57   // Constructor.
58   // Input slot #0 works with an Ntuple
59   DefineInput(0, TChain::Class());
60   // Output slot #0 writes into a TH1 container
61   DefineOutput(0,  TObjArray::Class()) ; 
62 }
63
64 //____________________________________________________________________________
65 AliPHOSQATask::AliPHOSQATask(const AliPHOSQATask& ta) :
66   AliAnalysisTask(ta.GetName(),""),  
67   fChain(ta.fChain),
68   fESD(ta.fESD), 
69   fOutputContainer(ta.fOutputContainer), 
70   fhPHOSPos(ta.fhPHOSPos),
71   fhPHOS(ta.fhPHOS),
72   fhPHOSEnergy(ta.fhPHOSEnergy),
73   fhPHOSDigits(ta.fhPHOSDigits),
74   fhPHOSRecParticles(ta.fhPHOSRecParticles),
75   fhPHOSPhotons(ta.fhPHOSPhotons),
76   fhPHOSInvariantMass(ta.fhPHOSInvariantMass),
77   fhPHOSDigitsEvent(ta.fhPHOSDigitsEvent)
78
79   // cpy ctor
80 }
81
82 //_____________________________________________________________________________
83 AliPHOSQATask& AliPHOSQATask::operator = (const AliPHOSQATask& ap)
84 {
85 // assignment operator
86
87   this->~AliPHOSQATask();
88   new(this) AliPHOSQATask(ap);
89   return *this;
90 }
91
92 //______________________________________________________________________________
93 AliPHOSQATask::~AliPHOSQATask()
94 {
95   // dtor
96   fOutputContainer->Clear() ; 
97   delete fOutputContainer ;
98
99   delete fhPHOSPos ;
100   delete fhPHOS ;
101   delete fhPHOSEnergy ;
102   delete fhPHOSDigits ;
103   delete fhPHOSRecParticles ;
104   delete fhPHOSPhotons ;
105   delete fhPHOSInvariantMass ;
106   delete fhPHOSDigitsEvent ;
107 }
108
109
110 //______________________________________________________________________________
111 void AliPHOSQATask::ConnectInputData(const Option_t*)
112 {
113   // Initialisation of branch container and histograms 
114     
115   AliInfo(Form("*** Initialization of %s", GetName())) ; 
116   
117   // Get input data
118   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
119   if (!fChain) {
120     AliError(Form("Input 0 for %s not found\n", GetName()));
121     return ;
122   }
123   
124   // One should first check if the branch address was taken by some other task
125   char ** address = (char **)GetBranchAddress(0, "ESD");
126   if (address) {
127     fESD = (AliESD*)(*address);
128   } else {
129     fESD = new AliESD();
130     SetBranchAddress(0, "ESD", &fESD);
131   }
132 }
133
134 //________________________________________________________________________
135 void AliPHOSQATask::CreateOutputObjects()
136 {  
137   // create histograms 
138   
139   OpenFile(0) ; 
140
141   fhPHOSPos            = new TNtuple("PHOSPos"         , "Position in PHOS"  , "x:y:z");
142   fhPHOS               = new TNtuple("PHOS"            , "PHOS"  , "event:digits:clusters:photons");
143   fhPHOSEnergy         = new TH1D("PHOSEnergy"         , "PHOSEnergy"        , 1000, 0., 10. ) ;
144   fhPHOSDigits         = new TH1I("PHOSDigitsCluster"  , "PHOSDigits"        , 20 , 0 , 20  ) ;
145   fhPHOSRecParticles   = new TH1D("PHOSRecParticles"   , "PHOSRecParticles" , 20 , 0., 20. ) ;
146   fhPHOSPhotons        = new TH1I("PHOSPhotons"        , "PHOSPhotons"       , 20 , 0 , 20  ) ;
147   fhPHOSInvariantMass  = new TH1D("PHOSInvariantMass"  , "PHOSInvariantMass" , 400, 0., 400.) ;
148   fhPHOSDigitsEvent    = new TH1I("PHOSDigitsEvent"    , "PHOSDigitsEvent"   , 30 , 0 , 30  ) ;
149   
150   // create output container
151   
152   fOutputContainer = new TObjArray(8) ; 
153   fOutputContainer->SetName(GetName()) ; 
154
155   fOutputContainer->AddAt(fhPHOSPos,             0) ; 
156   fOutputContainer->AddAt(fhPHOS,                1) ; 
157   fOutputContainer->AddAt(fhPHOSEnergy,          2) ; 
158   fOutputContainer->AddAt(fhPHOSDigits,          3) ; 
159   fOutputContainer->AddAt(fhPHOSRecParticles,    4) ; 
160   fOutputContainer->AddAt(fhPHOSPhotons,         5) ; 
161   fOutputContainer->AddAt(fhPHOSInvariantMass,   6) ; 
162   fOutputContainer->AddAt(fhPHOSDigitsEvent,     7) ; 
163 }
164
165 //______________________________________________________________________________
166 void AliPHOSQATask::Exec(Option_t *) 
167 {
168   // Processing of one event
169     
170   Long64_t entry = fChain->GetReadEntry() ;
171   
172   if (!fESD) {
173     AliError("fESD is not connected to the input!") ; 
174     return ; 
175   }
176   
177   if ( !((entry-1)%100) ) 
178     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
179   
180   //************************  PHOS *************************************
181       
182   Int_t       firstPhosCluster       = fESD->GetFirstPHOSCluster() ;
183   const Int_t kNumberOfPhosClusters   = fESD->GetNumberOfPHOSClusters() ;
184   
185   TVector3 ** phosVector       = new TVector3*[kNumberOfPhosClusters] ;
186   Float_t  * phosPhotonsEnergy = new Float_t[kNumberOfPhosClusters] ;
187   Int_t      phosCluster ; 
188   Int_t      numberOfPhotonsInPhos  = 0 ;
189   Int_t      numberOfDigitsInPhos   = 0 ;
190   
191   // loop over the PHOS Cluster
192   for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + kNumberOfPhosClusters ; phosCluster++) {
193     AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
194     if (caloCluster) {
195       Float_t pos[3] ;
196       caloCluster->GetPosition( pos ) ;
197       fhPHOSEnergy->Fill( caloCluster->E() ) ;
198       fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
199       fhPHOSDigits->Fill(entry, caloCluster->GetNumberOfDigits() ) ;
200       numberOfDigitsInPhos += caloCluster->GetNumberOfDigits() ;
201       Double_t * pid = caloCluster->GetPid() ;
202       if(pid[AliPID::kPhoton] > 0.9) {
203         phosVector[numberOfPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
204         phosPhotonsEnergy[numberOfPhotonsInPhos]=caloCluster->E() ;
205         numberOfPhotonsInPhos++;
206       }
207     }
208   } //PHOS clusters
209     
210   fhPHOSRecParticles->Fill(kNumberOfPhosClusters);
211   fhPHOSPhotons->Fill(numberOfPhotonsInPhos);
212   fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
213   fhPHOS->Fill(entry, numberOfDigitsInPhos, kNumberOfPhosClusters, numberOfPhotonsInPhos) ; 
214
215   // invariant Mass
216   if (numberOfPhotonsInPhos > 1 ) {
217     Int_t phosPhoton1, phosPhoton2 ; 
218     for(phosPhoton1 = 0 ; phosPhoton1 < numberOfPhotonsInPhos ; phosPhoton1++) {
219       for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < numberOfPhotonsInPhos ; phosPhoton2++) {      
220         Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
221                                         ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) ) 
222                                         );
223         fhPHOSInvariantMass->Fill(tempMass*1000.);
224       }
225     }    
226   }
227   
228   PostData(0, fOutputContainer);
229
230   delete [] phosVector ; 
231   delete [] phosPhotonsEnergy ; 
232   
233 }
234
235 //______________________________________________________________________________
236 void AliPHOSQATask::Terminate(Option_t *)
237 {
238   // Processing when the event loop is ended
239
240   fOutputContainer = (TObjArray*)GetOutputData(0);  
241   fhPHOSPos            = (TNtuple*)fOutputContainer->At(0);
242   fhPHOS               = (TNtuple*)fOutputContainer->At(1);
243   fhPHOSEnergy         = (TH1D*)fOutputContainer->At(2);
244   fhPHOSDigits         = (TH1I*)fOutputContainer->At(3);
245   fhPHOSRecParticles   = (TH1D*)fOutputContainer->At(4);
246   fhPHOSPhotons        = (TH1I*)fOutputContainer->At(5);
247   fhPHOSInvariantMass  = (TH1D*)fOutputContainer->At(6);
248   fhPHOSDigitsEvent    = (TH1I*)fOutputContainer->At(7);
249
250   Bool_t problem = kFALSE ; 
251   AliInfo(Form(" *** %s Report:", GetName())) ; 
252   printf("        PHOSEnergy Mean         : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(),         fhPHOSEnergy->GetRMS()         ) ;
253   printf("        PHOSDigits Mean         : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(),         fhPHOSDigits->GetRMS()         ) ;
254   printf("        PHOSRecParticles Mean   : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(),   fhPHOSRecParticles->GetRMS()   ) ;
255   printf("        PHOSPhotons Mean        : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(),        fhPHOSPhotons->GetRMS()        ) ;
256   printf("        PHOSInvariantMass Mean  : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(),  fhPHOSInvariantMass->GetRMS()  ) ;
257   printf("        PHOSDigitsEvent Mean    : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(),    fhPHOSDigitsEvent->GetRMS()    ) ;
258
259   TCanvas  * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
260   cPHOS->Divide(3, 2);
261
262   cPHOS->cd(1) ; 
263   if ( fhPHOSEnergy->GetMaximum() > 0. ) 
264     gPad->SetLogy();
265   fhPHOSEnergy->SetAxisRange(0, 25.);
266   fhPHOSEnergy->SetLineColor(2);
267   fhPHOSEnergy->Draw();
268
269   cPHOS->cd(2) ; 
270   fhPHOSDigits->SetAxisRange(0,25.);
271   fhPHOSDigits->SetLineColor(2);
272   fhPHOSDigits->Draw();
273
274   cPHOS->cd(3) ; 
275   if ( fhPHOSRecParticles->GetMaximum() > 0. ) 
276     gPad->SetLogy();
277   fhPHOSRecParticles->SetAxisRange(0, 25.);
278   fhPHOSRecParticles->SetLineColor(2);
279   fhPHOSRecParticles->Draw();
280
281   cPHOS->cd(4) ; 
282   if ( fhPHOSPhotons->GetMaximum() > 0. ) 
283     gPad->SetLogy();
284   fhPHOSPhotons->SetAxisRange(0,25.);
285   fhPHOSPhotons->SetLineColor(2);
286   fhPHOSPhotons->Draw();
287
288   cPHOS->cd(5) ; 
289   fhPHOSInvariantMass->SetLineColor(2);
290   fhPHOSInvariantMass->Draw();
291  
292   cPHOS->cd(6) ; 
293   if ( fhPHOSDigitsEvent->GetMaximum() > 0. ) 
294     gPad->SetLogy();
295   fhPHOSDigitsEvent->SetAxisRange(0,40.);
296   fhPHOSDigitsEvent->SetLineColor(2);
297   fhPHOSDigitsEvent->Draw();
298  
299   cPHOS->Print("PHOS.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 }