]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ANALYSIS/AliTriggerAnalysis.cxx
Adding a reminder for coders
[u/mrichter/AliRoot.git] / ANALYSIS / AliTriggerAnalysis.cxx
index 53e1b1d3f694d0df60e1c192e43ae750c82a4d09..f6876641f29dc8434711cd54e278ff166d9dc2fc 100644 (file)
@@ -54,6 +54,14 @@ AliTriggerAnalysis::AliTriggerAnalysis() :
   fV0HwAdcThr(2.5),
   fV0HwWinLow(61.5),
   fV0HwWinHigh(86.5),
+  fZDCCutRefSum(-568.5),
+  fZDCCutRefDelta(-2.1),
+  fZDCCutSigmaSum(3.25),
+  fZDCCutSigmaDelta(2.25),
+  fZDCCutRefSumCorr(-65.5),
+  fZDCCutRefDeltaCorr(-2.1),
+  fZDCCutSigmaSumCorr(6.0),
+  fZDCCutSigmaDeltaCorr(1.2),
   fDoFMD(kTRUE),
   fFMDLowCut(0.2),
   fFMDHitCut(0.5),
@@ -63,6 +71,8 @@ AliTriggerAnalysis::AliTriggerAnalysis() :
   fHistV0C(0),
   fHistZDC(0),    
   fHistTDCZDC(0),    
+  fHistTimeZDC(0),    
+  fHistTimeCorrZDC(0),    
   fHistFMDA(0),    
   fHistFMDC(0),   
   fHistFMDSingle(0),
@@ -109,8 +119,18 @@ AliTriggerAnalysis::~AliTriggerAnalysis()
   }
   if (fHistTDCZDC)
   {
-    delete fHistZDC;
-    fHistZDC = 0;
+    delete fHistTDCZDC;
+    fHistTDCZDC = 0;
+  }
+  if (fHistTimeZDC)
+  {
+    delete fHistTimeZDC;
+    fHistTimeZDC = 0;
+  }
+  if (fHistTimeCorrZDC)
+  {
+    delete fHistTimeCorrZDC;
+    fHistTimeCorrZDC = 0;
   }
 
   if (fHistFMDA)
@@ -164,6 +184,8 @@ void AliTriggerAnalysis::EnableHistograms()
   fHistV0C = new TH1F("fHistV0C", "V0C;leading time (ns);events", 400, -100, 100);
   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
   fHistTDCZDC = new TH1F("fHistTDCZDC", "ZDC;TDC bits;events", 32, -0.5, 32-0.5);
+  fHistTimeZDC = new TH2F("fHistTimeZDC", "ZDC;TDC timing A+C vs C-A; events", 120,-30,30,120,-600,-540);
+  fHistTimeCorrZDC = new TH2F("fHistTimeCorrZDC", "ZDC;Corrected TDC timing A+C vs C-A; events", 120,-30,30,120,-100,-40);
   
   // TODO check limits
   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
@@ -195,6 +217,8 @@ const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger)
     case kMB3 : str = "MB3"; break;
     case kSPDGFO : str = "SPD GFO"; break;
     case kSPDGFOBits : str = "SPD GFO Bits"; break;
+    case kSPDGFOL0 : str = "SPD GFO L0 (first layer)"; break;
+    case kSPDGFOL1 : str = "SPD GFO L1 (second layer)"; break;
     case kV0A : str = "V0 A BB"; break;
     case kV0C : str = "V0 C BB"; break;
     case kV0OR : str = "V0 OR BB"; break;
@@ -209,6 +233,9 @@ const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger)
     case kFPANY : str = "SPD GFO | V0 | ZDC | FMD"; break;
     case kNSD1 : str = "NSD1"; break;
     case kMB1Prime: str = "MB1prime"; break;
+    case kZDCTDCA : str = "ZDC TDC A"; break;
+    case kZDCTDCC : str = "ZDC TDC C"; break;
+    case kZDCTime : str = "ZDC Time Cut"; break;
     default: str = ""; break;
   }
    
@@ -298,6 +325,176 @@ Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t
   ULong64_t trigmask = aEsd->GetTriggerMask();
   return (trigmask & (1ull << (tclass-1)));
 }
