CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
CSCSkim Class Reference

#include <RecoLocalMuon/CSCSkim/src/CSCSkim.cc>

Inheritance diagram for CSCSkim:
edm::one::EDFilter<> edm::one::EDFilterBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void beginJob () override
 
 CSCSkim (const edm::ParameterSet &pset)
 
void endJob () override
 
bool filter (edm::Event &event, const edm::EventSetup &eventSetup) override
 
 ~CSCSkim () override
 
- Public Member Functions inherited from edm::one::EDFilter<>
 EDFilter ()=default
 
 EDFilter (const EDFilter &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDFilteroperator= (const EDFilter &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
- Public Member Functions inherited from edm::one::EDFilterBase
 EDFilterBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDFilterBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
std::vector< bool > const & recordProvenanceList () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
TypeLabelList const & typeLabelList () const
 used by the fwk to register the list of products of this module More...
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

int chamberSerial (int kE, int kS, int kR, int kCh)
 
bool doBFieldStudySelection (edm::Handle< reco::TrackCollection > saTracks, edm::Handle< reco::TrackCollection > Tracks, edm::Handle< reco::MuonCollection > gMuons)
 
bool doCertainChamberSelection (edm::Handle< CSCWireDigiCollection > wires, edm::Handle< CSCStripDigiCollection > strips)
 
bool doCSCSkimming (edm::Handle< CSCRecHit2DCollection > cscRecHits, edm::Handle< CSCSegmentCollection > cscSegments)
 
bool doDTOverlap (edm::Handle< CSCSegmentCollection > cscSegments)
 
bool doHaloLike (edm::Handle< CSCSegmentCollection > cscSegments)
 
bool doLongSATrack (edm::Handle< reco::TrackCollection > saTracks)
 
bool doMessyEventSkimming (edm::Handle< CSCRecHit2DCollection > cscRecHits, edm::Handle< CSCSegmentCollection > cscSegments)
 
bool doOverlapSkimming (edm::Handle< CSCSegmentCollection > cscSegments)
 

Private Attributes

bool demandChambersBothSides
 
edm::EDGetTokenT< reco::MuonCollectionglm_token
 
std::string histogramFileName
 
TH1F * hxnHitChambers
 
TH1F * hxnRecHits
 
TH1F * hxnRecHitsSel
 
TH1F * hxnSegments
 
int iEvent
 
int iRun
 
bool isSimulation
 
const edm::ESGetToken< CSCGeometry, MuonGeometryRecordm_CSCGeomToken
 
bool makeHistograms
 
bool makeHistogramsForMessyEvents
 
TH1F * mevnChambers0
 
TH1F * mevnChambers1
 
TH1F * mevnRecHits0
 
TH1F * mevnRecHits1
 
TH1F * mevnSegments0
 
TH1F * mevnSegments1
 
int minimumHitChambers
 
int minimumSegments
 
int nCSCHitsMin
 
int nEventsAnalyzed
 
int nEventsCertainChamber
 
int nEventsChambersBothSides
 
int nEventsDTOverlap
 
int nEventsForBFieldStudies
 
int nEventsHaloLike
 
int nEventsLongSATrack
 
int nEventsMessy
 
int nEventsOverlappingChambers
 
int nEventsSelected
 
int nLayersWithHitsMinimum
 
int nTrHitsMin
 
int nValidHitsMin
 
std::string outputFileName
 
float pMin
 
float redChiSqMax
 
float rExtMax
 
edm::EDGetTokenT< CSCRecHit2DCollectionrh_token
 
edm::EDGetTokenT< reco::TrackCollectionsam_token
 
edm::EDGetTokenT< CSCStripDigiCollectionsdr_token
 
edm::EDGetTokenT< CSCStripDigiCollectionsds_token
 
edm::EDGetTokenT< CSCSegmentCollectionseg_token
 
TFile * theHistogramFile
 
edm::EDGetTokenT< reco::TrackCollectiontrk_token
 
int typeOfSkim
 
edm::EDGetTokenT< CSCWireDigiCollectionwdr_token
 
edm::EDGetTokenT< CSCWireDigiCollectionwds_token
 
int whichChamber
 
int whichEndcap
 
int whichRing
 
int whichStation
 
TH1F * xxnCSCHits
 
TH1F * xxnTrackerHits
 
TH1F * xxnValidHits
 
TH1F * xxP
 
TH1F * xxredChiSq
 
float zInnerMax
 
float zLengthMin
 
float zLengthTrMin
 

Additional Inherited Members

- Public Types inherited from edm::one::EDFilterBase
typedef EDFilterBase ModuleType
 
- Public Types inherited from edm::ProducerBase
template<typename T >
using BranchAliasSetterT = ProductRegistryHelper::BranchAliasSetterT< T >
 
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex > >
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDFilterBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
template<Transition Tr = Transition::Event>
auto produces (std::string instanceName) noexcept
 declare what type of product will make and with which optional label More...
 
template<Transition B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<BranchType B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces ()
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces ()
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces ()
 
template<Transition Tr = Transition::Event>
auto produces () noexcept
 
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

This simple program selects minimal CSC events for output.

Michael Schmitt, Northwestern University, July 2008

Description: Offline skim module for CSC cosmic ray data

Implementation: <Notes on="" implementation>="">

Definition at line 71 of file CSCSkim.h.

Constructor & Destructor Documentation

◆ CSCSkim()

CSCSkim::CSCSkim ( const edm::ParameterSet pset)
explicit

Definition at line 51 of file CSCSkim.cc.

References demandChambersBothSides, glm_token, histogramFileName, ProducerED_cfi::InputTag, makeHistograms, makeHistogramsForMessyEvents, minimumHitChambers, minimumSegments, nCSCHitsMin, nLayersWithHitsMinimum, nTrHitsMin, nValidHitsMin, outputFileName, pMin, muonDTDigis_cfi::pset, redChiSqMax, rExtMax, rh_token, sam_token, sdr_token, sds_token, seg_token, AlCaHLTBitMon_QueryRunRegistry::string, trk_token, typeOfSkim, wdr_token, wds_token, whichChamber, whichEndcap, whichRing, whichStation, zInnerMax, zLengthMin, and zLengthTrMin.

52  // tokens from tags
53 
54  // Really should define the wire and digi tags in config, but for now, to avoid having to modify
55  // multiple config files, just hard-code those tags, to be equivalent to the pre-consumes code
56 
57  // wds_token = consumes<CSCWireDigiCollection>(pset.getParameter<InputTag>("simWireDigiTag"));
58  // sds_token = consumes<CSCStripDigiCollection>(pset.getParameter<InputTag>("simStripDigiTag"));
59  // wdr_token = consumes<CSCWireDigiCollection>(pset.getParameter<InputTag>("wireDigiTag"));
60  // sdr_token = consumes<CSCStripDigiCollection>(pset.getParameter<InputTag>("stripDigiTag"));
61 
62  wds_token = consumes<CSCWireDigiCollection>(edm::InputTag("simMuonCSCDigis", "MuonCSCWireDigi"));
63  sds_token = consumes<CSCStripDigiCollection>(edm::InputTag("simMuonCSCDigis", "MuonCSCStripDigi"));
64  wdr_token = consumes<CSCWireDigiCollection>(edm::InputTag("muonCSCDigis", "MuonCSCWireDigi"));
65  sdr_token = consumes<CSCStripDigiCollection>(edm::InputTag("muonCSCDigis", "MuonCSCStripDigi"));
66 
67  rh_token = consumes<CSCRecHit2DCollection>(pset.getParameter<InputTag>("cscRecHitTag"));
68  seg_token = consumes<CSCSegmentCollection>(pset.getParameter<InputTag>("cscSegmentTag"));
69 
70  sam_token = consumes<reco::TrackCollection>(pset.getParameter<InputTag>("SAMuonTag"));
71  trk_token = consumes<reco::TrackCollection>(pset.getParameter<InputTag>("trackTag"));
72  glm_token = consumes<reco::MuonCollection>(pset.getParameter<InputTag>("GLBMuonTag"));
73 
74  // Get the various input parameters
75  outputFileName = pset.getUntrackedParameter<std::string>("outputFileName", "outputSkim.root");
76  histogramFileName = pset.getUntrackedParameter<std::string>("histogramFileName", "histos.root");
77  typeOfSkim = pset.getUntrackedParameter<int>("typeOfSkim", 1);
78  nLayersWithHitsMinimum = pset.getUntrackedParameter<int>("nLayersWithHitsMinimum", 3);
79  minimumHitChambers = pset.getUntrackedParameter<int>("minimumHitChambers", 1);
80  minimumSegments = pset.getUntrackedParameter<int>("minimumSegments", 3);
81  demandChambersBothSides = pset.getUntrackedParameter<bool>("demandChambersBothSides", false);
82  makeHistograms = pset.getUntrackedParameter<bool>("makeHistograms", false);
83  makeHistogramsForMessyEvents = pset.getUntrackedParameter<bool>("makeHistogramsForMessyEvebts", false);
84  whichEndcap = pset.getUntrackedParameter<int>("whichEndcap", 2);
85  whichStation = pset.getUntrackedParameter<int>("whichStation", 3);
86  whichRing = pset.getUntrackedParameter<int>("whichRing", 2);
87  whichChamber = pset.getUntrackedParameter<int>("whichChamber", 24);
88 
89  // for BStudy selection (skim type 9)
90  pMin = pset.getUntrackedParameter<double>("pMin", 3.);
91  zLengthMin = pset.getUntrackedParameter<double>("zLengthMin", 200.);
92  nCSCHitsMin = pset.getUntrackedParameter<int>("nCSCHitsMin", 9);
93  zInnerMax = pset.getUntrackedParameter<double>("zInnerMax", 9000.);
94  nTrHitsMin = pset.getUntrackedParameter<int>("nTrHitsMin", 8);
95  zLengthTrMin = pset.getUntrackedParameter<double>("zLengthTrMin", 180.);
96  rExtMax = pset.getUntrackedParameter<double>("rExtMax", 3000.);
97  redChiSqMax = pset.getUntrackedParameter<double>("redChiSqMax", 20.);
98  nValidHitsMin = pset.getUntrackedParameter<int>("nValidHitsMin", 8);
99 
100  LogInfo("[CSCSkim] Setup") << "\n\t===== CSCSkim =====\n"
101  << "\t\ttype of skim ...............................\t" << typeOfSkim
102  << "\t\tminimum number of layers with hits .........\t" << nLayersWithHitsMinimum
103  << "\n\t\tminimum number of chambers w/ hit layers..\t" << minimumHitChambers
104  << "\n\t\tminimum number of segments ...............\t" << minimumSegments
105  << "\n\t\tdemand chambers on both sides.............\t" << demandChambersBothSides
106  << "\n\t\tmake histograms...........................\t" << makeHistograms
107  << "\n\t\t..for messy events........................\t" << makeHistogramsForMessyEvents
108  << "\n\t===================\n\n";
109 }
edm::EDGetTokenT< CSCSegmentCollection > seg_token
Definition: CSCSkim.h:149
edm::EDGetTokenT< reco::TrackCollection > sam_token
Definition: CSCSkim.h:150
float zInnerMax
Definition: CSCSkim.h:171
float rExtMax
Definition: CSCSkim.h:174
bool makeHistogramsForMessyEvents
Definition: CSCSkim.h:162
int whichEndcap
Definition: CSCSkim.h:163
int nCSCHitsMin
Definition: CSCSkim.h:170
int nLayersWithHitsMinimum
Definition: CSCSkim.h:157
float zLengthMin
Definition: CSCSkim.h:169
edm::EDGetTokenT< CSCWireDigiCollection > wds_token
Definition: CSCSkim.h:143
std::string outputFileName
Definition: CSCSkim.h:136
float pMin
Definition: CSCSkim.h:168
int minimumHitChambers
Definition: CSCSkim.h:158
int typeOfSkim
Definition: CSCSkim.h:156
edm::EDGetTokenT< reco::TrackCollection > trk_token
Definition: CSCSkim.h:151
float zLengthTrMin
Definition: CSCSkim.h:173
int whichChamber
Definition: CSCSkim.h:166
int minimumSegments
Definition: CSCSkim.h:159
edm::EDGetTokenT< reco::MuonCollection > glm_token
Definition: CSCSkim.h:152
std::string histogramFileName
Definition: CSCSkim.h:137
edm::EDGetTokenT< CSCStripDigiCollection > sds_token
Definition: CSCSkim.h:144
int nValidHitsMin
Definition: CSCSkim.h:176
edm::EDGetTokenT< CSCRecHit2DCollection > rh_token
Definition: CSCSkim.h:148
bool makeHistograms
Definition: CSCSkim.h:161
int nTrHitsMin
Definition: CSCSkim.h:172
Log< level::Info, false > LogInfo
int whichStation
Definition: CSCSkim.h:164
int whichRing
Definition: CSCSkim.h:165
const edm::ESGetToken< CSCGeometry, MuonGeometryRecord > m_CSCGeomToken
Definition: CSCSkim.h:140
edm::EDGetTokenT< CSCWireDigiCollection > wdr_token
Definition: CSCSkim.h:145
edm::EDGetTokenT< CSCStripDigiCollection > sdr_token
Definition: CSCSkim.h:146
bool demandChambersBothSides
Definition: CSCSkim.h:160
float redChiSqMax
Definition: CSCSkim.h:175

◆ ~CSCSkim()

CSCSkim::~CSCSkim ( )
override

Definition at line 114 of file CSCSkim.cc.

114 {}

Member Function Documentation

◆ beginJob()

void CSCSkim::beginJob ( )
overridevirtual

Reimplemented from edm::one::EDFilterBase.

Definition at line 119 of file CSCSkim.cc.

References histogramFileName, hxnHitChambers, hxnRecHits, hxnRecHitsSel, hxnSegments, iEvent, iRun, makeHistograms, makeHistogramsForMessyEvents, mevnChambers0, mevnChambers1, mevnRecHits0, mevnRecHits1, mevnSegments0, mevnSegments1, nEventsAnalyzed, nEventsCertainChamber, nEventsChambersBothSides, nEventsDTOverlap, nEventsForBFieldStudies, nEventsHaloLike, nEventsLongSATrack, nEventsMessy, nEventsOverlappingChambers, nEventsSelected, theHistogramFile, xxnCSCHits, xxnTrackerHits, xxnValidHits, xxP, and xxredChiSq.

119  {
120  // set counters to zero
121  nEventsAnalyzed = 0;
122  nEventsSelected = 0;
125  nEventsMessy = 0;
127  nEventsDTOverlap = 0;
128  nEventsHaloLike = 0;
129  nEventsLongSATrack = 0;
131  iRun = 0;
132  iEvent = 0;
133 
135  // Create the root file for the histograms
136  theHistogramFile = new TFile(histogramFileName.c_str(), "RECREATE");
137  theHistogramFile->cd();
138 
139  if (makeHistograms) {
140  // book histograms for the skimming module
141  hxnRecHits = new TH1F("hxnRecHits", "n RecHits", 61, -0.5, 60.5);
142  hxnSegments = new TH1F("hxnSegments", "n Segments", 11, -0.5, 10.5);
143  hxnHitChambers = new TH1F("hxnHitsChambers", "n chambers with hits", 11, -0.5, 10.5);
144  hxnRecHitsSel = new TH1F("hxnRecHitsSel", "n RecHits selected", 61, -0.5, 60.5);
145 
146  xxP = new TH1F("xxP", "P global", 100, 0., 200.);
147  xxnValidHits = new TH1F("xxnValidHits", "n valid hits global", 61, -0.5, 60.5);
148  xxnTrackerHits = new TH1F("xxnTrackerHits", "n tracker hits global", 61, -0.5, 60.5);
149  xxnCSCHits = new TH1F("xxnCSCHits", "n CSC hits global", 41, -0.5, 40.5);
150  xxredChiSq = new TH1F("xxredChiSq", "red chisq global", 100, 0., 100.);
151  }
153  // book histograms for the messy event skimming module
154  mevnRecHits0 = new TH1F("mevnRecHits0", "n RecHits", 121, -0.5, 120.5);
155  mevnChambers0 = new TH1F("mevnChambers0", "n chambers with hits", 21, -0.5, 20.5);
156  mevnSegments0 = new TH1F("mevnSegments0", "n Segments", 21, -0.5, 20.5);
157  mevnRecHits1 = new TH1F("mevnRecHits1", "n RecHits", 100, 0., 300.);
158  mevnChambers1 = new TH1F("mevnChambers1", "n chambers with hits", 50, 0., 50.);
159  mevnSegments1 = new TH1F("mevnSegments1", "n Segments", 30, 0., 30.);
160  }
161  }
162 }
TH1F * hxnRecHitsSel
Definition: CSCSkim.h:182
bool makeHistogramsForMessyEvents
Definition: CSCSkim.h:162
int nEventsForBFieldStudies
Definition: CSCSkim.h:126
int nEventsChambersBothSides
Definition: CSCSkim.h:119
TH1F * mevnRecHits0
Definition: CSCSkim.h:183
TH1F * mevnRecHits1
Definition: CSCSkim.h:186
TH1F * mevnSegments0
Definition: CSCSkim.h:185
TH1F * mevnSegments1
Definition: CSCSkim.h:188
TH1F * hxnRecHits
Definition: CSCSkim.h:179
TH1F * xxnCSCHits
Definition: CSCSkim.h:190
int nEventsDTOverlap
Definition: CSCSkim.h:123
TH1F * xxnTrackerHits
Definition: CSCSkim.h:190
int nEventsOverlappingChambers
Definition: CSCSkim.h:120
int nEventsSelected
Definition: CSCSkim.h:118
TH1F * hxnSegments
Definition: CSCSkim.h:180
std::string histogramFileName
Definition: CSCSkim.h:137
bool makeHistograms
Definition: CSCSkim.h:161
int nEventsLongSATrack
Definition: CSCSkim.h:125
int nEventsAnalyzed
Definition: CSCSkim.h:117
int nEventsCertainChamber
Definition: CSCSkim.h:122
TH1F * mevnChambers1
Definition: CSCSkim.h:187
TFile * theHistogramFile
Definition: CSCSkim.h:133
int nEventsMessy
Definition: CSCSkim.h:121
TH1F * hxnHitChambers
Definition: CSCSkim.h:181
int iEvent
Definition: CSCSkim.h:130
int nEventsHaloLike
Definition: CSCSkim.h:124
int iRun
Definition: CSCSkim.h:129
TH1F * xxnValidHits
Definition: CSCSkim.h:190
TH1F * mevnChambers0
Definition: CSCSkim.h:184
TH1F * xxredChiSq
Definition: CSCSkim.h:190
TH1F * xxP
Definition: CSCSkim.h:190

