Skip to content

Commit be7612c

Browse files
author
Ivan
committed
added switch in Approx2Approx and fixes as per PR request
1 parent 65c6831 commit be7612c

File tree

5 files changed

+48
-40
lines changed

5 files changed

+48
-40
lines changed

DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h

Lines changed: 0 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -15,15 +15,12 @@ class SiStripApproximateCluster {
1515
explicit SiStripApproximateCluster(float barycenter, uint8_t width, float avgCharge) {
1616
barycenter_ = barycenter;
1717
width_ = width;
18-
//if(width_>0x3F) width_=0x3F;
1918
avgCharge_ = avgCharge;
2019
}
2120

2221
explicit SiStripApproximateCluster(const SiStripCluster& cluster) {
23-
//barycenter_=std::round(cluster.barycenter());
2422
barycenter_ = cluster.barycenter();
2523
width_ = cluster.size();
26-
//if(width_>0x3F) width_=0x3F;
2724
avgCharge_ = cluster.charge() / cluster.size();
2825
}
2926

DataFormats/SiStripCluster/src/classes_def.xml

Lines changed: 4 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -23,7 +23,10 @@
2323
<class name="edm::Wrapper<edmNew::DetSetVector<edm::Ref<edmNew::DetSetVector<SiStripCluster>,SiStripCluster,edmNew::DetSetVector<SiStripCluster>::FindForDetSetVector> > >" />
2424

2525

26-
<class name="SiStripApproximateCluster"/>
26+
<class name="SiStripApproximateCluster" ClassVersion="3">
27+
<version ClassVersion="3" checksum="2041370183"/>
28+
</class>
29+
2730
<class name="edmNew::DetSetVector<SiStripApproximateCluster>"/>
2831
<class name="edm::Wrapper<edmNew::DetSetVector<SiStripApproximateCluster>>"/>
2932

RecoLocalTracker/SiStripClusterizer/BuildFile.xml

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,8 +1,10 @@
11
<use name="DataFormats/Common"/>
22
<use name="FWCore/Framework"/>
33
<use name="FWCore/ParameterSet"/>
4+
<use name="FWCore/Utilities"/>
45
<use name="DataFormats/SiStripDigi"/>
56
<use name="DataFormats/SiStripCluster"/>
7+
<use name="DataFormats/Common"/>
68
<use name="CondFormats/DataRecord"/>
79
<use name="CondFormats/SiStripObjects"/>
810
<use name="CalibFormats/SiStripObjects"/>
@@ -13,4 +15,4 @@
1315
<use name="root"/>
1416
<export>
1517
<lib name="1"/>
16-
</export>
18+
</export>

RecoLocalTracker/SiStripClusterizer/plugins/SiStripApprox2ApproxClusters.cc

Lines changed: 37 additions & 30 deletions
Original file line numberDiff line numberDiff line change
@@ -4,9 +4,9 @@
44
#include "FWCore/Framework/interface/Frameworkfwd.h"
55
#include "FWCore/Framework/interface/stream/EDProducer.h"
66
#include "FWCore/ParameterSet/interface/ParameterSet.h"
7-
#include "FWCore/Utilities/interface/InputTag.h"
87
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
98
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
9+
#include "FWCore/Utilities/interface/InputTag.h"
1010
#include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"
1111
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
1212
#include "DataFormats/Common/interface/DetSetVectorNew.h"
@@ -24,58 +24,65 @@ class SiStripApprox2ApproxClusters : public edm::stream::EDProducer<> {
2424

2525
private:
2626
edm::InputTag inputApproxClusters;
27-
std::string approxVersion;
27+
uint8_t approxVersion;
28+
std::string approxVersionS;
2829
edm::EDGetTokenT<edmNew::DetSetVector<SiStripApproximateCluster>> clusterToken;
2930
};
3031

3132
SiStripApprox2ApproxClusters::SiStripApprox2ApproxClusters(const edm::ParameterSet& conf) {
3233
inputApproxClusters = conf.getParameter<edm::InputTag>("inputApproxClusters");
33-
approxVersion = conf.getParameter<std::string>("approxVersion");
34-
34+
approxVersionS = conf.getParameter<std::string>("approxVersion");
35+
36+
approxVersion=-1;
37+
38+
if(approxVersionS=="ORIGINAL") approxVersion=0;
39+
else if(approxVersionS=="FULL_WIDTH") approxVersion=1;
40+
else if(approxVersionS=="BARY_RES_0.1") approxVersion=2;
41+
else if(approxVersionS=="BARY_CHARGE_RES_0.1") approxVersion=3;
42+
3543
clusterToken = consumes<edmNew::DetSetVector<SiStripApproximateCluster>>(inputApproxClusters);
3644
produces<edmNew::DetSetVector<SiStripApproximateCluster>>();
3745
}
3846

3947
void SiStripApprox2ApproxClusters::produce(edm::Event& event, edm::EventSetup const&) {
4048
auto result = std::make_unique<edmNew::DetSetVector<SiStripApproximateCluster>>();
41-
const auto& clusters = event.get(clusterToken);
49+
const auto& clusterCollection = event.get(clusterToken);
4250

43-
for (const auto& detClusters : *clusterCollection) {
51+
for (const auto& detClusters : clusterCollection) {
4452
edmNew::DetSetVector<SiStripApproximateCluster>::FastFiller ff{*result, detClusters.id()};
4553

4654
for (const auto& cluster : detClusters) {
4755
float barycenter = cluster.barycenter();
4856
uint8_t width = cluster.width();
4957
float avgCharge = cluster.avgCharge();
50-
51-
if (approxVersion == "ORIGINAL") {
52-
barycenter = std::round(barycenter);
53-
if (width > 0x3F)
54-
width = 0x3F;
55-
avgCharge = std::round(avgCharge);
56-
57-
} else if (approxVersion == "FULL_WIDTH") {
58-
barycenter = std::round(barycenter);
59-
avgCharge = std::round(avgCharge);
60-
61-
} else if (approxVersion == "BARY_RES_0.1") {
62-
barycenter = std::round(barycenter * 10) / 10;
63-
if (width > 0x3F)
64-
width = 0x3F;
65-
avgCharge = std::round(avgCharge);
66-
67-
} else if (approxVersion == "BARY_CHARGE_RES_0.1") {
68-
barycenter = std::round(barycenter * 10) / 10;
69-
if (width > 0x3F)
70-
width = 0x3F;
71-
avgCharge = std::round(avgCharge * 10) / 10;
58+
59+
switch(approxVersion){
60+
case 0: //ORIGINAL
61+
barycenter = std::round(barycenter);
62+
if (width > 0x3F) width = 0x3F;
63+
avgCharge = std::round(avgCharge);
64+
break;
65+
case 1: //FULL_WIDTH
66+
barycenter = std::round(barycenter);
67+
avgCharge = std::round(avgCharge);
68+
break;
69+
case 2: //BARY_RES_0.1
70+
barycenter = std::round(barycenter * 10) / 10;
71+
if (width > 0x3F) width = 0x3F;
72+
avgCharge = std::round(avgCharge);
73+
break;
74+
case 3: //BARY_CHARGE_RES_0.1
75+
barycenter = std::round(barycenter * 10) / 10;
76+
if (width > 0x3F) width = 0x3F;
77+
avgCharge = std::round(avgCharge * 10) / 10;
78+
break;
7279
}
73-
80+
7481
ff.push_back(SiStripApproximateCluster(barycenter, width, avgCharge));
7582
}
7683
}
7784

78-
e.put(std::move(result));
85+
event.put(std::move(result));
7986
}
8087

8188
void SiStripApprox2ApproxClusters::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {

RecoLocalTracker/SiStripClusterizer/plugins/SiStripClusters2ApproxClusters.cc

Lines changed: 4 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -41,20 +41,19 @@ SiStripClusters2ApproxClusters::SiStripClusters2ApproxClusters(const edm::Parame
4141

4242
}
4343

44-
void SiStripClusters2ApproxClusters::produce(edm::Event& e, edm::EventSetup const&){
44+
void SiStripClusters2ApproxClusters::produce(edm::Event& event, edm::EventSetup const&){
4545
auto result = std::make_unique<edmNew::DetSetVector< SiStripApproximateCluster > >();
46-
edm::Handle<edmNew::DetSetVector< SiStripCluster >> clusterCollection = e.getHandle(clusterToken);
46+
const auto& clusterCollection = event.get(clusterToken);
4747

4848

49-
for ( const auto& detClusters : *clusterCollection ) {
50-
std::vector< SiStripApproximateCluster > tempVec;
49+
for ( const auto& detClusters : clusterCollection ) {
5150
edmNew::DetSetVector<SiStripApproximateCluster>::FastFiller ff{*result, detClusters.id()};
5251

5352
for ( const auto& cluster : detClusters ) ff.push_back(SiStripApproximateCluster(cluster));
5453

5554
}
5655

57-
e.put(std::move(result));
56+
event.put(std::move(result));
5857
}
5958

6059
void SiStripClusters2ApproxClusters::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {

0 commit comments

Comments
 (0)