]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Skip Calibration events by default. If not skip do not clusterize, just fill the...
authorgconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 12 May 2010 13:16:44 +0000 (13:16 +0000)
committergconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 12 May 2010 13:16:44 +0000 (13:16 +0000)
Change print in RawUtils to debug mode

EMCAL/AliEMCALRawUtils.cxx
EMCAL/AliEMCALRecParam.cxx
EMCAL/AliEMCALRecParam.h
EMCAL/AliEMCALReconstructor.cxx

index ada2819f9794160694fe4b584b399881cd78ae95..149e7c1695e12dad0f7de02a8fe9d6b9d73ad455 100644 (file)
@@ -108,7 +108,7 @@ AliEMCALRawUtils::AliEMCALRawUtils(fitAlgorithm fitAlgo)
   if (rl && rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
   } else {
-    AliInfo(Form("Using default geometry in raw reco"));
+    AliDebug(1, Form("Using default geometry in raw reco"));
     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
   }
 
index 126cc97572621778749f7072a431f6b9c14e75e4..ce632feea40ea8e1488552546a20d6279c366deb 100644 (file)
@@ -62,7 +62,8 @@ AliEMCALRecParam::AliEMCALRecParam() :
   fNPedSamples(5), 
   fRemoveBadChannels(kFALSE),
   fFittingAlgorithm(0), 
-  fUseFALTRO(kTRUE)//raw signal
+  fUseFALTRO(kTRUE), 
+  fFitLEDEvents(kFALSE)//raw signal
 {
   // default reco values
   
@@ -254,7 +255,8 @@ AliEMCALRecParam::AliEMCALRecParam(const AliEMCALRecParam& rp) :
   fNPedSamples(rp.fNPedSamples),       
   fRemoveBadChannels(rp.fRemoveBadChannels),
   fFittingAlgorithm(rp.fFittingAlgorithm),  
-  fUseFALTRO(fUseFALTRO) //raw signal
+  fUseFALTRO(rp.fUseFALTRO),
+  fFitLEDEvents(rp.fFitLEDEvents) //raw signal
 {
   //copy constructor
   
@@ -306,7 +308,8 @@ AliEMCALRecParam& AliEMCALRecParam::operator = (const AliEMCALRecParam& rp)
     fNPedSamples       = rp.fNPedSamples; 
     fRemoveBadChannels = rp.fRemoveBadChannels;
     fFittingAlgorithm  = rp.fFittingAlgorithm;
-    fUseFALTRO         = rp.fUseFALTRO;//raw signal
+    fUseFALTRO         = rp.fUseFALTRO;
+       fFitLEDEvents      = rp.fFitLEDEvents;//raw signal
          
     //PID values
     Int_t i, j;
index 9753dd469930866971524311d9f86bad75580c9f..48af05747856cfbe4c778755a9c704f746ae419e 100644 (file)
@@ -99,6 +99,8 @@ class AliEMCALRecParam : public AliDetectorRecoParam
   void SetRemoveBadChannels(Bool_t val)     {fRemoveBadChannels=val; }
   void SetFittingAlgorithm(Int_t val)       {fFittingAlgorithm=val; }
   void SetFALTROUsage(Bool_t val)           {fUseFALTRO=val; }
+  void SetLEDFit(Bool_t val)                {fFitLEDEvents=val; }
+
        
   /* raw signal getters */
   Double_t GetHighLowGainFactor() const {return fHighLowGainFactor;}
@@ -109,7 +111,8 @@ class AliEMCALRecParam : public AliDetectorRecoParam
   Bool_t   GetRemoveBadChannels() const {return fRemoveBadChannels;}
   Int_t    GetFittingAlgorithm()  const {return fFittingAlgorithm; }
   Bool_t   UseFALTRO()            const {return fUseFALTRO; }
-       
+  Bool_t   FitLEDEvents()         const {return fFitLEDEvents; }
+
        virtual void Print(Option_t * option="") const ;
   
   static AliEMCALRecParam* GetDefaultParameters();
@@ -162,10 +165,11 @@ class AliEMCALRecParam : public AliDetectorRecoParam
   Bool_t   fRemoveBadChannels;     // select if bad channels are removed before fitting
   Int_t    fFittingAlgorithm;      // select the fitting algorithm
   Bool_t   fUseFALTRO;             // get FALTRO (trigger) and put it on trigger digits.
-               
+  Bool_t   fFitLEDEvents;          // fit LED events or not
+       
   static TObjArray* fgkMaps;       // ALTRO mappings for RCU0..RCUX
   
-  ClassDef(AliEMCALRecParam,11)     // Reconstruction parameters
+  ClassDef(AliEMCALRecParam,12)     // Reconstruction parameters
     
     } ;
 
index 4a0e894832b2aff6b9351499d9c932abc1702fca..0e5a828b4e766839004ccde1be63b9cea3ee2364 100644 (file)
@@ -156,15 +156,15 @@ void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree)
   ReadDigitsArrayFromTree(digitsTree);
   fgClusterizer->InitParameters();
   fgClusterizer->SetOutput(clustersTree);
