Просмотр исходного кода

Added AcceptanceFilter to imitate acceptance inefficiencies in TPC & FHCal

PeterParfenov лет назад: 3
Родитель
Сommit
d64acf72e7

+ 1 - 0
MpdFlow.LinkDef.h

@@ -30,6 +30,7 @@
 #pragma link C++ namespace mfFlow;
 
 #pragma link C++ class mfPicoDstReader+;
+#pragma link C++ class mfAcceptanceFilter+;
 #pragma link C++ class mfCentralityInfo+;
 #pragma link C++ class mfQvector+;
 #pragma link C++ class mfEventInfo+;

+ 22 - 2
bin/get-dca.cpp

@@ -1,5 +1,6 @@
 #include "mfPicoDstReader.h"
 #include "mfNamespace.h"
+#include "mfAcceptanceFilter.h"
 
 #include <TStopwatch.h>
 #include <TH1F.h>
@@ -25,16 +26,18 @@ using mfdca::GetPtBin;
 int main(int argc, char **argv)
 {
   TString iFileName, oFileName;
+  Bool_t isAcceptanceFilterOn = false;
 
   if (argc < 5)
   {
-    std::cerr << "./get-dca -i INPUTFILE/LIST -o OUTPUT.root" << std::endl;
+    std::cerr << "./get-dca -i INPUTFILE/LIST -o OUTPUT.root [OPTIONAL: --bad-acceptance]" << std::endl;
     return 1;
   }
   for (int i = 1; i < argc; i++)
   {
     if (std::string(argv[i]) != "-i" &&
-        std::string(argv[i]) != "-o")
+        std::string(argv[i]) != "-o" &&
+        std::string(argv[i]) != "--bad-acceptance")
     {
       std::cerr << "\n[ERROR]: Unknown parameter " << i << ": " << argv[i] << std::endl;
       return 2;
@@ -61,6 +64,12 @@ int main(int argc, char **argv)
         std::cerr << "\n[ERROR]: Output file name was not specified " << std::endl;
         return 1;
       }
+      if (std::string(argv[i]) == "--bad-acceptance")
+      {
+        std::cout << "Message: Acceptance filter is ON." << std::endl;
+        isAcceptanceFilterOn = true;
+        continue;
+      }
     }
   }
 
@@ -106,6 +115,10 @@ int main(int argc, char **argv)
   inChain->SetBranchAddress("recotracks", &recoTracks);
   // inChain->SetBranchAddress("FHCalModules", &fhcalmodules);
 
+  mfAcceptanceFilter *AccFilter = new mfAcceptanceFilter();
+  if (isAcceptanceFilterOn)
+    AccFilter->Activate();
+
   // Start event loop
   int n_entries = inChain->GetEntries();
   for (int iEv = 0; iEv < n_entries; iEv++)
@@ -119,6 +132,13 @@ int main(int argc, char **argv)
     for (int iTr = 0; iTr < reco_mult; iTr++)
     {
       auto recoTrack = (PicoDstRecoTrack *)recoTracks->UncheckedAt(iTr);
+      
+      if (!recoTrack)
+        continue;
+
+      if (!AccFilter->isGoodTpc(recoTrack->GetPhi()))
+        continue;
+
       Int_t pt_bin = GetPtBin(TMath::Abs(recoTrack->GetPt()));
       Int_t eta_bin = GetEtaBin(recoTrack->GetEta());
       if ((eta_bin == -1) || (pt_bin == -1))

+ 15 - 4
bin/get-flattening.cpp

@@ -3,6 +3,7 @@
 #include "mfProcessEvent.h"
 #include "mfInputReader.h"
 #include "mfRecentering.h"
+#include "mfAcceptanceFilter.h"
 
 #include <map>
 
@@ -38,10 +39,11 @@ using TMath::Sin;
 int main(int argc, char **argv)
 {
   TString iFileName, oFileName, dcaFileName, recenteringFileName;
+  Bool_t isAcceptanceFilterOn = false;
 
   if (argc < 7)
   {
-    std::cerr << "./get-flattening -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root -recentering RECENTERING.root" << std::endl;
+    std::cerr << "./get-flattening -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root -recentering RECENTERING.root [OPTIONAL: --bad-acceptance]" << std::endl;
     return 1;
   }
   for (int i = 1; i < argc; i++)
@@ -49,7 +51,8 @@ int main(int argc, char **argv)
     if (std::string(argv[i]) != "-i" &&
         std::string(argv[i]) != "-o" &&
         std::string(argv[i]) != "-dca" &&
-        std::string(argv[i]) != "-recentering")
+        std::string(argv[i]) != "-recentering" &&
+        std::string(argv[i]) != "--bad-acceptance")
     {
       std::cerr << "\n[ERROR]: Unknown parameter " << i << ": " << argv[i] << std::endl;
       return 2;
@@ -96,6 +99,12 @@ int main(int argc, char **argv)
         std::cerr << "\n[ERROR]: Input file name was not specified " << std::endl;
         return 1;
       }
+      if (std::string(argv[i]) == "--bad-acceptance")
+      {
+        std::cout << "Message: Acceptance filter is ON." << std::endl;
+        isAcceptanceFilterOn = true;
+        continue;
+      }
     }
   }
 
