CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
MuonTrackingRegionBuilder Class Reference

#include <MuonTrackingRegionBuilder.h>

Public Member Functions

void init (const MuonServiceProxy *)
 
 MuonTrackingRegionBuilder (const edm::ParameterSet &)
 constructor More...
 
 MuonTrackingRegionBuilder (const edm::ParameterSet &par, const MuonServiceProxy *service)
 
RectangularEtaPhiTrackingRegionregion (const reco::TrackRef &) const
 define tracking region More...
 
RectangularEtaPhiTrackingRegionregion (const reco::Track &) const
 define tracking region More...
 
virtual void setEvent (const edm::Event &)
 pass the Event to the algo at each event More...
 
virtual ~MuonTrackingRegionBuilder ()
 destructor More...
 

Private Member Functions

void build (const edm::ParameterSet &)
 

Private Attributes

edm::InputTag theBeamSpotTag
 
double theDeltaR
 
double theEtaFixed
 
double theEtaMin
 
double theEtaRegionPar1
 
double theEtaRegionPar2
 
const edm::EventtheEvent
 
double theHalfZ
 
std::string theMeasurementTrackerName
 
double theNsigmaDz
 
double theNsigmaEta
 
double theNsigmaPhi
 
double theOnDemand
 
double thePhiFixed
 
double thePhiMin
 
double thePhiRegionPar1
 
double thePhiRegionPar2
 
const MuonServiceProxytheService
 
double theTkEscapePt
 
edm::InputTag theVertexCollTag
 
GlobalPoint theVertexPos
 
bool useFixedRegion
 
bool useVertex
 

Detailed Description

Build a TrackingRegion around a standalone muon

Date:
2010/09/06 18:41:59
Revision:
1.11
Author
A. Everett - Purdue University
A. Grelli - Purdue University, Pavia University

Definition at line 30 of file MuonTrackingRegionBuilder.h.

Constructor & Destructor Documentation

MuonTrackingRegionBuilder::MuonTrackingRegionBuilder ( const edm::ParameterSet par)

constructor

Definition at line 43 of file MuonTrackingRegionBuilder.cc.

References newFWLiteAna::build.

44 {
45  build(par);
46 }
void build(const edm::ParameterSet &)
MuonTrackingRegionBuilder::MuonTrackingRegionBuilder ( const edm::ParameterSet par,
const MuonServiceProxy service 
)
inline

Definition at line 36 of file MuonTrackingRegionBuilder.h.

References build(), and init().

37  { build(par);init(service);}
void build(const edm::ParameterSet &)
void init(const MuonServiceProxy *)
virtual MuonTrackingRegionBuilder::~MuonTrackingRegionBuilder ( )
inlinevirtual

destructor

Definition at line 41 of file MuonTrackingRegionBuilder.h.

41 {}

Member Function Documentation

void MuonTrackingRegionBuilder::build ( const edm::ParameterSet par)
private

Definition at line 47 of file MuonTrackingRegionBuilder.cc.

References edm::ParameterSet::getParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by MuonTrackingRegionBuilder().

