]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - TPC/AliTPCcalibTimeGain.cxx
Cosmic tracker as an intergral part of the central ALICE reconstruction.
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibTimeGain.cxx
... / ...
CommitLineData
1/**************************************************************************
2 * Copyright(c) 1998-1999, 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
18This class provides the calibration of the time dependence of the TPC gain due to pressure and temperature changes etc.
19
20
21//0. Libraries to load
22gSystem->Load("libANALYSIS");
23gSystem->Load("libSTAT");
24gSystem->Load("libTPCcalib");
25
26
27//1. Do calibration ...
28//
29// compare reference
30
31//
32//2. Visualize results
33//
34TFile fcalib("CalibObjects.root");
35TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
36AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
37TGraphErrors * gr = gain->GetGraphGainVsTime(0,1000)
38
39// gain->GetHistGainTime()->GetAxis(1)->SetRangeUser(1213.8e6,1214.3e6)
40TH2D * GainTime = gain->GetHistGainTime()->Projection(0,1)
41GainTime->GetXaxis()->SetTimeDisplay(kTRUE)
42GainTime->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}")
43GainTime->Draw("colz")
44
45//
46// MakeSlineFit example
47//
48AliSplineFit *fit = AliTPCcalibTimeGain::MakeSplineFit(gr)
49
50TGraph * grfit = fit.MakeGraph(gr->GetX()[0],gr->GetX()[gr->GetN()-1],50000,0);
51
52gr->SetMarkerStyle(25);
53gr->Draw("lp");
54grfit->SetLineColor(2);
55grfit->Draw("lu");
56
57//
58// QA - dE/dx resoultion as a function of time
59//TCa
60
61TGraph * grSigma = AliTPCcalibBase::FitSlicesSigma(gain->GetHistGainTime(),0,1,1800,5)
62
63TCanvas *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
67 pad1->Draw();
68 pad1->cd();
69
70GainTime->Draw("colz")
71gr->Draw("lp")
72
73
74
75 c1->cd();
76 Double_t ymin = -0.04;
77 Double_t ymax = 0.12;
78Double_t dy = (ymax-ymin)/0.8;
79Double_t xmin = GainTime->GetXaxis()->GetXmin()
80Double_t xmax = GainTime->GetXaxis()->GetXmax()
81Double_t dx = (xmax-xmin)/0.8; //10 per cent margins left and right
82pad2->Range(xmin-0.1*dx,ymin-0.1*dy,xmax+0.1*dx,ymax+0.1*dy);
83 pad2->Draw();
84 pad2->cd();
85grSigma->SetLineColor(2);
86grSigma->SetLineWidth(2);
87grSigma->Draw("lp")
88TGaxis *axis = new TGaxis(xmax,ymin,xmax,ymax,ymin,ymax,50510,"+L");
89 axis->SetLabelColor(kRed);
90 axis->SetTitle("dE/dx resolution #sigma_{dE/dx}");
91 axis->Draw();
92
93////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
94
95 ----> make Attachment study
96
97TFile fcalib("CalibObjects40366b.root");
98TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
99AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
100TGraphErrors * grAttach = gain->GetGraphAttachment(2000,4)
101
102TCanvas *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
106 pad1->Draw();
107 pad1->cd();
108
109gain->GetHistGainTime()->GetAxis(1)->SetRangeUser(1213.8e6,1214.3e6)
110TH2D * GainTime = gain->GetHistGainTime()->Projection(0,1)
111GainTime->GetXaxis()->SetTimeDisplay(kTRUE)
112GainTime->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}")
113GainTime->Draw("colz")
114//gr->Draw("lp")
115
116 c1->cd();
117 Double_t ymin = -0.001;
118 Double_t ymax = 0.001;
119Double_t dy = (ymax-ymin)/0.8;
120Double_t xmin = GainTime->GetXaxis()->GetXmin()
121Double_t xmax = GainTime->GetXaxis()->GetXmax()
122Double_t dx = (xmax-xmin)/0.8; //10 per cent margins left and right
123pad2->Range(xmin-0.1*dx,ymin-0.1*dy,xmax+0.1*dx,ymax+0.1*dy);
124 pad2->Draw();
125 pad2->cd();
126grAttach->SetLineColor(2);
127grAttach->SetLineWidth(2);
128grAttach->Draw("lp")
129TGaxis *axis = new TGaxis(xmax,ymin,xmax,ymax,ymin,ymax,50510,"+L");
130 axis->SetLabelColor(kRed);
131 axis->SetTitle("attachment coefficient b");
132 axis->Draw();
133
134
135*/
136
137
138#include "Riostream.h"
139#include "TChain.h"
140#include "TTree.h"
141#include "TH1F.h"
142#include "TH2F.h"
143#include "TH3F.h"
144#include "THnSparse.h"
145#include "TGraphErrors.h"
146#include "TList.h"
147#include "TMath.h"
148#include "TCanvas.h"
149#include "TFile.h"
150#include "TF1.h"
151#include "TVectorD.h"
152#include "TProfile.h"
153#include "TGraphErrors.h"
154#include "TCanvas.h"
155
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"
163
164#include "AliTracker.h"
165#include "AliMagF.h"
166#include "AliTPCCalROC.h"
167#include "AliESDv0.h"
168
169#include "AliLog.h"
170
171#include "AliTPCcalibTimeGain.h"
172
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"
180
181ClassImp(AliTPCcalibTimeGain)
182
183
184AliTPCcalibTimeGain::AliTPCcalibTimeGain()
185 :AliTPCcalibBase(),
186 fHistGainTime(0),
187 fGainVsTime(0),
188 fHistDeDxTotal(0),
189 fIntegrationTimeDeDx(0),
190 fMIP(0),
191 fUseMax(0),
192 fLowerTrunc(0),
193 fUpperTrunc(0),
194 fUseShapeNorm(0),
195 fUsePosNorm(0),
196 fUsePadNorm(0),
197 fUseCookAnalytical(0),
198 fIsCosmic(0),
199 fLowMemoryConsumption(0)
200{
201 //
202 // Default constructor
203 //
204 AliInfo("Default Constructor");
205}
206
207
208AliTPCcalibTimeGain::AliTPCcalibTimeGain(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeGain)
209 :AliTPCcalibBase(),
210 fHistGainTime(0),
211 fGainVsTime(0),
212 fHistDeDxTotal(0),
213 fIntegrationTimeDeDx(0),
214 fMIP(0),
215 fUseMax(0),
216 fLowerTrunc(0),
217 fUpperTrunc(0),
218 fUseShapeNorm(0),
219 fUsePosNorm(0),
220 fUsePadNorm(0),
221 fUseCookAnalytical(0),
222 fIsCosmic(0),
223 fLowMemoryConsumption(0)
224{
225 //
226 // No default constructor
227 //
228 SetName(name);
229 SetTitle(title);
230
231 AliInfo("Non Default Constructor");
232
233 fIntegrationTimeDeDx = deltaIntegrationTimeGain;
234
235 Double_t deltaTime = EndTime - StartTime;
236
237
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);
245 //
246 fHistDeDxTotal = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",250,0.01,100.,1000,0.,8);
247 BinLogX(fHistDeDxTotal);
248
249 // default values for dE/dx
250 fMIP = 50.;
251 fUseMax = kTRUE;
252 fLowerTrunc = 0.02;
253 fUpperTrunc = 0.6;
254 fUseShapeNorm = kTRUE;
255 fUsePosNorm = kFALSE;
256 fUsePadNorm = kFALSE;
257 fUseCookAnalytical = kFALSE;
258 //
259 fIsCosmic = kTRUE;
260 fLowMemoryConsumption = kTRUE;
261 //
262
263}
264
265
266
267AliTPCcalibTimeGain::~AliTPCcalibTimeGain(){
268 //
269 // Destructor
270 //
271 delete fHistGainTime;
272 delete fGainVsTime;
273 delete fHistDeDxTotal;
274}
275
276
277void AliTPCcalibTimeGain::Process(AliESDEvent *event) {
278 //
279 // main track loop
280 //
281 if (!event) {
282 Printf("ERROR: ESD not available");
283 return;
284 }
285
286 if (fIsCosmic) { // this should be removed at some point based on trigger mask !?
287 ProcessCosmicEvent(event);
288 } else {
289 ProcessBeamEvent(event);
290 }
291
292
293
294
295}
296
297
298void AliTPCcalibTimeGain::ProcessCosmicEvent(AliESDEvent *event) {
299 //
300 // Process in case of cosmic event
301 //
302 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
303 if (!esdFriend) {
304 Printf("ERROR: ESDfriend not available");
305 return;
306 }
307 //
308 UInt_t time = event->GetTimeStamp();
309 Int_t ntracks = event->GetNumberOfTracks();
310 Int_t runNumber = event->GetRunNumber();
311 //
312 // track loop
313 //
314 for (Int_t i=0;i<ntracks;++i) {
315
316 AliESDtrack *track = event->GetTrack(i);
317 if (!track) continue;
318 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
319 if (!friendTrack) continue;
320 const AliExternalTrackParam * trackIn = track->GetInnerParam();
321 const AliExternalTrackParam * trackOut = friendTrack->GetTPCOut();
322 if (!trackIn) continue;
323 if (!trackOut) continue;
324
325 // calculate necessary track parameters
326 Double_t meanP = trackIn->GetP();
327 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ());
328 Int_t nclsDeDx = track->GetTPCNcls();
329
330 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
331 if (nclsDeDx < 60) continue;
332 if (TMath::Abs(trackIn->GetTgl()) > 1) continue;
333 if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue;
334
335 // Get seeds
336 TObject *calibObject;
337 AliTPCseed *seed = 0;
338 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
339 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
340 }
341
342 if (seed) {
343 Double_t tpcSignal = GetTPCdEdx(seed);
344 if (nclsDeDx > 100) fHistDeDxTotal->Fill(meanP, tpcSignal);
345 //
346 if (fLowMemoryConsumption) {
347 if (meanP < 20) continue;
348 meanP = 30; // set all momenta to one in order to save memory
349 }
350 //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta
351 Double_t vec[6] = {tpcSignal,time,1,meanDrift,meanP,runNumber};
352 fHistGainTime->Fill(vec);
353 }
354
355 }
356
357}
358
359
360
361void AliTPCcalibTimeGain::ProcessBeamEvent(AliESDEvent *event) {
362 //
363 // Process in case of beam event
364 //
365 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
366 if (!esdFriend) {
367 Printf("ERROR: ESDfriend not available");
368 return;
369 }
370 //
371 UInt_t time = event->GetTimeStamp();
372 Int_t ntracks = event->GetNumberOfTracks();
373 Int_t runNumber = event->GetRunNumber();
374 //
375 // track loop
376 //
377 for (Int_t i=0;i<ntracks;++i) { // begin track loop
378
379 AliESDtrack *track = event->GetTrack(i);
380 if (!track) continue;
381 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
382 if (!friendTrack) continue;
383
384 const AliExternalTrackParam * trackIn = track->GetInnerParam();
385 const AliExternalTrackParam * trackOut = friendTrack->GetTPCOut();
386 if (!trackIn) continue;
387 if (!trackOut) continue;
388
389 // calculate necessary track parameters
390 Double_t meanP = trackIn->GetP();
391 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ());
392 Int_t nclsDeDx = track->GetTPCNcls();
393
394 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
395 if (nclsDeDx < 60) continue;
396 //if (TMath::Abs(trackIn->GetTgl()) > 1) continue;
397 //if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue;
398 if (TMath::Abs(trackIn->Eta()) > 1) continue;
399 UInt_t status = track->GetStatus();
400 if ((status&AliESDtrack::kTPCrefit)==0) continue;
401 //if (track->GetNcls(0) < 3) continue; // ITS clusters
402 Float_t dca[2], cov[3];
403 track->GetImpactParameters(dca,cov);
404 if (TMath::Abs(dca[0]) > 7 || TMath::Abs(dca[0]) < 0.0000001 || TMath::Abs(dca[1]) > 25 ) continue; // cut in xy
405 Double_t eta = trackIn->Eta();
406
407 // Get seeds
408 TObject *calibObject;
409 AliTPCseed *seed = 0;
410 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
411 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
412 }
413
414 if (seed) {
415 Int_t particleCase = 0;
416 if (meanP < 0.5 && meanP > 0.4) particleCase = 2; // MIP pions
417 if (meanP < 0.57 && meanP > 0.56) particleCase = 3; // protons 1
418 if (meanP < 0.66 && meanP > 0.65) particleCase = 4; // protons 2
419 //
420 if (fLowMemoryConsumption && particleCase == 0) continue;
421 //
422 Double_t tpcSignal = GetTPCdEdx(seed);
423 fHistDeDxTotal->Fill(meanP, tpcSignal);
424 //
425 //dE/dx, time, type (1-muon cosmic,2-pion beam data, 3&4 protons), momenta, runNumner, eta
426 Double_t vec[7] = {tpcSignal,time,particleCase,meanDrift,meanP,runNumber, eta};
427 fHistGainTime->Fill(vec);
428 }
429
430 } // end track loop
431 //
432 // V0 loop -- in beam events the cosmic part of the histogram is filled with GammaConversions
433 //
434 for(Int_t iv0 = 0; iv0 < event->GetNumberOfV0s(); iv0++) {
435 AliESDv0 * v0 = event->GetV0(iv0);
436 if (!v0->GetOnFlyStatus()) continue;
437 if (v0->GetEffMass(0,0) > 0.02) continue; // select low inv. mass
438 Double_t xyz[3];
439 v0->GetXYZ(xyz[0], xyz[1], xyz[2]);
440 if (TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]) < 3 || TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]) > 30) continue;
441 //
442 // "loop over daughters"
443 //
444 for(Int_t idaughter = 0; idaughter < 2; idaughter++) { // daughter loop
445 Int_t index = idaughter == 0 ? v0->GetPindex() : v0->GetNindex();
446 AliESDtrack * trackP = event->GetTrack(index);
447 AliESDfriendTrack *friendTrackP = esdFriend->GetTrack(index);
448 if (!friendTrackP) continue;
449 const AliExternalTrackParam * trackPIn = trackP->GetInnerParam();
450 const AliExternalTrackParam * trackPOut = friendTrackP->GetTPCOut();
451 if (!trackPIn) continue;
452 if (!trackPOut) continue;
453 // calculate necessary track parameters
454 Double_t meanP = trackPIn->GetP();
455 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackPIn->GetZ() + trackPOut->GetZ());
456 Int_t nclsDeDx = trackP->GetTPCNcls();
457 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
458 if (nclsDeDx < 60) continue;
459 if (TMath::Abs(trackPIn->GetTgl()) > 1) continue;
460 //
461 TObject *calibObject;
462 AliTPCseed *seed = 0;
463 for (Int_t l=0;(calibObject=friendTrackP->GetCalibObject(l));++l) {
464 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
465 }
466 if (seed) {
467 if (fLowMemoryConsumption) {
468 if (meanP > 0.5 || meanP < 0.4) continue;
469 meanP = 0.45; // set all momenta to one in order to save memory
470 }
471 Double_t tpcSignal = GetTPCdEdx(seed);
472 //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta
473 Double_t vec[6] = {tpcSignal,time,1,meanDrift,meanP,runNumber};
474 fHistGainTime->Fill(vec);
475 }
476 }
477
478 }
479
480}
481
482
483Float_t AliTPCcalibTimeGain::GetTPCdEdx(AliTPCseed * seed) {
484 //
485 // calculate tpc dEdx
486 //
487 Double_t signal = 0;
488 //
489 if (!fUseCookAnalytical) {
490 signal = (1/fMIP)*seed->CookdEdxNorm(fLowerTrunc,fUpperTrunc,fUseMax,0,159,fUseShapeNorm,fUsePosNorm,fUsePadNorm,0);
491 } else {
492 signal = (1/fMIP)*seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,fUseMax);
493 }
494 //
495 return signal;
496}
497
498
499void AliTPCcalibTimeGain::AnalyzeRun(Int_t minEntries) {
500 //
501 // Analyze results of calibration
502 //
503 if (fIsCosmic) {
504 fHistGainTime->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
505 fHistGainTime->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons
506 } else {
507 fHistGainTime->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
508 fHistGainTime->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
509 }
510 //
511 fGainVsTime = AliTPCcalibBase::FitSlices(fHistGainTime,0,1,minEntries,10);
512 //
513 return;
514}
515
516
517TGraphErrors * AliTPCcalibTimeGain::GetGraphGainVsTime(Int_t runNumber, Int_t minEntries) {
518 //
519 // Analyze results and get the graph
520 //
521 if (runNumber == 0) {
522 if (!fGainVsTime) {
523 AnalyzeRun(minEntries);
524 }
525 } else {
526 // 1st check if the current run was cosmic or beam event
527 fHistGainTime->GetAxis(5)->SetRangeUser(runNumber,runNumber);
528 AnalyzeRun(minEntries);
529 }
530 if (fGainVsTime->GetN() == 0) return 0;
531 return fGainVsTime;
532}
533
534Long64_t AliTPCcalibTimeGain::Merge(TCollection *li) {
535 //
536 // merge component
537 //
538 TIterator* iter = li->MakeIterator();
539 AliTPCcalibTimeGain* cal = 0;
540
541 while ((cal = (AliTPCcalibTimeGain*)iter->Next())) {
542 if (!cal->InheritsFrom(AliTPCcalibTimeGain::Class())) {
543 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
544 return -1;
545 }
546
547 // add histograms here...
548 if (cal->GetHistGainTime()) fHistGainTime->Add(cal->GetHistGainTime());
549 if (cal->GetHistDeDxTotal()) fHistDeDxTotal->Add(cal->GetHistDeDxTotal());
550
551 }
552 delete iter;
553 return 0;
554
555}
556
557
558AliSplineFit * AliTPCcalibTimeGain::MakeSplineFit(TGraphErrors * graph) {
559 //
560 // make spline fit of gain
561 //
562 AliSplineFit *fit = new AliSplineFit();
563 fit->SetGraph(graph);
564 fit->SetMinPoints(graph->GetN()+1);
565 fit->InitKnots(graph,2,0,0.001);
566 fit->SplineFit(0);
567 return fit;
568
569}
570
571
572
573TGraphErrors * AliTPCcalibTimeGain::GetGraphAttachment(Int_t minEntries, Int_t nmaxBin, Float_t /*fracLow*/, Float_t /*fracUp*/) {
574 //
575 // For each time bin the driftlength-dependence of the signal is fitted.
576 //
577 TH3D * hist = fHistGainTime->Projection(1, 0, 3);
578 Double_t *xvec = new Double_t[hist->GetNbinsX()];
579 Double_t *yvec = new Double_t[hist->GetNbinsX()];
580 Double_t *xerr = new Double_t[hist->GetNbinsX()];
581 Double_t *yerr = new Double_t[hist->GetNbinsX()];
582 Int_t counter = 0;
583 TH2D * projectionHist = 0x0;
584 //
585 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
586 Int_t nsum=0;
587 Int_t imin = i;
588 Int_t imax = i;
589 for (Int_t idelta=0; idelta<nmaxBin; idelta++){
590 //
591 imin = TMath::Max(i-idelta,1);
592 imax = TMath::Min(i+idelta,hist->GetNbinsX());
593 nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
594 if (nsum==0) break;
595 if (nsum>minEntries) break;
596 }
597 if (nsum<minEntries) continue;
598 //
599 hist->GetXaxis()->SetRange(imin,imax);
600 projectionHist = (TH2D*)hist->Project3D("yzNUFNOF");
601 //
602 TObjArray arr;
603 projectionHist->FitSlicesY(0,2, projectionHist->GetNbinsX()-2,0,"QNR",&arr);
604 TH1D * histAttach = (TH1D*)arr.At(1);
605 TF1 pol1("polynom1","pol1",10,240);
606 histAttach->Fit(&pol1,"QNR");
607 xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin)+hist->GetXaxis()->GetBinCenter(imax));
608 yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
609 xerr[counter] = 0;
610 yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
611 counter++;
612 //
613 delete projectionHist;
614 }
615
616 TGraphErrors * graphErrors = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
617 delete [] xvec;
618 delete [] yvec;
619 delete [] xerr;
620 delete [] yerr;
621 delete hist;
622 return graphErrors;
623
624}
625
626
627
628void AliTPCcalibTimeGain::BinLogX(THnSparse *h, Int_t axisDim) {
629 //
630 // Method for the correct logarithmic binning of histograms
631 //
632 TAxis *axis = h->GetAxis(axisDim);
633 int bins = axis->GetNbins();
634
635 Double_t from = axis->GetXmin();
636 Double_t to = axis->GetXmax();
637 Double_t *newBins = new Double_t[bins + 1];
638
639 newBins[0] = from;
640 Double_t factor = pow(to/from, 1./bins);
641
642 for (int i = 1; i <= bins; i++) {
643 newBins[i] = factor * newBins[i-1];
644 }
645 axis->Set(bins, newBins);
646 delete[] newBins;
647
648}
649
650
651void AliTPCcalibTimeGain::BinLogX(TH1 *h) {
652 //
653 // Method for the correct logarithmic binning of histograms
654 //
655 TAxis *axis = h->GetXaxis();
656 int bins = axis->GetNbins();
657
658 Double_t from = axis->GetXmin();
659 Double_t to = axis->GetXmax();
660 Double_t *newBins = new Double_t[bins + 1];
661
662 newBins[0] = from;
663 Double_t factor = pow(to/from, 1./bins);
664
665 for (int i = 1; i <= bins; i++) {
666 newBins[i] = factor * newBins[i-1];
667 }
668 axis->Set(bins, newBins);
669 delete[] newBins;
670
671}
672
673
674
675void AliTPCcalibTimeGain::CalculateBetheAlephParams(TH2F *hist, Double_t * ini) {
676 //
677 // Fit the bethe bloch params
678 //
679 //{0.0762*MIP,10.632,1.34e-05,1.863,1.948}
680 const Double_t sigma = 0.06;
681
682 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());
683 BinLogX(histBG);
684
685 TF1 *foElectron = new TF1("foElectron", "AliExternalTrackParam::BetheBlochAleph(x/0.000511,[0],[1],[2],[3],[4])",0.01,100);
686 TF1 *foMuon = new TF1("foMuon", "AliExternalTrackParam::BetheBlochAleph(x/0.1056,[0],[1],[2],[3],[4])",0.01,100);
687 TF1 *foPion = new TF1("foPion", "AliExternalTrackParam::BetheBlochAleph(x/0.138,[0],[1],[2],[3],[4])",0.01,100);
688 TF1 *foKaon = new TF1("foKaon", "AliExternalTrackParam::BetheBlochAleph(x/0.498,[0],[1],[2],[3],[4])",0.01,100);
689 TF1 *foProton = new TF1("foProton", "AliExternalTrackParam::BetheBlochAleph(x/0.938,[0],[1],[2],[3],[4])",0.01,100);
690 foElectron->SetParameters(ini);
691 foMuon->SetParameters(ini);
692 foPion->SetParameters(ini);
693 foKaon->SetParameters(ini);
694 foProton->SetParameters(ini);
695
696 TCanvas *canvCheck1 = new TCanvas();
697 hist->Draw("colz");
698 foElectron->Draw("same");
699 foMuon->Draw("same");
700 foPion->Draw("same");
701 foKaon->Draw("same");
702 foProton->Draw("same");
703
704 // Loop over all points of the input histogram
705
706 for(Int_t i=1; i < hist->GetNbinsX(); i++) {
707 Double_t x = hist->GetXaxis()->GetBinCenter(i);
708 for(Int_t j=1; j < hist->GetNbinsY(); j++) {
709 Long64_t n = hist->GetBin(i, j);
710 Double_t y = hist->GetYaxis()->GetBinCenter(j);
711 Double_t entries = hist->GetBinContent(n);
712 Double_t mass = 0;
713
714 // 1. identify protons
715 mass = 0.938;
716 if (TMath::Abs(y - foProton->Eval(x))/foProton->Eval(x) < 4*sigma && x < 0.65 ) {
717 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
718 }
719
720 // 2. identify electrons
721 mass = 0.000511;
722 if (fIsCosmic) {
723 if (TMath::Abs(y - foElectron->Eval(x))/foElectron->Eval(x) < 4*sigma && (x>0.25&&x<0.7) && fIsCosmic) {
724 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
725 }
726 } else {
727 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))) {
728 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
729 }
730 }
731
732 // 3. identify either muons or pions depending on cosmic or not
733 if (fIsCosmic) {
734 mass = 0.1056;
735 if (TMath::Abs(y - foMuon->Eval(x))/foMuon->Eval(x) < 4*sigma && x > 0.25 && x < 100 ) {
736 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
737 }
738 } else {
739 mass = 0.1396;
740 if (TMath::Abs(y - foPion->Eval(x))/foPion->Eval(x) < 4*sigma && x > 0.15 && x < 2) {
741 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
742 }
743 }
744
745 // 4. for pp also Kaons must be included
746 if (!fIsCosmic) {
747 mass = 0.4936;
748 if (TMath::Abs(y - foKaon->Eval(x))/foKaon->Eval(x) < 4*sigma && x > 0.2 && x < 0.3 ) {
749 for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y);
750 }
751 }
752 }
753 }
754
755 // Fit Aleph-Parameters to the obtained profile
756 TF1 * funcAlephD = new TF1("AlephParametrizationD", "AliExternalTrackParam::BetheBlochAleph(x,[0],[1],[2],[3],[4])",0.3,10000);
757 funcAlephD->SetParameters(ini);
758
759 TCanvas *canvCheck2 = new TCanvas();
760 histBG->Draw();
761
762 //FitSlices
763 TObjArray * arr = new TObjArray();
764 histBG->FitSlicesY(0,0,-1,0,"QN",arr);
765 TH1D * fitMean = (TH1D*) arr->At(1);
766 fitMean->Draw("same");
767
768 funcAlephD->SetParLimits(2,1e-3,1e-7);
769 funcAlephD->SetParLimits(3,0.5,3.5);
770 funcAlephD->SetParLimits(4,0.5,3.5);
771 fitMean->Fit(funcAlephD, "QNR");
772 funcAlephD->Draw("same");
773
774 for(Int_t i=0;i<5;i++) ini[i] = funcAlephD->GetParameter(i);
775
776 foElectron->SetParameters(ini);
777 foMuon->SetParameters(ini);
778 foPion->SetParameters(ini);
779 foKaon->SetParameters(ini);
780 foProton->SetParameters(ini);
781
782 TCanvas *canvCheck3 = new TCanvas();
783 hist->Draw("colz");
784 foElectron->Draw("same");
785 foMuon->Draw("same");
786 foPion->Draw("same");
787 foKaon->Draw("same");
788 foProton->Draw("same");
789
790 canvCheck1->Print();
791 canvCheck2->Print();
792 canvCheck3->Print();
793
794 return;
795
796
797}