@@ -159,9 +168,11 @@ int main(int argc, char **argv)
     if (cent == -1)
       continue;
     procEvent->SetInputReader(inputReader);
-    // std::cout << "cent " << cent << std::endl;
 
-    // std::cout << "\tNtracks = " << recoTracks->GetEntriesFast() << std::endl;
+    mfAcceptanceFilter *AccFilter = new mfAcceptanceFilter();
+    if (isAcceptanceFilterOn)
+      AccFilter->Activate();
+    procEvent->SetAcceptaceFilter(AccFilter);
 
     mfEventInfo evInfoTpcRaw = procEvent->ProcessRecoTpcEP(recoTracks, 2.);
     if (evInfoTpcRaw.GetQvectors().size() == 0)

+ 23 - 6
bin/get-flow.cpp

@@ -5,6 +5,7 @@
 #include "mfSelector.h"
 #include "mfRecentering.h"
 #include "mfFlattening.h"
+#include "mfAcceptanceFilter.h"
 
 #include <TStopwatch.h>
 #include <TProfile2D.h>
@@ -45,11 +46,12 @@ using TMath::Sqrt;
 
 int main(int argc, char **argv)
 {
+  Bool_t isAcceptanceFilterOn = false;
   TString iFileName, oFileName, dcaFileName, resFileName, recenteringFileName = "", flatteningFileName = "";
 
   if (argc < 9)
   {
-    std::cerr << "./get-flow -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root -res RES_FIT.root [OPTIONAL: -recentering RECENTERING.root -flattening FLATTENING.root]" << std::endl;
+    std::cerr << "./get-flow -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root -res RES_FIT.root [OPTIONAL: -recentering RECENTERING.root -flattening FLATTENING.root --bad-acceptance]" << std::endl;
     return 1;
   }
   for (int i = 1; i < argc; i++)
@@ -59,7 +61,8 @@ int main(int argc, char **argv)
         std::string(argv[i]) != "-dca" &&
         std::string(argv[i]) != "-res" &&
         std::string(argv[i]) != "-recentering" &&
-        std::string(argv[i]) != "-flattening")
+        std::string(argv[i]) != "-flattening" &&
+        std::string(argv[i]) != "--bad-acceptance")
     {
       std::cerr << "\n[ERROR]: Unknown parameter " << i << ": " << argv[i] << std::endl;
       return 2;
@@ -126,6 +129,12 @@ int main(int argc, char **argv)
         std::cerr << "\n[ERROR]: Input file name was not specified " << std::endl;
         return 1;
       }
+      if (std::string(argv[i]) == "--bad-acceptance")
+      {
+        std::cout << "Message: Acceptance filter is ON." << std::endl;
+        isAcceptanceFilterOn = true;
+        continue;
+      }
     }
   }
 
@@ -287,6 +296,11 @@ int main(int argc, char **argv)
       continue;
     procEvent->SetInputReader(inputReader);
 
