CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes | Static Private Attributes
BeamSpotCompatibilityChecker Class Reference

#include <RecoVertex/BeamSpotProducer/plugins/BeamSpotCompatibilityChecker.cc>

Inheritance diagram for BeamSpotCompatibilityChecker:
edm::global::EDAnalyzer<> edm::global::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

void analyze (edm::StreamID, edm::Event const &, edm::EventSetup const &) const final
 
 BeamSpotCompatibilityChecker (const edm::ParameterSet &)
 
 ~BeamSpotCompatibilityChecker () override=default
 
- Public Member Functions inherited from edm::global::EDAnalyzer<>
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static std::array< float, 6 > compareBS (const reco::BeamSpot &BSA, const reco::BeamSpot &BSB)
 
static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::global::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Attributes

edm::ESGetToken< BeamSpotObjects, BeamSpotObjectsRcdbsFromDBToken_
 
const edm::EDGetTokenT< reco::BeamSpotbsFromEventToken_
 
edm::EDGetTokenT< reco::BeamSpotdbBSfromEventToken_
 
const bool dbFromEvent_
 
const double throwingThreshold_
 
const bool verbose_
 
const double warningThreshold_
 

Static Private Attributes

static constexpr int cmToum = 10000
 

Additional Inherited Members

- Public Types inherited from edm::global::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Description: Class to check the compatibilty between the BeamSpot payload in the database and the one in the event

Implementation: Makes use of the Significance struct to establish how compatible are the data members of the two BeamSpots in input

Definition at line 78 of file BeamSpotCompatibilityChecker.cc.

Constructor & Destructor Documentation

◆ BeamSpotCompatibilityChecker()

BeamSpotCompatibilityChecker::BeamSpotCompatibilityChecker ( const edm::ParameterSet iConfig)
explicit

Definition at line 102 of file BeamSpotCompatibilityChecker.cc.

References bsFromDBToken_, dbBSfromEventToken_, dbFromEvent_, Exception, edm::ParameterSet::getParameter(), throwingThreshold_, and warningThreshold_.

103  : warningThreshold_(iConfig.getParameter<double>("warningThr")),
104  throwingThreshold_(iConfig.getParameter<double>("errorThr")),
105  verbose_(iConfig.getUntrackedParameter<bool>("verbose", false)),
106  dbFromEvent_(iConfig.getParameter<bool>("dbFromEvent")),
107  bsFromEventToken_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("bsFromEvent"))) {
108  //now do what ever initialization is needed
110  throw cms::Exception("ConfigurationError")
111  << __PRETTY_FUNCTION__ << "\n Warning threshold (" << warningThreshold_
112  << ") cannot be smaller than the throwing threshold (" << throwingThreshold_ << ")" << std::endl;
113  }
114  if (dbFromEvent_) {
115  edm::LogWarning("BeamSpotCompatibilityChecker")
116  << "!!!! Warning !!!\nThe Database Beam Spot is going to be taken from the Event via BeamSpotProducer!";
117  dbBSfromEventToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("bsFromDB"));
118  } else {
119  bsFromDBToken_ = esConsumes<BeamSpotObjects, BeamSpotObjectsRcd>();
120  }
121 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::EDGetTokenT< reco::BeamSpot > bsFromEventToken_
edm::ESGetToken< BeamSpotObjects, BeamSpotObjectsRcd > bsFromDBToken_
edm::EDGetTokenT< reco::BeamSpot > dbBSfromEventToken_
T getUntrackedParameter(std::string const &, T const &) const
Log< level::Warning, false > LogWarning

◆ ~BeamSpotCompatibilityChecker()

BeamSpotCompatibilityChecker::~BeamSpotCompatibilityChecker ( )
overridedefault

Member Function Documentation

◆ analyze()

void BeamSpotCompatibilityChecker::analyze ( edm::StreamID  sid,
edm::Event const &  iEvent,
edm::EventSetup const &  iSetup 
) const
finalvirtual

Implements edm::global::EDAnalyzerBase.

Definition at line 128 of file BeamSpotCompatibilityChecker.cc.

