00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "DTRecHitQuality.h"
00010 #include "Validation/DTRecHits/interface/DTHitQualityUtils.h"
00011
00012 #include "Geometry/DTGeometry/interface/DTLayer.h"
00013 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00014 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00015
00016 #include "FWCore/Framework/interface/MakerMacros.h"
00017 #include "FWCore/Framework/interface/Frameworkfwd.h"
00018 #include "FWCore/Framework/interface/Event.h"
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020
00021 #include "Histograms.h"
00022
00023 #include "TFile.h"
00024
00025 #include <iostream>
00026 #include <map>
00027
00028 using namespace std;
00029 using namespace edm;
00030
00031
00032
00033
00034
00035
00036 DTRecHitQuality::DTRecHitQuality(const ParameterSet& pset){
00037
00038 debug = pset.getUntrackedParameter<bool>("debug");
00039 rootFileName = pset.getUntrackedParameter<string>("rootFileName");
00040
00041 simHitLabel = pset.getUntrackedParameter<string>("simHitLabel", "SimG4Object");
00042
00043 recHitLabel = pset.getUntrackedParameter<string>("recHitLabel", "DTRecHit1DProducer");
00044
00045 segment2DLabel = pset.getUntrackedParameter<string>("segment2DLabel");
00046
00047 segment4DLabel = pset.getUntrackedParameter<string>("segment4DLabel");
00048
00049
00050 doStep1 = pset.getUntrackedParameter<bool>("doStep1", false);
00051 doStep2 = pset.getUntrackedParameter<bool>("doStep2", false);
00052 doStep3 = pset.getUntrackedParameter<bool>("doStep3", false);
00053
00054
00055 theFile = new TFile(rootFileName.c_str(), "RECREATE");
00056 theFile->cd();
00057
00058
00059 hRes_S1RPhi= new HRes1DHit("S1RPhi");
00060 hRes_S2RPhi= new HRes1DHit("S2RPhi");
00061 hRes_S3RPhi= new HRes1DHit("S3RPhi");
00062
00063 hRes_S1RZ= new HRes1DHit("S1RZ");
00064 hRes_S2RZ= new HRes1DHit("S2RZ");
00065 hRes_S3RZ= new HRes1DHit("S3RZ");
00066
00067 hRes_S1RZ_W0= new HRes1DHit("S1RZ_W0");
00068 hRes_S2RZ_W0= new HRes1DHit("S2RZ_W0");
00069 hRes_S3RZ_W0= new HRes1DHit("S3RZ_W0");
00070
00071 hRes_S1RZ_W1= new HRes1DHit("S1RZ_W1");
00072 hRes_S2RZ_W1= new HRes1DHit("S2RZ_W1");
00073 hRes_S3RZ_W1= new HRes1DHit("S3RZ_W1");
00074
00075 hRes_S1RZ_W2= new HRes1DHit("S1RZ_W2");
00076 hRes_S2RZ_W2= new HRes1DHit("S2RZ_W2");
00077 hRes_S3RZ_W2= new HRes1DHit("S3RZ_W2");
00078
00079 hEff_S1RPhi= new HEff1DHit("S1RPhi");
00080 hEff_S2RPhi= new HEff1DHit("S2RPhi");
00081 hEff_S3RPhi= new HEff1DHit("S3RPhi");
00082
00083 hEff_S1RZ= new HEff1DHit("S1RZ");
00084 hEff_S2RZ= new HEff1DHit("S2RZ");
00085 hEff_S3RZ= new HEff1DHit("S3RZ");
00086
00087 hEff_S1RZ_W0= new HEff1DHit("S1RZ_W0");
00088 hEff_S2RZ_W0= new HEff1DHit("S2RZ_W0");
00089 hEff_S3RZ_W0= new HEff1DHit("S3RZ_W0");
00090
00091 hEff_S1RZ_W1= new HEff1DHit("S1RZ_W1");
00092 hEff_S2RZ_W1= new HEff1DHit("S2RZ_W1");
00093 hEff_S3RZ_W1= new HEff1DHit("S3RZ_W1");
00094
00095 hEff_S1RZ_W2= new HEff1DHit("S1RZ_W2");
00096 hEff_S2RZ_W2= new HEff1DHit("S2RZ_W2");
00097 hEff_S3RZ_W2= new HEff1DHit("S3RZ_W2");
00098 }
00099
00100
00101
00102
00103 DTRecHitQuality::~DTRecHitQuality(){
00104 }
00105
00106
00107
00108 void DTRecHitQuality::endJob() {
00109
00110 theFile->cd();
00111
00112 hEff_S1RPhi->ComputeEfficiency();
00113 hEff_S2RPhi->ComputeEfficiency();
00114 hEff_S3RPhi->ComputeEfficiency();
00115
00116 hEff_S1RZ->ComputeEfficiency();
00117 hEff_S2RZ->ComputeEfficiency();
00118 hEff_S3RZ->ComputeEfficiency();
00119
00120 hEff_S1RZ_W0->ComputeEfficiency();
00121 hEff_S2RZ_W0->ComputeEfficiency();
00122 hEff_S3RZ_W0->ComputeEfficiency();
00123
00124 hEff_S1RZ_W1->ComputeEfficiency();
00125 hEff_S2RZ_W1->ComputeEfficiency();
00126 hEff_S3RZ_W1->ComputeEfficiency();
00127
00128 hEff_S1RZ_W2->ComputeEfficiency();
00129 hEff_S2RZ_W2->ComputeEfficiency();
00130 hEff_S3RZ_W2->ComputeEfficiency();
00131
00132
00133 hRes_S1RPhi->Write();
00134 hRes_S2RPhi->Write();
00135 hRes_S3RPhi->Write();
00136
00137 hRes_S1RZ->Write();
00138 hRes_S2RZ->Write();
00139 hRes_S3RZ->Write();
00140
00141 hRes_S1RZ_W0->Write();
00142 hRes_S2RZ_W0->Write();
00143 hRes_S3RZ_W0->Write();
00144
00145 hRes_S1RZ_W1->Write();
00146 hRes_S2RZ_W1->Write();
00147 hRes_S3RZ_W1->Write();
00148
00149 hRes_S1RZ_W2->Write();
00150 hRes_S2RZ_W2->Write();
00151 hRes_S3RZ_W2->Write();
00152
00153
00154 hEff_S1RPhi->Write();
00155 hEff_S2RPhi->Write();
00156 hEff_S3RPhi->Write();
00157
00158 hEff_S1RZ->Write();
00159 hEff_S2RZ->Write();
00160 hEff_S3RZ->Write();
00161
00162 hEff_S1RZ_W0->Write();
00163 hEff_S2RZ_W0->Write();
00164 hEff_S3RZ_W0->Write();
00165
00166 hEff_S1RZ_W1->Write();
00167 hEff_S2RZ_W1->Write();
00168 hEff_S3RZ_W1->Write();
00169
00170 hEff_S1RZ_W2->Write();
00171 hEff_S2RZ_W2->Write();
00172 hEff_S3RZ_W2->Write();
00173
00174 theFile->Close();
00175 }
00176
00177
00178 void DTRecHitQuality::analyze(const Event & event, const EventSetup& eventSetup){
00179 if(debug)
00180 cout << "--- [DTRecHitQuality] Analysing Event: #Run: " << event.id().run()
00181 << " #Event: " << event.id().event() << endl;
00182 theFile->cd();
00183
00184 ESHandle<DTGeometry> dtGeom;
00185 eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00186
00187
00188
00189 Handle<PSimHitContainer> simHits;
00190 event.getByLabel(simHitLabel, "MuonDTHits", simHits);
00191
00192
00193 map<DTWireId, PSimHitContainer > simHitsPerWire =
00194 DTHitQualityUtils::mapSimHitsPerWire(*(simHits.product()));
00195
00196
00197
00198
00199
00200
00201 if(doStep1) {
00202 if(debug)
00203 cout << " -- DTRecHit S1: begin analysis:" << endl;
00204
00205 Handle<DTRecHitCollection> dtRecHits;
00206 event.getByLabel(recHitLabel, dtRecHits);
00207
00208
00209 map<DTWireId,vector<DTRecHit1DPair> > recHitsPerWire =
00210 map1DRecHitsPerWire(dtRecHits.product());
00211
00212 compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 1);
00213 }
00214
00215
00216
00217
00218
00219 if(doStep2) {
00220 if(debug)
00221 cout << " -- DTRecHit S2: begin analysis:" << endl;
00222
00223
00224 Handle<DTRecSegment2DCollection> segment2Ds;
00225 event.getByLabel(segment2DLabel, segment2Ds);
00226
00227
00228 map<DTWireId,vector<DTRecHit1D> > recHitsPerWire =
00229 map1DRecHitsPerWire(segment2Ds.product());
00230
00231 compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 2);
00232 }
00233
00234
00235
00236
00237 if(doStep3) {
00238 if(debug)
00239 cout << " -- DTRecHit S3: begin analysis:" << endl;
00240
00241
00242 Handle<DTRecSegment4DCollection> segment4Ds;
00243 event.getByLabel(segment4DLabel, segment4Ds);
00244
00245
00246 map<DTWireId,vector<DTRecHit1D> > recHitsPerWire =
00247 map1DRecHitsPerWire(segment4Ds.product());
00248
00249 compute(dtGeom.product(), simHitsPerWire, recHitsPerWire, 3);
00250 }
00251
00252 }
00253
00254
00255
00256
00257 map<DTWireId, vector<DTRecHit1DPair> >
00258 DTRecHitQuality::map1DRecHitsPerWire(const DTRecHitCollection* dt1DRecHitPairs) {
00259 map<DTWireId, vector<DTRecHit1DPair> > ret;
00260
00261 for(DTRecHitCollection::const_iterator rechit = dt1DRecHitPairs->begin();
00262 rechit != dt1DRecHitPairs->end(); rechit++) {
00263 ret[(*rechit).wireId()].push_back(*rechit);
00264 }
00265
00266 return ret;
00267 }
00268
00269
00270
00271 map<DTWireId, vector<DTRecHit1D> >
00272 DTRecHitQuality::map1DRecHitsPerWire(const DTRecSegment2DCollection* segment2Ds) {
00273 map<DTWireId, vector<DTRecHit1D> > ret;
00274
00275
00276 for(DTRecSegment2DCollection::const_iterator segment = segment2Ds->begin();
00277 segment != segment2Ds->end();
00278 segment++) {
00279 vector<DTRecHit1D> component1DHits= (*segment).specificRecHits();
00280
00281 for(vector<DTRecHit1D>::const_iterator hit = component1DHits.begin();
00282 hit != component1DHits.end();
00283 hit++) {
00284 ret[(*hit).wireId()].push_back(*hit);
00285 }
00286 }
00287 return ret;
00288 }
00289
00290
00291
00292
00293 map<DTWireId, std::vector<DTRecHit1D> >
00294 DTRecHitQuality::map1DRecHitsPerWire(const DTRecSegment4DCollection* segment4Ds) {
00295 map<DTWireId, vector<DTRecHit1D> > ret;
00296
00297 for(DTRecSegment4DCollection::const_iterator segment = segment4Ds->begin();
00298 segment != segment4Ds->end();
00299 segment++) {
00300
00301 vector<const TrackingRecHit*> segment2Ds = (*segment).recHits();
00302
00303 for(vector<const TrackingRecHit*>::const_iterator segment2D = segment2Ds.begin();
00304 segment2D != segment2Ds.end();
00305 segment2D++) {
00306
00307 vector<const TrackingRecHit*> hits = (*segment2D)->recHits();
00308
00309 for(vector<const TrackingRecHit*>::const_iterator hit = hits.begin();
00310 hit != hits.end(); hit++) {
00311 const DTRecHit1D* hit1D = dynamic_cast<const DTRecHit1D*>(*hit);
00312 ret[hit1D->wireId()].push_back(*hit1D);
00313 }
00314 }
00315 }
00316
00317 return ret;
00318 }
00319
00320
00321 float DTRecHitQuality::simHitDistFromWire(const DTLayer* layer,
00322 DTWireId wireId,
00323 const PSimHit& hit) {
00324 float xwire = layer->specificTopology().wirePosition(wireId.wire());
00325 LocalPoint entryP = hit.entryPoint();
00326 LocalPoint exitP = hit.exitPoint();
00327 float xEntry = entryP.x()-xwire;
00328 float xExit = exitP.x()-xwire;
00329
00330 return fabs(xEntry - (entryP.z()*(xExit-xEntry))/(exitP.z()-entryP.z()));
00331 }
00332
00333
00334 float DTRecHitQuality::simHitImpactAngle(const DTLayer* layer,
00335 DTWireId wireId,
00336 const PSimHit& hit) {
00337 LocalPoint entryP = hit.entryPoint();
00338 LocalPoint exitP = hit.exitPoint();
00339 float theta=(exitP.x()-entryP.x())/(exitP.z()-entryP.z());
00340 return atan(theta);
00341 }
00342
00343
00344 float DTRecHitQuality::simHitDistFromFE(const DTLayer* layer,
00345 DTWireId wireId,
00346 const PSimHit& hit) {
00347 LocalPoint entryP = hit.entryPoint();
00348 LocalPoint exitP = hit.exitPoint();
00349 float wireLenght=layer->specificTopology().cellLenght();
00350 return (entryP.y()+exitP.y())/2.+wireLenght;
00351 }
00352
00353
00354
00355 template <typename type>
00356 const type*
00357 DTRecHitQuality::findBestRecHit(const DTLayer* layer,
00358 DTWireId wireId,
00359 const vector<type>& recHits,
00360 const float simHitDist) {
00361 float res = 99999;
00362 const type* theBestRecHit = 0;
00363
00364 for(typename vector<type>::const_iterator recHit = recHits.begin();
00365 recHit != recHits.end();
00366 recHit++) {
00367 float distTmp = recHitDistFromWire(*recHit, layer);
00368 if(fabs(distTmp-simHitDist) < res) {
00369 res = fabs(distTmp-simHitDist);
00370 theBestRecHit = &(*recHit);
00371 }
00372 }
00373
00374 return theBestRecHit;
00375 }
00376
00377
00378
00379 float
00380 DTRecHitQuality::recHitDistFromWire(const DTRecHit1DPair& hitPair, const DTLayer* layer) {
00381
00382 return fabs(hitPair.localPosition(DTEnums::Left).x() -
00383 hitPair.localPosition(DTEnums::Right).x())/2.;
00384 }
00385
00386
00387
00388
00389 float
00390 DTRecHitQuality::recHitDistFromWire(const DTRecHit1D& recHit, const DTLayer* layer) {
00391 return fabs(recHit.localPosition().x() - layer->specificTopology().wirePosition(recHit.wireId().wire()));
00392 }
00393
00394
00395 template <typename type>
00396 void DTRecHitQuality::compute(const DTGeometry *dtGeom,
00397 std::map<DTWireId, std::vector<PSimHit> > simHitsPerWire,
00398 std::map<DTWireId, std::vector<type> > recHitsPerWire,
00399 int step) {
00400
00401 for(map<DTWireId, vector<PSimHit> >::const_iterator wireAndSHits = simHitsPerWire.begin();
00402 wireAndSHits != simHitsPerWire.end();
00403 wireAndSHits++) {
00404 DTWireId wireId = (*wireAndSHits).first;
00405 vector<PSimHit> simHitsInCell = (*wireAndSHits).second;
00406
00407
00408 const DTLayer* layer = dtGeom->layer(wireId);
00409
00410
00411 const PSimHit* muSimHit = DTHitQualityUtils::findMuSimHit(simHitsInCell);
00412 if (muSimHit==0) {
00413 if (debug)
00414 cout << " No mu SimHit in channel: " << wireId << ", skipping! " << endl;
00415 continue;
00416 }
00417
00418
00419 float simHitWireDist = simHitDistFromWire(layer, wireId, *muSimHit);
00420
00421 if(simHitWireDist>2.1) {
00422 if(debug)
00423 cout << " [DTRecHitQuality]###Warning: The mu SimHit in out of the cell, skipping!" << endl;
00424 continue;
00425 }
00426 GlobalPoint simHitGlobalPos = layer->toGlobal(muSimHit->localPosition());
00427
00428
00429 float simHitTheta = simHitImpactAngle(layer, wireId, *muSimHit);
00430
00431
00432 float simHitFEDist = simHitDistFromFE(layer, wireId, *muSimHit);
00433
00434 bool recHitReconstructed = false;
00435
00436
00437 if(recHitsPerWire.find(wireId) == recHitsPerWire.end()) {
00438
00439 if(debug)
00440 cout << " No RecHit found at Step: " << step << " in cell: " << wireId << endl;
00441 } else {
00442 recHitReconstructed = true;
00443
00444 vector<type> recHits = recHitsPerWire[wireId];
00445 if(debug)
00446 cout << " " << recHits.size() << " RecHits, Step " << step << " in channel: " << wireId << endl;
00447
00448
00449 const type* theBestRecHit = findBestRecHit(layer, wireId, recHits, simHitWireDist);
00450
00451
00452 float recHitWireDist = recHitDistFromWire(*theBestRecHit, layer);
00453 if(debug)
00454 cout << " SimHit distance from wire: " << simHitWireDist << endl
00455 << " SimHit distance from FE: " << simHitFEDist << endl
00456 << " SimHit distance angle " << simHitTheta << endl
00457 << " RecHit distance from wire: " << recHitWireDist << endl;
00458 float recHitErr = recHitPositionError(*theBestRecHit);
00459
00460 HRes1DHit *hRes = 0;
00461 HRes1DHit *hResTot = 0;
00462
00463
00464
00465 if(step == 1) {
00466
00467 if(wireId.superLayer() != 2) {
00468 hRes = hRes_S1RPhi;
00469 } else {
00470 hResTot = hRes_S1RZ;
00471 if(wireId.wheel() == 0)
00472 hRes = hRes_S1RZ_W0;
00473 if(abs(wireId.wheel()) == 1)
00474 hRes = hRes_S1RZ_W1;
00475 if(abs(wireId.wheel()) == 2)
00476 hRes = hRes_S1RZ_W2;
00477 }
00478
00479 } else if(step == 2) {
00480
00481 if(wireId.superlayer() != 2) {
00482 hRes = hRes_S2RPhi;
00483 } else {
00484 hResTot = hRes_S2RZ;
00485 if(wireId.wheel() == 0)
00486 hRes = hRes_S2RZ_W0;
00487 if(abs(wireId.wheel()) == 1)
00488 hRes = hRes_S2RZ_W1;
00489 if(abs(wireId.wheel()) == 2)
00490 hRes = hRes_S2RZ_W2;
00491 }
00492
00493 } else if(step == 3) {
00494
00495 if(wireId.superlayer() != 2) {
00496 hRes = hRes_S3RPhi;
00497 } else {
00498 hResTot = hRes_S3RZ;
00499 if(wireId.wheel() == 0)
00500 hRes = hRes_S3RZ_W0;
00501 if(abs(wireId.wheel()) == 1)
00502 hRes = hRes_S3RZ_W1;
00503 if(abs(wireId.wheel()) == 2)
00504 hRes = hRes_S3RZ_W2;
00505 }
00506
00507 }
00508
00509 hRes->Fill(simHitWireDist, simHitTheta, simHitFEDist, recHitWireDist, simHitGlobalPos.eta(),
00510 simHitGlobalPos.phi(),recHitErr);
00511 if(hResTot != 0)
00512 hResTot->Fill(simHitWireDist, simHitTheta, simHitFEDist, recHitWireDist, simHitGlobalPos.eta(),
00513 simHitGlobalPos.phi(),recHitErr);
00514
00515 }
00516
00517
00518 HEff1DHit *hEff = 0;
00519 HEff1DHit *hEffTot = 0;
00520
00521 if(step == 1) {
00522
00523 if(wireId.superlayer() != 2) {
00524 hEff = hEff_S1RPhi;
00525 } else {
00526 hEffTot = hEff_S1RZ;
00527 if(wireId.wheel() == 0)
00528 hEff = hEff_S1RZ_W0;
00529 if(abs(wireId.wheel()) == 1)
00530 hEff = hEff_S1RZ_W1;
00531 if(abs(wireId.wheel()) == 2)
00532 hEff = hEff_S1RZ_W2;
00533 }
00534
00535 } else if(step == 2) {
00536
00537 if(wireId.superlayer() != 2) {
00538 hEff = hEff_S2RPhi;
00539 } else {
00540 hEffTot = hEff_S2RZ;
00541 if(wireId.wheel() == 0)
00542 hEff = hEff_S2RZ_W0;
00543 if(abs(wireId.wheel()) == 1)
00544 hEff = hEff_S2RZ_W1;
00545 if(abs(wireId.wheel()) == 2)
00546 hEff = hEff_S2RZ_W2;
00547 }
00548
00549 } else if(step == 3) {
00550
00551 if(wireId.superlayer() != 2) {
00552 hEff = hEff_S3RPhi;
00553 } else {
00554 hEffTot = hEff_S3RZ;
00555 if(wireId.wheel() == 0)
00556 hEff = hEff_S3RZ_W0;
00557 if(abs(wireId.wheel()) == 1)
00558 hEff = hEff_S3RZ_W1;
00559 if(abs(wireId.wheel()) == 2)
00560 hEff = hEff_S3RZ_W2;
00561 }
00562
00563 }
00564
00565 hEff->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
00566 if(hEffTot != 0)
00567 hEffTot->Fill(simHitWireDist, simHitGlobalPos.eta(), simHitGlobalPos.phi(), recHitReconstructed);
00568 }
00569 }
00570
00571
00572 float DTRecHitQuality::recHitPositionError(const DTRecHit1DPair& recHit) {
00573 return sqrt(recHit.localPositionError(DTEnums::Left).xx());
00574 }
00575
00576
00577 float DTRecHitQuality::recHitPositionError(const DTRecHit1D& recHit) {
00578 return sqrt(recHit.localPositionError().xx());
00579 }