◆ chamberSerial()

int CSCSkim::chamberSerial ( int  kE,
int  kS,
int  kR,
int  kCh 
)
private

Definition at line 1280 of file CSCSkim.cc.

Referenced by doCSCSkimming(), doMessyEventSkimming(), and doOverlapSkimming().

1280  {
1281  int kSerial = kChamber;
1282  if (kStation == 1 && kRing == 1) {
1283  kSerial = kChamber;
1284  }
1285  if (kStation == 1 && kRing == 2) {
1286  kSerial = kChamber + 36;
1287  }
1288  if (kStation == 1 && kRing == 3) {
1289  kSerial = kChamber + 72;
1290  }
1291  if (kStation == 1 && kRing == 4) {
1292  kSerial = kChamber;
1293  }
1294  if (kStation == 2 && kRing == 1) {
1295  kSerial = kChamber + 108;
1296  }
1297  if (kStation == 2 && kRing == 2) {
1298  kSerial = kChamber + 126;
1299  }
1300  if (kStation == 3 && kRing == 1) {
1301  kSerial = kChamber + 162;
1302  }
1303  if (kStation == 3 && kRing == 2) {
1304  kSerial = kChamber + 180;
1305  }
1306  if (kStation == 4 && kRing == 1) {
1307  kSerial = kChamber + 216;
1308  }
1309  if (kStation == 4 && kRing == 2) {
1310  kSerial = kChamber + 234;
1311  } // one day...
1312  if (kEndcap == 2) {
1313  kSerial = kSerial + 300;
1314  }
1315  return kSerial;
1316 }

◆ doBFieldStudySelection()

bool CSCSkim::doBFieldStudySelection ( edm::Handle< reco::TrackCollection saTracks,
edm::Handle< reco::TrackCollection Tracks,
edm::Handle< reco::MuonCollection gMuons 
)
private

Definition at line 1104 of file CSCSkim.cc.

References funct::abs(), MuonSubdetId::CSC, l1tTrackerHTMiss_cfi::deltaZ, hcalRecHitTable_cff::detId, spr::goodTrack(), trackingPlots::hp, DetId::Muon, HLT_2024v14_cff::muon, dqmiodumpmetadata::n, nCSCHitsMin, nTrHitsMin, nValidHitsMin, pMin, redChiSqMax, rExtMax, mathSSE::sqrt(), HLT_2024v14_cff::track, DiMuonV_cfg::tracks, zInnerMax, zLengthMin, and zLengthTrMin.

Referenced by filter().

