CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Attributes
AlignmentMonitorMuonVsCurvature Class Reference
Inheritance diagram for AlignmentMonitorMuonVsCurvature:
AlignmentMonitorBase

Public Member Functions

void afterAlignment (const edm::EventSetup &iSetup) override
 
 AlignmentMonitorMuonVsCurvature (const edm::ParameterSet &cfg)
 
void book () override
 Book or retrieve histograms; MUST be reimplemented. More...
 
void event (const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks) override
 Called for each event (by "run()"): may be reimplemented. More...
 
void processMuonResidualsFromTrack (MuonResidualsFromTrack &mrft, const Trajectory *traj=NULL)
 
virtual ~AlignmentMonitorMuonVsCurvature ()
 
- Public Member Functions inherited from AlignmentMonitorBase
 AlignmentMonitorBase (const edm::ParameterSet &cfg, std::string name)
 Constructor. More...
 
void beginOfJob (AlignableTracker *pTracker, AlignableMuon *pMuon, AlignmentParameterStore *pStore)
 Called at beginning of job: don't reimplement. More...
 
void duringLoop (const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks)
 Called for each event: don't reimplement. More...
 
void endOfJob ()
 Called at end of processing: don't implement. More...
 
void endOfLoop (const edm::EventSetup &iSetup)
 Called at end of loop: don't reimplement. More...
 
void startingNewLoop ()
 Called at beginning of loop: don't reimplement. More...
 
virtual ~AlignmentMonitorBase ()
 Destructor. More...
 

Private Types

enum  { kDeltaX = 0, kDeltaDxDz, kNumComponents }
 

Private Attributes

bool m_allowTIDTEC
 
bool m_doCSC
 
bool m_doDT
 
int m_layer
 
double m_maxDxy
 
double m_maxTrackerRedChi2
 
int m_minCSCHits
 
int m_minDT13Hits
 
int m_minDT2Hits
 
bool m_minNCrossedChambers
 
int m_minTrackerHits
 
double m_minTrackP
 
double m_minTrackPt
 
edm::InputTag m_muonCollectionTag
 
std::string m_propagator
 
TH1F * th1f_trackerRedChi2
 
TH1F * th1f_trackerRedChi2Diff
 
TH2F * th2f_st_ring_chamber [8][3][36][kNumComponents]
 
TH2F * th2f_wheel_st_sector [5][4][14][kNumComponents]
 
TProfile * tprofile_st_ring_chamber [8][3][36][kNumComponents]
 
TProfile * tprofile_wheel_st_sector [5][4][14][kNumComponents]
 

Additional Inherited Members

- Public Types inherited from AlignmentMonitorBase
typedef std::pair< const
Trajectory *, const
reco::Track * > 
ConstTrajTrackPair
 
typedef std::vector
< ConstTrajTrackPair
ConstTrajTrackPairCollection
 
