CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaSCCorrectionMaker.cc
Go to the documentation of this file.
4 
8 
16 
19 
20 #include <string>
21 
23 {
24 
25 
26  // the input producers
27  rHTag_ = ps.getParameter<edm::InputTag>("recHitProducer");
28  rHInputProducer_ = consumes<EcalRecHitCollection>(rHTag_);
29  sCInputProducer_ = consumes<reco::SuperClusterCollection>(ps.getParameter<edm::InputTag>("rawSuperClusterProducer"));
30  std::string sCAlgo_str = ps.getParameter<std::string>("superClusterAlgo");
31 
32  // determine which BasicCluster algo we are correcting for
33  //And obtain forrection parameters form cfg file
34  edm::ParameterSet fCorrPset;
35  if (sCAlgo_str=="Hybrid") {
37  fCorrPset = ps.getParameter<edm::ParameterSet>("hyb_fCorrPset");
38  } else if (sCAlgo_str=="Island") {
40  fCorrPset = ps.getParameter<edm::ParameterSet>("isl_fCorrPset");
41  } else if (sCAlgo_str=="DynamicHybrid") {
43  fCorrPset = ps.getParameter<edm::ParameterSet>("dyn_fCorrPset");
44  } else if (sCAlgo_str=="Multi5x5") {
46  fCorrPset = ps.getParameter<edm::ParameterSet>("fix_fCorrPset");
47  } else {
48  edm::LogError("EgammaSCCorrectionMakerError")
49  << "Error! SuperClusterAlgo in config file must be Hybrid or Island: "
50  << sCAlgo_str << " Using Hybrid by default";
52  }
53 
54  // set correction algo parameters
55  applyEnergyCorrection_ = ps.getParameter<bool>("applyEnergyCorrection");
56  applyCrackCorrection_ = ps.getParameter<bool>("applyCrackCorrection");
57  applyLocalContCorrection_= ps.existsAs<bool>("applyLocalContCorrection") ?
58  ps.getParameter<bool>("applyLocalContCorrection") : false;
59 
60  energyCorrectorName_ = ps.getParameter<std::string>("energyCorrectorName");
61  crackCorrectorName_ = ps.existsAs<std::string>("crackCorrectorName") ?
62  ps.getParameter<std::string>("crackCorrectorName") : std::string("EcalClusterCrackCorrection");
63  localContCorrectorName_= ps.existsAs<std::string>("localContCorrectorName") ?
64  ps.getParameter<std::string>("localContCorrectorName") : std::string("EcalBasicClusterLocalContCorrection") ;
65 
66  modeEB_ = ps.getParameter<int>("modeEB");
67  modeEE_ = ps.getParameter<int>("modeEE");
68 
69  sigmaElectronicNoise_ = ps.getParameter<double>("sigmaElectronicNoise");
70 
71  etThresh_ = ps.getParameter<double>("etThresh");
72 
73 
74 
75  // set the producer parameters
76  outputCollection_ = ps.getParameter<std::string>("corectedSuperClusterCollection");
77  produces<reco::SuperClusterCollection>(outputCollection_);
78 
79  // instanciate the correction algo object
80  energyCorrector_ = new EgammaSCEnergyCorrectionAlgo(sigmaElectronicNoise_, sCAlgo_, fCorrPset);
81 
82  // energy correction class
85  //energyCorrectionFunction_ = EcalClusterFunctionFactory::get()->create("EcalClusterEnergyCorrection", ps);
86  else
88 
89  if (applyCrackCorrection_ )
90  crackCorrectionFunction_ = EcalClusterFunctionFactory::get()->create(crackCorrectorName_, ps);
91  else
93 
94 
95  if (applyLocalContCorrection_ )
97  else
99 
100 }
101 
103 {
107 
108  delete energyCorrector_;
109 }
110 
111 void
113 {
114  using namespace edm;
115 
116  // initialize energy correction class
119 
120  // initialize energy correction class
123 
124 
125  // initialize containemnt correction class
128 
129  // get the collection geometry:
130  edm::ESHandle<CaloGeometry> geoHandle;
131  es.get<CaloGeometryRecord>().get(geoHandle);
132  const CaloGeometry& geometry = *geoHandle;
133  const CaloSubdetectorGeometry *geometry_p;
134 
135  std::string rHInputCollection = rHTag_.instance();
136  if(rHInputCollection == "EcalRecHitsEB") {
137  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
138  } else if(rHInputCollection == "EcalRecHitsEE") {
139  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalEndcap);
140  } else if(rHInputCollection == "EcalRecHitsPS") {
141  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
142  } else {
143  std::string str = "\n\nSCCorrectionMaker encountered invalied ecalhitcollection type: " + rHInputCollection + ".\n\n";
144  throw(std::runtime_error( str.c_str() ));
145  }
146 
147  // Get raw SuperClusters from the event
148  Handle<reco::SuperClusterCollection> pRawSuperClusters;
149  evt.getByToken(sCInputProducer_, pRawSuperClusters);
150 
151  // Get the RecHits from the event
153  evt.getByToken(rHInputProducer_, pRecHits);
154 
155  // Create a pointer to the RecHits and raw SuperClusters
156  const EcalRecHitCollection *hitCollection = pRecHits.product();
157  const reco::SuperClusterCollection *rawClusters = pRawSuperClusters.product();
158 
159  // Define a collection of corrected SuperClusters to put back into the event
160  std::auto_ptr<reco::SuperClusterCollection> corrClusters(new reco::SuperClusterCollection);
161 
162  // Loop over raw clusters and make corrected ones
163  reco::SuperClusterCollection::const_iterator aClus;
164  int i=0;
165  for(aClus = rawClusters->begin(); aClus != rawClusters->end(); aClus++)
166  {
167  reco::SuperCluster enecorrClus,crackcorrClus,localContCorrClus;
168 
169  i++;
170 
173  else
174  enecorrClus=*aClus;
175 
176 
179  else
180  crackcorrClus=enecorrClus;
181 
183  localContCorrClus =
185  else
186  localContCorrClus = crackcorrClus;
187 
188 
189  if(localContCorrClus.energy()*sin(localContCorrClus.position().theta())>etThresh_) {
190 
191  corrClusters->push_back(localContCorrClus);
192  }
193  }
194 
195  // Put collection of corrected SuperClusters into the event
196  evt.put(corrClusters, outputCollection_);
197 
198 }
199 
virtual void produce(edm::Event &, const edm::EventSetup &)
T getParameter(std::string const &) const
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:126
int i
Definition: DBlmapReader.cc:9
reco::SuperCluster applyCorrection(const reco::SuperCluster &cl, const EcalRecHitCollection &rhc, reco::CaloCluster::AlgoId theAlgo, const CaloSubdetectorGeometry *geometry, EcalClusterFunctionBaseClass *energyCorrectionFunction, std::string energyCorrectorName_, int modeEB_, int modeEE_)
EcalClusterFunctionBaseClass * localContCorrectionFunction_
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
EcalClusterFunctionBaseClass * energyCorrectionFunction_
edm::EDGetTokenT< EcalRecHitCollection > rHInputProducer_
EgammaSCCorrectionMaker(const edm::ParameterSet &)
EcalClusterFunctionBaseClass * crackCorrectionFunction_
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
reco::CaloCluster::AlgoId sCAlgo_
edm::EDGetTokenT< reco::SuperClusterCollection > sCInputProducer_
double energy() const
cluster energy
Definition: CaloCluster.h:121
reco::SuperCluster applyCrackCorrection(const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
reco::SuperCluster applyLocalContCorrection(const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *localContCorrectionFunction)
EgammaSCEnergyCorrectionAlgo * energyCorrector_
const T & get() const
Definition: EventSetup.h:56
ESHandle< TrackerGeometry > geometry
virtual void init(const edm::EventSetup &es)=0
std::string const & instance() const
Definition: InputTag.h:44
T get(const Candidate &c)
Definition: component.h:55