Go to the documentation of this file.00001 #include "Validation/RPCRecHits/interface/RPCPointVsRecHit.h"
00002 #include "FWCore/Framework/interface/MakerMacros.h"
00003
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006
00007 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
00008 #include "Geometry/RPCGeometry/interface/RPCRoll.h"
00009 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
00010 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00011
00012 using namespace std;
00013
00014 typedef MonitorElement* MEP;
00015
00016 RPCPointVsRecHit::RPCPointVsRecHit(const edm::ParameterSet& pset)
00017 {
00018 rootFileName_ = pset.getUntrackedParameter<string>("rootFileName", "");
00019 refHitLabel_ = pset.getParameter<edm::InputTag>("refHit");
00020 recHitLabel_ = pset.getParameter<edm::InputTag>("recHit");
00021
00022 isStandAloneMode_ = pset.getUntrackedParameter<bool>("standAloneMode", false);
00023
00024 dbe_ = edm::Service<DQMStore>().operator->();
00025 if ( !dbe_ )
00026 {
00027 edm::LogError("RPCPointVsRecHit") << "No DQMStore instance\n";
00028 return;
00029 }
00030
00031
00032 const std::string subDir = pset.getParameter<std::string>("subDir");
00033 h_.bookHistograms(dbe_, subDir);
00034 }
00035
00036 RPCPointVsRecHit::~RPCPointVsRecHit()
00037 {
00038 if ( dbe_ )
00039 {
00040 if ( !rootFileName_.empty() ) dbe_->save(rootFileName_);
00041 }
00042 }
00043
00044 void RPCPointVsRecHit::beginJob()
00045 {
00046 }
00047
00048 void RPCPointVsRecHit::endJob()
00049 {
00050 }
00051
00052 void RPCPointVsRecHit::analyze(const edm::Event& event, const edm::EventSetup& eventSetup)
00053 {
00054 if ( !dbe_ )
00055 {
00056 edm::LogError("RPCPointVsRecHit") << "No DQMStore instance\n";
00057 return;
00058 }
00059
00060
00061 edm::ESHandle<RPCGeometry> rpcGeom;
00062 eventSetup.get<MuonGeometryRecord>().get(rpcGeom);
00063
00064
00065 edm::Handle<RPCRecHitCollection> refHitHandle;
00066 if ( !event.getByLabel(refHitLabel_, refHitHandle) )
00067 {
00068 edm::LogInfo("RPCPointVsRecHit") << "Cannot find reference hit collection\n";
00069 return;
00070 }
00071
00072
00073 edm::Handle<RPCRecHitCollection> recHitHandle;
00074 if ( !event.getByLabel(recHitLabel_, recHitHandle) )
00075 {
00076 edm::LogInfo("RPCPointVsRecHit") << "Cannot find recHit collection\n";
00077 return;
00078 }
00079
00080 typedef RPCRecHitCollection::const_iterator RecHitIter;
00081
00082
00083 for ( RecHitIter refHitIter = refHitHandle->begin();
00084 refHitIter != refHitHandle->end(); ++refHitIter )
00085 {
00086 const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
00087 const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
00088 if ( !roll ) continue;
00089
00090 const int region = roll->id().region();
00091 const int ring = roll->id().ring();
00092
00093 const int station = roll->id().station();
00094
00095
00096
00097 if ( region == 0 )
00098 {
00099 h_.nRefHit_W->Fill(ring);
00100 h_.nRefHit_WvsR->Fill(ring, station);
00101 }
00102 else
00103 {
00104 h_.nRefHit_D->Fill(region*station);
00105 h_.nRefHit_DvsR->Fill(region*station, ring);
00106 }
00107
00108 const GlobalPoint pos = roll->toGlobal(refHitIter->localPosition());
00109
00110 }
00111
00112
00113 for ( RecHitIter recHitIter = recHitHandle->begin();
00114 recHitIter != recHitHandle->end(); ++recHitIter )
00115 {
00116 const RPCDetId detId = static_cast<const RPCDetId>(recHitIter->rpcId());
00117 const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
00118 if ( !roll ) continue;
00119
00120 const int region = roll->id().region();
00121 const int ring = roll->id().ring();
00122
00123 const int station = (roll->id().station());
00124
00125
00126
00127 h_.clusterSize->Fill(recHitIter->clusterSize());
00128
00129 if ( region == 0 )
00130 {
00131 h_.nRecHit_W->Fill(ring);
00132 h_.nRecHit_WvsR->Fill(ring, station);
00133 }
00134 else
00135 {
00136 h_.nRecHit_D->Fill(region*station);
00137 h_.nRecHit_DvsR->Fill(ring, region*station);
00138 }
00139
00140
00141 const GlobalPoint pos = roll->toGlobal(recHitIter->localPosition());
00142
00143 }
00144
00145
00146 typedef std::map<RecHitIter, RecHitIter> RecToRecHitMap;
00147 RecToRecHitMap refToRecHitMap;
00148
00149 for ( RecHitIter refHitIter = refHitHandle->begin();
00150 refHitIter != refHitHandle->end(); ++refHitIter )
00151 {
00152 const RPCDetId refDetId = static_cast<const RPCDetId>(refHitIter->rpcId());
00153 const RPCRoll* refRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(refDetId));
00154 if ( !refRoll ) continue;
00155
00156 const double refX = refHitIter->localPosition().x();
00157
00158 for ( RecHitIter recHitIter = recHitHandle->begin();
00159 recHitIter != recHitHandle->end(); ++recHitIter )
00160 {
00161 const RPCDetId recDetId = static_cast<const RPCDetId>(recHitIter->rpcId());
00162 const RPCRoll* recRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(recDetId));
00163 if ( !recRoll ) continue;
00164
00165 if ( refDetId != recDetId ) continue;
00166
00167 const double recX = recHitIter->localPosition().x();
00168 const double newDx = fabs(recX - refX);
00169
00170
00171 RecToRecHitMap::const_iterator prevRefToReco = refToRecHitMap.find(refHitIter);
00172 if ( prevRefToReco == refToRecHitMap.end() )
00173 {
00174 refToRecHitMap.insert(std::make_pair(refHitIter, recHitIter));
00175 }
00176 else
00177 {
00178 const double oldDx = fabs(prevRefToReco->second->localPosition().x() - refX);
00179
00180 if ( newDx < oldDx )
00181 {
00182 refToRecHitMap[refHitIter] = recHitIter;
00183 }
00184 }
00185 }
00186 }
00187
00188
00189
00190 for ( RecToRecHitMap::const_iterator match = refToRecHitMap.begin();
00191 match != refToRecHitMap.end(); ++match )
00192 {
00193 RecHitIter refHitIter = match->first;
00194 RecHitIter recHitIter = match->second;
00195
00196 const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
00197 const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
00198
00199 const int region = roll->id().region();
00200 const int ring = roll->id().ring();
00201
00202 const int station = roll->id().station();
00203
00204
00205
00206 const double refX = refHitIter->localPosition().x();
00207 const double recX = recHitIter->localPosition().x();
00208 const double errX = recHitIter->localPositionError().xx();
00209 const double dX = recX - refX;
00210 const double pull = errX == 0 ? -999 : dX/errX;
00211
00212 const GlobalPoint refPos = roll->toGlobal(refHitIter->localPosition());
00213 const GlobalPoint recPos = roll->toGlobal(recHitIter->localPosition());
00214
00215 if ( region == 0 )
00216 {
00217 h_.res_W->Fill(dX);
00218 h_.pull_W->Fill(pull);
00219 h_.nMatchedRefHit_W->Fill(ring);
00220 h_.nMatchedRefHit_WvsR->Fill(ring, station);
00221
00222 h_.res2_W->Fill(ring, dX);
00223 h_.pull2_W->Fill(ring, pull);
00224 }
00225 else
00226 {
00227 h_.res_D->Fill(dX);
00228 h_.pull_D->Fill(pull);
00229 h_.nMatchedRefHit_D->Fill(region*station);
00230 h_.nMatchedRefHit_DvsR->Fill(region*station, ring);
00231
00232 h_.res2_D->Fill(region*station, dX);
00233 h_.pull2_D->Fill(region*station, pull);
00234 }
00235 }
00236
00237
00238 for ( RecHitIter refHitIter = refHitHandle->begin();
00239 refHitIter != refHitHandle->end(); ++refHitIter )
00240 {
00241 const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
00242 const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
00243
00244 const int region = roll->id().region();
00245 const int ring = roll->id().ring();
00246
00247 const int station = roll->id().station();
00248
00249
00250
00251 bool matched = false;
00252 for ( RecToRecHitMap::const_iterator match = refToRecHitMap.begin();
00253 match != refToRecHitMap.end(); ++match )
00254 {
00255 if ( refHitIter == match->first )
00256 {
00257 matched = true;
00258 break;
00259 }
00260 }
00261
00262 if ( !matched )
00263 {
00264 if ( region == 0 )
00265 {
00266 h_.nUnMatchedRefHit_W->Fill(ring);
00267 h_.nUnMatchedRefHit_WvsR->Fill(ring, station);
00268 }
00269 else
00270 {
00271 h_.nUnMatchedRefHit_D->Fill(region*station);
00272 h_.nUnMatchedRefHit_DvsR->Fill(region*station, ring);
00273 }
00274 }
00275 }
00276
00277
00278 for ( RecHitIter recHitIter = recHitHandle->begin();
00279 recHitIter != recHitHandle->end(); ++recHitIter )
00280 {
00281 const RPCDetId detId = static_cast<const RPCDetId>(recHitIter->rpcId());
00282 const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
00283
00284 const int region = roll->id().region();
00285 const int ring = roll->id().ring();
00286
00287 const int station = roll->id().station();
00288
00289
00290
00291 bool matched = false;
00292 for ( RecToRecHitMap::const_iterator match = refToRecHitMap.begin();
00293 match != refToRecHitMap.end(); ++match )
00294 {
00295 if ( recHitIter == match->second )
00296 {
00297 matched = true;
00298 break;
00299 }
00300 }
00301
00302 if ( !matched )
00303 {
00304 if ( region == 0 )
00305 {
00306 h_.nUnMatchedRecHit_W->Fill(ring);
00307 h_.nUnMatchedRecHit_WvsR->Fill(ring, station);
00308 }
00309 else
00310 {
00311 h_.nUnMatchedRecHit_D->Fill(region*station);
00312 h_.nUnMatchedRecHit_DvsR->Fill(region*station, ring);
00313 }
00314
00315 const GlobalPoint pos = roll->toGlobal(recHitIter->localPosition());
00316
00317 }
00318 }
00319 }
00320
00321 DEFINE_FWK_MODULE(RPCPointVsRecHit);
00322