test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentMonitorTemplate.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CommonAlignmentProducer
4 // Class : AlignmentMonitorTemplate
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: Jim Pivarski
10 // Created: Thu Mar 29 13:59:56 CDT 2007
11 // $Id: AlignmentMonitorTemplate.cc,v 1.6 2010/02/25 00:27:56 wmtan Exp $
12 //
13 
14 // system include files
16 #include "TH1.h"
17 #include "TObject.h"
18 // #include "PluginManager/ModuleDef.h"
21 
22 // user include files
23 
24 //
25 // class definition
26 //
27 
29  public:
30  AlignmentMonitorTemplate(const edm::ParameterSet& cfg): AlignmentMonitorBase(cfg, "AlignmentMonitorTemplate") { };
32 
33  void book() override;
34  void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& iTrajTracks) override;
35  void afterAlignment(const edm::EventSetup &iSetup) override;
36 
37  private:
39  std::map<Alignable*, TH1F*> m_residuals;
40 };
41 
42 //
43 // constants, enums and typedefs
44 //
45 
46 //
47 // static data member definitions
48 //
49 
50 //
51 // member functions
52 //
53 
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 }
81 
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 }
122 
124  m_otherdir->Fill(iteration()); // this one will only get one fill per iteration, because it's called in afterAlignment()
125 }
126 
127 //
128 // constructors and destructor
129 //
130 
131 // AlignmentMonitorTemplate::AlignmentMonitorTemplate(const AlignmentMonitorTemplate& rhs)
132 // {
133 // // do actual copying here;
134 // }
135 
136 //
137 // assignment operators
138 //
139 // const AlignmentMonitorTemplate& AlignmentMonitorTemplate::operator=(const AlignmentMonitorTemplate& rhs)
140 // {
141 // //An exception safe implementation is
142 // AlignmentMonitorTemplate temp(rhs);
143 // swap(rhs);
144 //
145 // return *this;
146 // }
147 
148 //
149 // const member functions
150 //
151 
152 //
153 // static member functions
154 //
155 
156 //
157 // SEAL definitions
158 //
159 
160 //
161 // DEFINE_SEAL_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorTemplate, "AlignmentMonitorTemplate");
AlignableNavigator * pNavigator()
AlignmentParameterStore * pStore()
tuple cfg
Definition: looper.py:293
ConstRecHitPointer const & recHit() const
TSOS combine(const TSOS &pTsos1, const TSOS &pTsos2) const
std::vector< T >::const_iterator search(const cond::Time_t &val, const std::vector< T > &container)
Definition: IOVProxy.cc:282
std::map< Alignable *, TH1F * > m_residuals
AlignmentMonitorTemplate(const edm::ParameterSet &cfg)
DataContainer const & measurements() const
Definition: Trajectory.h:250
void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks) override
Called for each event (by &quot;run()&quot;): may be reimplemented.
int iEvent
Definition: GenABIO.cc:230
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
TH1F * book1D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
Definition: DetId.h:18
tuple tracks
Definition: testEve_cfg.py:39
bool isValid() const
void book() override
Book or retrieve histograms; MUST be reimplemented.
void afterAlignment(const edm::EventSetup &iSetup) override
#define DEFINE_EDM_PLUGIN(factory, type, name)
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:62
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:90
const align::Alignables & alignables(void) const
get all alignables
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)