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 This class provides the calibration of the time dependence of the TPC gain due to pressure and temperature changes etc.
21 //0. Libraries to load
22 gSystem->Load("libANALYSIS");
23 gSystem->Load("libSTAT");
24 gSystem->Load("libTPCcalib");
27 //1. Do calibration ...
32 //2. Visualize results
34 TFile fcalib("CalibObjects.root");
35 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
36 AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
37 TGraphErrors * gr = gain->GetGraphGainVsTime(0,1000)
39 // gain->GetHistGainTime()->GetAxis(1)->SetRangeUser(1213.8e6,1214.3e6)
40 TH2D * GainTime = gain->GetHistGainTime()->Projection(0,1)
41 GainTime->GetXaxis()->SetTimeDisplay(kTRUE)
42 GainTime->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}")
43 GainTime->Draw("colz")
46 // MakeSlineFit example
48 AliSplineFit *fit = AliTPCcalibTimeGain::MakeSplineFit(gr)
50 TGraph * grfit = fit.MakeGraph(gr->GetX()[0],gr->GetX()[gr->GetN()-1],50000,0);
52 gr->SetMarkerStyle(25);
54 grfit->SetLineColor(2);
58 // QA - dE/dx resoultion as a function of time
61 TGraph * grSigma = AliTPCcalibBase::FitSlicesSigma(gain->GetHistGainTime(),0,1,1800,5)
63 TCanvas *c1 = new TCanvas("c1","transparent pad",200,10,700,500);
64 TPad *pad1 = new TPad("pad1","",0,0,1,1);
65 TPad *pad2 = new TPad("pad2","",0,0,1,1);
66 pad2->SetFillStyle(4000); //will be transparent
70 GainTime->Draw("colz")
76 Double_t ymin = -0.04;
78 Double_t dy = (ymax-ymin)/0.8;
79 Double_t xmin = GainTime->GetXaxis()->GetXmin()
80 Double_t xmax = GainTime->GetXaxis()->GetXmax()
81 Double_t dx = (xmax-xmin)/0.8; //10 per cent margins left and right
82 pad2->Range(xmin-0.1*dx,ymin-0.1*dy,xmax+0.1*dx,ymax+0.1*dy);
85 grSigma->SetLineColor(2);
86 grSigma->SetLineWidth(2);
88 TGaxis *axis = new TGaxis(xmax,ymin,xmax,ymax,ymin,ymax,50510,"+L");
89 axis->SetLabelColor(kRed);
90 axis->SetTitle("dE/dx resolution #sigma_{dE/dx}");
93 ////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
95 ----> make Attachment study
97 TFile fcalib("CalibObjects40366b.root");
98 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
99 AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
100 TGraphErrors * grAttach = gain->GetGraphAttachment(2000,4)
102 TCanvas *c1 = new TCanvas("c1","transparent pad",200,10,700,500);
103 TPad *pad1 = new TPad("pad1","",0,0,1,1);
104 TPad *pad2 = new TPad("pad2","",0,0,1,1);
105 pad2->SetFillStyle(4000); //will be transparent
109 gain->GetHistGainTime()->GetAxis(1)->SetRangeUser(1213.8e6,1214.3e6)
110 TH2D * GainTime = gain->GetHistGainTime()->Projection(0,1)
111 GainTime->GetXaxis()->SetTimeDisplay(kTRUE)
112 GainTime->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}")
113 GainTime->Draw("colz")
117 Double_t ymin = -0.001;
118 Double_t ymax = 0.001;
119 Double_t dy = (ymax-ymin)/0.8;
120 Double_t xmin = GainTime->GetXaxis()->GetXmin()
121 Double_t xmax = GainTime->GetXaxis()->GetXmax()
122 Double_t dx = (xmax-xmin)/0.8; //10 per cent margins left and right
123 pad2->Range(xmin-0.1*dx,ymin-0.1*dy,xmax+0.1*dx,ymax+0.1*dy);
126 grAttach->SetLineColor(2);
127 grAttach->SetLineWidth(2);
129 TGaxis *axis = new TGaxis(xmax,ymin,xmax,ymax,ymin,ymax,50510,"+L");
130 axis->SetLabelColor(kRed);
131 axis->SetTitle("attachment coefficient b");
138 #include "Riostream.h"
144 #include "THnSparse.h"
145 #include "TGraphErrors.h"
151 #include "TVectorD.h"
152 #include "TProfile.h"
153 #include "TGraphErrors.h"
156 #include "AliTPCclusterMI.h"
157 #include "AliTPCseed.h"
158 #include "AliESDVertex.h"
159 #include "AliESDEvent.h"
160 #include "AliESDfriend.h"
161 #include "AliESDInputHandler.h"
162 #include "AliAnalysisManager.h"
164 #include "AliTracker.h"
166 #include "AliTPCCalROC.h"
167 #include "AliESDv0.h"
171 #include "AliTPCcalibTimeGain.h"
173 #include "TTreeStream.h"
174 #include "AliTPCTracklet.h"
175 #include "TTimeStamp.h"
176 #include "AliTPCcalibDB.h"
177 #include "AliTPCcalibLaser.h"
178 #include "AliDCSSensorArray.h"
179 #include "AliDCSSensor.h"
181 ClassImp(AliTPCcalibTimeGain)
184 AliTPCcalibTimeGain::AliTPCcalibTimeGain()
189 fIntegrationTimeDeDx(0),
197 fUseCookAnalytical(0),
199 fLowMemoryConsumption(0)
202 // Default constructor
204 AliDebug(5,"Default Constructor");
208 AliTPCcalibTimeGain::AliTPCcalibTimeGain(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeGain)
213 fIntegrationTimeDeDx(0),
221 fUseCookAnalytical(0),
223 fLowMemoryConsumption(0)
226 // No default constructor
231 AliDebug(5,"Non Default Constructor");
233 fIntegrationTimeDeDx = deltaIntegrationTimeGain;
235 Double_t deltaTime = EndTime - StartTime;
238 // main histogram for time dependence: dE/dx, time, type (1-muon cosmic,2-pion beam data, 3&4 - proton points at higher dE/dx), meanDriftlength, momenta (only filled if enough space is available), run number, eta
239 Int_t timeBins = TMath::Nint(deltaTime/deltaIntegrationTimeGain);
240 Int_t binsGainTime[7] = {150, timeBins, 4, 25, 200, 10000000, 20};
241 Double_t xminGainTime[7] = {0.5, StartTime, 0.5, 0, 0.1, -0.5, -1};
242 Double_t xmaxGainTime[7] = { 8, EndTime, 4.5, 250, 50, 9999999.5, 1};
243 fHistGainTime = new THnSparseF("HistGainTime","dEdx time dep.;dEdx,time,type,driftlength,momenta,run number, eta;dEdx",7,binsGainTime,xminGainTime,xmaxGainTime);
244 BinLogX(fHistGainTime, 4);
246 fHistDeDxTotal = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",250,0.01,100.,1000,0.,8);
247 BinLogX(fHistDeDxTotal);
249 // default values for dE/dx
254 fUseShapeNorm = kTRUE;
255 fUsePosNorm = kFALSE;
256 fUsePadNorm = kFALSE;
257 fUseCookAnalytical = kFALSE;
260 fLowMemoryConsumption = kTRUE;
267 AliTPCcalibTimeGain::~AliTPCcalibTimeGain(){
271 delete fHistGainTime;
273 delete fHistDeDxTotal;
277 void AliTPCcalibTimeGain::Process(AliESDEvent *event) {
282 Printf("ERROR: ESD not available");
285 AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
289 if (ESDfriend->TestSkipBit()) return;
291 if (fIsCosmic) { // this should be removed at some point based on trigger mask !?
292 ProcessCosmicEvent(event);
294 ProcessBeamEvent(event);
303 void AliTPCcalibTimeGain::ProcessCosmicEvent(AliESDEvent *event) {
305 // Process in case of cosmic event
307 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
309 Printf("ERROR: ESDfriend not available");
313 UInt_t time = event->GetTimeStamp();
314 Int_t nFriendTracks = esdFriend->GetNumberOfTracks();
315 Int_t runNumber = event->GetRunNumber();
319 for (Int_t i=0;i<nFriendTracks;++i) {
321 AliESDtrack *track = event->GetTrack(i);
322 if (!track) continue;
323 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
324 if (!friendTrack) continue;
325 const AliExternalTrackParam * trackIn = track->GetInnerParam();
326 const AliExternalTrackParam * trackOut = friendTrack->GetTPCOut();
327 if (!trackIn) continue;
328 if (!trackOut) continue;
330 // calculate necessary track parameters
331 Double_t meanP = trackIn->GetP();
332 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ());
333 Int_t nclsDeDx = track->GetTPCNcls();
335 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
336 if (nclsDeDx < 60) continue;
337 if (TMath::Abs(trackIn->GetTgl()) > 1) continue;
338 if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue;
341 TObject *calibObject;
342 AliTPCseed *seed = 0;
343 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
344 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
348 Double_t tpcSignal = GetTPCdEdx(seed);
349 if (nclsDeDx > 100) fHistDeDxTotal->Fill(meanP, tpcSignal);
351 if (fLowMemoryConsumption) {
352 if (meanP < 20) continue;
353 meanP = 30; // set all momenta to one in order to save memory
355 //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta
356 Double_t vec[6] = {tpcSignal,time,1,meanDrift,meanP,runNumber};
357 fHistGainTime->Fill(vec);
366 void AliTPCcalibTimeGain::ProcessBeamEvent(AliESDEvent *event) {
368 // Process in case of beam event
370 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
372 Printf("ERROR: ESDfriend not available");
376 UInt_t time = event->GetTimeStamp();
377 Int_t nFriendTracks = esdFriend->GetNumberOfTracks();
378 Int_t runNumber = event->GetRunNumber();
382 for (Int_t i=0;i<nFriendTracks;++i) { // begin track loop
384 AliESDtrack *track = event->GetTrack(i);
385 if (!track) continue;
386 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
387 if (!friendTrack) continue;
389 const AliExternalTrackParam * trackIn = track->GetInnerParam();
390 const AliExternalTrackParam * trackOut = friendTrack->GetTPCOut();
391 if (!trackIn) continue;
392 if (!trackOut) continue;
394 // calculate necessary track parameters
395 Double_t meanP = trackIn->GetP();
396 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ());
397 Int_t nclsDeDx = track->GetTPCNcls();
399 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
400 if (nclsDeDx < 60) continue;
401 //if (TMath::Abs(trackIn->GetTgl()) > 1) continue;
402 //if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue;
403 if (TMath::Abs(trackIn->Eta()) > 1) continue;
404 UInt_t status = track->GetStatus();
405 if ((status&AliESDtrack::kTPCrefit)==0) continue;
406 //if (track->GetNcls(0) < 3) continue; // ITS clusters
407 Float_t dca[2], cov[3];
408 track->GetImpactParameters(dca,cov);
409 if (TMath::Abs(dca[0]) > 7 || TMath::Abs(dca[0]) < 0.0000001 || TMath::Abs(dca[1]) > 25 ) continue; // cut in xy
410 Double_t eta = trackIn->Eta();
413 TObject *calibObject;
414 AliTPCseed *seed = 0;
415 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
416 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
420 Int_t particleCase = 0;
421 if (meanP < 0.5 && meanP > 0.4) particleCase = 2; // MIP pions
422 if (meanP < 0.57 && meanP > 0.56) particleCase = 3; // protons 1
423 if (meanP < 0.66 && meanP > 0.65) particleCase = 4; // protons 2
425 if (fLowMemoryConsumption && particleCase == 0) continue;
427 Double_t tpcSignal = GetTPCdEdx(seed);
428 fHistDeDxTotal->Fill(meanP, tpcSignal);
430 //dE/dx, time, type (1-muon cosmic,2-pion beam data, 3&4 protons), momenta, runNumner, eta
431 Double_t vec[7] = {tpcSignal,time,particleCase,meanDrift,meanP,runNumber, eta};
432 fHistGainTime->Fill(vec);
437 // V0 loop -- in beam events the cosmic part of the histogram is filled with GammaConversions
439 for(Int_t iv0 = 0; iv0 < event->GetNumberOfV0s(); iv0++) {
440 AliESDv0 * v0 = event->GetV0(iv0);
441 if (!v0->GetOnFlyStatus()) continue;
442 if (v0->GetEffMass(0,0) > 0.02) continue; // select low inv. mass
444 v0->GetXYZ(xyz[0], xyz[1], xyz[2]);
445 if (TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]) < 3 || TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]) > 30) continue;
447 // "loop over daughters"
449 for(Int_t idaughter = 0; idaughter < 2; idaughter++) { // daughter loop
450 Int_t index = idaughter == 0 ? v0->GetPindex() : v0->GetNindex();
451 AliESDtrack * trackP = event->GetTrack(index);
452 AliESDfriendTrack *friendTrackP = esdFriend->GetTrack(index);
453 if (!friendTrackP) continue;
454 const AliExternalTrackParam * trackPIn = trackP->GetInnerParam();
455 const AliExternalTrackParam * trackPOut = friendTrackP->GetTPCOut();
456 if (!trackPIn) continue;
457 if (!trackPOut) continue;
458 // calculate necessary track parameters
459 Double_t meanP = trackPIn->GetP();
460 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackPIn->GetZ() + trackPOut->GetZ());
461 Int_t nclsDeDx = trackP->GetTPCNcls();
462 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
463 if (nclsDeDx < 60) continue;
464 if (TMath::Abs(trackPIn->GetTgl()) > 1) continue;
466 TObject *calibObject;
467 AliTPCseed *seed = 0;
468 for (Int_t l=0;(calibObject=friendTrackP->GetCalibObject(l));++l) {
469 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
472 if (fLowMemoryConsumption) {
473 if (meanP > 0.5 || meanP < 0.4) continue;
474 meanP = 0.45; // set all momenta to one in order to save memory
476 Double_t tpcSignal = GetTPCdEdx(seed);
477 //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta
478 Double_t vec[6] = {tpcSignal,time,1,meanDrift,meanP,runNumber};
479 fHistGainTime->Fill(vec);
488 Float_t AliTPCcalibTimeGain::GetTPCdEdx(AliTPCseed * seed) {
490 // calculate tpc dEdx
494 if (!fUseCookAnalytical) {
495 signal = (1/fMIP)*seed->CookdEdxNorm(fLowerTrunc,fUpperTrunc,fUseMax,0,159,fUseShapeNorm,fUsePosNorm,fUsePadNorm,0);
497 signal = (1/fMIP)*seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,fUseMax);
504 void AliTPCcalibTimeGain::AnalyzeRun(Int_t minEntries) {
506 // Analyze results of calibration
509 fHistGainTime->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
510 fHistGainTime->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
512 fHistGainTime->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
513 fHistGainTime->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
516 fGainVsTime = AliTPCcalibBase::FitSlices(fHistGainTime,0,1,minEntries,10);
522 TGraphErrors * AliTPCcalibTimeGain::GetGraphGainVsTime(Int_t runNumber, Int_t minEntries) {
524 // Analyze results and get the graph
526 if (runNumber == 0) {
528 AnalyzeRun(minEntries);
531 // 1st check if the current run was cosmic or beam event
532 fHistGainTime->GetAxis(5)->SetRangeUser(runNumber,runNumber);
533 AnalyzeRun(minEntries);
535 if (fGainVsTime->GetN() == 0) return 0;
539 Long64_t AliTPCcalibTimeGain::Merge(TCollection *li) {
543 TIterator* iter = li->MakeIterator();
544 AliTPCcalibTimeGain* cal = 0;
546 while ((cal = (AliTPCcalibTimeGain*)iter->Next())) {
547 if (!cal->InheritsFrom(AliTPCcalibTimeGain::Class())) {
548 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
552 // add histograms here...
553 if (cal->GetHistGainTime()) fHistGainTime->Add(cal->GetHistGainTime());
554 if (cal->GetHistDeDxTotal()) fHistDeDxTotal->Add(cal->GetHistDeDxTotal());
563 AliSplineFit * AliTPCcalibTimeGain::MakeSplineFit(TGraphErrors * graph) {
565 // make spline fit of gain
567 AliSplineFit *fit = new AliSplineFit();
568 fit->SetGraph(graph);
569 fit->SetMinPoints(graph->GetN()+1);
570 fit->InitKnots(graph,2,0,0.001);
578 TGraphErrors * AliTPCcalibTimeGain::GetGraphAttachment(Int_t minEntries, Int_t nmaxBin, Float_t /*fracLow*/, Float_t /*fracUp*/) {
580 // For each time bin the driftlength-dependence of the signal is fitted.
582 TH3D * hist = fHistGainTime->Projection(1, 0, 3);
583 Double_t *xvec = new Double_t[hist->GetNbinsX()];
584 Double_t *yvec = new Double_t[hist->GetNbinsX()];
585 Double_t *xerr = new Double_t[hist->GetNbinsX()];
586 Double_t *yerr = new Double_t[hist->GetNbinsX()];
588 TH2D * projectionHist = 0x0;
590 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
594 for (Int_t idelta=0; idelta<nmaxBin; idelta++){
596 imin = TMath::Max(i-idelta,1);
597 imax = TMath::Min(i+idelta,hist->GetNbinsX());
598 nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
600 if (nsum>minEntries) break;
602 if (nsum<minEntries) continue;
604 hist->GetXaxis()->SetRange(imin,imax);
605 projectionHist = (TH2D*)hist->Project3D("yzNUFNOF");
608 projectionHist->FitSlicesY(0,2, projectionHist->GetNbinsX()-2,0,"QNR",&arr);
609 TH1D * histAttach = (TH1D*)arr.At(1);
610 TF1 pol1("polynom1","pol1",10,240);
611 histAttach->Fit(&pol1,"QNR");
612 xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin)+hist->GetXaxis()->GetBinCenter(imax));
613 yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
615 yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
618 delete projectionHist;
621 TGraphErrors * graphErrors = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
633 void AliTPCcalibTimeGain::BinLogX(THnSparse *h, Int_t axisDim) {
635 // Method for the correct logarithmic binning of histograms
637 TAxis *axis = h->GetAxis(axisDim);
638 int bins = axis->GetNbins();
640 Double_t from = axis->GetXmin();
641 Double_t to = axis->GetXmax();
642 Double_t *newBins = new Double_t[bins + 1];
645 Double_t factor = pow(to/from, 1./bins);
647 for (int i = 1; i <= bins; i++) {
648 newBins[i] = factor * newBins[i-1];
650 axis->Set(bins, newBins);
656 void AliTPCcalibTimeGain::BinLogX(TH1 *h) {
658 // Method for the correct logarithmic binning of histograms
660 TAxis *axis = h->GetXaxis();
661 int bins = axis->GetNbins();
663 Double_t from = axis->GetXmin();
664 Double_t to = axis->GetXmax();
665 Double_t *newBins = new Double_t[bins + 1];
668 Double_t factor = pow(to/from, 1./bins);
670 for (int i = 1; i <= bins; i++) {
671 newBins[i] = factor * newBins[i-1];
673 axis->Set(bins, newBins);
680 void AliTPCcalibTimeGain::CalculateBetheAlephParams(TH2F *hist, Double_t * ini) {
682 // Fit the bethe bloch params
684 //{0.0762*MIP,10.632,1.34e-05,1.863,1.948}
685 const Double_t sigma = 0.06;
687 TH2F *histBG = new TH2F("histBG","dEdxBg; #beta #gamma; TPC signal (a.u.)",TMath::Nint(0.5*hist->GetNbinsX()),0.1,5000.,hist->GetNbinsY(),hist->GetYaxis()->GetXmin(),hist->GetYaxis()->GetXmax());
690 TF1 *foElectron = new TF1("foElectron", "AliExternalTrackParam::BetheBlochAleph(x/0.000511,[0],[1],[2],[3],[4])",0.01,100);
691 TF1 *foMuon = new TF1("foMuon", "AliExternalTrackParam::BetheBlochAleph(x/0.1056,[0],[1],[2],[3],[4])",0.01,100);
692 TF1 *foPion = new TF1("foPion", "AliExternalTrackParam::BetheBlochAleph(x/0.138,[0],[1],[2],[3],[4])",0.01,100);
693 TF1 *foKaon = new TF1("foKaon", "AliExternalTrackParam::BetheBlochAleph(x/0.498,[0],[1],[2],[3],[4])",0.01,100);
694 TF1 *foProton = new TF1("foProton", "AliExternalTrackParam::BetheBlochAleph(x/0.938,[0],[1],[2],[3],[4])",0.01,100);
695 foElectron->SetParameters(ini);
696 foMuon->SetParameters(ini);
697 foPion->SetParameters(ini);
698 foKaon->SetParameters(ini);
699 foProton->SetParameters(ini);
701 TCanvas *canvCheck1 = new TCanvas();
703 foElectron->Draw("same");
704 foMuon->Draw("same");
705 foPion->Draw("same");
706 foKaon->Draw("same");
707 foProton->Draw("same");
709 // Loop over all points of the input histogram
711 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
712 Double_t x = hist->GetXaxis()->GetBinCenter(i);
713 for(Int_t j=1; j < hist->GetNbinsY(); j++) {
714 Long64_t n = hist->GetBin(i, j);
715 Double_t y = hist->GetYaxis()->GetBinCenter(j);
716 Double_t entries = hist->GetBinContent(n);
719 // 1. identify protons
721 if (TMath::Abs(y - foProton->Eval(x))/foProton->Eval(x) < 4*sigma && x < 0.65 ) {
722 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
725 // 2. identify electrons
728 if (TMath::Abs(y - foElectron->Eval(x))/foElectron->Eval(x) < 4*sigma && (x>0.25&&x<0.7) && fIsCosmic) {
729 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
732 if (TMath::Abs(y - foElectron->Eval(x))/foElectron->Eval(x) < 2*sigma && ((x>0.25&&x<0.35) || (x>1.5&&x<1.8) || (x>0.65&&x<0.7))) {
733 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
737 // 3. identify either muons or pions depending on cosmic or not
740 if (TMath::Abs(y - foMuon->Eval(x))/foMuon->Eval(x) < 4*sigma && x > 0.25 && x < 100 ) {
741 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
745 if (TMath::Abs(y - foPion->Eval(x))/foPion->Eval(x) < 4*sigma && x > 0.15 && x < 2) {
746 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
750 // 4. for pp also Kaons must be included
753 if (TMath::Abs(y - foKaon->Eval(x))/foKaon->Eval(x) < 4*sigma && x > 0.2 && x < 0.3 ) {
754 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
760 // Fit Aleph-Parameters to the obtained profile
761 TF1 * funcAlephD = new TF1("AlephParametrizationD", "AliExternalTrackParam::BetheBlochAleph(x,[0],[1],[2],[3],[4])",0.3,10000);
762 funcAlephD->SetParameters(ini);
764 TCanvas *canvCheck2 = new TCanvas();
768 TObjArray * arr = new TObjArray();
769 histBG->FitSlicesY(0,0,-1,0,"QN",arr);
770 TH1D * fitMean = (TH1D*) arr->At(1);
771 fitMean->Draw("same");
773 funcAlephD->SetParLimits(2,1e-3,1e-7);
774 funcAlephD->SetParLimits(3,0.5,3.5);
775 funcAlephD->SetParLimits(4,0.5,3.5);
776 fitMean->Fit(funcAlephD, "QNR");
777 funcAlephD->Draw("same");
779 for(Int_t i=0;i<5;i++) ini[i] = funcAlephD->GetParameter(i);
781 foElectron->SetParameters(ini);
782 foMuon->SetParameters(ini);
783 foPion->SetParameters(ini);
784 foKaon->SetParameters(ini);
785 foProton->SetParameters(ini);
787 TCanvas *canvCheck3 = new TCanvas();
789 foElectron->Draw("same");
790 foMuon->Draw("same");
791 foPion->Draw("same");
792 foKaon->Draw("same");
793 foProton->Draw("same");