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