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.)) {}
21 #ifdef __INTEL_COMPILER
22  virtual ~BTransitionAnalyzer() = default;
23 #endif
24  // implicit copy constructor
25  // implicit assignment operator
26  // implicit destructor
27  void beginJob() final{};
28  void beginRun(edm::Run const&, edm::EventSetup const&) final{};
29  void analyze(edm::Event const&, edm::EventSetup const&) final{};
30  void endRun(edm::Run const& run, edm::EventSetup const& eventSetup) final {
31  edm::ESHandle<RunInfo> runInfoHandle;
32  edm::ESHandle<T> payloadHandle, payloadRefHandle;
33  eventSetup.get<RunInfoRcd>().get(runInfoHandle);
34  double avg_current = (double)runInfoHandle->m_avg_current;
35  double current_default = -1;
36  std::string bOnLabel = std::string("38T");
37  std::string bOffLabel = std::string("0T");
38  std::string bFieldLabel = bOnLabel;
39  LogDebug("BTransitionAnalyzer") << "Comparing value of magnet current: " << avg_current
40  << " A for run: " << run.run()
41  << " with the corresponding threshold: " << m_currentThreshold << " A."
42  << std::endl;
43  if (avg_current != current_default && avg_current <= m_currentThreshold)
44  bFieldLabel = bOffLabel;
45  edm::LogInfo("BTransitionAnalyzer")
46  << "The magnet was " << (bFieldLabel == bOnLabel ? "ON" : "OFF") << " during run " << run.run()
47  << ".\nLoading the product for the corrisponding label " << bFieldLabel << std::endl;
48  eventSetup.get<R>().get(bFieldLabel, payloadHandle);
49  eventSetup.get<R>().get(payloadRefHandle);
51  if (mydbservice.isAvailable()) {
52  if (!equalPayloads(payloadHandle, payloadRefHandle)) {
53  edm::LogInfo("BTransitionAnalyzer")
54  << "Exporting payload corresponding to the calibrations for magnetic field "
55  << (bFieldLabel == bOnLabel ? "ON" : "OFF") << " starting from run number: " << run.run() << std::endl;
56  mydbservice->writeOne(payloadHandle.product(), run.run(), demangledName(typeid(R)));
57  } else {
58  edm::LogInfo("BTransitionAnalyzer") << "The payload corresponding to the calibrations for magnetic field "
59  << (bFieldLabel == bOnLabel ? "ON" : "OFF") << " is still valid for run "
60  << run.run() << ".\nNo transfer needed." << std::endl;
61  }
62  } else {
63  edm::LogError("BTransitionAnalyzer") << "PoolDBOutputService unavailable";
64  }
65  }
66  void endJob() final{};
67  virtual bool equalPayloads(edm::ESHandle<T> const& payloadHandle, edm::ESHandle<T> const& payloadRefHandle) = 0;
68 
69  private:
71  };
72 } //namespace cond
73 #endif //BTRANSITIONANALYZER_H
#define LogDebug(id)
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
void beginRun(edm::Run const &, edm::EventSetup const &) final
bool isAvailable() const
Definition: Service.h:40
void writeOne(T *payload, Time_t time, const std::string &recordName, bool withlogging=false)
float m_avg_current
Definition: RunInfo.h:28
Definition: plugin.cc:23
T const * product() const
Definition: ESHandle.h:86
Definition: Run.h:45