References funct::abs(), align::BeamSpot, BeamSpotObjects::beamWidthX(), reco::BeamSpot::BeamWidthX(), reco::BeamSpot::BeamWidthY(), bsFromDBToken_, bsFromEventToken_, cmToum, compareBS(), BeamSpotObjects::covariance(), dbBSfromEventToken_, dbFromEvent_, reco::BeamSpot::dimension, BeamSpotObjects::dxdz(), BeamSpotObjects::dydz(), Exception, edm::EventSetup::getHandle(), mps_fire::i, iEvent, dqmiolumiharvest::j, SummaryClient_cfi::labels, makeMuonMisalignmentScenario::matrix, mps_check::msg, edm::ESHandle< T >::product(), BeamSpotObjects::sigmaZ(), reco::BeamSpot::sigmaZ(), AlCaHLTBitMon_QueryRunRegistry::string, throwingThreshold_, to_string(), verbose_, warningThreshold_, BeamSpotObjects::x(), reco::BeamSpot::x0(), BeamSpotObjects::y(), reco::BeamSpot::y0(), BeamSpotObjects::z(), and reco::BeamSpot::z0().

130  {
131  using namespace edm;
132  reco::BeamSpot spotEvent, spotDB;
133 
134  edm::Handle<reco::BeamSpot> beamSpotFromEventHandle;
135  iEvent.getByToken(bsFromEventToken_, beamSpotFromEventHandle);
136  spotEvent = *beamSpotFromEventHandle;
137 
138  double evt_BSx0 = spotEvent.x0();
139  double evt_BSy0 = spotEvent.y0();
140  double evt_BSz0 = spotEvent.z0();
141  double evt_Beamsigmaz = spotEvent.sigmaZ();
142  double evt_BeamWidthX = spotEvent.BeamWidthX();
143  double evt_BeamWidthY = spotEvent.BeamWidthY();
144 
145  if (!dbFromEvent_) {
146  edm::ESHandle<BeamSpotObjects> beamhandle = iSetup.getHandle(bsFromDBToken_);
147  const BeamSpotObjects* aSpot = beamhandle.product();
148 
149  // translate from BeamSpotObjects to reco::BeamSpot
150  reco::BeamSpot::Point apoint(aSpot->x(), aSpot->y(), aSpot->z());
151 
153  for (int i = 0; i < reco::BeamSpot::dimension; ++i) {
154  for (int j = 0; j < reco::BeamSpot::dimension; ++j) {
155  matrix(i, j) = aSpot->covariance(i, j);
156  }
157  }
158 
159  // this assume beam width same in x and y
160  spotDB = reco::BeamSpot(apoint, aSpot->sigmaZ(), aSpot->dxdz(), aSpot->dydz(), aSpot->beamWidthX(), matrix);
161  } else {
162  // take the DB beamspot from the event (different label)
163  edm::Handle<reco::BeamSpot> beamSpotFromDBHandle;
164  iEvent.getByToken(dbBSfromEventToken_, beamSpotFromDBHandle);
165  spotDB = *beamSpotFromDBHandle;
166  }
167 
168  double db_BSx0 = spotDB.x0();
169  double db_BSy0 = spotDB.y0();
170  double db_BSz0 = spotDB.z0();
171  double db_Beamsigmaz = spotDB.sigmaZ();
172  double db_BeamWidthX = spotDB.BeamWidthX();
173  double db_BeamWidthY = spotDB.BeamWidthY();
174 
175  double deltaX0 = evt_BSx0 - db_BSx0;
176  double deltaY0 = evt_BSy0 - db_BSy0;
177  double deltaZ0 = evt_BSz0 - db_BSz0;
178  double deltaSigmaX = evt_BeamWidthX - db_BeamWidthX;
179  double deltaSigmaY = evt_BeamWidthY - db_BeamWidthY;
180  double deltaSigmaZ = evt_Beamsigmaz - db_Beamsigmaz;
181 
182  if (verbose_) {
183  edm::LogPrint("BeamSpotCompatibilityChecker") << "BS from Event: \n" << spotEvent << std::endl;
184  edm::LogPrint("BeamSpotCompatibilityChecker") << "BS from DB: \n" << spotDB << std::endl;
185  }
186 
187  auto significances = compareBS(spotDB, spotEvent);
188  std::vector<std::string> labels = {"x0", "y0", "z0", "sigmaX", "sigmaY", "sigmaZ"};
189 
190  std::string msg = " |delta X0|=" + std::to_string(std::abs(deltaX0) * cmToum) +
191  " um |delta Y0|=" + std::to_string(std::abs(deltaY0) * cmToum) +
192  " um |delta Z0|=" + std::to_string(std::abs(deltaZ0) * cmToum) +
193  " um |delta sigmaX|=" + std::to_string(std::abs(deltaSigmaX) * cmToum) +
194  " um |delta sigmaY|=" + std::to_string(std::abs(deltaSigmaY) * cmToum) +
195  " um |delta sigmaZ|=" + std::to_string(std::abs(deltaSigmaZ)) + "cm";
196  if (verbose_) {
197  edm::LogPrint("BeamSpotCompatibilityChecker") << msg.c_str() << std::endl;
198  }
199 
200  for (unsigned int i = 0; i < 3; i++) {
201  auto sig = significances.at(i);
202  if (sig > throwingThreshold_) {
203  edm::LogError("BeamSpotCompatibilityChecker") << msg.c_str() << std::endl;
204  throw cms::Exception("BeamSpotCompatibilityChecker")
205  << "[" << __PRETTY_FUNCTION__ << "] \n DB-Event BeamSpot " << labels.at(i) << " distance sigificance " << sig
206  << ", exceeds the threshold of " << throwingThreshold_ << "!" << std::endl;
207  } else if (sig > warningThreshold_) {
208  edm::LogWarning("BeamSpotCompatibilityChecker") << msg.c_str() << std::endl;
209  edm::LogWarning("BeamSpotCompatibilityChecker")
210  << "[" << __PRETTY_FUNCTION__ << "] \n DB-Event BeamSpot " << labels.at(i) << " distance significance "
211  << sig << ", exceeds the threshold of " << warningThreshold_ << "!" << std::endl;
212  }
213  }
214 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:29
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:82
double z() const
get Z beam position
double dydz() const
get dydz slope, crossing angle in YZ
const edm::EDGetTokenT< reco::BeamSpot > bsFromEventToken_
double covariance(int i, int j) const
get i,j element of the full covariance matrix 7x7
static std::array< float, 6 > compareBS(const reco::BeamSpot &BSA, const reco::BeamSpot &BSB)
edm::ESGetToken< BeamSpotObjects, BeamSpotObjectsRcd > bsFromDBToken_
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:27
double beamWidthX() const
get average transverse beam width
edm::EDGetTokenT< reco::BeamSpot > dbBSfromEventToken_
Log< level::Error, false > LogError
static std::string to_string(const XMLCh *ch)
double x0() const
x coordinate
Definition: BeamSpot.h:61
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:84
double x() const
get X beam position
double y0() const
y coordinate
Definition: BeamSpot.h:63
Log< level::Warning, true > LogPrint
double y() const
get Y beam position
double sigmaZ() const
sigma z
Definition: BeamSpot.h:76
double sigmaZ() const
get sigma Z, RMS bunch length
tuple msg
Definition: mps_check.py:286
double z0() const
z coordinate
Definition: BeamSpot.h:65
HLT enums.
Log< level::Warning, false > LogWarning
double dxdz() const
get dxdz slope, crossing angle in XZ

◆ compareBS()

std::array< float, 6 > BeamSpotCompatibilityChecker::compareBS ( const reco::BeamSpot BSA,
const reco::BeamSpot BSB 
)
static

Definition at line 216 of file BeamSpotCompatibilityChecker.cc.

References reco::BeamSpot::BeamWidthX(), reco::BeamSpot::BeamWidthXError(), reco::BeamSpot::BeamWidthY(), reco::BeamSpot::BeamWidthYError(), runTheMatrix::ret, reco::BeamSpot::sigmaZ(), reco::BeamSpot::sigmaZ0Error(), reco::BeamSpot::x0(), reco::BeamSpot::x0Error(), reco::BeamSpot::y0(), reco::BeamSpot::y0Error(), reco::BeamSpot::z0(), and reco::BeamSpot::z0Error().

Referenced by analyze().

217  {
218  auto xS = Significance(beamSpotA.x0(), beamSpotB.x0(), beamSpotA.x0Error(), beamSpotB.x0Error());
219  auto yS = Significance(beamSpotA.y0(), beamSpotB.y0(), beamSpotA.y0Error(), beamSpotB.y0Error());
220  auto zS = Significance(beamSpotA.z0(), beamSpotB.z0(), beamSpotA.z0Error(), beamSpotB.z0Error());
221  auto xWS = Significance(
222  beamSpotA.BeamWidthX(), beamSpotB.BeamWidthX(), beamSpotA.BeamWidthXError(), beamSpotB.BeamWidthXError());
223  auto yWS = Significance(
224  beamSpotA.BeamWidthY(), beamSpotB.BeamWidthY(), beamSpotA.BeamWidthYError(), beamSpotB.BeamWidthYError());
225  auto zWS = Significance(beamSpotA.sigmaZ(), beamSpotB.sigmaZ(), beamSpotA.sigmaZ0Error(), beamSpotB.sigmaZ0Error());
226 
227  std::array<float, 6> ret = {
228  {xS.getSig(false), yS.getSig(false), zS.getSig(false), xWS.getSig(false), yWS.getSig(false), zWS.getSig(false)}};
229  return ret;
230 }
ret
prodAgent to be discontinued

◆ fillDescriptions()

void BeamSpotCompatibilityChecker::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 233 of file BeamSpotCompatibilityChecker.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

233  {
235  desc.add<double>("warningThr", 1.)->setComment("Threshold on the signficances to emit a warning");
236  desc.add<double>("errorThr", 3.)->setComment("Threshold on the signficances to abort the job");
237  desc.addUntracked<bool>("verbose", false)->setComment("verbose output");
238  desc.add<edm::InputTag>("bsFromEvent", edm::InputTag(""))
239  ->setComment("edm::InputTag on the BeamSpot from the Event (Reference)");
240  desc.add<bool>("dbFromEvent", false)
241  ->setComment("Switch to take the (target) DB beamspot from the event instead of the EventSetup");
242  desc.add<edm::InputTag>("bsFromDB", edm::InputTag(""))
243  ->setComment("edm::InputTag on the BeamSpot from the Event (Target)\n To be used only if dbFromEvent is True");
244  descriptions.addWithDefaultLabel(desc);
245 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ bsFromDBToken_

edm::ESGetToken<BeamSpotObjects, BeamSpotObjectsRcd> BeamSpotCompatibilityChecker::bsFromDBToken_
private

Definition at line 95 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze(), and BeamSpotCompatibilityChecker().

◆ bsFromEventToken_

const edm::EDGetTokenT<reco::BeamSpot> BeamSpotCompatibilityChecker::bsFromEventToken_
private

Definition at line 94 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze().

◆ cmToum

constexpr int BeamSpotCompatibilityChecker::cmToum = 10000
staticprivate

Definition at line 89 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze().

◆ dbBSfromEventToken_

edm::EDGetTokenT<reco::BeamSpot> BeamSpotCompatibilityChecker::dbBSfromEventToken_
private

Definition at line 96 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze(), and BeamSpotCompatibilityChecker().

◆ dbFromEvent_

const bool BeamSpotCompatibilityChecker::dbFromEvent_
private

Definition at line 93 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze(), and BeamSpotCompatibilityChecker().

◆ throwingThreshold_

const double BeamSpotCompatibilityChecker::throwingThreshold_
private

Definition at line 91 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze(), and BeamSpotCompatibilityChecker().

◆ verbose_

const bool BeamSpotCompatibilityChecker::verbose_
private

Definition at line 92 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze().

◆ warningThreshold_

const double BeamSpotCompatibilityChecker::warningThreshold_
private

Definition at line 90 of file BeamSpotCompatibilityChecker.cc.

Referenced by analyze(), and BeamSpotCompatibilityChecker().