CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTSegment2DSLPhiQuality.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author S. Bolognesi and G. Cerminara - INFN Torino
5  */
6 
9 
14 
18 
25 
26 #include "Histograms.h"
27 #include "TStyle.h"
28 #include "TFile.h"
29 
30 #include <iostream>
31 #include <map>
32 //#include "utils.C"
33 
34 using namespace std;
35 using namespace edm;
36 //TStyle * mystyle;
37 
38 
39 // Constructor
41  // Get the debug parameter for verbose output
42  debug = pset.getUntrackedParameter<bool>("debug");
44 
45  // the name of the simhit collection
46  simHitLabel = pset.getUntrackedParameter<InputTag>("simHitLabel");
47  // the name of the 4D rec hit collection
48  segment4DLabel = pset.getUntrackedParameter<InputTag>("segment4DLabel");
49 
50  //sigma resolution on position
51  sigmaResPos = pset.getParameter<double>("sigmaResPos");
52  //sigma resolution on angle
53  sigmaResAngle = pset.getParameter<double>("sigmaResAngle");
54  doall = pset.getUntrackedParameter<bool>("doall", false);
55  local = pset.getUntrackedParameter<bool>("local", false);
56 
57  // Create the root file
58  //theFile = new TFile(rootFileName.c_str(), "RECREATE");
59  //theFile->cd();
60  // ----------------------
61  // get hold of back-end interface
62  dbe_ = 0;
63  dbe_ = Service<DQMStore>().operator->();
64  if ( dbe_ ) {
65  if (debug) {
66  dbe_->setVerbose(1);
67  } else {
68  dbe_->setVerbose(0);
69  }
70  }
71  if ( dbe_ ) {
72  if ( debug ) dbe_->showDirStructure();
73  }
74 
75  // Book the histos
76  h2DHitSuperPhi = new HRes2DHit ("SuperPhi",dbe_,doall,local);
77  if(doall) h2DHitEff_SuperPhi = new HEff2DHit ("SuperPhi",dbe_);
78 }
79 
80 // Destructor
82 }
83 
85  edm::EventSetup const& c){
86 }
87 
88 
89 
90 
92  // Write the histos to file
93  //theFile->cd();
94  //h2DHitSuperPhi->Write();
95 
96  if(doall) h2DHitEff_SuperPhi->ComputeEfficiency();
97  //h2DHitEff_SuperPhi->Write();
98 
99  //if ( rootFileName.size() != 0 && dbe_ ) dbe_->save(rootFileName);
100  //theFile->Close();
101 }
102 
103 // The real analysis
104 void DTSegment2DSLPhiQuality::analyze(const Event & event, const EventSetup& eventSetup){
105  //theFile->cd();
106 
107  // Get the DT Geometry
108  ESHandle<DTGeometry> dtGeom;
109  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
110 
111  // Get the SimHit collection from the event
113  event.getByLabel(simHitLabel, simHits); //FIXME: second string to be removed
114 
115  //Map simHits by chamber
116  map<DTChamberId, PSimHitContainer > simHitsPerCh;
117  for(PSimHitContainer::const_iterator simHit = simHits->begin();
118  simHit != simHits->end(); simHit++){
119  // Create the id of the chamber (the simHits in the DT known their wireId)
120  DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
121  // Fill the map
122  simHitsPerCh[chamberId].push_back(*simHit);
123  }
124 
125  // Get the 4D rechits from the event
127  event.getByLabel(segment4DLabel, segment4Ds);
128 
129  if(!segment4Ds.isValid()) {
130  if(debug) cout << "[DTSegment2DSLPhiQuality]**Warning: no 4D Segments with label: " << segment4DLabel
131  << " in this event, skipping!" << endl;
132  return;
133  }
134 
135  // Loop over all chambers containing a segment
136  DTRecSegment4DCollection::id_iterator chamberId;
137  for (chamberId = segment4Ds->id_begin();
138  chamberId != segment4Ds->id_end();
139  ++chamberId){
140 
141  //------------------------- simHits ---------------------------//
142  //Get simHits of each chamber
143  PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
144 
145  // Map simhits per wire
146  map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
147  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
148  int nMuSimHit = muSimHitPerWire.size();
149  if(nMuSimHit == 0 || nMuSimHit == 1) {
150  if(debug && nMuSimHit == 1)
151  cout << "[DTSegment2DSLPhiQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
152  continue; // If no or only one mu SimHit is found skip this chamber
153  }
154  if(debug)
155  cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
156 
157  //Find outer and inner mu SimHit to build a segment
158  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
159 
160  //Find direction and position of the sim Segment in Chamber RF
161  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
162  (*chamberId),&(*dtGeom));
163 
164  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
165  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
166  const DTChamber* chamber = dtGeom->chamber(*chamberId);
167  GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
168 
169  //Atan(x/z) angle and x position in SL RF
170  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
171  float posSimSeg = simSegmLocalPos.x();
172  //Position (in eta,phi coordinates) in lobal RF
173  float etaSimSeg = simSegmGlobalPos.eta();
174  float phiSimSeg = simSegmGlobalPos.phi();
175 
176  if(debug)
177  cout<<" Simulated segment: local direction "<<simSegmLocalDir<<endl
178  <<" local position "<<simSegmLocalPos<<endl
179  <<" angle "<<angleSimSeg<<endl;
180 
181  //---------------------------- recHits --------------------------//
182  // Get the range of rechit for the corresponding chamberId
183  bool recHitFound = false;
184  DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
185  int nsegm = distance(range.first, range.second);
186  if(debug)
187  cout << " Chamber: " << *chamberId << " has " << nsegm
188  << " 4D segments" << endl;
189 
190  if (nsegm!=0) {
191  // Find the best RecHit: look for the 4D RecHit with the phi angle closest
192  // to that of segment made of SimHits.
193  // RecHits must have delta alpha and delta position within 5 sigma of
194  // the residual distribution (we are looking for residuals of segments
195  // usefull to the track fit) for efficency purpose
196  const DTRecSegment2D* bestRecHit = 0;
197  bool bestRecHitFound = false;
198  double deltaAlpha = 99999;
199 
200  // Loop over the recHits of this chamberId
201  for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
202  segment4D!=range.second;
203  ++segment4D){
204  // Check the dimension
205  if((*segment4D).dimension() != 4) {
206  if(debug) cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
207  continue;
208  }
209 
210  //Get 2D superPhi segments from 4D segments
211  const DTChamberRecSegment2D* phiSegment2D = (*segment4D).phiSegment();
212  if((*phiSegment2D).dimension() != 2) {
213  if(debug) cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
214  abort();
215  }
216 
217  // Segment Local Direction and position (in Chamber RF)
218  LocalVector recSegDirection = (*phiSegment2D).localDirection();
219 
220  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
221  if(debug)
222  cout << " RecSegment direction: " << recSegDirection << endl
223  << " position : " << (*phiSegment2D).localPosition() << endl
224  << " alpha : " << recSegAlpha << endl;
225 
226  if(fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
227  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
228  bestRecHit = &(*phiSegment2D);
229  bestRecHitFound = true;
230  }
231  } // End of Loop over all 4D RecHits of this chambers
232 
233  if(bestRecHitFound) {
234  // Best rechit direction and position in Chamber RF
235  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
236  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
237 
238  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
239  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
240 
241  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
242  if(fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle &&
243  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos) {
244  recHitFound = true;
245  }
246 
247  // Fill Residual histos
248  h2DHitSuperPhi->Fill(angleSimSeg,
249  angleBestRHit,
250  posSimSeg,
251  bestRecHitLocalPos.x(),
252  etaSimSeg,
253  phiSimSeg,
254  sqrt(bestRecHitLocalPosErr.xx()),
255  sqrt(bestRecHitLocalDirErr.xx())
256  );
257  }
258  } //end of if(nsegm!=0)
259 
260  // Fill Efficiency plot
261  if(doall) {h2DHitEff_SuperPhi->Fill(etaSimSeg,
262  phiSimSeg,
263  posSimSeg,
264  angleSimSeg,
265  recHitFound);}
266  } // End of loop over chambers
267 }
268 
269 
270 // Fit a histogram in the range (minfit, maxfit) with a gaussian and
271 // draw it in the range (min, max)
T getParameter(std::string const &) const
virtual LocalError localPositionError() const
local position error in SL frame
DTSegment2DSLPhiQuality(const edm::ParameterSet &pset)
Constructor.
T getUntrackedParameter(std::string const &, T const &) const
float xx() const
Definition: LocalError.h:24
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup)
Perform the real analysis.
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
virtual LocalError localDirectionError() const
the local direction error (xx,xy,yy) in SL frame: only xx is not 0.
static std::map< DTWireId, const PSimHit * > mapMuSimHitsPerWire(const std::map< DTWireId, edm::PSimHitContainer > &simHitWireMap)
Create a map between the Mu SimHits and corresponding MuBarWireId ;.
static std::pair< const PSimHit *, const PSimHit * > findMuSimSegment(const std::map< DTWireId, const PSimHit * > &mapWireAndMuSimHit)
Find Innermost and outermost SimHit from Mu in a SL (they identify a simulated segment) ...
T sqrt(T t)
Definition: SSEVec.h:48
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:76
static std::pair< double, double > findSegmentAlphaAndBeta(const LocalVector &direction)
Find the angles from a segment direction:
DQMStore * dbe_
virtual LocalPoint localPosition() const
local position in SL frame
#define debug
Definition: HDRShower.cc:19
static std::pair< LocalVector, LocalPoint > findMuSimSegmentDirAndPos(const std::pair< const PSimHit *, const PSimHit * > &inAndOutSimHit, const DetId detId, const DTGeometry *muonGeom)
Find direction and position of a segment (in local RF) from outer and inner mu SimHit in the RF of ob...
const T & get() const
Definition: EventSetup.h:55
tuple simHits
Definition: trackerHits.py:16
T eta() const
Definition: PV3DBase.h:76
void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c)
virtual LocalVector localDirection() const
the local direction in SL frame
tuple cout
Definition: gather_cfg.py:121
std::vector< PSimHit > PSimHitContainer
static std::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
T x() const
Definition: PV3DBase.h:62
virtual ~DTSegment2DSLPhiQuality()
Destructor.