]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - JETAN/AliJetFinder.cxx
Some additional updates for running under proof.
[u/mrichter/AliRoot.git] / JETAN / AliJetFinder.cxx
index 7a4d5f56a0ac1c0572300c1b607bb2c75b42132a..21d09e46ad2d87df42996081399e6f989097a70a 100644 (file)
@@ -13,6 +13,8 @@
  * provided "as is" without express or implied warranty.                  *
  **************************************************************************/
    
+/* $Id$ */
+
 //---------------------------------------------------------------------
 // Jet finder base class
 // manages the search for jets 
 ClassImp(AliJetFinder)
 
 AliJetFinder::AliJetFinder():
+    fTreeJ(0),
     fPlotMode(kFALSE),
     fJets(0),
     fGenJets(0),
     fLeading(0),
     fReader(0x0),
+    fHeader(0x0),
     fPlots(0x0),
     fOut(0x0)
     
@@ -71,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");
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -96,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);
-  TTree* jetT = new TTree(hname,"AliJet");
-  jetT->Branch("FoundJet",&fJets,1000);
-  jetT->Branch("GenJet",&fGenJets,1000);
-  jetT->Branch("LeadingPart",&fLeading,1000);
-  jetT->Fill();
-  jetT->Write(hname);
-  delete jetT;
+    // 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;
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -157,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();
@@ -203,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();
@@ -219,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();
+    }
 }