/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /* This class provides the calibration of the time dependence of the TPC gain due to pressure and temperature changes etc. //0. Libraries to load gSystem->Load("libANALYSIS"); gSystem->Load("libSTAT"); gSystem->Load("libTPCcalib"); //1. Do calibration ... // // compare reference // //2. Visualize results // TFile fcalib("CalibObjects.root"); TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib"); AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain"); TGraphErrors * gr = gain->GetGraphGainVsTime(0,1000) // gain->GetHistGainTime()->GetAxis(1)->SetRangeUser(1213.8e6,1214.3e6) TH2D * GainTime = gain->GetHistGainTime()->Projection(0,1) GainTime->GetXaxis()->SetTimeDisplay(kTRUE) GainTime->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}") GainTime->Draw("colz") // // MakeSlineFit example // AliSplineFit *fit = AliTPCcalibTimeGain::MakeSplineFit(gr) TGraph * grfit = fit.MakeGraph(gr->GetX()[0],gr->GetX()[gr->GetN()-1],50000,0); gr->SetMarkerStyle(25); gr->Draw("lp"); grfit->SetLineColor(2); grfit->Draw("lu"); // // QA - dE/dx resoultion as a function of time //TCa TGraph * grSigma = AliTPCcalibBase::FitSlicesSigma(gain->GetHistGainTime(),0,1,1800,5) TCanvas *c1 = new TCanvas("c1","transparent pad",200,10,700,500); TPad *pad1 = new TPad("pad1","",0,0,1,1); TPad *pad2 = new TPad("pad2","",0,0,1,1); pad2->SetFillStyle(4000); //will be transparent pad1->Draw(); pad1->cd(); GainTime->Draw("colz") gr->Draw("lp") c1->cd(); Double_t ymin = -0.04; Double_t ymax = 0.12; Double_t dy = (ymax-ymin)/0.8; Double_t xmin = GainTime->GetXaxis()->GetXmin() Double_t xmax = GainTime->GetXaxis()->GetXmax() Double_t dx = (xmax-xmin)/0.8; //10 per cent margins left and right pad2->Range(xmin-0.1*dx,ymin-0.1*dy,xmax+0.1*dx,ymax+0.1*dy); pad2->Draw(); pad2->cd(); grSigma->SetLineColor(2); grSigma->SetLineWidth(2); grSigma->Draw("lp") TGaxis *axis = new TGaxis(xmax,ymin,xmax,ymax,ymin,ymax,50510,"+L"); axis->SetLabelColor(kRed); axis->SetTitle("dE/dx resolution #sigma_{dE/dx}"); axis->Draw(); //////////////////////////////////////////////////////////////////////////////////////////////////////////////////////// ----> make Attachment study TFile fcalib("CalibObjects40366b.root"); TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib"); AliTPCcalibTimeGain * gain = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain"); TGraphErrors * grAttach = gain->GetGraphAttachment(2000,4) TCanvas *c1 = new TCanvas("c1","transparent pad",200,10,700,500); TPad *pad1 = new TPad("pad1","",0,0,1,1); TPad *pad2 = new TPad("pad2","",0,0,1,1); pad2->SetFillStyle(4000); //will be transparent pad1->Draw(); pad1->cd(); gain->GetHistGainTime()->GetAxis(1)->SetRangeUser(1213.8e6,1214.3e6) TH2D * GainTime = gain->GetHistGainTime()->Projection(0,1) GainTime->GetXaxis()->SetTimeDisplay(kTRUE) GainTime->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}") GainTime->Draw("colz") //gr->Draw("lp") c1->cd(); Double_t ymin = -0.001; Double_t ymax = 0.001; Double_t dy = (ymax-ymin)/0.8; Double_t xmin = GainTime->GetXaxis()->GetXmin() Double_t xmax = GainTime->GetXaxis()->GetXmax() Double_t dx = (xmax-xmin)/0.8; //10 per cent margins left and right pad2->Range(xmin-0.1*dx,ymin-0.1*dy,xmax+0.1*dx,ymax+0.1*dy); pad2->Draw(); pad2->cd(); grAttach->SetLineColor(2); grAttach->SetLineWidth(2); grAttach->Draw("lp") TGaxis *axis = new TGaxis(xmax,ymin,xmax,ymax,ymin,ymax,50510,"+L"); axis->SetLabelColor(kRed); axis->SetTitle("attachment coefficient b"); axis->Draw(); */ #include "Riostream.h" #include "TChain.h" #include "TTree.h" #include "TH1F.h" #include "TH2F.h" #include "TH3F.h" #include "THnSparse.h" #include "TGraphErrors.h" #include "TList.h" #include "TMath.h" #include "TCanvas.h" #include "TFile.h" #include "TF1.h" #include "TVectorD.h" #include "TProfile.h" #include "TGraphErrors.h" #include "TCanvas.h" #include "AliTPCclusterMI.h" #include "AliTPCseed.h" #include "AliESDVertex.h" #include "AliESDEvent.h" #include "AliESDfriend.h" #include "AliESDInputHandler.h" #include "AliAnalysisManager.h" #include "AliTracker.h" #include "AliMagF.h" #include "AliTPCCalROC.h" #include "AliLog.h" #include "AliTPCcalibTimeGain.h" #include "TTreeStream.h" #include "AliTPCTracklet.h" #include "TTimeStamp.h" #include "AliTPCcalibDB.h" #include "AliTPCcalibLaser.h" #include "AliDCSSensorArray.h" #include "AliDCSSensor.h" ClassImp(AliTPCcalibTimeGain) AliTPCcalibTimeGain::AliTPCcalibTimeGain() :AliTPCcalibBase(), fHistGainTime(0), fGainVsTime(0), fHistDeDxTotal(0), fIntegrationTimeDeDx(0), fMIP(0), fUseMax(0), fLowerTrunc(0), fUpperTrunc(0), fUseShapeNorm(0), fUsePosNorm(0), fUsePadNorm(0), fUseCookAnalytical(0), fIsCosmic(0), fLowMemoryConsumption(0) { AliInfo("Default Constructor"); } AliTPCcalibTimeGain::AliTPCcalibTimeGain(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeGain) :AliTPCcalibBase(), fHistGainTime(0), fGainVsTime(0), fHistDeDxTotal(0), fIntegrationTimeDeDx(0), fMIP(0), fUseMax(0), fLowerTrunc(0), fUpperTrunc(0), fUseShapeNorm(0), fUsePosNorm(0), fUsePadNorm(0), fUseCookAnalytical(0), fIsCosmic(0), fLowMemoryConsumption(0) { SetName(name); SetTitle(title); AliInfo("Non Default Constructor"); fIntegrationTimeDeDx = deltaIntegrationTimeGain; Double_t deltaTime = EndTime - StartTime; // main histogram for time dependence: dE/dx, time, type (1-muon cosmic,2-pion beam data), meanDriftlength, momenta (only filled if enough space is available), run number Int_t timeBins = TMath::Nint(deltaTime/deltaIntegrationTimeGain); Int_t binsGainTime[6] = {150, timeBins, 2, 25, 200, 10000000}; Double_t xminGainTime[6] = {0.5, StartTime, 0.5, 0, 0.1, -0.5}; Double_t xmaxGainTime[6] = { 8, EndTime, 2.5, 250, 50, 9999999.5}; fHistGainTime = new THnSparseF("HistGainTime","dEdx time dep.;dEdx,time,type,driftlength,momenta,run number;dEdx",6,binsGainTime,xminGainTime,xmaxGainTime); BinLogX(fHistGainTime, 4); // fHistDeDxTotal = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",250,0.01,100.,1000,0.,8); BinLogX(fHistDeDxTotal); // default values for dE/dx fMIP = 50.; fUseMax = kTRUE; fLowerTrunc = 0.0; fUpperTrunc = 0.7; fUseShapeNorm = kTRUE; fUsePosNorm = kFALSE; fUsePadNorm = kFALSE; fUseCookAnalytical = kFALSE; // fIsCosmic = kTRUE; fLowMemoryConsumption = kTRUE; // } AliTPCcalibTimeGain::~AliTPCcalibTimeGain(){ // // // } void AliTPCcalibTimeGain::Process(AliESDEvent *event) { // // main track loop // if (!event) { Printf("ERROR: ESD not available"); return; } if (fIsCosmic) { // this should be removed at some point based on trigger mask !? ProcessCosmicEvent(event); } else { ProcessBeamEvent(event); } } void AliTPCcalibTimeGain::ProcessCosmicEvent(AliESDEvent *event) { AliESDfriend *ESDfriend=static_cast(event->FindListObject("AliESDfriend")); if (!ESDfriend) { Printf("ERROR: ESDfriend not available"); return; } // UInt_t time = event->GetTimeStamp(); Int_t ntracks = event->GetNumberOfTracks(); Int_t runNumber = event->GetRunNumber(); // // track loop // for (Int_t i=0;iGetTrack(i); if (!track) continue; const AliExternalTrackParam * trackIn = track->GetInnerParam(); const AliExternalTrackParam * trackOut = track->GetOuterParam(); if (!trackIn) continue; if (!trackOut) continue; // calculate necessary track parameters Double_t meanP = trackIn->GetP(); Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ()); Int_t NclsDeDx = track->GetTPCNcls(); // exclude tracks which do not look like primaries or are simply too short or on wrong sectors if (NclsDeDx < 60) continue; if (TMath::Abs(trackIn->GetTgl()) > 1) continue; if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue; // Get seeds AliESDfriendTrack *friendTrack = ESDfriend->GetTrack(i); if (!friendTrack) continue; TObject *calibObject; AliTPCseed *seed = 0; for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) { if ((seed=dynamic_cast(calibObject))) break; } if (seed) { Double_t TPCsignal = GetTPCdEdx(seed); if (NclsDeDx > 100) fHistDeDxTotal->Fill(meanP, TPCsignal); // if (fLowMemoryConsumption) { if (meanP < 20) continue; meanP = 30; // set all momenta to one in order to save memory } //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta Double_t vec[6] = {TPCsignal,time,1,meanDrift,meanP,runNumber}; fHistGainTime->Fill(vec); } } } void AliTPCcalibTimeGain::ProcessBeamEvent(AliESDEvent *event) { AliESDfriend *ESDfriend=static_cast(event->FindListObject("AliESDfriend")); if (!ESDfriend) { Printf("ERROR: ESDfriend not available"); return; } // UInt_t time = event->GetTimeStamp(); Int_t ntracks = event->GetNumberOfTracks(); Int_t runNumber = event->GetRunNumber(); // // track loop // for (Int_t i=0;iGetTrack(i); if (!track) continue; const AliExternalTrackParam * trackIn = track->GetInnerParam(); const AliExternalTrackParam * trackOut = track->GetOuterParam(); if (!trackIn) continue; if (!trackOut) continue; // calculate necessary track parameters Double_t meanP = trackIn->GetP(); Double_t meanDrift = 250 - 0.5*TMath::Abs(trackIn->GetZ() + trackOut->GetZ()); Int_t NclsDeDx = track->GetTPCNcls(); // exclude tracks which do not look like primaries or are simply too short or on wrong sectors if (NclsDeDx < 60) continue; if (TMath::Abs(trackIn->GetTgl()) > 1) continue; if (TMath::Abs(trackIn->GetSnp()) > 0.6) continue; // Get seeds AliESDfriendTrack *friendTrack = ESDfriend->GetTrack(i); if (!friendTrack) continue; TObject *calibObject; AliTPCseed *seed = 0; for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) { if ((seed=dynamic_cast(calibObject))) break; } if (seed) { Double_t TPCsignal = GetTPCdEdx(seed); fHistDeDxTotal->Fill(meanP, TPCsignal); // if (fLowMemoryConsumption) { if (meanP > 0.5 || meanP < 0.4) continue; meanP = 0.45; // set all momenta to one in order to save memory } //dE/dx, time, type (1-muon cosmic,2-pion beam data), momenta Double_t vec[6] = {TPCsignal,time,2,meanDrift,meanP,runNumber}; fHistGainTime->Fill(vec); } } } Float_t AliTPCcalibTimeGain::GetTPCdEdx(AliTPCseed * seed) { Double_t signal = 0; // if (!fUseCookAnalytical) { signal = (1/fMIP)*seed->CookdEdxNorm(fLowerTrunc,fUpperTrunc,fUseMax,0,159,fUseShapeNorm,fUsePosNorm,fUsePadNorm,0); } else { signal = (1/fMIP)*seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,fUseMax); } // return signal; } void AliTPCcalibTimeGain::AnalyzeRun(Int_t minEntries) { // // // if (fIsCosmic) { fHistGainTime->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics fHistGainTime->GetAxis(4)->SetRangeUser(20,100); // only Fermi-Plateau muons } else { fHistGainTime->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data fHistGainTime->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions } // fGainVsTime = AliTPCcalibBase::FitSlices(fHistGainTime,0,1,minEntries,10); // return; } TGraphErrors * AliTPCcalibTimeGain::GetGraphGainVsTime(Int_t runNumber, Int_t minEntries) { // // // if (runNumber == 0) { if (!fGainVsTime) { AnalyzeRun(minEntries); } } else { // 1st check if the current run was cosmic or beam event fHistGainTime->GetAxis(5)->SetRangeUser(runNumber,runNumber); AnalyzeRun(minEntries); } if (fGainVsTime->GetN() == 0) return 0; return fGainVsTime; } Long64_t AliTPCcalibTimeGain::Merge(TCollection *li) { TIterator* iter = li->MakeIterator(); AliTPCcalibTimeGain* cal = 0; while ((cal = (AliTPCcalibTimeGain*)iter->Next())) { if (!cal->InheritsFrom(AliTPCcalibTimeGain::Class())) { Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName()); return -1; } // add histograms here... if (cal->GetHistGainTime()) fHistGainTime->Add(cal->GetHistGainTime()); if (cal->GetHistDeDxTotal()) fHistDeDxTotal->Add(cal->GetHistDeDxTotal()); } return 0; } AliSplineFit * AliTPCcalibTimeGain::MakeSplineFit(TGraphErrors * graph) { // // // AliSplineFit *fit = new AliSplineFit(); fit->SetGraph(graph); fit->SetMinPoints(graph->GetN()+1); fit->InitKnots(graph,2,0,0.001); fit->SplineFit(0); return fit; } TGraphErrors * AliTPCcalibTimeGain::GetGraphAttachment(Int_t minEntries, Int_t nmaxBin, Float_t /*fracLow*/, Float_t /*fracUp*/) { // // For each time bin the driftlength-dependence of the signal is fitted. // TH3D * hist = fHistGainTime->Projection(1, 0, 3); Double_t *xvec = new Double_t[hist->GetNbinsX()]; Double_t *yvec = new Double_t[hist->GetNbinsX()]; Double_t *xerr = new Double_t[hist->GetNbinsX()]; Double_t *yerr = new Double_t[hist->GetNbinsX()]; Int_t counter = 0; TH2D * projectionHist = 0x0; // for(Int_t i=1; i < hist->GetNbinsX(); i++) { Int_t nsum=0; Int_t imin = i; Int_t imax = i; for (Int_t idelta=0; ideltaGetNbinsX()); nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1)); if (nsum==0) break; if (nsum>minEntries) break; } if (nsumGetXaxis()->SetRange(imin,imax); projectionHist = (TH2D*)hist->Project3D("yzNUFNOF"); // TObjArray arr; projectionHist->FitSlicesY(0,2, projectionHist->GetNbinsX()-2,0,"QNR",&arr); TH1D * histAttach = (TH1D*)arr.At(1); TF1 pol1("polynom1","pol1",10,240); histAttach->Fit(&pol1,"QNR"); xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin)+hist->GetXaxis()->GetBinCenter(imax)); yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0); xerr[counter] = 0; yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0); counter++; // delete projectionHist; } TGraphErrors * graphErrors = new TGraphErrors(counter, xvec, yvec, xerr, yerr); delete [] xvec; delete [] yvec; delete [] xerr; delete [] yerr; delete hist; return graphErrors; } void AliTPCcalibTimeGain::BinLogX(THnSparse *h, Int_t axisDim) { // Method for the correct logarithmic binning of histograms TAxis *axis = h->GetAxis(axisDim); int bins = axis->GetNbins(); Double_t from = axis->GetXmin(); Double_t to = axis->GetXmax(); Double_t *new_bins = new Double_t[bins + 1]; new_bins[0] = from; Double_t factor = pow(to/from, 1./bins); for (int i = 1; i <= bins; i++) { new_bins[i] = factor * new_bins[i-1]; } axis->Set(bins, new_bins); delete new_bins; } void AliTPCcalibTimeGain::BinLogX(TH1 *h) { // Method for the correct logarithmic binning of histograms TAxis *axis = h->GetXaxis(); int bins = axis->GetNbins(); Double_t from = axis->GetXmin(); Double_t to = axis->GetXmax(); Double_t *new_bins = new Double_t[bins + 1]; new_bins[0] = from; Double_t factor = pow(to/from, 1./bins); for (int i = 1; i <= bins; i++) { new_bins[i] = factor * new_bins[i-1]; } axis->Set(bins, new_bins); delete new_bins; } void AliTPCcalibTimeGain::CalculateBetheAlephParams(TH2F *hist, Double_t * ini) { //{0.0762*MIP,10.632,1.34e-05,1.863,1.948} const Double_t sigma = 0.06; 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()); BinLogX(histBG); TF1 *foElectron = new TF1("foElectron", "AliExternalTrackParam::BetheBlochAleph(x/0.000511,[0],[1],[2],[3],[4])",0.01,100); TF1 *foMuon = new TF1("foMuon", "AliExternalTrackParam::BetheBlochAleph(x/0.1056,[0],[1],[2],[3],[4])",0.01,100); TF1 *foPion = new TF1("foPion", "AliExternalTrackParam::BetheBlochAleph(x/0.138,[0],[1],[2],[3],[4])",0.01,100); TF1 *foKaon = new TF1("foKaon", "AliExternalTrackParam::BetheBlochAleph(x/0.498,[0],[1],[2],[3],[4])",0.01,100); TF1 *foProton = new TF1("foProton", "AliExternalTrackParam::BetheBlochAleph(x/0.938,[0],[1],[2],[3],[4])",0.01,100); foElectron->SetParameters(ini); foMuon->SetParameters(ini); foPion->SetParameters(ini); foKaon->SetParameters(ini); foProton->SetParameters(ini); TCanvas *CanvCheck1 = new TCanvas(); hist->Draw("colz"); foElectron->Draw("same"); foMuon->Draw("same"); foPion->Draw("same"); foKaon->Draw("same"); foProton->Draw("same"); // Loop over all points of the input histogram for(Int_t i=1; i < hist->GetNbinsX(); i++) { Double_t x = hist->GetXaxis()->GetBinCenter(i); for(Int_t j=1; j < hist->GetNbinsY(); j++) { Long64_t n = hist->GetBin(i, j); Double_t y = hist->GetYaxis()->GetBinCenter(j); Double_t entries = hist->GetBinContent(n); Double_t mass = 0; // 1. identify protons mass = 0.938; if (TMath::Abs(y - foProton->Eval(x))/foProton->Eval(x) < 4*sigma && x < 0.65 ) { for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y); } // 2. identify electrons mass = 0.000511; if (fIsCosmic) { if (TMath::Abs(y - foElectron->Eval(x))/foElectron->Eval(x) < 4*sigma && (x>0.25&&x<0.7) && fIsCosmic) { for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y); } } else { 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))) { for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y); } } // 3. identify either muons or pions depending on cosmic or not if (fIsCosmic) { mass = 0.1056; if (TMath::Abs(y - foMuon->Eval(x))/foMuon->Eval(x) < 4*sigma && x > 0.25 && x < 100 ) { for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y); } } else { mass = 0.1396; if (TMath::Abs(y - foPion->Eval(x))/foPion->Eval(x) < 4*sigma && x > 0.15 && x < 2) { for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y); } } // 4. for pp also Kaons must be included if (!fIsCosmic) { mass = 0.4936; if (TMath::Abs(y - foKaon->Eval(x))/foKaon->Eval(x) < 4*sigma && x > 0.2 && x < 0.3 ) { for(Int_t iEntries =0; iEntries < entries; iEntries++) histBG->Fill(x/mass, y); } } } } // Fit Aleph-Parameters to the obtained profile TF1 * funcAlephD = new TF1("AlephParametrizationD", "AliExternalTrackParam::BetheBlochAleph(x,[0],[1],[2],[3],[4])",0.3,10000); funcAlephD->SetParameters(ini); TCanvas *CanvCheck2 = new TCanvas(); histBG->Draw(); //FitSlices TObjArray * arr = new TObjArray(); histBG->FitSlicesY(0,0,-1,0,"QN",arr); TH1D * fitMean = (TH1D*) arr->At(1); fitMean->Draw("same"); funcAlephD->SetParLimits(2,1e-3,1e-7); funcAlephD->SetParLimits(3,0.5,3.5); funcAlephD->SetParLimits(4,0.5,3.5); fitMean->Fit(funcAlephD, "QNR"); funcAlephD->Draw("same"); for(Int_t i=0;i<5;i++) ini[i] = funcAlephD->GetParameter(i); foElectron->SetParameters(ini); foMuon->SetParameters(ini); foPion->SetParameters(ini); foKaon->SetParameters(ini); foProton->SetParameters(ini); TCanvas *CanvCheck3 = new TCanvas(); hist->Draw("colz"); foElectron->Draw("same"); foMuon->Draw("same"); foPion->Draw("same"); foKaon->Draw("same"); foProton->Draw("same"); CanvCheck1->Print(); CanvCheck2->Print(); CanvCheck3->Print(); return; }