#include #include #include #include #include #include #include #include #include #include #include #include "AliTRDcluster.h" #include "AliESDHeader.h" #include "AliESDRun.h" #include "AliTRDseedV1.h" #include "AliTRDtrackV1.h" #include "AliTrackReference.h" #include "TTreeStream.h" #include "AliTRDtrackInfo/AliTRDtrackInfo.h" #include "AliTRDtrackInfo/AliTRDeventInfo.h" #include "AliTRDcheckDetector.h" #include #include //////////////////////////////////////////////////////////////////////////// // // // Reconstruction QA // // // // Task doing basic checks for tracking and detector performance // // // // Authors: // // Anton Andronic // // Markus Fasel // // // //////////////////////////////////////////////////////////////////////////// //_______________________________________________________ AliTRDcheckDetector::AliTRDcheckDetector(): AliTRDrecoTask("DetChecker", "Basic Detector Checker") ,fEventInfo(0x0) ,fTriggerNames(0x0) { // // Default constructor // DefineInput(1,AliTRDeventInfo::Class()); } //_______________________________________________________ AliTRDcheckDetector::~AliTRDcheckDetector(){ // // Destructor // if(fEventInfo) delete fEventInfo; if(fTriggerNames) delete fTriggerNames; } //_______________________________________________________ void AliTRDcheckDetector::ConnectInputData(Option_t *opt){ // // Connect the Input data with the task // AliTRDrecoTask::ConnectInputData(opt); fEventInfo = dynamic_cast(GetInputData(1)); } //_______________________________________________________ void AliTRDcheckDetector::CreateOutputObjects(){ // // Create Output Objects // fContainer = new TObjArray(25); // Register Histograms fContainer->Add(new TH1F("hNtrks", "Number of Tracks per event", 100, 0, 100)); fContainer->Add(new TH1F("hEventsTriggerTracks", "Trigger Class (Tracks)", 100, 0, 100)); fContainer->Add(new TH1F("hNcls", "Nr. of clusters per track", 181, -0.5, 180.5)); fContainer->Add(new TH1F("hNtls", "Nr. tracklets per track", 7, -0.5, 6.5)); fContainer->Add(new TH1F("hNclTls","Mean Number of clusters per tracklet", 31, -0.5, 30.5)); fContainer->Add(new TH1F("hChi2", "Chi2", 200, 0, 20)); fContainer->Add(new TH1F("hChi2n", "Norm. Chi2 (tracklets)", 50, 0, 5)); fContainer->Add(new TH1F("hSM", "Track Counts in Supermodule", 18, -0.5, 17.5)); // Detector signal on Detector-by-Detector basis fContainer->Add(new TProfile("hPHdetector", "Average PH", 31, -0.5, 30.5)); fContainer->Add(new TH1F("hQclDetector", "Cluster charge", 200, 0, 1200)); fContainer->Add(new TH1F("hQTdetector", "Total Charge Deposit", 6000, 0, 6000)); fContainer->Add(new TH1F("hEventsTrigger", "Trigger Class", 100, 0, 100)); fTriggerNames = new TMap(); } //_______________________________________________________ void AliTRDcheckDetector::Exec(Option_t *){ // // Execution function // Filling TRD quality histos // if(!HasMCdata() && fEventInfo->GetEventHeader()->GetEventType() != 7) return; // For real data we select only physical events Int_t nTracks = 0; // Count the number of tracks per event Int_t triggermask = fEventInfo->GetEventHeader()->GetTriggerMask(); TString triggername = fEventInfo->GetRunInfo()->GetFiredTriggerClasses(triggermask); if(fDebugLevel > 6)printf("Trigger cluster: %d, Trigger class: %s\n", triggermask, triggername.Data()); dynamic_cast(fContainer->UncheckedAt(kNEventsTrigger))->Fill(triggermask); AliTRDtrackInfo *fTrackInfo = 0x0; AliTRDtrackV1 *fTRDtrack = 0x0; AliTRDseedV1 *fTracklet = 0x0; AliTRDcluster *fTRDcluster = 0x0; Float_t momentum = 0.; // momentum information needed for systematic studies Int_t pdg = 0; Float_t theta = 0., phi = 0.; for(Int_t iti = 0; iti < fTracks->GetEntriesFast(); iti++){ fTrackInfo = dynamic_cast(fTracks->UncheckedAt(iti)); if(!fTrackInfo || !(fTRDtrack = fTrackInfo->GetTrack())) continue; Int_t nclusters = fTRDtrack->GetNumberOfClusters(); Int_t ntracklets = fTRDtrack->GetNumberOfTracklets(); Float_t chi2 = fTRDtrack->GetChi2(); // Fill Histograms dynamic_cast(fContainer->UncheckedAt(kNclustersHist))->Fill(nclusters); dynamic_cast(fContainer->UncheckedAt(kNtrackletsHist))->Fill(ntracklets); dynamic_cast(fContainer->UncheckedAt(kChi2))->Fill(chi2); dynamic_cast(fContainer->UncheckedAt(kChi2Normalized))->Fill(chi2/static_cast(ntracklets)); // now loop over single tracklets momentum = 0.; pdg = 0; if(HasMCdata()){ AliTrackReference *fRef = 0x0; Int_t iti = 0; while(!(fRef = fTrackInfo->GetTrackRef(iti++)) && (iti <=12)); if(fRef) momentum = fRef->P(); pdg = fTrackInfo->GetPDG(); } if(ntracklets > 2){ Int_t sector = -1; for(Int_t ilayer = 0; ilayer < kNLayers; ilayer++){ if(!(fTracklet = fTRDtrack->GetTracklet(ilayer)) || !fTracklet->IsOK()) continue; Float_t nClustersTracklet = fTracklet->GetN(); dynamic_cast(fContainer->UncheckedAt(kNclusterTrackletHist))->Fill(nClustersTracklet); if(nClustersTracklet){ Float_t Qtot = 0; Int_t detector = -1; theta = TMath::ATan(fTracklet->GetZfit(1)); phi = TMath::ATan(fTracklet->GetYfit(1)); for(Int_t itb = 0; itb < kNTimebins; itb++){ if(!(fTRDcluster = fTracklet->GetClusters(itb))) continue; Int_t localtime = fTRDcluster->GetLocalTimeBin(); Double_t clusterCharge = fTRDcluster->GetQ(); detector = fTRDcluster->GetDetector(); sector = static_cast(detector/kNDetectorsSector); Double_t absolute_charge = TMath::Abs(clusterCharge); Qtot += absolute_charge; dynamic_cast(fContainer->UncheckedAt(kPulseHeight))->Fill(localtime, absolute_charge); dynamic_cast(fContainer->UncheckedAt(kClusterCharge))->Fill(absolute_charge); if(fDebugLevel > 2){ (*fDebugStream) << "PulseHeight" << "Detector=" << detector << "Sector=" << sector << "Timebin=" << localtime << "Charge=" << absolute_charge << "momentum=" << momentum << "pdg=" << pdg << "theta=" << theta << "phi=" << phi << "\n"; } } dynamic_cast(fContainer->UncheckedAt(kChargeDeposit))->Fill(Qtot); Int_t crossing = fTracklet->GetNChange(); if(fDebugLevel > 3){ (*fDebugStream) << "ChargeDeposit" << "Detector=" << detector << "Sector=" << sector << "nclusters=" << nClustersTracklet << "crossing=" << crossing << "QT=" << Qtot << "momentum=" << momentum << "pdg=" << pdg << "theta=" << theta << "phi=" << phi << "\n"; } } } dynamic_cast(fContainer->UncheckedAt(kNTracksSectorHist))->Fill(sector); } nTracks++; } if(nTracks){ dynamic_cast(fContainer->UncheckedAt(kNEventsTriggerTracks))->Fill(triggermask); dynamic_cast(fContainer->UncheckedAt(kNTracksEventHist))->Fill(nTracks); } if(triggermask <= 20 && !fTriggerNames->FindObject(Form("%d", triggermask))){ fTriggerNames->Add(new TObjString(Form("%d", triggermask)), new TObjString(triggername)); // also set the label for both histograms TH1 *histo = dynamic_cast(fContainer->UncheckedAt(kNEventsTriggerTracks)); histo->GetXaxis()->SetBinLabel(histo->FindBin(triggermask), triggername); histo = dynamic_cast(fContainer->UncheckedAt(kNEventsTrigger)); histo->GetXaxis()->SetBinLabel(histo->FindBin(triggermask), triggername); } PostData(0, fContainer); } //_______________________________________________________ void AliTRDcheckDetector::Terminate(Option_t *){ // // Terminate function // } //_______________________________________________________ Bool_t AliTRDcheckDetector::PostProcess(){ // // Do Postprocessing (for the moment set the number of Reference histograms) // TH1 * histo = 0x0; histo = dynamic_cast(fContainer->UncheckedAt(kNTracksEventHist)); histo->GetXaxis()->SetTitle("Number of Tracks"); histo->GetYaxis()->SetTitle("Events"); histo = dynamic_cast(fContainer->UncheckedAt(kNclustersHist)); histo->GetXaxis()->SetTitle("Number of Clusters"); histo->GetYaxis()->SetTitle("Entries"); histo = dynamic_cast(fContainer->UncheckedAt(kNtrackletsHist)); histo->GetXaxis()->SetTitle("Number of Tracklets"); histo->GetYaxis()->SetTitle("Entries"); histo = dynamic_cast(fContainer->UncheckedAt(kNclusterTrackletHist)); histo->GetXaxis()->SetTitle("Number of Clusters"); histo->GetYaxis()->SetTitle("Entries"); histo = dynamic_cast(fContainer->UncheckedAt(kChi2)); histo->GetXaxis()->SetTitle("#chi^2"); histo->GetYaxis()->SetTitle("Entries"); histo = dynamic_cast(fContainer->UncheckedAt(kNTracksSectorHist)); histo->GetXaxis()->SetTitle("Sector"); histo->GetYaxis()->SetTitle("Number of Tracks"); histo = dynamic_cast(fContainer->UncheckedAt(kPulseHeight)); histo->GetXaxis()->SetTitle("Time / 100ns"); histo->GetYaxis()->SetTitle("Average Pulse Height (a. u.)"); histo = dynamic_cast(fContainer->UncheckedAt(kClusterCharge)); histo->GetXaxis()->SetTitle("Cluster Charge (a.u.)"); histo->GetYaxis()->SetTitle("Entries"); histo = dynamic_cast(fContainer->UncheckedAt(kChargeDeposit)); histo->GetXaxis()->SetTitle("Charge Deposit (a.u.)"); histo->GetYaxis()->SetTitle("Entries"); // Calculate the purity of the trigger clusters histo = dynamic_cast(fContainer->UncheckedAt(kNEventsTrigger)); TH1F *histoTracks = dynamic_cast(fContainer->UncheckedAt(kNEventsTriggerTracks)); histoTracks->Divide(histo); Float_t purities[20], val = 0; TString triggernames[20]; Int_t nTriggerClasses = 0; for(Int_t ibin = 1; ibin <= histo->GetNbinsX(); ibin++){ if((val = histoTracks->GetBinContent(ibin))){ purities[nTriggerClasses] = val; triggernames[nTriggerClasses] = histoTracks->GetXaxis()->GetBinLabel(ibin); nTriggerClasses++; } } TH1F *hTriggerInf = new TH1F("fTriggerInf", "Trigger Information", TMath::Max(nTriggerClasses, 1), 0, TMath::Max(nTriggerClasses, 1)); for(Int_t ibin = 1; ibin <= nTriggerClasses; ibin++){ hTriggerInf->SetBinContent(ibin, purities[ibin-1]); hTriggerInf->GetXaxis()->SetBinLabel(ibin, triggernames[ibin-1].Data()); } hTriggerInf->GetXaxis()->SetTitle("Trigger Cluster"); hTriggerInf->GetYaxis()->SetTitle("Ratio"); hTriggerInf->GetYaxis()->SetRangeUser(0,1); // hTriggerInf->SetMarkerColor(kBlue); // hTriggerInf->SetMarkerStyle(22); fContainer->Add(hTriggerInf); fNRefFigures = 10; return kTRUE; } //_______________________________________________________ void AliTRDcheckDetector::GetRefFigure(Int_t ifig){ // // Setting Reference Figures // TH1 *h = 0x0; switch(ifig){ case 0: ((TH1F*)fContainer->At(kNTracksEventHist))->Draw("pl"); break; case 1: ((TH1F*)fContainer->At(kNclustersHist))->Draw("pl"); break; case 2: h = (TH1F*)fContainer->At(kNtrackletsHist); if(!h->GetEntries()) break; h->Scale(100./h->Integral()); h->GetXaxis()->SetRangeUser(.5, 6.5); h->SetFillColor(kGreen); h->SetBarOffset(.2); h->SetBarWidth(.6); h->Draw("bar1"); break; case 3: ((TH1F*)fContainer->At(kNclusterTrackletHist))->Draw("pc"); break; case 4: ((TH1F*)fContainer->At(kChi2))->Draw(""); break; case 5: h = (TH1F*)fContainer->At(kNTracksSectorHist); if(!h->GetEntries()) break; h->Scale(100./h->Integral()); h->SetFillColor(kGreen); h->SetBarOffset(.2); h->SetBarWidth(.6); h->Draw("bar1"); break; case 6: h = (TH1F*)fContainer->At(kPulseHeight); h->SetMarkerStyle(24); h->Draw("e1"); break; case 7: ((TH1F*)fContainer->At(kClusterCharge))->Draw("c"); break; case 8: ((TH1F*)fContainer->At(kChargeDeposit))->Draw("c"); break; case 9: h=(TH1F*)fContainer->At(kPurity); h->SetBarOffset(.2); h->SetBarWidth(.6); h->SetFillColor(kGreen); h->Draw("bar1"); break; default: ((TH1F*)fContainer->At(kNTracksEventHist))->Draw("pl"); break; } }