+    mfAcceptanceFilter *AccFilter = new mfAcceptanceFilter();
+    if (isAcceptanceFilterOn)
+      AccFilter->Activate();
+    procEvent->SetAcceptaceFilter(AccFilter);
+
     mfQvector qv;
 
     mfEventInfo evInfoTpcRaw = procEvent->ProcessRecoTpcEP(recoTracks, 2.);
@@ -362,6 +376,9 @@ int main(int argc, char **argv)
       if (!recoTrack)
         continue;
 
+      if (!AccFilter->isGoodTpc(recoTrack->GetPhi()))
+        continue;
+
       // Primary track selection
       if (!selector->isRecoPrimary(recoTrack->GetDCAx(), inputReader->GetDcaSigmaX(recoTrack)))
         continue;
@@ -441,7 +458,7 @@ int main(int argc, char **argv)
             pv2YRecoMcTpc[iGap][2]->Fill(RapidityMc, cent, v2TpcMc);
           if (selector->isPdgProton(mcTrack)) //if (selector->isProton(recoTrack) && recoTrack->GetCharge() > 0)
             pv2YRecoMcTpc[iGap][3]->Fill(RapidityMc, cent, v2TpcMc);
-          
+
           if (recoTrack->GetCharge() < 0)
             pv2PtTpc[iGap][0]->Fill(recoTrack->GetPt(), cent, v2Tpc);
           if (selector->isPion(recoTrack) && recoTrack->GetCharge() < 0)
@@ -510,7 +527,7 @@ int main(int argc, char **argv)
             pv2YRecoMcTpc[iGap][2]->Fill(RapidityMc, cent, v2TpcMc);
           if (selector->isPdgProton(mcTrack)) //if (selector->isProton(recoTrack) && recoTrack->GetCharge() > 0)
             pv2YRecoMcTpc[iGap][3]->Fill(RapidityMc, cent, v2TpcMc);
-          
+
           if (recoTrack->GetCharge() < 0)
             pv2PtTpc[iGap][4]->Fill(recoTrack->GetPt(), cent, v2Tpc);
           if (selector->isPion(recoTrack) && recoTrack->GetCharge() < 0)
@@ -824,7 +841,7 @@ int main(int argc, char **argv)
             pv1PtRp[2]->Fill(mcTrack->GetPt(), cent, v1Rp);
           if (selector->isPdgProton(mcTrack))
             pv1PtRp[3]->Fill(mcTrack->GetPt(), cent, v1Rp);
-          
+
           if (charge < 0.)
             pv1PtRp[4]->Fill(mcTrack->GetPt(), cent, v1Rp);
           if (selector->isPdgAntiPion(mcTrack))
@@ -864,7 +881,7 @@ int main(int argc, char **argv)
           pv2YRp[2]->Fill(Rapidity, cent, v2Rp);
         if (selector->isPdgProton(mcTrack))
           pv2YRp[3]->Fill(Rapidity, cent, v2Rp);
-        
+
         if (charge < 0.)
           pv2PtRp[4]->Fill(mcTrack->GetPt(), cent, v2Rp);
         if (selector->isPdgAntiPion(mcTrack))

+ 58 - 2
bin/get-qv.cpp

@@ -4,6 +4,7 @@
 #include "mfInputReader.h"
 #include "mfRecentering.h"
 #include "mfFlattening.h"
+#include "mfAcceptanceFilter.h"
 
 #include <vector>
 #include <map>
@@ -23,9 +24,11 @@
 #include <PicoDstFHCal.h>
 
 using mfdca::Ndim;
+using mfEp::BadFHCalModules;
 using mfEp::EpNames;
 using mfEp::NetaGaps;
 using mfEp::Nflattening;
+using mfEp::Nmodules;
 using mfEp::NsubEvents;
 using mfEp::NsubEventsTpc;
 using mfEp::RangeAngle;
