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