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
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
27 //////////////////////////////////////////////////////////////////////////////
38 #include "AliPHOSQATask.h"
42 //______________________________________________________________________________
43 AliPHOSQATask::AliPHOSQATask(const char *name) :
44 AliAnalysisTask(name,""),
52 fhPHOSRecParticles(0),
54 fhPHOSInvariantMass(0),
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()) ;
64 //____________________________________________________________________________
65 AliPHOSQATask::AliPHOSQATask(const AliPHOSQATask& ta) :
66 AliAnalysisTask(ta.GetName(),""),
69 fOutputContainer(ta.fOutputContainer),
70 fhPHOSPos(ta.fhPHOSPos),
72 fhPHOSEnergy(ta.fhPHOSEnergy),
73 fhPHOSDigits(ta.fhPHOSDigits),
74 fhPHOSRecParticles(ta.fhPHOSRecParticles),
75 fhPHOSPhotons(ta.fhPHOSPhotons),
76 fhPHOSInvariantMass(ta.fhPHOSInvariantMass),
77 fhPHOSDigitsEvent(ta.fhPHOSDigitsEvent)
82 //_____________________________________________________________________________
83 AliPHOSQATask& AliPHOSQATask::operator = (const AliPHOSQATask& ap)
85 // assignment operator
87 this->~AliPHOSQATask();
88 new(this) AliPHOSQATask(ap);
92 //______________________________________________________________________________
93 AliPHOSQATask::~AliPHOSQATask()
96 fOutputContainer->Clear() ;
97 delete fOutputContainer ;
101 delete fhPHOSEnergy ;
102 delete fhPHOSDigits ;
103 delete fhPHOSRecParticles ;
104 delete fhPHOSPhotons ;
105 delete fhPHOSInvariantMass ;
106 delete fhPHOSDigitsEvent ;
110 //______________________________________________________________________________
111 void AliPHOSQATask::ConnectInputData(const Option_t*)
113 // Initialisation of branch container and histograms
115 AliInfo(Form("*** Initialization of %s", GetName())) ;
118 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
120 AliError(Form("Input 0 for %s not found\n", GetName()));
124 // One should first check if the branch address was taken by some other task
125 char ** address = (char **)GetBranchAddress(0, "ESD");
127 fESD = (AliESD*)(*address);
130 SetBranchAddress(0, "ESD", &fESD);
134 //________________________________________________________________________
135 void AliPHOSQATask::CreateOutputObjects()
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 ) ;
150 // create output container
152 fOutputContainer = new TObjArray(8) ;
153 fOutputContainer->SetName(GetName()) ;
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) ;
165 //______________________________________________________________________________
166 void AliPHOSQATask::Exec(Option_t *)
168 // Processing of one event
170 Long64_t entry = fChain->GetReadEntry() ;
173 AliError("fESD is not connected to the input!") ;
177 if ( !((entry-1)%100) )
178 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
180 //************************ PHOS *************************************
182 Int_t firstPhosCluster = fESD->GetFirstPHOSCluster() ;
183 const Int_t kNumberOfPhosClusters = fESD->GetNumberOfPHOSClusters() ;
185 TVector3 ** phosVector = new TVector3*[kNumberOfPhosClusters] ;
186 Float_t * phosPhotonsEnergy = new Float_t[kNumberOfPhosClusters] ;
188 Int_t numberOfPhotonsInPhos = 0 ;
189 Int_t numberOfDigitsInPhos = 0 ;
191 // loop over the PHOS Cluster
192 for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + kNumberOfPhosClusters ; phosCluster++) {
193 AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
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++;
210 fhPHOSRecParticles->Fill(kNumberOfPhosClusters);
211 fhPHOSPhotons->Fill(numberOfPhotonsInPhos);
212 fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
213 fhPHOS->Fill(entry, numberOfDigitsInPhos, kNumberOfPhosClusters, numberOfPhotonsInPhos) ;
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])) )
223 fhPHOSInvariantMass->Fill(tempMass*1000.);
228 PostData(0, fOutputContainer);
230 delete [] phosVector ;
231 delete [] phosPhotonsEnergy ;
235 //______________________________________________________________________________
236 void AliPHOSQATask::Terminate(Option_t *)
238 // Processing when the event loop is ended
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);
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() ) ;
259 TCanvas * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
263 if ( fhPHOSEnergy->GetMaximum() > 0. )
265 fhPHOSEnergy->SetAxisRange(0, 25.);
266 fhPHOSEnergy->SetLineColor(2);
267 fhPHOSEnergy->Draw();
270 fhPHOSDigits->SetAxisRange(0,25.);
271 fhPHOSDigits->SetLineColor(2);
272 fhPHOSDigits->Draw();
275 if ( fhPHOSRecParticles->GetMaximum() > 0. )
277 fhPHOSRecParticles->SetAxisRange(0, 25.);
278 fhPHOSRecParticles->SetLineColor(2);
279 fhPHOSRecParticles->Draw();
282 if ( fhPHOSPhotons->GetMaximum() > 0. )
284 fhPHOSPhotons->SetAxisRange(0,25.);
285 fhPHOSPhotons->SetLineColor(2);
286 fhPHOSPhotons->Draw();
289 fhPHOSInvariantMass->SetLineColor(2);
290 fhPHOSInvariantMass->Draw();
293 if ( fhPHOSDigitsEvent->GetMaximum() > 0. )
295 fhPHOSDigitsEvent->SetAxisRange(0,40.);
296 fhPHOSDigitsEvent->SetLineColor(2);
297 fhPHOSDigitsEvent->Draw();
299 cPHOS->Print("PHOS.eps");
302 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
303 gROOT->ProcessLine(line);
304 sprintf(line, ".!rm -fR *.eps");
305 gROOT->ProcessLine(line);
307 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
311 report="Problems found, please check!!!";
315 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ;