]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Introduction of VZERO event-plane selection task that can be used in order to flatten...
authorcvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 31 Jan 2012 16:06:08 +0000 (16:06 +0000)
committercvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 31 Jan 2012 16:06:08 +0000 (16:06 +0000)
ANALYSIS/ANALYSISaliceLinkDef.h
ANALYSIS/AliVZEROEPSelectionTask.cxx [new file with mode: 0644]
ANALYSIS/AliVZEROEPSelectionTask.h [new file with mode: 0644]
ANALYSIS/CMakelibANALYSISalice.pkg
ANALYSIS/macros/AddTaskVZEROEPSelection.C [new file with mode: 0644]
PWGPP/VZERO/FillVZEROEPOADB.C [new file with mode: 0644]
STEER/STEERBase/AliEventplane.cxx
STEER/STEERBase/AliEventplane.h

index 683fbe9a61721efb0d693bc50a4eb24817b41f3b..3511a88c88bee573e485bdba41cab8040c3a9d80 100644 (file)
@@ -35,6 +35,7 @@
 #pragma link C++ class AliCollisionNormalizationTask+;
 #pragma link C++ class AliCentralitySelectionTask+;
 #pragma link C++ class AliEPSelectionTask+;
+#pragma link C++ class AliVZEROEPSelectionTask+;
 #pragma link C++ class AliAnalysisTaskStat+;
 #pragma link C++ class AliMultiInputEventHandler+;
 #pragma link C++ class AliAnalysisTaskPIDResponse+;
