CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MagneticFieldFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/CommonAlignment
4 // Class: MagneticFieldFilter
5 //
20 //
21 // Original Author: Gregor Mittag
22 // Created: Wed, 25 Nov 2015 12:59:02 GMT
23 //
24 //
25 
26 
27 // user include files
35 
36 //
37 // class declaration
38 //
39 
41 public:
42  explicit MagneticFieldFilter(const edm::ParameterSet&);
43  ~MagneticFieldFilter() = default;
44 
46 
47 private:
48  virtual bool filter(edm::Event&, const edm::EventSetup&) override;
49  virtual void beginRun(const edm::Run&, const edm::EventSetup&) override;
50 
52  float currentToField(const float& current) const;
53 
54  // ----------member data ---------------------------
56  static constexpr float linearCoeffCurrentToField_ = 2.084287e-04;
58  static constexpr float constantTermCurrentToField_ = 1.704418e-02;
59 
60  const int magneticField_;
62 };
63 
64 //
65 // static data member definitions
66 //
69 
70 
71 //
72 // constructor
73 //
75  magneticField_(iConfig.getUntrackedParameter<int>("magneticField")),
76  magneticFieldCurrentRun_(-10000) {
77 }
78 
79 
80 //
81 // member functions
82 //
83 
84 // ------------ method called on each new Event ------------
85 bool
88 }
89 
90 // ------------ method called when starting to processes a run ------------
91 
92 void
95  iSetup.get<RunInfoRcd>().get(sum);
96  auto summary = sum.product();
97  // convert from Tesla to kGauss (multiply with 10) and
98  // round off to whole kGauss (add 0.5 and cast to int) as is done in
99  // 'MagneticField::computeNominalValue()':
101  static_cast<int>(currentToField(summary->m_avg_current)*10.0 + 0.5);
102 }
103 
104 
105 float
108 }
109 
110 
111 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
112 void
115  desc.setComment("Filters events with a magnetic field of 'magneticField'.");
116  desc.addUntracked<int>("magneticField", 38)
117  ->setComment("In units of kGauss (= 0.1 Tesla).");
118  descriptions.add("magneticFieldFilter", desc);
119 }
120 
121 //define this as a plug-in
static float constantTermCurrentToField_
see: https://hypernews.cern.ch/HyperNews/CMS/get/magnetic-field/63/1/1/1.html
static float linearCoeffCurrentToField_
see: https://hypernews.cern.ch/HyperNews/CMS/get/magnetic-field/63/1/1/1.html
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
#define constexpr
int magneticFieldCurrentRun_
magnetic field that is filtered
void setComment(std::string const &value)
float currentToField(const float &current) const
convert Ampere (A) to Tesla (T)
~MagneticFieldFilter()=default
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual bool filter(edm::Event &, const edm::EventSetup &) override
static void fillDescriptions(edm::ConfigurationDescriptions &)
MagneticFieldFilter(const edm::ParameterSet &)
Definition: Run.h:43