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