]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveTRDTrackList.cxx
update extra library list needed for class
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveTRDTrackList.cxx
1 // Author: Benjamin Hess   25/09/2008
2
3 /*************************************************************************
4  * Copyright (C) 2008-2009, 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 <../PWG1/TRD/AliTRDrecoTask.h>
81 #include <../PWG1/TRD/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   Char_t *libs[] = {"libANALYSIS.so", "libANALYSISalice.so", "libTENDER.so", "libPWG1.so"};
272   Int_t nlibs = static_cast<Int_t>(sizeof(libs)/sizeof(Char_t *));
273   for(Int_t ilib=0; ilib<nlibs; ilib++){
274     if(gSystem->Load(libs[ilib]) >= 0) continue;
275     AliError(Form("Fail loading %s.", libs[ilib]));
276     return;
277   }
278
279   AliTRDrecoTask *task = 0x0;
280   TList *fPlots = 0x0;
281   for(Int_t it=2; it<NTRDQATASKS; it++){
282     TClass c(fgkTRDtaskClassName[it]);
283     task = (AliTRDrecoTask*)c.New();
284     task->SetMCdata(kFALSE);
285     if(!(fPlots = task->GetPlotFunctors())){
286       //AliWarning(Form("No Plot functors defined for task \"%s\"", fgkTRDtaskClassName[it]));
287       delete task;
288       continue;
289     }
290     if(!(task->Histos())){
291       //AliWarning(Form("No Ref Histograms defined for task \"%s\"", fgkTRDtaskClassName[it]));
292       delete task;
293       continue;
294     }
295
296     // export task to CINT and add functions
297     gROOT->ProcessLine(Form("%s* %s = (%s*)0x%lx;", fgkTRDtaskClassName[it], task->GetName(), fgkTRDtaskClassName[it], (void*)task));
298     TIter iter(fPlots); TMethodCall *m = 0x0;
299     while((m = dynamic_cast<TMethodCall*>(iter()))){
300       AddMacroFast("", Form("%s->%s", task->GetName(), m->GetMethodName()), kSingleTrackHisto);
301     }
302   }
303 }
304
305 //______________________________________________________
306 Bool_t AliEveTRDTrackList::ApplyProcessMacros(const TList* selIterator, const TList* procIterator)
307 {
308   // Uses the procIterator (for the selected process macros) to apply the selected macros to the data.
309   // Returns kTRUE on success, otherwise kFALSE. If there no process macros selected, kTRUE is returned 
310   // (this is no error!).
311   // The single track process macros are applied to all selected tracks.
312   // The selIterator (for the selected selection macros) will be used to apply the correlated tracks selection
313   // macros to all track pairs (whereby BOTH tracks have to be selected, otherwise they will be skipped).
314   // All track pairs that have been selected by ALL correlated tracks selection macros will be processed by
315   // the correlated tracks process macros.
316
317   // No process macros need to be processed
318   if (procIterator->GetEntries() <= 0)  return kTRUE;
319
320   // Clear root
321   // A.B. gROOT->Reset();
322   
323   // Clear old data and re-allocate
324   if (fDataTree == 0x0){ 
325     TDirectory *cwd = gDirectory;
326     fDataTree = new TTreeSRedirector(Form("/tmp/TRD.TrackListMacroData_%s.root", gSystem->Getenv("USER")));
327     cwd->cd();
328   }
329   if (!fDataTree){
330     Error("Apply process macros", Form("File \"/tmp/TRD.TrackListMacroData_%s.root\" could not be accessed properly!", gSystem->Getenv("USER")));
331     return kFALSE;
332   }
333   
334   if (fDataFromMacroList != 0) {
335     fDataFromMacroList->Delete();
336     delete fDataFromMacroList;
337   }
338   fDataFromMacroList = new TList();
339   fDataFromMacroList->TCollection::SetOwner(kTRUE);
340
341   fHistoDataSelected = 0;
342
343
344   TMacroData* macro = 0;
345
346   Char_t** procCmds = 0;
347   AliEveTRDTrackListMacroType* mProcType = 0;
348   if (procIterator->GetEntries() > 0) {
349     procCmds = new Char_t*[procIterator->GetEntries()];
350     mProcType = new AliEveTRDTrackListMacroType[procIterator->GetEntries()];
351   }
352
353   Char_t** selCmds  = 0;
354   AliEveTRDTrackListMacroType* mSelType = 0;
355   if (selIterator->GetEntries() > 0) {
356     selCmds = new Char_t*[selIterator->GetEntries()];
357     mSelType = new AliEveTRDTrackListMacroType[selIterator->GetEntries()];
358   }
359   
360   Bool_t selectedByCorrSelMacro = kFALSE;
361
362   AliEveTRDTrackListMacroType macroType = kUnknown;
363   Int_t numHistoMacros = 0;
364   TH1** histos = 0;
365
366   AliEveTRDTrack* track1 = 0;
367   AliEveTRDTrack* track2 = 0;
368
369   // Collect the commands for each process macro and add them to "data-from-list"
370   for (Int_t i = 0; i < procIterator->GetEntries(); i++){
371     procCmds[i] = new Char_t[(fkMaxMacroPathNameLength + fkMaxApplyCommandLength)];
372     memset(procCmds[i], '\0', sizeof(Char_t) * (fkMaxMacroNameLength + fkMaxApplyCommandLength));
373
374     macro = (TMacroData*)fMacroList->GetValue(procIterator->At(i)->GetTitle());
375
376     if (!macro){
377       Error("Apply process macros", 
378         Form("Macro list is corrupted: Macro \"%s\" is not registered!", 
379         procIterator->At(i)->GetTitle()));
380       continue;
381     }
382
383 #ifdef ALIEVETRDTRACKLIST_DEBUG
384     printf("AliEveTRDTrackList: Checking process macro: %s\n", macro->GetName());
385 #endif 
386            
387     // Find the type of the process macro
388     macroType = macro->GetType();
389     if (macroType == kSingleTrackHisto || macroType == kCorrelTrackHisto){
390       mProcType[i] = macroType;
391       numHistoMacros++;
392       // Create the command 
393       sprintf(procCmds[i], macro->GetCmd());
394
395       // Add to "data-from-list" -> Mark as a histo macro with the substring "(histo macro)"
396       fDataFromMacroList->Add(new TObjString(Form("%s (histo macro)", macro->GetName())));
397     } else if (macroType == kSingleTrackAnalyse || macroType == kCorrelTrackAnalyse) {
398       mProcType[i] = macroType;
399       // Create the command 
400       sprintf(procCmds[i], macro->GetCmd());
401
402       // Add to "data-from-list"
403       fDataFromMacroList->Add(new TObjString(macro->GetName()));
404     } else {
405       Error("Apply process macros", 
406         Form("Macro list corrupted: Macro \"%s/%s.C\" is not registered as a process macro!", 
407         macro->GetPath(), macro->GetName()));
408       mProcType[i] = kUnknown;
409     } 
410   }  
411
412   // Collect the commands for each selection macro and add them to "data-from-list"
413   for (Int_t i = 0; i < selIterator->GetEntries(); i++){
414     selCmds[i] = new Char_t[(fkMaxMacroPathNameLength + fkMaxApplyCommandLength)];
415     memset(selCmds[i], '\0', sizeof(Char_t) * (fkMaxMacroNameLength + fkMaxApplyCommandLength));
416
417     macro = (TMacroData*)fMacroList->GetValue(selIterator->At(i)->GetTitle());
418
419     if (!macro){
420       Error("Apply process macros", 
421         Form("Macro list is corrupted: Macro \"%s\" is not registered!", 
422         selIterator->At(i)->GetTitle()));
423       continue;
424     }
425
426 #ifdef ALIEVETRDTRACKLIST_DEBUG
427     printf("AliEveTRDTrackList: Checking selection macro: %s\n", macro->GetName());
428 #endif
429        
430     // Find the type of the process macro
431     macroType = macro->GetType();
432
433     // Single track select macro
434     if (macroType == kSingleTrackSelect) {
435       // Has already been processed by ApplySTSelectionMacros(...)
436       mSelType[i] = macroType;         
437     }
438     // Correlated tracks select macro
439     else if (macroType == kCorrelTrackSelect) {
440       mSelType[i] = macroType;  
441  
442       // Create the command
443       sprintf(selCmds[i], macro->GetCmd());
444     } else {
445       Error("Apply process macros", 
446         Form("Macro list corrupted: Macro \"%s/%s.C\" is not registered as a selection macro!", 
447         macro->GetPath(), macro->GetName()));
448       mSelType[i] = kUnknown;
449     } 
450   }  
451
452   // Allocate memory for the histograms
453   if (numHistoMacros > 0)  histos = new TH1*[numHistoMacros];
454   for (Int_t i = 0; i < numHistoMacros; i++)  histos[i] = 0x0;
455
456
457   //////////////////////////////////
458   // WALK THROUGH THE LIST OF TRACKS
459   //////////////////////////////////     
460   for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter){
461     if(!(track1 = dynamic_cast<AliEveTRDTrack*>(*iter))) continue;
462
463     // Skip tracks that have not been selected
464     if (!track1->GetRnrState())  continue;
465     
466     // Cast to AliTRDtrackV1
467     gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)0x%xl;", track1));
468     gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_1 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
469
470     // Collect data for each macro
471     for (Int_t i = 0, histoIndex = 0; i < procIterator->GetEntries(); i++){
472       // Single track histo
473       if (mProcType[i] == kSingleTrackHisto){
474         histos[histoIndex++] = (TH1*)gROOT->ProcessLineSync(procCmds[i]);
475        // Correlated tracks histo
476       } else if (mProcType[i] == kCorrelTrackHisto) {
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         // The correlated tracks process macro will applied to all pairs that will be additionally selected by
481         // all correlated tracks selection macros.
482         TEveElement::List_i iter2 = iter;
483         iter2++;
484         for ( ; iter2 != this->EndChildren(); ++iter2)
485         {
486           if(!(track2 = dynamic_cast<AliEveTRDTrack*>(*iter2))) continue;
487
488           // Skip tracks that have not been selected
489           if (!track2->GetRnrState())  continue;
490       
491           // Cast to AliTRDtrackV1
492           gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)0x%xl;", track2));
493           gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_2 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
494
495           // Select track by default (so it will be processed, if there are no correlated tracks selection macros!)
496           selectedByCorrSelMacro = kTRUE;
497           for (Int_t j = 0; j < selIterator->GetEntries(); j++){
498             if (mSelType[j] == kCorrelTrackSelect){
499               selectedByCorrSelMacro = (Bool_t)gROOT->ProcessLineSync(selCmds[j]);
500               if (!selectedByCorrSelMacro)  break;
501             }
502           }       
503
504           // If the pair has not been selected by the correlated tracks selection macros, skip it!
505           if (!selectedByCorrSelMacro) continue;
506           
507           histos[histoIndex] = (TH1*)gROOT->ProcessLineSync(procCmds[i]);
508         }
509         histoIndex++;
510       }
511       // Single track analyse
512       else if (mProcType[i] == kSingleTrackAnalyse) {
513         // Create data pointers in CINT, execute the macro and get the data
514         gROOT->ProcessLineSync("Double_t* results = 0;");
515         gROOT->ProcessLineSync("Int_t n = 0;");
516         gROOT->ProcessLineSync(procCmds[i]);
517         Double_t* results = (Double_t*)gROOT->ProcessLineSync("results;");
518         Int_t nResults = (Int_t)gROOT->ProcessLineSync("n;");
519         
520         if (results == 0) {
521           Error("Apply macros", Form("Error reading data from macro \"%s\"", procIterator->At(i)->GetTitle()));
522           continue;
523         }
524         for (Int_t resInd = 0; resInd < nResults; resInd++){
525           (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << results[resInd] << (Char_t*)"\n";   
526         }
527
528         delete results;
529         results = 0;
530       }
531       // Correlated tracks analyse
532       else if (mProcType[i] == kCorrelTrackAnalyse){
533         // Loop over all pairs behind the current one - together with the other loop this will be a loop
534         // over all pairs. We have a pair of tracks, if and only if both tracks of the pair are selected (Rnr-state)
535         // and are not equal.
536         // The correlated tracks process macro will applied to all pairs that will be additionally selected by
537         // all correlated tracks selection macros.
538         TEveElement::List_i iter2 = iter;
539         iter2++;
540
541         for ( ; iter2 != this->EndChildren(); ++iter2) {
542           if(!(track2 = dynamic_cast<AliEveTRDTrack*>(*iter2))) continue;
543  
544           // Skip tracks that have not been selected
545           if (!track2->GetRnrState())  continue;
546     
547           // Cast to AliTRDtrackV1
548           gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)0x%xl;", track2));
549           gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_2 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
550
551           // Select track by default (so it will be processed, if there are no correlated tracks selection macros!)
552           selectedByCorrSelMacro = kTRUE;
553           for (Int_t j = 0; j < selIterator->GetEntries(); j++) {
554             if (mSelType[j] == kCorrelTrackSelect) {
555               selectedByCorrSelMacro = (Bool_t)gROOT->ProcessLineSync(selCmds[j]);
556               if (!selectedByCorrSelMacro)  break;
557             }
558           }       
559
560           // If the pair has not been selected by the correlated tracks selection macros, skip it!
561           if (!selectedByCorrSelMacro) continue;
562           
563           // Create data pointers in CINT, execute the macro and get the data
564           gROOT->ProcessLineSync("Double_t* results = 0;");
565           gROOT->ProcessLineSync("Int_t n = 0;");
566           gROOT->ProcessLineSync(procCmds[i]);
567           Double_t* results = (Double_t*)gROOT->ProcessLineSync("results;");
568           Int_t nResults = (Int_t)gROOT->ProcessLineSync("n;");
569      
570           if (results == 0) {
571             Error("Apply macros", Form("Error reading data from macro \"%s\"", procIterator->At(i)->GetTitle()));
572             continue;
573           }
574           for (Int_t resInd = 0; resInd < nResults; resInd++) {
575             (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << results[resInd] << (Char_t*)"\n";   
576           }
577
578           delete results;
579           results = 0;
580         }
581       }
582     }
583   }    
584
585   for (Int_t i = 0, histoIndex = 0; i < procIterator->GetEntries() && histoIndex < numHistoMacros; i++) {
586     if (mProcType[i] == kSingleTrackHisto || mProcType[i] == kCorrelTrackHisto) {
587       // Might be empty (e.g. no tracks have been selected)!
588       if (histos[histoIndex]) {
589         (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << histos[histoIndex] << (Char_t*)"\n";
590       }
591       histoIndex++;
592     }
593   }
594
595   if (fDataTree != 0) delete fDataTree;
596   fDataTree = 0;
597
598   if (procCmds != 0)  delete [] procCmds;
599   procCmds = 0;
600   if (mProcType != 0)  delete mProcType;
601   mProcType = 0;
602
603   if (selCmds != 0)  delete [] selCmds;
604   selCmds = 0;
605   if (mSelType != 0)  delete mSelType;
606   mSelType = 0;
607
608   if (histos != 0)  delete [] histos;
609   histos = 0;
610
611   // Clear root
612   // A.B. gROOT->Reset();
613   
614   // If there is data, select the first data set
615   if (procIterator->GetEntries() > 0) SETBIT(fHistoDataSelected, 0);
616
617   // Now the data is stored in "/tmp/TRD.TrackListMacroData_$USER.root"
618   // The editor will access this file to display the data
619   return kTRUE;
620 }
621
622 //______________________________________________________
623 void AliEveTRDTrackList::ApplySTSelectionMacros(const TList* iterator)
624 {
625   // Uses the iterator (for the selected selection macros) to apply the selected macros to the data.
626   // The rnr-states of the tracks are set according to the result of the macro calls (kTRUE, if all
627   // macros return kTRUE for this track, otherwise: kFALSE).
628   // "ST" stands for "single track". This means that only single track selection macros are applied.
629   // Correlated tracks selection macros will be used inside the call of ApplyProcessMacros(...)!
630
631   TMacroData* macro = 0;
632   AliEveTRDTrackListMacroType macroType = kUnknown;
633   AliEveTRDTrack* track1 = 0;
634   Bool_t selectedByMacro = kFALSE;
635
636   // Clear root
637   // A.B. gROOT->Reset();
638
639   // Select all tracks at first. A track is then deselect, if at least one selection macro
640   // returns kFALSE for this track
641   // Enable all tracks (Note: EnableListElements(..) will call "ElementChanged", which will cause unforeseen behavior!)
642   for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter) ((TEveElement*)(*iter))->SetRnrState(kTRUE);
643   SetRnrState(kTRUE);
644   
645   for (Int_t i = 0; i < iterator->GetEntries(); i++){
646     macro = (TMacroData*)fMacroList->GetValue(iterator->At(i)->GetTitle());
647
648     if (!macro){
649       Error("Apply selection macros", 
650             Form("Macro list is corrupted: Macro \"%s\" is not registered!", iterator->At(i)->GetTitle()));
651       continue;
652     }
653
654 #ifdef ALIEVETRDTRACKLIST_DEBUG
655     printf("AliEveTRDTrackList: Applying selection macro: %s\n", macro->GetName());
656 #endif
657     
658     // Determine macro type
659     macroType = macro->GetType();
660
661     // Single track select macro
662     if (macroType == kSingleTrackSelect){
663       // Walk through the list of tracks
664       for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter)
665       {
666         track1 = dynamic_cast<AliEveTRDTrack*>(*iter);
667
668         if (!track1) continue;
669
670         // If the track has already been deselected, nothing is to do here
671         if (!track1->GetRnrState()) continue;
672
673         // Cast to AliTRDtrackV1
674         gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)0x%xl;", track1));
675         gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_1 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
676         selectedByMacro = (Bool_t)gROOT->ProcessLineSync(macro->GetCmd());
677         track1->SetRnrState(selectedByMacro && track1->GetRnrState());               
678       }
679     }
680     // Correlated tracks select macro
681     else if (macroType == kCorrelTrackSelect){
682       // Will be processed in ApplyProcessMacros(...)
683       continue;
684     } else {
685       Error("Apply selection macros", 
686         Form("Macro list corrupted: Macro \"%s/%s.C\" is not registered as a selection macro!", 
687         macro->GetPath(), macro->GetName()));
688     } 
689   }
690
691   // Clear root
692   // A.B. gROOT->Reset();  
693 }
694
695 //______________________________________________________
696 AliEveTRDTrackList::AliEveTRDTrackListMacroType AliEveTRDTrackList::GetMacroType(const Char_t* name, Bool_t UseList) const
697 {
698   // Returns the type of the corresponding macro. 
699   // If "UseList" is kTRUE, the type will be looked up in the internal list (very fast). But if this list
700   // does not exist, you have to use kFALSE for this parameter. Then the type will be determined by the
701   // prototype! NOTE: It is assumed that the macro has been compiled! If not, the return value is not
702   // predictable, but normally will be kUnknown.
703   // Note: AddMacro(Fast) will update the internal list and RemoveMacros respectively.
704
705   AliEveTRDTrackListMacroType type = kUnknown;
706
707   // Re-do the check of the macro type
708   if (!UseList){
709     // Single track select macro or single track histo macro?
710     TFunction* f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*", kTRUE);
711     if (f != 0x0)
712     {
713       // Some additional check (is the parameter EXACTLY of the desired type?)
714       if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0)
715       {
716         // Single track select macro?
717         if (!strcmp(f->GetReturnTypeName(), "Bool_t")) 
718         { 
719           type = kSingleTrackSelect;     
720         }
721         // single track histo macro?
722         else if (!strcmp(f->GetReturnTypeName(), "TH1*"))
723         {
724           type = kSingleTrackHisto;
725         }
726       }
727     }
728     // Single track analyse macro?
729     else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*, Double_t*&, Int_t&", kTRUE)) 
730              != 0x0)
731     {
732       if (!strcmp(f->GetReturnTypeName(), "void"))
733       {
734         // Some additional check (are the parameters EXACTLY of the desired type?)
735         if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
736             strstr(f->GetMangledName(), "cODouble_tmUaNsP") != 0x0 &&
737             strstr(f->GetMangledName(), "cOInt_taNsP") != 0x0)
738         {
739           type = kSingleTrackAnalyse;
740         }
741       }
742     }    
743     // Correlated tracks select macro or correlated tracks histo macro?
744     else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*, const AliTRDtrackV1*", kTRUE)) 
745              != 0x0)
746     {
747       // Some additional check (is the parameter EXACTLY of the desired type?)
748       if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
749           strstr(f->GetMangledName(), "cOconstsPAliTRDtrackV1mUsP") != 0x0)
750       {
751         // Single track select macro?
752         if (!strcmp(f->GetReturnTypeName(), "Bool_t")) 
753         { 
754           type = kCorrelTrackSelect;     
755         }
756         // single track histo macro?
757         else if (!strcmp(f->GetReturnTypeName(), "TH1*"))
758         {
759           type = kCorrelTrackHisto;
760         }
761       }
762     }    
763     // Correlated tracks analyse macro?
764     else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, 
765                               "const AliTRDtrackV1*, const AliTRDtrackV1*, Double_t*&, Int_t&", kTRUE)) 
766              != 0x0)
767     {
768       if (!strcmp(f->GetReturnTypeName(), "void"))
769       {
770         // Some additional check (is the parameter EXACTLY of the desired type?)
771         if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
772             strstr(f->GetMangledName(), "cOconstsPAliTRDtrackV1mUsP") != 0x0 &&
773             strstr(f->GetMangledName(), "cODouble_tmUaNsP") != 0x0 &&
774             strstr(f->GetMangledName(), "cOInt_taNsP") != 0x0)
775         {
776           type = kCorrelTrackAnalyse;
777         }
778       }
779     }    
780   }
781   // Use list to look up the macro type
782   else
783   {
784     TMacroData* macro = 0;
785     macro = (TMacroData*)fMacroList->GetValue(name);
786     if (macro == 0)  return kUnknown; 
787     
788     type = macro->GetType();
789     switch (type)
790     {
791       case kSingleTrackSelect:
792       case kSingleTrackAnalyse:
793       case kSingleTrackHisto:
794       case kCorrelTrackSelect:
795       case kCorrelTrackAnalyse:
796       case kCorrelTrackHisto:      
797         break;
798     default:
799       type = kUnknown;
800       break;
801     }
802   }
803
804   return type;
805 }
806
807 //______________________________________________________
808 void AliEveTRDTrackList::RemoveSelectedMacros(const TList* iterator) 
809 {
810   // Uses the iterator (for the selected macros) to remove the selected macros from 
811   // the corresponding list.
812    
813   TObject* key = 0;
814   TPair*   entry = 0;
815   for (Int_t i = 0; i < iterator->GetEntries(); i++)
816   {
817     entry = (TPair*)fMacroList->FindObject(iterator->At(i)->GetTitle());
818
819     if (entry == 0)
820     {
821       Error("AliEveTRDTrackList::RemoveSelectedMacros", Form("Macro \"%s\" not found in list!", 
822                                                              iterator->At(i)->GetTitle()));
823       continue;
824     }
825     key = entry->Key();
826
827     if (key == 0)   
828     {
829       Error("AliEveTRDTrackList::RemoveSelectedMacros", Form("Key for macro \"%s\" not found in list!", 
830                                                              iterator->At(i)->GetTitle()));
831       continue;
832     }
833
834     // Key and value will be deleted, too, since fMacroList is the owner of them
835     Bool_t rem = fMacroList->DeleteEntry(key);
836
837     if (rem)
838     {
839 #ifdef ALIEVETRDTRACKLIST_DEBUG
840     printf("AliEveTRDTrackList::RemoveSelectedMacros(): Removed macro: %s\n", iterator->At(i)->GetTitle());
841 #endif
842     }
843     else
844     {
845       Error("AliEveTRDTrackList::RemoveSelectedMacros", Form("Macro \"%s\" could not be removed from the list!", 
846                                                              iterator->At(i)->GetTitle()));
847     }
848   }
849 }
850
851 //______________________________________________________
852 void AliEveTRDTrackList::UpdateTrackStyle(AliEveTRDTrack::AliEveTRDTrackState s, UChar_t ss)
853 {
854   // Updates the track style and sets this style for each track.
855
856   switch(s)
857   {
858     case AliEveTRDTrack::kSource:
859       SETBIT(fSelectedStyle, AliEveTRDTrack::kSource);
860       break;  
861     case AliEveTRDTrack::kPID:
862       CLRBIT(fSelectedStyle, AliEveTRDTrack::kSource);
863       switch(ss)
864       {
865       case AliTRDpidUtil::kLQ:
866         CLRBIT(fSelectedStyle, AliEveTRDTrack::kPID);
867         break;
868       case AliTRDpidUtil::kNN:
869         SETBIT(fSelectedStyle, AliEveTRDTrack::kPID);
870         break;
871       }
872       break;  
873     case AliEveTRDTrack::kTrackCosmics:
874       SETBIT(fSelectedStyle, AliEveTRDTrack::kTrackCosmics);
875       break;  
876     case AliEveTRDTrack::kTrackModel:
877       CLRBIT(fSelectedStyle, AliEveTRDTrack::kTrackCosmics);
878       switch(ss)
879       {
880       case AliEveTRDTrack::kRieman:
881         CLRBIT(fSelectedStyle, AliEveTRDTrack::kTrackModel);
882         break;
883       case AliEveTRDTrack::kKalman:
884         //AliWarning("Kalman fit under testing for the moment.");
885         SETBIT(fSelectedStyle, AliEveTRDTrack::kTrackModel);
886         break;
887       }
888       break;  
889   }
890
891
892   // Walk through the list of tracks     
893   AliEveTRDTrack* track = 0x0;
894   for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter) 
895   {
896     if (!(track = dynamic_cast<AliEveTRDTrack*>(*iter)))  continue;
897
898     track->SetStatus(fSelectedStyle);
899   }
900 }