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 //___________________________________________________________________________
82 AliTagAnalysis::AddTagsFile(const char* alienUrl, Bool_t checkFile)
84 /// Add a single tags file to the chain
86 /// If checkFile=kTRUE (default) the file is opened to check
87 /// it can be and that it contains data.
88 /// It's safer but a lot longer...
90 if (!fChain) fChain = new TChain("T");
94 return ( fChain->AddFile(alienUrl,-1) > 0 );
98 return ( fChain->AddFile(alienUrl) > 0 );
103 //___________________________________________________________________________
104 void AliTagAnalysis::ChainLocalTags(const char *dirname) {
105 //Searches the entries of the provided direcory
106 //Chains the tags that are stored locally
107 fTagDirName = dirname;
108 TString fTagFilename;
110 if (! fChain) fChain = new TChain("T");
111 const char * tagPattern = 0x0;
112 if(fAnalysisType == "ESD") tagPattern = "ESD.tag.root";
113 else if(fAnalysisType == "AOD") tagPattern = "AOD.tag.root";
114 else AliFatal("Only ESD and AOD type is implemented!!!");
116 // Open the working directory
117 void * dirp = gSystem->OpenDirectory(fTagDirName);
118 const char * name = 0x0;
119 // Add all files matching *pattern* to the chain
120 while((name = gSystem->GetDirEntry(dirp))) {
121 if (strstr(name,tagPattern)) {
122 fTagFilename = fTagDirName;
124 fTagFilename += name;
126 fChain->Add(fTagFilename);
127 printf("Tag file %s\n", fTagFilename.Data());
131 //AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
132 // AliDebug(Form("Chained tag files: %d ",fChain->GetEntries()));
138 //___________________________________________________________________________
139 TChain * AliTagAnalysis::ChainGridTags(TGridResult *res) {
140 //Loops overs the entries of the TGridResult
141 //Chains the tags that are stored in the GRID
143 Int_t nEntries = ftagresult->GetEntries();
145 if (! fChain) fChain = new TChain("T");
147 TString gridname = "alien://";
150 for(Int_t i = 0; i < nEntries; i++) {
151 alienUrl = ftagresult->GetKey(i,"turl");
152 fChain->Add(alienUrl);
158 //___________________________________________________________________________
159 TChain *AliTagAnalysis::QueryTags(AliRunTagCuts *runTagCuts,
160 AliLHCTagCuts *lhcTagCuts,
161 AliDetectorTagCuts *detTagCuts,
162 AliEventTagCuts *evTagCuts) {
163 //Queries the tag chain using the defined
164 //event tag cuts from the AliEventTagCuts object
165 //and returns a TChain along with the associated TEventList
166 AliInfo(Form("Querying the tags........"));
169 if(fAnalysisType == "ESD") aliceFile = "esdTree";
170 else if(fAnalysisType == "AOD") aliceFile = "aodTree";
171 else AliFatal("Only ESD and AOD type is implemented!!!");
174 TChain *esdChain = new TChain(aliceFile.Data());
176 fGlobalList = new TEntryList();
178 //Defining tag objects
179 AliRunTag *tag = new AliRunTag;
180 AliEventTag *evTag = new AliEventTag;
181 fChain->SetBranchAddress("AliTAG",&tag);
187 TEntryList* localList = new TEntryList();
191 for(Int_t iEntry = 0; iEntry < fChain->GetEntries(); iEntry++) {
192 fChain->GetEntry(iEntry);
194 if(runTagCuts->IsAccepted(tag)) {
195 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
196 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
198 Int_t iEvents = tag->GetNEvents();
199 const TClonesArray *tagList = tag->GetEventTags();
200 for(Int_t i = 0; i < iEvents; i++) {
201 evTag = (AliEventTag *) tagList->At(i);
202 guid = evTag->GetGUID();
203 turl = evTag->GetTURL();
204 path = evTag->GetPath();
205 localList->SetTreeName(aliceFile.Data());
206 if(turl!="") localList->SetFileName(turl.Data());
207 else localList->SetFileName(path.Data());
209 if(evTagCuts->IsAccepted(evTag)) localList->Enter(i);
211 iAccepted += localList->GetN();
212 if(turl != "") esdChain->AddFile(turl);
213 else if(path != "") esdChain->AddFile(path);
214 fGlobalList->Add(localList);
220 AliInfo(Form("Accepted events: %d", iAccepted));
222 esdChain->SetEntryList(fGlobalList,"ne");
229 //___________________________________________________________________________
230 TChain *AliTagAnalysis::QueryTags(const char *fRunCut,
232 const char *fDetectorCut,
233 const char *fEventCut) {
234 //Queries the tag chain using the defined
235 //event tag cuts from the AliEventTagCuts object
236 //and returns a TChain along with the associated TEventList
237 AliInfo(Form("Querying the tags........"));
240 if(fAnalysisType == "ESD") aliceFile = "esdTree";
241 else if(fAnalysisType == "AOD") aliceFile = "aodTree";
242 else AliFatal("Only ESD and AOD type is implemented!!!");
246 TChain *esdChain = new TChain(aliceFile.Data());
248 fGlobalList = new TEntryList();
250 //Defining tag objects
251 AliRunTag *tag = new AliRunTag;
252 AliEventTag *evTag = new AliEventTag;
253 fChain->SetBranchAddress("AliTAG",&tag);
259 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
260 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
261 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
262 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
264 TEntryList* localList = new TEntryList();
269 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
270 fChain->GetEntry(iTagFiles);
271 if (current != fChain->GetTreeNumber()) {
272 fRunFormula->UpdateFormulaLeaves();
273 fLHCFormula->UpdateFormulaLeaves();
274 fDetectorFormula->UpdateFormulaLeaves();
275 fEventFormula->UpdateFormulaLeaves();
276 current = fChain->GetTreeNumber();
279 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
280 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
281 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
283 Int_t iEvents = fEventFormula->GetNdata();
284 const TClonesArray *tagList = tag->GetEventTags();
285 for(Int_t i = 0; i < iEvents; i++) {
286 evTag = (AliEventTag *) tagList->At(i);
287 guid = evTag->GetGUID();
288 turl = evTag->GetTURL();
289 path = evTag->GetPath();
290 localList->SetTreeName(aliceFile.Data());
291 localList->SetFileName(turl.Data());
292 if(fEventFormula->EvalInstance(i) == 1) localList->Enter(i);
295 if(path != "") esdChain->AddFile(path);
296 else if(turl != "") esdChain->AddFile(turl);
297 fGlobalList->Add(localList);
298 iAccepted += localList->GetN();
304 AliInfo(Form("Accepted events: %d", iAccepted));
305 esdChain->SetEntryList(fGlobalList,"ne");
312 //___________________________________________________________________________
314 AliTagAnalysis::CreateXMLCollection(const char* name,
315 AliRunTagCuts *runTagCuts,
316 AliLHCTagCuts *lhcTagCuts,
317 AliDetectorTagCuts *detTagCuts,
318 AliEventTagCuts *evTagCuts)
320 /// Queries the tag chain using the defined run, lhc, detector and event tag objects
321 /// and create a XML collection named "name.xml"
322 /// if any of the runTagCuts, lhcTagCuts, detTagCuts or evTagCuts is NULL
323 /// check on that object will be skipped.
325 AliInfo(Form("Creating the collection........"));
329 AliError("fChain is NULL. Cannot make a collection from that !");
333 AliXMLCollection collection;
334 collection.SetCollectionName(name);
335 collection.WriteHeader();
341 TEntryList localList;
344 Int_t iRejectedRun = 0;
345 Int_t iRejectedLHC = 0;
346 Int_t iRejectedDet = 0;
347 Int_t iRejectedEvt = 0;
349 Int_t iTotalEvents = 0;
351 Int_t iAcceptedEvtInFile = 0;
352 Int_t iRejectedEvtInFile = 0;
354 //Defining tag objects
355 AliRunTag* tag = new AliRunTag;
356 fChain->SetBranchAddress("AliTAG",&tag);
358 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetListOfFiles()->GetEntries(); ++iTagFiles)
360 fChain->GetEntry(iTagFiles);
362 iTotalEvents += tag->GetNEvents();
365 if ( !runTagCuts || ( runTagCuts && runTagCuts->IsAccepted(tag) ) )
367 if ( !lhcTagCuts || ( lhcTagCuts && lhcTagCuts->IsAccepted(tag->GetLHCTag())) )
369 if ( !detTagCuts || ( detTagCuts && detTagCuts->IsAccepted(tag->GetDetectorTags())) )
372 TIter next(tag->GetEventTags());
373 AliEventTag* evTag(0x0);
374 iRejectedEvtInFile = 0;
375 iAcceptedEvtInFile = 0;
376 while ( ( evTag = static_cast<AliEventTag*>(next()) ) )
378 guid = evTag->GetGUID();
379 turl = evTag->GetTURL();
380 lfn = turl(8,turl.Length());
381 if( !evTagCuts || ( evTagCuts && evTagCuts->IsAccepted(evTag)) )
384 iAcceptedEvtInFile++;
389 ++iRejectedEvtInFile;
393 iAccepted += localList.GetN();
394 collection.WriteBody(iTagFiles+1,guid,lfn,turl,&localList,iAcceptedEvtInFile,iRejectedEvtInFile);
397 iRejectedDet += tag->GetNEvents();
401 iRejectedLHC += tag->GetNEvents();
405 iRejectedRun += tag->GetNEvents();
410 collection.WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
416 //___________________________________________________________________________
417 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name,
420 const char *fDetectorCut,
421 const char *fEventCut) {
422 //Queries the tag chain using the defined
423 //event tag cuts from the AliEventTagCuts object
424 //and returns a XML collection
425 AliInfo(Form("Creating the collection........"));
428 AliXMLCollection *collection = new AliXMLCollection();
429 collection->SetCollectionName(name);
430 collection->WriteHeader();
435 TEntryList* localList = new TEntryList();
439 Int_t iRejectedRun = 0;
440 Int_t iRejectedLHC = 0;
441 Int_t iRejectedDet = 0;
442 Int_t iRejectedEvt = 0;
444 Int_t iTotalEvents = 0;
446 Int_t iAcceptedEvtInFile = 0;
447 Int_t iRejectedEvtInFile = 0;
449 //Defining tag objects
450 AliRunTag *tag = new AliRunTag;
451 AliEventTag *evTag = new AliEventTag;
452 fChain->SetBranchAddress("AliTAG",&tag);
454 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
455 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
456 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
457 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
461 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
463 fChain->GetEntry(iTagFiles);
464 if (current != fChain->GetTreeNumber()) {
465 fRunFormula->UpdateFormulaLeaves();
466 fLHCFormula->UpdateFormulaLeaves();
467 fDetectorFormula->UpdateFormulaLeaves();
468 fEventFormula->UpdateFormulaLeaves();
469 current = fChain->GetTreeNumber();
473 iTotalEvents += tag->GetNEvents();
475 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
476 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
477 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
478 Int_t iEvents = fEventFormula->GetNdata();
479 const TClonesArray *tagList = tag->GetEventTags();
480 iRejectedEvtInFile = 0;
481 iAcceptedEvtInFile = 0;
482 for(Int_t i = 0; i < iEvents; i++) {
483 evTag = (AliEventTag *) tagList->At(i);
484 guid = evTag->GetGUID();
485 turl = evTag->GetTURL();
486 lfn = turl(8,turl.Length());
487 if(fEventFormula->EvalInstance(i) == 1) {
489 iAcceptedEvtInFile++;
493 iRejectedEvtInFile++;
496 collection->WriteBody(iTagFiles+1,guid,lfn,turl,localList,iAcceptedEvtInFile, iRejectedEvtInFile);
497 iAccepted += localList->GetN();
500 iRejectedDet += tag->GetNEvents();
504 iRejectedLHC += tag->GetNEvents();
508 iRejectedRun += tag->GetNEvents();
511 collection->WriteSummary(iTotalEvents, iAccepted, iRejectedRun, iRejectedLHC, iRejectedDet, iRejectedEvt);
512 collection->Export();
518 //___________________________________________________________________________
519 TChain *AliTagAnalysis::GetInputChain(const char* system, const char *wn) {
520 //returns the chain+event list - used in batch sessions
521 // this function will be removed once the new root
522 // improvements are committed
523 TString fsystem = system;
526 TChain *fAnalysisChain = 0;
527 if(fAnalysisType == "ESD") fAnalysisChain = new TChain("esdTree");
528 else if(fAnalysisType == "AOD") fAnalysisChain = new TChain("aodTree");
529 else AliFatal("Only ESD and AOD type is implemented!!!");
532 TEventList *fEventList = new TEventList();
533 AliXMLCollection *collection = AliXMLCollection::Open(wn);
536 while (collection->Next()) {
537 AliInfo(Form("Adding: %s",collection->GetTURL("")));
538 fAnalysisChain->Add(collection->GetTURL(""));
539 TEntryList *list = (TEntryList *)collection->GetEventList("");
540 for(Int_t i = 0; i < list->GetN(); i++) fEventList->Enter(iAccepted+list->GetEntry(i));
542 if(fsystem == "pp") iAccepted += 100;
543 else if(fsystem == "PbPb") iAccepted += 1;
546 fAnalysisChain->SetEventList(fEventList);
548 AliInfo(Form("Number of selected events: %d",fEventList->GetN()));
550 return fAnalysisChain;
553 //___________________________________________________________________________
555 AliTagAnalysis::CreateChainFromCollection(const char* collectionname, const char* treename)
557 /// Build a TChain (with its TEntryList object attached) from an XML collection.
558 /// Returned chain must be deleted by the client.
560 TString streename(treename);
561 if ( streename != "esdTree" && streename != "aodTree" )
563 AliErrorClass("Only esdTree and aodTree implemented so far...");
567 TChain* chain = new TChain(streename.Data());
569 // create the event list for the chain. Will be attached to the chain
570 // which thus becomes the owner of it.
571 TEntryList* elist = new TEntryList;
573 AliXMLCollection* collection = AliXMLCollection::Open(collectionname);
575 // Tag selection summary per file
576 TMap* tagCutSummary = new TMap();
577 tagCutSummary->SetName("TagCutSumm");
583 while (collection->Next())
585 AliDebugClass(1,Form("Adding: %s",collection->GetTURL("")));
586 chain->Add(collection->GetTURL(""));
587 TEntryList *list = collection->GetEventList("");
588 list->SetTreeName(streename.Data());
589 list->SetFileName(collection->GetTURL(""));
591 iAccepted += list->GetN();
592 if (collection->GetCutSumm())
594 tagCutSummary->Add(new TObjString(collection->GetTURL("")), new TObjString(collection->GetCutSumm()));
598 chain->SetEntryList(elist,"ne"); // ne => do not expand tree name and/or file names
600 AliDebugClass(1,Form("Number of selected events: %d",iAccepted));
602 TList *aUserInfo = chain->GetUserInfo();
603 aUserInfo->Add(tagCutSummary);
612 collection->GetCollectionSummary(&iTotalEvents, &iAccEv, &iRejRun, &iRejLHC, &iRejDet, &iRejEvt);
616 sprintf(nstr, "TotalEvents=%i", iTotalEvents);
617 TObjString *iTotStr = new TObjString(nstr);
618 aUserInfo->Add(iTotStr);
620 sprintf(nstr, "AcceptedEvents=%i", iAccepted);
621 TObjString *iAccStr = new TObjString(nstr);
622 aUserInfo->Add(iAccStr);
624 sprintf(nstr, "RejectedRun=%i", iRejRun);
625 TObjString *iRejRunStr = new TObjString(nstr);
626 aUserInfo->Add(iRejRunStr);
628 sprintf(nstr, "RejectedLHC=%i", iRejLHC);
629 TObjString *iRejLHCStr = new TObjString(nstr);
630 aUserInfo->Add(iRejLHCStr);
632 sprintf(nstr, "RejectedDet=%i", iRejDet);
633 TObjString *iRejDetStr = new TObjString(nstr);
634 aUserInfo->Add(iRejDetStr);
636 sprintf(nstr, "RejectedEvt=%i", iRejEvt);
637 TObjString *iRejEvtStr = new TObjString(nstr);
638 aUserInfo->Add(iRejEvtStr);