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 EMCAL photon data in simulated data
22 //////////////////////////////////////////////////////////////////////////////
35 #include "AliEMCALQATask.h"
39 //______________________________________________________________________________
40 AliEMCALQATask::AliEMCALQATask(const char *name) :
41 AliAnalysisTask(name,""),
48 fhEMCALRecParticles(0),
50 fhEMCALInvariantMass(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 AliEMCALQATask::~AliEMCALQATask()
64 fOutputContainer->Clear() ;
65 delete fOutputContainer ;
69 delete fhEMCALEnergy ;
70 delete fhEMCALDigits ;
71 delete fhEMCALRecParticles ;
72 delete fhEMCALPhotons ;
73 delete fhEMCALInvariantMass ;
74 delete fhEMCALDigitsEvent ;
77 //______________________________________________________________________________
78 void AliEMCALQATask::ConnectInputData(const Option_t*)
80 // Initialisation of branch container and histograms
82 AliInfo(Form("*** Initialization of %s", GetName())) ;
85 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
87 AliError(Form("Input 0 for %s not found\n", GetName()));
91 // One should first check if the branch address was taken by some other task
92 char ** address = (char **)GetBranchAddress(0, "ESD");
94 fESD = (AliESD*)(*address);
97 SetBranchAddress(0, "ESD", &fESD);
101 //______________________________________________________________________________
102 void AliEMCALQATask::CreateOutputObjects()
108 fhEMCALPos = new TNtuple("EMCALPos" , "Position in EMCAL" , "x:y:z");
109 fhEMCAL = new TNtuple("EMCAL" , "EMCAL" , "event:digits:clusters:photons");
110 fhEMCALEnergy = new TH1D("EMCALEnergy" , "EMCALEnergy" , 1000, 0., 10. ) ;
111 fhEMCALDigits = new TH1I("EMCALDigitsCluster" , "EMCALDigits" , 20 , 0 , 20 ) ;
112 fhEMCALRecParticles = new TH1D("EMCALRecParticles" , "EMCALRecParticles", 20 , 0., 20. ) ;
113 fhEMCALPhotons = new TH1I("EMCALPhotons" , "EMCALPhotons" , 20 , 0 , 20 ) ;
114 fhEMCALInvariantMass = new TH1D("EMCALInvariantMass" , "EMCALInvariantMass", 400, 0., 400.) ;
115 fhEMCALDigitsEvent = new TH1I("EMCALDigitsEvent" , "EMCALDigitsEvent" , 30 , 0 , 30 ) ;
117 // create output container
119 fOutputContainer = new TObjArray(8) ;
120 fOutputContainer->SetName(GetName()) ;
122 fOutputContainer->AddAt(fhEMCALPos, 0) ;
123 fOutputContainer->AddAt(fhEMCAL, 1) ;
124 fOutputContainer->AddAt(fhEMCALEnergy, 2) ;
125 fOutputContainer->AddAt(fhEMCALDigits, 3) ;
126 fOutputContainer->AddAt(fhEMCALRecParticles, 4) ;
127 fOutputContainer->AddAt(fhEMCALPhotons, 5) ;
128 fOutputContainer->AddAt(fhEMCALInvariantMass, 6) ;
129 fOutputContainer->AddAt(fhEMCALDigitsEvent, 7) ;
132 //______________________________________________________________________________
133 void AliEMCALQATask::Exec(Option_t *)
135 // Processing of one event
137 Long64_t entry = fChain->GetReadEntry() ;
140 AliError("fESD is not connected to the input!") ;
144 if ( !((entry-1)%100) )
145 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
147 //************************ EMCAL *************************************
148 Int_t firstEmcalCluster = fESD->GetFirstEMCALCluster() ;
149 const Int_t mumberOfEmcalClusters = fESD->GetNumberOfEMCALClusters() ;
151 TVector3 ** emcalVector = new TVector3*[mumberOfEmcalClusters] ;
152 Float_t * emcalPhotonsEnergy = new Float_t[mumberOfEmcalClusters] ;
154 Int_t numberOfEmcalClustersv1 = 0 ;
155 Int_t numberOfPhotonsInEmcal = 0 ;
156 Int_t numberOfDigitsInEmcal = 0 ;
159 // loop over all the EMCAL Cluster
160 for(emcalCluster = firstEmcalCluster ; emcalCluster < firstEmcalCluster + mumberOfEmcalClusters ; emcalCluster++) {
161 AliESDCaloCluster * caloCluster = fESD->GetCaloCluster(emcalCluster) ;
164 if(caloCluster->GetClusterType() == AliESDCaloCluster::kClusterv1) {
165 caloCluster->GetGlobalPosition(pos) ;
166 fhEMCALPos->Fill(pos[0],pos[1],pos[2]) ;
167 fhEMCALEnergy->Fill(caloCluster->GetClusterEnergy()) ;
168 fhEMCALDigits->Fill(entry, caloCluster->GetNumberOfDigits()) ;
169 numberOfEmcalClustersv1++ ;
170 numberOfDigitsInEmcal += caloCluster->GetNumberOfDigits() ;
171 // Float_t * pid = clus->GetPid() ;
172 // if(pid[AliPID::kPhoton]>0.9){
173 emcalVector[numberOfPhotonsInEmcal] = new TVector3(pos[0],pos[1],pos[2]) ;
174 emcalPhotonsEnergy[numberOfPhotonsInEmcal] = caloCluster->GetClusterEnergy() ;
175 numberOfPhotonsInEmcal++ ;
178 } // EMCAL clusters loop
180 fhEMCALRecParticles->Fill(numberOfEmcalClustersv1);
181 fhEMCALPhotons->Fill(numberOfPhotonsInEmcal);
182 fhEMCALDigitsEvent->Fill(numberOfDigitsInEmcal);
183 fhEMCAL->Fill(entry, numberOfDigitsInEmcal, numberOfEmcalClustersv1, numberOfPhotonsInEmcal) ;
186 if (numberOfPhotonsInEmcal > 1 ) {
187 Int_t emcalPhoton1, emcalPhoton2 ;
188 for(emcalPhoton1 = 0 ; emcalPhoton1 < numberOfPhotonsInEmcal ; emcalPhoton1++) {
189 for(emcalPhoton2 = emcalPhoton1 + 1 ; emcalPhoton2 < numberOfPhotonsInEmcal ; emcalPhoton2++) {
190 Float_t tempMass = TMath::Sqrt( 2 * emcalPhotonsEnergy[emcalPhoton1] * emcalPhotonsEnergy[emcalPhoton2] *
191 ( 1 - TMath::Cos( emcalVector[emcalPhoton1]->Angle(*emcalVector[emcalPhoton2]))
194 fhEMCALInvariantMass->Fill(tempMass*1000.);
199 PostData(0, fOutputContainer);
201 delete [] emcalVector ;
202 delete [] emcalPhotonsEnergy ;
205 //______________________________________________________________________________
206 void AliEMCALQATask::Terminate(Option_t *)
208 // Processing when the event loop is ended
209 fOutputContainer = (TObjArray*)GetOutputData(0);
210 fhEMCALEnergy = (TH1D*)fOutputContainer->At(2);
211 fhEMCALDigits = (TH1I*)fOutputContainer->At(3);
212 fhEMCALRecParticles = (TH1D*)fOutputContainer->At(4);
213 fhEMCALPhotons = (TH1I*)fOutputContainer->At(5);
214 fhEMCALInvariantMass = (TH1D*)fOutputContainer->At(6);
215 fhEMCALDigitsEvent = (TH1I*)fOutputContainer->At(7);
217 Bool_t problem = kFALSE ;
218 AliInfo(Form(" *** %s Report:", GetName())) ;
219 printf(" EMCALEnergy Mean : %5.3f , RMS : %5.3f \n", fhEMCALEnergy->GetMean(), fhEMCALEnergy->GetRMS() ) ;
220 printf(" EMCALDigits Mean : %5.3f , RMS : %5.3f \n", fhEMCALDigits->GetMean(), fhEMCALDigits->GetRMS() ) ;
221 printf(" EMCALRecParticles Mean : %5.3f , RMS : %5.3f \n", fhEMCALRecParticles->GetMean(), fhEMCALRecParticles->GetRMS() ) ;
222 printf(" EMCALPhotons Mean : %5.3f , RMS : %5.3f \n", fhEMCALPhotons->GetMean(), fhEMCALPhotons->GetRMS() ) ;
223 printf(" EMCALInvariantMass Mean : %5.3f , RMS : %5.3f \n", fhEMCALInvariantMass->GetMean(), fhEMCALInvariantMass->GetRMS() ) ;
224 printf(" EMCALDigitsEvent Mean : %5.3f , RMS : %5.3f \n", fhEMCALDigitsEvent->GetMean(), fhEMCALDigitsEvent->GetRMS() ) ;
226 TCanvas * cEMCAL = new TCanvas("EMCAL", "EMCAL ESD Test", 400, 10, 600, 700);
227 cEMCAL->Divide(3, 2) ;
230 if ( fhEMCALEnergy->GetMaximum() > 0. )
232 fhEMCALEnergy->SetAxisRange(0, 25.);
233 fhEMCALEnergy->SetXTitle("Energy (GeV)");
234 fhEMCALEnergy->Draw();
237 if ( fhEMCALDigits->GetMaximum() > 0. )
239 fhEMCALDigits->SetAxisRange(0, 25.);
240 fhEMCALDigits->SetXTitle("DigitsPerCluster");
241 fhEMCALDigits->Draw();
244 if ( fhEMCALRecParticles->GetMaximum() > 0. )
246 fhEMCALRecParticles->SetAxisRange(0, 25.);
247 fhEMCALRecParticles->SetXTitle("RecParticles");
248 fhEMCALRecParticles->Draw();
251 if ( fhEMCALPhotons->GetMaximum() > 0. )
253 fhEMCALPhotons->SetAxisRange(0, 25.);
254 fhEMCALPhotons->SetXTitle("Photons");
255 fhEMCALPhotons->Draw();
258 fhEMCALInvariantMass->SetXTitle("InvariantMass (MeV/c²)");
259 fhEMCALInvariantMass->Draw();
262 if ( fhEMCALDigitsEvent->GetMaximum() > 0. )
264 fhEMCALDigitsEvent->SetAxisRange(0, 40.);
265 fhEMCALDigitsEvent->SetXTitle("DigitsPerEvent");
266 fhEMCALDigitsEvent->Draw();
268 cEMCAL->Print("EMCAL.eps");
271 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
272 gROOT->ProcessLine(line);
273 sprintf(line, ".!rm -fR *.eps");
274 gROOT->ProcessLine(line);
276 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
280 report="Problems found, please check!!!";
284 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report)) ;