CMS 3D CMS Logo

HiSpikeCleaner.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //111
3 // Package: HiSpikeCleaner
4 // Class: HiSpikeCleaner
5 //
13 //
14 // Original Author: Yong Kim,32 4-A08,+41227673039,
15 // Created: Mon Nov 1 18:22:21 CET 2010
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
28 
31 
34 
37 
42 
43 
44 
45 //
46 // class declaration
47 //
48 
50 public:
51  explicit HiSpikeCleaner(const edm::ParameterSet&);
52  ~HiSpikeCleaner() override;
53 
54 private:
55 
56  void produce(edm::Event&, const edm::EventSetup&) override;
57 
58  // ----------member data ---------------------------
59 
63 
65  double TimingCut_;
66  double swissCutThr_;
67  double etCut_;
68 };
69 
71 {
72  //register your products
73 /* Examples
74  produces<ExampleData2>();
75 
76  //if do put with a label
77  produces<ExampleData2>("label");
78 */
79  //now do what ever other initialization is needed
80 
81  rHInputProducerBToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitProducerBarrel"));
82  rHInputProducerEToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitProducerEndcap"));
83 
84  sCInputProducerToken_ = consumes<reco::SuperClusterCollection>(iConfig.getParameter<edm::InputTag>("originalSuperClusterProducer"));
85  TimingCut_ = iConfig.getUntrackedParameter<double> ("TimingCut",4.0);
86  swissCutThr_ = iConfig.getUntrackedParameter<double>("swissCutThr",0.95);
87  etCut_ = iConfig.getParameter<double>("etCut");
88 
89  outputCollection_ = iConfig.getParameter<std::string>("outputColl");
90  produces<reco::SuperClusterCollection>(outputCollection_);
91 }
92 
93 
95 {
96  // do anything here that needs to be done at desctruction time
97  // (e.g. close files, deallocate resources etc.)
98 }
99 
100 
101 //
102 // member functions
103 //
104 
105 // ------------ method called to produce the data ------------
106 void
108 {
109  using namespace edm;
110 
111  // Get raw SuperClusters from the event
112  Handle<reco::SuperClusterCollection> pRawSuperClusters;
113  try {
114  iEvent.getByToken(sCInputProducerToken_, pRawSuperClusters);
115  } catch ( cms::Exception& ex ) {
116  edm::LogError("EgammaSCCorrectionMakerError")
117  << "Error! can't get the rawSuperClusters ";
118  }
119 
120  // Get the RecHits from the event
122  try {
123  iEvent.getByToken(rHInputProducerBToken_, pRecHitsB);
124  } catch ( cms::Exception& ex ) {
125  edm::LogError("EgammaSCCorrectionMakerError")
126  << "Error! can't get the RecHits ";
127  }
128 
129  // Get the RecHits from the event
131  try {
132  iEvent.getByToken(rHInputProducerEToken_, pRecHitsE);
133  } catch ( cms::Exception& ex ) {
134  edm::LogError("EgammaSCCorrectionMakerError")
135  << "Error! can't get the RecHits ";
136  }
137 
138 
139  // get the channel status from the DB
140  // edm::ESHandle<EcalChannelStatus> chStatus;
141  // iSetup.get<EcalChannelStatusRcd>().get(chStatus);
142 
143  edm::ESHandle<EcalSeverityLevelAlgo> ecalSevLvlAlgoHndl;
144  iSetup.get<EcalSeverityLevelAlgoRcd>().get(ecalSevLvlAlgoHndl);
145 
146 
147  // Create a pointer to the RecHits and raw SuperClusters
148  const reco::SuperClusterCollection *rawClusters = pRawSuperClusters.product();
149 
150 
152 
153  // Define a collection of corrected SuperClusters to put back into the event
154  auto corrClusters = std::make_unique<reco::SuperClusterCollection>();
155 
156  // Loop over raw clusters and make corrected ones
157  reco::SuperClusterCollection::const_iterator aClus;
158  for(aClus = rawClusters->begin(); aClus != rawClusters->end(); aClus++)
159  {
160  double theEt = aClus->energy()/cosh( aClus->eta() ) ;
161  // std::cout << " et of SC = " << theEt << std::endl;
162 
163  if ( theEt < etCut_ ) continue; // cut off low pT superclusters
164 
165  bool flagS = true;
166  float swissCrx(0);
167 
168  const reco::CaloClusterPtr seed = aClus->seed();
169  DetId id = lazyTool.getMaximum(*seed).first;
170  const EcalRecHitCollection & rechits = *pRecHitsB;
172 
173  if( it != rechits.end() ) {
174  ecalSevLvlAlgoHndl->severityLevel(id, rechits);
175  swissCrx = EcalTools::swissCross (id, rechits, 0.,true);
176  // std::cout << "swissCross = " << swissCrx <<std::endl;
177  // std::cout << " timing = " << it->time() << std::endl;
178  }
179 
180  if ( fabs(it->time()) > TimingCut_ ) {
181  flagS = false;
182  // std::cout << " timing = " << it->time() << std::endl;
183  // std::cout << " timing is bad........" << std::endl;
184  }
185  if ( swissCrx > (float)swissCutThr_ ) {
186  flagS = false ; // swissCross cut
187  // std::cout << "swissCross = " << swissCrx <<std::endl;
188  // std::cout << " removed by swiss cross cut" << std::endl;
189  }
190  // - kGood --> good channel
191  // - kProblematic --> problematic (e.g. noisy)
192  // - kRecovered --> recovered (e.g. an originally dead or saturated)
193  // - kTime --> the channel is out of time (e.g. spike)
194  // - kWeird --> weird (e.g. spike)
195  // - kBad --> bad, not suitable to be used in the reconstruction
196  // enum EcalSeverityLevel { kGood=0, kProblematic, kRecovered, kTime, kWeird, kBad };
197 
198 
199  reco::SuperCluster newClus;
200  if ( flagS == true)
201  newClus=*aClus;
202  else
203  continue;
204  corrClusters->push_back(newClus);
205  }
206 
207  // Put collection of corrected SuperClusters into the event
208  iEvent.put(std::move(corrClusters), outputCollection_);
209 
210 }
211 
212 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id) const
Evaluate status from id use channelStatus from DB.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
edm::EDGetTokenT< reco::SuperClusterCollection > sCInputProducerToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< EcalRecHit >::const_iterator const_iterator
int iEvent
Definition: GenABIO.cc:230
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
std::string outputCollection_
HiSpikeCleaner(const edm::ParameterSet &)
Definition: DetId.h:18
void produce(edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:58
static float swissCross(const DetId &id, const EcalRecHitCollection &recHits, float recHitThreshold, bool avoidIeta85=true)
the good old 1-e4/e1. Ignore hits below recHitThreshold
Definition: EcalTools.cc:11
HLT enums.
edm::EDGetTokenT< EcalRecHitCollection > rHInputProducerBToken_
~HiSpikeCleaner() override
edm::EDGetTokenT< EcalRecHitCollection > rHInputProducerEToken_
def move(src, dest)
Definition: eostools.py:510