Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:06

0001 // -*- C++ -*-
0002 //
0003 // Package:    Alignment/CommonAlignment
0004 // Class:      MagneticFieldFilter
0005 //
0006 /**\class MagneticFieldFilter MagneticFieldFilter.cc Alignment/CommonAlignment/plugins/MagneticFieldFilter.cc
0007 
0008  Description: Plugin to filter events based on the magnetic field value
0009 
0010  Implementation:
0011      Takes the magnet current from the RunInfoRcd and translates it into a
0012      magnetic field value using the parameterization given here:
0013 
0014      https://hypernews.cern.ch/HyperNews/CMS/get/magnetic-field/63/1/1/1.html
0015 
0016      Agrees within an accuracy of ~20 mT with results of:
0017      https://cmswbm.web.cern.ch/cmswbm/cmsdb/servlet/RunSummary
0018 
0019 */
0020 //
0021 // Original Author:  Gregor Mittag
0022 //         Created:  Wed, 25 Nov 2015 12:59:02 GMT
0023 //
0024 //
0025 
0026 // user include files
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/Framework/interface/EventSetup.h"
0030 #include "FWCore/Framework/interface/stream/EDFilter.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0032 #include "CondFormats/RunInfo/interface/RunInfo.h"
0033 #include "CondFormats/DataRecord/interface/RunSummaryRcd.h"
0034 
0035 //
0036 // class declaration
0037 //
0038 
0039 class MagneticFieldFilter : public edm::stream::EDFilter<> {
0040 public:
0041   explicit MagneticFieldFilter(const edm::ParameterSet&);
0042   ~MagneticFieldFilter() override = default;
0043 
0044   static void fillDescriptions(edm::ConfigurationDescriptions&);
0045 
0046 private:
0047   bool filter(edm::Event&, const edm::EventSetup&) override;
0048   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0049 
0050   /// convert Ampere (A) to Tesla (T)
0051   float currentToField(const float& current) const;
0052 
0053   // ----------member data ---------------------------
0054 
0055   // esConsumes
0056   const edm::ESGetToken<RunInfo, RunInfoRcd> runInfoToken_;
0057 
0058   /// see: https://hypernews.cern.ch/HyperNews/CMS/get/magnetic-field/63/1/1/1.html
0059   static constexpr float linearCoeffCurrentToField_ = 2.084287e-04;
0060   /// see: https://hypernews.cern.ch/HyperNews/CMS/get/magnetic-field/63/1/1/1.html
0061   static constexpr float constantTermCurrentToField_ = 1.704418e-02;
0062 
0063   const int magneticField_;      /// magnetic field that is filtered
0064   int magneticFieldCurrentRun_;  /// magnetic field estimate of the current run
0065 };
0066 
0067 //
0068 // static data member definitions
0069 //
0070 constexpr float MagneticFieldFilter::linearCoeffCurrentToField_;
0071 constexpr float MagneticFieldFilter::constantTermCurrentToField_;
0072 
0073 //
0074 // constructor
0075 //
0076 MagneticFieldFilter::MagneticFieldFilter(const edm::ParameterSet& iConfig)
0077     : runInfoToken_(esConsumes<edm::Transition::BeginRun>()),
0078       magneticField_(iConfig.getUntrackedParameter<int>("magneticField")),
0079       magneticFieldCurrentRun_(-10000) {}
0080 
0081 //
0082 // member functions
0083 //
0084 
0085 // ------------ method called on each new Event  ------------
0086 bool MagneticFieldFilter::filter(edm::Event&, const edm::EventSetup&) {
0087   return magneticField_ == magneticFieldCurrentRun_;
0088 }
0089 
0090 // ------------ method called when starting to processes a run  ------------
0091 
0092 void MagneticFieldFilter::beginRun(const edm::Run&, const edm::EventSetup& iSetup) {
0093   const auto& summary = &iSetup.getData(runInfoToken_);
0094   // convert from Tesla to kGauss (multiply with 10) and
0095   // round off to whole kGauss (add 0.5 and cast to int) as is done in
0096   // 'MagneticField::computeNominalValue()':
0097   magneticFieldCurrentRun_ = static_cast<int>(currentToField(summary->m_avg_current) * 10.0 + 0.5);
0098 }
0099 
0100 float MagneticFieldFilter::currentToField(const float& current) const {
0101   return linearCoeffCurrentToField_ * current + constantTermCurrentToField_;
0102 }
0103 
0104 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0105 void MagneticFieldFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0106   edm::ParameterSetDescription desc;
0107   desc.setComment("Filters events with a magnetic field of 'magneticField'.");
0108   desc.addUntracked<int>("magneticField", 38)->setComment("In units of kGauss (= 0.1 Tesla).");
0109   descriptions.add("magneticFieldFilter", desc);
0110 }
0111 
0112 //define this as a plug-in
0113 DEFINE_FWK_MODULE(MagneticFieldFilter);