1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // Here some comments on what this task does //
25 #include <TInterpreter.h>
32 #include "AliFMDAnalysisTaskCollector.h"
33 #include "AliAnalysisManager.h"
34 #include "AliESDFMD.h"
35 #include "AliESDEvent.h"
36 #include "AliAODEvent.h"
37 #include "AliAODHandler.h"
38 #include "AliMCEventHandler.h"
40 #include "AliMultiplicity.h"
41 #include "AliESDVertex.h"
42 #include "AliFMDAnaParameters.h"
43 //#include "AliFMDGeometry.h"
46 ClassImp(AliFMDAnalysisTaskCollector)
48 //____________________________________________________________________
49 Double_t AliFMDAnalysisTaskCollector::TripleLandau(Double_t *x, Double_t *par) {
51 Double_t energy = x[0];
52 Double_t constant = par[0];
53 Double_t mpv = par[1];
54 Double_t sigma = par[2];
55 Double_t alpha = par[3];
56 Double_t beta = par[4];
58 Double_t f = constant*(TMath::Landau(energy,mpv,sigma,kTRUE)+
59 alpha*TMath::Landau(energy,2*mpv+2*sigma*TMath::Log(2),2*sigma,kTRUE)+
60 beta*TMath::Landau(energy,3*mpv+3*sigma*TMath::Log(3),3*sigma,kTRUE) );
64 //____________________________________________________________________
66 AliFMDAnalysisTaskCollector::AliFMDAnalysisTaskCollector()
74 // Default constructor
78 //____________________________________________________________________
79 AliFMDAnalysisTaskCollector::AliFMDAnalysisTaskCollector(const char* name):
80 AliAnalysisTaskSE(name),
88 // Default constructor
90 DefineOutput(1, TList::Class());
92 //____________________________________________________________________
93 void AliFMDAnalysisTaskCollector::UserCreateOutputObjects()
95 // Create the output container
96 printf("AnalysisTaskFMD::CreateOutPutData() \n");
98 fOutputList = new TList();//(TList*)GetOutputData(0);
99 AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
101 fArray = new TObjArray();
102 fArray->SetName("FMD");
105 TH1F* hEmptyEdist = 0;
106 TH1F* hRingEdist = 0;
107 for(Int_t nEta = 0; nEta < pars->GetNetaBins(); nEta++) {
108 TObjArray* etaArray = new TObjArray();
109 fArray->AddAtAndExpand(etaArray,nEta);
110 for(Int_t det =1; det<=3;det++)
112 TObjArray* detArray = new TObjArray();
113 detArray->SetName(Form("FMD%d",det));
114 etaArray->AddAtAndExpand(detArray,det);
115 Int_t nRings = (det==1 ? 1 : 2);
116 for(Int_t ring = 0;ring<nRings;ring++)
118 Char_t ringChar = (ring == 0 ? 'I' : 'O');
119 hEdist = new TH1F(Form("FMD%d%c_etabin%d",det,ringChar,nEta),Form("FMD%d%c_etabin%d",det,ringChar,nEta),200,0,6);
120 hEdist->SetXTitle("#Delta E / E_{MIP}");
121 fOutputList->Add(hEdist);
122 detArray->AddAtAndExpand(hEdist,ring);
129 for(Int_t det =1; det<=3;det++)
131 Int_t nRings = (det==1 ? 1 : 2);
132 for(Int_t ring = 0;ring<nRings;ring++)
134 Char_t ringChar = (ring == 0 ? 'I' : 'O');
135 hRingEdist = new TH1F(Form("ringFMD%d%c",det,ringChar),Form("ringFMD%d%c",det,ringChar),200,0,6);
136 hRingEdist->SetXTitle("#Delta E / E_{MIP}");
137 fOutputList->Add(hRingEdist);
138 hEmptyEdist = new TH1F(Form("emptyFMD%d%c",det,ringChar),Form("emptyFMD%d%c",det,ringChar),200,0,6);
139 hEmptyEdist->SetXTitle("#Delta E / E_{MIP}");
140 fOutputList->Add(hEmptyEdist);
144 fZvtxDist = new TH1F("ZvtxDist","Vertex distribution",100,-30,30);
145 fZvtxDist->SetXTitle("z vertex");
146 fOutputList->Add(fZvtxDist);
149 //____________________________________________________________________
150 void AliFMDAnalysisTaskCollector::UserExec(Option_t */*option*/)
153 AliESDEvent* esd = (AliESDEvent*)InputEvent();
154 AliESD* old = esd->GetAliESDOld();
156 esd->CopyFromOldESD();
158 AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
159 TString triggers = esd->GetFiredTriggerClasses();
160 //if(!triggers.IsNull()) return;
161 //Bool_t trigger = pars->IsEventTriggered(esd);
163 Bool_t physics = pars->IsEventTriggered(esd);
164 Bool_t empty = pars->IsEventTriggered(esd,AliFMDAnaParameters::kEMPTY);
165 //std::cout<<physics<<" "<<empty<<std::endl;
170 Bool_t vtxStatus = pars->GetVertex(esd,vertex);
175 fZvtxDist->Fill(vertex[2]);
176 if(TMath::Abs(vertex[2]) > pars->GetVtxCutZ())
179 std::cout<<"Bananer "<<vtxStatus<<" "<<physics<<std::endl;
180 AliESDFMD* fmd = esd->GetFMDData();
187 if(!physics && !empty)
190 TH1F* emptyEdist = 0;
192 for(UShort_t det=1;det<=3;det++) {
194 Int_t nRings = (det==1 ? 1 : 2);
195 for (UShort_t ir = 0; ir < nRings; ir++) {
196 Char_t ring = (ir == 0 ? 'I' : 'O');
197 emptyEdist = (TH1F*)fOutputList->FindObject(Form("emptyFMD%d%c",det,ring));
198 ringEdist = (TH1F*)fOutputList->FindObject(Form("ringFMD%d%c",det,ring));
199 UShort_t nsec = (ir == 0 ? 20 : 40);
200 UShort_t nstr = (ir == 0 ? 512 : 256);
202 for(UShort_t sec =0; sec < nsec; sec++) {
203 for(UShort_t strip = 0; strip < nstr; strip++) {
206 Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
207 if(mult == AliESDFMD::kInvalidMult || mult == 0) continue;
209 Float_t eta = pars->GetEtaFromStrip(det,ring,sec,strip,vertex[2]);
211 Int_t nEta = pars->GetEtaBin(eta);
212 // std::cout<<det<<" "<<ring<<" "<<sec<<" "<<strip<<" "<<vertex[2]<<" "<<eta<<" "<<nEta<<std::endl;
214 Edist = (TH1F*)fOutputList->FindObject(Form("FMD%d%c_etabin%d",det,ring,nEta));
216 ringEdist->Fill(mult);
219 emptyEdist->Fill(mult);
223 AliWarning("Something is wrong - wrong trigger");
232 PostData(1, fOutputList);
236 //____________________________________________________________________
238 void AliFMDAnalysisTaskCollector::Terminate(Option_t */*option*/) {
239 std::cout<<"Analysed "<<fEvents<<" events and "<<fEmptyEvents<<" empty"<<std::endl;
240 AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
243 for(Int_t det = 1; det<=3; det++) {
244 Int_t nRings = (det == 1 ? 1 : 2);
245 for(Int_t ring = 0;ring<nRings; ring++) {
247 Char_t ringChar = (ring == 0 ? 'I' : 'O');
248 TH1F* hRingEdist = (TH1F*)fOutputList->FindObject(Form("ringFMD%d%c",det,ringChar));
249 TH1F* hEmptyEdist = (TH1F*)fOutputList->FindObject(Form("emptyFMD%d%c",det,ringChar));
251 hEmptyEdist->Scale(1./(Float_t)fEmptyEvents);
254 hRingEdist->Scale(1./(Float_t)fEvents);
255 for(Int_t nEta = 0; nEta < pars->GetNetaBins(); nEta++) {
256 TH1F* hEdist = (TH1F*)fOutputList->FindObject(Form("FMD%d%c_etabin%d",det,ringChar,nEta));
258 hEdist->Scale(1./(Float_t)fEvents);
266 //____________________________________________________________________
267 void AliFMDAnalysisTaskCollector::ReadFromFile(const Char_t* filename, Bool_t store, Int_t speciesOption) {
271 //1: Pb+Pb triple landau convolution fit
273 AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
274 //pars->Init(kTRUE,AliFMDAnaParameters::kBackgroundCorrection);
276 TFile fin(filename,"UPDATE");
278 TList* list = (TList*)fin.Get("energyDist");
280 AliFMDAnaCalibEnergyDistribution* EnergyDist = new AliFMDAnaCalibEnergyDistribution();
282 EnergyDist->SetNetaBins(pars->GetNetaBins());
283 EnergyDist->SetEtaLimits(pars->GetEtaMin(),pars->GetEtaMax());
287 for(Int_t det = 1; det<=3; det++) {
288 Int_t nRings = (det == 1 ? 1 : 2);
289 for(Int_t ring = 0;ring<nRings; ring++) {
290 Char_t ringChar = (ring == 0 ? 'I' : 'O');
291 TH1F* hEmptyEdist = (TH1F*)list->FindObject(Form("emptyFMD%d%c",det,ringChar));
292 TH1F* hRingEdist = (TH1F*)list->FindObject(Form("ringFMD%d%c",det,ringChar));
293 fitFunc = FitEnergyDistribution(hEmptyEdist, speciesOption) ;
295 fitFunc->Write(Form("emptyFMD%d%c_fitfunc",det,ringChar),TObject::kWriteDelete);
296 fitFunc = FitEnergyDistribution(hRingEdist, speciesOption) ;
298 fitFunc->Write(Form("FMD%d%c_fitfunc",det,ringChar),TObject::kWriteDelete);
301 EnergyDist->SetEmptyEnergyDistribution(det,ringChar,hEmptyEdist);
302 EnergyDist->SetRingEnergyDistribution(det,ringChar,hRingEdist);
303 for(Int_t nEta = 0; nEta < pars->GetNetaBins(); nEta++) {
304 TH1F* hEdist = (TH1F*)list->FindObject(Form("FMD%d%c_etabin%d",det,ringChar,nEta));
306 fitFunc = FitEnergyDistribution(hEdist, speciesOption) ;
308 fitFunc->Write(Form("FMD%d%c_etabin%d_fitfunc",det,ringChar,nEta),TObject::kWriteDelete);
309 EnergyDist->SetEnergyDistribution(det,ringChar,nEta,hEdist);
321 TFile fcalib(pars->GetPath(pars->GetEdistID() ),"RECREATE");
322 EnergyDist->Write(AliFMDAnaParameters::GetEdistID());
329 //____________________________________________________________________
330 TF1* AliFMDAnalysisTaskCollector::FitEnergyDistribution(TH1F* hEnergy, Int_t /*speciesOption*/) {
331 AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
333 if(hEnergy->GetEntries() != 0) {
335 hEnergy->GetXaxis()->SetRangeUser(0.2,hEnergy->GetXaxis()->GetXmax());
337 if(pars->GetCollisionSystem() == 0)
338 fitFunc = new TF1("FMDfitFunc","landau",hEnergy->GetBinCenter(hEnergy->GetMaximumBin())-0.2,3);
339 if(pars->GetCollisionSystem() == 1) {
340 fitFunc = new TF1("FMDfitFunc",TripleLandau,hEnergy->GetBinCenter(hEnergy->GetMaximumBin())-0.2,5,5);
341 fitFunc->SetParNames("constant","MPV","sigma","2-Mip weight","3-Mip weight");
342 fitFunc->SetParameters(10,0.8,0.1,0.05,0.01);
343 fitFunc->SetParLimits(1,0.6,1.2);
344 fitFunc->SetParLimits(3,0,1);
345 fitFunc->SetParLimits(4,0,1);
348 hEnergy->Fit(fitFunc,"","",hEnergy->GetBinCenter(hEnergy->GetMaximumBin())-0.2,3);
354 //____________________________________________________________________