]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliAnalysisTaskMultiparticleCorrelations.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWG / FLOW / Tasks / AliAnalysisTaskMultiparticleCorrelations.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2008, 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  * analysis task for flow analysis with *
18  *     multi-particle correlations      * 
19  *                                      * 
20  * author: Ante Bilandzic               *
21  *         (abilandzic@gmail.com)       * 
22  ***************************************/
23   
24 #include "Riostream.h"
25 #include "AliFlowEventSimple.h"
26 #include "AliAnalysisTaskMultiparticleCorrelations.h"
27 #include "AliFlowAnalysisWithMultiparticleCorrelations.h"
28 #include "AliLog.h"
29
30 using std::cout;
31 using std::endl;
32
33 ClassImp(AliAnalysisTaskMultiparticleCorrelations)
34
35 //================================================================================================================
36
37 AliAnalysisTaskMultiparticleCorrelations::AliAnalysisTaskMultiparticleCorrelations(const char *name, Bool_t useParticleWeights): 
38  AliAnalysisTaskSE(name), 
39  fEvent(NULL),
40  fMPC(NULL), 
41  fHistList(NULL),
42  fUseInternalFlags(kFALSE),
43  fMinNoRPs(-44),
44  fMaxNoRPs(-44),
45  fExactNoRPs(-44),
46  fAnalysisTag(""),
47  fDumpThePoints(kFALSE),
48  fMaxNoEventsPerFile(100),
49  fFillControlHistograms(kFALSE),
50  fFillKinematicsHist(kFALSE),
51  fFillMultDistributionsHist(kFALSE),
52  fFillMultCorrelationsHist(kFALSE),
53  fCalculateQvector(kFALSE),
54  fCalculateDiffQvectors(kFALSE),
55  fProduction(""),
56  fCalculateCorrelations(kFALSE),
57  fCalculateIsotropic(kFALSE),
58  fCalculateSame(kFALSE),
59  fSkipZeroHarmonics(kFALSE),
60  fCalculateSameIsotropic(kFALSE),
61  fCalculateAll(kFALSE),
62  fDontGoBeyond(0),
63  fCalculateOnlyForHarmonicQC(kFALSE),
64  fCalculateOnlyForSC(kFALSE),
65  fCalculateOnlyCos(kFALSE),
66  fCalculateOnlySin(kFALSE),
67  fCalculateEbECumulants(kFALSE),
68  fCrossCheckWithNestedLoops(kFALSE),
69  fCrossCheckDiffWithNestedLoops(kFALSE),
70  fCalculateStandardCandles(kFALSE),
71  fPropagateErrorSC(kTRUE),
72  fCalculateQcumulants(kFALSE),
73  fHarmonicQC(2),
74  fPropagateErrorQC(kTRUE),
75  fCalculateDiffCorrelations(kFALSE),
76  fCalculateDiffCos(kTRUE),
77  fCalculateDiffSin(kFALSE),
78  fCalculateDiffCorrelationsVsPt(kTRUE),
79  fUseDefaultBinning(kTRUE),
80  fnDiffBins(-44),
81  fRangesDiffBins(NULL)
82  {
83   // Constructor.
84  
85   AliDebug(2,"AliAnalysisTaskMultiparticleCorrelations::AliAnalysisTaskMultiparticleCorrelations(const char *name, Bool_t useParticleWeights)");
86
87   // Define input and output slots here
88   // Input slot #0 works with an AliFlowEventSimple
89   DefineInput(0, AliFlowEventSimple::Class());  
90   // Input slot #1 is needed for the weights input file:
91   if(useParticleWeights)
92   {
93    DefineInput(1, TList::Class());   
94   }  
95   // Output slot #0 is reserved              
96   // Output slot #1 writes into a TList container
97   DefineOutput(1, TList::Class());  
98
99   // Initialize all arrays:
100   for(Int_t rp=0;rp<2;rp++) // [RP,POI]
101   {
102    for(Int_t ppe=0;ppe<3;ppe++) // [phi,pt,eta]
103    {
104     fUseWeights[rp][ppe] = kFALSE;
105     fWeightsHist[rp][ppe] = NULL; 
106    }
107   }
108   // For nested loops arrays:
109   fCrossCheckDiffCSCOBN[0] = 0; // cos/sin
110   fCrossCheckDiffCSCOBN[1] = 2; // correlator order
111   fCrossCheckDiffCSCOBN[2] = 4; // bin number 
112   // Initialize default vaulues for fDontFill[3]:
113   fDontFill[0] = kFALSE;
114   fDontFill[1] = kFALSE;
115   fDontFill[2] = kFALSE;
116   // Initialize default binning values for fKinematicsHist[2][3]:
117   // nBins:
118   fnBins[0][0] = 360;  // [RP][phi]
119   fnBins[0][1] = 1000; // [RP][pt]
120   fnBins[0][2] = 1000; // [RP][eta]
121   fnBins[1][0] = 360;  // [POI][phi]
122   fnBins[1][1] = 1000; // [POI][pt]
123   fnBins[1][2] = 1000; // [POI][eta]
124   // Min:
125   fMin[0][0] = 0.;  // [RP][phi]
126   fMin[0][1] = 0.;  // [RP][pt]
127   fMin[0][2] = -1.; // [RP][eta]
128   fMin[1][0] = 0.;  // [POI][phi]
129   fMin[1][1] = 0.;  // [POI][pt]
130   fMin[1][2] = -1.; // [POI][eta]
131   // Max:
132   fMax[0][0] = TMath::TwoPi(); // [RP][phi]
133   fMax[0][1] = 10.;            // [RP][pt]
134   fMax[0][2] = 1.;             // [RP][eta]
135   fMax[1][0] = TMath::TwoPi(); // [POI][phi]
136   fMax[1][1] = 10.;            // [POI][pt]
137   fMax[1][2] = 1.;             // [POI][eta]
138   // Initialize default binning values for fMultCorrelationsHist[3]:
139   // nBins:
140   fnBinsMult[0] = 3000; // [RP]
141   fnBinsMult[1] = 3000; // [POI]
142   fnBinsMult[2] = 3000; // [REF]
143   // Min:
144   fMinMult[0] = 0.; // [RP]
145   fMinMult[1] = 0.; // [POI]
146   fMinMult[2] = 0.; // [REF]
147   // Max:
148   fMaxMult[0] = 3000.; // [RP]
149   fMaxMult[1] = 3000.; // [POI]
150   fMaxMult[2] = 3000.; // [REF]
151
152 } // AliAnalysisTaskMultiparticleCorrelations::AliAnalysisTaskMultiparticleCorrelations(const char *name, Bool_t useParticleWeights): 
153
154 //================================================================================================================
155
156 AliAnalysisTaskMultiparticleCorrelations::AliAnalysisTaskMultiparticleCorrelations(): 
157  AliAnalysisTaskSE(),
158  fEvent(NULL),
159  fMPC(NULL),
160  fHistList(NULL),
161  fUseInternalFlags(kFALSE),
162  fMinNoRPs(-44),
163  fMaxNoRPs(-44),
164  fExactNoRPs(-44),
165  fAnalysisTag(""),
166  fDumpThePoints(kFALSE),
167  fMaxNoEventsPerFile(0),
168  fFillControlHistograms(kFALSE),
169  fFillKinematicsHist(kFALSE),
170  fFillMultDistributionsHist(kFALSE),
171  fFillMultCorrelationsHist(kFALSE),
172  fCalculateQvector(kFALSE),
173  fCalculateDiffQvectors(kFALSE),
174  fProduction(""),
175  fCalculateCorrelations(kFALSE),
176  fCalculateIsotropic(kFALSE),
177  fCalculateSame(kFALSE),
178  fSkipZeroHarmonics(kFALSE),
179  fCalculateSameIsotropic(kFALSE),
180  fCalculateAll(kFALSE),
181  fDontGoBeyond(0),
182  fCalculateOnlyForHarmonicQC(kFALSE),
183  fCalculateOnlyForSC(kFALSE),
184  fCalculateOnlyCos(kFALSE),
185  fCalculateOnlySin(kFALSE),
186  fCalculateEbECumulants(kFALSE),
187  fCrossCheckWithNestedLoops(kFALSE),
188  fCrossCheckDiffWithNestedLoops(kFALSE),
189  fCalculateStandardCandles(kFALSE),
190  fPropagateErrorSC(kFALSE),
191  fCalculateQcumulants(kFALSE),
192  fHarmonicQC(0),
193  fPropagateErrorQC(kFALSE),
194  fCalculateDiffCorrelations(kFALSE),
195  fCalculateDiffCos(kTRUE),
196  fCalculateDiffSin(kFALSE),
197  fCalculateDiffCorrelationsVsPt(kTRUE),
198  fUseDefaultBinning(kTRUE),
199  fnDiffBins(-44),
200  fRangesDiffBins(NULL)
201  {
202   // Dummy constructor.
203  
204   AliDebug(2,"AliAnalysisTaskMultiparticleCorrelations::AliAnalysisTaskMultiparticleCorrelations()");
205
206   // Initialize all arrays:
207   for(Int_t rp=0;rp<2;rp++) // [RP,POI]
208   {
209    for(Int_t ppe=0;ppe<3;ppe++) // [phi,pt,eta]
210    {
211     fUseWeights[rp][ppe] = kFALSE;
212     fWeightsHist[rp][ppe] = NULL; 
213    }
214   }
215   // For nested loops arrays:
216   fCrossCheckDiffCSCOBN[0] = 0; // cos/sin
217   fCrossCheckDiffCSCOBN[1] = 2; // correlator order
218   fCrossCheckDiffCSCOBN[2] = 4; // bin number 
219   // Initialize default vaulues for fDontFill[3]:
220   fDontFill[0] = kFALSE;
221   fDontFill[1] = kFALSE;
222   fDontFill[2] = kFALSE;
223   // Initialize default binning values for fKinematicsHist[2][3]:
224   // nBins:
225   fnBins[0][0] = 360;  // [RP][phi]
226   fnBins[0][1] = 1000; // [RP][pt]
227   fnBins[0][2] = 1000; // [RP][eta]
228   fnBins[1][0] = 360;  // [POI][phi]
229   fnBins[1][1] = 1000; // [POI][pt]
230   fnBins[1][2] = 1000; // [POI][eta]
231   // Min:
232   fMin[0][0] = 0.;  // [RP][phi]
233   fMin[0][1] = 0.;  // [RP][pt]
234   fMin[0][2] = -1.; // [RP][eta]
235   fMin[1][0] = 0.;  // [POI][phi]
236   fMin[1][1] = 0.;  // [POI][pt]
237   fMin[1][2] = -1.; // [POI][eta]
238   // Max:
239   fMax[0][0] = TMath::TwoPi(); // [RP][phi]
240   fMax[0][1] = 10.;            // [RP][pt]
241   fMax[0][2] = 1.;             // [RP][eta]
242   fMax[1][0] = TMath::TwoPi(); // [POI][phi]
243   fMax[1][1] = 10.;            // [POI][pt]
244   fMax[1][2] = 1.;             // [POI][eta]
245   // Initialize default binning values for fMultCorrelationsHist[3]:
246   // nBins:
247   fnBinsMult[0] = 3000; // [RP]
248   fnBinsMult[1] = 3000; // [POI]
249   fnBinsMult[2] = 3000; // [REF]
250   // Min:
251   fMinMult[0] = 0.; // [RP]
252   fMinMult[1] = 0.; // [POI]
253   fMinMult[2] = 0.; // [REF]
254   // Max:
255   fMaxMult[0] = 3000.; // [RP]
256   fMaxMult[1] = 3000.; // [POI]
257   fMaxMult[2] = 3000.; // [REF]
258
259 } // AliAnalysisTaskMultiparticleCorrelations::AliAnalysisTaskMultiparticleCorrelations():
260
261 //================================================================================================================
262
263 void AliAnalysisTaskMultiparticleCorrelations::UserCreateOutputObjects() 
264 {
265  // Called at every worker node to initialize.
266   
267  AliDebug(2,"AliAnalysisTaskMultiparticleCorrelations::UserCreateOutputObjects()");
268  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::UserCreateOutputObjects()";
269
270  // Analyser:
271  fMPC = new AliFlowAnalysisWithMultiparticleCorrelations();
272  
273  // Setters:
274  if(fUseInternalFlags){fMPC->SetMinNoRPs(fMinNoRPs);}
275  if(fUseInternalFlags){fMPC->SetMaxNoRPs(fMaxNoRPs);}
276  if(fUseInternalFlags){fMPC->SetExactNoRPs(fExactNoRPs);}
277  fMPC->SetAnalysisTag(fAnalysisTag.Data());
278  fMPC->SetDumpThePoints(fDumpThePoints,fMaxNoEventsPerFile);
279  fMPC->SetFillControlHistograms(fFillControlHistograms);
280  if(fDontFill[0]){fMPC->SetDontFill("RP");}
281  if(fDontFill[1]){fMPC->SetDontFill("POI");}
282  if(fDontFill[2]){fMPC->SetDontFill("REF");}
283  fMPC->SetFillKinematicsHist(fFillKinematicsHist);
284  fMPC->SetFillMultDistributionsHist(fFillMultDistributionsHist);
285  fMPC->SetFillMultCorrelationsHist(fFillMultCorrelationsHist);
286  fMPC->SetCalculateQvector(fCalculateQvector);
287  fMPC->SetCalculateDiffQvectors(fCalculateDiffQvectors);
288  fMPC->SetCalculateCorrelations(fCalculateCorrelations);
289  fMPC->SetCalculateIsotropic(fCalculateIsotropic);
290  fMPC->SetCalculateSame(fCalculateSame);
291  fMPC->SetSkipZeroHarmonics(fSkipZeroHarmonics);
292  fMPC->SetCalculateSameIsotropic(fCalculateSameIsotropic);
293  fMPC->SetCalculateAll(fCalculateAll);
294  fMPC->SetDontGoBeyond(fDontGoBeyond);
295  fMPC->SetCalculateOnlyForHarmonicQC(fCalculateOnlyForHarmonicQC);
296  fMPC->SetCalculateOnlyForSC(fCalculateOnlyForSC);
297  fMPC->SetCalculateOnlyCos(fCalculateOnlyCos);
298  fMPC->SetCalculateOnlySin(fCalculateOnlySin);
299  fMPC->SetCalculateEbECumulants(fCalculateEbECumulants);
300  fMPC->SetCrossCheckWithNestedLoops(fCrossCheckWithNestedLoops);
301  fMPC->SetCrossCheckDiffWithNestedLoops(fCrossCheckDiffWithNestedLoops);
302  fMPC->SetCrossCheckDiffCSCOBN(fCrossCheckDiffCSCOBN[0],fCrossCheckDiffCSCOBN[1],fCrossCheckDiffCSCOBN[2]);  
303  fMPC->SetCalculateStandardCandles(fCalculateStandardCandles);
304  fMPC->SetPropagateErrorSC(fPropagateErrorSC);
305  fMPC->SetCalculateQcumulants(fCalculateQcumulants);
306  fMPC->SetHarmonicQC(fHarmonicQC);
307  fMPC->SetPropagateErrorQC(fPropagateErrorQC);
308  fMPC->SetCalculateDiffCorrelations(fCalculateDiffCorrelations);
309  fMPC->SetCalculateDiffCos(fCalculateDiffCos);
310  fMPC->SetCalculateDiffSin(fCalculateDiffSin);
311  fMPC->SetCalculateDiffCorrelationsVsPt(fCalculateDiffCorrelationsVsPt);
312  fMPC->SetUseDefaultBinning(fUseDefaultBinning);
313  fMPC->SetnDiffBins(fnDiffBins);
314  fMPC->SetRangesDiffBins(fRangesDiffBins);
315
316  // Weights:
317  TString type[2] = {"RP","POI"};
318  TString variable[3] = {"phi","pt","eta"}; 
319  for(Int_t rp=0;rp<2;rp++) // [RP,POI]
320  {
321   for(Int_t ppe=0;ppe<3;ppe++) // [phi,pt,eta]
322   {
323    if(fUseWeights[rp][ppe])
324    {
325     if(!fWeightsHist[rp][ppe])
326     {
327      if(fProduction.EqualTo("")){Fatal(sMethodName.Data(),"fProduction is empty, for one reason or another...");}
328      fWeightsHist[rp][ppe] = GetHistogramWithWeights(TString(Form("%s/%s",gSystem->pwd(),"weights.root")).Data(),TString(this->fName).Data(),type[rp].Data(),variable[ppe].Data(),fProduction.Data());
329     }
330     if(!fWeightsHist[rp][ppe])
331     {
332      Fatal(sMethodName.Data(),"fWeightsHist[%d][%d]",rp,ppe);
333     } else{fMPC->SetWeightsHist(fWeightsHist[rp][ppe],type[rp].Data(),variable[ppe].Data());}
334    } // if(fUseWeights[rp][ppe])
335   } // for(Int_t ppe=0;ppe<3;ppe++) // [phi,pt,eta]
336  } // for(Int_t rp=0;rp<2;rp++) // [RP,POI]
337
338  // Control histos:
339  // Kinematics:
340  //TString typeKine[2] = {"RP","POI"};
341  //TString variable[3] = {"phi","pt","eta"};
342  for(Int_t rp=0;rp<2;rp++) // [RP,POI]
343  {
344   for(Int_t ppe=0;ppe<3;ppe++) // [phi,pt,eta]
345   {
346    fMPC->SetnBins(type[rp].Data(),variable[ppe].Data(),fnBins[rp][ppe]);
347    fMPC->SetMin(type[rp].Data(),variable[ppe].Data(),fMin[rp][ppe]);
348    fMPC->SetMax(type[rp].Data(),variable[ppe].Data(),fMax[rp][ppe]);
349   } // for(Int_t ppe=0;ppe<3;ppe++) // [phi,pt,eta]
350  } // for(Int_t rp=0;rp<2;rp++) // [RP,POI]
351  // Multiplicites:
352  TString typeMult[3] = {"RP","POI","REF"};
353  for(Int_t rpr=0;rpr<3;rpr++) // [RP,POI,REF]
354  {  
355   fMPC->SetnBinsMult(typeMult[rpr].Data(),fnBinsMult[rpr]);
356   fMPC->SetMinMult(typeMult[rpr].Data(),fMinMult[rpr]);
357   fMPC->SetMaxMult(typeMult[rpr].Data(),fMaxMult[rpr]);
358  } // for(Int_t rpr=0;rpr<3;rpr++) // [RP,POI,REF]
359
360  // Initialize:
361  fMPC->Init();
362  if(fMPC->GetHistList()) 
363  {
364   fHistList = fMPC->GetHistList();
365   // fHistList->Print();
366  } else 
367    {
368     Printf("ERROR: Could not retrieve histogram list (MPC, Task::UserCreateOutputObjects()) !!!!"); 
369    }
370  
371  PostData(1,fHistList);
372   
373 } // void AliAnalysisTaskMultiparticleCorrelations::UserCreateOutputObjects() 
374
375 //================================================================================================================
376
377 void AliAnalysisTaskMultiparticleCorrelations::UserExec(Option_t *) 
378 {
379  // Main loop (called for each event).
380
381  fEvent = dynamic_cast<AliFlowEventSimple*>(GetInputData(0));
382
383  // It's time for multi-particle correlations:
384  if(fEvent) 
385  {
386   fMPC->Make(fEvent);
387  } else 
388    {
389     cout<<" WARNING: No input data (MPC, Task::UserExec()) !!!!"<<endl;
390     cout<<endl;
391    }
392   
393  PostData(1,fHistList);
394
395 } // void AliAnalysisTaskMultiparticleCorrelations::UserExec(Option_t *) 
396
397 //================================================================================================================
398
399 void AliAnalysisTaskMultiparticleCorrelations::Terminate(Option_t *) 
400 {
401  // Accessing the merged output list. 
402
403  fHistList = (TList*)GetOutputData(1);
404  
405  fMPC = new AliFlowAnalysisWithMultiparticleCorrelations(); 
406  
407  if(fHistList) 
408  {
409   fMPC->GetOutputHistograms(fHistList);
410   fMPC->Finish();
411   PostData(1,fHistList);
412  } else
413    {
414     cout<<" WARNING: fHistList is NULL (MPC, Task::Terminate()) !!!!"<<endl;
415     cout<<endl;
416    }
417     
418 } // end of void AliAnalysisTaskMultiparticleCorrelations::Terminate(Option_t *)
419
420 //================================================================================================================
421
422 void AliAnalysisTaskMultiparticleCorrelations::SetWeightsHist(TH1D* const hist, const char *type, const char *variable)
423 {
424  // Pass histogram holding weights from an external file to the corresponding data member. 
425  
426  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::SetWeightsHist(TH1D* const hist, const char *type, const char *variable)";
427  
428  // Basic protection:
429  if(!hist){Fatal(sMethodName.Data(),"hist");}
430  if(!(TString(type).EqualTo("RP") || TString(type).EqualTo("POI"))){Fatal(sMethodName.Data(),"!(TString(type).EqualTo... type = %s ",type);}
431  if(!(TString(variable).EqualTo("phi") || TString(variable).EqualTo("pt") || TString(variable).EqualTo("eta"))){Fatal(sMethodName.Data(),"!(TString(variable).EqualTo... variable = %s ",variable);}
432
433  Int_t rp = 0; // [RP,POI]
434  if(TString(type).EqualTo("POI")){rp=1;} 
435
436  Int_t ppe = 0; // [phi,pt,eta]
437  if(TString(variable).EqualTo("pt")){ppe=1;} 
438  if(TString(variable).EqualTo("eta")){ppe=2;} 
439
440  // Finally:
441  hist->SetDirectory(0);
442  fWeightsHist[rp][ppe] = (TH1D*)hist->Clone();
443  if(!fWeightsHist[rp][ppe]){Fatal(sMethodName.Data(),"fWeightsHist[%d][%d]",rp,ppe);}
444
445  fUseWeights[rp][ppe] = kTRUE; 
446
447 } // void AliAnalysisTaskMultiparticleCorrelations::SetWeightsHist(TH1D* const hwh, const char *type, const char *variable)
448
449 //================================================================================================================
450
451 void AliAnalysisTaskMultiparticleCorrelations::SetnBins(const char *type, const char *variable, Int_t nBins)
452 {
453  // Set number of bins for histograms fKinematicsHist[2][3].
454
455  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::SetnBins(const char *type, const char *variable, Int_t nBins)";
456  
457  // Basic protection:
458  if(!(TString(type).EqualTo("RP") || TString(type).EqualTo("POI")))
459  {
460   cout<<"Well, it would be better for you to use RP or POI here..."<<endl;
461   Fatal(sMethodName.Data(),"!(TString(type).EqualTo... type = %s ",type);
462  }
463  if(!(TString(variable).EqualTo("phi") || TString(variable).EqualTo("pt") || TString(variable).EqualTo("eta")))
464  {
465   cout<<"phi, pt or eta, please!"<<endl;
466   Fatal(sMethodName.Data(),"!(TString(variable).EqualTo... variable = %s ",variable);
467  }
468
469  Int_t rp = 0; // [RP,POI]
470  if(TString(type).EqualTo("POI")){rp=1;} 
471
472  Int_t ppe = 0; // [phi,pt,eta]
473  if(TString(variable).EqualTo("pt")){ppe=1;} 
474  if(TString(variable).EqualTo("eta")){ppe=2;} 
475
476  fnBins[rp][ppe] = nBins;
477
478 } // void AliAnalysisTaskMultiparticleCorrelations::SetnBins(const char *type, const char *variable, Int_t nBins)
479
480 //=======================================================================================================================
481
482 void AliAnalysisTaskMultiparticleCorrelations::SetMin(const char *type, const char *variable, Double_t min)
483 {
484  // Set min bin range for histograms fKinematicsHist[2][3].
485
486  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::SetMin(const char *type, const char *variable, Double_t min)";
487  
488  // Basic protection:
489  if(!(TString(type).EqualTo("RP") || TString(type).EqualTo("POI")))
490  {
491   cout<<"Well, it would be better for you to use RP or POI here..."<<endl;
492   Fatal(sMethodName.Data(),"!(TString(type).EqualTo... type = %s ",type);
493  }
494  if(!(TString(variable).EqualTo("phi") || TString(variable).EqualTo("pt") || TString(variable).EqualTo("eta")))
495  {
496   cout<<"phi, pt or eta, please!"<<endl;
497   Fatal(sMethodName.Data(),"!(TString(variable).EqualTo... variable = %s ",variable);
498  }
499
500  Int_t rp = 0; // [RP,POI]
501  if(TString(type).EqualTo("POI")){rp=1;} 
502
503  Int_t ppe = 0; // [phi,pt,eta]
504  if(TString(variable).EqualTo("pt")){ppe=1;} 
505  if(TString(variable).EqualTo("eta")){ppe=2;} 
506
507  fMin[rp][ppe] = min;
508
509 } // void AliAnalysisTaskMultiparticleCorrelations::SetMin(const char *type, const char *variable, Double_t min)
510
511 //=======================================================================================================================
512
513 void AliAnalysisTaskMultiparticleCorrelations::SetMax(const char *type, const char *variable, Double_t max)
514 {
515  // Set max bin range for histograms fKinematicsHist[2][3].
516
517  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::SetMax(const char *type, const char *variable, Double_t max)";
518  
519  // Basic protection:
520  if(!(TString(type).EqualTo("RP") || TString(type).EqualTo("POI")))
521  {
522   cout<<"Well, it would be better for you to use RP or POI here..."<<endl;
523   Fatal(sMethodName.Data(),"!(TString(type).EqualTo... type = %s ",type);
524  }
525  if(!(TString(variable).EqualTo("phi") || TString(variable).EqualTo("pt") || TString(variable).EqualTo("eta")))
526  {
527   cout<<"phi, pt or eta, please!"<<endl;
528   Fatal(sMethodName.Data(),"!(TString(variable).EqualTo... variable = %s ",variable);
529  }
530
531  Int_t rp = 0; // [RP,POI]
532  if(TString(type).EqualTo("POI")){rp=1;} 
533
534  Int_t ppe = 0; // [phi,pt,eta]
535  if(TString(variable).EqualTo("pt")){ppe=1;} 
536  if(TString(variable).EqualTo("eta")){ppe=2;} 
537
538  fMax[rp][ppe] = max;
539
540 } // void AliAnalysisTaskMultiparticleCorrelations::SetMax(const char *type, const char *variable, Double_t min)
541
542 //=======================================================================================================================
543
544 void AliAnalysisTaskMultiparticleCorrelations::SetnBinsMult(const char *type, Int_t nBinsMult)
545 {
546  // Set number of bins for histograms fMultDistributionsHist[3].
547
548  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::SetnBinsMult(const char *type, Int_t nBinsMult)";
549  
550  // Basic protection:
551  if(!(TString(type).EqualTo("RP") || TString(type).EqualTo("POI") || TString(type).EqualTo("REF")))
552  {
553   cout<<"Well, it would be better for you to use RP, POI or REF here..."<<endl;
554   Fatal(sMethodName.Data(),"!(TString(type).EqualTo... type = %s ",type);
555  }
556
557  Int_t rpr = 0; // [RP,POI,REF]
558  if(TString(type).EqualTo("POI")){rpr=1;} 
559  else if(TString(type).EqualTo("REF")){rpr=2;} 
560
561  fnBinsMult[rpr] = nBinsMult;
562
563 } // void AliAnalysisTaskMultiparticleCorrelations::SetnBinsMult(const char *type, Int_t nBinsMult)
564
565 //=======================================================================================================================
566
567 void AliAnalysisTaskMultiparticleCorrelations::SetMinMult(const char *type, Double_t minMult)
568 {
569  // Set min bin range for histograms fMultDistributionsHist[3].
570
571  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::SetMinMult(const char *type, Double_t minMult)";
572  
573  // Basic protection:
574  if(!(TString(type).EqualTo("RP") || TString(type).EqualTo("POI") || TString(type).EqualTo("REF")))
575  {
576   cout<<"Well, it would be better for you to use RP, POI or REF here..."<<endl;
577   Fatal(sMethodName.Data(),"!(TString(type).EqualTo... type = %s ",type);
578  }
579
580  Int_t rpr = 0; // [RP,POI,REF]
581  if(TString(type).EqualTo("POI")){rpr=1;} 
582  else if(TString(type).EqualTo("REF")){rpr=2;} 
583
584  fMinMult[rpr] = minMult;
585
586 } // void AliAnalysisTaskMultiparticleCorrelations::SetMinMult(const char *type, Double_t minMult)
587
588 //=======================================================================================================================
589
590 void AliAnalysisTaskMultiparticleCorrelations::SetMaxMult(const char *type, Double_t maxMult)
591 {
592  // Set max bin range for histograms fMultDistributionsHist[3].
593
594  TString sMethodName = "void AliAnalysisTaskMultiparticleCorrelations::SetMaxMult(const char *type, Double_t maxMult)";
595  
596  // Basic protection:
597  if(!(TString(type).EqualTo("RP") || TString(type).EqualTo("POI") || TString(type).EqualTo("REF")))
598  {
599   cout<<"Well, it would be better for you to use RP, POI or REF here..."<<endl;
600   Fatal(sMethodName.Data(),"!(TString(type).EqualTo... type = %s ",type);
601  }
602
603  Int_t rpr = 0; // [RP,POI,REF]
604  if(TString(type).EqualTo("POI")){rpr=1;} 
605  else if(TString(type).EqualTo("REF")){rpr=2;} 
606
607  fMaxMult[rpr] = maxMult;
608
609 } // void AliAnalysisTaskMultiparticleCorrelations::SetMaxMult(const char *type, Double_t minMult)
610
611 //=======================================================================================================================
612
613
614
615
616