47  {
48  // vertex Collection and Beam Spot
49  theBeamSpotTag = par.getParameter<edm::InputTag>("beamSpot");
50  theVertexCollTag = par.getParameter<edm::InputTag>("vertexCollection");
51 
52  // parmeters
53  theNsigmaEta = par.getParameter<double>("Rescale_eta");
54  theNsigmaPhi = par.getParameter<double>("Rescale_phi");
55  theNsigmaDz = par.getParameter<double>("Rescale_Dz");
56  theTkEscapePt = par.getParameter<double>("EscapePt");
57 
58  // upper limits parameters
59  theEtaRegionPar1 = par.getParameter<double>("EtaR_UpperLimit_Par1");
60  theEtaRegionPar2 = par.getParameter<double>("EtaR_UpperLimit_Par2");
61  thePhiRegionPar1 = par.getParameter<double>("PhiR_UpperLimit_Par1");
62  thePhiRegionPar2 = par.getParameter<double>("PhiR_UpperLimit_Par2");
63 
64  useVertex = par.getParameter<bool>("UseVertex");
65  useFixedRegion = par.getParameter<bool>("UseFixedRegion");
66 
67  // fixed limits
68  thePhiFixed = par.getParameter<double>("Phi_fixed");
69  theEtaFixed = par.getParameter<double>("Eta_fixed");
70 
71  thePhiMin = par.getParameter<double>("Phi_min");
72  theEtaMin = par.getParameter<double>("Eta_min");
73  theHalfZ = par.getParameter<double>("DeltaZ_Region");
74  theDeltaR = par.getParameter<double>("DeltaR");
75 
76  // perigee reference point
77 
78  theOnDemand = par.getParameter<double>("OnDemand");
79  theMeasurementTrackerName = par.getParameter<std::string>("MeasurementTrackerName");
80 }
T getParameter(std::string const &) const
void MuonTrackingRegionBuilder::init ( const MuonServiceProxy service)
RectangularEtaPhiTrackingRegion * MuonTrackingRegionBuilder::region ( const reco::TrackRef track) const

define tracking region

Definition at line 87 of file MuonTrackingRegionBuilder.cc.

Referenced by GlobalTrajectoryBuilderBase::defineRegionOfInterest(), TSGFromL2Muon::produce(), FastTSGFromL2Muon::produce(), and HIMuonTrackingRegionProducer::regions().

87  {
88 
89  return region(*track);
90 
91 }
RectangularEtaPhiTrackingRegion * region(const reco::TrackRef &) const
define tracking region
RectangularEtaPhiTrackingRegion * MuonTrackingRegionBuilder::region ( const reco::Track staTrack) const

define tracking region

Definition at line 108 of file MuonTrackingRegionBuilder.cc.

References deltaR(), eta(), PerigeeConversions::ftsToPerigeeError(), trajectoryStateTransform::initialFreeState(), reco::Track::innerMomentum(), TrajectoryStateClosestToBeamLine::isValid(), LogDebug, max(), min, gsfElectronCkfTrackCandidateMaker_cff::minPt, FreeTrajectoryState::momentum(), reco::TrackBase::p(), phi, PerigeeTrajectoryError::phiError(), funct::sin(), theta(), PerigeeTrajectoryError::thetaError(), TrajectoryStateClosestToBeamLine::trackStateAtPCA(), GoodVertex_cfg::vertexCollection, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().

