4 #include "TEventList.h"
9 #include "AliESDtrack.h"
10 #include "AliTrackReference.h"
12 #include "AliTRDReconstructor.h"
13 #include "AliTRDtrackV1.h"
14 #include "AliTRDseedV1.h"
15 #include "AliTRDpidRefMaker.h"
16 #include "info/AliTRDv0Info.h"
19 // Defines and implements basic functionality for building reference data for TRD PID.
21 // Here is the list of functionality provided by this class
26 // Alex Bercuci <A.Bercuci@gsi.de>
27 // Alex Wilk <wilka@uni-muenster.de>
28 // Markus Fasel <mfasel@gsi.de>
29 // Markus Heide <mheide@uni-muenster.de>
32 ClassImp(AliTRDpidRefMaker)
33 ClassImp(AliTRDpidRefMaker::AliTRDpidRefData)
34 ClassImp(AliTRDpidRefMaker::AliTRDpidRefDataArray)
36 //________________________________________________________________________
37 AliTRDpidRefMaker::AliTRDpidRefDataArray::AliTRDpidRefDataArray() :
41 // Constructor of data array
42 fData = new AliTRDpidRefData[AliTRDgeometry::kNlayer];
45 //________________________________________________________________________
46 AliTRDpidRefMaker::AliTRDpidRefDataArray::~AliTRDpidRefDataArray()
52 //________________________________________________________________________
53 void AliTRDpidRefMaker::AliTRDpidRefDataArray::PushBack(Int_t ly, Int_t p, Float_t *dedx)
55 // Add PID data to the end of the array
56 fData[fNtracklets].fPLbin= (ly<<4) | (p&0xf);
57 memcpy(fData[fNtracklets].fdEdx, dedx, 8*sizeof(Float_t));
61 //________________________________________________________________________
62 void AliTRDpidRefMaker::AliTRDpidRefDataArray::Reset()
66 if(!fNtracklets) return;
68 fData[fNtracklets].fPLbin = 0xff;
69 memset(fData[fNtracklets].fdEdx, 0, 8*sizeof(Float_t));
74 //________________________________________________________________________
75 AliTRDpidRefMaker::AliTRDpidRefMaker(const char *name, const char *title)
76 :AliTRDrecoTask(name, title)
89 // Default constructor
92 fReconstructor = new AliTRDReconstructor();
93 fReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
94 memset(fdEdx, 0, 10*sizeof(Float_t));
95 memset(fPID, 0, AliPID::kSPECIES*sizeof(Float_t));
97 DefineInput(1, TObjArray::Class());
98 DefineOutput(1, TTree::Class());
102 //________________________________________________________________________
103 AliTRDpidRefMaker::~AliTRDpidRefMaker()
105 if(fPIDdataArray) delete fPIDdataArray;
106 if(fReconstructor) delete fReconstructor;
109 //________________________________________________________________________
110 Bool_t AliTRDpidRefMaker::CheckQuality(AliTRDseedV1* /*trklt*/)
112 // Place holder for checking tracklet quality for PID.
117 //________________________________________________________________________
118 Float_t* AliTRDpidRefMaker::CookdEdx(AliTRDseedV1 *trklt)
120 trklt->CookdEdx(AliTRDpidUtil::kNNslices);
121 memcpy(fdEdx, trklt->GetdEdx(), AliTRDpidUtil::kNNslices*sizeof(Float_t));
126 //________________________________________________________________________
127 void AliTRDpidRefMaker::ConnectInputData(Option_t *opt)
129 AliTRDrecoTask::ConnectInputData(opt);
130 fV0s = dynamic_cast<TObjArray*>(GetInputData(1));
133 //________________________________________________________________________
134 void AliTRDpidRefMaker::CreateOutputObjects()
139 OpenFile(0, "RECREATE");
140 fContainer = new TObjArray();
141 fContainer->SetName(Form("Moni%s", GetName()));
143 TH2 *h2 = new TH2I("hPDG","Particle abundance", AliPID::kSPECIES, -0.5, 4.5, AliTRDCalPID::kNMom, -0.5, AliTRDCalPID::kNMom-0.5);
144 TAxis *ax = h2->GetXaxis();
145 ax->SetNdivisions(505);
146 ax->SetTitle("Particle species");
147 for(Int_t is=AliPID::kSPECIES; is--;) ax->SetBinLabel(is+1, AliPID::ParticleShortName(is));
148 h2->GetYaxis()->SetTitle("P bins");
149 h2->GetYaxis()->SetNdivisions(511);
151 fContainer->AddAt(h2, 0);
153 fData = new TTree(GetName(), Form("Reference data for %s", GetName()));
154 fData->Branch("s", &fPIDbin, "s/b");
155 fData->Branch("data", &fPIDdataArray);
158 //________________________________________________________________________
159 void AliTRDpidRefMaker::Exec(Option_t *)
162 // Called for each event
164 AliTRDtrackInfo *track = 0x0;
165 AliTRDtrackV1 *trackTRD = 0x0;
166 AliTrackReference *ref = 0x0;
167 const AliTRDtrackInfo::AliESDinfo *infoESD = 0x0;
168 for(Int_t itrk=0; itrk<fTracks->GetEntriesFast(); itrk++){
169 track = (AliTRDtrackInfo*)fTracks->UncheckedAt(itrk);
170 if(!track->HasESDtrack()) continue;
171 trackTRD = track->GetTrack();
172 infoESD = track->GetESDinfo();
173 Double32_t *infoPID = infoESD->GetSliceIter();
174 Int_t n = infoESD->GetNSlices() - AliTRDgeometry::kNlayer;
175 Double32_t *p = &infoPID[n];
178 ULong_t status = track->GetStatus();
179 if(!(status&AliESDtrack::kTPCout)) continue;
181 // fill the pid information
182 SetRefPID(fRefPID, track, fPID);
184 // prepare PID data array
186 fPIDdataArray = new AliTRDpidRefDataArray();
187 } else fPIDdataArray->Reset();
189 // fill PID information
190 for(Int_t ily = 0; ily < AliTRDgeometry::kNlayer; ily++){
192 // fill P & dE/dx information
193 if(HasFriends()){ // from TRD track
194 if(!trackTRD) continue;
195 AliTRDseedV1 *trackletTRD(0x0);
196 trackTRD -> SetReconstructor(fReconstructor);
197 if(!(trackletTRD = trackTRD -> GetTracklet(ily))) continue;
198 if(!CheckQuality(trackletTRD)) continue;
199 if(!CookdEdx(trackletTRD)) continue;
201 // fill momentum information
205 if(!(ref = track->GetTrackRef(trackletTRD))) continue;
209 fP = trackletTRD->GetMomentum();
213 } else { // from ESD track
214 // fill momentum information
217 if(!(ref = track->GetTrackRef(ily))) continue;
225 Double32_t *it = &infoPID[ily*AliTRDCalPID::kNSlicesNN];
226 for(Int_t is=AliTRDCalPID::kNSlicesNN; is--; it++) fdEdx[is] = (*it);
229 // momentum threshold
230 if(fP < fPthreshold) continue;
233 fPIDdataArray->PushBack(ily, AliTRDpidUtil::GetMomentumBin(fP), fdEdx);
239 PostData(0, fContainer);
244 //________________________________________________________________________
245 void AliTRDpidRefMaker::Fill()
249 if(!fPIDdataArray->fNtracklets) return;
250 fPIDbin = TMath::LocMax(AliPID::kSPECIES, fPID); // get particle type
256 for(Int_t ily=fPIDdataArray->fNtracklets; ily--;){
257 Int_t pBin = fPIDdataArray->fData[ily].fPLbin & 0xf;
258 ((TH2*)fContainer->At(0))->Fill(fPIDbin, pBin);
262 //________________________________________________________________________
263 void AliTRDpidRefMaker::LinkPIDdata()
265 // Link data tree to data members
266 fData->SetBranchAddress("s", &fPIDbin);
267 fData->SetBranchAddress("data", &fPIDdataArray);
270 //________________________________________________________________________
271 void AliTRDpidRefMaker::SetRefPID(ETRDpidRefMakerSource select, void *source, Float_t *pid)
273 // Fill the reference PID values "pid" from "source" object
274 // according to the option "select". Possible options are
275 // - kV0 - v0 based PID
276 // - kMC - MC truth [default]
277 // - kRec - outside detectors
280 memset(fPID, 0, AliPID::kSPECIES*sizeof(Float_t));
284 AliTRDtrackInfo *track = static_cast<AliTRDtrackInfo*>(source);
286 AliError("No trackInfo found");
290 //Get V0 PID decisions from the AliTRDv0Info for all particle species (implemented so far : electrons from conversions, pions from K0s and protons from Lambdas) :
292 for(Int_t is=AliPID::kSPECIES; is--;) fPID[is] = v0.GetV0PID(is, track);
297 AliError("Could not retrive reference PID from MC");
301 AliTRDtrackInfo *track = static_cast<AliTRDtrackInfo*>(source);
302 switch(track->GetPDG()){
305 fPID[AliPID::kElectron] = 1.;
309 fPID[AliPID::kMuon] = 1.;
313 fPID[AliPID::kPion] = 1.;
317 fPID[AliPID::kKaon] = 1.;
321 fPID[AliPID::kProton] = 1.;
328 AliTRDtrackInfo *track = static_cast<AliTRDtrackInfo*>(source);
329 AliTRDtrackV1 *trackTRD = track->GetTrack();
330 trackTRD -> SetReconstructor(fReconstructor);
331 //fReconstructor -> SetOption("nn");
332 trackTRD -> CookPID();
333 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
334 pid[iPart] = trackTRD -> GetPID(iPart);
335 AliDebug(4, Form("PDG is (in V0info) %d %f", iPart, pid[iPart]));
340 AliWarning("PID reference source not implemented");
343 AliDebug(4, Form("Ref PID [%] : %s[%5.2f] %s[%5.2f] %s[%5.2f] %s[%5.2f] %s[%5.2f]"
344 ,AliPID::ParticleShortName(0), 1.e2*fPID[0]
345 ,AliPID::ParticleShortName(1), 1.e2*fPID[1]
346 ,AliPID::ParticleShortName(2), 1.e2*fPID[2]
347 ,AliPID::ParticleShortName(3), 1.e2*fPID[3]
348 ,AliPID::ParticleShortName(4), 1.e2*fPID[4]
352 //________________________________________________________________________
353 void AliTRDpidRefMaker::SetAbundance(Float_t train)
355 // Split data sample between trainning and testing
357 if(train<0. || train >1.){
358 AliWarning("The input data should be in the interval [0, 1]");