CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Alignment/CommonAlignmentMonitor/plugins/AlignmentMonitorGeneric.cc

Go to the documentation of this file.
00001 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00002 #include "Alignment/CommonAlignment/interface/AlignableObjectId.h"
00003 
00004 #include "Alignment/CommonAlignmentMonitor/plugins/AlignmentMonitorGeneric.h"
00005 #include <DataFormats/GeometrySurface/interface/LocalError.h> 
00006 #include "TObject.h" 
00007 
00008 #include <TString.h>
00009 
00010 AlignmentMonitorGeneric::AlignmentMonitorGeneric(const edm::ParameterSet& cfg):
00011   AlignmentMonitorBase(cfg, "AlignmentMonitorGeneric")
00012 {
00013 }
00014 
00015 void AlignmentMonitorGeneric::book()
00016 {
00017   static AlignableObjectId idMap;
00018 
00019   std::vector<std::string> residNames; // names of residual histograms
00020 
00021   residNames.push_back("x hit residuals pos track");
00022   residNames.push_back("x hit residuals neg track");
00023   residNames.push_back("y hit residuals pos track");
00024   residNames.push_back("y hit residuals neg track");
00025 
00026   const std::vector<Alignable*>& alignables = pStore()->alignables();
00027 
00028   unsigned int nAlignable = alignables.size();
00029   unsigned int nResidName = residNames.size();
00030 
00031   for (unsigned int i = 0; i < nAlignable; ++i)
00032   {
00033     const Alignable* ali = alignables[i];
00034 
00035     Hist1Ds& hists = m_resHists[ali];
00036 
00037     hists.resize(nResidName, 0);
00038 
00039     align::ID id = ali->id();
00040     align::StructureType type = ali->alignableObjectId();
00041 
00042     for (unsigned int n = 0; n < nResidName; ++n)
00043     {
00044       const std::string& name = residNames[n];
00045 
00046       TString histName(name.c_str());
00047       histName += Form("_%s_%d", idMap.typeToName(type).c_str(), id);
00048       histName.ReplaceAll(" ", "");
00049 
00050       TString histTitle(name.c_str());
00051       histTitle += Form(" for %s with ID %d (subdet %d)",
00052                         idMap.typeToName(type).c_str(),
00053                         id, DetId(id).subdetId());
00054 
00055       hists[n] = book1D(std::string("/iterN/") + std::string(name) + std::string("/"), std::string(histName), std::string(histTitle), nBin_, -5., 5.);
00056     }
00057   }
00058 
00059   m_trkHists.resize(6, 0);
00060   
00061   m_trkHists[0] = book1D("/iterN/", "pt"  , "track p_{t} (GeV)" , nBin_,   0.0,100.0);
00062   m_trkHists[1] = book1D("/iterN/", "eta" , "track #eta"        , nBin_, - 3.0,  3.0);
00063   m_trkHists[2] = book1D("/iterN/", "phi" , "track #phi"        , nBin_, -M_PI, M_PI);
00064   m_trkHists[3] = book1D("/iterN/", "d0"  , "track d0 (cm)"     , nBin_, -0.02, 0.02);
00065   m_trkHists[4] = book1D("/iterN/", "dz"  , "track dz (cm)"     , nBin_, -20.0, 20.0);
00066   m_trkHists[5] = book1D("/iterN/", "chi2", "track #chi^{2}/dof", nBin_,   0.0, 20.0);
00067 
00068 }
00069 
00070 void AlignmentMonitorGeneric::event(const edm::Event &iEvent,
00071                                     const edm::EventSetup&,
00072                                     const ConstTrajTrackPairCollection& tracks)
00073 {
00074   static TrajectoryStateCombiner tsoscomb;
00075 
00076   for (unsigned int t = 0; t < tracks.size(); ++t)
00077   {
00078     const reco::Track* track = tracks[t].second;
00079 
00080     float charge = tracks[t].second->charge();
00081 
00082     const std::vector<TrajectoryMeasurement>& meass
00083       = tracks[t].first->measurements();
00084 
00085     for (unsigned int m = 0; m < meass.size(); ++m)
00086     {
00087       const TrajectoryMeasurement& meas = meass[m];
00088       const TransientTrackingRecHit& hit = *meas.recHit();
00089 
00090       if ( hit.isValid() )
00091       {
00092         const Alignable* ali = pNavigator()->alignableFromDetId( hit.geographicalId() );
00093 
00094         while (ali) {
00095           std::map<const Alignable*, Hist1Ds>::iterator h = m_resHists.find(ali);
00096           if ( h != m_resHists.end() )
00097             {
00098               TrajectoryStateOnSurface tsos = tsoscomb( meas.forwardPredictedState(), meas.backwardPredictedState() );
00099               
00100               align::LocalVector res = tsos.localPosition() - hit.localPosition();
00101               LocalError err1 = tsos.localError().positionError();
00102               LocalError err2 = hit.localPositionError();
00103               
00104               float errX = std::sqrt( err1.xx() + err2.xx() );
00105               float errY = std::sqrt( err1.yy() + err2.yy() );
00106               
00107               h->second[charge > 0 ? 0 : 1]->Fill(res.x() / errX);
00108               h->second[charge > 0 ? 2 : 3]->Fill(res.y() / errY);
00109             }
00110           ali = ali->mother();
00111         }
00112       }
00113     }
00114 
00115     m_trkHists[0]->Fill( track->pt() );
00116     m_trkHists[1]->Fill( track->eta() );
00117     m_trkHists[2]->Fill( track->phi() );
00118     m_trkHists[3]->Fill( track->d0() );
00119     m_trkHists[4]->Fill( track->dz() );
00120     m_trkHists[5]->Fill( track->normalizedChi2() );
00121   }
00122 }
00123 
00124 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorPluginFactory.h"
00125 DEFINE_EDM_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorGeneric, "AlignmentMonitorGeneric");
00126