/************************************************************************** * 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. * **************************************************************************/ /* $Id: AliTRDqaEnergyDeposit.cxx $ */ // // This class is a part of a package of high level QA monitoring for TRD. // In this class the dEdX is analyzed as a function of the particle type, // momentum and chamber. // // S. Radomski // radomski@physi.uni-heidelberg.de // March 2008 // #include "AliTRDqaEnergyDeposit.h" #include "TMath.h" #include "TH1D.h" #include "TH2D.h" #include "TFile.h" #include "TChain.h" #include "AliESDEvent.h" #include "AliESDtrack.h" //______________________________________________________________________________ AliTRDqaEnergyDeposit::AliTRDqaEnergyDeposit() : AliAnalysisTask("",""), fChain(0), fESD(0), fOutputContainer(0) { // // Default constructor // for (Int_t i = 0; i < 2; i++) { fSignalPtSum[i] = 0x0; } for (Int_t j = 0; j < 2*5; j++) { fSignalPtType[j] = 0x0; fSignalPtPure[j] = 0x0; fProb[j] = 0x0; } } //______________________________________________________________________________ AliTRDqaEnergyDeposit:: AliTRDqaEnergyDeposit(const AliTRDqaEnergyDeposit & /*trd*/) : AliAnalysisTask("",""), fChain(0), fESD(0), fOutputContainer(0) { // // Copy constructor // for (Int_t i = 0; i < 2; i++) { fSignalPtSum[i] = 0x0; } for (Int_t j = 0; j < 2*5; j++) { fSignalPtType[j] = 0x0; fSignalPtPure[j] = 0x0; fProb[j] = 0x0; } } //______________________________________________________________________________ AliTRDqaEnergyDeposit::AliTRDqaEnergyDeposit(const char *name) : AliAnalysisTask(name,""), fChain(0), fESD(0), fOutputContainer(0) { // // Constructor. // // Input slot #0 works with an Ntuple DefineInput(0, TChain::Class()); // Output slot #0 writes into a TH1 container DefineOutput(0, TObjArray::Class()) ; for (Int_t i = 0; i < 2; i++) { fSignalPtSum[i] = 0x0; } for (Int_t j = 0; j < 2*5; j++) { fSignalPtType[j] = 0x0; fSignalPtPure[j] = 0x0; fProb[j] = 0x0; } } //______________________________________________________________________________ void AliTRDqaEnergyDeposit::ConnectInputData(const Option_t *) { // Initialisation of branch container and histograms //AliInfo(Form("*** Initialization of %s", GetName())) ; fChain = (TChain*)GetInputData(0); fESD = new AliESDEvent(); fESD->ReadFromTree(fChain); } //________________________________________________________________________ void AliTRDqaEnergyDeposit::CreateOutputObjects() { // build histograms // prepare the scale from 0.5 to 10 GeV const Int_t knbinsx = 50; Double_t scalex[knbinsx+1]; Double_t dd = (TMath::Log(10) - TMath::Log(0.5)) / knbinsx; for(Int_t ix=0; ixAddAt(fSignalPtSum[i], c++); for(Int_t j=0; jAddAt(fSignalPtType[idx], c++); // fSignalPtPure[idx] = new TH2D(Form("ptSigPure%s%d", charge[i], j), title, knbinsx, scalex, knbinsy, scaley); fOutputContainer->AddAt(fSignalPtPure[idx], c++); fProb[idx] = new TH1D(Form("prob%s%d", charge[i], j), ";LQ", 100, 0, 1); fOutputContainer->AddAt(fProb[idx], c++); } } printf("n hist = %d\n", c); } //______________________________________________________________________________ void AliTRDqaEnergyDeposit::Exec(Option_t *) { // // Process one event // // Long64_t entry = fChain->GetReadEntry() ; // Processing of one event if (!fESD) { //AliError("fESD is not connected to the input!") ; return ; } Int_t nTracks = fESD->GetNumberOfTracks(); //fNTracks->Fill(nTracks); // track loop for(Int_t i=0; iGetTrack(i); const AliExternalTrackParam *paramOut = track->GetOuterParam(); const AliExternalTrackParam *paramIn = track->GetInnerParam(); // long track .. if (!paramIn) continue; if (!paramOut) continue; UInt_t status = track->GetStatus(); if (!(status & AliESDtrack::kTRDrefit)) continue; if (!(status & AliESDtrack::kTRDpid)) continue; if (track->GetTRDntrackletsPID() < 6) continue; // standard selection Int_t idx = (track->GetSign() > 0) ? 0 : 1; Double_t pt = paramOut->Pt(); Double_t signal = 0; for(Int_t k=0; k<6; ++k) signal += track->GetTRDslice(k, -1); signal /= 6; fSignalPtSum[idx]->Fill(pt, signal); for(Int_t k=0; kGetTRDpid(k); fProb[AliPID::kSPECIES*idx+k]->Fill(lq); if (lq > 0.8) fSignalPtType[AliPID::kSPECIES*idx+k]->Fill(pt, signal); } } PostData(0, fOutputContainer); } //______________________________________________________________________________ void AliTRDqaEnergyDeposit::FillElectrons() const { // // Fill the electrons // } //______________________________________________________________________________ void AliTRDqaEnergyDeposit::Terminate(Option_t *) { // // Retrieve histograms // fOutputContainer = (TObjArray*)GetOutputData(0); // post processing // save the results TFile *file = new TFile("outEnergyDeposit.root", "RECREATE"); fOutputContainer->Write(); file->Flush(); file->Close(); delete file; //for(Int_t i=0; iGetEntries(); i++) { // TObject *obj = fOu // } } //______________________________________________________________________________