1106  {
1107  bool acceptThisEvent = false;
1108 
1109  //-----------------------------------
1110  // examine the stand-alone tracks
1111  //-----------------------------------
1112  int nGoodSAMuons = 0;
1113  for (reco::TrackCollection::const_iterator muon = saMuons->begin(); muon != saMuons->end(); ++muon) {
1114  float preco = muon->p();
1115 
1116  math::XYZPoint innerPo = muon->innerPosition();
1117  GlobalPoint iPnt(innerPo.x(), innerPo.y(), innerPo.z());
1118  math::XYZPoint outerPo = muon->outerPosition();
1119  GlobalPoint oPnt(outerPo.x(), outerPo.y(), outerPo.z());
1120  float zLength = abs(iPnt.z() - oPnt.z());
1121 
1122  math::XYZVector innerMom = muon->innerMomentum();
1123  GlobalVector iP(innerMom.x(), innerMom.y(), innerMom.z());
1124  math::XYZVector outerMom = muon->outerMomentum();
1125  GlobalVector oP(outerMom.x(), outerMom.y(), outerMom.z());
1126 
1127  const float zRef = 300.;
1128  float xExt = 10000.;
1129  float yExt = 10000.;
1130  if (abs(oPnt.z()) < abs(iPnt.z())) {
1131  float deltaZ = 0.;
1132  if (oPnt.z() > 0) {
1133  deltaZ = zRef - oPnt.z();
1134  } else {
1135  deltaZ = -zRef - oPnt.z();
1136  }
1137  xExt = oPnt.x() + deltaZ * oP.x() / oP.z();
1138  yExt = oPnt.y() + deltaZ * oP.y() / oP.z();
1139  } else {
1140  float deltaZ = 0.;
1141  if (iPnt.z() > 0) {
1142  deltaZ = zRef - iPnt.z();
1143  } else {
1144  deltaZ = -zRef - iPnt.z();
1145  }
1146  xExt = iPnt.x() + deltaZ * iP.x() / iP.z();
1147  yExt = iPnt.y() + deltaZ * iP.y() / iP.z();
1148  }
1149  float rExt = sqrt(xExt * xExt + yExt * yExt);
1150 
1151  int nCSCHits = 0;
1152  for (trackingRecHit_iterator hit = muon->recHitsBegin(); hit != muon->recHitsEnd(); ++hit) {
1153  const DetId detId((*hit)->geographicalId());
1154  if (detId.det() == DetId::Muon) {
1155  if (detId.subdetId() == MuonSubdetId::CSC) {
1156  nCSCHits++;
1157  }
1158  }
1159  } // end loop over hits
1160 
1161  float zInner = -1.;
1162  if (nCSCHits >= nCSCHitsMin) {
1163  if (abs(iPnt.z()) < abs(iPnt.z())) {
1164  zInner = iPnt.z();
1165  } else {
1166  zInner = oPnt.z();
1167  }
1168  }
1169 
1170  bool goodSAMuon = (preco > pMin) && (zLength > zLengthMin) && (nCSCHits >= nCSCHitsMin) && (zInner < zInnerMax) &&
1171  (rExt < rExtMax);
1172 
1173  if (goodSAMuon) {
1174  nGoodSAMuons++;
1175  }
1176 
1177  } // end loop over stand-alone muon collection
1178 
1179  //-----------------------------------
1180  // examine the tracker tracks
1181  //-----------------------------------
1182  int nGoodTracks = 0;
1183  for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
1184  float preco = track->p();
1185  int n = track->recHitsSize();
1186 
1187  math::XYZPoint innerPo = track->innerPosition();
1188  GlobalPoint iPnt(innerPo.x(), innerPo.y(), innerPo.z());
1189  math::XYZPoint outerPo = track->outerPosition();
1190  GlobalPoint oPnt(outerPo.x(), outerPo.y(), outerPo.z());
1191  float zLength = abs(iPnt.z() - oPnt.z());
1192 
1193  math::XYZVector innerMom = track->innerMomentum();
1194  GlobalVector iP(innerMom.x(), innerMom.y(), innerMom.z());
1195  math::XYZVector outerMom = track->outerMomentum();
1196  GlobalVector oP(outerMom.x(), outerMom.y(), outerMom.z());
1197 
1198  const float zRef = 300.;
1199  float xExt = 10000.;
1200  float yExt = 10000.;
1201  if (abs(oPnt.z()) > abs(iPnt.z())) {
1202  float deltaZ = 0.;
1203  if (oPnt.z() > 0) {
1204  deltaZ = zRef - oPnt.z();
1205  } else {
1206  deltaZ = -zRef - oPnt.z();
1207  }
1208  xExt = oPnt.x() + deltaZ * oP.x() / oP.z();
1209  yExt = oPnt.y() + deltaZ * oP.y() / oP.z();
1210  } else {
1211  float deltaZ = 0.;
1212  if (iPnt.z() > 0) {
1213  deltaZ = zRef - iPnt.z();
1214  } else {
1215  deltaZ = -zRef - iPnt.z();
1216  }
1217  xExt = iPnt.x() + deltaZ * iP.x() / iP.z();
1218  yExt = iPnt.y() + deltaZ * iP.y() / iP.z();
1219  }
1220  float rExt = sqrt(xExt * xExt + yExt * yExt);
1221 
1222  bool goodTrack = (preco > pMin) && (n >= nTrHitsMin) && (zLength > zLengthTrMin) && (rExt < rExtMax);
1223 
1224  if (goodTrack) {
1225  nGoodTracks++;
1226  }
1227 
1228  } // end loop over tracker tracks
1229 
1230  //-----------------------------------
1231  // examine the global muons
1232  //-----------------------------------
1233  int nGoodGlobalMuons = 0;
1234  for (reco::MuonCollection::const_iterator global = gMuons->begin(); global != gMuons->end(); ++global) {
1235  if (global->isGlobalMuon()) {
1236  float pDef = global->p();
1237  float redChiSq = global->globalTrack()->normalizedChi2();
1238  const reco::HitPattern& hp = (global->globalTrack())->hitPattern();
1239  // int nTotalHits = hp.numberOfHits();
1240  // int nValidHits = hp.numberOfValidHits();
1241  int nTrackerHits = hp.numberOfValidTrackerHits();
1242  // int nPixelHits = hp.numberOfValidPixelHits();
1243  // int nStripHits = hp.numberOfValidStripHits();
1244 
1245  int nCSCHits = 0;
1246  for (trackingRecHit_iterator hit = (global->globalTrack())->recHitsBegin();
1247  hit != (global->globalTrack())->recHitsEnd();
1248  ++hit) {
1249  const DetId detId((*hit)->geographicalId());
1250  if (detId.det() == DetId::Muon) {
1251  if (detId.subdetId() == MuonSubdetId::CSC) {
1252  nCSCHits++;
1253  }
1254  }
1255  } // end loop over hits
1256 
1257  bool goodGlobalMuon =
1258  (pDef > pMin) && (nTrackerHits >= nValidHitsMin) && (nCSCHits >= nCSCHitsMin) && (redChiSq < redChiSqMax);
1259 
1260  if (goodGlobalMuon) {
1261  nGoodGlobalMuons++;
1262  }
1263 
1264  } // this is a global muon
1265  } // end loop over stand-alone muon collection
1266 
1267  //-----------------------------------
1268  // do we accept this event?
1269  //-----------------------------------
1270 
1271  acceptThisEvent = ((nGoodSAMuons > 0) && (nGoodTracks > 0)) || (nGoodGlobalMuons > 0);
1272 
1273  return acceptThisEvent;
1274 }
float zInnerMax
Definition: CSCSkim.h:171
float rExtMax
Definition: CSCSkim.h:174
int nCSCHitsMin
Definition: CSCSkim.h:170
float zLengthMin
Definition: CSCSkim.h:169
float pMin
Definition: CSCSkim.h:168
float zLengthTrMin
Definition: CSCSkim.h:173
bool goodTrack(const reco::Track *pTrack, math::XYZPoint leadPV, trackSelectionParameters parameters, bool debug=false)
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int nValidHitsMin
Definition: CSCSkim.h:176
int nTrHitsMin
Definition: CSCSkim.h:172
Definition: DetId.h:17
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
static constexpr int CSC
Definition: MuonSubdetId.h:12
float redChiSqMax
Definition: CSCSkim.h:175

◆ doCertainChamberSelection()

bool CSCSkim::doCertainChamberSelection ( edm::Handle< CSCWireDigiCollection wires,
edm::Handle< CSCStripDigiCollection strips 
)
private

Definition at line 728 of file CSCSkim.cc.

References CSCDetId::endcap(), DigiDM_cff::strips, whichChamber, whichEndcap, whichRing, whichStation, and DigiDM_cff::wires.

Referenced by filter().

729  {
730  // Loop through the wire DIGIs, looking for a match
731  bool certainChamberIsPresentInWires = false;
732  for (CSCWireDigiCollection::DigiRangeIterator jw = wires->begin(); jw != wires->end(); jw++) {
733  CSCDetId id = (CSCDetId)(*jw).first;
734  int kEndcap = id.endcap();
735  int kRing = id.ring();
736  int kStation = id.station();
737  int kChamber = id.chamber();
738  if ((kEndcap == whichEndcap) && (kStation == whichStation) && (kRing == whichRing) && (kChamber == whichChamber)) {
739  certainChamberIsPresentInWires = true;
740  }
741  } // end wire loop
742 
743  // Loop through the strip DIGIs, looking for a match
744  bool certainChamberIsPresentInStrips = false;
745  for (CSCStripDigiCollection::DigiRangeIterator js = strips->begin(); js != strips->end(); js++) {
746  CSCDetId id = (CSCDetId)(*js).first;
747  int kEndcap = id.endcap();
748  int kRing = id.ring();
749  int kStation = id.station();
750  int kChamber = id.chamber();
751  if ((kEndcap == whichEndcap) && (kStation == whichStation) && (kRing == whichRing) && (kChamber == whichChamber)) {
752  certainChamberIsPresentInStrips = true;
753  }
754  }
755 
756  bool certainChamberIsPresent = certainChamberIsPresentInWires || certainChamberIsPresentInStrips;
757 
758  return certainChamberIsPresent;
759 }
int whichEndcap
Definition: CSCSkim.h:163
int whichChamber
Definition: CSCSkim.h:166
int endcap() const
Definition: CSCDetId.h:85
int whichStation
Definition: CSCSkim.h:164
int whichRing
Definition: CSCSkim.h:165
strips
#turn off noise in all subdetectors simHcalUnsuppressedDigis.doNoise = False mix.digitizers.hcal.doNoise = False simEcalUnsuppressedDigis.doNoise = False mix.digitizers.ecal.doNoise = False simEcalUnsuppressedDigis.doESNoise = False simSiPixelDigis.AddNoise = False mix.digitizers.pixel.AddNoise = False simSiStripDigis.Noise = False mix.digitizers.strip.AddNoise = False
Definition: DigiDM_cff.py:32