108  {
109 
110  // get the free trajectory state of the muon updated at vertex
111  TSCPBuilderNoMaterial tscpBuilder;
112 
114 
115  LogDebug("MuonTrackingRegionBuilder")<<"from state: "<<muFTS;
116 
117  // get track direction at vertex
118  GlobalVector dirVector(muFTS.momentum());
119 
120  // get track momentum
121  const math::XYZVector& mo = staTrack.innerMomentum();
122  GlobalVector mom(mo.x(),mo.y(),mo.z());
123  if ( staTrack.p() > 1.5 ) {
124  mom = dirVector;
125  }
126 
127  // initial vertex position - in the following it is replaced with beam spot/vertexing
128  GlobalPoint vertexPos(0.0,0.0,0.0);
129  double deltaZatVTX = 0.0;
130 
131  // retrieve beam spot information
133  bool bsHandleFlag = theEvent->getByLabel(theBeamSpotTag, bsHandle);
134  // check the validity, otherwise vertexing
135  // inizialization of BS
136 
137  if ( bsHandleFlag && !useVertex ) {
138  const reco::BeamSpot& bs = *bsHandle;
139  vertexPos = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
140  } else {
141  // get originZPos from list of reconstructed vertices (first or all)
143  bool vtxHandleFlag = theEvent->getByLabel(theVertexCollTag,vertexCollection);
144  // check if there exists at least one reconstructed vertex
145  if ( vtxHandleFlag && !vertexCollection->empty() ) {
146  // use the first vertex in the collection and assume it is the primary event vertex
147  reco::VertexCollection::const_iterator vtx = vertexCollection->begin();
148  vertexPos = GlobalPoint(vtx->x(),vtx->y(),vtx->z());
149  // delta Z from vertex error
150  deltaZatVTX = vtx->zError() * theNsigmaDz;
151  }
152  }
153 
154 
155  // inizialize to the maximum possible value to avoit
156  // problems with TSCBL
157 
158  double deta = 0.4;
159  double dphi = 0.6;
160 
161  // take into account the correct beanspot rotation
162  if ( bsHandleFlag ) {
163 
164  const reco::BeamSpot& bs = *bsHandle;
165 
166  TSCBLBuilderNoMaterial tscblBuilder;
167  TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(muFTS,bs);
168 
169 
170  // evaluate the dynamical region if possible
171  if(tscbl.isValid()){
172 
174  GlobalVector pTrack = tscbl.trackStateAtPCA().momentum();
175 
176  // calculate deltaEta from deltaTheta
177  double deltaTheta = trackPerigeeErrors.thetaError();
178  double theta = pTrack.theta();
179  double sin_theta = sin(theta);
180 
181  // get dEta and dPhi
182  deta = theNsigmaEta*(1/fabs(sin_theta))*deltaTheta;
183  dphi = theNsigmaPhi*(trackPerigeeErrors.phiError());
184 
185  }
186  }
187 
188  /* Region_Parametrizations to take into account possible
189  L2 error matrix inconsistencies. Detailed Explanation in TWIKI
190  page.
191  */
192  double region_dEta = 0;
193  double region_dPhi = 0;
194  double eta=0; double phi=0;
195 
196  // eta, pt parametrization from MC study
197  float pt = fabs(mom.perp());
198 
199  if ( pt <= 10. ) {
200  // angular coefficients
201  float acoeff_Phi = (thePhiRegionPar2 - thePhiRegionPar1)/5;
202  float acoeff_Eta = (theEtaRegionPar2 - theEtaRegionPar1)/5;
203 
204  eta = theEtaRegionPar1 + (acoeff_Eta)*(mom.perp()-5.);
205  phi = thePhiRegionPar1 + (acoeff_Phi)*(mom.perp()-5.) ;
206  }
207  // parametrization 2nd bin in pt from MC study
208  if ( pt > 10. && pt < 100. ) {
209  eta = theEtaRegionPar2;
210  phi = thePhiRegionPar2;
211  }
212  // parametrization 3rd bin in pt from MC study
213  if ( pt >= 100. ) {
214  // angular coefficients
215  float acoeff_Phi = (thePhiRegionPar1 - thePhiRegionPar2)/900;
216  float acoeff_Eta = (theEtaRegionPar1 - theEtaRegionPar2)/900;
217 
218  eta = theEtaRegionPar2 + (acoeff_Eta)*(mom.perp()-100.);
219  phi = thePhiRegionPar2 + (acoeff_Phi)*(mom.perp()-100.);
220  }
221 
222  // decide to use either a parametrization or a dynamical region
223  double region_dPhi1 = min(phi,dphi);
224  double region_dEta1 = min(eta,deta);
225 
226  // minimum size
227  region_dPhi = max(thePhiMin,region_dPhi1);
228  region_dEta = max(theEtaMin,region_dEta1);
229 
230  float deltaZ = 0.0;
231  // standard 15.9 is useVertex than region from vertexing
232  if ( useVertex ) {
233  deltaZ = max(theHalfZ,deltaZatVTX);
234  } else {
235  deltaZ = theHalfZ;
236  }
237 
238  float deltaR = theDeltaR;
239  double minPt = max(theTkEscapePt,mom.perp()*0.6);
240 
242 
243  if (useFixedRegion) {
244  region_dEta = theEtaFixed;
245  region_dPhi = thePhiFixed;
246  }
247 
248  region = new RectangularEtaPhiTrackingRegion(dirVector, vertexPos,
249  minPt, deltaR,
250  deltaZ, region_dEta, region_dPhi,
251  theOnDemand,
252  true,/*default in the header*/
254 
255  LogDebug("MuonTrackingRegionBuilder")<<"the region parameters are:\n"
256  <<"\n dirVector: "<<dirVector
257  <<"\n vertexPos: "<<vertexPos
258  <<"\n minPt: "<<minPt
259  <<"\n deltaR:"<<deltaR
260  <<"\n deltaZ:"<<deltaZ
261  <<"\n region_dEta:"<<region_dEta
262  <<"\n region_dPhi:"<<region_dPhi
263  <<"\n on demand parameter: "<<theOnDemand;
264 
265 
266  return region;
267 
268 }
#define LogDebug(id)
RectangularEtaPhiTrackingRegion * region(const reco::TrackRef &) const
define tracking region
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
double z0() const
z coordinate
Definition: BeamSpot.h:69
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Geom::Theta< T > theta() const
#define min(a, b)
Definition: mlp_lapack.h:161
T eta() const
const MuonServiceProxy * theService
tuple vertexCollection
const T & max(const T &a, const T &b)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
GlobalVector momentum() const
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field)
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:45
double y0() const
y coordinate
Definition: BeamSpot.h:67
PerigeeTrajectoryError ftsToPerigeeError(const FTS &originalFTS)
double x0() const
x coordinate
Definition: BeamSpot.h:65
Definition: DDAxes.h:10
void MuonTrackingRegionBuilder::setEvent ( const edm::Event event)
virtual

