Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:34:28

0001 // -*- C++ -*-
0002 //
0003 // Package:     CommonAlignmentProducer
0004 // Class  :     AlignmentMonitorTemplate
0005 //
0006 // Implementation:
0007 //     <Notes on implementation>
0008 //
0009 // Original Author:  Jim Pivarski
0010 //         Created:  Thu Mar 29 13:59:56 CDT 2007
0011 // $Id: AlignmentMonitorTemplate.cc,v 1.6 2010/02/25 00:27:56 wmtan Exp $
0012 //
0013 
0014 // system include files
0015 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorPluginFactory.h"
0016 #include "TH1.h"
0017 #include "TObject.h"
0018 // #include "PluginManager/ModuleDef.h"
0019 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorBase.h"
0020 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
0021 
0022 // user include files
0023 
0024 //
0025 // class definition
0026 //
0027 
0028 class AlignmentMonitorTemplate : public AlignmentMonitorBase {
0029 public:
0030   AlignmentMonitorTemplate(const edm::ParameterSet& cfg, edm::ConsumesCollector iC)
0031       : AlignmentMonitorBase(cfg, iC, "AlignmentMonitorTemplate") {}
0032   ~AlignmentMonitorTemplate() override {}
0033 
0034   void book() override;
0035   void event(const edm::Event& iEvent,
0036              const edm::EventSetup& iSetup,
0037              const ConstTrajTrackPairCollection& iTrajTracks) override;
0038   void afterAlignment() override;
0039 
0040 private:
0041   TH1F *m_hist, *m_ihist, *m_otherdir, *m_otherdir2;
0042   std::map<Alignable*, TH1F*> m_residuals;
0043 };
0044 
0045 //
0046 // constants, enums and typedefs
0047 //
0048 
0049 //
0050 // static data member definitions
0051 //
0052 
0053 //
0054 // member functions
0055 //
0056 
0057 void AlignmentMonitorTemplate::book() {
0058   m_hist = book1D("/", "hist", "hist", 10, 0.5, 10.5);      // there's only one of these per job
0059   m_ihist = book1D("/iterN/", "ihist", "ihist", 10, 0, 1);  // there's a new one of these for each iteration
0060   // "/iterN/" is a special directory name, in which the "N" gets replaced by the current iteration number.
0061 
0062   m_otherdir = book1D("/otherdir/", "hist", "this is a histogram in another directory", 10, 0.5, 10.5);
0063   m_otherdir2 =
0064       book1D("/iterN/otherdir/", "hist", "here's one in another directory inside the iterN directories", 10, 0.5, 10.5);
0065 
0066   // This is a procedure that makes one histogram for each selected alignable, and puts them in the iterN directory.
0067   // This is not a constant-time lookup.  If you need something faster, see AlignmentMonitorMuonHIP, which has a
0068   // dynamically-allocated array of TH1F*s.
0069   const auto& alignables = pStore()->alignables();
0070   for (const auto& it : alignables) {
0071     char name[256], title[256];
0072     snprintf(name, sizeof(name), "xresid%d", it->geomDetId().rawId());
0073     snprintf(title, sizeof(title), "x track-hit for DetId %d", it->geomDetId().rawId());
0074 
0075     m_residuals[it] = book1D("/iterN/", name, title, 100, -5., 5.);
0076   }
0077 
0078   // Important: you create TObject pointers with the "new" operator, but YOU don't delete them.  They're deleted by the
0079   // base class, which knows when they are no longer needed (based on whether they are in the /iterN/ directory or not).
0080 
0081   // For more detail, see the twiki page:
0082   // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentMonitors    for creating a plugin like this one
0083   // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentAlgorithms#Monitoring    for configuring it
0084 }
0085 
0086 void AlignmentMonitorTemplate::event(const edm::Event& iEvent,
0087                                      const edm::EventSetup& iSetup,
0088                                      const ConstTrajTrackPairCollection& tracks) {
0089   m_hist->Fill(iteration());  // get the number of events per iteration
0090   m_ihist->Fill(0.5);         // get the number of events in this iteration in the central bin
0091 
0092   TrajectoryStateCombiner tsoscomb;
0093 
0094   // This is a procedure that loops over tracks/hits, calculates residuals, and fills the appropriate histogram.
0095   for (ConstTrajTrackPairCollection::const_iterator it = tracks.begin(); it != tracks.end(); ++it) {
0096     const Trajectory* traj = it->first;
0097     //      const reco::Track* track = it->second;
0098     // If your tracks are refit using the producer in RecoTracker, you'll get updated reco::Track objects with
0099     // each iteration, and it makes sense to make plots using these.
0100     // If your tracks are refit using TrackingTools/TrackRefitter, only the Trajectories will be updated.
0101     // We're working on that.  I'll try to remember to change this comment when the update is ready.
0102 
0103     std::vector<TrajectoryMeasurement> measurements = traj->measurements();
0104     for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end(); ++im) {
0105       const TrajectoryMeasurement& meas = *im;
0106       const TransientTrackingRecHit* hit = &(*meas.recHit());
0107       const DetId id = hit->geographicalId();
0108 
0109       if (hit->isValid() && pNavigator()->detAndSubdetInMap(id)) {
0110         // Combine the forward-propagated state with the backward-propagated state to get a minimally-biased residual.
0111         // This is the same procedure that is used to calculate residuals in the HIP algorithm
0112         TrajectoryStateOnSurface tsosc = tsoscomb.combine(meas.forwardPredictedState(), meas.backwardPredictedState());
0113 
0114         // Search for our histogram using the Alignable* -> TH1F* map
0115         // The "alignable = alignable->mother()" part ascends the alignable tree, because hits are on the lowest-level
0116         // while our histograms may be associated with higher-level Alignables.
0117         Alignable* alignable = pNavigator()->alignableFromDetId(id);
0118         std::map<Alignable*, TH1F*>::const_iterator search = m_residuals.find(alignable);
0119         while (search == m_residuals.end() && (alignable = alignable->mother()))
0120           search = m_residuals.find(alignable);
0121 
0122         if (search != m_residuals.end()) {
0123           search->second->Fill(tsosc.localPosition().x() - hit->localPosition().x());
0124         }
0125       }  // end if hit is valid
0126     }  // end loop over hits
0127   }  // end loop over tracks/trajectories
0128 }
0129 
0130 void AlignmentMonitorTemplate::afterAlignment() {
0131   m_otherdir->Fill(
0132       iteration());  // this one will only get one fill per iteration, because it's called in afterAlignment()
0133 }
0134 
0135 //
0136 // constructors and destructor
0137 //
0138 
0139 // AlignmentMonitorTemplate::AlignmentMonitorTemplate(const AlignmentMonitorTemplate& rhs)
0140 // {
0141 //    // do actual copying here;
0142 // }
0143 
0144 //
0145 // assignment operators
0146 //
0147 // const AlignmentMonitorTemplate& AlignmentMonitorTemplate::operator=(const AlignmentMonitorTemplate& rhs)
0148 // {
0149 //   //An exception safe implementation is
0150 //   AlignmentMonitorTemplate temp(rhs);
0151 //   swap(rhs);
0152 //
0153 //   return *this;
0154 // }
0155 
0156 //
0157 // const member functions
0158 //
0159 
0160 //
0161 // static member functions
0162 //
0163 
0164 //
0165 // SEAL definitions
0166 //
0167 
0168 //
0169 // DEFINE_SEAL_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorTemplate, "AlignmentMonitorTemplate");
0170 DEFINE_EDM_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorTemplate, "AlignmentMonitorTemplate");