CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GamIsoDetIdCollectionProducer.cc
Go to the documentation of this file.
2 
6 
8 
19 
22 
24 
27 
28 
30  recHitsLabel_(iConfig.getParameter< edm::InputTag > ("recHitsLabel")),
31  emObjectLabel_(iConfig.getParameter< edm::InputTag > ("emObjectLabel")),
32  energyCut_(iConfig.getParameter<double>("energyCut")),
33  etCut_(iConfig.getParameter<double>("etCut")),
34  etCandCut_(iConfig.getParameter<double> ("etCandCut")),
35  outerRadius_(iConfig.getParameter<double>("outerRadius")),
36  innerRadius_(iConfig.getParameter<double>("innerRadius")),
37  interestingDetIdCollection_(iConfig.getParameter<std::string>("interestingDetIdCollection")),
38  severityLevelCut_(iConfig.getParameter<int>("severityLevelCut")),
39  severityRecHitThreshold_(iConfig.getParameter<double>("severityRecHitThreshold")),
40  spIdString_(iConfig.getParameter<std::string>("spikeIdString")),
41  spIdThreshold_(iConfig.getParameter<double>("spikeIdThreshold")),
42  v_chstatus_(iConfig.getParameter<std::vector<int> >("recHitFlagsToBeExcluded")) {
43 
44  if ( !spIdString_.compare("kE1OverE9") ) spId_ = EcalSeverityLevelAlgo::kE1OverE9;
45  else if( !spIdString_.compare("kSwissCross") ) spId_ = EcalSeverityLevelAlgo::kSwissCross;
46  else if( !spIdString_.compare("kSwissCrossBordersIncluded") ) spId_ = EcalSeverityLevelAlgo::kSwissCrossBordersIncluded;
48 
49  //register your products
50  produces< DetIdCollection > (interestingDetIdCollection_) ;
51 }
52 
54 
56 {}
57 
58 // ------------ method called to produce the data ------------
59  void
61  const edm::EventSetup& iSetup)
62 {
63  using namespace edm;
64  using namespace std;
65 
66  //Get EM Object
68  iEvent.getByLabel(emObjectLabel_,emObjectH);
69 
70  // take EcalRecHits
72  iEvent.getByLabel(recHitsLabel_,recHitsH);
73  std::auto_ptr<CaloRecHitMetaCollectionV> recHits_(0);
74  recHits_ = std::auto_ptr<CaloRecHitMetaCollectionV>(new EcalRecHitMetaCollection(*recHitsH));
75 
77  iSetup.get<CaloGeometryRecord>().get(pG);
78  const CaloGeometry* caloGeom = pG.product();
79 
80  //Get the channel status from the db
82  iSetup.get<EcalChannelStatusRcd>().get(chStatus);
83 
84  CaloDualConeSelector *doubleConeSel_ = 0;
85  if(recHitsLabel_.instance() == "EcalRecHitsEB")
87  else if(recHitsLabel_.instance() == "EcalRecHitsEE")
89 
90  //Create empty output collections
91  std::auto_ptr< DetIdCollection > detIdCollection (new DetIdCollection() ) ;
92 
93  reco::PhotonCollection::const_iterator emObj;
94  if(doubleConeSel_) { //if cone selector was created
95  for (emObj = emObjectH->begin(); emObj != emObjectH->end(); emObj++) { //Loop over candidates
96 
97  if(emObj->et() < etCandCut_) continue;
98 
99  GlobalPoint pclu (emObj->caloPosition().x(),emObj->caloPosition().y(),emObj->caloPosition().z());
100  std::auto_ptr<CaloRecHitMetaCollectionV> chosen = doubleConeSel_->select(pclu,*recHits_);
101 
103  for (recIt = chosen->begin(); recIt!= chosen->end () ; ++recIt) { // Select RecHits
104 
105  if ( fabs(recIt->energy()) < energyCut_) continue; //dont fill if below E noise value
106 
107  double et = recIt->energy() *
108  caloGeom->getPosition(recIt->detid()).perp() /
109  caloGeom->getPosition(recIt->detid()).mag();
110 
111  if ( fabs(et) < etCut_) continue; //dont fill if below ET noise value
112 
113  if(recHitsLabel_.instance() == "EcalRecHitsEB" && //make sure we have a barrel rechit
114  EcalSeverityLevelAlgo::severityLevel( //call the severity level method
115  EBDetId(recIt->detid()), //passing the EBDetId
116  *recHitsH, //the rechit collection in order to calculate the swiss crss
117  *chStatus, //and the EcalChannelRecHitRcd
118  severityRecHitThreshold_, //only consider rechits with ET >
119  spId_, //the SpikeId method (currently kE1OverE9 or kSwissCross)
120  spIdThreshold_ //cut value for above
121  ) >= severityLevelCut_) continue; //then if the severity level is too high, we continue to the next rechit
122 
123  //Check based on flags to protect from recovered channels from non-read towers
124  //Assumption is that v_chstatus_ is empty unless doFlagChecks() has been called
125  std::vector<int>::const_iterator vit = std::find( v_chstatus_.begin(), v_chstatus_.end(), ((EcalRecHit*)(&*recIt))->recoFlag() );
126  if ( vit != v_chstatus_.end() ) continue; // the recHit has to be excluded from the iso sum
127 
128  if(std::find(detIdCollection->begin(),detIdCollection->end(),recIt->detid()) == detIdCollection->end())
129  detIdCollection->push_back(recIt->detid());
130  } //end rechits
131 
132  } //end candidates
133 
134  delete doubleConeSel_;
135  } //end if cone selector was created
136 
137  iEvent.put( detIdCollection, interestingDetIdCollection_ );
138 }
const DetId & detid() const
Definition: CaloRecHit.h:21
T perp() const
Magnitude of transverse component.
EcalSeverityLevelAlgo::SpikeId spId_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
int iEvent
Definition: GenABIO.cc:243
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
float energy() const
Definition: CaloRecHit.h:19
virtual void produce(edm::Event &, const edm::EventSetup &)
producer
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
GamIsoDetIdCollectionProducer(const edm::ParameterSet &)
ctor
const T & get() const
Definition: EventSetup.h:55
SimpleCaloRecHitMetaCollection< EcalRecHitCollection > EcalRecHitMetaCollection
static int severityLevel(const DetId, const EcalRecHitCollection &, const EcalChannelStatus &, float recHitEtThreshold=5., SpikeId spId=kSwissCross, float spIdThreshold=0.95, float recHitEnergyThresholdForTiming=2., float recHitEnergyThresholdForEE=15, float spIdThresholdIEta85=0.999)
edm::EDCollection< DetId > DetIdCollection
std::string const & instance() const
Definition: InputTag.h:26