CMS 3D CMS Logo

EcalTPSkimmer.cc
Go to the documentation of this file.
1 
9 
12 
14 
18 
21 
23 
25 {
26  skipModule_ = ps.getParameter<bool>("skipModule");
27 
28  doBarrel_ = ps.getParameter<bool>("doBarrel");
29  doEndcap_ = ps.getParameter<bool>("doEndcap");
30 
31  chStatusToSelectTP_ = ps.getParameter<std::vector<uint32_t> >("chStatusToSelectTP");
32 
33  tpOutputCollection_ = ps.getParameter<std::string>("tpOutputCollection");
34  tpInputToken_ = consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("tpInputCollection"));
35 
36  produces< EcalTrigPrimDigiCollection >(tpOutputCollection_);
37 }
38 
40 {
41 }
42 
43 void
45 {
46  insertedTP_.clear();
47 
48  using namespace edm;
49 
50  es.get<IdealGeometryRecord>().get(ttMap_);
51 
52  // collection of rechits to put in the event
53  auto tpOut = std::make_unique<EcalTrigPrimDigiCollection>();
54 
55  if ( skipModule_ ) {
56  evt.put(std::move(tpOut), tpOutputCollection_);
57  return;
58  }
59 
61  es.get<EcalChannelStatusRcd>().get(chStatus);
62 
64  evt.getByToken(tpInputToken_, tpIn);
65 
66 
67  if ( doBarrel_ ) {
69  uint16_t code = 0;
70  for ( int i = 0; i < EBDetId::kSizeForDenseIndexing; ++i )
71  {
72  if ( ! EBDetId::validDenseIndex( i ) ) continue;
74  chit = chStatus->find( id );
75  // check if the channel status means TP to be kept
76  if ( chit != chStatus->end() ) {
77  code = (*chit).getStatusCode();
78  if ( std::find( chStatusToSelectTP_.begin(), chStatusToSelectTP_.end(), code ) != chStatusToSelectTP_.end() ) {
79  // retrieve the TP DetId
80  EcalTrigTowerDetId ttDetId( ((EBDetId)id).tower() );
81  // insert the TP if not done already
82  if ( ! alreadyInserted( ttDetId ) ) insertTP( ttDetId, tpIn, *tpOut );
83  }
84  } else {
85  edm::LogError("EcalDetIdToBeRecoveredProducer") << "No channel status found for xtal "
86  << id.rawId()
87  << "! something wrong with EcalChannelStatus in your DB? ";
88  }
89  }
90  }
91 
92  if ( doEndcap_ ) {
94  uint16_t code = 0;
95  for ( int i = 0; i < EEDetId::kSizeForDenseIndexing; ++i )
96  {
97  if ( ! EEDetId::validDenseIndex( i ) ) continue;
99  chit = chStatus->find( id );
100  // check if the channel status means TP to be kept
101  if ( chit != chStatus->end() ) {
102  code = (*chit).getStatusCode() ;
103  if ( std::find( chStatusToSelectTP_.begin(), chStatusToSelectTP_.end(), code ) != chStatusToSelectTP_.end() ) {
104  // retrieve the TP DetId
105  EcalTrigTowerDetId ttDetId = ttMap_->towerOf( id );
106  // insert the TP if not done already
107  if ( ! alreadyInserted( ttDetId ) ) insertTP( ttDetId, tpIn, *tpOut );
108  }
109  } else {
110  edm::LogError("EcalDetIdToBeRecoveredProducer") << "No channel status found for xtal "
111  << id.rawId()
112  << "! something wrong with EcalChannelStatus in your DB? ";
113  }
114  }
115  }
116 
117  // put the collection of reconstructed hits in the event
118  LogInfo("EcalTPSkimmer") << "total # of TP inserted: " << tpOut->size();
119 
120  evt.put(std::move(tpOut), tpOutputCollection_);
121 }
122 
123 
125 {
126  return ( insertedTP_.find( ttId ) != insertedTP_.end() );
127 }
128 
129 
131 {
133  if ( tpIt != tpIn->end() ) {
134  tpOut.push_back( *tpIt );
135  insertedTP_.insert( ttId );
136  }
137 }
138 
139 
142 
static EEDetId detIdFromDenseIndex(uint32_t din)
Definition: EEDetId.h:220
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
std::vector< uint32_t > chStatusToSelectTP_
Definition: EcalTPSkimmer.h:40
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
EcalTPSkimmer(const edm::ParameterSet &ps)
std::vector< EcalTriggerPrimitiveDigi >::const_iterator const_iterator
void push_back(T const &t)
std::set< EcalTrigTowerDetId > insertedTP_
Definition: EcalTPSkimmer.h:43
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
unsigned ttId(DetId const &)
EcalTrigTowerDetId towerOf(const DetId &id) const
Get the tower id for this det id (or null if not known)
static EBDetId detIdFromDenseIndex(uint32_t di)
Definition: EBDetId.h:111
edm::ESHandle< EcalTrigTowerConstituentsMap > ttMap_
Definition: EcalTPSkimmer.h:41
std::string tpOutputCollection_
Definition: EcalTPSkimmer.h:47
const_iterator end() const
const T & get() const
Definition: EventSetup.h:56
std::vector< Item >::const_iterator const_iterator
bool alreadyInserted(EcalTrigTowerDetId ttId)
void insertTP(EcalTrigTowerDetId ttId, edm::Handle< EcalTrigPrimDigiCollection > &in, EcalTrigPrimDigiCollection &out)
static bool validDenseIndex(uint32_t din)
Definition: EEDetId.h:213
iterator find(key_type k)
HLT enums.
edm::EDGetTokenT< EcalTrigPrimDigiCollection > tpInputToken_
Definition: EcalTPSkimmer.h:45
const_iterator find(uint32_t rawId) const
const_iterator end() const
virtual void produce(edm::Event &evt, const edm::EventSetup &es)
def move(src, dest)
Definition: eostools.py:510
static bool validDenseIndex(uint32_t din)
Definition: EBDetId.h:109