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);
124 AliInfo(Form("Chained tag files: %d ",fChain->GetEntries()));
128 //___________________________________________________________________________
129 TChain * AliTagAnalysis::ChainGridTags(TGridResult *res) {
130 //Loops overs the entries of the TGridResult
131 //Chains the tags that are stored in the GRID
133 Int_t nEntries = ftagresult->GetEntries();
135 if (! fChain) fChain = new TChain("T");
137 TString gridname = "alien://";
140 for(Int_t i = 0; i < nEntries; i++) {
141 alienUrl = ftagresult->GetKey(i,"turl");
142 fChain->Add(alienUrl);
148 //___________________________________________________________________________
149 TChain *AliTagAnalysis::QueryTags(AliRunTagCuts *runTagCuts,
150 AliLHCTagCuts *lhcTagCuts,
151 AliDetectorTagCuts *detTagCuts,
152 AliEventTagCuts *evTagCuts) {
153 //Queries the tag chain using the defined
154 //event tag cuts from the AliEventTagCuts object
155 //and returns a TChain along with the associated TEventList
156 AliInfo(Form("Querying the tags........"));
159 if(fAnalysisType == "ESD") fAliceFile = "esdTree";
160 else if(fAnalysisType == "AOD") fAliceFile = "aodTree";
161 else AliFatal("Only ESD and AOD type is implemented!!!");
164 TChain *fESDchain = new TChain(fAliceFile.Data());
166 fGlobalList = new TEntryList();
168 //Defining tag objects
169 AliRunTag *tag = new AliRunTag;
170 AliEventTag *evTag = new AliEventTag;
171 fChain->SetBranchAddress("AliTAG",&tag);
178 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
179 fChain->GetEntry(iTagFiles);
180 if(runTagCuts->IsAccepted(tag)) {
181 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
182 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
183 TEntryList *fLocalList = new TEntryList();
184 Int_t iEvents = tag->GetNEvents();
185 const TClonesArray *tagList = tag->GetEventTags();
186 for(Int_t i = 0; i < iEvents; i++) {
187 evTag = (AliEventTag *) tagList->At(i);
188 guid = evTag->GetGUID();
189 turl = evTag->GetTURL();
190 path = evTag->GetPath();
191 fLocalList->SetTreeName(fAliceFile.Data());
192 if(turl!="") fLocalList->SetFileName(turl.Data());
193 else fLocalList->SetFileName(path.Data());
194 if(evTagCuts->IsAccepted(evTag)) fLocalList->Enter(i);
196 if(path != "") fESDchain->AddFile(path);
197 else if(turl != "") fESDchain->AddFile(turl);
198 fGlobalList->Add(fLocalList);
199 iAccepted += fLocalList->GetN();
205 AliInfo(Form("Accepted events: %d",iAccepted));
206 fESDchain->SetEntryList(fGlobalList,"ne");
211 //___________________________________________________________________________
212 TChain *AliTagAnalysis::QueryTags(const char *fRunCut,
214 const char *fDetectorCut,
215 const char *fEventCut) {
216 //Queries the tag chain using the defined
217 //event tag cuts from the AliEventTagCuts object
218 //and returns a TChain along with the associated TEventList
219 AliInfo(Form("Querying the tags........"));
222 if(fAnalysisType == "ESD") fAliceFile = "esdTree";
223 else if(fAnalysisType == "AOD") fAliceFile = "aodTree";
224 else AliFatal("Only ESD and AOD type is implemented!!!");
227 TChain *fESDchain = new TChain(fAliceFile.Data());
229 fGlobalList = new TEntryList();
231 //Defining tag objects
232 AliRunTag *tag = new AliRunTag;
233 AliEventTag *evTag = new AliEventTag;
234 fChain->SetBranchAddress("AliTAG",&tag);
240 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
241 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
242 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
243 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
247 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
248 fChain->GetEntry(iTagFiles);
249 if (current != fChain->GetTreeNumber()) {
250 fRunFormula->UpdateFormulaLeaves();
251 fLHCFormula->UpdateFormulaLeaves();
252 fDetectorFormula->UpdateFormulaLeaves();
253 fEventFormula->UpdateFormulaLeaves();
254 current = fChain->GetTreeNumber();
256 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
257 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
258 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
259 TEntryList *fLocalList = new TEntryList();
260 Int_t iEvents = fEventFormula->GetNdata();
261 const TClonesArray *tagList = tag->GetEventTags();
262 for(Int_t i = 0; i < iEvents; i++) {
263 evTag = (AliEventTag *) tagList->At(i);
264 guid = evTag->GetGUID();
265 turl = evTag->GetTURL();
266 path = evTag->GetPath();
267 fLocalList->SetTreeName(fAliceFile.Data());
268 fLocalList->SetFileName(turl.Data());
269 if(fEventFormula->EvalInstance(i) == 1) fLocalList->Enter(i);
271 iAccepted += fLocalList->GetN();
273 if(path != "") fESDchain->AddFile(path);
274 else if(turl != "") fESDchain->AddFile(turl);
275 fGlobalList->Add(fLocalList);
276 iAccepted += fLocalList->GetN();
281 AliInfo(Form("Accepted events: %d",iAccepted));
282 fESDchain->SetEntryList(fGlobalList,"ne");
287 //___________________________________________________________________________
288 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name,
289 AliRunTagCuts *runTagCuts,
290 AliLHCTagCuts *lhcTagCuts,
291 AliDetectorTagCuts *detTagCuts,
292 AliEventTagCuts *evTagCuts) {
293 //Queries the tag chain using the defined
294 //event tag cuts from the AliEventTagCuts object
295 //and returns a XML collection
296 AliInfo(Form("Creating the collection........"));
298 AliXMLCollection *collection = new AliXMLCollection();
299 collection->SetCollectionName(name);
300 collection->WriteHeader();
306 //Defining tag objects
307 AliRunTag *tag = new AliRunTag;
308 AliEventTag *evTag = new AliEventTag;
309 fChain->SetBranchAddress("AliTAG",&tag);
311 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
313 TEntryList *fList = new TEntryList();
314 fChain->GetEntry(iTagFiles);
315 if(runTagCuts->IsAccepted(tag)) {
316 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
317 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
318 Int_t iEvents = tag->GetNEvents();
319 const TClonesArray *tagList = tag->GetEventTags();
320 for(Int_t i = 0; i < iEvents; i++) {
321 evTag = (AliEventTag *) tagList->At(i);
322 guid = evTag->GetGUID();
323 turl = evTag->GetTURL();
324 lfn = turl(8,turl.Length());
325 if(evTagCuts->IsAccepted(evTag)) fList->Enter(i);
327 collection->WriteBody(iTagFiles+1,guid,lfn,turl,fList);
333 collection->Export();
338 //___________________________________________________________________________
339 Bool_t AliTagAnalysis::CreateXMLCollection(const char* name,
342 const char *fDetectorCut,
343 const char *fEventCut) {
344 //Queries the tag chain using the defined
345 //event tag cuts from the AliEventTagCuts object
346 //and returns a XML collection
347 AliInfo(Form("Creating the collection........"));
349 AliXMLCollection *collection = new AliXMLCollection();
350 collection->SetCollectionName(name);
351 collection->WriteHeader();
357 //Defining tag objects
358 AliRunTag *tag = new AliRunTag;
359 AliEventTag *evTag = new AliEventTag;
360 fChain->SetBranchAddress("AliTAG",&tag);
362 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
363 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
364 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
365 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
368 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
370 TEntryList *fList = new TEntryList();
371 fChain->GetEntry(iTagFiles);
372 if (current != fChain->GetTreeNumber()) {
373 fRunFormula->UpdateFormulaLeaves();
374 fLHCFormula->UpdateFormulaLeaves();
375 fDetectorFormula->UpdateFormulaLeaves();
376 fEventFormula->UpdateFormulaLeaves();
377 current = fChain->GetTreeNumber();
379 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
380 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
381 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
382 Int_t iEvents = fEventFormula->GetNdata();
383 const TClonesArray *tagList = tag->GetEventTags();
384 for(Int_t i = 0; i < iEvents; i++) {
385 evTag = (AliEventTag *) tagList->At(i);
386 guid = evTag->GetGUID();
387 turl = evTag->GetTURL();
388 lfn = turl(8,turl.Length());
389 if(fEventFormula->EvalInstance(i) == 1) fList->Enter(i);
391 collection->WriteBody(iTagFiles+1,guid,lfn,turl,fList);
396 collection->Export();
401 //___________________________________________________________________________
402 Bool_t AliTagAnalysis::CreateAsciiCollection(const char* name,
403 AliRunTagCuts *runTagCuts,
404 AliLHCTagCuts *lhcTagCuts,
405 AliDetectorTagCuts *detTagCuts,
406 AliEventTagCuts *evTagCuts) {
407 //Queries the tag chain using the defined
408 //event tag cuts from the AliEventTagCuts object
409 //and returns a XML collection
410 AliInfo(Form("Creating the collection........"));
421 //Defining tag objects
422 AliRunTag *tag = new AliRunTag;
423 AliEventTag *evTag = new AliEventTag;
424 fChain->SetBranchAddress("AliTAG",&tag);
426 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
428 TEntryList *fList = new TEntryList();
429 fChain->GetEntry(iTagFiles);
430 if(runTagCuts->IsAccepted(tag)) {
431 if(lhcTagCuts->IsAccepted(tag->GetLHCTag())) {
432 if(detTagCuts->IsAccepted(tag->GetDetectorTags())) {
433 Int_t iEvents = tag->GetNEvents();
434 const TClonesArray *tagList = tag->GetEventTags();
435 for(Int_t i = 0; i < iEvents; i++) {
436 evTag = (AliEventTag *) tagList->At(i);
437 guid = evTag->GetGUID();
438 turl = evTag->GetTURL();
439 lfn = turl(8,turl.Length());
440 if(evTagCuts->IsAccepted(evTag)) fList->Enter(i);
442 line0 = guid; line0 += " "; line0 += turl; line0 += " ";
443 for(Int_t i = 0; i < fList->GetN(); i++) {
444 line0 += fList->GetEntry(i);
459 //___________________________________________________________________________
460 Bool_t AliTagAnalysis::CreateAsciiCollection(const char* name,
463 const char *fDetectorCut,
464 const char *fEventCut) {
465 //Queries the tag chain using the defined
466 //event tag cuts from the AliEventTagCuts object
467 //and returns a XML collection
468 AliInfo(Form("Creating the collection........"));
479 //Defining tag objects
480 AliRunTag *tag = new AliRunTag;
481 AliEventTag *evTag = new AliEventTag;
482 fChain->SetBranchAddress("AliTAG",&tag);
484 TTreeFormula *fRunFormula = new TTreeFormula("fRun",fRunCut,fChain);
485 TTreeFormula *fLHCFormula = new TTreeFormula("fLHC",fLHCCut,fChain);
486 TTreeFormula *fDetectorFormula = new TTreeFormula("fDetector",fDetectorCut,fChain);
487 TTreeFormula *fEventFormula = new TTreeFormula("fEvent",fEventCut,fChain);
490 for(Int_t iTagFiles = 0; iTagFiles < fChain->GetEntries(); iTagFiles++) {
492 TEntryList *fList = new TEntryList();
493 fChain->GetEntry(iTagFiles);
494 if (current != fChain->GetTreeNumber()) {
495 fRunFormula->UpdateFormulaLeaves();
496 fLHCFormula->UpdateFormulaLeaves();
497 fDetectorFormula->UpdateFormulaLeaves();
498 fEventFormula->UpdateFormulaLeaves();
499 current = fChain->GetTreeNumber();
501 if(fRunFormula->EvalInstance(iTagFiles) == 1) {
502 if(fLHCFormula->EvalInstance(iTagFiles) == 1) {
503 if(fDetectorFormula->EvalInstance(iTagFiles) == 1) {
504 Int_t iEvents = fEventFormula->GetNdata();
505 const TClonesArray *tagList = tag->GetEventTags();
506 for(Int_t i = 0; i < iEvents; i++) {
507 evTag = (AliEventTag *) tagList->At(i);
508 guid = evTag->GetGUID();
509 turl = evTag->GetTURL();
510 lfn = turl(8,turl.Length());
511 if(fEventFormula->EvalInstance(i) == 1) fList->Enter(i);
513 line0 = guid; line0 += " "; line0 += turl; line0 += " ";
514 for(Int_t i = 0; i < fList->GetN(); i++) {
515 line0 += fList->GetEntry(i);
529 //___________________________________________________________________________
530 TChain *AliTagAnalysis::GetInputChain(const char* system, const char *wn) {
531 //returns the chain+event list - used in batch sessions
532 // this function will be removed once the new root
533 // improvements are committed
534 TString fsystem = system;
537 TChain *fAnalysisChain = 0;
538 if(fAnalysisType == "ESD") fAnalysisChain = new TChain("esdTree");
539 else if(fAnalysisType == "AOD") fAnalysisChain = new TChain("aodTree");
540 else AliFatal("Only ESD and AOD type is implemented!!!");
543 TEventList *fEventList = new TEventList();
544 AliXMLCollection *collection = AliXMLCollection::Open(wn);
547 while (collection->Next()) {
548 AliInfo(Form("Adding: %s",collection->GetTURL("")));
549 fAnalysisChain->Add(collection->GetTURL(""));
550 TEntryList *list = (TEntryList *)collection->GetEventList("");
551 for(Int_t i = 0; i < list->GetN(); i++) fEventList->Enter(iAccepted+list->GetEntry(i));
553 if(fsystem == "pp") iAccepted += 100;
554 else if(fsystem == "PbPb") iAccepted += 1;
557 fAnalysisChain->SetEventList(fEventList);
559 AliInfo(Form("Number of selected events: %d",fEventList->GetN()));
561 return fAnalysisChain;
564 //___________________________________________________________________________
565 TChain *AliTagAnalysis::GetChainFromCollection(const char* collectionname,
566 const char* treename) {
567 //returns the TChain+TEntryList object- used in batch sessions
568 TString fAliceFile = treename;
570 TChain *fAnalysisChain = 0;
571 if(fAliceFile == "esdTree") fAnalysisChain = new TChain("esdTree");
572 else if(fAliceFile == "aodTree") fAnalysisChain = new TChain("aodTree");
573 else AliFatal("Inconsistent tree name - use esdTree or aodTree!");
576 fGlobalList = new TEntryList();
577 AliXMLCollection *collection = AliXMLCollection::Open(collectionname);
580 while (collection->Next()) {
581 AliInfo(Form("Adding: %s",collection->GetTURL("")));
582 fAnalysisChain->Add(collection->GetTURL(""));
583 TEntryList *list = (TEntryList *)collection->GetEventList("");
584 list->SetTreeName(fAliceFile.Data());
585 list->SetFileName(collection->GetTURL(""));
586 fGlobalList->Add(list);
587 iAccepted += list->GetN();
590 fAnalysisChain->SetEntryList(fGlobalList,"ne");
592 AliInfo(Form("Number of selected events: %d",iAccepted));
594 return fAnalysisChain;