2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
20 Basic calibration and QA class for the TPC gain calibration based on tracks from BEAM events.
23 Send comments etc. to: A.Kalweit@gsi.de, marian.ivanov@cern.ch
27 #include "Riostream.h"
40 #include "AliTPCcalibDB.h"
41 #include "AliTPCclusterMI.h"
42 #include "AliTPCClusterParam.h"
43 #include "AliTPCseed.h"
44 #include "AliESDVertex.h"
45 #include "AliESDEvent.h"
46 #include "AliESDfriend.h"
47 #include "AliESDInputHandler.h"
48 #include "AliAnalysisManager.h"
49 #include "AliTPCParam.h"
51 #include "AliComplexCluster.h"
52 #include "AliTPCclusterMI.h"
56 #include "AliTPCcalibGainMult.h"
58 #include "TTreeStream.h"
61 ClassImp(AliTPCcalibGainMult)
64 AliTPCcalibGainMult::AliTPCcalibGainMult()
78 // Empty default cosntructor
80 AliInfo("Default Constructor");
84 AliTPCcalibGainMult::AliTPCcalibGainMult(const Text_t *name, const Text_t *title)
104 fLowerTrunc = 0.02; // IMPORTANT CHANGE --> REMOVE HARDWIRED TRUNCATION FROM TRACKER
106 fUseMax = kTRUE; // IMPORTANT CHANGE FOR PbPb; standard: kFALSE;
108 fHistNTracks = new TH1F("ntracks","Number of Tracks per Event; number of tracks per event; number of tracks",1001,-0.5,1000.5);
109 fHistClusterShape = new TH1F("fHistClusterShape","cluster shape; rms meas. / rms exp.;",300,0,3);
110 fHistQA = new TH3F("fHistQA","dEdx; momentum p (GeV); TPC signal (a.u.); pad",500,0.1,20.,500,0.,500,6,-0.5,5.5);
114 // MIP, sect, pad (short,med,long,full,oroc), run, ncl
115 Int_t binsGainSec[5] = { 145, 72, 4, 10000000, 65};
116 Double_t xminGainSec[5] = { 10., -0.5, -0.5, -0.5, -0.5};
117 Double_t xmaxGainSec[5] = {300., 71.5, 3.5, 9999999.5, 64.5};
118 TString axisNameSec[5]={"Q","sector","pad type","run", "ncl"};
119 TString axisTitleSec[5]={"Q (a.u)","sector","pad type","run","ncl"};
121 fHistGainSector = new THnSparseF("fHistGainSector","0:MIP, 1:sect, 2:pad, 3:run, 4:ncl", 5, binsGainSec, xminGainSec, xmaxGainSec);
122 for (Int_t iaxis=0; iaxis<5;iaxis++){
123 fHistGainSector->GetAxis(iaxis)->SetName(axisNameSec[iaxis]);
124 fHistGainSector->GetAxis(iaxis)->SetTitle(axisTitleSec[iaxis]);
129 Int_t binsPadEqual[6] = { 200, 200, 4, 20, 50, 100};
130 Double_t xminPadEqual[6] = { 0.5, 0.5, -0.5, 0, -250, 0};
131 Double_t xmaxPadEqual[6] = { 1.5, 1.5, 3.5, 4000, 250, 3};
132 TString axisNamePadEqual[6] = {"dEdxRatioMax","dEdxRatioTot","padType","mult","driftlength", "1_pt"};
133 TString axisTitlePadEqual[6] = {"dEdx_padRegion/mean_dEdx Qmax", "dEdx_padRegion/mean_dEdx Qtot","padType","mult","driftlength", "1/pt"};
135 fHistPadEqual = new THnSparseF("fHistPadEqual","0:dEdx_pad/dEdx_mean, 1:pad, 2:mult, 3:drift, 4:1/pt", 6, binsPadEqual, xminPadEqual, xmaxPadEqual);
136 for (Int_t iaxis=0; iaxis<6;iaxis++){
137 fHistPadEqual->GetAxis(iaxis)->SetName(axisNamePadEqual[iaxis]);
138 fHistPadEqual->GetAxis(iaxis)->SetTitle(axisTitlePadEqual[iaxis]);
142 // MIP Qmax, MIP Qtot, z, pad, vtx. contribut., ncl
143 Int_t binsGainMult[6] = { 145, 145, 25, 4, 100, 80};
144 Double_t xminGainMult[6] = { 10., 10., 0, -0.5, 0, -0.5};
145 Double_t xmaxGainMult[6] = {300., 300., 250, 3.5, 5000, 159.5};
146 TString axisNameMult[6]={"Qmax","Qtot","drift","padtype""multiplicity","ncl"};
147 TString axisTitleMult[6]={"Qmax (a.u)","Qtot (a.u.)","driftlenght l (cm)","Pad Type","multiplicity","ncl"};
149 fHistGainMult = new THnSparseF("fHistGainMult","MIP Qmax, MIP Qtot, z, type, vtx. contribut.", 6, binsGainMult, xminGainMult, xmaxGainMult);
150 for (Int_t iaxis=0; iaxis<6;iaxis++){
151 fHistGainMult->GetAxis(iaxis)->SetName(axisNameMult[iaxis]);
152 fHistGainMult->GetAxis(iaxis)->SetTitle(axisTitleMult[iaxis]);
155 AliInfo("Non Default Constructor");
160 AliTPCcalibGainMult::~AliTPCcalibGainMult(){
164 //delete fHistNTracks; // histogram showing number of ESD tracks per event
165 //delete fHistGainSector; // histogram showing the number of clusters per track
172 void AliTPCcalibGainMult::Process(AliESDEvent *event) {
177 Printf("ERROR: ESD not available");
180 Int_t ntracks=event->GetNumberOfTracks();
181 fHistNTracks->Fill(ntracks);
183 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
185 Printf("ERROR: esdFriend not available");
188 UInt_t runNumber = event->GetRunNumber();
189 Int_t nContributors = 0;
190 if (event->GetPrimaryVertexTPC()) nContributors = event->GetPrimaryVertexTPC()->GetNContributors();
194 for (Int_t i=0;i<ntracks;++i) {
196 AliESDtrack *track = event->GetTrack(i);
197 if (!track) continue;
199 AliExternalTrackParam * trackIn = (AliExternalTrackParam *)track->GetInnerParam();
200 if (!trackIn) continue;
202 // calculate necessary track parameters
203 Double_t meanP = trackIn->GetP();
204 Int_t ncls = track->GetTPCNcls();
206 if (ncls < 80) continue;
208 // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
210 if (TMath::Abs(trackIn->Eta()) > 0.8) continue;
211 UInt_t status = track->GetStatus();
212 if ((status&AliESDtrack::kTPCrefit)==0) continue;
213 //if (track->GetNcls(0) < 3) continue; // ITS clusters
214 Float_t dca[2], cov[3];
215 track->GetImpactParameters(dca,cov);
216 Float_t primVtxDCA = TMath::Sqrt(dca[0]*dca[0]);
217 if (primVtxDCA > 10 || primVtxDCA < 0.00001) continue;
218 if (TMath::Abs(dca[1]) > 5) continue;
220 // require that the track does not cross any dead area
222 //if (track->GetTPCNclsF() < 158) continue;
224 //if (seed->CookShape(1) > 1) continue;
225 //if (TMath::Abs(trackIn->GetY()) > 20) continue;
226 //if (TMath::Abs(d)>20) continue; // distance to the 0,0; select only tracks which cross chambers under proper angle
227 //if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue;
228 if (primVtxDCA < 3 && track->GetNcls(0) > 3 && track->GetKinkIndex(0) == 0 && ncls > 100) fHistQA->Fill(meanP, track->GetTPCsignal(), 5);
231 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
232 if (!friendTrack) continue;
233 TObject *calibObject;
234 AliTPCseed *seed = 0;
235 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
236 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
241 const AliExternalTrackParam * trackOut = friendTrack->GetTPCOut();
242 if (!trackIn) continue;
243 if (!trackOut) continue;
244 Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ());
246 for (Int_t irow =0; irow<160;irow++) {
247 AliTPCTrackerPoint * point = seed->GetTrackPoint(irow);
248 if (point==0) continue;
249 AliTPCclusterMI * cl = seed->GetClusterPointer(irow);
252 Float_t rsigmay = TMath::Sqrt(point->GetSigmaY());
253 fHistClusterShape->Fill(rsigmay);
259 Double_t signalShortMax = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,1,0,62);
260 Double_t signalMedMax = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,1,63,126);
261 Double_t signalLongMax = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,1,127,159);
262 Double_t signalMax = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,1,row0,row1);
263 Double_t signalArrayMax[4] = {signalShortMax, signalMedMax, signalLongMax, signalMax};
265 Double_t signalShortTot = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,0,0,62);
266 Double_t signalMedTot = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,0,63,126);
267 Double_t signalLongTot = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,0,127,159);
268 Double_t signalTot = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,0,row0,row1);
269 Double_t signalArrayTot[4] = {signalShortTot, signalMedTot, signalLongTot, signalTot};
271 Double_t mipSignalShort = fUseMax ? signalShortMax : signalShortTot;
272 Double_t mipSignalMed = fUseMax ? signalMedMax : signalMedTot;
273 Double_t mipSignalLong = fUseMax ? signalLongMax : signalLongTot;
274 Double_t mipSignalOroc = seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,fUseMax,63,159);
275 Double_t signal = fUseMax ? signalMax : signalTot;
277 fHistQA->Fill(meanP, mipSignalShort, 0);
278 fHistQA->Fill(meanP, mipSignalMed, 1);
279 fHistQA->Fill(meanP, mipSignalLong, 2);
280 fHistQA->Fill(meanP, signal, 3);
281 fHistQA->Fill(meanP, mipSignalOroc, 4);
283 // "dEdxRatioMax","dEdxRatioTot","padType","mult","driftlength", "1_pt"
284 Float_t meanMax = (1/3.)*(signalArrayMax[0] + signalArrayMax[1] + signalArrayMax[2]);
285 Float_t meanTot = (1/3.)*(signalArrayTot[0] + signalArrayTot[1] + signalArrayTot[2]);
286 if (meanMax || meanTot < 1e-5) continue;
287 for(Int_t ipad = 0; ipad < 4; ipad ++) {
288 Double_t vecPadEqual[6] = {signalArrayMax[ipad]/meanMax, signalArrayTot[ipad]/meanTot, ipad, nContributors, meanDrift, track->OneOverPt()};
289 fHistPadEqual->Fill(vecPadEqual);
292 if (meanP > 0.4 && meanP < 0.55) {
293 Double_t vecMult[6] = {seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,1,row0,row1),
294 seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,0,row0,row1),
300 fHistGainMult->Fill(vecMult);
301 vecMult[0]=mipSignalShort; vecMult[1]=mipSignalShort; vecMult[3]=0;
302 fHistGainMult->Fill(vecMult);
303 vecMult[0]=mipSignalMed; vecMult[1]=mipSignalMed; vecMult[3]=1;
304 fHistGainMult->Fill(vecMult);
305 vecMult[0]=mipSignalLong; vecMult[1]=mipSignalLong; vecMult[3]=2;
306 fHistGainMult->Fill(vecMult);
311 if (meanP > 0.5 || meanP < 0.4) continue; // only MIP pions
313 // for each track, we look at the three different pad regions, split it into tracklets, check if the sector does not change and fill the histogram
315 Bool_t isNotSplit[3] = {kTRUE, kTRUE, kTRUE}; // short, medium, long (true if the track is not split between two chambers)
317 Double_t sector[4] = {-1, -1, -1, -1}; // sector number short, medium, long, all
318 Int_t ncl[3] = {0,0,0};
320 for (Int_t irow=0; irow < 159; irow++){
322 if (irow > 62) padRegion = 1;
323 if (irow > 126) padRegion = 2;
325 AliTPCclusterMI* cluster = seed->GetClusterPointer(irow);
326 if (!cluster) continue;
327 if (sector[padRegion] == -1) {
328 sector[padRegion] = cluster->GetDetector();
331 if (sector[padRegion] != -1 && sector[padRegion] != cluster->GetDetector()) isNotSplit[padRegion] = kFALSE;
335 // MIP, sect, pad, run
337 Double_t vecMip[5] = {mipSignalShort, mipSignalMed, mipSignalLong, signal, mipSignalOroc};
339 for(Int_t ipad = 0; ipad < 3; ipad++) {
341 Double_t vecGainSec[5] = {vecMip[ipad], sector[ipad], ipad, runNumber, ncl[ipad]};
342 if (isNotSplit[ipad]) fHistGainSector->Fill(vecGainSec);
352 void AliTPCcalibGainMult::Analyze() {
360 Long64_t AliTPCcalibGainMult::Merge(TCollection *li) {
362 TIterator* iter = li->MakeIterator();
363 AliTPCcalibGainMult* cal = 0;
365 while ((cal = (AliTPCcalibGainMult*)iter->Next())) {
366 if (!cal->InheritsFrom(AliTPCcalibGainMult::Class())) {
367 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
371 if (cal->GetHistNTracks()) fHistNTracks->Add(cal->GetHistNTracks());
372 if (cal->GetHistClusterShape()) fHistClusterShape->Add(cal->GetHistClusterShape());
373 if (cal->GetHistQA()) fHistQA->Add(cal->GetHistQA());
374 if (cal->GetHistGainSector()) fHistGainSector->Add(cal->GetHistGainSector());
375 if (cal->GetHistPadEqual()) fHistPadEqual->Add(cal->GetHistPadEqual());
376 if (cal->GetHistGainMult()) fHistGainMult->Add(cal->GetHistGainMult());
386 void AliTPCcalibGainMult::BinLogX(const TH1 *h) {
388 // Method for the correct logarithmic binning of histograms
390 TAxis *axis = h->GetXaxis();
391 int bins = axis->GetNbins();
393 Double_t from = axis->GetXmin();
394 Double_t to = axis->GetXmax();
395 Double_t *newBins = new Double_t[bins + 1];
398 Double_t factor = pow(to/from, 1./bins);
400 for (int i = 1; i <= bins; i++) {
401 newBins[i] = factor * newBins[i-1];
403 axis->Set(bins, newBins);
410 void AliTPCcalibGainMult::UpdateGainMap() {
412 // read in the old gain map and scale it appropriately...
415 gSystem->Load("libANALYSIS");
416 gSystem->Load("libTPCcalib");
418 TFile jj("Run0_999999999_v1_s0.root");
419 AliTPCCalPad * pad = AliCDBEntry->GetObject()->Clone();
420 TFile hh("output.root");
421 AliTPCcalibGainMult * gain = calibTracksGain;
422 TH2D * histGainSec = gain->GetHistGainSector()->Projection(0,1);
425 histGainSec->FitSlicesY(0, 0, -1, 0, "QNR", &arr);
426 TH1D * meanGainSec = arr->At(1);
427 Double_t gainsIROC[36];
428 Double_t gainsOROC[36];
431 for(Int_t isec = 1; isec < meanGainSec->GetNbinsX() + 1; isec++) {
432 cout << isec << " " << meanGainSec->GetXaxis()->GetBinCenter(isec) << " " <<meanGainSec->GetBinContent(isec) << endl;
433 gains[isec-1] = meanGainSec->GetBinContent(isec);
435 gainsIROC[isec-1] = meanGainSec->GetBinContent(isec);
437 gainsOROC[isec - 36 -1] = meanGainSec->GetBinContent(isec);
440 Double_t meanIroc = TMath::Mean(36, gainsIROC);
441 Double_t meanOroc = TMath::Mean(36, gainsIROC);
442 for(Int_t i = 0; i < 36; i++) gains[i] /= meanIroc;
443 for(Int_t i = 36; i < 72; i++) gains[i] /= meanOroc;
445 for(Int_t i = 0; i < 72; i++) {
446 AliTPCCalROC * chamber = pad->GetCalROC(i);
447 chamber->Multiply(gains[i]);
448 cout << i << " "<< chamber->GetMean() << endl;
453 AliCDBMetaData *metaData= new AliCDBMetaData();
454 metaData->SetObjectClassName("AliTPCCalPad");
455 metaData->SetResponsible("Alexander Kalweit");
456 metaData->SetBeamPeriod(1);
457 metaData->SetAliRootVersion("04-19-05"); //root version
458 metaData->SetComment("New gain map for 1600V OROC gain increase and equalization. Valid for runs starting after technical stop beginning of September.");
459 AliCDBId id1("TPC/Calib/GainFactorDedx", 131541, AliCDBRunRange::Infinity()); // important: new gain runs here..
460 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage("local:///d/alice05/akalweit/projects/OCDBupdate/HighGain_2010-09-03/OCDB/");
461 gStorage->Put(pad, id1, metaData);
466 void AliTPCcalibGainMult::UpdateClusterParam() {
471 gSystem->Load("libANALYSIS");
472 gSystem->Load("libTPCcalib");
474 TFile ff("OldClsParam.root");
475 AliTPCClusterParam * param = AliCDBEntry->GetObject()->Clone();
477 TFile hh("output.root");
478 AliTPCcalibGainMult * gain = calibGainMult;
479 TH2D * histGainSec = gain->GetHistGainSector()->Projection(0,2);
481 histGainSec->FitSlicesY(0, 0, -1, 0, "QNR", &arr);
482 histGainSec->Draw("colz");
483 TH1D * fitVal = arr.At(1);
484 fitVal->Draw("same");
485 param->GetQnormCorrMatrix()->Print();
486 param->GetQnormCorrMatrix()(0,5) *= fitVal->GetBinContent(1)/fitVal->GetBinContent(1); // short pads Qtot
487 param->GetQnormCorrMatrix()(1,5) *= fitVal->GetBinContent(2)/fitVal->GetBinContent(1); // med pads Qtot
488 param->GetQnormCorrMatrix()(2,5) *= fitVal->GetBinContent(3)/fitVal->GetBinContent(1); // long pads Qtot
490 param->GetQnormCorrMatrix()(3,5) *= fitVal->GetBinContent(1)/fitVal->GetBinContent(1); // short pads Qmax -> scaling assumed
491 param->GetQnormCorrMatrix()(4,5) *= fitVal->GetBinContent(2)/fitVal->GetBinContent(1); // med pads Qmax -> scaling assumed
492 param->GetQnormCorrMatrix()(5,5) *= fitVal->GetBinContent(3)/fitVal->GetBinContent(1); // long pads Qmax -> scaling assumed
494 TFile jj("newClusterParam.root","RECREATE");
496 param->GetQnormCorrMatrix()->Print();
500 AliCDBMetaData *metaData= new AliCDBMetaData();
501 metaData->SetObjectClassName("AliTPCClusterParam");
502 metaData->SetResponsible("Alexander Kalweit");
503 metaData->SetBeamPeriod(1);
504 metaData->SetAliRootVersion("04-19-04"); //root version
505 metaData->SetComment("1600V OROC / hard thres. / new algorithm");
506 AliCDBId id1("TPC/Calib/ClusterParam", 0, AliCDBRunRange::Infinity()); // important: new gain runs here..
507 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage("local:///lustre/alice/akalweit/baseline/CalibrationEntries/OldThres_NewAlgo_PP");
508 gStorage->Put(param, id1, metaData);