Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Improve messages due VertexException in SequentialPrimaryVertexFitterAdapter and clean-up of BeamSpotOnlineProducer #46893

Merged
merged 2 commits into from
Dec 9, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
280 changes: 150 additions & 130 deletions RecoVertex/BeamSpotProducer/plugins/BeamSpotOnlineProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -41,6 +41,16 @@ class BeamSpotOnlineProducer : public edm::stream::EDProducer<> {
static void fillDescriptions(edm::ConfigurationDescriptions& iDesc);

private:
// helper methods
bool shouldShout(const edm::Event& iEvent) const;
bool processTransientRecord(const edm::EventSetup& iSetup, reco::BeamSpot& result, bool shoutMODE);
void createBeamSpotFromTransient(const BeamSpotObjects& spotDB, reco::BeamSpot& result) const;
bool processScalers(const edm::Event& iEvent, reco::BeamSpot& result, bool shoutMODE);
void createBeamSpotFromScaler(const BeamSpotOnline& spotOnline, reco::BeamSpot& result) const;
bool isInvalidScaler(const BeamSpotOnline& spotOnline, bool shoutMODE) const;
void createBeamSpotFromDB(const edm::EventSetup& iSetup, reco::BeamSpot& result, bool shoutMODE) const;

// data members
const bool changeFrame_;
const double theMaxZ, theSetSigmaZ;
double theMaxR2;
Expand Down Expand Up @@ -88,151 +98,161 @@ void BeamSpotOnlineProducer::fillDescriptions(edm::ConfigurationDescriptions& iD
iDesc.addWithDefaultLabel(ps);
}

void BeamSpotOnlineProducer::produce(Event& iEvent, const EventSetup& iSetup) {
// product is a reco::BeamSpot object
void BeamSpotOnlineProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
auto result = std::make_unique<reco::BeamSpot>();
reco::BeamSpot aSpot;
//shout MODE only in stable beam
bool shoutMODE = false;

// Determine if we should "shout" based on the beam mode
bool shoutMODE = shouldShout(iEvent);
bool fallBackToDB{false};

// Use transient record if enabled
if (useTransientRecord_) {
fallBackToDB = processTransientRecord(iSetup, *result, shoutMODE);
} else {
// Process online beam spot scalers
fallBackToDB = processScalers(iEvent, *result, shoutMODE);
}

// Fallback to DB if necessary
if (fallBackToDB) {
createBeamSpotFromDB(iSetup, *result, shoutMODE);
}

iEvent.put(std::move(result));
}

bool BeamSpotOnlineProducer::shouldShout(const edm::Event& iEvent) const {
edm::Handle<L1GlobalTriggerEvmReadoutRecord> gtEvmReadoutRecord;
if (iEvent.getByToken(l1GtEvmReadoutRecordToken_, gtEvmReadoutRecord)) {
if (gtEvmReadoutRecord->gtfeWord().beamMode() == theBeamShoutMode)
shoutMODE = true;
} else {
shoutMODE = true;
return gtEvmReadoutRecord->gtfeWord().beamMode() == theBeamShoutMode;
}
bool fallBackToDB = false;
if (useTransientRecord_) {
auto const& spotDB = iSetup.getData(beamTransientToken_);
if (spotDB.beamType() != 2) {
if (shoutMODE && beamTransientRcdESWatcher_.check(iSetup)) {
edm::LogWarning("BeamSpotOnlineProducer")
<< "Online Beam Spot producer falls back to DB value because the ESProducer returned a fake beamspot ";
}
fallBackToDB = true;
} else {
// translate from BeamSpotObjects to reco::BeamSpot
// in case we need to switch to LHC reference frame
// ignore for the moment rotations, and translations
double f = 1.;
if (changeFrame_)
f = -1.;
reco::BeamSpot::Point apoint(f * spotDB.x(), f * spotDB.y(), f * spotDB.z());

reco::BeamSpot::CovarianceMatrix matrix;
for (int i = 0; i < reco::BeamSpot::dimension; ++i) {
for (int j = 0; j < reco::BeamSpot::dimension; ++j) {
matrix(i, j) = spotDB.covariance(i, j);
}
}
double sigmaZ = spotDB.sigmaZ();
if (theSetSigmaZ > 0)
sigmaZ = theSetSigmaZ;

// this assume beam width same in x and y
aSpot = reco::BeamSpot(apoint, sigmaZ, spotDB.dxdz(), spotDB.dydz(), spotDB.beamWidthX(), matrix);
aSpot.setBeamWidthY(spotDB.beamWidthY());
aSpot.setEmittanceX(spotDB.emittanceX());
aSpot.setEmittanceY(spotDB.emittanceY());
aSpot.setbetaStar(spotDB.betaStar());
aSpot.setType(reco::BeamSpot::Tracker);
return true; // Default to "shout" if the record is missing
}

bool BeamSpotOnlineProducer::processTransientRecord(const edm::EventSetup& iSetup,
reco::BeamSpot& result,
bool shoutMODE) {
auto const& spotDB = iSetup.getData(beamTransientToken_);

if (spotDB.beamType() != 2) {
if (shoutMODE && beamTransientRcdESWatcher_.check(iSetup)) {
edm::LogWarning("BeamSpotOnlineProducer")
<< "Online Beam Spot producer falls back to DB value due to fake beamspot in transient record.";
}
} else {
// get scalar collection
Handle<BeamSpotOnlineCollection> handleScaler;
iEvent.getByToken(scalerToken_, handleScaler);

// beam spot scalar object
BeamSpotOnline spotOnline;

// product is a reco::BeamSpot object
auto result = std::make_unique<reco::BeamSpot>();

if (!handleScaler->empty()) {
// get one element
spotOnline = *(handleScaler->begin());

// in case we need to switch to LHC reference frame
// ignore for the moment rotations, and translations
double f = 1.;
if (changeFrame_)
f = -1.;

reco::BeamSpot::Point apoint(f * spotOnline.x(), spotOnline.y(), f * spotOnline.z());

reco::BeamSpot::CovarianceMatrix matrix;
matrix(0, 0) = spotOnline.err_x() * spotOnline.err_x();
matrix(1, 1) = spotOnline.err_y() * spotOnline.err_y();
matrix(2, 2) = spotOnline.err_z() * spotOnline.err_z();
matrix(3, 3) = spotOnline.err_sigma_z() * spotOnline.err_sigma_z();

double sigmaZ = spotOnline.sigma_z();
if (theSetSigmaZ > 0)
sigmaZ = theSetSigmaZ;

aSpot = reco::BeamSpot(apoint, sigmaZ, spotOnline.dxdz(), f * spotOnline.dydz(), spotOnline.width_x(), matrix);

aSpot.setBeamWidthY(spotOnline.width_y());
aSpot.setEmittanceX(0.);
aSpot.setEmittanceY(0.);
aSpot.setbetaStar(0.);
aSpot.setType(reco::BeamSpot::LHC); // flag value from scalars

// check if we have a valid beam spot fit result from online DQM
if (spotOnline.x() == 0 && spotOnline.y() == 0 && spotOnline.z() == 0 && spotOnline.width_x() == 0 &&
spotOnline.width_y() == 0) {
if (shoutMODE) {
edm::LogWarning("BeamSpotOnlineProducer")
<< "Online Beam Spot producer falls back to DB value because the scaler values are zero ";
}
fallBackToDB = true;
}
double r2 = spotOnline.x() * spotOnline.x() + spotOnline.y() * spotOnline.y();
if (std::abs(spotOnline.z()) >= theMaxZ || r2 >= theMaxR2) {
if (shoutMODE) {
edm::LogError("BeamSpotOnlineProducer")
<< "Online Beam Spot producer falls back to DB value because the scaler values are too big to be true :"
<< spotOnline.x() << " " << spotOnline.y() << " " << spotOnline.z();
}
fallBackToDB = true;
}
} else {
//empty online beamspot collection: FED data was empty
//the error should probably have been send at unpacker level
fallBackToDB = true;
return true; // Trigger fallback to DB
}

// Create BeamSpot from transient record
createBeamSpotFromTransient(spotDB, result);
return false; // No fallback needed
}

void BeamSpotOnlineProducer::createBeamSpotFromTransient(const BeamSpotObjects& spotDB, reco::BeamSpot& result) const {
double f = changeFrame_ ? -1.0 : 1.0;
reco::BeamSpot::Point apoint(f * spotDB.x(), f * spotDB.y(), f * spotDB.z());

reco::BeamSpot::CovarianceMatrix matrix;
for (int i = 0; i < reco::BeamSpot::dimension; ++i) {
for (int j = 0; j < reco::BeamSpot::dimension; ++j) {
matrix(i, j) = spotDB.covariance(i, j);
}
}
if (fallBackToDB) {
edm::ESHandle<BeamSpotObjects> beamhandle = iSetup.getHandle(beamToken_);
const BeamSpotObjects* spotDB = beamhandle.product();

// translate from BeamSpotObjects to reco::BeamSpot
reco::BeamSpot::Point apoint(spotDB->x(), spotDB->y(), spotDB->z());
double sigmaZ = (theSetSigmaZ > 0) ? theSetSigmaZ : spotDB.sigmaZ();
result = reco::BeamSpot(apoint, sigmaZ, spotDB.dxdz(), spotDB.dydz(), spotDB.beamWidthX(), matrix);
result.setBeamWidthY(spotDB.beamWidthY());
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Not related to this PR (if not because of the part of the original code which is touched by it), but why not also setBeamWidthX?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

but why not also setBeamWidthX

the x-width is already set in the constructor at L162:

  result = reco::BeamSpot(apoint, sigmaZ, spotDB.dxdz(), spotDB.dydz(), spotDB.beamWidthX(), matrix);

as per method signature:

BeamSpot(const Point& point,
double sigmaZ,
double dxdz,
double dydz,
double BeamWidthX,
const CovarianceMatrix& error,
BeamType type = Unknown) {

as both widths are constructed with the same value (in the hypothesis of round-optics)

BeamWidthX_ = BeamWidthX;
BeamWidthY_ = BeamWidthX;

setting explicitly the y component is necessary if the two differ (and in principle they will differ in 2025, if the flat optics is used).

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

setting explicitly the y component is necessary if the two differ (and in principle they will differ in 2025, if the flat optics is used).

Ah, ok: that was indeed my worry. Thank you for confirming that the two directions can be set independently,

result.setEmittanceX(spotDB.emittanceX());
result.setEmittanceY(spotDB.emittanceY());
result.setbetaStar(spotDB.betaStar());
result.setType(reco::BeamSpot::Tracker);
}

bool BeamSpotOnlineProducer::processScalers(const edm::Event& iEvent, reco::BeamSpot& result, bool shoutMODE) {
edm::Handle<BeamSpotOnlineCollection> handleScaler;
iEvent.getByToken(scalerToken_, handleScaler);

if (handleScaler->empty()) {
return true; // Fallback to DB if scaler collection is empty
}

reco::BeamSpot::CovarianceMatrix matrix;
for (int i = 0; i < reco::BeamSpot::dimension; ++i) {
for (int j = 0; j < reco::BeamSpot::dimension; ++j) {
matrix(i, j) = spotDB->covariance(i, j);
}
// Extract data from scaler
BeamSpotOnline spotOnline = *(handleScaler->begin());
createBeamSpotFromScaler(spotOnline, result);

// Validate the scaler data
if (isInvalidScaler(spotOnline, shoutMODE)) {
return true; // Trigger fallback to DB
}
return false; // No fallback needed
}

void BeamSpotOnlineProducer::createBeamSpotFromScaler(const BeamSpotOnline& spotOnline, reco::BeamSpot& result) const {
double f = changeFrame_ ? -1.0 : 1.0;
reco::BeamSpot::Point apoint(f * spotOnline.x(), spotOnline.y(), f * spotOnline.z());

reco::BeamSpot::CovarianceMatrix matrix;
matrix(0, 0) = spotOnline.err_x() * spotOnline.err_x();
matrix(1, 1) = spotOnline.err_y() * spotOnline.err_y();
matrix(2, 2) = spotOnline.err_z() * spotOnline.err_z();
matrix(3, 3) = spotOnline.err_sigma_z() * spotOnline.err_sigma_z();

double sigmaZ = (theSetSigmaZ > 0) ? theSetSigmaZ : spotOnline.sigma_z();
result = reco::BeamSpot(apoint, sigmaZ, spotOnline.dxdz(), f * spotOnline.dydz(), spotOnline.width_x(), matrix);
result.setBeamWidthY(spotOnline.width_y());
result.setEmittanceX(0.0);
result.setEmittanceY(0.0);
result.setbetaStar(0.0);
result.setType(reco::BeamSpot::LHC);
}

bool BeamSpotOnlineProducer::isInvalidScaler(const BeamSpotOnline& spotOnline, bool shoutMODE) const {
if (spotOnline.x() == 0 && spotOnline.y() == 0 && spotOnline.z() == 0 && spotOnline.width_x() == 0 &&
spotOnline.width_y() == 0) {
if (shoutMODE) {
edm::LogWarning("BeamSpotOnlineProducer")
<< "Online Beam Spot producer falls back to DB due to zero scaler values.";
}
return true;
}

// this assume beam width same in x and y
aSpot = reco::BeamSpot(apoint, spotDB->sigmaZ(), spotDB->dxdz(), spotDB->dydz(), spotDB->beamWidthX(), matrix);
aSpot.setBeamWidthY(spotDB->beamWidthY());
aSpot.setEmittanceX(spotDB->emittanceX());
aSpot.setEmittanceY(spotDB->emittanceY());
aSpot.setbetaStar(spotDB->betaStar());
aSpot.setType(reco::BeamSpot::Tracker);

GlobalError bse(aSpot.rotatedCovariance3D());
if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
edm::LogError("UnusableBeamSpot") << "Beamspot from fallback to DB with invalid errors: " << aSpot.covariance();
double r2 = spotOnline.x() * spotOnline.x() + spotOnline.y() * spotOnline.y();
if (std::abs(spotOnline.z()) >= theMaxZ || r2 >= theMaxR2) {
if (shoutMODE) {
edm::LogError("BeamSpotOnlineProducer")
<< "Online Beam Spot producer falls back to DB due to out-of-range scaler values: " << spotOnline.x() << ", "
<< spotOnline.y() << ", " << spotOnline.z();
}
return true;
}
return false;
}

*result = aSpot;
void BeamSpotOnlineProducer::createBeamSpotFromDB(const edm::EventSetup& iSetup,
reco::BeamSpot& result,
bool shoutMODE) const {
edm::ESHandle<BeamSpotObjects> beamhandle = iSetup.getHandle(beamToken_);
const BeamSpotObjects* spotDB = beamhandle.product();

iEvent.put(std::move(result));
reco::BeamSpot::Point apoint(spotDB->x(), spotDB->y(), spotDB->z());

reco::BeamSpot::CovarianceMatrix matrix;
for (int i = 0; i < reco::BeamSpot::dimension; ++i) {
for (int j = 0; j < reco::BeamSpot::dimension; ++j) {
matrix(i, j) = spotDB->covariance(i, j);
}
}

result = reco::BeamSpot(apoint, spotDB->sigmaZ(), spotDB->dxdz(), spotDB->dydz(), spotDB->beamWidthX(), matrix);
result.setBeamWidthY(spotDB->beamWidthY());
result.setEmittanceX(spotDB->emittanceX());
result.setEmittanceY(spotDB->emittanceY());
result.setbetaStar(spotDB->betaStar());
result.setType(reco::BeamSpot::Tracker);

GlobalError bse(result.rotatedCovariance3D());
if ((bse.cxx() <= 0.0) || (bse.cyy() <= 0.0) || (bse.czz() <= 0.0)) {
edm::LogError("UnusableBeamSpot") << "Beamspot from DB fallback has invalid errors: " << result.covariance();
}
}

#include "FWCore/Framework/interface/MakerMacros.h"
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,8 @@

*/

#include <sstream>

#include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
#include "RecoVertex/PrimaryVertexProducer/interface/PrimaryVertexFitterBase.h"
Expand All @@ -27,7 +29,15 @@ class SequentialPrimaryVertexFitterAdapter : public PrimaryVertexFitterBase {
const std::vector<reco::TransientTrack>& tracklist = cluster.originalTracks();
TransientVertex v;
if (useBeamConstraint && (tracklist.size() > 1)) {
v = fitter->vertex(tracklist, beamspot);
try {
v = fitter->vertex(tracklist, beamspot);
} catch (VertexException& ex) {
std::ostringstream beamspotInfo;
beamspotInfo << "While processing SequentialPrimaryVertexFitterAdapter::fit() with BeamSpot parameters: \n"
<< beamspot;
ex.addContext(beamspotInfo.str());
throw; // rethrow the exception
}
} else if (!(useBeamConstraint) && (tracklist.size() > 1)) {
v = fitter->vertex(tracklist);
} // else: no fit ==> v.isValid()=False
Expand Down