File indexing completed on 2024-04-06 11:56:06
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
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
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
0051 float currentToField(const float& current) const;
0052
0053
0054
0055
0056 const edm::ESGetToken<RunInfo, RunInfoRcd> runInfoToken_;
0057
0058
0059 static constexpr float linearCoeffCurrentToField_ = 2.084287e-04;
0060
0061 static constexpr float constantTermCurrentToField_ = 1.704418e-02;
0062
0063 const int magneticField_;
0064 int magneticFieldCurrentRun_;
0065 };
0066
0067
0068
0069
0070 constexpr float MagneticFieldFilter::linearCoeffCurrentToField_;
0071 constexpr float MagneticFieldFilter::constantTermCurrentToField_;
0072
0073
0074
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
0083
0084
0085
0086 bool MagneticFieldFilter::filter(edm::Event&, const edm::EventSetup&) {
0087 return magneticField_ == magneticFieldCurrentRun_;
0088 }
0089
0090
0091
0092 void MagneticFieldFilter::beginRun(const edm::Run&, const edm::EventSetup& iSetup) {
0093 const auto& summary = &iSetup.getData(runInfoToken_);
0094
0095
0096
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
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
0113 DEFINE_FWK_MODULE(MagneticFieldFilter);