Go to the documentation of this file.00001
00002
00003
00004
00005
00006
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
00038 DTSegment2DSLPhiQuality::DTSegment2DSLPhiQuality(const ParameterSet& pset) {
00039
00040 debug = pset.getUntrackedParameter<bool>("debug");
00041 DTHitQualityUtils::debug = debug;
00042
00043 rootFileName = pset.getUntrackedParameter<string>("rootFileName");
00044
00045 simHitLabel = pset.getUntrackedParameter<string>("simHitLabel", "SimG4Object");
00046
00047 segment4DLabel = pset.getUntrackedParameter<string>("segment4DLabel");
00048
00049
00050 sigmaResPos = pset.getParameter<double>("sigmaResPos");
00051
00052 sigmaResAngle = pset.getParameter<double>("sigmaResAngle");
00053
00054
00055 theFile = new TFile(rootFileName.c_str(), "RECREATE");
00056 theFile->cd();
00057
00058
00059 h2DHitSuperPhi = new HRes2DHit ("SuperPhi");
00060 h2DHitEff_SuperPhi = new HEff2DHit ("SuperPhi");
00061 }
00062
00063
00064 DTSegment2DSLPhiQuality::~DTSegment2DSLPhiQuality(){
00065 }
00066
00067 void DTSegment2DSLPhiQuality::endJob() {
00068
00069 theFile->cd();
00070
00071 h2DHitSuperPhi->Write();
00072
00073 h2DHitEff_SuperPhi->ComputeEfficiency();
00074 h2DHitEff_SuperPhi->Write();
00075
00076 theFile->Close();
00077 }
00078
00079
00080 void DTSegment2DSLPhiQuality::analyze(const Event & event, const EventSetup& eventSetup){
00081 theFile->cd();
00082
00083
00084 ESHandle<DTGeometry> dtGeom;
00085 eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00086
00087
00088 edm::Handle<PSimHitContainer> simHits;
00089 event.getByLabel(simHitLabel, "MuonDTHits", simHits);
00090
00091
00092 map<DTChamberId, PSimHitContainer > simHitsPerCh;
00093 for(PSimHitContainer::const_iterator simHit = simHits->begin();
00094 simHit != simHits->end(); simHit++){
00095
00096 DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
00097
00098 simHitsPerCh[chamberId].push_back(*simHit);
00099 }
00100
00101
00102 Handle<DTRecSegment4DCollection> segment4Ds;
00103 event.getByLabel(segment4DLabel, segment4Ds);
00104
00105
00106 DTRecSegment4DCollection::id_iterator chamberId;
00107 for (chamberId = segment4Ds->id_begin();
00108 chamberId != segment4Ds->id_end();
00109 ++chamberId){
00110
00111
00112
00113 PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
00114
00115
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;
00123 }
00124 if(debug)
00125 cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
00126
00127
00128 pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
00129
00130
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
00140 float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
00141 float posSimSeg = simSegmLocalPos.x();
00142
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
00152
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
00162
00163
00164
00165
00166 const DTRecSegment2D* bestRecHit = 0;
00167 bool bestRecHitFound = false;
00168 double deltaAlpha = 99999;
00169
00170
00171 for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
00172 segment4D!=range.second;
00173 ++segment4D){
00174
00175 if((*segment4D).dimension() != 4) {
00176 if(debug) cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
00177 continue;
00178 }
00179
00180
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
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 }
00202
00203 if(bestRecHitFound) {
00204
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
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 }
00229
00230
00231 h2DHitEff_SuperPhi->Fill(etaSimSeg,
00232 phiSimSeg,
00233 posSimSeg,
00234 angleSimSeg,
00235 recHitFound);
00236 }
00237 }
00238
00239
00240
00241
00242
00243