]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Some additional updates for running under proof.
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 24 Apr 2007 14:05:15 +0000 (14:05 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 24 Apr 2007 14:05:15 +0000 (14:05 +0000)
JETAN/AliAnalysisTaskJets.cxx
JETAN/AliAnalysisTaskJets.h
JETAN/AliJetESDReader.h
JETAN/AliJetFinder.cxx
JETAN/AliJetFinder.h
JETAN/AliJetHadronCorrection.cxx
JETAN/AliUA1JetFinderV1.h

index cb0e5b7a89302ad93c89fc41f1cfd3d8d6dfcf28..5ec6b6ad51db8becb339c5415a7592058726ae11 100644 (file)
@@ -31,7 +31,8 @@ AliAnalysisTaskJets::AliAnalysisTaskJets():
     fDebug(0),
     fJetFinder(0x0),
     fChain(0x0),
-    fESD(0x0)
+    fESD(0x0),
+    fTreeJ(0x0)
 {
   // Default constructor
 }
@@ -41,17 +42,25 @@ AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
     fDebug(0),
     fJetFinder(0x0),
     fChain(0x0),
-    fESD(0x0)
+    fESD(0x0),
+    fTreeJ(0x0)
 {
   // Default constructor
     DefineInput (0, TChain::Class());
-//    DefineOutput(0, TTree::Class());
+    DefineOutput(0, TTree::Class());
+}
+
+void AliAnalysisTaskJets::CreateOutputObjects()
+{
+// Create the output container
+    fTreeJ = fJetFinder->MakeTreeJ("TreeJ");
 }
 
 void AliAnalysisTaskJets::Init()
 {
     // Initialization
     if (fDebug > 1) printf("AnalysisTaskJets::Init() \n");
+
     // Call configuration file
     gROOT->LoadMacro("ConfigJetAnalysis.C");
     fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
@@ -64,6 +73,7 @@ void AliAnalysisTaskJets::ConnectInputData(Option_t */*option*/)
 // Connect the input data
 //
     if (fDebug > 1) printf("AnalysisTaskJets::ConnectInputData() \n");
+    
     fChain = (TChain*)GetInputData(0);
     fJetFinder->ConnectTree(fChain);
     fJetFinder->WriteHeaders();
@@ -74,9 +84,9 @@ void AliAnalysisTaskJets::Exec(Option_t */*option*/)
 // Execute analysis for current event
 //
     Long64_t ientry = fChain->GetReadEntry();
-     if (fDebug > 1) printf("Analysing event # %5d \n", (Int_t) ientry);
-    
+    if (fDebug > 1) printf("Analysing event # %5d \n", (Int_t) ientry);
     fJetFinder->ProcessEvent(ientry);
+    PostData(0, fTreeJ);
 }
 
 void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
@@ -84,6 +94,7 @@ void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
 // Terminate analysis
 //
     if (fDebug > 1) printf("AnalysisJets: Terminate() \n");
+   
     if (fJetFinder) fJetFinder->FinishRun();
 }
 
index 840810df1692dac898cab3289281f8f6ad89ccf7..02dba0470bf8b97fea796c821b66d3bbda800b3f 100644 (file)
@@ -18,6 +18,7 @@ class AliAnalysisTaskJets : public AliAnalysisTask
     virtual ~AliAnalysisTaskJets() {;}
     // Implementation of interface methods
     virtual void ConnectInputData(Option_t *option = "");
+    virtual void CreateOutputObjects();
     virtual void Init();
     virtual void Exec(Option_t *option);
     virtual void Terminate(Option_t *option);
@@ -28,6 +29,7 @@ class AliAnalysisTaskJets : public AliAnalysisTask
     AliJetFinder* fJetFinder; //  Pointer to the jet finder 
     TChain*       fChain;     //! chained files
     AliESD*       fESD;       //! ESD
+    TTree*        fTreeJ;     //  tree of reconstructed jets
     
     ClassDef(AliAnalysisTaskJets, 1); // Analysis task for standard jet analysis
 };
