CMS 3D CMS Logo

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:
31  : AlignmentMonitorBase(cfg, iC, "AlignmentMonitorTemplate"){};
33 
34  void book() override;
35  void event(const edm::Event& iEvent,
36  const edm::EventSetup& iSetup,
37  const ConstTrajTrackPairCollection& iTrajTracks) override;
38  void afterAlignment() override;
39 
40 private:
42  std::map<Alignable*, TH1F*> m_residuals;
43 };
44 
45 //
46 // constants, enums and typedefs
47 //
48 
49 //
50 // static data member definitions
51 //
52 
53 //
54 // member functions
55 //
56 
58  m_hist = book1D("/", "hist", "hist", 10, 0.5, 10.5); // there's only one of these per job
59  m_ihist = book1D("/iterN/", "ihist", "ihist", 10, 0, 1); // there's a new one of these for each iteration
60  // "/iterN/" is a special directory name, in which the "N" gets replaced by the current iteration number.
61 
62  m_otherdir = book1D("/otherdir/", "hist", "this is a histogram in another directory", 10, 0.5, 10.5);
63  m_otherdir2 =
64  book1D("/iterN/otherdir/", "hist", "here's one in another directory inside the iterN directories", 10, 0.5, 10.5);
65 
66  // This is a procedure that makes one histogram for each selected alignable, and puts them in the iterN directory.
67  // This is not a constant-time lookup. If you need something faster, see AlignmentMonitorMuonHIP, which has a
68  // dynamically-allocated array of TH1F*s.
69  const auto& alignables = pStore()->alignables();
70  for (const auto& it : alignables) {
71  char name[256], title[256];
72  snprintf(name, sizeof(name), "xresid%d", it->geomDetId().rawId());
73  snprintf(title, sizeof(title), "x track-hit for DetId %d", it->geomDetId().rawId());
74 
75  m_residuals[it] = book1D("/iterN/", name, title, 100, -5., 5.);
76  }
77 
78  // Important: you create TObject pointers with the "new" operator, but YOU don't delete them. They're deleted by the
79  // base class, which knows when they are no longer needed (based on whether they are in the /iterN/ directory or not).
80 
81  // For more detail, see the twiki page:
82  // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentMonitors for creating a plugin like this one
83  // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentAlgorithms#Monitoring for configuring it
84 }
85 
87  const edm::EventSetup& iSetup,
89  m_hist->Fill(iteration()); // get the number of events per iteration
90  m_ihist->Fill(0.5); // get the number of events in this iteration in the central bin
91 
92  TrajectoryStateCombiner tsoscomb;
93 
94  // This is a procedure that loops over tracks/hits, calculates residuals, and fills the appropriate histogram.
95  for (ConstTrajTrackPairCollection::const_iterator it = tracks.begin(); it != tracks.end(); ++it) {
96  const Trajectory* traj = it->first;
97  // const reco::Track* track = it->second;
98  // If your tracks are refit using the producer in RecoTracker, you'll get updated reco::Track objects with
99  // each iteration, and it makes sense to make plots using these.
100  // If your tracks are refit using TrackingTools/TrackRefitter, only the Trajectories will be updated.
101  // We're working on that. I'll try to remember to change this comment when the update is ready.
102 
103  std::vector<TrajectoryMeasurement> measurements = traj->measurements();
104  for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end(); ++im) {
105  const TrajectoryMeasurement meas = *im;
106  const TransientTrackingRecHit* hit = &(*meas.recHit());
107  const DetId id = hit->geographicalId();
108 
109  if (hit->isValid() && pNavigator()->detAndSubdetInMap(id)) {
110  // Combine the forward-propagated state with the backward-propagated state to get a minimally-biased residual.
111  // This is the same procedure that is used to calculate residuals in the HIP algorithm
113 
114  // Search for our histogram using the Alignable* -> TH1F* map
115  // The "alignable = alignable->mother()" part ascends the alignable tree, because hits are on the lowest-level
116  // while our histograms may be associated with higher-level Alignables.
117  Alignable* alignable = pNavigator()->alignableFromDetId(id);
118  std::map<Alignable*, TH1F*>::const_iterator search = m_residuals.find(alignable);
119  while (search == m_residuals.end() && (alignable = alignable->mother()))
120  search = m_residuals.find(alignable);
121 
122  if (search != m_residuals.end()) {
123  search->second->Fill(tsosc.localPosition().x() - hit->localPosition().x());
124  }
125  } // end if hit is valid
126  } // end loop over hits
127  } // end loop over tracks/trajectories
128 }
129 
131  m_otherdir->Fill(
132  iteration()); // this one will only get one fill per iteration, because it's called in afterAlignment()
133 }
134 
135 //
136 // constructors and destructor
137 //
138 
139 // AlignmentMonitorTemplate::AlignmentMonitorTemplate(const AlignmentMonitorTemplate& rhs)
140 // {
141 // // do actual copying here;
142 // }
143 
144 //
145 // assignment operators
146 //
147 // const AlignmentMonitorTemplate& AlignmentMonitorTemplate::operator=(const AlignmentMonitorTemplate& rhs)
148 // {
149 // //An exception safe implementation is
150 // AlignmentMonitorTemplate temp(rhs);
151 // swap(rhs);
152 //
153 // return *this;
154 // }
155 
156 //
157 // const member functions
158 //
159 
160 //
161 // static member functions
162 //
163 
164 //
165 // SEAL definitions
166 //
167 
168 //
169 // DEFINE_SEAL_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorTemplate, "AlignmentMonitorTemplate");
TrajectoryStateCombiner.h
AlignableNavigator::detAndSubdetInMap
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...
Definition: AlignableNavigator.cc:193
AlignmentMonitorTemplate
Definition: AlignmentMonitorTemplate.cc:28
runGCPTkAlMap.title
string title
Definition: runGCPTkAlMap.py:94
AlignmentMonitorTemplate::AlignmentMonitorTemplate
AlignmentMonitorTemplate(const edm::ParameterSet &cfg, edm::ConsumesCollector iC)
Definition: AlignmentMonitorTemplate.cc:30
AlignmentParameterStore::alignables
const align::Alignables & alignables(void) const
get all alignables
Definition: AlignmentParameterStore.h:47
AlignmentMonitorTemplate::m_residuals
std::map< Alignable *, TH1F * > m_residuals
Definition: AlignmentMonitorTemplate.cc:42
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
Alignable
Definition: Alignable.h:27
AlignmentMonitorBase::pStore
AlignmentParameterStore * pStore()
Definition: AlignmentMonitorBase.h:114
AlignmentMonitorTemplate::m_hist
TH1F * m_hist
Definition: AlignmentMonitorTemplate.cc:41
AlignmentMonitorTemplate::book
void book() override
Book or retrieve histograms; MUST be reimplemented.
Definition: AlignmentMonitorTemplate.cc:57
AlignmentMonitorPluginFactory
AlignmentMonitorTemplate::m_otherdir
TH1F * m_otherdir
Definition: AlignmentMonitorTemplate.cc:41
AlignmentMonitorBase::pNavigator
AlignableNavigator * pNavigator()
Definition: AlignmentMonitorBase.h:115
hit::x
double x
Definition: SiStripHitEffFromCalibTree.cc:89
AlignmentMonitorTemplate::event
void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks) override
Called for each event (by "run()"): may be reimplemented.
Definition: AlignmentMonitorTemplate.cc:86
DetId
Definition: DetId.h:17
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
TrajectoryMeasurement::backwardPredictedState
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)
Definition: TrajectoryMeasurement.h:179
TrajectoryMeasurement::forwardPredictedState
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
Definition: TrajectoryMeasurement.h:177
cond::persistency::search
std::vector< T >::const_iterator search(const cond::Time_t &val, const std::vector< T > &container)
Definition: IOVProxy.cc:21
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
AlignmentMonitorPluginFactory.h
TrajectoryStateCombiner
Definition: TrajectoryStateCombiner.h:13
TrajectoryStateOnSurface::localPosition
LocalPoint localPosition() const
Definition: TrajectoryStateOnSurface.h:74
edm::ParameterSet
Definition: ParameterSet.h:47
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:176
AlignmentMonitorBase.h
AlignmentMonitorBase
Definition: AlignmentMonitorBase.h:43
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup
Definition: EventSetup.h:58
AlignmentMonitorTemplate::m_ihist
TH1F * m_ihist
Definition: AlignmentMonitorTemplate.cc:41
Trajectory::measurements
DataContainer const & measurements() const
Definition: Trajectory.h:178
looper.cfg
cfg
Definition: looper.py:296
TrackingRecHit
Definition: TrackingRecHit.h:21
AlignmentMonitorBase::book1D
TH1F * book1D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
Definition: AlignmentMonitorBase.cc:109
TrajectoryMeasurement::recHit
ConstRecHitPointer const & recHit() const
Definition: TrajectoryMeasurement.h:190
AlignmentMonitorBase::iteration
int iteration()
Definition: AlignmentMonitorBase.h:111
AlignableNavigator::alignableFromDetId
AlignableDetOrUnitPtr alignableFromDetId(const DetId &detid)
Returns AlignableDetOrUnitPtr corresponding to given DetId.
Definition: AlignableNavigator.cc:91
Trajectory
Definition: Trajectory.h:38
TrajectoryStateCombiner::combine
TSOS combine(const TSOS &pTsos1, const TSOS &pTsos2) const
Definition: TrajectoryStateCombiner.cc:6
AlignmentMonitorTemplate::~AlignmentMonitorTemplate
~AlignmentMonitorTemplate() override
Definition: AlignmentMonitorTemplate.cc:32
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
AlignmentMonitorBase::ConstTrajTrackPairCollection
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
Definition: AlignmentMonitorBase.h:46
AlignmentMonitorTemplate::afterAlignment
void afterAlignment() override
Definition: AlignmentMonitorTemplate.cc:130
AlignmentMonitorTemplate::m_otherdir2
TH1F * m_otherdir2
Definition: AlignmentMonitorTemplate.cc:41
edm::Event
Definition: Event.h:73
TrajectoryMeasurement
Definition: TrajectoryMeasurement.h:25
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
hit
Definition: SiStripHitEffFromCalibTree.cc:88
Alignable::mother
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:91