CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TSOSHistogramMaker.cc
Go to the documentation of this file.
2 #include <iostream>
3 #include "TH1F.h"
4 #include "TH2F.h"
10 
12  : m_2dhistos(false),
13  m_detsels(),
14  m_selnames(),
15  m_seltitles(),
16  m_histocluslenangle(),
17  m_tsosy(),
18  m_tsosx(),
19  m_tsosxy(),
20  m_ttrhy(),
21  m_ttrhx(),
22  m_ttrhxy(),
23  m_tsosdy(),
24  m_tsosdx(),
25  m_tsosdxdy() {}
26 
28  : m_2dhistos(iConfig.getUntrackedParameter<bool>("wanted2DHistos", false)),
29  m_detsels(),
30  m_selnames(),
31  m_seltitles(),
32  m_histocluslenangle(),
33  m_tsosy(),
34  m_tsosx(),
35  m_tsosxy(),
36  m_ttrhy(),
37  m_ttrhx(),
38  m_ttrhxy(),
39  m_tsosdy(),
40  m_tsosdx(),
41  m_tsosdxdy() {
43 
44  std::vector<edm::ParameterSet> wantedsubds(iConfig.getParameter<std::vector<edm::ParameterSet> >("wantedSubDets"));
45 
46  std::cout << "selections found: " << wantedsubds.size() << std::endl;
47 
48  for (std::vector<edm::ParameterSet>::iterator ps = wantedsubds.begin(); ps != wantedsubds.end(); ++ps) {
49  m_selnames.push_back(ps->getParameter<std::string>("detLabel"));
50  if (ps->existsAs<std::string>("title")) {
51  m_seltitles.push_back(ps->getParameter<std::string>("title"));
52  } else {
53  m_seltitles.push_back(ps->getParameter<std::string>("detLabel"));
54  }
55  m_detsels.push_back(DetIdSelector(ps->getUntrackedParameter<std::vector<std::string> >("selection")));
56  }
57 
58  for (unsigned int isel = 0; isel < m_detsels.size(); ++isel) {
59  TFileDirectory subdir = tfserv->mkdir(m_selnames[isel]);
60 
61  std::string name = "tsosy_" + m_selnames[isel];
62  std::string title = "TSOS y " + m_seltitles[isel];
63  m_tsosy.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 200, -20., 20.));
64  name = "tsosx_" + m_selnames[isel];
65  title = "TSOS x " + m_seltitles[isel];
66  m_tsosx.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 200, -20., 20.));
67  if (m_2dhistos) {
68  name = "tsosxy_" + m_selnames[isel];
69  title = "TSOS y vs x " + m_seltitles[isel];
70  m_tsosxy.push_back(subdir.make<TH2F>(name.c_str(), title.c_str(), 200, -20., 20., 200, -20., 20.));
71  }
72 
73  name = "tsosprojx_" + m_selnames[isel];
74  title = "TSOS x projection " + m_seltitles[isel];
75  m_tsosprojx.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 400, -2., 2.));
76  name = "tsosprojy_" + m_selnames[isel];
77  title = "TSOS y projection " + m_seltitles[isel];
78  m_tsosprojy.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 400, -2., 2.));
79 
80  name = "ttrhy_" + m_selnames[isel];
81  title = "TT RecHit y " + m_seltitles[isel];
82  m_ttrhy.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 200, -20., 20.));
83  name = "ttrhx_" + m_selnames[isel];
84  title = "TT RecHit x " + m_seltitles[isel];
85  m_ttrhx.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 200, -20., 20.));
86  if (m_2dhistos) {
87  name = "ttrhxy_" + m_selnames[isel];
88  title = "TT RecHit y vs x " + m_seltitles[isel];
89  m_ttrhxy.push_back(subdir.make<TH2F>(name.c_str(), title.c_str(), 200, -20., 20., 200, -20., 20.));
90  }
91 
92  name = "tsosdy_" + m_selnames[isel];
93  title = "TSOS-TTRH y " + m_seltitles[isel];
94  m_tsosdy.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 200, -5., 5.));
95  name = "tsosdx_" + m_selnames[isel];
96  title = "TSOS-TTRH x " + m_seltitles[isel];
97  m_tsosdx.push_back(subdir.make<TH1F>(name.c_str(), title.c_str(), 200, -0.1, 0.1));
98  if (m_2dhistos) {
99  name = "tsosdxdy_" + m_selnames[isel];
100  title = "TSOS-TTRH dy vs dy " + m_seltitles[isel];
101  m_tsosdxdy.push_back(subdir.make<TH2F>(name.c_str(), title.c_str(), 200, -0.1, 0.1, 200, -5., 5.));
102  }
103 
104  name = "cluslenangle_" + m_selnames[isel];
105  title = "Cluster Length vs Track Angle " + m_seltitles[isel];
106  m_histocluslenangle.push_back(subdir.make<TH2F>(name.c_str(), title.c_str(), 200, -1., 1., 40, -0.5, 39.5));
107  }
108 }
109 
112  if (hit == nullptr || !hit->isValid())
113  return;
114 
115  for (unsigned int i = 0; i < m_detsels.size(); ++i) {
116  if (m_detsels[i].isSelected(hit->geographicalId())) {
117  m_tsosy[i]->Fill(tsos.localPosition().y());
118  m_tsosx[i]->Fill(tsos.localPosition().x());
119  if (m_2dhistos)
120  m_tsosxy[i]->Fill(tsos.localPosition().x(), tsos.localPosition().y());
121  m_ttrhy[i]->Fill(hit->localPosition().y());
122  m_ttrhx[i]->Fill(hit->localPosition().x());
123  if (m_2dhistos)
124  m_ttrhxy[i]->Fill(hit->localPosition().x(), hit->localPosition().y());
125  m_tsosdy[i]->Fill(tsos.localPosition().y() - hit->localPosition().y());
126  m_tsosdx[i]->Fill(tsos.localPosition().x() - hit->localPosition().x());
127  if (m_2dhistos)
128  m_tsosdxdy[i]->Fill(tsos.localPosition().x() - hit->localPosition().x(),
129  tsos.localPosition().y() - hit->localPosition().y());
130 
131  if (tsos.localDirection().z() != 0) {
132  m_tsosprojx[i]->Fill(tsos.localDirection().x() / tsos.localDirection().z());
133  m_tsosprojy[i]->Fill(tsos.localDirection().y() / tsos.localDirection().z());
134  }
135  }
136  }
137 }
std::vector< std::string > m_selnames
std::vector< TH1F * > m_tsosprojy
std::vector< TH1F * > m_tsosdx
std::vector< TH2F * > m_tsosdxdy
LocalVector localDirection() const
void fill(const TrajectoryStateOnSurface &tsos, TransientTrackingRecHit::ConstRecHitPointer hit) const
std::vector< TH1F * > m_tsosdy
T y() const
Definition: PV3DBase.h:60
std::vector< TH1F * > m_tsosprojx
std::vector< TH1F * > m_tsosy
std::vector< std::string > m_seltitles
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T z() const
Definition: PV3DBase.h:61
T * make(const Args &...args) const
make new ROOT object
std::vector< TH1F * > m_tsosx
std::vector< TH1F * > m_ttrhx
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
std::vector< DetIdSelector > m_detsels
std::vector< TH2F * > m_tsosxy
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool isSelected(const std::vector< L1HPSPFTauQualityCut > &qualityCuts, const l1t::PFCandidate &pfCand, float_t primaryVertexZ)
std::vector< TH2F * > m_histocluslenangle
std::vector< TH2F * > m_ttrhxy
std::vector< TH1F * > m_ttrhy
tuple cout
Definition: gather_cfg.py:144
T x() const
Definition: PV3DBase.h:59