CMS 3D CMS Logo

BTransitionAnalyzer.h
Go to the documentation of this file.
1 #ifndef BTRANSITIONANALYZER_H
2 #define BTRANSITIONANALYZER_H
3 
14 
15 namespace cond {
16  template <class T, class R>
17  class BTransitionAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
18  public:
20  : m_currentThreshold(pset.getUntrackedParameter<double>("currentThreshold", 18000.)),
22  m_ESToken(esConsumes<edm::Transition::EndRun>()),
23  m_ESTokenB0T(esConsumes<edm::Transition::EndRun>(edm::ESInputTag("", "0T"))),
24  m_ESTokenB38T(esConsumes<edm::Transition::EndRun>(edm::ESInputTag("", "38T"))) {}
25 #ifdef __INTEL_COMPILER
26  virtual ~BTransitionAnalyzer() = default;
27 #endif
28  // implicit copy constructor
29  // implicit assignment operator
30  // implicit destructor
31  void beginJob() final{};
32  void beginRun(edm::Run const&, edm::EventSetup const&) final{};
33  void analyze(edm::Event const&, edm::EventSetup const&) final{};
34  void endRun(edm::Run const& run, edm::EventSetup const& eventSetup) final {
35  edm::ESHandle<RunInfo> runInfoHandle = eventSetup.getHandle(m_RunInfoToken);
36  edm::ESHandle<T> payloadHandle, payloadRefHandle;
37  double avg_current = (double)runInfoHandle->m_avg_current;
38  double current_default = -1;
39  std::string bOnLabel = std::string("38T");
40  std::string bOffLabel = std::string("0T");
41  std::string bFieldLabel = bOnLabel;
42  LogDebug("BTransitionAnalyzer") << "Comparing value of magnet current: " << avg_current
43  << " A for run: " << run.run()
44  << " with the corresponding threshold: " << m_currentThreshold << " A."
45  << std::endl;
46  if (avg_current != current_default && avg_current <= m_currentThreshold) {
47  bFieldLabel = bOffLabel;
48  payloadHandle = eventSetup.getHandle(m_ESTokenB0T);
49  } else {
50  payloadHandle = eventSetup.getHandle(m_ESTokenB38T);
51  }
52  edm::LogInfo("BTransitionAnalyzer")
53  << "The magnet was " << (bFieldLabel == bOnLabel ? "ON" : "OFF") << " during run " << run.run()
54  << ".\nLoading the product for the corrisponding label " << bFieldLabel << std::endl;
55  payloadRefHandle = eventSetup.getHandle(m_ESToken);
57  if (mydbservice.isAvailable()) {
58  if (!equalPayloads(payloadHandle, payloadRefHandle)) {
59  edm::LogInfo("BTransitionAnalyzer")
60  << "Exporting payload corresponding to the calibrations for magnetic field "
61  << (bFieldLabel == bOnLabel ? "ON" : "OFF") << " starting from run number: " << run.run() << std::endl;
62  mydbservice->writeOneIOV(*payloadHandle.product(), run.run(), demangledName(typeid(R)));
63  } else {
64  edm::LogInfo("BTransitionAnalyzer") << "The payload corresponding to the calibrations for magnetic field "
65  << (bFieldLabel == bOnLabel ? "ON" : "OFF") << " is still valid for run "
66  << run.run() << ".\nNo transfer needed." << std::endl;
67  }
68  } else {
69  edm::LogError("BTransitionAnalyzer") << "PoolDBOutputService unavailable";
70  }
71  }
72  void endJob() final{};
73  virtual bool equalPayloads(edm::ESHandle<T> const& payloadHandle, edm::ESHandle<T> const& payloadRefHandle) = 0;
74 
75  private:
81  };
82 } //namespace cond
83 #endif //BTRANSITIONANALYZER_H
void analyze(edm::Event const &, edm::EventSetup const &) final
void endRun(edm::Run const &run, edm::EventSetup const &eventSetup) final
BTransitionAnalyzer(const edm::ParameterSet &pset)
virtual bool equalPayloads(edm::ESHandle< T > const &payloadHandle, edm::ESHandle< T > const &payloadRefHandle)=0
Log< level::Error, false > LogError
T const * product() const
Definition: ESHandle.h:86
void beginRun(edm::Run const &, edm::EventSetup const &) final
Hash writeOneIOV(const T &payload, Time_t time, const std::string &recordName)
Transition
Definition: Transition.h:12
const edm::ESGetToken< T, R > m_ESTokenB38T
Log< level::Info, false > LogInfo
const edm::ESGetToken< RunInfo, RunInfoRcd > m_RunInfoToken
const edm::ESGetToken< T, R > m_ESToken
float m_avg_current
Definition: RunInfo.h:28
Definition: plugin.cc:23
HLT enums.
const edm::ESGetToken< T, R > m_ESTokenB0T
bool isAvailable() const
Definition: Service.h:40
Definition: Run.h:45
#define LogDebug(id)