]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Calls to LoadXXX() modified to try to avoid wrong 0 event (Gustavo)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 12 Apr 2006 16:43:28 +0000 (16:43 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 12 Apr 2006 16:43:28 +0000 (16:43 +0000)
EMCAL/AliEMCALLoader.cxx
EMCAL/AliEMCALLoader.h

index d2b55a4ca2dd8ce9713c243632b8be32bd011368..3d6f818865e639a0aae8ffad904b3abe80cfab16 100644 (file)
@@ -149,42 +149,105 @@ Int_t AliEMCALLoader::CalibrateRaw(Double_t energy, Int_t module,
   return amp;
 }
 
-
 //____________________________________________________________________________ 
-Int_t AliEMCALLoader::GetEvent() {
-  AliLoader::GetEvent();  // First call AliLoader to do all the groundwork
-
-  // Now connect and fill TClonesArray
-
-  // Hits
+Int_t AliEMCALLoader::LoadHits(Option_t* opt) {
+  
+  Int_t status = AliLoader::LoadHits(opt);  // First call AliLoader to do all the groundwork
+  
   TTree *treeH = TreeH();
+  
   if (treeH) {
     treeH->SetBranchAddress(fDetectorName,&fHits);
     if (treeH->GetEntries() > 1)
       AliWarning("Multiple arrays in treeH no longer supported");
     treeH->GetEvent(0);
   }
+  return status;
+}
 
-  // SDigits
+//____________________________________________________________________________ 
+Int_t AliEMCALLoader::LoadSDigits(Option_t* opt) {
+  
+  Int_t status = AliLoader::LoadSDigits(opt);  // First call AliLoader to do all the groundwork
+  
   TTree *treeS = TreeS();
+  
   if (treeS) {
     treeS->SetBranchAddress(fDetectorName,&fSDigits);
     treeS->GetEvent(0);
   }
+  return status;
+}
 
-  // Digits
+//____________________________________________________________________________ 
+Int_t AliEMCALLoader::LoadDigits(Option_t* opt) {
+  
+  Int_t status = AliLoader::LoadDigits(opt);  // First call AliLoader to do all the groundwork
+  
   TTree *treeD = TreeD();
+  
   if (treeD) {
     treeD->SetBranchAddress(fDetectorName,&fDigits);
     treeD->GetEvent(0);
   }
+  return status;
+}
 
-  // RecPoints
+//____________________________________________________________________________ 
+Int_t AliEMCALLoader::LoadRecPoints(Option_t* opt) {
+  
+  Int_t status = AliLoader::LoadRecPoints(opt);  // First call AliLoader to do all the groundwork
+  
   TTree *treeR = TreeR();
   if (treeR) {
     treeR->SetBranchAddress(fgkECARecPointsBranchName,&fRecPoints);
     treeR->GetEvent(0);
   }
+  
+  return status;
+}
 
-  return 0;
+//____________________________________________________________________________ 
+Int_t AliEMCALLoader::GetEvent() {
+  
+  AliLoader::GetEvent();  // First call AliLoader to do all the groundwork
+  
+  // Now connect and fill TClonesArray
+
+  // Hits
+   TTree *treeH = TreeH();
+   
+   if (treeH) {
+     treeH->SetBranchAddress(fDetectorName,&fHits);
+     if (treeH->GetEntries() > 1)
+       AliWarning("Multiple arrays in treeH no longer supported");
+     treeH->GetEvent(0);
+   }
+   
+   
+   // SDigits
+   TTree *treeS = TreeS();
+   if (treeS) {
+     treeS->SetBranchAddress(fDetectorName,&fSDigits);
+     treeS->GetEvent(0);
+   }
+   
+   // Digits
+   TTree *treeD = TreeD();
+   if (treeD) {
+     treeD->SetBranchAddress(fDetectorName,&fDigits);
+     treeD->GetEvent(0);
+   }
+
+   // RecPoints
+   TTree *treeR = TreeR();
+   if (treeR) {
+     treeR->SetBranchAddress(fgkECARecPointsBranchName,&fRecPoints);
+     treeR->GetEvent(0);
+   }
+   
+   return 0;
 }
index 9bfe567a69f110f9751a3efb0b795f36e06168d4..97053fef23f64fd6d675c2c5dfa349d689bedd46 100644 (file)
@@ -50,6 +50,10 @@ class AliEMCALLoader : public AliLoader {
   const AliEMCALLoader & operator = (const AliEMCALLoader & ) {return *this;}
 
   virtual Int_t GetEvent();  // Overload to fill TClonesArray
+  virtual Int_t LoadHits(Option_t* opt);  // Overload to fill TClonesArray
+  virtual Int_t LoadSDigits(Option_t* opt); // Overload to fill TClonesArray
+  virtual Int_t LoadDigits(Option_t* opt); // Overload to fill TClonesArray
+  virtual Int_t LoadRecPoints(Option_t* opt); // Overload to fill TClonesArray
 
   virtual void    CleanHits() const 
     { if (fHits) fHits->Clear(); AliLoader::CleanHits(); }