]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/Base/AliTwoPlusOneContainer.cxx
Merge branch 'master' into TRDdev
[u/mrichter/AliRoot.git] / PWGCF / Correlations / Base / AliTwoPlusOneContainer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 //
17 //
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
24 //
25 //
26 // Author: Markus Zimmermann
27 // 
28
29
30 #include "AliTwoPlusOneContainer.h"
31
32 #include "TList.h"
33
34 #include "AliCFContainer.h"
35 #include "AliVParticle.h"
36
37 #include "TH1F.h"
38 #include "TMath.h"
39
40 ClassImp(AliTwoPlusOneContainer)
41
42 AliTwoPlusOneContainer::AliTwoPlusOneContainer(const char* name, const char* binning, Double_t alpha) : 
43   TNamed(name, name),
44   fTwoPlusOne(0),
45   fTriggerPt1Min(0),
46   fTriggerPt1Max(0),
47   fTriggerPt2Min(0),
48   fTriggerPt2Max(0),
49   fPtAssocMin(0),
50   fPtAssocMax(0),
51   fAlpha(alpha),
52   fMergeCount(0)
53 {
54   // Constructor
55   //
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: -2.0, -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, 2.0\n";
67
68   // combine customBinning with defaultBinningStr -> use customBinning where available and otherwise defaultBinningStr
69   TString binningStr = AliUEHist::CombineBinning(defaultBinningStr, TString(binning));
70
71   fTwoPlusOne = new AliUEHist("TwoPlusOne", binningStr);
72
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();
80
81   TH1::AddDirectory();
82 }
83
84 //_____________________________________________________________________________
85 AliTwoPlusOneContainer::AliTwoPlusOneContainer(const AliTwoPlusOneContainer &c) :
86   TNamed(fName, fTitle),
87   fTwoPlusOne(0),
88   fTriggerPt1Min(0),
89   fTriggerPt1Max(0),
90   fTriggerPt2Min(0),
91   fTriggerPt2Max(0),
92   fPtAssocMin(0),
93   fPtAssocMax(0),
94   fAlpha(0.2),
95   fMergeCount(0)
96 {
97   //
98   // AliTwoPlusOneContainer copy constructor
99   //
100
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;
110 }
111
112 //____________________________________________________________________
113 AliTwoPlusOneContainer::~AliTwoPlusOneContainer()
114 {
115   // Destructor
116   
117   DeleteContainers();
118 }
119
120 //____________________________________________________________________
121 void AliTwoPlusOneContainer::DeleteContainers()
122 {
123   if (fTwoPlusOne)
124   {
125     delete fTwoPlusOne;
126     fTwoPlusOne = 0;
127   }
128 }
129
130
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)
133 {
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.
137
138   AliCFContainer* track_hist = fTwoPlusOne->GetTrackHist(AliUEHist::kToward);
139   AliCFContainer* event_hist = fTwoPlusOne->GetEventHist();
140   AliUEHist::CFStep stepUEHist = static_cast<AliUEHist::CFStep>(step);
141
142   for (Int_t i=0; i<triggerNear->GetEntriesFast(); i++){
143     AliVParticle* part = (AliVParticle*) triggerNear->UncheckedAt(i);
144     
145     Double_t part_pt = part->Pt();
146     if(part_pt<fTriggerPt1Min || part_pt>fTriggerPt1Max)
147       continue;
148
149     Double_t part_eta = part->Eta();
150     Double_t part_phi = part->Phi();
151
152     //search for second trigger particle
153     //entries is the number of particles in the event, so there can only be less particles within the trigger 2 conditions
154     Int_t ind_found = 0;
155     Int_t ind_max_found_pt = -1;
156     Int_t triggerAway_entries = triggerAway->GetEntriesFast();
157     AliVParticle* found_particle[triggerAway_entries];
158
159     for (Int_t j=0; j<triggerAway_entries; j++){
160       AliVParticle* part2 = (AliVParticle*) triggerAway->UncheckedAt(j);
161
162       Double_t part2_pt = part2->Pt();
163       //check if pT of trigger 2 is within the trigger range
164       //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)
165       if(part2_pt<fTriggerPt2Min || part2_pt>fTriggerPt2Max || part2_pt>part_pt)
166         continue;
167       
168        // don't use the same particle (is in any case impossible because the Delta phi angle will be 0)
169       if(part==part2){
170         continue;
171       }
172         
173       Double_t dphi_triggers = part_phi-part2->Phi();
174
175       if(dphi_triggers>1.5*TMath::Pi()) dphi_triggers -= TMath::TwoPi();
176       else if(dphi_triggers<-0.5*TMath::Pi()) dphi_triggers += TMath::TwoPi();
177         
178       dphi_triggers -= TMath::Pi();
179       if(TMath::Abs(dphi_triggers)>fAlpha)
180         continue;
181
182       found_particle[ind_found] = part2;
183       if(ind_max_found_pt==-1 || part2_pt>found_particle[ind_max_found_pt]->Pt()) ind_max_found_pt = ind_found;
184       ind_found++;
185     }//end loop to search for the second trigger particle
186
187     //if no second trigger particle was found continue to search for the next first trigger particle
188     if(ind_found==0)
189       continue;
190
191     //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
192     // 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
193     //fill trigger particles
194     if(ind_found>0){
195         Double_t vars[4];
196         vars[0] = part_pt;
197         vars[1] = centrality;
198         vars[2] = zVtx;
199         vars[3] = found_particle[ind_max_found_pt]->Pt();
200
201         event_hist->Fill(vars, stepUEHist, weight);//near side
202         
203         for(Int_t k=0; k< ind_found; k++){
204           vars[3] = found_particle[k]->Pt();
205           event_hist->Fill(vars, stepUEHist+1, weight);//away side
206         }
207     }
208
209     //add correlated particles on the near side
210     for (Int_t k=0; k<assocNear->GetEntriesFast(); k++){
211       AliVParticle* part3 = (AliVParticle*) assocNear->UncheckedAt(k);
212
213       Double_t part3_pt = part3->Pt();
214       if(part3_pt<fPtAssocMin || part3_pt>fPtAssocMax)
215         continue;
216
217       //do not add the trigger 1 particle
218       if(part==part3)
219         continue;
220
221       //do not add the trigger 2 particle
222       if(found_particle[ind_max_found_pt]==part3)
223         continue;
224
225       Double_t dphi_near = part_phi-part3->Phi(); 
226       if(dphi_near>1.5*TMath::Pi()) dphi_near -= TMath::TwoPi();
227       else if(dphi_near<-0.5*TMath::Pi()) dphi_near += TMath::TwoPi();
228
229       Double_t deta_near = part_eta-part3->Eta();
230       
231       Double_t vars[7];
232       vars[0] = deta_near;
233       vars[1] = part3_pt;
234       vars[2] = part_pt;
235       vars[3] = centrality;
236       vars[4] = dphi_near;
237       vars[5] = zVtx;
238       vars[6] = found_particle[ind_max_found_pt]->Pt();
239
240       track_hist->Fill(vars, stepUEHist, weight);
241     }
242
243     //add correlated particles on the away side
244     for (Int_t k=0; k<assocAway->GetEntriesFast(); k++){
245       AliVParticle* part3 = (AliVParticle*) assocAway->UncheckedAt(k);
246         
247       Double_t part3_pt = part3->Pt();
248       if(part3_pt<fPtAssocMin || part3_pt>fPtAssocMax)
249         continue;
250
251       //do not add the trigger 1 particle
252       if(part==part3)
253         continue;
254
255       for(int l=0; l<ind_found; l++){
256         //do not add the trigger 2 particle
257         if(found_particle[l]==part3)
258           continue;
259
260         Double_t dphi_away = found_particle[l]->Phi()-part3->Phi();
261         if(dphi_away>1.5*TMath::Pi()) dphi_away -= TMath::TwoPi();
262         else if(dphi_away<-0.5*TMath::Pi()) dphi_away += TMath::TwoPi();
263         
264         Double_t deta_away = found_particle[l]->Eta()-part3->Eta();
265       
266         Double_t vars[7];
267         vars[0] = deta_away;
268         vars[1] = part3_pt;
269         vars[2] = part_pt;
270         vars[3] = centrality;
271         vars[4] = dphi_away;
272         vars[5] = zVtx;
273         vars[6] = found_particle[l]->Pt();
274         
275         track_hist->Fill(vars, stepUEHist+1, weight);//step +1 is the AS to the NS plot of step
276       }
277     }
278   }//end loop to search for the first trigger particle
279 }
280
281
282
283 //____________________________________________________________________
284 AliTwoPlusOneContainer &AliTwoPlusOneContainer::operator=(const AliTwoPlusOneContainer &c)
285 {
286   // assigment operator
287
288   DeleteContainers();
289
290   if (this != &c)
291     ((AliTwoPlusOneContainer &) c).Copy(*this);
292
293   return *this;
294 }
295
296 //____________________________________________________________________
297 void AliTwoPlusOneContainer::Copy(TObject& c) const
298 {
299   // copy function
300
301   AliTwoPlusOneContainer& target = (AliTwoPlusOneContainer &) c;
302
303   if (fTwoPlusOne)
304     target.fTwoPlusOne = dynamic_cast<AliUEHist*> (fTwoPlusOne->Clone());
305
306 }
307
308 //____________________________________________________________________
309 Long64_t AliTwoPlusOneContainer::Merge(TCollection* list)
310 {
311   // Merge a list of AliTwoPlusOneContainer objects with this. 
312   // Returns the number of merged objects (including this).
313
314   if (!list)
315     return 0;
316
317   if (list->IsEmpty())
318     return 1;
319
320   TIterator* iter = list->MakeIterator();
321   TObject* obj;
322
323   // collections of objects
324   TList* lists = new TList;
325
326   Int_t count = 0;
327   while ((obj = iter->Next())) {
328     
329     AliTwoPlusOneContainer* entry = dynamic_cast<AliTwoPlusOneContainer*> (obj);
330     if (entry == 0) 
331       continue;
332
333     lists->Add(entry->fTwoPlusOne);
334
335     fMergeCount += entry->fMergeCount;
336     count++;
337   }
338   
339   fTwoPlusOne->Merge(lists);
340
341   delete lists;
342   return count+1;
343
344 }