@@ -40,11 +43,12 @@ using TMath::Sin;
 
 int main(int argc, char **argv)
 {
+  Bool_t isAcceptanceFilterOn = false;
   TString iFileName, oFileName, dcaFileName, recenteringFileName, flatteningFileName;
 
   if (argc < 11)
   {
-    std::cerr << "./get-qv -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root -recentering RECENTERING.root -flattening FLATTENING.root" << std::endl;
+    std::cerr << "./get-qv -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root -recentering RECENTERING.root -flattening FLATTENING.root [OPTIONAL: --bad-acceptance]" << std::endl;
     return 1;
   }
   for (int i = 1; i < argc; i++)
@@ -53,7 +57,8 @@ int main(int argc, char **argv)
         std::string(argv[i]) != "-o" &&
         std::string(argv[i]) != "-dca" &&
         std::string(argv[i]) != "-recentering" &&
-        std::string(argv[i]) != "-flattening")
+        std::string(argv[i]) != "-flattening" &&
+        std::string(argv[i]) != "--bad-acceptance")
     {
       std::cerr << "\n[ERROR]: Unknown parameter " << i << ": " << argv[i] << std::endl;
       return 2;
@@ -110,6 +115,12 @@ int main(int argc, char **argv)
         std::cerr << "\n[ERROR]: Input file name was not specified " << std::endl;
         return 1;
       }
+      if (std::string(argv[i]) == "--bad-acceptance")
+      {
+        std::cout << "Message: Acceptance filter is ON." << std::endl;
+        isAcceptanceFilterOn = true;
+        continue;
+      }
     }
   }
 
@@ -133,7 +144,14 @@ int main(int argc, char **argv)
   mfFlattening *flattening = new mfFlattening();
   flattening->ReadFlatteningFile(flatteningFileName);
 
+  // Imitate bad acceptance
+  mfAcceptanceFilter *AccFilter = new mfAcceptanceFilter();
+  if (isAcceptanceFilterOn)
+    AccFilter->Activate();
+
   // Init output profiles
+  TH2F *hTpcAcceptancePhi[2];
+  TH2F *hFHCalAcceptance[2];
   TH2F *hQxVsCent[3][NsubEvents];
   TH2F *hQyVsCent[3][NsubEvents];
   TH2F *hPsiEpVsCent[3][NsubEvents];
@@ -142,6 +160,12 @@ int main(int argc, char **argv)
   TH2F *hPsiEpRLVsCent[3][NetaGaps + 1];
   TH2F *hPsiRpRVsCent[3][NetaGaps + 1];
   TH2F *hPsiRpLVsCent[3][NetaGaps + 1];
+
+  for (int i = 0; i < 2; i++)
+  {
+    hTpcAcceptancePhi[i] = new TH2F(Form("hTpcAcceptancePhi_%i", i), Form("hTpcAcceptancePhi_%i", i), 360, -Pi(), Pi(), NresBin, 0., 100.);
+    hFHCalAcceptance[i] = new TH2F(Form("hFHCalAcceptance_%i", i), Form("hFHCalAcceptance_%i", i), 2 * Nmodules, 0, 2 * Nmodules, NresBin, 0., 100.);
+  }
   for (int iType = 0; iType < 3; iType++)
   {
     for (int i = 0; i < NsubEvents; i++)
@@ -203,6 +227,7 @@ int main(int argc, char **argv)
     if (cent == -1)
       continue;
     procEvent->SetInputReader(inputReader);
+    procEvent->SetAcceptaceFilter(AccFilter);
 
     // std::cout << "\tNtracks = " << recoTracks->GetEntriesFast() << std::endl;
 
@@ -232,6 +257,18 @@ int main(int argc, char **argv)
       psi[2].push_back(qv[2].at(iSub).GetPhi());
     }
 
+    for (int iTr = 0; iTr < recoTracks->GetEntriesFast(); iTr++)
+    {
+      auto recoTrack = (PicoDstRecoTrack *)recoTracks->UncheckedAt(iTr);
+      if (!recoTrack)
+        continue;
+
+      hTpcAcceptancePhi[0]->Fill(recoTrack->GetPhi(), cent);
+      if (!AccFilter->isGoodTpc(recoTrack->GetPhi()))
+        continue;
+      hTpcAcceptancePhi[1]->Fill(recoTrack->GetPhi(), cent);
+    }
+
     mfEventInfo evInfoFHCalRaw = procEvent->ProcessRecoFHCalEP(fhcalmodules, 1.);
     mfEventInfo evInfoFHCalRec = recentering->ProcessRecentering(evInfoFHCalRaw, cent);
     mfEventInfo evInfoFHCalFlat = flattening->ProcessFlattening(evInfoFHCalRec, cent);
