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