1 // Author: Benjamin Hess 29/01/2010
3 /*************************************************************************
4 * Copyright (C) 2009-2010, Alexandru Bercuci, Benjamin Hess. *
5 * All rights reserved. *
6 *************************************************************************/
9 //////////////////////////////////////////////////////////////////////////
11 // AliEveTRDTrackList //
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. //
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*); //
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*); //
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 //////////////////////////////////////////////////////////////////////////
60 // Uncomment to display debugging infos
61 //#define ALIEVETRDTRACKLIST_DEBUG
64 #include <TFunction.h>
68 #include <TObjString.h>
72 #include <TTreeStream.h>
73 #include <TMethodCall.h>
75 #include <AliTRDReconstructor.h>
77 #include <EveDet/AliEveTRDTrackList.h>
78 #include <EveDet/AliEveTRDTrackListEditor.h>
80 #include <../PWG1/TRD/AliTRDrecoTask.h>
81 #include <../PWG1/TRD/AliTRDpwg1Helper.h>
83 ClassImp(AliEveTRDTrackList)
85 ///////////////////////////////////////////////////////////
86 ///////////// AliEveTRDTrackList ////////////////////////
87 ///////////////////////////////////////////////////////////
88 AliEveTRDTrackList::AliEveTRDTrackList(const Text_t* n, const Text_t* t, Bool_t doColor):
89 TEveElementList(n, t, doColor),
91 fDataFromMacroList(0x0),
94 fHistoDataSelected(0),
95 fMacroListSelected(0),
96 fSelectedTab(1), // Standard tab: "Apply macros" (index 1)
99 // Creates the AliEveTRDTrackList.
101 // Only accept childs of type AliEveTRDTrack
102 SetChildClass(AliEveTRDTrack::Class());
104 // Allocate memory for the lists and declare them as owners of their contents
105 fDataFromMacroList = new TList();
106 fDataFromMacroList->TCollection::SetOwner(kTRUE);
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);
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")));
116 AddStandardContent();
119 //______________________________________________________
120 AliEveTRDTrackList::~AliEveTRDTrackList()
122 // Frees allocated memory (lists etc.).
124 // Let the editor know that the list will be destroyed -> The editor will save the data
127 fEditor->SaveMacroList(fMacroList);
131 if (fDataFromMacroList != 0)
133 fDataFromMacroList->Delete();
134 delete fDataFromMacroList;
135 fDataFromMacroList = 0;
144 fMacroList->DeleteAll();
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")));
153 //______________________________________________________
154 Int_t AliEveTRDTrackList::AddMacro(const Char_t* path, const Char_t* nameC, Bool_t forceReload)
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:
161 // Bool_t YourMacro(const AliTRDtrackV1*)
162 // Bool_t YourMacro(const AliTRDtrackV1*, const AliTRDtrackV1*)
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*)
171 Char_t pathname[fkMaxMacroPathNameLength];
172 memset(pathname, '\0', sizeof(Char_t) * fkMaxMacroPathNameLength);
174 // Expand the path and create the pathname
175 Char_t* systemPath = gSystem->ExpandPathName(path);
176 snprintf(pathname, fkMaxMacroPathNameLength, "%s/%s", systemPath, nameC);
180 // Delete ".C" from filename
181 Char_t name[fkMaxMacroNameLength];
182 memset(name, '\0', sizeof(Char_t) * fkMaxMacroNameLength);
184 for (UInt_t ind = 0; ind < fkMaxMacroNameLength && ind < strlen(nameC) - 2; ind++) name[ind] = nameC[ind];
186 // Check, if files exists
188 if((fp = fopen(pathname, "rb"))){
191 } else return NOT_EXIST_ERROR;
193 // Clean up root, load the desired macro and then check the type of the macro
194 // A.B. gROOT->Reset();
196 gROOT->ProcessLineSync(Form(".L %s+%c", pathname, forceReload ? '+' : ' '));
198 // We need this line... otherwise, in some cases, there will be problems concerning ACLIC
199 gROOT->ProcessLineSync(Form(".L %s", pathname));
201 AliEveTRDTrackListMacroType type = GetMacroType(name, kFALSE);
204 // A.B. gROOT->Reset();
206 // Has not the correct signature!
207 if (type == kUnknown) return SIGNATURE_ERROR;
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;
217 //______________________________________________________
218 Bool_t AliEveTRDTrackList::AddMacroFast(const Char_t* path, const Char_t* name, AliEveTRDTrackListMacroType type)
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;
225 Bool_t success = kFALSE;
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));
237 // We do not know, where the element has been inserted - deselect this list
238 fMacroListSelected = 0;
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);
250 // Error will always be displayed
251 printf("AliEveTRDTrackList::AddMacroFast: ERROR: Could not add macro \"%s/%s\" to the corresponding list\n",
262 //______________________________________________________
263 void AliEveTRDTrackList::AddStandardContent()
265 // Adds standard macros to the macro list.
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(...)))
274 const Char_t *libs[] = {"libANALYSIS.so", "libANALYSISalice.so", "libTENDER.so", "libPWG1.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]));
282 const Char_t *fgkTRDPWG1taskClassName[AliTRDpwg1Helper::kNTRDQATASKS] = {
291 AliTRDrecoTask *task(NULL);
293 for(Int_t it=2; it<AliTRDpwg1Helper::kNTRDQATASKS; it++){
294 TClass c(fgkTRDPWG1taskClassName[it]);
295 task = (AliTRDrecoTask*)c.New();
296 task->SetMCdata(kFALSE);
297 if(!(fPlots = task->GetPlotFunctors())){
298 //AliWarning(Form("No Plot functors defined for task \"%s\"", fgkTRDtaskClassName[it]));
302 if(!(task->Histos())){
303 //AliWarning(Form("No Ref Histograms defined for task \"%s\"", fgkTRDtaskClassName[it]));
308 // export task to CINT and add functions
309 gROOT->ProcessLine(Form("%s* %s = (%s*)%p;", fgkTRDPWG1taskClassName[it], task->GetName(), fgkTRDPWG1taskClassName[it], (void*)task));
310 TIter iter(fPlots); TMethodCall *m = 0x0;
311 while((m = dynamic_cast<TMethodCall*>(iter()))){
312 AddMacroFast("", Form("%s->%s", task->GetName(), m->GetMethodName()), kSingleTrackHisto);
318 //______________________________________________________
319 Bool_t AliEveTRDTrackList::ApplyProcessMacros(const TList* selIterator, const TList* procIterator)
321 // Uses the procIterator (for the selected process macros) to apply the selected macros to the data.
322 // Returns kTRUE on success, otherwise kFALSE. If there no process macros selected, kTRUE is returned
323 // (this is no error!).
324 // The single track process macros are applied to all selected tracks.
325 // The selIterator (for the selected selection macros) will be used to apply the correlated tracks selection
326 // macros to all track pairs (whereby BOTH tracks have to be selected, otherwise they will be skipped).
327 // All track pairs that have been selected by ALL correlated tracks selection macros will be processed by
328 // the correlated tracks process macros.
330 // No process macros need to be processed
331 if (procIterator->GetEntries() <= 0) return kTRUE;
334 // A.B. gROOT->Reset();
336 // Clear old data and re-allocate
338 TDirectory *cwd = gDirectory;
339 fDataTree = new TTreeSRedirector(Form("/tmp/TRD.TrackListMacroData_%s.root", gSystem->Getenv("USER")));
343 Error("Apply process macros", "File \"/tmp/TRD.TrackListMacroData_%s.root\" could not be accessed properly!", gSystem->Getenv("USER"));
347 if (fDataFromMacroList != 0) {
348 fDataFromMacroList->Delete();
349 delete fDataFromMacroList;
351 fDataFromMacroList = new TList();
352 fDataFromMacroList->TCollection::SetOwner(kTRUE);
354 fHistoDataSelected = 0;
357 TMacroData* macro = 0;
359 TString* procCmds = new TString[procIterator->GetEntries()];
360 AliEveTRDTrackListMacroType* mProcType = new AliEveTRDTrackListMacroType[procIterator->GetEntries()];
362 TString* selCmds(NULL);
363 AliEveTRDTrackListMacroType* mSelType(NULL);
364 if (selIterator->GetEntries() > 0) {
365 selCmds = new TString[selIterator->GetEntries()];
366 mSelType = new AliEveTRDTrackListMacroType[selIterator->GetEntries()];
369 Bool_t selectedByCorrSelMacro = kFALSE;
371 AliEveTRDTrackListMacroType macroType = kUnknown;
372 Int_t numHistoMacros = 0;
375 AliEveTRDTrack* track1(NULL);
376 AliEveTRDTrack* track2(NULL);
378 // Collect the commands for each process macro and add them to "data-from-list"
379 for (Int_t i = 0; i < procIterator->GetEntries(); i++){
380 macro = (TMacroData*)fMacroList->GetValue(procIterator->At(i)->GetTitle());
383 Error("Apply process macros",
384 "Macro list is corrupted: Macro \"%s\" is not registered!",
385 procIterator->At(i)->GetTitle());
389 #ifdef ALIEVETRDTRACKLIST_DEBUG
390 printf("AliEveTRDTrackList: Checking process macro: %s\n", macro->GetName());
393 // Find the type of the process macro
394 macroType = macro->GetType();
395 if (macroType == kSingleTrackHisto || macroType == kCorrelTrackHisto){
396 mProcType[i] = macroType;
398 // Create the command
399 procCmds[i] = macro->GetCmd();
401 // Add to "data-from-list" -> Mark as a histo macro with the substring "(histo macro)"
402 fDataFromMacroList->Add(new TObjString(Form("%s (histo macro)", macro->GetName())));
403 } else if (macroType == kSingleTrackAnalyse || macroType == kCorrelTrackAnalyse) {
404 mProcType[i] = macroType;
405 // Create the command
406 procCmds[i] = macro->GetCmd();
408 // Add to "data-from-list"
409 fDataFromMacroList->Add(new TObjString(macro->GetName()));
411 Error("Apply process macros",
412 "Macro list corrupted: Macro \"%s/%s.C\" is not registered as a process macro!",
413 macro->GetPath(), macro->GetName());
414 mProcType[i] = kUnknown;
418 // Collect the commands for each selection macro and add them to "data-from-list"
419 for (Int_t i = 0; i < selIterator->GetEntries(); i++){
420 macro = (TMacroData*)fMacroList->GetValue(selIterator->At(i)->GetTitle());
423 Error("Apply process macros",
424 "Macro list is corrupted: Macro \"%s\" is not registered!",
425 selIterator->At(i)->GetTitle());
429 #ifdef ALIEVETRDTRACKLIST_DEBUG
430 printf("AliEveTRDTrackList: Checking selection macro: %s\n", macro->GetName());
433 // Find the type of the process macro
434 macroType = macro->GetType();
436 // Single track select macro
437 if (macroType == kSingleTrackSelect) {
438 // Has already been processed by ApplySTSelectionMacros(...)
439 if(mSelType) mSelType[i] = macroType;
441 // Correlated tracks select macro
442 else if (macroType == kCorrelTrackSelect) {
443 if(mSelType) mSelType[i] = macroType;
445 // Create the command
446 if(selCmds) selCmds[i] = macro->GetCmd();
448 Error("Apply process macros",
449 "Macro list corrupted: Macro \"%s/%s.C\" is not registered as a selection macro!",
450 macro->GetPath(), macro->GetName());
451 if(mSelType) mSelType[i] = kUnknown;
455 // Allocate memory for the histograms
456 if (numHistoMacros > 0){
457 histos = new TH1*[numHistoMacros];
458 memset(histos, 0, numHistoMacros*sizeof(TH1*));
461 //////////////////////////////////
462 // WALK THROUGH THE LIST OF TRACKS
463 //////////////////////////////////
464 for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter){
465 if(!(track1 = dynamic_cast<AliEveTRDTrack*>(*iter))) continue;
467 // Skip tracks that have not been selected
468 if (!track1->GetRnrState()) continue;
470 // Cast to AliTRDtrackV1
471 gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)%p;", (void*)track1));
472 gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_1 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
474 // Collect data for each macro
475 for (Int_t i = 0, histoIndex = 0; i < procIterator->GetEntries(); i++){
476 // Single track histo
477 if (mProcType[i] == kSingleTrackHisto){
478 if(histos) histos[histoIndex++] = (TH1*)gROOT->ProcessLineSync(procCmds[i]);
479 // Correlated tracks histo
480 } else if (mProcType[i] == kCorrelTrackHisto) {
481 // Loop over all pairs behind the current one - together with the other loop this will be a loop
482 // over all pairs. We have a pair of tracks, if and only if both tracks of the pair are selected (Rnr-state)
483 // and are not equal.
484 // The correlated tracks process macro will be applied to all pairs that will be additionally selected by
485 // all correlated tracks selection macros.
486 TEveElement::List_i iter2 = iter;
488 for ( ; iter2 != this->EndChildren(); ++iter2)
490 if(!(track2 = dynamic_cast<AliEveTRDTrack*>(*iter2))) continue;
492 // Skip tracks that have not been selected
493 if (!track2->GetRnrState()) continue;
495 // Cast to AliTRDtrackV1
496 gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)%p;", (void*)track2));
497 gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_2 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
499 // Select track by default (so it will be processed, if there are no correlated tracks selection macros!)
500 selectedByCorrSelMacro = kTRUE;
501 for (Int_t j = 0; j < selIterator->GetEntries(); j++){
502 if (mSelType && mSelType[j] == kCorrelTrackSelect){
503 selectedByCorrSelMacro = (Bool_t)gROOT->ProcessLineSync(selCmds[j]);
504 if (!selectedByCorrSelMacro) break;
508 // If the pair has not been selected by the correlated tracks selection macros, skip it!
509 if (!selectedByCorrSelMacro) continue;
511 if(histos) histos[histoIndex] = (TH1*)gROOT->ProcessLineSync(procCmds[i]);
515 // Single track analyse
516 else if (mProcType[i] == kSingleTrackAnalyse) {
517 // Create data pointers in CINT, execute the macro and get the data
518 gROOT->ProcessLineSync("Double_t* results = 0;");
519 gROOT->ProcessLineSync("Int_t n = 0;");
520 gROOT->ProcessLineSync(procCmds[i]);
521 Double_t* results = (Double_t*)gROOT->ProcessLineSync("results;");
522 Int_t nResults = (Int_t)gROOT->ProcessLineSync("n;");
525 Error("Apply macros", "Error reading data from macro \"%s\"", procIterator->At(i)->GetTitle());
528 for (Int_t resInd = 0; resInd < nResults; resInd++){
529 (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << results[resInd] << (Char_t*)"\n";
535 // Correlated tracks analyse
536 else if (mProcType[i] == kCorrelTrackAnalyse){
537 // Loop over all pairs behind the current one - together with the other loop this will be a loop
538 // over all pairs. We have a pair of tracks, if and only if both tracks of the pair are selected (Rnr-state)
539 // and are not equal.
540 // The correlated tracks process macro will be applied to all pairs that will be additionally selected by
541 // all correlated tracks selection macros.
542 TEveElement::List_i iter2 = iter;
545 for ( ; iter2 != this->EndChildren(); ++iter2) {
546 if(!(track2 = dynamic_cast<AliEveTRDTrack*>(*iter2))) continue;
548 // Skip tracks that have not been selected
549 if (!track2->GetRnrState()) continue;
551 // Cast to AliTRDtrackV1
552 gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)%p;", (void*)track2));
553 gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_2 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
555 // Select track by default (so it will be processed, if there are no correlated tracks selection macros!)
556 selectedByCorrSelMacro = kTRUE;
557 for (Int_t j = 0; j < selIterator->GetEntries(); j++) {
558 if (mSelType && mSelType[j] == kCorrelTrackSelect) {
559 selectedByCorrSelMacro = (Bool_t)gROOT->ProcessLineSync(selCmds[j]);
560 if (!selectedByCorrSelMacro) break;
564 // If the pair has not been selected by the correlated tracks selection macros, skip it!
565 if (!selectedByCorrSelMacro) continue;
567 // Create data pointers in CINT, execute the macro and get the data
568 gROOT->ProcessLineSync("Double_t* results = 0;");
569 gROOT->ProcessLineSync("Int_t n = 0;");
570 gROOT->ProcessLineSync(procCmds[i]);
571 Double_t* results = (Double_t*)gROOT->ProcessLineSync("results;");
572 Int_t nResults = (Int_t)gROOT->ProcessLineSync("n;");
575 Error("Apply macros", "Error reading data from macro \"%s\"", procIterator->At(i)->GetTitle());
578 for (Int_t resInd = 0; resInd < nResults; resInd++) {
579 (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << results[resInd] << (Char_t*)"\n";
589 for (Int_t i = 0, histoIndex = 0; i < procIterator->GetEntries() && histoIndex < numHistoMacros; i++) {
590 if (mProcType[i] == kSingleTrackHisto || mProcType[i] == kCorrelTrackHisto) {
591 // Might be empty (e.g. no tracks have been selected)!
592 if (histos[histoIndex]) {
593 (*fDataTree) << Form("TrackData%d", i) << Form("Macro%d=", i) << histos[histoIndex] << (Char_t*)"\n";
599 if (fDataTree) delete fDataTree;
602 if (procCmds) delete [] procCmds;
604 if (mProcType) delete [] mProcType;
607 if (selCmds) delete [] selCmds;
609 if (mSelType) delete [] mSelType;
612 if (histos) delete [] histos;
616 // A.B. gROOT->Reset();
618 // If there is data, select the first data set
619 if (procIterator->GetEntries() > 0) SETBIT(fHistoDataSelected, 0);
621 // Now the data is stored in "/tmp/TRD.TrackListMacroData_$USER.root"
622 // The editor will access this file to display the data
626 //______________________________________________________
627 void AliEveTRDTrackList::ApplySTSelectionMacros(const TList* iterator)
629 // Uses the iterator (for the selected selection macros) to apply the selected macros to the data.
630 // The rnr-states of the tracks are set according to the result of the macro calls (kTRUE, if all
631 // macros return kTRUE for this track, otherwise: kFALSE).
632 // "ST" stands for "single track". This means that only single track selection macros are applied.
633 // Correlated tracks selection macros will be used inside the call of ApplyProcessMacros(...)!
635 TMacroData* macro = 0;
636 AliEveTRDTrackListMacroType macroType = kUnknown;
637 AliEveTRDTrack* track1 = 0;
638 Bool_t selectedByMacro = kFALSE;
641 // A.B. gROOT->Reset();
643 // Select all tracks at first. A track is then deselected, if at least one selection macro
644 // returns kFALSE for this track.
645 // Enable all tracks (Note: EnableListElements(..) will call "ElementChanged", which will cause unforeseen behaviour!)
646 for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter) ((TEveElement*)(*iter))->SetRnrState(kTRUE);
649 for (Int_t i = 0; i < iterator->GetEntries(); i++){
650 macro = (TMacroData*)fMacroList->GetValue(iterator->At(i)->GetTitle());
653 Error("Apply selection macros",
654 "Macro list is corrupted: Macro \"%s\" is not registered!", iterator->At(i)->GetTitle());
658 #ifdef ALIEVETRDTRACKLIST_DEBUG
659 printf("AliEveTRDTrackList: Applying selection macro: %s\n", macro->GetName());
662 // Determine macro type
663 macroType = macro->GetType();
665 // Single track select macro
666 if (macroType == kSingleTrackSelect){
667 // Walk through the list of tracks
668 for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter)
670 track1 = dynamic_cast<AliEveTRDTrack*>(*iter);
672 if (!track1) continue;
674 // If the track has already been deselected, nothing is to do here
675 if (!track1->GetRnrState()) continue;
677 // Cast to AliTRDtrackV1
678 gROOT->ProcessLineSync(Form("AliEveTRDTrack *automaticTrack = (AliEveTRDTrack*)%p;", (void*)track1));
679 gROOT->ProcessLineSync("AliTRDtrackV1* automaticTrackV1_1 = (AliTRDtrackV1*)automaticTrack->GetUserData();");
680 selectedByMacro = (Bool_t)gROOT->ProcessLineSync(macro->GetCmd());
681 track1->SetRnrState(selectedByMacro && track1->GetRnrState());
684 // Correlated tracks select macro
685 else if (macroType == kCorrelTrackSelect){
686 // Will be processed in ApplyProcessMacros(...)
689 Error("Apply selection macros",
690 "Macro list corrupted: Macro \"%s/%s.C\" is not registered as a selection macro!",
691 macro->GetPath(), macro->GetName());
696 // A.B. gROOT->Reset();
699 //______________________________________________________
700 AliEveTRDTrackList::AliEveTRDTrackListMacroType AliEveTRDTrackList::GetMacroType(const Char_t* name, Bool_t UseList) const
702 // Returns the type of the corresponding macro.
703 // If "UseList" is kTRUE, the type will be looked up in the internal list (very fast). But if this list
704 // does not exist, you have to use kFALSE for this parameter. Then the type will be determined by the
705 // prototype! NOTE: It is assumed that the macro has been compiled! If not, the return value is not
706 // predictable, but normally will be kUnknown.
707 // Note: AddMacro(Fast) will update the internal list and RemoveMacros respectively.
709 AliEveTRDTrackListMacroType type = kUnknown;
711 // Re-do the check of the macro type
713 // Single track select macro or single track histo macro?
714 TFunction* f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*", kTRUE);
717 // Some additional check (is the parameter EXACTLY of the desired type?)
718 if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0)
720 // Single track select macro?
721 if (!strcmp(f->GetReturnTypeName(), "Bool_t"))
723 type = kSingleTrackSelect;
725 // single track histo macro?
726 else if (!strcmp(f->GetReturnTypeName(), "TH1*"))
728 type = kSingleTrackHisto;
732 // Single track analyse macro?
733 else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*, Double_t*&, Int_t&", kTRUE))
736 if (!strcmp(f->GetReturnTypeName(), "void"))
738 // Some additional check (are the parameters EXACTLY of the desired type?)
739 if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
740 strstr(f->GetMangledName(), "cODouble_tmUaNsP") != 0x0 &&
741 strstr(f->GetMangledName(), "cOInt_taNsP") != 0x0)
743 type = kSingleTrackAnalyse;
747 // Correlated tracks select macro or correlated tracks histo macro?
748 else if ((f = gROOT->GetGlobalFunctionWithPrototype(name, "const AliTRDtrackV1*, const AliTRDtrackV1*", kTRUE))
751 // Some additional check (is the parameter EXACTLY of the desired type?)
752 if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
753 strstr(f->GetMangledName(), "cOconstsPAliTRDtrackV1mUsP") != 0x0)
755 // Correlated track select macro?
756 if (!strcmp(f->GetReturnTypeName(), "Bool_t"))
758 type = kCorrelTrackSelect;
760 // Correlated track histo macro?
761 else if (!strcmp(f->GetReturnTypeName(), "TH1*"))
763 type = kCorrelTrackHisto;
767 // Correlated tracks analyse macro?
768 else if ((f = gROOT->GetGlobalFunctionWithPrototype(name,
769 "const AliTRDtrackV1*, const AliTRDtrackV1*, Double_t*&, Int_t&", kTRUE))
772 if (!strcmp(f->GetReturnTypeName(), "void"))
774 // Some additional check (is the parameter EXACTLY of the desired type?)
775 if (strstr(f->GetMangledName(), "oPconstsPAliTRDtrackV1mUsP") != 0x0 &&
776 strstr(f->GetMangledName(), "cOconstsPAliTRDtrackV1mUsP") != 0x0 &&
777 strstr(f->GetMangledName(), "cODouble_tmUaNsP") != 0x0 &&
778 strstr(f->GetMangledName(), "cOInt_taNsP") != 0x0)
780 type = kCorrelTrackAnalyse;
785 // Use list to look up the macro type
788 TMacroData* macro = 0;
789 macro = (TMacroData*)fMacroList->GetValue(name);
790 if (macro == 0) return kUnknown;
792 type = macro->GetType();
795 case kSingleTrackSelect:
796 case kSingleTrackAnalyse:
797 case kSingleTrackHisto:
798 case kCorrelTrackSelect:
799 case kCorrelTrackAnalyse:
800 case kCorrelTrackHisto:
811 //______________________________________________________
812 void AliEveTRDTrackList::RemoveSelectedMacros(const TList* iterator)
814 // Uses the iterator (for the selected macros) to remove the selected macros from
815 // the corresponding list.
819 for (Int_t i = 0; i < iterator->GetEntries(); i++)
821 entry = (TPair*)fMacroList->FindObject(iterator->At(i)->GetTitle());
825 Error("AliEveTRDTrackList::RemoveSelectedMacros", "Macro \"%s\" not found in list!",
826 iterator->At(i)->GetTitle());
833 Error("AliEveTRDTrackList::RemoveSelectedMacros", "Key for macro \"%s\" not found in list!",
834 iterator->At(i)->GetTitle());
838 // Key and value will be deleted, too, since fMacroList is the owner of them
839 Bool_t rem = fMacroList->DeleteEntry(key);
843 #ifdef ALIEVETRDTRACKLIST_DEBUG
844 printf("AliEveTRDTrackList::RemoveSelectedMacros(): Removed macro: %s\n", iterator->At(i)->GetTitle());
849 Error("AliEveTRDTrackList::RemoveSelectedMacros", "Macro \"%s\" could not be removed from the list!",
850 iterator->At(i)->GetTitle());
855 //______________________________________________________
856 void AliEveTRDTrackList::UpdateTrackStyle(AliEveTRDTrack::AliEveTRDTrackState s, UChar_t ss)
858 // Updates the track style and sets this style for each track.
862 case AliEveTRDTrack::kSource:
863 SETBIT(fSelectedStyle, AliEveTRDTrack::kSource);
865 case AliEveTRDTrack::kPID:
866 CLRBIT(fSelectedStyle, AliEveTRDTrack::kSource);
869 case AliTRDpidUtil::kLQ:
870 CLRBIT(fSelectedStyle, AliEveTRDTrack::kPID);
872 case AliTRDpidUtil::kNN:
873 SETBIT(fSelectedStyle, AliEveTRDTrack::kPID);
877 case AliEveTRDTrack::kTrackCosmics:
878 SETBIT(fSelectedStyle, AliEveTRDTrack::kTrackCosmics);
880 case AliEveTRDTrack::kTrackModel:
881 CLRBIT(fSelectedStyle, AliEveTRDTrack::kTrackCosmics);
884 case AliEveTRDTrack::kRieman:
885 CLRBIT(fSelectedStyle, AliEveTRDTrack::kTrackModel);
887 case AliEveTRDTrack::kKalman:
888 //AliWarning("Kalman fit under testing for the moment.");
889 SETBIT(fSelectedStyle, AliEveTRDTrack::kTrackModel);
896 // Walk through the list of tracks
897 AliEveTRDTrack* track = 0x0;
898 for (TEveElement::List_i iter = this->BeginChildren(); iter != this->EndChildren(); ++iter)
900 if (!(track = dynamic_cast<AliEveTRDTrack*>(*iter))) continue;
902 track->SetStatus(fSelectedStyle);