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)
183 Double_t AliTPCcalibTimeGain::fgMergeEntriesCut=10000000.;
185 AliTPCcalibTimeGain::AliTPCcalibTimeGain()
190 fIntegrationTimeDeDx(0),
194 fCutRequireITSrefit(0),
206 fUseCookAnalytical(0),
208 fLowMemoryConsumption(0)
211 // Default constructor
213 AliDebug(5,"Default Constructor");
217 AliTPCcalibTimeGain::AliTPCcalibTimeGain(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeGain)
222 fIntegrationTimeDeDx(0),
226 fCutRequireITSrefit(0),
238 fUseCookAnalytical(0),
240 fLowMemoryConsumption(0)
243 // No default constructor
248 AliDebug(5,"Non Default Constructor");
250 fIntegrationTimeDeDx = deltaIntegrationTimeGain;
252 Double_t deltaTime = EndTime - StartTime;
255 // 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
256 Int_t timeBins = TMath::Nint(deltaTime/deltaIntegrationTimeGain);
257 Int_t binsGainTime[7] = {150, static_cast<Int_t>(timeBins), 4, 25, 200, 10000000, 20};
258 Double_t xminGainTime[7] = {0.5, static_cast<Double_t>(StartTime), 0.5, 0, 0.1, -0.5, -1};
259 Double_t xmaxGainTime[7] = { 8, static_cast<Double_t>(EndTime), 4.5, 250, 50, 9999999.5, 1};
260 fHistGainTime = new THnSparseF("HistGainTime","dEdx time dep.;dEdx,time,type,driftlength,momenta,run number, eta;dEdx",7,binsGainTime,xminGainTime,xmaxGainTime);
261 BinLogX(fHistGainTime, 4);
263 fHistDeDxTotal = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",250,0.01,100.,1000,0.,8);
264 BinLogX(fHistDeDxTotal);
267 // default track selection cuts
270 fCutRequireITSrefit = kFALSE;
274 // default values for MIP window selection
275 fMinMomentumMIP = 0.4;
276 fMaxMomentumMIP = 0.6;
277 fAlephParameters[0] = 0.07657; // the following parameters work for most of the periods and are therefore default
278 fAlephParameters[1] = 10.6654; // but they can be overwritten in the train setup of cpass0
279 fAlephParameters[2] = 2.51466e-14;
280 fAlephParameters[3] = 2.05379;
281 fAlephParameters[4] = 1.84288;
283 // default values for dE/dx
288 fUseShapeNorm = kTRUE;
289 fUsePosNorm = kFALSE;
290 fUsePadNorm = kFALSE;
291 fUseCookAnalytical = kFALSE;
294 fLowMemoryConsumption = kTRUE;
301 AliTPCcalibTimeGain::~AliTPCcalibTimeGain(){
305 delete fHistGainTime;
307 delete fHistDeDxTotal;
311 void AliTPCcalibTimeGain::Process(AliESDEvent *event) {
316 Printf("ERROR: ESD not available");
319 AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
323 if (ESDfriend->TestSkipBit()) return;
325 if (fIsCosmic) { // this should be removed at some point based on trigger mask !?
326 ProcessCosmicEvent(event);
328 ProcessBeamEvent(event);
337 void AliTPCcalibTimeGain::ProcessCosmicEvent(AliESDEvent *event) {
339 // Process in case of cosmic event
341 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
343 Printf("ERROR: ESDfriend not available");
347 UInt_t time = event->GetTimeStamp();
348 Int_t nFriendTracks = esdFriend->GetNumberOfTracks();
349 Int_t runNumber = event->GetRunNumber();
353 for (Int_t i=0;i<nFriendTracks;++i) {
355 AliESDtrack *track = event->GetTrack(i);
356 if (!track) continue;
357 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
358 if (!friendTrack) continue;
359 const AliExternalTrackParam * trackIn = track->GetInnerParam();
360 const AliExternalTrackParam * trackOut = friendTrack->GetTPCOut();
361 if (!trackIn) continue;
362 if (!trackOut) continue;
364 // calculate necessary track parameters
365 Double_t meanP = trackIn->GetP();
366 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ());
367 Int_t nclsDeDx = track->GetTPCNcls();
369 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
370 if (nclsDeDx < 60) continue;
371 if (TMath::Abs(trackIn->GetTgl()) > 1) continue;
372 if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue;
375 TObject *calibObject;
376 AliTPCseed *seed = 0;
377 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
378 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
382 Double_t tpcSignal = GetTPCdEdx(seed);
383 if (nclsDeDx > 100) fHistDeDxTotal->Fill(meanP, tpcSignal);
385 if (fLowMemoryConsumption) {
386 if (meanP < 20) continue;
387 meanP = 30; // set all momenta to one in order to save memory
389 //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta
390 Double_t vec[6] = {tpcSignal,static_cast<Double_t>(time),1,meanDrift,meanP,static_cast<Double_t>(runNumber)};
391 fHistGainTime->Fill(vec);
400 void AliTPCcalibTimeGain::ProcessBeamEvent(AliESDEvent *event) {
402 // Process in case of beam event
404 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
406 Printf("ERROR: ESDfriend not available");
410 UInt_t time = event->GetTimeStamp();
411 Int_t nFriendTracks = esdFriend->GetNumberOfTracks();
412 Int_t runNumber = event->GetRunNumber();
416 for (Int_t i=0;i<nFriendTracks;++i) { // begin track loop
418 AliESDtrack *track = event->GetTrack(i);
419 if (!track) continue;
420 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
421 if (!friendTrack) continue;
423 const AliExternalTrackParam * trackIn = track->GetInnerParam();
424 const AliExternalTrackParam * trackOut = friendTrack->GetTPCOut();
425 if (!trackIn) continue;
426 if (!trackOut) continue;
428 // calculate necessary track parameters
429 Double_t meanP = trackIn->GetP();
430 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ());
431 Int_t nCrossedRows = track->GetTPCCrossedRows();
434 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
435 // fCutCrossRows = 80, fCutEtaWindow = 0.8, fCutRequireITSrefit, fCutMaxDcaXY = 3.5, fCutMaxDcaZ = 25
437 if (nCrossedRows < fCutCrossRows) continue;
438 if (TMath::Abs(trackIn->Eta()) > fCutEtaWindow) continue;
440 UInt_t status = track->GetStatus();
441 if ((status&AliESDtrack::kTPCrefit)==0) continue;
442 if ((status&AliESDtrack::kITSrefit)==0 && fCutRequireITSrefit) continue; // ITS cluster
444 Float_t dca[2], cov[3];
445 track->GetImpactParameters(dca,cov);
446 if (TMath::Abs(dca[0]) > fCutMaxDcaXY || TMath::Abs(dca[0]) < 0.0000001) continue; // cut in xy
447 if (((status&AliESDtrack::kITSrefit) == 1 && TMath::Abs(dca[1]) > 3.) || TMath::Abs(dca[1]) > fCutMaxDcaZ ) continue;
449 Double_t eta = trackIn->Eta();
452 TObject *calibObject;
453 AliTPCseed *seed = 0;
454 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
455 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
459 Int_t particleCase = 0;
460 if (meanP < fMaxMomentumMIP && meanP > fMinMomentumMIP) particleCase = 2; // MIP pions
461 if (meanP < 0.57 && meanP > 0.56) particleCase = 3; // protons 1
462 if (meanP < 0.66 && meanP > 0.65) particleCase = 4; // protons 2
464 if (fLowMemoryConsumption && particleCase == 0) continue;
466 Double_t tpcSignal = GetTPCdEdx(seed);
468 // flattens signal in MIP window
470 if (particleCase == 2) {
471 Float_t corrFactor = AliExternalTrackParam::BetheBlochAleph(meanP/0.13957,
476 fAlephParameters[4]);
477 tpcSignal /= corrFactor;
479 fHistDeDxTotal->Fill(meanP, tpcSignal);
481 //dE/dx, time, type (1-muon cosmic,2-pion beam data, 3&4 protons), momenta, runNumner, eta
482 Double_t vec[7] = {tpcSignal,static_cast<Double_t>(time),static_cast<Double_t>(particleCase),meanDrift,meanP,static_cast<Double_t>(runNumber), eta};
483 fHistGainTime->Fill(vec);
488 // V0 loop -- in beam events the cosmic part of the histogram is filled with GammaConversions
490 for(Int_t iv0 = 0; iv0 < event->GetNumberOfV0s(); iv0++) {
491 AliESDv0 * v0 = event->GetV0(iv0);
492 if (!v0->GetOnFlyStatus()) continue;
493 if (v0->GetEffMass(0,0) > 0.02) continue; // select low inv. mass
495 v0->GetXYZ(xyz[0], xyz[1], xyz[2]);
496 if (TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]) < 3 || TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]) > 30) continue;
498 // "loop over daughters"
500 for(Int_t idaughter = 0; idaughter < 2; idaughter++) { // daughter loop
501 Int_t index = idaughter == 0 ? v0->GetPindex() : v0->GetNindex();
502 AliESDtrack * trackP = event->GetTrack(index);
503 AliESDfriendTrack *friendTrackP = esdFriend->GetTrack(index);
504 if (!friendTrackP) continue;
505 const AliExternalTrackParam * trackPIn = trackP->GetInnerParam();
506 const AliExternalTrackParam * trackPOut = friendTrackP->GetTPCOut();
507 if (!trackPIn) continue;
508 if (!trackPOut) continue;
509 // calculate necessary track parameters
510 Double_t meanP = trackPIn->GetP();
511 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackPIn->GetZ() + trackPOut->GetZ());
512 Int_t nclsDeDx = trackP->GetTPCNcls();
513 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
514 if (nclsDeDx < 60) continue;
515 if (TMath::Abs(trackPIn->GetTgl()) > 1) continue;
517 TObject *calibObject;
518 AliTPCseed *seed = 0;
519 for (Int_t l=0;(calibObject=friendTrackP->GetCalibObject(l));++l) {
520 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
523 if (fLowMemoryConsumption) {
524 if (meanP > fMaxMomentumMIP || meanP < fMinMomentumMIP) continue;
525 meanP = 0.45; // set all momenta to one in order to save memory
527 Double_t tpcSignal = GetTPCdEdx(seed);
528 //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta
529 Double_t vec[6] = {tpcSignal,static_cast<Double_t>(time),1,meanDrift,meanP,static_cast<Double_t>(runNumber)};
530 fHistGainTime->Fill(vec);
539 Float_t AliTPCcalibTimeGain::GetTPCdEdx(AliTPCseed * seed) {
541 // calculate tpc dEdx
545 if (!fUseCookAnalytical) {
546 signal = (1/fMIP)*seed->CookdEdxNorm(fLowerTrunc,fUpperTrunc,fUseMax,0,159,fUseShapeNorm,fUsePosNorm,fUsePadNorm,0);
548 signal = (1/fMIP)*seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,fUseMax);
555 void AliTPCcalibTimeGain::AnalyzeRun(Int_t minEntries) {
557 // Analyze results of calibration
560 fHistGainTime->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
561 fHistGainTime->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
563 fHistGainTime->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
564 fHistGainTime->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
567 fGainVsTime = AliTPCcalibBase::FitSlices(fHistGainTime,0,1,minEntries,10);
573 TGraphErrors * AliTPCcalibTimeGain::GetGraphGainVsTime(Int_t runNumber, Int_t minEntries) {
575 // Analyze results and get the graph
577 if (runNumber == 0) {
579 AnalyzeRun(minEntries);
582 // 1st check if the current run was cosmic or beam event
583 fHistGainTime->GetAxis(5)->SetRangeUser(runNumber,runNumber);
584 AnalyzeRun(minEntries);
586 if (fGainVsTime->GetN() == 0) return 0;
590 Long64_t AliTPCcalibTimeGain::Merge(TCollection *li) {
594 TIterator* iter = li->MakeIterator();
595 AliTPCcalibTimeGain* cal = 0;
597 while ((cal = (AliTPCcalibTimeGain*)iter->Next())) {
598 if (!cal->InheritsFrom(AliTPCcalibTimeGain::Class())) {
599 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
603 // add histograms here...
604 if (cal->GetHistGainTime() && fHistGainTime )
606 if ((fHistGainTime->GetEntries() + cal->GetHistGainTime()->GetEntries()) < fgMergeEntriesCut)
608 //AliInfo(Form("fHistGainTime has %.0f tracks, going to add %.0f\n",fHistGainTime->GetEntries(),cal->GetHistGainTime()->GetEntries()));
609 fHistGainTime->Add(cal->GetHistGainTime());
613 AliInfo(Form("fHistGainTime full (has %.0f merged tracks, max allowed: %.0f)",fHistGainTime->GetEntries(),fgMergeEntriesCut));
617 if (cal->GetHistDeDxTotal()) fHistDeDxTotal->Add(cal->GetHistDeDxTotal());
626 AliSplineFit * AliTPCcalibTimeGain::MakeSplineFit(TGraphErrors * graph) {
628 // make spline fit of gain
630 AliSplineFit *fit = new AliSplineFit();
631 fit->SetGraph(graph);
632 fit->SetMinPoints(graph->GetN()+1);
633 fit->InitKnots(graph,2,0,0.001);
641 TGraphErrors * AliTPCcalibTimeGain::GetGraphAttachment(Int_t minEntries, Int_t nmaxBin, Float_t /*fracLow*/, Float_t /*fracUp*/) {
643 // For each time bin the driftlength-dependence of the signal is fitted.
645 TH3D * hist = fHistGainTime->Projection(1, 0, 3);
646 Double_t *xvec = new Double_t[hist->GetNbinsX()];
647 Double_t *yvec = new Double_t[hist->GetNbinsX()];
648 Double_t *xerr = new Double_t[hist->GetNbinsX()];
649 Double_t *yerr = new Double_t[hist->GetNbinsX()];
651 TH2D * projectionHist = 0x0;
653 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
657 for (Int_t idelta=0; idelta<nmaxBin; idelta++){
659 imin = TMath::Max(i-idelta,1);
660 imax = TMath::Min(i+idelta,hist->GetNbinsX());
661 nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
663 if (nsum>minEntries) break;
665 if (nsum<minEntries) continue;
667 hist->GetXaxis()->SetRange(imin,imax);
668 projectionHist = (TH2D*)hist->Project3D("yzNUFNOF");
671 projectionHist->FitSlicesY(0,2, projectionHist->GetNbinsX()-2,0,"QNR",&arr);
672 TH1D * histAttach = (TH1D*)arr.At(1);
673 TF1 pol1("polynom1","pol1",10,240);
674 histAttach->Fit(&pol1,"QNR");
675 xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin)+hist->GetXaxis()->GetBinCenter(imax));
676 yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
678 yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
681 delete projectionHist;
684 TGraphErrors * graphErrors = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
696 void AliTPCcalibTimeGain::BinLogX(THnSparse *h, Int_t axisDim) {
698 // Method for the correct logarithmic binning of histograms
700 TAxis *axis = h->GetAxis(axisDim);
701 int bins = axis->GetNbins();
703 Double_t from = axis->GetXmin();
704 Double_t to = axis->GetXmax();
705 Double_t *newBins = new Double_t[bins + 1];
708 Double_t factor = pow(to/from, 1./bins);
710 for (int i = 1; i <= bins; i++) {
711 newBins[i] = factor * newBins[i-1];
713 axis->Set(bins, newBins);
719 void AliTPCcalibTimeGain::BinLogX(TH1 *h) {
721 // Method for the correct logarithmic binning of histograms
723 TAxis *axis = h->GetXaxis();
724 int bins = axis->GetNbins();
726 Double_t from = axis->GetXmin();
727 Double_t to = axis->GetXmax();
728 Double_t *newBins = new Double_t[bins + 1];
731 Double_t factor = pow(to/from, 1./bins);
733 for (int i = 1; i <= bins; i++) {
734 newBins[i] = factor * newBins[i-1];
736 axis->Set(bins, newBins);
743 void AliTPCcalibTimeGain::CalculateBetheAlephParams(TH2F *hist, Double_t * ini) {
745 // Fit the bethe bloch params
747 //{0.0762*MIP,10.632,1.34e-05,1.863,1.948}
748 const Double_t sigma = 0.06;
750 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());
753 TF1 *foElectron = new TF1("foElectron", "AliExternalTrackParam::BetheBlochAleph(x/0.000511,[0],[1],[2],[3],[4])",0.01,100);
754 TF1 *foMuon = new TF1("foMuon", "AliExternalTrackParam::BetheBlochAleph(x/0.1056,[0],[1],[2],[3],[4])",0.01,100);
755 TF1 *foPion = new TF1("foPion", "AliExternalTrackParam::BetheBlochAleph(x/0.138,[0],[1],[2],[3],[4])",0.01,100);
756 TF1 *foKaon = new TF1("foKaon", "AliExternalTrackParam::BetheBlochAleph(x/0.498,[0],[1],[2],[3],[4])",0.01,100);
757 TF1 *foProton = new TF1("foProton", "AliExternalTrackParam::BetheBlochAleph(x/0.938,[0],[1],[2],[3],[4])",0.01,100);
758 foElectron->SetParameters(ini);
759 foMuon->SetParameters(ini);
760 foPion->SetParameters(ini);
761 foKaon->SetParameters(ini);
762 foProton->SetParameters(ini);
764 TCanvas *canvCheck1 = new TCanvas();
766 foElectron->Draw("same");
767 foMuon->Draw("same");
768 foPion->Draw("same");
769 foKaon->Draw("same");
770 foProton->Draw("same");
772 // Loop over all points of the input histogram
774 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
775 Double_t x = hist->GetXaxis()->GetBinCenter(i);
776 for(Int_t j=1; j < hist->GetNbinsY(); j++) {
777 Long64_t n = hist->GetBin(i, j);
778 Double_t y = hist->GetYaxis()->GetBinCenter(j);
779 Double_t entries = hist->GetBinContent(n);
782 // 1. identify protons
784 if (TMath::Abs(y - foProton->Eval(x))/foProton->Eval(x) < 4*sigma && x < 0.65 ) {
785 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
788 // 2. identify electrons
791 if (TMath::Abs(y - foElectron->Eval(x))/foElectron->Eval(x) < 4*sigma && (x>0.25&&x<0.7) && fIsCosmic) {
792 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
795 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))) {
796 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
800 // 3. identify either muons or pions depending on cosmic or not
803 if (TMath::Abs(y - foMuon->Eval(x))/foMuon->Eval(x) < 4*sigma && x > 0.25 && x < 100 ) {
804 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
808 if (TMath::Abs(y - foPion->Eval(x))/foPion->Eval(x) < 4*sigma && x > 0.15 && x < 2) {
809 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
813 // 4. for pp also Kaons must be included
816 if (TMath::Abs(y - foKaon->Eval(x))/foKaon->Eval(x) < 4*sigma && x > 0.2 && x < 0.3 ) {
817 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
823 // Fit Aleph-Parameters to the obtained profile
824 TF1 * funcAlephD = new TF1("AlephParametrizationD", "AliExternalTrackParam::BetheBlochAleph(x,[0],[1],[2],[3],[4])",0.3,10000);
825 funcAlephD->SetParameters(ini);
827 TCanvas *canvCheck2 = new TCanvas();
831 TObjArray * arr = new TObjArray();
832 histBG->FitSlicesY(0,0,-1,0,"QN",arr);
833 TH1D * fitMean = (TH1D*) arr->At(1);
834 fitMean->Draw("same");
836 funcAlephD->SetParLimits(2,1e-3,1e-7);
837 funcAlephD->SetParLimits(3,0.5,3.5);
838 funcAlephD->SetParLimits(4,0.5,3.5);
839 fitMean->Fit(funcAlephD, "QNR");
840 funcAlephD->Draw("same");
842 for(Int_t i=0;i<5;i++) ini[i] = funcAlephD->GetParameter(i);
844 foElectron->SetParameters(ini);
845 foMuon->SetParameters(ini);
846 foPion->SetParameters(ini);
847 foKaon->SetParameters(ini);
848 foProton->SetParameters(ini);
850 TCanvas *canvCheck3 = new TCanvas();
852 foElectron->Draw("same");
853 foMuon->Draw("same");
854 foPion->Draw("same");
855 foKaon->Draw("same");
856 foProton->Draw("same");