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 **************************************************************************/
16 //_________________________________________________________________________
17 // An analysis task to check the trigger data in ESD
18 // Creates an ntuple for 2x2 and NxN triggers
19 // Each ntuple connects the maximum trigger amplitudes
20 // and its positions with reconstructed clusters
22 //*-- Yves Schutz (CERN) & Gustavo Conesa Balbastre (INFN-LNF)
23 //////////////////////////////////////////////////////////////////////////////
30 #include "AliAnaCaloTrigger.h"
31 #include "AliAnalysisManager.h"
32 #include "AliESDEvent.h"
34 #include "AliESDCaloCluster.h"
36 //______________________________________________________________________________
37 AliAnaCaloTrigger::AliAnaCaloTrigger() :
46 // Default Constructor.
50 //______________________________________________________________________________
51 AliAnaCaloTrigger::AliAnaCaloTrigger(const char *name) :
52 AliAnalysisTask(name,"AnaCaloTrigger"),
62 // Input slot #0 works with an Ntuple
63 DefineInput(0, TChain::Class());
64 // Output slot #0 writes into a TH1 container
65 DefineOutput(0, TObjArray::Class()) ;
67 //____________________________________________________________________________
68 AliAnaCaloTrigger::AliAnaCaloTrigger(const AliAnaCaloTrigger & ct) :
69 AliAnalysisTask(ct), fChain(ct.fChain), fESD(ct.fESD),
70 fOutputContainer(ct.fOutputContainer), fCalorimeter(ct. fCalorimeter),
71 fNtTrigger22(ct.fNtTrigger22), fNtTriggerNN(ct.fNtTriggerNN)
75 SetName (ct.GetName()) ;
76 SetTitle(ct.GetTitle()) ;
80 //_________________________________________________________________________
81 AliAnaCaloTrigger & AliAnaCaloTrigger::operator = (const AliAnaCaloTrigger & source)
83 // assignment operator
85 if(&source == this) return *this;
87 fChain = source.fChain ;
89 fOutputContainer = source.fOutputContainer ;
90 fCalorimeter = source. fCalorimeter ;
91 fNtTrigger22 = source.fNtTrigger22 ;
92 fNtTriggerNN = source.fNtTriggerNN ;
98 //______________________________________________________________________________
99 AliAnaCaloTrigger::~AliAnaCaloTrigger()
102 //fOutputContainer->Clear() ;
103 //delete fOutputContainer ;
108 //______________________________________________________________________________
109 void AliAnaCaloTrigger::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 fESD = new AliESDEvent();
123 fESD->ReadFromTree(fChain);
127 //________________________________________________________________________
129 void AliAnaCaloTrigger::CreateOutputObjects()
132 // Create the outputs containers
136 fNtTrigger22 = new TNtuple(fCalorimeter+"trigger22", "Trigger data 2x2 patch", "a22:a220:enMax:phEnMax:eta22:phi22:etaMax:phiMax:phEtaMax:phPhiMax");
137 fNtTriggerNN = new TNtuple(fCalorimeter+"triggerNN", "Trigger data NxN patch", "aNN:aNN0:enMax:phEnMax:etaNN:phiNN:etaMax:phiMax:phEtaMax:phPhiMax");
139 // create output container
141 fOutputContainer = new TObjArray(2) ;
142 fOutputContainer->SetName(GetName()) ;
144 fOutputContainer->AddAt(fNtTrigger22, 0) ;
145 fOutputContainer->AddAt(fNtTriggerNN, 1) ;
149 //______________________________________________________________________________
150 void AliAnaCaloTrigger::Exec(Option_t *)
152 // Processing of one event
154 Long64_t entry = fChain->GetReadEntry() ;
157 AliError("fESD is not connected to the input!") ;
161 if ( !((entry-1)%100) )
162 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
164 // Get trigger information of fCalorimeter
165 TArrayF * triggerAmplitudes = 0x0 ;
166 TArrayF * triggerPosition = 0x0 ;
167 Int_t numberOfCaloClusters = fESD->GetNumberOfCaloClusters() ;
169 if(fCalorimeter == "PHOS"){
170 triggerAmplitudes = fESD->GetPHOSTriggerAmplitudes();
171 triggerPosition = fESD->GetPHOSTriggerPosition();
173 else if(fCalorimeter == "EMCAL"){
174 triggerAmplitudes = fESD->GetEMCALTriggerAmplitudes();
175 triggerPosition = fESD->GetEMCALTriggerPosition();
178 if( triggerAmplitudes && triggerPosition ){
179 // trigger amplitudes
180 const Float_t ka22 = static_cast<Float_t>(triggerAmplitudes->At(0)) ;
181 const Float_t ka22O = static_cast<Float_t>(triggerAmplitudes->At(1)) ;
182 const Float_t kaNN = static_cast<Float_t>(triggerAmplitudes->At(2)) ;
183 const Float_t kaNNO = static_cast<Float_t>(triggerAmplitudes->At(3)) ;
186 const Float_t kx22 = static_cast<Float_t>(triggerPosition->At(0)) ;
187 const Float_t ky22 = static_cast<Float_t>(triggerPosition->At(1)) ;
188 const Float_t kz22 = static_cast<Float_t>(triggerPosition->At(2)) ;
189 const Float_t kxNN = static_cast<Float_t>(triggerPosition->At(3)) ;
190 const Float_t kyNN = static_cast<Float_t>(triggerPosition->At(4)) ;
191 const Float_t kzNN = static_cast<Float_t>(triggerPosition->At(5)) ;
194 Float_t phEnMax = 0. ;
195 Float_t etaMax = 0.5 ;
196 Float_t phiMax = 0. ;
197 Float_t phEtaMax = 0.5 ;
198 Float_t phPhiMax = 0. ;
200 TVector3 vpos22(kx22, ky22, kz22) ;
201 TVector3 vposNN(kxNN, kyNN, kzNN) ;
202 Float_t eta22 = vpos22.Eta() ;
203 Float_t phi22 = vpos22.Phi() * TMath::RadToDeg() + 360. ;
204 Float_t etaNN = vposNN.Eta() ;
205 Float_t phiNN = vposNN.Phi() * TMath::RadToDeg() + 360. ;
209 // loop over the Calorimeters Clusters
211 for(icaloCluster = 0 ; icaloCluster < numberOfCaloClusters ; icaloCluster++) {
213 AliESDCaloCluster * cluster = fESD->GetCaloCluster(icaloCluster) ;
215 if (cluster && ( (fCalorimeter == "PHOS" && cluster->IsPHOS()) ||
216 (fCalorimeter == "EMCAL" && cluster->IsEMCAL()))) {
218 Float_t cluEnergy = cluster->E() ;
222 cluster->GetPosition( pos ) ;
224 if ( cluEnergy > enMax) {
226 vpos.SetXYZ(pos[0], pos[1], pos[2]) ;
227 etaMax = vpos.Eta() ;
228 phiMax = vpos.Phi() ;
231 Double_t * pid = cluster->GetPid() ;
233 if(pid[AliPID::kPhoton] > 0.9) {
234 if ( cluEnergy > phEnMax) {
235 phEnMax = cluEnergy ;
236 vpos.SetXYZ(pos[0], pos[1], pos[2]) ;
237 phEtaMax = vpos.Eta() ;
238 phPhiMax = vpos.Phi() ;
243 fNtTrigger22->Fill(ka22, ka22O, enMax, phEnMax, eta22, phi22, etaMax, phiMax * TMath::RadToDeg() + 360., phEtaMax, phPhiMax * TMath::RadToDeg() + 360.);
244 fNtTriggerNN->Fill(kaNN, kaNNO, enMax, phEnMax, etaNN, phiNN, etaMax, phiMax * TMath::RadToDeg() + 360., phEtaMax, phPhiMax * TMath::RadToDeg() + 360.);
247 }//If trigger arrays filled
249 PostData(0, fOutputContainer);
253 //______________________________________________________________________________
254 void AliAnaCaloTrigger::Terminate(Option_t *)
256 // Processing when the event loop is ended