00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorPluginFactory.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/Utilities/interface/InputTag.h"
00017 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00018 #include "DataFormats/GeometrySurface/interface/Surface.h"
00019 #include "TH1.h"
00020 #include "TObject.h"
00021 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorBase.h"
00022
00023 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
00024 #include "RecoMuon/TrackingTools/interface/MuonUpdatorAtVertex.h"
00025 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00026
00027 #include <fstream>
00028
00029
00030
00031
00032
00033
00034
00035 class AlignmentMonitorTracksFromTrajectories: public AlignmentMonitorBase {
00036 public:
00037 AlignmentMonitorTracksFromTrajectories(const edm::ParameterSet& cfg);
00038 ~AlignmentMonitorTracksFromTrajectories() {};
00039
00040 void book();
00041 void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& iTrajTracks);
00042 void afterAlignment(const edm::EventSetup &iSetup);
00043
00044 private:
00045 MuonServiceProxy *theMuonServiceProxy;
00046 MuonUpdatorAtVertex *theMuonUpdatorAtVertex;
00047 bool m_vertexConstraint;
00048 edm::InputTag m_beamSpot;
00049
00050 TH1F *m_diMuon_Z;
00051 TH1F *m_diMuon_Zforward;
00052 TH1F *m_diMuon_Zbarrel;
00053 TH1F *m_diMuon_Zbackward;
00054 TH1F *m_diMuon_Ups;
00055 TH1F *m_diMuon_Jpsi;
00056 TH1F *m_diMuon_log;
00057 TH1F *m_chi2_100;
00058 TH1F *m_chi2_10000;
00059 TH1F *m_chi2_1000000;
00060 TH1F *m_chi2_log;
00061 TH1F *m_chi2DOF_5;
00062 TH1F *m_chi2DOF_100;
00063 TH1F *m_chi2DOF_1000;
00064 TH1F *m_chi2DOF_log;
00065 TH1F *m_chi2_improvement;
00066 TH1F *m_chi2DOF_improvement;
00067 TH1F *m_pt[36];
00068 };
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099 AlignmentMonitorTracksFromTrajectories::AlignmentMonitorTracksFromTrajectories(const edm::ParameterSet& cfg)
00100 : AlignmentMonitorBase(cfg, "AlignmentMonitorTracksFromTrajectories")
00101 , m_vertexConstraint(cfg.getParameter<bool>("vertexConstraint"))
00102 , m_beamSpot(cfg.getParameter<edm::InputTag>("beamSpot"))
00103 {
00104 theMuonServiceProxy = new MuonServiceProxy(cfg.getParameter<edm::ParameterSet>("ServiceParameters"));
00105 theMuonUpdatorAtVertex = new MuonUpdatorAtVertex(cfg.getParameter<edm::ParameterSet>("MuonUpdatorAtVertexParameters"), theMuonServiceProxy);
00106 }
00107
00108
00109
00110
00111
00113
00115
00116 void AlignmentMonitorTracksFromTrajectories::book() {
00117 m_diMuon_Z = book1D("/iterN/", "diMuon_Z", "Di-muon mass (GeV)", 200, 90. - 50., 90. + 50.);
00118 m_diMuon_Zforward = book1D("/iterN/", "diMuon_Zforward", "Di-muon mass (GeV) eta > 1.4", 200, 90. - 50., 90. + 50.);
00119 m_diMuon_Zbarrel = book1D("/iterN/", "diMuon_Zbarrel", "Di-muon mass (GeV) -1.4 < eta < 1.4", 200, 90. - 50., 90. + 50.);
00120 m_diMuon_Zbackward = book1D("/iterN/", "diMuon_Zbackward", "Di-muon mass (GeV) eta < -1.4", 200, 90. - 50., 90. + 50.);
00121 m_diMuon_Ups = book1D("/iterN/", "diMuon_Ups", "Di-muon mass (GeV)", 200, 9. - 3., 9. + 3.);
00122 m_diMuon_Jpsi = book1D("/iterN/", "diMuon_Jpsi", "Di-muon mass (GeV)", 200, 3. - 1., 3. + 1.);
00123 m_diMuon_log = book1D("/iterN/", "diMuon_log", "Di-muon mass (log GeV)", 300, 0, 3);
00124 m_chi2_100 = book1D("/iterN/", "m_chi2_100", "Track chi^2", 100, 0, 100);
00125 m_chi2_10000 = book1D("/iterN/", "m_chi2_10000", "Track chi^2", 100, 0, 10000);
00126 m_chi2_1000000 = book1D("/iterN/", "m_chi2_1000000", "Track chi^2", 100, 1, 1000000);
00127 m_chi2_log = book1D("/iterN/", "m_chi2_log", "Log track chi^2", 100, -3, 7);
00128 m_chi2DOF_5 = book1D("/iterN/", "m_chi2DOF_5", "Track chi^2/nDOF", 100, 0, 5);
00129 m_chi2DOF_100 = book1D("/iterN/", "m_chi2DOF_100", "Track chi^2/nDOF", 100, 0, 100);
00130 m_chi2DOF_1000 = book1D("/iterN/", "m_chi2DOF_1000", "Track chi^2/nDOF", 100, 0, 1000);
00131 m_chi2DOF_log = book1D("/iterN/", "m_chi2DOF_log", "Log track chi^2/nDOF", 100, -3, 7);
00132 m_chi2_improvement = book1D("/iterN/", "m_chi2_improvement", "Track-by-track chi^2/original chi^2", 100, 0., 10.);
00133 m_chi2DOF_improvement = book1D("/iterN/", "m_chi2DOF_improvement", "Track-by-track (chi^2/DOF)/(original chi^2/original DOF)", 100, 0., 10.);
00134 for (int i = 0; i < 36; i++) {
00135 char name[100], title[100];
00136 snprintf(name, sizeof(name), "m_pt_phi%d", i);
00137 snprintf(title, sizeof(title), "Track pt (GeV) in phi bin %d/36", i);
00138 m_pt[i] = book1D("/iterN/", name, title, 100, 0, 100);
00139 }
00140 }
00141
00143
00145
00146 void AlignmentMonitorTracksFromTrajectories::event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& tracks) {
00147 theMuonServiceProxy->update(iSetup);
00148
00149 edm::Handle<reco::BeamSpot> beamSpot;
00150 iEvent.getByLabel(m_beamSpot, beamSpot);
00151
00152 GlobalVector p1, p2;
00153 double e1 = 0.;
00154 double e2 = 0.;
00155
00156 for (ConstTrajTrackPairCollection::const_iterator it = tracks.begin(); it != tracks.end(); ++it) {
00157 const Trajectory *traj = it->first;
00158 const reco::Track *track = it->second;
00159
00160 int nDOF = 0;
00161 TrajectoryStateOnSurface closestTSOS;
00162 double closest = 10000.;
00163
00164 std::vector<TrajectoryMeasurement> measurements = traj->measurements();
00165 for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end(); ++im) {
00166 const TrajectoryMeasurement meas = *im;
00167 const TransientTrackingRecHit* hit = &(*meas.recHit());
00168
00169
00170 nDOF += hit->dimension();
00171
00172 TrajectoryStateOnSurface TSOS = meas.backwardPredictedState();
00173 GlobalPoint where = TSOS.surface().toGlobal(LocalPoint(0,0,0));
00174
00175 if (where.mag() < closest) {
00176 closest = where.mag();
00177 closestTSOS = TSOS;
00178 }
00179
00180 }
00181 nDOF -= 5;
00182
00183 if (closest != 10000.) {
00184 std::pair<bool, FreeTrajectoryState> state;
00185
00186 if (m_vertexConstraint) {
00187 state = theMuonUpdatorAtVertex->propagateWithUpdate(closestTSOS, *beamSpot);
00188
00189 }
00190 else {
00191 state = theMuonUpdatorAtVertex->propagate(closestTSOS, *beamSpot);
00192 }
00193
00194 if (state.first) {
00195 double chi2 = traj->chiSquared();
00196 double chi2DOF = chi2 / double(nDOF);
00197
00198 m_chi2_100->Fill(chi2);
00199 m_chi2_10000->Fill(chi2);
00200 m_chi2_1000000->Fill(chi2);
00201 m_chi2_log->Fill(log10(chi2));
00202
00203 m_chi2DOF_5->Fill(chi2DOF);
00204 m_chi2DOF_100->Fill(chi2DOF);
00205 m_chi2DOF_1000->Fill(chi2DOF);
00206 m_chi2DOF_log->Fill(log10(chi2DOF));
00207 m_chi2_improvement->Fill(chi2 / track->chi2());
00208 m_chi2DOF_improvement->Fill(chi2DOF / track->normalizedChi2());
00209
00210 GlobalVector momentum = state.second.momentum();
00211 double energy = momentum.mag();
00212
00213 if (energy > e1) {
00214 e2 = e1;
00215 e1 = energy;
00216 p2 = p1;
00217 p1 = momentum;
00218 }
00219 else if (energy > e2) {
00220 e2 = energy;
00221 p2 = momentum;
00222 }
00223
00224 double pt = momentum.perp();
00225 double phi = momentum.phi();
00226 while (phi < -M_PI) phi += 2.* M_PI;
00227 while (phi > M_PI) phi -= 2.* M_PI;
00228
00229 double phibin = (phi + M_PI) / (2. * M_PI) * 36.;
00230
00231 for (int i = 0; i < 36; i++) {
00232 if (phibin < i+1) {
00233 m_pt[i]->Fill(pt);
00234 break;
00235 }
00236 }
00237
00238 }
00239 }
00240 }
00241
00242 if (e1 > 0. && e2 > 0.) {
00243 double energy_tot = e1 + e2;
00244 GlobalVector momentum_tot = p1 + p2;
00245 double mass = sqrt(energy_tot*energy_tot - momentum_tot.mag2());
00246 double eta = momentum_tot.eta();
00247
00248 m_diMuon_Z->Fill(mass);
00249 if (eta > 1.4) m_diMuon_Zforward->Fill(mass);
00250 else if (eta > -1.4) m_diMuon_Zbarrel->Fill(mass);
00251 else m_diMuon_Zbackward->Fill(mass);
00252
00253 m_diMuon_Ups->Fill(mass);
00254 m_diMuon_Jpsi->Fill(mass);
00255 m_diMuon_log->Fill(log10(mass));
00256 }
00257 }
00258
00260
00262
00263 void AlignmentMonitorTracksFromTrajectories::afterAlignment(const edm::EventSetup &iSetup) {
00264 }
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278 DEFINE_EDM_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorTracksFromTrajectories, "AlignmentMonitorTracksFromTrajectories");