]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliVZEROEPSelectionTask.cxx
end-of-line normalization
[u/mrichter/AliRoot.git] / ANALYSIS / AliVZEROEPSelectionTask.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 //*****************************************************
17 //   Class AliVZEROEPSelectionTask
18 //   author: Cvetan Cheshkov
19 //   30/01/2012
20 //   This analysis task reads the OADB and
21 //   provides the parameters needed to flatten
22 //   the VZERO event plane in AliEventplane
23 //*****************************************************
24
25 #include "AliVZEROEPSelectionTask.h"
26
27 #include <TList.h>
28 #include <TProfile.h>
29 #include <TFile.h>
30 #include <TString.h>
31 #include <TDirectory.h>
32
33 #include "AliLog.h"
34 #include "AliVEvent.h"
35 #include "AliAnalysisManager.h"
36 #include "AliOADBContainer.h"
37 #include "AliEventplane.h"
38 #include "AliCentrality.h"
39
40 ClassImp(AliVZEROEPSelectionTask)
41
42 //________________________________________________________________________
43 AliVZEROEPSelectionTask::AliVZEROEPSelectionTask():
44 AliAnalysisTaskSE(),
45   fRunNumber(-1),
46   fUserParams(kFALSE),
47   fUseVZEROCentrality(kFALSE),
48   fVZEROEPContainer(0)
49 {   
50   // Default constructor
51   // Initialize pointers
52   AliInfo("VZERO Event Plane Selection enabled.");
53   for(Int_t i = 0; i < 11; ++i) fX2In[i] = fY2In[i] = fX2Y2In[i] = fCos8PsiIn[i] = NULL;
54 }   
55
56 //________________________________________________________________________
57 AliVZEROEPSelectionTask::AliVZEROEPSelectionTask(const char *name):
58   AliAnalysisTaskSE(name),
59   fRunNumber(-1),
60   fUserParams(kFALSE),
61   fUseVZEROCentrality(kFALSE),
62   fVZEROEPContainer(0)
63 {
64   // Default constructor
65   // Initialize pointers
66   AliInfo("Event Plane Selection enabled.");
67   for(Int_t i = 0; i < 11; ++i) fX2In[i] = fY2In[i] = fX2Y2In[i] = fCos8PsiIn[i] = NULL;
68 }
69  
70 //________________________________________________________________________
71 AliVZEROEPSelectionTask::~AliVZEROEPSelectionTask()
72 {
73   // Destructor
74   // ...
75   if (fUserParams) {
76     for(Int_t i = 0; i < 11; ++i) {
77       delete fX2In[i];
78       fX2In[i] = NULL;
79       delete fY2In[i];
80       fY2In[i] = NULL;
81       delete fX2Y2In[i];
82       fX2Y2In[i] = NULL;
83       delete fCos8PsiIn[i];
84       fCos8PsiIn[i] = NULL;
85     }
86   }
87   if (fVZEROEPContainer){
88     delete fVZEROEPContainer;
89     fVZEROEPContainer = NULL;
90   }
91 }  
92
93 //________________________________________________________________________
94 void AliVZEROEPSelectionTask::UserCreateOutputObjects()
95 {  
96   // Create the output containers (none in this case)
97   // Open the OADB file
98   
99   if(!fUserParams) {
100     TString oadbFileName = Form("%s/COMMON/EVENTPLANE/data/vzero.root", AliAnalysisManager::GetOADBPath());
101     TFile *fOADB = TFile::Open(oadbFileName); 
102     if(!fOADB->IsOpen()) AliFatal(Form("Cannot open OADB file %s", oadbFileName.Data()));
103
104     AliInfo("Using Standard OADB");
105     AliOADBContainer *cont = (AliOADBContainer*)fOADB->Get("vzeroEP");
106     if (!cont) AliFatal("Cannot fetch OADB container for VZERO EP selection");
107     fVZEROEPContainer = new AliOADBContainer(*cont);
108     fOADB->Close();
109     delete fOADB;
110   }
111 }
112
113 //________________________________________________________________________
114 void AliVZEROEPSelectionTask::UserExec(Option_t */*option*/)
115
116   // Execute analysis for current event:
117   // Fill the flatenning parameters in
118   // AliEventplane object
119
120   AliVEvent* event = InputEvent();
121   if (!(fRunNumber == event->GetRunNumber())) {
122     fRunNumber = event->GetRunNumber();
123     SetParamsFromOADB();
124   }
125
126   AliCentrality *centrality = event->GetCentrality();
127   Float_t percentile = (fUseVZEROCentrality) ? centrality->GetCentralityPercentile("V0M") : centrality->GetCentralityPercentile("CL1");
128   AliEventplane *esdEP = event->GetEventplane();
129   if(esdEP) SetEventplaneParams(esdEP,percentile);
130 }
131
132 //________________________________________________________________________
133 void AliVZEROEPSelectionTask::Terminate(Option_t */*option*/)
134 {
135   // Terminate analysis
136   // Nothing here
137 }
138
139 //________________________________________________________________________
140 void AliVZEROEPSelectionTask::SetEventplaneParams(AliEventplane *esdEP,Float_t percentile)
141 {
142   // Read the OADB histograms and
143   // prepare parameters used in order to
144   // flatten the event-plane
145   if(!esdEP)
146     AliFatal("No event plane received");
147
148   if (percentile < 0 || percentile > 100) {
149     for(Int_t ring = 0; ring < 11; ++ring) esdEP->SetVZEROEPParams(ring,0.,0.,1.,1.,0.,0.,0.);
150     return;
151   }
152
153   for(Int_t ring = 0; ring < 11; ++ring) {
154     Int_t ibin = fX2In[ring]->FindBin(percentile);
155     if (fX2In[ring]->GetBinEntries(ibin) == 0) {
156       esdEP->SetVZEROEPParams(ring,0.,0.,1.,1.,0.,0.,0.);
157       continue;
158     }
159     Double_t meanX2 = fX2In[ring]->GetBinContent(ibin);
160     Double_t meanY2 = fY2In[ring]->GetBinContent(ibin);
161     Double_t sigmaX2 = fX2In[ring]->GetBinError(ibin);
162     Double_t sigmaY2 = fY2In[ring]->GetBinError(ibin);
163     Double_t rho = (fX2Y2In[ring]->GetBinContent(ibin)-meanX2*meanY2)/sigmaX2/sigmaY2;
164   
165     Double_t b = rho*sigmaX2*sigmaY2*
166       TMath::Sqrt(2.*(sigmaX2*sigmaX2+sigmaY2*sigmaY2-2.*sigmaX2*sigmaY2*TMath::Sqrt(1.-rho*rho))/
167                   ((sigmaX2*sigmaX2-sigmaY2*sigmaY2)*(sigmaX2*sigmaX2-sigmaY2*sigmaY2)+
168                    4.*sigmaX2*sigmaX2*sigmaY2*sigmaY2*rho*rho));
169     Double_t aPlus = TMath::Sqrt(2.*sigmaX2*sigmaX2-b*b);
170     Double_t aMinus= TMath::Sqrt(2.*sigmaY2*sigmaY2-b*b);
171
172     Double_t lambdaPlus = b/aPlus;
173     Double_t lambdaMinus = b/aMinus;
174
175     Double_t cos8Psi = fCos8PsiIn[ring]->GetBinContent(ibin);
176     esdEP->SetVZEROEPParams(ring,meanX2,meanY2,aPlus,aMinus,lambdaPlus,lambdaMinus,cos8Psi);
177   }
178 }
179
180 //__________________________________________________________________________
181 void AliVZEROEPSelectionTask::SetParamsFromOADB() 
182 {
183   if(!fUserParams) {
184     TList *list = (TList*)fVZEROEPContainer->GetObject(fRunNumber, "Default");
185     if (!list) AliFatal(Form("Cannot find VZERO OADB list for run %d", fRunNumber));
186     SetHistograms(list);
187   }
188   else
189     AliInfo("Using custom VZERO event-plane params");
190 }
191
192 //__________________________________________________________________________
193 void AliVZEROEPSelectionTask::SetUserParams(const char* inFileName, const char* listName)
194 {
195   
196   fUserParams = kTRUE;
197   
198   TFile f(inFileName);
199   TList* list = (TList*)f.Get(listName);
200   if (!list) AliFatal(Form("Cannot find list %s in file %s", listName, inFileName));
201   SetHistograms(list);
202   f.Close();
203
204
205 //__________________________________________________________________________
206 void AliVZEROEPSelectionTask::SetHistograms(TList *list)
207 {
208   // Set the flatenning parameters
209   // histograms from a given list
210
211   for(Int_t i = 0; i < 11; ++i) {
212     if (fX2In[i]) delete fX2In[i];
213     fX2In[i] = (TProfile*)list->FindObject(Form("fX2_%d",i))->Clone(Form("fX2In_%d",i));
214     fX2In[i]->SetDirectory(0);
215     if (fY2In[i]) delete fY2In[i];
216     fY2In[i] = (TProfile*)list->FindObject(Form("fY2_%d",i))->Clone(Form("fY2In_%d",i));
217     fY2In[i]->SetDirectory(0);
218     if (fX2Y2In[i]) delete fX2Y2In[i];
219     fX2Y2In[i] = (TProfile*)list->FindObject(Form("fX2Y2_%d",i))->Clone(Form("fX2Y2In_%d",i));
220     fX2Y2In[i]->SetDirectory(0);
221     if (fCos8PsiIn[i]) delete fCos8PsiIn[i];
222     fCos8PsiIn[i] = (TProfile*)list->FindObject(Form("fCos8Psi_%d",i))->Clone(Form("fCos8PsiIn_%d",i));
223     fCos8PsiIn[i]->SetDirectory(0);
224   }
225 }