CMS 3D CMS Logo

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_ = std::make_unique<EgammaSCEnergyCorrectionAlgo>(sigmaElectronicNoise_, sCAlgo_, fCorrPset);
81 
82  // energy correction class
84  energyCorrectionFunction_ = std::unique_ptr<EcalClusterFunctionBaseClass>(EcalClusterFunctionFactory::get()->create(energyCorrectorName_, ps));
85  //energyCorrectionFunction_ = EcalClusterFunctionFactory::get()->create("EcalClusterEnergyCorrection", ps);
86 
87  if (applyCrackCorrection_ )
88  crackCorrectionFunction_ = std::unique_ptr<EcalClusterFunctionBaseClass>(EcalClusterFunctionFactory::get()->create(crackCorrectorName_, ps));
89 
90 
91  if (applyLocalContCorrection_ )
92  localContCorrectionFunction_ = std::unique_ptr<EcalClusterFunctionBaseClass>(EcalClusterFunctionFactory::get()->create(localContCorrectorName_, ps));
93 
94 }
95 
97 
98 void
100 {
101  using namespace edm;
102 
103  // initialize energy correction class
105  energyCorrectionFunction_->init(es);
106 
107  // initialize energy correction class
109  crackCorrectionFunction_->init(es);
110 
111 
112  // initialize containemnt correction class
115 
116  // get the collection geometry:
117  edm::ESHandle<CaloGeometry> geoHandle;
118  es.get<CaloGeometryRecord>().get(geoHandle);
119  const CaloGeometry& geometry = *geoHandle;
120  const CaloSubdetectorGeometry *geometry_p;
121 
122  std::string rHInputCollection = rHTag_.instance();
123  if(rHInputCollection == "EcalRecHitsEB") {
124  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
125  } else if(rHInputCollection == "EcalRecHitsEE") {
126  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalEndcap);
127  } else if(rHInputCollection == "EcalRecHitsPS") {
128  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
129  } else {
130  std::string str = "\n\nSCCorrectionMaker encountered invalied ecalhitcollection type: " + rHInputCollection + ".\n\n";
131  throw(std::runtime_error( str.c_str() ));
132  }
133 
134  // Get raw SuperClusters from the event
135  Handle<reco::SuperClusterCollection> pRawSuperClusters;
136  evt.getByToken(sCInputProducer_, pRawSuperClusters);
137 
138  // Get the RecHits from the event
140  evt.getByToken(rHInputProducer_, pRecHits);
141 
142  // Create a pointer to the RecHits and raw SuperClusters
143  const EcalRecHitCollection *hitCollection = pRecHits.product();
144  const reco::SuperClusterCollection *rawClusters = pRawSuperClusters.product();
145 
146  // Define a collection of corrected SuperClusters to put back into the event
147  auto corrClusters = std::make_unique<reco::SuperClusterCollection>();
148 
149  // Loop over raw clusters and make corrected ones
150  reco::SuperClusterCollection::const_iterator aClus;
151  int i=0;
152  for(aClus = rawClusters->begin(); aClus != rawClusters->end(); aClus++)
153  {
154  reco::SuperCluster enecorrClus,crackcorrClus,localContCorrClus;
155 
156  i++;
157 
159  enecorrClus = energyCorrector_->applyCorrection(*aClus, *hitCollection, sCAlgo_, geometry_p, energyCorrectionFunction_.get(), energyCorrectorName_, modeEB_, modeEE_);
160  else
161  enecorrClus=*aClus;
162 
163 
165  crackcorrClus=energyCorrector_->applyCrackCorrection(enecorrClus,crackCorrectionFunction_.get());
166  else
167  crackcorrClus=enecorrClus;
168 
170  localContCorrClus =
171  energyCorrector_->applyLocalContCorrection(crackcorrClus,localContCorrectionFunction_.get());
172  else
173  localContCorrClus = crackcorrClus;
174 
175 
176  if(localContCorrClus.energy()*sin(localContCorrClus.position().theta())>etThresh_) {
177 
178  corrClusters->push_back(localContCorrClus);
179  }
180  }
181 
182  // Put collection of corrected SuperClusters into the event
183  evt.put(std::move(corrClusters), outputCollection_);
184 
185 }
186 
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:131
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:161
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
~EgammaSCCorrectionMaker() override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
edm::EDGetTokenT< EcalRecHitCollection > rHInputProducer_
EgammaSCCorrectionMaker(const edm::ParameterSet &)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
reco::CaloCluster::AlgoId sCAlgo_
edm::EDGetTokenT< reco::SuperClusterCollection > sCInputProducer_
double energy() const
cluster energy
Definition: CaloCluster.h:126
void produce(edm::Event &, const edm::EventSetup &) override
std::unique_ptr< EcalClusterFunctionBaseClass > energyCorrectionFunction_
T const * product() const
Definition: Handle.h:74
std::unique_ptr< EcalClusterFunctionBaseClass > crackCorrectionFunction_
std::unique_ptr< EcalClusterFunctionBaseClass > localContCorrectionFunction_
HLT enums.
T get() const
Definition: EventSetup.h:71
#define str(s)
std::string const & instance() const
Definition: InputTag.h:37
std::unique_ptr< EgammaSCEnergyCorrectionAlgo > energyCorrector_
def move(src, dest)
Definition: eostools.py:511
T get(const Candidate &c)
Definition: component.h:55