CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/src/DQM/DTMonitorModule/src/DTCalibValidation.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2011/11/12 09:18:41 $
00006  *  $Revision: 1.15 $
00007  *  \author G. Mila - INFN Torino
00008  */
00009 
00010 #include "DQM/DTMonitorModule/interface/DTCalibValidation.h"
00011 
00012 // Framework
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/ServiceRegistry/interface/Service.h"
00017 #include "DQMServices/Core/interface/DQMStore.h"
00018 #include "DQMServices/Core/interface/MonitorElement.h"
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 
00021 //Geometry
00022 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00023 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00024 
00025 //RecHit
00026 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00027 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
00028 
00029 
00030 #include <iterator>
00031 
00032 using namespace edm;
00033 using namespace std;
00034 
00035 
00036 DTCalibValidation::DTCalibValidation(const ParameterSet& pset) {
00037 
00038   //debug = pset.getUntrackedParameter<bool>("debug",false);
00039 
00040   // Get the DQM needed services
00041 
00042   // To remove into CMSSW versions before 20X
00043   theDbe = edm::Service<DQMStore>().operator->();
00044   // To add into CMSSW versions before 20X
00045   /*theDbe = edm::Service<DaqMonitorBEInterface>().operator->();
00046   theDbe->setVerbose(1);
00047   edm::Service<MonitorDaemon>().operator->();*/
00048 
00049   theDbe->setCurrentFolder("DT/DTCalibValidation");
00050 
00051   parameters = pset;
00052 }
00053 
00054 
00055 DTCalibValidation::~DTCalibValidation(){
00056 }
00057 
00058 
00059 void DTCalibValidation::beginJob(){
00060 
00061   // the name of the rechits collection at step 1
00062   recHits1DLabel = parameters.getUntrackedParameter<string>("recHits1DLabel");
00063   // the name of the 2D segments
00064   segment2DLabel = parameters.getUntrackedParameter<string>("segment2DLabel");
00065   // the name of the 4D segments
00066   segment4DLabel = parameters.getUntrackedParameter<string>("segment4DLabel");
00067   // the counter of segments not used to compute residuals
00068   wrongSegment = 0;
00069   // the counter of segments used to compute residuals
00070   rightSegment = 0;
00071   // the analysis type
00072   detailedAnalysis = parameters.getUntrackedParameter<bool>("detailedAnalysis",false);
00073 
00074   nevent=0;
00075 
00076 }
00077 
00078 void DTCalibValidation::beginRun(const Run& run, const EventSetup& setup) {
00079   
00080   // get the geometry
00081   setup.get<MuonGeometryRecord>().get(dtGeom);
00082 
00083   // Loop over all the chambers          
00084   vector<DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();         
00085   vector<DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();          
00086   for (; ch_it != ch_end; ++ch_it) {     
00087     vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();         
00088     vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();          
00089     // Loop over the SLs         
00090     for(; sl_it != sl_end; ++sl_it) { 
00091       DTSuperLayerId slId = (*sl_it)->id();
00092       if(detailedAnalysis){
00093          // Loop over the 3 steps
00094          for(int step = 1; step <= 3; ++step) bookHistos(slId,step);
00095       } else {
00096          // Only step 3
00097          bookHistos(slId,3);
00098       }
00099     }
00100   }
00101 
00102 }
00103 
00104 
00105 void DTCalibValidation::endJob(){
00106 
00107  LogVerbatim("DTCalibValidation") << "Segments used to compute residuals: " << rightSegment;
00108  LogVerbatim("DTCalibValidation") << "Segments not used to compute residuals: " << wrongSegment;
00109 
00110  //theDbe->showDirStructure();
00111  bool outputMEsInRootFile = parameters.getParameter<bool>("OutputMEsInRootFile");
00112  std::string outputFileName = parameters.getParameter<std::string>("OutputFileName");
00113  if(outputMEsInRootFile){
00114    theDbe->showDirStructure();
00115    theDbe->save(outputFileName);
00116  }
00117  
00118  theDbe->rmdir("DT/DTCalibValidation");
00119 }
00120   
00121 
00122 
00123 void DTCalibValidation::analyze(const edm::Event& event, const edm::EventSetup& setup) {
00124 
00125   ++nevent;
00126   LogTrace("DTCalibValidation") << "[DTCalibValidation] Analyze #Run: " << event.id().run()
00127                                 << " #Event: " << nevent;
00128 
00129   // RecHit mapping at Step 1 -------------------------------
00130   map<DTWireId,vector<DTRecHit1DPair> > recHitsPerWire_1S;
00131 
00132   // RecHit mapping at Step 2 ------------------------------
00133   map<DTWireId,vector<DTRecHit1D> > recHitsPerWire_2S;
00134 
00135   if(detailedAnalysis){
00136      LogTrace("DTCalibValidation") << "  -- DTRecHit S1: begin analysis:";
00137      // Get the rechit collection from the event
00138      Handle<DTRecHitCollection> dtRecHits;
00139      event.getByLabel(recHits1DLabel, dtRecHits);
00140      recHitsPerWire_1S = map1DRecHitsPerWire(dtRecHits.product());
00141  
00142      LogTrace("DTCalibValidation") << "  -- DTRecHit S2: begin analysis:";
00143      // Get the 2D rechits from the event
00144      Handle<DTRecSegment2DCollection> segment2Ds;
00145      event.getByLabel(segment2DLabel, segment2Ds);
00146      recHitsPerWire_2S =  map1DRecHitsPerWire(segment2Ds.product());
00147   }
00148 
00149   // RecHit mapping at Step 3 ---------------------------------
00150   LogTrace("DTCalibValidation") << "  -- DTRecHit S3: begin analysis:";
00151   // Get the 4D rechits from the event
00152   Handle<DTRecSegment4DCollection> segment4Ds;
00153   event.getByLabel(segment4DLabel, segment4Ds);
00154   map<DTWireId,vector<DTRecHit1D> > recHitsPerWire_3S =  map1DRecHitsPerWire(segment4Ds.product());
00155  
00156   
00157   // Loop over all 4D segments
00158   for(DTRecSegment4DCollection::const_iterator segment = segment4Ds->begin();
00159       segment != segment4Ds->end();
00160       ++segment) {
00161 
00162     if(detailedAnalysis){
00163        LogTrace("DTCalibValidation") << "Anlysis on recHit at step 1";
00164        compute(dtGeom.product(), (*segment), recHitsPerWire_1S, 1);
00165 
00166        LogTrace("DTCalibValidation") << "Anlysis on recHit at step 2";
00167        compute(dtGeom.product(), (*segment), recHitsPerWire_2S, 2);
00168     }
00169 
00170     LogTrace("DTCalibValidation") << "Anlysis on recHit at step 3";
00171     compute(dtGeom.product(), (*segment), recHitsPerWire_3S, 3);
00172   }
00173 
00174 }
00175 
00176 
00177 // Return a map between DTRecHit1DPair and wireId
00178 map<DTWireId, vector<DTRecHit1DPair> > 
00179 DTCalibValidation::map1DRecHitsPerWire(const DTRecHitCollection* dt1DRecHitPairs) {
00180   map<DTWireId, vector<DTRecHit1DPair> > ret;
00181 
00182   for(DTRecHitCollection::const_iterator rechit = dt1DRecHitPairs->begin();
00183       rechit != dt1DRecHitPairs->end(); ++rechit) {
00184     ret[(*rechit).wireId()].push_back(*rechit);
00185   }
00186 
00187   return ret;
00188 }
00189 
00190         
00191 // Return a map between DTRecHit1D at S2 and wireId
00192 map<DTWireId, vector<DTRecHit1D> >
00193 DTCalibValidation::map1DRecHitsPerWire(const DTRecSegment2DCollection* segment2Ds) {
00194   map<DTWireId, vector<DTRecHit1D> > ret;
00195 
00196   // Loop over all 2D segments
00197   for(DTRecSegment2DCollection::const_iterator segment = segment2Ds->begin();
00198       segment != segment2Ds->end();
00199       ++segment) {
00200     vector<DTRecHit1D> component1DHits= (*segment).specificRecHits();
00201     // Loop over all component 1D hits
00202     for(vector<DTRecHit1D>::const_iterator hit = component1DHits.begin();
00203         hit != component1DHits.end(); ++hit) {
00204       ret[(*hit).wireId()].push_back(*hit);
00205     }
00206   }
00207   return ret;
00208 }
00209 
00210 // Return a map between DTRecHit1D at S3 and wireId
00211 map<DTWireId, std::vector<DTRecHit1D> >
00212 DTCalibValidation::map1DRecHitsPerWire(const DTRecSegment4DCollection* segment4Ds) {
00213   map<DTWireId, vector<DTRecHit1D> > ret;
00214   // Loop over all 4D segments
00215   for(DTRecSegment4DCollection::const_iterator segment = segment4Ds->begin();
00216       segment != segment4Ds->end();
00217       ++segment) {
00218     // Get component 2D segments
00219     vector<const TrackingRecHit*> segment2Ds = (*segment).recHits();
00220     // Loop over 2D segments:
00221     for(vector<const TrackingRecHit*>::const_iterator segment2D = segment2Ds.begin();
00222         segment2D != segment2Ds.end();
00223         ++segment2D) {
00224       // Get 1D component rechits
00225       vector<const TrackingRecHit*> hits = (*segment2D)->recHits();
00226       // Loop over them
00227       for(vector<const TrackingRecHit*>::const_iterator hit = hits.begin();
00228           hit != hits.end(); ++hit) {
00229         const DTRecHit1D* hit1D = dynamic_cast<const DTRecHit1D*>(*hit);
00230         ret[hit1D->wireId()].push_back(*hit1D);
00231       }
00232     }
00233   }
00234 
00235   return ret;
00236 }
00237 
00238     
00239 
00240 // Find the RecHit closest to the segment4D
00241 template  <typename type>
00242 const type* 
00243 DTCalibValidation::findBestRecHit(const DTLayer* layer,
00244                                 DTWireId wireId,
00245                                 const vector<type>& recHits,
00246                                 const float segmDist) {
00247   float res = 99999;
00248   const type* theBestRecHit = 0;
00249   // Loop over RecHits within the cell
00250   for(typename vector<type>::const_iterator recHit = recHits.begin();
00251       recHit != recHits.end();
00252       ++recHit) {
00253     float distTmp = recHitDistFromWire(*recHit, layer);
00254     if(fabs(distTmp-segmDist) < res) {
00255       res = fabs(distTmp-segmDist);
00256       theBestRecHit = &(*recHit);
00257     }
00258   } // End of loop over RecHits within the cell
00259 
00260   return theBestRecHit;
00261 }
00262 
00263 
00264 // Compute the distance from wire (cm) of a hits in a DTRecHit1DPair
00265 float 
00266 DTCalibValidation::recHitDistFromWire(const DTRecHit1DPair& hitPair, const DTLayer* layer) {
00267   return fabs(hitPair.localPosition(DTEnums::Left).x() -
00268               hitPair.localPosition(DTEnums::Right).x())/2.;
00269 }
00270 
00271 
00272 // Compute the distance from wire (cm) of a hits in a DTRecHit1D
00273 float 
00274 DTCalibValidation::recHitDistFromWire(const DTRecHit1D& recHit, const DTLayer* layer) {
00275   return fabs(recHit.localPosition().x() - layer->specificTopology().wirePosition(recHit.wireId().wire()));
00276   
00277   //to check the compatibility position / distance
00278   /* // Get the recHit and the wire position
00279   const DTChamber* chamber = (*layer).chamber();
00280   GlobalPoint recHitPosGlob = layer->toGlobal(recHit.localPosition());
00281   LocalPoint recHitPosInChamber = chamber->toLocal(recHitPosGlob);
00282   float wireX = layer->specificTopology().wirePosition(recHit.wireId().wire());
00283   LocalPoint wirePosInLay(wireX,recHit.localPosition().y(),0);
00284   GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
00285   LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
00286   
00287   float recHitDist = -1;
00288   if(recHit.wireId().layerId().superlayerId().superlayer() != 2)
00289     recHitDist = fabs(recHitPosInChamber.x()-wirePosInChamber.x());
00290   else
00291     recHitDist = fabs(recHitPosInChamber.y()-wirePosInChamber.y());
00292   
00293     return recHitDist; */
00294   
00295 }
00296 
00297 
00298 
00299 
00300 // Compute the position (cm) of a hits in a DTRecHit1DPair
00301 float 
00302 DTCalibValidation::recHitPosition(const DTRecHit1DPair& hitPair, const DTLayer* layer, const DTChamber* chamber, float segmentPos, int sl) {
00303   
00304   // Get the layer and the wire position
00305   GlobalPoint hitPosGlob_right = layer->toGlobal(hitPair.localPosition(DTEnums::Right));
00306   LocalPoint hitPosInChamber_right = chamber->toLocal(hitPosGlob_right);
00307   GlobalPoint hitPosGlob_left = layer->toGlobal(hitPair.localPosition(DTEnums::Left));
00308   LocalPoint hitPosInChamber_left = chamber->toLocal(hitPosGlob_left);
00309   
00310   float recHitPos=-1;
00311   if(sl != 2){
00312     if(fabs(hitPosInChamber_left.x()-segmentPos)<fabs(hitPosInChamber_right.x()-segmentPos))
00313       recHitPos = hitPosInChamber_left.x();
00314     else
00315       recHitPos = hitPosInChamber_right.x();
00316   }
00317   else{
00318     if(fabs(hitPosInChamber_left.y()-segmentPos)<fabs(hitPosInChamber_right.y()-segmentPos))
00319       recHitPos = hitPosInChamber_left.y();
00320     else
00321       recHitPos = hitPosInChamber_right.y();
00322   }
00323   
00324   return recHitPos;
00325 }
00326 
00327 
00328 // Compute the position (cm) of a hits in a  DTRecHit1D
00329 float 
00330 DTCalibValidation::recHitPosition(const DTRecHit1D& recHit, const DTLayer* layer, const DTChamber* chamber, float segmentPos, int sl) {
00331   
00332   // Get the layer and the wire position
00333   GlobalPoint recHitPosGlob = layer->toGlobal(recHit.localPosition());
00334   LocalPoint recHitPosInChamber = chamber->toLocal(recHitPosGlob);
00335 
00336   float recHitPos = -1;
00337   if(sl != 2)
00338     recHitPos = recHitPosInChamber.x();
00339   else
00340     recHitPos = recHitPosInChamber.y();
00341 
00342   return recHitPos;
00343 
00344 }
00345 
00346 
00347 
00348 // Compute the residuals    
00349 template  <typename type>
00350 void DTCalibValidation::compute(const DTGeometry *dtGeom,
00351                               const DTRecSegment4D& segment,
00352                               std::map<DTWireId, std::vector<type> > recHitsPerWire,
00353                               int step) {
00354   bool computeResidual = true;
00355   
00356   // Get all 1D RecHits at step 3 within the 4D segment
00357   vector<DTRecHit1D> recHits1D_S3;
00358   
00359   // Get 1D RecHits at Step 3 and select only events with
00360   // 8 hits in phi and 4 hits in theta (if any)
00361   const DTChamberRecSegment2D* phiSeg = segment.phiSegment();
00362   if(phiSeg){
00363     vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
00364     if(phiRecHits.size() != 8) {
00365       LogTrace("DTCalibValidation") << "[DTCalibValidation] Phi segments has: " << phiRecHits.size()
00366                                     << " hits, skipping"; // FIXME: info output
00367       computeResidual = false;
00368     }
00369     copy(phiRecHits.begin(), phiRecHits.end(), back_inserter(recHits1D_S3));
00370   }
00371   if(!phiSeg){
00372     LogTrace("DTCalibValidation") << " [DTCalibValidation] 4D segment has not the phi segment! ";
00373     computeResidual = false;
00374   }
00375   
00376   if(segment.dimension() == 4) {
00377     const DTSLRecSegment2D* zSeg = segment.zSegment();
00378     if(zSeg){
00379       vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
00380       if(zRecHits.size() != 4) {
00381         LogTrace("DTCalibValidation") << "[DTCalibValidation] Theta segments has: " << zRecHits.size()
00382                                       << " hits, skipping"; // FIXME: info output
00383         computeResidual = false;
00384       }
00385       copy(zRecHits.begin(), zRecHits.end(), back_inserter(recHits1D_S3));
00386     }
00387     if(!zSeg){
00388       LogTrace("DTCalibValidation") << " [DTCalibValidation] 4D segment has not the z segment! ";
00389       computeResidual = false;
00390     }
00391   }
00392 
00393   if(!computeResidual)
00394     ++wrongSegment;
00395   if(computeResidual){
00396     ++rightSegment;
00397     // Loop over 1D RecHit inside 4D segment
00398     for(vector<DTRecHit1D>::const_iterator recHit1D = recHits1D_S3.begin();
00399         recHit1D != recHits1D_S3.end();
00400         ++recHit1D) {
00401       const DTWireId wireId = (*recHit1D).wireId();
00402       
00403       // Get the layer and the wire position
00404       const DTLayer* layer = dtGeom->layer(wireId);
00405       float wireX = layer->specificTopology().wirePosition(wireId.wire());
00406       
00407       // Extrapolate the segment to the z of the wire
00408       // Get wire position in chamber RF
00409       // (y and z must be those of the hit to be coherent in the transf. of RF in case of rotations of the layer alignment)
00410       LocalPoint wirePosInLay(wireX,(*recHit1D).localPosition().y(),(*recHit1D).localPosition().z());
00411       GlobalPoint wirePosGlob = layer->toGlobal(wirePosInLay);
00412       const DTChamber* chamber = dtGeom->chamber((*recHit1D).wireId().layerId().chamberId());
00413       LocalPoint wirePosInChamber = chamber->toLocal(wirePosGlob);
00414       
00415       // Segment position at Wire z in chamber local frame
00416       LocalPoint segPosAtZWire = segment.localPosition()
00417         + segment.localDirection()*wirePosInChamber.z()/cos(segment.localDirection().theta());
00418       
00419       // Compute the distance of the segment from the wire
00420       int sl = wireId.superlayer();
00421       float SegmDistance = -1;
00422       if(sl == 1 || sl == 3) {
00423         // RPhi SL
00424         SegmDistance = fabs(wirePosInChamber.x() - segPosAtZWire.x());
00425         LogTrace("DTCalibValidation") << "SegmDistance: " << SegmDistance;
00426       } else if(sl == 2) {
00427         // RZ SL
00428         SegmDistance =  fabs(segPosAtZWire.y() - wirePosInChamber.y());
00429         LogTrace("DTCalibValidation") << "SegmDistance: " << SegmDistance;
00430       }
00431       if(SegmDistance > 2.1)
00432         LogTrace("DTCalibValidation") << "  Warning: dist segment-wire: " << SegmDistance;
00433 
00434       // Look for RecHits in the same cell
00435       if(recHitsPerWire.find(wireId) == recHitsPerWire.end()) {
00436         LogTrace("DTCalibValidation") << "   No RecHit found at Step: " << step << " in cell: " << wireId;
00437       } else {
00438         vector<type> recHits = recHitsPerWire[wireId];
00439         LogTrace("DTCalibValidation") << "   " << recHits.size() << " RecHits, Step " << step << " in channel: " << wireId;
00440         
00441         // Get the layer
00442         const DTLayer* layer = dtGeom->layer(wireId);
00443         // Find the best RecHits
00444         const type* theBestRecHit = findBestRecHit(layer, wireId, recHits, SegmDistance);
00445         // Compute the distance of the recHit from the wire
00446         float recHitWireDist =  recHitDistFromWire(*theBestRecHit, layer);
00447         LogTrace("DTCalibValidation") << "recHitWireDist: " << recHitWireDist;
00448 
00449         // Compute the residuals
00450         float residualOnDistance = recHitWireDist - SegmDistance;
00451         LogTrace("DTCalibValidation") << "WireId: " << wireId << "  ResidualOnDistance: " << residualOnDistance;
00452         float residualOnPosition = -1;
00453         float recHitPos = -1;
00454         if(sl == 1 || sl == 3) {
00455           recHitPos = recHitPosition(*theBestRecHit, layer, chamber,  segPosAtZWire.x(), sl);
00456           residualOnPosition = recHitPos - segPosAtZWire.x();
00457         }
00458         else{
00459           recHitPos = recHitPosition(*theBestRecHit, layer, chamber,  segPosAtZWire.y(), sl);
00460           residualOnPosition = recHitPos - segPosAtZWire.y();
00461         }
00462         LogTrace("DTCalibValidation") << "WireId: " << wireId << "  ResidualOnPosition: " << residualOnPosition;
00463 
00464         // Fill the histos
00465         if(sl == 1 || sl == 3)
00466           fillHistos(wireId.superlayerId(), SegmDistance, residualOnDistance, (wirePosInChamber.x() - segPosAtZWire.x()), residualOnPosition, step);
00467         else
00468           fillHistos(wireId.superlayerId(), SegmDistance, residualOnDistance, (wirePosInChamber.y() - segPosAtZWire.y()), residualOnPosition, step);
00469           
00470         
00471       }
00472     }
00473   }
00474 
00475 }
00476 
00477 
00478 
00479 // Book a set of histograms for a given SL
00480 void DTCalibValidation::bookHistos(DTSuperLayerId slId, int step) {
00481   LogTrace("DTCalibValidation") << "   Booking histos for SL: " << slId;
00482 
00483   // Compose the chamber name
00484   stringstream wheel; wheel << slId.wheel();    
00485   stringstream station; station << slId.station();      
00486   stringstream sector; sector << slId.sector(); 
00487   stringstream superLayer; superLayer << slId.superlayer();
00488   // Define the step
00489   stringstream Step; Step << step;
00490 
00491   
00492   string slHistoName =
00493     "_STEP" + Step.str() +
00494     "_W" + wheel.str() +
00495     "_St" + station.str() +
00496     "_Sec" + sector.str() +
00497     "_SL" + superLayer.str();
00498   
00499   theDbe->setCurrentFolder("DT/DTCalibValidation/Wheel" + wheel.str() +
00500                            "/Station" + station.str() +
00501                            "/Sector" + sector.str());
00502   // Create the monitor elements
00503   vector<MonitorElement *> histos;
00504   // Note hte order matters
00505   histos.push_back(theDbe->book1D("hResDist"+slHistoName,
00506                                   "Residuals on the distance from wire (rec_hit - segm_extr) (cm)",
00507                                   200, -0.4, 0.4));
00508   histos.push_back(theDbe->book2D("hResDistVsDist"+slHistoName,
00509                                   "Residuals on the distance (cm) from wire (rec_hit - segm_extr) vs distance  (cm)",
00510                                   100, 0, 2.5, 200, -0.4, 0.4));
00511   if(detailedAnalysis){
00512     histos.push_back(theDbe->book1D("hResPos"+slHistoName,
00513                                     "Residuals on the position from wire (rec_hit - segm_extr) (cm)",
00514                                     200, -0.4, 0.4));
00515     histos.push_back(theDbe->book2D("hResPosVsPos"+slHistoName,
00516                                     "Residuals on the position (cm) from wire (rec_hit - segm_extr) vs distance  (cm)",
00517                                     200, -2.5, 2.5, 200, -0.4, 0.4));
00518   }
00519 
00520   histosPerSL[make_pair(slId, step)] = histos;
00521 }
00522 
00523 
00524 
00525 // Fill a set of histograms for a given SL 
00526 void DTCalibValidation::fillHistos(DTSuperLayerId slId,
00527                                      float distance,
00528                                      float residualOnDistance,
00529                                      float position,    
00530                                      float residualOnPosition,
00531                                      int step) {
00532   // FIXME: optimization of the number of searches
00533   vector<MonitorElement *> histos =  histosPerSL[make_pair(slId,step)];                          
00534   histos[0]->Fill(residualOnDistance);
00535   histos[1]->Fill(distance, residualOnDistance);
00536   if(detailedAnalysis){
00537     histos[2]->Fill(residualOnPosition);
00538     histos[3]->Fill(position, residualOnPosition);
00539   }
00540 }
00541