CMS 3D CMS Logo

DiMuonMassBiasClient.h
Go to the documentation of this file.
1 #ifndef DQMOffline_Alignment_DiMuonMassBiasClient_h
2 #define DQMOffline_Alignment_DiMuonMassBiasClient_h
3 // -*- C++ -*-
4 //
5 // Package: DQMOffline/Alignment
6 // Class : DiMuonMassBiasClient
7 //
8 // DQM class to plot di-muon mass bias in different kinematics bins
9 
10 // system includes
11 #include <string>
12 
13 // user includes
23 
24 namespace diMuonMassBias {
25 
26  struct fitOutputs {
27  public:
28  fitOutputs(const Measurement1D& bias, const Measurement1D& width) : m_bias(bias), m_width(width) {}
29 
30  // getters
31  const Measurement1D getBias() { return m_bias; }
32  const Measurement1D getWidth() { return m_width; }
33 
34  private:
37  };
38 
39  // helper functions to fill arrays from vectors
40  inline void fillArrayF(float* x, const edm::ParameterSet& cfg, const char* name) {
41  auto v = cfg.getParameter<std::vector<double>>(name);
42  assert(v.size() == 3);
43  std::copy(std::begin(v), std::end(v), x);
44  }
45 
46  inline void fillArrayI(int* x, const edm::ParameterSet& cfg, const char* name) {
47  auto v = cfg.getParameter<std::vector<int>>(name);
48  assert(v.size() == 3);
49  std::copy(std::begin(v), std::end(v), x);
50  }
51 
52  static constexpr int minimumHits = 10;
53 
54 } // namespace diMuonMassBias
55 
57 public:
60 
62  ~DiMuonMassBiasClient() override;
63 
64  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
65 
66 protected:
68  void beginJob(void) override;
69 
71  void beginRun(edm::Run const& run, edm::EventSetup const& eSetup) override;
72 
74  void dqmEndJob(DQMStore::IBooker& ibooker_, DQMStore::IGetter& igetter_) override;
75 
76 private:
78  void bookMEs(DQMStore::IBooker& ibooker);
79  void getMEsToHarvest(DQMStore::IGetter& igetter);
80  diMuonMassBias::fitOutputs fitVoigt(TH1* hist, const bool& fitBackground = false) const;
81 
82  // data members
84  const bool fitBackground_;
85  const bool debugMode_;
86 
87  float meanConfig_[3]; /* parmaeters for the fit: mean */
88  float widthConfig_[3]; /* parameters for the fit: width */
89  float sigmaConfig_[3]; /* parameters for the fit: sigma */
90 
91  // list of histograms to harvest
92  std::vector<std::string> MEtoHarvest_;
93 
94  // the histograms to be filled
95  std::map<std::string, MonitorElement*> meanProfiles_;
96  std::map<std::string, MonitorElement*> widthProfiles_;
97 
98  // the histograms than need to be fit and displays
99  std::map<std::string, MonitorElement*> harvestTargets_;
100 };
101 #endif
void getMEsToHarvest(DQMStore::IGetter &igetter)
~DiMuonMassBiasClient() override
Destructor.
const std::string TopFolder_
void beginRun(edm::Run const &run, edm::EventSetup const &eSetup) override
BeginRun.
assert(be >=bs)
std::map< std::string, MonitorElement * > harvestTargets_
std::map< std::string, MonitorElement * > meanProfiles_
std::map< std::string, MonitorElement * > widthProfiles_
const Measurement1D getWidth()
const Measurement1D getBias()
std::vector< std::string > MEtoHarvest_
void fillArrayI(int *x, const edm::ParameterSet &cfg, const char *name)
void dqmEndJob(DQMStore::IBooker &ibooker_, DQMStore::IGetter &igetter_) override
EndJob.
__shared__ Hist hist
void beginJob(void) override
BeginJob.
fitOutputs(const Measurement1D &bias, const Measurement1D &width)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static constexpr int minimumHits
DiMuonMassBiasClient(const edm::ParameterSet &ps)
Constructor.
void bookMEs(DQMStore::IBooker &ibooker)
book MEs
float x
void fillArrayF(float *x, const edm::ParameterSet &cfg, const char *name)
diMuonMassBias::fitOutputs fitVoigt(TH1 *hist, const bool &fitBackground=false) const
Definition: Run.h:45