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