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 Attributes
AlignmentMonitorTemplate Class Reference
Inheritance diagram for AlignmentMonitorTemplate:
AlignmentMonitorBase

Public Member Functions

void afterAlignment (const edm::EventSetup &iSetup)
 
 AlignmentMonitorTemplate (const edm::ParameterSet &cfg)
 
void book ()
 Book or retrieve histograms; MUST be reimplemented. More...
 
void event (const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks)
 Called for each event (by "run()"): may be reimplemented. More...
 
 ~AlignmentMonitorTemplate ()
 
- 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 Attributes

TH1F * m_hist
 
TH1F * m_ihist
 
TH1F * m_otherdir
 
TH1F * m_otherdir2
 
std::map< Alignable *, TH1F * > m_residuals
 

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 28 of file AlignmentMonitorTemplate.cc.

Constructor & Destructor Documentation

AlignmentMonitorTemplate::AlignmentMonitorTemplate ( const edm::ParameterSet cfg)
inline

Definition at line 30 of file AlignmentMonitorTemplate.cc.

30 : AlignmentMonitorBase(cfg, "AlignmentMonitorTemplate") { };
AlignmentMonitorBase(const edm::ParameterSet &cfg, std::string name)
Constructor.
AlignmentMonitorTemplate::~AlignmentMonitorTemplate ( )
inline

Definition at line 31 of file AlignmentMonitorTemplate.cc.

31 {};

Member Function Documentation

void AlignmentMonitorTemplate::afterAlignment ( const edm::EventSetup iSetup)
virtual

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

Reimplemented from AlignmentMonitorBase.

Definition at line 123 of file AlignmentMonitorTemplate.cc.

References AlignmentMonitorBase::iteration(), and m_otherdir.

123  {
124  m_otherdir->Fill(iteration()); // this one will only get one fill per iteration, because it's called in afterAlignment()
125 }
void AlignmentMonitorTemplate::book ( )
virtual

Book or retrieve histograms; MUST be reimplemented.

Implements AlignmentMonitorBase.

Definition at line 54 of file AlignmentMonitorTemplate.cc.

References AlignmentParameterStore::alignables(), AlignmentMonitorBase::book1D(), m_hist, m_ihist, m_otherdir, m_otherdir2, m_residuals, mergeVDriftHistosByStation::name, AlignmentMonitorBase::pStore(), and indexGen::title.

54  {
55  m_hist = book1D("/", "hist", "hist", 10, 0.5, 10.5); // there's only one of these per job
56  m_ihist = book1D("/iterN/", "ihist", "ihist", 10, 0, 1); // there's a new one of these for each iteration
57  // "/iterN/" is a special directory name, in which the "N" gets replaced by the current iteration number.
58 
59  m_otherdir = book1D("/otherdir/", "hist", "this is a histogram in another directory", 10, 0.5, 10.5);
60  m_otherdir2 = book1D("/iterN/otherdir/", "hist", "here's one in another directory inside the iterN directories", 10, 0.5, 10.5);
61 
62  // This is a procedure that makes one histogram for each selected alignable, and puts them in the iterN directory.
63  // This is not a constant-time lookup. If you need something faster, see AlignmentMonitorMuonHIP, which has a
64  // dynamically-allocated array of TH1F*s.
65  std::vector<Alignable*> alignables = pStore()->alignables();
66  for (std::vector<Alignable*>::const_iterator it = alignables.begin(); it != alignables.end(); ++it) {
67  char name[256], title[256];
68  snprintf(name, sizeof(name), "xresid%d", (*it)->geomDetId().rawId());
69  snprintf(title, sizeof(title), "x track-hit for DetId %d", (*it)->geomDetId().rawId());
70 
71  m_residuals[*it] = book1D("/iterN/", name, title, 100, -5., 5.);
72  }
73 
74  // Important: you create TObject pointers with the "new" operator, but YOU don't delete them. They're deleted by the
75  // base class, which knows when they are no longer needed (based on whether they are in the /iterN/ directory or not).
76 
77  // For more detail, see the twiki page:
78  // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentMonitors for creating a plugin like this one
79  // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentAlgorithms#Monitoring for configuring it
80 }
AlignmentParameterStore * pStore()
std::map< Alignable *, TH1F * > m_residuals
TH1F * book1D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
const align::Alignables & alignables(void) const
get all alignables
void AlignmentMonitorTemplate::event ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const ConstTrajTrackPairCollection iTrajTracks 
)
virtual

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

Reimplemented from AlignmentMonitorBase.

Definition at line 82 of file AlignmentMonitorTemplate.cc.

