1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
19 // An analysis task to check the PHOS photon data in simulated data
22 //////////////////////////////////////////////////////////////////////////////
36 #include "AliPHOSESDQA.h"
40 //______________________________________________________________________________
41 AliPHOSESDQA::AliPHOSESDQA(const char *name) :
42 AliAnalysisTask(name,""),
48 fhPHOSRecParticles(0),
50 fhPHOSInvariantMass(0),
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()) ;
60 //______________________________________________________________________________
61 AliPHOSESDQA::~AliPHOSESDQA()
64 fOutputContainer->Clear() ;
65 delete fOutputContainer ;
71 delete fhPHOSRecParticles ;
72 delete fhPHOSPhotons ;
73 delete fhPHOSInvariantMass ;
74 delete fhPHOSDigitsEvent ;
78 //______________________________________________________________________________
79 void AliPHOSESDQA::ConnectInputData(const Option_t*)
81 // Initialisation of branch container and histograms
83 AliInfo(Form("*** Initialization of %s", GetName())) ;
86 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
88 AliError(Form("Input 0 for %s not found\n", GetName()));
92 // One should first check if the branch address was taken by some other task
93 char ** address = (char **)GetBranchAddress(0, "ESD");
95 fESD = (AliESD*)(*address);
98 SetBranchAddress(0, "ESD", &fESD);
102 //________________________________________________________________________
103 void AliPHOSESDQA::CreateOutputObjects()
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 ) ;
118 // create output container
120 fOutputContainer = new TObjArray(8) ;
121 fOutputContainer->SetName(GetName()) ;
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) ;
133 //______________________________________________________________________________
134 void AliPHOSESDQA::Exec(Option_t *)
136 // Processing of one event
138 Long64_t entry = fChain->GetReadEntry() ;
141 AliError("fESD is not connected to the input!") ;
145 if ( !((entry-1)%100) )
146 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
148 //************************ PHOS *************************************
150 Int_t firstPhosCluster = fESD->GetFirstPHOSCluster() ;
151 const Int_t numberOfPhosClusters = fESD->GetNumberOfPHOSClusters() ;
153 TVector3 ** phosVector = new TVector3*[numberOfPhosClusters] ;
154 Float_t * phosPhotonsEnergy = new Float_t[numberOfPhosClusters] ;
156 Int_t numberOfPhotonsInPhos = 0 ;
157 Int_t numberOfDigitsInPhos = 0 ;
159 // loop over the PHOS Cluster
160 for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + numberOfPhosClusters ; phosCluster++) {
161 AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
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++;
178 fhPHOSRecParticles->Fill(numberOfPhosClusters);
179 fhPHOSPhotons->Fill(numberOfPhotonsInPhos);
180 fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
181 fhPHOS->Fill(entry, numberOfDigitsInPhos, numberOfPhosClusters, numberOfPhotonsInPhos) ;
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])) )
191 fhPHOSInvariantMass->Fill(tempMass*1000.);
196 PostData(0, fOutputContainer);
198 delete [] phosVector ;
199 delete [] phosPhotonsEnergy ;
203 //______________________________________________________________________________
204 void AliPHOSESDQA::Terminate(Option_t *)
206 // Processing when the event loop is ended
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);
218 Bool_t problem = kFALSE ;
219 AliDebug(AliQAv1::GetQADebugLevel(), Form(" *** %s Report:", GetName())) ;
220 if ( AliDebugLevel() == AliQAv1::GetQADebugLevel() ) {
221 printf(" PHOSEnergy Mean : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(), fhPHOSEnergy->GetRMS() ) ;
222 printf(" PHOSDigits Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(), fhPHOSDigits->GetRMS() ) ;
223 printf(" PHOSRecParticles Mean : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(), fhPHOSRecParticles->GetRMS() ) ;
224 printf(" PHOSPhotons Mean : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(), fhPHOSPhotons->GetRMS() ) ;
225 printf(" PHOSInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(), fhPHOSInvariantMass->GetRMS() ) ;
226 printf(" PHOSDigitsEvent Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(), fhPHOSDigitsEvent->GetRMS() ) ;
229 TCanvas * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
233 if ( fhPHOSEnergy->GetMaximum() > 0. )
235 fhPHOSEnergy->SetAxisRange(0, 25.);
236 fhPHOSEnergy->SetLineColor(2);
237 fhPHOSEnergy->Draw();
240 fhPHOSDigits->SetAxisRange(0,25.);
241 fhPHOSDigits->SetLineColor(2);
242 fhPHOSDigits->Draw();
245 if ( fhPHOSRecParticles->GetMaximum() > 0. )
247 fhPHOSRecParticles->SetAxisRange(0, 25.);
248 fhPHOSRecParticles->SetLineColor(2);
249 fhPHOSRecParticles->Draw();
252 if ( fhPHOSPhotons->GetMaximum() > 0. )
254 fhPHOSPhotons->SetAxisRange(0,25.);
255 fhPHOSPhotons->SetLineColor(2);
256 fhPHOSPhotons->Draw();
259 fhPHOSInvariantMass->SetLineColor(2);
260 fhPHOSInvariantMass->Draw();
263 if ( fhPHOSDigitsEvent->GetMaximum() > 0. )
265 fhPHOSDigitsEvent->SetAxisRange(0,40.);
266 fhPHOSDigitsEvent->SetLineColor(2);
267 fhPHOSDigitsEvent->Draw();
269 cPHOS->Print("PHOS.eps");
272 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
273 gROOT->ProcessLine(line);
274 sprintf(line, ".!rm -fR *.eps");
275 gROOT->ProcessLine(line);
277 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
281 report="Problems found, please check!!!";
285 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ;