CMS 3D CMS Logo

DTSegment2DSLPhiQuality Class Reference

Basic analyzer class which accesses 2D DTSegments reconstructed with both SL Phi and plot resolution comparing reconstructed and simulated quantities. More...

#include <Validation/DTRecHits/plugins/DTSegment2DSLPhiQuality.h>

Inheritance diagram for DTSegment2DSLPhiQuality:

edm::EDAnalyzer

List of all members.

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 Perform the real analysis.
 DTSegment2DSLPhiQuality (const edm::ParameterSet &pset)
 Constructor.
void endJob ()
virtual ~DTSegment2DSLPhiQuality ()
 Destructor.

Private Attributes

bool debug
HEff2DHith2DHitEff_SuperPhi
HRes2DHith2DHitSuperPhi
std::string rootFileName
std::string segment4DLabel
double sigmaResAngle
double sigmaResPos
std::string simHitLabel
TFile * theFile


Detailed Description

Basic analyzer class which accesses 2D DTSegments reconstructed with both SL Phi and plot resolution comparing reconstructed and simulated quantities.

Date
2007/10/25 11:58:37
Revision
1.1
Author:
S. Bolognesi and G. Cerminara - INFN Torino

Definition at line 28 of file DTSegment2DSLPhiQuality.h.


Constructor & Destructor Documentation

DTSegment2DSLPhiQuality::DTSegment2DSLPhiQuality ( const edm::ParameterSet pset  ) 

Constructor.

Definition at line 38 of file DTSegment2DSLPhiQuality.cc.

References GenMuonPlsPt100GeV_cfg::cout, DTHitQualityUtils::debug, debug, lat::endl(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), h2DHitEff_SuperPhi, h2DHitSuperPhi, rootFileName, segment4DLabel, sigmaResAngle, sigmaResPos, simHitLabel, and theFile.

00038                                                                           {
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   if(debug)
00059     cout << "[DTSegment2DSLPhiQuality] Constructor called" << endl;
00060 
00061   // Book the histos
00062   h2DHitSuperPhi = new HRes2DHit ("SuperPhi");
00063   h2DHitEff_SuperPhi = new HEff2DHit ("SuperPhi");
00064 }

DTSegment2DSLPhiQuality::~DTSegment2DSLPhiQuality (  )  [virtual]

Destructor.

Definition at line 67 of file DTSegment2DSLPhiQuality.cc.

References GenMuonPlsPt100GeV_cfg::cout, debug, and lat::endl().

00067                                                  {
00068 
00069   if(debug)
00070     cout << "[DTSegment2DSLPhiQuality] Destructor called" << endl;
00071 }


Member Function Documentation

void DTSegment2DSLPhiQuality::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
) [virtual]

Perform the real analysis.

Implements edm::EDAnalyzer.

Definition at line 86 of file DTSegment2DSLPhiQuality.cc.

References GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), PV3DBase< T, PVType, FrameType >::eta(), HEff2DHit::Fill(), HRes2DHit::Fill(), DTHitQualityUtils::findMuSimSegment(), DTHitQualityUtils::findMuSimSegmentDirAndPos(), DTHitQualityUtils::findSegmentAlphaAndBeta(), edm::EventSetup::get(), h2DHitEff_SuperPhi, h2DHitSuperPhi, DTRecSegment2D::localDirection(), DTRecSegment2D::localDirectionError(), DTRecSegment2D::localPosition(), DTRecSegment2D::localPositionError(), DTHitQualityUtils::mapMuSimHitsPerWire(), DTHitQualityUtils::mapSimHitsPerWire(), PV3DBase< T, PVType, FrameType >::phi(), range, segment4DLabel, sigmaResAngle, sigmaResPos, simHitLabel, trackerHits::simHits, funct::sqrt(), theFile, GeomDet::toGlobal(), PV3DBase< T, PVType, FrameType >::x(), and LocalError::xx().

