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 // Author : R. Vernet, Consorzio Cometa - Catania (it)
18 //-----------------------------------------------------------------------
21 #include "AliCFRsnTask.h"
24 #include "TParticle.h"
27 #include "AliMCEvent.h"
28 #include "AliESDEvent.h"
29 #include "AliCFManager.h"
30 #include "AliCFCutBase.h"
31 #include "AliCFContainer.h"
32 #include "AliESDtrack.h"
34 #include "AliRsnDaughter.h"
35 #include "AliCFPair.h"
36 #include "AliRsnParticle.h"
38 //__________________________________________________________________________
39 AliCFRsnTask::AliCFRsnTask() :
43 fHistEventsProcessed(0x0)
49 //___________________________________________________________________________
50 AliCFRsnTask::AliCFRsnTask(const Char_t* name) :
51 AliAnalysisTaskSE(name),
54 fHistEventsProcessed(0x0)
57 // Constructor. Initialization of Inputs and Outputs
59 Info("AliCFRsnTask","Calling Constructor");
61 DefineInput(0) and DefineOutput(0)
62 are taken care of by AliAnalysisTaskSE constructor
64 DefineOutput(1,TH1I::Class());
65 DefineOutput(2,AliCFContainer::Class());
66 // DefineOutput(3,TList::Class());
69 //___________________________________________________________________________
70 AliCFRsnTask& AliCFRsnTask::operator=(const AliCFRsnTask& c)
73 // Assignment operator
76 AliAnalysisTaskSE::operator=(c) ;
78 fCFManager = c.fCFManager;
79 fHistEventsProcessed = c.fHistEventsProcessed;
84 //___________________________________________________________________________
85 AliCFRsnTask::AliCFRsnTask(const AliCFRsnTask& c) :
88 fCFManager(c.fCFManager),
89 fHistEventsProcessed(c.fHistEventsProcessed)
96 //___________________________________________________________________________
97 AliCFRsnTask::~AliCFRsnTask() {
101 Info("~AliCFRsnTask","Calling Destructor");
102 if (fCFManager) delete fCFManager ;
103 if (fHistEventsProcessed) delete fHistEventsProcessed ;
106 //_________________________________________________
107 void AliCFRsnTask::UserExec(Option_t *)
110 // Main loop function
112 Info("UserExec","") ;
114 AliESDEvent* fESD = dynamic_cast<AliESDEvent*>(fInputEvent);
116 Error("UserExec","NO ESD FOUND!");
120 if (!fMCEvent) Error("UserExec","NO MC INFO FOUND!");
121 fCFManager->SetEventInfo(fMCEvent);
123 AliStack* stack = fMCEvent->Stack();
125 // MC-event selection
126 Double_t containerInput[2] ;
128 //loop on the MC event
129 Info("UserExec","Looping on MC event");
130 for (Int_t ipart=0; ipart<stack->GetNprimary(); ipart++) {
131 AliMCParticle *mcPart = fMCEvent->GetTrack(ipart);
133 //check the MC-level cuts
134 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcPart)) continue;
135 containerInput[0] = mcPart->Pt();
136 containerInput[1] = mcPart->Eta() ;
137 //fill the container for Gen-level selection
138 fCFManager->GetParticleContainer()->Fill(containerInput,kStepGenerated);
140 //check the Acceptance-level cuts
141 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartAccCuts,mcPart)) continue;
142 //fill the container for Acceptance-level selection
143 fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructible);
147 //Now go to rec level
148 Info("UserExec","Looping on ESD event");
150 //SET THE ESD AS EVENT INFO IN RECONSTRUCTION CUTS
151 TObjArray* fCutsReco = fCFManager->GetParticleCutsList(AliCFManager::kPartRecCuts);
152 TObjArray* fCutsPID = fCFManager->GetParticleCutsList(AliCFManager::kPartSelCuts);
153 TObjArrayIter iter1(fCutsReco);
154 TObjArrayIter iter2(fCutsPID);
155 AliCFCutBase *cut = 0;
156 while ( (cut = (AliCFCutBase*)iter1.Next()) ) {
157 cut->SetEvtInfo(fESD);
159 while ( (cut = (AliCFCutBase*)iter2.Next()) ) {
160 cut->SetEvtInfo(fESD);
163 for (Int_t iTrack1 = 0; iTrack1<fESD->GetNumberOfTracks(); iTrack1++) {
164 AliESDtrack* esdTrack1 = fESD->GetTrack(iTrack1);
166 if (esdTrack1->Charge()>=0) continue;
167 Int_t esdLabel1 = esdTrack1->GetLabel();
168 if (esdLabel1<0) continue;
170 for (Int_t iTrack2 = 0; iTrack2<fESD->GetNumberOfTracks(); iTrack2++) {
171 AliESDtrack* esdTrack2 = fESD->GetTrack(iTrack2);
173 if (esdTrack2->Charge()<=0) continue;
174 Int_t esdLabel2 = esdTrack2->GetLabel();
175 if (esdLabel2<0) continue;
177 AliRsnDaughter* tmp1 = new AliRsnDaughter(esdTrack1);
178 AliRsnDaughter* tmp2 = new AliRsnDaughter(esdTrack2);
179 AliRsnDaughter track1(*tmp1);
180 AliRsnDaughter track2(*tmp2);
184 TParticle *part1 = stack->Particle(esdLabel1);
185 track1.InitParticle(part1);
186 track1.GetParticle()->SetPDG(part1->GetPdgCode());
188 Int_t mother1 = part1->GetFirstMother();
189 track1.GetParticle()->SetMother(mother1);
191 TParticle *mum = stack->Particle(mother1);
192 track1.GetParticle()->SetMotherPDG(mum->GetPdgCode());
195 TParticle *part2 = stack->Particle(esdLabel2);
196 track2.InitParticle(part2);
197 track2.GetParticle()->SetPDG(part2->GetPdgCode());
199 Int_t mother2 = part2->GetFirstMother();
200 track2.GetParticle()->SetMother(mother2);
202 TParticle *mum = stack->Particle(mother2);
203 track2.GetParticle()->SetMotherPDG(mum->GetPdgCode());
206 //make a mother resonance from the 2 candidate daughters
207 AliRsnDaughter rsn = AliRsnDaughter::Sum(track1,track2) ;
208 AliCFPair pair(esdTrack1,esdTrack2);
210 //check if true resonance
211 if (rsn.GetParticle()->MotherPDG() != fRsnPDG) continue;
212 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,&pair)) continue;
214 //check if associated MC resonance passes the cuts
215 Int_t motherLabel=rsn.Label();
216 if (motherLabel<0) continue ;
217 AliMCParticle* mcRsn = fMCEvent->GetTrack(motherLabel);
218 if (!mcRsn) continue;
219 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,mcRsn)) continue;
222 containerInput[0] = rsn.Pt() ;
223 containerInput[1] = GetRapidity(rsn.E(),rsn.Pz());
224 fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructed) ;
226 if (!fCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,&pair)) continue ;
227 fCFManager->GetParticleContainer()->Fill(containerInput,kStepSelected);
231 fHistEventsProcessed->Fill(0);
232 /* PostData(0) is taken care of by AliAnalysisTaskSE */
233 PostData(1,fHistEventsProcessed) ;
234 PostData(2,fCFManager->GetParticleContainer()) ;
236 // TList * list = new TList();
237 // fCFManager->AddQAHistosToList(list);
238 // PostData(2,list) ;
242 //___________________________________________________________________________
243 void AliCFRsnTask::Terminate(Option_t*)
245 // The Terminate() function is the last function to be called during
246 // a query. It always runs on the client, it can be used to present
247 // the results graphically or save the results to file.
249 Info("Terminate","");
251 Double_t max1 = fCFManager->GetParticleContainer()->ShowProjection(0,0)->GetMaximum();
252 Double_t max2 = fCFManager->GetParticleContainer()->ShowProjection(1,0)->GetMaximum();
254 fCFManager->GetParticleContainer()->ShowProjection(0,0)->GetYaxis()->SetRangeUser(0,max1*1.2);
255 fCFManager->GetParticleContainer()->ShowProjection(0,1)->GetYaxis()->SetRangeUser(0,max1*1.2);
256 fCFManager->GetParticleContainer()->ShowProjection(0,2)->GetYaxis()->SetRangeUser(0,max1*1.2);
257 fCFManager->GetParticleContainer()->ShowProjection(0,3)->GetYaxis()->SetRangeUser(0,max1*1.2);
259 fCFManager->GetParticleContainer()->ShowProjection(1,0)->GetYaxis()->SetRangeUser(0,max2*1.2);
260 fCFManager->GetParticleContainer()->ShowProjection(1,1)->GetYaxis()->SetRangeUser(0,max2*1.2);
261 fCFManager->GetParticleContainer()->ShowProjection(1,2)->GetYaxis()->SetRangeUser(0,max2*1.2);
262 fCFManager->GetParticleContainer()->ShowProjection(1,3)->GetYaxis()->SetRangeUser(0,max2*1.2);
264 fCFManager->GetParticleContainer()->ShowProjection(0,0)->SetMarkerStyle(23) ;
265 fCFManager->GetParticleContainer()->ShowProjection(0,1)->SetMarkerStyle(24) ;
266 fCFManager->GetParticleContainer()->ShowProjection(0,2)->SetMarkerStyle(25) ;
267 fCFManager->GetParticleContainer()->ShowProjection(0,3)->SetMarkerStyle(26) ;
269 fCFManager->GetParticleContainer()->ShowProjection(1,0)->SetMarkerStyle(23) ;
270 fCFManager->GetParticleContainer()->ShowProjection(1,1)->SetMarkerStyle(24) ;
271 fCFManager->GetParticleContainer()->ShowProjection(1,2)->SetMarkerStyle(25) ;
272 fCFManager->GetParticleContainer()->ShowProjection(1,3)->SetMarkerStyle(26) ;
274 TCanvas * c =new TCanvas("c","",1400,800);
278 fCFManager->GetParticleContainer()->ShowProjection(0,0)->Draw("p");
280 fCFManager->GetParticleContainer()->ShowProjection(0,1)->Draw("p");
282 fCFManager->GetParticleContainer()->ShowProjection(0,2)->Draw("p");
284 fCFManager->GetParticleContainer()->ShowProjection(0,3)->Draw("p");
286 fCFManager->GetParticleContainer()->ShowProjection(1,0)->Draw("p");
288 fCFManager->GetParticleContainer()->ShowProjection(1,1)->Draw("p");
290 fCFManager->GetParticleContainer()->ShowProjection(1,2)->Draw("p");
292 fCFManager->GetParticleContainer()->ShowProjection(1,3)->Draw("p");
294 c->SaveAs("plots.eps");
296 delete fHistEventsProcessed ;
299 //___________________________________________________________________________
300 void AliCFRsnTask::UserCreateOutputObjects() {
301 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
302 //TO BE SET BEFORE THE EXECUTION OF THE TASK
304 Info("UserCreateOutputObjects","CreateOutputObjects of task %s\n", GetName());
308 fHistEventsProcessed = new TH1I("fHistEventsProcessed","",1,0,1) ;
311 //___________________________________________________________________________
312 Double_t AliCFRsnTask::GetRapidity(Double_t energy, Double_t pz) {
313 if (energy == pz || energy == -pz) {
314 printf("GetRapidity : ERROR : rapidity for 4-vector with E = Pz -- infinite result");
318 printf("GetRapidity : ERROR : rapidity for 4-vector with E = Pz -- infinite result");
321 Double_t y = 0.5 * log((energy + pz) / (energy - pz));