CMS 3D CMS Logo

DTRecHitQuality.cc

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

Generated on Tue Jun 9 17:49:02 2009 for CMSSW by  doxygen 1.5.4