00086                                                                                       {
00087   if(debug)
00088     cout << "--- [DTSegment2DSLPhiQuality] Analysing Event: #Run: " << event.id().run()
00089          << " #Event: " << event.id().event() << endl;
00090   theFile->cd();
00091 
00092   // Get the DT Geometry
00093   ESHandle<DTGeometry> dtGeom;
00094   eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00095 
00096   // Get the SimHit collection from the event
00097   edm::Handle<PSimHitContainer> simHits;
00098   event.getByLabel(simHitLabel, "MuonDTHits", simHits); //FIXME: second string to be removed
00099 
00100   //Map simHits by chamber
00101   map<DTChamberId, PSimHitContainer > simHitsPerCh;
00102   for(PSimHitContainer::const_iterator simHit = simHits->begin();
00103       simHit != simHits->end(); simHit++){
00104     // Create the id of the chamber (the simHits in the DT known their wireId)
00105     DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
00106     // Fill the map
00107     simHitsPerCh[chamberId].push_back(*simHit);
00108   }
00109 
00110   // Get the 4D rechits from the event
00111   Handle<DTRecSegment4DCollection> segment4Ds;
00112   event.getByLabel(segment4DLabel, segment4Ds);
00113     
00114   // Loop over all chambers containing a segment
00115   DTRecSegment4DCollection::id_iterator chamberId;
00116   for (chamberId = segment4Ds->id_begin();
00117        chamberId != segment4Ds->id_end();
00118        ++chamberId){
00119     
00120     //------------------------- simHits ---------------------------//
00121     //Get simHits of each chamber
00122     PSimHitContainer simHits =  simHitsPerCh[(*chamberId)];
00123        
00124     // Map simhits per wire
00125     map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
00126     map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
00127     int nMuSimHit = muSimHitPerWire.size();
00128     if(nMuSimHit == 0 || nMuSimHit == 1) {
00129       if(debug && nMuSimHit == 1)
00130         cout << "[DTSegment2DSLPhiQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
00131       continue; // If no or only one mu SimHit is found skip this chamber
00132     } 
00133     if(debug)
00134       cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
00135   
00136     //Find outer and inner mu SimHit to build a segment
00137     pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire); 
00138 
00139     //Find direction and position of the sim Segment in Chamber RF
00140     pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
00141                                                                                                   (*chamberId),&(*dtGeom));
00142 
00143     LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
00144     LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
00145     const DTChamber* chamber = dtGeom->chamber(*chamberId);
00146     GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
00147 
00148     //Atan(x/z) angle and x position in SL RF
00149     float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
00150     float posSimSeg = simSegmLocalPos.x(); 
00151     //Position (in eta,phi coordinates) in lobal RF
00152     float etaSimSeg = simSegmGlobalPos.eta(); 
00153     float phiSimSeg = simSegmGlobalPos.phi();
00154 
00155     if(debug)
00156       cout<<"  Simulated segment:  local direction "<<simSegmLocalDir<<endl
00157           <<"                      local position  "<<simSegmLocalPos<<endl
00158           <<"                      angle           "<<angleSimSeg<<endl;
00159     
00160     //---------------------------- recHits --------------------------//
00161     // Get the range of rechit for the corresponding chamberId
00162     bool recHitFound = false;
00163     DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
00164     int nsegm = distance(range.first, range.second);
00165     if(debug)
00166       cout << "   Chamber: " << *chamberId << " has " << nsegm
00167            << " 4D segments" << endl;
00168 
00169     if (nsegm!=0) {
00170       // Find the best RecHit: look for the 4D RecHit with the phi angle closest
00171       // to that of segment made of SimHits. 
00172       // RecHits must have delta alpha and delta position within 5 sigma of
00173       // the residual distribution (we are looking for residuals of segments
00174       // usefull to the track fit) for efficency purpose
00175       const DTRecSegment2D* bestRecHit = 0;
00176       bool bestRecHitFound = false;
00177       double deltaAlpha = 99999;
00178 
00179       // Loop over the recHits of this chamberId
00180       for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
00181            segment4D!=range.second;
00182            ++segment4D){
00183         // Check the dimension
00184         if((*segment4D).dimension() != 4) {
00185           if(debug) cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
00186           continue;
00187         }
00188 
00189         //Get 2D superPhi segments from 4D segments
00190         const DTChamberRecSegment2D* phiSegment2D = (*segment4D).phiSegment();
00191         if((*phiSegment2D).dimension() != 2) {
00192           if(debug) cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
00193           abort();
00194         }
00195 
00196         // Segment Local Direction and position (in Chamber RF)
00197         LocalVector recSegDirection = (*phiSegment2D).localDirection();
00198 
00199         float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
00200         if(debug)
00201           cout << "  RecSegment direction: " << recSegDirection << endl
00202             << "             position : " <<  (*phiSegment2D).localPosition() << endl
00203             << "             alpha    : " << recSegAlpha << endl;
00204 
00205         if(fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
00206           deltaAlpha = fabs(recSegAlpha - angleSimSeg);
00207           bestRecHit = &(*phiSegment2D);
00208           bestRecHitFound = true;
00209         }
00210       }  // End of Loop over all 4D RecHits of this chambers
00211 
00212       if(bestRecHitFound) {
00213         // Best rechit direction and position in Chamber RF
00214         LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
00215         LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
00216 
00217         LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
00218         LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
00219 
00220         float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
00221         if(fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle &&
00222            fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos) {
00223           recHitFound = true;
00224         }
00225 
00226         // Fill Residual histos
00227         h2DHitSuperPhi->Fill(angleSimSeg,
00228                             angleBestRHit,
00229                             posSimSeg,
00230                             bestRecHitLocalPos.x(),
00231                             etaSimSeg,
00232                             phiSimSeg,
00233                             sqrt(bestRecHitLocalPosErr.xx()),
00234                             sqrt(bestRecHitLocalDirErr.xx())
00235                            );
00236       }
00237     } //end of if(nsegm!=0)
00238 
00239       // Fill Efficiency plot
00240     h2DHitEff_SuperPhi->Fill(etaSimSeg,
00241                             phiSimSeg,
00242                             posSimSeg,
00243                             angleSimSeg,
00244                             recHitFound);
00245   } // End of loop over chambers
00246 }