@@ -253,6 +290,20 @@ int main(int argc, char **argv)
       psi[2].push_back(qv[2].at(iSub).GetPhi());
     }
 
+    for (int iMod = 0; iMod < fhcalmodules->GetEntriesFast(); iMod++)
+    {
+      auto fhcalModule = (PicoDstFHCal *)fhcalmodules->UncheckedAt(iMod);
+      if (!fhcalModule)
+        continue;
+      if ((iMod == BadFHCalModules.at(0)) || (iMod == BadFHCalModules.at(1)))
+        continue;
+
+      hFHCalAcceptance[0]->Fill(iMod, cent, fhcalModule->GetEnergy());
+      if (!AccFilter->isGoodFHCal(iMod))
+        continue;
+      hFHCalAcceptance[1]->Fill(iMod, cent, fhcalModule->GetEnergy());
+    }
+
     for (int i = 0; i < 3; i++)
     {
       for (int iGap = 0; iGap < NetaGaps; iGap++)
@@ -303,6 +354,11 @@ int main(int argc, char **argv)
   // Init output file
   TFile *fo = new TFile(oFileName.Data(), "recreate");
   fo->cd();
+  for (int i = 0; i < 2; i++)
+  {
+    hTpcAcceptancePhi[i]->Write();
+    hFHCalAcceptance[i]->Write();
+  }
   fo->mkdir("Raw");
   fo->mkdir("Recentering");
   fo->mkdir("Flattening");

+ 17 - 2
bin/get-recentering.cpp

@@ -2,6 +2,7 @@
 #include "mfNamespace.h"
 #include "mfProcessEvent.h"
 #include "mfInputReader.h"
+#include "mfAcceptanceFilter.h"
 
 #include <TStopwatch.h>
 #include <TProfile.h>
@@ -34,17 +35,19 @@ using TMath::Sin;
 int main(int argc, char **argv)
 {
   TString iFileName, oFileName, dcaFileName;
+  Bool_t isAcceptanceFilterOn = false;
 
   if (argc < 7)
   {
-    std::cerr << "./get-recentering -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root" << std::endl;
+    std::cerr << "./get-recentering -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root [OPTIONAL: --bad-acceptance]" << std::endl;
     return 1;
   }
   for (int i = 1; i < argc; i++)
   {
     if (std::string(argv[i]) != "-i" &&
         std::string(argv[i]) != "-o" &&
-        std::string(argv[i]) != "-dca")
+        std::string(argv[i]) != "-dca" &&
+        std::string(argv[i]) != "--bad-acceptance")
     {
       std::cerr << "\n[ERROR]: Unknown parameter " << i << ": " << argv[i] << std::endl;
       return 2;
@@ -81,6 +84,12 @@ int main(int argc, char **argv)
         std::cerr << "\n[ERROR]: Input file name was not specified " << std::endl;
         return 1;
       }
+      if (std::string(argv[i]) == "--bad-acceptance")
+      {
+        std::cout << "Message: Acceptance filter is ON." << std::endl;
+        isAcceptanceFilterOn = true;
+        continue;
+      }
     }
   }
 
@@ -138,6 +147,12 @@ int main(int argc, char **argv)
       continue;
     procEvent->SetInputReader(inputReader);
 
+    mfAcceptanceFilter *AccFilter = new mfAcceptanceFilter();
+    if (isAcceptanceFilterOn)
+      AccFilter->Activate();
+
+    procEvent->SetAcceptaceFilter(AccFilter);
+
     // std::cout << "\tNtracks = " << recoTracks->GetEntriesFast() << std::endl;
 
     mfEventInfo evInfoTpc = procEvent->ProcessRecoTpcEP(recoTracks, 2.);

+ 16 - 2
bin/get-res.cpp

@@ -5,6 +5,7 @@
 #include "mfSelector.h"
 #include "mfRecentering.h"
 #include "mfFlattening.h"
+#include "mfAcceptanceFilter.h"
 
 #include <TStopwatch.h>
 #include <TProfile.h>