+  
+Int_t AliTriggerAnalysis::EvaluateTrigger(const AliESDEvent* aEsd, Trigger trigger)
+{
+  // evaluates a given trigger "offline"
+  // trigger combinations are not supported, for that see IsOfflineTriggerFired
+  
+  UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
+  Bool_t offline = kFALSE;
+  if (trigger & kOfflineFlag)
+    offline = kTRUE;
+  
+  Int_t decision = 0;
+  switch (triggerNoFlags)
+  {
+    case kSPDGFO:
+    {
+      decision = SPDFiredChips(aEsd, (offline) ? 0 : 1);
+      break;
+    }
+    case kSPDGFOL0:
+    {
+      decision = SPDFiredChips(aEsd, (offline) ? 0 : 1, kFALSE, 1);
+      break;
+    }
+    case kSPDGFOL1:
+    {
+      decision = SPDFiredChips(aEsd, (offline) ? 0 : 1, kFALSE, 2);
+      break;
+    }
+    case kV0A:
+    {
+      if (V0Trigger(aEsd, kASide, !offline) == kV0BB)
+        decision = 1;
+      break;
+    }
+    case kV0C:
+    {
+      if (V0Trigger(aEsd, kCSide, !offline) == kV0BB)
+        decision = 1;
+      break;
+    }
+    case kV0ABG:
+    {
+      if (V0Trigger(aEsd, kASide, !offline) == kV0BG)
+        decision = 1;
+      break;
+    }
+    case kV0CBG:
+    {
+      if (V0Trigger(aEsd, kCSide, !offline) == kV0BG)
+        decision = 1;
+      break;
+    }
+    case kZDCA:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      if (ZDCTrigger(aEsd, kASide))
+        decision = 1;
+      break;
+    }
+    case kZDCC:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      if (ZDCTrigger(aEsd, kCSide))
+        decision = 1;
+      break;
+    }
+    case kZDCTDCA:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      if (ZDCTDCTrigger(aEsd, kASide))
+        decision = 1;
+      break;
+    }
+    case kZDCTDCC:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      if (ZDCTDCTrigger(aEsd, kCSide))
+        decision = 1;
+      break;
+    }
+    case kZDCTime:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      if (ZDCTimeTrigger(aEsd))
+        decision = 1;
+      break;
+    }
+    case kFMDA:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      if (FMDTrigger(aEsd, kASide))
+        decision = 1;
+      break;
+    }
+    case kFMDC:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      if (FMDTrigger(aEsd, kCSide))
+        decision = 1;
+      break;
+    }
+    case kCTPV0A:
+    {
+      if (offline)
+        AliFatal(Form("Offline trigger not available for trigger %d", triggerNoFlags));
+      if (IsL0InputFired(aEsd, 2))
+        decision = 1;
+      break;
+    }
+    case kCTPV0C:
+    {
+      if (offline)
+        AliFatal(Form("Offline trigger not available for trigger %d", triggerNoFlags));
+      if (IsL0InputFired(aEsd, 3))
+        decision = 1;
+      break;
+    }
+    case kTPCLaserWarmUp:
+    {
+      if (!offline)
+        AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
+      return IsLaserWarmUpTPCEvent(aEsd);
+    }
+    default:
+    {
+      AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
+    }
+  }  
+  
+  return decision;
+}
+
+Bool_t AliTriggerAnalysis::IsLaserWarmUpTPCEvent(const AliESDEvent* esd)
+{
+  //
+  // This function flags noisy TPC events which can happen during laser warm-up.
+  //
+  
+  Int_t trackCounter = 0;
+  for (Int_t i=0; i<esd->GetNumberOfTracks(); ++i) 
+  {
+    AliESDtrack *track = esd->GetTrack(i);
+    if (!track) 
+      continue;
+      
+    if (track->GetTPCNcls() < 30) continue;
+    if (TMath::Abs(track->Eta()) > 0.005) continue;
+    if (track->Pt() < 4) continue;
+    if (track->GetKinkIndex(0) > 0) continue;
+    
+    UInt_t status = track->GetStatus();
+    if ((status&AliESDtrack::kITSrefit)==AliESDtrack::kITSrefit) continue; // explicitly ask for tracks without ITS refit
+    if ((status&AliESDtrack::kTPCrefit)!=AliESDtrack::kTPCrefit) continue;
+    
+    if (track->GetTPCsignal() > 10) continue;          // explicitly ask for tracks without dE/dx
+    
+    trackCounter++;
+  }
+  if (trackCounter > 15) 
+    return kTRUE;
+  return kFALSE;
+}
 
 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
 {
@@ -529,11 +726,15 @@ Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Ch
     else if (op == '|') tcl =tcl || inpnext;
     else {
        AliError(Form("Syntax error in %s", tclass));
-       tcltokens->Delete();
+       delete tcltokens;
+       tcltokens = 0;
+       //      tcltokens->Delete();
        return kFALSE;
     }
   }
