Added protection against occasional non-LED STANDALONE runs
authorpolicheh <policheh@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 28 Aug 2008 15:39:42 +0000 (15:39 +0000)
committerpolicheh <policheh@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 28 Aug 2008 15:39:42 +0000 (15:39 +0000)
PHOS/AliPHOSDA2.cxx
PHOS/AliPHOSDA2.h
PHOS/AliPHOSPreprocessor.cxx
PHOS/PHOSBCMda.cxx
PHOS/PHOSLEDda.cxx

index 70e499d51a4b5f5a3f0ff807a69bbe8deb80e6f3..c4a86f50bc659da0b8b6e12c4d90e6c3af58d4cb 100644 (file)
@@ -5,7 +5,7 @@ ClassImp(AliPHOSDA2)
 
 //----------------------------------------------------------------
 AliPHOSDA2::AliPHOSDA2(int module) : TNamed(),
- fHistoFile(0),fMod(module)
+ fHistoFile(0),fFiredCells(0),fMod(module)
 
 {
   // Create AliPHOSDA2 ("Bad channels finder") object.
@@ -35,11 +35,13 @@ AliPHOSDA2::AliPHOSDA2(int module) : TNamed(),
   fMaps[0]=0;
   fMaps[1]=0;
 
+  fFiredCells = new TH1I("fFiredCells","Number of fired cells per event",100,0,1000);
+
 }
 
 //-------------------------------------------------------------------
 AliPHOSDA2::AliPHOSDA2(const AliPHOSDA2& da) : TNamed(da),
-  fHistoFile(0),fMod(da.fMod)
+  fHistoFile(0),fFiredCells(0),fMod(da.fMod)
 {
   // Copy constructor.
 
@@ -70,6 +72,7 @@ AliPHOSDA2::AliPHOSDA2(const AliPHOSDA2& da) : TNamed(da),
     fMaps[1] = 0;
   
   fHistoFile = new TFile(da.GetName(),"recreate");
+  fFiredCells = new TH1I(*da.fFiredCells);
   
 }
 
@@ -112,6 +115,11 @@ AliPHOSDA2& AliPHOSDA2::operator= (const AliPHOSDA2& da)
       fMaps[1] = da.fMaps[1];
     } 
     
+    if(fFiredCells) {
+      delete fFiredCells;
+      fFiredCells = da.fFiredCells;
+    }
+    
   }
   
   return *this;
@@ -162,6 +170,12 @@ void AliPHOSDA2::FillQualityHistograms(Float_t quality[64][56][2])
 
 }
 
+//-------------------------------------------------------------------
+void  AliPHOSDA2::FillFiredCellsHistogram(Int_t nCells)
+{
+  fFiredCells->Fill(nCells);
+}
+
 //-------------------------------------------------------------------
 void AliPHOSDA2::UpdateHistoFile()
 {
@@ -206,5 +220,7 @@ void AliPHOSDA2::UpdateHistoFile()
   fMaps[0]->Write(fMaps[0]->GetName(),TObject::kWriteDelete);
   fMaps[1]->Write(fMaps[1]->GetName(),TObject::kWriteDelete);
 
+  fFiredCells->Write();
+
 }
 
index a0d4ff560b7bf0a59d309cc7d05dc23d28391262..86b7237d05918b7c777231d08b29f08fd40c5d2f 100644 (file)
@@ -16,16 +16,20 @@ class AliPHOSDA2 : public TNamed {
   ~AliPHOSDA2();
   
   void  FillQualityHistograms(Float_t quality[64][56][2]);
+  void  FillFiredCellsHistogram(Int_t nCells);
   Int_t GetModule() { return fMod; }
   void  UpdateHistoFile();
 
   const TH1F* GetQualityHistogram(Int_t X, Int_t Z, Int_t gain) const
   { return fHQuality[X][Z][gain]; }
+
+  const TH1I* GetFiredCellsHistogram() { return fFiredCells; }
   
  private:
 
   TFile* fHistoFile;            // root file to store histograms in
   TH1F* fHQuality[64][56][2];   // "quality" for high and low gains
+  TH1I* fFiredCells;            // Number of fired cells pre event.
   Int_t fMod;                   // PHOS module number (0..4)
   TH2F* fMaps[2];               // 2D quality map for low and high gains.
   
index fd3baf7d9fa773bfaef8a4fe1ad5ed080ca80c63..9a0e628c465f1549daec60f21a423209a793c5fb 100644 (file)
@@ -130,6 +130,7 @@ Bool_t AliPHOSPreprocessor::ProcessLEDRun()
   
   TIter iter(list);
   TObjString *source;
+  Bool_t firedOK = kFALSE;
   
   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
     
@@ -145,6 +146,25 @@ Bool_t AliPHOSPreprocessor::ProcessLEDRun()
       return kFALSE;
     }
 