@@ -36,11 +37,12 @@ using TMath::Sin;
 
 int main(int argc, char **argv)
 {
+  Bool_t isAcceptanceFilterOn = false;
   TString iFileName, oFileName, dcaFileName, recenteringFileName = "", flatteningFileName = "";
 
   if (argc < 7)
   {
-    std::cerr << "./get-res -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root [OPTIONAL: -recentering RECENTERING.root -flattening FLATTENING.root]" << std::endl;
+    std::cerr << "./get-res -i INPUTFILE/LIST -o OUTPUT.root -dca DCA_FIT.root [OPTIONAL: -recentering RECENTERING.root -flattening FLATTENING.root --bad-acceptance]" << std::endl;
     return 1;
   }
   for (int i = 1; i < argc; i++)
@@ -49,7 +51,8 @@ int main(int argc, char **argv)
         std::string(argv[i]) != "-o" &&
         std::string(argv[i]) != "-dca" &&
         std::string(argv[i]) != "-recentering" &&
-        std::string(argv[i]) != "-flattening")
+        std::string(argv[i]) != "-flattening" &&
+        std::string(argv[i]) != "--bad-acceptance")
     {
       std::cerr << "\n[ERROR]: Unknown parameter " << i << ": " << argv[i] << std::endl;
       return 2;
@@ -106,6 +109,12 @@ int main(int argc, char **argv)
         std::cerr << "\n[ERROR]: Input file name was not specified " << std::endl;
         return 1;
       }
+      if (std::string(argv[i]) == "--bad-acceptance")
+      {
+        std::cout << "Message: Acceptance filter is ON." << std::endl;
+        isAcceptanceFilterOn = true;
+        continue;
+      }
     }
   }
 
@@ -224,6 +233,11 @@ int main(int argc, char **argv)
     if (cent == -1)
       continue;
     procEvent->SetInputReader(inputReader);
+    
+    mfAcceptanceFilter *AccFilter = new mfAcceptanceFilter();
+    if (isAcceptanceFilterOn)
+      AccFilter->Activate();
+    procEvent->SetAcceptaceFilter(AccFilter);
 
     // std::cout << "\tNtracks = " << recoTracks->GetEntriesFast() << std::endl;
 

+ 55 - 0
corrections/mfAcceptanceFilter.cxx

@@ -0,0 +1,55 @@
+#include "mfAcceptanceFilter.h"
+#include "mfNamespace.h"
+
+using mfAcceptance::AcceptanceFHCalModules;
+using mfAcceptance::AcceptanceTpcPhiMax;
+using mfAcceptance::AcceptanceTpcPhiMin;
+using mfEp::RangeAngle_02Pi;
+
+ClassImp(mfAcceptanceFilter);
+
+mfAcceptanceFilter::mfAcceptanceFilter(/* args */) : isActive(false),
+                                                     fTpcPhiMin(AcceptanceTpcPhiMin),
+                                                     fTpcPhiMax(AcceptanceTpcPhiMax),
+                                                     fFHCalModules(AcceptanceFHCalModules)
+{
+}
+
+mfAcceptanceFilter::~mfAcceptanceFilter()
+{
+  fFHCalModules.clear();
+}
+
+void mfAcceptanceFilter::SetTpcPhi(double phi_min, double phi_max)
+{
+  fTpcPhiMin = phi_min;
+  fTpcPhiMax = phi_max;
+}
+
+Bool_t mfAcceptanceFilter::isGoodTpc(double phi)
+{
+  if (!isActive)
+    return true;
+
+  double phi_02Pi = RangeAngle_02Pi(phi, 1.);
+
+  if (phi_02Pi > fTpcPhiMin && phi_02Pi < fTpcPhiMax)
+    return false;
+  else
+    return true;
+}
+
+Bool_t mfAcceptanceFilter::isGoodFHCal(int mod)
+{
+  if (!isActive)
+    return true;
+
+  if (fFHCalModules.empty())
+    return true;
+  for (int i = 0; i < (int)fFHCalModules.size(); i++)
+  {
+    if (fFHCalModules.at(i) == mod)
+      return false;
+  }
+  return true;
+}