index 3598eed168a8a5abd08452207e6215ba6f593a97..015f34c51131e96a3262e9ad32078c10879f7a02 100755 (executable)
@@ -45,8 +45,8 @@ class AliJetESDReader : public AliJetReader
   void InitParameters();
  protected:
   AliJetDummyGeo             *fGeom;             //! EMCAL Geometry 
-  TChain                     *fChain;            //!chain for reconstructed tracks
-  AliESD                     *fESD;              //!pointer to esd
+  TChain                     *fChain;            //! chain for reconstructed tracks
+  AliESD                     *fESD;              //! pointer to esd
   AliJetHadronCorrectionv1   *fHadCorr;          //! Pointer to Hadron Correction Object 
   AliJetGrid                 *fTpcGrid;          //! Pointer to grid object
   AliJetGrid                 *fEmcalGrid;        //! Pointer to grid object
index 31f886b1fb560b2a6f414c0174d3063b69bae01d..21d09e46ad2d87df42996081399e6f989097a70a 100644 (file)
@@ -34,7 +34,7 @@
 ClassImp(AliJetFinder)
 
 AliJetFinder::AliJetFinder():
-    fJetT(0),
+    fTreeJ(0),
     fPlotMode(kFALSE),
     fJets(0),
     fGenJets(0),
@@ -75,8 +75,8 @@ AliJetFinder::~AliJetFinder()
 
 void AliJetFinder::SetOutputFile(const char *name)
 {
-  // opens output file 
-  fOut = new TFile(name,"recreate");
+  //  opens output file 
+  //  fOut = new TFile(name,"recreate");
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -100,20 +100,16 @@ void AliJetFinder::SetPlotMode(Bool_t b)
 }
 
 ////////////////////////////////////////////////////////////////////////
-
-void AliJetFinder::WriteJetsToFile(Int_t i)
+TTree* AliJetFinder::MakeTreeJ(char* name)
 {
-  // Writes the jets to file
-  fOut->cd();
-  char hname[30];
-  sprintf(hname,"TreeJ%d",i);
-  fJetT = new TTree(hname,"AliJet");
-  fJetT->Branch("FoundJet",&fJets,1000);
-  fJetT->Branch("GenJet",&fGenJets,1000);
-  fJetT->Branch("LeadingPart",&fLeading,1000);
-  fJetT->Fill();
-  fJetT->Write(hname);
-  delete fJetT;
+    // Create the tree for reconstructed jets
+    fOut = new TFile("jets.root","recreate");
+    fOut->cd();
+    fTreeJ = new TTree(name, "AliJet");
+    fTreeJ->Branch("FoundJet",   &fJets,   1000);
+    fTreeJ->Branch("GenJet",     &fGenJets,1000);
+    fTreeJ->Branch("LeadingPart",&fLeading,1000);
+    return fTreeJ;
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -161,7 +157,11 @@ void AliJetFinder::Run()
           if(debug > 1) printf("In FindJets() routine: find jets with fUnitArray !!!\n");
           FindJets();
       }
-      if (fOut) WriteJetsToFile(i);
+      if (fOut) {
+         fOut->cd();
+         fTreeJ->Fill();
+      }
+      
       if (fPlots) fPlots->FillHistos(fJets);
       fLeading->Reset();
       fGenJets->ClearJets();
@@ -207,12 +207,17 @@ Bool_t AliJetFinder::ProcessEvent(Long64_t entry)
 //
 // Process one event
 //
-    printf("<<<<< Processing Event %5d >>>>> \n", (Int_t) entry);
+    Int_t debug  = fReader->GetReaderHeader()->GetDebug();
+    if (debug > 0) printf("<<<<< Processing Event %5d >>>>> \n", (Int_t) entry);
     Bool_t ok = fReader->FillMomentumArray(entry);
     if (!ok) return kFALSE;
     fLeading->FindLeading(fReader);
     FindJets();
