]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWGGA/GammaConv/AliConversionMesonCuts.cxx
fix coverity 24652
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliConversionMesonCuts.cxx
index cd42639db07d594baef6bb2bc75071a899565cfa..033e98baa382adbe7a06e8f70780677f1ce5e60f 100644 (file)
@@ -58,7 +58,7 @@ const char* AliConversionMesonCuts::fgkCutNames[AliConversionMesonCuts::kNCuts]
    "RapidityMesonCut", //4
    "RCut", //5
    "AlphaMesonCut", //6
-   "Chi2MesonCut", //7
+   "SelectionWindow", //7
    "SharedElectronCuts", //8
    "RejectToCloseV0s", //9
    "UseMCPSmearing", //10
@@ -74,7 +74,8 @@ AliConversionMesonCuts::AliConversionMesonCuts(const char *name,const char *titl
    fHistograms(NULL),
    fMesonKind(0),
    fMaxR(200),
-   fChi2CutMeson(1000),
+   fSelectionLow(0.08),
+   fSelectionHigh(0.145),
    fAlphaMinCutMeson(0),
    fAlphaCutMeson(1),
    fRapidityCutMeson(1),
@@ -101,6 +102,9 @@ AliConversionMesonCuts::AliConversionMesonCuts(const char *name,const char *titl
    fDCAGammaGammaCut(1000),
    fDCAZMesonPrimVtxCut(1000),
    fDCARMesonPrimVtxCut(1000),
+   fDCAGammaGammaCutOn(kFALSE),
+   fDCAZMesonPrimVtxCutOn(kFALSE),
+   fDCARMesonPrimVtxCutOn(kFALSE),
    fBackgroundHandler(0),
    fCutString(NULL),
    hMesonCuts(NULL),
@@ -131,7 +135,8 @@ AliConversionMesonCuts::AliConversionMesonCuts(const AliConversionMesonCuts &ref
    fHistograms(NULL),
    fMesonKind(ref.fMesonKind),
    fMaxR(ref.fMaxR),
-   fChi2CutMeson(ref.fChi2CutMeson),
+   fSelectionLow(ref.fSelectionLow),
+   fSelectionHigh(ref.fSelectionHigh),
    fAlphaMinCutMeson(ref.fAlphaMinCutMeson),
    fAlphaCutMeson(ref.fAlphaCutMeson),
    fRapidityCutMeson(ref.fRapidityCutMeson),
@@ -158,6 +163,9 @@ AliConversionMesonCuts::AliConversionMesonCuts(const AliConversionMesonCuts &ref
    fDCAGammaGammaCut(ref.fDCAGammaGammaCut),
    fDCAZMesonPrimVtxCut(ref.fDCAZMesonPrimVtxCut),
    fDCARMesonPrimVtxCut(ref.fDCARMesonPrimVtxCut),
+   fDCAGammaGammaCutOn(ref.fDCAGammaGammaCutOn),
+   fDCAZMesonPrimVtxCutOn(ref.fDCAZMesonPrimVtxCutOn),
+   fDCARMesonPrimVtxCutOn(ref.fDCARMesonPrimVtxCutOn),
    fBackgroundHandler(ref.fBackgroundHandler),
    fCutString(NULL),
    hMesonCuts(NULL),
@@ -449,6 +457,63 @@ Bool_t AliConversionMesonCuts::MesonIsSelectedMCEtaPiPlPiMiGamma(TParticle *fMCM
        return kFALSE;
 }
 
+//________________________________________________________________________
+Bool_t AliConversionMesonCuts::MesonIsSelectedMCPiPlPiMiPiZero(TParticle *fMCMother,AliStack *fMCStack, Int_t &labelNegPion, Int_t &labelPosPion, Int_t &labelNeutPion, Double_t fRapidityShift){
+
+       // Returns true for all pions within acceptance cuts for decay into 2 photons
+       // If bMCDaughtersInAcceptance is selected, it requires in addition that both daughter photons are within acceptance cuts
+
+       if( !fMCStack )return kFALSE;
+
+       if( !(fMCMother->GetPdgCode() == 221 || fMCMother->GetPdgCode() == 223) ) return kFALSE;
+       
+       if( fMCMother->R()>fMaxR ) return kFALSE; // cuts on distance from collision point
+
+       Double_t rapidity = 10.;
+
+       if( fMCMother->Energy() - fMCMother->Pz() == 0 || fMCMother->Energy() + fMCMother->Pz() == 0 ){
+               rapidity=8.-fRapidityShift;
+       }
+       else{
+               rapidity = 0.5*(TMath::Log((fMCMother->Energy()+fMCMother->Pz()) / (fMCMother->Energy()-fMCMother->Pz())))-fRapidityShift;
+       }
+
+       // Rapidity Cut
+       if( abs(rapidity) > fRapidityCutMeson )return kFALSE;
+
+       // Select only -> pi+ pi- pi0
+       if( fMCMother->GetNDaughters() != 3 )return kFALSE;
+
+       TParticle *posPion = 0x0;
+       TParticle *negPion = 0x0;
+       TParticle *neutPion = 0x0;
+
+//     cout << "\n"<< fMCMother->GetPdgCode() << "\n" << endl;
+       for(Int_t index= fMCMother->GetFirstDaughter();index<= fMCMother->GetLastDaughter();index++){
+               
+               TParticle* temp = (TParticle*)fMCStack->Particle( index );
+//             cout << temp->GetPdgCode() << endl;
+               switch( temp->GetPdgCode() ) {
+               case 211:
+                       posPion      =  temp;
+                       labelPosPion = index;
+                       break;
+               case -211:
+                       negPion      =  temp;
+                       labelNegPion = index;
+                       break;
+               case 111:
+                       neutPion         =  temp;
+                       labelNeutPion    = index;
+                       break;
+               }
+       }
+
+       if( posPion && negPion && neutPion ) return kTRUE;
+       return kFALSE;
+}
+
+
 //________________________________________________________________________
 Bool_t AliConversionMesonCuts::MesonIsSelectedMCChiC(TParticle *fMCMother,AliStack *fMCStack,Int_t & labelelectronChiC, Int_t & labelpositronChiC, Int_t & labelgammaChiC, Double_t fRapidityShift){
    // Returns true for all ChiC within acceptance cuts for decay into JPsi + gamma -> e+ + e- + gamma
@@ -521,108 +586,112 @@ Bool_t AliConversionMesonCuts::MesonIsSelectedMCChiC(TParticle *fMCMother,AliSta
    return kFALSE;
 }
 
-
-
 ///________________________________________________________________________
 Bool_t AliConversionMesonCuts::MesonIsSelected(AliAODConversionMother *pi0,Bool_t IsSignal, Double_t fRapidityShift)
 {
 
-   // Selection of reconstructed Meson candidates
-   // Use flag IsSignal in order to fill Fill different
-   // histograms for Signal and Background
-   TH1 *hist=0x0;
-
-   if(IsSignal){hist=hMesonCuts;}
-   else{hist=hMesonBGCuts;}
+       // Selection of reconstructed Meson candidates
+       // Use flag IsSignal in order to fill Fill different
+       // histograms for Signal and Background
+       TH1 *hist=0x0;
 
-   Int_t cutIndex=0;
+       if(IsSignal){hist=hMesonCuts;}
+       else{hist=hMesonBGCuts;}
 
-   if(hist)hist->Fill(cutIndex);
-   cutIndex++;
-
-   // Undefined Rapidity -> Floating Point exception
-   if((pi0->E()+pi0->Pz())/(pi0->E()-pi0->Pz())<=0){
-      if(hist)hist->Fill(cutIndex);
-      cutIndex++;
-//       cout << "undefined rapidity" << endl;
-      return kFALSE;
-   }
-   else{
-      // PseudoRapidity Cut --> But we cut on Rapidity !!!
-      cutIndex++;
-      if(abs(pi0->Rapidity()-fRapidityShift)>fRapidityCutMeson){
-         if(hist)hist->Fill(cutIndex);
-//              cout << abs(pi0->Rapidity()-fRapidityShift) << ">" << fRapidityCutMeson << endl;
-         return kFALSE;
-      }
-   }
-   cutIndex++;
+       Int_t cutIndex=0;
 
-   // Opening Angle Cut
-   //fOpeningAngle=2*TMath::ATan(0.134/pi0->P());// physical minimum opening angle
-   if( pi0->GetOpeningAngle() < fOpeningAngle){
-//       cout << pi0->GetOpeningAngle() << "<" << fOpeningAngle << endl; 
-      if(hist)hist->Fill(cutIndex);
-      return kFALSE;
-   }
-   cutIndex++;
+       if(hist)hist->Fill(cutIndex);
+       cutIndex++;
 
-   if ( fAlphaPtDepCut == kTRUE ) {
-       fAlphaCutMeson = fFAlphaCut->Eval( pi0->Pt() );
-   }
-   
-   
-   // Alpha Max Cut
-   if(pi0->GetAlpha()>fAlphaCutMeson){
-//        cout << pi0->GetAlpha() << ">" << fAlphaCutMeson << endl; 
-      if(hist)hist->Fill(cutIndex);
-      return kFALSE;
-   }
-   cutIndex++;
+       // Undefined Rapidity -> Floating Point exception
+       if((pi0->E()+pi0->Pz())/(pi0->E()-pi0->Pz())<=0){
+               if(hist)hist->Fill(cutIndex);
+               cutIndex++;
+               if (!IsSignal)cout << "undefined rapidity" << endl;
+               return kFALSE;
+       }
+       else{
+               // PseudoRapidity Cut --> But we cut on Rapidity !!!
+               cutIndex++;
+               if(abs(pi0->Rapidity()-fRapidityShift)>fRapidityCutMeson){
+                       if(hist)hist->Fill(cutIndex);
+//                     if (!IsSignal)   cout << abs(pi0->Rapidity()-fRapidityShift) << ">" << fRapidityCutMeson << endl;
+                       return kFALSE;
+               }
+       }
+       cutIndex++;
+
+       // Opening Angle Cut
+       //fOpeningAngle=2*TMath::ATan(0.134/pi0->P());// physical minimum opening angle
+       if( pi0->GetOpeningAngle() < fOpeningAngle){
+//             if (!IsSignal) cout << pi0->GetOpeningAngle() << "<" << fOpeningAngle << endl; 
+               if(hist)hist->Fill(cutIndex);
+               return kFALSE;
+       }
+       cutIndex++;
 
-   // Alpha Min Cut
-   if(pi0->GetAlpha()<fAlphaMinCutMeson){
-//       cout << pi0->GetAlpha() << "<" << fAlphaMinCutMeson << endl; 
-      if(hist)hist->Fill(cutIndex);
-      return kFALSE;
-   }
-   cutIndex++;
+       if ( fAlphaPtDepCut == kTRUE ) {
+       
+               fAlphaCutMeson = fFAlphaCut->Eval( pi0->Pt() );
+       }
+       
+       
+       // Alpha Max Cut
+       if(pi0->GetAlpha()>fAlphaCutMeson){
+//             if (!IsSignal) cout << pi0->GetAlpha() << ">" << fAlphaCutMeson << endl; 
+               if(hist)hist->Fill(cutIndex);
+               return kFALSE;
+       }
+       cutIndex++;
 
-   if (hDCAGammaGammaMesonBefore)hDCAGammaGammaMesonBefore->Fill(pi0->GetDCABetweenPhotons());
-   if (hDCARMesonPrimVtxBefore)hDCARMesonPrimVtxBefore->Fill(pi0->GetDCARMotherPrimVtx());
+       // Alpha Min Cut
+       if(pi0->GetAlpha()<fAlphaMinCutMeson){
+//             if (!IsSignal)cout << pi0->GetAlpha() << "<" << fAlphaMinCutMeson << endl; 
+               if(hist)hist->Fill(cutIndex);
+               return kFALSE;
+       }
+       cutIndex++;
 
-   if (pi0->GetDCABetweenPhotons() > fDCAGammaGammaCut){
-//       cout << pi0->GetDCABetweenPhotons() << ">" << fDCAGammaGammaCut << endl; 
-      if(hist)hist->Fill(cutIndex);
-      return kFALSE;
-   }
-   cutIndex++;
+       if (hDCAGammaGammaMesonBefore)hDCAGammaGammaMesonBefore->Fill(pi0->GetDCABetweenPhotons());
+       if (hDCARMesonPrimVtxBefore)hDCARMesonPrimVtxBefore->Fill(pi0->GetDCARMotherPrimVtx());
 
-   if (pi0->GetDCARMotherPrimVtx() > fDCARMesonPrimVtxCut){
-//        cout << pi0->GetDCARMotherPrimVtx() << ">" << fDCARMesonPrimVtxCut << endl; 
-      if(hist)hist->Fill(cutIndex);
-      return kFALSE;
-   }
-   cutIndex++;
+       if (fDCAGammaGammaCutOn){
+               if (pi0->GetDCABetweenPhotons() > fDCAGammaGammaCut){
+//                     if (!IsSignal)cout << pi0->GetDCABetweenPhotons() << ">" << fDCAGammaGammaCut << endl; 
+                       if(hist)hist->Fill(cutIndex);
+                       return kFALSE;
+               }
+       }       
+       cutIndex++;
+
+       if (fDCARMesonPrimVtxCutOn){
+               if (pi0->GetDCARMotherPrimVtx() > fDCARMesonPrimVtxCut){
+//                     if (!IsSignal) cout << pi0->GetDCARMotherPrimVtx() << ">" << fDCARMesonPrimVtxCut << endl; 
+                       if(hist)hist->Fill(cutIndex);
+                       return kFALSE;
+               }
+       }       
+       cutIndex++;
 
 
-   if (hDCAZMesonPrimVtxBefore)hDCAZMesonPrimVtxBefore->Fill(pi0->GetDCAZMotherPrimVtx());
+       if (hDCAZMesonPrimVtxBefore)hDCAZMesonPrimVtxBefore->Fill(pi0->GetDCAZMotherPrimVtx());
 
-   if (abs(pi0->GetDCAZMotherPrimVtx()) > fDCAZMesonPrimVtxCut){
-//       cout << pi0->GetDCAZMotherPrimVtx() << ">" << fDCAZMesonPrimVtxCut << endl; 
-      if(hist)hist->Fill(cutIndex);
-      return kFALSE;
-   }
-   cutIndex++;
+       if (fDCAZMesonPrimVtxCutOn){
+               if (abs(pi0->GetDCAZMotherPrimVtx()) > fDCAZMesonPrimVtxCut){
+//                     if (!IsSignal) cout << pi0->GetDCAZMotherPrimVtx() << ">" << fDCAZMesonPrimVtxCut << endl; 
+                       if(hist)hist->Fill(cutIndex);
+                       return kFALSE;
+               }
+       }
+       cutIndex++;
 
 
-   if (hDCAGammaGammaMesonAfter)hDCAGammaGammaMesonAfter->Fill(pi0->GetDCABetweenPhotons());
-   if (hDCARMesonPrimVtxAfter)hDCARMesonPrimVtxAfter->Fill(pi0->GetDCARMotherPrimVtx());
-   if (hDCAZMesonPrimVtxAfter)hDCAZMesonPrimVtxAfter->Fill(pi0->M(),pi0->GetDCAZMotherPrimVtx());
+       if (hDCAGammaGammaMesonAfter)hDCAGammaGammaMesonAfter->Fill(pi0->GetDCABetweenPhotons());
+       if (hDCARMesonPrimVtxAfter)hDCARMesonPrimVtxAfter->Fill(pi0->GetDCARMotherPrimVtx());
+       if (hDCAZMesonPrimVtxAfter)hDCAZMesonPrimVtxAfter->Fill(pi0->M(),pi0->GetDCAZMotherPrimVtx());
 
-   if(hist)hist->Fill(cutIndex);
-   return kTRUE;
+       if(hist)hist->Fill(cutIndex);
+       return kTRUE;
 }
 
 
@@ -664,7 +733,7 @@ Bool_t AliConversionMesonCuts::InitializeCutsFromCutString(const TString analysi
       if(!SetCut(cutIds(ii),fCuts[ii]))return kFALSE;
    }
 
-   //PrintCuts();
+   PrintCutsWithValues();
    return kTRUE;
 }
 ///________________________________________________________________________
@@ -679,9 +748,9 @@ Bool_t AliConversionMesonCuts::SetCut(cutIds cutID, const Int_t value) {
          UpdateCutString();
          return kTRUE;
       } else return kFALSE;
-   case kchi2MesonCut:
-      if( SetChi2MesonCut(value)) {
-         fCuts[kchi2MesonCut] = value;
+   case kSelectionCut:
+      if( SetSelectionWindowCut(value)) {
+         fCuts[kSelectionCut] = value;
          UpdateCutString();
          return kTRUE;
       } else return kFALSE;
@@ -782,6 +851,38 @@ void AliConversionMesonCuts::PrintCuts() {
    }
 }
 
+///________________________________________________________________________
+void AliConversionMesonCuts::PrintCutsWithValues() {
+   // Print out current Cut Selection with values
+       printf("\nMeson cutnumber \n");
+       for(Int_t ic = 0; ic < kNCuts; ic++) {
+               printf("%d",fCuts[ic]);
+       }
+       printf("\n\n");
+       
+       printf("Meson cuts \n");
+       printf("\t |y| < %3.2f \n", fRapidityCutMeson);
+       printf("\t theta_{open} < %3.2f\n", fOpeningAngle);
+       if (!fAlphaPtDepCut) printf("\t %3.2f < alpha < %3.2f\n", fAlphaMinCutMeson, fAlphaCutMeson);
+       if (fDCAGammaGammaCutOn)printf("\t dca_{gamma,gamma} > %3.2f\n", fDCAGammaGammaCut);
+       if (fDCARMesonPrimVtxCutOn)printf("\t dca_{R, prim Vtx} > %3.2f\n", fDCARMesonPrimVtxCut); 
+       if (fDCAZMesonPrimVtxCutOn)printf("\t dca_{Z, prim Vtx} > %3.2f\n\n", fDCAZMesonPrimVtxCut); 
+       printf("\t Meson selection window for further analysis %3.3f > M_{gamma,gamma} > %3.3f\n\n", fSelectionLow, fSelectionHigh); 
+       
+        printf("Meson BG settings \n");
+        if (!fDoBG){
+               if (!fUseRotationMethodInBG  & !fUseTrackMultiplicityForBG) printf("\t BG scheme: mixing V0 mult \n");
+               if (!fUseRotationMethodInBG  & fUseTrackMultiplicityForBG) printf("\t BG scheme: mixing track mult \n");
+               if (fUseRotationMethodInBG )printf("\t BG scheme: rotation \n");
+               if (fdoBGProbability) printf("\t -> use BG probability \n");
+               if (fBackgroundHandler) printf("\t -> use new BG handler \n");
+               printf("\t depth of pool: %d\n", fNumberOfBGEvents);
+               if (fUseRotationMethodInBG )printf("\t degree's for BG rotation: %d\n", fnDegreeRotationPMForBG);
+        }
+        
+}
+
+
 ///________________________________________________________________________
 Bool_t AliConversionMesonCuts::SetMesonKind(Int_t mesonKind){
    // Set Cut
@@ -833,29 +934,40 @@ Bool_t AliConversionMesonCuts::SetRCut(Int_t rCut){
 }
 
 ///________________________________________________________________________
-Bool_t AliConversionMesonCuts::SetChi2MesonCut(Int_t chi2MesonCut){
+Bool_t AliConversionMesonCuts::SetSelectionWindowCut(Int_t selectionCut){
    // Set Cut
-   switch(chi2MesonCut){
-   case 0:  // 100.
-      fChi2CutMeson = 100.;
-      break;
-   case 1:  // 50.
-      fChi2CutMeson = 50.;
-      break;
-   case 2:  // 30.
-      fChi2CutMeson = 30.;
-      break;
+   switch(selectionCut){
+   case 0:  
+      fSelectionLow = 0.08;
+         fSelectionHigh = 0.145;
+      break;
+   case 1:  
+      fSelectionLow = 0.1;
+         fSelectionHigh = 0.145;
+      break;
+   case 2:  
+         fSelectionLow = 0.11;
+         fSelectionHigh = 0.145;
+         break;
    case 3:
-      fChi2CutMeson = 200.;
+         fSelectionLow = 0.12;
+         fSelectionHigh = 0.145;
       break;
    case 4:
-      fChi2CutMeson = 500.;
+         fSelectionLow = 0.1;
+         fSelectionHigh = 0.15;
       break;
    case 5:
-      fChi2CutMeson = 1000.;
+         fSelectionLow = 0.11;
+         fSelectionHigh = 0.15;
+      break;
+   case 6:
+         fSelectionLow = 0.12;
+         fSelectionHigh = 0.15;
       break;
+         
    default:
-      cout<<"Warning: Chi2MesonCut not defined "<<chi2MesonCut<<endl;
+      cout<<"Warning: SelectionCut not defined "<<selectionCut<<endl;
       return kFALSE;
    }
    return kTRUE;
@@ -1213,125 +1325,155 @@ Bool_t AliConversionMesonCuts::SetMCPSmearing(Int_t useMCPSmearing)
 
 ///________________________________________________________________________
 Bool_t AliConversionMesonCuts::SetDCAGammaGammaCut(Int_t DCAGammaGamma){
-   // Set Cut
-   switch(DCAGammaGamma){
-   case 0:  //
-      fDCAGammaGammaCut   = 1000;
-      break;
-   case 1:  //
-      fDCAGammaGammaCut   = 10;
-      break;
-   case 2:  //
-      fDCAGammaGammaCut   = 5;
-      break;
-   case 3:  //
-      fDCAGammaGammaCut   = 4;
-      break;
-   case 4:  //
-      fDCAGammaGammaCut   = 3;
-      break;
-   case 5:  //
-      fDCAGammaGammaCut   = 2.5;
-      break;
-   case 6:  //
-      fDCAGammaGammaCut   = 2;
-      break;
-   case 7:  //
-      fDCAGammaGammaCut   = 1.5;
-      break;
-   case 8:  //
-      fDCAGammaGammaCut   = 1;
-      break;
-   case 9:  //
-      fDCAGammaGammaCut   = 0.5;
-      break;
-   default:
-      cout<<"Warning: DCAGammaGamma not defined "<<DCAGammaGamma<<endl;
-      return kFALSE;
-   }
-   return kTRUE;
+       // Set Cut
+       switch(DCAGammaGamma){
+       case 0:  //
+               fDCAGammaGammaCutOn = kFALSE;
+               fDCAGammaGammaCut   = 1000;
+               break;
+       case 1:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 10;
+               break;
+       case 2:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 5;
+               break;
+       case 3:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 4;
+               break;
+       case 4:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 3;
+               break;
+       case 5:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 2.5;
+               break;
+       case 6:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 2;
+               break;
+       case 7:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 1.5;
+               break;
+       case 8:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 1;
+               break;
+       case 9:  //
+               fDCAGammaGammaCutOn = kTRUE;
+               fDCAGammaGammaCut   = 0.5;
+               break;
+       default:
+               cout<<"Warning: DCAGammaGamma not defined "<<DCAGammaGamma<<endl;
+               return kFALSE;
+       }
+       return kTRUE;
 }
 
 ///________________________________________________________________________
 Bool_t AliConversionMesonCuts::SetDCAZMesonPrimVtxCut(Int_t DCAZMesonPrimVtx){
-   // Set Cut
-   switch(DCAZMesonPrimVtx){
-   case 0:  //
-      fDCAZMesonPrimVtxCut   = 1000;
-      break;
-   case 1:  //
-      fDCAZMesonPrimVtxCut   = 10;
-      break;
-   case 2:  //
-      fDCAZMesonPrimVtxCut   = 5;
-      break;
-   case 3:  //
-      fDCAZMesonPrimVtxCut   = 4;
-      break;
-   case 4:  //
-      fDCAZMesonPrimVtxCut   = 3;
-      break;
-   case 5:  //
-      fDCAZMesonPrimVtxCut   = 2.5;
-      break;
-   case 6:  //
-      fDCAZMesonPrimVtxCut   = 2;
-      break;
-   case 7:  //
-      fDCAZMesonPrimVtxCut   = 1.5;
-      break;
-   case 8:  //
-      fDCAZMesonPrimVtxCut   = 1;
-      break;
-   case 9:  //
-      fDCAZMesonPrimVtxCut   = 0.5;
-      break;
-   default:
-      cout<<"Warning: DCAZMesonPrimVtx not defined "<<DCAZMesonPrimVtx<<endl;
-      return kFALSE;
-   }
-   return kTRUE;
+       // Set Cut
+       switch(DCAZMesonPrimVtx){
+       case 0:  //
+               fDCAZMesonPrimVtxCutOn = kFALSE;
+               fDCAZMesonPrimVtxCut   = 1000;
+               break;
+       case 1:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 10;
+               break;
+       case 2:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 5;
+               break;
+       case 3:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 4;
+               break;
+       case 4:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 3;
+               break;
+       case 5:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 2.5;
+               break;
+       case 6:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 2;
+               break;
+       case 7:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 1.5;
+               break;
+       case 8:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 1;
+               break;
+       case 9:  //
+               fDCAZMesonPrimVtxCutOn = kTRUE;
+               fDCAZMesonPrimVtxCut   = 0.5;
+               break;
+       default:
+               cout<<"Warning: DCAZMesonPrimVtx not defined "<<DCAZMesonPrimVtx<<endl;
+               return kFALSE;
+       }
+       return kTRUE;
 }
 
 ///________________________________________________________________________
 Bool_t AliConversionMesonCuts::SetDCARMesonPrimVtxCut(Int_t DCARMesonPrimVtx){
-   // Set Cut
-   switch(DCARMesonPrimVtx){
-   case 0:  //
-      fDCARMesonPrimVtxCut   = 1000;
-      break;
-   case 1:  //
-      fDCARMesonPrimVtxCut   = 10;
-      break;
-   case 2:  //
-      fDCARMesonPrimVtxCut   = 5;
-      break;
-   case 3:  //
-      fDCARMesonPrimVtxCut   = 4;
-      break;
-   case 4:  //
-      fDCARMesonPrimVtxCut   = 3;
-      break;
-   case 5:  //
-      fDCARMesonPrimVtxCut   = 2.5;
-      break;
-   case 6:  //
-      fDCARMesonPrimVtxCut   = 2;
-      break;
-   case 7:  //
-      fDCARMesonPrimVtxCut   = 1.5;
-      break;
-   case 8:  //
-      fDCARMesonPrimVtxCut   = 1;
-      break;
-   case 9:  //
-      fDCARMesonPrimVtxCut   = 0.5;
-      break;
-   default:
-      cout<<"Warning: DCARMesonPrimVtx not defined "<<DCARMesonPrimVtx<<endl;
-      return kFALSE;
-   }
-   return kTRUE;
+       // Set Cut
+       switch(DCARMesonPrimVtx){
+       case 0:  //
+               fDCARMesonPrimVtxCutOn = kFALSE;
+               fDCARMesonPrimVtxCut   = 1000;
+               break;
+       case 1:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 10;
+               break;
+       case 2:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 5;
+               break;
+       case 3:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 4;
+               break;
+       case 4:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 3;
+               break;
+       case 5:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 2.5;
+               break;
+       case 6:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 2;
+               break;
+       case 7:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 1.5;
+               break;
+       case 8:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 1;
+               break;
+       case 9:  //
+               fDCARMesonPrimVtxCutOn = kTRUE;
+               fDCARMesonPrimVtxCut   = 0.5;
+               break;
+       default:
+               cout<<"Warning: DCARMesonPrimVtx not defined "<<DCARMesonPrimVtx<<endl;
+               return kFALSE;
+       }
+       return kTRUE;
 }
 
 
@@ -1431,3 +1573,89 @@ void AliConversionMesonCuts::SmearParticle(AliAODConversionPhoton* photon)
    photon->SetPz(facPBrem* (1+facPSig)* P*cos(theta)) ;
    photon->SetE(photon->P());
 }
+///________________________________________________________________________
+void AliConversionMesonCuts::SmearVirtualPhoton(AliAODConversionPhoton* photon)
+{
+
+   if (photon==NULL) return;
+   Double_t facPBrem = 1.;
+   Double_t facPSig = 0.;
+
+   Double_t phi=0.;
+   Double_t theta=0.;
+   Double_t P=0.;
+
+
+   P=photon->P();
+   phi=photon->Phi();
+   if( photon->P()!=0){
+      theta=acos( photon->Pz()/ photon->P());
+   }
+
+   if( fPSigSmearing != 0. || fPSigSmearingCte!=0. ){
+      facPSig = TMath::Sqrt(fPSigSmearingCte*fPSigSmearingCte+fPSigSmearing*fPSigSmearing*P*P)*fRandom.Gaus(0.,1.);
+   }
+
+   if( fPBremSmearing != 1.){
+      if(fBrem!=NULL){
+         facPBrem = fBrem->GetRandom();
+      }
+   }
+
+   photon->SetPx(facPBrem* (1+facPSig)* P*sin(theta)*cos(phi)) ;
+   photon->SetPy(facPBrem* (1+facPSig)* P*sin(theta)*sin(phi)) ;
+   photon->SetPz(facPBrem* (1+facPSig)* P*cos(theta)) ;
+   
+}
+///________________________________________________________________________
+TLorentzVector AliConversionMesonCuts::SmearElectron(TLorentzVector particle)
+{
+
+   //if (particle==0) return;
+   Double_t facPBrem = 1.;
+   Double_t facPSig = 0.;
+
+   Double_t phi=0.;
+   Double_t theta=0.;
+   Double_t P=0.;
+
+
+   P=particle.P();
+   
+   
+   phi=particle.Phi();
+   if (phi < 0.) phi += 2. * TMath::Pi();
+   
+   if( particle.P()!=0){
+      theta=acos( particle.Pz()/ particle.P());
+   }
+
+   
+   Double_t fPSigSmearingHalf    =  fPSigSmearing  / 2.0;  //The parameter was set for gammas with 2 particles and here we have just one electron
+   Double_t sqrtfPSigSmearingCteHalf =  fPSigSmearingCte / 2.0 ;  //The parameter was set for gammas with 2 particles and here we have just one electron
+
+   
+   
+   if( fPSigSmearingHalf != 0. || sqrtfPSigSmearingCteHalf!=0. ){
+      facPSig = TMath::Sqrt(sqrtfPSigSmearingCteHalf*sqrtfPSigSmearingCteHalf+fPSigSmearingHalf*fPSigSmearingHalf*P*P)*fRandom.Gaus(0.,1.);
+   }
+
+   if( fPBremSmearing != 1.){
+      if(fBrem!=NULL){
+         facPBrem = fBrem->GetRandom();
+      }
+   }
+   
+   TLorentzVector SmearedParticle;
+   
+   SmearedParticle.SetXYZM( facPBrem* (1+facPSig)* P*sin(theta)*cos(phi) , facPBrem* (1+facPSig)* P*sin(theta)*sin(phi)  , 
+                         facPBrem* (1+facPSig)* P*cos(theta) , TDatabasePDG::Instance()->GetParticle(  ::kElectron   )->Mass()) ;
+   
+   //particle.SetPx(facPBrem* (1+facPSig)* P*sin(theta)*cos(phi)) ;
+   //particle.SetPy(facPBrem* (1+facPSig)* P*sin(theta)*sin(phi)) ;
+   //particle.SetPz(facPBrem* (1+facPSig)* P*cos(theta)) ;
+   
+   return SmearedParticle;
+   
+}
+