]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliDielectronCF.cxx
deleting this class
[u/mrichter/AliRoot.git] / PWG3 / dielectron / AliDielectronCF.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2009, 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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////
19 //       Dielectron Correction framework manager                         //
20 //                                                                       //
21 /*
22
23
24
25
26
27
28
29
30
31 */
32 //                                                                       //
33 ///////////////////////////////////////////////////////////////////////////
34
35 #include <TList.h>
36 #include <TObjArray.h>
37 #include <TVectorD.h>
38 #include <TString.h>
39 #include <TObjString.h>
40
41 #include <AliCFContainer.h>
42 #include <AliAnalysisFilter.h>
43 #include <AliAnalysisCuts.h>
44 #include <AliVParticle.h>
45 #include <AliLog.h>
46
47 #include "AliDielectronCF.h"
48 #include "AliDielectronMC.h"
49 #include "AliDielectronPair.h"
50
51 ClassImp(AliDielectronCF)
52
53 AliDielectronCF::AliDielectronCF() :
54   TNamed("DielectronCF","DielectronCF"),
55   fNSteps(0),
56   fNVars(0),
57   fVarBinLimits(0x0),
58   fNVarsLeg(0),
59   fVarBinLimitsLeg(0x0),
60   fNCuts(0),
61   fValues(0x0),
62   fStepForMCtruth(kFALSE),
63   fStepForNoCutsMCmotherPid(kFALSE),
64   fStepForAfterAllCuts(kTRUE),
65   fStepForPreFilter(kFALSE),
66   fStepsForEachCut(kFALSE),
67   fStepsForCutsIncreasing(kFALSE),
68   fStepsForSignal(kTRUE),
69   fStepsForBackground(kFALSE),
70   fNStepMasks(0),
71   fPdgMother(-1),
72   fCfContainer(0x0),
73   fHasMC(kFALSE),
74   fNAddSteps(0)
75 {
76   //
77   // Default constructor
78   //
79   for (Int_t i=0; i<AliDielectronVarManager::kNMaxValues; ++i){
80     fVariables[i]=0;
81     fVariablesLeg[i]=0;
82   }
83
84   for (Int_t i=0; i<kNmaxAddSteps; ++i){
85     fStepMasks[i]=0xFFFFFF;
86   }
87 }
88
89 //________________________________________________________________
90 AliDielectronCF::AliDielectronCF(const char* name, const char* title) :
91   TNamed(name, title),
92   fNSteps(0),
93   fNVars(0),
94   fVarBinLimits(0x0),
95   fNVarsLeg(0),
96   fVarBinLimitsLeg(0x0),
97   fNCuts(0),
98   fValues(0x0),
99   fStepForMCtruth(kFALSE),
100   fStepForNoCutsMCmotherPid(kFALSE),
101   fStepForAfterAllCuts(kTRUE),
102   fStepForPreFilter(kFALSE),
103   fStepsForEachCut(kFALSE),
104   fStepsForCutsIncreasing(kFALSE),
105   fStepsForSignal(kTRUE),
106   fStepsForBackground(kFALSE),
107   fNStepMasks(0),
108   fPdgMother(-1),
109   fCfContainer(0x0),
110   fHasMC(kFALSE),
111   fNAddSteps(0)
112 {
113   //
114   // Named constructor
115   //
116   for (Int_t i=0; i<AliDielectronVarManager::kNMaxValues; ++i){
117     fVariables[i]=0;
118     fVariablesLeg[i]=0;
119   }
120   
121   for (Int_t i=0; i<kNmaxAddSteps; ++i){
122     fStepMasks[i]=0xFFFFFF;
123   }
124 }
125
126 //________________________________________________________________
127 AliDielectronCF::~AliDielectronCF()
128 {
129   //
130   // Destructor
131   //
132   if (fValues) delete [] fValues;
133   if (fVarBinLimits) delete fVarBinLimits;
134   if (fVarBinLimitsLeg) delete fVarBinLimitsLeg;
135 }
136
137 //________________________________________________________________
138 void AliDielectronCF::AddVariable(AliDielectronVarManager::ValueTypes type, Int_t nbins,
139                                   Double_t min, Double_t max, Bool_t leg, Bool_t log)
140 {
141   //
142   // Add a variable to the CF configuration
143   // if leg is true it will add the variables of the leg
144   // if log is true log binning will be created
145   //
146
147   TVectorD *binLimits=0x0;
148   if (!log) binLimits=MakeLinBinning(nbins,min,max);
149   else binLimits=MakeLogBinning(nbins,min,max);
150   AddVariable(type,binLimits,leg);
151 }
152
153 //________________________________________________________________
154 void AliDielectronCF::AddVariable(AliDielectronVarManager::ValueTypes type, const char* binLimitStr, Bool_t leg/*=kFALSE*/)
155 {
156   //
157   // Add a variable to the CF configuration
158   // specify arbitrary binning in a string.
159   // Bin limits need to be separated by a ","
160   //
161   TString limits(binLimitStr);
162   if (limits.IsNull()){
163     AliError(Form("Bin Limit string is empty, cannot add the variable '%s'",AliDielectronVarManager::GetValueName(type)));
164     return;
165   }
166
167   TObjArray *arr=limits.Tokenize(",");
168   Int_t nLimits=arr->GetEntries();
169   if (nLimits<2){
170     AliError(Form("Need at leas 2 bin limits, cannot add the variable '%s'",AliDielectronVarManager::GetValueName(type)));
171     delete arr;
172     return;
173   }
174
175   TVectorD *binLimits=new TVectorD(nLimits);
176   for (Int_t iLim=0; iLim<nLimits; ++iLim){
177     (*binLimits)[iLim]=(static_cast<TObjString*>(arr->At(iLim)))->GetString().Atof();
178   }
179
180   delete arr;
181   AddVariable(type,binLimits,leg);
182 }
183
184 //________________________________________________________________
185 void AliDielectronCF::AddVariable(AliDielectronVarManager::ValueTypes type, TVectorD *binLimits, Bool_t leg/*=kFALSE*/)
186 {
187   //
188   // Add variable with the binning given in the TVectorD
189   //
190   if (!leg){
191     if (!fVarBinLimits){
192       fVarBinLimits=new TObjArray;
193       fVarBinLimits->SetOwner();
194     }
195     fVarBinLimits->Add(binLimits);
196     fVariables[fNVars]  = (UInt_t)type;
197     ++fNVars;
198   } else {
199     if (!fVarBinLimitsLeg){
200       fVarBinLimitsLeg=new TObjArray;
201       fVarBinLimitsLeg->SetOwner();
202     }
203     fVarBinLimitsLeg->Add(binLimits);
204     fVariablesLeg[fNVarsLeg]  = (UInt_t)type;
205     ++fNVarsLeg;
206   }
207 }
208
209 //________________________________________________________________
210 void AliDielectronCF::InitialiseContainer(const AliAnalysisFilter& filter)
211 {
212   //
213   // Initialise container based on the cuts in the analysis filter
214   //
215
216   fNCuts=filter.GetCuts()->GetEntries();
217
218   fHasMC=AliDielectronMC::Instance()->HasMC();
219   fNAddSteps=1;
220   if (fHasMC){
221     if (fStepsForSignal) ++fNAddSteps;
222     if (fStepsForBackground) ++fNAddSteps;
223   } else {
224     //if 
225     fStepForMCtruth=kFALSE;
226     fStepForNoCutsMCmotherPid=kFALSE;
227     fStepsForSignal=kFALSE;
228     fStepsForBackground=kFALSE;
229   }
230     
231   fNSteps=0;
232   if (fStepForMCtruth) ++fNSteps;
233   if (fStepForNoCutsMCmotherPid) ++fNSteps;
234   if (fStepForAfterAllCuts) fNSteps+=fNAddSteps;
235
236   if (fStepsForEachCut&&fNCuts>1)        fNSteps+=(fNAddSteps*fNCuts);     //one step for each cut + Signal (MC)
237   if (fStepsForCutsIncreasing&&fNCuts>2) fNSteps+=(fNAddSteps*(fNCuts-2)); //one step for the increasing cuts + Signal (MC)
238                                                       // e.g. cut2&cut3, cut2&cut3&cut4
239   fNSteps+=(fNAddSteps*fNStepMasks);                            // cuts for the additional cut masks
240
241   if (fStepForPreFilter) fNSteps+=fNAddSteps; //Add at the end for Prefilter (maxcutmask+1)
242
243   // create the container
244   
245   Int_t *nbins=new Int_t[fNVars+2*fNVarsLeg];
246   for (Int_t i=0;i<fNVars;++i) {
247     Int_t nBins=(static_cast<TVectorD*>(fVarBinLimits->At(i)))->GetNrows()-1;
248     nbins[i]=nBins;
249   }
250   for (Int_t i=0;i<fNVarsLeg;++i){
251     Int_t nBins=(static_cast<TVectorD*>(fVarBinLimitsLeg->At(i)))->GetNrows()-1;
252     nbins[i+fNVars]=nBins;
253     nbins[i+fNVars+fNVarsLeg]=nBins;
254   }
255   
256   fCfContainer = new AliCFContainer(GetName(), GetTitle(), fNSteps, fNVars+2*fNVarsLeg, nbins);
257   delete [] nbins;
258   
259   // initialize the variables and their bin limits
260   for (Int_t iVar=0; iVar<fNVars; iVar++) {
261     UInt_t type=fVariables[iVar];
262     Double_t *binLim = (static_cast<TVectorD*>(fVarBinLimits->At(iVar)))->GetMatrixArray();
263
264     fCfContainer->SetBinLimits(iVar, binLim);
265     fCfContainer->SetVarTitle(iVar, AliDielectronVarManager::GetValueName(type));
266   }
267   
268   // initialize the variables and their bin limits for the Legs
269   for (Int_t iVar=0; iVar<fNVarsLeg; iVar++) {
270     UInt_t type=fVariablesLeg[iVar];
271     Double_t *binLim=(static_cast<TVectorD*>(fVarBinLimitsLeg->At(iVar)))->GetMatrixArray();
272
273     //Leg1
274     fCfContainer->SetBinLimits(iVar+fNVars, binLim);
275     fCfContainer->SetVarTitle(iVar+fNVars, Form("Leg1_%s",AliDielectronVarManager::GetValueName(type)));
276     
277     //Leg2
278     fCfContainer->SetBinLimits(iVar+fNVars+fNVarsLeg, binLim);
279     fCfContainer->SetVarTitle(iVar+fNVars+fNVarsLeg, Form("Leg2_%s",AliDielectronVarManager::GetValueName(type)));
280   }
281
282   // array for storing values
283   fValues = new Double_t[fNVars+2*fNVarsLeg];
284   
285   //=================//
286   // Set step titles //
287   //=================//
288   Int_t step=0;
289
290   //Pure MC truth
291   if (fStepForMCtruth){
292     fCfContainer->SetStepTitle(step++,"MC truth");
293   }
294
295   //before cuts (MC truth)
296   if (fStepForNoCutsMCmotherPid){
297     fCfContainer->SetStepTitle(step++,"No cuts (Signal)");
298   }
299   
300   TString cutName;
301   //Steps for each of the cuts
302   if (fStepsForEachCut&&fNCuts>1){
303     for (Int_t iCut=0; iCut<fNCuts;++iCut) {
304       cutName=filter.GetCuts()->At(iCut)->GetName(); //TODO: User GetTitle???
305       
306       fCfContainer->SetStepTitle(step++, cutName.Data()); //Step for the cut
307       
308       if (fHasMC){
309         if (fStepsForSignal)
310           fCfContainer->SetStepTitle(step++, (cutName+" (Signal)").Data()); //Step for the cut with MC truth
311         if (fStepsForBackground)
312           fCfContainer->SetStepTitle(step++, (cutName+" (Background)").Data()); //Step for the cut with MC truth
313       }
314     }
315   }
316
317   //Steps for increasing cut match
318   if (fStepsForCutsIncreasing&&fNCuts>2){
319     cutName=filter.GetCuts()->At(0)->GetName(); //TODO: User GetTitle???
320     for (Int_t iCut=1; iCut<fNCuts-1;++iCut) {
321       cutName+="&";
322       cutName+=filter.GetCuts()->At(iCut)->GetName();
323       
324       fCfContainer->SetStepTitle(step++, cutName.Data()); //Step for the cut
325       
326       if (fHasMC){
327         if (fStepsForSignal)
328           fCfContainer->SetStepTitle(step++, (cutName+" (Signal)").Data()); //Step for the cut with MC truth
329         if (fStepsForBackground)
330           fCfContainer->SetStepTitle(step++, (cutName+" (Background)").Data()); //Step for the cut with MC truth
331       }
332     }
333   }
334
335   //Steps of user defined cut combinations
336   for (UInt_t iComb=0; iComb<fNStepMasks; ++iComb){
337     cutName="";
338     UInt_t mask=fStepMasks[iComb];
339     for (Int_t iCut=0; iCut<fNCuts;++iCut) {
340       if (mask&(1<<iCut)){
341         if (cutName.IsNull()){
342           cutName=filter.GetCuts()->At(iCut)->GetName();
343         }else{
344           cutName+="&";
345           cutName+=filter.GetCuts()->At(iCut)->GetName();
346         }
347       }
348     }
349     
350     fCfContainer->SetStepTitle(step++, cutName.Data()); //Step for the cut
351     
352     if (fHasMC){
353       if (fStepsForSignal)
354         fCfContainer->SetStepTitle(step++, (cutName+" (Signal)").Data()); //Step for the cut with MC truth
355       if (fStepsForBackground)
356         fCfContainer->SetStepTitle(step++, (cutName+" (Background)").Data()); //Step for the cut with MC truth
357     }
358   }
359
360   //After All cuts
361   if (fStepForAfterAllCuts){
362     cutName="No pair cuts";
363     if (filter.GetCuts()->At(0)){
364       cutName=filter.GetCuts()->At(0)->GetName(); //TODO: User GetTitle???
365       for (Int_t iCut=1; iCut<fNCuts;++iCut) {
366         cutName+="&";
367         cutName+=filter.GetCuts()->At(iCut)->GetName();
368       }
369     }
370     fCfContainer->SetStepTitle(step++, cutName.Data()); //Step for the cut
371     if (fHasMC){
372       if (fStepsForSignal)
373         fCfContainer->SetStepTitle(step++, (cutName+" (Signal)").Data()); //Step for the cut with MC truth
374       if (fStepsForBackground)
375         fCfContainer->SetStepTitle(step++, (cutName+" (Background)").Data()); //Step for the cut with MC truth
376     }
377   }
378
379   //Additional Step for result after PreFilter
380   if (fStepForPreFilter){
381     cutName="PreFilter";
382     fCfContainer->SetStepTitle(step++, cutName.Data()); //Step for the cut
383     if (fHasMC){
384       if (fStepsForSignal)
385         fCfContainer->SetStepTitle(step++, (cutName+" (Signal)").Data()); //Step for the cut with MC truth
386       if (fStepsForBackground)
387         fCfContainer->SetStepTitle(step++, (cutName+" (Background)").Data()); //Step for the cut with MC truth
388     }
389   }
390
391
392
393   if (step!=fNSteps) {
394     AliError(Form("Something went wrong in the naming of the steps!!! (%d != %d)",step,fNSteps));
395   }
396 }
397
398 //________________________________________________________________
399 void AliDielectronCF::Fill(UInt_t mask, const AliDielectronPair *particle)
400 {
401   //
402   // Fill the containers
403   //
404
405   Bool_t isMCTruth=kFALSE;
406   if (fHasMC) isMCTruth=AliDielectronMC::Instance()->IsMotherPdg(particle,fPdgMother);
407
408   Double_t valuesPair[AliDielectronVarManager::kNMaxValues];
409   AliDielectronVarManager::Fill(particle,valuesPair);
410
411   for (Int_t iVar=0; iVar<fNVars; ++iVar){
412     Int_t var=fVariables[iVar];
413     fValues[iVar]=valuesPair[var];
414   }
415
416   if (fNVarsLeg>0){
417     Double_t valuesLeg1[AliDielectronVarManager::kNMaxValues];
418     AliDielectronVarManager::Fill(particle->GetFirstDaughter(),valuesLeg1);
419     Double_t valuesLeg2[AliDielectronVarManager::kNMaxValues];
420     AliDielectronVarManager::Fill(particle->GetSecondDaughter(),valuesLeg2);
421
422     for (Int_t iVar=0; iVar<fNVarsLeg; ++iVar){
423       Int_t var=fVariablesLeg[iVar];
424       fValues[iVar+fNVars]=valuesLeg1[var];
425       fValues[iVar+fNVars+fNVarsLeg]=valuesLeg2[var];
426     }
427   }
428   
429   UInt_t selectedMask=(1<<fNCuts)-1;
430
431   //============//
432   // Fill steps //
433   //============//
434   // step 0 would be full MC truth and is handled in FillMC
435   Int_t step=0;
436   if (fStepForMCtruth) ++step;
437   
438   //No cuts (MC truth)
439   if (fStepForNoCutsMCmotherPid){
440     if (isMCTruth) fCfContainer->Fill(fValues,step);
441     ++step;
442   }
443   
444   //Steps for each of the cuts
445   if (fStepsForEachCut&&fNCuts>1){
446     for (Int_t iCut=0; iCut<fNCuts;++iCut) {
447       if (mask&(1<<iCut)) {
448         fCfContainer->Fill(fValues,step);
449         ++step;
450
451         if (fHasMC){
452           if ( fStepsForSignal){
453             if (isMCTruth) fCfContainer->Fill(fValues,step);
454             ++step;
455           }
456           if ( fStepsForBackground ){
457             if (!isMCTruth) fCfContainer->Fill(fValues,step);
458             ++step;
459           }
460         }
461       } else {
462         step+=fNAddSteps;
463       }
464     }
465   }
466
467
468   //Steps for increasing cut match
469   if (fStepsForCutsIncreasing&&fNCuts>2){
470     for (Int_t iCut=1; iCut<fNCuts-1;++iCut) {
471       if (mask&(1<<((iCut+1)-1))) {
472         fCfContainer->Fill(fValues,step);
473         ++step;
474         
475         if (fHasMC){
476           if ( fStepsForSignal){
477             if (isMCTruth) fCfContainer->Fill(fValues,step);
478             ++step;
479           }
480           if ( fStepsForBackground ){
481             if (!isMCTruth) fCfContainer->Fill(fValues,step);
482             ++step;
483           }
484         }
485       } else {
486         step+=fNAddSteps;
487       }
488     }
489   }
490
491   //Steps of user defined cut combinations
492   for (UInt_t iComb=0; iComb<fNStepMasks; ++iComb){
493     UInt_t userMask=fStepMasks[iComb];
494     if (mask&userMask) {
495       fCfContainer->Fill(fValues,step);
496       ++step;
497       
498       if (fHasMC){
499         if ( fStepsForSignal){
500           if (isMCTruth) fCfContainer->Fill(fValues,step);
501           ++step;
502         }
503         if ( fStepsForBackground ){
504           if (!isMCTruth) fCfContainer->Fill(fValues,step);
505           ++step;
506         }
507       }
508     } else {
509       step+=fNAddSteps;
510     }
511   }
512   
513   //All cuts
514   if (fStepForAfterAllCuts){
515     if (mask == selectedMask){
516       fCfContainer->Fill(fValues,step);
517       ++step;
518       
519       if (fHasMC){
520         if ( fStepsForSignal){
521           if (isMCTruth) fCfContainer->Fill(fValues,step);
522           ++step;
523         }
524         if ( fStepsForBackground ){
525           if (!isMCTruth) fCfContainer->Fill(fValues,step);
526           ++step;
527         }
528       }
529     } else {
530       step+=fNAddSteps;
531     }
532   }
533   if (fStepForPreFilter) {
534         if (mask&(1<<fNCuts)) {
535           fCfContainer->Fill(fValues,step);
536           ++step;
537
538           if (fHasMC){
539                 if ( fStepsForSignal){
540                   if (isMCTruth) fCfContainer->Fill(fValues,step);
541                   ++step;
542                 }
543                 if ( fStepsForBackground ){
544                   if (!isMCTruth) fCfContainer->Fill(fValues,step);
545                   ++step;
546                 }
547           }
548         }
549         else {
550           step+=fNAddSteps;
551         }
552   } 
553   if (step!=fNSteps) {
554     AliError("Something went wrong in the step filling!!!");
555   }
556   
557 }
558
559 //________________________________________________________________
560 void AliDielectronCF::FillMC(const TObject *particle)
561 {
562   //
563   // fill MC part of the Container
564   //
565   if (!fStepForMCtruth) return;
566   
567   Double_t valuesPair[AliDielectronVarManager::kNMaxValues];
568   AliDielectronVarManager::Fill(particle,valuesPair);
569
570   AliVParticle *d1=0x0;
571   AliVParticle *d2=0x0;
572   AliDielectronMC::Instance()->GetDaughters(particle,d1,d2);
573   
574   valuesPair[AliDielectronVarManager::kThetaHE]=AliDielectronPair::ThetaPhiCM(d1,d2,kTRUE,kTRUE);
575   valuesPair[AliDielectronVarManager::kPhiHE]=AliDielectronPair::ThetaPhiCM(d1,d2,kTRUE,kFALSE);
576   valuesPair[AliDielectronVarManager::kThetaCS]=AliDielectronPair::ThetaPhiCM(d1,d2,kFALSE,kTRUE);
577   valuesPair[AliDielectronVarManager::kPhiCS]=AliDielectronPair::ThetaPhiCM(d1,d2,kFALSE,kFALSE);
578   
579   //TODO: temporary solution, set manually the pair type to 1: unlikesign SE
580   valuesPair[AliDielectronVarManager::kPairType]=1;
581   
582   for (Int_t iVar=0; iVar<fNVars; ++iVar){
583     Int_t var=fVariables[iVar];
584     fValues[iVar]=valuesPair[var];
585   }
586   
587   if (fNVarsLeg>0){
588     Double_t valuesLeg1[AliDielectronVarManager::kNMaxValues];
589     Double_t valuesLeg2[AliDielectronVarManager::kNMaxValues];
590     if (d1->Pt()>d2->Pt()){
591       AliDielectronVarManager::Fill(d1,valuesLeg1);
592       AliDielectronVarManager::Fill(d2,valuesLeg2);
593     } else {
594       AliDielectronVarManager::Fill(d2,valuesLeg1);
595       AliDielectronVarManager::Fill(d1,valuesLeg2);
596     }
597     
598     for (Int_t iVar=0; iVar<fNVarsLeg; ++iVar){
599       Int_t var=fVariablesLeg[iVar];
600       fValues[iVar+fNVars]=valuesLeg1[var];
601       fValues[iVar+fNVars+fNVarsLeg]=valuesLeg2[var];
602     }
603   }
604   
605   fCfContainer->Fill(fValues,0);
606 }
607
608 //_____________________________________________________________________________
609 TVectorD* AliDielectronCF::MakeLogBinning(Int_t nbinsX, Double_t xmin, Double_t xmax) const
610 {
611   //
612   // Make logarithmic binning
613   // the user has to delete the array afterwards!!!
614   //
615   
616   //check limits
617   if (xmin<1e-20 || xmax<1e-20){
618     AliError("For Log binning xmin and xmax must be > 1e-20. Using linear binning instead!");
619     return MakeLinBinning(nbinsX, xmin, xmax);
620   }
621   if (xmax<xmin){
622     Double_t tmp=xmin;
623     xmin=xmax;
624     xmax=tmp;
625   }
626   TVectorD *binLim=new TVectorD(nbinsX+1);
627   Double_t first=xmin;
628   Double_t last=xmax;
629   Double_t expMax=TMath::Log(last/first);
630   for (Int_t i=0; i<nbinsX+1; ++i){
631     (*binLim)[i]=first*TMath::Exp(expMax/nbinsX*(Double_t)i);
632   }
633   return binLim;
634 }
635
636 //_____________________________________________________________________________
637 TVectorD* AliDielectronCF::MakeLinBinning(Int_t nbinsX, Double_t xmin, Double_t xmax) const
638 {
639   //
640   // Make logarithmic binning
641   // the user has to delete the array afterwards!!!
642   //
643   if (xmax<xmin){
644     Double_t tmp=xmin;
645     xmin=xmax;
646     xmax=tmp;
647   }
648   TVectorD *binLim=new TVectorD(nbinsX+1);
649   Double_t first=xmin;
650   Double_t last=xmax;
651   Double_t binWidth=(last-first)/nbinsX;
652   for (Int_t i=0; i<nbinsX+1; ++i){
653     (*binLim)[i]=first+binWidth*(Double_t)i;
654   }
655   return binLim;
656 }
657