-  tcltokens->Delete();
+  delete tcltokens;
+  tcltokens = 0;
+       //  tcltokens->Delete();
   return tcl;
 }
 
@@ -585,7 +786,7 @@ void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd)
   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
   ZDCTDCTrigger(aEsd,kASide,kFALSE,kFALSE,kTRUE);
-
+  ZDCTimeTrigger(aEsd,kTRUE);
 
   AliESDZDC* zdcData = aEsd->GetESDZDC();
   if (zdcData)
@@ -675,7 +876,7 @@ Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, B
     if(layer == 2)
       firstChip = 400;
 
-    for (Int_t i=0; i<1200; i++)
+    for (Int_t i=firstChip; i<lastChip; i++)
     {
       if (mult->TestFastOrFiredChips(i) == kTRUE)
       {
@@ -923,28 +1124,105 @@ Bool_t AliTriggerAnalysis::ZDCTDCTrigger(const AliESDEvent* aEsd, AliceSide side
   
   AliESDZDC *esdZDC = aEsd->GetESDZDC();
 
-  Bool_t tdc[32] = {kFALSE};
-  for(Int_t itdc=0; itdc<32; itdc++){
-    for(Int_t i=0; i<4; i++){
-      if (0.025*esdZDC->GetZDCTDCData(itdc, i) != 0){
-       tdc[itdc] = kTRUE;
+  Bool_t zdcNA = kFALSE;
+  Bool_t zdcNC = kFALSE;
+  Bool_t zdcPA = kFALSE;
+  Bool_t zdcPC = kFALSE;
+
+  if (fMC) {
+        // If it's MC, we use the energy
+    Double_t minEnergy = 0;
+    Double_t  zNCEnergy = esdZDC->GetZDCN1Energy();
+    Double_t  zPCEnergy = esdZDC->GetZDCP1Energy();
+    Double_t  zNAEnergy = esdZDC->GetZDCN2Energy();
+    Double_t  zPAEnergy = esdZDC->GetZDCP2Energy();
+    zdcNA = (zNAEnergy>minEnergy);
+    zdcNC = (zNCEnergy>minEnergy);
+    zdcPA = (zPAEnergy>minEnergy);
+    zdcPC = (zPCEnergy>minEnergy);
+
+  }
+  else {
+
+    Bool_t tdc[32] = {kFALSE};
+    for(Int_t itdc=0; itdc<32; itdc++){
+      for(Int_t i=0; i<4; i++){
+       if (esdZDC->GetZDCTDCData(itdc, i) != 0){
+         tdc[itdc] = kTRUE;
+       }
       }
+      if(fillHists && tdc[itdc]) {
+       fHistTDCZDC->Fill(itdc);
+      }    
     }
-    if(fillHists && tdc[itdc]) {
-      fHistTDCZDC->Fill(itdc);
-    }    
+    zdcNA = tdc[12];
+    zdcNC = tdc[10];
+    zdcPA = tdc[13];
+    zdcPC = tdc[11];
   }
-  Bool_t zdcNA = tdc[12];
-  Bool_t zdcNC = tdc[10];
-  Bool_t zdcPA = tdc[13];
-  Bool_t zdcPC = tdc[11];
-
 
   if (side == kASide) return ((useZP&&zdcPA) || (useZN&&zdcNA)); 
   if (side == kCSide) return ((useZP&&zdcPC) || (useZN&&zdcNC)); 
   return kFALSE;
 }
 
+Bool_t AliTriggerAnalysis::ZDCTimeTrigger(const AliESDEvent *aEsd, Bool_t fillHists) const
+{
+  // This method implements a selection
+  // based on the timing in both sides of zdcN
+  // It can be used in order to eliminate
+  // parasitic collisions
+  Bool_t zdcAccept = kFALSE;
+  AliESDZDC *esdZDC = aEsd->GetESDZDC();
+
+  if(fMC) {
+     UInt_t esdFlag =  esdZDC->GetESDQuality();
+
+     Bool_t znaFired=kFALSE, zpaFired=kFALSE;
+     Bool_t zem1Fired=kFALSE, zem2Fired=kFALSE;
+     Bool_t zncFired=kFALSE, zpcFired=kFALSE;
+     //
+     // **** Trigger patterns
+     if((esdFlag & 0x00000001) == 0x00000001) znaFired=kTRUE;
+     if((esdFlag & 0x00000002) == 0x00000002) zpaFired=kTRUE;
+     if((esdFlag & 0x00000004) == 0x00000004) zem1Fired=kTRUE;
+     if((esdFlag & 0x00000008) == 0x00000008) zem2Fired=kTRUE;
+     if((esdFlag & 0x00000010) == 0x00000010) zncFired=kTRUE;
+     if((esdFlag & 0x00000020) == 0x00000020) zpcFired=kTRUE;
+     zdcAccept = (znaFired | zncFired);
+  }
+  else {
+    for(Int_t i = 0; i < 4; ++i) {
+      if (esdZDC->GetZDCTDCData(10,i) != 0) {
+       Float_t tdcC = 0.025*(esdZDC->GetZDCTDCData(10,i)-esdZDC->GetZDCTDCData(14,i)); 
+       Float_t tdcCcorr = esdZDC->GetZDCTDCCorrected(10,i); 
+       for(Int_t j = 0; j < 4; ++j) {
+         if (esdZDC->GetZDCTDCData(12,j) != 0) {
+           Float_t tdcA = 0.025*(esdZDC->GetZDCTDCData(12,j)-esdZDC->GetZDCTDCData(14,j));
+
+           Float_t tdcAcorr = esdZDC->GetZDCTDCCorrected(12,j);
+           if(fillHists) {
+             fHistTimeZDC->Fill(tdcC-tdcA,tdcC+tdcA);
+             fHistTimeCorrZDC->Fill(tdcCcorr-tdcAcorr,tdcCcorr+tdcAcorr);
+           }
+           if (esdZDC->TestBit(AliESDZDC::kCorrectedTDCFilled)) {
+             if (((tdcCcorr-tdcAcorr-fZDCCutRefDeltaCorr)*(tdcCcorr-tdcAcorr-fZDCCutRefDeltaCorr)/(fZDCCutSigmaDeltaCorr*fZDCCutSigmaDeltaCorr) +
+                  (tdcCcorr+tdcAcorr-fZDCCutRefSumCorr)*(tdcCcorr+tdcAcorr-fZDCCutRefSumCorr)/(fZDCCutSigmaSumCorr*fZDCCutSigmaSumCorr))< 1.0)
+               zdcAccept = kTRUE;
+           }
+           else {
+             if (((tdcC-tdcA-fZDCCutRefDelta)*(tdcC-tdcA-fZDCCutRefDelta)/(fZDCCutSigmaDelta*fZDCCutSigmaDelta) +
+                  (tdcC+tdcA-fZDCCutRefSum)*(tdcC+tdcA-fZDCCutRefSum)/(fZDCCutSigmaSum*fZDCCutSigmaSum))< 1.0)
+               zdcAccept = kTRUE;
+           }
+         }
+       }
+      }
+    }
+  }
+  return zdcAccept;
+}
+
 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
 {
   // Returns if ZDC triggered
@@ -1062,7 +1340,7 @@ Long64_t AliTriggerAnalysis::Merge(TCollection* list)
   TObject* obj;
 
   // collections of all histograms
-  const Int_t nHists = 10;
+  const Int_t nHists = 12;
   TList collections[nHists];
 
   Int_t count = 0;
@@ -1077,6 +1355,8 @@ Long64_t AliTriggerAnalysis::Merge(TCollection* list)
     collections[n++].Add(entry->fHistV0C);
     collections[n++].Add(entry->fHistZDC);
     collections[n++].Add(entry->fHistTDCZDC);
+    collections[n++].Add(entry->fHistTimeZDC);
+    collections[n++].Add(entry->fHistTimeCorrZDC);
     collections[n++].Add(entry->fHistFMDA);
     collections[n++].Add(entry->fHistFMDC);
     collections[n++].Add(entry->fHistFMDSingle);
@@ -1089,7 +1369,7 @@ Long64_t AliTriggerAnalysis::Merge(TCollection* list)
     TObjString* obj2 = 0;
     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
     {
-      TParameter<Long64_t>* param2 = dynamic_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
+      TParameter<Long64_t>* param2 = static_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
       
       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
       if (param1)
@@ -1113,6 +1393,14 @@ Long64_t AliTriggerAnalysis::Merge(TCollection* list)
   fHistV0C->Merge(&collections[n++]);
   fHistZDC->Merge(&collections[n++]);
   fHistTDCZDC->Merge(&collections[n++]);
+  if (fHistTimeZDC)
+    fHistTimeZDC->Merge(&collections[n++]);
+  else
+    n++;
+  if (fHistTimeCorrZDC)
+    fHistTimeCorrZDC->Merge(&collections[n++]);
+  else
+    n++;
   fHistFMDA->Merge(&collections[n++]);
   fHistFMDC->Merge(&collections[n++]);
   fHistFMDSingle->Merge(&collections[n++]);
@@ -1132,21 +1420,36 @@ void AliTriggerAnalysis::SaveHistograms() const
   if (!fHistBitsSPD)
     return;
     
-  fHistBitsSPD->Write();
-  fHistBitsSPD->ProjectionX();
-  fHistBitsSPD->ProjectionY();
-  fHistFiredBitsSPD->Write();
-  fHistV0A->Write();
-  fHistV0C->Write();
-  fHistZDC->Write();
-  fHistTDCZDC->Write();
-  fHistFMDA->Write();
-  fHistFMDC->Write();
-  fHistFMDSingle->Write();
-  fHistFMDSum->Write();
-  
-  if (fSPDGFOEfficiency)
-    fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
+  if (fHistBitsSPD) {
+    fHistBitsSPD->Write();
+    fHistBitsSPD->ProjectionX();
+    fHistBitsSPD->ProjectionY();
+  }
+  else Printf("Cannot save fHistBitsSPD");
+  if (fHistFiredBitsSPD) fHistFiredBitsSPD->Write();
+  else Printf("Cannot save fHistFiredBitsSPD");
+  if (fHistV0A) fHistV0A->Write();
+  else Printf("Cannot save fHistV0A");
+  if (fHistV0C) fHistV0C->Write();
+  else Printf("Cannot save fHistV0C");
+  if (fHistZDC) fHistZDC->Write();
+  else Printf("Cannot save fHistZDC");
+  if (fHistTDCZDC) fHistTDCZDC->Write();
+  else Printf("Cannot save fHistTDCZDC");
+  if (fHistTimeZDC) fHistTimeZDC->Write();
+  else Printf("Cannot save fHistTimeZDC");
+  if (fHistTimeCorrZDC) fHistTimeCorrZDC->Write();
+  else Printf("Cannot save fHistTimeCorrZDC");
+  if (fHistFMDA) fHistFMDA->Write();
+  else Printf("Cannot save fHistFMDA");
+  if (fHistFMDC) fHistFMDC->Write();
+  else Printf("Cannot save fHistFMDC");
+  if (fHistFMDSingle) fHistFMDSingle->Write();
+  else Printf("Cannot save fHistFMDSingle");
+  if (fHistFMDSum) fHistFMDSum->Write();
+  else Printf("Cannot save fHistFMDSum");
+  if (fSPDGFOEfficiency) fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
+  //  else Printf("Cannot save fSPDGFOEfficiency");
   
   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
 }
@@ -1199,3 +1502,5 @@ void AliTriggerAnalysis::PrintTriggerClasses() const
   
   singleTrigger.DeleteAll();
 }
+
+