]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliDielectronHF.cxx
merging trunk to TPCdev
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliDielectronHF.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 ///////////////////////////////////////////////////////////////////////////
17 //                Dielectron HF                                  //
18 //                                                                       //
19 //                                                                       //
20 /*
21 Detailed description
22
23
24 */
25 //                                                                       //
26 ///////////////////////////////////////////////////////////////////////////
27
28 #include <TVectorD.h>
29 #include <TH1.h>
30 #include <TH1F.h>
31 #include <TH2.h>
32 #include <TH3.h>
33 #include <TProfile.h>
34 #include <TProfile2D.h>
35 #include <TProfile3D.h>
36 #include <TAxis.h>
37 #include <TString.h>
38 #include <TObjString.h>
39 #include <TObjArray.h>
40
41 #include <AliVParticle.h>
42 #include <AliLog.h>
43
44 #include "AliDielectron.h"
45 #include "AliDielectronHelper.h"
46 #include "AliDielectronMC.h"
47 #include "AliDielectronPair.h"
48 #include "AliDielectronSignalMC.h"
49
50 #include "AliDielectronHistos.h"
51 #include "AliDielectronHF.h"
52
53 ClassImp(AliDielectronHF)
54
55 AliDielectronHF::AliDielectronHF() :
56   TNamed(),
57   fArrPairType(0x0),
58   fPairType(kOSonly),
59   fSignalsMC(0x0),
60   fAxes(kMaxCuts),
61   fHasMC(kFALSE),
62   fStepGenerated(kFALSE),
63   fRefObj(1)
64 {
65   //
66   // Default Constructor
67   //
68   for (Int_t i=0; i<kMaxCuts; ++i){
69     fVarCuts[i]=0;
70     fVarCutType[i]=0;
71     fBinType[i]=kStdBin;
72   }
73   fAxes.SetOwner(kTRUE);
74   fRefObj.SetOwner(kTRUE);
75 }
76
77 //______________________________________________
78 AliDielectronHF::AliDielectronHF(const char* name, const char* title) :
79   TNamed(name, title),
80   fArrPairType(0x0),
81   fPairType(kOSonly),
82   fSignalsMC(0x0),
83   fAxes(kMaxCuts),
84   fHasMC(kFALSE),
85   fStepGenerated(kFALSE),
86   fRefObj(1)
87 {
88   //
89   // Named Constructor
90   //
91   for (Int_t i=0; i<kMaxCuts; ++i){
92     fVarCuts[i]=0;
93     fVarCutType[i]=0;
94     fBinType[i]=kStdBin;
95   }
96   fAxes.SetOwner(kTRUE);
97   fRefObj.SetOwner(kTRUE);
98 }
99
100 //______________________________________________
101 AliDielectronHF::~AliDielectronHF()
102 {
103   //
104   // Default Destructor
105   //
106   fAxes.Delete();
107   fRefObj.Delete();
108   fArrPairType.Delete();
109 }
110
111 //_____________________________________________________________________________
112 void AliDielectronHF::UserProfile(const char* histClass, UInt_t valTypeP,
113                                       const TVectorD * const binsX,
114                                       UInt_t valTypeX, TString option)
115 {
116   //
117   // Histogram creation 1D case with arbitraty binning X
118   // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
119   //
120
121   TH1 *hist=0x0;
122   if(valTypeP==999)
123     hist=new TH1F("","",binsX->GetNrows()-1,binsX->GetMatrixArray());
124   else {
125     TString opt=""; Double_t pmin=0., pmax=0.;
126     if(!option.IsNull()) {
127       TObjArray *arr=option.Tokenize(";");
128       arr->SetOwner();
129       opt=((TObjString*)arr->At(0))->GetString();
130       if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
131       if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
132       delete arr;
133     }
134     hist=new TProfile("","",binsX->GetNrows()-1,binsX->GetMatrixArray());
135     ((TProfile*)hist)->BuildOptions(pmin,pmax,opt.Data());
136     //      printf(" name %s PROFILE options: pmin %.1f pmax %.1f err %s \n",name,((TProfile*)hist)->GetYmin(),((TProfile*)hist)->GetYmax(),((TProfile*)hist)->GetErrorOption() );
137   }
138
139   // store variales in axes
140   UInt_t valType[4] = {0};
141   valType[0]=valTypeX;     valType[1]=valTypeP;
142   AliDielectronHistos::StoreVariables(hist, valType);
143
144   // adapt the name and title of the histogram in case they are empty
145   AliDielectronHistos::AdaptNameTitle(hist, histClass);
146   hist->SetName(Form("HF_%s",hist->GetName()));
147
148   fRefObj.AddLast(hist);
149   delete binsX;
150 }
151
152 //_____________________________________________________________________________
153 void AliDielectronHF::UserProfile(const char* histClass, UInt_t valTypeP,
154                                       const TVectorD * const binsX, const TVectorD * const binsY,
155                                       UInt_t valTypeX, UInt_t valTypeY, TString option)
156 {
157   //
158   // Histogram creation 2D case with arbitraty binning X and Y
159   // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
160   //
161
162   TH1 *hist=0x0;
163   if(valTypeP==999) {
164     hist=new TH2F("","",
165                   binsX->GetNrows()-1,binsX->GetMatrixArray(),
166                   binsY->GetNrows()-1,binsY->GetMatrixArray()); 
167   }
168   else  {
169     TString opt=""; Double_t pmin=0., pmax=0.;
170     if(!option.IsNull()) {
171       TObjArray *arr=option.Tokenize(";");
172       arr->SetOwner();
173       opt=((TObjString*)arr->At(0))->GetString();
174       if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
175       if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
176       delete arr;
177     }
178     hist=new TProfile2D("","",
179                         binsX->GetNrows()-1,binsX->GetMatrixArray(),
180                         binsY->GetNrows()-1,binsY->GetMatrixArray());
181     ((TProfile2D*)hist)->BuildOptions(pmin,pmax,opt.Data());
182   }
183
184   // store variales in axes
185   UInt_t valType[4] = {0};
186   valType[0]=valTypeX;     valType[1]=valTypeY; valType[3]=valTypeP;
187   AliDielectronHistos::StoreVariables(hist, valType);
188
189   // adapt the name and title of the histogram in case they are empty
190   AliDielectronHistos::AdaptNameTitle(hist, histClass);
191   hist->SetName(Form("HF_%s",hist->GetName()));
192
193   fRefObj.AddLast(hist);
194   delete binsX;
195   delete binsY;
196 }
197
198 //_____________________________________________________________________________
199 void AliDielectronHF::UserProfile(const char* histClass, UInt_t valTypeP,
200                                       const TVectorD * const binsX, const TVectorD * const binsY, const TVectorD * const binsZ,
201                                       UInt_t valTypeX, UInt_t valTypeY, UInt_t valTypeZ, TString option)
202 {
203   //
204   // Histogram creation 3D case with arbitraty binning X, Y, Z
205   // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
206   //
207   TH1 *hist=0x0;
208   if(valTypeP==999) {
209     hist=new TH3F("","",
210                   binsX->GetNrows()-1,binsX->GetMatrixArray(),
211                   binsY->GetNrows()-1,binsY->GetMatrixArray(),
212                   binsZ->GetNrows()-1,binsZ->GetMatrixArray());
213   }
214   else {
215     TString opt=""; Double_t pmin=0., pmax=0.;
216     if(!option.IsNull()) {
217       TObjArray *arr=option.Tokenize(";");
218       arr->SetOwner();
219       opt=((TObjString*)arr->At(0))->GetString();
220       if(arr->GetEntriesFast()>1) pmin=(((TObjString*)arr->At(1))->GetString()).Atof();
221       if(arr->GetEntriesFast()>2) pmax=(((TObjString*)arr->At(2))->GetString()).Atof();
222       delete arr;
223     }
224     hist=new TProfile3D("","",
225                         binsX->GetNrows()-1,binsX->GetMatrixArray(),
226                         binsY->GetNrows()-1,binsY->GetMatrixArray(),
227                         binsZ->GetNrows()-1,binsZ->GetMatrixArray());
228     ((TProfile3D*)hist)->BuildOptions(pmin,pmax,opt.Data());
229   }
230
231   // store variales in axes
232   UInt_t valType[4] = {0};
233   valType[0]=valTypeX;     valType[1]=valTypeY;     valType[2]=valTypeZ;     valType[3]=valTypeP;
234   AliDielectronHistos::StoreVariables(hist, valType);
235
236   // adapt the name and title of the histogram in case they are empty
237   AliDielectronHistos::AdaptNameTitle(hist, histClass);
238   hist->SetName(Form("HF_%s",hist->GetName()));
239
240   fRefObj.AddLast(hist);
241   delete binsX;
242   delete binsY;
243   delete binsZ;
244 }
245
246 //________________________________________________________________
247 void AliDielectronHF::AddCutVariable(AliDielectronVarManager::ValueTypes type,
248                                      Int_t nbins, Double_t min, Double_t max, Bool_t log, Bool_t leg, EBinType btype)
249 {
250   //
251   // Add a variable to the mixing handler
252   //
253
254   // limit number of variables to kMaxCuts
255   if (fAxes.GetEntriesFast()>=kMaxCuts) return;
256   
257   TVectorD *binLimits=0x0;
258   if (!log) binLimits=AliDielectronHelper::MakeLinBinning(nbins,min,max);
259   else binLimits=AliDielectronHelper::MakeLogBinning(nbins,min,max);
260   if (!binLimits) return;
261
262   Int_t size=fAxes.GetEntriesFast();
263   fVarCuts[size]=(UShort_t)type;
264   fVarCutType[size]=leg;
265   fAxes.Add(binLimits->Clone());
266   fBinType[size]=btype;
267 }
268
269 //________________________________________________________________
270 void AliDielectronHF::AddCutVariable(AliDielectronVarManager::ValueTypes type,
271                                              const char* binLimitStr, Bool_t leg, EBinType btype)
272 {
273   //
274   // Add a variable to the mixing handler with arbitrary binning
275   //
276
277   // limit number of variables to kMaxCuts
278   if (fAxes.GetEntriesFast()>=kMaxCuts) return;
279   
280   TVectorD *binLimits=AliDielectronHelper::MakeArbitraryBinning(binLimitStr);
281   if (!binLimits) return;
282   
283   Int_t size=fAxes.GetEntriesFast();
284   fVarCuts[size]=(UShort_t)type;
285   fVarCutType[size]=leg;
286   fAxes.Add(binLimits);
287   fBinType[size]=btype;
288 }
289
290 //________________________________________________________________
291 void AliDielectronHF::AddCutVariable(AliDielectronVarManager::ValueTypes type,
292                                              TVectorD * binLimits, Bool_t leg, EBinType btype)
293 {
294   //
295   // Add a variable to the mixing handler with a vector
296   // the TVectorD is assumed to be surplus after the creation and will be deleted!!!
297   //
298
299   // limit number of variables to kMaxCuts
300   if (fAxes.GetEntriesFast()>=kMaxCuts) return;
301   
302   if (!binLimits) return;
303   
304   Int_t size=fAxes.GetEntriesFast();
305   fVarCuts[size]=(UShort_t)type;
306   fVarCutType[size]=leg;
307   fAxes.Add(binLimits);
308   fBinType[size]=btype;
309 }
310
311 //______________________________________________
312 void AliDielectronHF::Fill(Int_t label1, Int_t label2, Int_t nSignal) 
313 {
314   //
315   // fill the pure MC part of the container starting from a pair of 2 particles (part1 and part2 are legs)
316   //
317   // fill only if we have asked for these steps
318   if(!fStepGenerated) return;
319
320   AliVParticle* part1 = AliDielectronMC::Instance()->GetMCTrackFromMCEvent(label1);
321   AliVParticle* part2 = AliDielectronMC::Instance()->GetMCTrackFromMCEvent(label2);
322   if(!part1 || !part2) return;
323
324   AliDielectronMC* dieMC = AliDielectronMC::Instance();
325   
326   Int_t mLabel1 = dieMC->GetMothersLabel(label1);    // should work for both ESD and AOD
327   Int_t mLabel2 = dieMC->GetMothersLabel(label2);
328
329   // check the same mother option
330   AliDielectronSignalMC* sigMC = (AliDielectronSignalMC*)fSignalsMC->At(nSignal);
331   if(sigMC->GetMothersRelation()==AliDielectronSignalMC::kSame && mLabel1!=mLabel2) return;
332   if(sigMC->GetMothersRelation()==AliDielectronSignalMC::kDifferent && mLabel1==mLabel2) return;
333     
334   // fill the leg variables
335   Double_t valuesLeg1[AliDielectronVarManager::kNMaxValues];
336   Double_t valuesLeg2[AliDielectronVarManager::kNMaxValues];
337   AliDielectronVarManager::Fill(part1,valuesLeg1);
338   AliDielectronVarManager::Fill(part2,valuesLeg2);
339     
340   // fill the pair and event variables
341   Double_t valuesPair[AliDielectronVarManager::kNMaxValues];
342   AliDielectronVarManager::Fill(dieMC->GetMCEvent(), valuesPair);
343   AliDielectronVarManager::FillVarMCParticle2(part1,part2,valuesPair);
344
345   if(part1->Charge()*part2->Charge()<0) {
346     //valuesPair[AliDielectronVarManager::kPairType]=1;
347 //     for(Int_t i=0; i<fAxes.GetEntriesFast(); i++) {
348 //       printf("[D]: %s %f %f %f \n",
349 //           AliDielectronVarManager::GetValueName(fVarCuts[i]),
350 //           valuesLeg1[fVarCuts[i]], valuesLeg2[fVarCuts[i]], valuesPair[fVarCuts[i]]); 
351 //     }
352     Fill(nSignal+fSignalsMC->GetEntries(), valuesPair,  valuesLeg1, valuesLeg2);
353   }
354   // only OS at the moment
355   // else if(part1->Charge()>0)
356   //   valuesPair[AliDielectronVarManager::kPairType]=0;
357   // else
358   //   valuesPair[AliDielectronVarManager::kPairType]=2; // if one of the two particles is neutral, the pair will go here
359
360   return;
361 }
362 //______________________________________________
363 void AliDielectronHF::Fill(Int_t pairIndex, const AliDielectronPair *particle)
364 {
365   //
366   // fill histograms for event, pair and daughter cuts and pair types
367   //
368   
369   // only OS pairs in case of MC
370   if(fHasMC && pairIndex!=AliDielectron::kEv1PM) return;
371
372   // only selected pair types in case of data
373   if(!IsPairTypeSelected(pairIndex)) return;
374
375   // get event and pair variables
376   Double_t valuesPair[AliDielectronVarManager::kNMaxValues];
377   AliDielectronVarManager::Fill(particle,valuesPair);
378
379   // get leg variables
380   Double_t valuesLeg1[AliDielectronVarManager::kNMaxValues]={0};
381   AliDielectronVarManager::Fill(particle->GetFirstDaughter(),valuesLeg1);
382   Double_t valuesLeg2[AliDielectronVarManager::kNMaxValues]={0};
383   AliDielectronVarManager::Fill(particle->GetSecondDaughter(),valuesLeg2);
384
385   // fill
386   if(!fHasMC) { Fill(pairIndex, valuesPair,  valuesLeg1, valuesLeg2); }
387   if(fHasMC && fSignalsMC) {
388     for(Int_t i=0; i<fSignalsMC->GetEntries(); i++) {
389       if(AliDielectronMC::Instance()->IsMCTruth(particle, (AliDielectronSignalMC*)fSignalsMC->At(i))) 
390         Fill(i, valuesPair,  valuesLeg1, valuesLeg2);
391     }
392   }
393   
394   return;
395   
396 }
397
398 //______________________________________________
399 void AliDielectronHF::Fill(Int_t Index, Double_t * const valuesPair, Double_t * const valuesLeg1, Double_t * const valuesLeg2)
400 {
401   //
402   // main fill function using index and values as input
403   //
404
405   TObjArray *histArr = static_cast<TObjArray*>(fArrPairType.At(Index));
406   if(!histArr) return;
407
408   Int_t size  = GetNumberOfBins();
409   // loop over all histograms
410   for(Int_t ihist=0; ihist<size; ihist++) {
411
412     Int_t sizeAdd   = 1;
413     Bool_t selected = kTRUE;
414     
415     // loop over all cut variables
416     Int_t nvars = fAxes.GetEntriesFast();
417     for(Int_t ivar=0; ivar<nvars; ivar++) {
418       
419       // get bin limits
420       TVectorD *bins = static_cast<TVectorD*>(fAxes.At(ivar));
421       Int_t nbins    = bins->GetNrows()-1;
422
423       // bin limits for current ivar bin
424       Int_t ibin   = (ihist/sizeAdd)%nbins;
425       Double_t lowEdge = (*bins)[ibin];
426       Double_t upEdge  = (*bins)[ibin+1];
427       switch(fBinType[ivar]) {
428       case kStdBin:     upEdge=(*bins)[ibin+1];     break;
429       case kBinToMax:   upEdge=(*bins)[nbins];      break;
430       case kBinFromMin: lowEdge=(*bins)[0];         break;
431       case kSymBin:     upEdge=(*bins)[nbins-ibin];
432         if(ibin>=((Double_t)(nbins+1))/2) upEdge=(*bins)[nbins]; // to avoid low>up
433         break;
434       }
435
436       // leg variable
437       if(fVarCutType[ivar]) {
438         if( (valuesLeg1[fVarCuts[ivar]] < lowEdge || valuesLeg1[fVarCuts[ivar]] >= upEdge) ||
439             (valuesLeg2[fVarCuts[ivar]] < lowEdge || valuesLeg2[fVarCuts[ivar]] >= upEdge) ) {
440           selected=kFALSE;
441           break;
442         }
443       }
444       else { // pair and event variables
445         if( (valuesPair[fVarCuts[ivar]] < lowEdge || valuesPair[fVarCuts[ivar]] >= upEdge) ) {
446           selected=kFALSE;
447           break;
448         }
449       }
450
451       sizeAdd*=nbins;
452     } //end of var cut loop
453
454     // do not fill the histogram
455     if(!selected) continue;
456
457     // fill the object with Pair and event values (TODO: this needs to be changed)
458     TObjArray *tmp = (TObjArray*) histArr->At(ihist);
459     TString title = tmp->GetName();
460     AliDebug(10,title.Data());
461     for(Int_t i=0; i<tmp->GetEntriesFast(); i++) {
462       AliDielectronHistos::FillValues(tmp->At(i), valuesPair);
463     }
464     //    AliDebug(10,Form("Fill var %d %s value %f in %s \n",fVar,AliDielectronVarManager::GetValueName(fVar),valuesPair[fVar],tmp->GetName()));
465   } //end of hist loop
466
467 }
468
469 //______________________________________________
470 void AliDielectronHF::Init()
471 {
472   //
473   // initialise event buffers
474   //
475
476   // has MC signals
477   fHasMC=AliDielectronMC::Instance()->HasMC();
478   Int_t steps = 0;
479   if(fHasMC) steps=fSignalsMC->GetEntries();
480   if(fStepGenerated) steps*=2; 
481
482   // init pair type array
483   fArrPairType.SetName(Form("%s_HF",GetName()));
484   if(fHasMC) fArrPairType.Expand(steps);
485   else fArrPairType.Expand(AliDielectron::kEv1PMRot+1);
486
487   Int_t size  = GetNumberOfBins();
488   AliDebug(10,Form("Creating a histo array with size %d \n",size));
489
490   Int_t sizeAdd  = 1; 
491
492   // fill object array with the array of bin cells
493   TObjArray *histArr = new TObjArray();
494   histArr->Expand(size);
495
496   //  printf("fRefObj %p \n",fRefObj);
497   // array of histograms to each bin cell
498   for(Int_t ihist=0; ihist<size; ihist++) {
499     histArr->AddAt(fRefObj.Clone(""), ihist);
500     //histArr->AddAt(fRefObj.Clone(Form("h%04d",ihist)), ihist);
501   }
502
503   // loop over all cut variables and do the naming according to it bin cell
504   Int_t nvars = fAxes.GetEntriesFast();
505   for(Int_t ivar=0; ivar<nvars; ivar++) {
506     
507     // get bin limits
508     TVectorD *bins = static_cast<TVectorD*>(fAxes.At(ivar));
509     Int_t nbins    = bins->GetNrows()-1;
510
511
512     // loop over all bin cells an set unique titles
513     for(Int_t ihist=0; ihist<size; ihist++) {
514
515       // get the lower limit for current ivar bin
516       Int_t ibin   = (ihist/sizeAdd)%nbins; 
517       Double_t lowEdge = (*bins)[ibin];
518       Double_t upEdge  = (*bins)[ibin+1];
519       switch(fBinType[ivar]) {
520       case kStdBin:     upEdge=(*bins)[ibin+1];     break;
521       case kBinToMax:   upEdge=(*bins)[nbins];      break;
522       case kBinFromMin: lowEdge=(*bins)[0];         break;
523       case kSymBin:     upEdge=(*bins)[nbins-ibin];
524         if(ibin>=((Double_t)(nbins+1))/2) upEdge=(*bins)[nbins]; // to avoid low>up
525         break;
526       }
527
528       TObjArray *tmp= (TObjArray*) histArr->At(ihist);
529       TString title = tmp->GetName();
530       if(!ivar)             title ="";
531       if( ivar)             title+=":";
532       if(fVarCutType[ivar]) title+="Leg";
533       title+=AliDielectronVarManager::GetValueName(fVarCuts[ivar]);
534       title+=Form("#%.2f#%.2f",lowEdge,upEdge);
535       tmp->SetName(title.Data());
536       AliDebug(10,title.Data());
537     } // end: array of bin cell
538     sizeAdd*=nbins;
539   } //end: cut loop
540
541   // copy array to the selected pair types/ MC sources
542   if(fHasMC) {
543     for(Int_t i=0; i<fSignalsMC->GetEntries(); i++) {
544         TString title = Form("(Signal: %s)",fSignalsMC->At(i)->GetTitle());
545         fArrPairType[i]=(TObjArray*)histArr->Clone(title.Data());
546         if(fStepGenerated)  {
547           title+=" MC truth";
548           fArrPairType[i+fSignalsMC->GetEntries()]=(TObjArray*)histArr->Clone(title.Data());
549         }
550     } // end: loop over sources
551
552    }
553    else {
554     for(Int_t i=0; i<AliDielectron::kEv1PMRot+1; i++) {
555       fArrPairType[i]=(TObjArray*)histArr->Clone(Form("%s",AliDielectron::PairClassName(i)));
556       if(!IsPairTypeSelected(i))  ((TObjArray*)fArrPairType[i])->Delete();
557     } //end: loop over pair types
558
559   }
560   
561   // clean up
562   if(histArr) {
563     delete histArr;
564     histArr=0;
565   }
566
567 }
568
569 //______________________________________________
570 Int_t AliDielectronHF::GetNumberOfBins() const
571 {
572   //
573   // return the number of bins this mixing handler has
574   //
575   Int_t size=1;
576   for (Int_t i=0; i<fAxes.GetEntriesFast(); ++i)
577     size*=((static_cast<TVectorD*>(fAxes.At(i)))->GetNrows()-1);
578   return size;
579 }
580
581 //______________________________________________
582 Bool_t AliDielectronHF::IsPairTypeSelected(Int_t itype)
583 {
584   //
585   // check whether a pair type was selected
586   //
587   
588   Bool_t selected = kFALSE;
589
590   // fill all
591   if(fPairType==kAll) selected=kTRUE;
592
593   switch(itype) {
594   case AliDielectron::kEv1PP: 
595   case AliDielectron::kEv1MM:
596     if(fPairType==kOSandLS)   selected = kTRUE;
597     break;
598   case AliDielectron::kEv1PM: selected = kTRUE;
599     break;
600   case AliDielectron::kEv1PEv2P:
601   case AliDielectron::kEv1MEv2P:
602   case AliDielectron::kEv1PEv2M:
603   case AliDielectron::kEv1MEv2M:
604     if(fPairType==kOSandMIX)  selected = kTRUE;
605     break;
606   case AliDielectron::kEv2PP:
607   case AliDielectron::kEv2PM: 
608   case AliDielectron::kEv2MM:
609     selected = kFALSE;
610     break;
611   case AliDielectron::kEv1PMRot:
612     if(fPairType==kOSandROT)  selected = kTRUE;
613     break;
614   }
615   
616   return selected;
617
618 }
619
620