CMS 3D CMS Logo

AlignmentMonitorGeneric Class Reference

#include <Alignment/CommonAlignmentMonitor/plugins/AlignmentMonitorGeneric.h>

Inheritance diagram for AlignmentMonitorGeneric:

AlignmentMonitorBase

List of all members.

Public Member Functions

virtual void afterAlignment (const edm::EventSetup &)
 Called after updating AlignableTracker and AlignableMuon (by "endOfLoop()"): may be reimplemented.
 AlignmentMonitorGeneric (const edm::ParameterSet &)
virtual void book ()
 Book or retrieve histograms; MUST be reimplemented.
virtual void event (const edm::EventSetup &, const ConstTrajTrackPairCollection &)

Private Types

typedef std::vector< TH1F * > Hist1Ds

Private Attributes

std::map< const Alignable *,
Hist1Ds
m_resHists
Hist1Ds m_trkHists

Static Private Attributes

static const unsigned int nBin_ = 50


Detailed Description

Definition at line 28 of file AlignmentMonitorGeneric.h.


Member Typedef Documentation

typedef std::vector<TH1F*> AlignmentMonitorGeneric::Hist1Ds [private]

Definition at line 31 of file AlignmentMonitorGeneric.h.


Constructor & Destructor Documentation

AlignmentMonitorGeneric::AlignmentMonitorGeneric ( const edm::ParameterSet cfg  ) 

Definition at line 10 of file AlignmentMonitorGeneric.cc.

00010                                                                           :
00011   AlignmentMonitorBase(cfg, "AlignmentMonitorGeneric")
00012 {
00013 }


Member Function Documentation

virtual void AlignmentMonitorGeneric::afterAlignment ( const edm::EventSetup iSetup  )  [inline, virtual]

Called after updating AlignableTracker and AlignableMuon (by "endOfLoop()"): may be reimplemented.

Reimplemented from AlignmentMonitorBase.

Definition at line 46 of file AlignmentMonitorGeneric.h.

00048                                 {}

void AlignmentMonitorGeneric::book (  )  [virtual]

Book or retrieve histograms; MUST be reimplemented.

Implements AlignmentMonitorBase.

Definition at line 15 of file AlignmentMonitorGeneric.cc.

References Alignable::alignableObjectId(), AlignmentParameterStore::alignables(), AlignmentMonitorBase::book1D(), i, Alignable::id(), m_resHists, m_trkHists, n, name, nBin_, AlignmentMonitorBase::pStore(), and AlignableObjectId::typeToName().

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 }

void AlignmentMonitorGeneric::event ( const edm::EventSetup ,
const ConstTrajTrackPairCollection tracks 
) [virtual]

Definition at line 70 of file AlignmentMonitorGeneric.cc.

References AlignableNavigator::alignableFromDetId(), TrajectoryMeasurement::backwardPredictedState(), reco::TrackBase::d0(), reco::TrackBase::dz(), err1, err2, reco::TrackBase::eta(), TrajectoryMeasurement::forwardPredictedState(), TrackingRecHit::geographicalId(), h, TrackingRecHit::isValid(), TrajectoryStateOnSurface::localError(), TrackingRecHit::localPosition(), TrajectoryStateOnSurface::localPosition(), TrackingRecHit::localPositionError(), m, m_resHists, m_trkHists, Alignable::mother(), reco::TrackBase::normalizedChi2(), reco::TrackBase::phi(), AlignmentMonitorBase::pNavigator(), LocalTrajectoryError::positionError(), reco::TrackBase::pt(), TrajectoryMeasurement::recHit(), res, funct::sqrt(), t, track, LocalError::xx(), and LocalError::yy().

00072 {
00073   static TrajectoryStateCombiner tsoscomb;
00074 
00075   for (unsigned int t = 0; t < tracks.size(); ++t)
00076   {
00077     const reco::Track* track = tracks[t].second;
00078 
00079     float charge = tracks[t].second->charge();
00080 
00081     const std::vector<TrajectoryMeasurement>& meass
00082       = tracks[t].first->measurements();
00083 
00084     for (unsigned int m = 0; m < meass.size(); ++m)
00085     {
00086       const TrajectoryMeasurement& meas = meass[m];
00087       const TransientTrackingRecHit& hit = *meas.recHit();
00088 
00089       if ( hit.isValid() )
00090       {
00091         const Alignable* ali = pNavigator()->alignableFromDetId( hit.geographicalId() );
00092 
00093         while (ali) {
00094           std::map<const Alignable*, Hist1Ds>::iterator h = m_resHists.find(ali);
00095           if ( h != m_resHists.end() )
00096             {
00097               TrajectoryStateOnSurface tsos = tsoscomb( meas.forwardPredictedState(), meas.backwardPredictedState() );
00098               
00099               align::LocalVector res = tsos.localPosition() - hit.localPosition();
00100               LocalError err1 = tsos.localError().positionError();
00101               LocalError err2 = hit.localPositionError();
00102               
00103               float errX = std::sqrt( err1.xx() + err2.xx() );
00104               float errY = std::sqrt( err1.yy() + err2.yy() );
00105               
00106               h->second[charge > 0 ? 0 : 1]->Fill(res.x() / errX);
00107               h->second[charge > 0 ? 2 : 3]->Fill(res.y() / errY);
00108             }
00109           ali = ali->mother();
00110         }
00111       }
00112     }
00113 
00114     m_trkHists[0]->Fill( track->pt() );
00115     m_trkHists[1]->Fill( track->eta() );
00116     m_trkHists[2]->Fill( track->phi() );
00117     m_trkHists[3]->Fill( track->d0() );
00118     m_trkHists[4]->Fill( track->dz() );
00119     m_trkHists[5]->Fill( track->normalizedChi2() );
00120   }
00121 }


Member Data Documentation

std::map<const Alignable*, Hist1Ds> AlignmentMonitorGeneric::m_resHists [private]

Definition at line 56 of file AlignmentMonitorGeneric.h.

Referenced by book(), and event().

Hist1Ds AlignmentMonitorGeneric::m_trkHists [private]

Definition at line 54 of file AlignmentMonitorGeneric.h.

Referenced by book(), and event().

const unsigned int AlignmentMonitorGeneric::nBin_ = 50 [static, private]

Definition at line 52 of file AlignmentMonitorGeneric.h.

Referenced by book().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:14:28 2009 for CMSSW by  doxygen 1.5.4