Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:10

0001 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
0002 #include "Alignment/CommonAlignment/interface/AlignableObjectId.h"
0003 
0004 #include "Alignment/CommonAlignmentMonitor/plugins/AlignmentMonitorGeneric.h"
0005 #include <DataFormats/GeometrySurface/interface/LocalError.h>
0006 #include "Geometry/CommonDetUnit/interface/TrackerGeomDet.h"
0007 #include "TObject.h"
0008 
0009 #include <TString.h>
0010 
0011 AlignmentMonitorGeneric::AlignmentMonitorGeneric(const edm::ParameterSet& cfg, edm::ConsumesCollector iC)
0012     : AlignmentMonitorBase(cfg, iC, "AlignmentMonitorGeneric") {}
0013 
0014 void AlignmentMonitorGeneric::book() {
0015   std::vector<std::string> residNames;  // names of residual histograms
0016 
0017   residNames.push_back("x hit residuals pos track");
0018   residNames.push_back("x hit residuals neg track");
0019   residNames.push_back("y hit residuals pos track");
0020   residNames.push_back("y hit residuals neg track");
0021 
0022   auto alignableObjectId = AlignableObjectId::commonObjectIdProvider(pTracker(), pMuon());
0023 
0024   const auto& alignables = pStore()->alignables();
0025 
0026   unsigned int nAlignable = alignables.size();
0027   unsigned int nResidName = residNames.size();
0028 
0029   for (unsigned int i = 0; i < nAlignable; ++i) {
0030     const Alignable* ali = alignables[i];
0031 
0032     Hist1Ds& hists = m_resHists[ali];
0033 
0034     hists.resize(nResidName, nullptr);
0035 
0036     align::ID id = ali->id();
0037     align::StructureType type = ali->alignableObjectId();
0038 
0039     for (unsigned int n = 0; n < nResidName; ++n) {
0040       const std::string& name = residNames[n];
0041 
0042       TString histName(name.c_str());
0043       histName += Form("_%s_%d", alignableObjectId.idToString(type), id);
0044       histName.ReplaceAll(" ", "");
0045 
0046       TString histTitle(name.c_str());
0047       histTitle += Form(" for %s with ID %d (subdet %d)", alignableObjectId.idToString(type), id, DetId(id).subdetId());
0048 
0049       hists[n] = book1D(std::string("/iterN/") + std::string(name) + std::string("/"),
0050                         std::string(histName.Data()),
0051                         std::string(histTitle.Data()),
0052                         nBin_,
0053                         -5.,
0054                         5.);
0055     }
0056   }
0057 
0058   m_trkHists.resize(6, nullptr);
0059 
0060   m_trkHists[0] = book1D("/iterN/", "pt", "track p_{t} (GeV)", nBin_, 0.0, 100.0);
0061   m_trkHists[1] = book1D("/iterN/", "eta", "track #eta", nBin_, -3.0, 3.0);
0062   m_trkHists[2] = book1D("/iterN/", "phi", "track #phi", nBin_, -M_PI, M_PI);
0063   m_trkHists[3] = book1D("/iterN/", "d0", "track d0 (cm)", nBin_, -0.02, 0.02);
0064   m_trkHists[4] = book1D("/iterN/", "dz", "track dz (cm)", nBin_, -20.0, 20.0);
0065   m_trkHists[5] = book1D("/iterN/", "chi2", "track #chi^{2}/dof", nBin_, 0.0, 20.0);
0066 }
0067 
0068 void AlignmentMonitorGeneric::event(const edm::Event& iEvent,
0069                                     const edm::EventSetup&,
0070                                     const ConstTrajTrackPairCollection& tracks) {
0071   TrajectoryStateCombiner tsoscomb;
0072 
0073   for (unsigned int t = 0; t < tracks.size(); ++t) {
0074     const reco::Track* track = tracks[t].second;
0075 
0076     float charge = tracks[t].second->charge();
0077 
0078     const std::vector<TrajectoryMeasurement>& meass = tracks[t].first->measurements();
0079 
0080     for (unsigned int m = 0; m < meass.size(); ++m) {
0081       const TrajectoryMeasurement& meas = meass[m];
0082       const TransientTrackingRecHit& hit = *meas.recHit();
0083 
0084       if (hit.isValid()) {
0085         const Alignable* ali = pNavigator()->alignableFromDetId(hit.geographicalId());
0086 
0087         while (ali) {
0088           std::map<const Alignable*, Hist1Ds>::iterator h = m_resHists.find(ali);
0089           if (h != m_resHists.end()) {
0090             TrajectoryStateOnSurface tsos = tsoscomb(meas.forwardPredictedState(), meas.backwardPredictedState());
0091 
0092             align::LocalVector res = tsos.localPosition() - hit.localPosition();
0093             LocalError err1 = tsos.localError().positionError();
0094             LocalError err2 = hit.localPositionError();  // CPE+APE
0095 
0096             // subtract APEs from err2 (if existing) from covariance matrix
0097             auto det = static_cast<const TrackerGeomDet*>(hit.det());
0098             const auto localAPE = det->localAlignmentError();
0099             if (localAPE.valid()) {
0100               err2 = LocalError(err2.xx() - localAPE.xx(), err2.xy() - localAPE.xy(), err2.yy() - localAPE.yy());
0101             }
0102 
0103             float errX = std::sqrt(err1.xx() + err2.xx());
0104             float errY = std::sqrt(err1.yy() + err2.yy());
0105 
0106             h->second[charge > 0 ? 0 : 1]->Fill(res.x() / errX);
0107             h->second[charge > 0 ? 2 : 3]->Fill(res.y() / errY);
0108           }
0109           ali = ali->mother();
0110         }
0111       }
0112     }
0113 
0114     m_trkHists[0]->Fill(track->pt());
0115     m_trkHists[1]->Fill(track->eta());
0116     m_trkHists[2]->Fill(track->phi());
0117     m_trkHists[3]->Fill(track->d0());
0118     m_trkHists[4]->Fill(track->dz());
0119     m_trkHists[5]->Fill(track->normalizedChi2());
0120   }
0121 }
0122 
0123 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorPluginFactory.h"
0124 DEFINE_EDM_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorGeneric, "AlignmentMonitorGeneric");