CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CSCStubResolutionValidation.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
7  : CSCBaseValidation(pset) {
8  const auto& simVertex = pset.getParameter<edm::ParameterSet>("simVertex");
9  simVertexInput_ = iC.consumes<edm::SimVertexContainer>(simVertex.getParameter<edm::InputTag>("inputTag"));
10  const auto& simTrack = pset.getParameter<edm::ParameterSet>("simTrack");
11  simTrackInput_ = iC.consumes<edm::SimTrackContainer>(simTrack.getParameter<edm::InputTag>("inputTag"));
12 
13  // Initialize stub matcher
14  cscStubMatcher_ = std::make_unique<CSCStubMatcher>(pset, std::move(iC));
15 }
16 
18 
19 //create folder for resolution histograms and book them
21  for (int i = 1; i <= 10; ++i) {
22  int j = i - 1;
24 
25  //Position resolution; CLCT
26  std::string t1 = "CLCTPosRes_hs_" + cn;
27  std::string t2 = "CLCTPosRes_qs_" + cn;
28  std::string t3 = "CLCTPosRes_es_" + cn;
29 
30  iBooker.setCurrentFolder("MuonCSCDigisV/CSCDigiTask/CLCT/Resolution/");
31  posresCLCT_hs[j] = iBooker.book1D(
32  t1, cn + " CLCT Position Resolution (1/2-strip prec.); Strip_{L1T} - Strip_{SIM}; Entries", 50, -1, 1);
33  posresCLCT_qs[j] = iBooker.book1D(
34  t2, cn + " CLCT Position Resolution (1/4-strip prec.); Strip_{L1T} - Strip_{SIM}; Entries", 50, -1, 1);
35  posresCLCT_es[j] = iBooker.book1D(
36  t3, cn + " CLCT Position Resolution (1/8-strip prec.); Strip_{L1T} - Strip_{SIM}; Entries", 50, -1, 1);
37 
38  //Slope resolution; CLCT
39  std::string t4 = "CLCTBendRes_" + cn;
40 
41  bendresCLCT[j] =
42  iBooker.book1D(t4, cn + " CLCT Bend Resolution; Slope_{L1T} - Slope_{SIM}; Entries", 50, -0.5, 0.5);
43  }
44 }
45 
47  // Define handles
50 
51  // Use token to retreive event information
52  e.getByToken(simTrackInput_, sim_tracks);
53  e.getByToken(simVertexInput_, sim_vertices);
54 
55  // Initialize StubMatcher
56  cscStubMatcher_->init(e, eventSetup);
57 
58  const edm::SimTrackContainer& sim_track = *sim_tracks.product();
59  const edm::SimVertexContainer& sim_vert = *sim_vertices.product();
60 
61  // select simtracks for true muons
62  edm::SimTrackContainer sim_track_selected;
63  for (const auto& t : sim_track) {
64  if (!isSimTrackGood(t))
65  continue;
66  sim_track_selected.push_back(t);
67  }
68 
69  // Skip events with no selected simtracks
70  if (sim_track_selected.empty())
71  return;
72 
73  // Loop through good tracks, use corresponding vertex to match stubs, then fill hists of chambers where the stub appears.
74  for (const auto& t : sim_track_selected) {
75  std::vector<bool> hitCLCT(10);
76 
77  std::vector<float> delta_fhs_clct(10);
78  std::vector<float> delta_fqs_clct(10);
79  std::vector<float> delta_fes_clct(10);
80 
81  std::vector<float> dslope_clct(10);
82 
83  // Match track to stubs with appropriate vertex
84  cscStubMatcher_->match(t, sim_vert[t.vertIndex()]);
85 
86  // Store matched stubs.
87  // Key: ChamberID, Value : CSCStubDigiContainer
88  const auto& clcts = cscStubMatcher_->clcts();
89 
90  // CLCTs
91  for (auto& [id, container] : clcts) {
92  const CSCDetId cscId(id);
93  const unsigned chamberType(cscId.iChamberType());
94 
95  // get the best clct in chamber
96  const auto& clct = cscStubMatcher_->bestClctInChamber(id);
97  if (!clct.isValid())
98  continue;
99 
100  hitCLCT[chamberType - 1] = true;
101 
102  // calculate deltastrip for ME1/a. Basically, we need to subtract 64 from the CLCT key strip to
103  // compare with key strip as obtained through the fit to simhits positions.
104  int deltaStrip = 0;
105  if (cscId.station() == 1 and cscId.ring() == 4 and clct.getKeyStrip() > CSCConstants::MAX_HALF_STRIP_ME1B)
106  deltaStrip = CSCConstants::NUM_STRIPS_ME1B;
107 
108  // fractional strip
109  const float fhs_clct = clct.getFractionalStrip(2);
110  const float fqs_clct = clct.getFractionalStrip(4);
111  const float fes_clct = clct.getFractionalStrip(8);
112 
113  // in half-strips per layer
114  const float slopeHalfStrip(clct.getFractionalSlope());
115  const float slopeStrip(slopeHalfStrip / 2.);
116 
117  // get the fit hits in chamber for true value
118  float stripIntercept, stripSlope;
119  cscStubMatcher_->cscDigiMatcher()->muonSimHitMatcher()->fitHitsInChamber(id, stripIntercept, stripSlope);
120 
121  // add offset of +0.25 strips for non-ME1/1 chambers
122  const bool isME11(cscId.station() == 1 and (cscId.ring() == 4 or cscId.ring() == 1));
123  if (!isME11) {
124  stripIntercept -= 0.25;
125  }
126 
127  const float strip_csc_sh = stripIntercept;
128  const float bend_csc_sh = stripSlope;
129 
130  delta_fhs_clct[chamberType - 1] = fhs_clct - deltaStrip - strip_csc_sh;
131  delta_fqs_clct[chamberType - 1] = fqs_clct - deltaStrip - strip_csc_sh;
132  delta_fes_clct[chamberType - 1] = fes_clct - deltaStrip - strip_csc_sh;
133 
134  dslope_clct[chamberType - 1] = slopeStrip - bend_csc_sh;
135  }
136 
137  for (int i = 0; i < 10; i++) {
138  if (hitCLCT[i]) {
139  //fill histograms
140  posresCLCT_hs[i]->Fill(delta_fhs_clct[i]);
141  posresCLCT_qs[i]->Fill(delta_fqs_clct[i]);
142  posresCLCT_es[i]->Fill(delta_fes_clct[i]);
143 
144  bendresCLCT[i]->Fill(dslope_clct[i]);
145  }
146  }
147  }
148 }
CSCStubResolutionValidation(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
void bookHistograms(DQMStore::IBooker &)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
bool isSimTrackGood(const SimTrack &t) const
void Fill(long long x)
edm::EDGetTokenT< edm::SimVertexContainer > simVertexInput_
def move
Definition: eostools.py:511
unsigned short iChamberType() const
Definition: CSCDetId.h:96
void analyze(const edm::Event &, const edm::EventSetup &) override
int ring() const
Definition: CSCDetId.h:68
std::string chamberName() const
Definition: CSCDetId.cc:92
edm::EDGetTokenT< edm::SimTrackContainer > simTrackInput_
T const * product() const
Definition: Handle.h:70
std::vector< SimVertex > SimVertexContainer
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< CSCStubMatcher > cscStubMatcher_
int station() const
Definition: CSCDetId.h:79
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
std::vector< SimTrack > SimTrackContainer