diff --git a/ANALYSIS/AliVZEROEPSelectionTask.cxx b/ANALYSIS/AliVZEROEPSelectionTask.cxx
new file mode 100644 (file)
index 0000000..8ee3a43
--- /dev/null
@@ -0,0 +1,221 @@
+/**************************************************************************
+ * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+//*****************************************************
+//   Class AliVZEROEPSelectionTask
+//   author: Cvetan Cheshkov
+//   30/01/2012
+//   This analysis task reads the OADB and
+//   provides the parameters needed to flatten
+//   the VZERO event plane in AliEventplane
+//*****************************************************
+
+#include "AliVZEROEPSelectionTask.h"
+
+#include <TList.h>
+#include <TProfile.h>
+#include <TFile.h>
+#include <TString.h>
+#include <TDirectory.h>
+
+#include "AliLog.h"
+#include "AliVEvent.h"
+#include "AliAnalysisManager.h"
+#include "AliOADBContainer.h"
+#include "AliEventplane.h"
+#include "AliCentrality.h"
+
+ClassImp(AliVZEROEPSelectionTask)
+
+//________________________________________________________________________
+AliVZEROEPSelectionTask::AliVZEROEPSelectionTask():
+AliAnalysisTaskSE(),
+  fRunNumber(-1),
+  fUserParams(kFALSE),
+  fUseVZEROCentrality(kFALSE),
+  fVZEROEPContainer(0)
+{   
+  // Default constructor
+  // Initialize pointers
+  AliInfo("VZERO Event Plane Selection enabled.");
+  for(Int_t i = 0; i < 8; ++i) fX2In[i] = fY2In[i] = fX2Y2In[i] = fCos8PsiIn[i] = NULL;
+}   
+
+//________________________________________________________________________
+AliVZEROEPSelectionTask::AliVZEROEPSelectionTask(const char *name):
+  AliAnalysisTaskSE(name),
+  fRunNumber(-1),
+  fUserParams(kFALSE),
+  fUseVZEROCentrality(kFALSE),
+  fVZEROEPContainer(0)
+{
+  // Default constructor
+  // Initialize pointers
+  AliInfo("Event Plane Selection enabled.");
+  for(Int_t i = 0; i < 8; ++i) fX2In[i] = fY2In[i] = fX2Y2In[i] = fCos8PsiIn[i] = NULL;
+}
+//________________________________________________________________________
+AliVZEROEPSelectionTask::~AliVZEROEPSelectionTask()
+{
+  // Destructor
+  // ...
+  if (fUserParams) {
+    for(Int_t i = 0; i < 8; ++i) {
+      delete fX2In[i];
+      fX2In[i] = NULL;
+      delete fY2In[i];
+      fY2In[i] = NULL;
+      delete fX2Y2In[i];
+      fX2Y2In[i] = NULL;
+      delete fCos8PsiIn[i];
+      fCos8PsiIn[i] = NULL;
+    }
+  }
+  if (fVZEROEPContainer){
+    delete fVZEROEPContainer;
+    fVZEROEPContainer = NULL;
+  }
+}  
+
+//________________________________________________________________________
+void AliVZEROEPSelectionTask::UserCreateOutputObjects()
+{  
+  // Create the output containers (none in this case)
+  // Open the OADB file
+  
+  if(!fUserParams) {
+    TString oadbFileName = Form("%s/COMMON/EVENTPLANE/data/vzero.root", AliAnalysisManager::GetOADBPath());
+    TFile *fOADB = TFile::Open(oadbFileName); 
+    if(!fOADB->IsOpen()) AliFatal(Form("Cannot open OADB file %s", oadbFileName.Data()));
+
+    AliInfo("Using Standard OADB");
+    AliOADBContainer *cont = (AliOADBContainer*)fOADB->Get("vzeroEP");
+    if (!cont) AliFatal("Cannot fetch OADB container for VZERO EP selection");
+    fVZEROEPContainer = new AliOADBContainer(*cont);
+    fOADB->Close();
+    delete fOADB;
+  }
+}
+
+//________________________________________________________________________
+void AliVZEROEPSelectionTask::UserExec(Option_t */*option*/)
+{ 
+  // Execute analysis for current event:
+  // Fill the flatenning parameters in
+  // AliEventplane object
+
+  AliVEvent* event = InputEvent();
+  if (!(fRunNumber == event->GetRunNumber())) {
+    fRunNumber = event->GetRunNumber();
+    SetParamsFromOADB();
+  }
+
+  AliCentrality *centrality = event->GetCentrality();
+  Float_t percentile = (fUseVZEROCentrality) ? centrality->GetCentralityPercentile("V0M") : centrality->GetCentralityPercentile("CL1");
+  AliEventplane *esdEP = event->GetEventplane();
+  SetEventplaneParams(esdEP,percentile);
+}
+
+//________________________________________________________________________
+void AliVZEROEPSelectionTask::Terminate(Option_t */*option*/)
+{
+  // Terminate analysis
+  // Nothing here
+}
+
+//________________________________________________________________________
+void AliVZEROEPSelectionTask::SetEventplaneParams(AliEventplane *esdEP,Float_t percentile)
+{
+  // Read the OADB histograms and
+  // prepare parameters used in order to
+  // flatten the event-plane
+  if(!esdEP)
+    AliFatal("No event plane received");
+
+  if (percentile < 0 || percentile > 100) {
+    for(Int_t ring = 0; ring < 8; ++ring) esdEP->SetVZEROEPParams(ring,0.,0.,1.,1.,0.,0.,0.);
+    return;
+  }
+
+  for(Int_t ring = 0; ring < 8; ++ring) {
+    Int_t ibin = fX2In[ring]->FindBin(percentile);
+    if (fX2In[ring]->GetBinEntries(ibin) == 0) {
+      esdEP->SetVZEROEPParams(ring,0.,0.,1.,1.,0.,0.,0.);
+      continue;
+    }
+    Double_t meanX2 = fX2In[ring]->GetBinContent(ibin);
+    Double_t meanY2 = fY2In[ring]->GetBinContent(ibin);
+    Double_t sigmaX2 = fX2In[ring]->GetBinError(ibin);
+    Double_t sigmaY2 = fY2In[ring]->GetBinError(ibin);
+    Double_t rho = (fX2Y2In[ring]->GetBinContent(ibin)-meanX2*meanY2)/sigmaX2/sigmaY2;
+  
+    Double_t b = rho*sigmaX2*sigmaY2*
+      TMath::Sqrt(2.*(sigmaX2*sigmaX2+sigmaY2*sigmaY2-2.*sigmaX2*sigmaY2*TMath::Sqrt(1.-rho*rho))/
+                 ((sigmaX2*sigmaX2-sigmaY2*sigmaY2)*(sigmaX2*sigmaX2-sigmaY2*sigmaY2)+
+                  4.*sigmaX2*sigmaX2*sigmaY2*sigmaY2*rho*rho));
+    Double_t aPlus = TMath::Sqrt(2.*sigmaX2*sigmaX2-b*b);
+    Double_t aMinus= TMath::Sqrt(2.*sigmaY2*sigmaY2-b*b);
+
+    Double_t lambdaPlus = b/aPlus;
+    Double_t lambdaMinus = b/aMinus;
+
+    Double_t cos8Psi = fCos8PsiIn[ring]->GetBinContent(ibin);
+    esdEP->SetVZEROEPParams(ring,meanX2,meanY2,aPlus,aMinus,lambdaPlus,lambdaMinus,cos8Psi);
+  }
+}
+
+//__________________________________________________________________________
+void AliVZEROEPSelectionTask::SetParamsFromOADB() 
+{
+  if(!fUserParams) {
+    TList *list = (TList*)fVZEROEPContainer->GetObject(fRunNumber, "Default");
+    if (!list) AliFatal(Form("Cannot find VZERO OADB list for run %d", fRunNumber));
+    SetHistograms(list);
+  }
+  else
+    AliInfo("Using custom VZERO event-plane params");
+}
+
+//__________________________________________________________________________
+void AliVZEROEPSelectionTask::SetUserParams(const char* inFileName, const char* listName)
+{
+  
+  fUserParams = kTRUE;
+  
+  TFile f(inFileName);
+  TList* list = (TList*)f.Get(listName);
+  if (!list) AliFatal(Form("Cannot find list %s in file %s", listName, inFileName));
+  SetHistograms(list);
+  f.Close();
+} 
+
+//__________________________________________________________________________
+void AliVZEROEPSelectionTask::SetHistograms(TList *list)
+{
+  // Set the flatenning parameters
+  // histograms from a given list
+
+  for(Int_t i = 0; i < 8; ++i) {
+    fX2In[i] = (TProfile*)list->FindObject(Form("fX2_%d",i))->Clone(Form("fX2In_%d",i));
+    fX2In[i]->SetDirectory(0);
+    fY2In[i] = (TProfile*)list->FindObject(Form("fY2_%d",i))->Clone(Form("fY2In_%d",i));
+    fY2In[i]->SetDirectory(0);
+    fX2Y2In[i] = (TProfile*)list->FindObject(Form("fX2Y2_%d",i))->Clone(Form("fX2Y2In_%d",i));
+    fX2Y2In[i]->SetDirectory(0);
+    fCos8PsiIn[i] = (TProfile*)list->FindObject(Form("fCos8Psi_%d",i))->Clone(Form("fCos8PsiIn_%d",i));
+    fCos8PsiIn[i]->SetDirectory(0);
+  }
+}
diff --git a/ANALYSIS/AliVZEROEPSelectionTask.h b/ANALYSIS/AliVZEROEPSelectionTask.h
new file mode 100644 (file)
index 0000000..2825e85
--- /dev/null
@@ -0,0 +1,60 @@
+#ifndef ALIVZEROEPSELECTIONTASK_H
+#define ALIVZEROEPSELECTIONTASK_H
+
+/* Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//*****************************************************
+//   Class AliVZEROEPSelectionTask
+//   author: Cvetan Cheshkov
+//   30/01/2012
+//   This analysis task reads the OADB and
+//   provides the parameters needed to flatten
+//   the VZERO event plane in AliEventplane
+//*****************************************************
+
+#include "AliAnalysisTaskSE.h"
+
+class TProfile;
+
+class AliOADBContainer;
+class AliEventplane;
+
+class AliVZEROEPSelectionTask : public AliAnalysisTaskSE {
+
+ public:
+  AliVZEROEPSelectionTask();
+  AliVZEROEPSelectionTask(const char *name);
+  virtual ~AliVZEROEPSelectionTask();
+
+  // Implementation of interface methods
+  virtual void UserCreateOutputObjects();
+  virtual void UserExec(Option_t *option);
+  virtual void Terminate(Option_t *option);
+  
+  void SetUserParams(const char* inFileName, const char* listName);
+  void UseVZEROCentrality() {fUseVZEROCentrality = kTRUE;}
+  
+ private:
+  void SetEventplaneParams(AliEventplane *esdEP,Float_t percentile);
+  void SetHistograms(TList *list);
+  void SetParamsFromOADB();
+   
+  AliVZEROEPSelectionTask(const AliVZEROEPSelectionTask& ep);
+  AliVZEROEPSelectionTask& operator= (const AliVZEROEPSelectionTask& ep); 
+
+  Int_t    fRunNumber;                 // runnumber
+  Bool_t   fUserParams;                        // in case one wants to use custom flatenning params
+  Bool_t   fUseVZEROCentrality;         // use VZERO centrality estimator instead of SPD
+  AliOADBContainer* fVZEROEPContainer; // VZERO event-plane OADB Container
+
+  TProfile *fX2In[8];                   // Profile histogram for Q^2_x (read from input file)
+  TProfile *fY2In[8];                   // Profile histogram for Q^2_y (read from input file)
+  TProfile *fX2Y2In[8];                 // Profile histogram for Q^2_x*Q^2_y (read from input file)
+  TProfile *fCos8PsiIn[8];              // Profile histogram for Cos(8*Psi) (read from input file)
+
+  ClassDef(AliVZEROEPSelectionTask,1) 
+};
+
+#endif
+
index 0a07085866c81f6745988bafec5ca335ca499dc8..fdc2703667b58db86a34252e84c25835fdc0bb0d 100644 (file)
@@ -50,6 +50,7 @@ set ( SRCS
     AliCollisionNormalizationTask.cxx 
     AliCentralitySelectionTask.cxx 
     AliEPSelectionTask.cxx 
+    AliVZEROEPSelectionTask.cxx 
     AliAnalysisTaskStat.cxx 
     AliMultiInputEventHandler.cxx 
     AliESDv0KineCuts.cxx 
diff --git a/ANALYSIS/macros/AddTaskVZEROEPSelection.C b/ANALYSIS/macros/AddTaskVZEROEPSelection.C
new file mode 100644 (file)
index 0000000..ac69e60
--- /dev/null
@@ -0,0 +1,22 @@
+AliVZEROEPSelectionTask *AddTaskVZEROEPSelection()
+{
+  AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+  if (!mgr) {
+    ::Error("AddTaskVZEROEPSelection", "No analysis manager to connect to.");
+    return NULL;
+  }    
+  // Check the analysis type using the event handlers connected to the analysis manager.
+  //==============================================================================
+  if (!mgr->GetInputEventHandler()) {
+    ::Error("AddTaskVZEROEPSelection", "This task requires an input event handler");
+    return NULL;
+  }
+
+  AliVZEROEPSelectionTask *task = new AliVZEROEPSelectionTask("AliVZEROEPSelectionTask");
+  //  task->UseVZEROCentrality(); // Optional line to swith from SPD to VZERO centrality estimator
+  mgr->AddTask(task);
+  
+  mgr->ConnectInput(task, 0, mgr->GetCommonInputContainer());
+
+  return task;
+}   
diff --git a/PWGPP/VZERO/FillVZEROEPOADB.C b/PWGPP/VZERO/FillVZEROEPOADB.C
new file mode 100644 (file)
index 0000000..433234b
--- /dev/null
@@ -0,0 +1,75 @@
+#if !defined(__CINT__) || defined(__MAKECINT__)
+#include "TFile.h"
+#include "TString.h"
+#include "TProfile.h"
+#include "TList.h"
+#include "TSystem.h"
+
+#include "AliOADBContainer.h"
+#include "AliAnalysisManager.h"
+#endif
+
+void FillVZEROEPOADB()
+{
+  gSystem->Load("libCore.so");  
+  gSystem->Load("libTree.so");
+  gSystem->Load("libGeom.so");
+  gSystem->Load("libVMC.so");
+  gSystem->Load("libPhysics.so");
+  gSystem->Load("libMinuit");
+  gSystem->Load("libSTEERBase");
+  gSystem->Load("libESD");
+  gSystem->Load("libAOD");
+  gSystem->Load("libANALYSIS");
+  gSystem->Load("libANALYSISalice");   
+  gSystem->Load("libOADB");
+
+  AliOADBContainer * oadbCont = new AliOADBContainer("vzeroEP");
+
+  TList *defaultList = new TList;
+  defaultList->SetName("Default");
+  TList *inputList = NULL;
+  TProfile *profHisto = NULL;
+  TFile fInputDefault("VZERO.EPFlatenning.PS.LHC11h_000170162_p1_muon_.root");
+  inputList = (TList*)fInputDefault.Get("coutput");
+  for(Int_t i = 0; i < 8; ++i) {
+    profHisto = (TProfile*)inputList->FindObject(Form("fX2_%d",i))->Clone(Form("fX2_%d",i));
+    profHisto->SetDirectory(0);
+    defaultList->Add(profHisto);
+    profHisto = (TProfile*)inputList->FindObject(Form("fY2_%d",i))->Clone(Form("fY2_%d",i));
+    profHisto->SetDirectory(0);
+    defaultList->Add(profHisto);
+    profHisto = (TProfile*)inputList->FindObject(Form("fX2Y2_%d",i))->Clone(Form("fX2Y2_%d",i));
+    profHisto->SetDirectory(0);
+    defaultList->Add(profHisto);
+    profHisto = (TProfile*)inputList->FindObject(Form("fCos8Psi_%d",i))->Clone(Form("fCos8Psi_%d",i));
+    profHisto->SetDirectory(0);
+    defaultList->Add(profHisto);
+  }
+  fInputDefault.Close();
+  oadbCont->AddDefaultObject(defaultList);
+
+
+  TList *list1 = new TList;
+  TFile fInput1("VZERO.EPFlatenning.PS.LHC11h_000169683_p1_muon_ESDs.root");
+  inputList = (TList*)fInput1.Get("coutput");
+  for(Int_t i = 0; i < 8; ++i) {
+    profHisto = (TProfile*)inputList->FindObject(Form("fX2_%d",i))->Clone(Form("fX2_%d",i));
+    profHisto->SetDirectory(0);
+    list1->Add(profHisto);
+    profHisto = (TProfile*)inputList->FindObject(Form("fY2_%d",i))->Clone(Form("fY2_%d",i));
+    profHisto->SetDirectory(0);
+    list1->Add(profHisto);
+    profHisto = (TProfile*)inputList->FindObject(Form("fX2Y2_%d",i))->Clone(Form("fX2Y2_%d",i));
+    profHisto->SetDirectory(0);
+    list1->Add(profHisto);
+    profHisto = (TProfile*)inputList->FindObject(Form("fCos8Psi_%d",i))->Clone(Form("fCos8Psi_%d",i));
+    profHisto->SetDirectory(0);
+    list1->Add(profHisto);
+  }
+  oadbCont->AppendObject(list1, 169683, 169683);
+
+
+  TString oadbFileName = Form("%s/COMMON/EVENTPLANE/data/vzero.root", AliAnalysisManager::GetOADBPath());
+  oadbCont->WriteToFile(oadbFileName.Data());
+}
index 7216873fb6b96e6721e0677a758136921ebbfa2c..236baef5a265047f24c231f2f382c1788a12a7d3 100644 (file)
@@ -50,6 +50,12 @@ AliEventplane::AliEventplane() : TNamed("Eventplane", "Eventplane"),
   fQContributionYsub1 = new TArrayF(0);
   fQContributionXsub2 = new TArrayF(0);
   fQContributionYsub2 = new TArrayF(0);
+  for(Int_t i = 0; i < 8; ++i) {
+    fMeanX2[i]=fMeanY2[i]=0.;
+    fAPlus[i]=fAMinus[i]=1.;
+    fLambdaPlus[i]=fLambdaMinus[i]=0.;
+    fCos8Psi[i]=0.;
+  }
 }
 
 AliEventplane::AliEventplane(const AliEventplane& ep) : 
@@ -105,6 +111,16 @@ void AliEventplane::CopyEP(AliEventplane& ep) const
       target.fQsub2 = dynamic_cast<TVector2*> (fQsub2->Clone());
   if (fQsubRes)
       target.fQsubRes = fQsubRes;
+
+  for(Int_t i = 0; i < 8; ++i) {
+    target.fMeanX2[i]=fMeanX2[i];
+    target.fMeanY2[i]=fMeanY2[i];
+    target.fAPlus[i]=fAPlus[i];
+    target.fAMinus[i]=fAMinus[i];
+    target.fLambdaPlus[i]=fLambdaPlus[i];
+    target.fLambdaMinus[i]=fLambdaMinus[i];
+    target.fCos8Psi[i]=fCos8Psi[i];
+  }
 }
 
 AliEventplane::~AliEventplane()
@@ -180,18 +196,79 @@ Double_t AliEventplane::CalculateVZEROEventPlane(const AliVEvent *  event, Int_t
     return -1000.;
   }
 
-  Double_t qx=0., qy=0.;
-  for(Int_t iCh = firstRing*8; iCh < (lastRing+1)*8; ++iCh) {
-    Double_t phi = TMath::Pi()/8. + (iCh%8) * TMath::Pi()/4.;
+  Double_t qxTot=0., qyTot=0.;
+  Double_t qx, qy;
+  Double_t totMult = 0.;
+  for(Int_t ring = firstRing; ring <=lastRing; ++ring) {
+    qx=qy=0.;
+    Double_t trans[2][2];
+    trans[0][0] = 1./(1.-fLambdaPlus[ring]*fLambdaMinus[ring]);
+    trans[0][1] = -fLambdaMinus[ring]/(1.-fLambdaPlus[ring]*fLambdaMinus[ring]);
+    trans[1][0] = -fLambdaPlus[ring]/(1.-fLambdaPlus[ring]*fLambdaMinus[ring]);
+    trans[1][1] = 1./(1.-fLambdaPlus[ring]*fLambdaMinus[ring]);
+    Double_t multRing = 0.;
+    for(Int_t iCh = ring*8; iCh < (ring+1)*8; ++iCh) {
+      Double_t phi = TMath::Pi()/8. + TMath::Pi()/4.*(iCh%8);
+      Double_t mult = event->GetVZEROEqMultiplicity(iCh);
+      qx += mult*TMath::Cos(harmonic*phi);
+      qy += mult*TMath::Sin(harmonic*phi);
+      multRing += mult;
+    }
+
+    if (multRing < 1e-6) continue;
+    totMult += multRing;
+
+    if (harmonic == 2) {
+      // Recentering
+      Double_t qxPrime = qx - fMeanX2[ring];
+      Double_t qyPrime = qy - fMeanY2[ring];
+      // Twist
+      Double_t qxSeconde = trans[0][0]*qxPrime + trans[0][1]*qyPrime;
+      Double_t qySeconde = trans[1][0]*qxPrime + trans[1][1]*qyPrime;
+      // Rescaling
+      Double_t qxTierce = qxSeconde/fAPlus[ring];
+      Double_t qyTierce = qySeconde/fAMinus[ring];
+      qxTot += qxTierce;
+      qyTot += qyTierce;
+    }
+    else {
+      qxTot += qx;
+      qyTot += qy;
+    }
+  }
 
-    Double_t mult = event->GetVZEROEqMultiplicity(iCh);
+  // In case of no hits return an invalid event-plane angle
+  if (totMult<1e-6) return -999;
 
-    qx += mult*TMath::Cos(harmonic*phi);
-    qy += mult*TMath::Sin(harmonic*phi);
+  // 4th Fourier momenta in order to flatten the EP within a sector
+  Double_t psi = TMath::ATan2(qyTot,qxTot)/harmonic;
+  if ((harmonic == 2) && (firstRing == lastRing)) {
+    psi += (fCos8Psi[firstRing]*TMath::Sin(2.*4.*psi))/2.;
+    if (psi > TMath::Pi()/2) psi -= TMath::Pi();
+    if (psi <-TMath::Pi()/2) psi += TMath::Pi();
   }
-  return (TMath::ATan2(qy,qx)/harmonic);
+
+  return psi;
 }
 
+void AliEventplane::SetVZEROEPParams(Int_t ring,
+                                    Double_t meanX2, Double_t meanY2,
+                                    Double_t aPlus, Double_t aMinus,
+                                    Double_t lambdaPlus, Double_t lambdaMinus,
+                                    Double_t cos8Psi)
+{
+  // Set the VZERO event-plane
+  // flatenning parameters
+  if (ring < 0 || ring >=8) AliFatal(Form("Bad ring index (%d)",ring));
+
+  fMeanX2[ring] = meanX2;
+  fMeanY2[ring] = meanY2;
+  fAPlus[ring] = aPlus;
+  fAMinus[ring] = aMinus;
+  fLambdaPlus[ring] = lambdaPlus;
+  fLambdaMinus[ring] = lambdaMinus;
+  fCos8Psi[ring] = cos8Psi;
+}
 
 TVector2* AliEventplane::GetQsub1()
 {
index d19017157ad392a12c2edf18a33af84402db7384..b0ffbfaa1263436e09670456ce35591459e00eb1 100644 (file)
@@ -53,6 +53,11 @@ class AliEventplane : public TNamed
   Double_t  GetQsubRes();
   Bool_t    IsEventInEventplaneClass(Double_t a, Double_t b, const char *method);
   Double_t  CalculateVZEROEventPlane(const AliVEvent *event, Int_t firstRing, Int_t lastRing, Int_t harmonic) const;
+  void      SetVZEROEPParams(Int_t ring,
+                            Double_t meanX2, Double_t meanY2,
+                            Double_t aPlus, Double_t aMinus,
+                            Double_t lambdaPlus, Double_t lambdaMinus,
+                            Double_t cos8Psi);
 
   void Reset();
 
@@ -68,7 +73,14 @@ class AliEventplane : public TNamed
    TVector2* fQsub1;            // Q-Vector of subevent 1
    TVector2* fQsub2;            // Q-Vector of subevent 2
    Double_t fQsubRes;           // Difference of EP angles of subevents
-  ClassDef(AliEventplane, 3)
+   Double_t fMeanX2[8];          // Mean Q^2_X for VZERO EP
+   Double_t fMeanY2[8];          // Mean Q^2_Y for VZERO EP
+   Double_t fAPlus[8];           // Q^2_X rescaling parameter for VZERO EP
+   Double_t fAMinus[8];          // Q^2_Y rescaling parameter for VZERO EP
+   Double_t fLambdaPlus[8];      // Q^2_X twisting parameter for VZERO EP
+   Double_t fLambdaMinus[8];     // Q^2_Y twisting parameter for VZERO EP 
+   Double_t fCos8Psi[8];         // 4th Fourier momenta used to flatten VZERO EP within a sector 
+
+  ClassDef(AliEventplane, 4)
 };
 #endif //ALIEVENTPLANE_H