◆ doCSCSkimming()

bool CSCSkim::doCSCSkimming ( edm::Handle< CSCRecHit2DCollection cscRecHits,
edm::Handle< CSCSegmentCollection cscSegments 
)
private

Definition at line 375 of file CSCSkim.cc.

References CSCDetId::chamber(), chamberSerial(), dtChamberEfficiency_cfi::cscSegments, CSCDetId::endcap(), nano_mu_digi_cff::float, hxnHitChambers, hxnRecHits, hxnRecHitsSel, hxnSegments, mps_fire::i, createfilelist::int, dqmiolumiharvest::j, kLayer(), CSCDetId::layer(), LogDebug, makeHistograms, minimumHitChambers, minimumSegments, nEventsChambersBothSides, nLayersWithHitsMinimum, funct::pow(), CSCDetId::ring(), CSCDetId::station(), and typeOfSkim.

Referenced by filter().

376  {
377  // how many RecHits in the collection?
378  int nRecHits = cscRecHits->size();
379 
380  // zero the recHit counter
381  int cntRecHit[600];
382  for (int i = 0; i < 600; i++) {
383  cntRecHit[i] = 0;
384  }
385 
386  // ---------------------
387  // Loop over rechits
388  // ---------------------
389 
391  for (recIt = cscRecHits->begin(); recIt != cscRecHits->end(); recIt++) {
392  // which chamber is it?
393  CSCDetId idrec = (CSCDetId)(*recIt).cscDetId();
394  int kEndcap = idrec.endcap();
395  int kRing = idrec.ring();
396  int kStation = idrec.station();
397  int kChamber = idrec.chamber();
398  int kLayer = idrec.layer();
399 
400  // compute chamber serial number
401  int kSerial = chamberSerial(kEndcap, kStation, kRing, kChamber);
402 
403  // increment recHit counter
404  // (each layer is represented by a different power of 10)
405  int kDigit = (int)pow((float)10., (float)(kLayer - 1));
406  cntRecHit[kSerial] += kDigit;
407 
408  } //end rechit loop
409 
410  // ------------------------------------------------------
411  // Are there chambers with the minimum number of hits?
412  // ------------------------------------------------------
413 
414  int nChambersWithMinimalHits = 0;
415  int nChambersWithMinimalHitsPOS = 0;
416  int nChambersWithMinimalHitsNEG = 0;
417  if (nRecHits > 0) {
418  for (int i = 0; i < 600; i++) {
419  if (cntRecHit[i] > 0) {
420  int nLayersWithHits = 0;
421  float dummy = (float)cntRecHit[i];
422  for (int j = 5; j > -1; j--) {
423  float digit = dummy / pow((float)10., (float)j);
424  int kCount = (int)digit;
425  if (kCount > 0)
426  nLayersWithHits++;
427  dummy = dummy - ((float)kCount) * pow((float)10., (float)j);
428  }
429  if (nLayersWithHits > nLayersWithHitsMinimum) {
430  if (i < 300) {
431  nChambersWithMinimalHitsPOS++;
432  } else {
433  nChambersWithMinimalHitsNEG++;
434  }
435  }
436  }
437  }
438  nChambersWithMinimalHits = nChambersWithMinimalHitsPOS + nChambersWithMinimalHitsNEG;
439  }
440 
441  // how many Segments?
442  int nSegments = cscSegments->size();
443 
444  // ----------------------
445  // fill histograms
446  // ----------------------
447 
448  if (makeHistograms) {
449  hxnRecHits->Fill(nRecHits);
450  if (nRecHits > 0) {
451  hxnSegments->Fill(nSegments);
452  hxnHitChambers->Fill(nChambersWithMinimalHits);
453  }
454  if (nChambersWithMinimalHits > 0) {
455  hxnRecHitsSel->Fill(nRecHits);
456  }
457  }
458 
459  // ----------------------
460  // set the filter flag
461  // ----------------------
462  bool basicEvent = (nChambersWithMinimalHits >= minimumHitChambers) && (nSegments >= minimumSegments);
463 
464  bool chambersOnBothSides =
465  ((nChambersWithMinimalHitsPOS >= minimumHitChambers) && (nChambersWithMinimalHitsNEG >= minimumHitChambers));
466 
467  if (chambersOnBothSides) {
469  }
470 
471  bool selectEvent = false;
472  if (typeOfSkim == 1) {
473  selectEvent = basicEvent;
474  }
475  if (typeOfSkim == 2) {
476  selectEvent = chambersOnBothSides;
477  }
478 
479  // debug
480  LogDebug("[CSCSkim]") << "----- nRecHits = " << nRecHits
481  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits << "\tnSegments = " << nSegments
482  << "\tselect? " << selectEvent << std::endl;
483 
484  /*
485  if ((nChambersWithMinimalHitsPOS >= minimumHitChambers) && (nChambersWithMinimalHitsNEG >= minimumHitChambers)) {
486  std::cout << "\n==========================================================================\n"
487  << "\tinteresting event - chambers hit on both sides\n"
488  << "\t " << nEventsAnalyzed
489  << "\trun " << iRun << "\tevent " << iEvent << std::endl;
490  std::cout << "----- nRecHits = " << nRecHits
491  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits
492  << "\tnSegments = " << nSegments
493  << "\tselect? " << selectEvent << std::endl;
494  for (int i = 0; i < 600; i++) {
495  if (cntRecHit[i] > 0) {
496  cout << "\t\t" << i << "\tcntRecHit= " << cntRecHit[i] << std::endl;
497  }
498  }
499  std::cout << "==========================================================================\n\n" ;
500  }
501  */
502 
503  return selectEvent;
504 }
TH1F * hxnRecHitsSel
Definition: CSCSkim.h:182
int nLayersWithHitsMinimum
Definition: CSCSkim.h:157
int nEventsChambersBothSides
Definition: CSCSkim.h:119
TH1F * hxnRecHits
Definition: CSCSkim.h:179
int layer() const
Definition: CSCDetId.h:56
int minimumHitChambers
Definition: CSCSkim.h:158
int typeOfSkim
Definition: CSCSkim.h:156
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int minimumSegments
Definition: CSCSkim.h:159
TH1F * hxnSegments
Definition: CSCSkim.h:180
int chamber() const
Definition: CSCDetId.h:62
bool makeHistograms
Definition: CSCSkim.h:161
static const std::string kLayer("layer")
int station() const
Definition: CSCDetId.h:79
int endcap() const
Definition: CSCDetId.h:85
TH1F * hxnHitChambers
Definition: CSCSkim.h:181
int chamberSerial(int kE, int kS, int kR, int kCh)
Definition: CSCSkim.cc:1280
int ring() const
Definition: CSCDetId.h:68
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
#define LogDebug(id)

◆ doDTOverlap()

bool CSCSkim::doDTOverlap ( edm::Handle< CSCSegmentCollection cscSegments)
private

Definition at line 766 of file CSCSkim.cc.

References dtChamberEfficiency_cfi::cscSegments, CSCDetId::endcap(), mps_fire::i, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, and nhits.

Referenced by filter().

766  {
767  const float chisqMax = 100.;
768  const int nhitsMin = 5;
769  const int maxNSegments = 3;
770 
771  // initialize
772  bool DTOverlapCandidate = false;
773  int cntMEP13[36];
774  int cntMEN13[36];
775  int cntMEP22[36];
776  int cntMEN22[36];
777  int cntMEP32[36];
778  int cntMEN32[36];
779  for (int i = 0; i < 36; ++i) {
780  cntMEP13[i] = 0;
781  cntMEN13[i] = 0;
782  cntMEP22[i] = 0;
783  cntMEN22[i] = 0;
784  cntMEP32[i] = 0;
785  cntMEN32[i] = 0;
786  }
787 
788  // -----------------------
789  // loop over segments
790  // -----------------------
791 
792  int nSegments = cscSegments->size();
793  if (nSegments < 2)
794  return DTOverlapCandidate;
795 
796  for (CSCSegmentCollection::const_iterator it = cscSegments->begin(); it != cscSegments->end(); it++) {
797  // which chamber?
798  CSCDetId id = (CSCDetId)(*it).cscDetId();
799  int kEndcap = id.endcap();
800  int kStation = id.station();
801  int kRing = id.ring();
802  int kChamber = id.chamber();
803  // segment information
804  float chisq = (*it).chi2();
805  int nhits = (*it).nRecHits();
806  bool goodSegment = (chisq < chisqMax) && (nhits >= nhitsMin);
807  if (goodSegment) {
808  if ((kStation == 1) && (kRing == 3)) {
809  if (kEndcap == 1) {
810  cntMEP13[kChamber - 1]++;
811  }
812  if (kEndcap == 2) {
813  cntMEN13[kChamber - 1]++;
814  }
815  }
816  if ((kStation == 2) && (kRing == 2)) {
817  if (kEndcap == 1) {
818  cntMEP22[kChamber - 1]++;
819  }
820  if (kEndcap == 2) {
821  cntMEN22[kChamber - 1]++;
822  }
823  }
824  if ((kStation == 3) && (kRing == 2)) {
825  if (kEndcap == 1) {
826  cntMEP32[kChamber - 1]++;
827  }
828  if (kEndcap == 2) {
829  cntMEN32[kChamber - 1]++;
830  }
831  }
832  } // this is a good segment
833  } // end loop over segments
834 
835  // ---------------------------------------------
836  // veto messy events
837  // ---------------------------------------------
838  bool tooManySegments = false;
839  for (int i = 0; i < 36; ++i) {
840  if ((cntMEP13[i] > maxNSegments) || (cntMEN13[i] > maxNSegments) || (cntMEP22[i] > maxNSegments) ||
841  (cntMEN22[i] > maxNSegments) || (cntMEP32[i] > maxNSegments) || (cntMEN32[i] > maxNSegments))
842  tooManySegments = true;
843  }
844  if (tooManySegments) {
845  return DTOverlapCandidate;
846  }
847 
848  // ---------------------------------------------
849  // check for relevant matchup of segments
850  // ---------------------------------------------
851  bool matchup = false;
852  for (int i = 0; i < 36; ++i) {
853  if ((cntMEP13[i] > 0) && (cntMEP22[i] + cntMEP32[i] > 0)) {
854  matchup = true;
855  }
856  if ((cntMEN13[i] > 0) && (cntMEN22[i] + cntMEN32[i] > 0)) {
857  matchup = true;
858  }
859  }
860  /*
861  if (matchup) {
862  std::cout << "\tYYY looks like a good event. Select!\n";
863  std::cout << "-- pos endcap --\n"
864  << "ME1/3: ";
865  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP13[k];}
866  std::cout << "\nME2/2: ";
867  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP22[k];}
868  std::cout << "\nME3/2: ";
869  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP32[k];}
870  std::cout << std::endl;
871  }
872  */
873 
874  // set the selection flag
875  DTOverlapCandidate = matchup;
876  return DTOverlapCandidate;
877 }
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int endcap() const
Definition: CSCDetId.h:85

