17 #include "TLorentzVector.h" 50 histname =
"InvariantMass_";
61 histname =
"TrackPtPositive_";
63 histname +
AlgoName, histname +
AlgoName, TrackPtPositiveBin, TrackPtPositiveMin, TrackPtPositiveMax);
70 histname =
"TrackPtNegative_";
72 histname +
AlgoName, histname +
AlgoName, TrackPtNegativeBin, TrackPtNegativeMin, TrackPtNegativeMax);
75 histname =
"TrackQuality_";
88 histname =
"SumCharge_";
96 histname =
"TrackCurvature_";
114 histname =
"MinJetDeltaR_";
127 histname =
"MinTrackDeltaR_";
136 histname =
"AlCaRecoTrackEfficiency_";
150 histname =
"Hits_ZvsR_";
157 histname =
"Hits_XvsY_";
161 histname =
"Hits_perDetId_";
189 edm::LogError(
"Alignment") <<
"invalid trackcollection encountered!";
196 edm::LogError(
"Alignment") <<
"invalid reference track-collection encountered!";
202 edm::LogError(
"Alignment") <<
"invalid geometry found in event setup!";
207 edm::LogError(
"Alignment") <<
"invalid magnetic field configuration encountered!";
214 if (!
jets.isValid()) {
215 edm::LogError(
"Alignment") <<
"no jet collection found in event!";
224 double sumOfCharges = 0;
225 for (reco::TrackCollection::const_iterator
track = (*trackCollection).begin();
track < (*trackCollection).end();
230 for (reco::CaloJetCollection::const_iterator itJet =
jets->begin(); itJet !=
jets->end(); ++itJet) {
242 double minTrackDeltaR = 10;
243 for (reco::TrackCollection::const_iterator track2 = (*trackCollection).begin(); track2 < (*trackCollection).end();
246 if (dR < minTrackDeltaR && dR > 1
e-6)
256 GlobalPoint gPoint((*track).vx(), (*track).vy(), (*track).vz());
258 double curv = -(*track).charge() * 0.002998 *
B / (*track).pt();
261 if ((*track).charge() > 0)
263 if ((*track).charge() < 0)
268 sumOfCharges += (*track).charge();
274 if ((*trackCollection).size() == 2) {
275 TLorentzVector track0(
276 (*trackCollection).at(0).px(),
277 (*trackCollection).at(0).py(),
278 (*trackCollection).at(0).pz(),
280 TLorentzVector track1(
281 (*trackCollection).at(1).px(),
282 (*trackCollection).at(1).py(),
283 (*trackCollection).at(1).pz(),
285 TLorentzVector mother = track0 + track1;
289 edm::LogInfo(
"Alignment") <<
"wrong number of tracks trackcollection encountered: " << (*trackCollection).size();
296 if ((*iHit)->isValid()) {
298 const DetId geoId(
hit->geographicalId());
305 double r =
sqrt(globP.x() * globP.x() + globP.y() * globP.y());
307 r *= globP.y() / fabs(globP.y());
319 std::vector<int> sortedRawIds;
320 for (std::vector<DetId>::const_iterator iDetId =
geometry.detUnitIds().begin(); iDetId !=
geometry.detUnitIds().end();
322 sortedRawIds.push_back((*iDetId).rawId());
324 std::sort(sortedRawIds.begin(), sortedRawIds.end());
327 for (std::vector<int>::iterator iRawId = sortedRawIds.begin(); iRawId != sortedRawIds.end(); ++iRawId) {
std::map< int, int > binByRawId_
T getParameter(std::string const &) const
static std::string qualityName(TrackQuality)
void analyze(const edm::Event &, const edm::EventSetup &) override
virtual void setCurrentFolder(std::string const &fullpath)
void fillRawIdMap(const TrackerGeometry &geometry)
TrackQuality
track quality
#define DEFINE_FWK_MODULE(type)
MonitorElement * TrackQuality_
MonitorElement * minJetDeltaR_
Log< level::Error, false > LogError
MonitorElement * Hits_XvsY_
MonitorElement * Hits_ZvsR_
MonitorElement * AlCaRecoTrackEfficiency_
MonitorElement * sumCharge_
edm::EDGetTokenT< reco::TrackCollection > trackProducer_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > mfToken_
MonitorElement * invariantMass_
MonitorElement * TrackPtNegative_
void fillHitmaps(const reco::Track &track, const TrackerGeometry &geometry)
edm::EDGetTokenT< reco::CaloJetCollection > jetCollection_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Log< level::Info, false > LogInfo
MonitorElement * Hits_perDetId_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * TrackCurvature_
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
MonitorElement * minTrackDeltaR_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
edm::EDGetTokenT< reco::TrackCollection > referenceTrackProducer_
MonitorElement * TrackPtPositive_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
TkAlCaRecoMonitor(const edm::ParameterSet &)
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)