+    TH1I* fFiredCells = (TH1I*)f.Get("fFiredCells");
+
+    if(!fFiredCells) {
+      Log(Form("fFiredCells histogram not found in %s, skip this source.",fileName.Data()));
+      firedOK=kFALSE;
+    }
+    else {
+      const Double_t nFiredCells = fFiredCells->GetMean();
+      if(nFiredCells>100 && nFiredCells<600) { 
+       firedOK = kTRUE;
+       Log(Form("Number of fired cells per event is %d.",nFiredCells));
+      }
+      else {
+       Log(Form("Number of fired cells per event is %d, possibly not a LED run!",nFiredCells));
+       firedOK = kFALSE;
+      }
+      
+    }
+
     const Int_t nMod=5; // 1:5 modules
     const Int_t nCol=56; //1:56 columns in each module
     const Int_t nRow=64; //1:64 rows in each module
@@ -153,12 +173,14 @@ Bool_t AliPHOSPreprocessor::ProcessLEDRun()
       for(Int_t col=0; col<nCol; col++) {
        for(Int_t row=0; row<nRow; row++) {
 
-         //High Gain to Low Gain ratio 
-          Float_t ratio = HG2LG(mod,row,col,&f);
-          calibData.SetHighLowRatioEmc(mod+1,col+1,row+1,ratio);
-         if(ratio != 16.)
-           AliInfo(Form("mod %d iX %d iZ %d  ratio %.3f\n",mod,row,col,ratio));
-         
+         //High Gain to Low Gain ratio
+         if(firedOK) {
+           Float_t ratio = HG2LG(mod,row,col,&f);
+           calibData.SetHighLowRatioEmc(mod+1,col+1,row+1,ratio);
+           if(ratio != 16.)
+             AliInfo(Form("mod %d iX %d iZ %d  ratio %.3f\n",mod,row,col,ratio));
+         }       
+
          if(clb) {
            Double_t coeff = clb->GetADCchannelEmc(mod+1,col+1,row+1);
            calibData.SetADCchannelEmc(mod+1,col+1,row+1,coeff);
@@ -295,8 +317,22 @@ Bool_t AliPHOSPreprocessor::DoFindBadChannelsEmc(Int_t system, TList* list, AliP
       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
       return kFALSE;
     }
+
+    TH1I* fFiredCells = (TH1I*)f.Get("fFiredCells");
+
+    if(!fFiredCells) {
+      Log(Form("fFiredCells histogram not found in %s, skip this source.",fileName.Data()));
+      continue;
+    }
+
+    const Double_t nFiredCells = fFiredCells->GetMean();
+
+    if(nFiredCells<100. || nFiredCells>600. ) {
+      Log(Form("Number of fired cells per event is %d, possibly not a LED run!",nFiredCells));
+      continue; // not a LED run!
+    }
     
-    Log(Form("Begin check for bad channels."));
+    Log(Form("Number of fired cells per event %d. Begin check for bad channels.",nFiredCells));
     
     for(Int_t mod=0; mod<5; mod++) {
       for(Int_t iX=0; iX<64; iX++) {
index 4e0b2e676a102a7b3ec0fef8ae277c611ee39e23..1d56b54fdfe26cc84e1f3026b4829275cab46885 100644 (file)
@@ -106,13 +106,14 @@ int main(int argc, char **argv) {
 
   AliRawReader *rawReader = NULL;
 
-  AliPHOSDA2 da2(2); // DA2 ("Checking for bad channels") for module2
+  AliPHOSDA2* da2 = new AliPHOSDA2(2); // DA2 ("Checking for bad channels") for module2
   
   Float_t q[64][56][2];
 
   Int_t gain = -1;
   Int_t X = -1;
   Int_t Z = -1;
+  Int_t nFired = -1;
 
   /* main loop (infinite) */
   for(;;) {
@@ -153,6 +154,8 @@ int main(int argc, char **argv) {
        }
       }
 
+      nFired = 0;
+
       rawReader = new AliRawReaderDate((void*)event);
       AliPHOSRawDecoderv1 dc(rawReader,mapping);
       dc.SubtractPedestals(kTRUE);
@@ -167,10 +170,14 @@ int main(int argc, char **argv) {
          gain = 1;
        
        q[X][Z][gain] = dc.GetSampleQuality();
+
+       if(gain && dc.GetEnergy()>40)
+         nFired++;
        
       }
       
-      da2.FillQualityHistograms(q);       
+      da2->FillQualityHistograms(q);
+      da2->FillFiredCellsHistogram(nFired);
       //da1.UpdateHistoFile();
       
       delete rawReader;     
@@ -190,14 +197,12 @@ int main(int argc, char **argv) {
   }
   
   for(Int_t i = 0; i < 4; i++) delete mapping[i];  
+
+  /* Be sure that all histograms are saved */
+  delete da2;
   
   /* Store output files to the File Exchange Server */
-  char localfile[128];
-
-  for(Int_t iMod=0; iMod<5; iMod++) {
-    sprintf(localfile,"PHOS_Module%d_BCM.root",iMod);
-    daqDA_FES_storeFile(localfile,"BAD_CHANNELS");
-  }
+  daqDA_FES_storeFile("PHOS_Module2_BCM.root","BAD_CHANNELS");
 
   return status;
 }
index bdabd7f50a95e2719b9c01c88e63eb905f4d7b5c..3e652f194d77bb5a0ff8f7b1a13c2db45ae852f1 100644 (file)
@@ -109,6 +109,9 @@ int main(int argc, char **argv) {
   Int_t gain = -1;
   Int_t X = -1;
   Int_t Z = -1;
+  Int_t nFired = -1;
+
+  TH1I fFiredCells("fFiredCells","Number of fired cells per event",100,0,1000);
 
   /* main loop (infinite) */
   for(;;) {
@@ -150,6 +153,8 @@ int main(int argc, char **argv) {
        }
       }
 
+      nFired = 0;
+
       rawReader = new AliRawReaderDate((void*)event);
 //       AliPHOSRawDecoderv1 dc(rawReader,mapping);
       AliPHOSRawDecoder dc(rawReader,mapping);
@@ -167,10 +172,14 @@ int main(int argc, char **argv) {
        e[X][Z][gain] = dc.GetEnergy();
        t[X][Z][gain] = dc.GetTime();
        
+       if(gain && dc.GetEnergy()>40)
+         nFired++;
+
       }
       
       da1.FillHistograms(e,t);
-      
+      fFiredCells.Fill(nFired);
+
       delete rawReader;     
       nevents_physics++;
     }
@@ -221,14 +230,11 @@ int main(int argc, char **argv) {
     }
   }
   
+  fFiredCells.Write();
   f->Close();
   
   /* Store output files to the File Exchange Server */
-  
-  for(Int_t iMod=0; iMod<5; iMod++) {
-    sprintf(localfile,"PHOS_Module%d_LED.root",iMod);
-    daqDA_FES_storeFile(localfile,"LED");
-  }
+  daqDA_FES_storeFile("PHOS_Module2_LED.root","LED");
   
   printf("%d physics events of %d total processed.\n",nevents_physics,nevents_total);
   printf("%d histograms has >10 entries in maximum, max. is %d entries ",nGood,nMax);