◆ doHaloLike()

bool CSCSkim::doHaloLike ( edm::Handle< CSCSegmentCollection cscSegments)
private

Definition at line 885 of file CSCSkim.cc.

References dtChamberEfficiency_cfi::cscSegments, CSCDetId::endcap(), mps_fire::i, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, and nhits.

Referenced by filter().

885  {
886  const float chisqMax = 100.;
887  const int nhitsMin = 5; // on a segment
888  const int maxNSegments = 3; // in a chamber
889 
890  // initialize
891  bool HaloLike = false;
892  int cntMEP11[36];
893  int cntMEN11[36];
894  int cntMEP12[36];
895  int cntMEN12[36];
896  int cntMEP21[36];
897  int cntMEN21[36];
898  int cntMEP31[36];
899  int cntMEN31[36];
900  int cntMEP41[36];
901  int cntMEN41[36];
902  for (int i = 0; i < 36; ++i) {
903  cntMEP11[i] = 0;
904  cntMEN11[i] = 0;
905  cntMEP12[i] = 0;
906  cntMEN12[i] = 0;
907  cntMEP21[i] = 0;
908  cntMEN21[i] = 0;
909  cntMEP31[i] = 0;
910  cntMEN31[i] = 0;
911  cntMEP41[i] = 0;
912  cntMEN41[i] = 0;
913  }
914 
915  // -----------------------
916  // loop over segments
917  // -----------------------
918  int nSegments = cscSegments->size();
919  if (nSegments < 4)
920  return HaloLike;
921 
922  for (CSCSegmentCollection::const_iterator it = cscSegments->begin(); it != cscSegments->end(); it++) {
923  // which chamber?
924  CSCDetId id = (CSCDetId)(*it).cscDetId();
925  int kEndcap = id.endcap();
926  int kStation = id.station();
927  int kRing = id.ring();
928  int kChamber = id.chamber();
929  // segment information
930  float chisq = (*it).chi2();
931  int nhits = (*it).nRecHits();
932  bool goodSegment = (chisq < chisqMax) && (nhits >= nhitsMin);
933  if (goodSegment) {
934  if ((kStation == 1) && (kRing == 1)) {
935  if (kEndcap == 1) {
936  cntMEP11[kChamber - 1]++;
937  }
938  if (kEndcap == 2) {
939  cntMEN11[kChamber - 1]++;
940  }
941  }
942  if ((kStation == 1) && (kRing == 2)) {
943  if (kEndcap == 1) {
944  cntMEP12[kChamber - 1]++;
945  }
946  if (kEndcap == 2) {
947  cntMEN12[kChamber - 1]++;
948  }
949  }
950  if ((kStation == 2) && (kRing == 1)) {
951  if (kEndcap == 1) {
952  cntMEP21[kChamber - 1]++;
953  }
954  if (kEndcap == 2) {
955  cntMEN21[kChamber - 1]++;
956  }
957  }
958  if ((kStation == 3) && (kRing == 1)) {
959  if (kEndcap == 1) {
960  cntMEP31[kChamber - 1]++;
961  }
962  if (kEndcap == 2) {
963  cntMEN31[kChamber - 1]++;
964  }
965  }
966  if ((kStation == 4) && (kRing == 1)) {
967  if (kEndcap == 1) {
968  cntMEP41[kChamber - 1]++;
969  }
970  if (kEndcap == 2) {
971  cntMEN41[kChamber - 1]++;
972  }
973  }
974  } // this is a good segment
975  } // end loop over segments
976 
977  // ---------------------------------------------
978  // veto messy events
979  // ---------------------------------------------
980  bool tooManySegments = false;
981  for (int i = 0; i < 36; ++i) {
982  if ((cntMEP11[i] > 3 * maxNSegments) || (cntMEN11[i] > 3 * maxNSegments) || (cntMEP12[i] > maxNSegments) ||
983  (cntMEN12[i] > maxNSegments) || (cntMEP21[i] > maxNSegments) || (cntMEN21[i] > maxNSegments) ||
984  (cntMEP31[i] > maxNSegments) || (cntMEN31[i] > maxNSegments) || (cntMEP41[i] > maxNSegments) ||
985  (cntMEN41[i] > maxNSegments))
986  tooManySegments = true;
987  }
988  if (tooManySegments) {
989  return HaloLike;
990  }
991 
992  // ---------------------------------------------
993  // check for relevant matchup of segments
994  // ---------------------------------------------
995  bool matchup = false;
996  for (int i = 0; i < 36; ++i) {
997  if ((cntMEP11[i] + cntMEP12[i] > 0) && (cntMEP21[i] > 0) && (cntMEP31[i] > 0) && (cntMEP41[i] > 0)) {
998  matchup = true;
999  }
1000  if ((cntMEN11[i] + cntMEN12[i] > 0) && (cntMEN21[i] > 0) && (cntMEN31[i] > 0) && (cntMEN41[i] > 0)) {
1001  matchup = true;
1002  }
1003  }
1004  /*
1005  if (matchup) {
1006  std::cout << "\tYYY looks like a good event. Select!\n";
1007  std::cout << "-- pos endcap --\n"
1008  << "ME1/1: ";
1009  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP11[k];}
1010  std::cout << "\nME1/2: ";
1011  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP12[k];}
1012  std::cout << "\nME2/1: ";
1013  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP21[k];}
1014  std::cout << "\nME3/1: ";
1015  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP31[k];}
1016  std::cout << "\nME4/1: ";
1017  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEP41[k];}
1018  std::cout << std::endl;
1019  std::cout << "-- neg endcap --\n"
1020  << "ME1/1: ";
1021  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN11[k];}
1022  std::cout << "\nME1/2: ";
1023  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN12[k];}
1024  std::cout << "\nME2/1: ";
1025  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN21[k];}
1026  std::cout << "\nME3/1: ";
1027  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN31[k];}
1028  std::cout << "\nME4/1: ";
1029  for (int k=0; k<36; ++k) {std::cout << " " << setw(3) << cntMEN41[k];}
1030  std::cout << std::endl;
1031  std::cout << "\tn Analyzed = " << nEventsAnalyzed << "\tn Halo-like = " << nEventsHaloLike << std::endl;
1032  }
1033  */
1034 
1035  // set the selection flag
1036  HaloLike = matchup;
1037  return HaloLike;
1038 }
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int endcap() const
Definition: CSCDetId.h:85

◆ doLongSATrack()

bool CSCSkim::doLongSATrack ( edm::Handle< reco::TrackCollection saTracks)
private

Definition at line 1043 of file CSCSkim.cc.

References MuonSubdetId::CSC, hcalRecHitTable_cff::detId, SiStripPI::min, DetId::Muon, HLT_2024v14_cff::muon, nCSCHitsMin, findAndChange::op, singleTopDQM_cfi::select, and zInnerMax.

Referenced by filter().

1043  {
1044  const float zDistanceMax = 2500.;
1045  const float zDistanceMin = 700.;
1046  const int nCSCHitsMin = 25;
1047  const int nCSCHitsMax = 50;
1048  const float zInnerMax = 80000.;
1049 
1050  const int nNiceMuonsMin = 1;
1051 
1052  //
1053  // Loop through the track collection and test each one
1054  //
1055 
1056  int nNiceMuons = 0;
1057 
1058  for (reco::TrackCollection::const_iterator muon = saMuons->begin(); muon != saMuons->end(); ++muon) {
1059  // basic information
1060  math::XYZVector innerMo = muon->innerMomentum();
1061  GlobalVector im(innerMo.x(), innerMo.y(), innerMo.z());
1062  math::XYZPoint innerPo = muon->innerPosition();
1063  GlobalPoint ip(innerPo.x(), innerPo.y(), innerPo.z());
1064  math::XYZPoint outerPo = muon->outerPosition();
1065  GlobalPoint op(outerPo.x(), outerPo.y(), outerPo.z());
1066  float zInner = ip.z();
1067  float zOuter = op.z();
1068  float zDistance = fabs(zOuter - zInner);
1069 
1070  // loop over hits
1071  int nCSCHits = 0;
1072  for (trackingRecHit_iterator hit = muon->recHitsBegin(); hit != muon->recHitsEnd(); ++hit) {
1073  const DetId detId((*hit)->geographicalId());
1074  if (detId.det() == DetId::Muon) {
1075  if (detId.subdetId() == MuonSubdetId::CSC) {
1076  //CSCDetId cscId(detId.rawId());
1077  //int chamberId = cscId.chamber();
1078  nCSCHits++;
1079  }
1080  }
1081  }
1082 
1083  // is this a nice muon?
1084  if ((zDistance < zDistanceMax) && (zDistance > zDistanceMin) && (nCSCHits > nCSCHitsMin) &&
1085  (nCSCHits < nCSCHitsMax) && (min(fabs(zInner), fabs(zOuter)) < zInnerMax) &&
1086  (fabs(innerMo.z()) > 0.000000001)) {
1087  nNiceMuons++;
1088  }
1089  }
1090 
1091  bool select = (nNiceMuons >= nNiceMuonsMin);
1092 
1093  return select;
1094 }
float zInnerMax
Definition: CSCSkim.h:171
int nCSCHitsMin
Definition: CSCSkim.h:170
Definition: DetId.h:17
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
static constexpr int CSC
Definition: MuonSubdetId.h:12

◆ doMessyEventSkimming()

bool CSCSkim::doMessyEventSkimming ( edm::Handle< CSCRecHit2DCollection cscRecHits,
edm::Handle< CSCSegmentCollection cscSegments 
)
private

Definition at line 592 of file CSCSkim.cc.

References CSCDetId::chamber(), chamberSerial(), dtChamberEfficiency_cfi::cscSegments, CSCDetId::endcap(), nano_mu_digi_cff::float, mps_fire::i, createfilelist::int, dqmiolumiharvest::j, kLayer(), CSCDetId::layer(), LogDebug, makeHistogramsForMessyEvents, mevnChambers0, mevnChambers1, mevnRecHits0, mevnRecHits1, mevnSegments0, mevnSegments1, nLayersWithHitsMinimum, funct::pow(), CSCDetId::ring(), and CSCDetId::station().

Referenced by filter().

