]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Changes needed to run analysis with proof.
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 5 Apr 2007 16:28:52 +0000 (16:28 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 5 Apr 2007 16:28:52 +0000 (16:28 +0000)
JETAN/AliAnalysisTaskJets.cxx
JETAN/AliAnalysisTaskJets.h
JETAN/AliJetESDReader.cxx
JETAN/AliJetESDReader.h
JETAN/AliJetFinder.cxx
JETAN/AliJetFinder.h
JETAN/AliUA1JetFinder.h

index 7c12a94fff177c67e3e022f8684191f5d92ac071..cb0e5b7a89302ad93c89fc41f1cfd3d8d6dfcf28 100644 (file)
@@ -27,6 +27,15 @@ ClassImp(AliAnalysisTaskJets)
 
 ////////////////////////////////////////////////////////////////////////
 
+AliAnalysisTaskJets::AliAnalysisTaskJets():
+    fDebug(0),
+    fJetFinder(0x0),
+    fChain(0x0),
+    fESD(0x0)
+{
+  // Default constructor
+}
+
 AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
     AliAnalysisTask(name, "AnalysisTaskJets"),
     fDebug(0),
@@ -36,20 +45,25 @@ AliAnalysisTaskJets::AliAnalysisTaskJets(const char* name):
 {
   // Default constructor
     DefineInput (0, TChain::Class());
-    DefineOutput(0, TH1::Class());
+//    DefineOutput(0, TTree::Class());
 }
 
-void AliAnalysisTaskJets::ConnectInputData(Option_t */*option*/)
+void AliAnalysisTaskJets::Init()
 {
-// Initialisation
-//
-    if (fDebug > 1) printf("AnalysisJets::Init() \n");
-
+    // Initialization
+    if (fDebug > 1) printf("AnalysisTaskJets::Init() \n");
     // Call configuration file
     gROOT->LoadMacro("ConfigJetAnalysis.C");
     fJetFinder = (AliJetFinder*) gInterpreter->ProcessLine("ConfigJetAnalysis()");
     // Initialise Jet Analysis
     fJetFinder->Init();
+}
+
+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();
@@ -69,7 +83,7 @@ void AliAnalysisTaskJets::Terminate(Option_t */*option*/)
 {
 // Terminate analysis
 //
-    printf("AnalysisJets: Terminate() \n");
+    if (fDebug > 1) printf("AnalysisJets: Terminate() \n");
     if (fJetFinder) fJetFinder->FinishRun();
 }
 
index 2705581acba91c9fcc4b5179c8935ccc18895419..840810df1692dac898cab3289281f8f6ad89ccf7 100644 (file)
@@ -12,22 +12,24 @@ class TChain;
 
 class AliAnalysisTaskJets : public AliAnalysisTask
 {
- public: 
-  AliAnalysisTaskJets(const char* name);
-  virtual ~AliAnalysisTaskJets() {;}
-  // Implementation of interface methods
-  virtual void ConnectInputData(Option_t *option = "");
-  virtual void Exec(Option_t *option);
-  virtual void Terminate(Option_t *option);
-  virtual void SetDebugLevel(Int_t level) {fDebug = level;}
-  
+ public:
+    AliAnalysisTaskJets();
+    AliAnalysisTaskJets(const char* name);
+    virtual ~AliAnalysisTaskJets() {;}
+    // Implementation of interface methods
+    virtual void ConnectInputData(Option_t *option = "");
+    virtual void Init();
+    virtual void Exec(Option_t *option);
+    virtual void Terminate(Option_t *option);
+    virtual void SetDebugLevel(Int_t level) {fDebug = level;}
+    
  private:
-  Int_t         fDebug;     //  Debug flag
-  AliJetFinder* fJetFinder; //  Pointer to the jet finder 
-  TChain*       fChain;     //! chained files
-  AliESD*       fESD;       //! ESD
-  
-  ClassDef(AliAnalysisTaskJets, 1); // Analysis task for standard jet analysis
+    Int_t         fDebug;     //  Debug flag
+    AliJetFinder* fJetFinder; //  Pointer to the jet finder 
+    TChain*       fChain;     //! chained files
+    AliESD*       fESD;       //! ESD
+    
+    ClassDef(AliAnalysisTaskJets, 1); // Analysis task for standard jet analysis
 };
  
 #endif
index af92943b993c14e348f27354ea53047b61ca19d4..5e8e817cf9e79067968982811aa31949aed66d4f 100755 (executable)
@@ -106,7 +106,7 @@ void AliJetESDReader::OpenInputFiles()
   
   int nMax = fChain->GetEntries(); 
 
-  printf("\nTotal number of events in chain= %d \n",nMax);
+  printf("\n AliJetESDReader: Total number of events in chain= %d \n",nMax);
   
   // set number of events in header
   if (fReaderHeader->GetLastEvent() == -1)
@@ -123,7 +123,7 @@ void AliJetESDReader::ConnectTree(TTree* tree) {
      
      fChain->SetBranchAddress("ESD",    &fESD);
      Int_t nMax = fChain->GetEntries(); 
-     printf("\nTotal number of events in chain= %5d \n", nMax);
+     printf("\n AliJetESDReader: Total number of events in chain= %5d \n", nMax);
      // set number of events in header
      if (fReaderHeader->GetLastEvent() == -1)
         fReaderHeader->SetLastEvent(nMax);
index 242384b387c981e61a7cf2326c5dfd15486ebcda..3598eed168a8a5abd08452207e6215ba6f593a97 100755 (executable)
@@ -44,12 +44,12 @@ class AliJetESDReader : public AliJetReader
   void SetEMCALGeometry();
   void InitParameters();
  protected:
-  AliJetDummyGeo             *fGeom;             //!EMCAL Geometry 
-  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
+  AliJetDummyGeo             *fGeom;             //! EMCAL Geometry 
+  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
   Float_t                     fPtCut;            // Pt cut for tracks to minimise background contribution
   Int_t                       fHCorrection;      // Hadron correction flag
   Int_t                       fNumUnits;         // Number of units in the unit object array
index 1e9769be3314e6eb42ed4271bbfe29b771ec7d72..31f886b1fb560b2a6f414c0174d3063b69bae01d 100644 (file)
 ClassImp(AliJetFinder)
 
 AliJetFinder::AliJetFinder():
+    fJetT(0),
     fPlotMode(kFALSE),
     fJets(0),
     fGenJets(0),
     fLeading(0),
     fReader(0x0),
+    fHeader(0x0),
     fPlots(0x0),
     fOut(0x0)
     
@@ -105,13 +107,13 @@ void AliJetFinder::WriteJetsToFile(Int_t i)
   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;
+  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;
 }
 
 ////////////////////////////////////////////////////////////////////////
index 3f30aa60903fe3696162aca9fb23bd137044eb2b..bd068ef3be7290c4d6b41cee0bebc621afde4045 100755 (executable)
@@ -57,15 +57,15 @@ class AliJetFinder : public TObject
  protected:
   AliJetFinder(const AliJetFinder& rJetFinder);
   AliJetFinder& operator = (const AliJetFinder& rhsf);
-
-  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 
-  AliJetReader* fReader;         // pointer to reader
-  AliJetHeader* fHeader;         // pointer to header
-  AliJetControlPlots* fPlots;    // pointer to control plots
-  TFile* fOut;                   // output file
+  TTree* fJetT;                  //! 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 
+  AliJetReader* fReader;         //  pointer to reader
+  AliJetHeader* fHeader;         //  pointer to header
+  AliJetControlPlots* fPlots;    //! pointer to control plots
+  TFile* fOut;                   //! output file
   ClassDef(AliJetFinder,2)
 };
 
index afae577d061d73849d520dc3ad65a6e39a2c1915..3006224d6ba19fb4aa2eb59998af9bda2214628b 100755 (executable)
@@ -40,7 +40,7 @@ class AliUA1JetFinder : public AliJetFinder
   AliUA1JetFinder& operator = (const AliUA1JetFinder& ruaf);
 
   AliUA1JetHeader* fHeader;         // pointer to jet header
-  TH2F           * fLego;           //! Lego Histo
+  TH2F           * fLego;           // Lego Histo
 
   ClassDef(AliUA1JetFinder,1)
 };