Skip to content

Commit 9fa71b4

Browse files
authored
Merge pull request cms-sw#40744 from dzuolo/from-CMSSW_13_0_X_2023-02-10-1100
Adding a check on the BS transverse widths in OnlineBeamSpotESProducer
2 parents 3a272a9 + 49ad628 commit 9fa71b4

File tree

1 file changed

+24
-11
lines changed

1 file changed

+24
-11
lines changed

RecoVertex/BeamSpotProducer/plugins/OnlineBeamSpotESProducer.cc

Lines changed: 24 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -30,6 +30,7 @@ class OnlineBeamSpotESProducer : public edm::ESProducer {
3030
private:
3131
const BeamSpotOnlineObjects* compareBS(const BeamSpotOnlineObjects* bs1, const BeamSpotOnlineObjects* bs2);
3232
const BeamSpotOnlineObjects* checkSingleBS(const BeamSpotOnlineObjects* bs1);
33+
bool isGoodBS(const BeamSpotOnlineObjects* bs1) const;
3334

3435
edm::ESGetToken<BeamSpotObjects, BeamSpotTransientObjectsRcd> const bsToken_;
3536
edm::ESGetToken<BeamSpotOnlineObjects, BeamSpotOnlineHLTObjectsRcd> bsHLTToken_;
@@ -38,12 +39,14 @@ class OnlineBeamSpotESProducer : public edm::ESProducer {
3839
BeamSpotObjects fakeBS_;
3940
const int timeThreshold_;
4041
const double sigmaZThreshold_;
42+
const double sigmaXYThreshold_;
4143
};
4244

4345
OnlineBeamSpotESProducer::OnlineBeamSpotESProducer(const edm::ParameterSet& p)
4446
// get parameters
4547
: timeThreshold_(p.getParameter<int>("timeThreshold")),
46-
sigmaZThreshold_(p.getParameter<double>("sigmaZThreshold")) {
48+
sigmaZThreshold_(p.getParameter<double>("sigmaZThreshold")),
49+
sigmaXYThreshold_(p.getParameter<double>("sigmaXYThreshold") * 1E-4) {
4750
auto cc = setWhatProduced(this);
4851

4952
fakeBS_.setBeamWidthX(0.1);
@@ -60,6 +63,7 @@ void OnlineBeamSpotESProducer::fillDescriptions(edm::ConfigurationDescriptions&
6063
edm::ParameterSetDescription dsc;
6164
dsc.add<int>("timeThreshold", 48)->setComment("hours");
6265
dsc.add<double>("sigmaZThreshold", 2.)->setComment("cm");
66+
dsc.add<double>("sigmaXYThreshold", 4.)->setComment("um");
6367
desc.addWithDefaultLabel(dsc);
6468
}
6569

@@ -80,37 +84,37 @@ const BeamSpotOnlineObjects* OnlineBeamSpotESProducer::compareBS(const BeamSpotO
8084

8185
// Logic to choose between the two BeamSpots:
8286
// 1. If both BS are older than limitTime retun fake BS
83-
// 2. If only one BS is newer than limitTime return it only if it has
84-
// sigmaZ larger than sigmaZthreshold_ and the fit converged (BeamType 2)
85-
// 3. If both are newer than the limit threshold return
86-
// the BS that converged and has larger sigmaZ
87+
// 2. If only one BS is newer than limitTime return it only if
88+
// it passes isGoodBS (checks on sigmaZ, sigmaXY and fit convergence)
89+
// 3. If both are newer than the limit threshold return the BS that
90+
// passes isGoodBS and has larger sigmaZ
8791
if (diffBStime1 > limitTime && diffBStime2 > limitTime) {
8892
edm::LogInfo("OnlineBeamSpotESProducer") << "Defaulting to fake because both payloads are too old.";
8993
return nullptr;
9094
} else if (diffBStime2 > limitTime) {
91-
if (bs1->sigmaZ() > sigmaZThreshold_ && bs1->beamType() == 2) {
95+
if (isGoodBS(bs1)) {
9296
return bs1;
9397
} else {
9498
edm::LogInfo("OnlineBeamSpotESProducer")
9599
<< "Defaulting to fake because the legacy Beam Spot is not suitable and HLT one is too old.";
96100
return nullptr;
97101
}
98102
} else if (diffBStime1 > limitTime) {
99-
if (bs2->sigmaZ() > sigmaZThreshold_ && bs2->beamType() == 2) {
103+
if (isGoodBS(bs2)) {
100104
return bs2;
101105
} else {
102106
edm::LogInfo("OnlineBeamSpotESProducer")
103107
<< "Defaulting to fake because the HLT Beam Spot is not suitable and the legacy one too old.";
104108
return nullptr;
105109
}
106110
} else {
107-
if (bs1->sigmaZ() > bs2->sigmaZ() && bs1->beamType() == 2) {
111+
if (bs1->sigmaZ() > bs2->sigmaZ() && isGoodBS(bs1)) {
108112
return bs1;
109-
} else if (bs2->sigmaZ() >= bs1->sigmaZ() && bs2->beamType() == 2) {
113+
} else if (bs2->sigmaZ() >= bs1->sigmaZ() && isGoodBS(bs2)) {
110114
return bs2;
111115
} else {
112116
edm::LogInfo("OnlineBeamSpotESProducer")
113-
<< "Defaulting to fake because despite both payloads are young enough, none has the right BeamType.";
117+
<< "Defaulting to fake because despite both payloads are young enough, none has passed the fit sanity checks";
114118
return nullptr;
115119
}
116120
}
@@ -129,13 +133,22 @@ const BeamSpotOnlineObjects* OnlineBeamSpotESProducer::checkSingleBS(const BeamS
129133
auto limitTime = std::chrono::microseconds((std::chrono::hours)timeThreshold_).count();
130134

131135
// Check that the BS is within the timeThreshold, converges and passes the sigmaZthreshold
132-
if (diffBStime1 < limitTime && bs1->sigmaZ() > sigmaZThreshold_ && bs1->beamType() == 2) {
136+
if (diffBStime1 < limitTime && isGoodBS(bs1)) {
133137
return bs1;
134138
} else {
135139
return nullptr;
136140
}
137141
}
138142

143+
// This method is used to check the quality of the beamspot fit
144+
bool OnlineBeamSpotESProducer::isGoodBS(const BeamSpotOnlineObjects* bs1) const {
145+
if (bs1->sigmaZ() > sigmaZThreshold_ && bs1->beamType() == reco::BeamSpot::Tracker &&
146+
bs1->beamWidthX() > sigmaXYThreshold_ && bs1->beamWidthY() > sigmaXYThreshold_)
147+
return true;
148+
else
149+
return false;
150+
}
151+
139152
std::shared_ptr<const BeamSpotObjects> OnlineBeamSpotESProducer::produce(const BeamSpotTransientObjectsRcd& iRecord) {
140153
auto legacyRec = iRecord.tryToGetRecord<BeamSpotOnlineLegacyObjectsRcd>();
141154
auto hltRec = iRecord.tryToGetRecord<BeamSpotOnlineHLTObjectsRcd>();

0 commit comments

Comments
 (0)