593  {
594  // how many RecHits in the collection?
595  int nRecHits = cscRecHits->size();
596 
597  // zero the recHit counter
598  int cntRecHit[600];
599  for (int i = 0; i < 600; i++) {
600  cntRecHit[i] = 0;
601  }
602 
603  // ---------------------
604  // Loop over rechits
605  // ---------------------
606 
608  for (recIt = cscRecHits->begin(); recIt != cscRecHits->end(); recIt++) {
609  // which chamber is it?
610  CSCDetId idrec = (CSCDetId)(*recIt).cscDetId();
611  int kEndcap = idrec.endcap();
612  int kRing = idrec.ring();
613  int kStation = idrec.station();
614  int kChamber = idrec.chamber();
615  int kLayer = idrec.layer();
616 
617  // compute chamber serial number
618  int kSerial = chamberSerial(kEndcap, kStation, kRing, kChamber);
619 
620  // increment recHit counter
621  // (each layer is represented by a different power of 10)
622  int kDigit = (int)pow((float)10., (float)(kLayer - 1));
623  cntRecHit[kSerial] += kDigit;
624 
625  } //end rechit loop
626 
627  // ------------------------------------------------------
628  // Are there chambers with the minimum number of hits?
629  // ------------------------------------------------------
630 
631  int nChambersWithMinimalHits = 0;
632  int nChambersWithMinimalHitsPOS = 0;
633  int nChambersWithMinimalHitsNEG = 0;
634  if (nRecHits > 0) {
635  for (int i = 0; i < 600; i++) {
636  if (cntRecHit[i] > 0) {
637  int nLayersWithHits = 0;
638  float dummy = (float)cntRecHit[i];
639  for (int j = 5; j > -1; j--) {
640  float digit = dummy / pow((float)10., (float)j);
641  int kCount = (int)digit;
642  if (kCount > 0)
643  nLayersWithHits++;
644  dummy = dummy - ((float)kCount) * pow((float)10., (float)j);
645  }
646  if (nLayersWithHits > nLayersWithHitsMinimum) {
647  if (i < 300) {
648  nChambersWithMinimalHitsPOS++;
649  } else {
650  nChambersWithMinimalHitsNEG++;
651  }
652  }
653  }
654  }
655  nChambersWithMinimalHits = nChambersWithMinimalHitsPOS + nChambersWithMinimalHitsNEG;
656  }
657 
658  // how many Segments?
659  int nSegments = cscSegments->size();
660 
661  // ----------------------
662  // fill histograms
663  // ----------------------
664 
666  if (nRecHits > 8) {
667  mevnRecHits0->Fill(nRecHits);
668  mevnChambers0->Fill(nChambersWithMinimalHits);
669  mevnSegments0->Fill(nSegments);
670  }
671  if (nRecHits > 54) {
672  double dummy = (double)nRecHits;
673  if (dummy > 299.9)
674  dummy = 299.9;
675  mevnRecHits1->Fill(dummy);
676  dummy = (double)nChambersWithMinimalHits;
677  if (dummy > 49.9)
678  dummy = 49.9;
679  mevnChambers1->Fill(dummy);
680  dummy = (double)nSegments;
681  if (dummy > 29.9)
682  dummy = 29.9;
683  mevnSegments1->Fill(dummy);
684  }
685  }
686 
687  // ----------------------
688  // set the filter flag
689  // ----------------------
690 
691  bool selectEvent = false;
692  if ((nRecHits > 54) && (nChambersWithMinimalHits > 5)) {
693  selectEvent = true;
694  }
695 
696  // debug
697  LogDebug("[CSCSkim]") << "----- nRecHits = " << nRecHits
698  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits << "\tnSegments = " << nSegments
699  << "\tselect? " << selectEvent << std::endl;
700 
701  /*
702  if (selectEvent) {
703  std::cout << "\n==========================================================================\n"
704  << "\tmessy event!\n"
705  << "\t " << nEventsAnalyzed
706  << "\trun " << iRun << "\tevent " << iEvent << std::endl;
707  std::cout << "----- nRecHits = " << nRecHits
708  << "\tnChambersWithMinimalHits = " << nChambersWithMinimalHits
709  << "\tnSegments = " << nSegments
710  << "\tselect? " << selectEvent << std::endl;
711  for (int i = 0; i < 600; i++) {
712  if (cntRecHit[i] > 0) {
713  cout << "\t\t" << i << "\tcntRecHit= " << cntRecHit[i] << std::endl;
714  }
715  }
716  std::cout << "==========================================================================\n\n" ;
717  }
718  */
719 
720  return selectEvent;
721 }
bool makeHistogramsForMessyEvents
Definition: CSCSkim.h:162
int nLayersWithHitsMinimum
Definition: CSCSkim.h:157
TH1F * mevnRecHits0
Definition: CSCSkim.h:183
TH1F * mevnRecHits1
Definition: CSCSkim.h:186
TH1F * mevnSegments0
Definition: CSCSkim.h:185
TH1F * mevnSegments1
Definition: CSCSkim.h:188
int layer() const
Definition: CSCDetId.h:56
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int chamber() const
Definition: CSCDetId.h:62
static const std::string kLayer("layer")
int station() const
Definition: CSCDetId.h:79
int endcap() const
Definition: CSCDetId.h:85
TH1F * mevnChambers1
Definition: CSCSkim.h:187
int chamberSerial(int kE, int kS, int kR, int kCh)
Definition: CSCSkim.cc:1280
int ring() const
Definition: CSCDetId.h:68
TH1F * mevnChambers0
Definition: CSCSkim.h:184
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
#define LogDebug(id)

◆ doOverlapSkimming()

bool CSCSkim::doOverlapSkimming ( edm::Handle< CSCSegmentCollection cscSegments)
private

Definition at line 510 of file CSCSkim.cc.

References chamberSerial(), dtChamberEfficiency_cfi::cscSegments, CSCDetId::endcap(), mps_fire::i, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, and nhits.

Referenced by filter().

510  {
511  const int nhitsMinimum = 4;
512  const float chisqMaximum = 100.;
513  const int nAllMaximum = 3;
514 
515  // how many Segments?
516  // int nSegments = cscSegments->size();
517 
518  // zero arrays
519  int nAll[600];
520  int nGood[600];
521  for (int i = 0; i < 600; i++) {
522  nAll[i] = 0;
523  nGood[i] = 0;
524  }
525 
526  // -----------------------
527  // loop over segments
528  // -----------------------
529  for (CSCSegmentCollection::const_iterator it = cscSegments->begin(); it != cscSegments->end(); it++) {
530  // which chamber?
531  CSCDetId id = (CSCDetId)(*it).cscDetId();
532  int kEndcap = id.endcap();
533  int kStation = id.station();
534  int kRing = id.ring();
535  int kChamber = id.chamber();
536  int kSerial = chamberSerial(kEndcap, kStation, kRing, kChamber);
537 
538  // segment information
539  float chisq = (*it).chi2();
540  int nhits = (*it).nRecHits();
541 
542  // is this a good segment?
543  bool goodSegment = (nhits >= nhitsMinimum) && (chisq < chisqMaximum);
544 
545  /*
546  LocalPoint localPos = (*it).localPosition();
547  float segX = localPos.x();
548  float segY = localPos.y();
549  std::cout << "E/S/R/Ch: " << kEndcap << "/" << kStation << "/" << kRing << "/" << kChamber
550  << "\tnhits/chisq: " << nhits << "/" << chisq
551  << "\tX/Y: " << segX << "/" << segY
552  << "\tgood? " << goodSegment << std::endl;
553  */
554 
555  // count
556  nAll[kSerial - 1]++;
557  if (goodSegment)
558  nGood[kSerial]++;
559 
560  } // end loop over segments
561 
562  //----------------------
563  // select the event
564  //----------------------
565 
566  // does any chamber have too many segments?
567  bool messyChamber = false;
568  for (int i = 0; i < 600; i++) {
569  if (nAll[i] > nAllMaximum)
570  messyChamber = true;
571  }
572 
573  // are there consecutive chambers with good segments
574  // (This is a little sloppy but is probably fine for skimming...)
575  bool consecutiveChambers = false;
576  for (int i = 0; i < 599; i++) {
577  if ((nGood[i] > 0) && (nGood[i + 1] > 0))
578  consecutiveChambers = true;
579  }
580 
581  bool selectThisEvent = !messyChamber && consecutiveChambers;
582 
583  return selectThisEvent;
584 }
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int endcap() const
Definition: CSCDetId.h:85
int chamberSerial(int kE, int kS, int kR, int kCh)
Definition: CSCSkim.cc:1280

◆ endJob()

void CSCSkim::endJob ( void  )
overridevirtual

Reimplemented from edm::one::EDFilterBase.

Definition at line 167 of file CSCSkim.cc.

References nano_mu_digi_cff::float, HLT_2024v14_cff::fraction, hxnHitChambers, hxnRecHits, hxnRecHitsSel, hxnSegments, LogDebug, makeHistograms, makeHistogramsForMessyEvents, mevnChambers0, mevnChambers1, mevnRecHits0, mevnRecHits1, mevnSegments0, mevnSegments1, nEventsAnalyzed, nEventsCertainChamber, nEventsChambersBothSides, nEventsDTOverlap, nEventsForBFieldStudies, nEventsHaloLike, nEventsLongSATrack, nEventsMessy, nEventsOverlappingChambers, nEventsSelected, theHistogramFile, and typeOfSkim.

