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