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 // encapsulates several AliUEHist objects for a full UE analysis plus additional control histograms
19 // the method FillCorrelations is the only method which fills the container
20 // associated particles are filled once for each near side and once for each away side
21 // the pT 2 for the near side particles is set to the highest pT of all accepted away side pT 2
22 // as consequence the pT 2 axis can only be cut for the minimum
23 // if the maximum of the axis should be changed the full analysis needs to rerun
26 // Author: Markus Zimmermann
30 #include "AliTwoPlusOneContainer.h"
34 #include "AliCFContainer.h"
35 #include "AliVParticle.h"
40 ClassImp(AliTwoPlusOneContainer)
42 AliTwoPlusOneContainer::AliTwoPlusOneContainer(const char* name, const char* binning, Double_t alpha) :
56 TString defaultBinningStr;
57 defaultBinningStr = "eta: -1.0, 1.0\n"//needs to be defined for AliUEHist
58 "p_t_assoc: 0.5, 1.0, 2.0, 3.0, 4.0, 6.0, 8.0\n"
59 "p_t_leading: 6.0, 8.0, 10.0, 12.0, 14.0\n"
60 "p_t_leading_course: 4.0, 6.0, 8.0, 10.0\n"
61 "p_t_eff: 0.5, 8.0\n"//needs to be defined for AliUEHist
62 "vertex_eff: -10, 10\n"//needs to be defined for AliUEHist
63 "multiplicity: 0, 1, 2, 3, 4, 5, 10, 20, 30, 40, 50, 60, 70, 80, 90, 100.1\n"
64 "vertex: -7, -5, -3, -1, 1, 3, 5, 7\n"
65 "delta_phi: -1.570796, -1.396263, -1.221730, -1.047197, -0.872664, -0.698131, -0.523598, -0.349065, -0.174532, 0.0, 0.174532, 0.349065, 0.523598, 0.698131, 0.872664, 1.047197, 1.221730, 1.396263, 1.570796, 1.745329, 1.919862, 2.094395, 2.268928, 2.443460, 2.617993, 2.792526, 2.967059, 3.141592, 3.316125, 3.490658, 3.665191, 3.839724, 4.014257, 4.188790, 4.363323, 4.537856, 4.712388\n" // this binning starts at -pi/2 and it ends at +3pi/2, it has 36 bins
66 "delta_eta: -1.8, -1.6, -1.4, -1.2, -1.0, -0.8, -0.6, -0.4, -0.2, 0, 0.2, 0.4, 0.6, 0.8, 1.0, 1.2, 1.4, 1.6, 1.8\n";
68 // combine customBinning with defaultBinningStr -> use customBinning where available and otherwise defaultBinningStr
69 TString binningStr = AliUEHist::CombineBinning(defaultBinningStr, TString(binning));
71 fTwoPlusOne = new AliUEHist("TwoPlusOne", binningStr);
73 //set minimum and maximum trigger pt values
74 fTriggerPt1Min = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(2)->GetXmin();
75 fTriggerPt1Max = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(2)->GetXmax();
76 fTriggerPt2Min = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(6)->GetXmin();
77 fTriggerPt2Max = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(6)->GetXmax();
78 fPtAssocMin = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(1)->GetXmin();
79 fPtAssocMax = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(1)->GetXmax();
84 //_____________________________________________________________________________
85 AliTwoPlusOneContainer::AliTwoPlusOneContainer(const AliTwoPlusOneContainer &c) :
86 TNamed(fName, fTitle),
98 // AliTwoPlusOneContainer copy constructor
101 ((AliTwoPlusOneContainer &) c).Copy(*this);
102 fTriggerPt1Min = ((AliTwoPlusOneContainer &) c).getTriggerPt1Min();
103 fTriggerPt1Max = ((AliTwoPlusOneContainer &) c).getTriggerPt1Max();
104 fTriggerPt2Min = ((AliTwoPlusOneContainer &) c).getTriggerPt2Min();
105 fTriggerPt2Max = ((AliTwoPlusOneContainer &) c).getTriggerPt2Max();
106 fPtAssocMin = ((AliTwoPlusOneContainer &) c).getPtAssocMin();
107 fPtAssocMax = ((AliTwoPlusOneContainer &) c).getPtAssocMax();
108 fAlpha = ((AliTwoPlusOneContainer &) c).fAlpha;
109 fMergeCount = ((AliTwoPlusOneContainer &) c).fMergeCount;
112 //____________________________________________________________________
113 AliTwoPlusOneContainer::~AliTwoPlusOneContainer()
120 //____________________________________________________________________
121 void AliTwoPlusOneContainer::DeleteContainers()
131 //____________________________________________________________________
132 void AliTwoPlusOneContainer::FillCorrelations(Double_t centrality, Float_t zVtx, AliTwoPlusOneContainer::PlotKind step, TObjArray* triggerNear, TObjArray* triggerAway, TObjArray* assocNear, TObjArray* assocAway, Double_t weight, Bool_t is1plus1, Bool_t isBackgroundSame)
134 //Fill Correlations fills the UEHist fTwoPlusOne with the 2+1 correlation
135 //the input variables centrality and zVtx are the centrality and the z vertex of the event
136 //the 4 lists triggerNear, triggerAway, assocNear and assocAway are four lists of particles. For the same event analysis all for lists are identical. For the mixed event analysis assocNear and assocAway are from different events as the trigger lists and for the mixed combinatorics triggerNear and assocNear are from one event and triggerAway and assocAway are from another event.
138 AliCFContainer* track_hist = fTwoPlusOne->GetTrackHist(AliUEHist::kToward);
139 AliCFContainer* event_hist = fTwoPlusOne->GetEventHist();
140 AliUEHist::CFStep stepUEHist = static_cast<AliUEHist::CFStep>(step);
142 //in case of the computation of the background in the same event there are two possible positions: delta phi = +/- pi/2
143 //both positions are used so the results could only be weighted with 0.5*weight
147 for (Int_t i=0; i<triggerNear->GetEntriesFast(); i++){
148 AliVParticle* part = (AliVParticle*) triggerNear->UncheckedAt(i);
150 Double_t part_pt = part->Pt();
151 if(part_pt<fTriggerPt1Min || part_pt>fTriggerPt1Max)
154 Double_t part_eta = part->Eta();
155 Double_t part_phi = part->Phi();
157 //search for second trigger particle
158 //entries is the number of particles in the event, so there can only be less particles within the trigger 2 conditions
160 Int_t ind_max_found_pt = -1;
161 Int_t triggerAway_entries = triggerAway->GetEntriesFast();
162 AliVParticle* found_particle[triggerAway_entries];
164 //have to fake the away side triggers for the 1+1 analysis
166 found_particle[ind_found] = part;//in 1plus1 use first trigger particle also as pseudo second trigger particle
167 ind_max_found_pt = ind_found;
170 //normal 2+1 analysis
171 for (Int_t j=0; j<triggerAway_entries; j++){
172 AliVParticle* part2 = (AliVParticle*) triggerAway->UncheckedAt(j);
174 Double_t part2_pt = part2->Pt();
175 //check if pT of trigger 2 is within the trigger range
176 //also pt of trigger 2 needs to be smaller than the pt of trigger 1 (to have an ordering if both pt are close to each other)
177 if(part2_pt<fTriggerPt2Min || part2_pt>fTriggerPt2Max || part2_pt>part_pt)
180 // don't use the same particle (is in any case impossible because the Delta phi angle will be 0)
185 Double_t dphi_triggers = part_phi-part2->Phi();
187 if(dphi_triggers>1.5*TMath::Pi()) dphi_triggers -= TMath::TwoPi();
188 else if(dphi_triggers<-0.5*TMath::Pi()) dphi_triggers += TMath::TwoPi();
190 //if 2+1 analysis check if trigger particles have a delta phi = pi +/- alpha
191 if(!isBackgroundSame)
192 dphi_triggers -= TMath::Pi();
193 else if(isBackgroundSame){
194 //shift defined area of delta phi
195 if(dphi_triggers>TMath::Pi()) dphi_triggers -= TMath::TwoPi();
197 //look at delta phi = +/- pi/2
199 dphi_triggers += 0.5*TMath::Pi();
200 else if(dphi_triggers>0)
201 dphi_triggers -= 0.5*TMath::Pi();
203 if(TMath::Abs(dphi_triggers)>fAlpha)
206 found_particle[ind_found] = part2;
207 if(ind_max_found_pt==-1 || part2_pt>found_particle[ind_max_found_pt]->Pt()) ind_max_found_pt = ind_found;
209 }//end loop to search for the second trigger particle
212 //if no second trigger particle was found continue to search for the next first trigger particle
216 //the energy of the second trigger particle is set for the near side to the maximum energy of all trigger 2 particles on the away side
217 // this leads to the fact that the number of accepted trigger combinations can be artificial smaller than the real number if there is a cut on the pT 2 energy from the top; cutting away the smallest energy of pT 2 is still save; this is the reason why it is not allowed to use a cut on the top pt of trigger particle 2
218 //fill trigger particles
222 vars[1] = centrality;
224 vars[3] = found_particle[ind_max_found_pt]->Pt();
226 vars[3] = (fTriggerPt2Max+fTriggerPt2Min)/2;
228 event_hist->Fill(vars, stepUEHist, weight);//near side
231 for(Int_t k=0; k< ind_found; k++){
232 vars[3] = found_particle[k]->Pt();
233 event_hist->Fill(vars, stepUEHist+1, weight);//away side
237 //add correlated particles on the near side
238 for (Int_t k=0; k<assocNear->GetEntriesFast(); k++){
239 AliVParticle* part3 = (AliVParticle*) assocNear->UncheckedAt(k);
241 Double_t part3_pt = part3->Pt();
242 if(part3_pt<fPtAssocMin || part3_pt>fPtAssocMax)
245 //do not add the trigger 1 particle
249 //do not add the trigger 2 particle
250 if(found_particle[ind_max_found_pt]==part3)
253 Double_t dphi_near = part_phi-part3->Phi();
254 if(dphi_near>1.5*TMath::Pi()) dphi_near -= TMath::TwoPi();
255 else if(dphi_near<-0.5*TMath::Pi()) dphi_near += TMath::TwoPi();
257 Double_t deta_near = part_eta-part3->Eta();
263 vars[3] = centrality;
266 vars[6] = found_particle[ind_max_found_pt]->Pt();
268 track_hist->Fill(vars, stepUEHist, weight);
271 //search only for the distribution of the 2nd trigger particle
275 //add correlated particles on the away side
276 for (Int_t k=0; k<assocAway->GetEntriesFast(); k++){
277 AliVParticle* part3 = (AliVParticle*) assocAway->UncheckedAt(k);
279 Double_t part3_pt = part3->Pt();
280 if(part3_pt<fPtAssocMin || part3_pt>fPtAssocMax)
283 //do not add the trigger 1 particle
287 for(int l=0; l<ind_found; l++){
288 //do not add the trigger 2 particle
289 if(found_particle[l]==part3)
292 Double_t dphi_away = found_particle[l]->Phi()-part3->Phi();
293 if(dphi_away>1.5*TMath::Pi()) dphi_away -= TMath::TwoPi();
294 else if(dphi_away<-0.5*TMath::Pi()) dphi_away += TMath::TwoPi();
296 Double_t deta_away = found_particle[l]->Eta()-part3->Eta();
302 vars[3] = centrality;
305 vars[6] = found_particle[l]->Pt();
307 track_hist->Fill(vars, stepUEHist+1, weight);//step +1 is the AS to the NS plot of step
310 }//end loop to search for the first trigger particle
315 //____________________________________________________________________
316 AliTwoPlusOneContainer &AliTwoPlusOneContainer::operator=(const AliTwoPlusOneContainer &c)
318 // assigment operator
323 ((AliTwoPlusOneContainer &) c).Copy(*this);
328 //____________________________________________________________________
329 void AliTwoPlusOneContainer::Copy(TObject& c) const
333 AliTwoPlusOneContainer& target = (AliTwoPlusOneContainer &) c;
336 target.fTwoPlusOne = dynamic_cast<AliUEHist*> (fTwoPlusOne->Clone());
340 //____________________________________________________________________
341 Long64_t AliTwoPlusOneContainer::Merge(TCollection* list)
343 // Merge a list of AliTwoPlusOneContainer objects with this.
344 // Returns the number of merged objects (including this).
352 TIterator* iter = list->MakeIterator();
355 // collections of objects
356 TList* lists = new TList;
359 while ((obj = iter->Next())) {
361 AliTwoPlusOneContainer* entry = dynamic_cast<AliTwoPlusOneContainer*> (obj);
365 lists->Add(entry->fTwoPlusOne);
367 fMergeCount += entry->fMergeCount;
371 fTwoPlusOne->Merge(lists);