-    if (fOut)   WriteJetsToFile(entry);
+    if (fOut) {
+       fOut->cd();
+       fTreeJ->Fill();
+    }
+    
     if (fPlots) fPlots->FillHistos(fJets);
     fLeading->Reset();
     fGenJets->ClearJets();
@@ -223,16 +228,18 @@ Bool_t AliJetFinder::ProcessEvent(Long64_t entry)
 void AliJetFinder::FinishRun()
 {
     // Finish a run
-  if (fPlots) {
-      fPlots->Normalize();
-      fPlots->PlotHistos();
-      if (fOut) {
-         fOut->cd();
-         fPlots->Write();
-         fOut->Close();
-      }   
-  } else {
-      if (fOut) fOut->Close();
-  }
+    if (fPlots) {
+       fPlots->Normalize();
+       fPlots->PlotHistos();
+    }
+    
+    if (fOut) {
+        fOut->cd();
+        fTreeJ->Write();
+        if (fPlots) {
+            fPlots->Write();
+        }
+        fOut->Close();
+    }
 }
 
index bd068ef3be7290c4d6b41cee0bebc621afde4045..a1822c847ce64d2e95f7b506f4e3ff551b9f21db 100755 (executable)
@@ -36,11 +36,9 @@ class AliJetFinder : public TObject
   virtual void SetOutputFile(const char *name="jets.root");
   virtual void SetJetReader(AliJetReader* r) {fReader=r;}
   virtual void SetJetHeader(AliJetHeader* h) {fHeader=h;}
-
   // others
   virtual void   PrintJets();
   virtual void   Run();
-  virtual void   WriteJetsToFile(Int_t i);
   virtual void   WriteRHeaderToFile();  
   // the following have to be implemented for each specific finder
   virtual void Init() {}
@@ -52,13 +50,15 @@ class AliJetFinder : public TObject
   virtual Bool_t ProcessEvent(Long64_t entry);
   virtual void   FinishRun();
   virtual void   ConnectTree(TTree* tree);
+  virtual TTree* MakeTreeJ(char* name);
   virtual void   WriteHeaders();
-
+  virtual void   WriteJetsToFile() {;}
+  virtual void   TestJet() {;}
  protected:
   AliJetFinder(const AliJetFinder& rJetFinder);
   AliJetFinder& operator = (const AliJetFinder& rhsf);
-  TTree* fJetT;                  //! pointer to jet tree
-  Bool_t fPlotMode;              //!  do you want control plots?
+  TTree* fTreeJ;                 //! pointer to jet tree
+  Bool_t fPlotMode;              //! do you want control plots?
   AliJet* fJets;                 //! pointer to jet class
   AliJet* fGenJets;              //! pointer to generated jets
   AliLeading*   fLeading;        //! pointer to leading particle data 
index efd5a20365b0256d51f6d16a6465f3c2dfbfa754..d81a37d9228b70c2cd717dc583863cfeb12026e2 100644 (file)
@@ -31,5 +31,7 @@
 
 ClassImp(AliJetHadronCorrection)
 
+AliJetHadronCorrection::AliJetHadronCorrection() {}
+
 AliJetHadronCorrection::AliJetHadronCorrection(const char *name,const char *title) 
 :TNamed(name,title) { }
index 342b02273a72d460bb97717cbbc659fcf8cb8cca..b8e346542f6e37444e1b28326298e718245a982e 100644 (file)
@@ -55,8 +55,7 @@ class AliUA1JetFinderV1 : public AliJetFinder
  protected:
   AliUA1JetFinderV1(const AliUA1JetFinderV1& rJetF1);
   AliUA1JetFinderV1& operator = (const AliUA1JetFinderV1& rhsf);
-  TH2F           * fLego;           // Lego Histo
-
+  TH2F           * fLego;           //Lego Histo
   ClassDef(AliUA1JetFinderV1,1)
 };