CMS 3D CMS Logo

OnlineBeamSpotESProducer.cc
Go to the documentation of this file.
19 #include <memory>
20 #include <string>
21 
22 using namespace edm;
23 
25 public:
27  std::shared_ptr<const BeamSpotObjects> produce(const BeamSpotTransientObjectsRcd&);
29 
30 private:
31  const BeamSpotOnlineObjects* compareBS(const BeamSpotOnlineObjects* bs1, const BeamSpotOnlineObjects* bs2);
32  const BeamSpotOnlineObjects* checkSingleBS(const BeamSpotOnlineObjects* bs1);
33  bool isGoodBS(const BeamSpotOnlineObjects* bs1) const;
34 
38 
40  const int timeThreshold_;
41  const double sigmaZThreshold_;
42  const double sigmaXYThreshold_;
43 };
44 
46  // get parameters
47  : timeThreshold_(p.getParameter<int>("timeThreshold")),
48  sigmaZThreshold_(p.getParameter<double>("sigmaZThreshold")),
49  sigmaXYThreshold_(p.getParameter<double>("sigmaXYThreshold") * 1E-4) {
50  auto cc = setWhatProduced(this);
51 
54  fakeBS_.setSigmaZ(15.);
55  fakeBS_.setPosition(0.0001, 0.0001, 0.0001);
56  fakeBS_.setType(-1);
57  // Set diagonal covariance, i.e. errors on the parameters
58  fakeBS_.setCovariance(0, 0, 5e-10);
59  fakeBS_.setCovariance(1, 1, 5e-10);
60  fakeBS_.setCovariance(2, 2, 0.002);
61  fakeBS_.setCovariance(3, 3, 0.002);
62  fakeBS_.setCovariance(4, 4, 5e-11);
63  fakeBS_.setCovariance(5, 5, 5e-11);
64  fakeBS_.setCovariance(6, 6, 1e-09);
65 
68 }
69 
72  dsc.add<int>("timeThreshold", 48)->setComment("hours");
73  dsc.add<double>("sigmaZThreshold", 2.)->setComment("cm");
74  dsc.add<double>("sigmaXYThreshold", 4.)->setComment("um");
75  desc.addWithDefaultLabel(dsc);
76 }
77 
79  const BeamSpotOnlineObjects* bs2) {
80  // Current time to be compared with the one read from the payload
81  auto currentTime =
82  std::chrono::duration_cast<std::chrono::microseconds>(std::chrono::system_clock::now().time_since_epoch());
83 
84  // Get two beamspot creation times and compute the time difference wrt currentTime
85  auto bs1time = std::chrono::microseconds(bs1->creationTime());
86  auto diffBStime1 = (currentTime - bs1time).count();
87  auto bs2time = std::chrono::microseconds(bs2->creationTime());
88  auto diffBStime2 = (currentTime - bs2time).count();
89 
90  // Convert timeThreshold_ from hours to microseconds for comparison
91  auto limitTime = std::chrono::microseconds((std::chrono::hours)timeThreshold_).count();
92 
93  // Logic to choose between the two BeamSpots:
94  // 1. If both BS are older than limitTime retun fake BS
95  // 2. If only one BS is newer than limitTime return it only if
96  // it passes isGoodBS (checks on sigmaZ, sigmaXY and fit convergence)
97  // 3. If both are newer than the limit threshold return the BS that
98  // passes isGoodBS and has larger sigmaZ
99  if (diffBStime1 > limitTime && diffBStime2 > limitTime) {
100  edm::LogWarning("OnlineBeamSpotESProducer")
101  << "Defaulting to fake (fallback to PCL) because both payloads are too old.";
102  return nullptr;
103  } else if (diffBStime2 > limitTime) {
104  if (isGoodBS(bs1)) {
105  return bs1;
106  } else {
107  edm::LogWarning("OnlineBeamSpotESProducer") << "Defaulting to fake (fallback to PCL) because the legacy Beam "
108  "Spot is not suitable and HLT one is too old.";
109  return nullptr;
110  }
111  } else if (diffBStime1 > limitTime) {
112  if (isGoodBS(bs2)) {
113  return bs2;
114  } else {
115  edm::LogWarning("OnlineBeamSpotESProducer") << "Defaulting to fake (fallback to PCL) because the HLT Beam Spot "
116  "is not suitable and the legacy one too old.";
117  return nullptr;
118  }
119  } else {
120  if (bs1->sigmaZ() > bs2->sigmaZ() && isGoodBS(bs1)) {
121  return bs1;
122  } else if (bs2->sigmaZ() >= bs1->sigmaZ() && isGoodBS(bs2)) {
123  return bs2;
124  } else {
125  edm::LogWarning("OnlineBeamSpotESProducer") << "Defaulting to fake (fallback to PCL) because despite both "
126  "payloads are young enough, none has passed the fit sanity checks";
127  return nullptr;
128  }
129  }
130 }
131 
133  // Current time to be compared with the one read from the payload
134  auto currentTime =
135  std::chrono::duration_cast<std::chrono::microseconds>(std::chrono::system_clock::now().time_since_epoch());
136 
137  // Get the beamspot creation time and compute the time difference wrt currentTime
138  auto bs1time = std::chrono::microseconds(bs1->creationTime());
139  auto diffBStime1 = (currentTime - bs1time).count();
140 
141  // Convert timeThreshold_ from hours to microseconds for comparison
142  auto limitTime = std::chrono::microseconds((std::chrono::hours)timeThreshold_).count();
143 
144  // Check that the BS is within the timeThreshold, converges and passes the sigmaZthreshold
145  if (diffBStime1 < limitTime && isGoodBS(bs1)) {
146  return bs1;
147  } else {
148  return nullptr;
149  }
150 }
151 
152 // This method is used to check the quality of the beamspot fit
154  if (bs1->sigmaZ() > sigmaZThreshold_ && bs1->beamType() == reco::BeamSpot::Tracker &&
156  return true;
157  else
158  return false;
159 }
160 
161 std::shared_ptr<const BeamSpotObjects> OnlineBeamSpotESProducer::produce(const BeamSpotTransientObjectsRcd& iRecord) {
162  auto legacyRec = iRecord.tryToGetRecord<BeamSpotOnlineLegacyObjectsRcd>();
163  auto hltRec = iRecord.tryToGetRecord<BeamSpotOnlineHLTObjectsRcd>();
164  if (not legacyRec and not hltRec) {
165  edm::LogWarning("OnlineBeamSpotESProducer")
166  << "None of the Beam Spots in ES are available! \n returning a fake one (fallback to PCL).";
167  return std::shared_ptr<const BeamSpotObjects>(&fakeBS_, edm::do_nothing_deleter());
168  }
169 
170  const BeamSpotOnlineObjects* best;
171  if (legacyRec and hltRec) {
172  best = compareBS(&legacyRec->get(bsLegacyToken_), &hltRec->get(bsHLTToken_));
173  } else if (legacyRec) {
174  best = checkSingleBS(&legacyRec->get(bsLegacyToken_));
175  } else {
176  best = checkSingleBS(&hltRec->get(bsHLTToken_));
177  }
178  if (best) {
179  return std::shared_ptr<const BeamSpotObjects>(best, edm::do_nothing_deleter());
180  } else {
181  return std::shared_ptr<const BeamSpotObjects>(&fakeBS_, edm::do_nothing_deleter());
182  edm::LogWarning("OnlineBeamSpotESProducer")
183  << "None of the Online BeamSpots in the ES is suitable, \n returning a fake one(fallback to PCL).";
184  }
185 };
186 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:166
void setComment(std::string const &value)
cond::Time_t creationTime() const
void setCovariance(int i, int j, double val)
set i,j element of the full covariance matrix 7x7
const BeamSpotOnlineObjects * checkSingleBS(const BeamSpotOnlineObjects *bs1)
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
edm::ESGetToken< BeamSpotOnlineObjects, BeamSpotOnlineHLTObjectsRcd > bsHLTToken_
double beamWidthX() const
get average transverse beam width
int beamType() const
get beam type
void setType(int type)
set beam type
OnlineBeamSpotESProducer(const edm::ParameterSet &p)
static void fillDescriptions(edm::ConfigurationDescriptions &desc)
const BeamSpotOnlineObjects * compareBS(const BeamSpotOnlineObjects *bs1, const BeamSpotOnlineObjects *bs2)
bool isGoodBS(const BeamSpotOnlineObjects *bs1) const
edm::ESGetToken< BeamSpotOnlineObjects, BeamSpotOnlineLegacyObjectsRcd > bsLegacyToken_
edm::ESGetToken< BeamSpotObjects, BeamSpotTransientObjectsRcd > const bsToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:61
double beamWidthY() const
get average transverse beam width
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::shared_ptr< const BeamSpotObjects > produce(const BeamSpotTransientObjectsRcd &)
void setSigmaZ(double val)
set sigma Z, RMS bunch length
double sigmaZ() const
get sigma Z, RMS bunch length
HLT enums.
void setPosition(double x, double y, double z)
set XYZ position
void setBeamWidthX(double val)
set average transverse beam width X
Log< level::Warning, false > LogWarning
void setBeamWidthY(double val)
set average transverse beam width Y