CMS 3D CMS Logo

TrackFoldedOccupancyClient.cc
Go to the documentation of this file.
5 //-----------------------------------------------------------------------------------
7 //-----------------------------------------------------------------------------------
8 {
9  edm::LogInfo("TrackFoldedOccupancyClient") << "TrackFoldedOccupancyClient::Deleting TrackFoldedOccupancyClient ";
10  TopFolder_ = iConfig.getParameter<std::string>("FolderName");
11  quality_ = iConfig.getParameter<std::string>("TrackQuality");
12  algoName_ = iConfig.getParameter<std::string>("AlgoName");
13  state_ = iConfig.getParameter<std::string>("MeasurementState");
14  histTag_ = (state_ == "default") ? algoName_ : state_ + "_" + algoName_;
15  conf_ = iConfig;
16 }
17 
18 //-----------------------------------------------------------------------------------
20 //-----------------------------------------------------------------------------------
21 {
22  edm::LogInfo("TrackFoldedOccupancyClient") << "TrackFoldedOccupancyClient::Deleting TrackFoldedOccupancyClient ";
23 }
24 
25 //-----------------------------------------------------------------------------------
27 //-----------------------------------------------------------------------------------
28 {
29  edm::LogInfo("TrackFoldedOccupancyClient") << "TrackFoldedOccupancyClient::beginJob done";
30 }
31 
32 //-----------------------------------------------------------------------------------
34 //-----------------------------------------------------------------------------------
35 {
36  edm::LogInfo("TrackFoldedOccupancyClient") << "TrackFoldedOccupancyClient:: Begining of Run";
37 }
38 
39 //-----------------------------------------------------------------------------------
41 //-----------------------------------------------------------------------------------
42 {
43  ibooker.setCurrentFolder(TopFolder_ + "/" + quality_ + "/GeneralProperties/");
44  int Phi2DBin = conf_.getParameter<int>("Phi2DBin");
45  int Eta2DBin = conf_.getParameter<int>("Eta2DBin");
46  double EtaMin = conf_.getParameter<double>("EtaMin");
47  double EtaMax = conf_.getParameter<double>("EtaMax");
48  double PhiMin = conf_.getParameter<double>("PhiMin");
49  double PhiMax = conf_.getParameter<double>("PhiMax");
50 
51  // use the AlgoName and Quality Name
52  std::string histname = "TkEtaPhi_RelativeDifference_byFoldingmap_" + histTag_;
54  ibooker.book2D(histname, histname, Eta2DBin, EtaMin, EtaMax, Phi2DBin, PhiMin, PhiMax);
57 
58  histname = "TkEtaPhi_RelativeDifference_byFoldingmap_op_" + histTag_;
60  ibooker.book2D(histname, histname, Eta2DBin, EtaMin, EtaMax, Phi2DBin, PhiMin, PhiMax);
63 
64  histname = "TkEtaPhi_Ratio_byFoldingmap_" + histTag_;
65  TkEtaPhi_Ratio_byFoldingmap = ibooker.book2D(histname, histname, Eta2DBin, EtaMin, EtaMax, Phi2DBin, PhiMin, PhiMax);
66  TkEtaPhi_Ratio_byFoldingmap->setAxisTitle("Track #eta", 1);
67  TkEtaPhi_Ratio_byFoldingmap->setAxisTitle("Track #phi", 2);
68 
69  histname = "TkEtaPhi_Ratio_byFoldingmap_op_" + histTag_;
71  ibooker.book2D(histname, histname, Eta2DBin, EtaMin, EtaMax, Phi2DBin, PhiMin, PhiMax);
74 }
75 
76 //-----------------------------------------------------------------------------------
78 //-----------------------------------------------------------------------------------
79 {
80  edm::LogInfo("TrackFoldedOccupancyClient") << "TrackFoldedOccupancyClient::endLuminosityBlock";
81 
82  bookMEs(ibooker);
83  std::string inFolder = TopFolder_ + "/" + quality_ + "/GeneralProperties/";
84 
85  std::string hname;
86  hname = "TrackEtaPhi_";
87  MonitorElement* TrackEtaPhi = igetter.get(inFolder + hname + histTag_);
88 
89  hname = "TrackEtaPhiInverted_";
90  MonitorElement* TrackEtaPhiInverted = igetter.get(inFolder + hname + histTag_);
91 
92  hname = "TrackEtaPhiInvertedoutofphase_";
93  MonitorElement* TrackEtaPhiInvertedoutofphase = igetter.get(inFolder + hname + histTag_);
94 
95  if (TrackEtaPhi == nullptr || TrackEtaPhiInverted == nullptr || TrackEtaPhiInvertedoutofphase == nullptr) {
96  edm::LogWarning("TrackFoldedOccupancyClient") << "MEs needed for this module not found. Skipping.";
97  return;
98  }
99 
100  TkEtaPhi_Ratio_byFoldingmap->divide(TrackEtaPhi, TrackEtaPhiInverted, 1., 1., "");
101  TkEtaPhi_Ratio_byFoldingmap_op->divide(TrackEtaPhi, TrackEtaPhiInvertedoutofphase, 1., 1., "");
102 
103  int nx = TrackEtaPhi->getNbinsX();
104  int ny = TrackEtaPhi->getNbinsY();
105 
106  for (int ii = 1; ii <= nx; ii++) {
107  for (int jj = 1; jj <= ny; jj++) {
108  double Sum1 = TrackEtaPhi->getBinContent(ii, jj) + TrackEtaPhiInverted->getBinContent(ii, jj);
109  double Sum2 = TrackEtaPhi->getBinContent(ii, jj) + TrackEtaPhiInvertedoutofphase->getBinContent(ii, jj);
110 
111  double Sub1 = TrackEtaPhi->getBinContent(ii, jj) - TrackEtaPhiInverted->getBinContent(ii, jj);
112  double Sub2 = TrackEtaPhi->getBinContent(ii, jj) - TrackEtaPhiInvertedoutofphase->getBinContent(ii, jj);
113 
114  if (Sum1 == 0 || Sum2 == 0) {
117  } else {
118  double ratio1 = Sub1 / Sum1;
119  double ratio2 = Sub2 / Sum2;
122  }
123  }
124  }
125 }
126 
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
void beginJob(void) override
BeginJob.
TrackFoldedOccupancyClient(const edm::ParameterSet &ps)
Constructor.
MonitorElement * TkEtaPhi_RelativeDifference_byFoldingmap
void bookMEs(DQMStore::IBooker &ibooker_)
book MEs
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void dqmEndJob(DQMStore::IBooker &ibooker_, DQMStore::IGetter &igetter_) override
EndJob.
virtual int getNbinsY() const
get # of bins in Y-axis
ii
Definition: cuy.py:589
~TrackFoldedOccupancyClient() override
Destructor.
Log< level::Info, false > LogInfo
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
MonitorElement * TkEtaPhi_RelativeDifference_byFoldingmap_op
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:712
virtual int getNbinsX() const
get # of bins in X-axis
void beginRun(edm::Run const &run, edm::EventSetup const &eSetup) override
BeginRun.
Log< level::Warning, false > LogWarning
virtual void divide(const MonitorElement *, const MonitorElement *, double, double, const char *)
Replace entries with results of dividing num by denom.
Definition: Run.h:45
virtual double getBinContent(int binx) const
get content of bin (1-D)
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)