CMS 3D CMS Logo

EcalTPSkimmer.cc
Go to the documentation of this file.
1 
9 
12 
14 
18 
20  skipModule_ = ps.getParameter<bool>("skipModule");
21 
22  doBarrel_ = ps.getParameter<bool>("doBarrel");
23  doEndcap_ = ps.getParameter<bool>("doEndcap");
24 
25  chStatusToSelectTP_ = ps.getParameter<std::vector<uint32_t> >("chStatusToSelectTP");
26 
27  tpOutputCollection_ = ps.getParameter<std::string>("tpOutputCollection");
28  tpInputToken_ = consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("tpInputCollection"));
29  ttMapToken_ = esConsumes<EcalTrigTowerConstituentsMap, IdealGeometryRecord>();
30  if (not skipModule_) {
31  chStatusToken_ = esConsumes<EcalChannelStatus, EcalChannelStatusRcd>();
32  }
33  produces<EcalTrigPrimDigiCollection>(tpOutputCollection_);
34 }
35 
37 
39  insertedTP_.clear();
40 
41  using namespace edm;
42 
44 
45  // collection of rechits to put in the event
46  auto tpOut = std::make_unique<EcalTrigPrimDigiCollection>();
47 
48  if (skipModule_) {
49  evt.put(std::move(tpOut), tpOutputCollection_);
50  return;
51  }
52 
54 
56  evt.getByToken(tpInputToken_, tpIn);
57 
58  if (doBarrel_) {
60  uint16_t code = 0;
61  for (int i = 0; i < EBDetId::kSizeForDenseIndexing; ++i) {
63  continue;
65  chit = chStatus->find(id);
66  // check if the channel status means TP to be kept
67  if (chit != chStatus->end()) {
68  code = (*chit).getStatusCode();
69  if (std::find(chStatusToSelectTP_.begin(), chStatusToSelectTP_.end(), code) != chStatusToSelectTP_.end()) {
70  // retrieve the TP DetId
71  EcalTrigTowerDetId ttDetId(((EBDetId)id).tower());
72  // insert the TP if not done already
73  if (!alreadyInserted(ttDetId))
74  insertTP(ttDetId, tpIn, *tpOut);
75  }
76  } else {
77  edm::LogError("EcalDetIdToBeRecoveredProducer") << "No channel status found for xtal " << id.rawId()
78  << "! something wrong with EcalChannelStatus in your DB? ";
79  }
80  }
81  }
82 
83  if (doEndcap_) {
85  uint16_t code = 0;
86  for (int i = 0; i < EEDetId::kSizeForDenseIndexing; ++i) {
88  continue;
90  chit = chStatus->find(id);
91  // check if the channel status means TP to be kept
92  if (chit != chStatus->end()) {
93  code = (*chit).getStatusCode();
94  if (std::find(chStatusToSelectTP_.begin(), chStatusToSelectTP_.end(), code) != chStatusToSelectTP_.end()) {
95  // retrieve the TP DetId
96  EcalTrigTowerDetId ttDetId = ttMap_->towerOf(id);
97  // insert the TP if not done already
98  if (!alreadyInserted(ttDetId))
99  insertTP(ttDetId, tpIn, *tpOut);
100  }
101  } else {
102  edm::LogError("EcalDetIdToBeRecoveredProducer") << "No channel status found for xtal " << id.rawId()
103  << "! something wrong with EcalChannelStatus in your DB? ";
104  }
105  }
106  }
107 
108  // put the collection of reconstructed hits in the event
109  LogInfo("EcalTPSkimmer") << "total # of TP inserted: " << tpOut->size();
110 
111  evt.put(std::move(tpOut), tpOutputCollection_);
112 }
113 
115 
120  if (tpIt != tpIn->end()) {
121  tpOut.push_back(*tpIt);
122  insertedTP_.insert(ttId);
123  }
124 }
125 
static EEDetId detIdFromDenseIndex(uint32_t din)
Definition: EEDetId.h:220
edm::ESGetToken< EcalChannelStatus, EcalChannelStatusRcd > chStatusToken_
Definition: EcalTPSkimmer.h:45
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::vector< uint32_t > chStatusToSelectTP_
Definition: EcalTPSkimmer.h:42
EcalTPSkimmer(const edm::ParameterSet &ps)
void produce(edm::Event &evt, const edm::EventSetup &es) override
std::vector< T >::const_iterator const_iterator
void push_back(T const &t)
std::set< EcalTrigTowerDetId > insertedTP_
Definition: EcalTPSkimmer.h:47
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
unsigned ttId(DetId const &, EcalElectronicsMapping 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:107
edm::ESGetToken< EcalTrigTowerConstituentsMap, IdealGeometryRecord > ttMapToken_
Definition: EcalTPSkimmer.h:44
edm::ESHandle< EcalTrigTowerConstituentsMap > ttMap_
Definition: EcalTPSkimmer.h:43
std::string tpOutputCollection_
Definition: EcalTPSkimmer.h:51
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
const_iterator find(uint32_t rawId) const
const_iterator end() const
Log< level::Info, false > LogInfo
std::vector< Item >::const_iterator const_iterator
~EcalTPSkimmer() override
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:49
const_iterator end() const
def move(src, dest)
Definition: eostools.py:511
static bool validDenseIndex(uint32_t din)
Definition: EBDetId.h:105