CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/Validation/DTRecHits/plugins/DTSegment2DSLPhiQuality.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2008/10/21 10:52:20 $
00005  *  $Revision: 1.2 $
00006  *  \author S. Bolognesi and G. Cerminara - INFN Torino
00007  */
00008 
00009 #include "DTSegment2DSLPhiQuality.h"
00010 #include "Validation/DTRecHits/interface/DTHitQualityUtils.h"
00011 
00012 #include "DataFormats/MuonDetId/interface/DTWireId.h"
00013 #include "Geometry/DTGeometry/interface/DTChamber.h"
00014 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00015 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00016 
00017 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00018 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
00019 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00020 
00021 #include "FWCore/Framework/interface/MakerMacros.h"
00022 #include "FWCore/Framework/interface/Frameworkfwd.h"
00023 #include "FWCore/Framework/interface/Event.h"
00024 #include "FWCore/Framework/interface/ESHandle.h"
00025 
00026 #include "Histograms.h"
00027 
00028 #include "TFile.h"
00029 
00030 #include <iostream>
00031 #include <map>
00032 
00033 
00034 using namespace std;
00035 using namespace edm;
00036 
00037 // Constructor
00038 DTSegment2DSLPhiQuality::DTSegment2DSLPhiQuality(const ParameterSet& pset)  {
00039    // Get the debug parameter for verbose output
00040   debug = pset.getUntrackedParameter<bool>("debug");
00041   DTHitQualityUtils::debug = debug;
00042  
00043   rootFileName = pset.getUntrackedParameter<string>("rootFileName");
00044   // the name of the simhit collection
00045   simHitLabel = pset.getUntrackedParameter<string>("simHitLabel", "SimG4Object");
00046   // the name of the 4D rec hit collection
00047   segment4DLabel = pset.getUntrackedParameter<string>("segment4DLabel");
00048 
00049   //sigma resolution on position
00050   sigmaResPos = pset.getParameter<double>("sigmaResPos");
00051   //sigma resolution on angle
00052   sigmaResAngle = pset.getParameter<double>("sigmaResAngle");
00053 
00054   // Create the root file
00055   theFile = new TFile(rootFileName.c_str(), "RECREATE");
00056   theFile->cd();
00057 
00058   // Book the histos
00059   h2DHitSuperPhi = new HRes2DHit ("SuperPhi");
00060   h2DHitEff_SuperPhi = new HEff2DHit ("SuperPhi");
00061 }
00062 
00063 // Destructor
00064 DTSegment2DSLPhiQuality::~DTSegment2DSLPhiQuality(){
00065 }
00066 
00067 void DTSegment2DSLPhiQuality::endJob() {
00068   // Write the histos to file
00069   theFile->cd();
00070 
00071   h2DHitSuperPhi->Write();
00072 
00073   h2DHitEff_SuperPhi->ComputeEfficiency();
00074   h2DHitEff_SuperPhi->Write();
00075 
00076   theFile->Close();
00077 } 
00078 
00079 // The real analysis
00080 void DTSegment2DSLPhiQuality::analyze(const Event & event, const EventSetup& eventSetup){
00081   theFile->cd();
00082 
00083   // Get the DT Geometry
00084   ESHandle<DTGeometry> dtGeom;
00085   eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00086 
00087   // Get the SimHit collection from the event
00088   edm::Handle<PSimHitContainer> simHits;
00089   event.getByLabel(simHitLabel, "MuonDTHits", simHits); //FIXME: second string to be removed
00090 
00091   //Map simHits by chamber
00092   map<DTChamberId, PSimHitContainer > simHitsPerCh;
00093   for(PSimHitContainer::const_iterator simHit = simHits->begin();
00094       simHit != simHits->end(); simHit++){
00095     // Create the id of the chamber (the simHits in the DT known their wireId)
00096     DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
00097     // Fill the map
00098     simHitsPerCh[chamberId].push_back(*simHit);
00099   }
00100 
00101   // Get the 4D rechits from the event
00102   Handle<DTRecSegment4DCollection> segment4Ds;
00103   event.getByLabel(segment4DLabel, segment4Ds);
00104     
00105   // Loop over all chambers containing a segment
00106   DTRecSegment4DCollection::id_iterator chamberId;
00107   for (chamberId = segment4Ds->id_begin();
00108        chamberId != segment4Ds->id_end();
00109        ++chamberId){
00110     
00111     //------------------------- simHits ---------------------------//
00112     //Get simHits of each chamber
00113     PSimHitContainer simHits =  simHitsPerCh[(*chamberId)];
00114        
00115     // Map simhits per wire
00116     map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
00117     map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
00118     int nMuSimHit = muSimHitPerWire.size();
00119     if(nMuSimHit == 0 || nMuSimHit == 1) {
00120       if(debug && nMuSimHit == 1)
00121         cout << "[DTSegment2DSLPhiQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
00122       continue; // If no or only one mu SimHit is found skip this chamber
00123     } 
00124     if(debug)
00125       cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
00126   
00127     //Find outer and inner mu SimHit to build a segment
00128     pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire); 
00129 
00130     //Find direction and position of the sim Segment in Chamber RF
00131     pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
00132                                                                                                   (*chamberId),&(*dtGeom));
00133 
00134     LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
00135     LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
00136     const DTChamber* chamber = dtGeom->chamber(*chamberId);
00137     GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
00138 
00139     //Atan(x/z) angle and x position in SL RF
00140     float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
00141     float posSimSeg = simSegmLocalPos.x(); 
00142     //Position (in eta,phi coordinates) in lobal RF
00143     float etaSimSeg = simSegmGlobalPos.eta(); 
00144     float phiSimSeg = simSegmGlobalPos.phi();
00145 
00146     if(debug)
00147       cout<<"  Simulated segment:  local direction "<<simSegmLocalDir<<endl
00148           <<"                      local position  "<<simSegmLocalPos<<endl
00149           <<"                      angle           "<<angleSimSeg<<endl;
00150     
00151     //---------------------------- recHits --------------------------//
00152     // Get the range of rechit for the corresponding chamberId
00153     bool recHitFound = false;
00154     DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
00155     int nsegm = distance(range.first, range.second);
00156     if(debug)
00157       cout << "   Chamber: " << *chamberId << " has " << nsegm
00158            << " 4D segments" << endl;
00159 
00160     if (nsegm!=0) {
00161       // Find the best RecHit: look for the 4D RecHit with the phi angle closest
00162       // to that of segment made of SimHits. 
00163       // RecHits must have delta alpha and delta position within 5 sigma of
00164       // the residual distribution (we are looking for residuals of segments
00165       // usefull to the track fit) for efficency purpose
00166       const DTRecSegment2D* bestRecHit = 0;
00167       bool bestRecHitFound = false;
00168       double deltaAlpha = 99999;
00169 
00170       // Loop over the recHits of this chamberId
00171       for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
00172            segment4D!=range.second;
00173            ++segment4D){
00174         // Check the dimension
00175         if((*segment4D).dimension() != 4) {
00176           if(debug) cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
00177           continue;
00178         }
00179 
00180         //Get 2D superPhi segments from 4D segments
00181         const DTChamberRecSegment2D* phiSegment2D = (*segment4D).phiSegment();
00182         if((*phiSegment2D).dimension() != 2) {
00183           if(debug) cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
00184           abort();
00185         }
00186 
00187         // Segment Local Direction and position (in Chamber RF)
00188         LocalVector recSegDirection = (*phiSegment2D).localDirection();
00189 
00190         float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
00191         if(debug)
00192           cout << "  RecSegment direction: " << recSegDirection << endl
00193             << "             position : " <<  (*phiSegment2D).localPosition() << endl
00194             << "             alpha    : " << recSegAlpha << endl;
00195 
00196         if(fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
00197           deltaAlpha = fabs(recSegAlpha - angleSimSeg);
00198           bestRecHit = &(*phiSegment2D);
00199           bestRecHitFound = true;
00200         }
00201       }  // End of Loop over all 4D RecHits of this chambers
00202 
00203       if(bestRecHitFound) {
00204         // Best rechit direction and position in Chamber RF
00205         LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
00206         LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
00207 
00208         LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
00209         LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
00210 
00211         float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
00212         if(fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle &&
00213            fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos) {
00214           recHitFound = true;
00215         }
00216 
00217         // Fill Residual histos
00218         h2DHitSuperPhi->Fill(angleSimSeg,
00219                             angleBestRHit,
00220                             posSimSeg,
00221                             bestRecHitLocalPos.x(),
00222                             etaSimSeg,
00223                             phiSimSeg,
00224                             sqrt(bestRecHitLocalPosErr.xx()),
00225                             sqrt(bestRecHitLocalDirErr.xx())
00226                            );
00227       }
00228     } //end of if(nsegm!=0)
00229 
00230       // Fill Efficiency plot
00231     h2DHitEff_SuperPhi->Fill(etaSimSeg,
00232                             phiSimSeg,
00233                             posSimSeg,
00234                             angleSimSeg,
00235                             recHitFound);
00236   } // End of loop over chambers
00237 }
00238 
00239 
00240 
00241 
00242 
00243