]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Possibility of using an eta gap between subevents + VZERO plane updates (Robert,...
authorprino <prino@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 30 Sep 2011 22:16:03 +0000 (22:16 +0000)
committerprino <prino@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 30 Sep 2011 22:16:03 +0000 (22:16 +0000)
PWG3/vertexingHF/charmFlow/AddTaskHFv2.C
PWG3/vertexingHF/charmFlow/AliAnalysisTaskSEHFv2.cxx
PWG3/vertexingHF/charmFlow/AliAnalysisTaskSEHFv2.h

index 4d75e5d81ce6bfee304e048688c52fdc95b1523c..1205c7d797fdaef8f499bcd3e79568a6ffef7d33 100644 (file)
@@ -1,4 +1,4 @@
-AliAnalysisTaskSEHFv2 *AddTaskHFv2(TString filename="DplustoKpipiCuts.root", AliAnalysisTaskSEHFv2::DecChannel decCh=AliAnalysisTaskSEHFv2::kD0toKpi,Bool_t readMC=kFALSE,TString name="",Int_t flagep=0 /*0=tracks,1=V0*/)
+AliAnalysisTaskSEHFv2 *AddTaskHFv2(TString filename="DplustoKpipiCuts.root", AliAnalysisTaskSEHFv2::DecChannel decCh=AliAnalysisTaskSEHFv2::kD0toKpi,Bool_t readMC=kFALSE,TString name="",Int_t flagep=1 /*0=tracks,1=V0*/)
 {
   //
   // Test macro for the AliAnalysisTaskSE for  D 
@@ -71,22 +71,24 @@ AliAnalysisTaskSEHFv2 *AddTaskHFv2(TString filename="DplustoKpipiCuts.root", Ali
   Float_t pi=TMath::Pi();
   Float_t philimits[nphibinlimits]={0., pi/4.,pi/2., 3./4.*pi, pi};
 
-  //histogram for V0
-  TFile *fpar = TFile::Open("VZEROParHist.root");
-  TH2D *hh[6];
-  for(Int_t i=0;i<6;i++){
-    TString hhname;hhname.Form("parhist%d_%d",(i+2)*10,(i+3)*10);
-    hh[i]=(TH2D*)fpar->Get(hhname.Data());
-  }
-
   // Analysis task    
-  AliAnalysisTaskSEHFv2 *v2Task = new AliAnalysisTaskSEHFv2("HFv2Analysis",analysiscuts,decCh,nphibinlimits,philimits,hh);
+  AliAnalysisTaskSEHFv2 *v2Task = new AliAnalysisTaskSEHFv2("HFv2Analysis",analysiscuts,decCh,nphibinlimits,philimits);
   v2Task->SetReadMC(readMC);
 
+  v2Task->SetEtaGapFeatureForEventplaneFromTracks(kTRUE);
+  
   v2Task->SetDebugLevel(0);
   
-  v2Task->SetUseV0EP(flagep);
-
+  if(flagep){
+    //histogram for V0
+    TFile *fpar = TFile::Open("VZEROParHist.root");
+    TH2D *hh[6];
+    for(Int_t i=0;i<6;i++){
+      TString hhname;hhname.Form("parhist%d_%d",(i+2)*10,(i+3)*10);
+      hh[i]=(TH2D*)fpar->Get(hhname.Data());
+    } 
+    v2Task->SetVZEROParHist(hh);
+  }
   mgr->AddTask(v2Task);
 
   // Create containers for input/output
@@ -105,9 +107,6 @@ AliAnalysisTaskSEHFv2 *AddTaskHFv2(TString filename="DplustoKpipiCuts.root", Ali
 
   contname=Form("cutobj%s",suffix.Data());
   AliAnalysisDataContainer *cutobj = mgr->CreateContainer(contname.Data(),AliRDHFCuts::Class(),AliAnalysisManager::kOutputContainer,outputfile.Data());
-  
-  contname=Form("coutputVZEROpar%s",suffix.Data());
-  AliAnalysisDataContainer *coutputpar = mgr->CreateContainer(contname.Data(),TList::Class(),AliAnalysisManager::kOutputContainer,outputfile.Data());
 
   mgr->ConnectInput(v2Task,0,mgr->GetCommonInputContainer());
   
@@ -117,7 +116,11 @@ AliAnalysisTaskSEHFv2 *AddTaskHFv2(TString filename="DplustoKpipiCuts.root", Ali
 
   mgr->ConnectOutput(v2Task,3,cutobj);
  
-  mgr->ConnectOutput(v2Task,4,coutputpar);
+  if(flagep){
+    contname=Form("coutputVZEROpar%s",suffix.Data());
+    AliAnalysisDataContainer *coutputpar = mgr->CreateContainer(contname.Data(),TList::Class(),AliAnalysisManager::kOutputContainer,outputfile.Data());
+    mgr->ConnectOutput(v2Task,4,coutputpar);
+  }
 
   return v2Task;
 }
index 97d21396280afce8a3c0b9c3aaffd005c0ac3af7..b913d9132fd47a2d85c65a83b7e89f53903edb72 100644 (file)
@@ -70,7 +70,7 @@
 #include "AliFlowTrack.h"
 #include "AliFlowVector.h"
 #include "AliFlowTrackCuts.h"
-#include "AliFlowEvent.h"
+#include "AliFlowEvent.h"  
 
 #include "AliAnalysisTaskSEHFv2.h"
 
@@ -92,20 +92,21 @@ AliAnalysisTaskSE(),
   fCentLowLimit(0),
   fCentUpLimit(100),
   fNMassBins(200),
-  fReadMC(kFALSE),
+  fReadMC(kFALSE),    
   fUseAfterBurner(kFALSE),
   fDecChannel(0),
   fAfterBurner(0),
   fUseV0EP(kFALSE),
   fV0EPorder(2),
   fMinCentr(10),
-  fMaxCentr(80)
+  fMaxCentr(80),
+  fEtaGap(kFALSE)
 {
   // Default constructor
 }
 
 //________________________________________________________________________
-AliAnalysisTaskSEHFv2::AliAnalysisTaskSEHFv2(const char *name,AliRDHFCuts *rdCuts,Int_t decaychannel, Int_t nlimsphibin, Float_t *phibinlimits,TH2D **histPar):
+AliAnalysisTaskSEHFv2::AliAnalysisTaskSEHFv2(const char *name,AliRDHFCuts *rdCuts,Int_t decaychannel, Int_t nlimsphibin, Float_t *phibinlimits):
   AliAnalysisTaskSE(name),
   fhEventsInfo(0),
   fOutput(0),
@@ -122,10 +123,12 @@ AliAnalysisTaskSEHFv2::AliAnalysisTaskSEHFv2(const char *name,AliRDHFCuts *rdCut
   fReadMC(kFALSE),
   fUseAfterBurner(kFALSE),
   fDecChannel(decaychannel),
+  fAfterBurner(0),
   fUseV0EP(kFALSE),
   fV0EPorder(2),
   fMinCentr(10),
-  fMaxCentr(80)
+  fMaxCentr(80),
+  fEtaGap(kFALSE)
 {
 
   Int_t pdg=421;
@@ -150,8 +153,6 @@ AliAnalysisTaskSEHFv2::AliAnalysisTaskSEHFv2(const char *name,AliRDHFCuts *rdCut
     break;
   }
   fAfterBurner = new AliHFAfterBurner(fDecChannel);
-  for(Int_t i=0;i<6;i++)fHistvzero[i]=(TH2D*)histPar[i]->Clone();
-  for(Int_t i=0;i<6;i++)if(!fHistvzero[i])printf("No VZERO histograms!\n");
   if(pdg==413) SetMassLimits((Float_t)0.135,(Float_t)0.165);
   else SetMassLimits((Float_t)0.2,pdg); //check range
   fNPtBins=fRDCuts->GetNPtBins();
@@ -179,7 +180,7 @@ AliAnalysisTaskSEHFv2::AliAnalysisTaskSEHFv2(const char *name,AliRDHFCuts *rdCut
     break;
   }
   //DefineOutput(4,AliFlowEventSimple::Class());
-  DefineOutput(4,TList::Class());
+  //DefineOutput(4,TList::Class());
 }
 
 //________________________________________________________________________
@@ -216,7 +217,19 @@ AliAnalysisTaskSEHFv2::~AliAnalysisTaskSEHFv2()
     fAfterBurner=0;  
   }
 }
-
+//_________________________________________________________________
+void  AliAnalysisTaskSEHFv2::SetVZEROParHist(TH2D** histPar){
+  for(Int_t i=0;i<6;i++)fHistvzero[i]=(TH2D*)histPar[i]->Clone();
+  for(Int_t i=0;i<6;i++){
+    if(!fHistvzero[i]){
+      printf("No VZERO histograms!\n");
+      fUseV0EP=kFALSE;
+      return;
+    }
+  }
+  DefineOutput(4,TList::Class());
+  fUseV0EP=kTRUE;
+}
 //_________________________________________________________________
 void  AliAnalysisTaskSEHFv2::SetMassLimits(Float_t range, Int_t pdg){
   Float_t mass=0;
@@ -366,8 +379,11 @@ void AliAnalysisTaskSEHFv2::UserCreateOutputObjects()
     TH2F* hMphi=new TH2F(Form("hMphi%s",centrname.Data()),Form("Mass vs #Delta#phi %s;#Delta#phi;M (GeV/c^{2})",centrname.Data()),96,0,TMath::Pi(),fNMassBins,fLowmasslimit,fUpmasslimit);
     fOutput->Add(hMphi);
 
-    TH1F* hEvPlane=new TH1F(Form("hEvPlane%s",centrname.Data()),Form("Event plane angle %s;#phi Ev Plane;Entries",centrname.Data()),200,0.,TMath::Pi());
-    fOutput->Add(hEvPlane);
+    TH1F* hEvPlaneneg=new TH1F(Form("hEvPlaneneg%s",centrname.Data()),Form("Event plane angle %s;#phi Ev Plane;Entries",centrname.Data()),200,0.,TMath::Pi());
+    fOutput->Add(hEvPlaneneg);
+       
+       TH1F* hEvPlanepos=new TH1F(Form("hEvPlanepos%s",centrname.Data()),Form("Event plane angle %s;#phi Ev Plane;Entries",centrname.Data()),200,0.,TMath::Pi());
+    fOutput->Add(hEvPlanepos);
 
     //TH1F* hEvPlaneCheck=new TH1F(Form("hEvPlaneCheck%s",centrname.Data()),Form("Event plane angle - Event plane angle per candidate %s;(#phi(Ev Plane) - #phi(Ev Plane Candidate))/#phi(EvPlane);Entries",centrname.Data()),200,-0.2,0.2);
     //fOutput->Add(hEvPlaneCheck);
@@ -382,16 +398,17 @@ void AliAnalysisTaskSEHFv2::UserCreateOutputObjects()
   TH1F* hPhiBins=new TH1F("hPhiBins","Bins in #Delta#phi used in this analysis;#phi bin;n jobs",fNPhiBinLims-1,fPhiBins);
   fOutput->Add(hPhiBins);
   for(Int_t k=0;k<fNPhiBinLims-1;k++)hPhiBins->SetBinContent(k+1,1);
-  
   PostData(1,fhEventsInfo);
   PostData(2,fOutput);
-  fParHist = new TList();
-  fParHist->SetOwner();
-  fParHist->SetName("VZEROcorr");
-  for(Int_t i=0;i<6;i++){
-    fParHist->Add((TH2D*)fHistvzero[i]);
+  if(fUseV0EP){
+    fParHist = new TList();
+    fParHist->SetOwner();
+    fParHist->SetName("VZEROcorr");
+    for(Int_t i=0;i<6;i++){
+      fParHist->Add((TH2D*)fHistvzero[i]);
+    }
+    PostData(4,fParHist);
   }
-  PostData(4,fParHist);
   return;
 }
 
@@ -400,7 +417,6 @@ void AliAnalysisTaskSEHFv2::UserExec(Option_t */*option*/)
 {
   // Execute analysis for current event:
   // heavy flavor candidates association to MC truth
-   
   AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
   if(fDebug>2) printf("Analysing decay %d\n",fDecChannel);
   // Post the data already here
@@ -499,11 +515,16 @@ void AliAnalysisTaskSEHFv2::UserExec(Option_t */*option*/)
       fhEventsInfo->Fill(4);
     return;
   }
-
   AliEventplane *pl=0x0;
   TVector2* q=0x0;
   Double_t rpangleevent=0;
+  Double_t rpangleeventneg=0;
+  Double_t rpangleeventpos=0;
   Double_t eventplane=0;
+  TVector2 *qsub1=0x0;
+  TVector2 *qsub2=0x0;
+  
   //determine centrality bin
   Float_t centr=fRDCuts->GetCentrality(aod);
   Int_t icentr=0;
@@ -520,11 +541,13 @@ void AliAnalysisTaskSEHFv2::UserExec(Option_t */*option*/)
     rpangleevent=g->Rndm()*TMath::Pi();
     delete g;g=0x0;
     eventplane=rpangleevent;
+       ((TH1F*)fOutput->FindObject(Form("hEvPlanepos%s",centrbinname.Data())))->Fill(rpangleevent);
     if(fUseAfterBurner)fAfterBurner->SetEventPlane((Double_t)rpangleevent);
   }else{
     if(fUseV0EP){
       rpangleevent=GetEventPlaneFromV0(aod);
       eventplane=rpangleevent;
+         ((TH1F*)fOutput->FindObject(Form("hEvPlanepos%s",centrbinname.Data())))->Fill(rpangleevent);
     }else{
       // event plane and resolution 
       //--------------------------------------------------------------------------
@@ -534,8 +557,20 @@ void AliAnalysisTaskSEHFv2::UserExec(Option_t */*option*/)
        AliError("AliAnalysisTaskSEHFv2::UserExec:no eventplane! v2 analysis without eventplane not possible!\n");
        return;
       }
-      q = pl->GetQVector();
-      rpangleevent = pl->GetEventplane("Q"); // reaction plane angle without autocorrelations removal
+      if(fEtaGap){
+        qsub1 = pl->GetQsub1();
+       qsub2 = pl->GetQsub2();
+       rpangleeventpos = qsub1->Phi()/2.;
+       rpangleeventneg = qsub2->Phi()/2.;
+       ((TH1F*)fOutput->FindObject(Form("hEvPlanepos%s",centrbinname.Data())))->Fill(rpangleeventpos);
+       ((TH1F*)fOutput->FindObject(Form("hEvPlaneneg%s",centrbinname.Data())))->Fill(rpangleeventneg);
+       }
+      else if(!fEtaGap){
+       q = pl->GetQVector();
+        rpangleevent = pl->GetEventplane("Q");
+       ((TH1F*)fOutput->FindObject(Form("hEvPlanepos%s",centrbinname.Data())))->Fill(rpangleevent); // reaction plane angle without autocorrelations removal
+        } 
+       
       Double_t deltaPsi = pl->GetQsubRes();
       if(TMath::Abs(deltaPsi)>TMath::Pi()/2.){
        if(deltaPsi>0.) deltaPsi-=TMath::Pi();
@@ -543,10 +578,10 @@ void AliAnalysisTaskSEHFv2::UserExec(Option_t */*option*/)
       } // difference of subevents reaction plane angle cannot be bigger than phi/2
       Double_t planereso = TMath::Cos(2.*deltaPsi); // reaction plane resolution
       //--------------------------------------------------------------------------
-      ((TH1F*)fOutput->FindObject(Form("hEvPlaneReso%s",centrbinname.Data())))->Fill(planereso);
+      ((TH1F*)fOutput->FindObject(Form("hEvPlaneReso%s",centrbinname.Data())))->Fill(planereso);    
     }
   }
-  ((TH1F*)fOutput->FindObject(Form("hEvPlane%s",centrbinname.Data())))->Fill(rpangleevent);
+  
 
   for (Int_t iCand = 0; iCand < nCand; iCand++) {
     
@@ -573,7 +608,7 @@ void AliAnalysisTaskSEHFv2::UserExec(Option_t */*option*/)
       }
 
       if(!fUseV0EP) {
-       eventplane = GetEventPlaneForCandidate(d,q,pl); // remove autocorrelations
+       eventplane = GetEventPlaneForCandidate(d,q,pl,qsub1,qsub2); // remove autocorrelations
        ((TH1F*)fOutput->FindObject(Form("hEvPlaneCand%s",centrbinname.Data())))->Fill(rpangleevent-eventplane);
        //((TH1F*)fOutput->FindObject(Form("hEvPlaneCheck%s",centrbinname.Data())))->Fill((rpangleevent-eventplane)/100.*rpangleevent);
       }
@@ -829,13 +864,32 @@ Float_t AliAnalysisTaskSEHFv2::GetPhi0Pi(Float_t phi){
 }
 
 //________________________________________________________________________
-Float_t AliAnalysisTaskSEHFv2::GetEventPlaneForCandidate(AliAODRecoDecayHF* d, TVector2* q,AliEventplane *pl){
+Float_t AliAnalysisTaskSEHFv2::GetEventPlaneForCandidate(AliAODRecoDecayHF* d, TVector2* q,AliEventplane *pl,TVector2* qsub1,TVector2* qsub2){
   // remove autocorrelations 
  
-  TArrayF* qx = pl->GetQContributionXArray();
-  TArrayF* qy = pl->GetQContributionYArray();
-  TVector2 qcopy = *q;
+  TArrayF* qx = 0x0;
+  TArrayF* qy = 0x0;
+  TVector2 qcopy; 
+  if(!fEtaGap){
+    qx = pl->GetQContributionXArray();
+    qy = pl->GetQContributionYArray();
+    qcopy = *q;
+    }
+  else {
+    if(d->Eta()>0.){
+      qx = pl->GetQContributionXArraysub1();
+      qy = pl->GetQContributionYArraysub1();
+      qcopy = *qsub1;
+    }
+    else{
+      qx = pl->GetQContributionXArraysub2();
+      qy = pl->GetQContributionYArraysub2();
+      qcopy = *qsub2;
+    }
+  }
+  
   
   if(fDecChannel==2){
     //D* -- Yifei, Alessandro,Robert
     AliAODRecoDecayHF2Prong* theD0particle = ((AliAODRecoCascadeHF*)d)->Get2Prong();
index bc4268139503c0c01749cba6ffb909308b82be33..e71600fe4db2d9ef68ba58c230b033ed35654cfa 100644 (file)
@@ -21,6 +21,7 @@
 #include "AliAnalysisVertexingHF.h"
 #include "AliHFAfterBurner.h"
 
+
 class TH1F;
 class TH2D;
 class AliMultiDimVector;
@@ -35,7 +36,7 @@ class AliAnalysisTaskSEHFv2 : public AliAnalysisTaskSE
   enum DecChannel{kDplustoKpipi,kD0toKpi,kDstartoKpipi}; //more particles can be added
 
   AliAnalysisTaskSEHFv2();
-  AliAnalysisTaskSEHFv2(const char *name, AliRDHFCuts *rdCuts, Int_t decaychannel,Int_t nbinsphi, Float_t *phibinlimits,TH2D** histPar);
+  AliAnalysisTaskSEHFv2(const char *name, AliRDHFCuts *rdCuts, Int_t decaychannel,Int_t nbinsphi, Float_t *phibinlimits);
  
   virtual ~AliAnalysisTaskSEHFv2();
 
@@ -45,12 +46,14 @@ class AliAnalysisTaskSEHFv2 : public AliAnalysisTaskSE
   void SetNMassBins(Int_t nbins){fNMassBins=nbins;}
   void SetUpperCentLimit(Float_t lim){fCentUpLimit = lim;}
   void SetLowerCentLimit(Float_t lim){fCentLowLimit = lim;}
-  void SetUseV0EP(Bool_t flagV0EP){fUseV0EP=flagV0EP;}
+  //void SetUseV0EP(Bool_t flagV0EP){fUseV0EP=flagV0EP;}
   void SetV0EventPlaneOrder(Int_t n){fV0EPorder=n;}
   void SetMinCentrality(Int_t mincentr){fMinCentr=mincentr;}
   void SetMaxCentrality(Int_t maxcentr){fMaxCentr=maxcentr;}
   void SetUseAfterBurner(Bool_t ab){fUseAfterBurner=ab;}
   void SetAfterBurner(AliHFAfterBurner *ab){fAfterBurner=ab;}
+  void SetEtaGapFeatureForEventplaneFromTracks (Bool_t etaGap) {fEtaGap = etaGap;}
+  void SetVZEROParHist(TH2D** h);
 
   Float_t GetUpperMassLimit()const {return fUpmasslimit;}
   Float_t GetLowerMassLimit()const {return fLowmasslimit;}
@@ -77,7 +80,7 @@ class AliAnalysisTaskSEHFv2 : public AliAnalysisTaskSE
   void FillDplus(AliAODRecoDecayHF* d,TClonesArray *arrayMC,Int_t ptbin, Float_t dphi,Float_t* masses,Int_t isSel,Int_t icentr);
   void FillD02p(AliAODRecoDecayHF* d,TClonesArray *arrayMC,Int_t ptbin, Float_t dphi,Float_t* masses, Int_t isSel,Int_t icentr);
   void FillDstar(AliAODRecoDecayHF* d,TClonesArray *arrayMC,Int_t ptbin, Float_t dphi,Float_t* masses,Int_t isSel,Int_t icentr);
-  Float_t GetEventPlaneForCandidate(AliAODRecoDecayHF* d, TVector2* q,AliEventplane *pl);
+  Float_t GetEventPlaneForCandidate(AliAODRecoDecayHF* d, TVector2* q,AliEventplane *pl,TVector2* qsub1,TVector2* qsub2);
   Float_t GetEventPlaneFromV0(AliAODEvent *aodEvent);
 
 
@@ -102,6 +105,7 @@ class AliAnalysisTaskSEHFv2 : public AliAnalysisTaskSE
   Int_t  fV0EPorder;            //harmonic for VZERO event plane
   Int_t fMinCentr;              //minimum centrality
   Int_t fMaxCentr;              //maximum centrality
+  Bool_t fEtaGap;               // Eta gap feature for Eventplane from tracks; be careful that you do the correct settings in AddTaskEventPlane.C !!!!
 
   ClassDef(AliAnalysisTaskSEHFv2,1); // AliAnalysisTaskSE for the HF v2 analysis
 };