]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Extended Notify() functions to read cross section and trials also from the histogram...
authorkleinb <kleinb@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 3 Sep 2009 09:27:04 +0000 (09:27 +0000)
committerkleinb <kleinb@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 3 Sep 2009 09:27:04 +0000 (09:27 +0000)
PWG4/JetTasks/AliAnalysisTaskJFSystematics.cxx
PWG4/JetTasks/AliAnalysisTaskJetSpectrum.cxx

index dc47ee3ece027be17e1ee012fbe6d7772e683690..820273249aca942ef3b57e9628c17a3f64c29640 100644 (file)
@@ -34,6 +34,7 @@
 #include <TH3F.h>
 #include <TProfile.h>
 #include <TList.h>
+#include <TKey.h>
 #include <TLorentzVector.h>
 #include <TClonesArray.h>
 #include  "TDatabasePDG.h"
@@ -166,14 +167,14 @@ AliAnalysisTaskJFSystematics::AliAnalysisTaskJFSystematics(const char* name):
 
 Bool_t AliAnalysisTaskJFSystematics::Notify()
 {
-
-  //
+//
   // Implemented Notify() to read the cross sections
   // and number of trials from pyxsec.root
   // 
   TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
   Double_t xsection = 0;
   UInt_t   ntrials  = 0;
+  Float_t   ftrials  = 0;
   if(tree){
     TFile *curfile = tree->GetCurrentFile();
     if (!curfile) {
@@ -187,31 +188,53 @@ Bool_t AliAnalysisTaskJFSystematics::Notify()
 
     TString fileName(curfile->GetName());
     if(fileName.Contains("AliESDs.root")){
-        fileName.ReplaceAll("AliESDs.root", "pyxsec.root");
+        fileName.ReplaceAll("AliESDs.root", "");
     }
-    else if(fileName.Contains("AliAOD.root")){
-        fileName.ReplaceAll("AliAOD.root", "pyxsec.root");
+    else if(fileName.Contains("AliAODs.root")){
+        fileName.ReplaceAll("AliAODs.root", "");
     }
     else if(fileName.Contains("galice.root")){
         // for running with galice and kinematics alone...                      
-        fileName.ReplaceAll("galice.root", "pyxsec.root");
+        fileName.ReplaceAll("galice.root", "");
     }
-
-    TFile *fxsec = TFile::Open(fileName.Data());
+    TFile *fxsec = TFile::Open(Form("%s%s",fileName.Data(),"pyxsec.root"));
     if(!fxsec){
-      Printf("%s:%d %s not found in the Input",(char*)__FILE__,__LINE__,fileName.Data());
-      // no a severe condition, will happen for minbias data
-      return kTRUE;
+      if(fDebug>0)Printf("%s:%d %s not found in the Input",(char*)__FILE__,__LINE__,Form("%s%s",fileName.Data(),"pyxsec.root"));
+      // next trial fetch the histgram file
+      fxsec = TFile::Open(Form("%s%s",fileName.Data(),"pyxsec_hists.root"));
+      if(!fxsec){
+       // not a severe condition
+       if(fDebug>0)Printf("%s:%d %s not found in the Input",(char*)__FILE__,__LINE__,Form("%s%s",fileName.Data(),"pyxsec_hists.root"));        
+       return kTRUE;
+      }
+      else{
+       // find the tlist we want to be independtent of the name so use the Tkey
+       TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0); 
+       if(!key){
+         if(fDebug>0)Printf("%s:%d key not found in the file",(char*)__FILE__,__LINE__);       
+         return kTRUE;
+       }
+       TList *list = dynamic_cast<TList*>(key->ReadObj());
+       if(!list){
+         if(fDebug>0)Printf("%s:%d key is not a tlist",(char*)__FILE__,__LINE__);      
+         return kTRUE;
+       }
+       xsection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
+       ftrials  = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
+      }
     }
-    TTree *xtree = (TTree*)fxsec->Get("Xsection");
-    if(!xtree){
-      Printf("%s:%d tree not found in the pyxsec.root",(char*)__FILE__,__LINE__);
+    else{
+      TTree *xtree = (TTree*)fxsec->Get("Xsection");
+      if(!xtree){
+       Printf("%s:%d tree not found in the pyxsec.root",(char*)__FILE__,__LINE__);
+      }
+      xtree->SetBranchAddress("xsection",&xsection);
+      xtree->SetBranchAddress("ntrials",&ntrials);
+      ftrials = ntrials;
+      xtree->GetEntry(0);
     }
-    xtree->SetBranchAddress("xsection",&xsection);
-    xtree->SetBranchAddress("ntrials",&ntrials);
-    xtree->GetEntry(0);
     fh1Xsec->Fill("<#sigma>",xsection);
-    fh1Trials->Fill("#sum{ntrials}",ntrials);
+    fh1Trials->Fill("#sum{ntrials}",ftrials);
   }
   return kTRUE;
 }
index ea5379aa0e37bd5a9821f4526abbbd404127225c..e97dabff202a21a4a803d45779c39e0349f91aa7 100644 (file)
@@ -26,6 +26,7 @@
 #include <TInterpreter.h>
 #include <TChain.h>
 #include <TFile.h>
+#include <TKey.h>
 #include <TH1F.h>
 #include <TH2F.h>
 #include <TH3F.h>
@@ -157,6 +158,7 @@ Bool_t AliAnalysisTaskJetSpectrum::Notify()
   TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
   Double_t xsection = 0;
   UInt_t   ntrials  = 0;
+  Float_t   ftrials  = 0;
   if(tree){
     TFile *curfile = tree->GetCurrentFile();
     if (!curfile) {
@@ -170,30 +172,53 @@ Bool_t AliAnalysisTaskJetSpectrum::Notify()
 
     TString fileName(curfile->GetName());
     if(fileName.Contains("AliESDs.root")){
-        fileName.ReplaceAll("AliESDs.root", "pyxsec.root");
+        fileName.ReplaceAll("AliESDs.root", "");
     }
-    else if(fileName.Contains("AliAOD.root")){
-        fileName.ReplaceAll("AliAOD.root", "pyxsec.root");
+    else if(fileName.Contains("AliAODs.root")){
+        fileName.ReplaceAll("AliAODs.root", "");
     }
     else if(fileName.Contains("galice.root")){
         // for running with galice and kinematics alone...                      
-        fileName.ReplaceAll("galice.root", "pyxsec.root");
+        fileName.ReplaceAll("galice.root", "");
     }
-    TFile *fxsec = TFile::Open(fileName.Data());
+    TFile *fxsec = TFile::Open(Form("%s%s",fileName.Data(),"pyxsec.root"));
     if(!fxsec){
-      Printf("%s:%d %s not found in the Input",(char*)__FILE__,__LINE__,fileName.Data());
-      // no a severe condition
-      return kTRUE;
+      if(fDebug>0)Printf("%s:%d %s not found in the Input",(char*)__FILE__,__LINE__,Form("%s%s",fileName.Data(),"pyxsec.root"));
+      // next trial fetch the histgram file
+      fxsec = TFile::Open(Form("%s%s",fileName.Data(),"pyxsec_hists.root"));
+      if(!fxsec){
+       // not a severe condition
+       if(fDebug>0)Printf("%s:%d %s not found in the Input",(char*)__FILE__,__LINE__,Form("%s%s",fileName.Data(),"pyxsec_hists.root"));        
+       return kTRUE;
+      }
+      else{
+       // find the tlist we want to be independtent of the name so use the Tkey
+       TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0); 
+       if(!key){
+         if(fDebug>0)Printf("%s:%d key not found in the file",(char*)__FILE__,__LINE__);       
+         return kTRUE;
+       }
+       TList *list = dynamic_cast<TList*>(key->ReadObj());
+       if(!list){
+         if(fDebug>0)Printf("%s:%d key is not a tlist",(char*)__FILE__,__LINE__);      
+         return kTRUE;
+       }
+       xsection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
+       ftrials  = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
+      }
     }
-    TTree *xtree = (TTree*)fxsec->Get("Xsection");
-    if(!xtree){
-      Printf("%s:%d tree not found in the pyxsec.root",(char*)__FILE__,__LINE__);
+    else{
+      TTree *xtree = (TTree*)fxsec->Get("Xsection");
+      if(!xtree){
+       Printf("%s:%d tree not found in the pyxsec.root",(char*)__FILE__,__LINE__);
+      }
+      xtree->SetBranchAddress("xsection",&xsection);
+      xtree->SetBranchAddress("ntrials",&ntrials);
+      ftrials = ntrials;
+      xtree->GetEntry(0);
     }
-    xtree->SetBranchAddress("xsection",&xsection);
-    xtree->SetBranchAddress("ntrials",&ntrials);
-    xtree->GetEntry(0);
     fh1Xsec->Fill("<#sigma>",xsection);
-    fh1Trials->Fill("#sum{ntrials}",ntrials);
+    fh1Trials->Fill("#sum{ntrials}",ftrials);
   }
   
   return kTRUE;