CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalTPSkimmer.cc
Go to the documentation of this file.
1 
12 
15 
17 
21 
24 
26 
28 {
29  skipModule_ = ps.getParameter<bool>("skipModule");
30 
31  doBarrel_ = ps.getParameter<bool>("doBarrel");
32  doEndcap_ = ps.getParameter<bool>("doEndcap");
33 
34  chStatusToSelectTP_ = ps.getParameter<std::vector<uint32_t> >("chStatusToSelectTP");
35 
36  tpOutputCollection_ = ps.getParameter<std::string>("tpOutputCollection");
37  tpInputCollection_ = ps.getParameter<edm::InputTag>("tpInputCollection");
38 
39  produces< EcalTrigPrimDigiCollection >(tpOutputCollection_);
40 }
41 
43 {
44 }
45 
46 void
48 {
49  insertedTP_.clear();
50 
51  using namespace edm;
52 
53  es.get<IdealGeometryRecord>().get(ttMap_);
54 
55  // collection of rechits to put in the event
56  std::auto_ptr< EcalTrigPrimDigiCollection > tpOut( new EcalTrigPrimDigiCollection );
57 
58  if ( skipModule_ ) {
59  evt.put( tpOut, tpOutputCollection_ );
60  return;
61  }
62 
64  es.get<EcalChannelStatusRcd>().get(chStatus);
65 
67  evt.getByLabel(tpInputCollection_, tpIn);
68 
69  if ( ! tpIn.isValid() ) {
70  edm::LogError("EcalTPSkimmer") << "Can't get the product " << tpInputCollection_.instance()
71  << " with label " << tpInputCollection_.label();
72  return;
73  }
74 
75  if ( doBarrel_ ) {
77  uint16_t code = 0;
78  for ( int i = 0; i < EBDetId::kSizeForDenseIndexing; ++i )
79  {
80  if ( ! EBDetId::validDenseIndex( i ) ) continue;
82  chit = chStatus->find( id );
83  // check if the channel status means TP to be kept
84  if ( chit != chStatus->end() ) {
85  code = (*chit).getStatusCode() & 0x001F;
86  if ( std::find( chStatusToSelectTP_.begin(), chStatusToSelectTP_.end(), code ) != chStatusToSelectTP_.end() ) {
87  // retrieve the TP DetId
88  EcalTrigTowerDetId ttDetId( ((EBDetId)id).tower() );
89  // insert the TP if not done already
90  if ( ! alreadyInserted( ttDetId ) ) insertTP( ttDetId, tpIn, *tpOut );
91  }
92  } else {
93  edm::LogError("EcalDetIdToBeRecoveredProducer") << "No channel status found for xtal "
94  << id.rawId()
95  << "! something wrong with EcalChannelStatus in your DB? ";
96  }
97  }
98  }
99 
100  if ( doEndcap_ ) {
102  uint16_t code = 0;
103  for ( int i = 0; i < EEDetId::kSizeForDenseIndexing; ++i )
104  {
105  if ( ! EEDetId::validDenseIndex( i ) ) continue;
107  chit = chStatus->find( id );
108  // check if the channel status means TP to be kept
109  if ( chit != chStatus->end() ) {
110  code = (*chit).getStatusCode() & 0x001F;
111  if ( std::find( chStatusToSelectTP_.begin(), chStatusToSelectTP_.end(), code ) != chStatusToSelectTP_.end() ) {
112  // retrieve the TP DetId
113  EcalTrigTowerDetId ttDetId = ttMap_->towerOf( id );
114  // insert the TP if not done already
115  if ( ! alreadyInserted( ttDetId ) ) insertTP( ttDetId, tpIn, *tpOut );
116  }
117  } else {
118  edm::LogError("EcalDetIdToBeRecoveredProducer") << "No channel status found for xtal "
119  << id.rawId()
120  << "! something wrong with EcalChannelStatus in your DB? ";
121  }
122  }
123  }
124 
125  // put the collection of reconstructed hits in the event
126  LogInfo("EcalTPSkimmer") << "total # of TP inserted: " << tpOut->size();
127 
128  evt.put( tpOut, tpOutputCollection_ );
129 }
130 
131 
133 {
134  return ( insertedTP_.find( ttId ) != insertedTP_.end() );
135 }
136 
137 
139 {
140  EcalTrigPrimDigiCollection::const_iterator tpIt = tpIn->find( ttId );
141  if ( tpIt != tpIn->end() ) {
142  tpOut.push_back( *tpIt );
143  insertedTP_.insert( ttId );
144  }
145 }
146 
147 
150 
static EEDetId detIdFromDenseIndex(uint32_t din)
Definition: EEDetId.h:215
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< uint32_t > chStatusToSelectTP_
Definition: EcalTPSkimmer.h:43
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
EcalTPSkimmer(const edm::ParameterSet &ps)
edm::InputTag tpInputCollection_
Definition: EcalTPSkimmer.h:48
std::vector< T >::const_iterator const_iterator
void push_back(T const &t)
std::set< EcalTrigTowerDetId > insertedTP_
Definition: EcalTPSkimmer.h:46
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
static EBDetId detIdFromDenseIndex(uint32_t di)
Definition: EBDetId.h:98
edm::ESHandle< EcalTrigTowerConstituentsMap > ttMap_
Definition: EcalTPSkimmer.h:44
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
std::string tpOutputCollection_
Definition: EcalTPSkimmer.h:50
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
const T & get() const
Definition: EventSetup.h:55
std::vector< Item >::const_iterator const_iterator
bool alreadyInserted(EcalTrigTowerDetId ttId)
std::string const & label() const
Definition: InputTag.h:25
void insertTP(EcalTrigTowerDetId ttId, edm::Handle< EcalTrigPrimDigiCollection > &in, EcalTrigPrimDigiCollection &out)
unsigned ttId(const DetId &)
static bool validDenseIndex(uint32_t din)
Definition: EEDetId.h:208
virtual void produce(edm::Event &evt, const edm::EventSetup &es)
std::string const & instance() const
Definition: InputTag.h:26
static bool validDenseIndex(uint32_t din)
Definition: EBDetId.h:96