1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /* $Id: AliTRDinfoGen.cxx 27496 2008-07-22 08:35:45Z cblume $ */
18 ////////////////////////////////////////////////////////////////////////////
20 // Tender wagon for TRD performance/calibration train
22 // In this wagon the information from
24 // - Friends [if available]
25 // - MC [if available]
26 // are grouped into AliTRDtrackInfo objects and fed to worker tasks
29 // Markus Fasel <M.Fasel@gsi.de>
30 // Alexandru Bercuci <A.Bercuci@gsi.de>
32 ////////////////////////////////////////////////////////////////////////////
34 #include <TClonesArray.h>
35 #include <TObjArray.h>
44 #include <TParticle.h>
47 #include "AliAnalysisManager.h"
48 #include "AliCDBManager.h"
49 #include "AliESDEvent.h"
50 #include "AliMCEvent.h"
51 #include "AliESDInputHandler.h"
52 #include "AliMCEventHandler.h"
54 #include "AliESDfriend.h"
55 #include "AliESDfriendTrack.h"
56 #include "AliESDHeader.h"
57 #include "AliESDtrack.h"
59 #include "AliESDtrackCuts.h"
60 #include "AliMCParticle.h"
63 #include "AliTrackReference.h"
64 #include "TTreeStream.h"
72 #include "AliTRDcalibDB.h"
73 #include "AliTRDtrackerV1.h"
74 #include "AliTRDgeometry.h"
75 #include "AliTRDtrackV1.h"
76 #include "AliTRDseedV1.h"
77 #include "AliTRDcluster.h"
78 #include "AliTRDinfoGen.h"
79 #include "info/AliTRDtrackInfo.h"
80 #include "info/AliTRDeventInfo.h"
81 #include "info/AliTRDv0Info.h"
82 #include "info/AliTRDeventCuts.h"
83 #include "macros/AliTRDperformanceTrain.h"
85 ClassImp(AliTRDinfoGen)
87 const Float_t AliTRDinfoGen::fgkITS = 100.; // to be checked
88 const Float_t AliTRDinfoGen::fgkTPC = 290.;
89 const Float_t AliTRDinfoGen::fgkTRD = 365.;
91 const Float_t AliTRDinfoGen::fgkEvVertexZ = 15.;
92 const Int_t AliTRDinfoGen::fgkEvVertexN = 1;
94 const Float_t AliTRDinfoGen::fgkTrkDCAxy = 3.;
95 const Float_t AliTRDinfoGen::fgkTrkDCAz = 10.;
96 const Int_t AliTRDinfoGen::fgkNclTPC = 70;
97 const Float_t AliTRDinfoGen::fgkPt = 0.2;
98 const Float_t AliTRDinfoGen::fgkEta = 0.9;
100 //____________________________________________________________________
101 AliTRDinfoGen::AliTRDinfoGen()
109 ,fOCDB("local://$ALICE_ROOT/OCDB")
121 // Default constructor
123 SetNameTitle("TRDinfoGen", "MC-REC TRD-track list generator");
126 //____________________________________________________________________
127 AliTRDinfoGen::AliTRDinfoGen(char* name)
128 :AliAnalysisTaskSE(name)
135 ,fOCDB("local://$ALICE_ROOT/OCDB")
147 // Default constructor
149 SetTitle("MC-REC TRD-track list generator");
150 DefineOutput(kTracksBarrel, TObjArray::Class());
151 DefineOutput(kTracksSA, TObjArray::Class());
152 DefineOutput(kTracksKink, TObjArray::Class());
153 DefineOutput(kEventInfo, AliTRDeventInfo::Class());
154 DefineOutput(kV0List, TObjArray::Class());
155 DefineOutput(kMonitor, TObjArray::Class()); // histogram list
158 //____________________________________________________________________
159 AliTRDinfoGen::~AliTRDinfoGen()
163 if(fDebugStream) delete fDebugStream;
164 if(fEvTrigger) delete fEvTrigger;
165 if(fV0Cut) delete fV0Cut;
166 if(fTrackCut) delete fTrackCut;
167 if(fEventCut) delete fEventCut;
168 if(fTrackInfo) delete fTrackInfo; fTrackInfo = NULL;
169 if(fEventInfo) delete fEventInfo; fEventInfo = NULL;
170 if(fV0Info) delete fV0Info; fV0Info = NULL;
172 fTracksBarrel->Delete(); delete fTracksBarrel;
173 fTracksBarrel = NULL;
176 fTracksSA->Delete(); delete fTracksSA;
180 fTracksKink->Delete(); delete fTracksKink;
189 fContainer->Delete();
195 //____________________________________________________________________
196 Bool_t AliTRDinfoGen::GetRefFigure(Int_t)
199 AliWarning("Please provide a canvas to draw results.");
202 fContainer->At(0)->Draw("bar");
206 //____________________________________________________________________
207 void AliTRDinfoGen::UserCreateOutputObjects()
210 // Create Output Containers (TObjectArray containing 1D histograms)
213 fTrackInfo = new AliTRDtrackInfo();
214 fEventInfo = new AliTRDeventInfo();
215 fV0Info = new AliTRDv0Info();
216 fTracksBarrel = new TObjArray(200); fTracksBarrel->SetOwner(kTRUE);
217 fTracksSA = new TObjArray(20); fTracksSA->SetOwner(kTRUE);
218 fTracksKink = new TObjArray(20); fTracksKink->SetOwner(kTRUE);
219 fV0List = new TObjArray(10); fV0List->SetOwner(kTRUE);
221 // define general monitor
222 fContainer = new TObjArray(1); fContainer->SetOwner(kTRUE);
223 TH1 *h=new TH1S("hStat", "Run statistics;Observable;Entries", 12, -0.5, 11.5);
224 TAxis *ax(h->GetXaxis());
225 ax->SetBinLabel( 1, "ESD");
226 ax->SetBinLabel( 2, "MC");
227 ax->SetBinLabel( 3, "V0");
228 ax->SetBinLabel( 4, "TPC");
229 ax->SetBinLabel( 5, "TRDin");
230 ax->SetBinLabel( 6, "TRDout");
231 ax->SetBinLabel( 7, "Barrel");
232 ax->SetBinLabel( 8, "BarrelMC");
233 ax->SetBinLabel( 9, "SA");
234 ax->SetBinLabel(10, "SAMC");
235 ax->SetBinLabel(11, "Kink");
236 ax->SetBinLabel(12, "KinkMC");
237 fContainer->AddAt(h, 0);
238 PostData(kMonitor, fContainer);
241 //____________________________________________________________________
242 Bool_t AliTRDinfoGen::Load(const Char_t *file, const Char_t *dir, const Char_t *name)
244 // Load data from performance file
246 if(!TFile::Open(file)){
247 AliWarning(Form("Couldn't open file %s.", file));
252 AliWarning(Form("Couldn't cd to %s in %s.", dir, file));
257 const Char_t *tn=(name ? name : GetName());
258 if(!(o = (TObjArray*)gDirectory->Get(tn))){
259 AliWarning(Form("Missing histogram container %s.", tn));
262 fContainer = (TObjArray*)o->Clone(GetName());
267 //____________________________________________________________________
268 void AliTRDinfoGen::UserExec(Option_t *){
273 fTracksBarrel->Delete();
275 fTracksKink->Delete();
277 fEventInfo->Delete("");
279 fESDev = dynamic_cast<AliESDEvent*>(InputEvent());
281 AliError("Failed retrieving ESD event");
285 AliInfo("Initializing OCDB ...");
286 // prepare OCDB access
287 AliCDBManager* ocdb = AliCDBManager::Instance();
288 ocdb->SetDefaultStorage(fOCDB.Data());
289 ocdb->SetRun(fESDev->GetRunNumber());
290 AliTRDtrackerV1::SetNTimeBins(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
291 AliInfo(Form("OCDB : Loc[%s] Run[%d] TB[%d]", fOCDB.Data(), ocdb->GetRun(), AliTRDtrackerV1::GetNTimeBins()));
295 // link MC if available
298 // event selection : trigger cut
299 if(UseLocalEvSelection() && fEvTrigger){
300 Bool_t kTRIGGERED(kFALSE);
301 std::auto_ptr<TObjArray> trig(fEvTrigger->Tokenize(" "));
302 for(Int_t itrig=trig->GetEntriesFast(); itrig--;){
303 const Char_t *trigClass(((TObjString*)(*trig)[itrig])->GetName());
304 if(fESDev->IsTriggerClassFired(trigClass)) {
305 AliDebug(2, Form("Ev[%4d] Trigger[%s]", fESDev->GetEventNumberInFile(), trigClass));
311 AliDebug(2, Form("Reject Ev[%4d] Trigger", fESDev->GetEventNumberInFile()));
314 // select only physical events
315 if(fESDev->GetEventType() != 7){
316 AliDebug(2, Form("Reject Ev[%4d] EvType[%d]", fESDev->GetEventNumberInFile(), fESDev->GetEventType()));
321 // if the required trigger is a collision trigger then apply event vertex cut
322 if(UseLocalEvSelection() && IsCollision()){
323 const AliESDVertex *vertex = fESDev->GetPrimaryVertex();
324 if(TMath::Abs(vertex->GetZv())<1.e-10 ||
325 TMath::Abs(vertex->GetZv())>fgkEvVertexZ ||
326 vertex->GetNContributors()<fgkEvVertexN) {
327 AliDebug(2, Form("Reject Ev[%4d] Vertex Zv[%f] Nv[%d]", fESDev->GetEventNumberInFile(), TMath::Abs(vertex->GetZv()), vertex->GetNContributors()));
332 if(fEventCut && !fEventCut->IsSelected(fESDev, IsCollision())) return;
335 AliError("Failed retrieving ESD friend event");
338 if(HasMCdata() && !fMCev){
339 AliError("Failed retrieving MC event");
343 new(fEventInfo)AliTRDeventInfo(fESDev->GetHeader(), const_cast<AliESDRun *>(fESDev->GetESDRun()));
345 Bool_t *trackMap(NULL);
346 AliStack * mStack(NULL);
347 Int_t nTracksMC = HasMCdata() ? fMCev->GetNumberOfTracks() : 0, nTracksESD = fESDev->GetNumberOfTracks();
349 mStack = fMCev->Stack();
351 AliError("Failed retrieving MC Stack");
354 trackMap = new Bool_t[nTracksMC];
355 memset(trackMap, 0, sizeof(Bool_t) * nTracksMC);
358 Double32_t dedx[100]; Int_t nSlices(0);
359 Int_t nTRDout(0), nTRDin(0), nTPC(0)
361 ,nBarrel(0), nSA(0), nKink(0)
362 ,nBarrelMC(0), nSAMC(0), nKinkMC(0);
363 AliESDtrack *esdTrack = NULL;
364 AliESDfriendTrack *esdFriendTrack = NULL;
365 TObject *calObject = NULL;
366 AliTRDtrackV1 *track = NULL;
367 AliTRDseedV1 *tracklet = NULL;
368 AliTRDcluster *cl = NULL;
371 // LOOP 0 - over ESD v0s
372 Float_t bField(fESDev->GetMagneticField());
374 Int_t v0pid[AliPID::kSPECIES];
375 for(Int_t iv0(0); iv0<fESDev->GetNumberOfV0s(); iv0++){
376 if(!(v0 = fESDev->GetV0(iv0))) continue;
378 if(fV0Cut) new(fV0Info) AliTRDv0Info(*fV0Cut);
379 else new(fV0Info) AliTRDv0Info();
380 fV0Info->SetMagField(bField);
381 fV0Info->SetV0tracks(fESDev->GetTrack(v0->GetPindex()), fESDev->GetTrack(v0->GetNindex()));
382 fV0Info->SetV0Info(v0);
383 fV0List->Add(new AliTRDv0Info(*fV0Info));// kFOUND=kFALSE;
387 // LOOP 1 - over ESD tracks
388 AliTRDv0Info *v0info=NULL;
389 for(Int_t itrk = 0; itrk < nTracksESD; itrk++){
390 new(fTrackInfo) AliTRDtrackInfo();
391 esdTrack = fESDev->GetTrack(itrk);
392 AliDebug(3, Form("\n%3d ITS[%d] TPC[%d] TRD[%d]\n", itrk, esdTrack->GetNcls(0), esdTrack->GetNcls(1), esdTrack->GetNcls(2)));
394 if(esdTrack->GetStatus()&AliESDtrack::kTPCout) nTPC++;
395 if(esdTrack->GetStatus()&AliESDtrack::kTRDout) nTRDout++;
396 if(esdTrack->GetStatus()&AliESDtrack::kTRDin) nTRDin++;
398 // look at external track param
399 const AliExternalTrackParam *op = esdTrack->GetOuterParam();
403 op->Global2LocalPosition(xyz, op->GetAlpha());
404 AliDebug(3, Form("op @ X[%7.3f]\n", xyz[0]));
409 Int_t label = -1; UInt_t alab=UINT_MAX;
411 label = esdTrack->GetLabel();
412 alab = TMath::Abs(label);
413 // register the track
414 if(alab < UInt_t(nTracksMC)){
415 trackMap[alab] = kTRUE;
417 AliError(Form("MC label[%d] outside scope for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
420 AliMCParticle *mcParticle = NULL;
421 if(!(mcParticle = (AliMCParticle*) fMCev->GetTrack(alab))){
422 AliError(Form("MC particle label[%d] missing for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
425 fPdg = mcParticle->Particle()->GetPdgCode();
426 Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
427 Int_t iref = 0; AliTrackReference *ref = NULL;
429 ref = mcParticle->GetTrackReference(iref);
430 if(ref->LocalX() > fgkTPC) break;
435 fTrackInfo->SetPDG(fPdg);
436 fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
437 fTrackInfo->SetLabel(label);
438 Int_t jref = iref;//, kref = 0;
440 ref = mcParticle->GetTrackReference(jref);
441 if(ref->LocalX() > fgkTRD) break;
442 AliDebug(4, Form(" trackRef[%2d (%2d)] @ %7.3f OK", jref-iref, jref, ref->LocalX()));
443 fTrackInfo->AddTrackRef(ref);
446 AliDebug(3, Form("NtrackRefs[%d(%d)]", fTrackInfo->GetNTrackRefs(), nRefs));
449 // copy some relevant info to TRD track info
450 fTrackInfo->SetStatus(esdTrack->GetStatus());
451 fTrackInfo->SetTrackId(esdTrack->GetID());
452 Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);
453 fTrackInfo->SetESDpid(p);
454 fTrackInfo->SetESDpidQuality(esdTrack->GetTRDntrackletsPID());
455 if(!nSlices) nSlices = esdTrack->GetNumberOfTRDslices();
456 memset(dedx, 0, 100*sizeof(Double32_t));
458 for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++)
459 for(Int_t is=0; is<nSlices; is++)
460 dedx[in++]=esdTrack->GetTRDslice(il, is);
461 for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++) dedx[in++]=esdTrack->GetTRDmomentum(il);
462 fTrackInfo->SetSlices(in, dedx);
463 fTrackInfo->SetNumberOfClustersRefit(esdTrack->GetNcls(2));
464 // some other Informations which we may wish to store in order to find problematic cases
465 fTrackInfo->SetKinkIndex(esdTrack->GetKinkIndex(0));
466 fTrackInfo->SetTPCncls(static_cast<UShort_t>(esdTrack->GetNcls(1)));
470 for(Int_t iv(0); iv<fV0List->GetEntriesFast(); iv++){
471 if(!(v0info = (AliTRDv0Info*)fV0List->At(iv))) continue;
472 if(!v0info->fTrackP && !v0info->fTrackN) continue;
473 if(!v0info->HasTrack(fTrackInfo)) continue;
474 memset(v0pid, 0, AliPID::kSPECIES*sizeof(Int_t));
476 for(Int_t is=AliPID::kSPECIES; is--;){v0pid[is] = v0info->GetPID(is, fTrackInfo);}
477 fTrackInfo->SetV0pid(v0pid);
479 //const AliTRDtrackInfo::AliESDinfo *ei = fTrackInfo->GetESDinfo();
484 esdFriendTrack = fESDfriend->GetTrack(itrk);
487 while((calObject = esdFriendTrack->GetCalibObject(icalib++))){
488 if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
489 if(!(track = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
490 AliDebug(4, Form("TRD track OK"));
491 // Set the clusters to unused
492 for(Int_t ipl = 0; ipl < AliTRDgeometry::kNlayer; ipl++){
493 if(!(tracklet = track->GetTracklet(ipl))) continue;
494 tracklet->ResetClusterIter();
495 while((cl = tracklet->NextCluster())) cl->Use(0);
497 fTrackInfo->SetTrack(track);
500 AliDebug(3, Form("Ntracklets[%d]\n", fTrackInfo->GetNTracklets()));
501 } else AliDebug(3, "No ESD friends");
502 if(op) fTrackInfo->SetOuterParam(op);
504 if(DebugLevel() >= 1){
505 AliTRDtrackInfo info(*fTrackInfo);
506 (*DebugStream()) << "trackInfo"
507 << "TrackInfo.=" << &info
512 ULong_t status(esdTrack->GetStatus());
513 if((status&AliESDtrack::kTPCout)){
514 if(!esdTrack->GetKinkIndex(0)){ // Barrel Track Selection
515 Bool_t selected(kTRUE);
516 if(UseLocalTrkSelection()){
517 if(esdTrack->Pt() < fgkPt){
518 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Pt[%5.2f]", fESDev->GetEventNumberInFile(), itrk, esdTrack->Pt()));
521 if(selected && TMath::Abs(esdTrack->Eta()) > fgkEta){
522 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Eta[%5.2f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(esdTrack->Eta())));
525 if(selected && esdTrack->GetTPCNcls() < fgkNclTPC){
526 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] NclTPC[%d]", fESDev->GetEventNumberInFile(), itrk, esdTrack->GetTPCNcls()));
529 Float_t par[2], cov[3];
530 esdTrack->GetImpactParameters(par, cov);
531 if(IsCollision()){ // cuts on DCA
532 if(selected && TMath::Abs(par[0]) > fgkTrkDCAxy){
533 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAxy[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[0])));
536 if(selected && TMath::Abs(par[1]) > fgkTrkDCAz){
537 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAz[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[1])));
540 } else if(selected && fMCev && !fMCev->IsPhysicalPrimary(alab)){;
541 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Primary", fESDev->GetEventNumberInFile(), itrk));
545 if(fTrackCut && !fTrackCut->IsSelected(esdTrack)) selected = kFALSE;
547 fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
551 fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
554 } else if((status&AliESDtrack::kTRDout) && !(status&AliESDtrack::kTRDin)){
555 fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
558 fTrackInfo->Delete("");
561 // LOOP 2 - over MC tracks which are passing TRD where the track is not reconstructed
563 AliDebug(10, "Output of the MC track map:");
564 for(Int_t itk = 0; itk < nTracksMC; itk++) AliDebug(10, Form("trackMap[%d] = %s", itk, trackMap[itk] == kTRUE ? "TRUE" : "kFALSE"));
566 for(Int_t itk = 0; itk < nTracksMC; itk++){
567 if(trackMap[itk]) continue;
568 AliMCParticle *mcParticle = (AliMCParticle*) fMCev->GetTrack(TMath::Abs(itk));
569 Int_t fPdg = mcParticle->Particle()->GetPdgCode();
570 Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
571 Int_t iref = 0; AliTrackReference *ref = NULL;
573 new(fTrackInfo) AliTRDtrackInfo();
575 fTrackInfo->SetPDG(fPdg);
576 while(iref<nRefs){ // count TRD TR
578 ref = mcParticle->GetTrackReference(iref);
579 if(ref->LocalX() > fgkTPC && ref->LocalX() < fgkTRD){
580 fTrackInfo->AddTrackRef(ref);
581 nRefsTRD++;kIN=kTRUE;
583 AliDebug(4, Form(" trackRef[%2d] @ x[%7.3f] %s", iref, ref->LocalX(), kIN?"IN":"OUT"));
587 // In this stage we at least require 1 hit inside TRD. What will be done with this tracks is a task for the
589 fTrackInfo->Delete("");
592 fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
593 fTrackInfo->SetLabel(itk);
594 if(DebugLevel() >= 1){
595 AliTRDtrackInfo info(*fTrackInfo);
596 (*DebugStream()) << "trackInfo"
597 << "TrackInfo.=" << &info
601 AliDebug(3, Form("Add MC track @ label[%d] nTRDrefs[%d].", itk, nRefsTRD));
602 // check where the track starts
603 ref = mcParticle->GetTrackReference(0);
604 if(ref->LocalX() < fgkITS){
605 fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
607 } else if(ref->LocalX() < fgkTPC) {
608 fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
610 } else if(nRefsTRD>6){
611 fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
614 fTrackInfo->Delete("");
619 "\nEv[%3d] Tracks: ESD[%d] MC[%d] V0[%d]\n"
620 " TPCout[%d] TRDin[%d] TRDout[%d]\n"
621 " Barrel[%3d+%3d=%3d] SA[%2d+%2d=%2d] Kink[%2d+%2d=%2d]"
622 ,(Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTracksESD, nTracksMC, fV0List->GetEntries()
623 , nTPC, nTRDin, nTRDout
624 ,nBarrel, nBarrelMC, fTracksBarrel->GetEntries()
625 ,nSA, nSAMC, fTracksSA->GetEntries()
626 ,nKink, nKinkMC, fTracksKink->GetEntries()
628 // save track statistics
629 TH1 *h((TH1S*)fContainer->At(0));
630 h->Fill( 0., nTracksESD);
631 h->Fill( 1., nTracksMC);
632 h->Fill( 2., fV0List->GetEntries());
634 h->Fill( 4., nTRDin);
635 h->Fill( 5., nTRDout);
636 h->Fill( 6., nBarrel);
637 h->Fill( 7., nBarrelMC);
641 h->Fill(11., nKinkMC);
643 PostData(kTracksBarrel, fTracksBarrel);
644 PostData(kTracksSA, fTracksSA);
645 PostData(kTracksKink, fTracksKink);
646 PostData(kEventInfo, fEventInfo);
647 PostData(kV0List, fV0List);
650 //____________________________________________________________________
651 void AliTRDinfoGen::SetLocalV0Selection(AliTRDv0Info *v0)
653 // Set V0 cuts from outside
655 if(!fV0Cut) fV0Cut = new AliTRDv0Info(*v0);
656 else new(fV0Cut) AliTRDv0Info(*v0);
659 //____________________________________________________________________
660 void AliTRDinfoGen::SetTrigger(const Char_t *trigger)
662 if(!fEvTrigger) fEvTrigger = new TString(trigger);
663 else (*fEvTrigger) = trigger;
666 //____________________________________________________________________
667 TTreeSRedirector* AliTRDinfoGen::DebugStream()
670 TDirectory *savedir = gDirectory;
671 fDebugStream = new TTreeSRedirector("TRD.DebugInfoGen.root");