CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/Alignment/CommonAlignmentMonitor/plugins/AlignmentMonitorTemplate.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     CommonAlignmentProducer
00004 // Class  :     AlignmentMonitorTemplate
00005 // 
00006 // Implementation:
00007 //     <Notes on implementation>
00008 //
00009 // Original Author:  Jim Pivarski
00010 //         Created:  Thu Mar 29 13:59:56 CDT 2007
00011 // $Id: AlignmentMonitorTemplate.cc,v 1.7 2011/09/04 17:08:14 mussgill Exp $
00012 //
00013 
00014 // system include files
00015 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorPluginFactory.h"
00016 #include "TH1.h" 
00017 #include "TObject.h" 
00018 // #include "PluginManager/ModuleDef.h"
00019 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorBase.h"
00020 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00021 
00022 // user include files
00023 
00024 // 
00025 // class definition
00026 // 
00027 
00028 class AlignmentMonitorTemplate: public AlignmentMonitorBase {
00029    public:
00030       AlignmentMonitorTemplate(const edm::ParameterSet& cfg): AlignmentMonitorBase(cfg, "AlignmentMonitorTemplate") { };
00031       ~AlignmentMonitorTemplate() {};
00032 
00033       void book();
00034       void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& iTrajTracks);
00035       void afterAlignment(const edm::EventSetup &iSetup);
00036 
00037    private:
00038       TH1F *m_hist, *m_ihist, *m_otherdir, *m_otherdir2;
00039       std::map<Alignable*, TH1F*> m_residuals;
00040 };
00041 
00042 //
00043 // constants, enums and typedefs
00044 //
00045 
00046 //
00047 // static data member definitions
00048 //
00049 
00050 //
00051 // member functions
00052 //
00053 
00054 void AlignmentMonitorTemplate::book() {
00055    m_hist = book1D("/", "hist", "hist", 10, 0.5, 10.5);      // there's only one of these per job
00056    m_ihist = book1D("/iterN/", "ihist", "ihist", 10, 0, 1);  // there's a new one of these for each iteration
00057    // "/iterN/" is a special directory name, in which the "N" gets replaced by the current iteration number.
00058 
00059    m_otherdir = book1D("/otherdir/", "hist", "this is a histogram in another directory", 10, 0.5, 10.5);
00060    m_otherdir2 = book1D("/iterN/otherdir/", "hist", "here's one in another directory inside the iterN directories", 10, 0.5, 10.5);
00061 
00062    // This is a procedure that makes one histogram for each selected alignable, and puts them in the iterN directory.
00063    // This is not a constant-time lookup.  If you need something faster, see AlignmentMonitorMuonHIP, which has a
00064    // dynamically-allocated array of TH1F*s.
00065    std::vector<Alignable*> alignables = pStore()->alignables();
00066    for (std::vector<Alignable*>::const_iterator it = alignables.begin();  it != alignables.end();  ++it) {
00067       char name[256], title[256];
00068       snprintf(name, sizeof(name), "xresid%d", (*it)->geomDetId().rawId());
00069       snprintf(title, sizeof(title), "x track-hit for DetId %d", (*it)->geomDetId().rawId());
00070 
00071       m_residuals[*it] = book1D("/iterN/", name, title, 100, -5., 5.);
00072    }
00073 
00074    // Important: you create TObject pointers with the "new" operator, but YOU don't delete them.  They're deleted by the
00075    // base class, which knows when they are no longer needed (based on whether they are in the /iterN/ directory or not).
00076 
00077    // For more detail, see the twiki page:
00078    // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentMonitors    for creating a plugin like this one
00079    // https://twiki.cern.ch/twiki/bin/view/CMS/SWGuideAlignmentAlgorithms#Monitoring    for configuring it
00080 }
00081 
00082 void AlignmentMonitorTemplate::event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& tracks) {
00083    m_hist->Fill(iteration());  // get the number of events per iteration
00084    m_ihist->Fill(0.5);         // get the number of events in this iteration in the central bin
00085 
00086   TrajectoryStateCombiner tsoscomb;
00087 
00088   // This is a procedure that loops over tracks/hits, calculates residuals, and fills the appropriate histogram.
00089    for (ConstTrajTrackPairCollection::const_iterator it = tracks.begin();  it != tracks.end();  ++it) {
00090       const Trajectory* traj = it->first;
00091 //      const reco::Track* track = it->second;
00092       // If your tracks are refit using the producer in RecoTracker, you'll get updated reco::Track objects with
00093       // each iteration, and it makes sense to make plots using these.
00094       // If your tracks are refit using TrackingTools/TrackRefitter, only the Trajectories will be updated.
00095       // We're working on that.  I'll try to remember to change this comment when the update is ready.
00096 
00097       std::vector<TrajectoryMeasurement> measurements = traj->measurements();
00098       for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin();  im != measurements.end();  ++im) {
00099          const TrajectoryMeasurement meas = *im;
00100          const TransientTrackingRecHit* hit = &(*meas.recHit());
00101          const DetId id = hit->geographicalId();
00102 
00103          if (hit->isValid()  &&  pNavigator()->detAndSubdetInMap(id)) {
00104             // Combine the forward-propagated state with the backward-propagated state to get a minimally-biased residual.
00105             // This is the same procedure that is used to calculate residuals in the HIP algorithm
00106             TrajectoryStateOnSurface tsosc = tsoscomb.combine(meas.forwardPredictedState(), meas.backwardPredictedState());
00107 
00108             // Search for our histogram using the Alignable* -> TH1F* map
00109             // The "alignable = alignable->mother()" part ascends the alignable tree, because hits are on the lowest-level
00110             // while our histograms may be associated with higher-level Alignables.
00111             Alignable *alignable = pNavigator()->alignableFromDetId(id);
00112             std::map<Alignable*, TH1F*>::const_iterator search = m_residuals.find(alignable);
00113             while (search == m_residuals.end()  &&  (alignable = alignable->mother())) search = m_residuals.find(alignable);
00114 
00115             if (search != m_residuals.end()) {
00116                search->second->Fill(tsosc.localPosition().x() - hit->localPosition().x());
00117             }
00118          } // end if hit is valid
00119       } // end loop over hits
00120    } // end loop over tracks/trajectories
00121 }
00122 
00123 void AlignmentMonitorTemplate::afterAlignment(const edm::EventSetup &iSetup) {
00124    m_otherdir->Fill(iteration());  // this one will only get one fill per iteration, because it's called in afterAlignment()
00125 }
00126 
00127 //
00128 // constructors and destructor
00129 //
00130 
00131 // AlignmentMonitorTemplate::AlignmentMonitorTemplate(const AlignmentMonitorTemplate& rhs)
00132 // {
00133 //    // do actual copying here;
00134 // }
00135 
00136 //
00137 // assignment operators
00138 //
00139 // const AlignmentMonitorTemplate& AlignmentMonitorTemplate::operator=(const AlignmentMonitorTemplate& rhs)
00140 // {
00141 //   //An exception safe implementation is
00142 //   AlignmentMonitorTemplate temp(rhs);
00143 //   swap(rhs);
00144 //
00145 //   return *this;
00146 // }
00147 
00148 //
00149 // const member functions
00150 //
00151 
00152 //
00153 // static member functions
00154 //
00155 
00156 //
00157 // SEAL definitions
00158 //
00159 
00160 // 
00161 // DEFINE_SEAL_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorTemplate, "AlignmentMonitorTemplate");
00162 DEFINE_EDM_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorTemplate, "AlignmentMonitorTemplate");