167  {
168  // Write out results
169 
170  float fraction = 0.;
171  if (nEventsAnalyzed > 0) {
173  }
174 
175  LogInfo("[CSCSkim] Summary") << "\n\n\t====== CSCSkim ==========================================================\n"
176  << "\t\ttype of skim ...............................\t" << typeOfSkim << "\n"
177  << "\t\tevents analyzed ..............\t" << nEventsAnalyzed << "\n"
178  << "\t\tevents selected ..............\t" << nEventsSelected
179  << "\tfraction= " << fraction << std::endl
180  << "\t\tevents chambers both sides ...\t" << nEventsChambersBothSides << "\n"
181  << "\t\tevents w/ overlaps .......... \t" << nEventsOverlappingChambers << "\n"
182  << "\t\tevents lots of hit chambers . \t" << nEventsMessy << "\n"
183  << "\t\tevents from certain chamber . \t" << nEventsCertainChamber << "\n"
184  << "\t\tevents in DT-CSC overlap .... \t" << nEventsDTOverlap << "\n"
185  << "\t\tevents halo-like ............ \t" << nEventsHaloLike << "\n"
186  << "\t\tevents w/ long SA track ..... \t" << nEventsLongSATrack << "\n"
187  << "\t\tevents good for BField ..... \t" << nEventsForBFieldStudies << "\n"
188  << "\t=========================================================================\n\n";
189 
191  // Write the histos to file
192  LogDebug("[CSCSkim]") << "======= write out my histograms ====\n";
193  theHistogramFile->cd();
194  if (makeHistograms) {
195  hxnRecHits->Write();
196  hxnSegments->Write();
197  hxnHitChambers->Write();
198  hxnRecHitsSel->Write();
199  }
201  mevnRecHits0->Write();
202  mevnChambers0->Write();
203  mevnSegments0->Write();
204  mevnRecHits1->Write();
205  mevnChambers1->Write();
206  mevnSegments1->Write();
207  }
208  theHistogramFile->Close();
209  }
210 }
TH1F * hxnRecHitsSel
Definition: CSCSkim.h:182
bool makeHistogramsForMessyEvents
Definition: CSCSkim.h:162
int nEventsForBFieldStudies
Definition: CSCSkim.h:126
int nEventsChambersBothSides
Definition: CSCSkim.h:119
TH1F * mevnRecHits0
Definition: CSCSkim.h:183
TH1F * mevnRecHits1
Definition: CSCSkim.h:186
TH1F * mevnSegments0
Definition: CSCSkim.h:185
TH1F * mevnSegments1
Definition: CSCSkim.h:188
TH1F * hxnRecHits
Definition: CSCSkim.h:179
int nEventsDTOverlap
Definition: CSCSkim.h:123
int nEventsOverlappingChambers
Definition: CSCSkim.h:120
int typeOfSkim
Definition: CSCSkim.h:156
int nEventsSelected
Definition: CSCSkim.h:118
TH1F * hxnSegments
Definition: CSCSkim.h:180
bool makeHistograms
Definition: CSCSkim.h:161
int nEventsLongSATrack
Definition: CSCSkim.h:125
int nEventsAnalyzed
Definition: CSCSkim.h:117
Log< level::Info, false > LogInfo
int nEventsCertainChamber
Definition: CSCSkim.h:122
TH1F * mevnChambers1
Definition: CSCSkim.h:187
TFile * theHistogramFile
Definition: CSCSkim.h:133
int nEventsMessy
Definition: CSCSkim.h:121
TH1F * hxnHitChambers
Definition: CSCSkim.h:181
int nEventsHaloLike
Definition: CSCSkim.h:124
TH1F * mevnChambers0
Definition: CSCSkim.h:184
#define LogDebug(id)

◆ filter()

bool CSCSkim::filter ( edm::Event event,
const edm::EventSetup eventSetup 
)
overridevirtual

Implements edm::one::EDFilterBase.

Definition at line 215 of file CSCSkim.cc.

References dtChamberEfficiency_cfi::cscSegments, doBFieldStudySelection(), doCertainChamberSelection(), doCSCSkimming(), doDTOverlap(), doHaloLike(), doLongSATrack(), doMessyEventSkimming(), doOverlapSkimming(), options_cfi::eventSetup, glm_token, iEvent, iRun, LogDebug, m_CSCGeomToken, nEventsAnalyzed, nEventsCertainChamber, nEventsDTOverlap, nEventsForBFieldStudies, nEventsHaloLike, nEventsLongSATrack, nEventsMessy, nEventsOverlappingChambers, nEventsSelected, rh_token, sam_token, sdr_token, sds_token, seg_token, DigiDM_cff::strips, DiMuonV_cfg::tracks, trk_token, typeOfSkim, wdr_token, wds_token, and DigiDM_cff::wires.

