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 if(debug)
00059 cout << "[DTSegment2DSLPhiQuality] Constructor called" << endl;
00060
00061
00062 h2DHitSuperPhi = new HRes2DHit ("SuperPhi");
00063 h2DHitEff_SuperPhi = new HEff2DHit ("SuperPhi");
00064 }
00065
00066
00067 DTSegment2DSLPhiQuality::~DTSegment2DSLPhiQuality(){
00068
00069 if(debug)
00070 cout << "[DTSegment2DSLPhiQuality] Destructor called" << endl;
00071 }
00072
00073 void DTSegment2DSLPhiQuality::endJob() {
00074
00075 theFile->cd();
00076
00077 h2DHitSuperPhi->Write();
00078
00079 h2DHitEff_SuperPhi->ComputeEfficiency();
00080 h2DHitEff_SuperPhi->Write();
00081
00082 theFile->Close();
00083 }
00084
00085
00086 void DTSegment2DSLPhiQuality::analyze(const Event & event, const EventSetup& eventSetup){
00087 if(debug)
00088 cout << "--- [DTSegment2DSLPhiQuality] Analysing Event: #Run: " << event.id().run()
00089 << " #Event: " << event.id().event() << endl;
00090 theFile->cd();
00091
00092
00093 ESHandle<DTGeometry> dtGeom;
00094 eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00095
00096
00097 edm::Handle<PSimHitContainer> simHits;
00098 event.getByLabel(simHitLabel, "MuonDTHits", simHits);
00099
00100
00101 map<DTChamberId, PSimHitContainer > simHitsPerCh;
00102 for(PSimHitContainer::const_iterator simHit = simHits->begin();
00103 simHit != simHits->end(); simHit++){
00104
00105 DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
00106
00107 simHitsPerCh[chamberId].push_back(*simHit);
00108 }
00109
00110
00111 Handle<DTRecSegment4DCollection> segment4Ds;
00112 event.getByLabel(segment4DLabel, segment4Ds);
00113
00114
00115 DTRecSegment4DCollection::id_iterator chamberId;
00116 for (chamberId = segment4Ds->id_begin();
00117 chamberId != segment4Ds->id_end();
00118 ++chamberId){
00119
00120
00121
00122 PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
00123
00124
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;
00132 }
00133 if(debug)
00134 cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
00135
00136
00137 pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
00138
00139
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
00149 float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
00150 float posSimSeg = simSegmLocalPos.x();
00151
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
00161
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
00171
00172
00173
00174
00175 const DTRecSegment2D* bestRecHit = 0;
00176 bool bestRecHitFound = false;
00177 double deltaAlpha = 99999;
00178
00179
00180 for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
00181 segment4D!=range.second;
00182 ++segment4D){
00183
00184 if((*segment4D).dimension() != 4) {
00185 if(debug) cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
00186 continue;
00187 }
00188
00189
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
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 }
00211
00212 if(bestRecHitFound) {
00213
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
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 }
00238
00239
00240 h2DHitEff_SuperPhi->Fill(etaSimSeg,
00241 phiSimSeg,
00242 posSimSeg,
00243 angleSimSeg,
00244 recHitFound);
00245 }
00246 }
00247
00248
00249
00250
00251
00252