+ 35 - 0
corrections/mfAcceptanceFilter.h

@@ -0,0 +1,35 @@
+#ifndef MPD_FLOW_ACCEPTANCE_FILTER_H
+#define MPD_FLOW_ACCEPTANCE_FILTER_H
+
+#include <iostream>
+#include <vector>
+#include <Rtypes.h>
+
+class mfAcceptanceFilter
+{
+private:
+  Bool_t isActive;
+  double fTpcPhiMin;
+  double fTpcPhiMax;
+  std::vector<int> fFHCalModules;
+
+public:
+  mfAcceptanceFilter(/* args */);
+  virtual ~mfAcceptanceFilter();
+
+  // Activate acceptance filter
+  void Activate() { isActive = true; }
+  void Deactivate() { isActive = false; }
+
+  // Set bad acceptance
+  void SetTpcPhi(double phi_min, double phi_max);
+  void SetFHCal(std::vector<int> _a) { fFHCalModules = _a; }
+  void AddFHCalModule(int _a) { fFHCalModules.push_back(_a); }
+
+  Bool_t isGoodTpc(double phi);
+  Bool_t isGoodFHCal(int mod);
+
+  ClassDef(mfAcceptanceFilter, 0);
+};
+
+#endif

+ 11 - 0
eventanalysis/mfProcessEvent.cxx

@@ -57,6 +57,11 @@ mfEventInfo mfProcessEvent::ProcessRecoTpcEP(TClonesArray *const &recoTracks, Do
   if (!fInputReader)
     return fEventInfo;
 
+  if (!fAccFilter)
+  {
+    fAccFilter = new mfAcceptanceFilter();
+  }
+
   for (int i = 0; i < NsubEvents; i++)
   {
     Qx.push_back(0.);
@@ -72,6 +77,9 @@ mfEventInfo mfProcessEvent::ProcessRecoTpcEP(TClonesArray *const &recoTracks, Do
     if (!recoTrack)
       continue;
 
+    if (!fAccFilter->isGoodTpc(recoTrack->GetPhi()))
+      continue;
+
     // Get refMult
     if (selector->isRefMult(recoTrack))
     {
@@ -165,6 +173,9 @@ mfEventInfo mfProcessEvent::ProcessRecoFHCalEP(TClonesArray *const &fhcalModules
       continue;
     if ((iMod == BadFHCalModules.at(0)) || (iMod == BadFHCalModules.at(1)))
       continue;
+    if (!fAccFilter->isGoodFHCal(iMod))
+      continue;
+
     double phi = GetFHCalPhi(iMod);
 
     // std::cout << "\tModule " << iMod

+ 3 - 0
eventanalysis/mfProcessEvent.h

@@ -14,16 +14,19 @@
 #include "mfEventInfo.h"
 #include "mfNamespace.h"
 #include "mfInputReader.h"
+#include "mfAcceptanceFilter.h"
 
 class mfProcessEvent
 {
 private:
   mfInputReader *fInputReader;
+  mfAcceptanceFilter *fAccFilter;
 public:
   mfProcessEvent(/* args */);
   virtual ~mfProcessEvent();
 
   void SetInputReader(mfInputReader *const& _a) { fInputReader = _a; }
+  void SetAcceptaceFilter(mfAcceptanceFilter *const& _filter) { fAccFilter = _filter; }
 
   Float_t GetCentB(PicoDstMCEvent *const& mcEvent);
   mfEventInfo ProcessRecoTpcEP(TClonesArray *const& recoTracks, Double_t harmonic = 2.);

+ 10 - 0
utils/mfNamespace.h

@@ -92,4 +92,14 @@ namespace mfFlow
   const float MinV1PtEta = 0.2;
 }; // namespace mfFlow
 
+namespace mfAcceptance
+{
+  // Bad TPC sector (in 0 to 2Pi range)
+  const float AcceptanceTpcPhiMin = 0.26179939; // 15 deg
+  const float AcceptanceTpcPhiMax = 0.78539816; // 45 deg
+
+  // Bad FHCal modules (not 22, 67)
+  const std::vector<int> AcceptanceFHCalModules = {13, 50};
+}
+
 #endif