References AlignableNavigator::alignableFromDetId(), TrajectoryMeasurement::backwardPredictedState(), TrajectoryStateCombiner::combine(), AlignableNavigator::detAndSubdetInMap(), TrajectoryMeasurement::forwardPredictedState(), TrackingRecHit::geographicalId(), TrackingRecHit::isValid(), AlignmentMonitorBase::iteration(), TrackingRecHit::localPosition(), TrajectoryStateOnSurface::localPosition(), m_hist, m_ihist, m_residuals, Trajectory::measurements(), Alignable::mother(), AlignmentMonitorBase::pNavigator(), TrajectoryMeasurement::recHit(), dataformats::search(), and PV3DBase< T, PVType, FrameType >::x().

82  {
83  m_hist->Fill(iteration()); // get the number of events per iteration
84  m_ihist->Fill(0.5); // get the number of events in this iteration in the central bin
85 
86  TrajectoryStateCombiner tsoscomb;
87 
88  // This is a procedure that loops over tracks/hits, calculates residuals, and fills the appropriate histogram.
89  for (ConstTrajTrackPairCollection::const_iterator it = tracks.begin(); it != tracks.end(); ++it) {
90  const Trajectory* traj = it->first;
91 // const reco::Track* track = it->second;
92  // If your tracks are refit using the producer in RecoTracker, you'll get updated reco::Track objects with
93  // each iteration, and it makes sense to make plots using these.
94  // If your tracks are refit using TrackingTools/TrackRefitter, only the Trajectories will be updated.
95  // We're working on that. I'll try to remember to change this comment when the update is ready.
96 
97  std::vector<TrajectoryMeasurement> measurements = traj->measurements();
98  for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end(); ++im) {
99  const TrajectoryMeasurement meas = *im;
100  const TransientTrackingRecHit* hit = &(*meas.recHit());
101  const DetId id = hit->geographicalId();
102 
103  if (hit->isValid() && pNavigator()->detAndSubdetInMap(id)) {
104  // Combine the forward-propagated state with the backward-propagated state to get a minimally-biased residual.
105  // This is the same procedure that is used to calculate residuals in the HIP algorithm
107 
108  // Search for our histogram using the Alignable* -> TH1F* map
109  // The "alignable = alignable->mother()" part ascends the alignable tree, because hits are on the lowest-level
110  // while our histograms may be associated with higher-level Alignables.
111  Alignable *alignable = pNavigator()->alignableFromDetId(id);
112  std::map<Alignable*, TH1F*>::const_iterator search = m_residuals.find(alignable);
113  while (search == m_residuals.end() && (alignable = alignable->mother())) search = m_residuals.find(alignable);
114 
115  if (search != m_residuals.end()) {
116  search->second->Fill(tsosc.localPosition().x() - hit->localPosition().x());
117  }
118  } // end if hit is valid
119  } // end loop over hits
120  } // end loop over tracks/trajectories
121 }
AlignableNavigator * pNavigator()
TSOS combine(const TSOS &pTsos1, const TSOS &pTsos2) const
TrajectoryStateOnSurface forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
std::map< Alignable *, TH1F * > m_residuals
ConstRecHitPointer recHit() const
DataContainer const & measurements() const
Definition: Trajectory.h:203
Definition: DetId.h:20
tuple tracks
Definition: testEve_cfg.py:39
bool isValid() const
DetId geographicalId() const
bool detAndSubdetInMap(const DetId &detid) const
Given a DetId, returns true if DetIds with this detector and subdetector id are in the map (not neces...
T x() const
Definition: PV3DBase.h:61
virtual LocalPoint localPosition() const =0
AlignableDetOrUnitPtr alignableFromDetId(const DetId &detid)
Returns AlignableDetOrUnitPtr corresponding to given DetId.
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:85
TrajectoryStateOnSurface backwardPredictedState() const
Access to backward predicted state (from smoother)

Member Data Documentation

TH1F* AlignmentMonitorTemplate::m_hist
private

Definition at line 38 of file AlignmentMonitorTemplate.cc.

Referenced by book(), and event().

TH1F * AlignmentMonitorTemplate::m_ihist
private

Definition at line 38 of file AlignmentMonitorTemplate.cc.

Referenced by book(), and event().

TH1F * AlignmentMonitorTemplate::m_otherdir
private

Definition at line 38 of file AlignmentMonitorTemplate.cc.

Referenced by afterAlignment(), and book().

TH1F * AlignmentMonitorTemplate::m_otherdir2
private

Definition at line 38 of file AlignmentMonitorTemplate.cc.

Referenced by book().

std::map<Alignable*, TH1F*> AlignmentMonitorTemplate::m_residuals
private

Definition at line 39 of file AlignmentMonitorTemplate.cc.

Referenced by book(), and event().