void DTSegment2DSLPhiQuality::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 73 of file DTSegment2DSLPhiQuality.cc.

References HEff2DHit::ComputeEfficiency(), h2DHitEff_SuperPhi, h2DHitSuperPhi, theFile, HRes2DHit::Write(), and HEff2DHit::Write().

00073                                      {
00074   // Write the histos to file
00075   theFile->cd();
00076 
00077   h2DHitSuperPhi->Write();
00078 
00079   h2DHitEff_SuperPhi->ComputeEfficiency();
00080   h2DHitEff_SuperPhi->Write();
00081 
00082   theFile->Close();
00083 } 


Member Data Documentation

bool DTSegment2DSLPhiQuality::debug [private]

Definition at line 50 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), DTSegment2DSLPhiQuality(), and ~DTSegment2DSLPhiQuality().

HEff2DHit* DTSegment2DSLPhiQuality::h2DHitEff_SuperPhi [private]

Definition at line 62 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), DTSegment2DSLPhiQuality(), and endJob().

HRes2DHit* DTSegment2DSLPhiQuality::h2DHitSuperPhi [private]

Definition at line 61 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), DTSegment2DSLPhiQuality(), and endJob().

std::string DTSegment2DSLPhiQuality::rootFileName [private]

Definition at line 52 of file DTSegment2DSLPhiQuality.h.

Referenced by DTSegment2DSLPhiQuality().

std::string DTSegment2DSLPhiQuality::segment4DLabel [private]

Definition at line 55 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), and DTSegment2DSLPhiQuality().

double DTSegment2DSLPhiQuality::sigmaResAngle [private]

Definition at line 59 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), and DTSegment2DSLPhiQuality().

double DTSegment2DSLPhiQuality::sigmaResPos [private]

Definition at line 57 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), and DTSegment2DSLPhiQuality().

std::string DTSegment2DSLPhiQuality::simHitLabel [private]

Definition at line 54 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), and DTSegment2DSLPhiQuality().

TFile* DTSegment2DSLPhiQuality::theFile [private]

Definition at line 48 of file DTSegment2DSLPhiQuality.h.

Referenced by analyze(), DTSegment2DSLPhiQuality(), and endJob().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:19:04 2009 for CMSSW by  doxygen 1.5.4