CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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&);
53 
54 private:
55  virtual void beginJob() override ;
56  virtual void produce(edm::Event&, const edm::EventSetup&) override;
57  virtual void endJob() override ;
58 
59  // ----------member data ---------------------------
60 
64 
66  double TimingCut_;
67  double swissCutThr_;
68  double etCut_;
69 };
70 
72 {
73  //register your products
74 /* Examples
75  produces<ExampleData2>();
76 
77  //if do put with a label
78  produces<ExampleData2>("label");
79 */
80  //now do what ever other initialization is needed
81 
82  rHInputProducerBToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitProducerBarrel"));
83  rHInputProducerEToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitProducerEndcap"));
84 
85  sCInputProducerToken_ = consumes<reco::SuperClusterCollection>(iConfig.getParameter<edm::InputTag>("originalSuperClusterProducer"));
86  TimingCut_ = iConfig.getUntrackedParameter<double> ("TimingCut",4.0);
87  swissCutThr_ = iConfig.getUntrackedParameter<double>("swissCutThr",0.95);
88  etCut_ = iConfig.getParameter<double>("etCut");
89 
90  outputCollection_ = iConfig.getParameter<std::string>("outputColl");
91  produces<reco::SuperClusterCollection>(outputCollection_);
92 }
93 
94 
96 {
97  // do anything here that needs to be done at desctruction time
98  // (e.g. close files, deallocate resources etc.)
99 }
100 
101 
102 //
103 // member functions
104 //
105 
106 // ------------ method called to produce the data ------------
107 void
109 {
110  using namespace edm;
111 
112  // Get raw SuperClusters from the event
113  Handle<reco::SuperClusterCollection> pRawSuperClusters;
114  try {
115  iEvent.getByToken(sCInputProducerToken_, pRawSuperClusters);
116  } catch ( cms::Exception& ex ) {
117  edm::LogError("EgammaSCCorrectionMakerError")
118  << "Error! can't get the rawSuperClusters ";
119  }
120 
121  // Get the RecHits from the event
123  try {
124  iEvent.getByToken(rHInputProducerBToken_, pRecHitsB);
125  } catch ( cms::Exception& ex ) {
126  edm::LogError("EgammaSCCorrectionMakerError")
127  << "Error! can't get the RecHits ";
128  }
129 
130  // Get the RecHits from the event
132  try {
133  iEvent.getByToken(rHInputProducerEToken_, pRecHitsE);
134  } catch ( cms::Exception& ex ) {
135  edm::LogError("EgammaSCCorrectionMakerError")
136  << "Error! can't get the RecHits ";
137  }
138 
139 
140  // get the channel status from the DB
141  // edm::ESHandle<EcalChannelStatus> chStatus;
142  // iSetup.get<EcalChannelStatusRcd>().get(chStatus);
143 
144  edm::ESHandle<EcalSeverityLevelAlgo> ecalSevLvlAlgoHndl;
145  iSetup.get<EcalSeverityLevelAlgoRcd>().get(ecalSevLvlAlgoHndl);
146 
147 
148  // Create a pointer to the RecHits and raw SuperClusters
149  const reco::SuperClusterCollection *rawClusters = pRawSuperClusters.product();
150 
151 
153 
154  // Define a collection of corrected SuperClusters to put back into the event
155  std::auto_ptr<reco::SuperClusterCollection> corrClusters(new reco::SuperClusterCollection);
156 
157  // Loop over raw clusters and make corrected ones
158  reco::SuperClusterCollection::const_iterator aClus;
159  for(aClus = rawClusters->begin(); aClus != rawClusters->end(); aClus++)
160  {
161  double theEt = aClus->energy()/cosh( aClus->eta() ) ;
162  // std::cout << " et of SC = " << theEt << std::endl;
163 
164  if ( theEt < etCut_ ) continue; // cut off low pT superclusters
165 
166  bool flagS = true;
167  float swissCrx(0);
168 
169  const reco::CaloClusterPtr seed = aClus->seed();
170  DetId id = lazyTool.getMaximum(*seed).first;
171  const EcalRecHitCollection & rechits = *pRecHitsB;
173 
174  if( it != rechits.end() ) {
175  ecalSevLvlAlgoHndl->severityLevel(id, rechits);
176  swissCrx = EcalTools::swissCross (id, rechits, 0.,true);
177  // std::cout << "swissCross = " << swissCrx <<std::endl;
178  // std::cout << " timing = " << it->time() << std::endl;
179  }
180 
181  if ( fabs(it->time()) > TimingCut_ ) {
182  flagS = false;
183  // std::cout << " timing = " << it->time() << std::endl;
184  // std::cout << " timing is bad........" << std::endl;
185  }
186  if ( swissCrx > (float)swissCutThr_ ) {
187  flagS = false ; // swissCross cut
188  // std::cout << "swissCross = " << swissCrx <<std::endl;
189  // std::cout << " removed by swiss cross cut" << std::endl;
190  }
191  // - kGood --> good channel
192  // - kProblematic --> problematic (e.g. noisy)
193  // - kRecovered --> recovered (e.g. an originally dead or saturated)
194  // - kTime --> the channel is out of time (e.g. spike)
195  // - kWeird --> weird (e.g. spike)
196  // - kBad --> bad, not suitable to be used in the reconstruction
197  // enum EcalSeverityLevel { kGood=0, kProblematic, kRecovered, kTime, kWeird, kBad };
198 
199 
200  reco::SuperCluster newClus;
201  if ( flagS == true)
202  newClus=*aClus;
203  else
204  continue;
205  corrClusters->push_back(newClus);
206  }
207 
208  // Put collection of corrected SuperClusters into the event
209  iEvent.put(corrClusters, outputCollection_);
210 
211 }
212 
213 // ------------ method called once each job just before starting event loop ------------
214 void
216 {
217 }
218 
219 // ------------ method called once each job just after ending the event loop ------------
220 void
222 }
223 
224 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
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
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
std::string outputCollection_
HiSpikeCleaner(const edm::ParameterSet &)
Definition: DetId.h:18
virtual void produce(edm::Event &, const edm::EventSetup &) override
virtual void beginJob() override
const T & get() const
Definition: EventSetup.h:55
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
virtual void endJob() override
edm::EDGetTokenT< EcalRecHitCollection > rHInputProducerBToken_
volatile std::atomic< bool > shutdown_flag false
edm::EDGetTokenT< EcalRecHitCollection > rHInputProducerEToken_