1 /**************************************************************************
2 * Author: Panos Christakoglou. *
3 * Contributors are mentioned in the code where appropriate. *
5 * Permission to use, copy, modify and distribute this software and its *
6 * documentation strictly for non-commercial purposes is hereby granted *
7 * without fee, provided that the above copyright notice appears in all *
8 * copies and that both the copyright notice and this permission notice *
9 * appear in the supporting documentation. The authors make no claims *
10 * about the suitability of this software for any purpose. It is *
11 * provided "as is" without express or implied warranty. *
12 **************************************************************************/
16 //-----------------------------------------------------------------
17 // AliTagAnalysis class
18 // This is the class to deal with the tag analysis
19 // Origin: Panos Christakoglou, UOA-CERN, Panos.Christakoglou@cern.ch
20 //-----------------------------------------------------------------
23 #include <Riostream.h>
27 #include <TEventList.h>
28 #include <TEntryList.h>
29 #include <TTreeFormula.h>
32 #include <TGridResult.h>
36 #include "AliRunTag.h"
37 #include "AliEventTag.h"
38 #include "AliTagAnalysis.h"
39 #include "AliEventTagCuts.h"
40 #include "AliDetectorTagCuts.h"
41 #include "AliLHCTagCuts.h"
42 #include "AliRunTagCuts.h"
43 #include "AliXMLCollection.h"
47 ClassImp(AliTagAnalysis)
49 //___________________________________________________________________________
50 AliTagAnalysis::AliTagAnalysis():
57 //Default constructor for a AliTagAnalysis
60 //___________________________________________________________________________
61 AliTagAnalysis::AliTagAnalysis(const char* type):
68 //constructor for a AliTagAnalysis
71 //___________________________________________________________________________
72 AliTagAnalysis::~AliTagAnalysis() {
73 //Default destructor for a AliTagAnalysis
74 if(ftagresult) delete ftagresult;
75 if(fChain) delete fChain;
76 if(fGlobalList) delete fGlobalList;
79 //___________________________________________________________________________
80 Bool_t AliTagAnalysis::AddTagsFile(const char *alienUrl) {
81 // Add a single tags file to the chain
85 if (! fChain) fChain = new TChain("T");
87 TFile *f = TFile::Open(alienUrl,"READ");
88 fChain->Add(alienUrl);
89 AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
92 if (fChain->GetEntries() == 0 )
98 //___________________________________________________________________________
99 void AliTagAnalysis::ChainLocalTags(const char *dirname) {
100 //Searches the entries of the provided direcory
101 //Chains the tags that are stored locally
102 fTagDirName = dirname;
103 TString fTagFilename;
105 if (! fChain) fChain = new TChain("T");
106 const char * tagPattern = 0x0;
107 if(fAnalysisType == "ESD") tagPattern = "ESD.tag.root";
108 else if(fAnalysisType == "AOD") tagPattern = "AOD.tag.root";
109 else AliFatal("Only ESD and AOD type is implemented!!!");
111 // Open the working directory
112 void * dirp = gSystem->OpenDirectory(fTagDirName);
113 const char * name = 0x0;
114 // Add all files matching *pattern* to the chain
115 while((name = gSystem->GetDirEntry(dirp))) {
116 if (strstr(name,tagPattern)) {
117 fTagFilename = fTagDirName;
119 fTagFilename += name;
121 fChain->Add(fTagFilename);
122 printf("Tag file %s\n", fTagFilename.Data());
126 AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
132 //___________________________________________________________________________
133 TChain * AliTagAnalysis::ChainGridTags(TGridResult *res) {
134 //Loops overs the entries of the TGridResult
135 //Chains the tags that are stored in the GRID
137 Int_t nEntries = ftagresult->GetEntries();
139 if (! fChain) fChain = new TChain("T");
141 TString gridname = "alien://";
144 for(Int_t i = 0; i < nEntries; i++) {
145 alienUrl = ftagresult->GetKey(i,"turl");
146 fChain->Add(alienUrl);
152 //___________________________________________________________________________
153 TChain *AliTagAnalysis::QueryTags(AliRunTagCuts *runTagCuts,
154 AliLHCTagCuts *lhcTagCuts,
155 AliDetectorTagCuts *detTagCuts,
156 AliEventTagCuts *evTagCuts) {
157 //Queries the tag chain using the defined
158 //event tag cuts from the AliEventTagCuts object
159 //and returns a TChain along with the associated TEventList
160 AliInfo(Form("Querying the tags........"));
164 if(fAnalysisType == "ESD") aliceFile = "esdTree";
165 else if(fAnalysisType == "AOD") {
166 aliceFile = "aodTree";
169 else AliFatal("Only ESD and AOD type is implemented!!!");
172 TChain *esdChain = new TChain(aliceFile.Data());
174 fGlobalList = new TEntryList();
176 //Defining tag objects
177 AliRunTag *tag = new AliRunTag;
178 AliEventTag *evTag = new AliEventTag;
179 fChain->SetBranchAddress("AliTAG",&tag);
186 TEntryList* localList = 0;
193 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
194 fChain->GetEntry(iTagFiles);
195 TTree* tree = fChain->GetTree();
197 // Fix for aod tags: for each tree in the chain, merge the entries
199 cEntries = tree->GetEntries();
204 if(runTagCuts->IsAccepted(tag)) {
205 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
206 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
207 if ((iev == 0) || !aod) localList = new TEntryList();
208 Int_t iEvents = tag->GetNEvents();
209 const TClonesArray *tagList = tag->GetEventTags();
210 for(Int_t i = 0; i < iEvents; i++) {
211 evTag = (AliEventTag *) tagList->At(i);
212 guid = evTag->GetGUID();
213 turl = evTag->GetTURL();
214 path = evTag->GetPath();
215 localList->SetTreeName(aliceFile.Data());
216 if(turl!="") localList->SetFileName(turl.Data());
217 else localList->SetFileName(path.Data());
219 if(evTagCuts->IsAccepted(evTag)) {
220 if(aod) localList->Enter(iev);
221 else localList->Enter(i);
225 if ((ientry == cEntries-1) || !aod) {
226 iAccepted += localList->GetN();
227 if(path != "") esdChain->AddFile(path);
228 else if(turl != "") esdChain->AddFile(turl);
229 fGlobalList->Add(localList);
237 AliInfo(Form("Accepted events: %d",iAccepted));
238 esdChain->SetEntryList(fGlobalList,"ne");
243 //___________________________________________________________________________
244 TChain *AliTagAnalysis::QueryTags(const char *fRunCut,
246 const char *fDetectorCut,
247 const char *fEventCut) {
248 //Queries the tag chain using the defined
249 //event tag cuts from the AliEventTagCuts object
250 //and returns a TChain along with the associated TEventList
251 AliInfo(Form("Querying the tags........"));
255 if(fAnalysisType == "ESD") aliceFile = "esdTree";
256 else if(fAnalysisType == "AOD") {
257 aliceFile = "aodTree";
260 else AliFatal("Only ESD and AOD type is implemented!!!");
263 TChain *esdChain = new TChain(aliceFile.Data());
265 fGlobalList = new TEntryList();
267 //Defining tag objects
268 AliRunTag *tag = new AliRunTag;
269 AliEventTag *evTag = new AliEventTag;
270 fChain->SetBranchAddress("AliTAG",&tag);
276 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
277 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
278 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
279 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
281 TEntryList* localList = 0;
289 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
290 fChain->GetEntry(iTagFiles);
291 if (current != fChain->GetTreeNumber()) {
292 fRunFormula->UpdateFormulaLeaves();
293 fLHCFormula->UpdateFormulaLeaves();
294 fDetectorFormula->UpdateFormulaLeaves();
295 fEventFormula->UpdateFormulaLeaves();
296 // Fix for aod tags: for each tree in the chain, merge the entries
297 cEntries = (fChain->GetTree())->GetEntries();
301 current = fChain->GetTreeNumber();
303 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
304 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
305 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
306 if ((iev == 0) || !aod) localList = new TEntryList();
307 Int_t iEvents = fEventFormula->GetNdata();
308 const TClonesArray *tagList = tag->GetEventTags();
309 for(Int_t i = 0; i < iEvents; i++) {
310 evTag = (AliEventTag *) tagList->At(i);
311 guid = evTag->GetGUID();
312 turl = evTag->GetTURL();
313 path = evTag->GetPath();
314 localList->SetTreeName(aliceFile.Data());
315 localList->SetFileName(turl.Data());
316 if(fEventFormula->EvalInstance(i) == 1) {
317 if(aod) localList->Enter(iev);
318 else localList->Enter(i);
323 if ((ientry == cEntries-1) || !aod) {
324 if(path != "") esdChain->AddFile(path);
325 else if(turl != "") esdChain->AddFile(turl);
326 fGlobalList->Add(localList);
327 iAccepted += localList->GetN();
335 AliInfo(Form("Accepted events: %d",iAccepted));
336 esdChain->SetEntryList(fGlobalList,"ne");
341 //___________________________________________________________________________
342 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name,
343 AliRunTagCuts *runTagCuts,
344 AliLHCTagCuts *lhcTagCuts,
345 AliDetectorTagCuts *detTagCuts,
346 AliEventTagCuts *evTagCuts) {
347 //Queries the tag chain using the defined
348 //event tag cuts from the AliEventTagCuts object
349 //and returns a XML collection
350 AliInfo(Form("Creating the collection........"));
353 if(fAnalysisType == "AOD") aod = kTRUE;
356 AliXMLCollection *collection = new AliXMLCollection();
357 collection->SetCollectionName(name);
358 collection->WriteHeader();
365 TEntryList* localList = 0;
371 //Defining tag objects
372 AliRunTag *tag = new AliRunTag;
373 AliEventTag *evTag = new AliEventTag;
374 fChain->SetBranchAddress("AliTAG",&tag);
376 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
378 fChain->GetEntry(iTagFiles);
379 TTree* tree = fChain->GetTree();
381 // Fix for aod tags: for each tree in the chain, merge the entries
383 cEntries = tree->GetEntries();
388 if ((iev == 0) || !aod) localList = new TEntryList();
389 if(runTagCuts->IsAccepted(tag)) {
390 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
391 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
392 Int_t iEvents = tag->GetNEvents();
393 const TClonesArray *tagList = tag->GetEventTags();
394 for(Int_t i = 0; i < iEvents; i++) {
395 evTag = (AliEventTag *) tagList->At(i);
396 guid = evTag->GetGUID();
397 turl = evTag->GetTURL();
398 lfn = turl(8,turl.Length());
399 if(evTagCuts->IsAccepted(evTag)) {
400 if(aod) localList->Enter(iev);
401 else localList->Enter(i);
405 if ((ientry == cEntries-1) || !aod) {
406 collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList);
414 collection->Export();
419 //___________________________________________________________________________
420 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name,
423 const char *fDetectorCut,
424 const char *fEventCut) {
425 //Queries the tag chain using the defined
426 //event tag cuts from the AliEventTagCuts object
427 //and returns a XML collection
428 AliInfo(Form("Creating the collection........"));
431 if(fAnalysisType == "AOD") aod = kTRUE;
433 AliXMLCollection *collection = new AliXMLCollection();
434 collection->SetCollectionName(name);
435 collection->WriteHeader();
440 TEntryList* localList = 0;
448 //Defining tag objects
449 AliRunTag *tag = new AliRunTag;
450 AliEventTag *evTag = new AliEventTag;
451 fChain->SetBranchAddress("AliTAG",&tag);
453 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
454 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
455 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
456 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
459 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
461 fChain->GetEntry(iTagFiles);
462 if (current != fChain->GetTreeNumber()) {
463 fRunFormula->UpdateFormulaLeaves();
464 fLHCFormula->UpdateFormulaLeaves();
465 fDetectorFormula->UpdateFormulaLeaves();
466 fEventFormula->UpdateFormulaLeaves();
467 // Fix for aod tags: for each tree in the chain, merge the entries
468 cEntries = (fChain->GetTree())->GetEntries();
472 current = fChain->GetTreeNumber();
475 if ((iev == 0) || !aod) localList = new TEntryList();
476 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
477 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
478 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
479 Int_t iEvents = fEventFormula->GetNdata();
480 const TClonesArray *tagList = tag->GetEventTags();
481 for(Int_t i = 0; i < iEvents; i++) {
482 evTag = (AliEventTag *) tagList->At(i);
483 guid = evTag->GetGUID();
484 turl = evTag->GetTURL();
485 lfn = turl(8,turl.Length());
486 if(fEventFormula->EvalInstance(i) == 1) {
487 if(aod) localList->Enter(iev);
488 else localList->Enter(i);
492 if ((ientry == cEntries-1) || !aod) {
493 collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList);
500 collection->Export();
504 //___________________________________________________________________________
505 TChain *AliTagAnalysis::GetInputChain(const char* system, const char *wn) {
506 //returns the chain+event list - used in batch sessions
507 // this function will be removed once the new root
508 // improvements are committed
509 TString fsystem = system;
512 TChain *fAnalysisChain = 0;
513 if(fAnalysisType == "ESD") fAnalysisChain = new TChain("esdTree");
514 else if(fAnalysisType == "AOD") fAnalysisChain = new TChain("aodTree");
515 else AliFatal("Only ESD and AOD type is implemented!!!");
518 TEventList *fEventList = new TEventList();
519 AliXMLCollection *collection = AliXMLCollection::Open(wn);
522 while (collection->Next()) {
523 AliInfo(Form("Adding: %s",collection->GetTURL("")));
524 fAnalysisChain->Add(collection->GetTURL(""));
525 TEntryList *list = (TEntryList *)collection->GetEventList("");
526 for(Int_t i = 0; i < list->GetN(); i++) fEventList->Enter(iAccepted+list->GetEntry(i));
528 if(fsystem == "pp") iAccepted += 100;
529 else if(fsystem == "PbPb") iAccepted += 1;
532 fAnalysisChain->SetEventList(fEventList);
534 AliInfo(Form("Number of selected events: %d",fEventList->GetN()));
536 return fAnalysisChain;
539 //___________________________________________________________________________
540 TChain *AliTagAnalysis::GetChainFromCollection(const char* collectionname,
541 const char* treename) {
542 //returns the TChain+TEntryList object- used in batch sessions
543 TString aliceFile = treename;
545 TChain *fAnalysisChain = 0;
546 if(aliceFile == "esdTree") fAnalysisChain = new TChain("esdTree");
547 else if(aliceFile == "aodTree") fAnalysisChain = new TChain("aodTree");
548 else AliFatal("Inconsistent tree name - use esdTree or aodTree!");
551 fGlobalList = new TEntryList();
552 AliXMLCollection *collection = AliXMLCollection::Open(collectionname);
555 while (collection->Next()) {
556 AliInfo(Form("Adding: %s",collection->GetTURL("")));
557 fAnalysisChain->Add(collection->GetTURL(""));
558 TEntryList *list = (TEntryList *)collection->GetEventList("");
559 list->SetTreeName(aliceFile.Data());
560 list->SetFileName(collection->GetTURL(""));
561 fGlobalList->Add(list);
562 iAccepted += list->GetN();
565 fAnalysisChain->SetEntryList(fGlobalList,"ne");
567 AliInfo(Form("Number of selected events: %d",iAccepted));
569 return fAnalysisChain;