CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MillePedeDQMModule.cc
Go to the documentation of this file.
1 
10 /*** header-file ***/
12 
13 /*** ROOT objects ***/
14 #include "TH1F.h"
15 
16 /*** Core framework functionality ***/
20 
21 /*** Geometry ***/
25 
26 /*** Alignment ***/
30 
31 
32 
35  mpReaderConfig_(
36  config.getParameter<edm::ParameterSet>("MillePedeFileReader")
37  ),
38 
39  sigCut_ (mpReaderConfig_.getParameter<double>("sigCut")),
40  Xcut_ (mpReaderConfig_.getParameter<double>("Xcut")),
41  tXcut_ (mpReaderConfig_.getParameter<double>("tXcut")),
42  Ycut_ (mpReaderConfig_.getParameter<double>("Ycut")),
43  tYcut_ (mpReaderConfig_.getParameter<double>("tYcut")),
44  Zcut_ (mpReaderConfig_.getParameter<double>("Zcut")),
45  tZcut_ (mpReaderConfig_.getParameter<double>("tZcut")),
46  maxMoveCut_ (mpReaderConfig_.getParameter<double>("maxMoveCut")),
47  maxErrorCut_ (mpReaderConfig_.getParameter<double>("maxErrorCut"))
48 {
49 }
50 
53 {
54 }
55 
56 //=============================================================================
57 //=== INTERFACE IMPLEMENTATION ===
58 //=============================================================================
59 
62 {
63  edm::LogInfo("MillePedeDQMModule") << "Booking histograms";
64 
65  booker.cd();
66  booker.setCurrentFolder("AlCaReco/SiPixelAli/");
67 
68  h_xPos = booker.book1D("Xpos", "#Delta X;;#mu m", 10, 0, 10.);
69  h_xRot = booker.book1D("Xrot", "#Delta #theta_{X};;#mu rad", 10, 0, 10.);
70  h_yPos = booker.book1D("Ypos", "#Delta Y;;#mu m", 10, 0., 10.);
71  h_yRot = booker.book1D("Yrot", "#Delta #theta_{Y};;#mu rad", 10, 0, 10.);
72  h_zPos = booker.book1D("Zpos", "#Delta Z;;#mu m", 10, 0., 10.);
73  h_zRot = booker.book1D("Zrot", "#Delta #theta_{Z};;#mu rad", 10, 0, 10.);
74 
75  booker.cd();
76 }
77 
78 
81 {
82  bookHistograms(booker);
83  if (mpReader_) {
84  mpReader_->read();
85  } else {
86  throw cms::Exception("LogicError")
87  << "@SUB=MillePedeDQMModule::dqmEndJob\n"
88  << "Try to read MillePede results before initializing MillePedeFileReader";
89  }
90  fillExpertHistos();
91 }
92 
93 
94 
95 //=============================================================================
96 //=== PRIVATE METHOD IMPLEMENTATION ===
97 //=============================================================================
98 
101 
102  if (!setupChanged(setup)) return;
103 
105  setup.get<TrackerTopologyRcd>().get(tTopo);
106  edm::ESHandle<GeometricDet> geometricDet;
107  setup.get<IdealGeometryRecord>().get(geometricDet);
109  setup.get<PTrackerParametersRcd>().get(ptp);
110 
112 
113  const auto trackerGeometry = builder.build(&(*geometricDet), *ptp, &(*tTopo));
114  tracker_ = std::make_unique<AlignableTracker>(trackerGeometry, &(*tTopo));
115 
116  const std::string labelerPlugin{"PedeLabeler"};
117  edm::ParameterSet labelerConfig{};
118  labelerConfig.addUntrackedParameter("plugin", labelerPlugin);
119  labelerConfig.addUntrackedParameter("RunRangeSelection", edm::VParameterSet{});
120 
121  std::shared_ptr<PedeLabelerBase> pedeLabeler{
123  ->create(labelerPlugin,
124  PedeLabelerBase::TopLevelAlignables(tracker_.get(), nullptr, nullptr),
125  labelerConfig)
126  };
127 
128  mpReader_ = std::make_unique<MillePedeFileReader>(mpReaderConfig_, pedeLabeler);
129 }
130 
131 
134 {
135 
136  fillExpertHisto(h_xPos, Xcut_, sigCut_, maxMoveCut_, maxErrorCut_, mpReader_->getXobs(), mpReader_->getXobsErr());
137  fillExpertHisto(h_xRot, tXcut_, sigCut_, maxMoveCut_, maxErrorCut_, mpReader_->getTXobs(), mpReader_->getTXobsErr());
138 
139  fillExpertHisto(h_yPos, Ycut_, sigCut_, maxMoveCut_, maxErrorCut_, mpReader_->getYobs(), mpReader_->getYobsErr());
140  fillExpertHisto(h_yRot, tYcut_, sigCut_, maxMoveCut_, maxErrorCut_, mpReader_->getTYobs(), mpReader_->getTYobsErr());
141 
142  fillExpertHisto(h_zPos, Zcut_, sigCut_, maxMoveCut_, maxErrorCut_, mpReader_->getZobs(), mpReader_->getZobsErr());
143  fillExpertHisto(h_zRot, tZcut_, sigCut_, maxMoveCut_, maxErrorCut_, mpReader_->getTZobs(), mpReader_->getTZobsErr());
144 
145 }
146 
148 ::fillExpertHisto(MonitorElement* histo, const double cut, const double sigCut, const double maxMoveCut, const double maxErrorCut,
149  std::array<double, 6> obs, std::array<double, 6> obsErr)
150 {
151  TH1F* histo_0 = histo->getTH1F();
152 
153  histo_0->SetMinimum(-(maxMoveCut_));
154  histo_0->SetMaximum( maxMoveCut_);
155 
156  for (size_t i = 0; i < obs.size(); ++i) {
157  histo_0->SetBinContent(i+1, obs[i]);
158  histo_0->SetBinError(i+1, obsErr[i]);
159  }
160  histo_0->SetBinContent(8,cut);
161  histo_0->SetBinContent(9,sigCut);
162  histo_0->SetBinContent(10,maxMoveCut);
163  histo_0->SetBinContent(11,maxErrorCut);
164 
165 }
166 
169 {
170  bool changed{false};
171 
172  if (watchIdealGeometryRcd_.check(setup)) changed = true;
173  if (watchTrackerTopologyRcd_.check(setup)) changed = true;
174  if (watchPTrackerParametersRcd_.check(setup)) changed = true;
175 
176  return changed;
177 }
int i
Definition: DBlmapReader.cc:9
static const char tracker_[]
void cd(void)
Definition: DQMStore.cc:269
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
bool setupChanged(const edm::EventSetup &)
virtual void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
DQM Plotter for PCL-Alignment.
TrackerGeometry * build(const GeometricDet *gd, const PTrackerParameters &ptp, const TrackerTopology *tTopo)
void bookHistograms(fwlite::EventContainer &eventCont)
void fillExpertHisto(MonitorElement *histo, const double cut, const double sigCut, const double maxMoveCut, const double maxErrorCut, std::array< double, 6 > obs, std::array< double, 6 > obsErr)
MillePedeDQMModule(const edm::ParameterSet &)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
void bookHistograms(DQMStore::IBooker &)
virtual void beginRun(const edm::Run &, const edm::EventSetup &) override
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
void addUntrackedParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:208
const T & get() const
Definition: EventSetup.h:56
TH1F * getTH1F(void) const
T get(const Candidate &c)
Definition: component.h:55
Definition: Run.h:42