RecoMuon
TrackerSeedGenerator
src
L1MuonRegionProducer.cc
Go to the documentation of this file.
1
#include "
RecoMuon/TrackerSeedGenerator/interface/L1MuonRegionProducer.h
"
2
3
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
4
#include "
RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h
"
5
#include "
RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h
"
6
#include "
DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTCand.h
"
7
#include "
RecoMuon/TrackerSeedGenerator/interface/L1MuonPixelTrackFitter.h
"
8
9
using namespace
std
;
10
11
L1MuonRegionProducer::L1MuonRegionProducer
(
const
edm::ParameterSet
&
cfg
) {
12
edm::ParameterSet
regionPSet =
cfg
.getParameter<
edm::ParameterSet
>(
"RegionPSet"
);
13
14
thePtMin = regionPSet.
getParameter
<
double
>(
"ptMin"
);
15
theOriginRadius = regionPSet.
getParameter
<
double
>(
"originRadius"
);
16
theOriginHalfLength = regionPSet.
getParameter
<
double
>(
"originHalfLength"
);
17
theOrigin =
GlobalPoint
(regionPSet.
getParameter
<
double
>(
"originXPos"
),
18
regionPSet.
getParameter
<
double
>(
"originYPos"
),
19
regionPSet.
getParameter
<
double
>(
"originZPos"
));
20
}
21
22
void
L1MuonRegionProducer::setL1Constraint
(
const
L1MuGMTCand
&
muon
) {
23
thePhiL1 =
muon
.phiValue() + 0.021817;
24
theEtaL1 =
muon
.etaValue();
25
theChargeL1 =
muon
.charge();
26
}
27
28
std::vector<std::unique_ptr<TrackingRegion> >
L1MuonRegionProducer::regions
()
const
{
29
double
dx
=
cos
(thePhiL1);
30
double
dy
=
sin
(thePhiL1);
31
double
dz
= sinh(theEtaL1);
32
GlobalVector
direction(
dx
,
dy
,
dz
);
// muon direction
33
34
std::vector<std::unique_ptr<TrackingRegion> >
result
;
35
double
bending =
L1MuonPixelTrackFitter::getBending
(1. / thePtMin, theEtaL1, theChargeL1);
36
bending = fabs(bending);
37
double
errBending =
L1MuonPixelTrackFitter::getBendingError
(1. / thePtMin, theEtaL1);
38
39
result
.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
40
direction, theOrigin, thePtMin, theOriginRadius, theOriginHalfLength, 0.15, bending + 3 * errBending));
41
42
return
result
;
43
}
Vector3DBase
Definition:
Vector3DBase.h:8
GlobalTrackingRegion.h
muon
Definition:
MuonCocktails.h:17
L1MuonRegionProducer::L1MuonRegionProducer
L1MuonRegionProducer(const edm::ParameterSet &cfg)
Definition:
L1MuonRegionProducer.cc:11
L1MuGMTCand
Definition:
L1MuGMTCand.h:39
funct::sin
Sin< T >::type sin(const T &t)
Definition:
Sin.h:22
funct::cos
Cos< T >::type cos(const T &t)
Definition:
Cos.h:22
L1MuonPixelTrackFitter::getBending
static double getBending(double invPt, double eta, int charge)
Definition:
L1MuonPixelTrackFitter.cc:755
L1MuonRegionProducer::regions
std::vector< std::unique_ptr< TrackingRegion > > regions() const
Definition:
L1MuonRegionProducer.cc:28
GlobalPoint
Global3DPoint GlobalPoint
Definition:
GlobalPoint.h:10
L1MuonPixelTrackFitter::getBendingError
static double getBendingError(double invPt, double eta)
Definition:
L1MuonPixelTrackFitter.cc:761
edm::ParameterSet
Definition:
ParameterSet.h:36
L1MuonPixelTrackFitter.h
RectangularEtaPhiTrackingRegion.h
PVValHelper::dy
Definition:
PVValidationHelpers.h:49
looper.cfg
cfg
Definition:
looper.py:297
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
L1MuGMTCand.h
std
Definition:
JetResolutionObject.h:76
PVValHelper::dz
Definition:
PVValidationHelpers.h:50
mps_fire.result
result
Definition:
mps_fire.py:303
ParameterSet.h
L1MuonRegionProducer.h
PVValHelper::dx
Definition:
PVValidationHelpers.h:48
L1MuonRegionProducer::setL1Constraint
void setL1Constraint(const L1MuGMTCand &muon)
Definition:
L1MuonRegionProducer.cc:22
Generated for CMSSW Reference Manual by
1.8.16