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>
33 #include <TGridResult.h>
37 #include "AliRunTag.h"
38 #include "AliEventTag.h"
39 #include "AliTagAnalysis.h"
40 #include "AliEventTagCuts.h"
41 #include "AliDetectorTagCuts.h"
42 #include "AliLHCTagCuts.h"
43 #include "AliRunTagCuts.h"
44 #include "AliXMLCollection.h"
48 ClassImp(AliTagAnalysis)
50 //___________________________________________________________________________
51 AliTagAnalysis::AliTagAnalysis():
58 //Default constructor for a AliTagAnalysis
61 //___________________________________________________________________________
62 AliTagAnalysis::AliTagAnalysis(const char* type):
69 //constructor for a AliTagAnalysis
72 //___________________________________________________________________________
73 AliTagAnalysis::~AliTagAnalysis() {
74 //Default destructor for a AliTagAnalysis
75 if(ftagresult) delete ftagresult;
76 if(fChain) delete fChain;
77 if(fGlobalList) delete fGlobalList;
80 //___________________________________________________________________________
81 Bool_t AliTagAnalysis::AddTagsFile(const char *alienUrl) {
82 // Add a single tags file to the chain
86 if (! fChain) fChain = new TChain("T");
88 TFile *f = TFile::Open(alienUrl,"READ");
89 fChain->Add(alienUrl);
90 AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
93 if (fChain->GetEntries() == 0 )
99 //___________________________________________________________________________
100 void AliTagAnalysis::ChainLocalTags(const char *dirname) {
101 //Searches the entries of the provided direcory
102 //Chains the tags that are stored locally
103 fTagDirName = dirname;
104 TString fTagFilename;
106 if (! fChain) fChain = new TChain("T");
107 const char * tagPattern = 0x0;
108 if(fAnalysisType == "ESD") tagPattern = "ESD.tag.root";
109 else if(fAnalysisType == "AOD") tagPattern = "AOD.tag.root";
110 else AliFatal("Only ESD and AOD type is implemented!!!");
112 // Open the working directory
113 void * dirp = gSystem->OpenDirectory(fTagDirName);
114 const char * name = 0x0;
115 // Add all files matching *pattern* to the chain
116 while((name = gSystem->GetDirEntry(dirp))) {
117 if (strstr(name,tagPattern)) {
118 fTagFilename = fTagDirName;
120 fTagFilename += name;
122 fChain->Add(fTagFilename);
123 printf("Tag file %s\n", fTagFilename.Data());
127 AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
133 //___________________________________________________________________________
134 TChain * AliTagAnalysis::ChainGridTags(TGridResult *res) {
135 //Loops overs the entries of the TGridResult
136 //Chains the tags that are stored in the GRID
138 Int_t nEntries = ftagresult->GetEntries();
140 if (! fChain) fChain = new TChain("T");
142 TString gridname = "alien://";
145 for(Int_t i = 0; i < nEntries; i++) {
146 alienUrl = ftagresult->GetKey(i,"turl");
147 fChain->Add(alienUrl);
153 //___________________________________________________________________________
154 TChain *AliTagAnalysis::QueryTags(AliRunTagCuts *runTagCuts,
155 AliLHCTagCuts *lhcTagCuts,
156 AliDetectorTagCuts *detTagCuts,
157 AliEventTagCuts *evTagCuts) {
158 //Queries the tag chain using the defined
159 //event tag cuts from the AliEventTagCuts object
160 //and returns a TChain along with the associated TEventList
161 AliInfo(Form("Querying the tags........"));
165 if(fAnalysisType == "ESD") aliceFile = "esdTree";
166 else if(fAnalysisType == "AOD") {
167 aliceFile = "aodTree";
170 else AliFatal("Only ESD and AOD type is implemented!!!");
173 TChain *esdChain = new TChain(aliceFile.Data());
175 fGlobalList = new TEntryList();
177 //Defining tag objects
178 AliRunTag *tag = new AliRunTag;
179 AliEventTag *evTag = new AliEventTag;
180 fChain->SetBranchAddress("AliTAG",&tag);
187 TEntryList* localList = new TEntryList();
194 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
195 fChain->GetEntry(iTagFiles);
196 TTree* tree = fChain->GetTree();
198 // Fix for aod tags: for each tree in the chain, merge the entries
200 cEntries = tree->GetEntries();
205 if(runTagCuts->IsAccepted(tag)) {
206 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
207 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
208 if ((iev == 0) || !aod) localList->Reset();
209 Int_t iEvents = tag->GetNEvents();
210 const TClonesArray *tagList = tag->GetEventTags();
211 for(Int_t i = 0; i < iEvents; i++) {
212 evTag = (AliEventTag *) tagList->At(i);
213 guid = evTag->GetGUID();
214 turl = evTag->GetTURL();
215 path = evTag->GetPath();
216 localList->SetTreeName(aliceFile.Data());
217 if(turl!="") localList->SetFileName(turl.Data());
218 else localList->SetFileName(path.Data());
220 if(evTagCuts->IsAccepted(evTag)) {
221 if(aod) localList->Enter(iev);
222 else localList->Enter(i);
226 if ((ientry == cEntries-1) || !aod) {
227 iAccepted += localList->GetN();
228 if(path != "") esdChain->AddFile(path);
229 else if(turl != "") esdChain->AddFile(turl);
230 fGlobalList->Add(localList);
238 AliInfo(Form("Accepted events: %d",iAccepted));
240 esdChain->SetEntryList(fGlobalList,"ne");
246 //___________________________________________________________________________
247 TChain *AliTagAnalysis::QueryTags(const char *fRunCut,
249 const char *fDetectorCut,
250 const char *fEventCut) {
251 //Queries the tag chain using the defined
252 //event tag cuts from the AliEventTagCuts object
253 //and returns a TChain along with the associated TEventList
254 AliInfo(Form("Querying the tags........"));
258 if(fAnalysisType == "ESD") aliceFile = "esdTree";
259 else if(fAnalysisType == "AOD") {
260 aliceFile = "aodTree";
263 else AliFatal("Only ESD and AOD type is implemented!!!");
266 TChain *esdChain = new TChain(aliceFile.Data());
268 fGlobalList = new TEntryList();
270 //Defining tag objects
271 AliRunTag *tag = new AliRunTag;
272 AliEventTag *evTag = new AliEventTag;
273 fChain->SetBranchAddress("AliTAG",&tag);
279 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
280 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
281 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
282 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
284 TEntryList* localList = new TEntryList();
292 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
293 fChain->GetEntry(iTagFiles);
294 if (current != fChain->GetTreeNumber()) {
295 fRunFormula->UpdateFormulaLeaves();
296 fLHCFormula->UpdateFormulaLeaves();
297 fDetectorFormula->UpdateFormulaLeaves();
298 fEventFormula->UpdateFormulaLeaves();
299 // Fix for aod tags: for each tree in the chain, merge the entries
300 cEntries = (fChain->GetTree())->GetEntries();
304 current = fChain->GetTreeNumber();
306 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
307 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
308 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
309 if ((iev == 0) || !aod) localList->Reset();
310 Int_t iEvents = fEventFormula->GetNdata();
311 const TClonesArray *tagList = tag->GetEventTags();
312 for(Int_t i = 0; i < iEvents; i++) {
313 evTag = (AliEventTag *) tagList->At(i);
314 guid = evTag->GetGUID();
315 turl = evTag->GetTURL();
316 path = evTag->GetPath();
317 localList->SetTreeName(aliceFile.Data());
318 localList->SetFileName(turl.Data());
319 if(fEventFormula->EvalInstance(i) == 1) {
320 if(aod) localList->Enter(iev);
321 else localList->Enter(i);
326 if ((ientry == cEntries-1) || !aod) {
327 if(path != "") esdChain->AddFile(path);
328 else if(turl != "") esdChain->AddFile(turl);
329 fGlobalList->Add(localList);
330 iAccepted += localList->GetN();
338 AliInfo(Form("Accepted events: %d",iAccepted));
339 esdChain->SetEntryList(fGlobalList,"ne");
345 //___________________________________________________________________________
346 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name,
347 AliRunTagCuts *runTagCuts,
348 AliLHCTagCuts *lhcTagCuts,
349 AliDetectorTagCuts *detTagCuts,
350 AliEventTagCuts *evTagCuts) {
351 //Queries the tag chain using the defined
352 //event tag cuts from the AliEventTagCuts object
353 //and returns a XML collection
354 AliInfo(Form("Creating the collection........"));
357 if(fAnalysisType == "AOD") aod = kTRUE;
360 AliXMLCollection *collection = new AliXMLCollection();
361 collection->SetCollectionName(name);
362 collection->WriteHeader();
369 TEntryList* localList = new TEntryList();
375 Int_t iRejectedRun = 0;
376 Int_t iRejectedLHC = 0;
377 Int_t iRejectedDet = 0;
378 Int_t iRejectedEvt = 0;
380 Int_t iTotalEvents = 0;
382 Int_t iAcceptedEvtInFile = 0;
383 Int_t iRejectedEvtInFile = 0;
385 //Defining tag objects
386 AliRunTag *tag = new AliRunTag;
387 AliEventTag *evTag = new AliEventTag;
388 fChain->SetBranchAddress("AliTAG",&tag);
390 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
392 fChain->GetEntry(iTagFiles);
393 TTree* tree = fChain->GetTree();
395 // Fix for aod tags: for each tree in the chain, merge the entries
397 cEntries = tree->GetEntries();
402 iTotalEvents += tag->GetNEvents();
403 if ((iev == 0) || !aod) localList->Reset();
404 if(runTagCuts->IsAccepted(tag)) {
405 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
406 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
407 Int_t iEvents = tag->GetNEvents();
408 const TClonesArray *tagList = tag->GetEventTags();
409 iRejectedEvtInFile = 0;
410 iAcceptedEvtInFile = 0;
411 for(Int_t i = 0; i < iEvents; i++) {
412 evTag = (AliEventTag *) tagList->At(i);
413 guid = evTag->GetGUID();
414 turl = evTag->GetTURL();
415 lfn = turl(8,turl.Length());
416 if(evTagCuts->IsAccepted(evTag)) {
417 if(aod) localList->Enter(iev);
418 else localList->Enter(i);
419 iAcceptedEvtInFile++;
423 iRejectedEvtInFile++;
427 if ((ientry == cEntries-1) || !aod) {
428 iAccepted += localList->GetN();
429 collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList,iAcceptedEvtInFile,iRejectedEvtInFile);
433 iRejectedDet += tag->GetNEvents();
437 iRejectedLHC += tag->GetNEvents();
441 iRejectedRun += tag->GetNEvents();
446 collection->WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
447 collection->Export();
453 //___________________________________________________________________________
454 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name,
457 const char *fDetectorCut,
458 const char *fEventCut) {
459 //Queries the tag chain using the defined
460 //event tag cuts from the AliEventTagCuts object
461 //and returns a XML collection
462 AliInfo(Form("Creating the collection........"));
465 if(fAnalysisType == "AOD") aod = kTRUE;
467 AliXMLCollection *collection = new AliXMLCollection();
468 collection->SetCollectionName(name);
469 collection->WriteHeader();
474 TEntryList* localList = new TEntryList();
481 Int_t iRejectedRun = 0;
482 Int_t iRejectedLHC = 0;
483 Int_t iRejectedDet = 0;
484 Int_t iRejectedEvt = 0;
486 Int_t iTotalEvents = 0;
488 Int_t iAcceptedEvtInFile = 0;
489 Int_t iRejectedEvtInFile = 0;
491 //Defining tag objects
492 AliRunTag *tag = new AliRunTag;
493 AliEventTag *evTag = new AliEventTag;
494 fChain->SetBranchAddress("AliTAG",&tag);
496 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
497 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
498 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
499 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
502 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
504 fChain->GetEntry(iTagFiles);
505 if (current != fChain->GetTreeNumber()) {
506 fRunFormula->UpdateFormulaLeaves();
507 fLHCFormula->UpdateFormulaLeaves();
508 fDetectorFormula->UpdateFormulaLeaves();
509 fEventFormula->UpdateFormulaLeaves();
510 // Fix for aod tags: for each tree in the chain, merge the entries
511 cEntries = (fChain->GetTree())->GetEntries();
515 current = fChain->GetTreeNumber();
518 iTotalEvents += tag->GetNEvents();
519 if ((iev == 0) || !aod) localList->Reset();
520 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
521 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
522 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
523 Int_t iEvents = fEventFormula->GetNdata();
524 const TClonesArray *tagList = tag->GetEventTags();
525 iRejectedEvtInFile = 0;
526 iAcceptedEvtInFile = 0;
527 for(Int_t i = 0; i < iEvents; i++) {
528 evTag = (AliEventTag *) tagList->At(i);
529 guid = evTag->GetGUID();
530 turl = evTag->GetTURL();
531 lfn = turl(8,turl.Length());
532 if(fEventFormula->EvalInstance(i) == 1) {
533 if(aod) localList->Enter(iev);
534 else localList->Enter(i);
535 iAcceptedEvtInFile++;
539 iRejectedEvtInFile++;
543 if ((ientry == cEntries-1) || !aod) {
544 collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList,iAcceptedEvtInFile, iRejectedEvtInFile);
545 iAccepted += localList->GetN();
549 iRejectedDet += tag->GetNEvents();
553 iRejectedLHC += tag->GetNEvents();
557 iRejectedRun += tag->GetNEvents();
561 collection->WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
562 collection->Export();
568 //___________________________________________________________________________
569 TChain *AliTagAnalysis::GetInputChain(const char* system, const char *wn) {
570 //returns the chain+event list - used in batch sessions
571 // this function will be removed once the new root
572 // improvements are committed
573 TString fsystem = system;
576 TChain *fAnalysisChain = 0;
577 if(fAnalysisType == "ESD") fAnalysisChain = new TChain("esdTree");
578 else if(fAnalysisType == "AOD") fAnalysisChain = new TChain("aodTree");
579 else AliFatal("Only ESD and AOD type is implemented!!!");
582 TEventList *fEventList = new TEventList();
583 AliXMLCollection *collection = AliXMLCollection::Open(wn);
586 while (collection->Next()) {
587 AliInfo(Form("Adding: %s",collection->GetTURL("")));
588 fAnalysisChain->Add(collection->GetTURL(""));
589 TEntryList *list = (TEntryList *)collection->GetEventList("");
590 for(Int_t i = 0; i < list->GetN(); i++) fEventList->Enter(iAccepted+list->GetEntry(i));
592 if(fsystem == "pp") iAccepted += 100;
593 else if(fsystem == "PbPb") iAccepted += 1;
596 fAnalysisChain->SetEventList(fEventList);
598 AliInfo(Form("Number of selected events: %d",fEventList->GetN()));
600 return fAnalysisChain;
603 //___________________________________________________________________________
604 TChain *AliTagAnalysis::GetChainFromCollection(const char* collectionname,
605 const char* treename) {
606 //returns the TChain+TEntryList object- used in batch sessions
607 TString aliceFile = treename;
609 TChain *fAnalysisChain = 0;
610 if(aliceFile == "esdTree") fAnalysisChain = new TChain("esdTree");
611 else if(aliceFile == "aodTree") fAnalysisChain = new TChain("aodTree");
612 else AliFatal("Inconsistent tree name - use esdTree or aodTree!");
615 fGlobalList = new TEntryList();
616 AliXMLCollection *collection = AliXMLCollection::Open(collectionname);
618 // Tag selection summary per file
619 TMap *tagCutSummary = new TMap();
620 tagCutSummary->SetName("TagCutSumm");
623 while (collection->Next()) {
624 AliInfo(Form("Adding: %s",collection->GetTURL("")));
625 fAnalysisChain->Add(collection->GetTURL(""));
626 TEntryList *list = (TEntryList *)collection->GetEventList("");
627 list->SetTreeName(aliceFile.Data());
628 list->SetFileName(collection->GetTURL(""));
629 fGlobalList->Add(list);
630 iAccepted += list->GetN();
631 if (collection->GetCutSumm())
632 tagCutSummary->Add(new TObjString(collection->GetTURL("")), new TObjString(collection->GetCutSumm()));
635 fAnalysisChain->SetEntryList(fGlobalList,"ne");
637 AliInfo(Form("Number of selected events: %d",iAccepted));
639 TList *aUserInfo = fAnalysisChain->GetUserInfo();
640 aUserInfo->Add(tagCutSummary);
649 collection->GetCollectionSummary(&iTotalEvents, &iAccEv, &iRejRun, &iRejLHC, &iRejDet, &iRejEvt);
653 sprintf(nstr, "TotalEvents=%i", iTotalEvents);
654 TObjString *iTotStr = new TObjString(nstr);
655 aUserInfo->Add(iTotStr);
657 sprintf(nstr, "AcceptedEvents=%i", iAccepted);
658 TObjString *iAccStr = new TObjString(nstr);
659 aUserInfo->Add(iAccStr);
661 sprintf(nstr, "RejectedRun=%i", iRejRun);
662 TObjString *iRejRunStr = new TObjString(nstr);
663 aUserInfo->Add(iRejRunStr);
665 sprintf(nstr, "RejectedLHC=%i", iRejLHC);
666 TObjString *iRejLHCStr = new TObjString(nstr);
667 aUserInfo->Add(iRejLHCStr);
669 sprintf(nstr, "RejectedDet=%i", iRejDet);
670 TObjString *iRejDetStr = new TObjString(nstr);
671 aUserInfo->Add(iRejDetStr);
673 sprintf(nstr, "RejectedEvt=%i", iRejEvt);
674 TObjString *iRejEvtStr = new TObjString(nstr);
675 aUserInfo->Add(iRejEvtStr);
677 return fAnalysisChain;