]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add method to merge only certain directories from on output file
authorkleinb <kleinb@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 15 Nov 2010 08:33:12 +0000 (08:33 +0000)
committerkleinb <kleinb@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 15 Nov 2010 08:33:12 +0000 (08:33 +0000)
PWG4/JetTasks/AliAnalysisHelperJetTasks.cxx
PWG4/JetTasks/AliAnalysisHelperJetTasks.h

index aaf3800799757f5f01b83a3e1ad1e9f5786639a1..8aedc2ba3c1d330477ca3506d7f0fad3b7405aad 100644 (file)
@@ -7,6 +7,8 @@
 #include "TH1F.h"
 #include "TProfile.h"
 #include "THnSparse.h"
+#include "TSeqCollection.h"
+#include "TMethodCall.h"
 #include "TFile.h"
 #include "TString.h"
 #include "AliMCEvent.h"
@@ -216,6 +218,90 @@ void AliAnalysisHelperJetTasks::GetClosestJets(AliAODJet *genJets,const Int_t &k
 }
 
 
+void  AliAnalysisHelperJetTasks::MergeOutputDirs(const char* cFiles,const char* cPattern,char *cOutFile){
+  // Routine to merge only directories containing the pattern
+  //
+  TString outFile(cOutFile);
+  if(outFile.Length()==0)outFile = Form("%s.root",cPattern);
+  ifstream in1;
+  in1.open(cFiles);
+  // open all files
+  TList fileList;
+  char cFile[240];
+  if(in1>>cFile){// only open the first file
+    fileList.Add(TFile::Open(cFile));
+  }
+
+  TFile *fOut = 0;
+  if(fileList.GetEntries()){// open the first file
+    TFile* fIn = dynamic_cast<TFile*>(fileList.At(0));
+    if(fIn){
+      Printf("Input File not Found");
+    }
+    // fetch the keys for the directories
+    TList *ldKeys = fIn->GetListOfKeys();
+    for(int id = 0;id < ldKeys->GetEntries();id++){
+      // check if it is a directory
+      TKey *dKey = (TKey*)ldKeys->At(id);
+      TDirectory *dir = dynamic_cast<TDirectory*>(dKey->ReadObj());
+      if(dir){
+       TString dName(dir->GetName());
+       if(dName.Contains(cPattern)){
+         // open new file if first match
+         if(!fOut){
+           fOut = new TFile(outFile.Data(),"RECREATE");
+         }
+         // merge all the stuff that is in this directory
+         TList *llKeys = dir->GetListOfKeys();
+         TList *tmplist;
+         TMethodCall callEnv;
+         for(int il = 0;il < llKeys->GetEntries();il++){
+           TKey *lKey = (TKey*)llKeys->At(id);
+           TObject *object = dynamic_cast<TObject*>(lKey->ReadObj());
+           //  from TSeqCollections::Merge
+           if(!object)continue;
+           // If current object is not mergeable just skip it
+           if (!object->IsA()) {
+             continue;
+           }
+           callEnv.InitWithPrototype(object->IsA(), "Merge", "TCollection*");
+           if (!callEnv.IsValid()) {
+             continue;
+           }
+           // Current object mergeable - get corresponding objects in input lists
+           tmplist = new TList();
+           for(int i = 1; i < fileList.GetEntries();i++){
+             TDirectory *fInTmp = dynamic_cast<TDirectory*>(fileList.At(i)); 
+             if(!fInTmp){
+               Printf("File %d not found",i);
+               continue;
+             }
+             TDirectory *dTmp = (TDirectory*)fInTmp->Get(dName.Data());
+             if(!dTmp){
+               Printf("Directory %s not found",dName.Data());
+               continue;
+             }
+             TObject* oTmp = (TObject*)dTmp->Get(object->GetName());
+             if(!oTmp){
+               Printf("Object %s not found",object->GetName());
+               continue;
+             }
+             tmplist->Add(oTmp);
+           }
+           callEnv.SetParam((Long_t) tmplist);
+           callEnv.Execute(object);
+           delete tmplist;tmplist = 0;
+           fOut->cd();
+           object->Write();
+         }
+       }
+      }
+    }
+    if(fOut){
+      fOut->Close();
+    }
+  }
+}
 
 void  AliAnalysisHelperJetTasks::MergeOutput(char* cFiles, char* cDir, char *cList,char *cOutFile,Bool_t bUpdate){
 
index c991187381cde87a0883daa7322d61dd791d4de6..09cf8b60c9552ed962f0e0966efb3ff57f9b2727 100644 (file)
@@ -55,6 +55,8 @@ class AliAnalysisHelperJetTasks : public TObject {
                             Int_t *iRecIndex,
                             Int_t iDebug, Float_t maxDist = 0.5);
 
+  static void MergeOutputDirs(const char* cFiles,const char* cPattern,char *cOutFile); // merges all directories containing the pattern
+
   static void MergeOutput(char* cFiles, char* cDir = "",char *cList = "",char* cOutFile ="allpt.root",Bool_t bUpdate = false); // Merges the files in the input text file  needs the two histograms fh1PtHard_Trials, fh1Xsec and the name of the input list
   static Bool_t PythiaInfoFromFile(const char* currFile,Float_t &fXsec,Float_t &fTrials);// get the cross section and the trails either from pyxsec.root or from pysec_hists.root
   static Bool_t PrintDirectorySize(const char* currFile); // print the size of the output on a given file