]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWGHF/vertexingHF/AliCFTaskVertexingHF.cxx
fixing casts
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / AliCFTaskVertexingHF.cxx
index 816bfd7666cd313693453a52e9ba6b03676dbabd..3e3e49691e5e3c6fe616adfcdc1778ce3622dfc4 100644 (file)
@@ -107,6 +107,7 @@ AliCFTaskVertexingHF::AliCFTaskVertexingHF() :
   fUseZWeight(kFALSE),
   fUseNchWeight(kFALSE),
   fUseTrackletsWeight(kFALSE),
+  fUseMultRatioAsWeight(kFALSE),
   fNvar(0),
   fPartName(""),
   fDauNames(""),
@@ -129,7 +130,11 @@ AliCFTaskVertexingHF::AliCFTaskVertexingHF() :
   fMultiplicityEstimator(kNtrk10),
   fRefMult(9.26),
   fZvtxCorrectedNtrkEstimator(kFALSE),
-  fIsPPData(kFALSE)
+  fIsPPData(kFALSE),
+  fIsPPbData(kFALSE),
+  fUseAdditionalCuts(kFALSE),
+  fUseCutsForTMVA(kFALSE),
+  fUseCascadeTaskForLctoV0bachelor(kFALSE)
 {
   //
   //Default ctor
@@ -164,6 +169,7 @@ AliCFTaskVertexingHF::AliCFTaskVertexingHF(const Char_t* name, AliRDHFCuts* cuts
   fUseZWeight(kFALSE),
   fUseNchWeight(kFALSE),
   fUseTrackletsWeight(kFALSE),
+  fUseMultRatioAsWeight(kFALSE),
   fNvar(0),
   fPartName(""),
   fDauNames(""),
@@ -186,7 +192,11 @@ AliCFTaskVertexingHF::AliCFTaskVertexingHF(const Char_t* name, AliRDHFCuts* cuts
   fMultiplicityEstimator(kNtrk10),
   fRefMult(9.26),
   fZvtxCorrectedNtrkEstimator(kFALSE),
-  fIsPPData(kFALSE)
+  fIsPPData(kFALSE),
+  fIsPPbData(kFALSE),
+  fUseAdditionalCuts(kFALSE),
+  fUseCutsForTMVA(kFALSE),
+  fUseCascadeTaskForLctoV0bachelor(kFALSE)
 {
   //
   // Constructor. Initialization of Inputs and Outputs
@@ -252,6 +262,7 @@ AliCFTaskVertexingHF::AliCFTaskVertexingHF(const AliCFTaskVertexingHF& c) :
   fUseZWeight(c.fUseZWeight),
   fUseNchWeight(c.fUseNchWeight),
   fUseTrackletsWeight(c.fUseTrackletsWeight),
+  fUseMultRatioAsWeight(c.fUseMultRatioAsWeight),
   fNvar(c.fNvar),
   fPartName(c.fPartName),
   fDauNames(c.fDauNames),
@@ -274,7 +285,11 @@ AliCFTaskVertexingHF::AliCFTaskVertexingHF(const AliCFTaskVertexingHF& c) :
   fMultiplicityEstimator(c.fMultiplicityEstimator),
   fRefMult(c.fRefMult),
   fZvtxCorrectedNtrkEstimator(c.fZvtxCorrectedNtrkEstimator),
-  fIsPPData(c.fIsPPData)
+  fIsPPData(c.fIsPPData),
+  fIsPPbData(c.fIsPPbData),
+  fUseAdditionalCuts(c.fUseAdditionalCuts),
+  fUseCutsForTMVA(c.fUseCutsForTMVA),
+  fUseCascadeTaskForLctoV0bachelor(c.fUseCascadeTaskForLctoV0bachelor)
 {
   //
   // Copy Constructor
@@ -310,7 +325,7 @@ void AliCFTaskVertexingHF::Init()
   if(fUseWeight && fUseZWeight) { AliFatal("Can not use at the same time pt and z-vtx weights, please choose"); return; }
   if(fUseWeight && fUseNchWeight) { AliInfo("Beware, using at the same time pt and Nch weights, please check"); }
   if(fUseNchWeight && !fHistoMCNch) { AliFatal("Need to pass the MC Nch distribution to use Nch weights"); return; }
-  if(fUseNchWeight) CreateMeasuredNchHisto();
+  if(fUseNchWeight && !fHistoMeasNch) CreateMeasuredNchHisto();
 
   AliRDHFCuts *copyfCuts = 0x0;
   if (!fCuts){
@@ -442,8 +457,18 @@ void AliCFTaskVertexingHF::Init()
 
   fListProfiles = new TList();
   fListProfiles->SetOwner();
-  TString period[4]={"LHC10b","LHC10c","LHC10d","LHC10e"};
-  for(Int_t i=0; i<4; i++){
+  TString period[4];
+  Int_t nProfiles=4;
+  
+  if (fIsPPbData) { //if pPb, use only two estimator histos
+     period[0] = "LHC13b"; period[1] = "LHC13c";
+     nProfiles = 2;
+  } else {        // else assume pp (four histos for LHC10)
+     period[0] = "LHC10b"; period[1] = "LHC10c"; period[2] = "LHC10d"; period[3] = "LHC10e";
+     nProfiles = 4;
+  }
+
+  for(Int_t i=0; i<nProfiles; i++){
     if(fMultEstimatorAvg[i]){
       TProfile* hprof=new TProfile(*fMultEstimatorAvg[i]);
       hprof->SetName(Form("ProfileTrkVsZvtx%s\n",period[i].Data()));
@@ -466,6 +491,8 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
   PostData(2,fCFManager->GetParticleContainer()) ;
   PostData(3,fCorrelation) ;   
 
+  AliDebug(3,Form("*** Processing event %d\n", fEvents));
+
   if (fFillFromGenerated){
     AliWarning("Flag to fill container with generated value ON ---> dca, d0pi, d0K, d0xd0, cosPointingAngle will be set as dummy!");
   }
@@ -549,8 +576,11 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
   }
 
   AliAODVertex *aodVtx = (AliAODVertex*)aodEvent->GetPrimaryVertex();
-  if (!aodVtx) return;
-       
+  if (!aodVtx) {
+    AliDebug(3, "The event was skipped due to missing vertex");
+    return;
+  }
+
   if (!arrayBranch) {
     AliError("Could not find array of HF vertices");
     return;
@@ -601,10 +631,31 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
   }
   case 21:{ 
     cfVtxHF = new AliCFVertexingHFCascade(mcArray, fOriginDselection);
+    ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGcascade(413);
+    ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGbachelor(211);
+    ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaugh(421);
+    ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaughForMC(421);
+    ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaughPositive(211);
+    ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaughNegative(321);
+    ((AliCFVertexingHFCascade*)cfVtxHF)->SetPrimaryVertex(aodVtx);
     break;
   }
   case 22:{
-    cfVtxHF = new AliCFVertexingHFLctoV0bachelor(mcArray, fOriginDselection,fGenLctoV0bachelorOption); // Lc -> K0S+proton
+    // Lc ->  K0S+proton
+    if (fUseCascadeTaskForLctoV0bachelor){
+      cfVtxHF = new AliCFVertexingHFCascade(mcArray, fOriginDselection);
+      ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGcascade(4122);
+      ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGbachelor(2212);
+      ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaugh(310);
+      ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaughForMC(311);
+      ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaughPositive(211);
+      ((AliCFVertexingHFCascade*)cfVtxHF)->SetPDGneutrDaughNegative(211);
+      ((AliCFVertexingHFCascade*)cfVtxHF)->SetPrimaryVertex(aodVtx);
+      if (fUseAdditionalCuts) ((AliCFVertexingHFCascade*)cfVtxHF)->SetUseCutsForTMVA(fUseCutsForTMVA);
+    }
+    else {
+      cfVtxHF = new AliCFVertexingHFLctoV0bachelor(mcArray, fOriginDselection,fGenLctoV0bachelorOption);
+    }
     break;
   }
   case 31:
@@ -639,8 +690,8 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
 
   // Multiplicity definition with tracklets
   Double_t nTracklets = 0;
-  Int_t nTrackletsEta10 = AliVertexingHFUtils::GetNumberOfTrackletsInEtaRange(aodEvent,-1.,1.);
-  Int_t nTrackletsEta16 = AliVertexingHFUtils::GetNumberOfTrackletsInEtaRange(aodEvent,-1.6,1.6);
+  Int_t nTrackletsEta10 = static_cast<Int_t>(AliVertexingHFUtils::GetNumberOfTrackletsInEtaRange(aodEvent,-1.,1.));
+  Int_t nTrackletsEta16 = static_cast<Int_t>(AliVertexingHFUtils::GetNumberOfTrackletsInEtaRange(aodEvent,-1.6,1.6));
   nTracklets = (Double_t)nTrackletsEta10;
   if(fMultiplicityEstimator==kNtrk10to16) { nTracklets = (Double_t)(nTrackletsEta16 - nTrackletsEta10); }
 
@@ -661,7 +712,7 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
   if(fUseNchWeight){
     Int_t nChargedMCPhysicalPrimary=AliVertexingHFUtils::GetGeneratedPhysicalPrimariesInEtaRange(mcArray,-1.0,1.0);
     if(!fUseTrackletsWeight) fWeight *= GetNchWeight(nChargedMCPhysicalPrimary);
-    else fWeight *= GetNchWeight(nTracklets);
+    else fWeight *= GetNchWeight(static_cast<Int_t>(nTracklets));
     AliDebug(2,Form("Using Nch weights, Mult=%d Weight=%f\n",nChargedMCPhysicalPrimary,fWeight));
   }
   Double_t eventWeight=fWeight;
@@ -670,6 +721,7 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
     AliDebug(3,Form("z coordinate of MC vertex = %f, it was required to be within [-%f, +%f], skipping event", zMCVertex, fCuts->GetMaxVtxZ(), fCuts->GetMaxVtxZ()));
     delete[] containerInput;
     delete[] containerInputMC;
+    delete cfVtxHF;
     return;
   }
 
@@ -678,6 +730,7 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
     AliDebug(3,"Event rejected because of null trigger mask");
     delete[] containerInput;
     delete[] containerInputMC;
+    delete cfVtxHF;
     return;
   }
 
@@ -717,10 +770,12 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
   if (fUseMCVertex) fCuts->SetUseMCVertex(); 
 
   if (fCentralitySelection){ // keep only the requested centrality
+
     if(fCuts->IsEventSelectedInCentrality(aodEvent)!=0) {
       delete[] containerInput;
       delete[] containerInputMC;
       delete [] trackCuts;
+      delete cfVtxHF;
       return;
     }    
   }  else { // keep all centralities
@@ -740,7 +795,6 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
   Double_t multiplicity = nTracklets; // set to the Ntracklet estimator
   if(fMultiplicityEstimator==kVZERO) { multiplicity = vzeroMult; }
 
-
   cfVtxHF->SetMultiplicity(multiplicity);
 
   //  printf("Multiplicity estimator %d, value %2.2f\n",fMultiplicityEstimator,multiplicity);
@@ -760,26 +814,29 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
       AliDebug(2,"Check the MC-level cuts - not desidered particle");
       continue;  // 0 stands for MC level
     }
+    else {
+      AliDebug(3, Form("\n\n---> COOL! we found a particle (particle %d)!!! with PDG code = %d \n\n", iPart, mcPart->GetPdgCode()));
+    }
     cfVtxHF->SetMCCandidateParam(iPart);
         
          
     if (!(cfVtxHF->SetLabelArray())){
-      AliDebug(2,Form("Impossible to set the label array (decaychannel = %d)",fDecayChannel));
+      AliDebug(2,Form("Impossible to set the label array for particle %d (decaychannel = %d)", iPart, fDecayChannel));
       continue;
     }             
 
     //check the candiate family at MC level
     if (!(cfVtxHF->CheckMCPartFamily(mcPart, mcArray))) {
-      AliDebug(2,Form("Check on the family wrong!!! (decaychannel = %d)",fDecayChannel));
+      AliDebug(2,Form("Check on the family wrong for particle %d!!! (decaychannel = %d)", iPart, fDecayChannel));
       continue;
     }
     else{
-      AliDebug(2,Form("Check on the family OK!!! (decaychannel = %d)",fDecayChannel));
+      AliDebug(2,Form("Check on the family OK for particle %d!!! (decaychannel = %d)", iPart, fDecayChannel));
     }
                
     //Fill the MC container
     Bool_t mcContainerFilled = cfVtxHF -> FillMCContainer(containerInputMC);
-    AliDebug(2,Form("mcContainerFilled = %d)",mcContainerFilled));
+    AliDebug(2, Form("particle = %d mcContainerFilled = %d",iPart, mcContainerFilled));
     if (mcContainerFilled) {
       if (fUseWeight){
        if (fFuncWeight){ // using user-defined function
@@ -792,17 +849,24 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
        }
        AliDebug(2,Form("pt = %f, weight = %f",containerInputMC[0], fWeight));
       }
-      if (!fCuts->IsInFiducialAcceptance(containerInputMC[0],containerInputMC[1])) continue;
+      if (!fCuts->IsInFiducialAcceptance(containerInputMC[0],containerInputMC[1])) {
+       AliDebug(3, Form("Not in limited acceptance, containerInputMC[0] = %f, containerInputMC[1] = %f", containerInputMC[0], containerInputMC[1]));
+       continue;
+      }
+      else{
+       AliDebug(3, Form("YES!! in limited acceptance, containerInputMC[0] = %f, containerInputMC[1] = %f", containerInputMC[0],containerInputMC[1]));
+      }
+
       //MC Limited Acceptance
       if (TMath::Abs(containerInputMC[1]) < 0.5) {
        fCFManager->GetParticleContainer()->Fill(containerInputMC,kStepGeneratedLimAcc, fWeight);
        AliDebug(3,"MC Lim Acc container filled\n");
-      }            
+      }
                        
       //MC 
       fCFManager->GetParticleContainer()->Fill(containerInputMC, kStepGenerated, fWeight);
       icountMC++;
-      AliDebug(3,"MC cointainer filled \n");
+      AliDebug(3,"MC container filled \n");
                        
       // MC in acceptance
       // check the MC-Acceptance level cuts
@@ -895,19 +959,21 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
                
     Bool_t signAssociation = cfVtxHF->SetRecoCandidateParam((AliAODRecoDecayHF*)charmCandidate);
     if (!signAssociation){
-      charmCandidate = 0x0;
+      if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
       continue;
     }
 
     Int_t isPartOrAntipart = cfVtxHF->CheckReflexion(fSign);
     if (isPartOrAntipart == 0){
       AliDebug(2, Form("The candidate pdg code doesn't match the requirement set in the task (fSign = %d)",fSign));
+      if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
       continue;
     }
 
     AliDebug(3,Form("iCandid=%d - signAssociation=%d, isPartOrAntipart=%d",iCandid, signAssociation, isPartOrAntipart));
 
     Bool_t recoContFilled = cfVtxHF->FillRecoContainer(containerInput);
+    AliDebug(3, Form("CF task: RecoContFilled for candidate %d is %d", iCandid, (Int_t)recoContFilled));
     if (recoContFilled){
 
       // weight according to pt
@@ -923,10 +989,13 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
        AliDebug(2, Form("pt = %f, weight = %f",containerInput[0], fWeight));
       }
 
-      if (!fCuts->IsInFiducialAcceptance(containerInput[0],containerInput[1])) continue;
-                       
+      if (!fCuts->IsInFiducialAcceptance(containerInput[0],containerInput[1])){
+       if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
+       continue;
+      }                
       //Reco Step
       Bool_t recoStep = cfVtxHF->RecoStep();
+      if (recoStep) AliDebug(2, Form("CF task: Reco step for candidate %d is %d", iCandid, (Int_t)recoStep));
       Bool_t vtxCheck = fCuts->IsEventSelected(aodEvent);
                        
 
@@ -939,8 +1008,10 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
       }else if(fDecayChannel==33){
        if(fUseSelectionBit && !charmCandidate->HasSelectionBit(AliRDHFCuts::kDsCuts)) isBitSelected = kFALSE;
       }
-      if(!isBitSelected) continue;
-
+      if(!isBitSelected){
+       if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
+       continue;
+      }
 
 
       if (recoStep && recoContFilled && vtxCheck){
@@ -980,9 +1051,13 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
            else if (fDecayChannel==22){ // Lc->V0+bachelor case, where more possibilities are considered
              Bool_t keepLctoV0bachelor=ProcessLctoV0Bachelor(recoAnalysisCuts);
              if (keepLctoV0bachelor) recoAnalysisCuts=3;
-           }
-
-                                                   
+             AliAnalysisManager *man = AliAnalysisManager::GetAnalysisManager();
+             AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
+             AliPIDResponse* pidResponse = inputHandler->GetPIDResponse();
+             if (fUseAdditionalCuts){
+               if (!((AliCFVertexingHFCascade*)cfVtxHF)->CheckAdditionalCuts(pidResponse)) recoAnalysisCuts = -1;
+             }
+           }                                       
 
            fCuts->SetUsePID(iscutsusingpid); //restoring usage of the PID from the cuts object 
            Bool_t tempAn=(recoAnalysisCuts == 3 || recoAnalysisCuts == isPartOrAntipart);
@@ -1009,7 +1084,34 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
              if (fDecayChannel == 32) tempPid=(recoPidSelection >0 || recoPidSelection == isPartOrAntipart);
 
              if (tempPid){
-               fCFManager->GetParticleContainer()->Fill(containerInput, kStepRecoPID, fWeight);
+               Double_t weigPID = 1.;
+               if (fDecayChannel == 2){ // D0 with Bayesian PID using weights
+                 if(((AliRDHFCutsD0toKpi*)fCuts)->GetCombPID() && (((AliRDHFCutsD0toKpi*)fCuts)->GetBayesianStrategy() == AliRDHFCutsD0toKpi::kBayesWeight || ((AliRDHFCutsD0toKpi*)fCuts)->GetBayesianStrategy() == AliRDHFCutsD0toKpi::kBayesWeightNoFilter)){
+                   if (isPartOrAntipart == 1){
+                     weigPID = ((AliRDHFCutsD0toKpi*)fCuts)->GetWeightsNegative()[AliPID::kKaon] * ((AliRDHFCutsD0toKpi*)fCuts)->GetWeightsPositive()[AliPID::kPion];
+                   }else if (isPartOrAntipart == 2){
+                     weigPID = ((AliRDHFCutsD0toKpi*)fCuts)->GetWeightsPositive()[AliPID::kKaon] * ((AliRDHFCutsD0toKpi*)fCuts)->GetWeightsNegative()[AliPID::kPion];
+                   }
+                   if ((weigPID  < 0) || (weigPID > 1)) weigPID = 0.;
+                 }
+               }else if (fDecayChannel == 33){ // Ds with Bayesian PID using weights
+                 if(((AliRDHFCutsDstoKKpi*)fCuts)->GetPidOption()==AliRDHFCutsDstoKKpi::kBayesianWeights){
+                   Int_t labDau0=((AliAODTrack*)charmCandidate->GetDaughter(0))->GetLabel();
+                   AliAODMCParticle* firstDau=(AliAODMCParticle*)mcArray->UncheckedAt(TMath::Abs(labDau0));
+                   if(firstDau){
+                     Int_t pdgCode0=TMath::Abs(firstDau->GetPdgCode());
+                     if(pdgCode0==321){
+                       weigPID=((AliRDHFCutsDstoKKpi*)fCuts)->GetWeightForKKpi();
+                     }else if(pdgCode0==211){
+                       weigPID=((AliRDHFCutsDstoKKpi*)fCuts)->GetWeightForpiKK();
+                     }
+                     if ((weigPID  < 0) || (weigPID > 1)) weigPID = 0.;
+                   }else{
+                     weigPID=0.;
+                   }
+                 }
+               }
+               fCFManager->GetParticleContainer()->Fill(containerInput, kStepRecoPID, fWeight*weigPID);
                icountRecoPID++;
                AliDebug(3,"Reco PID cuts passed and container filled \n");
                if(!fAcceptanceUnf){
@@ -1020,26 +1122,31 @@ void AliCFTaskVertexingHF::UserExec(Option_t *)
              }
              else {
                AliDebug(3, "Analysis Cuts step not passed \n");
+               if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
                continue;
              }
            }
            else {
              AliDebug(3, "PID selection not passed \n");
+             if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
              continue;
            }
          }
          else{
            AliDebug(3, "Number of ITS cluster step not passed\n");
+           if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
            continue;
          }
        }
        else{
          AliDebug(3, "Reco acceptance step not passed\n");
+         if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
          continue;
        }
       }
       else {
        AliDebug(3, "Reco step not passed\n");
+       if(unsetvtx) charmCandidate->UnsetOwnPrimaryVtx();
        continue;
       }
     }
@@ -1153,22 +1260,25 @@ void AliCFTaskVertexingHF::Terminate(Option_t*)
     } else if (fDecayChannel==22) {
       //nvarToPlot = 16;
       titles = new TString[nvarToPlot];
-      titles[0]="pT_Lc (GeV/c)";
-      titles[1]="rapidity";
-      titles[2]="phi (rad)";
-      titles[3]="cosPAV0";
-      titles[4]="onTheFlyStatusV0";
+      titles[0]="p_{T}(#Lambda_{c}) [GeV/c]";
+      titles[1]="y(#Lambda_{c})";
+      titles[2]="#varphi(#Lambda_{c}) [rad]";
+      titles[3]="onTheFlyStatusV0";
+      titles[4]="z_{vtx} [cm]";
       titles[5]="centrality";
       titles[6]="fake";
       titles[7]="multiplicity";
-      titles[8]="pT_bachelor (GeV/c)";
-      titles[9]="pT_V0pos (GeV/c)";
-      titles[10]="pT_V0neg (GeV/c)";
-      titles[11]="invMassV0 (GeV/c2)";
-      titles[12]="dcaV0 (nSigma)";
-      titles[13]="c#tauV0 (#mum)";
-      titles[14]="c#tau (#mum)";
-      titles[15]="cosPA";
+      //titles[8]="pT(bachelor) [GeV/c]";
+      titles[8]="p(bachelor) [GeV/c]";
+      titles[9]="p_{T}(V0) [GeV/c]";
+      titles[10]="y(V0)";
+      titles[11]="#varphi(V0) [rad]";
+      titles[12]="m_{inv}(#pi^{+}#pi^{+}) [GeV/c^{2}]";
+      titles[13]="dcaV0 (nSigma)";
+      titles[14]="cosine pointing angle (V0)";
+      titles[15]="cosine pointing angle (#Lambda_{c})";
+      //titles[16]="c#tauV0 (#mum)";
+      //titles[17]="c#tau (#mum)";
     } else {
       //nvarToPlot = 12;
       titles = new TString[nvarToPlot];
@@ -1190,11 +1300,11 @@ void AliCFTaskVertexingHF::Terminate(Option_t*)
     //nvarToPlot = 8;
     titles = new TString[nvarToPlot];
     if (fDecayChannel==22) {
-      titles[0]="pT_candidate (GeV/c)";
-      titles[1]="rapidity";
-      titles[2]="phi (rad)";
-      titles[3]="cosPAV0";
-      titles[4]="onTheFlyStatusV0";
+      titles[0]="p_{T}(#Lambda_{c}) [GeV/c]";
+      titles[1]="y(#Lambda_{c})";
+      titles[2]="#varphi(#Lambda_{c}) [rad]";
+      titles[3]="onTheFlyStatusV0";
+      titles[4]="z_{vtx} [cm]";
       titles[5]="centrality";
       titles[6]="fake";
       titles[7]="multiplicity";
@@ -1308,7 +1418,7 @@ void AliCFTaskVertexingHF::Terminate(Option_t*)
   file_projection->Close();
   for (Int_t ih = 0; ih<3; ih++) delete [] h[ih];
   delete [] h;
-
+  delete [] titles;
        
 }
 
@@ -1319,23 +1429,6 @@ void AliCFTaskVertexingHF::UserCreateOutputObjects()
   //TO BE SET BEFORE THE EXECUTION OF THE TASK
   //
   Info("UserCreateOutputObjects","CreateOutputObjects of task %s\n", GetName());
-       
-  AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
-  AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
-  AliPIDResponse *localPIDResponse = (AliPIDResponse*)inputHandler->GetPIDResponse();
-
-  if (fCuts->GetIsUsePID() && fDecayChannel==22) {
-         
-    fCuts->GetPidHF()->SetPidResponse(localPIDResponse);
-    fCuts->GetPidHF()->SetOldPid(kFALSE);
-    AliRDHFCutsLctoV0* lcv0Cuts=dynamic_cast<AliRDHFCutsLctoV0*>(fCuts);
-    if(lcv0Cuts){
-      lcv0Cuts->GetPidV0pos()->SetPidResponse(localPIDResponse);
-      lcv0Cuts->GetPidV0neg()->SetPidResponse(localPIDResponse);
-      lcv0Cuts->GetPidV0pos()->SetOldPid(kFALSE);
-      lcv0Cuts->GetPidV0neg()->SetOldPid(kFALSE);
-    }
-  }
 
   //slot #1
   OpenFile(1);
@@ -1399,7 +1492,7 @@ void AliCFTaskVertexingHF::SetPtWeightsFromFONLL276andBAMPSoverLHC12a17b(){
 
 //_________________________________________________________________________
 void AliCFTaskVertexingHF::SetPtWeightsFromFONLL5overLHC13d3(){
-  // weight function from the ratio of the LHC12a17b MC
+  // weight function from the ratio of the LHC13d3 MC
   // and FONLL calculations for pp data
   if(fFuncWeight) delete fFuncWeight;
   fFuncWeight=new TF1("funcWeight","([0]*x)/TMath::Power([2],(1+TMath::Power([3],x/[1])))+[4]*TMath::Exp([5]+[6]*x)+[7]*TMath::Exp([8]*x)",0.15,30.);
@@ -1407,6 +1500,66 @@ void AliCFTaskVertexingHF::SetPtWeightsFromFONLL5overLHC13d3(){
   fUseWeight=kTRUE;
 }
 
+//_________________________________________________________________________
+void AliCFTaskVertexingHF::SetPtWeightsFromFONLL7overLHC10f6a(){
+  // weight function from the ratio of the LHC10f6a MC
+  // and FONLL calculations for pp data
+  if(fFuncWeight) delete fFuncWeight;
+  fFuncWeight=new TF1("funcWeight","([0]*x)/TMath::Power([2],(1+TMath::Power([3],x/[1])))+[4]*TMath::Exp([5]+[6]*x)+[7]*TMath::Exp([8]*x)",0.15,40.);
+  fFuncWeight->SetParameters(2.41522e+01,4.92146e+00,6.72495e+00,2.5,6.15361e-03,4.78995e+00,-4.29135e-01,3.99421e-01,-1.57220e-02);
+  fUseWeight=kTRUE;
+}
+
+//_________________________________________________________________________
+void AliCFTaskVertexingHF::SetPtWeightsFromFONLL7overLHC12a12(){
+  // weight function from the ratio of the LHC12a12 MC
+  // and FONLL calculations for pp data
+  if(fFuncWeight) delete fFuncWeight;
+  fFuncWeight=new TF1("funcWeight","([0]*x)/TMath::Power([2],(1+TMath::Power([3],x/[1])))+[4]*TMath::Exp([5]+[6]*x)+[7]*TMath::Exp([8]*x+[9])",0.15,50.);
+  fFuncWeight->SetParameters(1.31497e+01,3.30503e+00,3.45594e+00,2.5,2.28642e-02,1.42372e+00,2.32892e-04,5.21986e-02,-2.14236e-01,3.86200e+00);
+  fUseWeight=kTRUE;
+}
+
+//_________________________________________________________________________
+void AliCFTaskVertexingHF::SetPtWeightsFromFONLL7overLHC12a12bis(){
+  // weight function from the ratio of the LHC12a12bis MC
+  // and FONLL calculations for pp data
+  if(fFuncWeight) delete fFuncWeight;
+  fFuncWeight=new TF1("funcWeight","([0]*x)/TMath::Power([2],(1+TMath::Power([3],x/[1])))+[4]*TMath::Exp([5]+[6]*x)+[7]*TMath::Exp([8]*x+[9])",0.15,50.);
+  fFuncWeight->SetParameters(6.54519e+00,2.74007e+00,2.48325e+00,2.5,1.61113e-01,-5.32546e-01,-3.75916e-04,2.38189e-01,-2.17561e-01,2.35975e+00);
+  fUseWeight=kTRUE;
+}
+
+//_________________________________________________________________________
+void AliCFTaskVertexingHF::SetPtWeightsFromFONLL7overLHC13e2fix(){
+  // weight function from the ratio of the LHC13e2fix MC
+  // and FONLL calculations for pp data
+  if(fFuncWeight) delete fFuncWeight;
+  fFuncWeight=new TF1("funcWeight","([0]*x)/TMath::Power([2],(1+TMath::Power([3],x/[1])))+[4]*TMath::Exp([5]+[6]*x)+[7]*TMath::Exp([8]*x+[9])",0.15,50.);
+  fFuncWeight->SetParameters(1.85862e+01,2.48171e+00,3.39356e+00,2.5,1.70426e-02,2.28453e+00,-4.57749e-02,5.84585e-02,-3.19719e-01,4.16789e+00);
+  fUseWeight=kTRUE;
+}
+
+//________________________________________________________________
+void AliCFTaskVertexingHF::SetPtWeightsFromFONLL5overLHC10f6a(){
+  // weight function from the ratio of the LHC10f6a MC
+  // and FONLL calculations for pp data
+  if(fFuncWeight) delete fFuncWeight;
+  fFuncWeight=new TF1("funcWeight","([0]*x)/TMath::Power([2],(1+TMath::Power([3],x/[1])))+[4]*TMath::Exp([5]+[6]*x)+[7]*TMath::Exp([8]*x)",0.15,40.);
+  fFuncWeight->SetParameters(2.77730e+01,4.78942e+00,7.45378e+00,2.5,9.86255e-02,2.30120e+00,-4.16435e-01,3.43770e-01,-2.29380e-02);
+  fUseWeight=kTRUE;
+}
+
+//________________________________________________________________
+void AliCFTaskVertexingHF::SetPtWeightsFromFONLL276overLHC10f6a(){
+  // weight function from the ratio of the LHC10f6a MC
+  // and FONLL calculations for pp data
+  if(fFuncWeight) delete fFuncWeight;
+  fFuncWeight=new TF1("funcWeight","([0]*x)/TMath::Power([2],(1+TMath::Power([3],x/[1])))+[4]*TMath::Exp([5]+[6]*x)+[7]*TMath::Exp([8]*x+[9])",0.15,40.);
+  fFuncWeight->SetParameters(1.34412e+01,3.20068e+00,5.14481e+00,2.5,7.59405e-04,7.51821e+00,-3.93811e-01,2.16849e-02,-3.37768e-02,2.40308e+00);
+  fUseWeight=kTRUE;
+}
+
 //_________________________________________________________________________
 Double_t AliCFTaskVertexingHF::GetWeight(Float_t pt)
 {
@@ -1488,7 +1641,7 @@ Double_t AliCFTaskVertexingHF::GetNchWeight(Int_t nch){
   Double_t pMeas=fHistoMeasNch->GetBinContent(fHistoMeasNch->FindBin(nch));
   Double_t pMC=fHistoMCNch->GetBinContent(fHistoMCNch->FindBin(nch));
   Double_t weight = pMC>0 ? pMeas/pMC : 0.;
-  if(fUseTrackletsWeight)  weight = pMC;
+  if(fUseMultRatioAsWeight)  weight = pMC;
   return weight;
 }
 //__________________________________________________________________________________________________
@@ -1584,12 +1737,20 @@ TProfile* AliCFTaskVertexingHF::GetEstimatorHistogram(const AliVEvent* event){
   //
 
   Int_t runNo  = event->GetRunNumber();
-  Int_t period = -1;   // 0-LHC10b, 1-LHC10c, 2-LHC10d, 3-LHC10e
-  if(runNo>114930 && runNo<117223) period = 0;
-  if(runNo>119158 && runNo<120830) period = 1;
-  if(runNo>122373 && runNo<126438) period = 2;
-  if(runNo>127711 && runNo<130841) period = 3;
-  if(period<0 || period>3) return 0;
+  Int_t period = -1;   // pp:  0-LHC10b, 1-LHC10c, 2-LHC10d, 3-LHC10e
+                       // pPb: 0-LHC13b, 1-LHC13c
+
+  if (fIsPPbData) {    // setting run numbers for LHC13 if pPb
+      if (runNo>195343 && runNo<195484) period = 0;
+      if (runNo>195528 && runNo<195678) period = 1;
+      if (period<0 || period>1) return 0;
+  } else {             //else assume pp 2010                 
+      if(runNo>114930 && runNo<117223) period = 0;
+      if(runNo>119158 && runNo<120830) period = 1;
+      if(runNo>122373 && runNo<126438) period = 2;
+      if(runNo>127711 && runNo<130841) period = 3;
+      if(period<0 || period>3) return 0;
+  }
 
   return fMultEstimatorAvg[period];
 }