215  {
216  // increment counter
217  nEventsAnalyzed++;
218 
219  iRun = event.id().run();
220  iEvent = event.id().event();
221 
222  LogDebug("[CSCSkim] EventInfo") << "Run: " << iRun << "\tEvent: " << iEvent << "\tn Analyzed: " << nEventsAnalyzed;
223 
224  // Get the CSC Geometry :
225  ESHandle<CSCGeometry> cscGeom = eventSetup.getHandle(m_CSCGeomToken);
226 
227  // Get the DIGI collections
230 
231  if (event.eventAuxiliary().isRealData()) {
232  event.getByToken(wdr_token, wires);
233  event.getByToken(sdr_token, strips);
234  } else {
235  event.getByToken(wds_token, wires);
236  event.getByToken(sds_token, strips);
237  }
238 
239  // Get the RecHits collection :
241  event.getByToken(rh_token, cscRecHits);
242 
243  // get CSC segment collection
245  event.getByToken(seg_token, cscSegments);
246 
247  // get the cosmic muons collection
249  if (typeOfSkim == 8) {
250  event.getByToken(sam_token, saMuons);
251  }
252 
253  // get the stand-alone muons collection
256  if (typeOfSkim == 9) {
257  event.getByToken(sam_token, saMuons);
258  event.getByToken(trk_token, tracks);
259  event.getByToken(glm_token, gMuons);
260  }
261 
262  //======================================
263  // evaluate the skimming routines
264  //======================================
265 
266  // basic skimming
267  bool basicEvent = false;
268  if (typeOfSkim == 1 || typeOfSkim == 2) {
269  basicEvent = doCSCSkimming(cscRecHits, cscSegments);
270  }
271 
272  // overlapping chamber skim
273  bool goodOverlapEvent = false;
274  if (typeOfSkim == 3) {
275  goodOverlapEvent = doOverlapSkimming(cscSegments);
276  if (goodOverlapEvent) {
278  }
279  }
280 
281  // messy events skim
282  bool messyEvent = false;
283  if (typeOfSkim == 4) {
284  messyEvent = doMessyEventSkimming(cscRecHits, cscSegments);
285  if (messyEvent) {
286  nEventsMessy++;
287  }
288  }
289 
290  // select events with DIGIs in a certain chamber
291  bool hasChamber = false;
292  if (typeOfSkim == 5) {
293  hasChamber = doCertainChamberSelection(wires, strips);
294  if (hasChamber) {
296  }
297  }
298 
299  // select events in the DT-CSC overlap region
300  bool DTOverlapCandidate = false;
301  if (typeOfSkim == 6) {
302  DTOverlapCandidate = doDTOverlap(cscSegments);
303  if (DTOverlapCandidate) {
305  }
306  }
307 
308  // select halo-like events
309  bool HaloLike = false;
310  if (typeOfSkim == 7) {
311  HaloLike = doHaloLike(cscSegments);
312  if (HaloLike) {
313  nEventsHaloLike++;
314  }
315  }
316 
317  // select long cosmic tracks
318  bool LongSATrack = false;
319  if (typeOfSkim == 8) {
320  LongSATrack = doLongSATrack(saMuons);
321  if (LongSATrack) {
323  }
324  }
325 
326  // select events suitable for a B-field study. They have tracks in the tracker.
327  bool GoodForBFieldStudy = false;
328  if (typeOfSkim == 9) {
329  GoodForBFieldStudy = doBFieldStudySelection(saMuons, tracks, gMuons);
330  if (GoodForBFieldStudy) {
332  }
333  }
334 
335  // set filter flag
336  bool selectThisEvent = false;
337  if (typeOfSkim == 1 || typeOfSkim == 2) {
338  selectThisEvent = basicEvent;
339  }
340  if (typeOfSkim == 3) {
341  selectThisEvent = goodOverlapEvent;
342  }
343  if (typeOfSkim == 4) {
344  selectThisEvent = messyEvent;
345  }
346  if (typeOfSkim == 5) {
347  selectThisEvent = hasChamber;
348  }
349  if (typeOfSkim == 6) {
350  selectThisEvent = DTOverlapCandidate;
351  }
352  if (typeOfSkim == 7) {
353  selectThisEvent = HaloLike;
354  }
355  if (typeOfSkim == 8) {
356  selectThisEvent = LongSATrack;
357  }
358  if (typeOfSkim == 9) {
359  selectThisEvent = GoodForBFieldStudy;
360  }
361 
362  if (selectThisEvent) {
363  nEventsSelected++;
364  }
365 
366  return selectThisEvent;
367 }
edm::EDGetTokenT< CSCSegmentCollection > seg_token
Definition: CSCSkim.h:149
edm::EDGetTokenT< reco::TrackCollection > sam_token
Definition: CSCSkim.h:150
int nEventsForBFieldStudies
Definition: CSCSkim.h:126
edm::EDGetTokenT< CSCWireDigiCollection > wds_token
Definition: CSCSkim.h:143
bool doLongSATrack(edm::Handle< reco::TrackCollection > saTracks)
Definition: CSCSkim.cc:1043
bool doOverlapSkimming(edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:510
int nEventsDTOverlap
Definition: CSCSkim.h:123
bool doDTOverlap(edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:766
int nEventsOverlappingChambers
Definition: CSCSkim.h:120
int typeOfSkim
Definition: CSCSkim.h:156
edm::EDGetTokenT< reco::TrackCollection > trk_token
Definition: CSCSkim.h:151
int nEventsSelected
Definition: CSCSkim.h:118
edm::EDGetTokenT< reco::MuonCollection > glm_token
Definition: CSCSkim.h:152
edm::EDGetTokenT< CSCStripDigiCollection > sds_token
Definition: CSCSkim.h:144
edm::EDGetTokenT< CSCRecHit2DCollection > rh_token
Definition: CSCSkim.h:148
int nEventsLongSATrack
Definition: CSCSkim.h:125
bool doCSCSkimming(edm::Handle< CSCRecHit2DCollection > cscRecHits, edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:375
bool doCertainChamberSelection(edm::Handle< CSCWireDigiCollection > wires, edm::Handle< CSCStripDigiCollection > strips)
Definition: CSCSkim.cc:728
bool doBFieldStudySelection(edm::Handle< reco::TrackCollection > saTracks, edm::Handle< reco::TrackCollection > Tracks, edm::Handle< reco::MuonCollection > gMuons)
Definition: CSCSkim.cc:1104
int nEventsAnalyzed
Definition: CSCSkim.h:117
int nEventsCertainChamber
Definition: CSCSkim.h:122
bool doMessyEventSkimming(edm::Handle< CSCRecHit2DCollection > cscRecHits, edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:592
int nEventsMessy
Definition: CSCSkim.h:121
const edm::ESGetToken< CSCGeometry, MuonGeometryRecord > m_CSCGeomToken
Definition: CSCSkim.h:140
strips
#turn off noise in all subdetectors simHcalUnsuppressedDigis.doNoise = False mix.digitizers.hcal.doNoise = False simEcalUnsuppressedDigis.doNoise = False mix.digitizers.ecal.doNoise = False simEcalUnsuppressedDigis.doESNoise = False simSiPixelDigis.AddNoise = False mix.digitizers.pixel.AddNoise = False simSiStripDigis.Noise = False mix.digitizers.strip.AddNoise = False
Definition: DigiDM_cff.py:32
int iEvent
Definition: CSCSkim.h:130
edm::EDGetTokenT< CSCWireDigiCollection > wdr_token
Definition: CSCSkim.h:145
int nEventsHaloLike
Definition: CSCSkim.h:124
edm::EDGetTokenT< CSCStripDigiCollection > sdr_token
Definition: CSCSkim.h:146
int iRun
Definition: CSCSkim.h:129
Definition: event.py:1
#define LogDebug(id)
bool doHaloLike(edm::Handle< CSCSegmentCollection > cscSegments)
Definition: CSCSkim.cc:885

Member Data Documentation

◆ demandChambersBothSides

bool CSCSkim::demandChambersBothSides
private

Definition at line 160 of file CSCSkim.h.

Referenced by CSCSkim().

◆ glm_token

edm::EDGetTokenT<reco::MuonCollection> CSCSkim::glm_token
private

Definition at line 152 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ histogramFileName

std::string CSCSkim::histogramFileName
private

Definition at line 137 of file CSCSkim.h.

Referenced by beginJob(), and CSCSkim().

◆ hxnHitChambers

TH1F* CSCSkim::hxnHitChambers
private

Definition at line 181 of file CSCSkim.h.

Referenced by beginJob(), doCSCSkimming(), and endJob().

◆ hxnRecHits

TH1F* CSCSkim::hxnRecHits
private

Definition at line 179 of file CSCSkim.h.

Referenced by beginJob(), doCSCSkimming(), and endJob().

◆ hxnRecHitsSel

TH1F* CSCSkim::hxnRecHitsSel
private

Definition at line 182 of file CSCSkim.h.

Referenced by beginJob(), doCSCSkimming(), and endJob().

◆ hxnSegments

TH1F* CSCSkim::hxnSegments
private

Definition at line 180 of file CSCSkim.h.

Referenced by beginJob(), doCSCSkimming(), and endJob().

◆ iEvent

int CSCSkim::iEvent
private

Definition at line 130 of file CSCSkim.h.

Referenced by beginJob(), and filter().

◆ iRun

int CSCSkim::iRun
private

Definition at line 129 of file CSCSkim.h.

Referenced by beginJob(), and filter().

◆ isSimulation

bool CSCSkim::isSimulation
private

Definition at line 155 of file CSCSkim.h.

◆ m_CSCGeomToken

const edm::ESGetToken<CSCGeometry, MuonGeometryRecord> CSCSkim::m_CSCGeomToken
private

Definition at line 140 of file CSCSkim.h.

Referenced by filter().

◆ makeHistograms

bool CSCSkim::makeHistograms
private

Definition at line 161 of file CSCSkim.h.

Referenced by beginJob(), CSCSkim(), doCSCSkimming(), and endJob().

◆ makeHistogramsForMessyEvents

bool CSCSkim::makeHistogramsForMessyEvents
private

Definition at line 162 of file CSCSkim.h.

Referenced by beginJob(), CSCSkim(), doMessyEventSkimming(), and endJob().

◆ mevnChambers0

TH1F* CSCSkim::mevnChambers0
private

Definition at line 184 of file CSCSkim.h.

Referenced by beginJob(), doMessyEventSkimming(), and endJob().

◆ mevnChambers1

TH1F* CSCSkim::mevnChambers1
private

Definition at line 187 of file CSCSkim.h.

Referenced by beginJob(), doMessyEventSkimming(), and endJob().

◆ mevnRecHits0

TH1F* CSCSkim::mevnRecHits0
private

Definition at line 183 of file CSCSkim.h.

Referenced by beginJob(), doMessyEventSkimming(), and endJob().

◆ mevnRecHits1

TH1F* CSCSkim::mevnRecHits1
private

Definition at line 186 of file CSCSkim.h.

Referenced by beginJob(), doMessyEventSkimming(), and endJob().

◆ mevnSegments0

TH1F* CSCSkim::mevnSegments0
private

Definition at line 185 of file CSCSkim.h.

Referenced by beginJob(), doMessyEventSkimming(), and endJob().

◆ mevnSegments1

TH1F* CSCSkim::mevnSegments1
private

Definition at line 188 of file CSCSkim.h.

Referenced by beginJob(), doMessyEventSkimming(), and endJob().

◆ minimumHitChambers

int CSCSkim::minimumHitChambers
private

Definition at line 158 of file CSCSkim.h.

Referenced by CSCSkim(), and doCSCSkimming().

◆ minimumSegments

int CSCSkim::minimumSegments
private

Definition at line 159 of file CSCSkim.h.

Referenced by CSCSkim(), and doCSCSkimming().

◆ nCSCHitsMin

int CSCSkim::nCSCHitsMin
private

Definition at line 170 of file CSCSkim.h.

Referenced by CSCSkim(), doBFieldStudySelection(), and doLongSATrack().

◆ nEventsAnalyzed

int CSCSkim::nEventsAnalyzed
private

Definition at line 117 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsCertainChamber

int CSCSkim::nEventsCertainChamber
private

Definition at line 122 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsChambersBothSides

int CSCSkim::nEventsChambersBothSides
private

Definition at line 119 of file CSCSkim.h.

Referenced by beginJob(), doCSCSkimming(), and endJob().

◆ nEventsDTOverlap

int CSCSkim::nEventsDTOverlap
private

Definition at line 123 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsForBFieldStudies

int CSCSkim::nEventsForBFieldStudies
private

Definition at line 126 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsHaloLike

int CSCSkim::nEventsHaloLike
private

Definition at line 124 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsLongSATrack

int CSCSkim::nEventsLongSATrack
private

Definition at line 125 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsMessy

int CSCSkim::nEventsMessy
private

Definition at line 121 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsOverlappingChambers

int CSCSkim::nEventsOverlappingChambers
private

Definition at line 120 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nEventsSelected

int CSCSkim::nEventsSelected
private

Definition at line 118 of file CSCSkim.h.

Referenced by beginJob(), endJob(), and filter().

◆ nLayersWithHitsMinimum

int CSCSkim::nLayersWithHitsMinimum
private

Definition at line 157 of file CSCSkim.h.

Referenced by CSCSkim(), doCSCSkimming(), and doMessyEventSkimming().

◆ nTrHitsMin

int CSCSkim::nTrHitsMin
private

Definition at line 172 of file CSCSkim.h.

Referenced by CSCSkim(), and doBFieldStudySelection().

◆ nValidHitsMin

int CSCSkim::nValidHitsMin
private

Definition at line 176 of file CSCSkim.h.

Referenced by CSCSkim(), and doBFieldStudySelection().

◆ outputFileName

std::string CSCSkim::outputFileName
private

◆ pMin

float CSCSkim::pMin
private

Definition at line 168 of file CSCSkim.h.

Referenced by CSCSkim(), and doBFieldStudySelection().

◆ redChiSqMax

float CSCSkim::redChiSqMax
private

Definition at line 175 of file CSCSkim.h.

Referenced by CSCSkim(), and doBFieldStudySelection().

◆ rExtMax

float CSCSkim::rExtMax
private

Definition at line 174 of file CSCSkim.h.

Referenced by CSCSkim(), and doBFieldStudySelection().

◆ rh_token

edm::EDGetTokenT<CSCRecHit2DCollection> CSCSkim::rh_token
private

Definition at line 148 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ sam_token

edm::EDGetTokenT<reco::TrackCollection> CSCSkim::sam_token
private

Definition at line 150 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ sdr_token

edm::EDGetTokenT<CSCStripDigiCollection> CSCSkim::sdr_token
private

Definition at line 146 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ sds_token

edm::EDGetTokenT<CSCStripDigiCollection> CSCSkim::sds_token
private

Definition at line 144 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ seg_token

edm::EDGetTokenT<CSCSegmentCollection> CSCSkim::seg_token
private

Definition at line 149 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ theHistogramFile

TFile* CSCSkim::theHistogramFile
private

Definition at line 133 of file CSCSkim.h.

Referenced by beginJob(), and endJob().

◆ trk_token

edm::EDGetTokenT<reco::TrackCollection> CSCSkim::trk_token
private

Definition at line 151 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ typeOfSkim

int CSCSkim::typeOfSkim
private

Definition at line 156 of file CSCSkim.h.

Referenced by CSCSkim(), doCSCSkimming(), endJob(), and filter().

◆ wdr_token

edm::EDGetTokenT<CSCWireDigiCollection> CSCSkim::wdr_token
private

Definition at line 145 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ wds_token

edm::EDGetTokenT<CSCWireDigiCollection> CSCSkim::wds_token
private

Definition at line 143 of file CSCSkim.h.

Referenced by CSCSkim(), and filter().

◆ whichChamber

int CSCSkim::whichChamber
private

Definition at line 166 of file CSCSkim.h.

Referenced by CSCSkim(), and doCertainChamberSelection().

◆ whichEndcap

int CSCSkim::whichEndcap
private

Definition at line 163 of file CSCSkim.h.

Referenced by CSCSkim(), and doCertainChamberSelection().

◆ whichRing

int CSCSkim::whichRing
private

Definition at line 165 of file CSCSkim.h.

Referenced by CSCSkim(), and doCertainChamberSelection().

◆ whichStation

int CSCSkim::whichStation
private

Definition at line 164 of file CSCSkim.h.

Referenced by CSCSkim(), and doCertainChamberSelection().

◆ xxnCSCHits

TH1F * CSCSkim::xxnCSCHits
private

Definition at line 190 of file CSCSkim.h.

Referenced by beginJob().

◆ xxnTrackerHits

TH1F * CSCSkim::xxnTrackerHits
private

Definition at line 190 of file CSCSkim.h.

Referenced by beginJob().

◆ xxnValidHits

TH1F * CSCSkim::xxnValidHits
private

Definition at line 190 of file CSCSkim.h.

Referenced by beginJob().

◆ xxP

TH1F* CSCSkim::xxP
private

Definition at line 190 of file CSCSkim.h.

Referenced by beginJob().

◆ xxredChiSq

TH1F * CSCSkim::xxredChiSq
private

Definition at line 190 of file CSCSkim.h.

Referenced by beginJob().

◆ zInnerMax

float CSCSkim::zInnerMax
private

Definition at line 171 of file CSCSkim.h.

Referenced by CSCSkim(), doBFieldStudySelection(), and doLongSATrack().

◆ zLengthMin

float CSCSkim::zLengthMin
private

Definition at line 169 of file CSCSkim.h.

Referenced by CSCSkim(), and doBFieldStudySelection().

◆ zLengthTrMin

float CSCSkim::zLengthTrMin
private

Definition at line 173 of file CSCSkim.h.

Referenced by CSCSkim(), and doBFieldStudySelection().