pass the Event to the algo at each event

Definition at line 97 of file MuonTrackingRegionBuilder.cc.

References event().

Referenced by TSGFromL2Muon::produce(), FastTSGFromL2Muon::produce(), HIMuonTrackingRegionProducer::regions(), and GlobalTrajectoryBuilderBase::setEvent().

97  {
98 
99  theEvent = &event;
100 
101 }
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

Member Data Documentation

edm::InputTag MuonTrackingRegionBuilder::theBeamSpotTag
private

Definition at line 55 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theDeltaR
private

Definition at line 74 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theEtaFixed
private

Definition at line 75 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theEtaMin
private

Definition at line 73 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theEtaRegionPar1
private

Definition at line 67 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theEtaRegionPar2
private

Definition at line 68 of file MuonTrackingRegionBuilder.h.

const edm::Event* MuonTrackingRegionBuilder::theEvent
private

Definition at line 58 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theHalfZ
private

Definition at line 74 of file MuonTrackingRegionBuilder.h.

std::string MuonTrackingRegionBuilder::theMeasurementTrackerName
private

Definition at line 80 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theNsigmaDz
private

Definition at line 65 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theNsigmaEta
private

Definition at line 65 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theNsigmaPhi
private

Definition at line 65 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theOnDemand
private

Definition at line 79 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::thePhiFixed
private

Definition at line 75 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::thePhiMin
private

Definition at line 72 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::thePhiRegionPar1
private

Definition at line 69 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::thePhiRegionPar2
private

Definition at line 70 of file MuonTrackingRegionBuilder.h.

const MuonServiceProxy* MuonTrackingRegionBuilder::theService
private

Definition at line 59 of file MuonTrackingRegionBuilder.h.

double MuonTrackingRegionBuilder::theTkEscapePt
private

Definition at line 64 of file MuonTrackingRegionBuilder.h.

edm::InputTag MuonTrackingRegionBuilder::theVertexCollTag
private

Definition at line 56 of file MuonTrackingRegionBuilder.h.

GlobalPoint MuonTrackingRegionBuilder::theVertexPos
private

Definition at line 77 of file MuonTrackingRegionBuilder.h.

bool MuonTrackingRegionBuilder::useFixedRegion
private

Definition at line 61 of file MuonTrackingRegionBuilder.h.

bool MuonTrackingRegionBuilder::useVertex
private

Definition at line 62 of file MuonTrackingRegionBuilder.h.