- Protected Member Functions inherited from AlignmentMonitorBase
TH1F * book1D (std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
 
TH2F * book2D (std::string dir, std::string name, std::string title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
 
TProfile * bookProfile (std::string dir, std::string name, std::string title, int nchX, double lowX, double highX, int nchY=1, double lowY=0., double highY=0., const char *option="s")
 
TFileDirectorydirectory (std::string dir)
 
int iteration ()
 
AlignableMuonpMuon ()
 
AlignableNavigatorpNavigator ()
 
AlignmentParameterStorepStore ()
 
AlignableTrackerpTracker ()
 
- Protected Attributes inherited from AlignmentMonitorBase
const edm::InputTag m_beamSpotTag
 

Detailed Description

Definition at line 34 of file AlignmentMonitorMuonVsCurvature.cc.

Member Enumeration Documentation

anonymous enum
private

Constructor & Destructor Documentation

AlignmentMonitorMuonVsCurvature::AlignmentMonitorMuonVsCurvature ( const edm::ParameterSet cfg)

Definition at line 84 of file AlignmentMonitorMuonVsCurvature.cc.

85 : AlignmentMonitorBase(cfg, "AlignmentMonitorMuonVsCurvature")
86 , m_muonCollectionTag(cfg.getParameter<edm::InputTag>("muonCollectionTag"))
87 , m_minTrackPt(cfg.getParameter<double>("minTrackPt"))
88 , m_minTrackP(cfg.getParameter<double>("minTrackP"))
89 , m_minTrackerHits(cfg.getParameter<int>("minTrackerHits"))
90 , m_maxTrackerRedChi2(cfg.getParameter<double>("maxTrackerRedChi2"))
91 , m_allowTIDTEC(cfg.getParameter<bool>("allowTIDTEC"))
92 , m_minNCrossedChambers(cfg.getParameter<int>("minNCrossedChambers"))
93 , m_maxDxy(cfg.getParameter<double>("maxDxy"))
94 , m_minDT13Hits(cfg.getParameter<int>("minDT13Hits"))
95 , m_minDT2Hits(cfg.getParameter<int>("minDT2Hits"))
96 , m_minCSCHits(cfg.getParameter<int>("minCSCHits"))
97 , m_layer(cfg.getParameter<int>("layer"))
98 , m_propagator(cfg.getParameter<std::string>("propagator"))
99 , m_doDT(cfg.getParameter<bool>("doDT"))
100 , m_doCSC(cfg.getParameter<bool>("doCSC"))
101 {}
T getParameter(std::string const &) const
AlignmentMonitorBase(const edm::ParameterSet &cfg, std::string name)
Constructor.
virtual AlignmentMonitorMuonVsCurvature::~AlignmentMonitorMuonVsCurvature ( )
inlinevirtual

Definition at line 38 of file AlignmentMonitorMuonVsCurvature.cc.

38 {}

Member Function Documentation

void AlignmentMonitorMuonVsCurvature::afterAlignment ( const edm::EventSetup iSetup)
inlineoverridevirtual

Called after updating AlignableTracker and AlignableMuon (by "endOfLoop()"): may be reimplemented

Reimplemented from AlignmentMonitorBase.

Definition at line 45 of file AlignmentMonitorMuonVsCurvature.cc.

45 {}
void AlignmentMonitorMuonVsCurvature::book ( )
overridevirtual

Book or retrieve histograms; MUST be reimplemented.

Implements AlignmentMonitorBase.

Definition at line 104 of file AlignmentMonitorMuonVsCurvature.cc.

References AlignmentMonitorBase::book1D(), AlignmentMonitorBase::book2D(), AlignmentMonitorBase::bookProfile(), kDeltaDxDz, kDeltaX, kNumComponents, m_doCSC, m_doDT, m_minTrackPt, relativeConstraints::ring, relativeConstraints::station, AlCaHLTBitMon_QueryRunRegistry::string, th1f_trackerRedChi2, th1f_trackerRedChi2Diff, th2f_st_ring_chamber, th2f_wheel_st_sector, tprofile_st_ring_chamber, and tprofile_wheel_st_sector.

105 {
106  // DT
107  std::string wheelname[5] = {"wheelm2_", "wheelm1_", "wheelz_", "wheelp1_", "wheelp2_"};
108  if (m_doDT)
109  for (int wheel = -2; wheel <=2 ; wheel++)
110  for (int station = 1; station <= 4; station++)
111  for (int sector = 1; sector <= 14; sector++)
112  {
113  if (station != 4 && sector > 12) continue;
114 
115  char stationname[20];
116  sprintf(stationname,"st%d_", station);
117 
118  char sectorname[20];
119  sprintf(sectorname,"sector%02d_", sector);
120 
121  for (int component = 0; component < kNumComponents; component++)
122  {
123  std::stringstream th2f_name, tprofile_name;
124  th2f_name << "th2f_" << wheelname[wheel+2] <<stationname<< sectorname;
125  tprofile_name << "tprofile_" << wheelname[wheel+2] <<stationname<< sectorname;
126 
127  double yminmax = 50., xminmax = 0.05;
128  if (m_minTrackPt>0.) xminmax = 1./m_minTrackPt;
129  int ynbins = 50;
130  if (component == kDeltaX) {
131  th2f_name << "deltax";
132  tprofile_name << "deltax";
133  }
134  else if (component == kDeltaDxDz) {
135  th2f_name << "deltadxdz";
136  tprofile_name << "deltadxdz";
137  }
138 
139  th2f_wheel_st_sector[wheel+2][station-1][sector-1][component] =
140  book2D("/iterN/", th2f_name.str().c_str(), "", 30, -xminmax, xminmax, ynbins, -yminmax, yminmax);
141  tprofile_wheel_st_sector[wheel+2][station-1][sector-1][component] =
142  bookProfile("/iterN/", tprofile_name.str().c_str(), "", 30, -xminmax, xminmax);
143  }
144  }
145 
146  // CSC
147  std::string stname[8] = {"Ep_S1_", "Ep_S2_", "Ep_S3_", "Ep_S4_", "Em_S1_", "Em_S2_", "Em_S3_", "Em_S4_"};
148  if (m_doCSC)
149  for (int station = 0; station < 8; station++)
150  for (int ring = 1; ring <= 3; ring++)
151  for (int chamber = 1; chamber <= 36; chamber++)
152  {
153  int st = station%4+1;
154  if (st > 1 && ring > 2) continue; // only station 1 has more then 2 rings
155  if (st > 1 && ring == 1 && chamber > 18) continue; // ring 1 stations 1,2,3 have 18 chambers
156 
157  char ringname[20];
158  sprintf(ringname,"R%d_", ring);
159 
160  char chname[20];
161  sprintf(chname,"C%02d_", chamber);
162 
163  for (int component = 0; component < kNumComponents; component++)
164  {
165  std::stringstream componentname;
166  double yminmax = 50., xminmax = 0.05;
167  if (m_minTrackPt>0.) xminmax = 1./m_minTrackPt;
168  if (ring == 1) xminmax *= 0.5;
169  if (component == kDeltaX) {
170  componentname << "deltax";
171  }
172  else if (component == kDeltaDxDz) {
173  componentname << "deltadxdz";
174  }
175 
176  std::stringstream th2f_name, tprofile_name;
177  th2f_name << "th2f_" << stname[station] << ringname << chname << componentname.str();
178  tprofile_name << "tprofile_" << stname[station] << ringname << chname << componentname.str();
179 
180  th2f_st_ring_chamber[station][ring-1][chamber-1][component] =
181  book2D("/iterN/", th2f_name.str().c_str(), "", 30, -xminmax, xminmax, 100, -yminmax, yminmax);
182  tprofile_st_ring_chamber[station][ring-1][chamber-1][component] =
183  bookProfile("/iterN/", tprofile_name.str().c_str(), "", 30, -xminmax, xminmax);
184  }
185  }
186 
187  th1f_trackerRedChi2 = book1D("/iterN/", "trackerRedChi2", "Refit tracker reduced chi^2", 100, 0., 30.);
188  th1f_trackerRedChi2Diff = book1D("/iterN/", "trackerRedChi2Diff", "Fit-minus-refit tracker reduced chi^2", 100, -5., 5.);
189 }
TH2F * book2D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
TH2F * th2f_st_ring_chamber[8][3][36][kNumComponents]
TProfile * tprofile_st_ring_chamber[8][3][36][kNumComponents]
TProfile * bookProfile(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX, int nchY=1, double lowY=0., double highY=0., const char *option="s")
TProfile * tprofile_wheel_st_sector[5][4][14][kNumComponents]
TH1F * book1D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
TH2F * th2f_wheel_st_sector[5][4][14][kNumComponents]
void AlignmentMonitorMuonVsCurvature::event ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const ConstTrajTrackPairCollection iTrajTracks 
)
overridevirtual

Called for each event (by "run()"): may be reimplemented.

Reimplemented from AlignmentMonitorBase.

Definition at line 192 of file AlignmentMonitorMuonVsCurvature.cc.

References SiPixelRawToDigiRegional_cfi::beamSpot, reco::TrackBase::dxy(), edm::EventSetup::get(), edm::Event::getByLabel(), edm::InputTag::label(), AlignmentMonitorBase::m_beamSpotTag, m_maxDxy, m_minTrackP, m_minTrackPt, m_muonCollectionTag, HLT_FULL_cff::magneticField, metsig::muon, patZpeak::muons, reco::TrackBase::p(), AlignmentMonitorBase::pNavigator(), processMuonResidualsFromTrack(), EnergyCorrector::pt, and reco::TrackBase::pt().

Referenced by Types.EventID::cppID().

193 {
195  iSetup.get<GlobalTrackingGeometryRecord>().get(globalGeometry);
196 
198  iEvent.getByLabel(m_beamSpotTag, beamSpot);
199 
200  edm::ESHandle<DetIdAssociator> muonDetIdAssociator_;
201  iSetup.get<DetIdAssociatorRecord>().get("MuonDetIdAssociator", muonDetIdAssociator_);
202 
204  iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny",prop);
205 
207  iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
208 
209  if (m_muonCollectionTag.label().empty()) // use trajectories
210  {
211  for (ConstTrajTrackPairCollection::const_iterator trajtrack = trajtracks.begin(); trajtrack != trajtracks.end(); ++trajtrack)
212  {
213  const Trajectory* traj = (*trajtrack).first;
214  const reco::Track* track = (*trajtrack).second;
215 
216  if (track->pt() > m_minTrackPt && track->p() > m_minTrackP && fabs(track->dxy(beamSpot->position())) < m_maxDxy )
217  {
218  MuonResidualsFromTrack muonResidualsFromTrack(iSetup, magneticField, globalGeometry, muonDetIdAssociator_, prop, traj, track, pNavigator(), 1000.);
219  processMuonResidualsFromTrack(muonResidualsFromTrack, traj );
220  } // end if track pT is within range
221  } // end loop over tracks
222  }
223  else
224  {
226  iEvent.getByLabel(m_muonCollectionTag, muons);
227 
228  for (reco::MuonCollection::const_iterator muon = muons->begin(); muon != muons->end(); ++muon)
229  {
230  if ( !(muon->isTrackerMuon() && muon->innerTrack().isNonnull() ) ) continue;
231 
232  if (m_minTrackPt < muon->pt() && m_minTrackP < muon->p() && fabs(muon->innerTrack()->dxy(beamSpot->position())) < m_maxDxy)
233  {
234  MuonResidualsFromTrack muonResidualsFromTrack(globalGeometry, &(*muon), pNavigator(), 100.);
235  processMuonResidualsFromTrack(muonResidualsFromTrack);
236  }
237  }
238  }
239 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:610
AlignableNavigator * pNavigator()
tuple magneticField
double pt() const
track transverse momentum
Definition: TrackBase.h:616
const edm::InputTag m_beamSpotTag
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:413
const T & get() const
Definition: EventSetup.h:56
std::string const & label() const
Definition: InputTag.h:36
tuple muons
Definition: patZpeak.py:38
void processMuonResidualsFromTrack(MuonResidualsFromTrack &mrft, const Trajectory *traj=NULL)
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:586
void AlignmentMonitorMuonVsCurvature::processMuonResidualsFromTrack ( MuonResidualsFromTrack mrft,
const Trajectory traj = NULL 
)

Definition at line 242 of file AlignmentMonitorMuonVsCurvature.cc.

References MuonResidualsFromTrack::chamberIds(), MuonResidualsFromTrack::chamberResidual(), reco::TrackBase::charge(), MuonResidualsFromTrack::contains_TIDTEC(), MuonSubdetId::CSC, MuonSubdetId::DT, MuonResidualsFromTrack::getTrack(), MuonChamberResidual::kCSC, kDeltaDxDz, kDeltaX, MuonChamberResidual::kDT13, m_allowTIDTEC, m_doCSC, m_doDT, m_maxTrackerRedChi2, m_minCSCHits, m_minDT13Hits, m_minNCrossedChambers, m_minTrackerHits, DetId::Muon, MuonResidualsFromTrack::normalizedChi2(), reco::TrackBase::normalizedChi2(), NULL, reco::TrackBase::pt(), reco::TrackBase::pz(), relativeConstraints::ring, relativeConstraints::station, th1f_trackerRedChi2, th1f_trackerRedChi2Diff, th2f_st_ring_chamber, th2f_wheel_st_sector, tprofile_st_ring_chamber, tprofile_wheel_st_sector, MuonResidualsFromTrack::trackerNumHits(), and MuonResidualsFromTrack::trackerRedChi2().

Referenced by event().

243 {
244  if (mrft.trackerNumHits() < m_minTrackerHits) return;
245  if (!m_allowTIDTEC && mrft.contains_TIDTEC()) return;
246 
247  int nMuChambers = 0;
248  std::vector<DetId> chamberIds = mrft.chamberIds();
249  for (unsigned ch=0; ch < chamberIds.size(); ch++) if (chamberIds[ch].det() == DetId::Muon) nMuChambers++;
250  if (nMuChambers < m_minNCrossedChambers ) return;
251 
252  th1f_trackerRedChi2->Fill(mrft.trackerRedChi2());
254 
255  if (mrft.normalizedChi2() > m_maxTrackerRedChi2) return;
256 
257  double qoverpt = mrft.getTrack()->charge() / mrft.getTrack()->pt();
258  double qoverpz = 0.;
259  if (fabs(mrft.getTrack()->pz()) > 0.01) qoverpz = mrft.getTrack()->charge() / fabs(mrft.getTrack()->pz());
260 
261  for (std::vector<DetId>::const_iterator chamberId = chamberIds.begin(); chamberId != chamberIds.end(); ++chamberId)
262  {
263  if (chamberId->det() != DetId::Muon ) continue;
264 
265  if (m_doDT && chamberId->subdetId() == MuonSubdetId::DT)
266  {
267  DTChamberId dtid(chamberId->rawId());
269 
270  if (dt13 != NULL && dt13->numHits() >= m_minDT13Hits)
271  {
272  int wheel = dtid.wheel() + 2;
273  int station = dtid.station() -1;
274  int sector = dtid.sector() - 1;
275 
276  double resid_x = 10. * dt13->global_residual();
277  double resid_dxdz = 1000. * dt13->global_resslope();
278 
279  if (fabs(resid_x) < 100. && fabs(resid_dxdz) < 100.)
280  {
281  th2f_wheel_st_sector[wheel][station][sector][kDeltaX]->Fill(qoverpt, resid_x);
282  tprofile_wheel_st_sector[wheel][station][sector][kDeltaX]->Fill(qoverpt, resid_x);
283  th2f_wheel_st_sector[wheel][station][sector][kDeltaDxDz]->Fill(qoverpt, resid_dxdz);
284  tprofile_wheel_st_sector[wheel][station][sector][kDeltaDxDz]->Fill(qoverpt, resid_dxdz);
285  }
286  } // if it's a good segment
287  } // if DT
288 
289  if (m_doCSC && chamberId->subdetId() == MuonSubdetId::CSC)
290  {
291  CSCDetId cscid(chamberId->rawId());
293 
294  if (csc != NULL && csc->numHits() >= m_minCSCHits)
295  {
296  int station = 4*cscid.endcap() + cscid.station() - 5;
297  int ring = cscid.ring() - 1;
298  if (cscid.station()==1 && cscid.ring()==4) ring = 0; // join ME1/a to ME1/b
299  int chamber = cscid.chamber() - 1;
300 
301  double resid_x = 10. * csc->global_residual();
302  double resid_dxdz = 1000. * csc->global_resslope();
303 
304  if (fabs(resid_x) < 100. && fabs(resid_dxdz) < 100.)
305  {
306  th2f_st_ring_chamber[station][ring][chamber][kDeltaX]->Fill(qoverpz, resid_x);
307  tprofile_st_ring_chamber[station][ring][chamber][kDeltaX]->Fill(qoverpz, resid_x);
308  th2f_st_ring_chamber[station][ring][chamber][kDeltaDxDz]->Fill(qoverpz, resid_dxdz);
309  tprofile_st_ring_chamber[station][ring][chamber][kDeltaDxDz]->Fill(qoverpz, resid_dxdz);
310  }
311  } // if it's a good segment
312  } // if CSC
313 
314  } // end loop over chamberIds
315 }
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:556
TH2F * th2f_st_ring_chamber[8][3][36][kNumComponents]
#define NULL
Definition: scimark2.h:8
TProfile * tprofile_st_ring_chamber[8][3][36][kNumComponents]
const std::vector< DetId > chamberIds() const
static const int CSC
Definition: MuonSubdetId.h:13
double pt() const
track transverse momentum
Definition: TrackBase.h:616
const reco::Track * getTrack()
TProfile * tprofile_wheel_st_sector[5][4][14][kNumComponents]
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:634
int charge() const
track electric charge
Definition: TrackBase.h:562
static const int DT
Definition: MuonSubdetId.h:12
TH2F * th2f_wheel_st_sector[5][4][14][kNumComponents]
MuonChamberResidual * chamberResidual(DetId chamberId, int type)

Member Data Documentation

bool AlignmentMonitorMuonVsCurvature::m_allowTIDTEC
private

Definition at line 54 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by processMuonResidualsFromTrack().

bool AlignmentMonitorMuonVsCurvature::m_doCSC
private

Definition at line 63 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().

bool AlignmentMonitorMuonVsCurvature::m_doDT
private

Definition at line 62 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().

int AlignmentMonitorMuonVsCurvature::m_layer
private

Definition at line 60 of file AlignmentMonitorMuonVsCurvature.cc.

double AlignmentMonitorMuonVsCurvature::m_maxDxy
private

Definition at line 56 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by event().

double AlignmentMonitorMuonVsCurvature::m_maxTrackerRedChi2
private

Definition at line 53 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by processMuonResidualsFromTrack().

int AlignmentMonitorMuonVsCurvature::m_minCSCHits
private

Definition at line 59 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by processMuonResidualsFromTrack().

int AlignmentMonitorMuonVsCurvature::m_minDT13Hits
private

Definition at line 57 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by processMuonResidualsFromTrack().

int AlignmentMonitorMuonVsCurvature::m_minDT2Hits
private

Definition at line 58 of file AlignmentMonitorMuonVsCurvature.cc.

bool AlignmentMonitorMuonVsCurvature::m_minNCrossedChambers
private

Definition at line 55 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by processMuonResidualsFromTrack().

int AlignmentMonitorMuonVsCurvature::m_minTrackerHits
private

Definition at line 52 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by processMuonResidualsFromTrack().

double AlignmentMonitorMuonVsCurvature::m_minTrackP
private

Definition at line 51 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by event().

double AlignmentMonitorMuonVsCurvature::m_minTrackPt
private

Definition at line 50 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and event().

edm::InputTag AlignmentMonitorMuonVsCurvature::m_muonCollectionTag
private

Definition at line 49 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by event().

std::string AlignmentMonitorMuonVsCurvature::m_propagator
private

Definition at line 61 of file AlignmentMonitorMuonVsCurvature.cc.

TH1F* AlignmentMonitorMuonVsCurvature::th1f_trackerRedChi2
private

Definition at line 79 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().

TH1F* AlignmentMonitorMuonVsCurvature::th1f_trackerRedChi2Diff
private

Definition at line 80 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().

TH2F* AlignmentMonitorMuonVsCurvature::th2f_st_ring_chamber[8][3][36][kNumComponents]
private

Definition at line 76 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().

TH2F* AlignmentMonitorMuonVsCurvature::th2f_wheel_st_sector[5][4][14][kNumComponents]
private

Definition at line 72 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().

TProfile* AlignmentMonitorMuonVsCurvature::tprofile_st_ring_chamber[8][3][36][kNumComponents]
private

Definition at line 77 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().

TProfile* AlignmentMonitorMuonVsCurvature::tprofile_wheel_st_sector[5][4][14][kNumComponents]
private

Definition at line 73 of file AlignmentMonitorMuonVsCurvature.cc.

Referenced by book(), and processMuonResidualsFromTrack().