]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
q vector calibration for generated data
authorlmilano <lmilano@cern.ch>
Wed, 1 Oct 2014 11:21:00 +0000 (13:21 +0200)
committerlmilano <lmilano@cern.ch>
Wed, 1 Oct 2014 11:21:26 +0000 (13:21 +0200)
PWGLF/SPECTRA/PiKaPr/TestAOD/AliAnalysisTaskSpectraAllChAOD.cxx
PWGLF/SPECTRA/PiKaPr/TestAOD/AliSpectraAODEventCuts.cxx
PWGLF/SPECTRA/PiKaPr/TestAOD/AliSpectraAODEventCuts.h

index be60e0032a15dcfa648028dcbd6bd217bab20b8a..3da5dfb65180e497e2c5e67bb8f1ac628eca64b8 100644 (file)
@@ -206,6 +206,7 @@ void AliAnalysisTaskSpectraAllChAOD::UserExec(Option_t *)
     if(fIsQvecCalibMode){
       QvecMC = fEventCuts->CalculateQVectorMC(fVZEROside);
     }
+  else QvecMC = fEventCuts->GetQvecPercentileMC(fVZEROside);
   }
   
   Double_t Cent=fEventCuts->GetCent();
index 4aba5cc3269727f4fb006cc95d79ed417a662624..58ec908d13b6d76f6e5f4d74304caedd8e62b7ad 100644 (file)
@@ -88,6 +88,9 @@ AliSpectraAODEventCuts::AliSpectraAODEventCuts(const char *name) :
   fQvecIntegral(0), 
   fSplineArrayV0A(0),
   fSplineArrayV0C(0),
+  fQgenIntegral(0), 
+  fSplineArrayV0Agen(0),
+  fSplineArrayV0Cgen(0),
   fNch(0),
   fQvecCalibType(0)
 {
@@ -123,6 +126,10 @@ AliSpectraAODEventCuts::AliSpectraAODEventCuts(const char *name) :
   fSplineArrayV0A->SetOwner();
   fSplineArrayV0C = new TObjArray();
   fSplineArrayV0C->SetOwner();
+  fSplineArrayV0Agen = new TObjArray();
+  fSplineArrayV0Agen->SetOwner();
+  fSplineArrayV0Cgen = new TObjArray();
+  fSplineArrayV0Cgen->SetOwner();
   
   fOutput->Add(fHistoCuts);
   fOutput->Add(fHistoVtxBefSel);
@@ -728,6 +735,58 @@ Double_t AliSpectraAODEventCuts::CalculateQVectorMC(Int_t v0side){
   
 }
 
+//______________________________________________________
+Double_t AliSpectraAODEventCuts::GetQvecPercentileMC(Int_t v0side){
+  //check Qvec and Centrality consistency
+  if(fCent>90.) return -999.;
+  
+  Double_t qvec = CalculateQVectorMC(v0side);
+  if(qvec==-999.) return -999.;
+  
+  fQgenIntegral = 0x0;
+
+  if(v0side==0/*V0A*/){ fQgenIntegral = (TH2D*)fQvecIntList->FindObject("VZEROAgen"); }
+  if(v0side==1/*V0C*/){ fQgenIntegral = (TH2D*)fQvecIntList->FindObject("VZEROCgen"); }
+  //FIXME you need a check on the TH2D, add AliFatal or a return.
+
+  Int_t ic = -999;
+  
+  if(fQvecCalibType==1){
+    if(fNch<0.) return -999.;
+    ic = GetNchBin();
+  } else ic = (Int_t)fCent; //fQvecIntegral: 1% centrality bin
+  
+  TH1D *h1D = (TH1D*)fQgenIntegral->ProjectionY("h1Dgen",ic+1,ic+1);
+  
+  TSpline *spline = 0x0;
+  
+  if(v0side==0/*V0A*/){
+    if( CheckSplineArray(fSplineArrayV0Agen, ic) ) {
+      spline = (TSpline*)fSplineArrayV0Agen->At(ic);
+      //cout<<"Qvec V0A - ic: "<<ic<<" - Found TSpline..."<<endl;
+    } else {
+      spline = new TSpline3(h1D,"sp3");
+      fSplineArrayV0Agen->AddAtAndExpand(spline,ic);
+    }
+  }
+  else if(v0side==1/*V0C*/){
+    if( CheckSplineArray(fSplineArrayV0Cgen, ic) ) {
+      spline = (TSpline*)fSplineArrayV0Cgen->At(ic);
+    } else {
+      spline = new TSpline3(h1D,"sp3");
+      fSplineArrayV0Cgen->AddAtAndExpand(spline,ic);
+    }
+  }
+  
+  Double_t percentile = 100*spline->Eval(qvec);
+  
+  if(percentile>100. || percentile<0.) return -999.;
+
+  return percentile;
+
+}
+
 //______________________________________________________
 Bool_t AliSpectraAODEventCuts::CheckSplineArray(TObjArray * splarr, Int_t n){
   //check TSpline array consistency
index 8710e33f7af1ab058cd87e285c0ba8354cf9e373..79094407a386cba09c9fbc8a1f9b1296c2936fc5 100644 (file)
@@ -71,6 +71,9 @@ class AliSpectraAODEventCuts : public TNamed
     fQvecIntegral(0), 
     fSplineArrayV0A(0),
     fSplineArrayV0C(0),
+    fQgenIntegral(0), 
+    fSplineArrayV0Agen(0),
+    fSplineArrayV0Cgen(0),
     fNch(0),
     fQvecCalibType(0)
       {
@@ -167,6 +170,7 @@ class AliSpectraAODEventCuts : public TNamed
   Int_t GetNchBin();
   
   Double_t CalculateQVectorMC(Int_t v0side);
+  Double_t GetQvecPercentileMC(Int_t v0side);
 
  private:
   
@@ -216,6 +220,9 @@ class AliSpectraAODEventCuts : public TNamed
   TH2D * fQvecIntegral;           // ! Integrated Qvec distribution
   TObjArray * fSplineArrayV0A;    // TSpline array for VZERO-A
   TObjArray * fSplineArrayV0C;    // TSpline array for VZERO-C
+  TH2D * fQgenIntegral;           // ! Integrated Qvec distribution for generated tracks
+  TObjArray * fSplineArrayV0Agen;    // TSpline array for VZERO-A for generated tracks
+  TObjArray * fSplineArrayV0Cgen;    // TSpline array for VZERO-C for generated tracks
   
   Int_t fNch;
   Int_t fQvecCalibType; //0. centrality - 1. Nch
@@ -223,7 +230,7 @@ class AliSpectraAODEventCuts : public TNamed
   AliSpectraAODEventCuts(const AliSpectraAODEventCuts&);
   AliSpectraAODEventCuts& operator=(const AliSpectraAODEventCuts&);
   
-  ClassDef(AliSpectraAODEventCuts, 7);
+  ClassDef(AliSpectraAODEventCuts, 8);
   
 };
 #endif