]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveTRDTrackList.cxx
first attempt to introduce analysis of correlations between tracks (Ben)
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveTRDTrackList.cxx
1 // Uncomment to display debugging infos
2 //#define ALIEVETRDTRACKLIST_DEBUG
3
4 #include "AliEveTRDTrackList.h"
5
6 #include <AliTRDReconstructor.h>
7 #include <TFile.h>
8 #include <TFunction.h>
9 #include <TH1.h>
10 #include <TList.h>
11 #include <TMap.h>
12 #include <TObjString.h>
13 #include <TROOT.h>
14 #include <TSystem.h>
15 #include <TTree.h>
16 #include <TTreeStream.h>
17
18 ClassImp(AliEveTRDTrackList)
19
20 ///////////////////////////////////////////////////////////
21 /////////////   AliEveTRDTrackList ////////////////////////
22 ///////////////////////////////////////////////////////////
23 AliEveTRDTrackList::AliEveTRDTrackList(const Text_t* n, const Text_t* t, Bool_t doColor):
24   TEveElementList(n, t, doColor),
25   fMacroList(0),
26   fMacroSelList(0),
27   fDataFromMacroList(0),
28   fMacroTypes(0),
29   fDataTree(0),
30   fHistoDataSelected(0),
31   fMacroListSelected(0),
32   fMacroSelListSelected(0),
33   fSelectedTab(1),                              // Standard tab: "Apply macros" (index 1)
34   fSelectedStyle(0)
35 {
36   // Only accept childs of type AliEveTRDTrack
37   SetChildClass(AliEveTRDTrack::Class());
38
39   // Allocate memory for the lists and declare them as owners of their contents
40   fMacroList = new TList();
41   fMacroList->TCollection::SetOwner(kTRUE);
42   fMacroSelList = new TList();
43   fMacroSelList->TCollection::SetOwner(kTRUE);
44   fDataFromMacroList = new TList();
45   fDataFromMacroList->TCollection::SetOwner(kTRUE);
46
47   fMacroTypes = new TMap();
48   // Set map to owner of it's objects to delete them, if they are removed from the map
49   fMacroTypes->SetOwnerKeyValue(kTRUE, kTRUE);
50
51   // Set the build directory for AClic
52   gSystem->SetBuildDir("$HOME/.trdQArec");
53
54   AddStandardMacros();
55 }
56
57 //______________________________________________________
58 AliEveTRDTrackList::~AliEveTRDTrackList()
59 {
60   if (fMacroList != 0)
61   {
62     fMacroList->Delete();
63     delete fMacroList;
64     fMacroList = 0;
65   }
66   if (fMacroSelList != 0)
67   {
68     fMacroSelList->Delete();
69     delete fMacroSelList;
70     fMacroSelList = 0;
71   } 
72   if (fDataFromMacroList != 0)
73   {
74     fDataFromMacroList->Delete();
75     delete fDataFromMacroList;
76     fDataFromMacroList = 0;
77   } 
78   if (fDataTree != 0)
79   {
80     delete fDataTree;
81     fDataTree = 0;
82   } 
83   if (fMacroTypes != 0)
84   {
85     fMacroTypes->DeleteAll();
86     delete fMacroTypes;
87     fMacroTypes = 0;
88   }
89   // Note: gSystem->AccessPathName(...) returns kTRUE, if the access FAILED!
90   if(!gSystem->AccessPathName(Form("/tmp/TRD.TrackListMacroData_%s.root", gSystem->Getenv("USER")))) 
91     gSystem->Exec(Form("rm /tmp/TRD.TrackListMacroData_%s.root", gSystem->Getenv("USER")));
92 }
93
94 //______________________________________________________
95 Int_t AliEveTRDTrackList::AddMacro(const Char_t* path, const Char_t* nameC, Bool_t forceReload)
96 {
97   // First check the type of the macro:
98   // If it has the signature of a selection macro:
99   // Bool_t MacroName(AliTRDtrackV1*)
100   // it is assumed to be a selection macro.
101   // If it has the signature of a process macro:
102   // void MacroName(AliTRDtrackV1*, Double_t*&, Int_t&)
103   // it is assumed to be a process macro.
104   // In all other cases: Macro is rejected
105   
106   Char_t* entryName = MakeMacroEntry(path, nameC);
107
108   Char_t pathname[fkMaxMacroPathNameLength];
109   memset(pathname, '\0', sizeof(Char_t) * fkMaxMacroPathNameLength);
110
111   // Expand the path and create the pathname
112   Char_t* systemPath = gSystem->ExpandPathName(path);
113   sprintf(pathname, "%s/%s", systemPath, nameC);
114   delete systemPath;
115   systemPath = 0;
116
117   // Delete ".C" from filename
118   Char_t name[fkMaxMacroNameLength];
119   memset(name, '\0', sizeof(Char_t) * fkMaxMacroNameLength);
120   
121   for (UInt_t ind = 0; ind < fkMaxMacroNameLength && ind < strlen(nameC) - 2; ind++)  name[ind] = nameC[ind];
122
123   // Check, if files exists
124   FILE* fp = 0;
125
126   fp = fopen(pathname, "rb");
127   if (fp != 0)
128   {
129     fclose(fp);
130     fp = 0;
131   }
132   else
133   {
134     if (entryName != 0)  delete entryName;
135     entryName = 0;
136
137     return NOT_EXIST_ERROR;
138   }
139
140   // Clean up root, load the desired macro and then check the type of the macro
141   gROOT->Reset();
142  
143   if (forceReload)  gROOT->ProcessLineSync(Form(".L %s++", pathname));
144   else              gROOT->ProcessLineSync(Form(".L %s+", pathname));
145
146   AliEveTRDTrackListMacroType type = GetMacroType(entryName, kFALSE);
147
148   // Clean up again
149   gROOT->Reset();
150   
151   // Has not the correct signature!
152   if (type == kUnknown) 
153   {
154     if (entryName != 0)  delete entryName;
155     entryName = 0;
156     return SIGNATURE_ERROR;
157   }
158
159   Int_t returnValue = WARNING;
160
161   // Only add macro, if it is not already in the list
162   if ((type == kSingleTrackAnalyse || type == kSingleTrackHisto 
163       || type == kCorrelTrackAnalyse || type == kCorrelTrackHisto) && fMacroList->FindObject(entryName) == 0)
164   {
165     fMacroList->Add(new TObjString(entryName));
166     fMacroList->Sort();
167
168     fMacroTypes->Add(new TObjString(entryName), new TObjString(Form("%d", type)));
169
170     // We do not know, where the element has been inserted - deselect this list
171     fMacroListSelected = 0;
172
173     returnValue = SUCCESS;
174   }
175   else if ((type == kSingleTrackSelect || type == kCorrelTrackSelect) && fMacroSelList->FindObject(entryName) == 0)
176   {
177     fMacroSelList->Add(new TObjString(entryName));
178     fMacroSelList->Sort();
179
180     fMacroTypes->Add(new TObjString(entryName), new TObjString(Form("%d", type)));
181   
182     // We do not know, where the element has been inserted - deselect this list
183     fMacroSelListSelected = 0;
184     
185     returnValue = SUCCESS;
186   }
187   else  returnValue = WARNING;
188
189   if (entryName != 0)  delete entryName;
190   entryName = 0;
191
192   return returnValue;
193 }
194
195 //______________________________________________________
196 void AliEveTRDTrackList::AddMacroFast(const Char_t* entry, AliEveTRDTrackListMacroType type)
197 {
198   switch (type)
199   {
200     case kSingleTrackSelect:
201     case kCorrelTrackSelect:
202       fMacroSelList->Add(new TObjString(entry));
203       fMacroSelList->Sort();
204
205       fMacroTypes->Add(new TObjString(entry), new TObjString(Form("%d", type)));
206
207       // We do not know, where the element has been inserted - deselect this list
208       fMacroSelListSelected = 0;
209
210       break;
211     case kSingleTrackAnalyse:
212     case kSingleTrackHisto:
213     case kCorrelTrackAnalyse:
214     case kCorrelTrackHisto:
215       fMacroList->Add(new TObjString(entry));
216       fMacroList->Sort();
217
218       fMacroTypes->Add(new TObjString(entry), new TObjString(Form("%d", type)));
219
220       // We do not know, where the element has been inserted - deselect this list
221       fMacroListSelected = 0;
222       break;
223     default:
224       Error("AliEveTRDTrackList::AddMacroFast", Form("Unknown macro type for entry \"%s\"!", entry));
225       break;
226   }
227 }
228
229 //______________________________________________________
230 void AliEveTRDTrackList::AddMacroFast(const Char_t* path, const Char_t* name, AliEveTRDTrackListMacroType type)
231 {
232   Char_t* entry = MakeMacroEntry(path, name);
233   if (entry != 0)
234   {
235     AddMacroFast(entry, type);
236
237 #ifdef ALIEVETRDTRACKLIST_DEBUG
238     // Successfull add will only be displayed in debug mode
239     printf("#AliEveTRDTrackList::AddMacroFast: Added macro \"%s/%s\" to the corresponding list\n", path, name);
240 #endif
241     
242     delete entry;
243     entry = 0;
244   }
245   else
246   {
247     // Error will always be displayed
248     printf("#AliEveTRDTrackList::AddMacroFast: ERROR: Could not add macro \"%s/%s\" to the corresponding list\n", 
249            path, name);
250   }
251 }
252
253 //______________________________________________________
254 void AliEveTRDTrackList::AddStandardMacros()
255 {
256   // Add your standard macros here, e.g.: 
257   // To add a macro use:
258   // AddMacro("$(ALICE_ROOT)/myFolder", "myMacroName.C");
259   // -> If the file does not exist, nothing happens. So if you want to handle this,
260   // use the return value of AddMacro (NOT_EXIST_ERROR is returned, if file does not exist)
261   // (-> You can also check for other return values (see AddMacro(...)))
262   AddMacro("$(ALICE_ROOT)/TRD/qaRec/macros", "clusterSelection.C");
263   AddMacro("$(ALICE_ROOT)/TRD/qaRec/macros", "chargeDistr.C");
264   AddMacro("$(ALICE_ROOT)/TRD/qaRec/macros", "clusterResiduals.C");
265   AddMacro("$(ALICE_ROOT)/TRD/qaRec/macros", "PH.C");
266 }
267
268 //______________________________________________________
269 Bool_t AliEveTRDTrackList::ApplyProcessMacros(TList* iterator)
270 {
271   // No process macros need to be processed
272   if (iterator->GetEntries() <= 0)  return kTRUE;
273
274   // Clear root
275   gROOT->Reset();
276   
277   // Clear old data and re-allocate
278   if (fDataTree == 0) fDataTree = new TTreeSRedirector(Form("/tmp/TRD.TrackListMacroData_%s.root", 
279                                                             gSystem->Getenv("USER")));
280   if (!fDataTree)
281   {
282     Error("Apply process macros", Form("File \"/tmp/TRD.TrackListMacroData_%s.root\" could not be accessed properly!", 
283                                        gSystem->Getenv("USER")));
284     return kFALSE;
285   }
286   
287   if (fDataFromMacroList != 0)
288   {
289     fDataFromMacroList->Delete();
290     delete fDataFromMacroList;
291   }
292   fDataFromMacroList = new TList();
293   fDataFromMacroList->TCollection::SetOwner(kTRUE);
294
295   fHistoDataSelected = 0;
296
297
298   Char_t name[fkMaxMacroNameLength];
299   Char_t** cmds = new Char_t*[iterator->GetEntries()];
300   AliEveTRDTrackListMacroType* mType = new AliEveTRDTrackListMacroType[iterator->GetEntries()];
301
302   AliEveTRDTrackListMacroType macroType = kUnknown;
303   Int_t numHistoMacros = 0;
304   TH1** histos = 0;
305
306   AliEveTRDTrack* track1 = 0;
307   AliEveTRDTrack* track2 = 0;
308   TH1* returnedHist = 0x0;
309
310   // Collect the commands for each macro and add them to "data-from-list"
311   for (Int_t i = 0; i < iterator->GetEntries(); i++)
312   {
313     memset(name, '\0', sizeof(Char_t) * fkMaxMacroNameLength);
314     
315     cmds[i] = new Char_t[(fkMaxMacroPathNameLength + fkMaxApplyCommandLength)];
316     memset(cmds[i], '\0', sizeof(Char_t) * (fkMaxMacroNameLength + fkMaxApplyCommandLength));
317
318 #ifdef ALIEVETRDTRACKLIST_DEBUG
319     printf("AliEveTRDTrackList: Applying process macro: %s\n", iterator->At(i)->GetTitle());
320 #endif
321  
322     // Extract the name
323     sscanf(iterator->At(i)->GetTitle(), "%s (Path: %*s)", name);
324    
325     // Delete ".C" at the end 
326     // -> Note: Physical address pointer, do NOT delete. / Changes "name" as well!
327     Char_t* dotC = (Char_t*)strrchr(name, '.');
328     if (dotC != 0)
329     {
330       *dotC = '\0';
331       dotC++;
332       *dotC = '\0';
333     }
334        
335     // Find the type of the process macro
336     macroType = GetMacroType(iterator->At(i)->GetTitle(), kTRUE);
337     if (macroType == kSingleTrackHisto)
338     {
339       mType[i] = macroType;
340       numHistoMacros++;
341       // Create the command 
342       sprintf(cmds[i], "%s(automaticTrackV1_1);", name);
343
344       // Add to "data-from-list" -> Mark as a histo macro with the substring "(histo macro)"
345       fDataFromMacroList->Add(new TObjString(Form("%s (histo macro)", name)));
346     }
347     else if (macroType == kSingleTrackAnalyse)
348     {
349       mType[i] = macroType;
350       // Create the command 
351       sprintf(cmds[i], "%s(automaticTrackV1_1, results, n);", name);
352
353       // Add to "data-from-list"
354       fDataFromMacroList->Add(new TObjString(name));
355     }
356     else if (macroType == kCorrelTrackHisto)
357     {
358       mType[i] = macroType;
359       numHistoMacros++;
360       // Create the command 
361       sprintf(cmds[i], "%s(automaticTrackV1_1, automaticTrackV1_2);", name);
362
363       // Add to "data-from-list" -> Mark as a histo macro with the substring "(histo macro)"
364       fDataFromMacroList->Add(new TObjString(Form("%s (histo macro)", name)));
365     }
366     else if (macroType == kCorrelTrackAnalyse)
367     {
368       mType[i] = macroType;
369       // Create the command 
370       sprintf(cmds[i], "%s(automaticTrackV1_1, automaticTrackV1_2, results, n);", name);
371
372       // Add to "data-from-list"
373       fDataFromMacroList->Add(new TObjString(name));
374     }
375     else
376     {
377       Error("Apply process macros", 
378             Form("Process macro list corrupted: Macro \"%s\" is not registered as a process macro!", name));
379       mType[i] = kUnknown;
380     } 
381   }  
382
383   // Allocate memory for the histograms
384   if (numHistoMacros > 0)  histos = new TH1*[numHistoMacros];
385   for (Int_t i = 0; i < numHistoMacros; i++)  histos[i] = 0;
386   
387   // Walk through the list of tracks     
388   for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter)
389   {
390     track1 = dynamic_cast<AliEveTRDTrack*>(*iter);
391
392     if (!track1)  continue;
393     
394     // Skip tracks that have not been selected
395     if (!track1->GetRnrState())  continue;
396       
397     track1->ExportToCINT((Text_t*)"automaticTrack");
398     // Cast to AliTRDtrackV1
399     gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_1 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
400
401     // Collect data for each macro
402     for (Int_t i = 0, histoIndex = 0; i < iterator->GetEntries(); i++)
403     {
404       // Single track histo
405       if (mType[i] == kSingleTrackHisto)
406       {
407         returnedHist = (TH1*)gROOT->ProcessLineSync(cmds[i]);
408         if (returnedHist != 0x0)
409         {
410           if (histos[histoIndex] == 0)  histos[histoIndex] = (TH1*)gROOT->ProcessLineSync(cmds[i]);
411           else  histos[histoIndex]->Add((const TH1*)gROOT->ProcessLineSync(cmds[i]));
412       
413           delete returnedHist;
414           returnedHist = 0;
415         }
416         histoIndex++;
417       }
418       // Correlated tracks histo
419       else if (mType[i] == kCorrelTrackHisto)
420       {
421         // Loop over all pairs behind the current one - together with the other loop this will be a loop
422         // over all pairs. We have a pair of tracks, if and only if both tracks of the pair are selected (Rnr-state)
423         // and are not equal.
424         TEveElement::List_i iter2 = iter;
425         iter2++;
426         for ( ; iter2 != this->EndChildren(); ++iter2)
427         {
428           track2 = dynamic_cast<AliEveTRDTrack*>(*iter2);
429
430           if (!track2)  continue;
431     
432           // Skip tracks that have not been selected
433           if (!track2->GetRnrState())  continue;
434       
435           track2->ExportToCINT((Text_t*)"automaticTrack");
436           // Cast to AliTRDtrackV1
437           gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_2 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
438
439           returnedHist = (TH1*)gROOT->ProcessLineSync(cmds[i]);
440           if (returnedHist != 0x0)
441           {
442             if (histos[histoIndex] == 0)  histos[histoIndex] = (TH1*)gROOT->ProcessLineSync(cmds[i]);
443             else  histos[histoIndex]->Add((const TH1*)gROOT->ProcessLineSync(cmds[i]));
444         
445             delete returnedHist;
446             returnedHist = 0;
447           }
448         }
449         histoIndex++;
450       }
451       // Single track analyse
452       else if (mType[i] == kSingleTrackAnalyse)
453       {
454         // Create data pointers in CINT, execute the macro and get the data
455         gROOT->ProcessLineSync("Double_t* results = 0;");
456         gROOT->ProcessLineSync("Int_t n = 0;");
457         gROOT->ProcessLineSync(cmds[i]);
458         Double_t* results = (Double_t*)gROOT->ProcessLineSync("results;");
459         Int_t nResults = (Int_t)gROOT->ProcessLineSync("n;");
460         
461         if (results == 0)
462         {
463           Error("Apply macros", Form("Error reading data from macro \"%s\"", iterator->At(i)->GetTitle()));
464           continue;
465         }
466         for (Int_t resInd = 0; resInd < nResults; resInd++)
467         {
468           (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << results[resInd] << (Char_t*)"\n";   
469         }
470
471         delete results;
472         results = 0;
473       }
474       // Correlated tracks analyse
475       else if (mType[i] == kCorrelTrackAnalyse)
476       {
477         // Loop over all pairs behind the current one - together with the other loop this will be a loop
478         // over all pairs. We have a pair of tracks, if and only if both tracks of the pair are selected (Rnr-state)
479         // and are not equal.
480         TEveElement::List_i iter2 = iter;
481         iter2++;
482         for ( ; iter2 != this->EndChildren(); ++iter2)
483         {
484           track2 = dynamic_cast<AliEveTRDTrack*>(*iter2);
485
486           if (!track2)  continue;
487     
488           // Skip tracks that have not been selected
489           if (!track2->GetRnrState())  continue;
490       
491           track2->ExportToCINT((Text_t*)"automaticTrack");
492           // Cast to AliTRDtrackV1
493           gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_2 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
494
495           // Create data pointers in CINT, execute the macro and get the data
496           gROOT->ProcessLineSync("Double_t* results = 0;");
497           gROOT->ProcessLineSync("Int_t n = 0;");
498           gROOT->ProcessLineSync(cmds[i]);
499           Double_t* results = (Double_t*)gROOT->ProcessLineSync("results;");
500           Int_t nResults = (Int_t)gROOT->ProcessLineSync("n;");
501           
502           if (results == 0)
503           {
504             Error("Apply macros", Form("Error reading data from macro \"%s\"", iterator->At(i)->GetTitle()));
505             continue;
506           }
507           for (Int_t resInd = 0; resInd < nResults; resInd++)
508           {
509             (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << results[resInd] << (Char_t*)"\n";   
510           }
511
512           delete results;
513           results = 0;
514         }
515       }
516     }
517   }    
518
519   for (Int_t i = 0, histoIndex = 0; i < iterator->GetEntries() && histoIndex < numHistoMacros; i++)
520   {
521     if (mType[i] == kSingleTrackHisto || mType[i] == kCorrelTrackHisto)
522     {
523       // Might be empty (e.g. no tracks have been selected)!
524       if (histos[histoIndex] != 0)
525       {
526         (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << histos[histoIndex] << (Char_t*)"\n";
527       }
528       histoIndex++;
529     }
530   }
531
532   if (fDataTree != 0) delete fDataTree;
533   fDataTree = 0;
534
535   if (cmds != 0)  delete [] cmds;
536   if (mType != 0)  delete mType;
537   mType = 0;
538
539   if (histos != 0)  delete [] histos;
540   histos = 0;
541
542   // Clear root
543   gROOT->Reset();
544   
545   // If there is data, select the first data set
546   if (iterator->GetEntries() > 0) SETBIT(fHistoDataSelected, 0);
547
548   // Now the data is stored in "/tmp/TRD.TrackListMacroData_$USER.root"
549   // The editor will access this file to display the data
550   return kTRUE;
551 }
552
553 //______________________________________________________
554 void AliEveTRDTrackList::ApplySelectionMacros(TList* iterator)
555 {
556   Char_t name[fkMaxMacroNameLength];
557   Char_t cmd[(fkMaxMacroNameLength + fkMaxApplyCommandLength)];
558
559   AliEveTRDTrackListMacroType macroType = kUnknown;
560   AliEveTRDTrack* track1 = 0;
561   AliEveTRDTrack* track2 = 0;
562   Bool_t selectedByMacro = kFALSE;
563
564   // Clear root
565   gROOT->Reset();
566
567   // Select all tracks at first. A track is then deselect, if at least one selection macro
568   // returns kFALSE for this track
569   // Enable all tracks (Note: EnableListElements(..) will call "ElementChanged", which will cause unforeseen behavior!)
570   for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter)
571   {
572     ((TEveElement*)(*iter))->SetRnrState(kTRUE);
573   }
574   SetRnrState(kTRUE);
575   
576   for (Int_t i = 0; i < iterator->GetEntries(); i++)
577   {
578
579     memset(name, '\0', sizeof(Char_t) * fkMaxMacroNameLength);
580     memset(cmd, '\0', sizeof(Char_t) * (fkMaxMacroNameLength + fkMaxApplyCommandLength));
581
582 #ifdef ALIEVETRDTRACKLIST_DEBUG
583     printf("AliEveTRDTrackList: Applying selection macro: %s\n", iterator->At(i)->GetTitle());
584 #endif
585     
586     // Extract the name
587     sscanf(iterator->At(i)->GetTitle(), "%s (Path: %*s)", name);
588     // Delete ".C" at the end 
589     // -> Note: Physical address pointer, do NOT delete. / Changes "name" as well!
590     Char_t* dotC = (Char_t*)strrchr(name, '.');
591     if (dotC != 0)
592     {
593       *dotC = '\0';
594       dotC++;
595       *dotC = '\0';
596     }
597
598     // Determine macro type
599     macroType = GetMacroType(iterator->At(i)->GetTitle(), kTRUE);
600
601     // Single track select macro
602     if (macroType == kSingleTrackSelect)
603     {
604       // Create the command
605       sprintf(cmd, "%s(automaticTrackV1);", name);
606
607       // Walk through the list of tracks
608       for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter)
609       {
610         track1 = dynamic_cast<AliEveTRDTrack*>(*iter);
611
612         if (!track1) continue;
613
614         // If the track has already been deselected, nothing is to do here
615         if (!track1->GetRnrState()) continue;
616
617         track1->ExportToCINT((Text_t*)"automaticTrack");
618         // Cast to AliTRDtrackV1
619         gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
620         selectedByMacro = (Bool_t)gROOT->ProcessLineSync(cmd);
621         track1->SetRnrState(selectedByMacro && track1->GetRnrState());               
622       }
623     }
624     // Correlated tracks select macro
625     else if (macroType == kCorrelTrackSelect)
626     {
627       // Create the command
628       sprintf(cmd, "%s(automaticTrackV1_1, automaticTrackV1_2);", name);
629
630       // Walk through the list of tracks and loop over all possible track pairs
631       for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter)
632       {
633         track1 = dynamic_cast<AliEveTRDTrack*>(*iter);
634
635         if (!track1) continue;
636
637         // If the track has already been deselected, nothing is to do here
638         if (!track1->GetRnrState()) continue;
639         
640         track1->ExportToCINT((Text_t*)"automaticTrack");
641         // Cast to AliTRDtrackV1
642         gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_1 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
643
644         // Correlate the track with each other track (except for itself!)
645         for (TEveElement::List_i iter2 = this->BeginChildren() ; iter2 != this->EndChildren(); ++iter2)
646         {
647           // Do not correlate track with itself
648           if (iter == iter2)  continue;
649
650           track2 = dynamic_cast<AliEveTRDTrack*>(*iter2);
651
652           if (!track2) continue;
653         
654           track2->ExportToCINT((Text_t*)"automaticTrack");
655           // Cast to AliTRDtrackV1
656           gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_2 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
657           selectedByMacro = (Bool_t)gROOT->ProcessLineSync(cmd);
658           track1->SetRnrState(selectedByMacro && track1->GetRnrState());               
659
660           // If track has not been selected by this call, it is switched off and we are done with this track
661           if (!selectedByMacro) break;
662         }
663       }
664     }
665     else
666     {
667       Error("Apply selection macros", 
668             Form("Selection macro list corrupted: Macro \"%s\" is not registered as a selection macro!", name));
669     } 
670   }
671
672   // Clear root
673   gROOT->Reset();  
674 }
675
676 //______________________________________________________
677 AliEveTRDTrackList::AliEveTRDTrackListMacroType AliEveTRDTrackList::GetMacroType(const Char_t* entry, Bool_t UseList)
678 {
679   AliEveTRDTrackListMacroType type = kUnknown;
680
681   // Re do the check of the macro type
682   if (!UseList)
683   {
684     Char_t name[fkMaxMacroNameLength];
685   
686     memset(name, '\0', sizeof(Char_t) * fkMaxMacroNameLength);
687
688     // Extract the name
689     sscanf(entry, "%s (Path: %*s)", name);
690    
691     // Delete ".C" at the end 
692     // -> Note: Physical address pointer, do NOT delete. / Changes "name" as well!
693     Char_t* dotC = (Char_t*)strrchr(name, '.');
694     if (dotC != 0)
695     {
696       *dotC = '\0';
697       dotC++;
698       *dotC = '\0';
699     }
700
701     // Single track select macro or single track histo macro?
702     TFunction* f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*", kTRUE);
703     if (f != 0x0)
704     {
705       // Some additional check (is the parameter EXACTLY of the desired type?)
706       if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0)
707       {
708         // Single track select macro?
709         if (!strcmp(f->GetReturnTypeName(), "Bool_t")) 
710         { 
711           type = kSingleTrackSelect;     
712         }
713         // single track histo macro?
714         else if (!strcmp(f->GetReturnTypeName(), "TH1*"))
715         {
716           type = kSingleTrackHisto;
717         }
718       }
719     }
720     // Single track analyse macro?
721     else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*, Double_t*&, Int_t&", kTRUE)) 
722              != 0x0)
723     {
724       if (!strcmp(f->GetReturnTypeName(), "void"))
725       {
726         // Some additional check (are the parameters EXACTLY of the desired type?)
727         if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
728             strstr(f->GetMangledName(), "cODouble_tmUaNsP") != 0x0 &&
729             strstr(f->GetMangledName(), "cOInt_taNsP") != 0x0)
730         {
731           type = kSingleTrackAnalyse;
732         }
733       }
734     }    
735     // Correlated tracks select macro or correlated tracks histo macro?
736     else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*, const AliTRDtrackV1*", kTRUE)) 
737              != 0x0)
738     {
739       // Some additional check (is the parameter EXACTLY of the desired type?)
740       if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
741           strstr(f->GetMangledName(), "cOconstsPAliTRDtrackV1mUsP") != 0x0)
742       {
743         // Single track select macro?
744         if (!strcmp(f->GetReturnTypeName(), "Bool_t")) 
745         { 
746           type = kCorrelTrackSelect;     
747         }
748         // single track histo macro?
749         else if (!strcmp(f->GetReturnTypeName(), "TH1*"))
750         {
751           type = kCorrelTrackHisto;
752         }
753       }
754     }    
755     // Correlated tracks analyse macro?
756     else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, 
757                               "const AliTRDtrackV1*, const AliTRDtrackV1*, Double_t*&, Int_t&", kTRUE)) 
758              != 0x0)
759     {
760       if (!strcmp(f->GetReturnTypeName(), "void"))
761       {
762         // Some additional check (is the parameter EXACTLY of the desired type?)
763         if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
764             strstr(f->GetMangledName(), "cOconstsPAliTRDtrackV1mUsP") != 0x0 &&
765             strstr(f->GetMangledName(), "cODouble_tmUaNsP") != 0x0 &&
766             strstr(f->GetMangledName(), "cOInt_taNsP") != 0x0)
767         {
768           type = kCorrelTrackAnalyse;
769         }
770       }
771     }    
772   }
773   // Use list to look up the macro type
774   else
775   {
776     TObjString* objEntry = 0;
777     objEntry = (TObjString*)fMacroTypes->GetValue(entry);
778     if (objEntry == 0)  return kUnknown; 
779     
780     type = (AliEveTRDTrackListMacroType)objEntry->GetString().Atoi();
781     switch (type)
782     {
783       case kSingleTrackSelect:
784       case kSingleTrackAnalyse:
785       case kSingleTrackHisto:
786       case kCorrelTrackSelect:
787       case kCorrelTrackAnalyse:
788       case kCorrelTrackHisto:      
789         break;
790     default:
791       type = kUnknown;
792       break;
793     }
794   }
795
796   return type;
797 }
798
799 //______________________________________________________
800 Char_t* AliEveTRDTrackList::MakeMacroEntry(const Char_t* path, const Char_t* name)
801 {
802   Char_t* entry = new Char_t[(fkMaxMacroPathNameLength + 30)];
803   memset(entry, '\0', sizeof(Char_t) * (fkMaxMacroPathNameLength + 30));
804
805   Char_t* systemPath = gSystem->ExpandPathName(path);
806   sprintf(entry, "%s (Path: %s)", name, systemPath);
807   delete systemPath;
808   systemPath = 0;
809
810   return entry;
811 }
812
813 //______________________________________________________
814 void AliEveTRDTrackList::RemoveProcessMacros(TList* iterator) 
815 {
816   TObjString* obj = 0;
817   for (Int_t i = 0; i < iterator->GetEntries(); i++)
818   {
819     // Key and value will be deleted, too, since fMacroTypes is the owner of them
820     fMacroTypes->DeleteEntry(fMacroTypes->FindObject(iterator->At(i)->GetTitle()));
821
822     obj = (TObjString*)fMacroList->Remove(fMacroList->FindObject(iterator->At(i)->GetTitle()));   
823     if (obj != 0) delete obj;
824   }
825   obj = 0;
826 }
827
828 //______________________________________________________
829 void AliEveTRDTrackList::RemoveSelectionMacros(TList* iterator) 
830 {
831   TObjString* obj = 0;
832   for (Int_t i = 0; i < iterator->GetEntries(); i++)
833   {
834     // Key and value will be deleted, too, since fMacroTypes is the owner of them
835     fMacroTypes->DeleteEntry(fMacroTypes->FindObject(iterator->At(i)->GetTitle()));
836
837     obj = (TObjString*)fMacroSelList->Remove(fMacroSelList->FindObject(iterator->At(i)->GetTitle()));
838     if (obj != 0) delete obj;
839   }
840   obj = 0;
841 }
842
843 //______________________________________________________
844 void AliEveTRDTrackList::UpdateTrackStyle(AliEveTRDTrack::AliEveTRDTrackState s, UChar_t ss)
845 {
846   switch(s)
847   {
848     case AliEveTRDTrack::kSource:
849       SETBIT(fSelectedStyle, AliEveTRDTrack::kSource);
850       break;  
851     case AliEveTRDTrack::kPID:
852       CLRBIT(fSelectedStyle, AliEveTRDTrack::kSource);
853       switch(ss)
854       {
855       case AliTRDReconstructor::kLQPID:
856         CLRBIT(fSelectedStyle, AliEveTRDTrack::kPID);
857         break;
858       case AliTRDReconstructor::kNNPID:
859         SETBIT(fSelectedStyle, AliEveTRDTrack::kPID);
860         break;
861       }
862       break;  
863     case AliEveTRDTrack::kTrackCosmics:
864       SETBIT(fSelectedStyle, AliEveTRDTrack::kTrackCosmics);
865       break;  
866     case AliEveTRDTrack::kTrackModel:
867       CLRBIT(fSelectedStyle, AliEveTRDTrack::kTrackCosmics);
868       switch(ss)
869       {
870       case AliEveTRDTrack::kRieman:
871         CLRBIT(fSelectedStyle, AliEveTRDTrack::kTrackModel);
872         break;
873       case AliEveTRDTrack::kKalman:
874         AliWarning("Kalman fit under testing for the moment.");
875         //SETBIT(fSelectedStyle, AliEveTRDTrack::kTrackModel);
876         break;
877       }
878       break;  
879   }
880
881
882   // Walk through the list of tracks     
883   AliEveTRDTrack* track = 0x0;
884   for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter) 
885   {
886     if (!(track = dynamic_cast<AliEveTRDTrack*>(*iter)))  continue;
887
888     track->SetStatus(fSelectedStyle);
889   }
890 }