CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
37 
39  const int timeThreshold_;
40  const double sigmaZThreshold_;
41 };
42 
44  // get parameters
45  : timeThreshold_(p.getParameter<int>("timeThreshold")),
46  sigmaZThreshold_(p.getParameter<double>("sigmaZThreshold")) {
47  auto cc = setWhatProduced(this);
48 
51  fakeBS_.SetSigmaZ(15.);
52  fakeBS_.SetPosition(0.0001, 0.0001, 0.0001);
53  fakeBS_.SetType(-1);
54 
57 }
58 
61  dsc.add<int>("timeThreshold", 48)->setComment("hours");
62  dsc.add<double>("sigmaZThreshold", 2.)->setComment("cm");
63  desc.addWithDefaultLabel(dsc);
64 }
65 
67  const BeamSpotOnlineObjects* bs2) {
68  // Current time to be compared with the one read from the payload
69  auto currentTime =
70  std::chrono::duration_cast<std::chrono::microseconds>(std::chrono::system_clock::now().time_since_epoch());
71 
72  // Get two beamspot creation times and compute the time difference wrt currentTime
73  auto bs1time = std::chrono::microseconds(bs1->GetCreationTime());
74  auto diffBStime1 = (currentTime - bs1time).count();
75  auto bs2time = std::chrono::microseconds(bs2->GetCreationTime());
76  auto diffBStime2 = (currentTime - bs2time).count();
77 
78  // Convert timeThreshold_ from hours to microseconds for comparison
79  auto limitTime = std::chrono::microseconds((std::chrono::hours)timeThreshold_).count();
80 
81  // Logic to choose between the two BeamSpots:
82  // 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  if (diffBStime1 > limitTime && diffBStime2 > limitTime) {
88  edm::LogInfo("OnlineBeamSpotESProducer") << "Defaulting to fake becuase both payloads are too old.";
89  return nullptr;
90  } else if (diffBStime2 > limitTime) {
91  if (bs1->GetSigmaZ() > sigmaZThreshold_ && bs1->GetBeamType() == 2) {
92  return bs1;
93  } else {
94  edm::LogInfo("OnlineBeamSpotESProducer")
95  << "Defaulting to fake because the legacy Beam Spot is not suitable and HLT one is too old.";
96  return nullptr;
97  }
98  } else if (diffBStime1 > limitTime) {
99  if (bs2->GetSigmaZ() > sigmaZThreshold_ && bs2->GetBeamType() == 2) {
100  return bs2;
101  } else {
102  edm::LogInfo("OnlineBeamSpotESProducer")
103  << "Defaulting to fake because the HLT Beam Spot is not suitable and the legacy one too old.";
104  return nullptr;
105  }
106  } else {
107  if (bs1->GetSigmaZ() > bs2->GetSigmaZ() && bs1->GetBeamType() == 2) {
108  return bs1;
109  } else if (bs2->GetSigmaZ() >= bs1->GetSigmaZ() && bs2->GetBeamType() == 2) {
110  return bs2;
111  } else {
112  edm::LogInfo("OnlineBeamSpotESProducer")
113  << "Defaulting to fake because despite both payloads are young enough, none has the right BeamType.";
114  return nullptr;
115  }
116  }
117 }
118 
120  // Current time to be compared with the one read from the payload
121  auto currentTime =
122  std::chrono::duration_cast<std::chrono::microseconds>(std::chrono::system_clock::now().time_since_epoch());
123 
124  // Get the beamspot creation time and compute the time difference wrt currentTime
125  auto bs1time = std::chrono::microseconds(bs1->GetCreationTime());
126  auto diffBStime1 = (currentTime - bs1time).count();
127 
128  // Convert timeThreshold_ from hours to microseconds for comparison
129  auto limitTime = std::chrono::microseconds((std::chrono::hours)timeThreshold_).count();
130 
131  // Check that the BS is within the timeThreshold, converges and passes the sigmaZthreshold
132  if (diffBStime1 < limitTime && bs1->GetSigmaZ() > sigmaZThreshold_ && bs1->GetBeamType() == 2) {
133  return bs1;
134  } else {
135  return nullptr;
136  }
137 }
138 
139 std::shared_ptr<const BeamSpotObjects> OnlineBeamSpotESProducer::produce(const BeamSpotTransientObjectsRcd& iRecord) {
140  auto legacyRec = iRecord.tryToGetRecord<BeamSpotOnlineLegacyObjectsRcd>();
141  auto hltRec = iRecord.tryToGetRecord<BeamSpotOnlineHLTObjectsRcd>();
142  if (not legacyRec and not hltRec) {
143  edm::LogInfo("OnlineBeamSpotESProducer") << "None of the Beam Spots in ES are available! \n returning a fake one.";
144  return std::shared_ptr<const BeamSpotObjects>(&fakeBS_, edm::do_nothing_deleter());
145  }
146 
147  const BeamSpotOnlineObjects* best;
148  if (legacyRec and hltRec) {
149  best = compareBS(&legacyRec->get(bsLegacyToken_), &hltRec->get(bsHLTToken_));
150  } else if (legacyRec) {
151  best = checkSingleBS(&legacyRec->get(bsLegacyToken_));
152  } else {
153  best = checkSingleBS(&hltRec->get(bsHLTToken_));
154  }
155  if (best) {
156  return std::shared_ptr<const BeamSpotObjects>(best, edm::do_nothing_deleter());
157  } else {
158  return std::shared_ptr<const BeamSpotObjects>(&fakeBS_, edm::do_nothing_deleter());
159  edm::LogInfo("OnlineBeamSpotESProducer")
160  << "None of the Online BeamSpots in the ES is suitable, \n returning a fake one. ";
161  }
162 };
163 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:163
void setComment(std::string const &value)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const BeamSpotOnlineObjects * checkSingleBS(const BeamSpotOnlineObjects *bs1)
cond::Time_t GetCreationTime() const
void SetSigmaZ(double val)
set sigma Z, RMS bunch length
double GetSigmaZ() const
get sigma Z, RMS bunch length
edm::ESGetToken< BeamSpotOnlineObjects, BeamSpotOnlineHLTObjectsRcd > bsHLTToken_
OnlineBeamSpotESProducer(const edm::ParameterSet &p)
static void fillDescriptions(edm::ConfigurationDescriptions &desc)
int GetBeamType() const
get beam type
const BeamSpotOnlineObjects * compareBS(const BeamSpotOnlineObjects *bs1, const BeamSpotOnlineObjects *bs2)
edm::ESGetToken< BeamSpotOnlineObjects, BeamSpotOnlineLegacyObjectsRcd > bsLegacyToken_
edm::ESGetToken< BeamSpotObjects, BeamSpotTransientObjectsRcd > const bsToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Info, false > LogInfo
void SetType(int type)
set beam type
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::shared_ptr< const BeamSpotObjects > produce(const BeamSpotTransientObjectsRcd &)
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
void SetBeamWidthX(double val)
set average transverse beam width X
void SetBeamWidthY(double val)
set average transverse beam width Y
void SetPosition(double x, double y, double z)
set XYZ position