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,""),
50 fhPHOSRecParticles(0),
52 fhPHOSInvariantMass(0),
56 // Input slot #0 works with an Ntuple
57 DefineInput(0, TChain::Class());
58 // Output slot #0 writes into a TH1 container
59 DefineOutput(0, TObjArray::Class()) ;
62 //____________________________________________________________________________
63 AliPHOSQATask::AliPHOSQATask(const AliPHOSQATask& ta) :
64 AliAnalysisTask(ta.GetName(),""),
67 fOutputContainer(ta.fOutputContainer),
68 fhPHOSPos(ta.fhPHOSPos),
70 fhPHOSEnergy(ta.fhPHOSEnergy),
71 fhPHOSDigits(ta.fhPHOSDigits),
72 fhPHOSRecParticles(ta.fhPHOSRecParticles),
73 fhPHOSPhotons(ta.fhPHOSPhotons),
74 fhPHOSInvariantMass(ta.fhPHOSInvariantMass),
75 fhPHOSDigitsEvent(ta.fhPHOSDigitsEvent)
80 //_____________________________________________________________________________
81 AliPHOSQATask& AliPHOSQATask::operator = (const AliPHOSQATask& ap)
83 // assignment operator
85 this->~AliPHOSQATask();
86 new(this) AliPHOSQATask(ap);
90 //______________________________________________________________________________
91 AliPHOSQATask::~AliPHOSQATask()
94 fOutputContainer->Clear() ;
95 delete fOutputContainer ;
100 delete fhPHOSDigits ;
101 delete fhPHOSRecParticles ;
102 delete fhPHOSPhotons ;
103 delete fhPHOSInvariantMass ;
104 delete fhPHOSDigitsEvent ;
108 //______________________________________________________________________________
109 void AliPHOSQATask::ConnectInputData(const Option_t*)
111 // Initialisation of branch container and histograms
113 AliInfo(Form("*** Initialization of %s", GetName())) ;
116 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
118 AliError(Form("Input 0 for %s not found\n", GetName()));
122 // One should first check if the branch address was taken by some other task
123 char ** address = (char **)GetBranchAddress(0, "ESD");
125 fESD = (AliESD*)(*address);
128 SetBranchAddress(0, "ESD", &fESD);
132 //________________________________________________________________________
133 void AliPHOSQATask::CreateOutputObjects()
139 fhPHOSPos = new TNtuple("PHOSPos" , "Position in PHOS" , "x:y:z");
140 fhPHOS = new TNtuple("PHOS" , "PHOS" , "event:digits:clusters:photons");
141 fhPHOSEnergy = new TH1D("PHOSEnergy" , "PHOSEnergy" , 1000, 0., 10. ) ;
142 fhPHOSDigits = new TH1I("PHOSDigitsCluster" , "PHOSDigits" , 20 , 0 , 20 ) ;
143 fhPHOSRecParticles = new TH1D("PHOSRecParticles" , "PHOSRecParticles" , 20 , 0., 20. ) ;
144 fhPHOSPhotons = new TH1I("PHOSPhotons" , "PHOSPhotons" , 20 , 0 , 20 ) ;
145 fhPHOSInvariantMass = new TH1D("PHOSInvariantMass" , "PHOSInvariantMass" , 400, 0., 400.) ;
146 fhPHOSDigitsEvent = new TH1I("PHOSDigitsEvent" , "PHOSDigitsEvent" , 30 , 0 , 30 ) ;
148 // create output container
150 fOutputContainer = new TObjArray(8) ;
151 fOutputContainer->SetName(GetName()) ;
153 fOutputContainer->AddAt(fhPHOSPos, 0) ;
154 fOutputContainer->AddAt(fhPHOS, 1) ;
155 fOutputContainer->AddAt(fhPHOSEnergy, 2) ;
156 fOutputContainer->AddAt(fhPHOSDigits, 3) ;
157 fOutputContainer->AddAt(fhPHOSRecParticles, 4) ;
158 fOutputContainer->AddAt(fhPHOSPhotons, 5) ;
159 fOutputContainer->AddAt(fhPHOSInvariantMass, 6) ;
160 fOutputContainer->AddAt(fhPHOSDigitsEvent, 7) ;
163 //______________________________________________________________________________
164 void AliPHOSQATask::Exec(Option_t *)
166 // Processing of one event
168 Long64_t entry = fChain->GetReadEntry() ;
171 AliError("fESD is not connected to the input!") ;
175 if ( !((entry-1)%100) )
176 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
178 //************************ PHOS *************************************
180 Int_t firstPhosCluster = fESD->GetFirstPHOSCluster() ;
181 const Int_t kNumberOfPhosClusters = fESD->GetNumberOfPHOSClusters() ;
183 TVector3 ** phosVector = new TVector3*[kNumberOfPhosClusters] ;
184 Float_t * phosPhotonsEnergy = new Float_t[kNumberOfPhosClusters] ;
186 Int_t numberOfPhotonsInPhos = 0 ;
187 Int_t numberOfDigitsInPhos = 0 ;
189 // loop over the PHOS Cluster
190 for(phosCluster = firstPhosCluster ; phosCluster < firstPhosCluster + kNumberOfPhosClusters ; phosCluster++) {
191 AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(phosCluster) ;
194 caloCluster->GetPosition( pos ) ;
195 fhPHOSEnergy->Fill( caloCluster->E() ) ;
196 fhPHOSPos->Fill( pos[0], pos[1], pos[2] ) ;
197 fhPHOSDigits->Fill(entry, caloCluster->GetNumberOfDigits() ) ;
198 numberOfDigitsInPhos += caloCluster->GetNumberOfDigits() ;
199 Double_t * pid = caloCluster->GetPid() ;
200 if(pid[AliPID::kPhoton] > 0.9) {
201 phosVector[numberOfPhotonsInPhos] = new TVector3(pos[0],pos[1],pos[2]) ;
202 phosPhotonsEnergy[numberOfPhotonsInPhos]=caloCluster->E() ;
203 numberOfPhotonsInPhos++;
208 fhPHOSRecParticles->Fill(kNumberOfPhosClusters);
209 fhPHOSPhotons->Fill(numberOfPhotonsInPhos);
210 fhPHOSDigitsEvent->Fill(numberOfDigitsInPhos);
211 fhPHOS->Fill(entry, numberOfDigitsInPhos, kNumberOfPhosClusters, numberOfPhotonsInPhos) ;
214 if (numberOfPhotonsInPhos > 1 ) {
215 Int_t phosPhoton1, phosPhoton2 ;
216 for(phosPhoton1 = 0 ; phosPhoton1 < numberOfPhotonsInPhos ; phosPhoton1++) {
217 for(phosPhoton2 = phosPhoton1 + 1 ; phosPhoton2 < numberOfPhotonsInPhos ; phosPhoton2++) {
218 Float_t tempMass = TMath::Sqrt( 2 * phosPhotonsEnergy[phosPhoton1] * phosPhotonsEnergy[phosPhoton2] *
219 ( 1 - TMath::Cos(phosVector[phosPhoton1]->Angle(*phosVector[phosPhoton2])) )
221 fhPHOSInvariantMass->Fill(tempMass*1000.);
226 PostData(0, fOutputContainer);
228 delete [] phosVector ;
229 delete [] phosPhotonsEnergy ;
233 //______________________________________________________________________________
234 void AliPHOSQATask::Terminate(Option_t *)
236 // Processing when the event loop is ended
238 fOutputContainer = (TObjArray*)GetOutputData(0);
239 fhPHOSPos = (TNtuple*)fOutputContainer->At(0);
240 fhPHOS = (TNtuple*)fOutputContainer->At(1);
241 fhPHOSEnergy = (TH1D*)fOutputContainer->At(2);
242 fhPHOSDigits = (TH1I*)fOutputContainer->At(3);
243 fhPHOSRecParticles = (TH1D*)fOutputContainer->At(4);
244 fhPHOSPhotons = (TH1I*)fOutputContainer->At(5);
245 fhPHOSInvariantMass = (TH1D*)fOutputContainer->At(6);
246 fhPHOSDigitsEvent = (TH1I*)fOutputContainer->At(7);
248 Bool_t problem = kFALSE ;
249 AliInfo(Form(" *** %s Report:", GetName())) ;
250 printf(" PHOSEnergy Mean : %5.3f , RMS : %5.3f \n", fhPHOSEnergy->GetMean(), fhPHOSEnergy->GetRMS() ) ;
251 printf(" PHOSDigits Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigits->GetMean(), fhPHOSDigits->GetRMS() ) ;
252 printf(" PHOSRecParticles Mean : %5.3f , RMS : %5.3f \n", fhPHOSRecParticles->GetMean(), fhPHOSRecParticles->GetRMS() ) ;
253 printf(" PHOSPhotons Mean : %5.3f , RMS : %5.3f \n", fhPHOSPhotons->GetMean(), fhPHOSPhotons->GetRMS() ) ;
254 printf(" PHOSInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhPHOSInvariantMass->GetMean(), fhPHOSInvariantMass->GetRMS() ) ;
255 printf(" PHOSDigitsEvent Mean : %5.3f , RMS : %5.3f \n", fhPHOSDigitsEvent->GetMean(), fhPHOSDigitsEvent->GetRMS() ) ;
257 TCanvas * cPHOS = new TCanvas("cPHOS", "PHOS ESD Test", 400, 10, 600, 700) ;
261 if ( fhPHOSEnergy->GetMaximum() > 0. )
263 fhPHOSEnergy->SetAxisRange(0, 25.);
264 fhPHOSEnergy->SetLineColor(2);
265 fhPHOSEnergy->Draw();
268 fhPHOSDigits->SetAxisRange(0,25.);
269 fhPHOSDigits->SetLineColor(2);
270 fhPHOSDigits->Draw();
273 if ( fhPHOSRecParticles->GetMaximum() > 0. )
275 fhPHOSRecParticles->SetAxisRange(0, 25.);
276 fhPHOSRecParticles->SetLineColor(2);
277 fhPHOSRecParticles->Draw();
280 if ( fhPHOSPhotons->GetMaximum() > 0. )
282 fhPHOSPhotons->SetAxisRange(0,25.);
283 fhPHOSPhotons->SetLineColor(2);
284 fhPHOSPhotons->Draw();
287 fhPHOSInvariantMass->SetLineColor(2);
288 fhPHOSInvariantMass->Draw();
291 if ( fhPHOSDigitsEvent->GetMaximum() > 0. )
293 fhPHOSDigitsEvent->SetAxisRange(0,40.);
294 fhPHOSDigitsEvent->SetLineColor(2);
295 fhPHOSDigitsEvent->Draw();
297 cPHOS->Print("PHOS.eps");
300 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
301 gROOT->ProcessLine(line);
302 sprintf(line, ".!rm -fR *.eps");
303 gROOT->ProcessLine(line);
305 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
309 report="Problems found, please check!!!";
313 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ;