Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:27

0001 /** \file RecoAnalyzerRecHits.cc
0002 *  plots for RecHits
0003   *
0004   *  $Date: 2012/12/26 20:35:50 $
0005   *  $Revision: 1.11 $
0006   *  \author Maarten Thomas
0007  */
0008 
0009 #include "Alignment/LaserAlignment/test/RecoAnalyzer.h"
0010 #include "FWCore/Framework/interface/Event.h" 
0011 #include "FWCore/Framework/interface/ESHandle.h" 
0012 #include "FWCore/Framework/interface/EventSetup.h" 
0013 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h" 
0014 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h" 
0015 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h" 
0016 #include "DataFormats/DetId/interface/DetId.h" 
0017 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h" 
0018 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0019 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0020 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h" 
0021 
0022   void RecoAnalyzer::trackerRecHits(edm::Event const& theEvent, edm::EventSetup const& theSetup)
0023 {
0024   //Retrieve tracker topology from geometry
0025   edm::ESHandle<TrackerTopology> tTopoHandle;
0026   theSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
0027   const TrackerTopology* const tTopo = tTopoHandle.product();
0028 
0029 
0030   // access the Tracker
0031   edm::ESHandle<TrackerGeometry> theTrackerGeometry;
0032   theSetup.get<TrackerDigiGeometryRecord>().get(theTrackerGeometry);
0033   const TrackerGeometry& theTracker(*theTrackerGeometry);
0034 
0035   edm::Handle<SiStripMatchedRecHit2DCollection> rechitsMatchedHandle;
0036   edm::Handle<SiStripRecHit2DCollection> rechitsRPhiHandle;
0037   edm::Handle<SiStripRecHit2DCollection> rechitsStereoHandle;
0038 
0039   // get the RecHits from the event
0040   theEvent.getByLabel(theRecHitProducer,"matchedRecHit", rechitsMatchedHandle);
0041   theEvent.getByLabel(theRecHitProducer,"rphiRecHit", rechitsRPhiHandle);
0042   theEvent.getByLabel(theRecHitProducer,"stereoRecHit", rechitsStereoHandle);
0043 
0044   // the RecHit collections
0045   const SiStripMatchedRecHit2DCollection * theMatchedRecHitCollection = rechitsMatchedHandle.product();
0046   const SiStripRecHit2DCollection * theRPhiRecHitCollection = rechitsRPhiHandle.product();
0047   const SiStripRecHit2DCollection * theStereoRecHitCollection = rechitsStereoHandle.product();
0048 
0049   // loop over the detIds for each RecHit Collection
0050   for ( SiStripMatchedRecHit2DCollection::const_iterator det_iter = theMatchedRecHitCollection->begin(), det_end = theMatchedRecHitCollection->end();
0051         det_iter != det_end; ++det_iter) {
0052     SiStripMatchedRecHit2DCollection::DetSet rechitRange = *det_iter;
0053     DetId detid(rechitRange.detId());
0054       // get the DetUnit
0055     const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>(theTracker.idToDet(detid));
0056 
0057       // some variables we need later on in the program
0058     int theBeam     = 0;
0059     int theRing     = 0;
0060     std::string thePart  = "";
0061     int theTIBLayer = 0;
0062     int theTOBLayer = 0;
0063     int theTECWheel = 0;
0064     int theTOBStereoDet = 0;
0065 
0066     switch (detid.subdetId())
0067     {
0068       case StripSubdetector::TIB:
0069       {
0070         
0071         thePart = "TIB";
0072         theTIBLayer = tTopo->tibLayer(detid.rawId);
0073         break;
0074       }
0075       case StripSubdetector::TOB:
0076       {
0077         
0078         thePart = "TOB";
0079         theTOBLayer = tTopo->tobLayer(detid.rawId);
0080         theTOBStereoDet = tTopo->tobStereo(detid.rawId);
0081         break;
0082       }
0083       case StripSubdetector::TEC:
0084       {
0085         
0086 
0087       // is this module in TEC+ or TEC-?
0088         if (tTopo->tecSide(detid.rawId) == 1) { thePart = "TEC-"; }
0089         else if (tTopo->tecSide(detid.rawId) == 2) { thePart = "TEC+"; }
0090 
0091       // in which ring is this module?
0092         if ( theStripDet->surface().position().perp() > 55.0 && theStripDet->surface().position().perp() < 59.0 )
0093           { theRing = 4; } // Ring 4
0094         else if ( theStripDet->surface().position().perp() > 81.0 && theStripDet->surface().position().perp() < 85.0 )
0095           { theRing = 6; } // Ring 6
0096         else
0097           { theRing = -1; } // probably not a Laser Hit!
0098 
0099       // on which disk is this module
0100         theTECWheel = tTopo->tecWheel(detid.rawId);
0101         break;
0102       }
0103     }
0104 
0105       // which beam belongs these digis to
0106     if ( thePart == "TIB" && theTIBLayer == 4 )
0107     {
0108       if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTIB) 
0109         && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTIB))          { theBeam = 0; } // beam 0 
0110 
0111       else if ( (theStripDet->surface().position().phi() > 1.29 - theSearchPhiTIB) 
0112         && (theStripDet->surface().position().phi() < 1.29 + theSearchPhiTIB))     { theBeam = 1; } // beam 1
0113 
0114       else if ( (theStripDet->surface().position().phi() > 1.85 - theSearchPhiTIB) 
0115         && (theStripDet->surface().position().phi() < 1.85 + theSearchPhiTIB))     { theBeam = 2; } // beam 2
0116 
0117       else if ( (theStripDet->surface().position().phi() > 2.75 - theSearchPhiTIB) 
0118         && (theStripDet->surface().position().phi() < 2.75 + theSearchPhiTIB))     { theBeam = 3; } // beam 3
0119 
0120       else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTIB) 
0121         && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTIB))    { theBeam = 4; } // beam 4
0122 
0123       else if ( (theStripDet->surface().position().phi() > -2.00 - theSearchPhiTIB) 
0124         && (theStripDet->surface().position().phi() < -2.00 + theSearchPhiTIB))    { theBeam = 5; } // beam 5
0125 
0126       else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTIB) 
0127         && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTIB))    { theBeam = 6; } // beam 6
0128 
0129       else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTIB) 
0130         && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTIB))    { theBeam = 7; } // beam 7
0131       else
0132         { theBeam = -1; } // probably not a Laser Hit!
0133     }
0134     else if ( thePart == "TOB" && theTOBLayer == 1 )
0135     {
0136       if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTOB) 
0137         && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTOB))          { theBeam = 0; } // beam 0 
0138 
0139       else if ( (theStripDet->surface().position().phi() > 1.29 - theSearchPhiTOB) 
0140         && (theStripDet->surface().position().phi() < 1.29 + theSearchPhiTOB))     { theBeam = 1; } // beam 1
0141 
0142       else if ( (theStripDet->surface().position().phi() > 1.85 - theSearchPhiTOB)
0143         && (theStripDet->surface().position().phi() < 1.85 + theSearchPhiTOB))     { theBeam = 2; } // beam 2
0144 
0145       else if ( (theStripDet->surface().position().phi() > 2.75 - theSearchPhiTOB)
0146         && (theStripDet->surface().position().phi() < 2.75 + theSearchPhiTOB))     { theBeam = 3; } // beam 3
0147 
0148       else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTOB)
0149         && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTOB))    { theBeam = 4; } // beam 4
0150 
0151       else if ( (theStripDet->surface().position().phi() > -2.00 - theSearchPhiTOB)
0152         && (theStripDet->surface().position().phi() < -2.00 + theSearchPhiTOB))    { theBeam = 5; } // beam 5
0153 
0154       else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTOB)
0155         && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTOB))    { theBeam = 6; } // beam 6
0156 
0157       else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTOB)
0158         && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTOB))    { theBeam = 7; } // beam 7
0159       else
0160         { theBeam = -1; } // probably not a Laser Hit!
0161     }
0162     else if ( thePart == "TEC+" || thePart == "TEC-" )
0163     {
0164       if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTEC)
0165         && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTEC))          { theBeam = 0; } // beam 0 
0166 
0167       else if ( (theStripDet->surface().position().phi() > 1.18 - theSearchPhiTEC)
0168         && (theStripDet->surface().position().phi() < 1.18 + theSearchPhiTEC))     { theBeam = 1; } // beam 1
0169 
0170       else if ( (theStripDet->surface().position().phi() > 1.96 - theSearchPhiTEC)
0171         && (theStripDet->surface().position().phi() < 1.96 + theSearchPhiTEC))     { theBeam = 2; } // beam 2
0172 
0173       else if ( (theStripDet->surface().position().phi() > 2.74 - theSearchPhiTEC)
0174         && (theStripDet->surface().position().phi() < 2.74 + theSearchPhiTEC))     { theBeam = 3; } // beam 3
0175 
0176       else if ( (theStripDet->surface().position().phi() > -2.74 - theSearchPhiTEC)
0177         && (theStripDet->surface().position().phi() < -2.74 + theSearchPhiTEC))    { theBeam = 4; } // beam 4
0178 
0179       else if ( (theStripDet->surface().position().phi() > -1.96 - theSearchPhiTEC)
0180         && (theStripDet->surface().position().phi() < -1.96 + theSearchPhiTEC))    { theBeam = 5; } // beam 5
0181 
0182       else if ( (theStripDet->surface().position().phi() > -1.18 - theSearchPhiTEC)
0183         && (theStripDet->surface().position().phi() < -1.18 + theSearchPhiTEC))    { theBeam = 6; } // beam 6
0184 
0185       else if ( (theStripDet->surface().position().phi() > -0.39 - theSearchPhiTEC)
0186         && (theStripDet->surface().position().phi() < -0.39 + theSearchPhiTEC))    { theBeam = 7; } // beam 7
0187 
0188       else if ( (theStripDet->surface().position().phi() > 1.28 - theSearchPhiTEC)
0189         && (theStripDet->surface().position().phi() < 1.28 + theSearchPhiTEC))     { theBeam = 21; } // beam 1 TEC2TEC
0190 
0191       else if ( (theStripDet->surface().position().phi() > 1.84 - theSearchPhiTEC)
0192         && (theStripDet->surface().position().phi() < 1.84 + theSearchPhiTEC))     { theBeam = 22; } // beam 2 TEC2TEC
0193 
0194       else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTEC)
0195         && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTEC))    { theBeam = 24; } // beam 4 TEC2TEC
0196 
0197       else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTEC)
0198         && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTEC))    { theBeam = 26; } // beam 6 TEC2TEC
0199 
0200       else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTEC)
0201         && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTEC))    { theBeam = 27; } // beam 7 TEC2TEC
0202       else 
0203         { theBeam = -1; } // probably not a Laser Hit!
0204     }
0205 
0206 
0207       // get the RecHits
0208     SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitRangeIteratorBegin = rechitRange.begin();
0209     SiStripMatchedRecHit2DCollection::DetSet::const_iterator rechitRangeIteratorEnd = rechitRange.end();
0210     SiStripMatchedRecHit2DCollection::DetSet::const_iterator iRecHit = rechitRangeIteratorBegin;
0211       // loop on the RecHits
0212     for (; iRecHit != rechitRangeIteratorEnd; iRecHit++)
0213     {
0214       SiStripMatchedRecHit2D const rechit = *iRecHit;
0215       theRecHitPositionsX->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).x());
0216       theRecHitPositionsY->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).y());
0217       theRecHitPositionsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z());
0218       theRecHitPositionsYvsX->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).x(),
0219         theStripDet->surface().toGlobal(rechit.localPosition()).y());
0220       theRecHitPositionsPhivsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z(),
0221         theStripDet->surface().toGlobal(rechit.localPosition()).phi());
0222       theRecHitPositionsRvsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z(),
0223         theStripDet->surface().toGlobal(rechit.localPosition()).perp());
0224 
0225       if (thePart == "TEC+" || thePart == "TEC-")   
0226       {             
0227         double r_ = sqrt(pow(theStripDet->surface().toGlobal(rechit.localPosition()).x(),2) + pow(theStripDet->surface().toGlobal(rechit.localPosition()).y(),2));
0228         fillLaserBeamPlots(r_,theTECWheel,thePart,theRing,theBeam);
0229       }
0230     }
0231   }
0232 
0233   for ( SiStripRecHit2DCollection::const_iterator det_iter = theRPhiRecHitCollection->begin(), det_end = theRPhiRecHitCollection->end();
0234         det_iter != det_end; ++det_iter) {
0235     SiStripRecHit2DCollection::DetSet rechitRange = *det_iter;
0236     DetId detid(rechitRange.detId());
0237       // get the DetUnit
0238     const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>(theTracker.idToDet(detid));
0239 
0240         // some variables we need later on in the program
0241       int theBeam     = 0;
0242       int theRing     = 0;
0243       std::string thePart  = "";
0244       int theTIBLayer = 0;
0245       int theTOBLayer = 0;
0246       int theTECWheel = 0;
0247       int theTOBStereoDet = 0;
0248 
0249       switch (detid.subdetId())
0250       {
0251         case StripSubdetector::TIB:
0252         {
0253           
0254           thePart = "TIB";
0255           theTIBLayer = tTopo->tibLayer(detid.rawId);
0256           break;
0257         }
0258         case StripSubdetector::TOB:
0259         {
0260           
0261           thePart = "TOB";
0262           theTOBLayer = tTopo->tobLayer(detid.rawId);
0263           theTOBStereoDet = tTopo->tobStereo(detid.rawId);
0264           break;
0265         }
0266         case StripSubdetector::TEC:
0267         {
0268           
0269 
0270         // is this module in TEC+ or TEC-?
0271           if (tTopo->tecSide(detid.rawId) == 1) { thePart = "TEC-"; }
0272           else if (tTopo->tecSide(detid.rawId) == 2) { thePart = "TEC+"; }
0273 
0274         // in which ring is this module?
0275           if ( theStripDet->surface().position().perp() > 55.0 && theStripDet->surface().position().perp() < 59.0 )
0276             { theRing = 4; } // Ring 4
0277           else if ( theStripDet->surface().position().perp() > 81.0 && theStripDet->surface().position().perp() < 85.0 )
0278             { theRing = 6; } // Ring 6
0279           else
0280             { theRing = -1; } // probably not a Laser Hit!
0281 
0282         // on which disk is this module
0283           theTECWheel = tTopo->tecWheel(detid.rawId);
0284           break;
0285         }
0286       }
0287 
0288         // which beam belongs these digis to
0289       if ( thePart == "TIB" && theTIBLayer == 4 )
0290       {
0291         if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTIB) 
0292           && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTIB))          { theBeam = 0; } // beam 0 
0293 
0294         else if ( (theStripDet->surface().position().phi() > 1.29 - theSearchPhiTIB) 
0295           && (theStripDet->surface().position().phi() < 1.29 + theSearchPhiTIB))     { theBeam = 1; } // beam 1
0296 
0297         else if ( (theStripDet->surface().position().phi() > 1.85 - theSearchPhiTIB) 
0298           && (theStripDet->surface().position().phi() < 1.85 + theSearchPhiTIB))     { theBeam = 2; } // beam 2
0299 
0300         else if ( (theStripDet->surface().position().phi() > 2.75 - theSearchPhiTIB) 
0301           && (theStripDet->surface().position().phi() < 2.75 + theSearchPhiTIB))     { theBeam = 3; } // beam 3
0302 
0303         else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTIB) 
0304           && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTIB))    { theBeam = 4; } // beam 4
0305 
0306         else if ( (theStripDet->surface().position().phi() > -2.00 - theSearchPhiTIB) 
0307           && (theStripDet->surface().position().phi() < -2.00 + theSearchPhiTIB))    { theBeam = 5; } // beam 5
0308 
0309         else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTIB) 
0310           && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTIB))    { theBeam = 6; } // beam 6
0311 
0312         else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTIB) 
0313           && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTIB))    { theBeam = 7; } // beam 7
0314         else
0315           { theBeam = -1; } // probably not a Laser Hit!
0316       }
0317       else if ( thePart == "TOB" && theTOBLayer == 1 )
0318       {
0319         if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTOB) 
0320           && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTOB))          { theBeam = 0; } // beam 0 
0321 
0322         else if ( (theStripDet->surface().position().phi() > 1.29 - theSearchPhiTOB) 
0323           && (theStripDet->surface().position().phi() < 1.29 + theSearchPhiTOB))     { theBeam = 1; } // beam 1
0324 
0325         else if ( (theStripDet->surface().position().phi() > 1.85 - theSearchPhiTOB)
0326           && (theStripDet->surface().position().phi() < 1.85 + theSearchPhiTOB))     { theBeam = 2; } // beam 2
0327 
0328         else if ( (theStripDet->surface().position().phi() > 2.75 - theSearchPhiTOB)
0329           && (theStripDet->surface().position().phi() < 2.75 + theSearchPhiTOB))     { theBeam = 3; } // beam 3
0330 
0331         else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTOB)
0332           && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTOB))    { theBeam = 4; } // beam 4
0333 
0334         else if ( (theStripDet->surface().position().phi() > -2.00 - theSearchPhiTOB)
0335           && (theStripDet->surface().position().phi() < -2.00 + theSearchPhiTOB))    { theBeam = 5; } // beam 5
0336 
0337         else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTOB)
0338           && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTOB))    { theBeam = 6; } // beam 6
0339 
0340         else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTOB)
0341           && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTOB))    { theBeam = 7; } // beam 7
0342         else
0343           { theBeam = -1; } // probably not a Laser Hit!
0344       }
0345       else if ( thePart == "TEC+" || thePart == "TEC-" )
0346       {
0347         if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTEC)
0348           && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTEC))          { theBeam = 0; } // beam 0 
0349 
0350         else if ( (theStripDet->surface().position().phi() > 1.18 - theSearchPhiTEC)
0351           && (theStripDet->surface().position().phi() < 1.18 + theSearchPhiTEC))     { theBeam = 1; } // beam 1
0352 
0353         else if ( (theStripDet->surface().position().phi() > 1.96 - theSearchPhiTEC)
0354           && (theStripDet->surface().position().phi() < 1.96 + theSearchPhiTEC))     { theBeam = 2; } // beam 2
0355 
0356         else if ( (theStripDet->surface().position().phi() > 2.74 - theSearchPhiTEC)
0357           && (theStripDet->surface().position().phi() < 2.74 + theSearchPhiTEC))     { theBeam = 3; } // beam 3
0358 
0359         else if ( (theStripDet->surface().position().phi() > -2.74 - theSearchPhiTEC)
0360           && (theStripDet->surface().position().phi() < -2.74 + theSearchPhiTEC))    { theBeam = 4; } // beam 4
0361 
0362         else if ( (theStripDet->surface().position().phi() > -1.96 - theSearchPhiTEC)
0363           && (theStripDet->surface().position().phi() < -1.96 + theSearchPhiTEC))    { theBeam = 5; } // beam 5
0364 
0365         else if ( (theStripDet->surface().position().phi() > -1.18 - theSearchPhiTEC)
0366           && (theStripDet->surface().position().phi() < -1.18 + theSearchPhiTEC))    { theBeam = 6; } // beam 6
0367 
0368         else if ( (theStripDet->surface().position().phi() > -0.39 - theSearchPhiTEC)
0369           && (theStripDet->surface().position().phi() < -0.39 + theSearchPhiTEC))    { theBeam = 7; } // beam 7
0370 
0371         else if ( (theStripDet->surface().position().phi() > 1.28 - theSearchPhiTEC)
0372           && (theStripDet->surface().position().phi() < 1.28 + theSearchPhiTEC))     { theBeam = 21; } // beam 1 TEC2TEC
0373 
0374         else if ( (theStripDet->surface().position().phi() > 1.84 - theSearchPhiTEC)
0375           && (theStripDet->surface().position().phi() < 1.84 + theSearchPhiTEC))     { theBeam = 22; } // beam 2 TEC2TEC
0376 
0377         else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTEC)
0378           && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTEC))    { theBeam = 24; } // beam 4 TEC2TEC
0379 
0380         else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTEC)
0381           && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTEC))    { theBeam = 26; } // beam 6 TEC2TEC
0382 
0383         else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTEC)
0384           && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTEC))    { theBeam = 27; } // beam 7 TEC2TEC
0385         else 
0386           { theBeam = -1; } // probably not a Laser Hit!
0387       }
0388 
0389       // get the RecHits
0390     SiStripRecHit2DCollection::DetSet::const_iterator rechitRangeIteratorBegin = rechitRange.begin();
0391     SiStripRecHit2DCollection::DetSet::const_iterator rechitRangeIteratorEnd = rechitRange.end();
0392     SiStripRecHit2DCollection::DetSet::const_iterator iRecHit = rechitRangeIteratorBegin;
0393       // loop on the RecHits
0394     for (; iRecHit != rechitRangeIteratorEnd; iRecHit++)
0395     {
0396       SiStripRecHit2D const rechit = *iRecHit;
0397       theRecHitPositionsX->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).x());
0398       theRecHitPositionsY->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).y());
0399       theRecHitPositionsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z());
0400       theRecHitPositionsYvsX->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).x(),
0401         theStripDet->surface().toGlobal(rechit.localPosition()).y());
0402       theRecHitPositionsPhivsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z(),
0403         theStripDet->surface().toGlobal(rechit.localPosition()).phi());
0404       theRecHitPositionsRvsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z(),
0405         theStripDet->surface().toGlobal(rechit.localPosition()).perp());
0406 
0407       if (thePart == "TEC+" || thePart == "TEC-")   
0408       {             
0409         double r_ = sqrt(pow(theStripDet->surface().toGlobal(rechit.localPosition()).x(),2) + pow(theStripDet->surface().toGlobal(rechit.localPosition()).y(),2));
0410         fillLaserBeamPlots(r_,theTECWheel,thePart,theRing,theBeam);
0411       }
0412 
0413     }
0414   }
0415 
0416   for ( SiStripRecHit2DCollection::const_iterator det_iter = theStereoRecHitCollection->begin(), det_end = theStereoRecHitCollection->end();
0417         det_iter != det_end; ++det_iter) {
0418     SiStripRecHit2DCollection::DetSet rechitRange = *det_iter;
0419     DetId detid(rechitRange.detId());
0420       // get the DetUnit
0421     const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>(theTracker.idToDet(detid));
0422 
0423         // some variables we need later on in the program
0424       int theBeam     = 0;
0425       int theRing     = 0;
0426       std::string thePart  = "";
0427       int theTIBLayer = 0;
0428       int theTOBLayer = 0;
0429       int theTECWheel = 0;
0430       int theTOBStereoDet = 0;
0431 
0432       switch (detid.subdetId())
0433       {
0434         case StripSubdetector::TIB:
0435         {
0436           
0437           thePart = "TIB";
0438           theTIBLayer = tTopo->tibLayer(detid.rawId);
0439           break;
0440         }
0441         case StripSubdetector::TOB:
0442         {
0443           
0444           thePart = "TOB";
0445           theTOBLayer = tTopo->tobLayer(detid.rawId);
0446           theTOBStereoDet = tTopo->tobStereo(detid.rawId);
0447           break;
0448         }
0449         case StripSubdetector::TEC:
0450         {
0451           
0452 
0453         // is this module in TEC+ or TEC-?
0454           if (tTopo->tecSide(detid.rawId) == 1) { thePart = "TEC-"; }
0455           else if (tTopo->tecSide(detid.rawId) == 2) { thePart = "TEC+"; }
0456 
0457         // in which ring is this module?
0458           if ( theStripDet->surface().position().perp() > 55.0 && theStripDet->surface().position().perp() < 59.0 )
0459             { theRing = 4; } // Ring 4
0460           else if ( theStripDet->surface().position().perp() > 81.0 && theStripDet->surface().position().perp() < 85.0 )
0461             { theRing = 6; } // Ring 6
0462           else
0463             { theRing = -1; } // probably not a Laser Hit!
0464 
0465         // on which disk is this module
0466           theTECWheel = tTopo->tecWheel(detid.rawId);
0467           break;
0468         }
0469       }
0470 
0471         // which beam belongs these digis to
0472       if ( thePart == "TIB" && theTIBLayer == 4 )
0473       {
0474         if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTIB) 
0475           && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTIB))          { theBeam = 0; } // beam 0 
0476 
0477         else if ( (theStripDet->surface().position().phi() > 1.29 - theSearchPhiTIB) 
0478           && (theStripDet->surface().position().phi() < 1.29 + theSearchPhiTIB))     { theBeam = 1; } // beam 1
0479 
0480         else if ( (theStripDet->surface().position().phi() > 1.85 - theSearchPhiTIB) 
0481           && (theStripDet->surface().position().phi() < 1.85 + theSearchPhiTIB))     { theBeam = 2; } // beam 2
0482 
0483         else if ( (theStripDet->surface().position().phi() > 2.75 - theSearchPhiTIB) 
0484           && (theStripDet->surface().position().phi() < 2.75 + theSearchPhiTIB))     { theBeam = 3; } // beam 3
0485 
0486         else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTIB) 
0487           && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTIB))    { theBeam = 4; } // beam 4
0488 
0489         else if ( (theStripDet->surface().position().phi() > -2.00 - theSearchPhiTIB) 
0490           && (theStripDet->surface().position().phi() < -2.00 + theSearchPhiTIB))    { theBeam = 5; } // beam 5
0491 
0492         else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTIB) 
0493           && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTIB))    { theBeam = 6; } // beam 6
0494 
0495         else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTIB) 
0496           && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTIB))    { theBeam = 7; } // beam 7
0497         else
0498           { theBeam = -1; } // probably not a Laser Hit!
0499       }
0500       else if ( thePart == "TOB" && theTOBLayer == 1 )
0501       {
0502         if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTOB) 
0503           && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTOB))          { theBeam = 0; } // beam 0 
0504 
0505         else if ( (theStripDet->surface().position().phi() > 1.29 - theSearchPhiTOB) 
0506           && (theStripDet->surface().position().phi() < 1.29 + theSearchPhiTOB))     { theBeam = 1; } // beam 1
0507 
0508         else if ( (theStripDet->surface().position().phi() > 1.85 - theSearchPhiTOB)
0509           && (theStripDet->surface().position().phi() < 1.85 + theSearchPhiTOB))     { theBeam = 2; } // beam 2
0510 
0511         else if ( (theStripDet->surface().position().phi() > 2.75 - theSearchPhiTOB)
0512           && (theStripDet->surface().position().phi() < 2.75 + theSearchPhiTOB))     { theBeam = 3; } // beam 3
0513 
0514         else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTOB)
0515           && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTOB))    { theBeam = 4; } // beam 4
0516 
0517         else if ( (theStripDet->surface().position().phi() > -2.00 - theSearchPhiTOB)
0518           && (theStripDet->surface().position().phi() < -2.00 + theSearchPhiTOB))    { theBeam = 5; } // beam 5
0519 
0520         else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTOB)
0521           && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTOB))    { theBeam = 6; } // beam 6
0522 
0523         else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTOB)
0524           && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTOB))    { theBeam = 7; } // beam 7
0525         else
0526           { theBeam = -1; } // probably not a Laser Hit!
0527       }
0528       else if ( thePart == "TEC+" || thePart == "TEC-" )
0529       {
0530         if ( (theStripDet->surface().position().phi() > 0.39 - theSearchPhiTEC)
0531           && (theStripDet->surface().position().phi() < 0.39 + theSearchPhiTEC))          { theBeam = 0; } // beam 0 
0532 
0533         else if ( (theStripDet->surface().position().phi() > 1.18 - theSearchPhiTEC)
0534           && (theStripDet->surface().position().phi() < 1.18 + theSearchPhiTEC))     { theBeam = 1; } // beam 1
0535 
0536         else if ( (theStripDet->surface().position().phi() > 1.96 - theSearchPhiTEC)
0537           && (theStripDet->surface().position().phi() < 1.96 + theSearchPhiTEC))     { theBeam = 2; } // beam 2
0538 
0539         else if ( (theStripDet->surface().position().phi() > 2.74 - theSearchPhiTEC)
0540           && (theStripDet->surface().position().phi() < 2.74 + theSearchPhiTEC))     { theBeam = 3; } // beam 3
0541 
0542         else if ( (theStripDet->surface().position().phi() > -2.74 - theSearchPhiTEC)
0543           && (theStripDet->surface().position().phi() < -2.74 + theSearchPhiTEC))    { theBeam = 4; } // beam 4
0544 
0545         else if ( (theStripDet->surface().position().phi() > -1.96 - theSearchPhiTEC)
0546           && (theStripDet->surface().position().phi() < -1.96 + theSearchPhiTEC))    { theBeam = 5; } // beam 5
0547 
0548         else if ( (theStripDet->surface().position().phi() > -1.18 - theSearchPhiTEC)
0549           && (theStripDet->surface().position().phi() < -1.18 + theSearchPhiTEC))    { theBeam = 6; } // beam 6
0550 
0551         else if ( (theStripDet->surface().position().phi() > -0.39 - theSearchPhiTEC)
0552           && (theStripDet->surface().position().phi() < -0.39 + theSearchPhiTEC))    { theBeam = 7; } // beam 7
0553 
0554         else if ( (theStripDet->surface().position().phi() > 1.28 - theSearchPhiTEC)
0555           && (theStripDet->surface().position().phi() < 1.28 + theSearchPhiTEC))     { theBeam = 21; } // beam 1 TEC2TEC
0556 
0557         else if ( (theStripDet->surface().position().phi() > 1.84 - theSearchPhiTEC)
0558           && (theStripDet->surface().position().phi() < 1.84 + theSearchPhiTEC))     { theBeam = 22; } // beam 2 TEC2TEC
0559 
0560         else if ( (theStripDet->surface().position().phi() > -2.59 - theSearchPhiTEC)
0561           && (theStripDet->surface().position().phi() < -2.59 + theSearchPhiTEC))    { theBeam = 24; } // beam 4 TEC2TEC
0562 
0563         else if ( (theStripDet->surface().position().phi() > -1.10 - theSearchPhiTEC)
0564           && (theStripDet->surface().position().phi() < -1.10 + theSearchPhiTEC))    { theBeam = 26; } // beam 6 TEC2TEC
0565 
0566         else if ( (theStripDet->surface().position().phi() > -0.50 - theSearchPhiTEC)
0567           && (theStripDet->surface().position().phi() < -0.50 + theSearchPhiTEC))    { theBeam = 27; } // beam 7 TEC2TEC
0568         else 
0569           { theBeam = -1; } // probably not a Laser Hit!
0570       }
0571 
0572       // get the RecHits
0573     SiStripRecHit2DCollection::DetSet::const_iterator rechitRangeIteratorBegin = rechitRange.begin();
0574     SiStripRecHit2DCollection::DetSet::const_iterator rechitRangeIteratorEnd = rechitRange.end();
0575     SiStripRecHit2DCollection::DetSet::const_iterator iRecHit = rechitRangeIteratorBegin;
0576       // loop on the RecHits
0577     for (; iRecHit != rechitRangeIteratorEnd; iRecHit++)
0578     {
0579       SiStripRecHit2D const rechit = *iRecHit;
0580       theRecHitPositionsX->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).x());
0581       theRecHitPositionsY->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).y());
0582       theRecHitPositionsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z());
0583       theRecHitPositionsYvsX->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).x(),
0584         theStripDet->surface().toGlobal(rechit.localPosition()).y());
0585       theRecHitPositionsPhivsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z(),
0586         theStripDet->surface().toGlobal(rechit.localPosition()).phi());
0587       theRecHitPositionsRvsZ->Fill(theStripDet->surface().toGlobal(rechit.localPosition()).z(),
0588         theStripDet->surface().toGlobal(rechit.localPosition()).phi());
0589 
0590       if (thePart == "TEC+" || thePart == "TEC-")   
0591       {             
0592         double r_ = sqrt(pow(theStripDet->surface().toGlobal(rechit.localPosition()).x(),2) + pow(theStripDet->surface().toGlobal(rechit.localPosition()).y(),2));
0593         fillLaserBeamPlots(r_,theTECWheel,thePart,theRing,theBeam);
0594       }
0595     }
0596   }
0597 }