#include <DTChamberEfficiency.h>
This class provides the histograms for the calculation of the efficiency of muons reconstruction in the DTs. It is applicable both in presence or absence of a magnetic field. Histos are 2D Sector vs Chamber plots for each wheel
Modification:
Definition at line 54 of file DTChamberEfficiency.h.
DTChamberEfficiency::DTChamberEfficiency | ( | const edm::ParameterSet & | pset | ) |
Definition at line 57 of file DTChamberEfficiency.cc.
References Chi2MeasurementEstimatorESProducer_cfi::Chi2MeasurementEstimator, debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), MuonServiceProxy_cff::MuonServiceProxy, cmsCodeRules::cppFunctionSkipper::operator, and createJobs::theNSigma.
{ // Get the debug parameter for verbose output debug = pSet.getUntrackedParameter<bool>("debug",false); LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: constructor called"; // Get the DQM needed services theDbe = Service<DQMStore>().operator->(); theDbe->setCurrentFolder("DT/05-ChamberEff/Task"); // service parameters ParameterSet serviceParameters = pSet.getParameter<ParameterSet>("ServiceParameters"); theService = new MuonServiceProxy(serviceParameters); theTracksLabel = pSet.getParameter<InputTag>("TrackCollection"); theMaxChi2 = static_cast<unsigned int>(pSet.getParameter<double>("theMaxChi2")); theNSigma = pSet.getParameter<double>("theNSigma"); theMinNrec = static_cast<int>(pSet.getParameter<double>("theMinNrec")); labelRPCRecHits = pSet.getParameter<InputTag>("theRPCRecHits"); thedt4DSegments = pSet.getParameter<InputTag>("dt4DSegments"); thecscSegments = pSet.getParameter<InputTag>("cscSegments"); theMeasurementExtractor = new MuonDetLayerMeasurements(thedt4DSegments,thecscSegments, labelRPCRecHits,true,false,false); theNavigationType = pSet.getParameter<string>("NavigationType"); theEstimator = new Chi2MeasurementEstimator(theMaxChi2,theNSigma); }
DTChamberEfficiency::~DTChamberEfficiency | ( | ) |
Definition at line 93 of file DTChamberEfficiency.cc.
References LogTrace.
{ LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: destructor called"; // free memory delete theService; delete theMeasurementExtractor; delete theEstimator; }
void DTChamberEfficiency::analyze | ( | const edm::Event & | event, |
const edm::EventSetup & | eventSetup | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 164 of file DTChamberEfficiency.cc.
References alongMomentum, debug, DetId::det(), GeomDetEnumerators::DT, DTChamberId, event(), HcalObjRepresent::Fill(), TrajectoryStateOnSurface::freeState(), mergeVDriftHistosByStation::histos, i, TrajectoryStateOnSurface::isValid(), LogTrace, autoMagneticFieldProducer_cfi::magfield, DetId::Muon, oppositeToMomentum, LargeD0_PixelPairStep_cff::propagator, DetId::rawId(), reco::TransientTrack::recHitsSize(), DTChamberId::sector(), DTChamberId::station(), DetId::subdetId(), TrajectoryStateOnSurface::surface(), testEve_cfg::tracks, and DTChamberId::wheel().
{ LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "--- [DTChamberEfficiency] Event analysed #Run: " << event.id().run() << " #Event: " << event.id().event() << endl; theService->update(eventSetup); theMeasurementExtractor->setEvent(event); // set navigation school NavigationSetter setter(*theService->muonNavigationSchool()); //Read tracks from event Handle<reco::TrackCollection> tracks; event.getByLabel(theTracksLabel, tracks); if(tracks.isValid()) { // check the validity of the collection //loop over the muons for(reco::TrackCollection::const_iterator track = tracks->begin(); track!=tracks->end(); ++track) { reco::TransientTrack trans_track(*track,magfield.product(),theTrackingGeometry); const int recHitsize = (int)trans_track.recHitsSize(); if(recHitsize < theMinNrec) continue; // printout the DT rechits used by the track if (debug) { LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "--- New track" << endl; set<DTChamberId> chAlrUsed; for(trackingRecHit_iterator rHit = trans_track.recHitsBegin(); rHit != trans_track.recHitsEnd(); ++rHit) { DetId rHitid = (*rHit)->geographicalId(); if(!(rHitid.det() == DetId::Muon && rHitid.subdetId() == MuonSubdetId::DT)) continue; DTChamberId wId(rHitid.rawId()); if(chAlrUsed.find(wId) != chAlrUsed.end()) continue; chAlrUsed.insert(wId); LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << " " << wId << endl; } } // Get the layer on which the seed relies DetId id = trans_track.track().innerDetId(); const DetLayer *initialLayer = theService->detLayerGeometry()->idToLayer(id); TrajectoryStateOnSurface init_fs = trans_track.innermostMeasurementState(); FreeTrajectoryState *init_fs_free = init_fs.freeState(); //get the list of compatible layers vector<const DetLayer*> layer_list = compatibleLayers(initialLayer,*init_fs_free,alongMomentum); vector<const DetLayer*> layer_list_2 = compatibleLayers(initialLayer,*init_fs_free,oppositeToMomentum); layer_list.insert(layer_list.end(),layer_list_2.begin(),layer_list_2.end()); set<DTChamberId> alreadyCheckedCh; //loop over the list of compatible layers for(int i=0;i< (int)layer_list.size();i++) { //propagate the track to the i-th layer TrajectoryStateOnSurface tsos = propagator()->propagate(init_fs,layer_list.at(i)->surface()); if(!tsos.isValid()) continue; //determine the chambers kinematically compatible with the track on the i-th layer vector<DetWithState> dss = layer_list.at(i)->compatibleDets(tsos, *propagator(), *theEstimator); if(dss.size() == 0) continue; // get the first det (it's the most compatible) const DetWithState detWithState = dss.front(); const DetId idDetLay = detWithState.first->geographicalId(); // check if this is a DT and the track has the needed quality if(!chamberSelection(idDetLay,trans_track)) continue; DTChamberId DTid = (DTChamberId) idDetLay; // check if the chamber has already been counted if(alreadyCheckedCh.find(DTid) != alreadyCheckedCh.end()) continue; alreadyCheckedCh.insert(DTid); // get the compatible measurements MeasurementContainer detMeasurements_initial = theMeasurementExtractor->measurements(layer_list.at(i), detWithState.first, detWithState.second, *theEstimator, event); LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << " chamber: " << DTid << " has: " << detMeasurements_initial.size() << " comp. meas." << endl; //we want to be more picky about the quality of the segments: //exclude the segments with less than 12 hits MeasurementContainer detMeasurements = segQualityCut(detMeasurements_initial); // get the histos for this chamber vector<MonitorElement *> histos = histosPerW[DTid.wheel()+2]; // fill them if (detMeasurements_initial.size() != 0) histos[0]->Fill(DTid.sector(),DTid.station(),1.); if (detMeasurements.size() != 0) histos[1]->Fill(DTid.sector(),DTid.station(),1.); histos[2]->Fill(DTid.sector(),DTid.station(),1.); } } } else { LogInfo("DTDQM|DTMonitorModule|DTChamberEfficiency") << "[DTChamberEfficiency] Collection: " << theTracksLabel << " is not valid!" << endl; } return; }
void DTChamberEfficiency::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 104 of file DTChamberEfficiency.cc.
References bookHistos(), and LogTrace.
{ LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: beginOfJob"; bookHistos(); return; }
void DTChamberEfficiency::beginRun | ( | const edm::Run & | run, |
const edm::EventSetup & | setup | ||
) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 114 of file DTChamberEfficiency.cc.
References edm::EventSetup::get(), and autoMagneticFieldProducer_cfi::magfield.
{ // Get the DT Geometry setup.get<MuonGeometryRecord>().get(dtGeom); setup.get<IdealMagneticFieldRecord>().get(magfield); setup.get<GlobalTrackingGeometryRecord>().get(theTrackingGeometry); return; }
void DTChamberEfficiency::bookHistos | ( | ) | [private] |
Definition at line 134 of file DTChamberEfficiency.cc.
References mergeVDriftHistosByStation::histos, and LogTrace.
{ LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: booking histos"; // Create the monitor elements theDbe->setCurrentFolder("DT/05-ChamberEff/Task"); for(int wheel=-2;wheel<=2;wheel++){ vector<MonitorElement *> histos; stringstream wheel_str; wheel_str << wheel; histos.push_back(theDbe->book2D("hCountSectVsChamb_All_W"+ wheel_str.str(), "Countings for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.)); histos.push_back(theDbe->book2D("hCountSectVsChamb_Qual_W"+ wheel_str.str(), "Countings for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.)); histos.push_back(theDbe->book2D("hExtrapSectVsChamb_W"+ wheel_str.str(), "Extrapolations for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.)); histosPerW.push_back(histos); } return; }
bool DTChamberEfficiency::chamberSelection | ( | const DetId & | idDetLay, |
reco::TransientTrack & | trans_track | ||
) | const [private] |
Definition at line 276 of file DTChamberEfficiency.cc.
References DetId::det(), GeomDetEnumerators::DT, DetId::Muon, reco::TransientTrack::recHit(), reco::TransientTrack::recHitsSize(), and DetId::subdetId().
{ //check that we have a muon and that is a DT detector if(!(idDetLay.det() == DetId::Muon && idDetLay.subdetId() == MuonSubdetId::DT)) return false; if(trans_track.recHitsSize() == 2) if(trans_track.recHit(0)->geographicalId() == idDetLay || trans_track.recHit(1)->geographicalId() == idDetLay) return false; return true; }
vector< const DetLayer * > DTChamberEfficiency::compatibleLayers | ( | const DetLayer * | initialLayer, |
const FreeTrajectoryState & | fts, | ||
PropagationDirection | propDir | ||
) | [private] |
Definition at line 318 of file DTChamberEfficiency.cc.
References DirectMuonNavigation::compatibleLayers(), and DetLayer::compatibleLayers().
{ vector<const DetLayer*> detLayers; if(theNavigationType == "Standard"){ // ask for compatible layers detLayers = initialLayer->compatibleLayers(fts,propDir); // I have to fit by hand the first layer until the seedTSOS is defined on the first rechit layer // In fact the first layer is not returned by initialLayer->compatibleLayers. detLayers.insert(detLayers.begin(),initialLayer); } else if (theNavigationType == "Direct"){ DirectMuonNavigation navigation(ESHandle<MuonDetLayerGeometry>(&*theService->detLayerGeometry())); detLayers = navigation.compatibleLayers(fts,propDir); } else LogError("DTDQM|DTMonitorModule|DTChamberEfficiency") << "No Properly Navigation Selected!!"<<endl; return detLayers; }
void DTChamberEfficiency::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 125 of file DTChamberEfficiency.cc.
References LogTrace.
{ LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: endOfJob"; return; }
ESHandle< Propagator > DTChamberEfficiency::propagator | ( | void | ) | const [inline, private] |
Definition at line 344 of file DTChamberEfficiency.cc.
{ return theService->propagator("SteppingHelixPropagatorAny"); }
MeasurementContainer DTChamberEfficiency::segQualityCut | ( | const MeasurementContainer | seg_list | ) | const [private] |
Definition at line 290 of file DTChamberEfficiency.cc.
References DTChamberId, query::result, and DTChamberId::station().
{ MeasurementContainer result; for(MeasurementContainer::const_iterator mescont_Itr = seg_list.begin(); mescont_Itr != seg_list.end(); ++mescont_Itr){ //get the rechits of the segment TransientTrackingRecHit::ConstRecHitContainer recHit_list = mescont_Itr->recHit()->transientHits(); //loop over the rechits and get the number of hits int nhit_seg(0); for(TransientTrackingRecHit::ConstRecHitContainer::const_iterator recList_Itr = recHit_list.begin(); recList_Itr != recHit_list.end(); ++recList_Itr){ nhit_seg += (int)(*recList_Itr)->transientHits().size(); } DTChamberId tmpId = (DTChamberId)mescont_Itr->recHit()->hit()->geographicalId(); if(tmpId.station() < 4 && nhit_seg >= 12) result.push_back(*mescont_Itr); if(tmpId.station() == 4 && nhit_seg >= 8) result.push_back(*mescont_Itr); } return result; }
bool DTChamberEfficiency::debug [private] |
Definition at line 83 of file DTChamberEfficiency.h.
edm::ESHandle<DTGeometry> DTChamberEfficiency::dtGeom [private] |
Definition at line 97 of file DTChamberEfficiency.h.
std::vector<std::vector<MonitorElement*> > DTChamberEfficiency::histosPerW [private] |
Definition at line 110 of file DTChamberEfficiency.h.
Definition at line 87 of file DTChamberEfficiency.h.
Definition at line 105 of file DTChamberEfficiency.h.
Definition at line 89 of file DTChamberEfficiency.h.
DQMStore* DTChamberEfficiency::theDbe [private] |
Definition at line 99 of file DTChamberEfficiency.h.
Definition at line 88 of file DTChamberEfficiency.h.
Definition at line 103 of file DTChamberEfficiency.h.
double DTChamberEfficiency::theMaxChi2 [private] |
Definition at line 91 of file DTChamberEfficiency.h.
Definition at line 102 of file DTChamberEfficiency.h.
int DTChamberEfficiency::theMinNrec [private] |
Definition at line 93 of file DTChamberEfficiency.h.
std::string DTChamberEfficiency::theNavigationType [private] |
Definition at line 95 of file DTChamberEfficiency.h.
double DTChamberEfficiency::theNSigma [private] |
Definition at line 92 of file DTChamberEfficiency.h.
MuonServiceProxy* DTChamberEfficiency::theService [private] |
Definition at line 101 of file DTChamberEfficiency.h.
Definition at line 106 of file DTChamberEfficiency.h.
Definition at line 85 of file DTChamberEfficiency.h.