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 
60 }
61 
64  dsc.add<int>("timeThreshold", 48)->setComment("hours");
65  dsc.add<double>("sigmaZThreshold", 2.)->setComment("cm");
66  dsc.add<double>("sigmaXYThreshold", 4.)->setComment("um");
67  desc.addWithDefaultLabel(dsc);
68 }
69 
71  const BeamSpotOnlineObjects* bs2) {
72  // Current time to be compared with the one read from the payload
73  auto currentTime =
74  std::chrono::duration_cast<std::chrono::microseconds>(std::chrono::system_clock::now().time_since_epoch());
75 
76  // Get two beamspot creation times and compute the time difference wrt currentTime
77  auto bs1time = std::chrono::microseconds(bs1->creationTime());
78  auto diffBStime1 = (currentTime - bs1time).count();
79  auto bs2time = std::chrono::microseconds(bs2->creationTime());
80  auto diffBStime2 = (currentTime - bs2time).count();
81 
82  // Convert timeThreshold_ from hours to microseconds for comparison
83  auto limitTime = std::chrono::microseconds((std::chrono::hours)timeThreshold_).count();
84 
85  // Logic to choose between the two BeamSpots:
86  // 1. If both BS are older than limitTime retun fake BS
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
91  if (diffBStime1 > limitTime && diffBStime2 > limitTime) {
92  edm::LogInfo("OnlineBeamSpotESProducer") << "Defaulting to fake because both payloads are too old.";
93  return nullptr;
94  } else if (diffBStime2 > limitTime) {
95  if (isGoodBS(bs1)) {
96  return bs1;
97  } else {
98  edm::LogInfo("OnlineBeamSpotESProducer")
99  << "Defaulting to fake because the legacy Beam Spot is not suitable and HLT one is too old.";
100  return nullptr;
101  }
102  } else if (diffBStime1 > limitTime) {
103  if (isGoodBS(bs2)) {
104  return bs2;
105  } else {
106  edm::LogInfo("OnlineBeamSpotESProducer")
107  << "Defaulting to fake because the HLT Beam Spot is not suitable and the legacy one too old.";
108  return nullptr;
109  }
110  } else {
111  if (bs1->sigmaZ() > bs2->sigmaZ() && isGoodBS(bs1)) {
112  return bs1;
113  } else if (bs2->sigmaZ() >= bs1->sigmaZ() && isGoodBS(bs2)) {
114  return bs2;
115  } else {
116  edm::LogInfo("OnlineBeamSpotESProducer")
117  << "Defaulting to fake because despite both payloads are young enough, none has passed the fit sanity checks";
118  return nullptr;
119  }
120  }
121 }
122 
124  // Current time to be compared with the one read from the payload
125  auto currentTime =
126  std::chrono::duration_cast<std::chrono::microseconds>(std::chrono::system_clock::now().time_since_epoch());
127 
128  // Get the beamspot creation time and compute the time difference wrt currentTime
129  auto bs1time = std::chrono::microseconds(bs1->creationTime());
130  auto diffBStime1 = (currentTime - bs1time).count();
131 
132  // Convert timeThreshold_ from hours to microseconds for comparison
133  auto limitTime = std::chrono::microseconds((std::chrono::hours)timeThreshold_).count();
134 
135  // Check that the BS is within the timeThreshold, converges and passes the sigmaZthreshold
136  if (diffBStime1 < limitTime && isGoodBS(bs1)) {
137  return bs1;
138  } else {
139  return nullptr;
140  }
141 }
142 
143 // This method is used to check the quality of the beamspot fit
145  if (bs1->sigmaZ() > sigmaZThreshold_ && bs1->beamType() == reco::BeamSpot::Tracker &&
147  return true;
148  else
149  return false;
150 }
151 
152 std::shared_ptr<const BeamSpotObjects> OnlineBeamSpotESProducer::produce(const BeamSpotTransientObjectsRcd& iRecord) {
153  auto legacyRec = iRecord.tryToGetRecord<BeamSpotOnlineLegacyObjectsRcd>();
154  auto hltRec = iRecord.tryToGetRecord<BeamSpotOnlineHLTObjectsRcd>();
155  if (not legacyRec and not hltRec) {
156  edm::LogInfo("OnlineBeamSpotESProducer") << "None of the Beam Spots in ES are available! \n returning a fake one.";
157  return std::shared_ptr<const BeamSpotObjects>(&fakeBS_, edm::do_nothing_deleter());
158  }
159 
160  const BeamSpotOnlineObjects* best;
161  if (legacyRec and hltRec) {
162  best = compareBS(&legacyRec->get(bsLegacyToken_), &hltRec->get(bsHLTToken_));
163  } else if (legacyRec) {
164  best = checkSingleBS(&legacyRec->get(bsLegacyToken_));
165  } else {
166  best = checkSingleBS(&hltRec->get(bsHLTToken_));
167  }
168  if (best) {
169  return std::shared_ptr<const BeamSpotObjects>(best, edm::do_nothing_deleter());
170  } else {
171  return std::shared_ptr<const BeamSpotObjects>(&fakeBS_, edm::do_nothing_deleter());
172  edm::LogInfo("OnlineBeamSpotESProducer")
173  << "None of the Online BeamSpots in the ES is suitable, \n returning a fake one. ";
174  }
175 };
176 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
void setComment(std::string const &value)
cond::Time_t creationTime() const
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
Log< level::Info, false > LogInfo
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
void setBeamWidthY(double val)
set average transverse beam width Y