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