CMS 3D CMS Logo

Public Member Functions | Private Types | Private Attributes | Static Private Attributes

AlignmentMonitorGeneric Class Reference

#include <AlignmentMonitorGeneric.h>

Inheritance diagram for AlignmentMonitorGeneric:
AlignmentMonitorBase

List of all members.

Public Member Functions

virtual void afterAlignment (const edm::EventSetup &)
 AlignmentMonitorGeneric (const edm::ParameterSet &)
virtual void book ()
 Book or retrieve histograms; MUST be reimplemented.
virtual void event (const edm::Event &, const edm::EventSetup &, const ConstTrajTrackPairCollection &)
 Called for each event (by "run()"): may be reimplemented.

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.

                                                                          :
  AlignmentMonitorBase(cfg, "AlignmentMonitorGeneric")
{
}

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.

                                {}
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, errorMatrix2Lands_multiChannel::id, Alignable::id(), M_PI, m_resHists, m_trkHists, n, mergeVDriftHistosByStation::name, nBin_, AlignmentMonitorBase::pStore(), DetId::subdetId(), and AlignableObjectId::typeToName().

{
  static AlignableObjectId idMap;

  std::vector<std::string> residNames; // names of residual histograms

  residNames.push_back("x hit residuals pos track");
  residNames.push_back("x hit residuals neg track");
  residNames.push_back("y hit residuals pos track");
  residNames.push_back("y hit residuals neg track");

  const std::vector<Alignable*>& alignables = pStore()->alignables();

  unsigned int nAlignable = alignables.size();
  unsigned int nResidName = residNames.size();

  for (unsigned int i = 0; i < nAlignable; ++i)
  {
    const Alignable* ali = alignables[i];

    Hist1Ds& hists = m_resHists[ali];

    hists.resize(nResidName, 0);

    align::ID id = ali->id();
    align::StructureType type = ali->alignableObjectId();

    for (unsigned int n = 0; n < nResidName; ++n)
    {
      const std::string& name = residNames[n];

      TString histName(name.c_str());
      histName += Form("_%s_%d", idMap.typeToName(type).c_str(), id);
      histName.ReplaceAll(" ", "");

      TString histTitle(name.c_str());
      histTitle += Form(" for %s with ID %d (subdet %d)",
                        idMap.typeToName(type).c_str(),
                        id, DetId(id).subdetId());

      hists[n] = book1D(std::string("/iterN/") + std::string(name) + std::string("/"), std::string(histName), std::string(histTitle), nBin_, -5., 5.);
    }
  }

  m_trkHists.resize(6, 0);
  
  m_trkHists[0] = book1D("/iterN/", "pt"  , "track p_{t} (GeV)" , nBin_,   0.0,100.0);
  m_trkHists[1] = book1D("/iterN/", "eta" , "track #eta"        , nBin_, - 3.0,  3.0);
  m_trkHists[2] = book1D("/iterN/", "phi" , "track #phi"        , nBin_, -M_PI, M_PI);
  m_trkHists[3] = book1D("/iterN/", "d0"  , "track d0 (cm)"     , nBin_, -0.02, 0.02);
  m_trkHists[4] = book1D("/iterN/", "dz"  , "track dz (cm)"     , nBin_, -20.0, 20.0);
  m_trkHists[5] = book1D("/iterN/", "chi2", "track #chi^{2}/dof", nBin_,   0.0, 20.0);

}
void AlignmentMonitorGeneric::event ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const ConstTrajTrackPairCollection iTrajTracks 
) [virtual]

Called for each event (by "run()"): may be reimplemented.

Reimplemented from AlignmentMonitorBase.

Definition at line 70 of file AlignmentMonitorGeneric.cc.

References AlignableNavigator::alignableFromDetId(), TrajectoryMeasurement::backwardPredictedState(), DeDxDiscriminatorTools::charge(), reco::TrackBase::d0(), reco::TrackBase::dz(), 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(), mathSSE::sqrt(), lumiQTWidget::t, LocalError::xx(), and LocalError::yy().

{
  static TrajectoryStateCombiner tsoscomb;

  for (unsigned int t = 0; t < tracks.size(); ++t)
  {
    const reco::Track* track = tracks[t].second;

    float charge = tracks[t].second->charge();

    const std::vector<TrajectoryMeasurement>& meass
      = tracks[t].first->measurements();

    for (unsigned int m = 0; m < meass.size(); ++m)
    {
      const TrajectoryMeasurement& meas = meass[m];
      const TransientTrackingRecHit& hit = *meas.recHit();

      if ( hit.isValid() )
      {
        const Alignable* ali = pNavigator()->alignableFromDetId( hit.geographicalId() );

        while (ali) {
          std::map<const Alignable*, Hist1Ds>::iterator h = m_resHists.find(ali);
          if ( h != m_resHists.end() )
            {
              TrajectoryStateOnSurface tsos = tsoscomb( meas.forwardPredictedState(), meas.backwardPredictedState() );
              
              align::LocalVector res = tsos.localPosition() - hit.localPosition();
              LocalError err1 = tsos.localError().positionError();
              LocalError err2 = hit.localPositionError();
              
              float errX = std::sqrt( err1.xx() + err2.xx() );
              float errY = std::sqrt( err1.yy() + err2.yy() );
              
              h->second[charge > 0 ? 0 : 1]->Fill(res.x() / errX);
              h->second[charge > 0 ? 2 : 3]->Fill(res.y() / errY);
            }
          ali = ali->mother();
        }
      }
    }

    m_trkHists[0]->Fill( track->pt() );
    m_trkHists[1]->Fill( track->eta() );
    m_trkHists[2]->Fill( track->phi() );
    m_trkHists[3]->Fill( track->d0() );
    m_trkHists[4]->Fill( track->dz() );
    m_trkHists[5]->Fill( track->normalizedChi2() );
  }
}

Member Data Documentation

Definition at line 56 of file AlignmentMonitorGeneric.h.

Referenced by book(), and event().

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().