CMS 3D CMS Logo

SiStripApprox2ApproxClusters.cc
Go to the documentation of this file.
1 
2 
14 
15 #include <vector>
16 #include <memory>
17 
19 public:
20  explicit SiStripApprox2ApproxClusters(const edm::ParameterSet& conf);
21  void produce(edm::Event&, const edm::EventSetup&) override;
22 
23  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24 
25 private:
27  uint8_t approxVersion;
30 };
31 
33  inputApproxClusters = conf.getParameter<edm::InputTag>("inputApproxClusters");
34  approxVersionS = conf.getParameter<std::string>("approxVersion");
35 
36  approxVersion = -1;
37 
38  if (approxVersionS == "ORIGINAL")
39  approxVersion = 0;
40  else if (approxVersionS == "FULL_WIDTH")
41  approxVersion = 1;
42  else if (approxVersionS == "BARY_RES_0.1")
43  approxVersion = 2;
44  else if (approxVersionS == "BARY_CHARGE_RES_0.1")
45  approxVersion = 3;
46 
47  clusterToken = consumes<edmNew::DetSetVector<SiStripApproximateCluster>>(inputApproxClusters);
48  produces<edmNew::DetSetVector<SiStripApproximateCluster>>();
49 }
50 
52  auto result = std::make_unique<edmNew::DetSetVector<SiStripApproximateCluster>>();
53  const auto& clusterCollection = event.get(clusterToken);
54 
55  for (const auto& detClusters : clusterCollection) {
57 
58  for (const auto& cluster : detClusters) {
59  float barycenter = cluster.barycenter();
60  uint8_t width = cluster.width();
61  float avgCharge = cluster.avgCharge();
62 
63  switch (approxVersion) {
64  case 0: //ORIGINAL
65  barycenter = std::round(barycenter);
66  if (width > 0x3F)
67  width = 0x3F;
68  avgCharge = std::round(avgCharge);
69  break;
70  case 1: //FULL_WIDTH
71  barycenter = std::round(barycenter);
72  avgCharge = std::round(avgCharge);
73  break;
74  case 2: //BARY_RES_0.1
75  barycenter = std::round(barycenter * 10) / 10;
76  if (width > 0x3F)
77  width = 0x3F;
78  avgCharge = std::round(avgCharge);
79  break;
80  case 3: //BARY_CHARGE_RES_0.1
81  barycenter = std::round(barycenter * 10) / 10;
82  if (width > 0x3F)
83  width = 0x3F;
84  avgCharge = std::round(avgCharge * 10) / 10;
85  break;
86  }
87 
88  ff.push_back(SiStripApproximateCluster(barycenter, width, avgCharge));
89  }
90  }
91 
92  event.put(std::move(result));
93 }
94 
97  desc.add<edm::InputTag>("inputApproxClusters", edm::InputTag("siStripClusters"));
98  desc.add<std::string>("approxVersion", std::string("ORIGINAL"));
99 
100  descriptions.add("SiStripApprox2ApproxClusters", desc);
101 }
102 
ConfigurationDescriptions.h
ApeEstimator_cff.width
width
Definition: ApeEstimator_cff.py:24
SiStripApprox2ApproxClusters::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: SiStripApprox2ApproxClusters.cc:95
SiStripApprox2ApproxClusters::clusterToken
edm::EDGetTokenT< edmNew::DetSetVector< SiStripApproximateCluster > > clusterToken
Definition: SiStripApprox2ApproxClusters.cc:29
SiStripApprox2ApproxClusters::inputApproxClusters
edm::InputTag inputApproxClusters
Definition: SiStripApprox2ApproxClusters.cc:26
edm::EDGetTokenT
Definition: EDGetToken.h:33
SiStripApproximateCluster.h
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
EDProducer.h
SiStripApprox2ApproxClusters::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: SiStripApprox2ApproxClusters.cc:51
SiStripApprox2ApproxClusters
Definition: SiStripApprox2ApproxClusters.cc:18
SiStripApprox2ApproxClusters::SiStripApprox2ApproxClusters
SiStripApprox2ApproxClusters(const edm::ParameterSet &conf)
Definition: SiStripApprox2ApproxClusters.cc:32
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
alignCSCRings.ff
ff
Definition: alignCSCRings.py:148
ParameterSetDescription.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
SiStripApprox2ApproxClusters::approxVersion
uint8_t approxVersion
Definition: SiStripApprox2ApproxClusters.cc:27
edm::ParameterSet
Definition: ParameterSet.h:47
SiStripCluster.h
SiStripApproximateCluster
Definition: SiStripApproximateCluster.h:11
edm::stream::EDProducer
Definition: EDProducer.h:36
edm::EventSetup
Definition: EventSetup.h:58
DetSetVector.h
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
InputTag.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
Frameworkfwd.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
SiStripApprox2ApproxClusters::approxVersionS
std::string approxVersionS
Definition: SiStripApprox2ApproxClusters.cc:28
edmNew::DetSetVector::FastFiller
Definition: DetSetVectorNew.h:202
mps_fire.result
result
Definition: mps_fire.py:311
ParameterSet.h
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
DetSetVectorNew.h
edm::InputTag
Definition: InputTag.h:15