]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDinfoGen.cxx
remove confusing setter for MC
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDinfoGen.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
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 **************************************************************************/
15
16 /* $Id: AliTRDinfoGen.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //
20 //  Tender wagon for TRD performance/calibration train
21 //
22 // In this wagon the information from
23 //   - ESD
24 //   - Friends [if available]
25 //   - MC [if available]
26 // are grouped into AliTRDtrackInfo objects and fed to worker tasks
27 //
28 //  Authors:
29 //    Markus Fasel <M.Fasel@gsi.de>
30 //    Alexandru Bercuci <A.Bercuci@gsi.de>
31 //
32 ////////////////////////////////////////////////////////////////////////////
33
34 #include <TClonesArray.h>
35 #include <TObjArray.h>
36 #include <TObject.h>
37 #include <TString.h>
38 #include <TH1F.h>
39 #include <TFile.h>
40 #include <TTree.h>
41 #include <TROOT.h>
42 #include <TChain.h>
43 #include <TParticle.h>
44
45 #include "AliLog.h"
46 #include "AliAnalysisManager.h"
47 #include "AliESDEvent.h"
48 #include "AliMCEvent.h"
49 #include "AliESDInputHandler.h"
50 #include "AliMCEventHandler.h"
51
52 #include "AliESDfriend.h"
53 #include "AliESDfriendTrack.h"
54 #include "AliESDHeader.h"
55 #include "AliESDtrack.h"
56 #include "AliESDtrackCuts.h"
57 #include "AliMCParticle.h"
58 #include "AliPID.h"
59 #include "AliStack.h"
60 #include "AliTRDtrackV1.h"
61 #include "AliTrackReference.h"
62 #include "AliTRDgeometry.h"
63 #include "AliTRDcluster.h"
64 #include "AliTRDseedV1.h"
65 #include "TTreeStream.h"
66
67 #include <cstdio>
68 #include <climits>
69 #include <cstring>
70 #include <iostream>
71
72 #include "AliTRDinfoGen.h"
73 #include "info/AliTRDtrackInfo.h"
74 #include "info/AliTRDeventInfo.h"
75 #include "info/AliTRDv0Info.h"
76 #include "info/AliTRDeventCuts.h"
77 #include "macros/AliTRDperformanceTrain.h"
78
79 ClassImp(AliTRDinfoGen)
80
81 const Float_t AliTRDinfoGen::fgkITS = 100.; // to be checked
82 const Float_t AliTRDinfoGen::fgkTPC = 290.;
83 const Float_t AliTRDinfoGen::fgkTRD = 365.;
84
85 const Float_t AliTRDinfoGen::fgkEvVertexZ = 15.;
86 const Int_t   AliTRDinfoGen::fgkEvVertexN = 1;
87
88 const Float_t AliTRDinfoGen::fgkTrkDCAxy  = 40.;
89 const Float_t AliTRDinfoGen::fgkTrkDCAz   = 15.;
90 const Int_t   AliTRDinfoGen::fgkNclTPC    = 100;
91 const Float_t AliTRDinfoGen::fgkPt        = 0.2;
92 const Float_t AliTRDinfoGen::fgkEta       = 0.9;
93
94 //____________________________________________________________________
95 AliTRDinfoGen::AliTRDinfoGen()
96   :AliAnalysisTaskSE()
97   ,fEvTrigger(NULL)
98   ,fESDev(NULL)
99   ,fMCev(NULL)
100   ,fEventCut(NULL)
101   ,fTrackCut(NULL)
102   ,fTrackInfo(NULL)
103   ,fEventInfo(NULL)
104   ,fV0Info(NULL)
105   ,fTracksBarrel(NULL)
106   ,fTracksSA(NULL)
107   ,fTracksKink(NULL)
108   ,fV0List(NULL)
109   ,fDebugStream(NULL)
110 {
111   //
112   // Default constructor
113   //
114   SetNameTitle("infoGen", "MC-REC TRD-track list generator");
115 }
116
117 //____________________________________________________________________
118 AliTRDinfoGen::AliTRDinfoGen(char* name)
119   :AliAnalysisTaskSE(name)
120   ,fEvTrigger(NULL)
121   ,fESDev(NULL)
122   ,fMCev(NULL)
123   ,fEventCut(NULL)
124   ,fTrackCut(NULL)
125   ,fTrackInfo(NULL)
126   ,fEventInfo(NULL)
127   ,fV0Info(NULL)
128   ,fTracksBarrel(NULL)
129   ,fTracksSA(NULL)
130   ,fTracksKink(NULL)
131   ,fV0List(NULL)
132   ,fDebugStream(NULL)
133 {
134   //
135   // Default constructor
136   //
137   SetTitle("MC-REC TRD-track list generator");
138   DefineOutput(kTracksBarrel, TObjArray::Class());
139   DefineOutput(kTracksSA, TObjArray::Class());
140   DefineOutput(kTracksKink, TObjArray::Class());
141   DefineOutput(kEventInfo, AliTRDeventInfo::Class());
142   DefineOutput(kV0List, TObjArray::Class());
143 }
144
145 //____________________________________________________________________
146 AliTRDinfoGen::~AliTRDinfoGen()
147 {
148 // Destructor
149   
150   if(fDebugStream) delete fDebugStream;
151   if(fEvTrigger) delete fEvTrigger;
152   if(fTrackCut) delete fTrackCut;
153   if(fEventCut) delete fEventCut;
154   if(fTrackInfo) delete fTrackInfo; fTrackInfo = NULL;
155   if(fEventInfo) delete fEventInfo; fEventInfo = NULL;
156   if(fV0Info) delete fV0Info; fV0Info = NULL;
157   if(fTracksBarrel){ 
158     fTracksBarrel->Delete(); delete fTracksBarrel;
159     fTracksBarrel = NULL;
160   }
161   if(fTracksSA){ 
162     fTracksSA->Delete(); delete fTracksSA;
163     fTracksSA = NULL;
164   }
165   if(fTracksKink){ 
166     fTracksKink->Delete(); delete fTracksKink;
167     fTracksKink = NULL;
168   }
169   if(fV0List){ 
170     fV0List->Delete(); delete fV0List;
171     fV0List = NULL;
172   }
173 }
174
175 //____________________________________________________________________
176 void AliTRDinfoGen::UserCreateOutputObjects()
177 {       
178   //
179   // Create Output Containers (TObjectArray containing 1D histograms)
180   //
181  
182   fTrackInfo = new AliTRDtrackInfo();
183   fEventInfo = new AliTRDeventInfo();
184   fV0Info    = new AliTRDv0Info();
185   fTracksBarrel = new TObjArray(200); fTracksBarrel->SetOwner(kTRUE);
186   fTracksSA = new TObjArray(20); fTracksSA->SetOwner(kTRUE);
187   fTracksKink = new TObjArray(20); fTracksKink->SetOwner(kTRUE);
188   fV0List = new TObjArray(10); fV0List->SetOwner(kTRUE);
189 }
190
191 //____________________________________________________________________
192 void AliTRDinfoGen::UserExec(Option_t *){
193   //
194   // Run the Analysis
195   //
196
197   fESDev = dynamic_cast<AliESDEvent*>(InputEvent());
198   fMCev = MCEvent();
199
200   if(!fESDev){
201     AliError("Failed retrieving ESD event");
202     return;
203   }
204
205   // event selection : trigger cut
206   if(UseLocalEvSelection() && fEvTrigger){ 
207     Bool_t kTRIGGERED(kFALSE);
208     const TObjArray *trig = fEvTrigger->Tokenize(" ");
209     for(Int_t itrig=trig->GetEntriesFast(); itrig--;){
210       const Char_t *trigClass(((TObjString*)(*trig)[itrig])->GetName());
211       if(fESDev->IsTriggerClassFired(trigClass)) {
212         AliDebug(2, Form("Ev[%4d] Trigger[%s]", fESDev->GetEventNumberInFile(), trigClass));
213         kTRIGGERED = kTRUE;
214         break; 
215       }
216     }
217     if(!kTRIGGERED){ 
218       AliDebug(2, Form("Reject Ev[%4d] Trigger", fESDev->GetEventNumberInFile()));
219       return;
220     }
221     // select only physical events
222     if(fESDev->GetEventType() != 7){ 
223       AliDebug(2, Form("Reject Ev[%4d] EvType[%d]", fESDev->GetEventNumberInFile(), fESDev->GetEventType()));
224       return;
225     }
226   }
227
228   // if the required trigger is a collision trigger then apply event vertex cut
229   if(UseLocalEvSelection() && IsCollision()){
230     const AliESDVertex *vertex = fESDev->GetPrimaryVertex();
231     if(TMath::Abs(vertex->GetZv())<1.e-10 || 
232        TMath::Abs(vertex->GetZv())>fgkEvVertexZ || 
233        vertex->GetNContributors()<fgkEvVertexN) {
234       AliDebug(2, Form("Reject Ev[%4d] Vertex Zv[%f] Nv[%d]", fESDev->GetEventNumberInFile(), TMath::Abs(vertex->GetZv()), vertex->GetNContributors()));
235       return;
236     }
237   }
238
239   if(fEventCut && !fEventCut->IsSelected(fESDev, IsCollision())) return;
240
241   if(!fESDfriend){
242     AliError("Failed retrieving ESD friend event");
243     return;
244   }
245   if(HasMCdata() && !fMCev){
246     AliError("Failed retrieving MC event");
247     return;
248   }
249
250   fTracksBarrel->Delete();
251   fTracksSA->Delete();
252   fTracksKink->Delete();
253   fV0List->Delete();
254   fEventInfo->Delete("");
255   new(fEventInfo)AliTRDeventInfo(fESDev->GetHeader(), const_cast<AliESDRun *>(fESDev->GetESDRun()));
256   
257   Bool_t *trackMap(NULL);
258   AliStack * mStack(NULL);
259   Int_t nTracksMC = HasMCdata() ? fMCev->GetNumberOfTracks() : 0, nTracksESD = fESDev->GetNumberOfTracks();
260   if(HasMCdata()){
261     mStack = fMCev->Stack();
262     if(!mStack){
263       AliError("Failed retrieving MC Stack");
264       return;
265     }
266     trackMap = new Bool_t[nTracksMC];
267     memset(trackMap, 0, sizeof(Bool_t) * nTracksMC);
268   }
269   
270   Double32_t dedx[100]; Int_t nSlices(0);
271   Int_t nTRDout(0), nTRDin(0), nTPC(0)
272        ,nclsTrklt
273        ,nBarrel(0), nSA(0), nKink(0)
274        ,nBarrelMC(0), nSAMC(0), nKinkMC(0);
275   AliESDtrack *esdTrack = NULL;
276   AliESDfriendTrack *esdFriendTrack = NULL;
277   TObject *calObject = NULL;
278   AliTRDtrackV1 *track = NULL;
279   AliTRDseedV1 *tracklet = NULL;
280   AliTRDcluster *cl = NULL;
281   for(Int_t itrk = 0; itrk < nTracksESD; itrk++){
282     new(fTrackInfo) AliTRDtrackInfo();
283     esdTrack = fESDev->GetTrack(itrk);
284     AliDebug(3, Form("\n%3d ITS[%d] TPC[%d] TRD[%d]\n", itrk, esdTrack->GetNcls(0), esdTrack->GetNcls(1), esdTrack->GetNcls(2)));
285
286     if(esdTrack->GetStatus()&AliESDtrack::kTPCout) nTPC++;
287     if(esdTrack->GetStatus()&AliESDtrack::kTRDout) nTRDout++;
288     if(esdTrack->GetStatus()&AliESDtrack::kTRDin) nTRDin++;
289
290     // look at external track param
291     const AliExternalTrackParam *op = esdTrack->GetOuterParam();
292     Double_t xyz[3];
293     if(op){
294       op->GetXYZ(xyz);
295       op->Global2LocalPosition(xyz, op->GetAlpha());
296       AliDebug(3, Form("op @ X[%7.3f]\n", xyz[0]));
297     }
298
299     // read MC info
300     Int_t fPdg = -1;
301     Int_t label = -1; UInt_t alab=UINT_MAX;
302     if(HasMCdata()){
303       label = esdTrack->GetLabel(); 
304       alab = TMath::Abs(label);
305       // register the track
306       if(alab < UInt_t(nTracksMC)){ 
307         trackMap[alab] = kTRUE; 
308       } else { 
309         AliError(Form("MC label[%d] outside scope for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
310         continue; 
311       }
312       AliMCParticle *mcParticle = NULL; 
313       if(!(mcParticle = (AliMCParticle*) fMCev->GetTrack(alab))){
314         AliError(Form("MC particle label[%d] missing for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
315         continue;
316       }
317       fPdg = mcParticle->Particle()->GetPdgCode();
318       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
319       Int_t iref = 0; AliTrackReference *ref = NULL; 
320       while(iref<nRefs){
321         ref = mcParticle->GetTrackReference(iref);
322         if(ref->LocalX() > fgkTPC) break;
323         iref++;
324       }
325
326       fTrackInfo->SetMC();
327       fTrackInfo->SetPDG(fPdg);
328       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
329       Int_t jref = iref;//, kref = 0;
330       while(jref<nRefs){
331         ref = mcParticle->GetTrackReference(jref);
332         if(ref->LocalX() > fgkTRD) break;
333         AliDebug(4, Form("  trackRef[%2d (%2d)] @ %7.3f OK", jref-iref, jref, ref->LocalX()));
334         fTrackInfo->AddTrackRef(ref);
335         jref++;
336       }
337       AliDebug(3, Form("NtrackRefs[%d(%d)]", fTrackInfo->GetNTrackRefs(), nRefs));
338     }
339
340     // copy some relevant info to TRD track info
341     fTrackInfo->SetStatus(esdTrack->GetStatus());
342     fTrackInfo->SetTrackId(esdTrack->GetID());
343     Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);
344     fTrackInfo->SetESDpid(p);
345     fTrackInfo->SetESDpidQuality(esdTrack->GetTRDntrackletsPID());
346     if(!nSlices) nSlices = esdTrack->GetNumberOfTRDslices();
347     memset(dedx, 0, 100*sizeof(Double32_t));
348     Int_t in(0);
349     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++)
350       for(Int_t is=0; is<nSlices; is++) 
351         dedx[in++]=esdTrack->GetTRDslice(il, is);
352     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++) dedx[in++]=esdTrack->GetTRDmomentum(il);
353     fTrackInfo->SetSlices(in, dedx);
354     fTrackInfo->SetLabel(label);
355     fTrackInfo->SetNumberOfClustersRefit(esdTrack->GetNcls(2));
356     // some other Informations which we may wish to store in order to find problematic cases
357     fTrackInfo->SetKinkIndex(esdTrack->GetKinkIndex(0));
358     fTrackInfo->SetTPCncls(static_cast<UShort_t>(esdTrack->GetNcls(1)));
359     nclsTrklt = 0;
360   
361
362     // read REC info
363     esdFriendTrack = fESDfriend->GetTrack(itrk);
364     if(esdFriendTrack){
365       Int_t icalib = 0;
366       while((calObject = esdFriendTrack->GetCalibObject(icalib++))){
367         if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
368         if(!(track = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
369         AliDebug(4, Form("TRD track OK"));
370         // Set the clusters to unused
371         for(Int_t ipl = 0; ipl < AliTRDgeometry::kNlayer; ipl++){
372           if(!(tracklet = track->GetTracklet(ipl))) continue;
373           tracklet->ResetClusterIter();
374           while((cl = tracklet->NextCluster())) cl->Use(0);
375         }
376         fTrackInfo->SetTrack(track);
377         break;
378       }
379       AliDebug(3, Form("Ntracklets[%d]\n", fTrackInfo->GetNTracklets()));
380     } else AliDebug(3, "No ESD friends");
381     if(op) fTrackInfo->SetOuterParam(op);
382
383     if(DebugLevel() >= 1){
384       AliTRDtrackInfo info(*fTrackInfo);
385       (*DebugStream()) << "trackInfo"
386       << "TrackInfo.=" << &info
387       << "\n";
388       info.Delete("");
389     }
390
391     ULong_t status(esdTrack->GetStatus());
392     if((status&AliESDtrack::kTPCout)){
393       if(!esdTrack->GetKinkIndex(0)){ // Barrel  Track Selection
394         Bool_t selected(kTRUE);
395         if(UseLocalTrkSelection()){
396           if(esdTrack->Pt() < fgkPt){ 
397             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Pt[%5.2f]", fESDev->GetEventNumberInFile(), itrk, esdTrack->Pt()));
398             selected = kFALSE;
399           }
400           if(selected && TMath::Abs(esdTrack->Eta()) > fgkEta){
401             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Eta[%5.2f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(esdTrack->Eta())));
402             selected = kFALSE;
403           }
404           if(selected && esdTrack->GetTPCNcls() < fgkNclTPC){ 
405             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] NclTPC[%d]", fESDev->GetEventNumberInFile(), itrk, esdTrack->GetTPCNcls()));
406             selected = kFALSE;
407           }
408           Float_t par[2], cov[3];
409           esdTrack->GetImpactParameters(par, cov);
410           if(IsCollision()){ // cuts on DCA
411             if(selected && TMath::Abs(par[0]) > fgkTrkDCAxy){ 
412               AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAxy[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[0])));
413               selected = kFALSE;
414             }
415             if(selected && TMath::Abs(par[1]) > fgkTrkDCAz){ 
416               AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAz[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[1])));
417               selected = kFALSE;
418             }
419           } else if(selected && fMCev && !fMCev->IsPhysicalPrimary(alab)){;
420             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Primary", fESDev->GetEventNumberInFile(), itrk));
421             selected = kFALSE;
422           }
423         }
424         if(fTrackCut && !fTrackCut->IsSelected(esdTrack)) selected = kFALSE;
425         if(selected){ 
426           fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
427           nBarrel++;
428         }
429       } else {
430         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
431         nKink++;
432       }
433     } else if((status&AliESDtrack::kTRDout) && !(status&AliESDtrack::kTRDin)){ 
434       fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
435       nSA++;
436     }
437     fTrackInfo->Delete("");
438   }
439
440 //   AliESDv0 *v0 = NULL;
441 //   for(Int_t iv0=0; iv0<fESD->GetNumberOfV0s(); iv0++){
442 //     if(!(v0 = fESD->GetV0(iv0))) continue;
443 //     fV0List->Add(new AliTRDv0Info(v0));
444 //   }
445
446   // Insert also MC tracks which are passing TRD where the track is not reconstructed
447   if(HasMCdata()){
448     AliDebug(10, "Output of the MC track map:");
449     for(Int_t itk = 0; itk < nTracksMC;  itk++) AliDebug(10, Form("trackMap[%d] = %s", itk, trackMap[itk] == kTRUE ? "TRUE" : "kFALSE"));
450   
451     for(Int_t itk = 0; itk < nTracksMC; itk++){
452       if(trackMap[itk]) continue;
453       AliMCParticle *mcParticle =  (AliMCParticle*) fMCev->GetTrack(TMath::Abs(itk));
454       Int_t fPdg = mcParticle->Particle()->GetPdgCode();
455       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
456       Int_t iref = 0; AliTrackReference *ref = NULL; 
457       Int_t nRefsTRD = 0;
458       new(fTrackInfo) AliTRDtrackInfo();
459       fTrackInfo->SetPDG(fPdg);
460       while(iref<nRefs){ // count TRD TR
461         Bool_t kIN(kFALSE);
462         ref = mcParticle->GetTrackReference(iref);
463         if(ref->LocalX() > fgkTPC && ref->LocalX() < fgkTRD){
464           fTrackInfo->AddTrackRef(ref);
465           nRefsTRD++;kIN=kTRUE;
466         }
467         AliDebug(4, Form("  trackRef[%2d] @ x[%7.3f] %s", iref, ref->LocalX(), kIN?"IN":"OUT"));
468         iref++;
469       }
470       if(!nRefsTRD){
471         // In this stage we at least require 1 hit inside TRD. What will be done with this tracks is a task for the 
472         // analysis job
473         fTrackInfo->Delete("");
474         continue;
475       }
476       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
477       fTrackInfo->SetLabel(itk);
478       if(DebugLevel() >= 1){
479         AliTRDtrackInfo info(*fTrackInfo);
480         (*DebugStream()) << "trackInfo"
481         << "TrackInfo.=" << &info
482         << "\n";
483         info.Delete("");
484       }
485       AliDebug(3, Form("Add MC track @ label[%d] nTRDrefs[%d].", itk, nRefsTRD));
486       // check where the track starts
487       ref = mcParticle->GetTrackReference(0);
488       if(ref->LocalX() < fgkITS){ 
489         fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
490         nBarrelMC++;
491       } else if(ref->LocalX() < fgkTPC) {
492         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
493         nKinkMC++;
494       } else if(nRefsTRD>6){
495         fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
496         nSAMC++;
497       }
498       fTrackInfo->Delete("");
499     }
500     delete[] trackMap;
501   }
502   AliDebug(2, Form(
503     "\nEv[%3d] Tracks: ESD[%d] MC[%d]\n"
504     "        TPCout[%d] TRDin[%d] TRDout[%d]\n"
505     "        Barrel[%3d+%3d=%3d] SA[%2d+%2d=%2d] Kink[%2d+%2d=%2d]"
506     ,(Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTracksESD, nTracksMC
507     , nTPC, nTRDin, nTRDout
508     ,nBarrel, nBarrelMC, fTracksBarrel->GetEntries()
509     ,nSA, nSAMC, fTracksSA->GetEntries()
510     ,nKink, nKinkMC, fTracksKink->GetEntries()
511   ));
512
513   PostData(kTracksBarrel, fTracksBarrel);
514   PostData(kTracksSA, fTracksSA);
515   PostData(kTracksKink, fTracksKink);
516   PostData(kEventInfo, fEventInfo);
517   PostData(kV0List, fV0List);
518 }
519
520 //____________________________________________________________________
521 void AliTRDinfoGen::SetTrigger(const Char_t *trigger)
522 {
523   if(!fEvTrigger) fEvTrigger = new TString(trigger);
524   else (*fEvTrigger) = trigger;
525 }
526
527 //____________________________________________________________________
528 TTreeSRedirector* AliTRDinfoGen::DebugStream()
529 {
530   if(!fDebugStream){
531     TDirectory *savedir = gDirectory;
532     fDebugStream = new TTreeSRedirector("TRD.DebugInfoGen.root");
533     savedir->cd();
534   }
535   return fDebugStream;
536 }
537
538