-
   AliEMCALTriggerData* trgData = new AliEMCALTriggerData();
        
   Int_t bufferSize = 32000;
 
   if (TBranch* triggerBranch = clustersTree->GetBranch("EMTRG"))
-         triggerBranch->SetAddress(&trgData);
-  else
-         clustersTree->Branch("EMTRG","AliEMCALTriggerData",&trgData,bufferSize);
+               triggerBranch->SetAddress(&trgData);
+       else
+               clustersTree->Branch("EMTRG","AliEMCALTriggerData",&trgData,bufferSize);
 
   TClonesArray *trgDigits = new TClonesArray("AliEMCALRawDigit",1000);
   TBranch *branchdig = digitsTree->GetBranch("EMTRG");
@@ -176,29 +176,33 @@ void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree)
 
   branchdig->SetAddress(&trgDigits);
   branchdig->GetEntry(0);
-
-  Int_t v0M[2] = {0,0};
-  fgTriggerProcessor->Digits2Trigger(trgDigits, v0M, trgData);
        
-  trgDigits->Delete();
-  delete trgDigits; trgDigits = 0x0;
+  //Skip clusterization of LED events
+  if (GetRecParam()->GetEventSpecie()!=AliRecoParam::kCalib){
 
-  if(fgDigitsArr && fgDigitsArr->GetEntries()) {
+         Int_t v0M[2] = {0,0};
+         fgTriggerProcessor->Digits2Trigger(trgDigits, v0M, trgData);
+       
+         if(fgDigitsArr && fgDigitsArr->GetEntries()) {
 
-    fgClusterizer->SetInput(digitsTree);
+                 fgClusterizer->SetInput(digitsTree);
     
-    if(Debug())
-      fgClusterizer->Digits2Clusters("deb all") ;
-    else
-      fgClusterizer->Digits2Clusters("");
+                 if(Debug())
+                         fgClusterizer->Digits2Clusters("deb all") ;
+                 else
+                         fgClusterizer->Digits2Clusters("");
     
-    fgClusterizer->Clear();
-
-  }
+                 fgClusterizer->Clear();
 
+         }//digits array exists and has somethind
+  }//not a LED event
+       
   clustersTree->Fill();        
+  trgDigits->Delete();
+  delete trgDigits; trgDigits = 0x0;
+  delete trgData;   trgData   = 0x0;
 
-  delete trgData; trgData = 0x0;
 }
 
 //____________________________________________________________________________
@@ -217,23 +221,32 @@ void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digits
   Int_t bufsize = 32000;
   digitsTree->Branch("EMCAL", &digitsArr, bufsize);
   digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
-
-  //must be done here because, in constructor, option is not yet known
-  fgRawUtils->SetOption(GetOption());
-
-  fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
-  fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
-  fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
-  fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
-  fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
-  fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
-  fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
-  fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
-  fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
-  fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
        
-  fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData,digitsTrg);
-
+  //Skip calibration events do the rest
+  Bool_t doFit = kTRUE;
+  if ( !(GetRecParam()->FitLEDEvents()) && GetRecParam()->GetEventSpecie()==AliRecoParam::kCalib) doFit = kFALSE;
+  if (doFit){
+         //must be done here because, in constructor, option is not yet known
+         fgRawUtils->SetOption(GetOption());
+
+         fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
+         fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
+         fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
+         fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
+         fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
+         fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
+         fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
+         fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
+         fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
+         fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
+       
+         fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData,digitsTrg);
+  }
+  else{
+       AliDebug(1," Calibration Event, skip!");
+       printf("**** AliEMCALReconstructor::ConvertDigits() Calibration Event, skip!!\n");
+  }
+       
   digitsTree->Fill();
   digitsArr->Delete();
   digitsTrg->Delete();