]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
ctor modified
authorbnandi <bnandi@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 20 Feb 2007 09:04:11 +0000 (09:04 +0000)
committerbnandi <bnandi@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 20 Feb 2007 09:04:11 +0000 (09:04 +0000)
PMD/AliPMDPreprocessor.cxx
PMD/AliPMDPreprocessor.h

index dd979bb96c5f3c72cc56e78b7f0f9da65379c558..b7c8660522be4194fdc8e7ea397b934f1814a7e6 100644 (file)
@@ -16,8 +16,8 @@
 ClassImp(AliPMDPreprocessor)
 
 //______________________________________________________________________________________________
-AliPMDPreprocessor::AliPMDPreprocessor(const char* detector, AliShuttleInterface* shuttle) :
-  AliPreprocessor(detector, shuttle)
+AliPMDPreprocessor::AliPMDPreprocessor(AliShuttleInterface* shuttle) :
+  AliPreprocessor("PMD", shuttle)
 {
   // constructor
 }
@@ -57,7 +57,7 @@ UInt_t AliPMDPreprocessor::Process(TMap* pdaqAliasMap)
         TList* filesources = GetFileSources(kDAQ, "PMDGAINS");
 
         if(!filesources) {
-                AliError(Form("No sources found for PMDGAINS for run %d !", fRun));
+                Log(Form("No sources found for PMDGAINS!"));
                 return 0;
         }
 
@@ -70,18 +70,29 @@ UInt_t AliPMDPreprocessor::Process(TMap* pdaqAliasMap)
         UInt_t result = 0;
        TString filename;
         while((source=dynamic_cast<TObjString*> (iter.Next()))){
-                AliInfo(Form("\n\n Getting file #%d\n",++i));
                 filename = GetFile(kDAQ, "PMDGAINS", source->GetName());
-                if(!filename.Data()) {
-                        AliError(Form("Error retrieving file from source %s failed!", source->GetName()));
+                if(filename.Length() == 0) {
+                        Log(Form("Error retrieving file from source %s failed!", source->GetName()));
                         delete filesources;
                         return 0;
                 }
 
+                Log(Form("File with id PMDGAINS got from %s", source->GetName()));
                Int_t DET,SM,ROW,COL;
                Float_t GAIN;
                TFile *f= new TFile(filename.Data());
-               TTree *tree = (TTree*)f->Get("ic");
+               if(!f || !f->IsOpen()) 
+               {
+                       Log(Form("Error opening file with Id PMDGAINS from source %s!", source->GetName()));
+                       return 0;
+               } 
+               TTree *tree = dynamic_cast<TTree *> (f->Get("ic"));
+               if (!tree) 
+               {
+                       Log("Could not find object \"ic\" in DAQ file!");
+                       return 0;
+               }
+               
                tree->SetBranchAddress("DET",       &DET);
                tree->SetBranchAddress("SM",        &SM);
                tree->SetBranchAddress("ROW",        &ROW);
index 5babe9e466095adb0c35442856ec88872b2fa74a..d5cb4fcd4afbce49751df88835f9c3eee8f08fd2 100644 (file)
@@ -9,7 +9,7 @@
 class AliPMDPreprocessor : public AliPreprocessor
 {
   public:
-    AliPMDPreprocessor(const char* detector, AliShuttleInterface* shuttle);
+    AliPMDPreprocessor(AliShuttleInterface* shuttle);
     virtual ~AliPMDPreprocessor();
 
   protected:
@@ -19,7 +19,7 @@ class AliPMDPreprocessor : public AliPreprocessor
   private:
 //    AliPMDDataDAQ *fData;    // CDB class that stores the data
 
-    ClassDef(AliPMDPreprocessor, 0);
+    ClassDef(AliPMDPreprocessor, 1);
 };
 
 #endif